25 #include "libmesh/dof_map.h"
26 #include "libmesh/elem.h"
27 #include "libmesh/equation_systems.h"
28 #include "libmesh/error_vector.h"
29 #include "libmesh/fe.h"
30 #include "libmesh/libmesh_common.h"
31 #include "libmesh/libmesh_logging.h"
32 #include "libmesh/mesh_base.h"
33 #include "libmesh/mesh_refinement.h"
34 #include "libmesh/numeric_vector.h"
35 #include "libmesh/quadrature.h"
36 #include "libmesh/system.h"
37 #include "libmesh/uniform_refinement_estimator.h"
38 #include "libmesh/partitioner.h"
39 #include "libmesh/tensor_tools.h"
40 #include "libmesh/enum_error_estimator_type.h"
41 #include "libmesh/enum_norm_type.h"
42 #include "libmesh/int_range.h"
43 #include "libmesh/auto_ptr.h"
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 auto error_norms = libmesh_make_unique<std::map<const System *, SystemNorm>>();
135 libmesh_assert_equal_to (&(_system->get_equation_systems()), _es);
140 system_list.push_back(const_cast<System *>(&(_es->
get_system(i))));
154 _error_norms = error_norms.get();
161 std::vector<Real> weights(
n_vars, 0.0);
162 for (
unsigned int v = 0; v !=
n_vars; ++v)
164 if (errors_per_cell->find(std::make_pair(&sys, v)) ==
165 errors_per_cell->end())
170 (*error_norms)[&sys] =
180 system_list.push_back(const_cast<System *>(_system));
184 _error_norms = error_norms.get();
196 const unsigned int dim =
mesh.mesh_dimension();
202 error_per_cell->clear();
203 error_per_cell->resize (
mesh.max_elem_id(), 0.);
208 for (
const auto & pr : *errors_per_cell)
212 e->resize(
mesh.max_elem_id(), 0.);
217 std::vector<std::map<std::string, std::unique_ptr<NumericVector<Number>>>> coarse_vectors(system_list.size());
218 std::vector<std::unique_ptr<NumericVector<Number>>> coarse_solutions(system_list.size());
219 std::vector<std::unique_ptr<NumericVector<Number>>> coarse_local_solutions(system_list.size());
221 std::vector<std::unique_ptr<NumericVector<Number>>> projected_solutions(system_list.size());
224 std::vector<bool> old_projection_settings(system_list.size());
227 std::unique_ptr<Partitioner> old_partitioner(
mesh.partitioner().release());
230 const bool old_renumbering_setting =
mesh.allow_renumbering();
231 mesh.allow_renumbering(
false);
235 System & system = *system_list[i];
239 _error_norms->end());
242 coarse_solutions[i] = system.
solution->clone();
250 const std::string & var_name = vec->first;
252 coarse_vectors[i][var_name] = vec->second->clone();
256 if (solution_vectors &&
257 solution_vectors->find(&system) != solution_vectors->end() &&
258 solution_vectors->find(&system)->second &&
259 solution_vectors->find(&system)->second != system.
solution.get())
263 (solution_vectors->find(&system)->second);
304 System & system = *system_list[i];
309 projected_solutions[i]->init(system.
solution->size(),
true,
SERIAL);
310 system.
solution->localize(*projected_solutions[i],
315 bool solve_adjoint =
false;
316 if (solution_vectors)
318 System * sys = system_list[0];
320 solution_vectors->end());
324 std::ostringstream adjoint_name;
325 adjoint_name <<
"adjoint_solution" << j;
329 solve_adjoint =
true;
345 if (!solution_vectors)
349 libmesh_assert_equal_to (solution_vectors->size(), es.
n_systems());
351 solution_vectors->end());
353 (solution_vectors->find(system_list[0])->second ==
354 system_list[0]->solution.get()) ||
355 !solution_vectors->find(system_list[0])->second);
358 for (
const auto & sys : system_list)
361 solution_vectors->end());
365 bool found_vec =
false;
368 std::ostringstream adjoint_name;
369 adjoint_name <<
"adjoint_solution" << j;
371 if (vec == sys->request_vector(adjoint_name.str()))
386 std::vector<unsigned int> adjs(system_list.size(),
391 System * sys = system_list[i];
393 solution_vectors->end());
397 std::ostringstream adjoint_name;
398 adjoint_name <<
"adjoint_solution" << j;
416 system_list[i]->get_adjoint_solution(adjs[i]).swap(*system_list[i]->solution);
417 system_list[i]->update();
426 System * sys = system_list[0];
433 if (!solution_vectors)
438 solution_vectors->end());
443 (solution_vectors->find(sys)->second ==
445 !solution_vectors->find(sys)->second);
450 for (
unsigned int j=0, n_qois = sys->
n_qois();
453 std::ostringstream adjoint_name;
454 adjoint_name <<
"adjoint_solution" << j;
480 System & system = *system_list[sysnum];
487 _error_norms->find(&system)->second;
492 for (
unsigned int var=0; var<
n_vars; var++)
500 (*errors_per_cell)[std::make_pair(&system,var)];
511 fe->attach_quadrature_rule (qrule.get());
513 const std::vector<Real> & JxW = fe->get_JxW();
514 const std::vector<std::vector<Real>> & phi = fe->get_phi();
515 const std::vector<std::vector<RealGradient>> & dphi =
517 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
518 const std::vector<std::vector<RealTensor>> & d2phi =
523 std::vector<dof_id_type> dof_indices;
527 for (
const auto & elem :
mesh.active_local_element_ptr_range())
530 const Elem * coarse = elem;
532 while (e_id >= max_coarse_elem_id)
535 coarse = coarse->
parent();
539 Real L2normsq = 0., H1seminormsq = 0.;
540 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
541 Real H2seminormsq = 0.;
552 const unsigned int n_qp = qrule->n_points();
555 const unsigned int n_sf =
556 cast_int<unsigned int>(dof_indices.size());
561 for (
unsigned int qp=0; qp<n_qp; qp++)
563 Number u_fine = 0., u_coarse = 0.;
565 Gradient grad_u_fine, grad_u_coarse;
566 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
567 Tensor grad2_u_fine, grad2_u_coarse;
574 for (
unsigned int i=0; i<n_sf; i++)
577 u_coarse += phi[i][qp]*(*projected_solution) (dof_indices[i]);
579 grad_u_coarse += dphi[i][qp]*(*projected_solution) (dof_indices[i]);
580 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
582 grad2_u_coarse += d2phi[i][qp]*(*projected_solution) (dof_indices[i]);
587 const Number val_error = u_fine - u_coarse;
590 if (system_i_norm.
type(var) ==
L2 ||
591 system_i_norm.
type(var) ==
H1 ||
592 system_i_norm.
type(var) ==
H2)
594 L2normsq += JxW[qp] * system_i_norm.
weight_sq(var) *
596 libmesh_assert_greater_equal (L2normsq, 0.);
602 if (system_i_norm.
type(var) ==
H1 ||
603 system_i_norm.
type(var) ==
H2 ||
606 Gradient grad_error = grad_u_fine - grad_u_coarse;
608 H1seminormsq += JxW[qp] * system_i_norm.
weight_sq(var) *
610 libmesh_assert_greater_equal (H1seminormsq, 0.);
615 if (system_i_norm.
type(var) ==
H2 ||
618 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
619 Tensor grad2_error = grad2_u_fine - grad2_u_coarse;
621 H2seminormsq += JxW[qp] * system_i_norm.
weight_sq(var) *
623 libmesh_assert_greater_equal (H2seminormsq, 0.);
626 (
"libMesh was not configured with --enable-second");
631 if (system_i_norm.
type(var) ==
L2 ||
632 system_i_norm.
type(var) ==
H1 ||
633 system_i_norm.
type(var) ==
H2)
635 static_cast<ErrorVectorReal>(L2normsq);
636 if (system_i_norm.
type(var) ==
H1 ||
637 system_i_norm.
type(var) ==
H2 ||
640 static_cast<ErrorVectorReal>(H1seminormsq);
642 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
643 if (system_i_norm.
type(var) ==
H2 ||
646 static_cast<ErrorVectorReal>(H2seminormsq);
674 libmesh_assert_equal_to (n_coarse_elem,
mesh.n_elem());
690 LOG_SCOPE(
"std::sqrt()",
"UniformRefinementEstimator");
691 for (
auto & val : *error_per_cell)
697 for (
const auto & pr : *errors_per_cell)
704 LOG_SCOPE(
"std::sqrt()",
"UniformRefinementEstimator");
714 System & system = *system_list[i];
719 *system.
solution = *coarse_solutions[i];
726 const std::string & var_name = vec->first;
728 system.
get_vector(var_name) = *coarse_vectors[i][var_name];
730 coarse_vectors[i][var_name]->
clear();
735 mesh.partitioner().reset(old_partitioner.release());
736 mesh.allow_renumbering(old_renumbering_setting);
741 #endif // #ifdef LIBMESH_ENABLE_AMR