590 LOG_SCOPE (
"operator()",
"GenericProjector");
594 std::unique_ptr<GFunctor> g;
605 DenseMatrix<Real> Ke;
606 DenseVector<FValue> Fe;
608 DenseVector<FValue> Ue;
611 FEMContext context(
system );
619 FEBase * side_fe =
nullptr;
620 FEBase * edge_fe =
nullptr;
622 const std::set<unsigned char> & elem_dims =
623 context.elem_dimensions();
625 for (
const auto & dim : elem_dims)
627 context.get_element_fe( var, fe, dim );
628 if (fe->get_fe_type().family ==
SCALAR)
631 context.get_side_fe( var, side_fe, dim );
633 context.get_edge_fe( var, edge_fe );
659 f.init_context(context);
661 g->init_context(context);
666 for (
const auto & elem : range)
668 unsigned char dim = cast_int<unsigned char>(elem->dim());
670 context.pre_fe_reinit(
system, elem);
674 #ifdef LIBMESH_ENABLE_AMR 679 if (!elem->old_dof_object &&
682 f.is_grid_projection())
684 #endif // LIBMESH_ENABLE_AMR 690 const Variable & variable = dof_map.variable(var);
692 const FEType & base_fe_type = variable.type();
694 FEType fe_type = base_fe_type;
700 if (fe_type.family ==
SCALAR)
705 if (!variable.active_on_subdomain(elem->subdomain_id()))
708 const std::vector<dof_id_type> & dof_indices =
709 context.get_dof_indices(var);
712 const unsigned int n_dofs =
713 cast_int<unsigned int>(dof_indices.size());
715 const unsigned int var_component =
719 Ue.resize (n_dofs); Ue.zero();
722 #ifdef LIBMESH_ENABLE_AMR 723 if (f.is_grid_projection() &&
733 elem->p_level() == 0 &&
740 LOG_SCOPE (
"copy_dofs",
"GenericProjector");
742 f.eval_old_dofs(context, var, Ue.get_values());
744 action.insert(context, var, Ue);
748 #endif // LIBMESH_ENABLE_AMR 751 FEBase * side_fe =
nullptr;
752 FEBase * edge_fe =
nullptr;
754 context.get_element_fe( var, fe, dim );
755 if (fe->get_fe_type().family ==
SCALAR)
758 context.get_side_fe( var, side_fe, dim );
760 context.get_edge_fe( var, edge_fe );
764 std::vector<unsigned int> side_dofs;
767 std::vector<char> dof_is_fixed(n_dofs,
false);
768 std::vector<int> free_dof(n_dofs, 0);
771 const ElemType elem_type = elem->type();
774 const unsigned int n_nodes = elem->n_nodes();
776 START_LOG (
"project_nodes",
"GenericProjector");
782 #ifdef LIBMESH_ENABLE_AMR 784 i_am_child = elem->parent()->which_child_am_i(elem);
798 unsigned int current_dof = 0;
799 for (
unsigned int n=0; n!=
n_nodes; ++n)
803 const unsigned int nc =
806 if (!elem->is_vertex(n) &&
815 libmesh_assert_equal_to (nc, 0);
821 libmesh_assert (!elem->is_vertex(n));
822 libmesh_assert_equal_to(fe_type.family,
LAGRANGE);
828 Ue(current_dof) = f.eval_at_node(context,
833 dof_is_fixed[current_dof] =
true;
836 else if (cont ==
C_ONE)
838 const bool is_parent_vertex = (i_am_child == -1) ||
839 elem->parent()->is_vertex_on_parent(i_am_child, n);
845 f.eval_at_node(context,
850 dof_is_fixed[current_dof] =
true;
852 VectorValue<FValue> grad =
854 g->eval_at_node(context,
859 g->eval_at_point(context,
864 Ue(current_dof) = grad(0);
865 dof_is_fixed[current_dof] =
true;
872 Point nxminus = elem->point(n),
873 nxplus = elem->point(n);
874 nxminus(0) -= delta_x;
875 nxplus(0) += delta_x;
876 VectorValue<FValue> gxminus =
877 g->eval_at_point(context,
881 VectorValue<FValue> gxplus =
882 g->eval_at_point(context,
887 Ue(current_dof) = grad(1);
888 dof_is_fixed[current_dof] =
true;
891 Ue(current_dof) = (gxplus(1) - gxminus(1))
893 dof_is_fixed[current_dof] =
true;
899 Ue(current_dof) = grad(2);
900 dof_is_fixed[current_dof] =
true;
903 Ue(current_dof) = (gxplus(2) - gxminus(2))
905 dof_is_fixed[current_dof] =
true;
908 Point nyminus = elem->point(n),
909 nyplus = elem->point(n);
910 nyminus(1) -= delta_x;
911 nyplus(1) += delta_x;
912 VectorValue<FValue> gyminus =
913 g->eval_at_point(context,
917 VectorValue<FValue> gyplus =
918 g->eval_at_point(context,
923 Ue(current_dof) = (gyplus(2) - gyminus(2))
925 dof_is_fixed[current_dof] =
true;
928 Point nxmym = elem->point(n),
929 nxmyp = elem->point(n),
930 nxpym = elem->point(n),
931 nxpyp = elem->point(n);
940 VectorValue<FValue> gxmym =
941 g->eval_at_point(context,
945 VectorValue<FValue> gxmyp =
946 g->eval_at_point(context,
950 VectorValue<FValue> gxpym =
951 g->eval_at_point(context,
955 VectorValue<FValue> gxpyp =
956 g->eval_at_point(context,
960 FValue gxzplus = (gxpyp(2) - gxmyp(2))
962 FValue gxzminus = (gxpym(2) - gxmym(2))
965 Ue(current_dof) = (gxzplus - gxzminus)
967 dof_is_fixed[current_dof] =
true;
977 libmesh_assert_equal_to (nc, (
unsigned int)(1 + dim));
978 Ue(current_dof) = f.eval_at_node(context,
983 dof_is_fixed[current_dof] =
true;
985 VectorValue<FValue> grad =
987 g->eval_at_node(context,
992 g->eval_at_point(context,
996 for (
unsigned int i=0; i!= dim; ++i)
998 Ue(current_dof) = grad(i);
999 dof_is_fixed[current_dof] =
true;
1005 libmesh_error_msg(
"Unknown continuity " << cont);
1008 STOP_LOG (
"project_nodes",
"GenericProjector");
1010 START_LOG (
"project_edges",
"GenericProjector");
1016 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1023 const std::vector<Point> & xyz_values =
1024 #ifdef LIBMESH_ENABLE_AMR 1029 const std::vector<Real> & JxW =
1030 #ifdef LIBMESH_ENABLE_AMR 1036 const std::vector<std::vector<Real>> & phi =
1037 #ifdef LIBMESH_ENABLE_AMR 1042 const std::vector<std::vector<RealGradient>> * dphi =
nullptr;
1045 #ifdef LIBMESH_ENABLE_AMR 1049 &(edge_fe->get_dphi());
1051 for (
auto e : elem->edge_index_range())
1053 context.edge = cast_int<unsigned char>(e);
1055 #ifdef LIBMESH_ENABLE_AMR 1058 std::vector<Point> fine_points;
1060 std::unique_ptr<FEBase> fine_fe
1063 std::unique_ptr<QBase> qrule
1064 (base_fe_type.default_quadrature_rule(1));
1065 fine_fe->attach_quadrature_rule(qrule.get());
1067 const std::vector<Point> & child_xyz =
1070 for (
unsigned int c = 0, nc = elem->n_children();
1073 if (!elem->is_child_on_edge(c, e))
1076 fine_fe->edge_reinit(elem->child_ptr(c), e);
1077 fine_points.insert(fine_points.end(),
1082 std::vector<Point> fine_qp;
1084 fine_points, fine_qp);
1086 context.elem_fe_reinit(&fine_qp);
1089 #endif // LIBMESH_ENABLE_AMR 1090 context.edge_fe_reinit();
1092 const unsigned int n_qp =
1093 cast_int<unsigned int>(xyz_values.size());
1098 const unsigned int n_side_dofs =
1099 cast_int<unsigned int>(side_dofs.size());
1103 unsigned int free_dofs = 0;
1104 for (
unsigned int i=0; i != n_side_dofs; ++i)
1105 if (!dof_is_fixed[side_dofs[i]])
1106 free_dof[free_dofs++] = i;
1112 Ke.resize (free_dofs, free_dofs); Ke.zero();
1113 Fe.resize (free_dofs); Fe.zero();
1115 DenseVector<FValue> Uedge(free_dofs);
1118 for (
unsigned int qp=0; qp<n_qp; qp++)
1121 FValue fineval = f.eval_at_point(context,
1126 VectorValue<FValue> finegrad;
1128 finegrad = g->eval_at_point(context,
1134 for (
unsigned int sidei=0, freei=0;
1135 sidei != n_side_dofs; ++sidei)
1137 unsigned int i = side_dofs[sidei];
1139 if (dof_is_fixed[i])
1141 for (
unsigned int sidej=0, freej=0;
1142 sidej != n_side_dofs; ++sidej)
1144 unsigned int j = side_dofs[sidej];
1145 if (dof_is_fixed[j])
1146 Fe(freei) -= phi[i][qp] * phi[j][qp] *
1149 Ke(freei,freej) += phi[i][qp] *
1150 phi[j][qp] * JxW[qp];
1153 if (dof_is_fixed[j])
1154 Fe(freei) -= ( (*dphi)[i][qp] *
1158 Ke(freei,freej) += ( (*dphi)[i][qp] *
1162 if (!dof_is_fixed[j])
1165 Fe(freei) += phi[i][qp] * fineval * JxW[qp];
1167 Fe(freei) += (finegrad * (*dphi)[i][qp] ) *
1173 Ke.cholesky_solve(Fe, Uedge);
1176 for (
unsigned int i=0; i != free_dofs; ++i)
1178 FValue & ui = Ue(side_dofs[free_dof[i]]);
1182 dof_is_fixed[side_dofs[free_dof[i]]] =
true;
1187 STOP_LOG (
"project_edges",
"GenericProjector");
1189 START_LOG (
"project_sides",
"GenericProjector");
1196 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1203 const std::vector<Point> & xyz_values =
1204 #ifdef LIBMESH_ENABLE_AMR 1207 #endif // LIBMESH_ENABLE_AMR 1209 const std::vector<Real> & JxW =
1210 #ifdef LIBMESH_ENABLE_AMR 1213 #endif // LIBMESH_ENABLE_AMR 1215 const std::vector<std::vector<Real>> & phi =
1216 #ifdef LIBMESH_ENABLE_AMR 1219 #endif // LIBMESH_ENABLE_AMR 1221 const std::vector<std::vector<RealGradient>> * dphi =
nullptr;
1224 #ifdef LIBMESH_ENABLE_AMR 1227 &(side_fe->get_dphi());
1229 for (
auto s : elem->side_index_range())
1234 unsigned int n_side_dofs =
1235 cast_int<unsigned int>(side_dofs.size());
1239 unsigned int free_dofs = 0;
1240 for (
unsigned int i=0; i != n_side_dofs; ++i)
1241 if (!dof_is_fixed[side_dofs[i]])
1242 free_dof[free_dofs++] = i;
1248 Ke.resize (free_dofs, free_dofs); Ke.zero();
1249 Fe.resize (free_dofs); Fe.zero();
1251 DenseVector<FValue> Uside(free_dofs);
1253 context.side = cast_int<unsigned char>(s);
1255 #ifdef LIBMESH_ENABLE_AMR 1258 std::vector<Point> fine_points;
1260 std::unique_ptr<FEBase> fine_fe
1263 std::unique_ptr<QBase> qrule
1264 (base_fe_type.default_quadrature_rule(dim-1));
1265 fine_fe->attach_quadrature_rule(qrule.get());
1267 const std::vector<Point> & child_xyz =
1270 for (
unsigned int c = 0, nc = elem->n_children();
1273 if (!elem->is_child_on_side(c, s))
1276 fine_fe->reinit(elem->child_ptr(c), s);
1277 fine_points.insert(fine_points.end(),
1282 std::vector<Point> fine_qp;
1284 fine_points, fine_qp);
1286 context.elem_fe_reinit(&fine_qp);
1289 #endif // LIBMESH_ENABLE_AMR 1290 context.side_fe_reinit();
1292 const unsigned int n_qp =
1293 cast_int<unsigned int>(xyz_values.size());
1296 for (
unsigned int qp=0; qp<n_qp; qp++)
1299 FValue fineval = f.eval_at_point(context,
1304 VectorValue<FValue> finegrad;
1306 finegrad = g->eval_at_point(context,
1312 for (
unsigned int sidei=0, freei=0;
1313 sidei != n_side_dofs; ++sidei)
1315 unsigned int i = side_dofs[sidei];
1317 if (dof_is_fixed[i])
1319 for (
unsigned int sidej=0, freej=0;
1320 sidej != n_side_dofs; ++sidej)
1322 unsigned int j = side_dofs[sidej];
1323 if (dof_is_fixed[j])
1324 Fe(freei) -= phi[i][qp] * phi[j][qp] *
1327 Ke(freei,freej) += phi[i][qp] *
1328 phi[j][qp] * JxW[qp];
1331 if (dof_is_fixed[j])
1332 Fe(freei) -= ( (*dphi)[i][qp] *
1336 Ke(freei,freej) += ( (*dphi)[i][qp] *
1340 if (!dof_is_fixed[j])
1343 Fe(freei) += (fineval * phi[i][qp]) * JxW[qp];
1345 Fe(freei) += (finegrad * (*dphi)[i][qp] ) *
1351 Ke.cholesky_solve(Fe, Uside);
1354 for (
unsigned int i=0; i != free_dofs; ++i)
1356 FValue & ui = Ue(side_dofs[free_dof[i]]);
1360 dof_is_fixed[side_dofs[free_dof[i]]] =
true;
1365 STOP_LOG (
"project_sides",
"GenericProjector");
1367 START_LOG (
"project_interior",
"GenericProjector");
1373 unsigned int free_dofs = 0;
1374 for (
unsigned int i=0; i != n_dofs; ++i)
1375 if (!dof_is_fixed[i])
1376 free_dof[free_dofs++] = i;
1382 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1387 const std::vector<Point> & xyz_values = fe->get_xyz();
1388 const std::vector<Real> & JxW = fe->get_JxW();
1390 const std::vector<std::vector<Real>> & phi = fe->get_phi();
1391 const std::vector<std::vector<RealGradient>> * dphi =
nullptr;
1393 dphi = &(fe->get_dphi());
1395 #ifdef LIBMESH_ENABLE_AMR 1398 std::vector<Point> fine_points;
1400 std::unique_ptr<FEBase> fine_fe
1403 std::unique_ptr<QBase> qrule
1404 (base_fe_type.default_quadrature_rule(dim));
1405 fine_fe->attach_quadrature_rule(qrule.get());
1407 const std::vector<Point> & child_xyz =
1410 for (
auto & child : elem->child_ref_range())
1412 fine_fe->reinit(&child);
1413 fine_points.insert(fine_points.end(),
1418 std::vector<Point> fine_qp;
1420 fine_points, fine_qp);
1422 context.elem_fe_reinit(&fine_qp);
1425 #endif // LIBMESH_ENABLE_AMR 1426 context.elem_fe_reinit();
1428 const unsigned int n_qp =
1429 cast_int<unsigned int>(xyz_values.size());
1431 Ke.resize (free_dofs, free_dofs); Ke.zero();
1432 Fe.resize (free_dofs); Fe.zero();
1434 DenseVector<FValue> Uint(free_dofs);
1437 for (
unsigned int qp=0; qp<n_qp; qp++)
1440 FValue fineval = f.eval_at_point(context,
1445 VectorValue<FValue> finegrad;
1447 finegrad = g->eval_at_point(context,
1453 for (
unsigned int i=0, freei=0; i != n_dofs; ++i)
1456 if (dof_is_fixed[i])
1458 for (
unsigned int j=0, freej=0; j != n_dofs; ++j)
1460 if (dof_is_fixed[j])
1461 Fe(freei) -= phi[i][qp] * phi[j][qp] * JxW[qp]
1464 Ke(freei,freej) += phi[i][qp] * phi[j][qp] *
1468 if (dof_is_fixed[j])
1469 Fe(freei) -= ( (*dphi)[i][qp] *
1470 (*dphi)[j][qp] ) * JxW[qp] *
1473 Ke(freei,freej) += ( (*dphi)[i][qp] *
1477 if (!dof_is_fixed[j])
1480 Fe(freei) += phi[i][qp] * fineval * JxW[qp];
1482 Fe(freei) += (finegrad * (*dphi)[i][qp] ) * JxW[qp];
1486 Ke.cholesky_solve(Fe, Uint);
1489 for (
unsigned int i=0; i != free_dofs; ++i)
1491 FValue & ui = Ue(free_dof[i]);
1495 dof_is_fixed[free_dof[i]] =
true;
1500 STOP_LOG (
"project_interior",
"GenericProjector");
1503 for (
unsigned int i=0; i != n_dofs; ++i)
1504 libmesh_assert(dof_is_fixed[i]);
1506 action.insert(context, var, Ue);
static void dofs_on_edge(const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int e, std::vector< unsigned int > &di)
unsigned int variable_scalar_number(const std::string &var, unsigned int component) const
static void dofs_on_side(const Elem *const elem, const unsigned int dim, const FEType &fe_t, unsigned int s, std::vector< unsigned int > &di)
static const Real TOLERANCE
const GFunctor * master_g
const dof_id_type n_nodes
const FFunctor & master_f
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
static Point inverse_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
const ProjectionAction & master_action
FEGenericBase< Real > FEBase
const std::vector< unsigned int > & variables
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const DofMap & get_dof_map() const