Go to the documentation of this file.
18 #ifndef LIBMESH_EXACT_SOLUTION_H
19 #define LIBMESH_EXACT_SOLUTION_H
23 #include "libmesh/libmesh_common.h"
25 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS
31 #include "libmesh/enum_norm_type.h"
46 class EquationSystems;
49 template <
typename Output>
class FunctionBase;
136 const std::string & sys_name,
137 const std::string & unknown_name);
159 const std::string & sys_name,
160 const std::string & unknown_name);
182 const std::string & sys_name,
183 const std::string & unknown_name);
201 const std::string & unknown_name);
211 const std::string & unknown_name);
221 const std::string & unknown_name);
236 const std::string & unknown_name);
246 const std::string & unknown_name);
259 const std::string & unknown_name);
272 const std::string & unknown_name);
282 const std::string & unknown_name);
295 const std::string & unknown_name,
305 template<
typename OutputShape>
307 const std::string & unknown_name,
308 std::vector<Real> & error_vals);
316 std::vector<Real> &
_check_inputs(
const std::string & sys_name,
317 const std::string & unknown_name);
353 std::map<std::string, SystemErrorMap>
_errors;
384 #endif // LIBMESH_EXACT_SOLUTION_H
int _extra_order
Extra order to use for quadrature rule.
This class handles the computation of the L2 and/or H1 error for the Systems in the EquationSystems o...
ExactSolution(const EquationSystems &es)
Constructor.
void attach_exact_derivs(const std::vector< FunctionBase< Gradient > * > &g)
Clone and attach arbitrary functors which compute the exact gradients of the EquationSystems' solutio...
std::set< subdomain_id_type > _excluded_subdomains
Elements in a subdomain from this set are skipped during the error computation.
void set_excluded_subdomains(const std::set< subdomain_id_type > &excluded)
The user can indicate that elements in certain subdomains should be excluded from the error calculati...
Real l1_error(const std::string &sys_name, const std::string &unknown_name)
Gradient(* GradientFunctionPointer)(const Point &p, const Parameters ¶meters, const std::string &sys_name, const std::string &unknown_name)
Attach an arbitrary function which computes the exact gradient of the solution at any point.
The libMesh namespace provides an interface to certain functionality in the library.
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
User-provided functors which compute the exact hessians of the solution for each system.
Real h1_error(const std::string &sys_name, const std::string &unknown_name)
Tensor(* HessianFunctionPointer)(const Point &p, const Parameters ¶meters, const std::string &sys_name, const std::string &unknown_name)
Attach an arbitrary function which computes the exact second derivatives of the solution at any point...
Real h2_error(const std::string &sys_name, const std::string &unknown_name)
Real hcurl_error(const std::string &sys_name, const std::string &unknown_name)
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
std::map< std::string, SystemErrorMap > _errors
A map of SystemErrorMaps, which contains entries for each system in the EquationSystems object.
void attach_exact_value(unsigned int sys_num, FunctionBase< Number > *f)
Clone and attach an arbitrary functor which computes the exact value of the system sys_num solution a...
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values
User-provided functors which compute the exact value of the solution for each system.
void extra_quadrature_order(const int extraorder)
Increases or decreases the order of the quadrature rule used for numerical integration.
void attach_reference_solution(const EquationSystems *es_fine)
Attach function similar to system.h which allows the user to attach a second EquationSystems object w...
void attach_exact_hessian(unsigned int sys_num, FunctionBase< Tensor > *h)
Clone and attach an arbitrary functor which computes the exact second derivatives of the system sys_n...
Real hdiv_error(const std::string &sys_name, const std::string &unknown_name)
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
A Point defines a location in LIBMESH_DIM dimensional Real space.
Real error_norm(const std::string &sys_name, const std::string &unknown_name, const FEMNormType &norm)
TensorValue< Number > NumberTensorValue
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
User-provided functors which compute the exact derivative of the solution for each system.
void attach_exact_hessians(std::vector< FunctionBase< Tensor > * > h)
Clone and attach arbitrary functors which compute the exact second derivatives of the EquationSystems...
const EquationSystems & _equation_systems
Constant reference to the EquationSystems object used for the simulation.
Gradient gptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Real l2_error(const std::string &sys_name, const std::string &unknown_name)
void compute_error(const std::string &sys_name, const std::string &unknown_name)
Computes and stores the error in the solution value e = u-u_h, the gradient grad(e) = grad(u) - grad(...
This is the EquationSystems class.
void attach_exact_values(const std::vector< FunctionBase< Number > * > &f)
Clone and attach arbitrary functors which compute the exact values of the EquationSystems' solutions ...
ExactSolution & operator=(const ExactSolution &)=delete
NumberVectorValue Gradient
Number(* ValueFunctionPointer)(const Point &p, const Parameters &Parameters, const std::string &sys_name, const std::string &unknown_name)
Attach an arbitrary function which computes the exact value of the solution at any point.
VectorValue< Number > NumberVectorValue
void _compute_error(const std::string &sys_name, const std::string &unknown_name, std::vector< Real > &error_vals)
This function computes the error (in the solution and its first derivative) for a single unknown in a...
MetaPhysicL::DualNumber< T, D > norm(const MetaPhysicL::DualNumber< T, D > &in)
std::vector< Real > & _check_inputs(const std::string &sys_name, const std::string &unknown_name)
This function is responsible for checking the validity of the sys_name and unknown_name inputs.
Real l_inf_error(const std::string &sys_name, const std::string &unknown_name)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
std::map< std::string, std::vector< Real > > SystemErrorMap
Data structure which stores the errors: ||e|| = ||u - u_h|| ||grad(e)|| = ||grad(u) - grad(u_h)|| for...
void attach_exact_deriv(unsigned int sys_num, FunctionBase< Gradient > *g)
Clone and attach an arbitrary functor which computes the exact gradient of the system sys_num solutio...
const EquationSystems * _equation_systems_fine
Constant pointer to the EquationSystems object containing the fine grid solution.
void ErrorVector unsigned int
This class provides the ability to map between arbitrary, user-defined strings and several data types...