NSF Postdoctoral Research
A set of C++ code developed by Andrew E. Slaughter
|
A base class including nodal data. More...
#include <fem/common/eq_data_base.h>
Public Member Functions | |
void | add_variable (const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=NULL) |
A short-cut for adding a variable (see libmesh) | |
void | add_variable (const std::string &var, const Order order=FIRST, const FEFamily family=MONOMIAL, const std::set< subdomain_id_type > *const active_subdomains=NULL) |
A short-cut for adding a variable (see libmesh) | |
void | update_solution (Real t=0) |
Projects the nodal data at the prescribed time. | |
const unsigned int | dimension () |
Returns the number of dimensions for the system. | |
Protected Member Functions | |
EqDataBase (EquationSystems &sys, string data_name) | |
Class constructor. | |
virtual void | value (DenseVector< Number > &output, const Point &p, const Real t=0)=0 |
A pure virtual function that will act as a function for libmesh. |
A base class including nodal data.
Provides mechanims for defining nodal data using libMesh a libmesh system. The data may very spatially and temporally. In the inherited class the pure virtual function value
dictates the behavior.
EqDataBase::EqDataBase | ( | EquationSystems & | sys, |
string | data_name | ||
) | [protected] |
Class constructor.
This class is meant to be inherited, as such the constructor is protected.
The class requires that an existing EquationSystem object be passed in, the system with the name given by data_name
is added to the EquationSystem object.
sys | libmesh::EquationSystems object that the new system for storing data will be added |
data_name | The name given the the system for storing data |
void EqDataBase::add_variable | ( | const std::string & | var, |
const FEType & | type, | ||
const std::set< subdomain_id_type > *const | active_subdomains = NULL |
||
) |
A short-cut for adding a variable (see libmesh)
Replicates the behavior of libMesh add_variable from the TransientExplicitSystem, see libmesh documentation for further details.
var | The name of the variable |
type | The type of FE variable |
active_subdomains | Limit the variable to a sub-domain |
void EqDataBase::add_variable | ( | const std::string & | var, |
const Order | order = FIRST , |
||
const FEFamily | family = MONOMIAL , |
||
const std::set< subdomain_id_type > *const | active_subdomains = NULL |
||
) |
A short-cut for adding a variable (see libmesh)
Replicates the behavior of libMesh add_variable from the TransientExplicitSystem, see libmesh documentation for further details.
var | The name of the variable |
order | The order of the data (assumed FIRST) |
family | The FEfamily, assumed to be MONOMIAL which is recommended for storing data |
active_subdomains | Limit the variable to a sub-domain |
const unsigned int SlaughterFEM::EqDataBase::dimension | ( | ) |
Returns the number of dimensions for the system.
void EqDataBase::update_solution | ( | Real | t = 0 | ) |
Projects the nodal data at the prescribed time.
This function creates a boost::function pointer to the value member of this class, packages that pointer into a libmesh acceptable format with MyAnalyticFunction, and then uses it to project the nodal data.
t | A libmesh::Real containing the current time, this time is used to redefine the TransientExplicitSystem time that contains the nodal data, thus it may be used for temporal data or it may just be ignored by letting it default to zero. |
virtual void SlaughterFEM::EqDataBase::value | ( | DenseVector< Number > & | output, |
const Point & | p, | ||
const Real | t = 0 |
||
) | [protected, pure virtual] |
A pure virtual function that will act as a function for libmesh.
The inheriting class must redefine this function, see example3::cpp .
output | A reference to the vector that stores the stores the computed flux |
p | A libMesh point |
t | The time, this time is not actually used libmesh automatically uses zero. |