libMesh
mesh_function.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 // Local Includes
21 #include "libmesh/mesh_function.h"
22 #include "libmesh/dense_vector.h"
23 #include "libmesh/equation_systems.h"
24 #include "libmesh/numeric_vector.h"
25 #include "libmesh/dof_map.h"
26 #include "libmesh/point_locator_tree.h"
27 #include "libmesh/fe_base.h"
28 #include "libmesh/fe_interface.h"
29 #include "libmesh/fe_compute_data.h"
30 #include "libmesh/mesh_base.h"
31 #include "libmesh/point.h"
32 #include "libmesh/elem.h"
33 #include "libmesh/int_range.h"
34 #include "libmesh/fe_map.h"
35 
36 namespace libMesh
37 {
38 
39 
40 //------------------------------------------------------------------
41 // MeshFunction methods
43  const NumericVector<Number> & vec,
44  const DofMap & dof_map,
45  const std::vector<unsigned int> & vars,
46  const FunctionBase<Number> * master) :
47  FunctionBase<Number> (master),
48  ParallelObject (eqn_systems),
49  _eqn_systems (eqn_systems),
50  _vector (vec),
51  _dof_map (dof_map),
52  _system_vars (vars),
53  _point_locator (nullptr),
54  _out_of_mesh_mode (false),
55  _out_of_mesh_value ()
56 {
57 }
58 
59 
60 
62  const NumericVector<Number> & vec,
63  const DofMap & dof_map,
64  const unsigned int var,
65  const FunctionBase<Number> * master) :
66  FunctionBase<Number> (master),
67  ParallelObject (eqn_systems),
68  _eqn_systems (eqn_systems),
69  _vector (vec),
70  _dof_map (dof_map),
71  _system_vars (1,var),
72  _point_locator (nullptr),
73  _out_of_mesh_mode (false),
74  _out_of_mesh_value ()
75 {
76  // std::vector<unsigned int> buf (1);
77  // buf[0] = var;
78  // _system_vars (buf);
79 }
80 
81 
82 
83 
84 
85 
86 
88 {
89  delete this->_point_locator;
90 }
91 
92 
93 
94 
95 void MeshFunction::init (const Trees::BuildType /*point_locator_build_type*/)
96 {
97  // are indices of the desired variable(s) provided?
98  libmesh_assert_greater (this->_system_vars.size(), 0);
99 
100  // Don't do twice...
101  if (this->_initialized)
102  {
104  return;
105  }
106 
107  /*
108  * set up the PointLocator: currently we always get this from the
109  * MeshBase we're associated with.
110  */
111  const MeshBase & mesh = this->_eqn_systems.get_mesh();
112 
113  this->_point_locator = mesh.sub_point_locator().release();
114 
115  // ready for use
116  this->_initialized = true;
117 }
118 
119 
120 void
122 {
123  // only delete the point locator when we are the master
124  if ((this->_point_locator != nullptr) && (this->_master == nullptr))
125  {
126  delete this->_point_locator;
127  this->_point_locator = nullptr;
128  }
129  this->_initialized = false;
130 }
131 
132 
133 
134 std::unique_ptr<FunctionBase<Number>> MeshFunction::clone () const
135 {
136  MeshFunction * mf_clone =
138 
139  if (this->initialized())
140  mf_clone->init();
141  mf_clone->set_point_locator_tolerance(
142  this->get_point_locator().get_close_to_point_tol());
143 
144  return std::unique_ptr<FunctionBase<Number>>(mf_clone);
145 }
146 
147 
148 
150  const Real time)
151 {
152  libmesh_assert (this->initialized());
153 
154  DenseVector<Number> buf (1);
155  this->operator() (p, time, buf);
156  return buf(0);
157 }
158 
159 std::map<const Elem *, Number> MeshFunction::discontinuous_value (const Point & p,
160  const Real time)
161 {
162  libmesh_assert (this->initialized());
163 
164  std::map<const Elem *, DenseVector<Number>> buffer;
165  this->discontinuous_value (p, time, buffer);
166  std::map<const Elem *, Number> return_value;
167  for (const auto & pr : buffer)
168  return_value[pr.first] = pr.second(0);
169  // NOTE: If no suitable element is found, then the map return_value is empty. This
170  // puts burden on the user of this function but I don't really see a better way.
171  return return_value;
172 }
173 
175  const Real time)
176 {
177  libmesh_assert (this->initialized());
178 
179  std::vector<Gradient> buf (1);
180  this->gradient(p, time, buf);
181  return buf[0];
182 }
183 
184 std::map<const Elem *, Gradient> MeshFunction::discontinuous_gradient (const Point & p,
185  const Real time)
186 {
187  libmesh_assert (this->initialized());
188 
189  std::map<const Elem *, std::vector<Gradient>> buffer;
190  this->discontinuous_gradient (p, time, buffer);
191  std::map<const Elem *, Gradient> return_value;
192  for (const auto & pr : buffer)
193  return_value[pr.first] = pr.second[0];
194  // NOTE: If no suitable element is found, then the map return_value is empty. This
195  // puts burden on the user of this function but I don't really see a better way.
196  return return_value;
197 }
198 
199 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
201  const Real time)
202 {
203  libmesh_assert (this->initialized());
204 
205  std::vector<Tensor> buf (1);
206  this->hessian(p, time, buf);
207  return buf[0];
208 }
209 #endif
210 
212  const Real time,
213  DenseVector<Number> & output)
214 {
215  this->operator() (p, time, output, nullptr);
216 }
217 
219  const Real,
220  DenseVector<Number> & output,
221  const std::set<subdomain_id_type> * subdomain_ids)
222 {
223  libmesh_assert (this->initialized());
224 
225  const Elem * element = this->find_element(p,subdomain_ids);
226 
227  if (!element)
228  {
229  // We'd better be in out_of_mesh_mode if we couldn't find an
230  // element in the mesh
232  output = _out_of_mesh_value;
233  }
234  else
235  {
236  // resize the output vector to the number of output values
237  // that the user told us
238  output.resize (cast_int<unsigned int>
239  (this->_system_vars.size()));
240 
241 
242  {
243  const unsigned int dim = element->dim();
244 
245 
246  /*
247  * Get local coordinates to feed these into compute_data().
248  * Note that the fe_type can safely be used from the 0-variable,
249  * since the inverse mapping is the same for all FEFamilies
250  */
251  const Point mapped_point (FEMap::inverse_map (dim, element,
252  p));
253 
254  // loop over all vars
255  for (auto index : index_range(this->_system_vars))
256  {
257  /*
258  * the data for this variable
259  */
260  const unsigned int var = _system_vars[index];
261 
262  if (var == libMesh::invalid_uint)
263  {
265  index < _out_of_mesh_value.size());
266  output(index) = _out_of_mesh_value(index);
267  continue;
268  }
269 
270  const FEType & fe_type = this->_dof_map.variable_type(var);
271 
276  {
277  FEComputeData data (this->_eqn_systems, mapped_point);
278 
279  FEInterface::compute_data (dim, fe_type, element, data);
280 
281  // where the solution values for the var-th variable are stored
282  std::vector<dof_id_type> dof_indices;
283  this->_dof_map.dof_indices (element, dof_indices, var);
284 
285  // interpolate the solution
286  {
287  Number value = 0.;
288 
289  for (auto i : index_range(dof_indices))
290  value += this->_vector(dof_indices[i]) * data.shape[i];
291 
292  output(index) = value;
293  }
294 
295  }
296 
297  // next variable
298  }
299  }
300  }
301 }
302 
303 
305  const Real time,
306  std::map<const Elem *, DenseVector<Number>> & output)
307 {
308  this->discontinuous_value (p, time, output, nullptr);
309 }
310 
311 
312 
314  const Real,
315  std::map<const Elem *, DenseVector<Number>> & output,
316  const std::set<subdomain_id_type> * subdomain_ids)
317 {
318  libmesh_assert (this->initialized());
319 
320  // clear the output map
321  output.clear();
322 
323  // get the candidate elements
324  std::set<const Elem *> candidate_element = this->find_elements(p,subdomain_ids);
325 
326  // loop through all candidates, if the set is empty this function will return an
327  // empty map
328  for (const auto & element : candidate_element)
329  {
330  const unsigned int dim = element->dim();
331 
332  // define a temporary vector to store all values
333  DenseVector<Number> temp_output (cast_int<unsigned int>(this->_system_vars.size()));
334 
335  /*
336  * Get local coordinates to feed these into compute_data().
337  * Note that the fe_type can safely be used from the 0-variable,
338  * since the inverse mapping is the same for all FEFamilies
339  */
340  const Point mapped_point (FEMap::inverse_map (dim, element, p));
341 
342  // loop over all vars
343  for (auto index : index_range(this->_system_vars))
344  {
345  /*
346  * the data for this variable
347  */
348  const unsigned int var = _system_vars[index];
349 
350  if (var == libMesh::invalid_uint)
351  {
353  index < _out_of_mesh_value.size());
354  temp_output(index) = _out_of_mesh_value(index);
355  continue;
356  }
357 
358  const FEType & fe_type = this->_dof_map.variable_type(var);
359 
364  {
365  FEComputeData data (this->_eqn_systems, mapped_point);
366 
367  FEInterface::compute_data (dim, fe_type, element, data);
368 
369  // where the solution values for the var-th variable are stored
370  std::vector<dof_id_type> dof_indices;
371  this->_dof_map.dof_indices (element, dof_indices, var);
372 
373  // interpolate the solution
374  {
375  Number value = 0.;
376 
377  for (auto i : index_range(dof_indices))
378  value += this->_vector(dof_indices[i]) * data.shape[i];
379 
380  temp_output(index) = value;
381  }
382 
383  }
384 
385  // next variable
386  }
387 
388  // Insert temp_output into output
389  output[element] = temp_output;
390  }
391 }
392 
393 
394 
396  const Real,
397  std::vector<Gradient> & output,
398  const std::set<subdomain_id_type> * subdomain_ids)
399 {
400  libmesh_assert (this->initialized());
401 
402  const Elem * element = this->find_element(p,subdomain_ids);
403 
404  if (!element)
405  {
406  output.resize(0);
407  return;
408  }
409  else
410  {
411  // resize the output vector to the number of output values
412  // that the user told us
413  output.resize (this->_system_vars.size());
414 
415 
416  {
417  const unsigned int dim = element->dim();
418 
419 
420  /*
421  * Get local coordinates to feed these into compute_data().
422  * Note that the fe_type can safely be used from the 0-variable,
423  * since the inverse mapping is the same for all FEFamilies
424  */
425  const Point mapped_point (FEMap::inverse_map (dim, element,
426  p));
427 
428  std::vector<Point> point_list (1, mapped_point);
429 
430  // loop over all vars
431  for (auto index : index_range(this->_system_vars))
432  {
433  /*
434  * the data for this variable
435  */
436  const unsigned int var = _system_vars[index];
437 
438  if (var == libMesh::invalid_uint)
439  {
441  index < _out_of_mesh_value.size());
442  output[index] = Gradient(_out_of_mesh_value(index));
443  continue;
444  }
445 
446  const FEType & fe_type = this->_dof_map.variable_type(var);
447 
448  // where the solution values for the var-th variable are stored
449  std::vector<dof_id_type> dof_indices;
450  this->_dof_map.dof_indices (element, dof_indices, var);
451 
452  // interpolate the solution
453  Gradient grad(0.);
454 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
455  //The other algorithm works in case of finite elements as well,
456  //but this one is faster.
457  if (!element->infinite())
458  {
459 #endif
460  std::unique_ptr<FEBase> point_fe (FEBase::build(dim, fe_type));
461  const std::vector<std::vector<RealGradient>> & dphi = point_fe->get_dphi();
462  point_fe->reinit(element, &point_list);
463 
464  for (auto i : index_range(dof_indices))
465  grad.add_scaled(dphi[i][0], this->_vector(dof_indices[i]));
466 
467 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
468  }
469  else
470  {
475  FEComputeData data (this->_eqn_systems, mapped_point);
476  data.enable_derivative();
477  FEInterface::compute_data (dim, fe_type, element, data);
478  //grad [x] = data.dshape[i](v) * dv/dx * dof_index [i]
479  // sum over all indices
480  for (auto i : index_range(dof_indices))
481  {
482  // local coordinates
483  for (std::size_t v=0; v<dim; v++)
484  for (std::size_t xyz=0; xyz<LIBMESH_DIM; xyz++)
485  {
486  // FIXME: this needs better syntax: It is matrix-vector multiplication.
487  grad(xyz) += data.local_transform[v][xyz]
488  * data.dshape[i](v)
489  * this->_vector(dof_indices[i]);
490  }
491  }
492  }
493 #endif
494  output[index] = grad;
495  }
496  }
497  }
498 }
499 
500 
502  const Real time,
503  std::map<const Elem *, std::vector<Gradient>> & output)
504 {
505  this->discontinuous_gradient (p, time, output, nullptr);
506 }
507 
508 
509 
511  const Real,
512  std::map<const Elem *, std::vector<Gradient>> & output,
513  const std::set<subdomain_id_type> * subdomain_ids)
514 {
515  libmesh_assert (this->initialized());
516 
517  // clear the output map
518  output.clear();
519 
520  // get the candidate elements
521  std::set<const Elem *> candidate_element = this->find_elements(p,subdomain_ids);
522 
523  // loop through all candidates, if the set is empty this function will return an
524  // empty map
525  for (const auto & element : candidate_element)
526  {
527  const unsigned int dim = element->dim();
528 
529  // define a temporary vector to store all values
530  std::vector<Gradient> temp_output (cast_int<unsigned int>(this->_system_vars.size()));
531 
532  /*
533  * Get local coordinates to feed these into compute_data().
534  * Note that the fe_type can safely be used from the 0-variable,
535  * since the inverse mapping is the same for all FEFamilies
536  */
537  const Point mapped_point (FEMap::inverse_map (dim, element, p));
538 
539 
540  // loop over all vars
541  std::vector<Point> point_list (1, mapped_point);
542  for (auto index : index_range(this->_system_vars))
543  {
544  /*
545  * the data for this variable
546  */
547  const unsigned int var = _system_vars[index];
548 
549  if (var == libMesh::invalid_uint)
550  {
552  index < _out_of_mesh_value.size());
553  temp_output[index] = Gradient(_out_of_mesh_value(index));
554  continue;
555  }
556 
557  const FEType & fe_type = this->_dof_map.variable_type(var);
558 
559  // where the solution values for the var-th variable are stored
560  std::vector<dof_id_type> dof_indices;
561  this->_dof_map.dof_indices (element, dof_indices, var);
562 
563  Gradient grad(0.);
564 
565  // for performance-reasons, we use different algorithms now.
566  // TODO: Check that both give the same result for finite elements.
567  // Otherwive it is wrong...
568 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
569  if (!element->infinite())
570  {
571 #endif
572  std::unique_ptr<FEBase> point_fe (FEBase::build(dim, fe_type));
573  const std::vector<std::vector<RealGradient>> & dphi = point_fe->get_dphi();
574  point_fe->reinit(element, & point_list);
575 
576  for (auto i : index_range(dof_indices))
577  grad.add_scaled(dphi[i][0], this->_vector(dof_indices[i]));
578 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
579  }
580  else
581  {
586  //TODO: enable this for a vector of points as well...
587  FEComputeData data (this->_eqn_systems, mapped_point);
588  data.enable_derivative();
589  FEInterface::compute_data (dim, fe_type, element, data);
590 
591  //grad [x] = data.dshape[i](v) * dv/dx * dof_index [i]
592  // sum over all indices
593  for (auto i : index_range(dof_indices))
594  {
595  // local coordinates.
596  for (std::size_t v=0; v<dim; v++)
597  for (std::size_t xyz=0; xyz<LIBMESH_DIM; xyz++)
598  {
599  // FIXME: this needs better syntax: It is matrix-vector multiplication.
600  grad(xyz) += data.local_transform[v][xyz]
601  * data.dshape[i](v)
602  * this->_vector(dof_indices[i]);
603  }
604  }
605  }
606 #endif
607  temp_output[index] = grad;
608 
609  // next variable
610  }
611 
612  // Insert temp_output into output
613  output[element] = temp_output;
614  }
615 }
616 
617 
618 
619 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
620 void MeshFunction::hessian (const Point & p,
621  const Real,
622  std::vector<Tensor> & output,
623  const std::set<subdomain_id_type> * subdomain_ids)
624 {
625  libmesh_assert (this->initialized());
626 
627  const Elem * element = this->find_element(p,subdomain_ids);
628 
629  if (!element)
630  {
631  output.resize(0);
632  }
633  else
634  {
635 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
636  if(element->infinite())
637  libmesh_warning("Warning: Requested the Hessian of an Infinite element."
638  << "Second derivatives for Infinite elements"
639  << " are not yet implemented!"
640  << std::endl);
641 #endif
642  // resize the output vector to the number of output values
643  // that the user told us
644  output.resize (this->_system_vars.size());
645 
646 
647  {
648  const unsigned int dim = element->dim();
649 
650 
651  /*
652  * Get local coordinates to feed these into compute_data().
653  * Note that the fe_type can safely be used from the 0-variable,
654  * since the inverse mapping is the same for all FEFamilies
655  */
656  const Point mapped_point (FEMap::inverse_map (dim, element,
657  p));
658 
659  std::vector<Point> point_list (1, mapped_point);
660 
661  // loop over all vars
662  for (auto index : index_range(this->_system_vars))
663  {
664  /*
665  * the data for this variable
666  */
667  const unsigned int var = _system_vars[index];
668 
669  if (var == libMesh::invalid_uint)
670  {
672  index < _out_of_mesh_value.size());
673  output[index] = Tensor(_out_of_mesh_value(index));
674  continue;
675  }
676  const FEType & fe_type = this->_dof_map.variable_type(var);
677 
678  std::unique_ptr<FEBase> point_fe (FEBase::build(dim, fe_type));
679  const std::vector<std::vector<RealTensor>> & d2phi =
680  point_fe->get_d2phi();
681  point_fe->reinit(element, &point_list);
682 
683  // where the solution values for the var-th variable are stored
684  std::vector<dof_id_type> dof_indices;
685  this->_dof_map.dof_indices (element, dof_indices, var);
686 
687  // interpolate the solution
688  Tensor hess;
689 
690  for (auto i : index_range(dof_indices))
691  hess.add_scaled(d2phi[i][0], this->_vector(dof_indices[i]));
692 
693  output[index] = hess;
694  }
695  }
696  }
697 }
698 #endif
699 
701  const std::set<subdomain_id_type> * subdomain_ids) const
702 {
703  /* Ensure that in the case of a master mesh function, the
704  out-of-mesh mode is enabled either for both or for none. This is
705  important because the out-of-mesh mode is also communicated to
706  the point locator. Since this is time consuming, enable it only
707  in debug mode. */
708 #ifdef DEBUG
709  if (this->_master != nullptr)
710  {
711  const MeshFunction * master =
712  cast_ptr<const MeshFunction *>(this->_master);
714  libmesh_error_msg("ERROR: If you use out-of-mesh-mode in connection with master mesh " \
715  << "functions, you must enable out-of-mesh mode for both the master and the slave mesh function.");
716  }
717 #endif
718 
719  // locate the point in the other mesh
720  const Elem * element = this->_point_locator->operator()(p,subdomain_ids);
721 
722  // If we have an element, but it's not a local element, then we
723  // either need to have a serialized vector or we need to find a
724  // local element sharing the same point.
725  if (element &&
726  (element->processor_id() != this->processor_id()) &&
727  _vector.type() != SERIAL)
728  {
729  // look for a local element containing the point
730  std::set<const Elem *> point_neighbors;
731  element->find_point_neighbors(p, point_neighbors);
732  element = nullptr;
733  for (const auto & elem : point_neighbors)
734  if (elem->processor_id() == this->processor_id())
735  {
736  element = elem;
737  break;
738  }
739  }
740 
741  return element;
742 }
743 
744 std::set<const Elem *> MeshFunction::find_elements(const Point & p,
745  const std::set<subdomain_id_type> * subdomain_ids) const
746 {
747  /* Ensure that in the case of a master mesh function, the
748  out-of-mesh mode is enabled either for both or for none. This is
749  important because the out-of-mesh mode is also communicated to
750  the point locator. Since this is time consuming, enable it only
751  in debug mode. */
752 #ifdef DEBUG
753  if (this->_master != nullptr)
754  {
755  const MeshFunction * master =
756  cast_ptr<const MeshFunction *>(this->_master);
758  libmesh_error_msg("ERROR: If you use out-of-mesh-mode in connection with master mesh " \
759  << "functions, you must enable out-of-mesh mode for both the master and the slave mesh function.");
760  }
761 #endif
762 
763  // locate the point in the other mesh
764  std::set<const Elem *> candidate_elements;
765  std::set<const Elem *> final_candidate_elements;
766  this->_point_locator->operator()(p,candidate_elements,subdomain_ids);
767  for (const auto & element : candidate_elements)
768  {
769  // If we have an element, but it's not a local element, then we
770  // either need to have a serialized vector or we need to find a
771  // local element sharing the same point.
772  if (element &&
773  (element->processor_id() != this->processor_id()) &&
774  _vector.type() != SERIAL)
775  {
776  // look for a local element containing the point
777  std::set<const Elem *> point_neighbors;
778  element->find_point_neighbors(p, point_neighbors);
779  for (const auto & elem : point_neighbors)
780  if (elem->processor_id() == this->processor_id())
781  {
782  final_candidate_elements.insert(elem);
783  break;
784  }
785  }
786  else
787  final_candidate_elements.insert(element);
788  }
789 
790  return final_candidate_elements;
791 }
792 
794 {
795  libmesh_assert (this->initialized());
796  return *_point_locator;
797 }
798 
800 {
801  libmesh_assert (this->initialized());
803  _out_of_mesh_mode = true;
805 }
806 
808 {
809  DenseVector<Number> v(1);
810  v(0) = value;
811  this->enable_out_of_mesh_mode(v);
812 }
813 
815 {
816  libmesh_assert (this->initialized());
818  _out_of_mesh_mode = false;
819 }
820 
822 {
824 }
825 
827 {
829 }
830 
831 } // namespace libMesh
libMesh::MeshFunction::find_element
const Elem * find_element(const Point &p, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Helper function to reduce code duplication.
Definition: mesh_function.C:700
libMesh::Elem::find_point_neighbors
void find_point_neighbors(const Point &p, std::set< const Elem * > &neighbor_set) const
This function finds all active elements (including this one) which are in the same manifold as this e...
Definition: elem.C:560
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::FunctionBase< Number >
libMesh::invalid_uint
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value.
Definition: libmesh.h:249
libMesh::MeshFunction::_out_of_mesh_value
DenseVector< Number > _out_of_mesh_value
Value to return outside the mesh if out-of-mesh mode is enabled.
Definition: mesh_function.h:363
libMesh::FEComputeData
class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEIn...
Definition: fe_compute_data.h:51
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::DofMap::dof_indices
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
Definition: dof_map.C:1967
libMesh::MeshFunction::init
virtual void init() override
Override the FunctionBase::init() member function by calling our own and specifying the Trees::NODES ...
Definition: mesh_function.h:110
libMesh::MeshFunction::operator()
Number operator()(const Point &p, const Real time=0.) override
Definition: mesh_function.C:149
libMesh::MeshFunction
This class provides function-like objects for data distributed over a mesh.
Definition: mesh_function.h:53
libMesh::MeshFunction::find_elements
std::set< const Elem * > find_elements(const Point &p, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Definition: mesh_function.C:744
libMesh::PointLocatorBase::set_close_to_point_tol
virtual void set_close_to_point_tol(Real close_to_point_tol)
Set a tolerance to use when determining if a point is contained within the mesh.
Definition: point_locator_base.C:92
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
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::MeshFunction::discontinuous_gradient
std::map< const Elem *, Gradient > discontinuous_gradient(const Point &p, const Real time=0.)
Definition: mesh_function.C:184
libMesh::Elem::infinite
virtual bool infinite() const =0
libMesh::MeshFunction::clone
virtual std::unique_ptr< FunctionBase< Number > > clone() const override
Definition: mesh_function.C:134
libMesh::MeshFunction::discontinuous_value
std::map< const Elem *, Number > discontinuous_value(const Point &p, const Real time=0.)
Definition: mesh_function.C:159
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::PointLocatorBase::disable_out_of_mesh_mode
virtual void disable_out_of_mesh_mode()=0
Disables out-of-mesh mode (default).
libMesh::Trees::BuildType
BuildType
enum defining how to build the tree.
Definition: tree_base.h:55
libMesh::MeshFunction::hessian
Tensor hessian(const Point &p, const Real time=0.)
Definition: mesh_function.C:200
libMesh::DofObject::processor_id
processor_id_type processor_id() const
Definition: dof_object.h:829
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::MeshFunction::get_point_locator
const PointLocatorBase & get_point_locator(void) const
Definition: mesh_function.C:793
libMesh::FunctionBase< Number >::initialized
bool initialized() const
Definition: function_base.h:205
libMesh::VectorValue< Number >
libMesh::NumericVector< Number >
libMesh::MeshFunction::set_point_locator_tolerance
void set_point_locator_tolerance(Real tol)
We may want to specify a tolerance for the PointLocator to use, since in some cases the point we want...
Definition: mesh_function.C:821
libMesh::PointLocatorBase::enable_out_of_mesh_mode
virtual void enable_out_of_mesh_mode()=0
Enables out-of-mesh mode.
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::MeshFunction::enable_out_of_mesh_mode
void enable_out_of_mesh_mode(const DenseVector< Number > &value)
Enables out-of-mesh mode.
Definition: mesh_function.C:799
libMesh::FEMap::inverse_map
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe_map.C:1622
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::FunctionBase< Number >::_initialized
bool _initialized
When init() was called so that everything is ready for calls to operator() (...), then this bool is t...
Definition: function_base.h:179
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libMesh::TensorValue
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
Definition: exact_solution.h:53
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::TypeVector::add_scaled
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
Definition: type_vector.h:665
libMesh::MeshFunction::_point_locator
PointLocatorBase * _point_locator
A point locator is needed to locate the points in the mesh.
Definition: mesh_function.h:351
libMesh::TypeTensor::add_scaled
void add_scaled(const TypeTensor< T2 > &, const T &)
Add a scaled tensor to this tensor without creating a temporary.
Definition: type_tensor.h:869
libMesh::MeshFunction::_eqn_systems
const EquationSystems & _eqn_systems
The equation systems handler, from which the data are gathered.
Definition: mesh_function.h:328
libMesh::DofMap::variable_type
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:1924
libMesh::DenseVector::size
virtual unsigned int size() const override
Definition: dense_vector.h:92
libMesh::MeshFunction::~MeshFunction
~MeshFunction()
Destructor.
Definition: mesh_function.C:87
libMesh::MeshFunction::unset_point_locator_tolerance
void unset_point_locator_tolerance()
Turn off the user-specified PointLocator tolerance.
Definition: mesh_function.C:826
libMesh::MeshFunction::_dof_map
const DofMap & _dof_map
Need access to the DofMap of the other system.
Definition: mesh_function.h:339
libMesh::FunctionBase< Number >::_master
const FunctionBase * _master
Const pointer to our master, initialized to nullptr.
Definition: function_base.h:173
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::MeshFunction::_vector
const NumericVector< Number > & _vector
A reference to the vector that holds the data that is to be interpolated.
Definition: mesh_function.h:334
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::MeshBase::sub_point_locator
std::unique_ptr< PointLocatorBase > sub_point_locator() const
Definition: mesh_base.C:672
libMesh::Gradient
NumberVectorValue Gradient
Definition: exact_solution.h:58
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
value
static const bool value
Definition: xdr_io.C:56
libMesh::MeshFunction::_out_of_mesh_mode
bool _out_of_mesh_mode
true if out-of-mesh mode is enabled.
Definition: mesh_function.h:357
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::MeshFunction::disable_out_of_mesh_mode
void disable_out_of_mesh_mode(void)
Disables out-of-mesh mode.
Definition: mesh_function.C:814
libMesh::MeshFunction::_system_vars
const std::vector< unsigned int > _system_vars
The indices of the variables within the other system for which data are to be gathered.
Definition: mesh_function.h:345
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
data
IterBase * data
Ideally this private member data should have protected access.
Definition: variant_filter_iterator.h:337
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::ParallelObject
An object whose state is distributed along a set of processors.
Definition: parallel_object.h:55
libMesh::Tensor
NumberTensorValue Tensor
Definition: exact_solution.h:56
libMesh::PointLocatorBase::unset_close_to_point_tol
virtual void unset_close_to_point_tol()
Specify that we do not want to use a user-specified tolerance to determine if a point is contained wi...
Definition: point_locator_base.C:99
libMesh::PointLocatorBase
This is the base class for point locators.
Definition: point_locator_base.h:62
libMesh::MeshFunction::clear
virtual void clear() override
Clears the function.
Definition: mesh_function.C:121
libMesh::MeshFunction::MeshFunction
MeshFunction(const EquationSystems &eqn_systems, const NumericVector< Number > &vec, const DofMap &dof_map, const std::vector< unsigned int > &vars, const FunctionBase< Number > *master=nullptr)
Constructor for mesh based functions with vectors as return value.
Definition: mesh_function.C:42
libMesh::NumericVector::type
ParallelType type() const
Definition: numeric_vector.h:160
libMesh::MeshFunction::gradient
Gradient gradient(const Point &p, const Real time=0.)
Definition: mesh_function.C:174
libMesh::DenseVector< Number >
libMesh::FEInterface::compute_data
static void compute_data(const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data)
Lets the appropriate child of FEBase compute the requested data for the input specified in data,...
Definition: fe_interface.C:1028