This provides a wrapper FieldFunction compatible class that interpolates the solution using libMesh's MeshFunction class.
Definition at line 44 of file mesh_field_function.h.
#include <mesh_field_function.h>
Classes | |
struct | SolFunc |
Public Member Functions | |
MeshFieldFunction (MAST::SystemInitialization &sys, const std::string &nm, libMesh::ParallelType p_type) | |
constructor More... | |
MeshFieldFunction (libMesh::System &sys, const std::string &nm, libMesh::ParallelType p_type) | |
constructor More... | |
virtual | ~MeshFieldFunction () |
destructor More... | |
void | clear () |
clear the solution More... | |
virtual void | clear_element_quadrature_point_solution () |
clears the quadrature point solution provided by the corresponding set method above. More... | |
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 . More... | |
virtual void | derivative_gradient (const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &v) const |
calculates the value of the function at the specified point, p , and time, t , and returns it in v . More... | |
libMesh::MeshFunction & | get_function () |
libMesh::MeshFunction & | get_perturbed_function () |
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 . More... | |
void | init (const libMesh::NumericVector< Real > &sol, bool reuse_vector) |
initializes the data structures to perform the interpolation function of sol . More... | |
void | init_sens (const MAST::FunctionBase &f, const libMesh::NumericVector< Real > &dsol, bool reuse_vector) |
initializes the the data structures for computation of sensitivity for the specified function. More... | |
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 . More... | |
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, t , and returns it in v . More... | |
virtual void | perturbation_gradient (const libMesh::Point &p, const Real t, RealMatrixX &v) const |
calculates the value of the function at the specified point, p , and time, t , and returns it in v . More... | |
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. More... | |
Public Member Functions inherited from MAST::FieldFunction< RealVectorX > | |
FieldFunction (const std::string &nm) | |
virtual void | derivative (const MAST::FunctionBase &f, RealVectorX &v) const |
calculates the value of the function derivative and returns it in v . More... | |
virtual void | operator() (RealVectorX &v) const |
calculates the value of the function and returns it in v . More... | |
virtual void | perturbation (RealVectorX &v) const |
calculates the perturbation and returns it in v . More... | |
Public Member Functions inherited from MAST::FunctionBase | |
FunctionBase (const std::string &nm, const bool is_field_func) | |
initializes the parameter to the given name More... | |
FunctionBase (const MAST::FunctionBase &f) | |
Copy constructor. More... | |
virtual | ~FunctionBase () |
virtual destructor More... | |
virtual bool | depends_on (const MAST::FunctionBase &f) const |
returns true if the function depends on the provided value More... | |
virtual bool | is_shape_parameter () const |
virtual bool | is_topology_parameter () const |
const std::string & | name () const |
returns the name of this function More... | |
virtual void | set_as_shape_parameter (bool f) |
virtual void | set_as_topology_parameter (bool f) |
Protected Member Functions | |
void | _init_sol_func (bool reuse_sol, const libMesh::NumericVector< Real > &sol, MAST::MeshFieldFunction::SolFunc &sol_func) |
Protected Attributes | |
SolFunc * | _function |
current solution that is going to be interpolated More... | |
std::map< const MAST::FunctionBase *, MAST::MeshFieldFunction::SolFunc * > | _function_sens |
solution sensitivity for specified value More... | |
libMesh::ParallelType | _p_type |
type of parallel vector required for this mesh function. More... | |
SolFunc * | _perturbed_function |
current perturbation solution that is going to be interpolated More... | |
RealVectorX | _qp_sol |
quadrature point solution of the element More... | |
libMesh::System * | _sys |
current system for which solution is to be interpolated More... | |
bool | _use_qp_sol |
flag is set to true when the quadrature point solution is provided by an element More... | |
Protected Attributes inherited from MAST::FunctionBase | |
std::set< const MAST::FunctionBase * > | _functions |
set of functions that this function depends on More... | |
bool | _is_field_func |
flag to store the nature of field function More... | |
bool | _is_shape_parameter |
bool | _is_topology_parameter |
std::string | _name |
name of this parameter More... | |
MAST::MeshFieldFunction::MeshFieldFunction | ( | MAST::SystemInitialization & | sys, |
const std::string & | nm, | ||
libMesh::ParallelType | p_type | ||
) |
constructor
Definition at line 30 of file mesh_field_function.cpp.
MAST::MeshFieldFunction::MeshFieldFunction | ( | libMesh::System & | sys, |
const std::string & | nm, | ||
libMesh::ParallelType | p_type | ||
) |
constructor
Definition at line 45 of file mesh_field_function.cpp.
|
virtual |
destructor
Definition at line 60 of file mesh_field_function.cpp.
|
protected |
Definition at line 366 of file mesh_field_function.cpp.
void MAST::MeshFieldFunction::clear | ( | ) |
clear the solution
Definition at line 315 of file mesh_field_function.cpp.
|
virtual |
clears the quadrature point solution provided by the corresponding set method above.
Definition at line 358 of file mesh_field_function.cpp.
|
virtual |
calculates the value of the function at the specified point, p
, and time, t
, and returns it in v
.
Reimplemented from MAST::FieldFunction< RealVectorX >.
Definition at line 205 of file mesh_field_function.cpp.
|
virtual |
calculates the value of the function at the specified point, p
, and time, t
, and returns it in v
.
Definition at line 242 of file mesh_field_function.cpp.
|
inline |
Definition at line 146 of file mesh_field_function.h.
|
inline |
Definition at line 156 of file mesh_field_function.h.
|
virtual |
calculates the gradient of value of the function at the specified point, p
, and time, t
, and returns it in g
.
g(i,j) = dv(i)/dx(j)
Definition at line 105 of file mesh_field_function.cpp.
void MAST::MeshFieldFunction::init | ( | const libMesh::NumericVector< Real > & | sol, |
bool | reuse_vector | ||
) |
initializes the data structures to perform the interpolation function of sol
.
If dsol
is provided, then it is used as the perturbation of sol
. if reuse_vector
is true
then instead of cloning sol
this vector will be used and it should be of type p_type
set in the constructor of this object. Note that this requires that sol
not be destroyed till this object needs it.
Definition at line 281 of file mesh_field_function.cpp.
void MAST::MeshFieldFunction::init_sens | ( | const MAST::FunctionBase & | f, |
const libMesh::NumericVector< Real > & | dsol, | ||
bool | reuse_vector | ||
) |
initializes the the data structures for computation of sensitivity for the specified function.
Definition at line 294 of file mesh_field_function.cpp.
|
virtual |
calculates the value of the function at the specified point, p
, and time, t
, and returns it in v
.
Reimplemented from MAST::FieldFunction< RealVectorX >.
Definition at line 72 of file mesh_field_function.cpp.
|
virtual |
calculates the value of perturbation in the function at the specified point, p
, and time, t
, and returns it in v
.
Reimplemented from MAST::FieldFunction< RealVectorX >.
Definition at line 139 of file mesh_field_function.cpp.
|
virtual |
calculates the value of the function at the specified point, p
, and time, t
, and returns it in v
.
Definition at line 171 of file mesh_field_function.cpp.
|
virtual |
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.
This method allows the element to provide the solution and in turn override the mesh function evaluation.
Definition at line 348 of file mesh_field_function.cpp.
|
protected |
current solution that is going to be interpolated
Definition at line 232 of file mesh_field_function.h.
|
protected |
solution sensitivity for specified value
Definition at line 243 of file mesh_field_function.h.
|
protected |
type of parallel vector required for this mesh function.
Definition at line 217 of file mesh_field_function.h.
|
protected |
current perturbation solution that is going to be interpolated
Definition at line 238 of file mesh_field_function.h.
|
protected |
quadrature point solution of the element
Definition at line 222 of file mesh_field_function.h.
|
protected |
current system for which solution is to be interpolated
Definition at line 227 of file mesh_field_function.h.
|
protected |
flag is set to true when the quadrature point solution is provided by an element
Definition at line 212 of file mesh_field_function.h.