26 #include "libmesh/libmesh_common.h"
27 #include "libmesh/patch_recovery_error_estimator.h"
28 #include "libmesh/dof_map.h"
29 #include "libmesh/fe_base.h"
30 #include "libmesh/dense_matrix.h"
31 #include "libmesh/dense_vector.h"
32 #include "libmesh/error_vector.h"
33 #include "libmesh/libmesh_logging.h"
34 #include "libmesh/elem.h"
35 #include "libmesh/patch.h"
36 #include "libmesh/quadrature_grid.h"
37 #include "libmesh/system.h"
38 #include "libmesh/mesh_base.h"
39 #include "libmesh/numeric_vector.h"
40 #include "libmesh/tensor_value.h"
41 #include "libmesh/threads.h"
42 #include "libmesh/tensor_tools.h"
43 #include "libmesh/enum_error_estimator_type.h"
44 #include "libmesh/enum_norm_type.h"
45 #include "libmesh/int_range.h"
67 target_patch_size(20),
68 patch_growth_strategy(&
Patch::add_local_face_neighbors),
79 const unsigned int matsize)
81 std::vector<Real> psi;
83 unsigned int npows = order+1;
84 std::vector<Real> xpow(npows,1.), ypow, zpow;
88 xpow[i] = xpow[i-1] * x;
93 ypow.resize(npows,1.);
95 ypow[i] = ypow[i-1] * y;
100 zpow.resize(npows,1.);
102 zpow[i] = zpow[i-1] * z;
107 for (
unsigned int poly_deg=0; poly_deg <= static_cast<unsigned int>(order) ; poly_deg++)
115 for (
int xexp=poly_deg; xexp >= 0; xexp--)
116 for (
int yexp=poly_deg-xexp; yexp >= 0; yexp--)
118 int zexp = poly_deg - xexp - yexp;
119 psi.push_back(xpow[xexp]*ypow[yexp]*zpow[zexp]);
127 for (
int xexp=poly_deg; xexp >= 0; xexp--)
129 int yexp = poly_deg - xexp;
130 psi.push_back(xpow[xexp]*ypow[yexp]);
139 psi.push_back(xpow[xexp]);
144 libmesh_error_msg(
"Invalid dimension dim " <<
dim);
158 LOG_SCOPE(
"estimate_error()",
"PatchRecoveryErrorEstimator");
165 error_per_cell.resize (
mesh.max_elem_id());
166 std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
170 if (solution_vector && solution_vector != system.
solution.get())
174 System & sys = const_cast<System &>(system);
175 newsol->
swap(*sys.solution);
183 mesh.active_local_elements_end(),
198 if (solution_vector && solution_vector != system.
solution.get())
202 System & sys = const_cast<System &>(system);
203 newsol->
swap(*sys.solution);
216 const unsigned int dim =
mesh.mesh_dimension();
226 for (
const auto & elem : range)
251 std::vector<Real> new_error_per_cell(1, 0.);
253 new_error_per_cell.resize(patch.size(), 0.);
257 for (
unsigned int var=0; var<
n_vars; var++)
259 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
261 bool is_valid_norm_type =
288 bool is_valid_norm_combo =
317 const Order element_order = static_cast<Order>
318 (fe_type.
order + elem->p_level());
327 fe->attach_quadrature_rule (qrule.get());
330 const std::vector<Real> & JxW = fe->get_JxW();
331 const std::vector<Point> & q_point = fe->get_xyz();
337 const std::vector<std::vector<Real>> * phi =
nullptr;
345 phi = &(fe->get_phi());
347 const std::vector<std::vector<RealGradient>> * dphi =
nullptr;
353 dphi = &(fe->get_dphi());
355 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
356 const std::vector<std::vector<RealTensor>> * d2phi =
nullptr;
359 d2phi = &(fe->get_d2phi());
363 std::vector<dof_id_type> dof_indices;
367 unsigned int matsize = element_order + 1;
370 matsize *= (element_order + 2);
375 matsize *= (element_order + 3);
406 libmesh_assert_greater (LIBMESH_DIM, 1);
411 libmesh_assert_greater (LIBMESH_DIM, 2);
430 for (
const auto & e_p : patch)
438 libmesh_assert_equal_to (dof_indices.size(), phi->size());
440 const unsigned int n_dofs =
441 cast_int<unsigned int>(dof_indices.size());
442 const unsigned int n_qp = qrule->n_points();
446 for (
unsigned int qp=0; qp<n_qp; qp++)
449 std::vector<Real> psi(
specpoly(
dim, element_order, q_point[qp], matsize));
451 const unsigned int psi_size = cast_int<unsigned int>(psi.size());
454 const unsigned int m = Kp.
m(), n = Kp.
n();
455 for (
unsigned int i=0; i<m; i++)
456 for (
unsigned int j=0; j<n; j++)
457 Kp(i,j) += JxW[qp]*psi[i]*psi[j];
466 for (
unsigned int i=0; i<n_dofs; i++)
470 for (
unsigned int i=0; i != psi_size; i++)
471 F(i) += JxW[qp]*u_h*psi[i];
481 for (
unsigned int i=0; i<n_dofs; i++)
486 for (
unsigned int i=0; i != psi_size; i++)
488 Fx(i) += JxW[qp]*grad_u_h(0)*psi[i];
490 Fy(i) += JxW[qp]*grad_u_h(1)*psi[i];
493 Fz(i) += JxW[qp]*grad_u_h(2)*psi[i];
503 for (
unsigned int i=0; i<n_dofs; i++)
508 for (
unsigned int i=0; i != psi_size; i++)
510 Fx(i) += JxW[qp]*grad_u_h(0)*psi[i];
520 for (
unsigned int i=0; i<n_dofs; i++)
525 for (
unsigned int i=0; i != psi_size; i++)
527 Fy(i) += JxW[qp]*grad_u_h(1)*psi[i];
530 #endif // LIBMESH_DIM > 1
538 for (
unsigned int i=0; i<n_dofs; i++)
543 for (
unsigned int i=0; i != psi_size; i++)
545 Fz(i) += JxW[qp]*grad_u_h(2)*psi[i];
548 #endif // LIBMESH_DIM > 2
552 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
557 for (
unsigned int i=0; i<n_dofs; i++)
562 for (
unsigned int i=0; i != psi_size; i++)
564 Fx(i) += JxW[qp]*hess_u_h(0,0)*psi[i];
566 Fy(i) += JxW[qp]*hess_u_h(1,1)*psi[i];
567 Fxy(i) += JxW[qp]*hess_u_h(0,1)*psi[i];
570 Fz(i) += JxW[qp]*hess_u_h(2,2)*psi[i];
571 Fxz(i) += JxW[qp]*hess_u_h(0,2)*psi[i];
572 Fyz(i) += JxW[qp]*hess_u_h(1,2)*psi[i];
576 libmesh_error_msg(
"ERROR: --enable-second-derivatives is required \nfor _sobolev_order == 2!");
580 libmesh_error_msg(
"Unsupported error norm type!");
637 Patch::const_iterator patch_re_it;
638 Patch::const_iterator patch_re_end;
646 patch_re_it = patch.begin();
647 patch_re_end = patch.end();
657 patch_re_it = patch_re.begin();
658 patch_re_end = patch_re.end();
668 for (
unsigned int e = 0 ; patch_re_it != patch_re_end; ++patch_re_it, ++e)
673 const Elem * e_p = *patch_re_it;
692 libmesh_assert_equal_to (dof_indices.size(), phi->size());
695 const unsigned int n_dofs =
696 cast_int<unsigned int>(dof_indices.size());
699 Real element_error = 0;
709 fe->attach_quadrature_rule (&samprule);
712 const unsigned int n_sp =
713 cast_int<unsigned int>(JxW.size());
716 for (
unsigned int sp=0; sp<n_sp; sp++)
720 std::vector<Number> temperr(6,0.0);
728 for (
unsigned int i=0; i<n_dofs; i++)
732 std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
733 for (
unsigned int i=0; i<matsize; i++)
735 temperr[0] += psi[i]*Pu_h(i);
746 for (
unsigned int i=0; i<n_dofs; i++)
751 std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
753 for (
unsigned int i=0; i<matsize; i++)
755 temperr[0] += psi[i]*Pu_x_h(i);
757 temperr[1] += psi[i]*Pu_y_h(i);
760 temperr[2] += psi[i]*Pu_z_h(i);
763 temperr[0] -= grad_u_h(0);
765 temperr[1] -= grad_u_h(1);
768 temperr[2] -= grad_u_h(2);
776 for (
unsigned int i=0; i<n_dofs; i++)
781 std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
782 for (
unsigned int i=0; i<matsize; i++)
784 temperr[0] += psi[i]*Pu_x_h(i);
787 temperr[0] -= grad_u_h(0);
795 for (
unsigned int i=0; i<n_dofs; i++)
800 std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
801 for (
unsigned int i=0; i<matsize; i++)
803 temperr[1] += psi[i]*Pu_y_h(i);
806 temperr[1] -= grad_u_h(1);
808 #endif // LIBMESH_DIM > 1
815 for (
unsigned int i=0; i<n_dofs; i++)
820 std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
821 for (
unsigned int i=0; i<matsize; i++)
823 temperr[2] += psi[i]*Pu_z_h(i);
826 temperr[2] -= grad_u_h(2);
828 #endif // LIBMESH_DIM > 2
832 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
836 for (
unsigned int i=0; i<n_dofs; i++)
841 std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
842 for (
unsigned int i=0; i<matsize; i++)
844 temperr[0] += psi[i]*Pu_x_h(i);
846 temperr[1] += psi[i]*Pu_y_h(i);
847 temperr[3] += psi[i]*Pu_xy_h(i);
850 temperr[2] += psi[i]*Pu_z_h(i);
851 temperr[4] += psi[i]*Pu_xz_h(i);
852 temperr[5] += psi[i]*Pu_yz_h(i);
856 temperr[0] -= hess_u_h(0,0);
858 temperr[1] -= hess_u_h(1,1);
859 temperr[3] -= hess_u_h(0,1);
862 temperr[2] -= hess_u_h(2,2);
863 temperr[4] -= hess_u_h(0,2);
864 temperr[5] -= hess_u_h(1,2);
867 libmesh_error_msg(
"ERROR: --enable-second-derivatives is required \nfor _sobolev_order == 2!");
875 element_error = std::max(element_error,
std::abs(temperr[0]));
877 for (
unsigned int i=0; i != LIBMESH_DIM; ++i)
878 element_error = std::max(element_error,
std::abs(temperr[i]));
880 for (
unsigned int i=0; i != 6; ++i)
881 element_error = std::max(element_error,
std::abs(temperr[i]));
885 for (
unsigned int i=0; i != LIBMESH_DIM; ++i)
895 for (
unsigned int i=0; i != LIBMESH_DIM; ++i)
898 for (
unsigned int i=3; i != 6; ++i)
916 libmesh_error_msg(
"Unsupported error norm type!");
927 Patch::const_iterator patch_re_it;
928 Patch::const_iterator patch_re_end;
931 Patch current_elem_patch(
mesh.processor_id());
936 patch_re_it = patch.begin();
937 patch_re_end = patch.end();
947 patch_re_it = current_elem_patch.begin();
948 patch_re_end = current_elem_patch.end();
952 for (
unsigned int i = 0 ; patch_re_it != patch_re_end; ++patch_re_it, ++i)
955 const Elem * e_p = *patch_re_it;
971 static_cast<ErrorVectorReal>(
std::sqrt(new_error_per_cell[i]));
981 static_cast<ErrorVectorReal>(new_error_per_cell[i]);