40 _equation_systems(es),
41 _equation_systems_fine(nullptr),
52 const std::string & sys_name = system.
name();
57 for (
unsigned int var=0; var<system.
n_vars(); ++var)
61 sem[var_name] = std::vector<Real>(5, 0.);
75 libmesh_assert(es_fine);
128 libmesh_assert(gptr);
169 libmesh_assert(hptr);
209 const std::string & unknown_name)
216 std::map<std::string, SystemErrorMap>::iterator sys_iter =
220 libmesh_error_msg(
"Sorry, couldn't find the requested system '" << sys_name <<
"'.");
223 SystemErrorMap::iterator var_iter = (*sys_iter).second.find(unknown_name);
225 if (var_iter == (*sys_iter).second.end())
226 libmesh_error_msg(
"Sorry, couldn't find the requested variable '" << unknown_name <<
"'.");
229 return (*var_iter).second;
235 const std::string & unknown_name)
239 std::vector<Real> & error_vals = this->
_check_inputs(sys_name,
245 libmesh_assert( sys.has_variable( unknown_name ) );
250 this->_compute_error<Real>(sys_name,
257 this->_compute_error<RealGradient>(sys_name,
263 libmesh_error_msg(
"Invalid variable type!");
272 const std::string & unknown_name,
277 std::vector<Real> & error_vals = this->
_check_inputs(sys_name,
287 return std::sqrt(error_vals[0]);
289 return std::sqrt(error_vals[0] + error_vals[1]);
291 return std::sqrt(error_vals[0] + error_vals[1] + error_vals[2]);
295 libmesh_error_msg(
"Cannot compute HCurl error norm of scalar-valued variables!");
297 return std::sqrt(error_vals[0] + error_vals[5]);
302 libmesh_error_msg(
"Cannot compute HDiv error norm of scalar-valued variables!");
304 return std::sqrt(error_vals[0] + error_vals[6]);
307 return std::sqrt(error_vals[1]);
309 return std::sqrt(error_vals[2]);
313 libmesh_error_msg(
"Cannot compute HCurl error seminorm of scalar-valued variables!");
315 return std::sqrt(error_vals[5]);
320 libmesh_error_msg(
"Cannot compute HDiv error seminorm of scalar-valued variables!");
322 return std::sqrt(error_vals[6]);
325 return error_vals[3];
327 return error_vals[4];
330 libmesh_error_msg(
"Currently only Sobolev norms/seminorms are supported!");
341 const std::string & unknown_name)
346 std::vector<Real> & error_vals = this->
_check_inputs(sys_name,
351 return std::sqrt(error_vals[0]);
361 const std::string & unknown_name)
366 std::vector<Real> & error_vals = this->
_check_inputs(sys_name,
371 return error_vals[3];
381 const std::string & unknown_name)
386 std::vector<Real> & error_vals = this->
_check_inputs(sys_name,
391 return error_vals[4];
401 const std::string & unknown_name)
409 std::vector<Real> & error_vals = this->
_check_inputs(sys_name,
413 return std::sqrt(error_vals[0] + error_vals[1]);
418 const std::string & unknown_name)
425 const std::string & unknown_name)
433 const std::string & unknown_name)
441 std::vector<Real> & error_vals = this->
_check_inputs(sys_name,
445 return std::sqrt(error_vals[0] + error_vals[1] + error_vals[2]);
454 template<
typename OutputShape>
456 const std::string & unknown_name,
457 std::vector<Real> & error_vals)
477 const unsigned int sys_num = computed_system.
number();
478 const unsigned int var = computed_system.
variable_number(unknown_name);
479 const unsigned int var_component =
483 std::unique_ptr<MeshFunction> coarse_values;
487 const System & comparison_system
490 std::vector<Number> global_soln;
492 comparison_soln->init(comparison_system.
solution->size(),
true,
SERIAL);
493 (*comparison_soln) = global_soln;
495 coarse_values = std::unique_ptr<MeshFunction>
500 coarse_values->init();
532 error_vals = std::vector<Real>(7, 0.);
543 libMesh::err <<
"Error calculation using reference solution not yet\n" 544 <<
"supported for vector-valued elements." 546 libmesh_not_implemented();
551 std::vector<std::unique_ptr<FEGenericBase<OutputShape>>> fe_ptrs(4);
552 std::vector<std::unique_ptr<QBase>> q_rules(4);
555 for (
const auto dim : elem_dims)
564 fe_ptrs[dim]->attach_quadrature_rule (q_rules[dim].
get());
569 std::vector<dof_id_type> dof_indices;
582 const unsigned int dim = elem->dim();
594 std::set<subdomain_id_type> subdomain_id;
595 subdomain_id.insert(elem_subid);
598 QBase * qrule = q_rules[dim].get();
600 libmesh_assert(qrule);
603 const std::vector<Real> & JxW = fe->
get_JxW();
607 const std::vector<std::vector<OutputShape>> & phi_values = fe->
get_phi();
610 const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient>> &
615 const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputShape>> * curl_values =
nullptr;
619 const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence>> * div_values =
nullptr;
627 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 629 const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor>> &
634 const std::vector<Point> & q_point = fe->
get_xyz();
641 computed_dof_map.
dof_indices (elem, dof_indices, var);
644 const unsigned int n_qp = qrule->
n_points();
647 const unsigned int n_sf =
648 cast_int<unsigned int>(dof_indices.size());
653 for (
unsigned int qp=0; qp<n_qp; qp++)
661 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 671 for (
unsigned int i=0; i<n_sf; i++)
675 grad_u_h += dphi_values[i][qp]*computed_system.
current_solution (dof_indices[i]);
676 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 677 grad2_u_h += d2phi_values[i][qp]*computed_system.
current_solution (dof_indices[i]);
681 curl_u_h += (*curl_values)[i][qp]*computed_system.
current_solution (dof_indices[i]);
682 div_u_h += (*div_values)[i][qp]*computed_system.
current_solution (dof_indices[i]);
691 for (
unsigned int c = 0; c < n_vec_dim; c++)
692 exact_val_accessor(c) =
694 component(var_component+c, q_point[qp], time);
700 (*coarse_values)(q_point[qp],time,output,&subdomain_id);
701 exact_val = output(0);
707 error_vals[0] += JxW[qp]*error_sq;
709 Real norm = sqrt(error_sq);
710 error_vals[3] += JxW[qp]*norm;
712 if (error_vals[4]<norm) { error_vals[4] = norm; }
720 for (
unsigned int c = 0; c < n_vec_dim; c++)
721 for (
unsigned int d = 0; d < LIBMESH_DIM; d++)
722 exact_grad_accessor(d + c*LIBMESH_DIM) =
724 component(var_component+c, q_point[qp], time)(d);
729 std::vector<Gradient> output(1);
730 coarse_values->gradient(q_point[qp],time,output,&subdomain_id);
731 exact_grad = output[0];
736 error_vals[1] += JxW[qp]*grad_error.
norm_sq();
776 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 787 libmesh_not_implemented();
789 for (
unsigned int c = 0; c < n_vec_dim; c++)
790 for (
unsigned int d = 0; d < dim; d++)
791 for (
unsigned int e =0; e < dim; e++)
792 exact_hess_accessor(d + e*dim + c*dim*dim) =
794 component(var_component+c, q_point[qp], time)(d,e);
799 std::vector<Tensor> output(1);
800 coarse_values->hessian(q_point[qp],time,output,&subdomain_id);
801 exact_hess = output[0];
807 error_vals[2] += JxW[qp]*grad2_error.
norm_sq();
815 Real l_infty_norm = error_vals[4];
818 error_vals[4] = l_infty_norm;
822 template void ExactSolution::_compute_error<Real>(
const std::string &,
const std::string &, std::vector<Real> &);
823 template void ExactSolution::_compute_error<RealGradient>(
const std::string &,
const std::string &, std::vector<Real> &);
Manages the family, order, etc. parameters for a given FE.
Wrap a libMesh-style function pointer into a FunctionBase object.
Manages multiples systems of equations.
void _compute_error(const std::string &sys_name, const std::string &unknown_name, std::vector< Real > &error_vals)
unsigned int variable_scalar_number(const std::string &var, unsigned int component) const
void update_global_solution(std::vector< Number > &global_soln) const
unsigned int n_systems() const
const Variable & variable(unsigned int var) const
const std::vector< Point > & get_xyz() const
TensorTools::DecrementRank< OutputNumber >::type OutputNumberDivergence
const EquationSystems & _equation_systems
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
ExactSolution(const EquationSystems &es)
void attach_exact_values(const std::vector< FunctionBase< Number > *> &f)
const T_sys & get_system(const std::string &name) const
const FEType & variable_type(const unsigned int c) const
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)
static FEFieldType field_type(const FEType &fe_type)
const Parallel::Communicator & comm() const
void attach_exact_derivs(const std::vector< FunctionBase< Gradient > *> &g)
const std::vector< std::vector< OutputDivergence > > & get_div_phi() const
Real error_norm(const std::string &sys_name, const std::string &unknown_name, const FEMNormType &norm)
const MeshBase & get_mesh() const
void compute_error(const std::string &sys_name, const std::string &unknown_name)
const EquationSystems * _equation_systems_fine
Number current_solution(const dof_id_type global_dof_number) const
std::map< std::string, std::vector< Real > > SystemErrorMap
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
Manages the degrees of freedom (DOFs) in a simulation.
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
const std::vector< std::vector< OutputGradient > > & get_dphi() const
unsigned int number() const
void attach_exact_hessians(std::vector< FunctionBase< Tensor > *> h)
Manages consistently variables, degrees of freedom, and coefficient vectors.
unsigned short int variable_number(const std::string &var) const
const std::vector< Real > & get_JxW() const
std::unique_ptr< NumericVector< Number > > solution
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr)=0
const std::string & variable_name(const unsigned int i) const
unsigned int n_points() const
const std::set< unsigned char > & elem_dimensions() const
OStreamProxy err(std::cerr)
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
bool active_on_subdomain(subdomain_id_type sid) const
TensorTools::MakeNumber< OutputShape >::type OutputNumber
void attach_exact_hessian(unsigned int sys_num, FunctionBase< Tensor > *h)
Real l1_error(const std::string &sys_name, const std::string &unknown_name)
static unsigned int n_vec_dim(const MeshBase &mesh, const FEType &fe_type)
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
void attach_exact_deriv(unsigned int sys_num, FunctionBase< Gradient > *g)
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
const std::vector< std::vector< OutputShape > > & get_curl_phi() const
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)
bool has_system(const std::string &name) const
virtual std::unique_ptr< FunctionBase< Output > > clone() const =0
void attach_exact_value(unsigned int sys_num, FunctionBase< Number > *f)
const std::string & name() const
void attach_reference_solution(const EquationSystems *es_fine)
Real l2_error(const std::string &sys_name, const std::string &unknown_name)
unsigned int n_vars() const
const DofMap & get_dof_map() const
Base class for all quadrature families and orders.
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values
std::vector< Real > & _check_inputs(const std::string &sys_name, const std::string &unknown_name)
const std::vector< std::vector< OutputShape > > & get_phi() const