45 #ifdef LIBMESH_ENABLE_AMR 55 number_h_refinements(1),
56 number_p_refinements(0)
72 bool estimate_parent_error)
74 LOG_SCOPE(
"estimate_error()",
"UniformRefinementEstimator");
75 std::map<const System *, const NumericVector<Number> *> solution_vectors;
76 solution_vectors[&_system] = solution_vector;
83 estimate_parent_error);
88 const std::map<const System *, SystemNorm> & error_norms,
90 bool estimate_parent_error)
92 LOG_SCOPE(
"estimate_errors()",
"UniformRefinementEstimator");
99 estimate_parent_error);
105 bool estimate_parent_error)
107 LOG_SCOPE(
"estimate_errors()",
"UniformRefinementEstimator");
114 estimate_parent_error);
121 const std::map<const System *, SystemNorm> * _error_norms,
127 std::vector<System *> system_list;
128 std::unique_ptr<std::map<const System *, SystemNorm>> error_norms =
129 std::unique_ptr<std::map<const System *, SystemNorm>>
130 (
new std::map<const System *, SystemNorm>);
134 libmesh_assert(!_system);
137 libmesh_assert_equal_to (&(_system->get_equation_systems()), _es);
140 for (
unsigned int i=0; i != _es->
n_systems(); ++i)
142 system_list.push_back(const_cast<System *>(&(_es->
get_system(i))));
148 libmesh_assert(!errors_per_cell);
154 libmesh_assert (errors_per_cell);
156 _error_norms = error_norms.get();
158 for (
unsigned int i=0; i!= _es->
n_systems(); ++i)
163 std::vector<Real> weights(
n_vars, 0.0);
164 for (
unsigned int v = 0; v !=
n_vars; ++v)
166 if (errors_per_cell->find(std::make_pair(&sys, v)) ==
167 errors_per_cell->end())
172 (*error_norms)[&sys] =
180 libmesh_assert(_system);
182 system_list.push_back(const_cast<System *>(_system));
184 libmesh_assert(!_error_norms);
186 _error_norms = error_norms.get();
198 const unsigned int dim =
mesh.mesh_dimension();
204 error_per_cell->clear();
205 error_per_cell->resize (
mesh.max_elem_id(), 0.);
209 libmesh_assert(errors_per_cell);
210 for (
const auto & pr : *errors_per_cell)
214 e->resize(
mesh.max_elem_id(), 0.);
219 std::vector<std::map<std::string, std::unique_ptr<NumericVector<Number>>>> coarse_vectors(system_list.size());
220 std::vector<std::unique_ptr<NumericVector<Number>>> coarse_solutions(system_list.size());
221 std::vector<std::unique_ptr<NumericVector<Number>>> coarse_local_solutions(system_list.size());
223 std::vector<std::unique_ptr<NumericVector<Number>>> projected_solutions(system_list.size());
226 std::vector<bool> old_projection_settings(system_list.size());
229 std::unique_ptr<Partitioner> old_partitioner(
mesh.partitioner().release());
232 const bool old_renumbering_setting =
mesh.allow_renumbering();
233 mesh.allow_renumbering(
false);
237 System & system = *system_list[i];
240 libmesh_assert (_error_norms->find(&system) !=
241 _error_norms->end());
244 coarse_solutions[i] = system.
solution->clone();
252 const std::string & var_name = vec->first;
254 coarse_vectors[i][var_name] = vec->second->clone();
258 if (solution_vectors &&
259 solution_vectors->find(&system) != solution_vectors->end() &&
260 solution_vectors->find(&system)->second &&
261 solution_vectors->find(&system)->second != system.
solution.get())
265 (solution_vectors->find(&system)->second);
306 System & system = *system_list[i];
311 projected_solutions[i]->init(system.
solution->size(),
true,
SERIAL);
312 system.
solution->localize(*projected_solutions[i],
317 bool solve_adjoint =
false;
318 if (solution_vectors)
320 System * sys = system_list[0];
321 libmesh_assert (solution_vectors->find(sys) !=
322 solution_vectors->end());
324 for (
unsigned int j=0; j != sys->
n_qois(); ++j)
326 std::ostringstream adjoint_name;
327 adjoint_name <<
"adjoint_solution" << j;
331 solve_adjoint =
true;
343 for (
unsigned int i=0; i != es.
n_systems(); ++i)
347 if (!solution_vectors)
351 libmesh_assert_equal_to (solution_vectors->size(), es.
n_systems());
352 libmesh_assert (solution_vectors->find(system_list[0]) !=
353 solution_vectors->end());
354 libmesh_assert(solve_adjoint ||
355 (solution_vectors->find(system_list[0])->second ==
356 system_list[0]->solution.get()) ||
357 !solution_vectors->find(system_list[0])->second);
360 for (
const auto & sys : system_list)
362 libmesh_assert (solution_vectors->find(sys) !=
363 solution_vectors->end());
367 bool found_vec =
false;
368 for (
unsigned int j=0; j != sys->n_qois(); ++j)
370 std::ostringstream adjoint_name;
371 adjoint_name <<
"adjoint_solution" << j;
373 if (vec == sys->request_vector(adjoint_name.str()))
379 libmesh_assert(found_vec);
382 libmesh_assert(vec == sys->solution.get() || !vec);
388 std::vector<unsigned int> adjs(system_list.size(),
393 System * sys = system_list[i];
394 libmesh_assert (solution_vectors->find(sys) !=
395 solution_vectors->end());
397 for (
unsigned int j=0, n_qois = sys->
n_qois();
400 std::ostringstream adjoint_name;
401 adjoint_name <<
"adjoint_solution" << j;
419 system_list[i]->get_adjoint_solution(adjs[i]).swap(*system_list[i]->solution);
420 system_list[i]->update();
429 System * sys = system_list[0];
436 if (!solution_vectors)
440 libmesh_assert (solution_vectors->find(sys) !=
441 solution_vectors->end());
445 libmesh_assert(solve_adjoint ||
446 (solution_vectors->find(sys)->second ==
448 !solution_vectors->find(sys)->second);
453 for (
unsigned int j=0, n_qois = sys->
n_qois();
456 std::ostringstream adjoint_name;
457 adjoint_name <<
"adjoint_solution" << j;
483 System & system = *system_list[sysnum];
490 _error_norms->find(&system)->second;
495 for (
unsigned int var=0; var<
n_vars; var++)
501 libmesh_assert(errors_per_cell);
503 (*errors_per_cell)[std::make_pair(&system,var)];
514 fe->attach_quadrature_rule (qrule.get());
516 const std::vector<Real> & JxW = fe->get_JxW();
517 const std::vector<std::vector<Real>> & phi = fe->get_phi();
518 const std::vector<std::vector<RealGradient>> & dphi =
520 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 521 const std::vector<std::vector<RealTensor>> & d2phi =
526 std::vector<dof_id_type> dof_indices;
530 for (
const auto & elem :
mesh.active_local_element_ptr_range())
533 const Elem * coarse = elem;
535 while (e_id >= max_coarse_elem_id)
537 libmesh_assert (coarse->
parent());
538 coarse = coarse->
parent();
542 Real L2normsq = 0., H1seminormsq = 0.;
543 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 544 Real H2seminormsq = 0.;
555 const unsigned int n_qp = qrule->n_points();
558 const unsigned int n_sf =
559 cast_int<unsigned int>(dof_indices.size());
564 for (
unsigned int qp=0; qp<n_qp; qp++)
566 Number u_fine = 0., u_coarse = 0.;
568 Gradient grad_u_fine, grad_u_coarse;
569 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 570 Tensor grad2_u_fine, grad2_u_coarse;
577 for (
unsigned int i=0; i<n_sf; i++)
580 u_coarse += phi[i][qp]*(*projected_solution) (dof_indices[i]);
582 grad_u_coarse += dphi[i][qp]*(*projected_solution) (dof_indices[i]);
583 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 585 grad2_u_coarse += d2phi[i][qp]*(*projected_solution) (dof_indices[i]);
590 const Number val_error = u_fine - u_coarse;
593 if (system_i_norm.
type(var) ==
L2 ||
594 system_i_norm.
type(var) ==
H1 ||
595 system_i_norm.
type(var) ==
H2)
597 L2normsq += JxW[qp] * system_i_norm.
weight_sq(var) *
599 libmesh_assert_greater_equal (L2normsq, 0.);
605 if (system_i_norm.
type(var) ==
H1 ||
606 system_i_norm.
type(var) ==
H2 ||
609 Gradient grad_error = grad_u_fine - grad_u_coarse;
611 H1seminormsq += JxW[qp] * system_i_norm.
weight_sq(var) *
613 libmesh_assert_greater_equal (H1seminormsq, 0.);
618 if (system_i_norm.
type(var) ==
H2 ||
621 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 622 Tensor grad2_error = grad2_u_fine - grad2_u_coarse;
624 H2seminormsq += JxW[qp] * system_i_norm.
weight_sq(var) *
626 libmesh_assert_greater_equal (H2seminormsq, 0.);
629 (
"libMesh was not configured with --enable-second");
634 if (system_i_norm.
type(var) ==
L2 ||
635 system_i_norm.
type(var) ==
H1 ||
636 system_i_norm.
type(var) ==
H2)
639 if (system_i_norm.
type(var) ==
H1 ||
640 system_i_norm.
type(var) ==
H2 ||
645 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 646 if (system_i_norm.
type(var) ==
H2 ||
677 libmesh_assert_equal_to (n_coarse_elem,
mesh.n_elem());
693 LOG_SCOPE(
"std::sqrt()",
"UniformRefinementEstimator");
694 for (
auto & val : *error_per_cell)
696 val = std::sqrt(val);
700 for (
const auto & pr : *errors_per_cell)
707 LOG_SCOPE(
"std::sqrt()",
"UniformRefinementEstimator");
710 val = std::sqrt(val);
717 System & system = *system_list[i];
722 *system.
solution = *coarse_solutions[i];
729 const std::string & var_name = vec->first;
731 system.
get_vector(var_name) = *coarse_vectors[i][var_name];
733 coarse_vectors[i][var_name]->
clear();
738 mesh.partitioner().reset(old_partitioner.release());
739 mesh.allow_renumbering(old_renumbering_setting);
744 #endif // #ifdef LIBMESH_ENABLE_AMR Manages the family, order, etc. parameters for a given FE.
vectors_iterator vectors_end()
void uniformly_p_refine(unsigned int n=1)
Manages multiples systems of equations.
const Elem * parent() const
unsigned int n_systems() const
const unsigned int invalid_uint
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
unsigned int n_qois() const
virtual void disable_cache()
const T_sys & get_system(const std::string &name) const
const FEType & variable_type(const unsigned int c) const
const EquationSystems & get_equation_systems() const
The base class for all geometric element types.
IntRange< std::size_t > index_range(const std::vector< T > &vec)
void uniformly_p_coarsen(unsigned int n=1)
const Parallel::Communicator & comm() const
virtual void adjoint_solve(const QoISet &qoi_indices=QoISet())
const unsigned int n_vars
vectors_iterator vectors_begin()
DIE A HORRIBLE DEATH HERE typedef float ErrorVectorReal
Number current_solution(const dof_id_type global_dof_number) const
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
Responsible for mesh refinement algorithms and data.
Manages consistently variables, degrees of freedom, and coefficient vectors.
FEMNormType type(unsigned int var) const
void uniformly_coarsen(unsigned int n=1)
const NumericVector< Number > & get_vector(const std::string &vec_name) const
std::unique_ptr< NumericVector< Number > > solution
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
virtual std::pair< unsigned int, Real > adjoint_solve(const QoISet &qoi_indices=QoISet())
std::map< std::pair< const System *, unsigned int >, ErrorVector * > ErrorMap
void reduce_error(std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm) const
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
std::map< std::string, NumericVector< Number > * >::iterator vectors_iterator
Real weight_sq(unsigned int var) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void swap(NumericVector< T > &v)
const NumericVector< Number > * request_vector(const std::string &vec_name) const
bool & project_solution_on_reinit(void)
const MeshBase & get_mesh() const
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
std::unique_ptr< NumericVector< Number > > current_local_solution
unsigned int n_vars() const
const DofMap & get_dof_map() const
const std::vector< dof_id_type > & get_send_list() const
void uniformly_refine(unsigned int n=1)