20 #ifndef GENERIC_PROJECTOR_H 
   21 #define GENERIC_PROJECTOR_H 
   27 #include "libmesh/dense_matrix.h" 
   28 #include "libmesh/dense_vector.h" 
   29 #include "libmesh/dof_map.h" 
   30 #include "libmesh/elem.h" 
   31 #include "libmesh/fe_base.h" 
   32 #include "libmesh/fe_interface.h" 
   33 #include "libmesh/int_range.h" 
   34 #include "libmesh/libmesh_logging.h" 
   35 #include "libmesh/mesh_tools.h" 
   36 #include "libmesh/numeric_vector.h" 
   37 #include "libmesh/quadrature.h" 
   38 #include "libmesh/system.h" 
   39 #include "libmesh/threads.h" 
   40 #include "libmesh/auto_ptr.h"  
   64 template <
typename SendT, 
typename T>
 
   66 { converted = received; }
 
   78 template <
typename FFunctor, 
typename GFunctor,
 
   79           typename FValue, 
typename ProjectionAction>
 
   92   std::unordered_map<dof_id_type, std::vector<dof_id_type>> * 
nodes_to_elem;
 
   97                     const FFunctor & f_in,
 
   98                     const GFunctor * g_in,
 
   99                     const ProjectionAction & act_in,
 
  100                     const std::vector<unsigned int> & variables_in,
 
  101                     std::unordered_map<
dof_id_type, std::vector<dof_id_type>> *
 
  102                       nodes_to_elem_in = 
nullptr) :
 
  115           std::unordered_map<dof_id_type, std::vector<dof_id_type>>;
 
  165     std::unordered_map<dof_id_type, std::pair<FValue, processor_id_type>> 
new_ids_to_push;
 
  174                             unsigned short node_num,
 
  183     void insert_ids(
const std::vector<dof_id_type> & ids,
 
  184                     const std::vector<FValue> & vals,
 
  207   typedef std::pair<const Node *, std::tuple<const Elem *, unsigned short, var_set>>
 
  227     std::unique_ptr<GFunctor> 
g;
 
  230       (
const std::vector<dof_id_type> & dof_indices_var,
 
  231        const std::vector<unsigned int> & involved_dofs,
 
  232        unsigned int var_component,
 
  260     typedef std::unordered_multimap<const Node *, std::tuple<const Elem *, unsigned short, var_set>> 
ves_multimap;
 
  349     (std::unordered_map<
dof_id_type, std::pair<FValue, processor_id_type>> & ids_to_push,
 
  350      ProjectionAction & action) 
const;
 
  363 template <
typename Val>
 
  384   void insert(
const std::vector<dof_id_type> & dof_indices,
 
  391     unsigned int size = Ue.
size();
 
  393     libmesh_assert_equal_to(size, dof_indices.size());
 
  399       for (
unsigned int i = 0; i != size; ++i)
 
  400         if ((dof_indices[i] >= first) && (dof_indices[i] <  last))
 
  416 template <
typename Output>
 
  423     _f(fw.
_f->clone()) {}
 
  433   { 
return _f->component(c, i, n, time); }
 
  439   { 
return _f->component(c, i, n, time); }
 
  446                       std::vector<dof_id_type> &,
 
  447                       std::vector<Output> &)
 
  454                       std::vector<dof_id_type> &,
 
  455                       std::vector<Output> &)
 
  459   std::unique_ptr<FEMFunctionBase<Output>> 
_f;
 
  473 #ifdef LIBMESH_ENABLE_AMR 
  474 template <
typename Output,
 
  491         const unsigned int vcomp = sys.variable_scalar_number(v,0);
 
  492         if (vcomp >= component_to_var.size())
 
  493           component_to_var.resize(vcomp+1, static_cast<unsigned int>(-1));
 
  494         component_to_var[vcomp] = v;
 
  502     component_to_var(in.component_to_var)
 
  506   static void get_shape_outputs(
FEBase & fe);
 
  518         const std::set<unsigned char> & elem_dims =
 
  519           old_context.elem_dimensions();
 
  521         for (
const auto & 
dim : elem_dims)
 
  523             old_context.get_element_fe(var, fe, 
dim);
 
  524             get_shape_outputs(*fe);
 
  534     LOG_SCOPE (
"check_old_context(c)", 
"OldSolutionBase");
 
  536     if (last_elem != &elem)
 
  540             old_context.pre_fe_reinit(sys, elem.
parent());
 
  553             old_context.pre_fe_reinit(sys, &elem);
 
  567     LOG_SCOPE (
"check_old_context(c,p)", 
"OldSolutionBase");
 
  569     if (last_elem != &elem)
 
  573             old_context.pre_fe_reinit(sys, elem.
parent());
 
  584             const Real master_tol = out_of_elem_tol / elem.
hmax() * 2;
 
  587               if (child.close_to_point(p, master_tol))
 
  589                   old_context.pre_fe_reinit(sys, &child);
 
  594               (old_context.get_elem().close_to_point(p, master_tol));
 
  601             old_context.pre_fe_reinit(sys, &elem);
 
  610         const Real master_tol = out_of_elem_tol / elem.
hmax() * 2;
 
  612         if (!old_context.get_elem().close_to_point(p, master_tol))
 
  614             libmesh_assert_equal_to
 
  618               if (child.close_to_point(p, master_tol))
 
  620                   old_context.pre_fe_reinit(sys, &child);
 
  625               (old_context.get_elem().close_to_point(p, master_tol));
 
  650 template <
typename Output,
 
  661     old_solution(old_sol)
 
  664     this->old_context.set_custom_solution(&old_solution);
 
  669     old_solution(in.old_solution)
 
  672     this->old_context.set_custom_solution(&old_solution);
 
  678                        unsigned int elem_dim,
 
  688     LOG_SCOPE (
"eval_at_point()", 
"OldSolutionValue");
 
  690     if (!this->check_old_context(c, p))
 
  694     libmesh_assert_less(i, this->component_to_var.size());
 
  695     unsigned int var = this->component_to_var[i];
 
  698     (this->old_context.*point_output)(var, p, n, this->out_of_elem_tol);
 
  703                       unsigned int node_num,
 
  704                       unsigned int var_num,
 
  705                       std::vector<dof_id_type> & indices,
 
  706                       std::vector<Output> & values)
 
  708     LOG_SCOPE (
"eval_old_dofs(node)", 
"OldSolutionValue");
 
  710     this->sys.get_dof_map().dof_indices(elem, node_num, indices, var_num);
 
  712     std::vector<dof_id_type> old_indices;
 
  714     this->sys.get_dof_map().old_dof_indices(elem, node_num, old_indices, var_num);
 
  716     libmesh_assert_equal_to (old_indices.size(), indices.size());
 
  718     values.resize(old_indices.size());
 
  720     old_solution.get(old_indices, values);
 
  726                       unsigned int sys_num,
 
  727                       unsigned int var_num,
 
  728                       std::vector<dof_id_type> & indices,
 
  729                       std::vector<Output> & values)
 
  731     LOG_SCOPE (
"eval_old_dofs(elem)", 
"OldSolutionValue");
 
  735     const Elem & old_elem =
 
  744     std::vector<dof_id_type> old_dof_indices(nc);
 
  754         const std::pair<unsigned int, unsigned int>
 
  756         const unsigned int vg = vg_and_offset.first;
 
  757         const unsigned int vig = vg_and_offset.second;
 
  759         const unsigned int n_comp = elem.
