20 #ifndef __mast__mesh_field_function__ 21 #define __mast__mesh_field_function__ 28 #include "libmesh/numeric_vector.h" 29 #include "libmesh/mesh_function.h" 30 #include "libmesh/system.h" 37 class SystemInitialization;
52 const std::string& nm);
71 virtual void operator() (
const libMesh::Point& p,
80 virtual void gradient(
const libMesh::Point& p,
99 const libMesh::Point& p,
109 void init(
const libMesh::NumericVector<Real>& sol,
110 const libMesh::NumericVector<Real>* dsol =
nullptr);
186 #endif // __mast__mesh_field_function__
MeshFieldFunction(MAST::SystemInitialization &sys, const std::string &nm)
constructor
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealVectorX &v) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
libMesh::NumericVector< Real > * _sol
current solution that is going to be interpolated
libMesh::NumericVector< Real > * _dsol
virtual void gradient(const libMesh::Point &p, const Real t, RealMatrixX &g) const
calculates the gradient of value of the function at the specified point, p, and time, t, and returns it in g.
libMesh::MeshFunction * _function
the MeshFunction object that performs the interpolation
virtual void perturbation(const libMesh::Point &p, const Real t, RealVectorX &v) const
calculates the value of perturbation in the function at the specified point, p, and time...
This provides a wrapper FieldFunction compatible class that interpolates the solution using libMesh's...
libMesh::MeshFunction * _perturbed_function
virtual void clear_element_quadrature_point_solution()
clears the quadrature point solution provided by the corresponding set method above.
virtual void operator()(const libMesh::Point &p, const Real t, RealVectorX &v) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
RealVectorX _qp_sol
quadrature point solution of the element
virtual void set_element_quadrature_point_solution(RealVectorX &sol)
When a mesh field function is attached to an assembly routine during system assembly, then the current solution can be provided by the element quadrature point update.
libMesh::MeshFunction & get_perturbed_function()
Matrix< Real, Dynamic, Dynamic > RealMatrixX
This creates the base class for functions that have a saptial and temporal dependence, and provide sensitivity operations with respect to the functions and parameters.
void clear()
clear the solution
Matrix< Real, Dynamic, 1 > RealVectorX
libMesh::MeshFunction & get_function()
bool _use_qp_sol
flag is set to true when the quadrature point solution is provided by an element
void init(const libMesh::NumericVector< Real > &sol, const libMesh::NumericVector< Real > *dsol=nullptr)
initializes the data structures to perform the interpolation function of sol.
libMesh::System * _sys
current system for which solution is to be interpolated
virtual ~MeshFieldFunction()
destructor