libMesh
generic_projector.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef GENERIC_PROJECTOR_H
21 #define GENERIC_PROJECTOR_H
22 
23 // C++ includes
24 #include <vector>
25 
26 // libMesh includes
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" // libmesh_make_unique
41 
42 // TIMPI includes
43 #include "timpi/parallel_sync.h"
44 
45 
46 namespace libMesh
47 {
48 
55 template <typename T>
56 struct TypeToSend {
57  typedef T type;
58 };
59 
60 template <typename T>
61 const typename TypeToSend<T>::type convert_to_send (const T& in)
62 { return in; }
63 
64 template <typename SendT, typename T>
65 void convert_from_receive (SendT & received, T & converted)
66 { converted = received; }
67 
78 template <typename FFunctor, typename GFunctor,
79  typename FValue, typename ProjectionAction>
81 {
82 private:
83  const System & system;
84 
85  // For TBB compatibility and thread safety we'll copy these in
86  // operator()
87  const FFunctor & master_f;
88  const GFunctor * master_g; // Needed for C1 type elements only
90  const ProjectionAction & master_action;
91  const std::vector<unsigned int> & variables;
92  std::unordered_map<dof_id_type, std::vector<dof_id_type>> * nodes_to_elem;
94 
95 public:
96  GenericProjector (const System & system_in,
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) :
103  system(system_in),
104  master_f(f_in),
105  master_g(g_in),
106  g_was_copied(false),
107  map_was_created(!nodes_to_elem_in),
108  master_action(act_in),
109  variables(variables_in),
110  nodes_to_elem(nodes_to_elem_in)
111  {
112  if (map_was_created) // past tense misnomer here
113  {
114  nodes_to_elem = new
115  std::unordered_map<dof_id_type, std::vector<dof_id_type>>;
117  }
118  }
119 
121  system(in.system),
122  master_f(in.master_f),
123  master_g(in.master_g ? new GFunctor(*in.master_g) : nullptr),
126  variables(in.variables),
128  {}
129 
131  {
132  if (g_was_copied)
133  delete master_g;
134  if (map_was_created)
135  delete nodes_to_elem;
136  }
137 
138  // The single "project" function handles all the threaded projection
139  // calculations and intervening MPI communications
140  void project(const ConstElemRange & range);
141 
142  // We generally need to hang on to every value we've calculated
143  // until we're all done, because later projection calculations
144  // depend on boundary data from earlier calculations.
145  std::unordered_map<dof_id_type, FValue> ids_to_save;
146 
147  // This needs to be sorted so we can get sorted dof indices cheaply
148  // later
149  typedef std::set<unsigned int> var_set;
150 
151  // We now need to divide our threadable work into stages, to allow
152  // for the potential necessity of MPI communication in between them.
153  struct SubFunctor {
155 
157 
158  // While we're working on nodes, also figure out which ghosted
159  // nodes will have data we might need to send to their owners
160  // instead of being acted on by ourselves.
161  //
162  // We keep track of which dof ids we might need to send, and what
163  // values those ids should get (along with a pprocessor_id to
164  // leave invalid in case *we* can't compute those values either).
165  std::unordered_map<dof_id_type, std::pair<FValue, processor_id_type>> new_ids_to_push;
166 
167  // We'll hang on to any new ids to save on a per-thread basis
168  // because we won't need them until subsequent phases
169  std::unordered_map<dof_id_type, FValue> new_ids_to_save;
170 
171  // Helper function for filling new_ids_to_push
172  void find_dofs_to_send (const Node & node,
173  const Elem & elem,
174  unsigned short node_num,
175  const var_set & vars);
176 
177  // When we have new data to act on, we may also need to save it
178  // and get ready to push it.
179  void insert_id(dof_id_type id,
180  const FValue & val,
181  processor_id_type pid);
182 
183  void insert_ids(const std::vector<dof_id_type> & ids,
184  const std::vector<FValue> & vals,
185  processor_id_type pid);
186 
187  // Our subclasses will need to merge saved ids
188  void join (const SubFunctor & other);
189 
190  protected:
191  // For thread safety and efficiency we'll want thread-local copies
192  // of various functors
193  ProjectionAction action;
194  FFunctor f;
195 
196  // Each stage of the work we do will require us to prepare
197  // FEMContext objects to assist in the work.
199 
200  // These get used often enough to cache them
201  std::vector<FEContinuity> conts;
202 
203  const System & system;
204  };
205 
206 
207  typedef std::pair<const Node *, std::tuple<const Elem *, unsigned short, var_set>>
209 
212 
213 
214  struct SubProjector : public SubFunctor {
216 
217  using SubFunctor::action;
218  using SubFunctor::f;
219  using SubFunctor::system;
220  using SubFunctor::context;
221  using SubFunctor::conts;
222  using SubFunctor::insert_id;
224 
225  protected:
226  // Projections of C1 elements require a gradient as well
227  std::unique_ptr<GFunctor> g;
228 
230  (const std::vector<dof_id_type> & dof_indices_var,
231  const std::vector<unsigned int> & involved_dofs,
232  unsigned int var_component,
233  const Node * node,
234  const FEBase & fe);
235  };
236 
237 
238  // First we'll copy DoFs where we can, while sorting remaining DoFs
239  // for interpolation and projection later.
240  struct SortAndCopy : public SubFunctor {
242 
244 
245  using SubFunctor::action;
246  using SubFunctor::f;
247  using SubFunctor::system;
248  using SubFunctor::context;
250 
251  void operator() (const ConstElemRange & range);
252 
253  // We need to merge saved multimaps when working from multiple
254  // threads
255  void join (const SortAndCopy & other);
256 
257  // For vertices, edges and sides, we need to know what element,
258  // what local vertex/edge/side number, and what set of variables
259  // to evaluate.
260  typedef std::unordered_multimap<const Node *, std::tuple<const Elem *, unsigned short, var_set>> ves_multimap;
261 
263  std::vector<const Elem *> interiors;
264  };
265 
266 
267  struct ProjectVertices : public SubProjector {
269 
271 
272  using SubProjector::action;
273  using SubProjector::f;
274  using SubProjector::g;
275  using SubProjector::system;
276  using SubProjector::context;
277  using SubProjector::conts;
278 
282 
283  void operator() (const node_range & range);
284  };
285 
286 
287  struct ProjectEdges : public SubProjector {
289 
291 
292  using SubProjector::action;
293  using SubProjector::f;
294  using SubProjector::g;
295  using SubProjector::system;
296  using SubProjector::context;
297  using SubProjector::conts;
298 
302 
303  void operator() (const node_range & range);
304  };
305 
306  struct ProjectSides : public SubProjector {
308 
310 
311  using SubProjector::action;
312  using SubProjector::f;
313  using SubProjector::g;
314  using SubProjector::system;
315  using SubProjector::context;
316  using SubProjector::conts;
317 
321 
322  void operator() (const node_range & range);
323  };
324 
325 
328 
329  struct ProjectInteriors : public SubProjector {
331 
333 
334  using SubProjector::action;
335  using SubProjector::f;
336  using SubProjector::g;
337  using SubProjector::system;
338  using SubProjector::context;
339  using SubProjector::conts;
340 
344 
345  void operator() (const interior_range & range);
346  };
347 
349  (std::unordered_map<dof_id_type, std::pair<FValue, processor_id_type>> & ids_to_push,
350  ProjectionAction & action) const;
351 };
352 
353 
363 template <typename Val>
365 {
366 private:
368 
369 public:
371  target_vector(target_vec) {}
372 
374  Val val)
375  {
376  // Lock the new vector since it is shared among threads.
377  {
378  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
379  target_vector.set(id, val);
380  }
381  }
382 
383 
384  void insert(const std::vector<dof_id_type> & dof_indices,
385  const DenseVector<Val> & Ue)
386  {
387  const numeric_index_type
390 
391  unsigned int size = Ue.size();
392 
393  libmesh_assert_equal_to(size, dof_indices.size());
394 
395  // Lock the new vector since it is shared among threads.
396  {
397  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
398 
399  for (unsigned int i = 0; i != size; ++i)
400  if ((dof_indices[i] >= first) && (dof_indices[i] < last))
401  target_vector.set(dof_indices[i], Ue(i));
402  }
403  }
404 
405 };
406 
407 
416 template <typename Output>
418 {
419 public:
420  FEMFunctionWrapper(const FEMFunctionBase<Output> & f) : _f(f.clone()) {}
421 
423  _f(fw._f->clone()) {}
424 
425  void init_context (FEMContext & c) { _f->init_context(c); }
426 
427  Output eval_at_node (const FEMContext & c,
428  unsigned int i,
429  unsigned int /*elem_dim*/,
430  const Node & n,
431  bool /*extra_hanging_dofs*/,
432  const Real time)
433  { return _f->component(c, i, n, time); }
434 
435  Output eval_at_point (const FEMContext & c,
436  unsigned int i,
437  const Point & n,
438  const Real time)
439  { return _f->component(c, i, n, time); }
440 
441  bool is_grid_projection() { return false; }
442 
443  void eval_old_dofs (const Elem &,
444  unsigned int,
445  unsigned int,
446  std::vector<dof_id_type> &,
447  std::vector<Output> &)
448  { libmesh_error(); }
449 
450  void eval_old_dofs (const Elem &,
451  const FEType &,
452  unsigned int,
453  unsigned int,
454  std::vector<dof_id_type> &,
455  std::vector<Output> &)
456  { libmesh_error(); }
457 
458 private:
459  std::unique_ptr<FEMFunctionBase<Output>> _f;
460 };
461 
462 
473 #ifdef LIBMESH_ENABLE_AMR
474 template <typename Output,
475  void (FEMContext::*point_output) (unsigned int,
476  const Point &,
477  Output &,
478  const Real) const>
480 {
481 public:
483  last_elem(nullptr),
484  sys(sys_in),
485  old_context(sys_in)
486  {
487  // We'll be queried for components but we'll typically be looking
488  // up data by variables, and those indices don't always match
489  for (auto v : IntRange<unsigned int>(0, sys.n_vars()))
490  {
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;
495  }
496  }
497 
499  last_elem(nullptr),
500  sys(in.sys),
501  old_context(sys),
502  component_to_var(in.component_to_var)
503  {
504  }
505 
506  static void get_shape_outputs(FEBase & fe);
507 
508  // Integrating on new mesh elements, we won't yet have an up to date
509  // current_local_solution.
511  {
513 
514  // Loop over variables, to prerequest
515  for (auto var : IntRange<unsigned int>(0, sys.n_vars()))
516  {
517  FEBase * fe = nullptr;
518  const std::set<unsigned char> & elem_dims =
519  old_context.elem_dimensions();
520 
521  for (const auto & dim : elem_dims)
522  {
523  old_context.get_element_fe(var, fe, dim);
524  get_shape_outputs(*fe);
525  }
526  }
527  }
528 
529  bool is_grid_projection() { return true; }
530 
531 protected:
532  void check_old_context (const FEMContext & c)
533  {
534  LOG_SCOPE ("check_old_context(c)", "OldSolutionBase");
535  const Elem & elem = c.get_elem();
536  if (last_elem != &elem)
537  {
538  if (elem.refinement_flag() == Elem::JUST_REFINED)
539  {
540  old_context.pre_fe_reinit(sys, elem.parent());
541  }
542  else if (elem.refinement_flag() == Elem::JUST_COARSENED)
543  {
544  libmesh_error();
545  }
546  else
547  {
548  if (!elem.old_dof_object)
549  {
550  libmesh_error();
551  }
552 
553  old_context.pre_fe_reinit(sys, &elem);
554  }
555 
556  last_elem = &elem;
557  }
558  else
559  {
560  libmesh_assert(old_context.has_elem());
561  }
562  }
563 
564 
565  bool check_old_context (const FEMContext & c, const Point & p)
566  {
567  LOG_SCOPE ("check_old_context(c,p)", "OldSolutionBase");
568  const Elem & elem = c.get_elem();
569  if (last_elem != &elem)
570  {
571  if (elem.refinement_flag() == Elem::JUST_REFINED)
572  {
573  old_context.pre_fe_reinit(sys, elem.parent());
574  }
575  else if (elem.refinement_flag() == Elem::JUST_COARSENED)
576  {
577  // Find the child with this point. Use out_of_elem_tol
578  // (in physical space, which may correspond to a large
579  // tolerance in master space!) to allow for out-of-element
580  // finite differencing of mixed gradient terms. Pray we
581  // have no quadrature locations which are within 1e-5 of
582  // the element subdivision boundary but are not exactly on
583  // that boundary.
584  const Real master_tol = out_of_elem_tol / elem.hmax() * 2;
585 
586  for (auto & child : elem.child_ref_range())
587  if (child.close_to_point(p, master_tol))
588  {
589  old_context.pre_fe_reinit(sys, &child);
590  break;
591  }
592 
594  (old_context.get_elem().close_to_point(p, master_tol));
595  }
596  else
597  {
598  if (!elem.old_dof_object)
599  return false;
600 
601  old_context.pre_fe_reinit(sys, &elem);
602  }
603 
604  last_elem = &elem;
605  }
606  else
607  {
608  libmesh_assert(old_context.has_elem());
609 
610  const Real master_tol = out_of_elem_tol / elem.hmax() * 2;
611 
612  if (!old_context.get_elem().close_to_point(p, master_tol))
613  {
614  libmesh_assert_equal_to
616 
617  for (auto & child : elem.child_ref_range())
618  if (child.close_to_point(p, master_tol))
619  {
620  old_context.pre_fe_reinit(sys, &child);
621  break;
622  }
623 
625  (old_context.get_elem().close_to_point(p, master_tol));
626  }
627  }
628 
629  return true;
630  }
631 
632 protected:
633  const Elem * last_elem;
634  const System & sys;
636  std::vector<unsigned int> component_to_var;
637 
638  static const Real out_of_elem_tol;
639 };
640 
641 
650 template <typename Output,
651  void (FEMContext::*point_output) (unsigned int,
652  const Point &,
653  Output &,
654  const Real) const>
656 {
657 public:
659  const NumericVector<Number> & old_sol) :
660  OldSolutionBase<Output, point_output>(sys_in),
661  old_solution(old_sol)
662  {
663  this->old_context.set_algebraic_type(FEMContext::OLD);
664  this->old_context.set_custom_solution(&old_solution);
665  }
666 
668  OldSolutionBase<Output, point_output>(in.sys),
669  old_solution(in.old_solution)
670  {
671  this->old_context.set_algebraic_type(FEMContext::OLD);
672  this->old_context.set_custom_solution(&old_solution);
673  }
674 
675 
676  Output eval_at_node (const FEMContext & c,
677  unsigned int i,
678  unsigned int elem_dim,
679  const Node & n,
680  bool /* extra_hanging_dofs */,
681  Real /* time */ =0.);
682 
683  Output eval_at_point(const FEMContext & c,
684  unsigned int i,
685  const Point & p,
686  Real /* time */ =0.)
687  {
688  LOG_SCOPE ("eval_at_point()", "OldSolutionValue");
689 
690  if (!this->check_old_context(c, p))
691  return 0;
692 
693  // Handle offset from non-scalar components in previous variables
694  libmesh_assert_less(i, this->component_to_var.size());
695  unsigned int var = this->component_to_var[i];
696 
697  Output n;
698  (this->old_context.*point_output)(var, p, n, this->out_of_elem_tol);
699  return n;
700  }
701 
702  void eval_old_dofs (const Elem & elem,
703  unsigned int node_num,
704  unsigned int var_num,
705  std::vector<dof_id_type> & indices,
706  std::vector<Output> & values)
707  {
708  LOG_SCOPE ("eval_old_dofs(node)", "OldSolutionValue");
709 
710  this->sys.get_dof_map().dof_indices(elem, node_num, indices, var_num);
711 
712  std::vector<dof_id_type> old_indices;
713 
714  this->sys.get_dof_map().old_dof_indices(elem, node_num, old_indices, var_num);
715 
716  libmesh_assert_equal_to (old_indices.size(), indices.size());
717 
718  values.resize(old_indices.size());
719 
720  old_solution.get(old_indices, values);
721  }
722 
723 
724  void eval_old_dofs (const Elem & elem,
725  const FEType & fe_type,
726  unsigned int sys_num,
727  unsigned int var_num,
728  std::vector<dof_id_type> & indices,
729  std::vector<Output> & values)
730  {
731  LOG_SCOPE ("eval_old_dofs(elem)", "OldSolutionValue");
732 
733  // We're only to be asked for old dofs on elements that can copy
734  // them through DO_NOTHING or through refinement.
735  const Elem & old_elem =
736  (elem.refinement_flag() == Elem::JUST_REFINED) ?
737  *elem.parent() : elem;
738 
739  // If there are any element-based DOF numbers, get them
740  const unsigned int nc = FEInterface::n_dofs_per_elem(elem.dim(),
741  fe_type,
742  elem.type());
743 
744  std::vector<dof_id_type> old_dof_indices(nc);
745  indices.resize(nc);
746 
747  // We should never have fewer dofs than necessary on an
748  // element unless we're getting indices on a parent element,
749  // and we should never need those indices
750  if (nc != 0)
751  {
752  libmesh_assert(old_elem.old_dof_object);
753 
754  const std::pair<unsigned int, unsigned int>
755  vg_and_offset = elem.var_to_vg_and_offset(sys_num,var_num);
756  const unsigned int vg = vg_and_offset.first;
757  const unsigned int vig = vg_and_offset.second;
758 
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);
762 
763  for (unsigned int i=0; i<nc; i++)
764  {
765  const dof_id_type d_old =
766  old_elem.old_dof_object->dof_number(sys_num, vg, vig, i, n_comp);
767  const dof_id_type d_new =
768  elem.dof_number(sys_num, vg, vig, i, n_comp);
769  libmesh_assert_not_equal_to (d_old, DofObject::invalid_id);
770  libmesh_assert_not_equal_to (d_new, DofObject::invalid_id);
771 
772  old_dof_indices[i] = d_old;
773  indices[i] = d_new;
774  }
775  }
776 
777  values.resize(nc);
778 
779  old_solution.get(old_dof_indices, values);
780  }
781 
782 private:
784 };
785 
786 
787 template<>
788 inline void
790 {
791  fe.get_phi();
792 }
793 
794 
795 template<>
796 inline void
798 {
799  fe.get_dphi();
800 }
801 
802 
803 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
804 template<>
805 inline void
807 {
808  fe.get_phi();
809 }
810 
811 
812 template<>
813 inline void
815 {
816  fe.get_dphi();
817 }
818 #endif // LIBMESH_USE_COMPLEX_NUMBERS
819 
820 
821 template<>
822 inline
823 Number
826  unsigned int i,
827  unsigned int /* elem_dim */,
828  const Node & n,
829  bool extra_hanging_dofs,
830  Real /* time */)
831 {
832  LOG_SCOPE ("Number eval_at_node()", "OldSolutionValue");
833 
834  // This should only be called on vertices
835  libmesh_assert_less(c.get_elem().get_node_index(&n),
836  c.get_elem().n_vertices());
837 
838  // Handle offset from non-scalar components in previous variables
839  libmesh_assert_less(i, this->component_to_var.size());
840  unsigned int var = this->component_to_var[i];
841 
842  // Optimize for the common case, where this node was part of the
843  // old solution.
844  //
845  // Be sure to handle cases where the variable wasn't defined on
846  // this node (due to changing subdomain support) or where the
847  // variable has no components on this node (due to Elem order
848  // exceeding FE order) or where the old_dof_object dofs might
849  // correspond to non-vertex dofs (due to extra_hanging_dofs and
850  // refinement)
851 
853 
854  if (n.old_dof_object &&
855  (!extra_hanging_dofs ||
856  flag == Elem::JUST_COARSENED ||
857  flag == Elem::DO_NOTHING) &&
858  n.old_dof_object->n_vars(sys.number()) &&
859  n.old_dof_object->n_comp(sys.number(), var))
860  {
861  const dof_id_type old_id =
862  n.old_dof_object->dof_number(sys.number(), var, 0);
863  return old_solution(old_id);
864  }
865 
866  return this->eval_at_point(c, i, n, 0);
867 }
868 
869 
870 
871 template<>
872 inline
873 Gradient
876  unsigned int i,
877  unsigned int elem_dim,
878  const Node & n,
879  bool extra_hanging_dofs,
880  Real /* time */)
881 {
882  LOG_SCOPE ("Gradient eval_at_node()", "OldSolutionValue");
883 
884  // This should only be called on vertices
885  libmesh_assert_less(c.get_elem().get_node_index(&n),
886  c.get_elem().n_vertices());
887 
888  // Handle offset from non-scalar components in previous variables
889  libmesh_assert_less(i, this->component_to_var.size());
890  unsigned int var = this->component_to_var[i];
891 
892  // Optimize for the common case, where this node was part of the
893  // old solution.
894  //
895  // Be sure to handle cases where the variable wasn't defined on
896  // this node (due to changing subdomain support) or where the
897  // variable has no components on this node (due to Elem order
898  // exceeding FE order) or where the old_dof_object dofs might
899  // correspond to non-vertex dofs (due to extra_hanging_dofs and
900  // refinement)
901 
903 
904  if (n.old_dof_object &&
905  (!extra_hanging_dofs ||
906  flag == Elem::JUST_COARSENED ||
907  flag == Elem::DO_NOTHING) &&
908  n.old_dof_object->n_vars(sys.number()) &&
909  n.old_dof_object->n_comp(sys.number(), var))
910  {
911  Gradient g;
912  for (unsigned int d = 0; d != elem_dim; ++d)
913  {
914  const dof_id_type old_id =
915  n.old_dof_object->dof_number(sys.number(), var, d+1);
916  g(d) = old_solution(old_id);
917  }
918  return g;
919  }
920 
921  return this->eval_at_point(c, i, n, 0);
922 }
923 
924 
925 
926 
927 
928 template <>
930 
931 template <>
933 
934 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
935 template <>
937 
938 template <>
940 #endif // LIBMESH_USE_COMPLEX_NUMBERS
941 
942 #endif // LIBMESH_ENABLE_AMR
943 
948 template <typename FFunctor, typename GFunctor,
949  typename FValue, typename ProjectionAction>
951  (const ConstElemRange & range)
952 {
953  LOG_SCOPE ("project", "GenericProjector");
954 
955  // Unless we split sort and copy into two passes we can't know for
956  // sure ahead of time whether we need to save the copied ids
957  done_saving_ids = false;
958 
959  SortAndCopy sort_work(*this);
960  Threads::parallel_reduce (range, sort_work);
961  ProjectionAction action(master_action);
962 
963  // Keep track of dof ids and values to send to other ranks
964  std::unordered_map<dof_id_type, std::pair<FValue, processor_id_type>> ids_to_push;
965 
966  ids_to_push.insert(sort_work.new_ids_to_push.begin(),
967  sort_work.new_ids_to_push.end());
968  ids_to_save.insert(sort_work.new_ids_to_save.begin(),
969  sort_work.new_ids_to_save.end());
970 
971  std::vector<node_projection> vertices(sort_work.vertices.begin(),
972  sort_work.vertices.end());
973 
974  done_saving_ids = sort_work.edges.empty() &&
975  sort_work.sides.empty() && sort_work.interiors.empty();
976  system.comm().max(done_saving_ids);
977 
978  {
979  ProjectVertices project_vertices(*this);
980  Threads::parallel_reduce (node_range(&vertices), project_vertices);
981  ids_to_push.insert(project_vertices.new_ids_to_push.begin(),
982  project_vertices.new_ids_to_push.end());
983  ids_to_save.insert(project_vertices.new_ids_to_save.begin(),
984  project_vertices.new_ids_to_save.end());
985  }
986 
987  done_saving_ids = sort_work.sides.empty() && sort_work.interiors.empty();
988  system.comm().max(done_saving_ids);
989 
990  this->send_and_insert_dof_values(ids_to_push, action);
991 
992  {
993  std::vector<node_projection> edges(sort_work.edges.begin(), sort_work.edges.end());
994  ProjectEdges project_edges(*this);
995  Threads::parallel_reduce (node_range(&edges), project_edges);
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());
1000  }
1001 
1002  done_saving_ids = sort_work.interiors.empty();
1003  system.comm().max(done_saving_ids);
1004 
1005  this->send_and_insert_dof_values(ids_to_push, action);
1006 
1007  {
1008  std::vector<node_projection> sides(sort_work.sides.begin(), sort_work.sides.end());
1009  ProjectSides project_sides(*this);
1010  Threads::parallel_reduce (node_range(&sides), project_sides);
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());
1015  }
1016 
1017  done_saving_ids = true;
1018  this->send_and_insert_dof_values(ids_to_push, action);
1019 
1020  // No ids to save or push this time, but we still use a reduce since
1021  // nominally ProjectInteriors still has non-const operator()
1022  ProjectInteriors project_interiors(*this);
1024  project_interiors);
1025 }
1026 
1027 
1028 template <typename FFunctor, typename GFunctor,
1029  typename FValue, typename ProjectionAction>
1032  projector(p), action(p.master_action), f(p.master_f),
1033  context(p.system), conts(p.system.n_vars()), system(p.system)
1034 {
1035  // Loop over all the variables we've been requested to project, to
1036  // pre-request
1037  for (const auto & var : this->projector.variables)
1038  {
1039  // FIXME: Need to generalize this to vector-valued elements. [PB]
1040  FEBase * fe = nullptr;
1041  FEBase * side_fe = nullptr;
1042  FEBase * edge_fe = nullptr;
1043 
1044  const std::set<unsigned char> & elem_dims =
1045  context.elem_dimensions();
1046 
1047  for (const auto & dim : elem_dims)
1048  {
1049  context.get_element_fe( var, fe, dim );
1050  if (fe->get_fe_type().family == SCALAR)
1051  continue;
1052  if (dim > 1)
1053  context.get_side_fe( var, side_fe, dim );
1054  if (dim > 2)
1055  context.get_edge_fe( var, edge_fe );
1056 
1057  fe->get_JxW();
1058  fe->get_xyz();
1059  fe->get_JxW();
1060 
1061  fe->get_phi();
1062  if (dim > 1)
1063  {
1064  side_fe->get_JxW();
1065  side_fe->get_xyz();
1066  side_fe->get_phi();
1067  }
1068  if (dim > 2)
1069  {
1070  edge_fe->get_JxW();
1071  edge_fe->get_xyz();
1072  edge_fe->get_phi();
1073  }
1074 
1075  const FEContinuity cont = fe->get_continuity();
1076  this->conts[var] = cont;
1077  if (cont == C_ONE)
1078  {
1079  fe->get_dphi();
1080  if (dim > 1)
1081  side_fe->get_dphi();
1082  if (dim > 2)
1083  edge_fe->get_dphi();
1084  }
1085  }
1086  }
1087 
1088  // Now initialize any user requested shape functions, xyz vals, etc
1089  f.init_context(context);
1090 }
1091 
1092 
1093 template <typename FFunctor, typename GFunctor,
1094  typename FValue, typename ProjectionAction>
1097  SubFunctor(p)
1098 {
1099  if (p.master_g)
1100  g = libmesh_make_unique<GFunctor>(*p.master_g);
1101 
1102 #ifndef NDEBUG
1103  // Our C1 elements need gradient information
1104  for (const auto & var : this->projector.variables)
1105  if (this->conts[var] == C_ONE)
1106  libmesh_assert(g);
1107 #endif
1108 
1109  if (g)
1110  g->init_context(context);
1111 }
1112 
1113 
1114 template <typename FFunctor, typename GFunctor,
1115  typename FValue, typename ProjectionAction>
1117  (dof_id_type id, const FValue & val, processor_id_type pid)
1118 {
1119  auto iter = new_ids_to_push.find(id);
1120  if (iter == new_ids_to_push.end())
1121  action.insert(id, val);
1122  else
1123  {
1125  iter->second = std::make_pair(val, pid);
1126  }
1127  if (!this->projector.done_saving_ids)
1128  {
1129  libmesh_assert(!new_ids_to_save.count(id));
1130  new_ids_to_save[id] = val;
1131  }
1132 }
1133 
1134 
1135 template <typename FFunctor, typename GFunctor,
1136  typename FValue, typename ProjectionAction>
1138  (const std::vector<dof_id_type> & ids, const std::vector<FValue> & vals, processor_id_type pid)
1139 {
1140  libmesh_assert_equal_to(ids.size(), vals.size());
1141  for (auto i : index_range(ids))
1142  {
1143  const dof_id_type id = ids[i];
1144  const FValue & val = vals[i];
1145 
1146  auto iter = new_ids_to_push.find(id);
1147  if (iter == new_ids_to_push.end())
1148  action.insert(id, val);
1149  else
1150  {
1152  iter->second = std::make_pair(val, pid);
1153  }
1154  if (!this->projector.done_saving_ids)
1155  {
1156  libmesh_assert(!new_ids_to_save.count(id));
1157  new_ids_to_save[id] = val;
1158  }
1159  }
1160 }
1161 
1162 
1163 template <typename FFunctor, typename GFunctor,
1164  typename FValue, typename ProjectionAction>
1167 {
1168  new_ids_to_push.insert(other.new_ids_to_push.begin(), other.new_ids_to_push.end());
1169  new_ids_to_save.insert(other.new_ids_to_save.begin(), other.new_ids_to_save.end());
1170 }
1171 
1172 
1173 template <typename FFunctor, typename GFunctor,
1174  typename FValue, typename ProjectionAction>
1176  (const ConstElemRange & range)
1177 {
1178  // Look at all the elements in the range. Directly copy values from
1179  // unchanged elements. For other elements, determine sets of
1180  // vertices, edge nodes, and side nodes to project.
1181  //
1182  // As per our other weird nomenclature, "sides" means faces in 3D
1183  // and edges in 2D, and "edges" gets skipped in 2D
1184  //
1185  // This gets tricky in the case of subdomain-restricted
1186  // variables, for which we might need to do the same projection
1187  // from different elements when evaluating different variables.
1188  // We'll keep track of which variables can be projected from which
1189  // elements.
1190 
1191  // If we copy DoFs on a node, keep track of it so we don't bother
1192  // with any redundant interpolations or projections later.
1193  //
1194  // We still need to keep track of *which* variables get copied,
1195  // since we may be able to copy elements which lack some
1196  // subdomain-restricted variables.
1197  //
1198  // For extra_hanging_dofs FE types, we'll keep track of which
1199  // variables were copied from vertices, and which from edges/sides,
1200  // because those types need copies made from *both* on hanging
1201  // nodes.
1202  std::unordered_map<const Node *, std::pair<var_set, var_set>> copied_nodes;
1203 
1204  const unsigned int sys_num = system.number();
1205 
1206  // At hanging nodes for variables with extra hanging dofs we'll need
1207  // to do projections *separately* from vertex elements and side/edge
1208  // elements.
1209  std::vector<unsigned short> extra_hanging_dofs;
1210  bool all_extra_hanging_dofs = true;
1211  for (auto v_num : this->projector.variables)
1212  {
1213  if (extra_hanging_dofs.size() <= v_num)
1214  extra_hanging_dofs.resize(v_num+1, false);
1215  extra_hanging_dofs[v_num] =
1216  FEInterface::extra_hanging_dofs(system.variable(v_num).type());
1217 
1218  if (!extra_hanging_dofs[v_num])
1219  all_extra_hanging_dofs = false;
1220  }
1221 
1222  for (const auto & elem : range)
1223  {
1224  // If we're doing AMR, we might be able to copy more DoFs than
1225  // we interpolate or project.
1226  bool copy_this_elem = false;
1227 
1228 #ifdef LIBMESH_ENABLE_AMR
1229  // If we're projecting from an old grid
1230  if (f.is_grid_projection())
1231  {
1232  // If this element doesn't have an old_dof_object, but it
1233  // wasn't just refined or just coarsened into activity, then
1234  // it must be newly added, so the user is responsible for
1235  // setting the new dofs on it during a grid projection.
1236  if (!elem->old_dof_object &&
1237  elem->refinement_flag() != Elem::JUST_REFINED &&
1238  elem->refinement_flag() != Elem::JUST_COARSENED)
1239  continue;
1240 
1241  // If this is an unchanged element, just copy everything
1242  if ((elem->refinement_flag() != Elem::JUST_REFINED &&
1243  elem->refinement_flag() != Elem::JUST_COARSENED &&
1244  elem->p_refinement_flag() != Elem::JUST_REFINED &&
1245  elem->p_refinement_flag() != Elem::JUST_COARSENED))
1246  copy_this_elem = true;
1247  else
1248  {
1249  bool reinitted = false;
1250 
1251  // If this element has a low order monomial which has
1252  // merely been h refined, copy it.
1253  for (auto v_num : this->projector.variables)
1254  {
1255  const Variable & var = system.variable(v_num);
1256  if (!var.active_on_subdomain(elem->subdomain_id()))
1257  continue;
1258  FEType fe_type = var.type();
1259 
1260  if (fe_type.family == MONOMIAL &&
1261  fe_type.order == CONSTANT &&
1262  elem->p_level() == 0 &&
1263  elem->refinement_flag() != Elem::JUST_COARSENED &&
1264  elem->p_refinement_flag() != Elem::JUST_COARSENED)
1265  {
1266  if (!reinitted)
1267  {
1268  reinitted = true;
1269  context.pre_fe_reinit(system, elem);
1270  }
1271 
1272  std::vector<FValue> Ue(1);
1273  std::vector<dof_id_type> elem_dof_ids(1);
1274 
1275  f.eval_old_dofs(*elem, fe_type, sys_num, v_num,
1276  elem_dof_ids, Ue);
1277 
1278  action.insert(elem_dof_ids[0], Ue[0]);
1279  }
1280  }
1281  }
1282  }
1283 #endif // LIBMESH_ENABLE_AMR
1284 
1285  const int dim = elem->dim();
1286 
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();
1290 
1291  // In 1-D we already handle our sides as vertices
1292  const unsigned int n_sides = (dim > 1) * elem->n_sides();
1293 
1294  // What variables are supported on each kind of node on this elem?
1295  var_set vertex_vars, edge_vars, side_vars;
1296 
1297  // If we have non-vertex nodes, the first is an edge node, but
1298  // if we're in 2D we'll call that a side node
1299  const bool has_edge_nodes = (n_nodes > n_vertices && dim > 2);
1300 
1301  // If we have even more nodes, the next is a side node.
1302  const bool has_side_nodes =
1303  (n_nodes > n_vertices + ((dim > 2) * n_edges));
1304 
1305  // We may be out of nodes at this point or we have interior
1306  // nodes which may have DoFs to project too
1307  const bool has_interior_nodes =
1308  (n_nodes > n_vertices + ((dim > 2) * n_edges) + n_sides);
1309 
1310  for (auto v_num : this->projector.variables)
1311  {
1312  const Variable & var = this->projector.system.variable(v_num);
1313  if (!var.active_on_subdomain(elem->subdomain_id()))
1314  continue;
1315  FEType fe_type = var.type();
1316  fe_type.order =
1317  libMesh::Order (fe_type.order + elem->p_level());
1318  const ElemType elem_type = elem->type();
1319 
1320  if (FEInterface::n_dofs_at_node(dim, fe_type, elem_type, 0))
1321  vertex_vars.insert(vertex_vars.end(), v_num);
1322 
1323  // The first non-vertex node is always an edge node if those
1324  // exist. All edge nodes have the same number of DoFs
1325  if (has_edge_nodes)
1326  if (FEInterface::n_dofs_at_node(dim, fe_type, elem_type, n_vertices))
1327  edge_vars.insert(edge_vars.end(), v_num);
1328 
1329  if (has_side_nodes)
1330  {
1331  if (dim != 3)
1332  {
1333  if (FEInterface::n_dofs_at_node(dim, fe_type, elem_type, n_vertices))
1334  side_vars.insert(side_vars.end(), v_num);
1335  }
1336  else
1337  // In 3D, not all face nodes always have the same number of
1338  // DoFs! We'll loop over all sides to be safe.
1339  for (unsigned int n = 0; n != n_nodes; ++n)
1340  if (elem->is_face(n))
1341  if (FEInterface::n_dofs_at_node(dim, fe_type,
1342  elem_type, n))
1343  {
1344  side_vars.insert(side_vars.end(), v_num);
1345  break;
1346  }
1347  }
1348 
1349  if (FEInterface::n_dofs_per_elem(dim, fe_type, elem_type) ||
1350  (has_interior_nodes &&
1351  FEInterface::n_dofs_at_node(dim, fe_type, elem_type, n_nodes-1)))
1352  {
1353 #ifdef LIBMESH_ENABLE_AMR
1354  // We may have already just copied constant monomials,
1355  // or we may be about to copy the whole element
1356  if ((f.is_grid_projection() &&
1357  fe_type.family == MONOMIAL &&
1358  fe_type.order == CONSTANT &&
1359  elem->p_level() == 0 &&
1360  elem->refinement_flag() != Elem::JUST_COARSENED &&
1361  elem->p_refinement_flag() != Elem::JUST_COARSENED)
1362  || copy_this_elem
1363  )
1364  continue;
1365 #endif // LIBMESH_ENABLE_AMR
1366 
1367  // We need to project any other variables
1368  if (interiors.empty() || interiors.back() != elem)
1369  interiors.push_back(elem);
1370  }
1371  }
1372 
1373  // We'll use a greedy algorithm in most cases: if another
1374  // element has already claimed some of our DoFs, we'll let it do
1375  // the work.
1376  auto erase_covered_vars = []
1377  (const Node * node, var_set & remaining, const ves_multimap & covered)
1378  {
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);
1384  };
1385 
1386  auto erase_nonhanging_vars = [&extra_hanging_dofs]
1387  (const Node * node, var_set & remaining, const ves_multimap & covered)
1388  {
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);
1395  };
1396 
1397  auto erase_copied_vars = [&copied_nodes, &extra_hanging_dofs]
1398  (const Node * node, bool is_vertex, var_set & remaining)
1399  {
1400  auto copying_range = copied_nodes.equal_range(node);
1401  for (const auto & v_ent : as_range(copying_range))
1402  {
1403  for (const unsigned int var_covered :
1404  v_ent.second.first)
1405  if (is_vertex || !extra_hanging_dofs[var_covered])
1406  remaining.erase(var_covered);
1407 
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);
1412  }
1413  };
1414 
1415  for (unsigned int v=0; v != n_vertices; ++v)
1416  {
1417  const Node * node = elem->node_ptr(v);
1418 
1419  auto remaining_vars = vertex_vars;
1420 
1421  erase_covered_vars(node, remaining_vars, vertices);
1422 
1423  if (remaining_vars.empty())
1424  continue;
1425 
1426  if (!all_extra_hanging_dofs)
1427  {
1428  erase_nonhanging_vars(node, remaining_vars, edges);
1429  if (remaining_vars.empty())
1430  continue;
1431 
1432  erase_nonhanging_vars(node, remaining_vars, sides);
1433  if (remaining_vars.empty())
1434  continue;
1435  }
1436 
1437  erase_copied_vars(node, true, remaining_vars);
1438  if (remaining_vars.empty())
1439  continue;
1440 
1441  if (copy_this_elem)
1442  {
1443  for (auto var : remaining_vars)
1444  {
1445  std::vector<dof_id_type> node_dof_ids;
1446  std::vector<FValue> values;
1447 
1448  f.eval_old_dofs(*elem, v, var, node_dof_ids, values);
1449 
1450  insert_ids(node_dof_ids, values, node->processor_id());
1451  }
1452  copied_nodes[node].first.insert(remaining_vars.begin(),
1453  remaining_vars.end());
1454  this->find_dofs_to_send(*node, *elem, v, remaining_vars);
1455  }
1456  else
1457  vertices.emplace
1458  (node, std::make_tuple(elem, v, std::move(remaining_vars)));
1459  }
1460 
1461  if (has_edge_nodes)
1462  {
1463  for (unsigned int e=0; e != n_edges; ++e)
1464  {
1465  const Node * node = elem->node_ptr(n_vertices+e);
1466 
1467  auto remaining_vars = edge_vars;
1468 
1469  erase_covered_vars(node, remaining_vars, edges);
1470  if (remaining_vars.empty())
1471  continue;
1472 
1473  erase_covered_vars(node, remaining_vars, sides);
1474  if (remaining_vars.empty())
1475  continue;
1476 
1477  if (!all_extra_hanging_dofs)
1478  {
1479  erase_nonhanging_vars(node, remaining_vars, vertices);
1480  if (remaining_vars.empty())
1481  continue;
1482  }
1483 
1484  erase_copied_vars(node, false, remaining_vars);
1485  if (remaining_vars.empty())
1486  continue;
1487 
1488  if (copy_this_elem)
1489  {
1490  for (auto var : remaining_vars)
1491  {
1492  std::vector<dof_id_type> edge_dof_ids;
1493  std::vector<FValue> values;
1494 
1495  f.eval_old_dofs(*elem, n_vertices+e, var, edge_dof_ids, values);
1496 
1497  insert_ids(edge_dof_ids, values, node->processor_id());
1498  }
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);
1502  }
1503  else
1504  edges.emplace
1505  (node, std::make_tuple(elem, e, std::move(remaining_vars)));
1506  }
1507  }
1508 
1509  if (has_side_nodes)
1510  {
1511  for (unsigned int side=0; side != n_sides; ++side)
1512  {
1513  const Node * node = nullptr;
1514  unsigned short node_num = n_vertices+(dim>2)*n_edges+side;
1515  if (dim != 3)
1516  node = elem->node_ptr(node_num);
1517  else
1518  {
1519  // In 3D only some sides may have nodes
1520  for (unsigned int n = 0; n != n_nodes; ++n)
1521  {
1522  if (!elem->is_face(n))
1523  continue;
1524 
1525  if (elem->is_node_on_side(n, side))
1526  {
1527  node_num = n;
1528  node = elem->node_ptr(node_num);
1529  break;
1530  }
1531  }
1532  }
1533 
1534  if (!node)
1535  continue;
1536 
1537  auto remaining_vars = side_vars;
1538 
1539  erase_covered_vars(node, remaining_vars, edges);
1540  if (remaining_vars.empty())
1541  continue;
1542 
1543  erase_covered_vars(node, remaining_vars, sides);
1544  if (remaining_vars.empty())
1545  continue;
1546 
1547  if (!all_extra_hanging_dofs)
1548  {
1549  erase_nonhanging_vars(node, remaining_vars, vertices);
1550  if (remaining_vars.empty())
1551  continue;
1552  }
1553 
1554  erase_copied_vars(node, false, remaining_vars);
1555  if (remaining_vars.empty())
1556  continue;
1557 
1558  if (copy_this_elem)
1559  {
1560  for (auto var : remaining_vars)
1561  {
1562  std::vector<dof_id_type> side_dof_ids;
1563  std::vector<FValue> values;
1564 
1565  f.eval_old_dofs(*elem, node_num, var, side_dof_ids, values);
1566 
1567  insert_ids(side_dof_ids, values, node->processor_id());
1568  }
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);
1572  }
1573  else
1574  sides.emplace
1575  (node, std::make_tuple(elem, side, std::move(remaining_vars)));
1576  }
1577  }
1578 
1579  // Elements with elemental dofs might need those copied too.
1580  if (copy_this_elem)
1581  {
1582  for (auto v_num : this->projector.variables)
1583  {
1584  const Variable & var = system.variable(v_num);
1585  if (!var.active_on_subdomain(elem->subdomain_id()))
1586  continue;
1587  FEType fe_type = var.type();
1588 
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,
1592  elem_dof_ids, Ue);
1593  action.insert(elem_dof_ids, Ue);
1594 
1595  if (has_interior_nodes)
1596  {
1597  std::vector<FValue> Un;
1598  std::vector<dof_id_type> node_dof_ids;
1599 
1600  f.eval_old_dofs(*elem, n_nodes-1, v_num, node_dof_ids, Un);
1601  action.insert(node_dof_ids, Un);
1602  }
1603  }
1604  }
1605  }
1606 }
1607 
1608 
1609 template <typename FFunctor, typename GFunctor,
1610  typename FValue, typename ProjectionAction>
1613 {
1614  auto merge_multimaps = [](ves_multimap & self, const ves_multimap & other_mm)
1615  {
1616  for (const auto & pair : other_mm)
1617  {
1618  const Node * node = pair.first;
1619  auto remaining_vars = std::get<2>(pair.second);
1620 
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);
1626 
1627  if (!remaining_vars.empty())
1628  self.emplace
1629  (node, std::make_tuple(std::get<0>(pair.second),
1630  std::get<1>(pair.second),
1631  std::move(remaining_vars)));
1632  }
1633  };
1634 
1635  merge_multimaps(vertices, other.vertices);
1636  merge_multimaps(edges, other.edges);
1637  merge_multimaps(sides, other.sides);
1638 
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);
1647 
1648  SubFunctor::join(other);
1649 }
1650 
1651 
1652 template <typename FFunctor, typename GFunctor,
1653  typename FValue, typename ProjectionAction>
1655  (const node_range & range)
1656 {
1657  LOG_SCOPE ("project_vertices","GenericProjector");
1658 
1659  const unsigned int sys_num = system.number();
1660 
1661  // Variables with extra hanging dofs can't safely use eval_at_node
1662  // in as many places as variables without can.
1663  std::vector<unsigned short> extra_hanging_dofs;
1664  for (auto v_num : this->projector.variables)
1665  {
1666  if (extra_hanging_dofs.size() <= v_num)
1667  extra_hanging_dofs.resize(v_num+1, false);
1668  extra_hanging_dofs[v_num] =
1669  FEInterface::extra_hanging_dofs(system.variable(v_num).type());
1670  }
1671 
1672  for (const auto & v_pair : range)
1673  {
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);
1678 
1679  context.pre_fe_reinit(system, &elem);
1680 
1681  this->find_dofs_to_send(vertex, elem, n, vertex_vars);
1682 
1683  // Look at all the variables we're supposed to interpolate from
1684  // this element on this vertex
1685  for (const auto & var : vertex_vars)
1686  {
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);
1691 
1692  if (base_fe_type.family == SCALAR)
1693  continue;
1694 
1695  const FEContinuity & cont = this->conts[var];
1696  if (cont == DISCONTINUOUS)
1697  {
1698  libmesh_assert_equal_to(vertex.n_comp(sys_num, var), 0);
1699  }
1700  else if (cont == C_ZERO)
1701  {
1702  libmesh_assert(vertex.n_comp(sys_num, var));
1703  const dof_id_type id = vertex.dof_number(sys_num, var, 0);
1704  // C_ZERO elements have a single nodal value DoF at vertices
1705  const FValue val = f.eval_at_node
1706  (context, var_component, /*dim=*/ 0, // Don't care w/C0
1707  vertex, extra_hanging_dofs[var], system.time);
1708  insert_id(id, val, vertex.processor_id());
1709  }
1710  else if (cont == C_ONE)
1711  {
1712  libmesh_assert(vertex.n_comp(sys_num, var));
1713  const dof_id_type first_id = vertex.dof_number(sys_num, var, 0);
1714 
1715  // C_ONE elements have a single nodal value and dim
1716  // gradient values at vertices, as well as cross
1717  // gradients for HERMITE. We need to have an element in
1718  // hand to figure out dim and to have in case this
1719  // vertex is a new vertex.
1720  const int dim = elem.dim();
1721 #ifndef NDEBUG
1722  // For now all C1 elements at a vertex had better have
1723  // the same dimension. If anyone hits these asserts let
1724  // me know; we could probably support a mixed-dimension
1725  // mesh IFF the 2D elements were all parallel to xy and
1726  // the 1D elements all parallel to x.
1727  for (const auto e_id : (*this->projector.nodes_to_elem)[vertex.id()])
1728  {
1729  const Elem & e = system.get_mesh().elem_ref(e_id);
1730  libmesh_assert_equal_to(dim, e.dim());
1731  }
1732 #endif
1733 #ifdef LIBMESH_ENABLE_AMR
1734  bool is_parent_vertex = false;
1735  if (elem.parent())
1736  {
1737  const int i_am_child =
1738  elem.parent()->which_child_am_i(&elem);
1739  is_parent_vertex =
1740  elem.parent()->is_vertex_on_parent(i_am_child, n);
1741  }
1742 #else
1743  const bool is_parent_vertex = false;
1744 #endif
1745 
1746  // The hermite element vertex shape functions are weird
1747  if (base_fe_type.family == HERMITE)
1748  {
1749  const FValue val =
1750  f.eval_at_node(context,
1751  var_component,
1752  dim,
1753  vertex,
1754  extra_hanging_dofs[var],
1755  system.time);
1756  insert_id(first_id, val, vertex.processor_id());
1757 
1758  VectorValue<FValue> grad =
1759  is_parent_vertex ?
1760  g->eval_at_node(context,
1761  var_component,
1762  dim,
1763  vertex,
1764  extra_hanging_dofs[var],
1765  system.time) :
1766  g->eval_at_point(context,
1767  var_component,
1768  vertex,
1769  system.time);
1770  // x derivative
1771  insert_id(first_id+1, grad(0),
1772  vertex.processor_id());
1773 #if LIBMESH_DIM > 1
1774  if (dim > 1)
1775  {
1776  // We'll finite difference mixed derivatives
1777  Real delta_x = TOLERANCE * elem.hmin();
1778 
1779  Point nxminus = elem.point(n),
1780  nxplus = elem.point(n);
1781  nxminus(0) -= delta_x;
1782  nxplus(0) += delta_x;
1783  VectorValue<FValue> gxminus =
1784  g->eval_at_point(context,
1785  var_component,
1786  nxminus,
1787  system.time);
1788  VectorValue<FValue> gxplus =
1789  g->eval_at_point(context,
1790  var_component,
1791  nxplus,
1792  system.time);
1793  // y derivative
1794  insert_id(first_id+2, grad(1),
1795  vertex.processor_id());
1796  // xy derivative
1797  insert_id(first_id+3,
1798  (gxplus(1) - gxminus(1)) / 2. / delta_x,
1799  vertex.processor_id());
1800 
1801 #if LIBMESH_DIM > 2
1802  if (dim > 2)
1803  {
1804  // z derivative
1805  insert_id(first_id+4, grad(2),
1806  vertex.processor_id());
1807  // xz derivative
1808  insert_id(first_id+5,
1809  (gxplus(2) - gxminus(2)) / 2. / delta_x,
1810  vertex.processor_id());
1811 
1812  // We need new points for yz
1813  Point nyminus = elem.point(n),
1814  nyplus = elem.point(n);
1815  nyminus(1) -= delta_x;
1816  nyplus(1) += delta_x;
1817  VectorValue<FValue> gyminus =
1818  g->eval_at_point(context,
1819  var_component,
1820  nyminus,
1821  system.time);
1822  VectorValue<FValue> gyplus =
1823  g->eval_at_point(context,
1824  var_component,
1825  nyplus,
1826  system.time);
1827  // yz derivative
1828  insert_id(first_id+6,
1829  (gyplus(2) - gyminus(2)) / 2. / delta_x,
1830  vertex.processor_id());
1831  // Getting a 2nd order xyz is more tedious
1832  Point nxmym = elem.point(n),
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;
1844  VectorValue<FValue> gxmym =
1845  g->eval_at_point(context,
1846  var_component,
1847  nxmym,
1848  system.time);
1849  VectorValue<FValue> gxmyp =
1850  g->eval_at_point(context,
1851  var_component,
1852  nxmyp,
1853  system.time);
1854  VectorValue<FValue> gxpym =
1855  g->eval_at_point(context,
1856  var_component,
1857  nxpym,
1858  system.time);
1859  VectorValue<FValue> gxpyp =
1860  g->eval_at_point(context,
1861  var_component,
1862  nxpyp,
1863  system.time);
1864  FValue gxzplus = (gxpyp(2) - gxmyp(2))
1865  / 2. / delta_x;
1866  FValue gxzminus = (gxpym(2) - gxmym(2))
1867  / 2. / delta_x;
1868  // xyz derivative
1869  insert_id(first_id+7,
1870  (gxzplus - gxzminus) / 2. / delta_x,
1871  vertex.processor_id());
1872  }
1873 #endif // LIBMESH_DIM > 2
1874  }
1875 #endif // LIBMESH_DIM > 1
1876  }
1877  else
1878  {
1879  // Currently other C_ONE elements have a single nodal
1880  // value shape function and nodal gradient component
1881  // shape functions
1882  libmesh_assert_equal_to
1884  (dim, base_fe_type, elem.type(),
1885  elem.get_node_index(&vertex)),
1886  (unsigned int)(1 + dim));
1887  const FValue val =
1888  f.eval_at_node(context, var_component, dim,
1889  vertex, extra_hanging_dofs[var],
1890  system.time);
1891  insert_id(first_id, val, vertex.processor_id());
1892  VectorValue<FValue> grad =
1893  is_parent_vertex ?
1894  g->eval_at_node(context, var_component, dim,
1895  vertex, extra_hanging_dofs[var],
1896  system.time) :
1897  g->eval_at_point(context, var_component, vertex,
1898  system.time);
1899  for (int i=0; i!= dim; ++i)
1900  insert_id(first_id + i + 1, grad(i),
1901  vertex.processor_id());
1902  }
1903  }
1904  else
1905  libmesh_error_msg("Unknown continuity " << cont);
1906  }
1907  }
1908 }
1909 
1910 
1911 template <typename FFunctor, typename GFunctor,
1912  typename FValue, typename ProjectionAction>
1914  (const node_range & range)
1915 {
1916  LOG_SCOPE ("project_edges","GenericProjector");
1917 
1918  const unsigned int sys_num = system.number();
1919 
1920  for (const auto & e_pair : range)
1921  {
1922  const Elem & elem = *std::get<0>(e_pair.second);
1923 
1924  // If this is an unchanged element then we already copied all
1925  // its dofs
1926 #ifdef LIBMESH_ENABLE_AMR
1927  if (f.is_grid_projection() &&
1928  (elem.refinement_flag() != Elem::JUST_REFINED &&
1932  continue;
1933 #endif // LIBMESH_ENABLE_AMR
1934 
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);
1938 
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);
1943 
1944  this->find_dofs_to_send(edge_node, elem, node_num, edge_vars);
1945 
1946  // Look at all the variables we're supposed to interpolate from
1947  // this element on this edge
1948  for (const auto & var : edge_vars)
1949  {
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);
1954 
1955  if (base_fe_type.family == SCALAR)
1956  continue;
1957 
1958  FEType fe_type = base_fe_type;
1959 
1960  // This may be a p refined element
1961  fe_type.order =
1962  libMesh::Order (fe_type.order + elem.p_level());
1963 
1964  // If this is a Lagrange element with DoFs on edges then by
1965  // convention we interpolate at the node rather than project
1966  // along the edge.
1967  if (fe_type.family == LAGRANGE)
1968  {
1969  if (fe_type.order > 1)
1970  {
1971  const dof_id_type dof_id =
1972  edge_node.dof_number(sys_num, var, 0);
1973  const processor_id_type pid =
1974  edge_node.processor_id();
1975  FValue fval = f.eval_at_point
1976  (context, var_component, edge_node, system.time);
1977  insert_id(dof_id, fval, pid);
1978  }
1979  continue;
1980  }
1981 
1982 #ifdef LIBMESH_ENABLE_AMR
1983  // If this is a low order monomial element which has merely
1984  // been h refined then we already copied all its dofs
1985  if (fe_type.family == MONOMIAL &&
1986  fe_type.order == CONSTANT &&
1989  continue;
1990 #endif // LIBMESH_ENABLE_AMR
1991 
1992  // FIXME: Need to generalize this to vector-valued elements. [PB]
1993  FEBase * fe = nullptr;
1994  context.get_element_fe( var, fe, dim );
1995  FEBase * edge_fe = nullptr;
1996  context.get_edge_fe( var, edge_fe );
1997 
1998  // If we're JUST_COARSENED we'll need a custom
1999  // evaluation, not just the standard edge FE
2000  const FEBase & proj_fe =
2001 #ifdef LIBMESH_ENABLE_AMR
2003  *fe :
2004 #endif
2005  *edge_fe;
2006 
2007 #ifdef LIBMESH_ENABLE_AMR
2008  if (elem.refinement_flag() == Elem::JUST_COARSENED)
2009  {
2010  std::vector<Point> fine_points;
2011 
2012  std::unique_ptr<FEBase> fine_fe
2013  (FEBase::build (dim, base_fe_type));
2014 
2015  std::unique_ptr<QBase> qrule
2016  (base_fe_type.default_quadrature_rule(1));
2017  fine_fe->attach_quadrature_rule(qrule.get());
2018 
2019  const std::vector<Point> & child_xyz =
2020  fine_fe->get_xyz();
2021 
2022  for (unsigned int c = 0, nc = elem.n_children();
2023  c != nc; ++c)
2024  {
2025  if (!elem.is_child_on_edge(c, context.edge))
2026  continue;
2027 
2028  fine_fe->edge_reinit(elem.child_ptr(c), context.edge);
2029  fine_points.insert(fine_points.end(),
2030  child_xyz.begin(),
2031  child_xyz.end());
2032  }
2033 
2034  std::vector<Point> fine_qp;
2035  FEMap::inverse_map (dim, &elem, fine_points, fine_qp);
2036 
2037  context.elem_fe_reinit(&fine_qp);
2038  }
2039  else
2040 #endif // LIBMESH_ENABLE_AMR
2041  context.edge_fe_reinit();
2042 
2043  const std::vector<dof_id_type> & dof_indices =
2044  context.get_dof_indices(var);
2045 
2046  std::vector<unsigned int> edge_dofs;
2047  FEInterface::dofs_on_edge(&elem, dim, base_fe_type,
2048  context.edge, edge_dofs);
2049 
2050  this->construct_projection
2051  (dof_indices, edge_dofs, var_component,
2052  &edge_node, proj_fe);
2053  }
2054  }
2055 }
2056 
2057 
2058 template <typename FFunctor, typename GFunctor,
2059  typename FValue, typename ProjectionAction>
2061  (const node_range & range)
2062 {
2063  LOG_SCOPE ("project_sides","GenericProjector");
2064 
2065  const unsigned int sys_num = system.number();
2066 
2067  for (const auto & s_pair : range)
2068  {
2069  const Elem & elem = *std::get<0>(s_pair.second);
2070 
2071  // If this is an unchanged element then we already copied all
2072  // its dofs
2073 #ifdef LIBMESH_ENABLE_AMR
2074  if (f.is_grid_projection() &&
2075  (elem.refinement_flag() != Elem::JUST_REFINED &&
2079  continue;
2080 #endif // LIBMESH_ENABLE_AMR
2081 
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);
2085 
2086  const unsigned int side_num = std::get<1>(s_pair.second);
2087  unsigned short node_num = elem.n_vertices()+side_num;
2088  // In 3D only some sides may have nodes
2089  if (dim == 3)
2090  for (auto n : IntRange<unsigned int>(0, elem.n_nodes()))
2091  {
2092  if (!elem.is_face(n))
2093  continue;
2094 
2095  if (elem.is_node_on_side(n, side_num))
2096  {
2097  node_num = n;
2098  break;
2099  }
2100  }
2101 
2102  context.side = side_num;
2103  context.pre_fe_reinit(system, &elem);
2104 
2105  this->find_dofs_to_send(side_node, elem, node_num, side_vars);
2106 
2107  // Look at all the variables we're supposed to interpolate from
2108  // this element on this side
2109  for (const auto & var : side_vars)
2110  {
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);
2115 
2116  if (base_fe_type.family == SCALAR)
2117  continue;
2118 
2119  FEType fe_type = base_fe_type;
2120 
2121  // This may be a p refined element
2122  fe_type.order =
2123  libMesh::Order (fe_type.order + elem.p_level());
2124 
2125  // If this is a Lagrange element with DoFs on sides then by
2126  // convention we interpolate at the node rather than project
2127  // along the side.
2128  if (fe_type.family == LAGRANGE)
2129  {
2130  if (fe_type.order > 1)
2131  {
2132  const dof_id_type dof_id =
2133  side_node.dof_number(sys_num, var, 0);
2134  const processor_id_type pid =
2135  side_node.processor_id();
2136  FValue fval = f.eval_at_point
2137  (context, var_component, side_node, system.time);
2138  insert_id(dof_id, fval, pid);
2139  }
2140  continue;
2141  }
2142 
2143 #ifdef LIBMESH_ENABLE_AMR
2144  // If this is a low order monomial element which has merely
2145  // been h refined then we already copied all its dofs
2146  if (fe_type.family == MONOMIAL &&
2147  fe_type.order == CONSTANT &&
2150  continue;
2151 #endif // LIBMESH_ENABLE_AMR
2152 
2153  // FIXME: Need to generalize this to vector-valued elements. [PB]
2154  FEBase * fe = nullptr;
2155  context.get_element_fe( var, fe, dim );
2156  FEBase * side_fe = nullptr;
2157  context.get_side_fe( var, side_fe );
2158 
2159  // If we're JUST_COARSENED we'll need a custom
2160  // evaluation, not just the standard side FE
2161  const FEBase & proj_fe =
2162 #ifdef LIBMESH_ENABLE_AMR
2164  *fe :
2165 #endif
2166  *side_fe;
2167 
2168 #ifdef LIBMESH_ENABLE_AMR
2169  if (elem.refinement_flag() == Elem::JUST_COARSENED)
2170  {
2171  std::vector<Point> fine_points;
2172 
2173  std::unique_ptr<FEBase> fine_fe
2174  (FEBase::build (dim, base_fe_type));
2175 
2176  std::unique_ptr<QBase> qrule
2177  (base_fe_type.default_quadrature_rule(1));
2178  fine_fe->attach_quadrature_rule(qrule.get());
2179 
2180  const std::vector<Point> & child_xyz =
2181  fine_fe->get_xyz();
2182 
2183  for (unsigned int c = 0, nc = elem.n_children();
2184  c != nc; ++c)
2185  {
2186  if (!elem.is_child_on_side(c, context.side))
2187  continue;
2188 
2189  fine_fe->reinit(elem.child_ptr(c), context.side);
2190  fine_points.insert(fine_points.end(),
2191  child_xyz.begin(),
2192  child_xyz.end());
2193  }
2194 
2195  std::vector<Point> fine_qp;
2196  FEMap::inverse_map (dim, &elem, fine_points, fine_qp);
2197 
2198  context.elem_fe_reinit(&fine_qp);
2199  }
2200  else
2201 #endif // LIBMESH_ENABLE_AMR
2202  context.side_fe_reinit();
2203 
2204  const std::vector<dof_id_type> & dof_indices =
2205  context.get_dof_indices(var);
2206 
2207  std::vector<unsigned int> side_dofs;
2208  FEInterface::dofs_on_side(&elem, dim, base_fe_type,
2209  context.side, side_dofs);
2210 
2211  this->construct_projection
2212  (dof_indices, side_dofs, var_component,
2213  &side_node, proj_fe);
2214  }
2215  }
2216 }
2217 
2218 
2219 template <typename FFunctor, typename GFunctor,
2220  typename FValue, typename ProjectionAction>
2222  (const interior_range & range)
2223 {
2224  LOG_SCOPE ("project_interiors","GenericProjector");
2225 
2226  const unsigned int sys_num = system.number();
2227 
2228  // Iterate over all dof-bearing element interiors in the range
2229  for (const auto & elem : range)
2230  {
2231  unsigned char dim = cast_int<unsigned char>(elem->dim());
2232 
2233  context.pre_fe_reinit(system, elem);
2234 
2235  // Loop over all the variables we've been requested to project, to
2236  // do the projection
2237  for (const auto & var : this->projector.variables)
2238  {
2239  const Variable & variable = system.variable(var);
2240 
2241  if (!variable.active_on_subdomain(elem->subdomain_id()))
2242  continue;
2243 
2244  const FEType & base_fe_type = variable.type();
2245 
2246  if (base_fe_type.family == SCALAR)
2247  continue;
2248 
2249  FEBase * fe = nullptr;
2250  context.get_element_fe( var, fe, dim );
2251 
2252  FEType fe_type = base_fe_type;
2253 
2254  // This may be a p refined element
2255  fe_type.order =
2256  libMesh::Order (fe_type.order + elem->p_level());
2257 
2258  const unsigned int var_component =
2259  system.variable_scalar_number(var, 0);
2260 
2261  // If this is a Lagrange element with interior DoFs then by
2262  // convention we interpolate at the interior node rather
2263  // than project along the interior.
2264  if (fe_type.family == LAGRANGE)
2265  {
2266  if (fe_type.order > 1)
2267  {
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();
2273 
2274  // < vs != is important here for HEX20, QUAD8!
2275  for (unsigned int n = first_interior_node; n < n_nodes; ++n)
2276  {
2277  const Node & interior_node = elem->node_ref(n);
2278  const dof_id_type dof_id =
2279  interior_node.dof_number(sys_num, var, 0);
2280  const processor_id_type pid =
2281  interior_node.processor_id();
2282  FValue fval = f.eval_at_point
2283  (context, var_component, interior_node, system.time);
2284  insert_id(dof_id, fval, pid);
2285  }
2286  }
2287  continue;
2288  }
2289 
2290 #ifdef LIBMESH_ENABLE_AMR
2291  if (elem->refinement_flag() == Elem::JUST_COARSENED)
2292  {
2293  std::vector<Point> fine_points;
2294 
2295  std::unique_ptr<FEBase> fine_fe
2296  (FEBase::build (dim, base_fe_type));
2297 
2298  std::unique_ptr<QBase> qrule
2299  (base_fe_type.default_quadrature_rule(dim));
2300  fine_fe->attach_quadrature_rule(qrule.get());
2301 
2302  const std::vector<Point> & child_xyz =
2303  fine_fe->get_xyz();
2304 
2305  for (auto & child : elem->child_ref_range())
2306  {
2307  fine_fe->reinit(&child);
2308  fine_points.insert(fine_points.end(),
2309  child_xyz.begin(),
2310  child_xyz.end());
2311  }
2312 
2313  std::vector<Point> fine_qp;
2314  FEMap::inverse_map (dim, elem, fine_points, fine_qp);
2315 
2316  context.elem_fe_reinit(&fine_qp);
2317  }
2318  else
2319 #endif // LIBMESH_ENABLE_AMR
2320  context.elem_fe_reinit();
2321 
2322  const std::vector<dof_id_type> & dof_indices =
2323  context.get_dof_indices(var);
2324 
2325  const unsigned int n_dofs =
2326  cast_int<unsigned int>(dof_indices.size());
2327 
2328  std::vector<unsigned int> all_dofs(n_dofs);
2329  std::iota(all_dofs.begin(), all_dofs.end(), 0);
2330 
2331  this->construct_projection
2332  (dof_indices, all_dofs, var_component,
2333  nullptr, *fe);
2334  } // end variables loop
2335  } // end elements loop
2336 }
2337 
2338 
2339 
2340 template <typename FFunctor, typename GFunctor,
2341  typename FValue, typename ProjectionAction>
2342 void
2343 GenericProjector<FFunctor, GFunctor, FValue,
2344  ProjectionAction>::SubFunctor::find_dofs_to_send
2345  (const Node & node, const Elem & elem, unsigned short node_num, const var_set & vars)
2346 {
2347  libmesh_assert (&node == elem.node_ptr(node_num));
2348 
2349  // Any ghosted node in our range might have an owner who needs our
2350  // data
2351  const processor_id_type owner = node.processor_id();
2352  if (owner != system.processor_id())
2353  {
2354  const MeshBase & mesh = system.get_mesh();
2355  const DofMap & dof_map = system.get_dof_map();
2356 
2357  // But let's check and see if we can be certain the owner can
2358  // compute any or all of its own dof coefficients on that node.
2359  std::vector<dof_id_type> node_dof_ids, patch_dof_ids;
2360  for (const auto & var : vars)
2361  {
2362  const Variable & variable = system.variable(var);
2363 
2364  if (!variable.active_on_subdomain(elem.subdomain_id()))
2365  continue;
2366 
2367  dof_map.dof_indices(elem, node_num, node_dof_ids, var);
2368  }
2369  libmesh_assert(std::is_sorted(node_dof_ids.begin(),
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)
2374  {
2375  const Elem & patch_elem = mesh.elem_ref(elem_id);
2376  if (!patch_elem.active() || owner != patch_elem.processor_id())
2377  continue;
2378  dof_map.dof_indices(&patch_elem, patch_dof_ids);
2379  std::sort(patch_dof_ids.begin(), patch_dof_ids.end());
2380 
2381  // Remove any node_dof_ids that we see can be calculated on
2382  // this element
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())
2389  break;
2390  }
2391 
2392  // Give ids_to_push default invalid pid: not yet computed
2393  for (auto id : node_dof_ids)
2394  new_ids_to_push[id].second = DofObject::invalid_processor_id;
2395  }
2396 }
2397 
2398 
2399 
2400 template <typename FFunctor, typename GFunctor,
2401  typename FValue, typename ProjectionAction>
2402 void
2403 GenericProjector<FFunctor, GFunctor, FValue,
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
2407 {
2408  // See if we calculated any ids that need to be pushed; get them
2409  // ready to push.
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)
2415  {
2416  processor_id_type pid = id_val_pid.second.second;
2418  {
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));
2421  }
2422  }
2423 
2424  // If and when we get ids pushed to us, act on them if we have
2425  // corresponding values or save them if not
2426  auto ids_action_functor =
2427  [&action, &received_dof_ids, &received_dof_values]
2428  (processor_id_type pid,
2429  const std::vector<dof_id_type> & data)
2430  {
2431  auto iter = received_dof_values.find(pid);
2432  if (iter == received_dof_values.end())
2433  {
2434  libmesh_assert(received_dof_ids.find(pid) ==
2435  received_dof_ids.end());
2436  received_dof_ids[pid] = data;
2437  }
2438  else
2439  {
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)
2444  {
2445  FValue val;
2446  convert_from_receive(vals[i], val);
2447  action.insert(data[i], val);
2448  }
2449  received_dof_values.erase(iter);
2450  }
2451  };
2452 
2453  // If and when we get values pushed to us, act on them if we have
2454  // corresponding ids or save them if not
2455  auto values_action_functor =
2456  [&action, &received_dof_ids, &received_dof_values]
2457  (processor_id_type pid,
2458  const std::vector<typename TypeToSend<FValue>::type> & data)
2459  {
2460  auto iter = received_dof_ids.find(pid);
2461  if (iter == received_dof_ids.end())
2462  {
2463  // We get no more than 1 values vector from anywhere
2464  libmesh_assert(received_dof_values.find(pid) ==
2465  received_dof_values.end());
2466  received_dof_values[pid] = data;
2467  }
2468  else
2469  {
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)
2474  {
2475  FValue val;
2476  convert_from_receive(data[i], val);
2477  action.insert(ids[i], val);
2478  }
2479  received_dof_ids.erase(iter);
2480  }
2481  };
2482 
2483  Parallel::push_parallel_vector_data
2484  (system.comm(), pushed_dof_ids, ids_action_functor);
2485 
2486  Parallel::push_parallel_vector_data
2487  (system.comm(), pushed_dof_values, values_action_functor);
2488 
2489  // At this point we shouldn't have any unprocessed data left
2490  libmesh_assert(received_dof_ids.empty());
2491  libmesh_assert(received_dof_values.empty());
2492 
2493 }
2494 
2495 
2496 
2497 template <typename FFunctor, typename GFunctor,
2498  typename FValue, typename ProjectionAction>
2499 void
2501  (const std::vector<dof_id_type> & dof_indices_var,
2502  const std::vector<unsigned int> & involved_dofs,
2503  unsigned int var_component,
2504  const Node * node,
2505  const FEBase & fe)
2506 {
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();
2511  const FEContinuity cont = fe.get_continuity();
2512  const std::unordered_map<dof_id_type, FValue> & ids_to_save =
2513  this->projector.ids_to_save;
2514 
2515  if (cont == C_ONE)
2516  dphi = &(fe.get_dphi());
2517 
2518  const unsigned int n_involved_dofs =
2519  cast_int<unsigned int>(involved_dofs.size());
2520 
2521  std::vector<dof_id_type> free_dof_ids;
2522  DenseVector<FValue> Uinvolved(n_involved_dofs);
2523  std::vector<char> dof_is_fixed(n_involved_dofs, false); // bools
2524 
2525  for (auto i : IntRange<unsigned int>(0, n_involved_dofs))
2526  {
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);
2531  else
2532  {
2533  dof_is_fixed[i] = true;
2534  Uinvolved(i) = iter->second;
2535  }
2536  }
2537 
2538  const unsigned int free_dofs = free_dof_ids.size();
2539 
2540  // There may be nothing to project
2541  if (!free_dofs)
2542  return;
2543 
2544  // The element matrix and RHS for projections.
2545  // Note that Ke is always real-valued, whereas
2546  // Fe may be complex valued if complex number
2547  // support is enabled
2548  DenseMatrix<Real> Ke(free_dofs, free_dofs);
2549  DenseVector<FValue> Fe(free_dofs);
2550  // The new degree of freedom coefficients to solve for
2551  DenseVector<FValue> Ufree(free_dofs);
2552 
2553  const unsigned int n_qp =
2554  cast_int<unsigned int>(xyz_values.size());
2555 
2556  // Loop over the quadrature points
2557  for (unsigned int qp=0; qp<n_qp; qp++)
2558  {
2559  // solution at the quadrature point
2560  FValue fineval = f.eval_at_point(context,
2561  var_component,
2562  xyz_values[qp],
2563  system.time);
2564  // solution grad at the quadrature point
2565  VectorValue<FValue> finegrad;
2566  if (cont == C_ONE)
2567  finegrad = g->eval_at_point(context,
2568  var_component,
2569  xyz_values[qp],
2570  system.time);
2571 
2572  // Form edge projection matrix
2573  for (unsigned int sidei=0, freei=0;
2574  sidei != n_involved_dofs; ++sidei)
2575  {
2576  unsigned int i = involved_dofs[sidei];
2577  // fixed DoFs aren't test functions
2578  if (dof_is_fixed[sidei])
2579  continue;
2580  for (unsigned int sidej=0, freej=0;
2581  sidej != n_involved_dofs; ++sidej)
2582  {
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);
2587  else
2588  Ke(freei,freej) += phi[i][qp] *
2589  phi[j][qp] * JxW[qp];
2590  if (cont == C_ONE)
2591  {
2592  if (dof_is_fixed[sidej])
2593  Fe(freei) -= ( (*dphi)[i][qp] *
2594  (*dphi)[j][qp] ) *
2595  JxW[qp] * Uinvolved(sidej);
2596  else
2597  Ke(freei,freej) += ( (*dphi)[i][qp] *
2598  (*dphi)[j][qp] )
2599  * JxW[qp];
2600  }
2601  if (!dof_is_fixed[sidej])
2602  freej++;
2603  }
2604  Fe(freei) += phi[i][qp] * fineval * JxW[qp];
2605  if (cont == C_ONE)
2606  Fe(freei) += (finegrad * (*dphi)[i][qp] ) *
2607  JxW[qp];
2608  freei++;
2609  }
2610  }
2611 
2612  Ke.cholesky_solve(Fe, Ufree);
2613 
2614  // Transfer new edge solutions to element
2615  const processor_id_type pid = node ?
2617  insert_ids(free_dof_ids, Ufree.get_values(), pid);
2618 }
2619 
2620 
2621 } // namespace libMesh
2622 
2623 #endif // GENERIC_PROJECTOR_H
libMesh::C_ONE
Definition: enum_fe_family.h:80
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::OldSolutionValue::old_solution
const NumericVector< Number > & old_solution
Definition: generic_projector.h:783
libMesh::GenericProjector::g_was_copied
bool g_was_copied
Definition: generic_projector.h:89
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::GenericProjector::SubFunctor::new_ids_to_save
std::unordered_map< dof_id_type, FValue > new_ids_to_save
Definition: generic_projector.h:169
libMesh::System::n_vars
unsigned int n_vars() const
Definition: system.h:2155
libMesh::Elem::child_ptr
const Elem * child_ptr(unsigned int i) const
Definition: elem.h:2567
libMesh::Elem::is_node_on_side
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const =0
libMesh::GenericProjector::ProjectVertices::ProjectVertices
ProjectVertices(ProjectVertices &p_v, Threads::split)
Definition: generic_projector.h:270
libMesh::OldSolutionValue
The OldSolutionValue input functor class can be used with GenericProjector to read values from a solu...
Definition: generic_projector.h:655
libMesh::OldSolutionBase::check_old_context
bool check_old_context(const FEMContext &c, const Point &p)
Definition: generic_projector.h:565
libMesh::GenericProjector::SortAndCopy::operator()
void operator()(const ConstElemRange &range)
Definition: generic_projector.h:1176
libMesh::GenericProjector::ProjectVertices::ProjectVertices
ProjectVertices(GenericProjector &p)
Definition: generic_projector.h:268
libMesh::GenericProjector::ProjectInteriors
Definition: generic_projector.h:329
libMesh::GenericProjector::map_was_created
bool map_was_created
Definition: generic_projector.h:89
libMesh::OldSolutionValue::eval_at_point
Output eval_at_point(const FEMContext &c, unsigned int i, const Point &p, Real=0.)
Definition: generic_projector.h:683
libMesh::FEInterface::extra_hanging_dofs
static bool extra_hanging_dofs(const FEType &fe_t)
Definition: fe_interface.C:1656
libMesh::FEAbstract::get_xyz
const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:237
libMesh::Variable::system
System * system() const
Definition: variable.h:92
libMesh::OldSolutionBase::is_grid_projection
bool is_grid_projection()
Definition: generic_projector.h:529
libMesh::GenericProjector::SortAndCopy::sides
ves_multimap sides
Definition: generic_projector.h:262
libMesh::DofObject::n_systems
unsigned int n_systems() const
Definition: dof_object.h:861
libMesh::GenericProjector::ProjectSides
Definition: generic_projector.h:306
libMesh::GenericProjector::interior_range
StoredRange< std::vector< const Elem * >::const_iterator, const Elem * > interior_range
Definition: generic_projector.h:327
libMesh::GenericProjector::GenericProjector
GenericProjector(const System &system_in, const FFunctor &f_in, const GFunctor *g_in, const ProjectionAction &act_in, const std::vector< unsigned int > &variables_in, std::unordered_map< dof_id_type, std::vector< dof_id_type >> *nodes_to_elem_in=nullptr)
Definition: generic_projector.h:96
libMesh::Elem::JUST_REFINED
Definition: elem.h:1172
libMesh::GenericProjector::ProjectSides::ProjectSides
ProjectSides(GenericProjector &p)
Definition: generic_projector.h:307
libMesh::GenericProjector::var_set
std::set< unsigned int > var_set
Definition: generic_projector.h:149
libMesh::FEType::family
FEFamily family
The type of finite element.
Definition: fe_type.h:203
libMesh::DofObject::old_dof_object
DofObject * old_dof_object
This object on the last mesh.
Definition: dof_object.h:81
libMesh::DofMap::dof_indices
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
Definition: dof_map.C:1967
libMesh::Elem::child_ref_range
SimpleRange< ChildRefIter > child_ref_range()
Returns a range with all children of a parent element, usable in range-based for loops.
Definition: elem.h:1839
libMesh::GenericProjector::ProjectSides::ProjectSides
ProjectSides(ProjectSides &p_s, Threads::split)
Definition: generic_projector.h:309
libMesh::Elem::n_nodes
virtual unsigned int n_nodes() const =0
libMesh::GenericProjector::SortAndCopy::vertices
ves_multimap vertices
Definition: generic_projector.h:262
libMesh::FEMFunctionBase
FEMFunctionBase is a base class from which users can derive in order to define "function-like" object...
Definition: dirichlet_boundaries.h:43
libMesh::GenericProjector::SubFunctor::SubFunctor
SubFunctor(GenericProjector &p)
Definition: generic_projector.h:1031
libMesh::GenericProjector::master_g
const GFunctor * master_g
Definition: generic_projector.h:88
libMesh::HERMITE
Definition: enum_fe_family.h:54
libMesh::GenericProjector::master_f
const FFunctor & master_f
Definition: generic_projector.h:87
libMesh::FEInterface::n_dofs_at_node
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
Definition: fe_interface.C:503
libMesh::GenericProjector::SubFunctor::insert_id
void insert_id(dof_id_type id, const FValue &val, processor_id_type pid)
Definition: generic_projector.h:1117
libMesh::VectorSetAction::target_vector
NumericVector< Val > & target_vector
Definition: generic_projector.h:367
libMesh::MeshBase::elem_ref
virtual const Elem & elem_ref(const dof_id_type i) const
Definition: mesh_base.h:521
libMesh::GenericProjector::SortAndCopy::ves_multimap
std::unordered_multimap< const Node *, std::tuple< const Elem *, unsigned short, var_set > > ves_multimap
Definition: generic_projector.h:260
libMesh::NumericVector::last_local_index
virtual numeric_index_type last_local_index() const =0
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::FEMFunctionWrapper::init_context
void init_context(FEMContext &c)
Definition: generic_projector.h:425
libMesh::Elem::dim
virtual unsigned short dim() const =0
libMesh::GenericProjector::master_action
const ProjectionAction & master_action
Definition: generic_projector.h:90
libMesh::FEMFunctionWrapper::_f
std::unique_ptr< FEMFunctionBase< Output > > _f
Definition: generic_projector.h:459
libMesh::FEGenericBase
This class forms the foundation from which generic finite elements may be derived.
Definition: exact_error_estimator.h:39
libMesh::GenericProjector::GenericProjector
GenericProjector(const GenericProjector &in)
Definition: generic_projector.h:120
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::FEAbstract::get_JxW
const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:244
libMesh::OldSolutionBase::out_of_elem_tol
static const Real out_of_elem_tol
Definition: generic_projector.h:638
libMesh::FEMContext::get_elem
const Elem & get_elem() const
Accessor for current Elem object.
Definition: fem_context.h:896
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::DenseMatrix< Real >
libMesh::Threads::split
Dummy "splitting object" used to distinguish splitting constructors from copy constructors.
Definition: threads_none.h:63
libMesh::TypeToSend
For ease of communication, we allow users to translate their own value types to a more easily computa...
Definition: generic_projector.h:56
libMesh::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::OldSolutionBase::init_context
void init_context(FEMContext &c)
Definition: generic_projector.h:510
libMesh::Elem::RefinementState
RefinementState
Enumeration of possible element refinement states.
Definition: elem.h:1169
libMesh::Variable::type
const FEType & type() const
Definition: variable.h:119
libMesh::GenericProjector::ProjectSides::operator()
void operator()(const node_range &range)
Definition: generic_projector.h:2061
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::GenericProjector::SubProjector
Definition: generic_projector.h:214
libMesh::MeshTools::build_nodes_to_elem_map
void build_nodes_to_elem_map(const MeshBase &mesh, std::vector< std::vector< dof_id_type >> &nodes_to_elem_map)
After calling this function the input vector nodes_to_elem_map will contain the node to element conne...
Definition: mesh_tools.C:248
libMesh::FEGenericBase::get_dphi
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:214
libMesh::GenericProjector::SubFunctor::new_ids_to_push
std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type > > new_ids_to_push
Definition: generic_projector.h:165
libMesh::System::variable
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Definition: system.h:2183
libMesh::DofObject::processor_id
processor_id_type processor_id() const
Definition: dof_object.h:829
libMesh::Elem::is_child_on_edge
virtual bool is_child_on_edge(const unsigned int c, const unsigned int e) const
Definition: elem.C:1646
libMesh::Elem::active
bool active() const
Definition: elem.h:2345
parallel_sync.h
libMesh::GenericProjector::SortAndCopy::SortAndCopy
SortAndCopy(GenericProjector &p)
Definition: generic_projector.h:241
libMesh::GenericProjector::SubProjector::SubProjector
SubProjector(GenericProjector &p)
Definition: generic_projector.h:1096
libMesh::OldSolutionValue::OldSolutionValue
OldSolutionValue(const libMesh::System &sys_in, const NumericVector< Number > &old_sol)
Definition: generic_projector.h:658
libMesh::GenericProjector::ProjectEdges::operator()
void operator()(const node_range &range)
Definition: generic_projector.h:1914
libMesh::Parallel::Utils::is_sorted
bool is_sorted(const std::vector< KeyType > &v)
Definition: parallel_conversion_utils.h:52
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::Elem::is_vertex_on_parent
virtual bool is_vertex_on_parent(unsigned int c, unsigned int n) const
Definition: elem.C:2671
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
libMesh::GenericProjector::project
void project(const ConstElemRange &range)
Function definitions.
Definition: generic_projector.h:951
libMesh::VectorValue
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
Definition: exact_solution.h:54
libMesh::OldSolutionBase::component_to_var
std::vector< unsigned int > component_to_var
Definition: generic_projector.h:636
libMesh::NumericVector< Val >
libMesh::FEAbstract::get_continuity
virtual FEContinuity get_continuity() const =0
libMesh::FEMFunctionWrapper::eval_old_dofs
void eval_old_dofs(const Elem &, unsigned int, unsigned int, std::vector< dof_id_type > &, std::vector< Output > &)
Definition: generic_projector.h:443
libMesh::FEMFunctionWrapper::FEMFunctionWrapper
FEMFunctionWrapper(const FEMFunctionBase< Output > &f)
Definition: generic_projector.h:420
libMesh::DofObject::dof_number
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:956
libMesh::GenericProjector::SubFunctor::conts
std::vector< FEContinuity > conts
Definition: generic_projector.h:201
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::OldSolutionValue::eval_old_dofs
void eval_old_dofs(const Elem &elem, unsigned int node_num, unsigned int var_num, std::vector< dof_id_type > &indices, std::vector< Output > &values)
Definition: generic_projector.h:702
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::FEMFunctionWrapper::eval_at_node
Output eval_at_node(const FEMContext &c, unsigned int i, unsigned int, const Node &n, bool, const Real time)
Definition: generic_projector.h:427
libMesh::OldSolutionValue::eval_at_node
Output eval_at_node(const FEMContext &c, unsigned int i, unsigned int elem_dim, const Node &n, bool, Real=0.)
libMesh::GenericProjector::SubFunctor::insert_ids
void insert_ids(const std::vector< dof_id_type > &ids, const std::vector< FValue > &vals, processor_id_type pid)
Definition: generic_projector.h:1138
libMesh::GenericProjector::SubFunctor::find_dofs_to_send
void find_dofs_to_send(const Node &node, const Elem &elem, unsigned short node_num, const var_set &vars)
Definition: generic_projector.h:2345
libMesh::OldSolutionBase::OldSolutionBase
OldSolutionBase(const libMesh::System &sys_in)
Definition: generic_projector.h:482
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::GenericProjector::done_saving_ids
bool done_saving_ids
Definition: generic_projector.h:93
libMesh::GenericProjector::ProjectEdges
Definition: generic_projector.h:287
libMesh::GenericProjector::ProjectInteriors::ProjectInteriors
ProjectInteriors(GenericProjector &p)
Definition: generic_projector.h:330
libMesh::FEMap::inverse_map
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe_map.C:1622
libMesh::GenericProjector::SubProjector::construct_projection
void construct_projection(const std::vector< dof_id_type > &dof_indices_var, const std::vector< unsigned int > &involved_dofs, unsigned int var_component, const Node *node, const FEBase &fe)
Definition: generic_projector.h:2501
libMesh::FEMFunctionWrapper::eval_old_dofs
void eval_old_dofs(const Elem &, const FEType &, unsigned int, unsigned int, std::vector< dof_id_type > &, std::vector< Output > &)
Definition: generic_projector.h:450
libMesh::FEInterface::n_dofs_per_elem
static unsigned int n_dofs_per_elem(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:534
libMesh::FEMFunctionWrapper
The FEMFunctionWrapper input functor class can be used with a GenericProjector to read values from an...
Definition: generic_projector.h:417
libMesh::GenericProjector::SubFunctor::context
FEMContext context
Definition: generic_projector.h:198
libMesh::OldSolutionBase::OldSolutionBase
OldSolutionBase(const OldSolutionBase &in)
Definition: generic_projector.h:498
libMesh::Elem::p_refinement_flag
RefinementState p_refinement_flag() const
Definition: elem.h:2630
libMesh::GenericProjector::SortAndCopy::edges
ves_multimap edges
Definition: generic_projector.h:262
libMesh::OldSolutionBase::check_old_context
void check_old_context(const FEMContext &c)
Definition: generic_projector.h:532
libMesh::DofObject::n_comp
unsigned int n_comp(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:926
libMesh::GenericProjector::SubFunctor
Definition: generic_projector.h:153
libMesh::GenericProjector::node_projection
std::pair< const Node *, std::tuple< const Elem *, unsigned short, var_set > > node_projection
Definition: generic_projector.h:208
libMesh::OldSolutionBase
The OldSolutionBase input functor abstract base class is the root of the OldSolutionValue and OldSolu...
Definition: generic_projector.h:479
libMesh::Utility::iota
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota is a duplication of the SGI STL extension std::iota.
Definition: utility.h:105
libMesh::Threads::spin_mtx
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:29
libMesh::System::get_mesh
const MeshBase & get_mesh() const
Definition: system.h:2083
libMesh::VectorSetAction
The VectorSetAction output functor class can be used with a GenericProjector to set projection values...
Definition: generic_projector.h:364
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libMesh::DofObject::invalid_id
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:421
libMesh::GenericProjector::SortAndCopy::SortAndCopy
SortAndCopy(SortAndCopy &other, Threads::split)
Definition: generic_projector.h:243
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::FEInterface::dofs_on_edge
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)
Fills the vector di with the local degree of freedom indices associated with edge e of element elem A...
Definition: fe_interface.C:566
libMesh::Variable
This class defines the notion of a variable in the system.
Definition: variable.h:49
libMesh::processor_id_type
uint8_t processor_id_type
Definition: id_types.h:104
libMesh::DISCONTINUOUS
Definition: enum_fe_family.h:78
libMesh::GenericProjector::ProjectInteriors::operator()
void operator()(const interior_range &range)
Definition: generic_projector.h:2222
libMesh::DofObject::var_to_vg_and_offset
std::pair< unsigned int, unsigned int > var_to_vg_and_offset(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:1112
libMesh::C_ZERO
Definition: enum_fe_family.h:79
libMesh::Elem::n_vertices
virtual unsigned int n_vertices() const =0
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::Threads::parallel_reduce
void parallel_reduce(const Range &range, Body &body)
Execute the provided reduction operation in parallel on the specified range.
Definition: threads_none.h:101
libMesh::CONSTANT
Definition: enum_order.h:41
libMesh::GenericProjector
The GenericProjector class implements the core of other projection operations, using two input functo...
Definition: generic_projector.h:80
libMesh::Elem::DO_NOTHING
Definition: elem.h:1170
libMesh::VectorSetAction::insert
void insert(dof_id_type id, Val val)
Definition: generic_projector.h:373
libMesh::as_range
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
Definition: simple_range.h:57
libMesh::FEType::default_quadrature_rule
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
Definition: fe_type.C:31
libMesh::OldSolutionValue::OldSolutionValue
OldSolutionValue(const OldSolutionValue &in)
Definition: generic_projector.h:667
libMesh::GenericProjector::nodes_to_elem
std::unordered_map< dof_id_type, std::vector< dof_id_type > > * nodes_to_elem
Definition: generic_projector.h:92
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::Elem::JUST_COARSENED
Definition: elem.h:1173
libMesh::Elem::is_child_on_side
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const =0
libMesh::DenseVector::size
virtual unsigned int size() const override
Definition: dense_vector.h:92
libMesh::FEAbstract::get_fe_type
FEType get_fe_type() const
Definition: fe_abstract.h:427
libMesh::DofObject::n_vars
unsigned int n_vars(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:891
libMesh::GenericProjector::SubFunctor::join
void join(const SubFunctor &other)
Definition: generic_projector.h:1166
libMesh::numeric_index_type
dof_id_type numeric_index_type
Definition: id_types.h:99
libMesh::MONOMIAL
Definition: enum_fe_family.h:39
libMesh::GenericProjector::SortAndCopy::interiors
std::vector< const Elem * > interiors
Definition: generic_projector.h:263
libMesh::GenericProjector::SortAndCopy
Definition: generic_projector.h:240
libMesh::Elem::parent
const Elem * parent() const
Definition: elem.h:2434
libMesh::FEInterface::dofs_on_side
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)
Fills the vector di with the local degree of freedom indices associated with side s of element elem A...
Definition: fe_interface.C:553
libMesh::DenseMatrix::cholesky_solve
void cholesky_solve(const DenseVector< T2 > &b, DenseVector< T2 > &x)
For symmetric positive definite (SPD) matrices.
Definition: dense_matrix_impl.h:947
libMesh::Elem::refinement_flag
RefinementState refinement_flag() const
Definition: elem.h:2614
libMesh::OldSolutionBase::sys
const System & sys
Definition: generic_projector.h:634
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::GenericProjector::ProjectVertices
Definition: generic_projector.h:267
libMesh::FEType::order
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:197
libMesh::Elem::hmax
virtual Real hmax() const
Definition: elem.C:379
libMesh::Elem::n_children
virtual unsigned int n_children() const =0
libMesh::NumericVector::set
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
libMesh::FEMContext::DOFS_ONLY
Definition: fem_context.h:959
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::FEMContext::set_algebraic_type
void set_algebraic_type(const AlgebraicType atype)
Setting which determines whether to initialize algebraic structures (elem_*) on each element and set ...
Definition: fem_context.h:973
libMesh::Elem::hmin
virtual Real hmin() const
Definition: elem.C:359
libMesh::Elem::subdomain_id
subdomain_id_type subdomain_id() const
Definition: elem.h:2069
libMesh::DofObject::id
dof_id_type id() const
Definition: dof_object.h:767
libMesh::DenseVector::get_values
std::vector< T > & get_values()
Definition: dense_vector.h:256
libMesh::OldSolutionBase::last_elem
const Elem * last_elem
Definition: generic_projector.h:633
libMesh::convert_from_receive
void convert_from_receive(SendT &received, T &converted)
Definition: generic_projector.h:65
libMesh::StoredRange
The StoredRange class defines a contiguous, divisible set of objects.
Definition: stored_range.h:52
libMesh::DofObject::invalid_processor_id
static const processor_id_type invalid_processor_id
An invalid processor_id to distinguish DoFs that have not been assigned to a processor.
Definition: dof_object.h:432
libMesh::GenericProjector::ProjectEdges::ProjectEdges
ProjectEdges(ProjectEdges &p_e, Threads::split)
Definition: generic_projector.h:290
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::GenericProjector::SortAndCopy::join
void join(const SortAndCopy &other)
Definition: generic_projector.h:1612
libMesh::GenericProjector::SubFunctor::action
ProjectionAction action
Definition: generic_projector.h:193
libMesh::Elem::is_face
virtual bool is_face(const unsigned int i) const =0
libMesh::GenericProjector::~GenericProjector
~GenericProjector()
Definition: generic_projector.h:130
libMesh::OldSolutionValue::eval_old_dofs
void eval_old_dofs(const Elem &elem, const FEType &fe_type, unsigned int sys_num, unsigned int var_num, std::vector< dof_id_type > &indices, std::vector< Output > &values)
Definition: generic_projector.h:724
libMesh::GenericProjector::SubFunctor::projector
GenericProjector & projector
Definition: generic_projector.h:154
data
IterBase * data
Ideally this private member data should have protected access.
Definition: variant_filter_iterator.h:337
libMesh::DofObject::n_comp_group
unsigned int n_comp_group(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:939
libMesh::GenericProjector::ids_to_save
std::unordered_map< dof_id_type, FValue > ids_to_save
Definition: generic_projector.h:145
libMesh::GenericProjector::system
const System & system
Definition: generic_projector.h:83
libMesh::Elem::get_node_index
unsigned int get_node_index(const Node *node_ptr) const
Definition: elem.h:2047
libMesh::GenericProjector::SubFunctor::f
FFunctor f
Definition: generic_projector.h:194
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::FEGenericBase::get_phi
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:206
libMesh::OldSolutionBase::old_context
FEMContext old_context
Definition: generic_projector.h:635
libMesh::GenericProjector::SubFunctor::system
const System & system
Definition: generic_projector.h:203
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::GenericProjector::variables
const std::vector< unsigned int > & variables
Definition: generic_projector.h:91
libMesh::Variable::active_on_subdomain
bool active_on_subdomain(subdomain_id_type sid) const
Definition: variable.h:136
libMesh::FEMFunctionWrapper::FEMFunctionWrapper
FEMFunctionWrapper(const FEMFunctionWrapper< Output > &fw)
Definition: generic_projector.h:422
libMesh::FEMContext::OLD
Definition: fem_context.h:962
libMesh::GenericProjector::send_and_insert_dof_values
void send_and_insert_dof_values(std::unordered_map< dof_id_type, std::pair< FValue, processor_id_type >> &ids_to_push, ProjectionAction &action) const
Definition: generic_projector.h:2405
libMesh::VectorSetAction::VectorSetAction
VectorSetAction(NumericVector< Val > &target_vec)
Definition: generic_projector.h:370
libMesh::FEContinuity
FEContinuity
Definition: enum_fe_family.h:77
libMesh::OldSolutionBase::get_shape_outputs
static void get_shape_outputs(FEBase &fe)
libMesh::GenericProjector::node_range
StoredRange< std::vector< node_projection >::const_iterator, node_projection > node_range
Definition: generic_projector.h:211
libMesh::Elem::node_ptr
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2009
libMesh::TypeToSend::type
T type
Definition: generic_projector.h:57
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::convert_to_send
const TypeToSend< T >::type convert_to_send(const T &in)
Definition: generic_projector.h:61
libMesh::FEMFunctionWrapper::eval_at_point
Output eval_at_point(const FEMContext &c, unsigned int i, const Point &n, const Real time)
Definition: generic_projector.h:435
libMesh::GenericProjector::ProjectInteriors::ProjectInteriors
ProjectInteriors(ProjectInteriors &p_i, Threads::split)
Definition: generic_projector.h:332
libMesh::Elem::which_child_am_i
unsigned int which_child_am_i(const Elem *e) const
Definition: elem.h:2596
int
void ErrorVector unsigned int
Definition: adjoints_ex3.C:360
libMesh::GenericProjector::SubProjector::g
std::unique_ptr< GFunctor > g
Definition: generic_projector.h:227
libMesh::DenseVector
Defines a dense vector for use in Finite Element-type computations.
Definition: meshless_interpolation_function.h:39
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::FEMFunctionWrapper::is_grid_projection
bool is_grid_projection()
Definition: generic_projector.h:441
libMesh::GenericProjector::ProjectVertices::operator()
void operator()(const node_range &range)
Definition: generic_projector.h:1655
libMesh::NumericVector::first_local_index
virtual numeric_index_type first_local_index() const =0
libMesh::GenericProjector::ProjectEdges::ProjectEdges
ProjectEdges(GenericProjector &p)
Definition: generic_projector.h:288
libMesh::VectorSetAction::insert
void insert(const std::vector< dof_id_type > &dof_indices, const DenseVector< Val > &Ue)
Definition: generic_projector.h:384
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62