42 int main(
int argc,
char *argv[])
45 Mpi::Init(argc, argv);
46 int num_procs = Mpi::WorldSize();
47 int myid = Mpi::WorldRank();
51 const char *mesh_file =
"../data/star.mesh";
53 bool visualization = 1;
56 args.
AddOption(&mesh_file,
"-m",
"--mesh",
59 "Finite element order (polynomial degree).");
60 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
62 "Enable or disable GLVis visualization.");
80 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
90 for (
int l = 0; l < ref_levels; l++)
102 int par_ref_levels = 1;
103 for (
int l = 0; l < par_ref_levels; l++)
117 unsigned int trial_order = order;
118 unsigned int trace_order = order - 1;
119 unsigned int test_order = order;
121 if (dim == 2 && (order%2 == 0 || (pmesh->
MeshGenerator() & 2 && order > 1)))
125 if (test_order < trial_order)
129 cerr <<
"Warning, test space not enriched enough to handle primal"
151 cout <<
"\nNumber of Unknowns:\n"
152 <<
" Trial space, X0 : " << glob_true_s0
153 <<
" (order " << trial_order <<
")\n"
154 <<
" Interface space, Xhat : " << glob_true_s1
155 <<
" (order " << trace_order <<
")\n"
156 <<
" Test space, Y : " << glob_true_s_test
157 <<
" (order " << test_order <<
")\n\n";
213 enum {x0_var, xhat_var, NVAR};
217 int true_s_test = test_space->
TrueVSize();
221 true_offsets[1] = true_s0;
222 true_offsets[2] = true_s0+true_s1;
225 true_offsets_test[0] = 0;
226 true_offsets_test[1] = true_s_test;
244 matSinv->
Mult(*trueF, SinvF);
263 if (dim == 2) { Shatinv =
new HypreAMS(*Shat, xhat_space); }
264 else { Shatinv =
new HypreADS(*Shat, xhat_space); }
285 matSinv->
Mult(LSres, tmp);
289 cout <<
"\n|| B0*x0 + Bhat*xhat - F ||_{S^-1} = " << res << endl;
298 ostringstream mesh_name, sol_name;
299 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
300 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
302 ofstream mesh_ofs(mesh_name.str().c_str());
303 mesh_ofs.precision(8);
304 pmesh->
Print(mesh_ofs);
306 ofstream sol_ofs(sol_name.str().c_str());
307 sol_ofs.precision(8);
317 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
318 sol_sock.precision(8);
319 sol_sock <<
"solution\n" << *pmesh << *x0 << flush;
Class for domain integration L(v) := (f, v)
Conjugate gradient method.
The Auxiliary-space Maxwell Solver in hypre.
Integrator defining a sum of multiple Integrators.
The Auxiliary-space Divergence Solver in hypre.
A class to handle Vectors in a block fashion.
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const
Determine the boundary degrees of freedom.
A coefficient that is constant across space and time.
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int TrueVSize() const
Obsolete, kept for backward compatibility.
HYPRE_BigInt GlobalTrueVSize() const
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator.
int GetNE() const
Returns number of elements.
virtual void Save(std::ostream &out) const
Abstract parallel finite element space.
void AddIntegrator(BilinearFormIntegrator *integ)
Integrator that inverts the matrix assembled by another integrator.
The BoomerAMG solver in hypre.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
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. ...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetPrintLevel(int print_level)
void SetMaxIter(int max_it)
int MeshGenerator()
Get the mesh generator/type.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
The operator x -> R*A*P*x constructed through the actions of R^T, A and P.
void PrintUsage(std::ostream &out) const
Print the usage message.
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Wrapper for hypre's parallel vector class.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void SetRelTol(double rtol)
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 Distribute(const Vector *tv)
FDualNumber< tbase > log(const FDualNumber< tbase > &f)
log([dual number])
double InnerProduct(HypreParVector *x, HypreParVector *y)
FDualNumber< tbase > sqrt(const FDualNumber< tbase > &f)
sqrt([dual number])
void PrintOptions(std::ostream &out) const
Print the options.
Abstract class for hypre's solvers and preconditioners.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the interface between mesh elem...
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A * x + beta * y.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
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.
A class to handle Block systems in a matrix-free implementation.
Class for parallel meshes.
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).
Arbitrary order "L2-conforming" discontinuous finite elements.
bool Good() const
Return true if the command line options were parsed successfully.