59int main(
int argc,
char *argv[])
62 const char *mesh_file =
"../data/star.mesh";
65 bool static_cond =
false;
66 bool hybridization =
false;
68 const char *device_config =
"cpu";
69 bool visualization = 1;
72 args.
AddOption(&mesh_file,
"-m",
"--mesh",
75 "Finite element order (polynomial degree).");
76 args.
AddOption(&set_bc,
"-bc",
"--impose-bc",
"-no-bc",
"--dont-impose-bc",
77 "Impose or not essential boundary conditions.");
78 args.
AddOption(&
freq,
"-f",
"--frequency",
"Set the frequency for the exact"
80 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
81 "--no-static-condensation",
"Enable static condensation.");
82 args.
AddOption(&hybridization,
"-hb",
"--hybridization",
"-no-hb",
83 "--no-hybridization",
"Enable hybridization.");
84 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
85 "--no-partial-assembly",
"Enable Partial Assembly.");
86 args.
AddOption(&device_config,
"-d",
"--device",
87 "Device configuration string, see Device::Configure().");
88 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
90 "Enable or disable GLVis visualization.");
102 Device device(device_config);
108 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
118 (int)floor(log(25000./mesh->
GetNE())/log(2.)/
dim);
119 for (
int l = 0; l < ref_levels; l++)
129 cout <<
"Number of finite element unknowns: "
140 ess_bdr = set_bc ? 1 : 0;
168 if (pa) {
a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
180 a->EnableStaticCondensation();
182 else if (hybridization)
193 a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
195 cout <<
"Size of linear system: " << A->
Height() << endl;
200#ifndef MFEM_USE_SUITESPARSE
203 PCG(*A, M, B, X, 1, 10000, 1e-20, 0.0);
207 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
209 umf_solver.
Mult(B, X);
217 PCG(*A, M, B, X, 1, 10000, 1e-20, 0.0);
221 CG(*A, B, X, 1, 10000, 1e-20, 0.0);
226 a->RecoverFEMSolution(X, *
b, x);
229 cout <<
"\n|| F_h - F ||_{L^2} = " << x.
ComputeL2Error(F) <<
'\n' << endl;
234 ofstream mesh_ofs(
"refined.mesh");
235 mesh_ofs.precision(8);
236 mesh->
Print(mesh_ofs);
237 ofstream sol_ofs(
"sol.gf");
238 sol_ofs.precision(8);
248 sol_sock.precision(8);
249 sol_sock <<
"solution\n" << *mesh << x << flush;
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
int Size() const
Return the logical size of the array.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
A coefficient that is constant across space and time.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
for Raviart-Thomas elements
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
int SpaceDimension() const
Dimension of the physical space containing the mesh.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Pointer to an Operator of a specified type.
Jacobi smoothing for a given bilinear form (no matrix necessary).
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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...
bool Good() const
Return true if the command line options were parsed successfully.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Direct sparse solver using UMFPACK.
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
real_t Control[UMFPACK_CONTROL]
virtual void Mult(const Vector &b, Vector &x) const
Direct solution of the linear system using UMFPACK.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
A general vector function coefficient.
int Size() const
Returns the size of the vector.
void f_exact(const Vector &, Vector &)
void F_exact(const Vector &, Vector &)
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
real_t p(const Vector &x, real_t t)