65 int main(
int argc,
char *argv[])
68 const char *mesh_file =
"../../data/star.mesh";
74 bool use_transfer =
false;
77 args.
AddOption(&mesh_file,
"-m",
"--mesh",
80 "Problem type (see the RHO_exact function).");
82 "Finite element order (polynomial degree) or -1 for"
83 " isoparametric space.");
84 args.
AddOption(&lref,
"-lref",
"--lor-ref-level",
"LOR refinement level.");
85 args.
AddOption(&lorder,
"-lo",
"--lor-order",
86 "LOR space order (polynomial degree, zero by default).");
87 args.
AddOption(&vis,
"-vis",
"--visualization",
"-no-vis",
89 "Enable or disable GLVis visualization.");
90 args.
AddOption(&useH1,
"-h1",
"--use-h1",
"-l2",
"--use-l2",
91 "Use H1 spaces instead of L2.");
92 args.
AddOption(&use_transfer,
"-t",
"--use-pointwise-transfer",
"-no-t",
93 "--dont-use-pointwise-transfer",
94 "Use pointwise transfer operators instead of L2 projection.");
104 Mesh mesh(mesh_file, 1, 1);
108 int basis_lor = BasisType::GaussLobatto;
109 Mesh mesh_lor(&mesh, lref, basis_lor);
119 cerr <<
"Switching the H1 LOR space order from 0 to 1\n";
148 double ho_mass =
compute_mass(&fespace, -1.0, HO_dc,
"HO ");
165 R.
Mult(rho, rho_lor);
172 P.
Mult(rho_lor, rho);
178 cout <<
"|HO - P(R(HO))|_∞ = " << rho_prev.
Normlinf() << endl << endl;
184 double lor_mass =
compute_mass(&fespace_lor, -1.0, LOR_dc,
"LOR ");
189 P.
Mult(rho_lor, rho);
196 R.
Mult(rho, rho_lor);
197 compute_mass(&fespace_lor, lor_mass, LOR_dc,
"R(P(LOR))");
200 rho_lor_prev -= rho_lor;
202 cout <<
"|LOR - R(P(LOR))|_∞ = " << rho_lor_prev.
Normlinf() << endl;
216 return x(1)+0.25*cos(2*M_PI*x.
Norml2());
218 return x(1)*x(1)*x(1) + 2*x(0)*x(1) + x(0);
220 return M_PI/2-atan(5*(2*x.
Norml2()-1));
222 return (x.
Norml2() < 0.1) ? 1 : 0;
233 char vishost[] =
"localhost";
237 sol_sockL2.precision(8);
239 <<
"window_geometry " << x <<
" " << y <<
" " << w <<
" " << h
240 <<
"plot_caption '" <<
space <<
" " << prefix <<
" Density'"
241 <<
"window_title '" <<
direction <<
"'" << flush;
258 cout <<
space <<
" " << prefix <<
" mass = " << newmass;
262 cout <<
" (" << fabs(newmass-massL2)*100/massL2 <<
"%)";
void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Class for grid function - Vector with associated FE space.
Subclass constant coefficient.
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
double Norml2() const
Returns the l2 norm of the vector.
virtual const Operator & BackwardOperator()=0
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
double Normlinf() const
Returns the l_infinity norm of the vector.
double RHO_exact(const Vector &x)
int main(int argc, char *argv[])
virtual const Operator & ForwardOperator()=0
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Data collection with VisIt I/O routines.
void PrintUsage(std::ostream &out) const
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
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)
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
void PrintOptions(std::ostream &out) const
Transfer data between a coarse mesh and an embedded refined mesh using L2 projection.
virtual void ProjectCoefficient(Coefficient &coeff)
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
class for C-function coefficient
Arbitrary order H1-conforming (continuous) finite elements.
double compute_mass(FiniteElementSpace *, double, VisItDataCollection &, string)
Arbitrary order "L2-conforming" discontinuous finite elements.