18 #ifndef LIBMESH_EXACT_SOLUTION_H    19 #define LIBMESH_EXACT_SOLUTION_H    23 #include "libmesh/libmesh_common.h"     37 class EquationSystems;
    40 template <
typename Output> 
class FunctionBase;
    41 template <
typename Output> 
class FEMFunctionBase;
   142                                          const std::string & sys_name,
   143                                          const std::string & unknown_name);
   178                                               const std::string & sys_name,
   179                                               const std::string & unknown_name);
   213                                            const std::string & sys_name,
   214                                            const std::string & unknown_name);
   235                      std::string_view unknown_name);
   245                 std::string_view unknown_name);
   255                 std::string_view unknown_name);
   270                    std::string_view unknown_name);
   280                 std::string_view unknown_name);
   293                    std::string_view unknown_name);
   306                   std::string_view unknown_name);
   316                 std::string_view unknown_name);
   329                   std::string_view unknown_name,
   339   template<
typename OutputShape>
   341                       std::string_view unknown_name,
   342                       std::vector<Real> & error_vals);
   351                                     std::string_view unknown_name);
   379   typedef std::map<std::string, std::vector<Real>, std::less<>>
   388   std::map<std::string, SystemErrorMap, std::less<>> 
_errors;
   419 #endif // LIBMESH_EXACT_SOLUTION_H std::map< std::string, std::vector< Real >, std::less<> > SystemErrorMap
Data structure which stores the errors: ||e|| = ||u - u_h|| ||grad(e)|| = ||grad(u) - grad(u_h)|| for...
This class handles the computation of the L2 and/or H1 error for the Systems in the EquationSystems o...
Real l_inf_error(std::string_view sys_name, std::string_view unknown_name)
This is the EquationSystems class. 
Real h2_error(std::string_view sys_name, std::string_view unknown_name)
Real l1_error(std::string_view sys_name, std::string_view unknown_name)
Real hdiv_error(std::string_view sys_name, std::string_view unknown_name)
This class provides the ability to map between arbitrary, user-defined strings and several data types...
const EquationSystems & _equation_systems
Constant reference to the EquationSystems object used for the simulation. 
std::map< std::string, SystemErrorMap, std::less<> > _errors
A map of SystemErrorMaps, which contains entries for each system in the EquationSystems object...
ExactSolution(const EquationSystems &es)
Constructor. 
void attach_exact_values(const std::vector< FunctionBase< Number > *> &f)
Clone and attach arbitrary functors which compute the exact values of the EquationSystems' solutions ...
int _extra_order
Extra order to use for quadrature rule. 
FEMNormType
defines an enum for norms defined on vectors of finite element coefficients 
Real hcurl_error(std::string_view sys_name, std::string_view unknown_name)
void extra_quadrature_order(const int extraorder)
Increases or decreases the order of the quadrature rule used for numerical integration. 
std::vector< Real > & _check_inputs(std::string_view sys_name, std::string_view unknown_name)
This function is responsible for checking the validity of the sys_name and unknown_name inputs...
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space. 
std::vector< std::unique_ptr< FEMFunctionBase< Number > > > _exact_values
User-provided functors which compute the exact value of the solution for each system. 
The libMesh namespace provides an interface to certain functionality in the library. 
void attach_exact_derivs(const std::vector< FunctionBase< Gradient > *> &g)
Clone and attach arbitrary functors which compute the exact gradients of the EquationSystems' solutio...
const EquationSystems * _equation_systems_fine
Constant pointer to the EquationSystems object containing the fine grid solution. ...
std::vector< std::unique_ptr< FEMFunctionBase< Gradient > > > _exact_derivs
User-provided functors which compute the exact derivative of the solution for each system...
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
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...
std::set< subdomain_id_type > _excluded_subdomains
Elements in a subdomain from this set are skipped during the error computation. 
TensorValue< Number > NumberTensorValue
void attach_exact_hessians(std::vector< FunctionBase< Tensor > *> h)
Clone and attach arbitrary functors which compute the exact second derivatives of the EquationSystems...
ExactSolution & operator=(const ExactSolution &)=delete
void compute_error(std::string_view sys_name, std::string_view unknown_name)
Computes and stores the error in the solution value e = u-u_h, the gradient grad(e) = grad(u) - grad(...
NumberVectorValue Gradient
Real l2_error(std::string_view sys_name, std::string_view unknown_name)
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...
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...
void _compute_error(std::string_view sys_name, std::string_view 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...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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...
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...
Gradient gptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
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...
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...
VectorValue< Number > NumberVectorValue
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...
Real h1_error(std::string_view sys_name, std::string_view unknown_name)
A Point defines a location in LIBMESH_DIM dimensional Real space. 
std::vector< std::unique_ptr< FEMFunctionBase< Tensor > > > _exact_hessians
User-provided functors which compute the exact hessians of the solution for each system. 
void ErrorVector unsigned int
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space. 
Real error_norm(std::string_view sys_name, std::string_view unknown_name, const FEMNormType &norm)