26 #include "libmesh/dense_matrix.h" 
   27 #include "libmesh/dense_vector.h" 
   28 #include "libmesh/dof_map.h" 
   29 #include "libmesh/elem.h" 
   30 #include "libmesh/error_vector.h" 
   31 #include "libmesh/fe_base.h" 
   32 #include "libmesh/fem_context.h" 
   33 #include "libmesh/libmesh_logging.h" 
   34 #include "libmesh/mesh_base.h" 
   35 #include "libmesh/numeric_vector.h" 
   36 #include "libmesh/quadrature_grid.h" 
   37 #include "libmesh/system.h" 
   38 #include "libmesh/tensor_value.h" 
   39 #include "libmesh/threads.h" 
   40 #include "libmesh/weighted_patch_recovery_error_estimator.h" 
   41 #include "libmesh/enum_error_estimator_type.h" 
   42 #include "libmesh/enum_order.h" 
   43 #include "libmesh/enum_norm_type.h" 
   60   LOG_SCOPE(
"estimate_error()", 
"WeightedPatchRecoveryErrorEstimator");
 
   67   error_per_cell.resize (
mesh.max_elem_id());
 
   68   std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
 
   72   if (solution_vector && solution_vector != system.
solution.get())
 
   76       System & sys = const_cast<System &>(system);
 
   77       newsol->
swap(*sys.solution);
 
   85                                         mesh.active_local_elements_end(),
 
  100   if (solution_vector && solution_vector != system.
solution.get())
 
  104       System & sys = const_cast<System &>(system);
 
  105       newsol->
swap(*sys.solution);
 
  118   const unsigned int dim = 
mesh.mesh_dimension();
 
  128   for (
const auto & elem : range)
 
  153       std::vector<Real> new_error_per_cell(1, 0.);
 
  155         new_error_per_cell.resize(patch.size(), 0.);
 
  159       for (
unsigned int var=0; var<
n_vars; var++)
 
  161 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  163           bool is_valid_norm_type =
 
  189               bool is_valid_norm_combo =
 
  218           const Order element_order  = static_cast<Order>
 
  219             (fe_type.
order + elem->p_level());
 
  228           fe->attach_quadrature_rule (qrule.get());
 
  231           const std::vector<Real> & JxW = fe->get_JxW();
 
  232           const std::vector<Point> & q_point = fe->get_xyz();
 
  238           const std::vector<std::vector<Real>> * phi = 
nullptr;
 
  246             phi = &(fe->get_phi());
 
  248           const std::vector<std::vector<RealGradient>> * dphi = 
nullptr;
 
  254             dphi = &(fe->get_dphi());
 
  256 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  257           const std::vector<std::vector<RealTensor>> * d2phi = 
nullptr;
 
  260             d2phi = &(fe->get_d2phi());
 
  264           std::vector<dof_id_type> dof_indices;
 
  268           unsigned int matsize = element_order + 1;
 
  271               matsize *= (element_order + 2);
 
  276               matsize *= (element_order + 3);
 
  332           for (
const auto & e_p : patch)
 
  342               const unsigned int n_dofs =
 
  343                 cast_int<unsigned int>(dof_indices.size());
 
  344               const unsigned int n_qp   = qrule->n_points();
 
  348               for (
unsigned int qp=0; qp<n_qp; qp++)
 
  351                   std::vector<Real> psi(
specpoly(
dim, element_order, q_point[qp], matsize));
 
  353                   const unsigned int psi_size = cast_int<unsigned int>(psi.size());
 
  356                   const unsigned int m = Kp.
