45 int main(
int argc,
char *argv[])
53 const char *mesh_file =
"../data/inline-quad.mesh";
54 int ser_ref_levels = 2;
55 int par_ref_levels = 1;
58 bool visualization = 1;
61 args.
AddOption(&mesh_file,
"-m",
"--mesh",
63 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
64 "Number of times to refine the mesh uniformly in serial.");
65 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
66 "Number of times to refine the mesh uniformly in parallel.");
68 "Finite element order (polynomial degree).");
69 args.
AddOption(&
freq,
"-f",
"--frequency",
"Set the frequency for the exact"
71 args.
AddOption(&use_ams,
"-ams",
"--hypre-ams",
"-slu",
72 "--superlu",
"Use AMS or SuperLU solver.");
73 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
75 "Enable or disable GLVis visualization.");
83 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
89 for (
int lev = 0; lev < ser_ref_levels; lev++)
98 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
100 for (
int lev = 0; lev < par_ref_levels; lev++)
122 if (mpi.
Root()) { cout <<
"Number of H(Curl) unknowns: " << size << endl; }
159 sigmaMat(0,0) = 2.0; sigmaMat(1,1) = 2.0; sigmaMat(2,2) = 2.0;
160 sigmaMat(0,2) = 0.0; sigmaMat(2,0) = 0.0;
161 sigmaMat(0,1) = M_SQRT1_2; sigmaMat(1,0) = M_SQRT1_2;
162 sigmaMat(1,2) = M_SQRT1_2; sigmaMat(2,1) = M_SQRT1_2;
186 cout <<
"Size of linear system: "
200 #ifdef MFEM_USE_SUPERLU
204 cout <<
"Size of linear system: "
215 if (mpi.
Root()) { cout <<
"No solvers available." << endl; }
229 cout <<
"\n|| E_h - E ||_{H(Curl)} = " << error <<
'\n' << endl;
237 ostringstream mesh_name, sol_name;
238 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
239 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
241 ofstream mesh_ofs(mesh_name.str().c_str());
242 mesh_ofs.precision(8);
243 pmesh.
Print(mesh_ofs);
245 ofstream sol_ofs(sol_name.str().c_str());
246 sol_ofs.precision(8);
269 dy_sock.precision(8);
270 dz_sock.precision(8);
272 Vector xVec(3); xVec = 0.0; xVec(0) = 1;
273 Vector yVec(3); yVec = 0.0; yVec(1) = 1;
274 Vector zVec(3); zVec = 0.0; zVec(2) = 1;
300 x_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
301 <<
"solution\n" << pmesh << xComp << flush
302 <<
"window_title 'X component'" << endl;
303 y_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
304 <<
"solution\n" << pmesh << yComp << flush
305 <<
"window_geometry 403 0 400 350 "
306 <<
"window_title 'Y component'" << endl;
307 z_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
308 <<
"solution\n" << pmesh << zComp << flush
309 <<
"window_geometry 806 0 400 350 "
310 <<
"window_title 'Z component'" << endl;
318 dy_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
319 <<
"solution\n" << pmesh << dyComp << flush
320 <<
"window_geometry 403 375 400 350 "
321 <<
"window_title 'Y component of Curl'" << endl;
322 dz_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
323 <<
"solution\n" << pmesh << dzComp << flush
324 <<
"window_geometry 806 375 400 350 "
325 <<
"window_title 'Z component of Curl'" << endl;
335 xyMat(0,0) = 1.0; xyMat(1,1) = 1.0;
337 Vector zVec(3); zVec = 0.0; zVec(2) = 1;
362 xy_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
363 xy_sock.precision(8);
364 xy_sock <<
"solution\n" << pmesh << xyComp
365 <<
"window_title 'XY components'\n" << flush;
366 z_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
367 <<
"solution\n" << pmesh << zComp << flush
368 <<
"window_geometry 403 0 400 350 "
369 <<
"window_title 'Z component'" << endl;
377 dxy_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
378 <<
"solution\n" << pmesh << dxyComp << flush
379 <<
"window_geometry 0 375 400 350 "
380 <<
"window_title 'XY components of Curl'" << endl;
381 dz_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
382 <<
"solution\n" << pmesh << dzComp << flush
383 <<
"window_geometry 403 375 400 350 "
384 <<
"window_title 'Z component of Curl'" << endl;
399 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
400 sol_sock.precision(8);
401 sol_sock <<
"solution\n" << pmesh << sol
402 <<
"window_title 'Solution'" << flush << endl;
403 dsol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
404 <<
"solution\n" << pmesh << dsol << flush
405 <<
"window_geometry 0 375 400 350 "
406 <<
"window_title 'Curl of solution'" << endl;
421 E(0) = 1.1 *
sin(
kappa * x(0) + 0.0 * M_PI);
422 E(1) = 1.2 *
sin(
kappa * x(0) + 0.4 * M_PI);
423 E(2) = 1.3 *
sin(
kappa * x(0) + 0.9 * M_PI);
427 E(0) = 1.1 *
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
428 E(1) = 1.2 *
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
429 E(2) = 1.3 *
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
433 E(0) = 1.1 *
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
434 E(1) = 1.2 *
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
435 E(2) = 1.3 *
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
444 double c4 =
cos(
kappa * x(0) + 0.4 * M_PI);
445 double c9 =
cos(
kappa * x(0) + 0.9 * M_PI);
454 double c0 =
cos(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
455 double c4 =
cos(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
456 double c9 =
cos(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
460 dE(2) = 1.2 * c4 - 1.1 * c0;
461 dE *=
kappa * M_SQRT1_2;
465 double s0 =
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
466 double c0 =
cos(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
467 double s4 =
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
468 double c4 =
cos(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
469 double c9 =
cos(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
473 dE(0) = 1.2 * s4 * sk + 1.3 * M_SQRT1_2 * c9 * ck;
474 dE(1) = -1.1 * s0 * sk - 1.3 * M_SQRT1_2 * c9 * ck;
475 dE(2) = -M_SQRT1_2 * (1.1 * c0 - 1.2 * c4) * ck;
484 double s0 =
sin(
kappa * x(0) + 0.0 * M_PI);
485 double s4 =
sin(
kappa * x(0) + 0.4 * M_PI);
486 double s9 =
sin(
kappa * x(0) + 0.9 * M_PI);
488 f(0) = 2.2 * s0 + 1.2 * M_SQRT1_2 * s4;
490 M_SQRT1_2 * (1.1 * s0 + 1.3 * s9);
491 f(2) = 1.3 * (2.0 +
kappa *
kappa) * s9 + 1.2 * M_SQRT1_2 * s4;
495 double s0 =
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
496 double s4 =
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
497 double s9 =
sin(
kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
501 f(1) = 0.55 * (M_SQRT2 - kappa *
kappa) * s0 +
502 0.6 * (4.0 + kappa * kappa) * s4 +
504 f(2) = 0.6 * M_SQRT2 * s4 + 1.3 * (2.0 + kappa *
kappa) * s9;
508 double s0 =
sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
509 double c0 =
cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
510 double s4 =
sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
511 double c4 =
cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
512 double s9 =
sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
513 double c9 =
cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
514 double sk =
sin(kappa * x(2));
515 double ck =
cos(kappa * x(2));
517 f(0) = 0.55 * (4.0 + 3.0 * kappa *
kappa) * s0 * ck +
518 0.6 * (M_SQRT2 - kappa * kappa) * s4 * ck -
519 0.65 * M_SQRT2 * kappa * kappa * c9 * sk;
521 f(1) = 0.55 * (M_SQRT2 - kappa *
kappa) * s0 * ck +
522 0.6 * (4.0 + 3.0 * kappa * kappa) * s4 * ck +
523 0.65 * M_SQRT2 * s9 * ck -
524 0.65 * M_SQRT2 * kappa * kappa * c9 * sk;
526 f(2) = 0.6 * M_SQRT2 * s4 * ck -
527 M_SQRT2 * kappa * kappa * (0.55 * c0 + 0.6 * c4) * sk
528 + 1.3 * (2.0 + kappa * kappa) * s9 * ck;
int WorldSize() const
Return MPI_COMM_WORLD's size.
A matrix coefficient that is constant in space and time.
int Size() const
Return the logical size of the array.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
The Auxiliary-space Maxwell Solver in hypre.
A coefficient that is constant across space and time.
Scalar coefficient defined as the inner product of two vector coefficients.
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 1D.
Integrator for (curl u, curl v) for Nedelec elements.
Vector coefficient that is constant in space and time.
Pointer to an Operator of a specified type.
HYPRE_BigInt GlobalTrueVSize() const
Data type dense matrix using column-major storage.
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
virtual void Save(std::ostream &out) const
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void SetPrintLevel(int print_lvl)
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
A simple convenience class based on the Mpi singleton class above. Preserved for backward compatibili...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Vector coefficient defined as the Curl of a vector GridFunction.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
bool Root() const
Return true if WorldRank() == 0.
void SetMaxIter(int max_iter)
FDualNumber< tbase > cos(const FDualNumber< tbase > &f)
cos([dual number])
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
A general vector function coefficient.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
int WorldRank() const
Return MPI_COMM_WORLD's rank.
FDualNumber< tbase > sin(const FDualNumber< tbase > &f)
sin([dual number])
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
void f_exact(const Vector &, Vector &)
Vector coefficient defined as a product of a matrix coefficient and a vector coefficient.
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void E_exact(const Vector &, Vector &)
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Arbitrary order H(curl)-conforming Nedelec finite elements.
Vector coefficient defined by a vector GridFunction.
void CurlE_exact(const Vector &, Vector &)
Arbitrary order H1-conforming (continuous) finite elements.
void Print(std::ostream &out=mfem::out) const override
Class for parallel grid function.
Wrapper for hypre's ParCSR matrix class.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
Class for parallel meshes.
void ParseCheck(std::ostream &out=mfem::out)
virtual double ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured H(curl)-norm for ND elements.
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)
double sigma(const Vector &x)