18 #ifndef LIBMESH_EXACT_SOLUTION_H 19 #define LIBMESH_EXACT_SOLUTION_H 25 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS 45 class EquationSystems;
48 template <
typename Output>
class FunctionBase;
127 const std::string & sys_name,
128 const std::string & unknown_name);
150 const std::string & sys_name,
151 const std::string & unknown_name);
173 const std::string & sys_name,
174 const std::string & unknown_name);
192 const std::string & unknown_name);
202 const std::string & unknown_name);
212 const std::string & unknown_name);
227 const std::string & unknown_name);
237 const std::string & unknown_name);
250 const std::string & unknown_name);
263 const std::string & unknown_name);
273 const std::string & unknown_name);
286 const std::string & unknown_name,
296 template<
typename OutputShape>
298 const std::string & unknown_name,
299 std::vector<Real> & error_vals);
307 std::vector<Real> &
_check_inputs(
const std::string & sys_name,
308 const std::string & unknown_name);
344 std::map<std::string, SystemErrorMap>
_errors;
369 #endif // LIBMESH_EXACT_SOLUTION_H
Manages multiples systems of equations.
void _compute_error(const std::string &sys_name, const std::string &unknown_name, std::vector< Real > &error_vals)
const EquationSystems & _equation_systems
ExactSolution(const EquationSystems &es)
void attach_exact_values(const std::vector< FunctionBase< Number > *> &f)
Real h2_error(const std::string &sys_name, const std::string &unknown_name)
Real hdiv_error(const std::string &sys_name, const std::string &unknown_name)
void extra_quadrature_order(const int extraorder)
void attach_exact_derivs(const std::vector< FunctionBase< Gradient > *> &g)
Real error_norm(const std::string &sys_name, const std::string &unknown_name, const FEMNormType &norm)
void compute_error(const std::string &sys_name, const std::string &unknown_name)
const EquationSystems * _equation_systems_fine
Number(* ValueFunctionPointer)(const Point &p, const Parameters &Parameters, const std::string &sys_name, const std::string &unknown_name)
std::map< std::string, std::vector< Real > > SystemErrorMap
TensorValue< Number > NumberTensorValue
void attach_exact_hessians(std::vector< FunctionBase< Tensor > *> h)
ExactSolution & operator=(const ExactSolution &)=delete
NumberVectorValue Gradient
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
void attach_exact_hessian(unsigned int sys_num, FunctionBase< Tensor > *h)
Real l1_error(const std::string &sys_name, const std::string &unknown_name)
void attach_exact_deriv(unsigned int sys_num, FunctionBase< Gradient > *g)
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
Real hcurl_error(const std::string &sys_name, const std::string &unknown_name)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::map< std::string, SystemErrorMap > _errors
Real h1_error(const std::string &sys_name, const std::string &unknown_name)
Real l_inf_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)
void attach_exact_value(unsigned int sys_num, FunctionBase< Number > *f)
void attach_reference_solution(const EquationSystems *es_fine)
Real l2_error(const std::string &sys_name, const std::string &unknown_name)
VectorValue< Number > NumberVectorValue
Gradient(* GradientFunctionPointer)(const Point &p, const Parameters ¶meters, const std::string &sys_name, const std::string &unknown_name)
A geometric point in (x,y,z) space.
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values
std::vector< Real > & _check_inputs(const std::string &sys_name, const std::string &unknown_name)