20 #ifndef __mast__physics_discipline_base__    21 #define __mast__physics_discipline_base__    30 #include "libmesh/equation_systems.h"    36     class BoundaryConditionBase;
    37     class DirichletBoundaryCondition;
    38     class PropertyCardBase;
    40     class FunctionSetBase;
    41     class ElementPropertyCardBase;
    42     class SystemInitialization;
    44     class PointLoadCondition;
    45     class NonlinearSystem;
    50     typedef std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>  
SideBCMapType;
    51     typedef std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*> 
VolumeBCMapType;
    52     typedef std::map<libMesh::subdomain_id_type, const MAST::ElementPropertyCardBase*>      
PropertyCardMapType;
    53     typedef std::map<libMesh::boundary_id_type, MAST::DirichletBoundaryCondition*>  
DirichletBCMapType;
   181                                               const unsigned int var);
   201                                           std::set<unsigned int>& dof_ids) 
const;
   268 #endif // __mast__physics_discipline_base__ void add_volume_load(libMesh::subdomain_id_type bid, MAST::BoundaryConditionBase &load)
adds the specified volume loads for the elements with subdomain tag s_id 
MAST::SideBCMapType & side_loads()
void clear_volume_load(libMesh::subdomain_id_type sid, MAST::BoundaryConditionBase &load)
clear the specified volume load from the applied loads 
std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > VolumeBCMapType
const MAST::ElementPropertyCardBase & get_property_card(const libMesh::Elem &elem) const
get property card for the specified element 
std::map< libMesh::boundary_id_type, MAST::DirichletBoundaryCondition * > DirichletBCMapType
This class implements a system for solution of nonlinear systems. 
MAST::VolumeBCMapType _vol_bc_map
volume boundary condition map of boundary id and load 
libMesh::EquationSystems & get_equation_systems()
returns a reference to the libMesh::System object 
std::map< libMesh::subdomain_id_type, const MAST::ElementPropertyCardBase * > PropertyCardMapType
virtual ~PhysicsDisciplineBase()
virtual destructor 
void init_system_dirichlet_bc(MAST::NonlinearSystem &sys) const
initializes the system for dirichlet boundary conditions 
void remove_volume_load(libMesh::subdomain_id_type bid, MAST::BoundaryConditionBase &load)
remove the specified volume loads for the elements with subdomain tag s_id 
const MAST::PointLoadSetType & point_loads() const
void add_dirichlet_bc(libMesh::boundary_id_type bid, MAST::DirichletBoundaryCondition &load)
adds the specified Dirichlet boundary condition for the boundary with tag b_id 
void add_side_load(libMesh::boundary_id_type bid, MAST::BoundaryConditionBase &load)
adds the specified side loads for the boudnary with tag b_id 
void clear_loads()
clear the loads and pointer to static solution system for this structural model 
std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > SideBCMapType
PhysicsDisciplineBase(libMesh::EquationSystems &eq_sys)
void remove_side_load(libMesh::boundary_id_type bid, MAST::BoundaryConditionBase &load)
remove the specified side loads for the boudnary with tag b_id 
std::map< libMesh::subdomain_id_type, std::vector< unsigned int > > _subdomain_var_constraint
variables constrained on subdomain 
This class allows for the specification of load associated with specified nodes in a user-provided se...
const MAST::SideBCMapType & side_loads() const
void get_system_dirichlet_bc_dofs(libMesh::System &sys, std::set< unsigned int > &dof_ids) const
Prepares a list of the constrained dofs for system sys and returns in dof_ids. 
MAST::PointLoadSetType & point_loads()
void clear_system_dirichlet_bc(MAST::NonlinearSystem &sys) const
clears the system dirichlet boundary conditions 
void constrain_subdomain_dofs_for_var(const libMesh::subdomain_id_type sid, const unsigned int var)
constrain dofs on a subdomain to zero 
libMesh::EquationSystems & _eq_systems
libMesh::System for which analysis is to be performed 
MAST::VolumeBCMapType & volume_loads()
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
void add_point_load(MAST::PointLoadCondition &load)
adds the specified point load 
const MAST::VolumeBCMapType & volume_loads() const
MAST::DirichletBCMapType _dirichlet_bc_map
Dirichlet boundary condition map of boundary id and load. 
void set_property_for_subdomain(const libMesh::subdomain_id_type sid, const MAST::ElementPropertyCardBase &prop)
sets the same property for all elements in the specified subdomain 
MAST::PointLoadSetType _point_loads
point loads 
std::set< MAST::PointLoadCondition * > PointLoadSetType
MAST::PropertyCardMapType _element_property
map of element property cards for each element 
MAST::SideBCMapType _side_bc_map
side boundary condition map of boundary id and load