68 int main(
int argc,
char *argv[])
71 const char *mesh_file =
"../../data/star.mesh";
77 bool use_pointwise_transfer =
false;
80 args.
AddOption(&mesh_file,
"-m",
"--mesh",
83 "Problem type (see the RHO_exact function).");
85 "Finite element order (polynomial degree) or -1 for" 86 " isoparametric space.");
87 args.
AddOption(&lref,
"-lref",
"--lor-ref-level",
"LOR refinement level.");
88 args.
AddOption(&lorder,
"-lo",
"--lor-order",
89 "LOR space order (polynomial degree, zero by default).");
90 args.
AddOption(&vis,
"-vis",
"--visualization",
"-no-vis",
92 "Enable or disable GLVis visualization.");
93 args.
AddOption(&useH1,
"-h1",
"--use-h1",
"-l2",
"--use-l2",
94 "Use H1 spaces instead of L2.");
95 args.
AddOption(&use_pointwise_transfer,
"-t",
"--use-pointwise-transfer",
96 "-no-t",
"--dont-use-pointwise-transfer",
97 "Use pointwise transfer operators instead of L2 projection.");
107 Mesh mesh(mesh_file, 1, 1);
111 int basis_lor = BasisType::GaussLobatto;
112 Mesh mesh_lor = Mesh::MakeRefined(mesh, lref, basis_lor);
122 cerr <<
"Switching the H1 LOR space order from 0 to 1\n";
160 double ho_mass =
compute_mass(&fespace, -1.0, HO_dc,
"HO ");
164 if (use_pointwise_transfer)
176 R.
Mult(rho, rho_lor);
186 P.
Mult(rho_lor, rho);
192 cout <<
"|HO - P(R(HO))|_∞ = " << rho_prev.
Normlinf() << endl;
199 LinearForm M_rho(&fespace), M_rho_lor(&fespace_lor);
203 M_ho.
Mult(rho, M_rho);
205 cout <<
"HO -> LOR dual field: " << fabs(M_rho(ones)-M_rho_lor(ones_lor))
213 double lor_mass =
compute_mass(&fespace_lor, -1.0, LOR_dc,
"LOR ");
221 P.
Mult(rho_lor, rho);
228 R.
Mult(rho, rho_lor);
229 compute_mass(&fespace_lor, lor_mass, LOR_dc,
"R(P(LOR))");
232 rho_lor_prev -= rho_lor;
234 cout <<
"|LOR - R(P(LOR))|_∞ = " << rho_lor_prev.
Normlinf() << endl;
238 if (!use_pointwise_transfer)
240 M_lor.
Mult(rho_lor, M_rho_lor);
242 cout <<
"LOR -> HO dual field: " << fabs(M_rho(ones)-M_rho_lor(ones_lor))
258 return x(1)+0.25*cos(2*M_PI*x.
Norml2());
260 return x(1)*x(1)*x(1) + 2*x(0)*x(1) + x(0);
262 return M_PI/2-atan(5*(2*x.
Norml2()-1));
264 return (x.
Norml2() < 0.1) ? 1 : 0;
279 sol_sockL2.precision(8);
281 <<
"window_geometry " << x <<
" " << y <<
" " << w <<
" " << h
282 <<
"plot_caption '" <<
space <<
" " << prefix <<
" Density'" 283 <<
"window_title '" <<
direction <<
"'" << flush;
300 cout <<
space <<
" " << prefix <<
" mass = " << newmass;
304 cout <<
" (" << fabs(newmass-massL2)*100/massL2 <<
"%)";
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 ...
Class for grid function - Vector with associated FE space.
int main(int argc, char *argv[])
A coefficient that is constant across space and time.
void PrintOptions(std::ostream &out) const
Print the options.
void PrintUsage(std::ostream &out) const
Print the usage message.
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
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).
bool Good() const
Return true if the command line options were parsed successfully.
double RHO_exact(const Vector &x)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
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.
virtual bool SupportsBackwardsOperator() 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.
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...
Transfer data in L2 and H1 finite element spaces between a coarse mesh and an embedded refined mesh u...
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
double Norml2() const
Returns the l2 norm of the vector.
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
A general function coefficient.
Arbitrary order H1-conforming (continuous) finite elements.
double compute_mass(FiniteElementSpace *, double, VisItDataCollection &, string)
double Normlinf() const
Returns the l_infinity norm of the vector.
void visualize(VisItDataCollection &, string, int, int)
Arbitrary order "L2-conforming" discontinuous finite elements.
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.