34 int main(
int argc,
char *argv[])
37 const char *mesh_file =
"../data/star.mesh";
39 bool visualization = 1;
42 args.
AddOption(&mesh_file,
"-m",
"--mesh",
45 "Finite element order (polynomial degree).");
46 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
48 "Enable or disable GLVis visualization.");
61 ifstream imesh(mesh_file);
64 cerr <<
"\nCan not open mesh file: " << mesh_file <<
'\n' << endl;
67 mesh =
new Mesh(imesh, 1, 1);
77 (int)floor(log(10000./mesh->
GetNE())/log(2.)/dim);
78 for (
int l = 0; l < ref_levels; l++)
89 int trial_order = order;
90 int trace_order = order - 1;
91 int test_order = order;
92 if (test_order < trial_order)
93 cerr <<
"Warning, test space not enriched enough to handle primal"
109 enum {x0_var, xhat_var, NVAR};
113 int s_test = test_space->
GetVSize();
122 offsets_test[1] = s_test;
124 std::cout <<
"\nNumber of Unknowns:\n"
125 <<
" Trial space, X0 : " << s0 <<
'\n'
126 <<
" Interface space, Xhat : " << s1 <<
'\n'
127 <<
" Test space, Y : " << s_test <<
"\n\n";
187 matSinv.
Mult(F,SinvF);
199 #ifndef MFEM_USE_SUITESPARSE
200 const double prec_rtol = 1e-3;
201 const int prec_maxit = 200;
228 PCG(A, P, b, x, 1, 200, 1e-12, 0.0);
235 cout <<
"\n|| B0*x0 + Bhat*xhat - F ||_{S^-1} = " << res << endl;
239 x0.
Update(x0_space, x.GetBlock(x0_var), 0);
244 ofstream mesh_ofs(
"refined.mesh");
245 mesh_ofs.precision(8);
246 mesh->
Print(mesh_ofs);
247 ofstream sol_ofs(
"sol.gf");
248 sol_ofs.precision(8);
255 char vishost[] =
"localhost";
258 sol_sock.precision(8);
259 sol_sock <<
"solution\n" << *mesh << x0 << flush;
Class for domain integration L(v) := (f, v)
Conjugate gradient method.
Class for grid function - Vector with associated FE space.
Integrator defining a sum of multiple Integrators.
Subclass constant coefficient.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator.
int GetNE() const
Returns number of elements.
void AddIntegrator(BilinearFormIntegrator *integ)
void SetBlock(int iRow, int iCol, Operator *op)
Add a block op in the block-entry (iblock, jblock).
bool iterative_mode
If true, use the second argument of Mult as an initial guess.
Direct sparse solver using UMFPACK.
Integrator that inverts the matrix assembled by another integrator.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
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 SetMaxIter(int max_it)
HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P)
Returns the matrix P^t * A * P.
The operator x -> R*A*P*x.
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
void PrintUsage(std::ostream &out) const
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
Array< int > bdr_attributes
int main(int argc, char *argv[])
void SetRelTol(double rtol)
Abstract finite element space.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
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 PrintOptions(std::ostream &out) const
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Arbitrary order H1-conforming (continuous) finite elements.
A class to handle Block systems in a matrix-free implementation.
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Arbitrary order "L2-conforming" discontinuous finite elements.