57 int main(
int argc,
char *argv[])
61 MPI_Init(&argc, &argv);
62 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
63 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
66 const char *mesh_file =
"../data/star.mesh";
69 bool static_cond =
false;
70 bool hybridization =
false;
72 const char *device_config =
"cpu";
73 bool visualization = 1;
76 args.
AddOption(&mesh_file,
"-m",
"--mesh",
79 "Finite element order (polynomial degree).");
80 args.
AddOption(&set_bc,
"-bc",
"--impose-bc",
"-no-bc",
"--dont-impose-bc",
81 "Impose or not essential boundary conditions.");
82 args.
AddOption(&
freq,
"-f",
"--frequency",
"Set the frequency for the exact"
84 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
85 "--no-static-condensation",
"Enable static condensation.");
86 args.
AddOption(&hybridization,
"-hb",
"--hybridization",
"-no-hb",
87 "--no-hybridization",
"Enable hybridization.");
88 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
89 "--no-partial-assembly",
"Enable Partial Assembly.");
90 args.
AddOption(&device_config,
"-d",
"--device",
91 "Device configuration string, see Device::Configure().");
92 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
94 "Enable or disable GLVis visualization.");
113 Device device(device_config);
114 if (myid == 0) { device.
Print(); }
119 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
129 (int)floor(log(1000./mesh->
GetNE())/log(2.)/
dim);
130 for (
int l = 0; l < ref_levels; l++)
144 int par_ref_levels = 2;
145 for (
int l = 0; l < par_ref_levels; l++)
159 cout <<
"Number of finite element unknowns: " << size << endl;
170 ess_bdr = set_bc ? 1 : 0;
213 else if (hybridization)
226 if (myid == 0 && !pa)
228 cout <<
"Size of linear system: "
263 cout <<
"\n|| F_h - F ||_{L^2} = " << err <<
'\n' << endl;
270 ostringstream mesh_name, sol_name;
271 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
272 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
274 ofstream mesh_ofs(mesh_name.str().c_str());
275 mesh_ofs.precision(8);
276 pmesh->
Print(mesh_ofs);
278 ofstream sol_ofs(sol_name.str().c_str());
279 sol_ofs.precision(8);
289 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
290 sol_sock.precision(8);
291 sol_sock <<
"solution\n" << *pmesh << x << flush;
341 f(0) = temp*cos(kappa*x)*sin(kappa*y);
342 f(1) = temp*cos(kappa*y)*sin(kappa*x);
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)
Conjugate gradient method.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
The Auxiliary-space Maxwell Solver in hypre.
The Auxiliary-space Divergence Solver in hypre.
A coefficient that is constant across space and time.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Pointer to an Operator of a specified type.
HYPRE_BigInt GlobalTrueVSize() const
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
virtual void Save(std::ostream &out) const
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
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...
(Q div u, div v) for RT elements
The BoomerAMG solver in hypre.
void SetPrintLevel(int print_lvl)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Jacobi smoothing for a given bilinear form (no matrix necessary).
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetMaxIter(int max_it)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void F_exact(const Vector &, Vector &)
virtual void Print(std::ostream &out=mfem::out) const
void PrintUsage(std::ostream &out) const
Print the usage message.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
A general vector function coefficient.
int SpaceDimension() const
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
double p(const Vector &x, double t)
void SetRelTol(double rtol)
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
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 &)
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
void PrintOptions(std::ostream &out) const
Print the options.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Class for parallel grid function.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.