53 int main(
int argc,
char *argv[])
58 const char *mesh_file =
"../data/star.mesh";
61 const char *device_config =
"cpu";
62 bool visualization = 1;
65 args.
AddOption(&mesh_file,
"-m",
"--mesh",
68 "Finite element order (polynomial degree).");
69 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
70 "--no-partial-assembly",
"Enable Partial Assembly.");
71 args.
AddOption(&device_config,
"-d",
"--device",
72 "Device configuration string, see Device::Configure().");
73 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
75 "Enable or disable GLVis visualization.");
86 Device device(device_config);
92 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
101 (int)floor(log(10000./mesh->
GetNE())/log(2.)/
dim);
102 for (
int l = 0; l < ref_levels; l++)
120 block_offsets[0] = 0;
121 block_offsets[1] = R_space->
GetVSize();
122 block_offsets[2] = W_space->
GetVSize();
125 std::cout <<
"***********************************************************\n";
126 std::cout <<
"dim(R) = " << block_offsets[1] - block_offsets[0] <<
"\n";
127 std::cout <<
"dim(W) = " << block_offsets[2] - block_offsets[1] <<
"\n";
128 std::cout <<
"dim(R+W) = " << block_offsets.
Last() <<
"\n";
129 std::cout <<
"***********************************************************\n";
147 BlockVector x(block_offsets, mt), rhs(block_offsets, mt);
225 auto Md_host = Md.HostRead();
227 for (
int i=0; i<mVarf->
Height(); ++i)
229 invMd(i) = 1.0 / Md_host[i];
249 for (
int i = 0; i < Md.Size(); i++)
251 MinvBt->ScaleRow(i, 1./Md(i));
254 S =
Mult(B, *MinvBt);
258 #ifndef MFEM_USE_SUITESPARSE
268 darcyPrec.SetDiagonalBlock(0, invM);
269 darcyPrec.SetDiagonalBlock(1, invS);
288 if (device.
IsEnabled()) { x.HostRead(); }
294 <<
" iterations with a residual norm of "
300 <<
" iterations. Residual norm is " << solver.
GetFinalNorm()
303 std::cout <<
"MINRES solver took " << chrono.
RealTime() <<
"s.\n";
307 u.
MakeRef(R_space, x.GetBlock(0), 0);
308 p.
MakeRef(W_space, x.GetBlock(1), 0);
310 int order_quad = max(2, 2*order+1);
312 for (
int i=0; i < Geometry::NumGeom; ++i)
322 std::cout <<
"|| u_h - u_ex || / || u_ex || = " << err_u / norm_u <<
"\n";
323 std::cout <<
"|| p_h - p_ex || / || p_ex || = " << err_p / norm_p <<
"\n";
329 ofstream mesh_ofs(
"ex5.mesh");
330 mesh_ofs.precision(8);
331 mesh->
Print(mesh_ofs);
333 ofstream u_ofs(
"sol_u.gf");
337 ofstream p_ofs(
"sol_p.gf");
367 u_sock <<
"solution\n" << *mesh << u <<
"window_title 'Velocity'" << endl;
370 p_sock <<
"solution\n" << *mesh << p <<
"window_title 'Pressure'" << endl;
403 u(0) = - exp(xi)*sin(yi)*cos(zi);
404 u(1) = - exp(xi)*cos(yi)*cos(zi);
408 u(2) = exp(xi)*sin(yi)*sin(zi);
424 return exp(xi)*sin(yi)*cos(zi);
Class for domain integration L(v) := (f, v)
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
virtual void Print(std::ostream &out=mfem::out) const
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
int GetNumIterations() const
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Data type for scaled Jacobi-type smoother of sparse matrix.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void fFun(const Vector &x, Vector &f)
A class to handle Vectors in a block fashion.
void SetDataFormat(VTKFormat fmt)
A coefficient that is constant across space and time.
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Helper class for ParaView visualization data.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
virtual void Save()
Save the collection and a VisIt root file.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
static bool IsEnabled()
Return true if any backend other than Backend::CPU is enabled.
double GetFinalNorm() const
Data type for Gauss-Seidel smoother of sparse matrix.
void Stop()
Stop the stopwatch.
Direct sparse solver using UMFPACK.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
Compute the Lp norm of a function f. .
void SetPrintLevel(int print_lvl)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Jacobi smoothing for a given bilinear form (no matrix necessary).
double f_natural(const Vector &x)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
double pFun_ex(const Vector &x)
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
void SetHighOrderOutput(bool high_order_output_)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void PrintUsage(std::ostream &out) const
Print the usage message.
void SetTime(double t)
Set physical time (for time-dependent simulations)
void Start()
Clear the elapsed time and start the stopwatch.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
A general vector function coefficient.
void SetAbsTol(double atol)
double p(const Vector &x, double t)
void SetRelTol(double rtol)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
MemoryType
Memory types supported by MFEM.
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 PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
The transpose of a given operator. Switches the roles of the methods Mult() and MultTranspose().
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
void GetDiag(Vector &d) const
Returns the Diagonal of A.
void uFun_ex(const Vector &x, Vector &u)
double gFun(const Vector &x)
void PrintOptions(std::ostream &out) const
Print the options.
T & Last()
Return the last element in the array.
void SetLevelsOfDetail(int levels_of_detail_)
for VectorFiniteElements (Nedelec, Raviart-Thomas)
A general function coefficient.
virtual void Save() override
double u(const Vector &xvec)
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
A class to handle Block systems in a matrix-free implementation.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
Vector & GetBlock(int i)
Get the i-th vector in the block.
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.