Go to the documentation of this file.
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"
45 const std::vector<unsigned int> & vars,
49 _eqn_systems (eqn_systems),
53 _point_locator (nullptr),
54 _out_of_mesh_mode (false),
64 const unsigned int var,
68 _eqn_systems (eqn_systems),
72 _point_locator (nullptr),
73 _out_of_mesh_mode (false),
144 return std::unique_ptr<FunctionBase<Number>>(mf_clone);
164 std::map<const Elem *, DenseVector<Number>> buffer;
166 std::map<const Elem *, Number> return_value;
167 for (
const auto & pr : buffer)
168 return_value[pr.first] = pr.second(0);
179 std::vector<Gradient> buf (1);
189 std::map<const Elem *, std::vector<Gradient>> buffer;
191 std::map<const Elem *, Gradient> return_value;
192 for (
const auto & pr : buffer)
193 return_value[pr.first] = pr.second[0];
199 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
205 std::vector<Tensor> buf (1);
221 const std::set<subdomain_id_type> * subdomain_ids)
238 output.
resize (cast_int<unsigned int>
243 const unsigned int dim = element->
dim();
282 std::vector<dof_id_type> dof_indices;
292 output(index) =
value;
316 const std::set<subdomain_id_type> * subdomain_ids)
324 std::set<const Elem *> candidate_element = this->
find_elements(p,subdomain_ids);
328 for (
const auto & element : candidate_element)
330 const unsigned int dim = element->dim();
370 std::vector<dof_id_type> dof_indices;
380 temp_output(index) =
value;
389 output[element] = temp_output;
397 std::vector<Gradient> & output,
398 const std::set<subdomain_id_type> * subdomain_ids)
417 const unsigned int dim = element->
dim();
428 std::vector<Point> point_list (1, mapped_point);
449 std::vector<dof_id_type> dof_indices;
454 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
461 const std::vector<std::vector<RealGradient>> & dphi = point_fe->get_dphi();
462 point_fe->reinit(element, &point_list);
467 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
476 data.enable_derivative();
483 for (std::size_t v=0; v<
dim; v++)
484 for (std::size_t xyz=0; xyz<LIBMESH_DIM; xyz++)
487 grad(xyz) +=
data.local_transform[v][xyz]
489 * this->
_vector(dof_indices[i]);
494 output[index] = grad;
503 std::map<
const Elem *, std::vector<Gradient>> & output)
512 std::map<
const Elem *, std::vector<Gradient>> & output,
513 const std::set<subdomain_id_type> * subdomain_ids)
521 std::set<const Elem *> candidate_element = this->
find_elements(p,subdomain_ids);
525 for (
const auto & element : candidate_element)
527 const unsigned int dim = element->dim();
530 std::vector<Gradient> temp_output (cast_int<unsigned int>(this->
_system_vars.size()));
541 std::vector<Point> point_list (1, mapped_point);
560 std::vector<dof_id_type> dof_indices;
568 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
569 if (!element->infinite())
573 const std::vector<std::vector<RealGradient>> & dphi = point_fe->get_dphi();
574 point_fe->reinit(element, & point_list);
578 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
588 data.enable_derivative();
596 for (std::size_t v=0; v<
dim; v++)
597 for (std::size_t xyz=0; xyz<LIBMESH_DIM; xyz++)
600 grad(xyz) +=
data.local_transform[v][xyz]
602 * this->
_vector(dof_indices[i]);
607 temp_output[index] = grad;
613 output[element] = temp_output;
619 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
622 std::vector<Tensor> & output,
623 const std::set<subdomain_id_type> * subdomain_ids)
635 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
637 libmesh_warning(
"Warning: Requested the Hessian of an Infinite element."
638 <<
"Second derivatives for Infinite elements"
639 <<
" are not yet implemented!"
648 const unsigned int dim = element->
dim();
659 std::vector<Point> point_list (1, mapped_point);
679 const std::vector<std::vector<RealTensor>> & d2phi =
680 point_fe->get_d2phi();
681 point_fe->reinit(element, &point_list);
684 std::vector<dof_id_type> dof_indices;
693 output[index] = hess;
701 const std::set<subdomain_id_type> * subdomain_ids)
const
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.");
730 std::set<const Elem *> point_neighbors;
733 for (
const auto & elem : point_neighbors)
745 const std::set<subdomain_id_type> * subdomain_ids)
const
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.");
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)
777 std::set<const Elem *> point_neighbors;
778 element->find_point_neighbors(p, point_neighbors);
779 for (
const auto & elem : point_neighbors)
782 final_candidate_elements.insert(elem);
787 final_candidate_elements.insert(element);
790 return final_candidate_elements;
const Elem * find_element(const Point &p, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
Helper function to reduce code duplication.
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...
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value.
DenseVector< Number > _out_of_mesh_value
Value to return outside the mesh if out-of-mesh mode is enabled.
class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEIn...
const MeshBase & get_mesh() const
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.
virtual void init() override
Override the FunctionBase::init() member function by calling our own and specifying the Trees::NODES ...
Number operator()(const Point &p, const Real time=0.) override
This class provides function-like objects for data distributed over a mesh.
std::set< const Elem * > find_elements(const Point &p, const std::set< subdomain_id_type > *subdomain_ids=nullptr) const
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.
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...
The libMesh namespace provides an interface to certain functionality in the library.
virtual unsigned short dim() const =0
std::map< const Elem *, Gradient > discontinuous_gradient(const Point &p, const Real time=0.)
virtual bool infinite() const =0
virtual std::unique_ptr< FunctionBase< Number > > clone() const override
std::map< const Elem *, Number > discontinuous_value(const Point &p, const Real time=0.)
virtual void disable_out_of_mesh_mode()=0
Disables out-of-mesh mode (default).
BuildType
enum defining how to build the tree.
Tensor hessian(const Point &p, const Real time=0.)
processor_id_type processor_id() const
const PointLocatorBase & get_point_locator(void) const
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...
virtual void enable_out_of_mesh_mode()=0
Enables out-of-mesh mode.
This is the MeshBase class.
void enable_out_of_mesh_mode(const DenseVector< Number > &value)
Enables out-of-mesh mode.
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
processor_id_type processor_id() const
bool _initialized
When init() was called so that everything is ready for calls to operator() (...), then this bool is t...
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
A Point defines a location in LIBMESH_DIM dimensional Real space.
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
PointLocatorBase * _point_locator
A point locator is needed to locate the points in the mesh.
void add_scaled(const TypeTensor< T2 > &, const T &)
Add a scaled tensor to this tensor without creating a temporary.
const EquationSystems & _eqn_systems
The equation systems handler, from which the data are gathered.
const FEType & variable_type(const unsigned int c) const
virtual unsigned int size() const override
~MeshFunction()
Destructor.
void unset_point_locator_tolerance()
Turn off the user-specified PointLocator tolerance.
const DofMap & _dof_map
Need access to the DofMap of the other system.
const FunctionBase * _master
Const pointer to our master, initialized to nullptr.
This is the EquationSystems class.
const NumericVector< Number > & _vector
A reference to the vector that holds the data that is to be interpolated.
void resize(const unsigned int n)
Resize the vector.
std::unique_ptr< PointLocatorBase > sub_point_locator() const
NumberVectorValue Gradient
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
bool _out_of_mesh_mode
true if out-of-mesh mode is enabled.
This class handles the numbering of degrees of freedom on a mesh.
void disable_out_of_mesh_mode(void)
Disables out-of-mesh mode.
const std::vector< unsigned int > _system_vars
The indices of the variables within the other system for which data are to be gathered.
This is the base class from which all geometric element types are derived.
IterBase * data
Ideally this private member data should have protected access.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
An object whose state is distributed along a set of processors.
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...
This is the base class for point locators.
virtual void clear() override
Clears the function.
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.
ParallelType type() const
Gradient gradient(const Point &p, const Real time=0.)
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,...