Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
system_projection.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 #include "libmesh/compare_types.h"
33 
34 // Template specialization declarations in here need to *precede* code
35 // using them.
36 #include "metaphysicl/dynamicsparsenumberarray_decl.h"
37 
38 using MetaPhysicL::DynamicSparseNumberArray;
39 
40 namespace libMesh
41 {
42 // From the perspective of libMesh gradient vectors, a DSNA is a
43 // scalar component
44 template <typename T, typename IndexType>
45 struct ScalarTraits<MetaPhysicL::DynamicSparseNumberArray<T,IndexType> >
46 {
47  static const bool value = true;
48 };
49 
50 // And although MetaPhysicL knows how to combine DSNA with something
51 // else, we need to teach libMesh too.
52 template <typename T, typename IndexType, typename T2>
53 struct CompareTypes<MetaPhysicL::DynamicSparseNumberArray<T,IndexType>, T2>
54 {
55  typedef typename
56  MetaPhysicL::DynamicSparseNumberArray
58 };
59 
60 template <typename T> struct TypeToSend;
61 
62 template <typename T, typename IndexType>
63 struct TypeToSend<MetaPhysicL::DynamicSparseNumberArray<T,IndexType>> {
64  typedef std::vector<std::pair<IndexType,T>> type;
65 };
66 
67 template <typename T, typename IndexType>
68 const std::vector<std::pair<IndexType,T>>
69 convert_to_send(MetaPhysicL::DynamicSparseNumberArray<T,IndexType> & in)
70 {
71  const std::size_t in_size = in.size();
72  std::vector<std::pair<IndexType,T>> returnval(in_size);
73 
74  for (std::size_t i=0; i != in_size; ++i)
75  {
76  returnval[i].first = in.raw_index(i);
77  returnval[i].second = in.raw_at(i);
78  }
79  return returnval;
80 }
81 
82 template <typename SendT, typename T, typename IndexType>
83 void convert_from_receive (SendT & received,
84  MetaPhysicL::DynamicSparseNumberArray<T,IndexType> & converted)
85 {
86  const std::size_t received_size = received.size();
87  converted.resize(received_size);
88  for (std::size_t i=0; i != received_size; ++i)
89  {
90  converted.raw_index(i) = received[i].first;
91  converted.raw_at(i) = received[i].second;
92  }
93 }
94 
95 }
96 
97 
98 #endif
99 
100 #include "libmesh/boundary_info.h"
101 #include "libmesh/dense_matrix.h"
102 #include "libmesh/dense_vector.h"
103 #include "libmesh/dof_map.h"
104 #include "libmesh/elem.h"
105 #include "libmesh/fe_base.h"
106 #include "libmesh/fe_interface.h"
107 #include "libmesh/generic_projector.h"
108 #include "libmesh/int_range.h"
109 #include "libmesh/libmesh_logging.h"
110 #include "libmesh/linear_solver.h"
111 #include "libmesh/mesh_base.h"
112 #include "libmesh/numeric_vector.h"
113 #include "libmesh/quadrature.h"
114 #include "libmesh/sparse_matrix.h"
115 #include "libmesh/system.h"
116 #include "libmesh/threads.h"
117 #include "libmesh/wrapped_function.h"
118 #include "libmesh/wrapped_functor.h"
119 #include "libmesh/fe_interface.h"
120 
121 
122 
123 #ifdef LIBMESH_HAVE_METAPHYSICL
124 // Include MetaPhysicL definitions finally
125 #include "metaphysicl/dynamicsparsenumberarray.h"
126 
127 // And make sure we instantiate the methods we'll need to use on them.
128 #include "libmesh/dense_matrix_impl.h"
129 
130 namespace libMesh {
131 typedef DynamicSparseNumberArray<Real, dof_id_type> DSNAN;
132 
133 template LIBMESH_EXPORT void
136 template LIBMESH_EXPORT void
138  DenseVector<DSNAN> &) const;
139 }
140 #endif
141 
142 
143 
144 namespace libMesh
145 {
146 
147 // ------------------------------------------------------------
148 // Helper class definitions
149 
150 #ifdef LIBMESH_ENABLE_AMR
151 
163 {
164 private:
165  const System & system;
166 
167 public:
168  BuildProjectionList (const System & system_in) :
169  system(system_in),
170  send_list()
171  {}
172 
174  system(other.system),
175  send_list()
176  {}
177 
178  void unique();
179  void operator()(const ConstElemRange & range);
180  void join (const BuildProjectionList & other);
181  std::vector<dof_id_type> send_list;
182 };
183 
184 #endif // LIBMESH_ENABLE_AMR
185 
186 
193 {
194 private:
195  const std::set<boundary_id_type> & b;
196  const std::vector<unsigned int> & variables;
197  const System & system;
198  std::unique_ptr<FunctionBase<Number>> f;
199  std::unique_ptr<FunctionBase<Gradient>> g;
202 
203 public:
204  BoundaryProjectSolution (const std::set<boundary_id_type> & b_in,
205  const std::vector<unsigned int> & variables_in,
206  const System & system_in,
207  FunctionBase<Number> * f_in,
208  FunctionBase<Gradient> * g_in,
209  const Parameters & parameters_in,
210  NumericVector<Number> & new_v_in) :
211  b(b_in),
212  variables(variables_in),
213  system(system_in),
214  f(f_in ? f_in->clone() : std::unique_ptr<FunctionBase<Number>>()),
215  g(g_in ? g_in->clone() : std::unique_ptr<FunctionBase<Gradient>>()),
216  parameters(parameters_in),
217  new_vector(new_v_in)
218  {
219  libmesh_assert(f.get());
220  f->init();
221  if (g.get())
222  g->init();
223  }
224 
226  b(in.b),
227  variables(in.variables),
228  system(in.system),
229  f(in.f.get() ? in.f->clone() : std::unique_ptr<FunctionBase<Number>>()),
230  g(in.g.get() ? in.g->clone() : std::unique_ptr<FunctionBase<Gradient>>()),
231  parameters(in.parameters),
232  new_vector(in.new_vector)
233  {
234  libmesh_assert(f.get());
235  f->init();
236  if (g.get())
237  g->init();
238  }
239 
240  void operator()(const ConstElemRange & range) const;
241 };
242 
243 
244 
245 // ------------------------------------------------------------
246 // System implementation
247 void System::project_vector (NumericVector<Number> & vector,
248  int is_adjoint) const
249 {
250  // Create a copy of the vector, which currently
251  // contains the old data.
252  std::unique_ptr<NumericVector<Number>>
253  old_vector (vector.clone());
254 
255  // Project the old vector to the new vector
256  this->project_vector (*old_vector, vector, is_adjoint);
257 }
258 
259 
265 void System::project_vector (const NumericVector<Number> & old_v,
266  NumericVector<Number> & new_v,
267  int is_adjoint) const
268 {
269  LOG_SCOPE ("project_vector(old,new)", "System");
270 
277  new_v.clear();
278 
279 #ifdef LIBMESH_ENABLE_AMR
280 
281  // Resize the new vector and get a serial version.
282  NumericVector<Number> * new_vector_ptr = nullptr;
283  std::unique_ptr<NumericVector<Number>> new_vector_built;
284  NumericVector<Number> * local_old_vector;
285  std::unique_ptr<NumericVector<Number>> local_old_vector_built;
286  const NumericVector<Number> * old_vector_ptr = nullptr;
287 
288  ConstElemRange active_local_elem_range
289  (this->get_mesh().active_local_elements_begin(),
290  this->get_mesh().active_local_elements_end());
291 
292  // If the old vector was uniprocessor, make the new
293  // vector uniprocessor
294  if (old_v.type() == SERIAL)
295  {
296  new_v.init (this->n_dofs(), false, SERIAL);
297  new_vector_ptr = &new_v;
298  old_vector_ptr = &old_v;
299  }
300 
301  // Otherwise it is a parallel, distributed vector, which
302  // we need to localize.
303  else if (old_v.type() == PARALLEL)
304  {
305  // Build a send list for efficient localization
306  BuildProjectionList projection_list(*this);
307  Threads::parallel_reduce (active_local_elem_range,
308  projection_list);
309 
310  // Create a sorted, unique send_list
311  projection_list.unique();
312 
313  new_v.init (this->n_dofs(), this->n_local_dofs(), false, PARALLEL);
314  new_vector_built = NumericVector<Number>::build(this->comm());
315  local_old_vector_built = NumericVector<Number>::build(this->comm());
316  new_vector_ptr = new_vector_built.get();
317  local_old_vector = local_old_vector_built.get();
318  new_vector_ptr->init(this->n_dofs(), this->n_local_dofs(),
319  this->get_dof_map().get_send_list(), false,
320  GHOSTED);
321  local_old_vector->init(old_v.size(), old_v.local_size(),
322  projection_list.send_list, false, GHOSTED);
323  old_v.localize(*local_old_vector, projection_list.send_list);
324  local_old_vector->close();
325  old_vector_ptr = local_old_vector;
326  }
327  else if (old_v.type() == GHOSTED)
328  {
329  // Build a send list for efficient localization
330  BuildProjectionList projection_list(*this);
331  Threads::parallel_reduce (active_local_elem_range,
332  projection_list);
333 
334  // Create a sorted, unique send_list
335  projection_list.unique();
336 
337  new_v.init (this->n_dofs(), this->n_local_dofs(),
338  this->get_dof_map().get_send_list(), false, GHOSTED);
339 
340  local_old_vector_built = NumericVector<Number>::build(this->comm());
341  new_vector_ptr = &new_v;
342  local_old_vector = local_old_vector_built.get();
343  local_old_vector->init(old_v.size(), old_v.local_size(),
344  projection_list.send_list, false, GHOSTED);
345  old_v.localize(*local_old_vector, projection_list.send_list);
346  local_old_vector->close();
347  old_vector_ptr = local_old_vector;
348  }
349  else // unknown old_v.type()
350  libmesh_error_msg("ERROR: Unknown old_v.type() == " << old_v.type());
351 
352  // Note that the above will have zeroed the new_vector.
353  // Just to be sure, assert that new_vector_ptr and old_vector_ptr
354  // were successfully set before trying to deref them.
355  libmesh_assert(new_vector_ptr);
356  libmesh_assert(old_vector_ptr);
357 
358  NumericVector<Number> & new_vector = *new_vector_ptr;
359  const NumericVector<Number> & old_vector = *old_vector_ptr;
360 
361  const unsigned int n_variables = this->n_vars();
362 
363  if (n_variables)
364  {
365  std::vector<unsigned int> vars(n_variables);
366  std::iota(vars.begin(), vars.end(), 0);
367  std::vector<unsigned int> regular_vars, vector_vars;
368  for (auto var : vars)
369  {
370  if (FEInterface::field_type(this->variable_type(var)) == TYPE_SCALAR)
371  regular_vars.push_back(var);
372  else
373  vector_vars.push_back(var);
374  }
375 
376  VectorSetAction<Number> setter(new_vector);
377 
378  if (!regular_vars.empty())
379  {
380  // Use a typedef to make the calling sequence for parallel_for() a bit more readable
381  typedef
384  Number, VectorSetAction<Number>> FEMProjector;
385 
387  f(*this, old_vector, &regular_vars);
389  g(*this, old_vector, &regular_vars);
390 
391  FEMProjector projector(*this, f, &g, setter, regular_vars);
392  projector.project(active_local_elem_range);
393  }
394 
395  if (!vector_vars.empty())
396  {
397  typedef
400  Gradient, VectorSetAction<Number>> FEMVectorProjector;
401 
402  OldSolutionValue<Gradient, &FEMContext::point_value> f_vector(*this, old_vector, &vector_vars);
403  OldSolutionValue<Tensor, &FEMContext::point_gradient> g_vector(*this, old_vector, &vector_vars);
404 
405  FEMVectorProjector vector_projector(*this, f_vector, &g_vector, setter, vector_vars);
406  vector_projector.project(active_local_elem_range);
407  }
408 
409  // Copy the SCALAR dofs from old_vector to new_vector
410  // Note: We assume that all SCALAR dofs are on the
411  // processor with highest ID
412  if (this->processor_id() == (this->n_processors()-1))
413  {
414  const DofMap & dof_map = this->get_dof_map();
415  for (auto var : make_range(this->n_vars()))
416  if (this->variable(var).type().family == SCALAR)
417  {
418  // We can just map SCALAR dofs directly across
419  std::vector<dof_id_type> new_SCALAR_indices, old_SCALAR_indices;
420  dof_map.SCALAR_dof_indices (new_SCALAR_indices, var, false);
421  dof_map.SCALAR_dof_indices (old_SCALAR_indices, var, true);
422  for (auto i : index_range(new_SCALAR_indices))
423  new_vector.set(new_SCALAR_indices[i], old_vector(old_SCALAR_indices[i]));
424  }
425  }
426  }
427 
428  new_vector.close();
429 
430  // If the old vector was serial, we probably need to send our values
431  // to other processors
432  //
433  // FIXME: I'm not sure how to make a NumericVector do that without
434  // creating a temporary parallel vector to use localize! - RHS
435  if (old_v.type() == SERIAL)
436  {
437  std::unique_ptr<NumericVector<Number>> dist_v = NumericVector<Number>::build(this->comm());
438  dist_v->init(this->n_dofs(), this->n_local_dofs(), false, PARALLEL);
439  dist_v->close();
440 
441  for (auto i : make_range(dist_v->size()))
442  if (new_vector(i) != 0.0)
443  dist_v->set(i, new_vector(i));
444 
445  dist_v->close();
446 
447  dist_v->localize (new_v, this->get_dof_map().get_send_list());
448  new_v.close();
449  }
450  // If the old vector was parallel, we need to update it
451  // and free the localized copies
452  else if (old_v.type() == PARALLEL)
453  {
454  // We may have to set dof values that this processor doesn't
455  // own in certain special cases, like LAGRANGE FIRST or
456  // HERMITE THIRD elements on second-order meshes?
457  new_v = new_vector;
458  new_v.close();
459  }
460 
461 
462  // Apply constraints only if we we are asked to
463  if(this->project_with_constraints)
464  {
465  if (is_adjoint == -1)
466  {
467  this->get_dof_map().enforce_constraints_exactly(*this, &new_v);
468  }
469  else if (is_adjoint >= 0)
470  {
471  this->get_dof_map().enforce_adjoint_constraints_exactly(new_v,
472  is_adjoint);
473  }
474  }
475 #else
476 
477  // AMR is disabled: simply copy the vector
478  new_v = old_v;
479 
480  libmesh_ignore(is_adjoint);
481 
482 #endif // #ifdef LIBMESH_ENABLE_AMR
483 }
484 
485 
486 #ifdef LIBMESH_ENABLE_AMR
487 #ifdef LIBMESH_HAVE_METAPHYSICL
488 
489 template <typename Output>
491 {
492 public:
493  typedef DynamicSparseNumberArray<Output, dof_id_type> type;
494 };
495 
496 template <typename InnerOutput>
497 class DSNAOutput<VectorValue<InnerOutput> >
498 {
499 public:
501 };
502 
512 template <typename Output,
513  void (FEMContext::*point_output) (unsigned int,
514  const Point &,
515  Output &,
516  const Real) const>
518 {
519 public:
520  typedef typename DSNAOutput<Output>::type DSNA;
523 
525  const std::vector<unsigned int> * vars) :
526  OldSolutionBase<Output, point_output>(sys_in, vars)
527  {
528  this->old_context.set_algebraic_type(FEMContext::OLD_DOFS_ONLY);
529  }
530 
532  OldSolutionBase<Output, point_output>(in.sys, in.old_context.active_vars())
533  {
534  this->old_context.set_algebraic_type(FEMContext::OLD_DOFS_ONLY);
535  }
536 
537  DSNA eval_at_node (const FEMContext & c,
538  unsigned int i,
539  unsigned int elem_dim,
540  const Node & n,
541  bool extra_hanging_dofs,
542  Real /* time */ = 0.);
543 
544  DSNA eval_at_point(const FEMContext & c,
545  unsigned int i,
546  const Point & p,
547  Real time,
548  bool skip_context_check);
549 
550  void eval_mixed_derivatives (const FEMContext & libmesh_dbg_var(c),
551  unsigned int i,
552  unsigned int dim,
553  const Node & n,
554  std::vector<DSNA> & derivs)
555  {
556  LOG_SCOPE ("eval_mixed_derivatives", "OldSolutionCoefs");
557 
558  // This should only be called on vertices
559  libmesh_assert_less(c.get_elem().get_node_index(&n),
560  c.get_elem().n_vertices());
561 
562  // Handle offset from non-scalar components in previous variables
563  libmesh_assert_less(i, this->component_to_var.size());
564  unsigned int var = this->component_to_var[i];
565 
566  // We have 1 mixed derivative in 2D, 4 in 3D
567  const unsigned int n_mixed = (dim-1) * (dim-1);
568  derivs.resize(n_mixed);
569 
570  // Be sure to handle cases where the variable wasn't defined on
571  // this node (e.g. due to changing subdomain support)
572  const DofObject * old_dof_object = n.get_old_dof_object();
573  if (old_dof_object &&
574  old_dof_object->n_vars(this->sys.number()) &&
575  old_dof_object->n_comp(this->sys.number(), var))
576  {
577  const dof_id_type first_old_id =
578  old_dof_object->dof_number(this->sys.number(), var, dim);
579  std::vector<dof_id_type> old_ids(n_mixed);
580  std::iota(old_ids.begin(), old_ids.end(), first_old_id);
581 
582  for (auto d_i : index_range(derivs))
583  {
584  derivs[d_i].resize(1);
585  derivs[d_i].raw_at(0) = 1;
586  derivs[d_i].raw_index(0) = old_ids[d_i];
587  }
588  }
589  else
590  {
591  std::fill(derivs.begin(), derivs.end(), 0);
592  }
593  }
594 
595 
596  void eval_old_dofs (const Elem & elem,
597  unsigned int node_num,
598  unsigned int var_num,
599  std::vector<dof_id_type> & indices,
600  std::vector<DSNA> & values)
601  {
602  LOG_SCOPE ("eval_old_dofs(node)", "OldSolutionCoefs");
603 
604  // We may be reusing a std::vector here, but the following
605  // dof_indices call appends without first clearing.
606  indices.clear();
607 
608  this->sys.get_dof_map().dof_indices(elem, node_num, indices, var_num);
609 
610  std::vector<dof_id_type> old_indices;
611 
612  this->sys.get_dof_map().old_dof_indices(elem, node_num, old_indices, var_num);
613 
614  libmesh_assert_equal_to (old_indices.size(), indices.size());
615 
616  values.resize(old_indices.size());
617 
618  for (auto i : index_range(values))
619  {
620  values[i].resize(1);
621  values[i].raw_at(0) = 1;
622  values[i].raw_index(0) = old_indices[i];
623  }
624  }
625 
626 
627  void eval_old_dofs (const Elem & elem,
628  const FEType & fe_type,
629  unsigned int sys_num,
630  unsigned int var_num,
631  std::vector<dof_id_type> & indices,
632  std::vector<DSNA> & values)
633  {
634  LOG_SCOPE ("eval_old_dofs(elem)", "OldSolutionCoefs");
635 
636  // We're only to be asked for old dofs on elements that can copy
637  // them through DO_NOTHING or through refinement.
638  const Elem & old_elem =
639  (elem.refinement_flag() == Elem::JUST_REFINED) ?
640  *elem.parent() : elem;
641 
642  // If there are any element-based DOF numbers, get them
643  const unsigned int nc =
644  FEInterface::n_dofs_per_elem(fe_type, &elem);
645 
646  std::vector<dof_id_type> old_dof_indices(nc);
647  indices.resize(nc);
648 
649  // We should never have fewer dofs than necessary on an
650  // element unless we're getting indices on a parent element,
651  // and we should never need those indices
652  if (nc != 0)
653  {
654  const DofObject & old_dof_object = old_elem.get_old_dof_object_ref();
655 
656  const auto [vg, vig] =
657  elem.var_to_vg_and_offset(sys_num,var_num);
658 
659  const unsigned int n_comp = elem.n_comp_group(sys_num,vg);
660  libmesh_assert_greater(elem.n_systems(), sys_num);
661  libmesh_assert_greater_equal(n_comp, nc);
662 
663  for (unsigned int i=0; i<nc; i++)
664  {
665  const dof_id_type d_old =
666  old_dof_object.dof_number(sys_num, vg, vig, i, n_comp);
667  const dof_id_type d_new =
668  elem.dof_number(sys_num, vg, vig, i, n_comp);
669  libmesh_assert_not_equal_to (d_old, DofObject::invalid_id);
670  libmesh_assert_not_equal_to (d_new, DofObject::invalid_id);
671 
672  old_dof_indices[i] = d_old;
673  indices[i] = d_new;
674  }
675  }
676 
677  values.resize(old_dof_indices.size());
678 
679  for (auto i : index_range(values))
680  {
681  values[i].resize(1);
682  values[i].raw_at(0) = 1;
683  values[i].raw_index(0) = old_dof_indices[i];
684  }
685  }
686 };
687 
688 
689 
690 template<>
691 inline
692 DynamicSparseNumberArray<Real, dof_id_type>
693 OldSolutionCoefs<Real, &FEMContext::point_value>::
694 eval_at_point(const FEMContext & c,
695  unsigned int i,
696  const Point & p,
697  Real /* time */,
698  bool skip_context_check)
699 {
700  LOG_SCOPE ("eval_at_point()", "OldSolutionCoefs");
701 
702  if (!skip_context_check)
703  if (!this->check_old_context(c, p))
704  return 0;
705 
706  // Get finite element object
707  FEGenericBase<Real> * fe = nullptr;
708  this->old_context.get_element_fe<Real>
709  (i, fe, this->old_context.get_elem_dim());
710 
711  // Build a FE for calculating phi(p)
712  FEGenericBase<Real> * fe_new =
713  this->old_context.build_new_fe(fe, p);
714 
715  // Get the values and global indices of the shape functions
716  const std::vector<std::vector<Real> > & phi = fe_new->get_phi();
717  const std::vector<dof_id_type> & dof_indices =
718  this->old_context.get_dof_indices(i);
719 
720  const std::size_t n_dofs = phi.size();
721  libmesh_assert_equal_to(n_dofs, dof_indices.size());
722 
723  DynamicSparseNumberArray<Real, dof_id_type> returnval;
724  returnval.resize(n_dofs);
725 
726  for (auto j : index_range(phi))
727  {
728  returnval.raw_at(j) = phi[j][0];
729  returnval.raw_index(j) = dof_indices[j];
730  }
731 
732  return returnval;
733 }
734 
735 
736 
737 template<>
738 inline
742  unsigned int i,
743  const Point & p,
744  Real /* time */,
745  bool skip_context_check)
746 {
747  LOG_SCOPE ("eval_at_point()", "OldSolutionCoefs");
748 
749  if (!skip_context_check)
750  if (!this->check_old_context(c, p))
751  return 0;
752 
753  // Get finite element object
754  FEGenericBase<Real> * fe = nullptr;
755  this->old_context.get_element_fe<Real>
756  (i, fe, this->old_context.get_elem_dim());
757 
758  // Build a FE for calculating phi(p)
759  FEGenericBase<Real> * fe_new =
760  this->old_context.build_new_fe(fe, p);
761 
762  // Get the values and global indices of the shape functions
763  const std::vector<std::vector<RealGradient> > & dphi = fe_new->get_dphi();
764  const std::vector<dof_id_type> & dof_indices =
765  this->old_context.get_dof_indices(i);
766 
767  const std::size_t n_dofs = dphi.size();
768  libmesh_assert_equal_to(n_dofs, dof_indices.size());
769 
771 
772  for (unsigned int d = 0; d != LIBMESH_DIM; ++d)
773  returnval(d).resize(n_dofs);
774 
775  for (auto j : index_range(dphi))
776  for (int d = 0; d != LIBMESH_DIM; ++d)
777  {
778  returnval(d).raw_at(j) = dphi[j][0](d);
779  returnval(d).raw_index(j) = dof_indices[j];
780  }
781 
782  return returnval;
783 }
784 
785 
786 template<>
787 inline
788 DynamicSparseNumberArray<Real, dof_id_type>
791  unsigned int i,
792  unsigned int /* elem_dim */,
793  const Node & n,
794  bool extra_hanging_dofs,
795  Real /* time */)
796 {
797  LOG_SCOPE ("Real eval_at_node()", "OldSolutionCoefs");
798 
799  // Optimize for the common case, where this node was part of the
800  // old solution.
801  //
802  // Be sure to handle cases where the variable wasn't defined on
803  // this node (due to changing subdomain support) or where the
804  // variable has no components on this node (due to Elem order
805  // exceeding FE order) or where the old_dof_object dofs might
806  // correspond to non-vertex dofs (due to extra_hanging_dofs and
807  // refinement)
808 
810 
811  const DofObject * old_dof_object = n.get_old_dof_object();
812  if (old_dof_object &&
813  (!extra_hanging_dofs ||
814  flag == Elem::JUST_COARSENED ||
815  flag == Elem::DO_NOTHING) &&
816  old_dof_object->n_vars(sys.number()) &&
817  old_dof_object->n_comp(sys.number(), i))
818  {
819  DynamicSparseNumberArray<Real, dof_id_type> returnval;
820  const dof_id_type old_id =
821  old_dof_object->dof_number(sys.number(), i, 0);
822  returnval.resize(1);
823  returnval.raw_at(0) = 1;
824  returnval.raw_index(0) = old_id;
825  return returnval;
826  }
827 
828  return this->eval_at_point(c, i, n, 0, false);
829 }
830 
831 
832 
833 template<>
834 inline
838  unsigned int i,
839  unsigned int elem_dim,
840  const Node & n,
841  bool extra_hanging_dofs,
842  Real /* time */)
843 {
844  LOG_SCOPE ("RealGradient eval_at_node()", "OldSolutionCoefs");
845 
846  // Optimize for the common case, where this node was part of the
847  // old solution.
848  //
849  // Be sure to handle cases where the variable wasn't defined on
850  // this node (due to changing subdomain support) or where the
851  // variable has no components on this node (due to Elem order
852  // exceeding FE order) or where the old_dof_object dofs might
853  // correspond to non-vertex dofs (due to extra_hanging_dofs and
854  // refinement)
855 
857 
858  const DofObject * old_dof_object = n.get_old_dof_object();
859  if (old_dof_object &&
860  (!extra_hanging_dofs ||
861  flag == Elem::JUST_COARSENED ||
862  flag == Elem::DO_NOTHING) &&
863  old_dof_object->n_vars(sys.number()) &&
864  old_dof_object->n_comp(sys.number(), i))
865  {
867  for (unsigned int d = 0; d != elem_dim; ++d)
868  {
869  const dof_id_type old_id =
870  old_dof_object->dof_number(sys.number(), i, d+1);
871  g(d).resize(1);
872  g(d).raw_at(0) = 1;
873  g(d).raw_index(0) = old_id;
874  }
875  return g;
876  }
877 
878  return this->eval_at_point(c, i, n, 0, false);
879 }
880 
881 
882 
891 template <typename ValIn, typename ValOut>
893 {
894 public:
895  typedef DynamicSparseNumberArray<ValIn, dof_id_type> InsertInput;
896 private:
898 
899 public:
901  target_matrix(target_mat) {}
902 
904  const DynamicSparseNumberArray<ValIn, dof_id_type> & val)
905  {
906  // Lock the target matrix since it is shared among threads.
907  {
908  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
909 
910  const std::size_t dnsa_size = val.size();
911  for (unsigned int j = 0; j != dnsa_size; ++j)
912  {
913  const dof_id_type dof_j = val.raw_index(j);
914  const ValIn dof_val = val.raw_at(j);
915  target_matrix.set(id, dof_j, dof_val);
916  }
917  }
918  }
919 
920 
921  void insert(const std::vector<dof_id_type> & dof_indices,
922  const std::vector<DynamicSparseNumberArray<ValIn, dof_id_type> > & Ue)
923  {
924  const numeric_index_type
925  begin_dof = target_matrix.row_start(),
926  end_dof = target_matrix.row_stop();
927 
928  unsigned int size = Ue.size();
929 
930  libmesh_assert_equal_to(size, dof_indices.size());
931 
932  // Lock the target matrix since it is shared among threads.
933  {
934  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
935 
936  for (unsigned int i = 0; i != size; ++i)
937  {
938  const dof_id_type dof_i = dof_indices[i];
939  if ((dof_i >= begin_dof) && (dof_i < end_dof))
940  {
941  const DynamicSparseNumberArray<ValIn,dof_id_type> & dnsa = Ue[i];
942  const std::size_t dnsa_size = dnsa.size();
943  for (unsigned int j = 0; j != dnsa_size; ++j)
944  {
945  const dof_id_type dof_j = dnsa.raw_index(j);
946  const ValIn dof_val = dnsa.raw_at(j);
947  target_matrix.set(dof_i, dof_j, dof_val);
948  }
949  }
950  }
951  }
952  }
953 };
954 
955 
956 
961 void System::projection_matrix (SparseMatrix<Number> & proj_mat) const
962 {
963  LOG_SCOPE ("projection_matrix()", "System");
964 
965  const unsigned int n_variables = this->n_vars();
966 
967  if (n_variables)
968  {
969  ConstElemRange active_local_elem_range
970  (this->get_mesh().active_local_elements_begin(),
971  this->get_mesh().active_local_elements_end());
972 
973  std::vector<unsigned int> vars(n_variables);
974  std::iota(vars.begin(), vars.end(), 0);
975 
976  // Use a typedef to make the calling sequence for parallel_for() a bit more readable
977  typedef OldSolutionCoefs<Real, &FEMContext::point_value> OldSolutionValueCoefs;
978  typedef OldSolutionCoefs<RealGradient, &FEMContext::point_gradient> OldSolutionGradientCoefs;
979 
980  typedef
981  GenericProjector<OldSolutionValueCoefs,
982  OldSolutionGradientCoefs,
983  DynamicSparseNumberArray<Real,dof_id_type>,
984  MatrixFillAction<Real, Number> > ProjMatFiller;
985 
986  OldSolutionValueCoefs f(*this, &vars);
987  OldSolutionGradientCoefs g(*this, &vars);
988  MatrixFillAction<Real, Number> setter(proj_mat);
989 
990  ProjMatFiller mat_filler(*this, f, &g, setter, vars);
991  mat_filler.project(active_local_elem_range);
992 
993  // Set the SCALAR dof transfer entries too.
994  // Note: We assume that all SCALAR dofs are on the
995  // processor with highest ID
996  if (this->processor_id() == (this->n_processors()-1))
997  {
998  const DofMap & dof_map = this->get_dof_map();
999  for (auto var : make_range(this->n_vars()))
1000  if (this->variable(var).type().family == SCALAR)
1001  {
1002  // We can just map SCALAR dofs directly across
1003  std::vector<dof_id_type> new_SCALAR_indices, old_SCALAR_indices;
1004  dof_map.SCALAR_dof_indices (new_SCALAR_indices, var, false);
1005  dof_map.SCALAR_dof_indices (old_SCALAR_indices, var, true);
1006  const unsigned int new_n_dofs =
1007  cast_int<unsigned int>(new_SCALAR_indices.size());
1008 
1009  for (unsigned int i=0; i<new_n_dofs; i++)
1010  {
1011  proj_mat.set( new_SCALAR_indices[i],
1012  old_SCALAR_indices[i], 1);
1013  }
1014  }
1015  }
1016  }
1017 }
1018 #endif // LIBMESH_HAVE_METAPHYSICL
1019 #endif // LIBMESH_ENABLE_AMR
1020 
1021 
1022 
1027 void System::project_solution (ValueFunctionPointer fptr,
1028  GradientFunctionPointer gptr,
1029  const Parameters & parameters) const
1030 {
1031  WrappedFunction<Number> f(*this, fptr, &parameters);
1032  WrappedFunction<Gradient> g(*this, gptr, &parameters);
1033  this->project_solution(&f, &g);
1034 }
1035 
1036 
1041 void System::project_solution (FunctionBase<Number> * f,
1042  FunctionBase<Gradient> * g) const
1043 {
1044  this->project_vector(*solution, f, g);
1045 
1046  solution->localize(*current_local_solution, _dof_map->get_send_list());
1047 }
1048 
1049 
1054 void System::project_solution (FEMFunctionBase<Number> * f,
1055  FEMFunctionBase<Gradient> * g) const
1056 {
1057  this->project_vector(*solution, f, g);
1058 
1059  solution->localize(*current_local_solution, _dof_map->get_send_list());
1060 }
1061 
1062 
1067 void System::project_vector (ValueFunctionPointer fptr,
1068  GradientFunctionPointer gptr,
1069  const Parameters & parameters,
1070  NumericVector<Number> & new_vector,
1071  int is_adjoint) const
1072 {
1073  WrappedFunction<Number> f(*this, fptr, &parameters);
1074  WrappedFunction<Gradient> g(*this, gptr, &parameters);
1075  this->project_vector(new_vector, &f, &g, is_adjoint);
1076 }
1077 
1082 void System::project_vector (NumericVector<Number> & new_vector,
1085  int is_adjoint) const
1086 {
1087  LOG_SCOPE ("project_vector(FunctionBase)", "System");
1088 
1089  libmesh_assert(f);
1090 
1091  WrappedFunctor<Number> f_fem(*f);
1092 
1093  if (g)
1094  {
1095  WrappedFunctor<Gradient> g_fem(*g);
1096 
1097  this->project_vector(new_vector, &f_fem, &g_fem, is_adjoint);
1098  }
1099  else
1100  this->project_vector(new_vector, &f_fem, nullptr, is_adjoint);
1101 }
1102 
1103 
1108 void System::project_vector (NumericVector<Number> & new_vector,
1111  int is_adjoint) const
1112 {
1113  LOG_SCOPE ("project_fem_vector()", "System");
1114 
1115  libmesh_assert (f);
1116 
1117  ConstElemRange active_local_range
1118  (this->get_mesh().active_local_elements_begin(),
1119  this->get_mesh().active_local_elements_end() );
1120 
1121  VectorSetAction<Number> setter(new_vector);
1122 
1123  const unsigned int n_variables = this->n_vars();
1124 
1125  std::vector<unsigned int> vars(n_variables);
1126  std::iota(vars.begin(), vars.end(), 0);
1127 
1128  // Use a typedef to make the calling sequence for parallel_for() a bit more readable
1129  typedef
1131  Number, VectorSetAction<Number>> FEMProjector;
1132 
1134 
1135  if (g)
1136  {
1138 
1139  FEMProjector projector(*this, fw, &gw, setter, vars);
1140  projector.project(active_local_range);
1141  }
1142  else
1143  {
1144  FEMProjector projector(*this, fw, nullptr, setter, vars);
1145  projector.project(active_local_range);
1146  }
1147 
1148  // Also, load values into the SCALAR dofs
1149  // Note: We assume that all SCALAR dofs are on the
1150  // processor with highest ID
1151  if (this->processor_id() == (this->n_processors()-1))
1152  {
1153  // FIXME: Do we want to first check for SCALAR vars before building this? [PB]
1154  FEMContext context( *this );
1155 
1156  const DofMap & dof_map = this->get_dof_map();
1157  for (auto var : make_range(this->n_vars()))
1158  if (this->variable(var).type().family == SCALAR)
1159  {
1160  // FIXME: We reinit with an arbitrary element in case the user
1161  // doesn't override FEMFunctionBase::component. Is there
1162  // any use case we're missing? [PB]
1163  context.pre_fe_reinit(*this, *(this->get_mesh().active_local_elements_begin()));
1164 
1165  std::vector<dof_id_type> SCALAR_indices;
1166  dof_map.SCALAR_dof_indices (SCALAR_indices, var);
1167  const unsigned int n_SCALAR_dofs =
1168  cast_int<unsigned int>(SCALAR_indices.size());
1169 
1170  for (unsigned int i=0; i<n_SCALAR_dofs; i++)
1171  {
1172  const dof_id_type global_index = SCALAR_indices[i];
1173  const unsigned int component_index =
1174  this->variable_scalar_number(var,i);
1175 
1176  new_vector.set(global_index, f->component(context, component_index, Point(), this->time));
1177  }
1178  }
1179  }
1180 
1181  new_vector.close();
1182 
1183  // Look for spline bases, in which case we need to backtrack
1184  // to calculate the spline DoF values.
1185  std::vector<const Variable *> rational_vars;
1186  for (auto varnum : vars)
1187  {
1188  const Variable & var = this->get_dof_map().variable(varnum);
1189  if (var.type().family == RATIONAL_BERNSTEIN)
1190  rational_vars.push_back(&var);
1191  }
1192 
1193  // Okay, but are we really using any *spline* bases, or just
1194  // unconstrained rational bases?
1195  bool using_spline_bases = false;
1196  if (!rational_vars.empty())
1197  {
1198  // Look for a spline node: a NodeElem with a rational variable
1199  // on it.
1200  for (auto & elem : active_local_range)
1201  if (elem->type() == NODEELEM)
1202  for (auto rational_var : rational_vars)
1203  if (rational_var->active_on_subdomain(elem->subdomain_id()))
1204  {
1205  using_spline_bases = true;
1206  goto checked_on_splines;
1207  }
1208  }
1209 
1210 checked_on_splines:
1211 
1212  // Not every processor may have a NodeElem, especially while
1213  // we're not partitioning them efficiently yet.
1214  this->comm().max(using_spline_bases);
1215 
1216  if (using_spline_bases)
1217  this->solve_for_unconstrained_dofs(new_vector, is_adjoint);
1218 
1219 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1220  if (is_adjoint == -1)
1221  this->get_dof_map().enforce_constraints_exactly(*this, &new_vector);
1222  else if (is_adjoint >= 0)
1223  this->get_dof_map().enforce_adjoint_constraints_exactly(new_vector,
1224  is_adjoint);
1225 #else
1226  libmesh_ignore(is_adjoint);
1227 #endif
1228 }
1229 
1230 
1236 void System::boundary_project_solution (const std::set<boundary_id_type> & b,
1237  const std::vector<unsigned int> & variables,
1238  ValueFunctionPointer fptr,
1239  GradientFunctionPointer gptr,
1240  const Parameters & parameters)
1241 {
1242  WrappedFunction<Number> f(*this, fptr, &parameters);
1243  WrappedFunction<Gradient> g(*this, gptr, &parameters);
1244  this->boundary_project_solution(b, variables, &f, &g);
1245 }
1246 
1247 
1253 void System::boundary_project_solution (const std::set<boundary_id_type> & b,
1254  const std::vector<unsigned int> & variables,
1257 {
1258  this->boundary_project_vector(b, variables, *solution, f, g);
1259 
1260  solution->localize(*current_local_solution);
1261 }
1262 
1263 
1264 
1265 
1266 
1271 void System::boundary_project_vector (const std::set<boundary_id_type> & b,
1272  const std::vector<unsigned int> & variables,
1273  ValueFunctionPointer fptr,
1274  GradientFunctionPointer gptr,
1275  const Parameters & parameters,
1276  NumericVector<Number> & new_vector,
1277  int is_adjoint) const
1278 {
1279  WrappedFunction<Number> f(*this, fptr, &parameters);
1280  WrappedFunction<Gradient> g(*this, gptr, &parameters);
1281  this->boundary_project_vector(b, variables, new_vector, &f, &g,
1282  is_adjoint);
1283 }
1284 
1289 void System::boundary_project_vector (const std::set<boundary_id_type> & b,
1290  const std::vector<unsigned int> & variables,
1291  NumericVector<Number> & new_vector,
1294  int is_adjoint) const
1295 {
1296  LOG_SCOPE ("boundary_project_vector()", "System");
1297 
1299  (ConstElemRange (this->get_mesh().active_local_elements_begin(),
1300  this->get_mesh().active_local_elements_end() ),
1301  BoundaryProjectSolution(b, variables, *this, f, g,
1302  this->get_equation_systems().parameters,
1303  new_vector)
1304  );
1305 
1306  // We don't do SCALAR dofs when just projecting the boundary, so
1307  // we're done here.
1308 
1309  new_vector.close();
1310 
1311 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1312  if (is_adjoint == -1)
1313  this->get_dof_map().enforce_constraints_exactly(*this, &new_vector);
1314  else if (is_adjoint >= 0)
1315  this->get_dof_map().enforce_adjoint_constraints_exactly(new_vector,
1316  is_adjoint);
1317 #else
1318  libmesh_ignore(is_adjoint);
1319 #endif
1320 }
1321 
1322 
1323 
1324 #ifdef LIBMESH_ENABLE_AMR
1325 void BuildProjectionList::unique()
1326 {
1327  // Sort the send list. After this duplicated
1328  // elements will be adjacent in the vector
1329  std::sort(this->send_list.begin(),
1330  this->send_list.end());
1331 
1332  // Now use std::unique to remove duplicate entries
1333  std::vector<dof_id_type>::iterator new_end =
1334  std::unique (this->send_list.begin(),
1335  this->send_list.end());
1336 
1337  // Remove the end of the send_list. Use the "swap trick"
1338  // from Effective STL
1339  std::vector<dof_id_type>
1340  (this->send_list.begin(), new_end).swap (this->send_list);
1341 }
1342 
1343 
1344 
1345 void BuildProjectionList::operator()(const ConstElemRange & range)
1346 {
1347  // The DofMap for this system
1348  const DofMap & dof_map = system.get_dof_map();
1349 
1350  const dof_id_type first_old_dof = dof_map.first_old_dof();
1351  const dof_id_type end_old_dof = dof_map.end_old_dof();
1352 
1353  // We can handle all the variables at once.
1354  // The old global DOF indices
1355  std::vector<dof_id_type> di;
1356 
1357  // Iterate over the elements in the range
1358  for (const auto & elem : range)
1359  {
1360  // If this element doesn't have an old_dof_object with dofs for the
1361  // current system, then it must be newly added, so the user
1362  // is responsible for setting the new dofs.
1363 
1364  // ... but we need a better way to test for that; the code
1365  // below breaks on any FE type for which the elem stores no
1366  // dofs.
1367  // if (!elem->get_old_dof_object() || !elem->get_old_dof_object()->has_dofs(system.number()))
1368  // continue;
1369 
1370  // Examining refinement flags instead should distinguish
1371  // between refinement-added and user-added elements lacking
1372  // old_dof_object
1373  const DofObject * old_dof_object = elem->get_old_dof_object();
1374  if (!old_dof_object &&
1375  elem->refinement_flag() != Elem::JUST_REFINED &&
1376  elem->refinement_flag() != Elem::JUST_COARSENED)
1377  continue;
1378 
1379  const Elem * parent = elem->parent();
1380 
1381  if (elem->refinement_flag() == Elem::JUST_REFINED)
1382  {
1383  libmesh_assert(parent);
1384 
1385  // We used to hack_p_level here, but that wasn't thread-safe
1386  // so now we take p refinement flags into account in
1387  // old_dof_indices
1388 
1389  dof_map.old_dof_indices (parent, di);
1390 
1391  for (auto & node : elem->node_ref_range())
1392  {
1393  const DofObject * old_dofs = node.get_old_dof_object();
1394 
1395  if (old_dofs)
1396  {
1397  const unsigned int sysnum = system.number();
1398  const unsigned int nvg = old_dofs->n_var_groups(sysnum);
1399 
1400  for (unsigned int vg=0; vg != nvg; ++vg)
1401  {
1402  const unsigned int nvig =
1403  old_dofs->n_vars(sysnum, vg);
1404  for (unsigned int vig=0; vig != nvig; ++vig)
1405  {
1406  const unsigned int n_comp =
1407  old_dofs->n_comp_group(sysnum, vg);
1408  for (unsigned int c=0; c != n_comp; ++c)
1409  {
1410  const dof_id_type old_id =
1411  old_dofs->dof_number(sysnum, vg, vig,
1412  c, n_comp);
1413 
1414  // We should either have no old id
1415  // (e.g. on a newly expanded subdomain)
1416  // or an id from the old system.
1417  libmesh_assert(old_id < dof_map.n_old_dofs() ||
1418  old_id == DofObject::invalid_id);
1419  di.push_back(old_id);
1420  }
1421  }
1422  }
1423  }
1424  }
1425 
1426  std::sort(di.begin(), di.end());
1427  std::vector<dof_id_type>::iterator new_end =
1428  std::unique(di.begin(), di.end());
1429  std::vector<dof_id_type>(di.begin(), new_end).swap(di);
1430  }
1431  else if (elem->refinement_flag() == Elem::JUST_COARSENED)
1432  {
1433  std::vector<dof_id_type> di_child;
1434  di.clear();
1435  for (auto & child : elem->child_ref_range())
1436  {
1437  dof_map.old_dof_indices (&child, di_child);
1438  di.insert(di.end(), di_child.begin(), di_child.end());
1439  }
1440  }
1441  else
1442  dof_map.old_dof_indices (elem, di);
1443 
1444  for (auto di_i : di)
1445  {
1446  // If we've just expanded a subdomain for a
1447  // subdomain-restricted variable, then we may have an
1448  // old_dof_object that doesn't have an old DoF for every
1449  // local index.
1450  if (di_i == DofObject::invalid_id)
1451  continue;
1452 
1453  libmesh_assert_less(di_i, dof_map.n_old_dofs());
1454  if (di_i < first_old_dof || di_i >= end_old_dof)
1455  this->send_list.push_back(di_i);
1456  }
1457  } // end elem loop
1458 }
1459 
1460 
1461 
1462 void BuildProjectionList::join(const BuildProjectionList & other)
1463 {
1464  // Joining simply requires I add the dof indices from the other object
1465  this->send_list.insert(this->send_list.end(),
1466  other.send_list.begin(),
1467  other.send_list.end());
1468 }
1469 #endif // LIBMESH_ENABLE_AMR
1470 
1471 
1472 
1473 void BoundaryProjectSolution::operator()(const ConstElemRange & range) const
1474 {
1475  // We need data to project
1476  libmesh_assert(f.get());
1477 
1485  // The dimensionality of the current mesh
1486  const unsigned int dim = system.get_mesh().mesh_dimension();
1487 
1488  // The DofMap for this system
1489  const DofMap & dof_map = system.get_dof_map();
1490 
1491  // Boundary info for the current mesh
1492  const BoundaryInfo & boundary_info =
1493  system.get_mesh().get_boundary_info();
1494 
1495  // The element matrix and RHS for projections.
1496  // Note that Ke is always real-valued, whereas
1497  // Fe may be complex valued if complex number
1498  // support is enabled
1499  DenseMatrix<Real> Ke;
1501  // The new element coefficients
1503 
1504 
1505  // Loop over all the variables we've been requested to project
1506  for (auto v : make_range(variables.size()))
1507  {
1508  const unsigned int var = variables[v];
1509 
1510  const Variable & variable = dof_map.variable(var);
1511 
1512  const FEType & fe_type = variable.type();
1513 
1514  if (fe_type.family == SCALAR)
1515  continue;
1516 
1517  const unsigned int var_component =
1518  system.variable_scalar_number(var, 0);
1519 
1520  // Get FE objects of the appropriate type
1521  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
1522 
1523  // Prepare variables for projection
1524  std::unique_ptr<QBase> qedgerule (fe_type.default_quadrature_rule(1));
1525  std::unique_ptr<QBase> qsiderule (fe_type.default_quadrature_rule(dim-1));
1526 
1527  // The values of the shape functions at the quadrature
1528  // points
1529  const std::vector<std::vector<Real>> & phi = fe->get_phi();
1530 
1531  // The gradients of the shape functions at the quadrature
1532  // points on the child element.
1533  const std::vector<std::vector<RealGradient>> * dphi = nullptr;
1534 
1535  const FEContinuity cont = fe->get_continuity();
1536 
1537  if (cont == C_ONE)
1538  {
1539  // We'll need gradient data for a C1 projection
1540  libmesh_assert(g.get());
1541 
1542  const std::vector<std::vector<RealGradient>> &
1543  ref_dphi = fe->get_dphi();
1544  dphi = &ref_dphi;
1545  }
1546 
1547  // The Jacobian * quadrature weight at the quadrature points
1548  const std::vector<Real> & JxW =
1549  fe->get_JxW();
1550 
1551  // The XYZ locations of the quadrature points
1552  const std::vector<Point> & xyz_values =
1553  fe->get_xyz();
1554 
1555  // The global DOF indices
1556  std::vector<dof_id_type> dof_indices;
1557  // Side/edge DOF indices
1558  std::vector<unsigned int> side_dofs;
1559 
1560  // Container to catch IDs passed back from BoundaryInfo.
1561  std::vector<boundary_id_type> bc_ids;
1562 
1563  // Iterate over all the elements in the range
1564  for (const auto & elem : range)
1565  {
1566  // Per-subdomain variables don't need to be projected on
1567  // elements where they're not active
1568  if (!variable.active_on_subdomain(elem->subdomain_id()))
1569  continue;
1570 
1571  const unsigned short n_nodes = elem->n_nodes();
1572  const unsigned short n_edges = elem->n_edges();
1573  const unsigned short n_sides = elem->n_sides();
1574 
1575  // Find out which nodes, edges and sides are on a requested
1576  // boundary:
1577  std::vector<bool> is_boundary_node(n_nodes, false),
1578  is_boundary_edge(n_edges, false),
1579  is_boundary_side(n_sides, false);
1580 
1581  // We also maintain a separate list of nodeset-based boundary nodes
1582  std::vector<bool> is_boundary_nodeset(n_nodes, false);
1583 
1584  for (unsigned char s=0; s != n_sides; ++s)
1585  {
1586  // First see if this side has been requested
1587  boundary_info.boundary_ids (elem, s, bc_ids);
1588  bool do_this_side = false;
1589  for (const auto & bc_id : bc_ids)
1590  if (b.count(bc_id))
1591  {
1592  do_this_side = true;
1593  break;
1594  }
1595  if (!do_this_side)
1596  continue;
1597 
1598  is_boundary_side[s] = true;
1599 
1600  // Then see what nodes and what edges are on it
1601  for (unsigned int n=0; n != n_nodes; ++n)
1602  if (elem->is_node_on_side(n,s))
1603  is_boundary_node[n] = true;
1604  for (unsigned int e=0; e != n_edges; ++e)
1605  if (elem->is_edge_on_side(e,s))
1606  is_boundary_edge[e] = true;
1607  }
1608 
1609  // We can also project on nodes, so we should also independently
1610  // check whether the nodes have been requested
1611  for (unsigned int n=0; n != n_nodes; ++n)
1612  {
1613  boundary_info.boundary_ids (elem->node_ptr(n), bc_ids);
1614 
1615  for (const auto & bc_id : bc_ids)
1616  if (b.count(bc_id))
1617  {
1618  is_boundary_node[n] = true;
1619  is_boundary_nodeset[n] = true;
1620  }
1621  }
1622 
1623  // We can also project on edges, so we should also independently
1624  // check whether the edges have been requested
1625  for (unsigned short e=0; e != n_edges; ++e)
1626  {
1627  boundary_info.edge_boundary_ids (elem, e, bc_ids);
1628 
1629  for (const auto & bc_id : bc_ids)
1630  if (b.count(bc_id))
1631  is_boundary_edge[e] = true;
1632  }
1633 
1634  // Update the DOF indices for this element based on
1635  // the current mesh
1636  dof_map.dof_indices (elem, dof_indices, var);
1637 
1638  // The number of DOFs on the element
1639  const unsigned int n_dofs =
1640  cast_int<unsigned int>(dof_indices.size());
1641 
1642  // Fixed vs. free DoFs on edge/face projections
1643  std::vector<char> dof_is_fixed(n_dofs, false); // bools
1644  std::vector<int> free_dof(n_dofs, 0);
1645 
1646  // Zero the interpolated values
1647  Ue.resize (n_dofs); Ue.zero();
1648 
1649  // In general, we need a series of
1650  // projections to ensure a unique and continuous
1651  // solution. We start by interpolating boundary nodes, then
1652  // hold those fixed and project boundary edges, then hold
1653  // those fixed and project boundary faces,
1654 
1655  // Interpolate node values first
1656  unsigned int current_dof = 0;
1657  for (unsigned short n = 0; n != n_nodes; ++n)
1658  {
1659  // FIXME: this should go through the DofMap,
1660  // not duplicate dof_indices code badly!
1661 
1662  // This call takes into account elem->p_level() internally.
1663  const unsigned int nc =
1664  FEInterface::n_dofs_at_node (fe_type, elem, n);
1665 
1666  if ((!elem->is_vertex(n) || !is_boundary_node[n]) &&
1667  !is_boundary_nodeset[n])
1668  {
1669  current_dof += nc;
1670  continue;
1671  }
1672  if (cont == DISCONTINUOUS)
1673  {
1674  libmesh_assert_equal_to (nc, 0);
1675  }
1676  // Assume that C_ZERO elements have a single nodal
1677  // value shape function
1678  else if (cont == C_ZERO)
1679  {
1680  libmesh_assert_equal_to (nc, 1);
1681  Ue(current_dof) = f->component(var_component,
1682  elem->point(n),
1683  system.time);
1684  dof_is_fixed[current_dof] = true;
1685  current_dof++;
1686  }
1687  // The hermite element vertex shape functions are weird
1688  else if (fe_type.family == HERMITE)
1689  {
1690  Ue(current_dof) = f->component(var_component,
1691  elem->point(n),
1692  system.time);
1693  dof_is_fixed[current_dof] = true;
1694  current_dof++;
1695  Gradient grad = g->component(var_component,
1696  elem->point(n),
1697  system.time);
1698  // x derivative
1699  Ue(current_dof) = grad(0);
1700  dof_is_fixed[current_dof] = true;
1701  current_dof++;
1702 #if LIBMESH_DIM > 1
1703  if (dim > 1)
1704  {
1705  // We'll finite difference mixed derivatives
1706  Point nxminus = elem->point(n),
1707  nxplus = elem->point(n);
1708  nxminus(0) -= TOLERANCE;
1709  nxplus(0) += TOLERANCE;
1710  Gradient gxminus = g->component(var_component,
1711  nxminus,
1712  system.time);
1713  Gradient gxplus = g->component(var_component,
1714  nxplus,
1715  system.time);
1716  // y derivative
1717  Ue(current_dof) = grad(1);
1718  dof_is_fixed[current_dof] = true;
1719  current_dof++;
1720  // xy derivative
1721  Ue(current_dof) = (gxplus(1) - gxminus(1))
1722  / 2. / TOLERANCE;
1723  dof_is_fixed[current_dof] = true;
1724  current_dof++;
1725 
1726 #if LIBMESH_DIM > 2
1727  if (dim > 2)
1728  {
1729  // z derivative
1730  Ue(current_dof) = grad(2);
1731  dof_is_fixed[current_dof] = true;
1732  current_dof++;
1733  // xz derivative
1734  Ue(current_dof) = (gxplus(2) - gxminus(2))
1735  / 2. / TOLERANCE;
1736  dof_is_fixed[current_dof] = true;
1737  current_dof++;
1738  // We need new points for yz
1739  Point nyminus = elem->point(n),
1740  nyplus = elem->point(n);
1741  nyminus(1) -= TOLERANCE;
1742  nyplus(1) += TOLERANCE;
1743  Gradient gyminus = g->component(var_component,
1744  nyminus,
1745  system.time);
1746  Gradient gyplus = g->component(var_component,
1747  nyplus,
1748  system.time);
1749  // xz derivative
1750  Ue(current_dof) = (gyplus(2) - gyminus(2))
1751  / 2. / TOLERANCE;
1752  dof_is_fixed[current_dof] = true;
1753  current_dof++;
1754  // Getting a 2nd order xyz is more tedious
1755  Point nxmym = elem->point(n),
1756  nxmyp = elem->point(n),
1757  nxpym = elem->point(n),
1758  nxpyp = elem->point(n);
1759  nxmym(0) -= TOLERANCE;
1760  nxmym(1) -= TOLERANCE;
1761  nxmyp(0) -= TOLERANCE;
1762  nxmyp(1) += TOLERANCE;
1763  nxpym(0) += TOLERANCE;
1764  nxpym(1) -= TOLERANCE;
1765  nxpyp(0) += TOLERANCE;
1766  nxpyp(1) += TOLERANCE;
1767  Gradient gxmym = g->component(var_component,
1768  nxmym,
1769  system.time);
1770  Gradient gxmyp = g->component(var_component,
1771  nxmyp,
1772  system.time);
1773  Gradient gxpym = g->component(var_component,
1774  nxpym,
1775  system.time);
1776  Gradient gxpyp = g->component(var_component,
1777  nxpyp,
1778  system.time);
1779  Number gxzplus = (gxpyp(2) - gxmyp(2))
1780  / 2. / TOLERANCE;
1781  Number gxzminus = (gxpym(2) - gxmym(2))
1782  / 2. / TOLERANCE;
1783  // xyz derivative
1784  Ue(current_dof) = (gxzplus - gxzminus)
1785  / 2. / TOLERANCE;
1786  dof_is_fixed[current_dof] = true;
1787  current_dof++;
1788  }
1789 #endif // LIBMESH_DIM > 2
1790  }
1791 #endif // LIBMESH_DIM > 1
1792  }
1793  // Assume that other C_ONE elements have a single nodal
1794  // value shape function and nodal gradient component
1795  // shape functions
1796  else if (cont == C_ONE)
1797  {
1798  libmesh_assert_equal_to (nc, 1 + dim);
1799  Ue(current_dof) = f->component(var_component,
1800  elem->point(n),
1801  system.time);
1802  dof_is_fixed[current_dof] = true;
1803  current_dof++;
1804  Gradient grad = g->component(var_component,
1805  elem->point(n),
1806  system.time);
1807  for (unsigned int i=0; i!= dim; ++i)
1808  {
1809  Ue(current_dof) = grad(i);
1810  dof_is_fixed[current_dof] = true;
1811  current_dof++;
1812  }
1813  }
1814  else
1815  libmesh_error_msg("Unknown continuity " << cont);
1816  }
1817 
1818  // In 3D, project any edge values next
1819  if (dim > 2 && cont != DISCONTINUOUS)
1820  for (unsigned short e = 0; e != n_edges; ++e)
1821  {
1822  if (!is_boundary_edge[e])
1823  continue;
1824 
1825  FEInterface::dofs_on_edge(elem, dim, fe_type, e,
1826  side_dofs);
1827 
1828  const unsigned int n_side_dofs =
1829  cast_int<unsigned int>(side_dofs.size());
1830 
1831  // Some edge dofs are on nodes and already
1832  // fixed, others are free to calculate
1833  unsigned int free_dofs = 0;
1834  for (auto i : make_range(n_side_dofs))
1835  if (!dof_is_fixed[side_dofs[i]])
1836  free_dof[free_dofs++] = i;
1837 
1838  // There may be nothing to project
1839  if (!free_dofs)
1840  continue;
1841 
1842  Ke.resize (free_dofs, free_dofs); Ke.zero();
1843  Fe.resize (free_dofs); Fe.zero();
1844  // The new edge coefficients
1845  DenseVector<Number> Uedge(free_dofs);
1846 
1847  // Initialize FE data on the edge
1848  fe->attach_quadrature_rule (qedgerule.get());
1849  fe->edge_reinit (elem, e);
1850  const unsigned int n_qp = qedgerule->n_points();
1851 
1852  // Loop over the quadrature points
1853  for (unsigned int qp=0; qp<n_qp; qp++)
1854  {
1855  // solution at the quadrature point
1856  Number fineval = f->component(var_component,
1857  xyz_values[qp],
1858  system.time);
1859  // solution grad at the quadrature point
1860  Gradient finegrad;
1861  if (cont == C_ONE)
1862  finegrad = g->component(var_component,
1863  xyz_values[qp],
1864  system.time);
1865 
1866  // Form edge projection matrix
1867  for (unsigned int sidei=0, freei=0;
1868  sidei != n_side_dofs; ++sidei)
1869  {
1870  unsigned int i = side_dofs[sidei];
1871  // fixed DoFs aren't test functions
1872  if (dof_is_fixed[i])
1873  continue;
1874  for (unsigned int sidej=0, freej=0;
1875  sidej != n_side_dofs; ++sidej)
1876  {
1877  unsigned int j = side_dofs[sidej];
1878  if (dof_is_fixed[j])
1879  Fe(freei) -= phi[i][qp] * phi[j][qp] *
1880  JxW[qp] * Ue(j);
1881  else
1882  Ke(freei,freej) += phi[i][qp] *
1883  phi[j][qp] * JxW[qp];
1884  if (cont == C_ONE)
1885  {
1886  if (dof_is_fixed[j])
1887  Fe(freei) -= ((*dphi)[i][qp] *
1888  (*dphi)[j][qp]) *
1889  JxW[qp] * Ue(j);
1890  else
1891  Ke(freei,freej) += ((*dphi)[i][qp] *
1892  (*dphi)[j][qp])
1893  * JxW[qp];
1894  }
1895  if (!dof_is_fixed[j])
1896  freej++;
1897  }
1898  Fe(freei) += phi[i][qp] * fineval * JxW[qp];
1899  if (cont == C_ONE)
1900  Fe(freei) += (finegrad * (*dphi)[i][qp]) *
1901  JxW[qp];
1902  freei++;
1903  }
1904  }
1905 
1906  Ke.cholesky_solve(Fe, Uedge);
1907 
1908  // Transfer new edge solutions to element
1909  for (unsigned int i=0; i != free_dofs; ++i)
1910  {
1911  Number & ui = Ue(side_dofs[free_dof[i]]);
1912  libmesh_assert(std::abs(ui) < TOLERANCE ||
1913  std::abs(ui - Uedge(i)) < TOLERANCE);
1914  ui = Uedge(i);
1915  dof_is_fixed[side_dofs[free_dof[i]]] = true;
1916  }
1917  }
1918 
1919  // Project any side values (edges in 2D, faces in 3D)
1920  if (dim > 1 && cont != DISCONTINUOUS)
1921  for (unsigned short s = 0; s != n_sides; ++s)
1922  {
1923  if (!is_boundary_side[s])
1924  continue;
1925 
1926  FEInterface::dofs_on_side(elem, dim, fe_type, s,
1927  side_dofs);
1928 
1929  // Some side dofs are on nodes/edges and already
1930  // fixed, others are free to calculate
1931  unsigned int free_dofs = 0;
1932  for (auto i : index_range(side_dofs))
1933  if (!dof_is_fixed[side_dofs[i]])
1934  free_dof[free_dofs++] = i;
1935 
1936  // There may be nothing to project
1937  if (!free_dofs)
1938  continue;
1939 
1940  Ke.resize (free_dofs, free_dofs); Ke.zero();
1941  Fe.resize (free_dofs); Fe.zero();
1942  // The new side coefficients
1943  DenseVector<Number> Uside(free_dofs);
1944 
1945  // Initialize FE data on the side
1946  fe->attach_quadrature_rule (qsiderule.get());
1947  fe->reinit (elem, s);
1948  const unsigned int n_qp = qsiderule->n_points();
1949 
1950  const unsigned int n_side_dofs =
1951  cast_int<unsigned int>(side_dofs.size());
1952 
1953  // Loop over the quadrature points
1954  for (unsigned int qp=0; qp<n_qp; qp++)
1955  {
1956  // solution at the quadrature point
1957  Number fineval = f->component(var_component,
1958  xyz_values[qp],
1959  system.time);
1960  // solution grad at the quadrature point
1961  Gradient finegrad;
1962  if (cont == C_ONE)
1963  finegrad = g->component(var_component,
1964  xyz_values[qp],
1965  system.time);
1966 
1967  // Form side projection matrix
1968  for (unsigned int sidei=0, freei=0;
1969  sidei != n_side_dofs; ++sidei)
1970  {
1971  unsigned int i = side_dofs[sidei];
1972  // fixed DoFs aren't test functions
1973  if (dof_is_fixed[i])
1974  continue;
1975  for (unsigned int sidej=0, freej=0;
1976  sidej != n_side_dofs; ++sidej)
1977  {
1978  unsigned int j = side_dofs[sidej];
1979  if (dof_is_fixed[j])
1980  Fe(freei) -= phi[i][qp] * phi[j][qp] *
1981  JxW[qp] * Ue(j);
1982  else
1983  Ke(freei,freej) += phi[i][qp] *
1984  phi[j][qp] * JxW[qp];
1985  if (cont == C_ONE)
1986  {
1987  if (dof_is_fixed[j])
1988  Fe(freei) -= ((*dphi)[i][qp] *
1989  (*dphi)[j][qp]) *
1990  JxW[qp] * Ue(j);
1991  else
1992  Ke(freei,freej) += ((*dphi)[i][qp] *
1993  (*dphi)[j][qp])
1994  * JxW[qp];
1995  }
1996  if (!dof_is_fixed[j])
1997  freej++;
1998  }
1999  Fe(freei) += (fineval * phi[i][qp]) * JxW[qp];
2000  if (cont == C_ONE)
2001  Fe(freei) += (finegrad * (*dphi)[i][qp]) *
2002  JxW[qp];
2003  freei++;
2004  }
2005  }
2006 
2007  Ke.cholesky_solve(Fe, Uside);
2008 
2009  // Transfer new side solutions to element
2010  for (unsigned int i=0; i != free_dofs; ++i)
2011  {
2012  Number & ui = Ue(side_dofs[free_dof[i]]);
2013  libmesh_assert(std::abs(ui) < TOLERANCE ||
2014  std::abs(ui - Uside(i)) < TOLERANCE);
2015  ui = Uside(i);
2016  dof_is_fixed[side_dofs[free_dof[i]]] = true;
2017  }
2018  }
2019 
2020  const dof_id_type
2021  first = new_vector.first_local_index(),
2022  last = new_vector.last_local_index();
2023 
2024  // Lock the new_vector since it is shared among threads.
2025  {
2026  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
2027 
2028  for (unsigned int i = 0; i < n_dofs; i++)
2029  if (dof_is_fixed[i] &&
2030  (dof_indices[i] >= first) &&
2031  (dof_indices[i] < last))
2032  {
2033  new_vector.set(dof_indices[i], Ue(i));
2034  }
2035  }
2036  } // end elem loop
2037  } // end variables loop
2038 }
2039 
2040 
2041 void System::solve_for_unconstrained_dofs(NumericVector<Number> & vec,
2042  int is_adjoint) const
2043 {
2044  const DofMap & dof_map = this->get_dof_map();
2045 
2046  std::unique_ptr<SparseMatrix<Number>> mat =
2047  SparseMatrix<Number>::build(this->comm());
2048 
2049  std::unique_ptr<SparsityPattern::Build> sp;
2050 
2051  if (dof_map.computed_sparsity_already())
2052  dof_map.update_sparsity_pattern(*mat);
2053  else
2054  {
2055  mat->attach_dof_map(dof_map);
2056  sp = dof_map.build_sparsity(this->get_mesh());
2057  mat->attach_sparsity_pattern(*sp);
2058  }
2059 
2060  mat->init();
2061 
2062  libmesh_assert_equal_to(vec.size(), dof_map.n_dofs());
2063  libmesh_assert_equal_to(vec.local_size(), dof_map.n_local_dofs());
2064 
2065  std::unique_ptr<NumericVector<Number>> rhs =
2066  NumericVector<Number>::build(this->comm());
2067 
2068  rhs->init(dof_map.n_dofs(), dof_map.n_local_dofs(), false,
2069  PARALLEL);
2070 
2071  // Here we start with the unconstrained (and indeterminate) linear
2072  // system, K*u = f, where K is the identity matrix for constrained
2073  // DoFs and 0 elsewhere, and f is the current solution values for
2074  // constrained DoFs and 0 elsewhere.
2075  // We then apply the usual heterogeneous constraint matrix C and
2076  // offset h, where u = C*x + h,
2077  // to get C^T*K*C*x = C^T*f - C^T*K*h
2078  // - a constrained and no-longer-singular system that finds the
2079  // closest approximation for the unconstrained degrees of freedom.
2080  //
2081  // Here, though "closest" is in an algebraic sense; we're
2082  // effectively using a pseudoinverse that optimizes in a
2083  // discretization-dependent norm. That only seems to give ~0.1%
2084  // excess error even in coarse unit test cases, but at some point it
2085  // might be reasonable to weight K and f properly.
2086 
2087  for (dof_id_type d : IntRange<dof_id_type>(dof_map.first_dof(),
2088  dof_map.end_dof()))
2089  {
2090  if (dof_map.is_constrained_dof(d))
2091  {
2092  DenseMatrix<Number> K(1,1);
2093  DenseVector<Number> F(1);
2094  std::vector<dof_id_type> dof_indices(1, d);
2095  K(0,0) = 1;
2096  F(0) = (*this->solution)(d);
2098  (K, F, dof_indices, false, is_adjoint);
2099  mat->add_matrix(K, dof_indices);
2100  rhs->add_vector(F, dof_indices);
2101  }
2102  }
2103 
2104  std::unique_ptr<LinearSolver<Number>> linear_solver =
2105  LinearSolver<Number>::build(this->comm());
2106 
2107  linear_solver->solve(*mat, vec, *rhs,
2108  double(this->get_equation_systems().parameters.get<Real>("linear solver tolerance")),
2109  this->get_equation_systems().parameters.get<unsigned int>("linear solver maximum iterations"));
2110 }
2111 
2112 
2113 } // namespace libMesh
The GenericProjector class implements the core of other projection operations, using two input functo...
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
FEFamily family
The type of finite element.
Definition: fe_type.h:221
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...
RefinementState refinement_flag() const
Definition: elem.h:3119
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)
Wrap a libMesh-style function pointer into a FunctionBase object.
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:1032
const Elem * parent() const
Definition: elem.h:2939
dof_id_type n_local_dofs() const
Definition: dof_map.h:694
OldSolutionCoefs(const libMesh::System &sys_in, const std::vector< unsigned int > *vars)
A Node is like a Point, but with more information.
Definition: node.h:52
virtual void zero() override final
Set every element in the vector to 0.
Definition: dense_vector.h:420
virtual void pre_fe_reinit(const System &, const Elem *e)
Reinitializes local data vectors/matrices on the current geometric element.
Definition: fem_context.C:1678
unsigned int n_comp(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:1002
BuildProjectionList(const System &system_in)
SparseMatrix< ValOut > & target_matrix
This class provides the ability to map between arbitrary, user-defined strings and several data types...
Definition: parameters.h:67
unsigned int n_var_groups(const unsigned int s) const
Definition: dof_object.h:957
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
Access multiple components at once.
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:2184
virtual void zero() override final
Sets all elements of the matrix to 0 and resets any decomposition flag which may have been previously...
Definition: dense_matrix.h:911
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
static constexpr Real TOLERANCE
const Elem & get_elem() const
Accessor for current Elem object.
Definition: fem_context.h:908
unsigned int dim
MetaPhysicL::DynamicSparseNumberArray< typename CompareTypes< T, T2 >::supertype, IndexType > supertype
std::unique_ptr< FunctionBase< Number > > f
The OldSolutionBase input functor abstract base class is the root of the OldSolutionValue and OldSolu...
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)
virtual numeric_index_type size() const =0
Dummy "splitting object" used to distinguish splitting constructors from copy constructors.
Definition: threads_none.h:63
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:396
DynamicSparseNumberArray< Output, dof_id_type > type
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
This class implements projecting an arbitrary boundary function to the current mesh.
RefinementState
Enumeration of possible element refinement states.
Definition: elem.h:1377
std::unique_ptr< SparsityPattern::Build > build_sparsity(const MeshBase &mesh, bool calculate_constrained=false) const
Builds a sparsity pattern for matrices using the current degree-of-freedom numbering and coupling...
Definition: dof_map.C:61
virtual std::unique_ptr< NumericVector< T > > clone() const =0
void boundary_ids(const Node *node, std::vector< boundary_id_type > &vec_to_fill) const
Fills a user-provided std::vector with the boundary ids associated with Node node.
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
The StoredRange class defines a contiguous, divisible set of objects.
Definition: stored_range.h:54
void eval_mixed_derivatives(const FEMContext &libmesh_dbg_var(c), unsigned int i, unsigned int dim, const Node &n, std::vector< DSNA > &derivs)
The libMesh namespace provides an interface to certain functionality in the library.
The OldSolutionValue input functor class can be used with GenericProjector to read values from a solu...
dof_id_type n_old_dofs() const
Definition: dof_map.h:1630
const std::vector< unsigned int > & variables
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.
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota was created back when std::iota was just an SGI STL extension.
Definition: utility.h:229
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:2552
virtual numeric_index_type row_stop() const =0
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value)=0
Set the element (i,j) to value.
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:80
This class provides a wrapper with which to evaluate a (libMesh-style) function pointer in a Function...
BuildProjectionList(BuildProjectionList &other, Threads::split)
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:178
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
Definition: fe_type.C:34
This base class can be inherited from to provide interfaces to linear solvers from different packages...
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:230
This class builds the send_list of old dof indices whose coefficients are needed to perform a project...
void libmesh_ignore(const Args &...)
const dof_id_type n_nodes
Definition: tecplot_io.C:67
dof_id_type n_dofs() const
Definition: dof_map.h:673
This class defines the notion of a variable in the system.
Definition: variable.h:49
const Variable & variable(const unsigned int c) const
Definition: dof_map.h:2284
dof_id_type numeric_index_type
Definition: id_types.h:99
unsigned int n_vars
The MatrixFillAction output functor class can be used with GenericProjector to write solution transfe...
NumericVector< Number > & new_vector
unsigned int n_vars(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:967
VectorValue< DynamicSparseNumberArray< InnerOutput, dof_id_type > > type
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
unsigned int n_systems() const
Definition: dof_object.h:937
NumberVectorValue Gradient
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
bool is_constrained_dof(const dof_id_type dof) const
Definition: dof_map.h:2352
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
void insert(dof_id_type id, const DynamicSparseNumberArray< ValIn, dof_id_type > &val)
dof_id_type end_old_dof(const processor_id_type proc) const
Definition: dof_map.h:779
DofObject * get_old_dof_object()
Pointer accessor for previously public old_dof_object.
Definition: dof_object.h:96
bool active_on_subdomain(subdomain_id_type sid) const
Definition: variable.h:157
bool computed_sparsity_already() const
Returns true iff a sparsity pattern has already been computed.
Definition: dof_map.C:289
void update_sparsity_pattern(SparseMatrix< Number > &matrix) const
Additional matrices may be be temporarily initialized by this DofMap.
Definition: dof_map.C:300
BoundaryProjectSolution(const BoundaryProjectSolution &in)
virtual void close()=0
Calls the NumericVector&#39;s internal assembly routines, ensuring that the values are consistent across ...
DofObject & get_old_dof_object_ref()
As above, but do not use in situations where the old_dof_object may be nullptr, since this function a...
Definition: dof_object.h:104
OldSolutionCoefs(const OldSolutionCoefs &in)
ParallelType type() const
void heterogenously_constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true, int qoi_index=-1) const
Definition: dof_map.h:1325
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void convert_from_receive(SendT &received, T &converted)
DynamicSparseNumberArray< Real, dof_id_type > DSNAN
std::vector< dof_id_type > send_list
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
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)
virtual numeric_index_type row_start() const =0
For ease of communication, we allow users to translate their own value types to a more easily computa...
virtual numeric_index_type local_size() const =0
The OldSolutionCoefs input functor class can be used with GenericProjector to read solution transfer ...
virtual Output component(const FEMContext &, unsigned int i, const Point &p, Real time=0.)
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
Definition: dense_matrix.h:895
Gradient gptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:95
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
The DofObject defines an abstract base class for objects that have degrees of freedom associated with...
Definition: dof_object.h:54
void parallel_reduce(const Range &range, Body &body)
Execute the provided reduction operation in parallel on the specified range.
Definition: threads_none.h:101
DynamicSparseNumberArray< ValIn, dof_id_type > InsertInput
MatrixFillAction(SparseMatrix< ValOut > &target_mat)
std::pair< unsigned int, unsigned int > var_to_vg_and_offset(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:1200
virtual void clear()
Restores the NumericVector<T> to a pristine state.
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
void edge_boundary_ids(const Elem *const elem, const unsigned short int edge, std::vector< boundary_id_type > &vec_to_fill) const
Defines a dense vector for use in Finite Element-type computations.
Definition: dof_map.h:73
void cholesky_solve(const DenseVector< T2 > &b, DenseVector< T2 > &x)
For symmetric positive definite (SPD) matrices.
dof_id_type first_dof(const processor_id_type proc) const
Definition: dof_map.h:732
The VectorSetAction output functor class can be used with a GenericProjector to set projection values...
dof_id_type end_dof(const processor_id_type proc) const
Definition: dof_map.h:756
const std::set< boundary_id_type > & b
std::unique_ptr< FunctionBase< Gradient > > g
const TypeToSend< T >::type convert_to_send(const T &in)
unsigned int n_comp_group(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:1015
static const bool value
Definition: compare_types.h:64
The FEMFunctionWrapper input functor class can be used with a GenericProjector to read values from an...
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
DSNAOutput< Output >::type DSNA
void ErrorVector unsigned int
Definition: adjoints_ex3.C:360
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
void insert(const std::vector< dof_id_type > &dof_indices, const std::vector< DynamicSparseNumberArray< ValIn, dof_id_type > > &Ue)
const Elem & get(const ElemType type_in)
This class forms the foundation from which generic finite elements may be derived.
dof_id_type first_old_dof(const processor_id_type proc) const
Definition: dof_map.h:742
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), for one variable vn.
Definition: dof_map.C:2423
uint8_t dof_id_type
Definition: id_types.h:67
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:207
const FEType & type() const
Definition: variable.h:140
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:30
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.