libMesh
system_projection.C
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 // C++ includes
21 #include <vector>
22 #include <numeric> // std::iota
23 
24 // Local includes
25 #include "libmesh/libmesh_config.h"
26 
27 #ifdef LIBMESH_HAVE_METAPHYSICL
28 
29 // With quad precision we need the shim function declarations to
30 // precede the MetaPhysicL use of them
31 #include "libmesh/libmesh_common.h"
32 
33 // FIXME - having to do this with MetaPhysicL brings me shame - RHS
34 #include "libmesh/ignore_warnings.h"
35 // Template specialization declarations in here need to *precede* code
36 // using them.
37 #include "metaphysicl/dynamicsparsenumberarray_decl.h"
38 #include "libmesh/restore_warnings.h"
39 
40 #include "libmesh/compare_types.h"
41 #include "libmesh/int_range.h"
42 
43 using MetaPhysicL::DynamicSparseNumberArray;
44 
45 namespace libMesh
46 {
47 // From the perspective of libMesh gradient vectors, a DSNA is a
48 // scalar component
49 template <typename T, typename IndexType>
50 struct ScalarTraits<MetaPhysicL::DynamicSparseNumberArray<T,IndexType> >
51 {
52  static const bool value = true;
53 };
54 
55 // And although MetaPhysicL knows how to combine DSNA with something
56 // else, we need to teach libMesh too.
57 template <typename T, typename IndexType, typename T2>
58 struct CompareTypes<MetaPhysicL::DynamicSparseNumberArray<T,IndexType>, T2>
59 {
60  typedef typename
61  MetaPhysicL::DynamicSparseNumberArray
63 };
64 
65 template <typename T> struct TypeToSend;
66 
67 template <typename T, typename IndexType>
68 struct TypeToSend<MetaPhysicL::DynamicSparseNumberArray<T,IndexType>> {
69  typedef std::vector<std::pair<IndexType,T>> type;
70 };
71 
72 template <typename T, typename IndexType>
73 const std::vector<std::pair<IndexType,T>>
74 convert_to_send(MetaPhysicL::DynamicSparseNumberArray<T,IndexType> & in)
75 {
76  const std::size_t in_size = in.size();
77  std::vector<std::pair<IndexType,T>> returnval(in_size);
78 
79  for (std::size_t i=0; i != in_size; ++i)
80  {
81  returnval[i].first = in.raw_index(i);
82  returnval[i].second = in.raw_at(i);
83  }
84  return returnval;
85 }
86 
87 template <typename SendT, typename T, typename IndexType>
88 void convert_from_receive (SendT & received,
89  MetaPhysicL::DynamicSparseNumberArray<T,IndexType> & converted)
90 {
91  const std::size_t received_size = received.size();
92  converted.resize(received_size);
93  for (std::size_t i=0; i != received_size; ++i)
94  {
95  converted.raw_index(i) = received[i].first;
96  converted.raw_at(i) = received[i].second;
97  }
98 }
99 
100 }
101 
102 
103 #endif
104 
105 #include "libmesh/boundary_info.h"
106 #include "libmesh/dense_matrix.h"
107 #include "libmesh/dense_vector.h"
108 #include "libmesh/dof_map.h"
109 #include "libmesh/elem.h"
110 #include "libmesh/fe_base.h"
111 #include "libmesh/fe_interface.h"
112 #include "libmesh/generic_projector.h"
113 #include "libmesh/libmesh_logging.h"
114 #include "libmesh/mesh_base.h"
115 #include "libmesh/numeric_vector.h"
116 #include "libmesh/quadrature.h"
117 #include "libmesh/sparse_matrix.h"
118 #include "libmesh/system.h"
119 #include "libmesh/threads.h"
120 #include "libmesh/wrapped_function.h"
121 #include "libmesh/wrapped_functor.h"
122 
123 
124 
125 #ifdef LIBMESH_HAVE_METAPHYSICL
126 // FIXME - wrapping MetaPhysicL is shameful - RHS
127 #include "libmesh/ignore_warnings.h"
128 // Include MetaPhysicL definitions finally
129 #include "metaphysicl/dynamicsparsenumberarray.h"
130 #include "libmesh/restore_warnings.h"
131 
132 // And make sure we instantiate the methods we'll need to use on them.
133 #include "libmesh/dense_matrix_impl.h"
134 
135 namespace libMesh {
136 typedef DynamicSparseNumberArray<Real, dof_id_type> DSNAN;
137 
138 template void
141 template void
143  DenseVector<DSNAN> &) const;
144 }
145 #endif
146 
147 
148 
149 namespace libMesh
150 {
151 
152 // ------------------------------------------------------------
153 // Helper class definitions
154 
155 #ifdef LIBMESH_ENABLE_AMR
156 
168 {
169 private:
170  const System & system;
171 
172 public:
173  BuildProjectionList (const System & system_in) :
174  system(system_in),
175  send_list()
176  {}
177 
179  system(other.system),
180  send_list()
181  {}
182 
183  void unique();
184  void operator()(const ConstElemRange & range);
185  void join (const BuildProjectionList & other);
186  std::vector<dof_id_type> send_list;
187 };
188 
189 #endif // LIBMESH_ENABLE_AMR
190 
191 
198 {
199 private:
200  const std::set<boundary_id_type> & b;
201  const std::vector<unsigned int> & variables;
202  const System & system;
203  std::unique_ptr<FunctionBase<Number>> f;
204  std::unique_ptr<FunctionBase<Gradient>> g;
207 
208 public:
209  BoundaryProjectSolution (const std::set<boundary_id_type> & b_in,
210  const std::vector<unsigned int> & variables_in,
211  const System & system_in,
212  FunctionBase<Number> * f_in,
213  FunctionBase<Gradient> * g_in,
214  const Parameters & parameters_in,
215  NumericVector<Number> & new_v_in) :
216  b(b_in),
217  variables(variables_in),
218  system(system_in),
219  f(f_in ? f_in->clone() : std::unique_ptr<FunctionBase<Number>>()),
220  g(g_in ? g_in->clone() : std::unique_ptr<FunctionBase<Gradient>>()),
221  parameters(parameters_in),
222  new_vector(new_v_in)
223  {
224  libmesh_assert(f.get());
225  f->init();
226  if (g.get())
227  g->init();
228  }
229 
231  b(in.b),
232  variables(in.variables),
233  system(in.system),
234  f(in.f.get() ? in.f->clone() : std::unique_ptr<FunctionBase<Number>>()),
235  g(in.g.get() ? in.g->clone() : std::unique_ptr<FunctionBase<Gradient>>()),
236  parameters(in.parameters),
237  new_vector(in.new_vector)
238  {
239  libmesh_assert(f.get());
240  f->init();
241  if (g.get())
242  g->init();
243  }
244 
245  void operator()(const ConstElemRange & range) const;
246 };
247 
248 
249 
250 // ------------------------------------------------------------
251 // System implementation
252 void System::project_vector (NumericVector<Number> & vector,
253  int is_adjoint) const
254 {
255  // Create a copy of the vector, which currently
256  // contains the old data.
257  std::unique_ptr<NumericVector<Number>>
258  old_vector (vector.clone());
259 
260  // Project the old vector to the new vector
261  this->project_vector (*old_vector, vector, is_adjoint);
262 }
263 
264 
270 void System::project_vector (const NumericVector<Number> & old_v,
271  NumericVector<Number> & new_v,
272  int is_adjoint) const
273 {
274  LOG_SCOPE ("project_vector(old,new)", "System");
275 
282  new_v.clear();
283 
284 #ifdef LIBMESH_ENABLE_AMR
285 
286  // Resize the new vector and get a serial version.
287  NumericVector<Number> * new_vector_ptr = nullptr;
288  std::unique_ptr<NumericVector<Number>> new_vector_built;
289  NumericVector<Number> * local_old_vector;
290  std::unique_ptr<NumericVector<Number>> local_old_vector_built;
291  const NumericVector<Number> * old_vector_ptr = nullptr;
292 
293  ConstElemRange active_local_elem_range
294  (this->get_mesh().active_local_elements_begin(),
295  this->get_mesh().active_local_elements_end());
296 
297  // If the old vector was uniprocessor, make the new
298  // vector uniprocessor
299  if (old_v.type() == SERIAL)
300  {
301  new_v.init (this->n_dofs(), false, SERIAL);
302  new_vector_ptr = &new_v;
303  old_vector_ptr = &old_v;
304  }
305 
306  // Otherwise it is a parallel, distributed vector, which
307  // we need to localize.
308  else if (old_v.type() == PARALLEL)
309  {
310  // Build a send list for efficient localization
311  BuildProjectionList projection_list(*this);
312  Threads::parallel_reduce (active_local_elem_range,
313  projection_list);
314 
315  // Create a sorted, unique send_list
316  projection_list.unique();
317 
318  new_v.init (this->n_dofs(), this->n_local_dofs(), false, PARALLEL);
319  new_vector_built = NumericVector<Number>::build(this->comm());
320  local_old_vector_built = NumericVector<Number>::build(this->comm());
321  new_vector_ptr = new_vector_built.get();
322  local_old_vector = local_old_vector_built.get();
323  new_vector_ptr->init(this->n_dofs(), false, SERIAL);
324  local_old_vector->init(old_v.size(), false, SERIAL);
325  old_v.localize(*local_old_vector, projection_list.send_list);
326  local_old_vector->close();
327  old_vector_ptr = local_old_vector;
328  }
329  else if (old_v.type() == GHOSTED)
330  {
331  // Build a send list for efficient localization
332  BuildProjectionList projection_list(*this);
333  Threads::parallel_reduce (active_local_elem_range,
334  projection_list);
335 
336  // Create a sorted, unique send_list
337  projection_list.unique();
338 
339  new_v.init (this->n_dofs(), this->n_local_dofs(),
340  this->get_dof_map().get_send_list(), false, GHOSTED);
341 
342  local_old_vector_built = NumericVector<Number>::build(this->comm());
343  new_vector_ptr = &new_v;
344  local_old_vector = local_old_vector_built.get();
345  local_old_vector->init(old_v.size(), old_v.local_size(),
346  projection_list.send_list, false, GHOSTED);
347  old_v.localize(*local_old_vector, projection_list.send_list);
348  local_old_vector->close();
349  old_vector_ptr = local_old_vector;
350  }
351  else // unknown old_v.type()
352  libmesh_error_msg("ERROR: Unknown old_v.type() == " << old_v.type());
353 
354  // Note that the above will have zeroed the new_vector.
355  // Just to be sure, assert that new_vector_ptr and old_vector_ptr
356  // were successfully set before trying to deref them.
357  libmesh_assert(new_vector_ptr);
358  libmesh_assert(old_vector_ptr);
359 
360  NumericVector<Number> & new_vector = *new_vector_ptr;
361  const NumericVector<Number> & old_vector = *old_vector_ptr;
362 
363  const unsigned int n_variables = this->n_vars();
364 
365  if (n_variables)
366  {
367  std::vector<unsigned int> vars(n_variables);
368  std::iota(vars.begin(), vars.end(), 0);
369 
370  // Use a typedef to make the calling sequence for parallel_for() a bit more readable
371  typedef
374  Number, VectorSetAction<Number>> FEMProjector;
375 
378  VectorSetAction<Number> setter(new_vector);
379 
380  FEMProjector projector(*this, f, &g, setter, vars);
381  projector.project(active_local_elem_range);
382 
383  // Copy the SCALAR dofs from old_vector to new_vector
384  // Note: We assume that all SCALAR dofs are on the
385  // processor with highest ID
386  if (this->processor_id() == (this->n_processors()-1))
387  {
388  const DofMap & dof_map = this->get_dof_map();
389  for (auto var : IntRange<unsigned int>(0, this->n_vars()))
390  if (this->variable(var).type().family == SCALAR)
391  {
392  // We can just map SCALAR dofs directly across
393  std::vector<dof_id_type> new_SCALAR_indices, old_SCALAR_indices;
394  dof_map.SCALAR_dof_indices (new_SCALAR_indices, var, false);
395  dof_map.SCALAR_dof_indices (old_SCALAR_indices, var, true);
396  for (auto i : index_range(new_SCALAR_indices))
397  new_vector.set(new_SCALAR_indices[i], old_vector(old_SCALAR_indices[i]));
398  }
399  }
400  }
401 
402  new_vector.close();
403 
404  // If the old vector was serial, we probably need to send our values
405  // to other processors
406  //
407  // FIXME: I'm not sure how to make a NumericVector do that without
408  // creating a temporary parallel vector to use localize! - RHS
409  if (old_v.type() == SERIAL)
410  {
411  std::unique_ptr<NumericVector<Number>> dist_v = NumericVector<Number>::build(this->comm());
412  dist_v->init(this->n_dofs(), this->n_local_dofs(), false, PARALLEL);
413  dist_v->close();
414 
415  for (auto i : IntRange<dof_id_type>(0, dist_v->size()))
416  if (new_vector(i) != 0.0)
417  dist_v->set(i, new_vector(i));
418 
419  dist_v->close();
420 
421  dist_v->localize (new_v, this->get_dof_map().get_send_list());
422  new_v.close();
423  }
424  // If the old vector was parallel, we need to update it
425  // and free the localized copies
426  else if (old_v.type() == PARALLEL)
427  {
428  // We may have to set dof values that this processor doesn't
429  // own in certain special cases, like LAGRANGE FIRST or
430  // HERMITE THIRD elements on second-order meshes
431  for (auto i : IntRange<dof_id_type>(0, new_v.size()))
432  if (new_vector(i) != 0.0)
433  new_v.set(i, new_vector(i));
434  new_v.close();
435  }
436 
437  if (is_adjoint == -1)
438  this->get_dof_map().enforce_constraints_exactly(*this, &new_v);
439  else if (is_adjoint >= 0)
440  this->get_dof_map().enforce_adjoint_constraints_exactly(new_v,
441  is_adjoint);
442 
443 #else
444 
445  // AMR is disabled: simply copy the vector
446  new_v = old_v;
447 
448  libmesh_ignore(is_adjoint);
449 
450 #endif // #ifdef LIBMESH_ENABLE_AMR
451 }
452 
453 
454 #ifdef LIBMESH_ENABLE_AMR
455 #ifdef LIBMESH_HAVE_METAPHYSICL
456 
457 template <typename Output>
459 {
460 public:
461  typedef DynamicSparseNumberArray<Output, dof_id_type> type;
462 };
463 
464 template <typename InnerOutput>
465 class DSNAOutput<VectorValue<InnerOutput> >
466 {
467 public:
469 };
470 
480 template <typename Output,
481  void (FEMContext::*point_output) (unsigned int,
482  const Point &,
483  Output &,
484  const Real) const>
486 {
487 public:
488  typedef typename DSNAOutput<Output>::type DSNA;
489 
491  OldSolutionBase<Output, point_output>(sys_in)
492  {
493  this->old_context.set_algebraic_type(FEMContext::OLD_DOFS_ONLY);
494  }
495 
497  OldSolutionBase<Output, point_output>(in.sys)
498  {
499  this->old_context.set_algebraic_type(FEMContext::OLD_DOFS_ONLY);
500  }
501 
502  DSNA eval_at_node (const FEMContext & c,
503  unsigned int i,
504  unsigned int elem_dim,
505  const Node & n,
506  bool extra_hanging_dofs,
507  Real /* time */ = 0.);
508 
509  DSNA eval_at_point(const FEMContext & c,
510  unsigned int i,
511  const Point & p,
512  Real /* time */ = 0.);
513 
514  void eval_old_dofs (const Elem & elem,
515  unsigned int node_num,
516  unsigned int var_num,
517  std::vector<dof_id_type> & indices,
518  std::vector<DSNA> & values)
519  {
520  LOG_SCOPE ("eval_old_dofs(node)", "OldSolutionCoefs");
521 
522  this->sys.get_dof_map().dof_indices(elem, node_num, indices, var_num);
523 
524  std::vector<dof_id_type> old_indices;
525 
526  this->sys.get_dof_map().old_dof_indices(elem, node_num, old_indices, var_num);
527 
528  libmesh_assert_equal_to (old_indices.size(), indices.size());
529 
530  values.resize(old_indices.size());
531 
532  for (auto i : index_range(values))
533  {
534  values[i].resize(1);
535  values[i].raw_at(0) = 1;
536  values[i].raw_index(0) = old_indices[i];
537  }
538  }
539 
540 
541  void eval_old_dofs (const Elem & elem,
542  const FEType & fe_type,
543  unsigned int sys_num,
544  unsigned int var_num,
545  std::vector<dof_id_type> & indices,
546  std::vector<DSNA> & values)
547  {
548  LOG_SCOPE ("eval_old_dofs(elem)", "OldSolutionCoefs");
549 
550  // We're only to be asked for old dofs on elements that can copy
551  // them through DO_NOTHING or through refinement.
552  const Elem & old_elem =
553  (elem.refinement_flag() == Elem::JUST_REFINED) ?
554  *elem.parent() : elem;
555 
556  // If there are any element-based DOF numbers, get them
557  const unsigned int nc = FEInterface::n_dofs_per_elem(elem.dim(),
558  fe_type,
559  elem.type());
560 
561  std::vector<dof_id_type> old_dof_indices(nc);
562  indices.resize(nc);
563 
564  // We should never have fewer dofs than necessary on an
565  // element unless we're getting indices on a parent element,
566  // and we should never need those indices
567  if (nc != 0)
568  {
569  libmesh_assert(old_elem.old_dof_object);
570 
571  const std::pair<unsigned int, unsigned int>
572  vg_and_offset = elem.var_to_vg_and_offset(sys_num,var_num);
573  const unsigned int vg = vg_and_offset.first;
574  const unsigned int vig = vg_and_offset.second;
575 
576  const unsigned int n_comp = elem.n_comp_group(sys_num,vg);
577  libmesh_assert_greater(elem.n_systems(), sys_num);
578  libmesh_assert_greater_equal(n_comp, nc);
579 
580  for (unsigned int i=0; i<nc; i++)
581  {
582  const dof_id_type d_old =
583  old_elem.old_dof_object->dof_number(sys_num, vg, vig, i, n_comp);
584  const dof_id_type d_new =
585  elem.dof_number(sys_num, vg, vig, i, n_comp);
586  libmesh_assert_not_equal_to (d_old, DofObject::invalid_id);
587  libmesh_assert_not_equal_to (d_new, DofObject::invalid_id);
588 
589  old_dof_indices[i] = d_old;
590  indices[i] = d_new;
591  }
592  }
593 
594  values.resize(old_dof_indices.size());
595 
596  for (auto i : index_range(values))
597  {
598  values[i].resize(1);
599  values[i].raw_at(0) = 1;
600  values[i].raw_index(0) = old_dof_indices[i];
601  }
602  }
603 };
604 
605 
606 
607 template<>
608 inline
609 DynamicSparseNumberArray<Real, dof_id_type>
610 OldSolutionCoefs<Real, &FEMContext::point_value>::
611 eval_at_point(const FEMContext & c,
612  unsigned int i,
613  const Point & p,
614  Real /* time */)
615 {
616  LOG_SCOPE ("eval_at_point()", "OldSolutionCoefs");
617 
618  if (!this->check_old_context(c, p))
619  return 0;
620 
621  // Get finite element object
622  FEGenericBase<Real> * fe = nullptr;
623  this->old_context.get_element_fe<Real>
624  (i, fe, this->old_context.get_elem_dim());
625 
626  // Build a FE for calculating phi(p)
627  FEGenericBase<Real> * fe_new =
628  this->old_context.build_new_fe(fe, p);
629 
630  // Get the values and global indices of the shape functions
631  const std::vector<std::vector<Real> > & phi = fe_new->get_phi();
632  const std::vector<dof_id_type> & dof_indices =
633  this->old_context.get_dof_indices(i);
634 
635  const std::size_t n_dofs = phi.size();
636  libmesh_assert_equal_to(n_dofs, dof_indices.size());
637 
638  DynamicSparseNumberArray<Real, dof_id_type> returnval;
639  returnval.resize(n_dofs);
640 
641  for (auto j : index_range(phi))
642  {
643  returnval.raw_at(j) = phi[j][0];
644  returnval.raw_index(j) = dof_indices[j];
645  }
646 
647  return returnval;
648 }
649 
650 
651 
652 template<>
653 inline
657  unsigned int i,
658  const Point & p,
659  Real /* time */)
660 {
661  LOG_SCOPE ("eval_at_point()", "OldSolutionCoefs");
662 
663  if (!this->check_old_context(c, p))
664  return 0;
665 
666  // Get finite element object
667  FEGenericBase<Real> * fe = nullptr;
668  this->old_context.get_element_fe<Real>
669  (i, fe, this->old_context.get_elem_dim());
670 
671  // Build a FE for calculating phi(p)
672  FEGenericBase<Real> * fe_new =
673  this->old_context.build_new_fe(fe, p);
674 
675  // Get the values and global indices of the shape functions
676  const std::vector<std::vector<RealGradient> > & dphi = fe_new->get_dphi();
677  const std::vector<dof_id_type> & dof_indices =
678  this->old_context.get_dof_indices(i);
679 
680  const std::size_t n_dofs = dphi.size();
681  libmesh_assert_equal_to(n_dofs, dof_indices.size());
682 
684 
685  for (unsigned int d = 0; d != LIBMESH_DIM; ++d)
686  returnval(d).resize(n_dofs);
687 
688  for (auto j : index_range(dphi))
689  for (int d = 0; d != LIBMESH_DIM; ++d)
690  {
691  returnval(d).raw_at(j) = dphi[j][0](d);
692  returnval(d).raw_index(j) = dof_indices[j];
693  }
694 
695  return returnval;
696 }
697 
698 
699 template<>
700 inline
701 DynamicSparseNumberArray<Real, dof_id_type>
704  unsigned int i,
705  unsigned int /* elem_dim */,
706  const Node & n,
707  bool extra_hanging_dofs,
708  Real /* time */)
709 {
710  LOG_SCOPE ("Real eval_at_node()", "OldSolutionCoefs");
711 
712  // Optimize for the common case, where this node was part of the
713  // old solution.
714  //
715  // Be sure to handle cases where the variable wasn't defined on
716  // this node (due to changing subdomain support) or where the
717  // variable has no components on this node (due to Elem order
718  // exceeding FE order) or where the old_dof_object dofs might
719  // correspond to non-vertex dofs (due to extra_hanging_dofs and
720  // refinement)
721 
723 
724  if (n.old_dof_object &&
725  (!extra_hanging_dofs ||
726  flag == Elem::JUST_COARSENED ||
727  flag == Elem::DO_NOTHING) &&
728  n.old_dof_object->n_vars(sys.number()) &&
729  n.old_dof_object->n_comp(sys.number(), i))
730  {
731  DynamicSparseNumberArray<Real, dof_id_type> returnval;
732  const dof_id_type old_id =
733  n.old_dof_object->dof_number(sys.number(), i, 0);
734  returnval.resize(1);
735  returnval.raw_at(0) = 1;
736  returnval.raw_index(0) = old_id;
737  return returnval;
738  }
739 
740  return this->eval_at_point(c, i, n, 0);
741 }
742 
743 
744 
745 template<>
746 inline
750  unsigned int i,
751  unsigned int elem_dim,
752  const Node & n,
753  bool extra_hanging_dofs,
754  Real /* time */)
755 {
756  LOG_SCOPE ("RealGradient eval_at_node()", "OldSolutionCoefs");
757 
758  // Optimize for the common case, where this node was part of the
759  // old solution.
760  //
761  // Be sure to handle cases where the variable wasn't defined on
762  // this node (due to changing subdomain support) or where the
763  // variable has no components on this node (due to Elem order
764  // exceeding FE order) or where the old_dof_object dofs might
765  // correspond to non-vertex dofs (due to extra_hanging_dofs and
766  // refinement)
767 
769 
770  if (n.old_dof_object &&
771  (!extra_hanging_dofs ||
772  flag == Elem::JUST_COARSENED ||
773  flag == Elem::DO_NOTHING) &&
774  n.old_dof_object->n_vars(sys.number()) &&
775  n.old_dof_object->n_comp(sys.number(), i))
776  {
778  for (unsigned int d = 0; d != elem_dim; ++d)
779  {
780  const dof_id_type old_id =
781  n.old_dof_object->dof_number(sys.number(), i, d+1);
782  g(d).resize(1);
783  g(d).raw_at(0) = 1;
784  g(d).raw_index(0) = old_id;
785  }
786  return g;
787  }
788 
789  return this->eval_at_point(c, i, n, 0);
790 }
791 
792 
793 
802 template <typename ValIn, typename ValOut>
804 {
805 private:
807 
808 public:
810  target_matrix(target_mat) {}
811 
813  const DynamicSparseNumberArray<ValIn, dof_id_type> & val)
814  {
815  // Lock the target matrix since it is shared among threads.
816  {
817  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
818 
819  const std::size_t dnsa_size = val.size();
820  for (unsigned int j = 0; j != dnsa_size; ++j)
821  {
822  const dof_id_type dof_j = val.raw_index(j);
823  const ValIn dof_val = val.raw_at(j);
824  target_matrix.set(id, dof_j, dof_val);
825  }
826  }
827  }
828 
829 
830  void insert(const std::vector<dof_id_type> & dof_indices,
831  const std::vector<DynamicSparseNumberArray<ValIn, dof_id_type> > & Ue)
832  {
833  const numeric_index_type
834  begin_dof = target_matrix.row_start(),
835  end_dof = target_matrix.row_stop();
836 
837  unsigned int size = Ue.size();
838 
839  libmesh_assert_equal_to(size, dof_indices.size());
840 
841  // Lock the target matrix since it is shared among threads.
842  {
843  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
844 
845  for (unsigned int i = 0; i != size; ++i)
846  {
847  const dof_id_type dof_i = dof_indices[i];
848  if ((dof_i >= begin_dof) && (dof_i < end_dof))
849  {
850  const DynamicSparseNumberArray<ValIn,dof_id_type> & dnsa = Ue[i];
851  const std::size_t dnsa_size = dnsa.size();
852  for (unsigned int j = 0; j != dnsa_size; ++j)
853  {
854  const dof_id_type dof_j = dnsa.raw_index(j);
855  const ValIn dof_val = dnsa.raw_at(j);
856  target_matrix.set(dof_i, dof_j, dof_val);
857  }
858  }
859  }
860  }
861  }
862 };
863 
864 
865 
870 void System::projection_matrix (SparseMatrix<Number> & proj_mat) const
871 {
872  LOG_SCOPE ("projection_matrix()", "System");
873 
874  const unsigned int n_variables = this->n_vars();
875 
876  if (n_variables)
877  {
878  ConstElemRange active_local_elem_range
879  (this->get_mesh().active_local_elements_begin(),
880  this->get_mesh().active_local_elements_end());
881 
882  std::vector<unsigned int> vars(n_variables);
883  std::iota(vars.begin(), vars.end(), 0);
884 
885  // Use a typedef to make the calling sequence for parallel_for() a bit more readable
886  typedef OldSolutionCoefs<Real, &FEMContext::point_value> OldSolutionValueCoefs;
887  typedef OldSolutionCoefs<RealGradient, &FEMContext::point_gradient> OldSolutionGradientCoefs;
888 
889  typedef
890  GenericProjector<OldSolutionValueCoefs,
891  OldSolutionGradientCoefs,
892  DynamicSparseNumberArray<Real,dof_id_type>,
893  MatrixFillAction<Real, Number> > ProjMatFiller;
894 
895  OldSolutionValueCoefs f(*this);
896  OldSolutionGradientCoefs g(*this);
897  MatrixFillAction<Real, Number> setter(proj_mat);
898 
899  ProjMatFiller mat_filler(*this, f, &g, setter, vars);
900  mat_filler.project(active_local_elem_range);
901 
902  // Set the SCALAR dof transfer entries too.
903  // Note: We assume that all SCALAR dofs are on the
904  // processor with highest ID
905  if (this->processor_id() == (this->n_processors()-1))
906  {
907  const DofMap & dof_map = this->get_dof_map();
908  for (auto var : IntRange<unsigned int>(0, this->n_vars()))
909  if (this->variable(var).type().family == SCALAR)
910  {
911  // We can just map SCALAR dofs directly across
912  std::vector<dof_id_type> new_SCALAR_indices, old_SCALAR_indices;
913  dof_map.SCALAR_dof_indices (new_SCALAR_indices, var, false);
914  dof_map.SCALAR_dof_indices (old_SCALAR_indices, var, true);
915  const unsigned int new_n_dofs =
916  cast_int<unsigned int>(new_SCALAR_indices.size());
917 
918  for (unsigned int i=0; i<new_n_dofs; i++)
919  {
920  proj_mat.set( new_SCALAR_indices[i],
921  old_SCALAR_indices[i], 1);
922  }
923  }
924  }
925  }
926 }
927 #endif // LIBMESH_HAVE_METAPHYSICL
928 #endif // LIBMESH_ENABLE_AMR
929 
930 
931 
936 void System::project_solution (ValueFunctionPointer fptr,
937  GradientFunctionPointer gptr,
938  const Parameters & parameters) const
939 {
940  WrappedFunction<Number> f(*this, fptr, &parameters);
941  WrappedFunction<Gradient> g(*this, gptr, &parameters);
942  this->project_solution(&f, &g);
943 }
944 
945 
950 void System::project_solution (FunctionBase<Number> * f,
951  FunctionBase<Gradient> * g) const
952 {
953  this->project_vector(*solution, f, g);
954 
955  solution->localize(*current_local_solution, _dof_map->get_send_list());
956 }
957 
958 
963 void System::project_solution (FEMFunctionBase<Number> * f,
964  FEMFunctionBase<Gradient> * g) const
965 {
966  this->project_vector(*solution, f, g);
967 
968  solution->localize(*current_local_solution, _dof_map->get_send_list());
969 }
970 
971 
976 void System::project_vector (ValueFunctionPointer fptr,
977  GradientFunctionPointer gptr,
978  const Parameters & parameters,
979  NumericVector<Number> & new_vector,
980  int is_adjoint) const
981 {
982  WrappedFunction<Number> f(*this, fptr, &parameters);
983  WrappedFunction<Gradient> g(*this, gptr, &parameters);
984  this->project_vector(new_vector, &f, &g, is_adjoint);
985 }
986 
991 void System::project_vector (NumericVector<Number> & new_vector,
994  int is_adjoint) const
995 {
996  LOG_SCOPE ("project_vector(FunctionBase)", "System");
997 
998  libmesh_assert(f);
999 
1000  WrappedFunctor<Number> f_fem(*f);
1001 
1002  if (g)
1003  {
1004  WrappedFunctor<Gradient> g_fem(*g);
1005 
1006  this->project_vector(new_vector, &f_fem, &g_fem, is_adjoint);
1007  }
1008  else
1009  this->project_vector(new_vector, &f_fem, nullptr, is_adjoint);
1010 }
1011 
1012 
1017 void System::project_vector (NumericVector<Number> & new_vector,
1020  int is_adjoint) const
1021 {
1022  LOG_SCOPE ("project_fem_vector()", "System");
1023 
1024  libmesh_assert (f);
1025 
1026  ConstElemRange active_local_range
1027  (this->get_mesh().active_local_elements_begin(),
1028  this->get_mesh().active_local_elements_end() );
1029 
1030  VectorSetAction<Number> setter(new_vector);
1031 
1032  const unsigned int n_variables = this->n_vars();
1033 
1034  std::vector<unsigned int> vars(n_variables);
1035  std::iota(vars.begin(), vars.end(), 0);
1036 
1037  // Use a typedef to make the calling sequence for parallel_for() a bit more readable
1038  typedef
1040  Number, VectorSetAction<Number>> FEMProjector;
1041 
1043 
1044  if (g)
1045  {
1047 
1048  FEMProjector projector(*this, fw, &gw, setter, vars);
1049  projector.project(active_local_range);
1050  }
1051  else
1052  {
1053  FEMProjector projector(*this, fw, nullptr, setter, vars);
1054  projector.project(active_local_range);
1055  }
1056 
1057  // Also, load values into the SCALAR dofs
1058  // Note: We assume that all SCALAR dofs are on the
1059  // processor with highest ID
1060  if (this->processor_id() == (this->n_processors()-1))
1061  {
1062  // FIXME: Do we want to first check for SCALAR vars before building this? [PB]
1063  FEMContext context( *this );
1064 
1065  const DofMap & dof_map = this->get_dof_map();
1066  for (auto var : IntRange<unsigned int>(0, this->n_vars()))
1067  if (this->variable(var).type().family == SCALAR)
1068  {
1069  // FIXME: We reinit with an arbitrary element in case the user
1070  // doesn't override FEMFunctionBase::component. Is there
1071  // any use case we're missing? [PB]
1072  context.pre_fe_reinit(*this, *(this->get_mesh().active_local_elements_begin()));
1073 
1074  std::vector<dof_id_type> SCALAR_indices;
1075  dof_map.SCALAR_dof_indices (SCALAR_indices, var);
1076  const unsigned int n_SCALAR_dofs =
1077  cast_int<unsigned int>(SCALAR_indices.size());
1078 
1079  for (unsigned int i=0; i<n_SCALAR_dofs; i++)
1080  {
1081  const dof_id_type global_index = SCALAR_indices[i];
1082  const unsigned int component_index =
1083  this->variable_scalar_number(var,i);
1084 
1085  new_vector.set(global_index, f->component(context, component_index, Point(), this->time));
1086  }
1087  }
1088  }
1089 
1090  new_vector.close();
1091 
1092 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1093  if (is_adjoint == -1)
1094  this->get_dof_map().enforce_constraints_exactly(*this, &new_vector);
1095  else if (is_adjoint >= 0)
1096  this->get_dof_map().enforce_adjoint_constraints_exactly(new_vector,
1097  is_adjoint);
1098 #else
1099  libmesh_ignore(is_adjoint);
1100 #endif
1101 }
1102 
1103 
1109 void System::boundary_project_solution (const std::set<boundary_id_type> & b,
1110  const std::vector<unsigned int> & variables,
1111  ValueFunctionPointer fptr,
1112  GradientFunctionPointer gptr,
1113  const Parameters & parameters)
1114 {
1115  WrappedFunction<Number> f(*this, fptr, &parameters);
1116  WrappedFunction<Gradient> g(*this, gptr, &parameters);
1117  this->boundary_project_solution(b, variables, &f, &g);
1118 }
1119 
1120 
1126 void System::boundary_project_solution (const std::set<boundary_id_type> & b,
1127  const std::vector<unsigned int> & variables,
1130 {
1131  this->boundary_project_vector(b, variables, *solution, f, g);
1132 
1133  solution->localize(*current_local_solution);
1134 }
1135 
1136 
1137 
1138 
1139 
1144 void System::boundary_project_vector (const std::set<boundary_id_type> & b,
1145  const std::vector<unsigned int> & variables,
1146  ValueFunctionPointer fptr,
1147  GradientFunctionPointer gptr,
1148  const Parameters & parameters,
1149  NumericVector<Number> & new_vector,
1150  int is_adjoint) const
1151 {
1152  WrappedFunction<Number> f(*this, fptr, &parameters);
1153  WrappedFunction<Gradient> g(*this, gptr, &parameters);
1154  this->boundary_project_vector(b, variables, new_vector, &f, &g,
1155  is_adjoint);
1156 }
1157 
1162 void System::boundary_project_vector (const std::set<boundary_id_type> & b,
1163  const std::vector<unsigned int> & variables,
1164  NumericVector<Number> & new_vector,
1167  int is_adjoint) const
1168 {
1169  LOG_SCOPE ("boundary_project_vector()", "System");
1170 
1172  (ConstElemRange (this->get_mesh().active_local_elements_begin(),
1173  this->get_mesh().active_local_elements_end() ),
1174  BoundaryProjectSolution(b, variables, *this, f, g,
1175  this->get_equation_systems().parameters,
1176  new_vector)
1177  );
1178 
1179  // We don't do SCALAR dofs when just projecting the boundary, so
1180  // we're done here.
1181 
1182  new_vector.close();
1183 
1184 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1185  if (is_adjoint == -1)
1186  this->get_dof_map().enforce_constraints_exactly(*this, &new_vector);
1187  else if (is_adjoint >= 0)
1188  this->get_dof_map().enforce_adjoint_constraints_exactly(new_vector,
1189  is_adjoint);
1190 #else
1191  libmesh_ignore(is_adjoint);
1192 #endif
1193 }
1194 
1195 
1196 
1197 #ifdef LIBMESH_ENABLE_AMR
1198 void BuildProjectionList::unique()
1199 {
1200  // Sort the send list. After this duplicated
1201  // elements will be adjacent in the vector
1202  std::sort(this->send_list.begin(),
1203  this->send_list.end());
1204 
1205  // Now use std::unique to remove duplicate entries
1206  std::vector<dof_id_type>::iterator new_end =
1207  std::unique (this->send_list.begin(),
1208  this->send_list.end());
1209 
1210  // Remove the end of the send_list. Use the "swap trick"
1211  // from Effective STL
1212  std::vector<dof_id_type>
1213  (this->send_list.begin(), new_end).swap (this->send_list);
1214 }
1215 
1216 
1217 
1218 void BuildProjectionList::operator()(const ConstElemRange & range)
1219 {
1220  // The DofMap for this system
1221  const DofMap & dof_map = system.get_dof_map();
1222 
1223  const dof_id_type first_old_dof = dof_map.first_old_dof();
1224  const dof_id_type end_old_dof = dof_map.end_old_dof();
1225 
1226  // We can handle all the variables at once.
1227  // The old global DOF indices
1228  std::vector<dof_id_type> di;
1229 
1230  // Iterate over the elements in the range
1231  for (const auto & elem : range)
1232  {
1233  // If this element doesn't have an old_dof_object with dofs for the
1234  // current system, then it must be newly added, so the user
1235  // is responsible for setting the new dofs.
1236 
1237  // ... but we need a better way to test for that; the code
1238  // below breaks on any FE type for which the elem stores no
1239  // dofs.
1240  // if (!elem->old_dof_object || !elem->old_dof_object->has_dofs(system.number()))
1241  // continue;
1242 
1243  // Examining refinement flags instead should distinguish
1244  // between refinement-added and user-added elements lacking
1245  // old_dof_object
1246  if (!elem->old_dof_object &&
1247  elem->refinement_flag() != Elem::JUST_REFINED &&
1248  elem->refinement_flag() != Elem::JUST_COARSENED)
1249  continue;
1250 
1251  const Elem * parent = elem->parent();
1252 
1253  if (elem->refinement_flag() == Elem::JUST_REFINED)
1254  {
1255  libmesh_assert(parent);
1256 
1257  // We used to hack_p_level here, but that wasn't thread-safe
1258  // so now we take p refinement flags into account in
1259  // old_dof_indices
1260 
1261  dof_map.old_dof_indices (parent, di);
1262 
1263  for (auto & node : elem->node_ref_range())
1264  {
1265  const DofObject * old_dofs = node.old_dof_object;
1266 
1267  if (old_dofs)
1268  {
1269  const unsigned int sysnum = system.number();
1270  const unsigned int nvg = old_dofs->n_var_groups(sysnum);
1271 
1272  for (unsigned int vg=0; vg != nvg; ++vg)
1273  {
1274  const unsigned int nvig =
1275  old_dofs->n_vars(sysnum, vg);
1276  for (unsigned int vig=0; vig != nvig; ++vig)
1277  {
1278  const unsigned int n_comp =
1279  old_dofs->n_comp_group(sysnum, vg);
1280  for (unsigned int c=0; c != n_comp; ++c)
1281  di.push_back
1282  (old_dofs->dof_number(sysnum, vg, vig, c, n_comp));
1283  }
1284  }
1285  }
1286  }
1287 
1288  std::sort(di.begin(), di.end());
1289  std::vector<dof_id_type>::iterator new_end =
1290  std::unique(di.begin(), di.end());
1291  std::vector<dof_id_type>(di.begin(), new_end).swap(di);
1292  }
1293  else if (elem->refinement_flag() == Elem::JUST_COARSENED)
1294  {
1295  std::vector<dof_id_type> di_child;
1296  di.clear();
1297  for (auto & child : elem->child_ref_range())
1298  {
1299  dof_map.old_dof_indices (&child, di_child);
1300  di.insert(di.end(), di_child.begin(), di_child.end());
1301  }
1302  }
1303  else
1304  dof_map.old_dof_indices (elem, di);
1305 
1306  for (auto di_i : di)
1307  if (di_i < first_old_dof || di_i >= end_old_dof)
1308  this->send_list.push_back(di_i);
1309  } // end elem loop
1310 }
1311 
1312 
1313 
1314 void BuildProjectionList::join(const BuildProjectionList & other)
1315 {
1316  // Joining simply requires I add the dof indices from the other object
1317  this->send_list.insert(this->send_list.end(),
1318  other.send_list.begin(),
1319  other.send_list.end());
1320 }
1321 #endif // LIBMESH_ENABLE_AMR
1322 
1323 
1324 
1325 void BoundaryProjectSolution::operator()(const ConstElemRange & range) const
1326 {
1327  // We need data to project
1328  libmesh_assert(f.get());
1329 
1337  // The dimensionality of the current mesh
1338  const unsigned int dim = system.get_mesh().mesh_dimension();
1339 
1340  // The DofMap for this system
1341  const DofMap & dof_map = system.get_dof_map();
1342 
1343  // Boundary info for the current mesh
1344  const BoundaryInfo & boundary_info =
1345  system.get_mesh().get_boundary_info();
1346 
1347  // The element matrix and RHS for projections.
1348  // Note that Ke is always real-valued, whereas
1349  // Fe may be complex valued if complex number
1350  // support is enabled
1351  DenseMatrix<Real> Ke;
1353  // The new element coefficients
1355 
1356 
1357  // Loop over all the variables we've been requested to project
1358  for (auto v : IntRange<std::size_t>(0, variables.size()))
1359  {
1360  const unsigned int var = variables[v];
1361 
1362  const Variable & variable = dof_map.variable(var);
1363 
1364  const FEType & fe_type = variable.type();
1365 
1366  if (fe_type.family == SCALAR)
1367  continue;
1368 
1369  const unsigned int var_component =
1370  system.variable_scalar_number(var, 0);
1371 
1372  // Get FE objects of the appropriate type
1373  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
1374 
1375  // Prepare variables for projection
1376  std::unique_ptr<QBase> qedgerule (fe_type.default_quadrature_rule(1));
1377  std::unique_ptr<QBase> qsiderule (fe_type.default_quadrature_rule(dim-1));
1378 
1379  // The values of the shape functions at the quadrature
1380  // points
1381  const std::vector<std::vector<Real>> & phi = fe->get_phi();
1382 
1383  // The gradients of the shape functions at the quadrature
1384  // points on the child element.
1385  const std::vector<std::vector<RealGradient>> * dphi = nullptr;
1386 
1387  const FEContinuity cont = fe->get_continuity();
1388 
1389  if (cont == C_ONE)
1390  {
1391  // We'll need gradient data for a C1 projection
1392  libmesh_assert(g.get());
1393 
1394  const std::vector<std::vector<RealGradient>> &
1395  ref_dphi = fe->get_dphi();
1396  dphi = &ref_dphi;
1397  }
1398 
1399  // The Jacobian * quadrature weight at the quadrature points
1400  const std::vector<Real> & JxW =
1401  fe->get_JxW();
1402 
1403  // The XYZ locations of the quadrature points
1404  const std::vector<Point> & xyz_values =
1405  fe->get_xyz();
1406 
1407  // The global DOF indices
1408  std::vector<dof_id_type> dof_indices;
1409  // Side/edge DOF indices
1410  std::vector<unsigned int> side_dofs;
1411 
1412  // Container to catch IDs passed back from BoundaryInfo.
1413  std::vector<boundary_id_type> bc_ids;
1414 
1415  // Iterate over all the elements in the range
1416  for (const auto & elem : range)
1417  {
1418  // Per-subdomain variables don't need to be projected on
1419  // elements where they're not active
1420  if (!variable.active_on_subdomain(elem->subdomain_id()))
1421  continue;
1422 
1423  const unsigned short n_nodes = elem->n_nodes();
1424  const unsigned short n_edges = elem->n_edges();
1425  const unsigned short n_sides = elem->n_sides();
1426 
1427  // Find out which nodes, edges and sides are on a requested
1428  // boundary:
1429  std::vector<bool> is_boundary_node(n_nodes, false),
1430  is_boundary_edge(n_edges, false),
1431  is_boundary_side(n_sides, false);
1432 
1433  // We also maintain a separate list of nodeset-based boundary nodes
1434  std::vector<bool> is_boundary_nodeset(n_nodes, false);
1435 
1436  for (unsigned char s=0; s != n_sides; ++s)
1437  {
1438  // First see if this side has been requested
1439  boundary_info.boundary_ids (elem, s, bc_ids);
1440  bool do_this_side = false;
1441  for (const auto & bc_id : bc_ids)
1442  if (b.count(bc_id))
1443  {
1444  do_this_side = true;
1445  break;
1446  }
1447  if (!do_this_side)
1448  continue;
1449 
1450  is_boundary_side[s] = true;
1451 
1452  // Then see what nodes and what edges are on it
1453  for (unsigned int n=0; n != n_nodes; ++n)
1454  if (elem->is_node_on_side(n,s))
1455  is_boundary_node[n] = true;
1456  for (unsigned int e=0; e != n_edges; ++e)
1457  if (elem->is_edge_on_side(e,s))
1458  is_boundary_edge[e] = true;
1459  }
1460 
1461  // We can also project on nodes, so we should also independently
1462  // check whether the nodes have been requested
1463  for (unsigned int n=0; n != n_nodes; ++n)
1464  {
1465  boundary_info.boundary_ids (elem->node_ptr(n), bc_ids);
1466 
1467  for (const auto & bc_id : bc_ids)
1468  if (b.count(bc_id))
1469  {
1470  is_boundary_node[n] = true;
1471  is_boundary_nodeset[n] = true;
1472  }
1473  }
1474 
1475  // We can also project on edges, so we should also independently
1476  // check whether the edges have been requested
1477  for (unsigned short e=0; e != n_edges; ++e)
1478  {
1479  boundary_info.edge_boundary_ids (elem, e, bc_ids);
1480 
1481  for (const auto & bc_id : bc_ids)
1482  if (b.count(bc_id))
1483  is_boundary_edge[e] = true;
1484  }
1485 
1486  // Update the DOF indices for this element based on
1487  // the current mesh
1488  dof_map.dof_indices (elem, dof_indices, var);
1489 
1490  // The number of DOFs on the element
1491  const unsigned int n_dofs =
1492  cast_int<unsigned int>(dof_indices.size());
1493 
1494  // Fixed vs. free DoFs on edge/face projections
1495  std::vector<char> dof_is_fixed(n_dofs, false); // bools
1496  std::vector<int> free_dof(n_dofs, 0);
1497 
1498  // The element type
1499  const ElemType elem_type = elem->type();
1500 
1501  // Zero the interpolated values
1502  Ue.resize (n_dofs); Ue.zero();
1503 
1504  // In general, we need a series of
1505  // projections to ensure a unique and continuous
1506  // solution. We start by interpolating boundary nodes, then
1507  // hold those fixed and project boundary edges, then hold
1508  // those fixed and project boundary faces,
1509 
1510  // Interpolate node values first
1511  unsigned int current_dof = 0;
1512  for (unsigned short n = 0; n != n_nodes; ++n)
1513  {
1514  // FIXME: this should go through the DofMap,
1515  // not duplicate dof_indices code badly!
1516  const unsigned int nc =
1517  FEInterface::n_dofs_at_node (dim, fe_type, elem_type,
1518  n);
1519  if ((!elem->is_vertex(n) || !is_boundary_node[n]) &&
1520  !is_boundary_nodeset[n])
1521  {
1522  current_dof += nc;
1523  continue;
1524  }
1525  if (cont == DISCONTINUOUS)
1526  {
1527  libmesh_assert_equal_to (nc, 0);
1528  }
1529  // Assume that C_ZERO elements have a single nodal
1530  // value shape function
1531  else if (cont == C_ZERO)
1532  {
1533  libmesh_assert_equal_to (nc, 1);
1534  Ue(current_dof) = f->component(var_component,
1535  elem->point(n),
1536  system.time);
1537  dof_is_fixed[current_dof] = true;
1538  current_dof++;
1539  }
1540  // The hermite element vertex shape functions are weird
1541  else if (fe_type.family == HERMITE)
1542  {
1543  Ue(current_dof) = f->component(var_component,
1544  elem->point(n),
1545  system.time);
1546  dof_is_fixed[current_dof] = true;
1547  current_dof++;
1548  Gradient grad = g->component(var_component,
1549  elem->point(n),
1550  system.time);
1551  // x derivative
1552  Ue(current_dof) = grad(0);
1553  dof_is_fixed[current_dof] = true;
1554  current_dof++;
1555 #if LIBMESH_DIM > 1
1556  if (dim > 1)
1557  {
1558  // We'll finite difference mixed derivatives
1559  Point nxminus = elem->point(n),
1560  nxplus = elem->point(n);
1561  nxminus(0) -= TOLERANCE;
1562  nxplus(0) += TOLERANCE;
1563  Gradient gxminus = g->component(var_component,
1564  nxminus,
1565  system.time);
1566  Gradient gxplus = g->component(var_component,
1567  nxplus,
1568  system.time);
1569  // y derivative
1570  Ue(current_dof) = grad(1);
1571  dof_is_fixed[current_dof] = true;
1572  current_dof++;
1573  // xy derivative
1574  Ue(current_dof) = (gxplus(1) - gxminus(1))
1575  / 2. / TOLERANCE;
1576  dof_is_fixed[current_dof] = true;
1577  current_dof++;
1578 
1579 #if LIBMESH_DIM > 2
1580  if (dim > 2)
1581  {
1582  // z derivative
1583  Ue(current_dof) = grad(2);
1584  dof_is_fixed[current_dof] = true;
1585  current_dof++;
1586  // xz derivative
1587  Ue(current_dof) = (gxplus(2) - gxminus(2))
1588  / 2. / TOLERANCE;
1589  dof_is_fixed[current_dof] = true;
1590  current_dof++;
1591  // We need new points for yz
1592  Point nyminus = elem->point(n),
1593  nyplus = elem->point(n);
1594  nyminus(1) -= TOLERANCE;
1595  nyplus(1) += TOLERANCE;
1596  Gradient gyminus = g->component(var_component,
1597  nyminus,
1598  system.time);
1599  Gradient gyplus = g->component(var_component,
1600  nyplus,
1601  system.time);
1602  // xz derivative
1603  Ue(current_dof) = (gyplus(2) - gyminus(2))
1604  / 2. / TOLERANCE;
1605  dof_is_fixed[current_dof] = true;
1606  current_dof++;
1607  // Getting a 2nd order xyz is more tedious
1608  Point nxmym = elem->point(n),
1609  nxmyp = elem->point(n),
1610  nxpym = elem->point(n),
1611  nxpyp = elem->point(n);
1612  nxmym(0) -= TOLERANCE;
1613  nxmym(1) -= TOLERANCE;
1614  nxmyp(0) -= TOLERANCE;
1615  nxmyp(1) += TOLERANCE;
1616  nxpym(0) += TOLERANCE;
1617  nxpym(1) -= TOLERANCE;
1618  nxpyp(0) += TOLERANCE;
1619  nxpyp(1) += TOLERANCE;
1620  Gradient gxmym = g->component(var_component,
1621  nxmym,
1622  system.time);
1623  Gradient gxmyp = g->component(var_component,
1624  nxmyp,
1625  system.time);
1626  Gradient gxpym = g->component(var_component,
1627  nxpym,
1628  system.time);
1629  Gradient gxpyp = g->component(var_component,
1630  nxpyp,
1631  system.time);
1632  Number gxzplus = (gxpyp(2) - gxmyp(2))
1633  / 2. / TOLERANCE;
1634  Number gxzminus = (gxpym(2) - gxmym(2))
1635  / 2. / TOLERANCE;
1636  // xyz derivative
1637  Ue(current_dof) = (gxzplus - gxzminus)
1638  / 2. / TOLERANCE;
1639  dof_is_fixed[current_dof] = true;
1640  current_dof++;
1641  }
1642 #endif // LIBMESH_DIM > 2
1643  }
1644 #endif // LIBMESH_DIM > 1
1645  }
1646  // Assume that other C_ONE elements have a single nodal
1647  // value shape function and nodal gradient component
1648  // shape functions
1649  else if (cont == C_ONE)
1650  {
1651  libmesh_assert_equal_to (nc, 1 + dim);
1652  Ue(current_dof) = f->component(var_component,
1653  elem->point(n),
1654  system.time);
1655  dof_is_fixed[current_dof] = true;
1656  current_dof++;
1657  Gradient grad = g->component(var_component,
1658  elem->point(n),
1659  system.time);
1660  for (unsigned int i=0; i!= dim; ++i)
1661  {
1662  Ue(current_dof) = grad(i);
1663  dof_is_fixed[current_dof] = true;
1664  current_dof++;
1665  }
1666  }
1667  else
1668  libmesh_error_msg("Unknown continuity " << cont);
1669  }
1670 
1671  // In 3D, project any edge values next
1672  if (dim > 2 && cont != DISCONTINUOUS)
1673  for (unsigned short e = 0; e != n_edges; ++e)
1674  {
1675  if (!is_boundary_edge[e])
1676  continue;
1677 
1678  FEInterface::dofs_on_edge(elem, dim, fe_type, e,
1679  side_dofs);
1680 
1681  const unsigned int n_side_dofs =
1682  cast_int<unsigned int>(side_dofs.size());
1683 
1684  // Some edge dofs are on nodes and already
1685  // fixed, others are free to calculate
1686  unsigned int free_dofs = 0;
1687  for (auto i : IntRange<unsigned int>(0, n_side_dofs))
1688  if (!dof_is_fixed[side_dofs[i]])
1689  free_dof[free_dofs++] = i;
1690 
1691  // There may be nothing to project
1692  if (!free_dofs)
1693  continue;
1694 
1695  Ke.resize (free_dofs, free_dofs); Ke.zero();
1696  Fe.resize (free_dofs); Fe.zero();
1697  // The new edge coefficients
1698  DenseVector<Number> Uedge(free_dofs);
1699 
1700  // Initialize FE data on the edge
1701  fe->attach_quadrature_rule (qedgerule.get());
1702  fe->edge_reinit (elem, e);
1703  const unsigned int n_qp = qedgerule->n_points();
1704 
1705  // Loop over the quadrature points
1706  for (unsigned int qp=0; qp<n_qp; qp++)
1707  {
1708  // solution at the quadrature point
1709  Number fineval = f->component(var_component,
1710  xyz_values[qp],
1711  system.time);
1712  // solution grad at the quadrature point
1713  Gradient finegrad;
1714  if (cont == C_ONE)
1715  finegrad = g->component(var_component,
1716  xyz_values[qp],
1717  system.time);
1718 
1719  // Form edge projection matrix
1720  for (unsigned int sidei=0, freei=0;
1721  sidei != n_side_dofs; ++sidei)
1722  {
1723  unsigned int i = side_dofs[sidei];
1724  // fixed DoFs aren't test functions
1725  if (dof_is_fixed[i])
1726  continue;
1727  for (unsigned int sidej=0, freej=0;
1728  sidej != n_side_dofs; ++sidej)
1729  {
1730  unsigned int j = side_dofs[sidej];
1731  if (dof_is_fixed[j])
1732  Fe(freei) -= phi[i][qp] * phi[j][qp] *
1733  JxW[qp] * Ue(j);
1734  else
1735  Ke(freei,freej) += phi[i][qp] *
1736  phi[j][qp] * JxW[qp];
1737  if (cont == C_ONE)
1738  {
1739  if (dof_is_fixed[j])
1740  Fe(freei) -= ((*dphi)[i][qp] *
1741  (*dphi)[j][qp]) *
1742  JxW[qp] * Ue(j);
1743  else
1744  Ke(freei,freej) += ((*dphi)[i][qp] *
1745  (*dphi)[j][qp])
1746  * JxW[qp];
1747  }
1748  if (!dof_is_fixed[j])
1749  freej++;
1750  }
1751  Fe(freei) += phi[i][qp] * fineval * JxW[qp];
1752  if (cont == C_ONE)
1753  Fe(freei) += (finegrad * (*dphi)[i][qp]) *
1754  JxW[qp];
1755  freei++;
1756  }
1757  }
1758 
1759  Ke.cholesky_solve(Fe, Uedge);
1760 
1761  // Transfer new edge solutions to element
1762  for (unsigned int i=0; i != free_dofs; ++i)
1763  {
1764  Number & ui = Ue(side_dofs[free_dof[i]]);
1766  std::abs(ui - Uedge(i)) < TOLERANCE);
1767  ui = Uedge(i);
1768  dof_is_fixed[side_dofs[free_dof[i]]] = true;
1769  }
1770  }
1771 
1772  // Project any side values (edges in 2D, faces in 3D)
1773  if (dim > 1 && cont != DISCONTINUOUS)
1774  for (unsigned short s = 0; s != n_sides; ++s)
1775  {
1776  if (!is_boundary_side[s])
1777  continue;
1778 
1779  FEInterface::dofs_on_side(elem, dim, fe_type, s,
1780  side_dofs);
1781 
1782  // Some side dofs are on nodes/edges and already
1783  // fixed, others are free to calculate
1784  unsigned int free_dofs = 0;
1785  for (auto i : index_range(side_dofs))
1786  if (!dof_is_fixed[side_dofs[i]])
1787  free_dof[free_dofs++] = i;
1788 
1789  // There may be nothing to project
1790  if (!free_dofs)
1791  continue;
1792 
1793  Ke.resize (free_dofs, free_dofs); Ke.zero();
1794  Fe.resize (free_dofs); Fe.zero();
1795  // The new side coefficients
1796  DenseVector<Number> Uside(free_dofs);
1797 
1798  // Initialize FE data on the side
1799  fe->attach_quadrature_rule (qsiderule.get());
1800  fe->reinit (elem, s);
1801  const unsigned int n_qp = qsiderule->n_points();
1802 
1803  const unsigned int n_side_dofs =
1804  cast_int<unsigned int>(side_dofs.size());
1805 
1806  // Loop over the quadrature points
1807  for (unsigned int qp=0; qp<n_qp; qp++)
1808  {
1809  // solution at the quadrature point
1810  Number fineval = f->component(var_component,
1811  xyz_values[qp],
1812  system.time);
1813  // solution grad at the quadrature point
1814  Gradient finegrad;
1815  if (cont == C_ONE)
1816  finegrad = g->component(var_component,
1817  xyz_values[qp],
1818  system.time);
1819 
1820  // Form side projection matrix
1821  for (unsigned int sidei=0, freei=0;
1822  sidei != n_side_dofs; ++sidei)
1823  {
1824  unsigned int i = side_dofs[sidei];
1825  // fixed DoFs aren't test functions
1826  if (dof_is_fixed[i])
1827  continue;
1828  for (unsigned int sidej=0, freej=0;
1829  sidej != n_side_dofs; ++sidej)
1830  {
1831  unsigned int j = side_dofs[sidej];
1832  if (dof_is_fixed[j])
1833  Fe(freei) -= phi[i][qp] * phi[j][qp] *
1834  JxW[qp] * Ue(j);
1835  else
1836  Ke(freei,freej) += phi[i][qp] *
1837  phi[j][qp] * JxW[qp];
1838  if (cont == C_ONE)
1839  {
1840  if (dof_is_fixed[j])
1841  Fe(freei) -= ((*dphi)[i][qp] *
1842  (*dphi)[j][qp]) *
1843  JxW[qp] * Ue(j);
1844  else
1845  Ke(freei,freej) += ((*dphi)[i][qp] *
1846  (*dphi)[j][qp])
1847  * JxW[qp];
1848  }
1849  if (!dof_is_fixed[j])
1850  freej++;
1851  }
1852  Fe(freei) += (fineval * phi[i][qp]) * JxW[qp];
1853  if (cont == C_ONE)
1854  Fe(freei) += (finegrad * (*dphi)[i][qp]) *
1855  JxW[qp];
1856  freei++;
1857  }
1858  }
1859 
1860  Ke.cholesky_solve(Fe, Uside);
1861 
1862  // Transfer new side solutions to element
1863  for (unsigned int i=0; i != free_dofs; ++i)
1864  {
1865  Number & ui = Ue(side_dofs[free_dof[i]]);
1867  std::abs(ui - Uside(i)) < TOLERANCE);
1868  ui = Uside(i);
1869  dof_is_fixed[side_dofs[free_dof[i]]] = true;
1870  }
1871  }
1872 
1873  const dof_id_type
1874  first = new_vector.first_local_index(),
1875  last = new_vector.last_local_index();
1876 
1877  // Lock the new_vector since it is shared among threads.
1878  {
1879  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
1880 
1881  for (unsigned int i = 0; i < n_dofs; i++)
1882  if (dof_is_fixed[i] &&
1883  (dof_indices[i] >= first) &&
1884  (dof_indices[i] < last))
1885  {
1886  new_vector.set(dof_indices[i], Ue(i));
1887  }
1888  }
1889  } // end elem loop
1890  } // end variables loop
1891 }
1892 
1893 
1894 } // namespace libMesh
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::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::WrappedFunctor
This class provides a wrapper with which to evaluate a (libMesh-style) function pointer in a Function...
Definition: wrapped_functor.h:45
libMesh::BoundaryInfo::boundary_ids
std::vector< boundary_id_type > boundary_ids(const Node *node) const
Definition: boundary_info.C:985
libMesh::FunctionBase< Number >
libMesh::BoundaryInfo
The BoundaryInfo class contains information relevant to boundary conditions including storing faces,...
Definition: boundary_info.h:57
libMesh::OldSolutionValue
The OldSolutionValue input functor class can be used with GenericProjector to read values from a solu...
Definition: generic_projector.h:655
libMesh::BoundaryProjectSolution::parameters
const Parameters & parameters
Definition: system_projection.C:205
libMesh::DofObject::n_systems
unsigned int n_systems() const
Definition: dof_object.h:861
libMesh::OldSolutionCoefs
The OldSolutionCoefs input functor class can be used with GenericProjector to read solution transfer ...
Definition: system_projection.C:485
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
n_vars
unsigned int n_vars
Definition: adaptivity_ex3.C:116
libMesh::MatrixFillAction
The MatrixFillAction output functor class can be used with GenericProjector to write solution transfe...
Definition: system_projection.C:803
libMesh::BoundaryProjectSolution::new_vector
NumericVector< Number > & new_vector
Definition: system_projection.C:206
libMesh::CompareTypes< MetaPhysicL::DynamicSparseNumberArray< T, IndexType >, T2 >::supertype
MetaPhysicL::DynamicSparseNumberArray< typename CompareTypes< T, T2 >::supertype, IndexType > supertype
Definition: system_projection.C:62
libMesh::PARALLEL
Definition: enum_parallel_type.h:36
libMesh::FEMFunctionBase< Number >
libMesh::DenseMatrix::_cholesky_back_substitute
void _cholesky_back_substitute(const DenseVector< T2 > &b, DenseVector< T2 > &x) const
Solves the equation Ax=b for the unknown value x and rhs b based on the Cholesky factorization of A.
Definition: dense_matrix_impl.h:1025
libMesh::BoundaryInfo::edge_boundary_ids
std::vector< boundary_id_type > edge_boundary_ids(const Elem *const elem, const unsigned short int edge) const
Definition: boundary_info.C:1018
libMesh::HERMITE
Definition: enum_fe_family.h:54
libMesh::SERIAL
Definition: enum_parallel_type.h:35
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::OldSolutionCoefs::OldSolutionCoefs
OldSolutionCoefs(const OldSolutionCoefs &in)
Definition: system_projection.C:496
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Elem::dim
virtual unsigned short dim() const =0
libMesh::DenseVector::zero
virtual void zero() override
Set every element in the vector to 0.
Definition: dense_vector.h:379
libMesh::FEMContext::pre_fe_reinit
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
Definition: fem_context.C:1642
libMesh::DSNAOutput< VectorValue< InnerOutput > >::type
VectorValue< DynamicSparseNumberArray< InnerOutput, dof_id_type > > type
Definition: system_projection.C:468
libMesh::NumericVector::close
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
libMesh::FEGenericBase
This class forms the foundation from which generic finite elements may be derived.
Definition: exact_error_estimator.h:39
libMesh::BoundaryProjectSolution::g
std::unique_ptr< FunctionBase< Gradient > > g
Definition: system_projection.C:204
libMesh::NumericVector::init
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
Change the dimension of the vector to n.
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::MatrixFillAction::insert
void insert(const std::vector< dof_id_type > &dof_indices, const std::vector< DynamicSparseNumberArray< ValIn, dof_id_type > > &Ue)
Definition: system_projection.C:830
libMesh::DenseMatrix< Real >
libMesh::Threads::split
Dummy "splitting object" used to distinguish splitting constructors from copy constructors.
Definition: threads_none.h:63
libMesh::DofObject::n_var_groups
unsigned int n_var_groups(const unsigned int s) const
Definition: dof_object.h:881
MetaPhysicL
Definition: dense_matrix.h:45
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::GHOSTED
Definition: enum_parallel_type.h:37
libMesh::NumericVector::clone
virtual std::unique_ptr< NumericVector< T > > clone() const =0
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::CompareTypes
Definition: compare_types.h:99
libMesh::FEMFunctionBase::component
virtual Output component(const FEMContext &, unsigned int i, const Point &p, Real time=0.)
Definition: fem_function_base.h:132
libMesh::DofMap::end_old_dof
dof_id_type end_old_dof(const processor_id_type proc) const
Definition: dof_map.h:715
libMesh::FEGenericBase::get_dphi
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:214
libMesh::MatrixFillAction::insert
void insert(dof_id_type id, const DynamicSparseNumberArray< ValIn, dof_id_type > &val)
Definition: system_projection.C:812
libMesh::DSNAOutput
Definition: system_projection.C:458
libMesh::BoundaryProjectSolution::BoundaryProjectSolution
BoundaryProjectSolution(const BoundaryProjectSolution &in)
Definition: system_projection.C:230
libMesh::DofObject
The DofObject defines an abstract base class for objects that have degrees of freedom associated with...
Definition: dof_object.h:53
libMesh::SparseMatrix< ValOut >
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::OldSolutionCoefs::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< DSNA > &values)
Definition: system_projection.C:541
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
libMesh::NumericVector::size
virtual numeric_index_type size() const =0
libMesh::VectorValue< Number >
libMesh::DofMap::old_dof_indices
void old_dof_indices(const Elem &elem, unsigned int n, std::vector< dof_id_type > &di, const unsigned int vn) const
Appends to the vector di the old global degree of freedom indices for elem.node_ref(n),...
Definition: dof_map.C:2258
libMesh::NumericVector< Number >
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::libmesh_assert
libmesh_assert(ctx)
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::BoundaryProjectSolution::system
const System & system
Definition: system_projection.C:202
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::Threads::parallel_for
void parallel_for(const Range &range, const Body &body)
Execute the provided function object in parallel on the specified range.
Definition: threads_none.h:73
libMesh::NumericVector::localize
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
libMesh::WrappedFunction
Wrap a libMesh-style function pointer into a FunctionBase object.
Definition: wrapped_function.h:51
libMesh::FEMFunctionWrapper
The FEMFunctionWrapper input functor class can be used with a GenericProjector to read values from an...
Definition: generic_projector.h:417
libMesh::DofObject::n_comp
unsigned int n_comp(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:926
libMesh::OldSolutionBase
The OldSolutionBase input functor abstract base class is the root of the OldSolutionValue and OldSolu...
Definition: generic_projector.h:479
libMesh::DofMap::variable
const Variable & variable(const unsigned int c) const
Definition: dof_map.h:1894
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::VectorSetAction
The VectorSetAction output functor class can be used with a GenericProjector to set projection values...
Definition: generic_projector.h:364
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::DofMap::SCALAR_dof_indices
void SCALAR_dof_indices(std::vector< dof_id_type > &di, const unsigned int vn, const bool old_dofs=false) const
Fills the vector di with the global degree of freedom indices corresponding to the SCALAR variable vn...
Definition: dof_map.C:2488
libMesh::Variable
This class defines the notion of a variable in the system.
Definition: variable.h:49
libMesh::DISCONTINUOUS
Definition: enum_fe_family.h:78
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::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::BoundaryProjectSolution::b
const std::set< boundary_id_type > & b
Definition: system_projection.C:200
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::ScalarTraits
Definition: compare_types.h:63
libMesh::GenericProjector
The GenericProjector class implements the core of other projection operations, using two input functo...
Definition: generic_projector.h:80
libMesh::NumericVector::local_size
virtual numeric_index_type local_size() const =0
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
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
gptr
Gradient gptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:95
libMesh::BoundaryProjectSolution::variables
const std::vector< unsigned int > & variables
Definition: system_projection.C:201
libMesh::CompareTypes::supertype
void supertype
Definition: compare_types.h:100
libMesh::NumericVector::clear
virtual void clear()
Restores the NumericVector<T> to a pristine state.
Definition: numeric_vector.h:811
libMesh::DofObject::n_vars
unsigned int n_vars(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:891
libMesh::MatrixFillAction::target_matrix
SparseMatrix< ValOut > & target_matrix
Definition: system_projection.C:806
libMesh::BuildProjectionList::unique
void unique()
Definition: system_projection.C:1198
libMesh::DenseMatrix::zero
virtual void zero() override
Set every element in the matrix to 0.
Definition: dense_matrix.h:838
libMesh::BuildProjectionList::send_list
std::vector< dof_id_type > send_list
Definition: system_projection.C:186
libMesh::numeric_index_type
dof_id_type numeric_index_type
Definition: id_types.h:99
libMesh::OldSolutionCoefs::DSNA
DSNAOutput< Output >::type DSNA
Definition: system_projection.C:488
libMesh::Elem::parent
const Elem * parent() const
Definition: elem.h:2434
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::BuildProjectionList
This class builds the send_list of old dof indices whose coefficients are needed to perform a project...
Definition: system_projection.C:167
libMesh::Elem::refinement_flag
RefinementState refinement_flag() const
Definition: elem.h:2614
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::BoundaryProjectSolution
This class implements projecting an arbitrary boundary function to the current mesh.
Definition: system_projection.C:197
libMesh::OldSolutionCoefs::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< DSNA > &values)
Definition: system_projection.C:514
libMesh::BuildProjectionList::BuildProjectionList
BuildProjectionList(const System &system_in)
Definition: system_projection.C:173
swap
void swap(Iterator &lhs, Iterator &rhs)
swap, used to implement op=
Definition: variant_filter_iterator.h:478
libMesh::ReferenceElem::get
const Elem & get(const ElemType type_in)
Definition: reference_elem.C:237
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::SparseMatrix::row_stop
virtual numeric_index_type row_stop() const =0
libMesh::NumericVector::set
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
libMesh::BuildProjectionList::BuildProjectionList
BuildProjectionList(BuildProjectionList &other, Threads::split)
Definition: system_projection.C:178
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
std
Definition: float128_shims.h:27
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::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::BuildProjectionList::system
const System & system
Definition: system_projection.C:170
libMesh::BoundaryProjectSolution::BoundaryProjectSolution
BoundaryProjectSolution(const std::set< boundary_id_type > &b_in, const std::vector< unsigned int > &variables_in, const System &system_in, FunctionBase< Number > *f_in, FunctionBase< Gradient > *g_in, const Parameters &parameters_in, NumericVector< Number > &new_v_in)
Definition: system_projection.C:209
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::NumericVector::get
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
Access multiple components at once.
Definition: numeric_vector.h:821
libMesh::SparseMatrix::set
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value)=0
Set the element (i,j) to value.
libMesh::ScalarTraits::value
static const bool value
Definition: compare_types.h:64
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::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Variable::active_on_subdomain
bool active_on_subdomain(subdomain_id_type sid) const
Definition: variable.h:136
libMesh::TypeToSend< MetaPhysicL::DynamicSparseNumberArray< T, IndexType > >::type
std::vector< std::pair< IndexType, T > > type
Definition: system_projection.C:69
libMesh::BoundaryProjectSolution::f
std::unique_ptr< FunctionBase< Number > > f
Definition: system_projection.C:203
libMesh::DSNAN
DynamicSparseNumberArray< Real, dof_id_type > DSNAN
Definition: system_projection.C:136
libMesh::FEContinuity
FEContinuity
Definition: enum_fe_family.h:77
libMesh::MatrixFillAction::MatrixFillAction
MatrixFillAction(SparseMatrix< ValOut > &target_mat)
Definition: system_projection.C:809
fptr
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:80
libMesh::DofMap::first_old_dof
dof_id_type first_old_dof(const processor_id_type proc) const
Definition: dof_map.h:660
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::NumericVector::type
ParallelType type() const
Definition: numeric_vector.h:160
int
void ErrorVector unsigned int
Definition: adjoints_ex3.C:360
libMesh::Parameters
This class provides the ability to map between arbitrary, user-defined strings and several data types...
Definition: parameters.h:59
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::OldSolutionCoefs::OldSolutionCoefs
OldSolutionCoefs(const libMesh::System &sys_in)
Definition: system_projection.C:490
libMesh::SparseMatrix::row_start
virtual numeric_index_type row_start() const =0
libMesh::DSNAOutput::type
DynamicSparseNumberArray< Output, dof_id_type > type
Definition: system_projection.C:461
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62