m(), n = Kp.
n();
 
  357                   for (
unsigned int i=0; i<m; i++)
 
  358                     for (
unsigned int j=0; j<n; j++)
 
  359                       Kp(i,j) += JxW[qp]*psi[i]*psi[j];
 
  368                       for (
unsigned int i=0; i<n_dofs; i++)
 
  372                       for (
unsigned int i=0; i != psi_size; i++)
 
  373                         F(i) = JxW[qp]*u_h*psi[i];
 
  383                       for (std::size_t i=0; i<n_dofs; i++)
 
  390                       for (
unsigned int i=0; i != psi_size; i++)
 
  392                           Fx(i) += JxW[qp]*grad_u_h(0)*psi[i];
 
  394                           Fy(i) += JxW[qp]*grad_u_h(1)*psi[i];
 
  397                           Fz(i) += JxW[qp]*grad_u_h(2)*psi[i];
 
  407                       for (
unsigned int i=0; i<n_dofs; i++)
 
  414                       for (
unsigned int i=0; i != psi_size; i++)
 
  416                           Fx(i) += JxW[qp]*grad_u_h(0)*psi[i];
 
  426                       for (
unsigned int i=0; i<n_dofs; i++)
 
  433                       for (
unsigned int i=0; i != psi_size; i++)
 
  435                           Fy(i) += JxW[qp]*grad_u_h(1)*psi[i];
 
  438 #endif // LIBMESH_DIM > 1 
  446                       for (
unsigned int i=0; i<n_dofs; i++)
 
  453                       for (
unsigned int i=0; i != psi_size; i++)
 
  455                           Fz(i) += JxW[qp]*grad_u_h(2)*psi[i];
 
  458 #endif // LIBMESH_DIM > 2 
  462 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  467                       for (
unsigned int i=0; i<n_dofs; i++)
 
  474                       for (
unsigned int i=0; i != psi_size; i++)
 
  476                           Fx(i)  += JxW[qp]*hess_u_h(0,0)*psi[i];
 
  478                           Fy(i)  += JxW[qp]*hess_u_h(1,1)*psi[i];
 
  479                           Fxy(i) += JxW[qp]*hess_u_h(0,1)*psi[i];
 
  482                           Fz(i)  += JxW[qp]*hess_u_h(2,2)*psi[i];
 
  483                           Fxz(i) += JxW[qp]*hess_u_h(0,2)*psi[i];
 
  484                           Fyz(i) += JxW[qp]*hess_u_h(1,2)*psi[i];
 
  488                       libmesh_error_msg(
"ERROR:  --enable-second-derivatives is required \nfor _sobolev_order == 2!");
 
  492                     libmesh_error_msg(
"Unsupported error norm type!");
 
  549           Patch::const_iterator patch_re_it;
 
  550           Patch::const_iterator patch_re_end;
 
  558               patch_re_it  = patch.begin();
 
  559               patch_re_end = patch.end();
 
  569               patch_re_it = patch_re.begin();
 
  570               patch_re_end = patch_re.end();
 
  585           for (
unsigned int e = 0 ; patch_re_it != patch_re_end; ++patch_re_it, ++e)
 
  590               const Elem * e_p = *patch_re_it;
 
  615               const unsigned int n_dofs =
 
  616                 cast_int<unsigned int>(dof_indices.size());
 
  619               Real element_error = 0;
 
  629                 fe->attach_quadrature_rule (&samprule);
 
  632               const unsigned int n_sp =
 
  633                 cast_int<unsigned int>(JxW.size());
 
  636               for (
unsigned int sp=0; sp<n_sp; sp++)
 
  640                   std::vector<Number> temperr(6,0.0); 
 
  648                       for (
unsigned int i=0; i<n_dofs; i++)
 
  652                       std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
 
  653                       for (
unsigned int i=0; i<matsize; i++)
 
  655                           temperr[0] += psi[i]*Pu_h(i);
 
  666                       for (
unsigned int i=0; i<n_dofs; i++)
 
  671                       std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
 
  673                       for (
unsigned int i=0; i<matsize; i++)
 
  675                           temperr[0] += psi[i]*Pu_x_h(i);
 
  677                           temperr[1] += psi[i]*Pu_y_h(i);
 
  680                           temperr[2] += psi[i]*Pu_z_h(i);
 
  683                       temperr[0] -= grad_u_h(0);
 
  685                       temperr[1] -= grad_u_h(1);
 
  688                       temperr[2] -= grad_u_h(2);
 
  696                       for (
unsigned int i=0; i<n_dofs; i++)
 
  701                       std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
 
  702                       for (
unsigned int i=0; i<matsize; i++)
 
  704                           temperr[0] += psi[i]*Pu_x_h(i);
 
  707                       temperr[0] -= grad_u_h(0);
 
  715                       for (
unsigned int i=0; i<n_dofs; i++)
 
  720                       std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
 
  721                       for (
unsigned int i=0; i<matsize; i++)
 
  723                           temperr[1] += psi[i]*Pu_y_h(i);
 
  726                       temperr[1] -= grad_u_h(1);
 
  728 #endif // LIBMESH_DIM > 1 
  735                       for (
unsigned int i=0; i<n_dofs; i++)
 
  740                       std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
 
  741                       for (
unsigned int i=0; i<matsize; i++)
 
  743                           temperr[2] += psi[i]*Pu_z_h(i);
 
  746                       temperr[2] -= grad_u_h(2);
 
  748 #endif // LIBMESH_DIM > 2 
  752 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  756                       for (
unsigned int i=0; i<n_dofs; i++)
 
  761                       std::vector<Real> psi(
specpoly(
dim, element_order, q_point[sp], matsize));
 
  762                       for (
unsigned int i=0; i<matsize; i++)
 
  764                           temperr[0] += psi[i]*Pu_x_h(i);
 
  766                           temperr[1] += psi[i]*Pu_y_h(i);
 
  767                           temperr[3] += psi[i]*Pu_xy_h(i);
 
  770                           temperr[2] += psi[i]*Pu_z_h(i);
 
  771                           temperr[4] += psi[i]*Pu_xz_h(i);
 
  772                           temperr[5] += psi[i]*Pu_yz_h(i);
 
  776                       temperr[0] -= hess_u_h(0,0);
 
  778                       temperr[1] -= hess_u_h(1,1);
 
  779                       temperr[3] -= hess_u_h(0,1);
 
  782                       temperr[2] -= hess_u_h(2,2);
 
  783                       temperr[4] -= hess_u_h(0,2);
 
  784                       temperr[5] -= hess_u_h(1,2);
 
  787                       libmesh_error_msg(
"ERROR:  --enable-second-derivatives is required \nfor _sobolev_order == 2!");
 
  799                     element_error = std::max(element_error, 
std::abs(
weight*temperr[0]));
 
  801                     for (
unsigned int i=0; i != LIBMESH_DIM; ++i)
 
  802                       element_error = std::max(element_error, 
std::abs(
weight*temperr[i]));
 
  804                     for (
unsigned int i=0; i != 6; ++i)
 
  805                       element_error = std::max(element_error, 
std::abs(
weight*temperr[i]));
 
  809                     for (
unsigned int i=0; i != LIBMESH_DIM; ++i)
 
  819                       for (
unsigned int i=0; i != LIBMESH_DIM; ++i)
 
  822                       for (
unsigned int i=3; i != 6; ++i)
 
  840                 libmesh_error_msg(
"Unsupported error norm type!");
 
  851       Patch::const_iterator patch_re_it;
 
  852       Patch::const_iterator patch_re_end;
 
  855       Patch current_elem_patch(
mesh.processor_id());
 
  860           patch_re_it  = patch.begin();
 
  861           patch_re_end = patch.end();
 
  871           patch_re_it = current_elem_patch.begin();
 
  872           patch_re_end = current_elem_patch.end();
 
  876       for (
unsigned int i = 0 ; patch_re_it != patch_re_end; ++patch_re_it, ++i)
 
  879           const Elem * e_p = *patch_re_it;
 
  905                   (new_error_per_cell[i]);