54 int main(
int argc,
char *argv[])
59 Mpi::Init(argc, argv);
60 int num_procs = Mpi::WorldSize();
61 int myid = Mpi::WorldRank();
63 bool verbose = (myid == 0);
66 const char *mesh_file =
"../data/star.mesh";
69 bool par_format =
false;
71 const char *device_config =
"cpu";
72 bool visualization = 1;
76 args.
AddOption(&mesh_file,
"-m",
"--mesh",
78 args.
AddOption(&ref_levels,
"-r",
"--refine",
79 "Number of times to refine the mesh uniformly.");
81 "Finite element order (polynomial degree).");
82 args.
AddOption(&par_format,
"-pf",
"--parallel-format",
"-sf",
84 "Format to use when saving the results for VisIt.");
85 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
86 "--no-partial-assembly",
"Enable Partial Assembly.");
87 args.
AddOption(&device_config,
"-d",
"--device",
88 "Device configuration string, see Device::Configure().");
89 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
91 "Enable or disable GLVis visualization.");
92 args.
AddOption(&adios2,
"-adios2",
"--adios2-streams",
"-no-adios2",
93 "--no-adios2-streams",
94 "Save data using adios2 streams.");
111 Device device(device_config);
112 if (myid == 0) { device.
Print(); }
117 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
125 if (ref_levels == -1)
127 ref_levels = (int)floor(log(10000./mesh->
GetNE())/log(2.)/
dim);
130 for (
int l = 0; l < ref_levels; l++)
142 int par_ref_levels = 2;
143 for (
int l = 0; l < par_ref_levels; l++)
162 std::cout <<
"***********************************************************\n";
163 std::cout <<
"dim(R) = " << dimR <<
"\n";
164 std::cout <<
"dim(W) = " << dimW <<
"\n";
165 std::cout <<
"dim(R+W) = " << dimR + dimW <<
"\n";
166 std::cout <<
"***********************************************************\n";
175 block_offsets[0] = 0;
176 block_offsets[1] = R_space->
GetVSize();
177 block_offsets[2] = W_space->
GetVSize();
181 block_trueOffsets[0] = 0;
182 block_trueOffsets[1] = R_space->
TrueVSize();
183 block_trueOffsets[2] = W_space->
TrueVSize();
199 BlockVector x(block_offsets, mt), rhs(block_offsets, mt);
200 BlockVector trueX(block_trueOffsets, mt), trueRhs(block_trueOffsets, mt);
291 for (
int i=0; i<Md_PA.
Size(); ++i)
293 invMd(i) = 1.0 / Md_host[i];
328 int maxIter(pa ? 1000 : 500);
342 solver.
Mult(trueRhs, trueX);
343 if (device.
IsEnabled()) { trueX.HostRead(); }
350 <<
" iterations with a residual norm of " << solver.
GetFinalNorm() <<
".\n";
353 <<
" iterations. Residual norm is " << solver.
GetFinalNorm() <<
".\n";
354 std::cout <<
"MINRES solver took " << chrono.
RealTime() <<
"s. \n";
362 u->MakeRef(R_space, x.GetBlock(0), 0);
363 p->MakeRef(W_space, x.GetBlock(1), 0);
364 u->Distribute(&(trueX.GetBlock(0)));
365 p->Distribute(&(trueX.GetBlock(1)));
367 int order_quad = max(2, 2*order+1);
369 for (
int i=0; i < Geometry::NumGeom; ++i)
374 double err_u =
u->ComputeL2Error(ucoeff, irs);
376 double err_p =
p->ComputeL2Error(pcoeff, irs);
381 std::cout <<
"|| u_h - u_ex || / || u_ex || = " << err_u / norm_u <<
"\n";
382 std::cout <<
"|| p_h - p_ex || / || p_ex || = " << err_p / norm_p <<
"\n";
388 ostringstream mesh_name, u_name, p_name;
389 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
390 u_name <<
"sol_u." << setfill(
'0') << setw(6) << myid;
391 p_name <<
"sol_p." << setfill(
'0') << setw(6) << myid;
393 ofstream mesh_ofs(mesh_name.str().c_str());
394 mesh_ofs.precision(8);
395 pmesh->
Print(mesh_ofs);
397 ofstream u_ofs(u_name.str().c_str());
401 ofstream p_ofs(p_name.str().c_str());
411 DataCollection::SERIAL_FORMAT :
412 DataCollection::PARALLEL_FORMAT);
429 #ifdef MFEM_USE_ADIOS2 432 std::string postfix(mesh_file);
433 postfix.erase(0, std::string(
"../data/").size() );
434 postfix +=
"_o" + std::to_string(order);
435 const std::string collection_name =
"ex5-p_" + postfix +
".bp";
453 u_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
455 u_sock <<
"solution\n" << *pmesh << *
u <<
"window_title 'Velocity'" 461 p_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
463 p_sock <<
"solution\n" << *pmesh << *
p <<
"window_title 'Pressure'" 504 u(0) = - exp(xi)*sin(yi)*cos(zi);
505 u(1) = - exp(xi)*cos(yi)*cos(zi);
509 u(2) = exp(xi)*sin(yi)*sin(zi);
525 return exp(xi)*sin(yi)*cos(zi);
void uFun_ex(const Vector &x, Vector &u)
Class for domain integration L(v) := (f, v)
Class for an integration rule - an Array of IntegrationPoint.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
A class to handle Vectors in a block fashion.
void SetDataFormat(VTKFormat fmt)
A coefficient that is constant across space and time.
void PrintOptions(std::ostream &out) const
Print the options.
int TrueVSize() const
Obsolete, kept for backward compatibility.
void SetSize(int s)
Resize the vector to size s.
Helper class for ParaView visualization data.
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
void PrintUsage(std::ostream &out) const
Print the usage message.
double pFun_ex(const Vector &x)
Pointer to an Operator of a specified type.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int Size() const
Returns the size of the vector.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
bool Good() const
Return true if the command line options were parsed successfully.
Abstract parallel finite element space.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
static bool IsEnabled()
Return true if any backend other than Backend::CPU is enabled.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
double GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult()...
void Stop()
Stop the stopwatch.
double f_natural(const Vector &x)
The BoomerAMG solver in hypre.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
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 main(int argc, char *argv[])
Jacobi preconditioner in hypre.
Jacobi smoothing for a given bilinear form (no matrix necessary).
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
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)
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
void SetHighOrderOutput(bool high_order_output_)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
void SetTime(double t)
Set physical time (for time-dependent simulations)
HYPRE_BigInt GlobalTrueVSize() const
void Start()
Start the stopwatch. The elapsed time is not cleared.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
A general vector function coefficient.
Wrapper for hypre's parallel vector class.
void SetAbsTol(double atol)
double p(const Vector &x, double t)
void SetRelTol(double rtol)
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
virtual void Save() override
Save the collection and a VisIt root file.
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...
bool GetConverged() const
Returns true if the last call to Mult() converged successfully.
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
The transpose of a given operator. Switches the roles of the methods Mult() and MultTranspose().
HypreParMatrix * Transpose() const
Returns the transpose of *this.
int GetNE() const
Returns number of elements.
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
double gFun(const Vector &x)
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
void fFun(const Vector &x, Vector &f)
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Operator * Ptr() const
Access the underlying Operator pointer.
void SetLevelsOfDetail(int levels_of_detail_)
for VectorFiniteElements (Nedelec, Raviart-Thomas)
A general function coefficient.
void SetLevelsOfDetail(const int levels_of_detail) noexcept
virtual void Save() override
void Print(std::ostream &out=mfem::out) const override
double u(const Vector &xvec)
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.
A class to handle Block systems in a matrix-free implementation.
Class for parallel meshes.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
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).
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Vector & GetBlock(int i)
Get the i-th vector in the block.
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
HYPRE_BigInt * GetRowStarts() const
Return the parallel row partitioning array.