n_comp_group(sys_num,vg);
 
  760         libmesh_assert_greater(elem.
n_systems(), sys_num);
 
  761         libmesh_assert_greater_equal(n_comp, nc);
 
  763         for (
unsigned int i=0; i<nc; i++)
 
  772             old_dof_indices[i] = d_old;
 
  779     old_solution.get(old_dof_indices, values);
 
  803 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 
  818 #endif // LIBMESH_USE_COMPLEX_NUMBERS 
  829              bool extra_hanging_dofs,
 
  832   LOG_SCOPE (
"Number eval_at_node()", 
"OldSolutionValue");
 
  839   libmesh_assert_less(i, this->component_to_var.size());
 
  840   unsigned int var = this->component_to_var[i];
 
  855       (!extra_hanging_dofs ||
 
  863       return old_solution(old_id);
 
  866   return this->eval_at_point(c, i, n, 0);
 
  877              unsigned int elem_dim,
 
  879              bool extra_hanging_dofs,
 
  882   LOG_SCOPE (
"Gradient eval_at_node()", 
"OldSolutionValue");
 
  889   libmesh_assert_less(i, this->component_to_var.size());
 
  890   unsigned int var = this->component_to_var[i];
 
  905       (!extra_hanging_dofs ||
 
  912       for (
unsigned int d = 0; d != elem_dim; ++d)
 
  916           g(d) = old_solution(old_id);
 
  921   return this->eval_at_point(c, i, n, 0);
 
  934 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 
  940 #endif // LIBMESH_USE_COMPLEX_NUMBERS 
  942 #endif // LIBMESH_ENABLE_AMR 
  948 template <
typename FFunctor, 
typename GFunctor,
 
  949           typename FValue, 
typename ProjectionAction>
 
  953   LOG_SCOPE (
"project", 
"GenericProjector");
 
  957   done_saving_ids = 
false;
 
  961   ProjectionAction action(master_action);
 
  964   std::unordered_map<dof_id_type, std::pair<FValue, processor_id_type>> ids_to_push;
 
  971   std::vector<node_projection> vertices(sort_work.
vertices.begin(),
 
  974   done_saving_ids = sort_work.
edges.empty() &&
 
  976   system.comm().max(done_saving_ids);
 
  987   done_saving_ids = sort_work.
sides.empty() && sort_work.
interiors.empty();
 
  988   system.comm().max(done_saving_ids);
 
  990   this->send_and_insert_dof_values(ids_to_push, action);
 
  993     std::vector<node_projection> edges(sort_work.
edges.begin(), sort_work.
edges.end());
 
  996     ids_to_push.insert(project_edges.new_ids_to_push.begin(),
 
  997                        project_edges.new_ids_to_push.end());
 
  998     ids_to_save.insert(project_edges.new_ids_to_save.begin(),
 
  999                        project_edges.new_ids_to_save.end());
 
 1002   done_saving_ids = sort_work.
interiors.empty();
 
 1003   system.comm().max(done_saving_ids);
 
 1005   this->send_and_insert_dof_values(ids_to_push, action);
 
 1008     std::vector<node_projection> sides(sort_work.
sides.begin(), sort_work.
sides.end());
 
 1011     ids_to_push.insert(project_sides.new_ids_to_push.begin(),
 
 1012                        project_sides.new_ids_to_push.end());
 
 1013     ids_to_save.insert(project_sides.new_ids_to_save.begin(),
 
 1014                        project_sides.new_ids_to_save.end());
 
 1017   done_saving_ids = 
true;
 
 1018   this->send_and_insert_dof_values(ids_to_push, action);
 
 1028 template <
typename FFunctor, 
typename GFunctor,
 
 1029           typename FValue, 
typename ProjectionAction>
 
 1037   for (
const auto & var : this->projector.variables)
 
 1041       FEBase * side_fe = 
nullptr;
 
 1042       FEBase * edge_fe = 
nullptr;
 
 1044       const std::set<unsigned char> & elem_dims =
 
 1045         context.elem_dimensions();
 
 1047       for (
const auto & 
dim : elem_dims)
 
 1049           context.get_element_fe( var, fe, 
dim );
 
 1053             context.get_side_fe( var, side_fe, 
dim );
 
 1055             context.get_edge_fe( var, edge_fe );
 
 1076           this->conts[var] = cont;
 
 1089   f.init_context(context);
 
 1093 template <
typename FFunctor, 
typename GFunctor,
 
 1094           typename FValue, 
typename ProjectionAction>
 
 1100     g = libmesh_make_unique<GFunctor>(*p.
master_g);
 
 1104   for (
const auto & var : this->projector.variables)
 
 1105     if (this->conts[var] == 
C_ONE)
 
 1110     g->init_context(context);
 
 1114 template <
typename FFunctor, 
typename GFunctor,
 
 1115           typename FValue, 
typename ProjectionAction>
 
 1119   auto iter = new_ids_to_push.find(
id);
 
 1120   if (iter == new_ids_to_push.end())
 
 1121     action.insert(
id, val);
 
 1125       iter->second = std::make_pair(val, pid);
 
 1127   if (!this->projector.done_saving_ids)
 
 1130       new_ids_to_save[id] = val;
 
 1135 template <
typename FFunctor, 
typename GFunctor,
 
 1136           typename FValue, 
typename ProjectionAction>
 
 1140   libmesh_assert_equal_to(ids.size(), vals.size());
 
 1144       const FValue & val = vals[i];
 
 1146       auto iter = new_ids_to_push.find(
id);
 
 1147       if (iter == new_ids_to_push.end())
 
 1148         action.insert(
id, val);
 
 1152           iter->second = std::make_pair(val, pid);
 
 1154       if (!this->projector.done_saving_ids)
 
 1157           new_ids_to_save[id] = val;
 
 1163 template <
typename FFunctor, 
typename GFunctor,
 
 1164           typename FValue, 
typename ProjectionAction>
 
 1173 template <
typename FFunctor, 
typename GFunctor,
 
 1174           typename FValue, 
typename ProjectionAction>
 
 1202   std::unordered_map<const Node *, std::pair<var_set, var_set>> copied_nodes;
 
 1204   const unsigned int sys_num = system.number();
 
 1209   std::vector<unsigned short> extra_hanging_dofs;
 
 1210   bool all_extra_hanging_dofs = 
true;
 
 1211   for (
auto v_num : this->projector.variables)
 
 1213       if (extra_hanging_dofs.size() <= v_num)
 
 1214         extra_hanging_dofs.resize(v_num+1, 
false);
 
 1215       extra_hanging_dofs[v_num] =
 
 1218       if (!extra_hanging_dofs[v_num])
 
 1219         all_extra_hanging_dofs = 
false;
 
 1222   for (
const auto & elem : range)
 
 1226       bool copy_this_elem = 
false;
 
 1228 #ifdef LIBMESH_ENABLE_AMR 
 1230       if (f.is_grid_projection())
 
 1236           if (!elem->old_dof_object &&
 
 1246             copy_this_elem = 
true;
 
 1249               bool reinitted = 
false;
 
 1253               for (
auto v_num : this->projector.variables)
 
 1255                   const Variable & var = system.variable(v_num);
 
 1262                       elem->p_level() == 0 &&
 
 1269                           context.pre_fe_reinit(system, elem);
 
 1272                       std::vector<FValue> Ue(1);
 
 1273                       std::vector<dof_id_type> elem_dof_ids(1);
 
 1275                       f.eval_old_dofs(*elem, fe_type, sys_num, v_num,
 
 1278                       action.insert(elem_dof_ids[0], Ue[0]);
 
 1283 #endif // LIBMESH_ENABLE_AMR 
 1285       const int dim = elem->dim();
 
 1287       const unsigned int n_vertices = elem->n_vertices();
 
 1288       const unsigned int n_edges = elem->n_edges();
 
 1289       const unsigned int n_nodes = elem->n_nodes();
 
 1292       const unsigned int n_sides = (
dim > 1) * elem->n_sides();
 
 1295       var_set vertex_vars, edge_vars, side_vars;
 
 1299       const bool has_edge_nodes = (
n_nodes > n_vertices && 
dim > 2);
 
 1302       const bool has_side_nodes =
 
 1303         (
n_nodes > n_vertices + ((
dim > 2) * n_edges));
 
 1307       const bool has_interior_nodes =
 
 1308         (
n_nodes > n_vertices + ((
dim > 2) * n_edges) + n_sides);
 
 1310       for (
auto v_num : this->projector.variables)
 
 1318           const ElemType elem_type = elem->type();
 
 1321             vertex_vars.insert(vertex_vars.end(), v_num);
 
 1327               edge_vars.insert(edge_vars.end(), v_num);
 
 1334                     side_vars.insert(side_vars.end(), v_num);
 
 1339                 for (
unsigned int n = 0; n != 
n_nodes; ++n)
 
 1340                   if (elem->is_face(n))
 
 1344                         side_vars.insert(side_vars.end(), v_num);
 
 1350               (has_interior_nodes &&
 
 1353 #ifdef LIBMESH_ENABLE_AMR 
 1356               if ((f.is_grid_projection() &&
 
 1359                    elem->p_level() == 0 &&
 
 1365 #endif // LIBMESH_ENABLE_AMR 
 1368               if (interiors.empty() || interiors.back() != elem)
 
 1369                 interiors.push_back(elem);
 
 1376       auto erase_covered_vars = []
 
 1379           auto covered_range = covered.equal_range(node);
 
 1380           for (
const auto & v_ent : 
as_range(covered_range))
 
 1381             for (
const unsigned int var_covered :
 
 1382                  std::get<2>(v_ent.second))
 
 1383               remaining.erase(var_covered);
 
 1386       auto erase_nonhanging_vars = [&extra_hanging_dofs]
 
 1389           auto covered_range = covered.equal_range(node);
 
 1390           for (
const auto & v_ent : 
as_range(covered_range))
 
 1391             for (
const unsigned int var_covered :
 
 1392                  std::get<2>(v_ent.second))
 
 1393               if (!extra_hanging_dofs[var_covered])
 
 1394                 remaining.erase(var_covered);
 
 1397       auto erase_copied_vars = [&copied_nodes, &extra_hanging_dofs]
 
 1398         (
const Node * node, 
bool is_vertex, 
var_set & remaining)
 
 1400           auto copying_range = copied_nodes.equal_range(node);
 
 1401           for (
const auto & v_ent : 
as_range(copying_range))
 
 1403               for (
const unsigned int var_covered :
 
 1405                 if (is_vertex || !extra_hanging_dofs[var_covered])
 
 1406                   remaining.erase(var_covered);
 
 1408               for (
const unsigned int var_covered :
 
 1409                    v_ent.second.second)
 
 1410                 if (!is_vertex || !extra_hanging_dofs[var_covered])
 
 1411                   remaining.erase(var_covered);
 
 1415       for (
unsigned int v=0; v != n_vertices; ++v)
 
 1417           const Node * node = elem->node_ptr(v);
 
 1419           auto remaining_vars = vertex_vars;
 
 1421           erase_covered_vars(node, remaining_vars, vertices);
 
 1423           if (remaining_vars.empty())
 
 1426           if (!all_extra_hanging_dofs)
 
 1428               erase_nonhanging_vars(node, remaining_vars, edges);
 
 1429               if (remaining_vars.empty())
 
 1432               erase_nonhanging_vars(node, remaining_vars, sides);
 
 1433               if (remaining_vars.empty())
 
 1437           erase_copied_vars(node, 
true, remaining_vars);
 
 1438           if (remaining_vars.empty())
 
 1443               for (
auto var : remaining_vars)
 
 1445                   std::vector<dof_id_type> node_dof_ids;
 
 1446                   std::vector<FValue> values;
 
 1448                   f.eval_old_dofs(*elem, v, var, node_dof_ids, values);
 
 1450                   insert_ids(node_dof_ids, values, node->
processor_id());
 
 1452               copied_nodes[node].first.insert(remaining_vars.begin(),
 
 1453                                               remaining_vars.end());
 
 1454               this->find_dofs_to_send(*node, *elem, v, remaining_vars);
 
 1458               (node, std::make_tuple(elem, v, std::move(remaining_vars)));
 
 1463           for (
unsigned int e=0; e != n_edges; ++e)
 
 1465               const Node * node = elem->node_ptr(n_vertices+e);
 
 1467               auto remaining_vars = edge_vars;
 
 1469               erase_covered_vars(node, remaining_vars, edges);
 
 1470               if (remaining_vars.empty())
 
 1473               erase_covered_vars(node, remaining_vars, sides);
 
 1474               if (remaining_vars.empty())
 
 1477               if (!all_extra_hanging_dofs)
 
 1479                   erase_nonhanging_vars(node, remaining_vars, vertices);
 
 1480                   if (remaining_vars.empty())
 
 1484               erase_copied_vars(node, 
false, remaining_vars);
 
 1485               if (remaining_vars.empty())
 
 1490                   for (
auto var : remaining_vars)
 
 1492                       std::vector<dof_id_type> edge_dof_ids;
 
 1493                       std::vector<FValue> values;
 
 1495                       f.eval_old_dofs(*elem, n_vertices+e, var, edge_dof_ids, values);
 
 1497                       insert_ids(edge_dof_ids, values, node->
processor_id());
 
 1499                   copied_nodes[node].second.insert(remaining_vars.begin(),
 
 1500                                                    remaining_vars.end());
 
 1501                   this->find_dofs_to_send(*node, *elem, n_vertices+e, remaining_vars);
 
 1505                   (node, std::make_tuple(elem, e, std::move(remaining_vars)));
 
 1511           for (
unsigned int side=0; side != n_sides; ++side)
 
 1513               const Node * node = 
nullptr;
 
 1514               unsigned short node_num = n_vertices+(
dim>2)*n_edges+side;
 
 1516                 node = elem->node_ptr(node_num);
 
 1520                   for (
unsigned int n = 0; n != 
n_nodes; ++n)
 
 1522                       if (!elem->is_face(n))
 
 1525                       if (elem->is_node_on_side(n, side))
 
 1528                           node = elem->node_ptr(node_num);
 
 1537               auto remaining_vars = side_vars;
 
 1539               erase_covered_vars(node, remaining_vars, edges);
 
 1540               if (remaining_vars.empty())
 
 1543               erase_covered_vars(node, remaining_vars, sides);
 
 1544               if (remaining_vars.empty())
 
 1547               if (!all_extra_hanging_dofs)
 
 1549                   erase_nonhanging_vars(node, remaining_vars, vertices);
 
 1550                   if (remaining_vars.empty())
 
 1554               erase_copied_vars(node, 
false, remaining_vars);
 
 1555               if (remaining_vars.empty())
 
 1560                   for (
auto var : remaining_vars)
 
 1562                       std::vector<dof_id_type> side_dof_ids;
 
 1563                       std::vector<FValue> values;
 
 1565                       f.eval_old_dofs(*elem, node_num, var, side_dof_ids, values);
 
 1567                       insert_ids(side_dof_ids, values, node->
processor_id());
 
 1569                   copied_nodes[node].second.insert(remaining_vars.begin(),
 
 1570                                                    remaining_vars.end());
 
 1571                   this->find_dofs_to_send(*node, *elem, node_num, remaining_vars);
 
 1575                   (node, std::make_tuple(elem, side, std::move(remaining_vars)));
 
 1582           for (
auto v_num : this->projector.variables)
 
 1584               const Variable & var = system.variable(v_num);
 
 1589               std::vector<FValue> Ue;
 
 1590               std::vector<dof_id_type> elem_dof_ids;
 
 1591               f.eval_old_dofs(*elem, fe_type, sys_num, v_num,
 
 1593               action.insert(elem_dof_ids, Ue);
 
 1595               if (has_interior_nodes)
 
 1597                   std::vector<FValue> Un;
 
 1598                   std::vector<dof_id_type> node_dof_ids;
 
 1600                   f.eval_old_dofs(*elem, 
n_nodes-1, v_num, node_dof_ids, Un);
 
 1601                   action.insert(node_dof_ids, Un);
 
 1609 template <
typename FFunctor, 
typename GFunctor,
 
 1610           typename FValue, 
typename ProjectionAction>
 
 1616       for (
const auto & pair : other_mm)
 
 1618           const Node * node = pair.first;
 
 1619           auto remaining_vars = std::get<2>(pair.second);
 
 1621           auto my_range = 
self.equal_range(node);
 
 1622           for (
const auto & v_ent : 
as_range(my_range))
 
 1623             for (
const unsigned int var_covered :
 
 1624                  std::get<2>(v_ent.second))
 
 1625               remaining_vars.erase(var_covered);
 
 1627           if (!remaining_vars.empty())
 
 1629               (node, std::make_tuple(std::get<0>(pair.second),
 
 1630                                      std::get<1>(pair.second),
 
 1631                                      std::move(remaining_vars)));
 
 1635   merge_multimaps(vertices, other.
vertices);
 
 1636   merge_multimaps(edges, other.
edges);
 
 1637   merge_multimaps(sides, other.
sides);
 
 1639   std::sort(interiors.begin(), interiors.end());
 
 1640   std::vector<const Elem *> other_interiors = other.
interiors;
 
 1641   std::sort(other_interiors.begin(), other_interiors.end());
 
 1642   std::vector<const Elem *> merged_interiors;
 
 1643   std::set_union(interiors.begin(), interiors.end(),
 
 1644                  other_interiors.begin(), other_interiors.end(),
 
 1645                  std::back_inserter(merged_interiors));
 
 1646   interiors.swap(merged_interiors);
 
 1648   SubFunctor::join(other);
 
 1652 template <
typename FFunctor, 
typename GFunctor,
 
 1653           typename FValue, 
typename ProjectionAction>
 
 1657   LOG_SCOPE (
"project_vertices",
"GenericProjector");
 
 1659   const unsigned int sys_num = system.number();
 
 1663   std::vector<unsigned short> extra_hanging_dofs;
 
 1664   for (
auto v_num : this->projector.variables)
 
 1666       if (extra_hanging_dofs.size() <= v_num)
 
 1667         extra_hanging_dofs.resize(v_num+1, 
false);
 
 1668       extra_hanging_dofs[v_num] =
 
 1672   for (
const auto & v_pair : range)
 
 1674       const Node & vertex = *v_pair.first;
 
 1675       const Elem & elem = *std::get<0>(v_pair.second);
 
 1676       const unsigned int n = std::get<1>(v_pair.second);
 
 1677       const var_set & vertex_vars = std::get<2>(v_pair.second);
 
 1679       context.pre_fe_reinit(system, &elem);
 
 1681       this->find_dofs_to_send(vertex, elem, n, vertex_vars);
 
 1685       for (
const auto & var : vertex_vars)
 
 1687           const Variable & variable = system.variable(var);
 
 1688           const FEType & base_fe_type = variable.
type();
 
 1689           const unsigned int var_component =
 
 1690             system.variable_scalar_number(var, 0);
 
 1698               libmesh_assert_equal_to(vertex.
n_comp(sys_num, var), 0);
 
 1705               const FValue val = f.eval_at_node
 
 1706                 (context, var_component,  0, 
 
 1707                  vertex, extra_hanging_dofs[var], system.time);
 
 1710           else if (cont == 
C_ONE)
 
 1720               const int dim = elem.
dim();
 
 1727               for (
const auto e_id : (*this->projector.nodes_to_elem)[vertex.
id()])
 
 1729                   const Elem & e = system.get_mesh().elem_ref(e_id);
 
 1730                   libmesh_assert_equal_to(
dim, e.
dim());
 
 1733 #ifdef LIBMESH_ENABLE_AMR 
 1734               bool is_parent_vertex = 
false;
 
 1737                   const int i_am_child =
 
 1743               const bool is_parent_vertex = 
false;
 
 1750                     f.eval_at_node(context,
 
 1754                                    extra_hanging_dofs[var],
 
 1760                     g->eval_at_node(context,
 
 1764                                     extra_hanging_dofs[var],
 
 1766                     g->eval_at_point(context,
 
 1771                   insert_id(first_id+1, grad(0),
 
 1780                             nxplus = elem.
point(n);
 
 1781                       nxminus(0) -= delta_x;
 
 1782                       nxplus(0) += delta_x;
 
 1784                         g->eval_at_point(context,
 
 1789                         g->eval_at_point(context,
 
 1794                       insert_id(first_id+2, grad(1),
 
 1797                       insert_id(first_id+3,
 
 1798                         (gxplus(1) - gxminus(1)) / 2. / delta_x,
 
 1805                           insert_id(first_id+4, grad(2),
 
 1808                           insert_id(first_id+5,
 
 1809                             (gxplus(2) - gxminus(2)) / 2. / delta_x,
 
 1814                             nyplus = elem.
point(n);
 
 1815                           nyminus(1) -= delta_x;
 
 1816                           nyplus(1) += delta_x;
 
 1818                             g->eval_at_point(context,
 
 1823                             g->eval_at_point(context,
 
 1828                           insert_id(first_id+6,
 
 1829                             (gyplus(2) - gyminus(2)) / 2. / delta_x,
 
 1833                             nxmyp = elem.
point(n),
 
 1834                             nxpym = elem.
point(n),
 
 1835                             nxpyp = elem.
point(n);
 
 1836                           nxmym(0) -= delta_x;
 
 1837                           nxmym(1) -= delta_x;
 
 1838                           nxmyp(0) -= delta_x;
 
 1839                           nxmyp(1) += delta_x;
 
 1840                           nxpym(0) += delta_x;
 
 1841                           nxpym(1) -= delta_x;
 
 1842                           nxpyp(0) += delta_x;
 
 1843                           nxpyp(1) += delta_x;
 
 1845                             g->eval_at_point(context,
 
 1850                             g->eval_at_point(context,
 
 1855                             g->eval_at_point(context,
 
 1860                             g->eval_at_point(context,
 
 1864                           FValue gxzplus = (gxpyp(2) - gxmyp(2))
 
 1866                           FValue gxzminus = (gxpym(2) - gxmym(2))
 
 1869                           insert_id(first_id+7,
 
 1870                             (gxzplus - gxzminus) / 2. / delta_x,
 
 1873 #endif // LIBMESH_DIM > 2 
 1875 #endif // LIBMESH_DIM > 1 
 1882                   libmesh_assert_equal_to
 
 1884                       (
dim, base_fe_type, elem.
type(),
 
 1886                     (
unsigned int)(1 + 
dim));
 
 1888                     f.eval_at_node(context, var_component, 
dim,
 
 1889                                    vertex, extra_hanging_dofs[var],
 
 1894                     g->eval_at_node(context, var_component, 
dim,
 
 1895                                     vertex, extra_hanging_dofs[var],
 
 1897                     g->eval_at_point(context, var_component, vertex,
 
 1899                   for (
int i=0; i!= 
dim; ++i)
 
 1900                     insert_id(first_id + i + 1, grad(i),
 
 1905             libmesh_error_msg(
"Unknown continuity " << cont);
 
 1911 template <
typename FFunctor, 
typename GFunctor,
 
 1912           typename FValue, 
typename ProjectionAction>
 
 1916   LOG_SCOPE (
"project_edges",
"GenericProjector");
 
 1918   const unsigned int sys_num = system.number();
 
 1920   for (
const auto & e_pair : range)
 
 1922       const Elem & elem = *std::get<0>(e_pair.second);
 
 1926 #ifdef LIBMESH_ENABLE_AMR 
 1927       if (f.is_grid_projection() &&
 
 1933 #endif // LIBMESH_ENABLE_AMR 
 1935       const Node & edge_node = *e_pair.first;
 
 1936       const int dim = elem.
dim();
 
 1937       const var_set & edge_vars = std::get<2>(e_pair.second);
 
 1939       const unsigned short edge_num = std::get<1>(e_pair.second);
 
 1940       const unsigned short node_num = elem.
n_vertices() + edge_num;
 
 1941       context.edge = edge_num;
 
 1942       context.pre_fe_reinit(system, &elem);
 
 1944       this->find_dofs_to_send(edge_node, elem, node_num, edge_vars);
 
 1948       for (
const auto & var : edge_vars)
 
 1950           const Variable & variable = system.variable(var);
 
 1951           const FEType & base_fe_type = variable.
type();
 
 1952           const unsigned int var_component =
 
 1953             system.variable_scalar_number(var, 0);
 
 1958           FEType fe_type = base_fe_type;
 
 1969               if (fe_type.
order > 1)
 
 1975                   FValue fval = f.eval_at_point
 
 1976                     (context, var_component, edge_node, system.time);
 
 1977                   insert_id(dof_id, fval, pid);
 
 1982 #ifdef LIBMESH_ENABLE_AMR 
 1990 #endif // LIBMESH_ENABLE_AMR 
 1994           context.get_element_fe( var, fe, 
dim );
 
 1995           FEBase * edge_fe = 
nullptr;
 
 1996           context.get_edge_fe( var, edge_fe );
 
 2001 #ifdef LIBMESH_ENABLE_AMR 
 2007 #ifdef LIBMESH_ENABLE_AMR 
 2010               std::vector<Point> fine_points;
 
 2012               std::unique_ptr<FEBase> fine_fe
 
 2015               std::unique_ptr<QBase> qrule
 
 2017               fine_fe->attach_quadrature_rule(qrule.get());
 
 2019               const std::vector<Point> & child_xyz =
 
 2022               for (
unsigned int c = 0, nc = elem.
n_children();
 
 2028                   fine_fe->edge_reinit(elem.
child_ptr(c), context.edge);
 
 2029                   fine_points.insert(fine_points.end(),
 
 2034               std::vector<Point> fine_qp;
 
 2037               context.elem_fe_reinit(&fine_qp);
 
 2040 #endif // LIBMESH_ENABLE_AMR 
 2041             context.edge_fe_reinit();
 
 2043           const std::vector<dof_id_type> & dof_indices =
 
 2044             context.get_dof_indices(var);
 
 2046           std::vector<unsigned int> edge_dofs;
 
 2048                                     context.edge, edge_dofs);
 
 2050           this->construct_projection
 
 2051             (dof_indices, edge_dofs, var_component,
 
 2052              &edge_node, proj_fe);
 
 2058 template <
typename FFunctor, 
typename GFunctor,
 
 2059           typename FValue, 
typename ProjectionAction>
 
 2063   LOG_SCOPE (
"project_sides",
"GenericProjector");
 
 2065   const unsigned int sys_num = system.number();
 
 2067   for (
const auto & s_pair : range)
 
 2069       const Elem & elem = *std::get<0>(s_pair.second);
 
 2073 #ifdef LIBMESH_ENABLE_AMR 
 2074       if (f.is_grid_projection() &&
 
 2080 #endif // LIBMESH_ENABLE_AMR 
 2082       const Node & side_node = *s_pair.first;
 
 2083       const int dim = elem.
dim();
 
 2084       const var_set & side_vars = std::get<2>(s_pair.second);
 
 2086       const unsigned int side_num = std::get<1>(s_pair.second);
 
 2087       unsigned short node_num = elem.
n_vertices()+side_num;
 
 2102       context.side = side_num;
 
 2103       context.pre_fe_reinit(system, &elem);
 
 2105       this->find_dofs_to_send(side_node, elem, node_num, side_vars);
 
 2109       for (
const auto & var : side_vars)
 
 2111           const Variable & variable = system.variable(var);
 
 2112           const FEType & base_fe_type = variable.
type();
 
 2113           const unsigned int var_component =
 
 2114             system.variable_scalar_number(var, 0);
 
 2119           FEType fe_type = base_fe_type;
 
 2130               if (fe_type.
order > 1)
 
 2136                   FValue fval = f.eval_at_point
 
 2137                     (context, var_component, side_node, system.time);
 
 2138                   insert_id(dof_id, fval, pid);
 
 2143 #ifdef LIBMESH_ENABLE_AMR 
 2151 #endif // LIBMESH_ENABLE_AMR 
 2155           context.get_element_fe( var, fe, 
dim );
 
 2156           FEBase * side_fe = 
nullptr;
 
 2157           context.get_side_fe( var, side_fe );
 
 2162 #ifdef LIBMESH_ENABLE_AMR 
 2168 #ifdef LIBMESH_ENABLE_AMR 
 2171               std::vector<Point> fine_points;
 
 2173               std::unique_ptr<FEBase> fine_fe
 
 2176               std::unique_ptr<QBase> qrule
 
 2178               fine_fe->attach_quadrature_rule(qrule.get());
 
 2180               const std::vector<Point> & child_xyz =
 
 2183               for (
unsigned int c = 0, nc = elem.
n_children();
 
 2189                   fine_fe->reinit(elem.
child_ptr(c), context.side);
 
 2190                   fine_points.insert(fine_points.end(),
 
 2195               std::vector<Point> fine_qp;
 
 2198               context.elem_fe_reinit(&fine_qp);
 
 2201 #endif // LIBMESH_ENABLE_AMR 
 2202             context.side_fe_reinit();
 
 2204           const std::vector<dof_id_type> & dof_indices =
 
 2205             context.get_dof_indices(var);
 
 2207           std::vector<unsigned int> side_dofs;
 
 2209                                     context.side, side_dofs);
 
 2211           this->construct_projection
 
 2212             (dof_indices, side_dofs, var_component,
 
 2213              &side_node, proj_fe);
 
 2219 template <
typename FFunctor, 
typename GFunctor,
 
 2220           typename FValue, 
typename ProjectionAction>
 
 2224   LOG_SCOPE (
"project_interiors",
"GenericProjector");
 
 2226   const unsigned int sys_num = system.number();
 
 2229   for (
const auto & elem : range)
 
 2231       unsigned char dim = cast_int<unsigned char>(elem->dim());
 
 2233       context.pre_fe_reinit(system, elem);
 
 2237       for (
const auto & var : this->projector.variables)
 
 2239           const Variable & variable = system.variable(var);
 
 2244           const FEType & base_fe_type = variable.
type();
 
 2250           context.get_element_fe( var, fe, 
dim );
 
 2252           FEType fe_type = base_fe_type;
 
 2258           const unsigned int var_component =
 
 2259             system.variable_scalar_number(var, 0);
 
 2266               if (fe_type.
order > 1)
 
 2268                   const unsigned int first_interior_node =
 
 2269                     (elem->n_vertices() +
 
 2270                      ((elem->dim() > 2) * elem->n_edges()) +
 
 2271                      ((elem->dim() > 1) * elem->n_sides()));
 
 2272                   const unsigned int n_nodes = elem->n_nodes();
 
 2275                   for (
unsigned int n = first_interior_node; n < 
n_nodes; ++n)
 
 2277                       const Node & interior_node = elem->node_ref(n);
 
 2282                       FValue fval = f.eval_at_point
 
 2283                         (context, var_component, interior_node, system.time);
 
 2284                       insert_id(dof_id, fval, pid);
 
 2290 #ifdef LIBMESH_ENABLE_AMR 
 2293               std::vector<Point> fine_points;
 
 2295               std::unique_ptr<FEBase> fine_fe
 
 2298               std::unique_ptr<QBase> qrule
 
 2300               fine_fe->attach_quadrature_rule(qrule.get());
 
 2302               const std::vector<Point> & child_xyz =
 
 2305               for (
auto & child : elem->child_ref_range())
 
 2307                   fine_fe->reinit(&child);
 
 2308                   fine_points.insert(fine_points.end(),
 
 2313               std::vector<Point> fine_qp;
 
 2316               context.elem_fe_reinit(&fine_qp);
 
 2319 #endif // LIBMESH_ENABLE_AMR 
 2320             context.elem_fe_reinit();
 
 2322           const std::vector<dof_id_type> & dof_indices =
 
 2323             context.get_dof_indices(var);
 
 2325           const unsigned int n_dofs =
 
 2326             cast_int<unsigned int>(dof_indices.size());
 
 2328           std::vector<unsigned int> all_dofs(n_dofs);
 
 2329           std::iota(all_dofs.begin(), all_dofs.end(), 0);
 
 2331           this->construct_projection
 
 2332             (dof_indices, all_dofs, var_component,
 
 2340 template <
typename FFunctor, 
typename GFunctor,
 
 2341           typename FValue, 
typename ProjectionAction>
 
 2344                  ProjectionAction>::SubFunctor::find_dofs_to_send
 
 2352   if (owner != system.processor_id())
 
 2355       const DofMap & dof_map = system.get_dof_map();
 
 2359       std::vector<dof_id_type> node_dof_ids, patch_dof_ids;
 
 2360       for (
const auto & var : vars)
 
 2362           const Variable & variable = system.variable(var);
 
 2367           dof_map.
dof_indices(elem, node_num, node_dof_ids, var);
 
 2370                                     node_dof_ids.end()));
 
 2371       const std::vector<dof_id_type> & patch =
 
 2372         (*this->projector.nodes_to_elem)[node.
id()];
 
 2373       for (
const auto & elem_id : patch)
 
 2379           std::sort(patch_dof_ids.begin(), patch_dof_ids.end());
 
 2383           std::vector<dof_id_type> diff_ids(node_dof_ids.size());
 
 2384           auto it = std::set_difference(node_dof_ids.begin(), node_dof_ids.end(),
 
 2385                                         patch_dof_ids.begin(), patch_dof_ids.end(), diff_ids.begin());
 
 2386           diff_ids.resize(it-diff_ids.begin());
 
 2387           node_dof_ids.swap(diff_ids);
 
 2388           if (node_dof_ids.empty())
 
 2393       for (
auto id : node_dof_ids)
 
 2400 template <
typename FFunctor, 
typename GFunctor,
 
 2401           typename FValue, 
typename ProjectionAction>
 
 2404                  ProjectionAction>::send_and_insert_dof_values
 
 2405   (std::unordered_map<
dof_id_type, std::pair<FValue, processor_id_type>> & ids_to_push,
 
 2406    ProjectionAction & action) 
const 
 2410   std::unordered_map<processor_id_type, std::vector<dof_id_type>>
 
 2411     pushed_dof_ids, received_dof_ids;
 
 2412   std::unordered_map<processor_id_type, std::vector<typename TypeToSend<FValue>::type>>
 
 2413     pushed_dof_values, received_dof_values;
 
 2414   for (
auto & id_val_pid : ids_to_push)
 
 2419           pushed_dof_ids[pid].push_back(id_val_pid.first);
 
 2420           pushed_dof_values[pid].push_back(
convert_to_send(id_val_pid.second.first));
 
 2426   auto ids_action_functor =
 
 2427     [&action, &received_dof_ids, &received_dof_values]
 
 2429      const std::vector<dof_id_type> & 
data)
 
 2431       auto iter = received_dof_values.find(pid);
 
 2432       if (iter == received_dof_values.end())
 
 2435                          received_dof_ids.end());
 
 2436           received_dof_ids[pid] = 
data;
 
 2440           auto & vals = iter->second;
 
 2441           std::size_t vals_size = vals.size();
 
 2442           libmesh_assert_equal_to(vals_size, 
data.size());
 
 2443           for (std::size_t i=0; i != vals_size; ++i)
 
 2447               action.insert(
data[i], val);
 
 2449           received_dof_values.erase(iter);
 
 2455   auto values_action_functor =
 
 2456     [&action, &received_dof_ids, &received_dof_values]
 
 2458      const std::vector<typename TypeToSend<FValue>::type> & 
data)
 
 2460       auto iter = received_dof_ids.find(pid);
 
 2461       if (iter == received_dof_ids.end())
 
 2465                          received_dof_values.end());
 
 2466           received_dof_values[pid] = 
data;
 
 2470           auto & ids = iter->second;
 
 2471           std::size_t ids_size = ids.size();
 
 2472           libmesh_assert_equal_to(ids_size, 
data.size());
 
 2473           for (std::size_t i=0; i != ids_size; ++i)
 
 2477               action.insert(ids[i], val);
 
 2479           received_dof_ids.erase(iter);
 
 2483   Parallel::push_parallel_vector_data
 
 2484     (system.comm(), pushed_dof_ids, ids_action_functor);
 
 2486   Parallel::push_parallel_vector_data
 
 2487     (system.comm(), pushed_dof_values, values_action_functor);
 
 2497 template <
typename FFunctor, 
typename GFunctor,
 
 2498           typename FValue, 
typename ProjectionAction>
 
 2501   (
const std::vector<dof_id_type> & dof_indices_var,
 
 2502    const std::vector<unsigned int> & involved_dofs,
 
 2503    unsigned int var_component,
 
 2507   const std::vector<Real> & JxW = fe.
get_JxW();
 
 2508   const std::vector<std::vector<Real>> & phi = fe.
get_phi();
 
 2509   const std::vector<std::vector<RealGradient>> * dphi = 
nullptr;
 
 2510   const std::vector<Point> & xyz_values = fe.
get_xyz();
 
 2512   const std::unordered_map<dof_id_type, FValue> & ids_to_save =
 
 2513     this->projector.ids_to_save;
 
 2518   const unsigned int n_involved_dofs =
 
 2519     cast_int<unsigned int>(involved_dofs.size());
 
 2521   std::vector<dof_id_type> free_dof_ids;
 
 2523   std::vector<char> dof_is_fixed(n_involved_dofs, 
false); 
 
 2527       const dof_id_type id = dof_indices_var[involved_dofs[i]];
 
 2528       auto iter = ids_to_save.find(
id);
 
 2529       if (iter == ids_to_save.end())
 
 2530         free_dof_ids.push_back(
id);
 
 2533           dof_is_fixed[i] = 
true;
 
 2534           Uinvolved(i) = iter->second;
 
 2538   const unsigned int free_dofs = free_dof_ids.
size();
 
 2553   const unsigned int n_qp =
 
 2554     cast_int<unsigned int>(xyz_values.size());
 
 2557   for (
unsigned int qp=0; qp<n_qp; qp++)
 
 2560       FValue fineval = f.eval_at_point(context,
 
 2567         finegrad = g->eval_at_point(context,
 
 2573       for (
unsigned int sidei=0, freei=0;
 
 2574            sidei != n_involved_dofs; ++sidei)
 
 2576           unsigned int i = involved_dofs[sidei];
 
 2578           if (dof_is_fixed[sidei])
 
 2580           for (
unsigned int sidej=0, freej=0;
 
 2581                sidej != n_involved_dofs; ++sidej)
 
 2583               unsigned int j = involved_dofs[sidej];
 
 2584               if (dof_is_fixed[sidej])
 
 2585                 Fe(freei) -= phi[i][qp] * phi[j][qp] *
 
 2586                   JxW[qp] * Uinvolved(sidej);
 
 2588                 Ke(freei,freej) += phi[i][qp] *
 
 2589                   phi[j][qp] * JxW[qp];
 
 2592                   if (dof_is_fixed[sidej])
 
 2593                     Fe(freei) -= ( (*dphi)[i][qp] *
 
 2595                       JxW[qp] * Uinvolved(sidej);
 
 2597                     Ke(freei,freej) += ( (*dphi)[i][qp] *
 
 2601               if (!dof_is_fixed[sidej])
 
 2604           Fe(freei) += phi[i][qp] * fineval * JxW[qp];
 
 2606             Fe(freei) += (finegrad * (*dphi)[i][qp] ) *
 
 2617   insert_ids(free_dof_ids, Ufree.
get_values(), pid);
 
 2623 #endif // GENERIC_PROJECTOR_H