40 int main(
int argc,
char *argv[])
44 MPI_Init(&argc, &argv);
45 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
46 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
49 const char *mesh_file =
"../data/star.mesh";
51 bool visualization = 1;
54 args.
AddOption(&mesh_file,
"-m",
"--mesh",
57 "Finite element order (polynomial degree).");
58 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
60 "Enable or disable GLVis visualization.");
80 ifstream imesh(mesh_file);
85 cerr <<
"\nCan not open mesh file: " << mesh_file <<
'\n' << endl;
90 mesh =
new Mesh(imesh, 1, 1);
100 (int)floor(log(10000./mesh->
GetNE())/log(2.)/
dim);
101 for (
int l = 0; l < ref_levels; l++)
113 int par_ref_levels = 1;
114 for (
int l = 0; l < par_ref_levels; l++)
129 int trial_order = order;
130 int trace_order = order - 1;
131 int test_order = order;
132 if (dim == 2 && (order%2 == 0 || (pmesh->
MeshGenerator() & 2 && order > 1)))
136 if (test_order < trial_order)
138 cerr <<
"Warning, test space not enriched enough to handle primal"
158 cout <<
"\nNumber of Unknowns:\n"
159 <<
" Trial space, X0 : " << glob_true_s0
160 <<
" (order " << trial_order <<
")\n"
161 <<
" Interface space, Xhat : " << glob_true_s1
162 <<
" (order " << trace_order <<
")\n"
163 <<
" Test space, Y : " << glob_true_s_test
164 <<
" (order " << test_order <<
")\n\n";
220 enum {x0_var, xhat_var, NVAR};
224 int true_s_test = test_space->
TrueVSize();
228 true_offsets[1] = true_s0;
229 true_offsets[2] = true_s0+true_s1;
232 true_offsets_test[0] = 0;
233 true_offsets_test[1] = true_s_test;
251 matSinv->
Mult(*trueF, SinvF);
270 if (dim == 2) { Shatinv =
new HypreAMS(*Shat, xhat_space); }
271 else { Shatinv =
new HypreADS(*Shat, xhat_space); }
292 matSinv->
Mult(LSres, tmp);
296 cout <<
"\n|| B0*x0 + Bhat*xhat - F ||_{S^-1} = " << res << endl;
305 ostringstream mesh_name, sol_name;
306 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
307 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
309 ofstream mesh_ofs(mesh_name.str().c_str());
310 mesh_ofs.precision(8);
311 pmesh->
Print(mesh_ofs);
313 ofstream sol_ofs(sol_name.str().c_str());
314 sol_ofs.precision(8);
321 char vishost[] =
"localhost";
324 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
325 sol_sock.precision(8);
326 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.
Subclass constant coefficient.
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
HYPRE_Int GlobalTrueVSize()
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)
void SetBlock(int iRow, int iCol, Operator *op)
Add a block op in the block-entry (iblock, jblock).
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
Integrator that inverts the matrix assembled by another integrator.
The BoomerAMG solver in hypre.
void SetPrintLevel(int print_lvl)
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)
HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P)
Returns the matrix P^t * A * P.
int MeshGenerator()
Truegrid or NetGen?
The operator x -> R*A*P*x.
void PrintUsage(std::ostream &out) const
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Wrapper for hypre's parallel vector class.
Array< int > bdr_attributes
int main(int argc, char *argv[])
void SetRelTol(double rtol)
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)
void Distribute(const Vector *tv)
double InnerProduct(HypreParVector *x, HypreParVector *y)
void PrintOptions(std::ostream &out) const
Abstract class for hypre's solvers and preconditioners.
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.
Arbitrary order H1-conforming (continuous) finite elements.
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs) const
Determine the boundary degrees of freedom.
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 SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Arbitrary order "L2-conforming" discontinuous finite elements.
virtual void Print(std::ostream &out=std::cout) const