68int main(
int argc,
char *argv[])
71 const char *mesh_file =
"../../data/star.mesh";
78 bool use_pointwise_transfer =
false;
79 const char *device_config =
"cpu";
83 args.
AddOption(&mesh_file,
"-m",
"--mesh",
86 "Problem type (see the RHO_exact function).");
88 "Finite element order (polynomial degree) or -1 for"
89 " isoparametric space.");
90 args.
AddOption(&lref,
"-lref",
"--lor-ref-level",
"LOR refinement level.");
91 args.
AddOption(&lorder,
"-lo",
"--lor-order",
92 "LOR space order (polynomial degree, zero by default).");
93 args.
AddOption(&vis,
"-vis",
"--visualization",
"-no-vis",
95 "Enable or disable GLVis visualization.");
96 args.
AddOption(&useH1,
"-h1",
"--use-h1",
"-l2",
"--use-l2",
97 "Use H1 spaces instead of L2.");
98 args.
AddOption(&use_pointwise_transfer,
"-t",
"--use-pointwise-transfer",
99 "-no-t",
"--dont-use-pointwise-transfer",
100 "Use pointwise transfer operators instead of L2 projection.");
101 args.
AddOption(&device_config,
"-d",
"--device",
102 "Device configuration string, see Device::Configure().");
103 args.
AddOption(&use_ea,
"-ea",
"--ea-version",
"-no-ea",
104 "--no-ea-version",
"Use element assembly version.");
108 Device device(device_config);
111 Mesh mesh(mesh_file, 1, 1);
126 cerr <<
"Switching the H1 LOR space order from 0 to 1\n";
172 if (use_pointwise_transfer)
188 R.
Mult(rho, rho_lor);
198 P.
Mult(rho_lor, rho);
204 cout <<
"|HO - P(R(HO))|_∞ = " << rho_prev.
Normlinf() << endl;
208 LinearForm M_rho(&fespace), M_rho_lor(&fespace_lor);
212 M_ho.
Mult(rho, M_rho);
214 cout <<
"HO -> LOR dual field: " << abs(M_rho.Sum()-M_rho_lor.
Sum()) <<
"\n\n";
229 P.
Mult(rho_lor, rho);
236 R.
Mult(rho, rho_lor);
237 compute_mass(&fespace_lor, lor_mass, LOR_dc,
"R(P(LOR))");
238 if (vis) {
visualize(LOR_dc,
"R(P(LOR))",
Wx,
Wy, visport); }
240 rho_lor_prev -= rho_lor;
242 cout <<
"|LOR - R(P(LOR))|_∞ = " << rho_lor_prev.
Normlinf() << endl;
246 if (!use_pointwise_transfer)
248 M_lor.
Mult(rho_lor, M_rho_lor);
250 cout <<
"LOR -> HO dual field: " << abs(M_rho.Sum() - M_rho_lor.
Sum()) <<
'\n';
266 return x(1)+0.25*cos(2*M_PI*x.
Norml2());
268 return x(1)*x(1)*x(1) + 2*x(0)*x(1) + x(0);
270 return M_PI/2-atan(5*(2*x.
Norml2()-1));
272 return (x.
Norml2() < 0.1) ? 1 : 0;
287 sol_sockL2.precision(8);
289 <<
"window_geometry " << x <<
" " << y <<
" " << w <<
" " << h
290 <<
"plot_caption '" <<
space <<
" " << prefix <<
" Density'"
291 <<
"window_title '" <<
direction <<
"'" << flush;
305 cout <<
space <<
" " << prefix <<
" mass = " << newmass;
309 cout <<
" (" << fabs(newmass-massL2)*100/massL2 <<
"%)";
@ GaussLobatto
Closed type.
A coefficient that is constant across space and time.
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Class for domain integration .
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
A general function coefficient.
Class for grid function - Vector with associated FE space.
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
virtual bool SupportsBackwardsOperator() const
virtual const Operator & ForwardOperator()=0
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
virtual const Operator & BackwardOperator()=0
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Arbitrary order H1-conforming (continuous) finite elements.
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Transfer data in L2 and H1 finite element spaces between a coarse mesh and an embedded refined mesh u...
Arbitrary order "L2-conforming" discontinuous finite elements.
int Dimension() const
Dimension of the reference space used within the elements.
static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type)
Create a refined (by any factor) version of orig_mesh.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void ParseCheck(std::ostream &out=mfem::out)
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...
real_t Normlinf() const
Returns the l_infinity norm of the vector.
real_t Norml2() const
Returns the l2 norm of the vector.
real_t Sum() const
Return the sum of the vector entries.
Data collection with VisIt I/O routines.
void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
void visualize(VisItDataCollection &, string, int, int, int visport=19916)
real_t RHO_exact(const Vector &x)
real_t compute_mass(FiniteElementSpace *, real_t, VisItDataCollection &, string)