OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
oofem::DofManager Class Reference

Base class for dof managers. More...

#include <dofmanager.h>

+ Inheritance diagram for oofem::DofManager:
+ Collaboration diagram for oofem::DofManager:

Public Member Functions

std::vector< Dof * >::iterator begin ()
 
std::vector< Dof * >::iterator end ()
 
std::vector< Dof * >::const_iterator begin () const
 
std::vector< Dof * >::const_iterator end () const
 
 DofManager (int n, Domain *aDomain)
 Constructor. More...
 
virtual ~DofManager ()
 Destructor. More...
 
virtual void printOutputAt (FILE *file, TimeStep *tStep)
 Prints output of receiver to stream, for given time step. More...
 
virtual void updateYourself (TimeStep *tStep)
 Updates receiver after equilibrium in time step has been reached. More...
 
bool isBoundary ()
 
void setBoundaryFlag (bool isBoundary)
 Sets the boundary flag. More...
 
virtual bool hasAnySlaveDofs ()
 
virtual bool giveMasterDofMans (IntArray &masters)
 Returns true if the receiver is linked (its slave DOFs depend on master values) to some other dof managers. More...
 
virtual IRResultType initializeFrom (InputRecord *ir)
 Initializes receiver according to object description stored in input record. More...
 
virtual void giveInputRecord (DynamicInputRecord &input)
 Setups the input record string of receiver. More...
 
virtual void printYourself ()
 Prints receiver state on stdout. Useful for debugging. More...
 
virtual contextIOResultType saveContext (DataStream &stream, ContextMode mode, void *obj=NULL)
 Stores receiver state to output stream. More...
 
virtual contextIOResultType restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL)
 Restores the receiver state previously written in stream. More...
 
virtual bool isDofTypeCompatible (dofType type) const
 Returns true if dof of given type is allowed to be associated to receiver. More...
 
virtual void postInitialize ()
 Performs post-initialization such like checking if there are any slave dofs etc. More...
 
virtual void updateLocalNumbering (EntityRenumberingFunctor &f)
 Local renumbering support. More...
 
void appendDof (Dof *dof)
 Adds the given Dof into the receiver. More...
 
void removeDof (DofIDItem id)
 Removes Dof with given id from dofArray. More...
 
bool hasDofID (DofIDItem id) const
 Checks if receiver contains dof with given ID. More...
 
virtual void drawYourself (oofegGraphicContext &gc, TimeStep *tStep)
 
int giveGlobalNumber () const
 
int giveLabel () const
 
void setGlobalNumber (int newNumber)
 Sets receiver global number. More...
 
dofManagerParallelMode giveParallelMode () const
 Return dofManagerParallelMode of receiver. More...
 
void setParallelMode (dofManagerParallelMode _mode)
 Sets parallel mode of receiver. More...
 
const IntArraygivePartitionList ()
 Returns partition list of receiver. More...
 
void setPartitionList (const IntArray *_p)
 Sets receiver's partition list. More...
 
void removePartitionFromList (int _part)
 Removes given partition from receiver list. More...
 
void mergePartitionList (IntArray &_p)
 Merges receiver partition list with given lists. More...
 
int givePartitionsConnectivitySize ()
 Returns number of partitions sharing given receiver (=number of shared partitions + local one). More...
 
bool isLocal ()
 Returns true if receiver is locally maintained. More...
 
bool isShared ()
 Returns true if receiver is shared. More...
 
bool isNull ()
 Returns true if receiver is shared. More...
 
Dof management methods
DofgiveDofWithID (int dofID) const
 Returns DOF with given dofID; issues error if not present. More...
 
int giveNumberOfDofs () const
 
void askNewEquationNumbers (TimeStep *tStep)
 Renumbers all contained DOFs. More...
 
int giveNumberOfPrimaryMasterDofs (const IntArray &dofIDArray) const
 Returns the number of primary dofs on which receiver dofs (given in dofArray) depend on. More...
 
void giveLocationArray (const IntArray &dofIDArry, IntArray &locationArray, const UnknownNumberingScheme &s) const
 Returns location array (array containing for each requested dof related equation number) for given numbering scheme. More...
 
void giveMasterDofIDArray (const IntArray &dofIDArry, IntArray &masterDofIDs) const
 Returns master dof ID array of receiver. More...
 
void giveCompleteLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s) const
 Returns full location array of receiver containing equation numbers of all dofs of receiver. More...
 
void giveCompleteMasterDofIDArray (IntArray &dofIDArray) const
 Returns the full dof ID array of receiver. More...
 
std::vector< Dof * >::const_iterator findDofWithDofId (DofIDItem dofID) const
 Finds index of DOF with required physical meaning of receiver. More...
 
void giveUnknownVector (FloatArray &answer, const IntArray &dofMask, ValueModeType mode, TimeStep *tStep, bool padding=false)
 Assembles the vector of unknowns in global c.s for given dofs of receiver. More...
 
void giveUnknownVector (FloatArray &answer, const IntArray &dofMask, PrimaryField &field, ValueModeType mode, TimeStep *tStep, bool padding=false)
 Assembles the vector of unknowns of given filed in global c.s for given dofs of receiver. More...
 
void giveCompleteUnknownVector (FloatArray &answer, ValueModeType mode, TimeStep *tStep)
 Assembles the complete unknown vector in node. More...
 
void giveUnknownVectorOfType (FloatArray &answer, UnknownType ut, ValueModeType mode, TimeStep *tStep)
 Constructs the requested vector by assembling e.g. More...
 
void givePrescribedUnknownVector (FloatArray &answer, const IntArray &dofMask, ValueModeType mode, TimeStep *tStep)
 Assembles the vector of prescribed unknowns in nodal c.s for given dofs of receiver. More...
 
Transformation functions

The governing equations can be assembled not only in global coordinate system, but also in user-defined local coordinate system of each DOF manager.

Methods in this section introduce necessary transformation methods, allowing receiver DOFs to be expressed in their own local c.s. or to be dependent on other DOFs on other DOF manager (to implement slave or rigid arm nodes etc.). The method for computing global c.s to receiver c.s transformation matrix is provided.

bool computeM2GTransformation (FloatMatrix &answer, const IntArray &dofIDArry)
 Computes receiver transformation matrix from global CS to dofManager specific coordinate system; $ u_g = R\cdot u_m $. More...
 
virtual bool computeL2GTransformation (FloatMatrix &answer, const IntArray &dofIDArry)
 Computes transformation matrix from global c.s. More...
 
virtual bool computeM2LTransformation (FloatMatrix &answer, const IntArray &dofIDArry)
 Computes transformation matrix from local DOFs to master DOFs; $ u_l = M\cdot u_m $. More...
 
virtual bool requiresTransformation ()
 Indicates, whether dofManager requires the transformation. More...
 
Load related functions
virtual void computeLoadVector (FloatArray &answer, Load *load, CharType type, TimeStep *tStep, ValueModeType mode)
 Computes the load vector for given load. More...
 
IntArraygiveLoadArray ()
 Returns the array containing applied loadings of the receiver. More...
 
void setLoadArray (IntArray &load)
 Sets the array of applied loadings of the receiver. More...
 
Position query functions
virtual bool hasCoordinates ()
 
virtual double giveCoordinate (int i)
 
virtual FloatArraygiveCoordinates ()
 
Functions necessary for dof creation. All optional.
const IntArraygiveForcedDofIDs ()
 Returns list of specific dofs that should be included in node. More...
 
std::map< int, int > * giveDofTypeMap ()
 Returns map from DofIDItem to dofType. More...
 
std::map< int, int > * giveMasterMap ()
 Returns map from DofIDItem to dofType. More...
 
std::map< int, int > * giveBcMap ()
 Returns map from DofIDItem to dofType. More...
 
std::map< int, int > * giveIcMap ()
 Returns map from DofIDItem to initial condition. More...
 
Advanced functions
void setNumberOfDofs (int _ndofs)
 Sets number of dofs of the receiver; Deallocates existing DOFs; Resizes the dofArray accordingly. More...
 
- Public Member Functions inherited from oofem::FEMComponent
 FEMComponent (int n, Domain *d)
 Regular constructor, creates component with given number and belonging to given domain. More...
 
virtual ~FEMComponent ()
 Virtual destructor. More...
 
virtual const char * giveClassName () const =0
 
virtual const char * giveInputRecordName () const =0
 
DomaingiveDomain () const
 
virtual void setDomain (Domain *d)
 Sets associated Domain. More...
 
int giveNumber () const
 
void setNumber (int num)
 Sets number of receiver. More...
 
virtual int checkConsistency ()
 Allows programmer to test some internal data, before computation begins. More...
 
virtual InterfacegiveInterface (InterfaceType t)
 Interface requesting service. More...
 
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros). More...
 

Protected Attributes

std::vector< Dof * > dofArray
 Array of DOFs. More...
 
IntArray loadArray
 List of applied loads. More...
 
bool isBoundaryFlag
 Indicates if dofManager is boundary (true boundary or on boundary between regions) or interior. More...
 
bool hasSlaveDofs
 Flag indicating whether receiver has slave DOFs. More...
 
int globalNumber
 In parallel mode, globalNumber contains globally unique DoFManager number. More...
 
dofManagerParallelMode parallel_mode
 
IntArray partitions
 List of partition sharing the shared dof manager or remote partition containing remote dofmanager counterpart. More...
 
IntArraydofidmask
 List of additional dof ids to include. More...
 
std::map< int, int > * dofTypemap
 Map from DofIDItem to dofType. More...
 
std::map< int, int > * dofMastermap
 Map from DofIDItem to master node. More...
 
std::map< int, int > * dofBCmap
 Map from DofIDItem to bc (to be removed). More...
 
std::map< int, int > * dofICmap
 Map from DofIDItem to ic (to be removed). More...
 
IntArray mBC
 
- Protected Attributes inherited from oofem::FEMComponent
int number
 Component number. More...
 
Domaindomain
 Link to domain object, useful for communicating with other FEM components. More...
 

Detailed Description

Base class for dof managers.

Dof manager is an abstraction for object possessing degrees of freedom. Dof managers (respectively derived classes like nodes or sides) are usually attributes of elements and are maintained by domain. Degrees of freedom belonging to dof manager are stored in 'dofArray'. Dof manager also maintain a list of applied loads it is subjected to. Number and physical meaning of dofs can be specified by user in input file (see input file description). If it is not specified, default values are obtained from domain, based on domain type of problem.

Tasks:

  • managing its degrees of freedom.
  • calculating its nodal vector.
  • printing and updating at end of step .
  • managing its swapping to and from disk.
  • managing the list of subjected loads (typically, concentrated forces and moments).
  • constructing the location array for assembling loads (also used by elements for assembling matrices).

Definition at line 113 of file dofmanager.h.

Constructor & Destructor Documentation

oofem::DofManager::DofManager ( int  n,
Domain aDomain 
)

Constructor.

Creates DofManager with given number belonging to domain aDomain.

Parameters
nDofManager's number in domain
aDomainreference to DofManager's domain

Definition at line 54 of file dofmanager.C.

References dofBCmap, dofICmap, dofidmask, oofem::DofManager_local, dofMastermap, dofTypemap, hasSlaveDofs, isBoundaryFlag, and parallel_mode.

oofem::DofManager::~DofManager ( )
virtual

Destructor.

Definition at line 68 of file dofmanager.C.

References dofArray, dofBCmap, dofICmap, dofidmask, dofMastermap, and dofTypemap.

Member Function Documentation

void oofem::DofManager::appendDof ( Dof dof)

Adds the given Dof into the receiver.

The dofID of scheduled DOF should not be present in receiver as multiple DOFs with same DofID are not allowed. The given DOF is appended at the end of the dofArray.

Parameters
dof

Definition at line 134 of file dofmanager.C.

References dofArray, end(), findDofWithDofId(), oofem::Dof::giveDofID(), oofem::FEMComponent::number, and OOFEM_ERROR.

Referenced by oofem::ContactDefinition::createContactDofs(), oofem::Domain::createDofs(), oofem::Subdivision::createMesh(), oofem::PrescribedGradientBCWeak::createTractionMesh(), oofem::PrescribedMean::initializeFrom(), oofem::SolutionbasedShapeFunction::initializeFrom(), oofem::Beam2d::initializeFrom(), oofem::Beam3d::initializeFrom(), and restoreContext().

void oofem::DofManager::askNewEquationNumbers ( TimeStep tStep)

Renumbers all contained DOFs.

Definition at line 295 of file dofmanager.C.

References oofem::Dof::askNewEquationNumber().

std::vector< Dof* >:: iterator oofem::DofManager::begin ( )
inline
std::vector< Dof* >:: const_iterator oofem::DofManager::begin ( ) const
inline

Definition at line 159 of file dofmanager.h.

bool oofem::DofManager::computeL2GTransformation ( FloatMatrix answer,
const IntArray dofIDArry 
)
virtual

Computes transformation matrix from global c.s.

to DOF-manager specific c.s; $ u_g = Q\cdot u_l $.

Parameters
answerComputed transformation matrix.
dofIDArryArray containing DofIDItem-type values for which transformation matrix is assembled. If dofIDArry is NULL, then all receiver DOFs are assumed.
Returns
True is transformation is needed, false otherwise.

Reimplemented in oofem::Node.

Definition at line 899 of file dofmanager.C.

Referenced by computeM2GTransformation(), givePrescribedUnknownVector(), giveUnknownVector(), and giveUnknownVectorOfType().

void oofem::DofManager::computeLoadVector ( FloatArray answer,
Load load,
CharType  type,
TimeStep tStep,
ValueModeType  mode 
)
virtual

Computes the load vector for given load.

Parameters
answerLoad vector for given load.
loadGiven load.
typeCharacteristic type of the vector.
tStepTime step when answer is computed.
modeDetermines response mode.

Reimplemented in oofem::Node.

Definition at line 97 of file dofmanager.C.

References oofem::FloatArray::assemble(), oofem::FloatArray::clear(), oofem::Load::computeComponentArrayAt(), oofem::GeneralBoundaryCondition::giveBCGeoType(), oofem::GeneralBoundaryCondition::giveDofIDs(), giveNumberOfDofs(), oofem::IntArray::giveSize(), oofem::NodalLoadBGT, OOFEM_ERROR, and oofem::FloatArray::resize().

Referenced by oofem::ExternalForceAssembler::vectorFromNodeLoad(), and oofem::ReferenceForceAssembler::vectorFromNodeLoad().

bool oofem::DofManager::computeM2GTransformation ( FloatMatrix answer,
const IntArray dofIDArry 
)

Computes receiver transformation matrix from global CS to dofManager specific coordinate system; $ u_g = R\cdot u_m $.

Parameters
answerComputed transformation matrix. It has generally dofIDArry.size rows and if loc is obtained using giveLocationArray(dofIDArry, loc) call, loc.giveSize() columns. This is because this transformation should generally include not only transformation to dof manager local coordinate system, but receiver dofs can be expressed using dofs of another dofManager (In this case, square answer is produced only if all dof transformation is required).
dofIDArryArray containing DofIDItem-type values (this is enumeration identifying physical meaning of particular DOF, see cltypes.h) for which transformation matrix is assembled. if dofIDArry is NULL, then all receiver dofs are assumed.
Returns
True if transformation is needed, false otherwise.

Definition at line 874 of file dofmanager.C.

References oofem::FloatMatrix::beProductOf(), oofem::FloatMatrix::clear(), computeL2GTransformation(), and computeM2LTransformation().

bool oofem::DofManager::computeM2LTransformation ( FloatMatrix answer,
const IntArray dofIDArry 
)
virtual

Computes transformation matrix from local DOFs to master DOFs; $ u_l = M\cdot u_m $.

Parameters
answerComputed transformation matrix.
dofIDArryArray containing DofIDItem-type values for which transformation matrix is assembled. If dofIDArry is NULL, then all receiver DOFs are assumed.
Returns
True is transformation is needed, false otherwise.
Todo:
I don't think this should be called like this since it relies on the order of the dofs.

Definition at line 905 of file dofmanager.C.

References oofem::IntArray::at(), oofem::Dof::computeDofTransformation(), oofem::FloatMatrix::copySubVectorRow(), giveDofWithID(), giveNumberOfDofs(), giveNumberOfPrimaryMasterDofs(), oofem::Dof::giveNumberOfPrimaryMasterDofs(), oofem::IntArray::giveSize(), oofem::FloatArray::giveSize(), hasAnySlaveDofs(), oofem::IntArray::isEmpty(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().

Referenced by oofem::EngngModel::assembleVectorFromBC(), oofem::EngngModel::assembleVectorFromDofManagers(), and computeM2GTransformation().

virtual void oofem::DofManager::drawYourself ( oofegGraphicContext gc,
TimeStep tStep 
)
inlinevirtual

Reimplemented in oofem::Node.

Definition at line 497 of file dofmanager.h.

std::vector< Dof* >:: const_iterator oofem::DofManager::end ( ) const
inline

Definition at line 160 of file dofmanager.h.

std::vector< Dof * >::const_iterator oofem::DofManager::findDofWithDofId ( DofIDItem  dofID) const

Finds index of DOF with required physical meaning of receiver.

This index can be different for different DOFManagers (user can alter dof order and type in input file).

Parameters
dofIDPhysical meaning of DOF.
Returns
Index of requested DOF. If such DOF with dofID does not exists, returns zero.

Definition at line 266 of file dofmanager.C.

References begin(), end(), and oofem::Dof::giveDofID().

Referenced by appendDof(), oofem::DofDistributedPrimaryField::applyInitialCondition(), oofem::Shell7BaseXFEM::computeOrderingArray(), oofem::PFEMParticle::drawScalar(), giveDofWithID(), oofem::SUPGElement2::giveInternalStateAtNode(), oofem::CBSElement::giveInternalStateAtNode(), oofem::PFEMElement::giveInternalStateAtNode(), oofem::TransportElement::giveInternalStateAtNode(), oofem::SUPGElement::giveInternalStateAtNode(), giveLocationArray(), giveMasterDofIDArray(), giveNumberOfPrimaryMasterDofs(), giveUnknownVector(), and oofem::EIPrimaryUnknownMapper::mapAndUpdate().

std :: map< int, int >* oofem::DofManager::giveBcMap ( )
inline

Returns map from DofIDItem to dofType.

Returns
NULL if no specific BCs are required, otherwise a map.
Deprecated:
This method of applying dirichlet b.c.s is soon to be deprecated.

Definition at line 408 of file dofmanager.h.

Referenced by oofem::Domain::createDofs().

void oofem::DofManager::giveCompleteLocationArray ( IntArray locationArray,
const UnknownNumberingScheme s 
) const
void oofem::DofManager::giveCompleteMasterDofIDArray ( IntArray dofIDArray) const

Returns the full dof ID array of receiver.

Mainly used at EngngModel level to assemble internal norms fronm DofManager contribution (typically load vector).

Parameters
dofIDArrayComplete dof ID array of receiver.

Definition at line 249 of file dofmanager.C.

References oofem::IntArray::followedBy(), oofem::Dof::giveDofID(), oofem::Dof::giveDofIDs(), giveNumberOfDofs(), hasSlaveDofs, and oofem::IntArray::resizeWithValues().

Referenced by oofem::EngngModel::assembleVectorFromDofManagers(), oofem::Beam2d::giveInternalDofManDofIDMask(), and oofem::Beam3d::giveInternalDofManDofIDMask().

void oofem::DofManager::giveCompleteUnknownVector ( FloatArray answer,
ValueModeType  mode,
TimeStep tStep 
)

Assembles the complete unknown vector in node.

Does not transform and local->global coordinate systems.

Parameters
answerComplete vector of all dof values in receiver.
modeMode of unknowns.
tStepTime step when unknown is requested.

Definition at line 737 of file dofmanager.C.

References oofem::FloatArray::at(), giveNumberOfDofs(), oofem::Dof::giveUnknown(), and oofem::FloatArray::resize().

Referenced by oofem::PrimaryField::__evaluateAt(), oofem::InteractionPFEMParticle::giveCoupledVelocities(), oofem::InteractionPFEMParticle::givePrescribedUnknownVector(), and oofem::GnuplotExportModule::outputNodeDisp().

virtual FloatArray* oofem::DofManager::giveCoordinates ( )
inlinevirtual
Returns
Pointer to node coordinate array.

Reimplemented in oofem::Node.

Definition at line 382 of file dofmanager.h.

Referenced by oofem::PrimaryField::__evaluateAt(), oofem::PrescribedGradientBCPeriodic::computeDofTransformation(), oofem::TransportGradientPeriodic::computeDofTransformation(), oofem::SolutionbasedShapeFunction::computeDofTransformation(), oofem::MixedGradientPressureDirichlet::computeDofTransformation(), oofem::PrescribedGradientBCWeak::computeDomainBoundingBox(), oofem::Node2NodeContact::computeGap(), oofem::IntElLine1IntPen::computeGlobalCoordinates(), oofem::StructuralInterfaceElement::computeGlobalCoordinates(), oofem::GeometryBasedEI::computeIntersectionPoints(), oofem::Line::computeIntersectionPoints(), oofem::Circle::computeIntersectionPoints(), oofem::PolygonLine::computeIntersectionPoints(), oofem::MaterialForceEvaluator::computeMaterialForce(), oofem::Quasicontinuum::computeStiffnessTensorOf1Link(), oofem::SolutionbasedShapeFunction::copyDofManagersToSurfaceData(), oofem::Quasicontinuum::createGlobalStiffnesMatrix(), oofem::VTKXMLExportModule::doOutput(), oofem::InternalVariableField::evaluateAt(), oofem::UniformGridField::evaluateAt(), oofem::UnstructuredGridField::evaluateAt(), oofem::Shell7BaseXFEM::EvaluateEnrFuncInDofMan(), oofem::PrescribedGradientBCWeak::findCrackBndIntersecCoord(), oofem::PrescribedGradientBCWeak::findHoleCoord(), oofem::QClinearStatic::findNearestParticle(), oofem::PrescribedGradient::give(), oofem::RotatingBoundary::give(), oofem::PrescribedGenStrainShell7::give(), oofem::UserDefDirichletBC::give(), oofem::TransportGradientDirichlet::give(), oofem::StokesFlowVelocityHomogenization::giveAreaOfRVE(), oofem::QTrPlaneStress2dXFEM::giveCompositeExportData(), oofem::PlaneStress2dXfem::giveCompositeExportData(), oofem::TrPlaneStress2dXFEM::giveCompositeExportData(), oofem::DelaunayTriangle::giveEdgeLength(), oofem::MacroLSpace::giveInternalForcesVector(), oofem::MicroMaterial::giveMacroStiffnessMatrix(), oofem::PrescribedGradientBCPeriodic::giveMasterDof(), oofem::XfemStructuralElementInterface::giveSubtriangulationCompositeExportData(), oofem::PrescribedGradientBCPeriodic::giveUnknown(), oofem::TransportGradientPeriodic::giveUnknown(), oofem::MixedGradientPressureDirichlet::giveUnknown(), oofem::IntElLine1::initializeFrom(), oofem::SolutionbasedShapeFunction::initializeSurfaceData(), oofem::Node2NodeContact::instanciateYourself(), oofem::Circle::intersects(), oofem::Circle::isInside(), oofem::SolutionbasedShapeFunction::loadProblem(), oofem::EIPrimaryUnknownMapper::mapAndUpdate(), oofem::EnrFrontLinearBranchFuncRadius::MarkNodesAsFront(), oofem::GnuplotExportModule::outputBoundaryCondition(), oofem::PLPrincipalStrain::propagateInterface(), oofem::PLHoopStressCirc::propagateInterface(), oofem::SolutionbasedShapeFunction::splitBoundaryNodeIDs(), oofem::GeometryBasedEI::updateNodeEnrMarker(), oofem::XfemElementInterface::XfemElementInterface_createEnrNmatrixAt(), and oofem::XfemElementInterface::XfemElementInterface_prepareNodesForDelaunay().

std :: map< int, int >* oofem::DofManager::giveDofTypeMap ( )
inline

Returns map from DofIDItem to dofType.

Returns
NULL if no specific dofTypes are required, otherwise a map.

Definition at line 396 of file dofmanager.h.

Referenced by oofem::Domain::createDofs(), oofem::qcNode::setAsHanging(), and oofem::qcNode::setAsRepnode().

Dof * oofem::DofManager::giveDofWithID ( int  dofID) const

Returns DOF with given dofID; issues error if not present.

Parameters
dofIDThe ID for the requested DOF.
Returns
The requested DOF.
See also
DofIDItem

Definition at line 119 of file dofmanager.C.

References end(), findDofWithDofId(), and OOFEM_ERROR.

Referenced by oofem::DofDistributedPrimaryField::applyBoundaryCondition(), oofem::PrimaryField::applyBoundaryCondition(), oofem::DofDistributedPrimaryField::applyInitialCondition(), oofem::MacroLSpace::changeMicroBoundaryConditions(), oofem::NodeErrorCheckingRule::check(), oofem::Node2NodeContactL::computeContactTractionAt(), oofem::Tr_Warp::computeEdgeLoadVectorAt(), oofem::PhaseFieldElement::computeLocationArrayOfDofIDs(), oofem::CoupledFieldsElement::computeLocationArrayOfDofIDs(), oofem::IntElLine1PF::computeLocationArrayOfDofIDs(), computeM2LTransformation(), oofem::RigidArmNode::computeMasterContribution(), oofem::InteractionLoad::computeValueAt(), oofem::CoupledFieldsElement::computeVectorOfDofIDs(), oofem::SparseNonLinearSystemNM::convertPertMap(), oofem::CohesiveSurface3d::drawDeformedGeometry(), oofem::CohesiveSurface3d::drawScalar(), oofem::tet21ghostsolid::EIPrimaryUnknownMI_computePrimaryUnknownVectorAtLocal(), oofem::Tet21Stokes::EIPrimaryUnknownMI_computePrimaryUnknownVectorAtLocal(), oofem::PrescribedGenStrainShell7::evaluateHigherOrderContribution(), oofem::VTKXMLExportModule::exportExternalForces(), oofem::VTKExportModule::getDofManPrimaryVariable(), oofem::VTKXMLExportModule::getNodalVariableFromPrimaryField(), oofem::RefinedElement::giveCompatibleBcDofArray(), oofem::Node2NodeContact::giveLocationArray(), oofem::Node2NodeContactL::giveLocationArray(), oofem::MicroMaterial::giveMacroStiffnessMatrix(), oofem::PrescribedGradientBCPeriodic::giveMasterDof(), oofem::TransportGradientPeriodic::giveMasterDof(), oofem::SimpleSlaveDof::giveMasterDof(), oofem::SolutionbasedShapeFunction::giveMasterDof(), oofem::InteractionPFEMParticle::givePrescribedUnknownVector(), givePrescribedUnknownVector(), oofem::tet21ghostsolid::giveRowTransformationMatrix(), oofem::PrescribedGradientBCPeriodic::giveUnknown(), oofem::TransportGradientPeriodic::giveUnknown(), oofem::NRSolver::initPrescribedEqs(), oofem::Hexa21Stokes::NodalAveragingRecoveryMI_computeNodalValue(), oofem::tet21ghostsolid::NodalAveragingRecoveryMI_computeNodalValue(), oofem::Tet21Stokes::NodalAveragingRecoveryMI_computeNodalValue(), oofem::Tr21Stokes::NodalAveragingRecoveryMI_computeNodalValue(), oofem::Shell7Base::NodalAveragingRecoveryMI_computeNodalValue(), oofem::GnuplotExportModule::outputBoundaryCondition(), oofem::GnuplotExportModule::outputReactionForces(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem1D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem2D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem3D(), oofem::HuertaErrorEstimator::solveRefinedElementProblem(), oofem::HuertaErrorEstimator::solveRefinedWholeProblem(), oofem::IncrementalLinearStatic::solveYourselfAt(), oofem::PrescribedGenStrainShell7::updateCoefficientMatrix(), oofem::PFEM::updateDofUnknownsDictionaryPressure(), oofem::SolutionbasedShapeFunction::updateModelWithFactors(), and oofem::DelaunayTriangulator::writeMesh().

const IntArray* oofem::DofManager::giveForcedDofIDs ( )
inline

Returns list of specific dofs that should be included in node.

Returns
NULL if no additional dofs are necessary, otherwise a list of DofIDItem's.

Definition at line 391 of file dofmanager.h.

int oofem::DofManager::giveGlobalNumber ( ) const
inline
Returns
Receivers globally unique number.

Definition at line 501 of file dofmanager.h.

Referenced by oofem::Element::addDofManager(), oofem::Domain::BuildDofManPlaceInArrayMap(), oofem::ParmetisLoadBalancer::calculateLoadTransfer(), oofem::XfemElementInterface::ComputeBOrBHMatrix(), oofem::Shell7BaseXFEM::computeEnrichedBmatrixAt(), oofem::Shell7BaseXFEM::computeEnrichedNmatrixAt(), oofem::GeometryBasedEI::computeIntersectionPoints(), oofem::XfemElementInterface::computeNCohesive(), oofem::Subdivision::createMesh(), oofem::PrescribedGradientBCWeak::createTractionMesh(), oofem::LoadBalancer::deleteRemoteDofManagers(), oofem::ProblemCommunicator::DofManCmp(), oofem::Node::drawYourself(), oofem::Shell7BaseXFEM::edgeComputeEnrichedBmatrixAt(), oofem::Shell7BaseXFEM::edgeComputeEnrichedNmatrixAt(), oofem::Shell7BaseXFEM::EvaluateEnrFuncInDofMan(), oofem::Subdivision::exchangeRemoteElements(), oofem::QTrPlaneStress2dXFEM::giveCompositeExportData(), oofem::PlaneStress2dXfem::giveCompositeExportData(), oofem::TrPlaneStress2dXFEM::giveCompositeExportData(), oofem::AuxVelocityNumberingScheme::giveDofEquationNumber(), oofem::PlaneStress2dXfem::giveDofManDofIDMask(), oofem::QTrPlaneStress2dXFEM::giveDofManDofIDMask(), oofem::TrPlaneStress2dXFEM::giveDofManDofIDMask(), oofem::Dof::giveDofManGlobalNumber(), oofem::MicroMaterial::giveMacroStiffnessMatrix(), oofem::EnrichmentItem::giveNumDofManEnrichments(), oofem::XfemStructuralElementInterface::giveSubtriangulationCompositeExportData(), oofem::Natural2GlobalOrdering::init(), oofem::NRSolver::initPrescribedEqs(), oofem::Node2NodeContact::instanciateYourself(), oofem::EnrichmentItem::isDofManEnriched(), oofem::QCFullsolveddomain::isNodeInside(), oofem::LoadBalancer::migrateLoad(), oofem::LoadBalancer::packMigratingData(), oofem::NonlocalMaterialWTP::packRemoteElements(), oofem::Subdivision::packRemoteElements(), oofem::ParmetisLoadBalancer::packSharedDmanPartitions(), oofem::qcNode::printOutputAt(), oofem::REGISTER_EnrichmentFront(), oofem::SimpleSlaveDof::saveContext(), oofem::Element::saveContext(), oofem::Domain::setDofManager(), oofem::MicroMaterial::setMacroProperties(), oofem::QClinearStatic::setRepNodesInVerticesOfInterpolationMesh(), oofem::FETICommunicator::setUpCommunicationMaps(), oofem::FETISolver::setUpCommunicationMaps(), oofem::NonlocalMaterialWTP::unpackRemoteElements(), oofem::GeometryBasedEI::updateNodeEnrMarker(), oofem::XfemElementInterface::XfemElementInterface_createEnrNmatrixAt(), and oofem::XfemElementInterface::XfemElementInterface_giveNumDofManEnrichments().

std :: map< int, int >* oofem::DofManager::giveIcMap ( )
inline

Returns map from DofIDItem to initial condition.

Returns
NULL if no specific ICs are required, otherwise a map.
Deprecated:
This method of applying i.c.s is soon to be deprecated.

Definition at line 414 of file dofmanager.h.

Referenced by oofem::Domain::createDofs().

void oofem::DofManager::giveLocationArray ( const IntArray dofIDArry,
IntArray locationArray,
const UnknownNumberingScheme s 
) const

Returns location array (array containing for each requested dof related equation number) for given numbering scheme.

Parameters
dofIDArryArray containing dof mask. This mask containing DofIDItem values (they describe physical meaning of dofs, see cltypes.h) is used to extract only required values. If dof with requested physical meaning does not exist in receiver, an error is generated and execution exits.
locationArrayReturn parameter containing required equation numbers.
sDetermines the equation numbering scheme.
See also
Element::giveDofManDofIDMask.

Definition at line 177 of file dofmanager.C.

References oofem::IntArray::at(), oofem::IntArray::clear(), end(), findDofWithDofId(), oofem::IntArray::followedBy(), oofem::UnknownNumberingScheme::giveDofEquationNumber(), oofem::IntArray::giveSize(), hasSlaveDofs, OOFEM_ERROR, and oofem::IntArray::resize().

Referenced by oofem::PrescribedMean::assemble(), oofem::PrescribedGradientBCWeak::assemble(), oofem::PrescribedGradientBCWeak::assembleExtraDisplock(), oofem::PrescribedGradientBCWeak::assembleVector(), oofem::EngngModel::assembleVectorFromBC(), oofem::EngngModel::assembleVectorFromDofManagers(), oofem::PrescribedMean::giveExternalForcesVector(), oofem::PrescribedMean::giveInternalForcesVector(), oofem::Element::giveLocationArray(), and oofem::PrescribedGradientBCWeak::giveTractionLocationArray().

void oofem::DofManager::giveMasterDofIDArray ( const IntArray dofIDArry,
IntArray masterDofIDs 
) const

Returns master dof ID array of receiver.

Parameters
dofIDArryArray containing dof mask. This mask containing DofIDItem values (they describe physical meaning of dofs, see cltypes.h) is used to extract only required values. If dof with requested physical meaning does not exist in receiver, an error is generated and execution exits.
masterDofIDsMaster dof ID array.

Definition at line 209 of file dofmanager.C.

References oofem::IntArray::clear(), end(), findDofWithDofId(), oofem::IntArray::followedBy(), hasSlaveDofs, and OOFEM_ERROR.

Referenced by oofem::Element::giveLocationArray().

bool oofem::DofManager::giveMasterDofMans ( IntArray masters)
virtual

Returns true if the receiver is linked (its slave DOFs depend on master values) to some other dof managers.

In this case, the masters array should contain the list of masters. In both serial and parallel modes, local numbers are be provided.

Parameters
mastersIndices of dof managers which receiver has slaves to.
Returns
If receiver contains only primary DOFs, false is returned.

Definition at line 842 of file dofmanager.C.

References oofem::IntArray::at(), oofem::IntArray::clear(), oofem::Dof::giveMasterDofManArray(), oofem::IntArray::giveSize(), oofem::IntArray::insertSortedOnce(), and oofem::Dof::isPrimaryDof().

Referenced by oofem::ParmetisLoadBalancer::handleMasterSlaveDofManLinks().

std :: map< int, int >* oofem::DofManager::giveMasterMap ( )
inline

Returns map from DofIDItem to dofType.

Returns
NULL if no specific BCs are required, otherwise a map.
Deprecated:
This method of applying dirichlet b.c.s is soon to be deprecated.

Definition at line 402 of file dofmanager.h.

Referenced by oofem::Domain::createDofs().

int oofem::DofManager::giveNumberOfDofs ( ) const
Returns
Total number of dofs managed by receiver.

Definition at line 279 of file dofmanager.C.

References dofArray.

Referenced by oofem::PetscSparseMtrx::buildInternalStructure(), oofem::RigidArmNode::checkConsistency(), oofem::Node::computeL2GTransformation(), computeLoadVector(), oofem::PhaseFieldElement::computeLocationArrayOfDofIDs(), oofem::CoupledFieldsElement::computeLocationArrayOfDofIDs(), oofem::StructuralInterfaceElementPhF::computeLocationArrayOfDofIDs(), oofem::IntElLine1PF::computeLocationArrayOfDofIDs(), computeM2LTransformation(), oofem::RigidArmNode::computeMasterContribution(), oofem::PlaneStress2dXfem::computeNumberOfDofs(), oofem::QTrPlaneStress2dXFEM::computeNumberOfDofs(), oofem::TrPlaneStress2dXFEM::computeNumberOfDofs(), oofem::StructuralInterfaceElementPhF::computeNumberOfDofs(), oofem::ElementSide::computeTransformation(), oofem::HuertaErrorEstimator::estimateError(), oofem::RefinedElement::giveBcDofArray1D(), oofem::RefinedElement::giveBcDofArray2D(), oofem::RefinedElement::giveBcDofArray3D(), giveCompleteLocationArray(), giveCompleteMasterDofIDArray(), giveCompleteUnknownVector(), oofem::MacroLSpace::giveInternalForcesVector(), oofem::Natural2GlobalOrdering::init(), oofem::LSPrimaryVariableMapper::mapPrimaryVariables(), oofem::RigidArmNode::postInitialize(), oofem::MicroMaterial::setMacroProperties(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem1D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem2D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem3D(), oofem::HuertaErrorEstimator::solveRefinedElementProblem(), oofem::HuertaErrorEstimator::solveRefinedPatchProblem(), and oofem::HuertaErrorEstimator::solveRefinedWholeProblem().

int oofem::DofManager::giveNumberOfPrimaryMasterDofs ( const IntArray dofIDArray) const

Returns the number of primary dofs on which receiver dofs (given in dofArray) depend on.

If receiver has only primary dofs, the answer is the size of dofArray.

Parameters
dofArrayArray with Dof IDs.
Returns
Number of primary dofs from given array.

Definition at line 303 of file dofmanager.C.

References end(), findDofWithDofId(), oofem::IntArray::giveSize(), hasSlaveDofs, and OOFEM_ERROR.

Referenced by computeM2LTransformation(), and oofem::Element::computeNumberOfPrimaryMasterDofs().

dofManagerParallelMode oofem::DofManager::giveParallelMode ( ) const
inline

Return dofManagerParallelMode of receiver.

Definition at line 512 of file dofmanager.h.

Referenced by oofem::OutputManager::_testDofManOutput(), oofem::ActiveDof::askNewEquationNumber(), oofem::MasterDof::askNewEquationNumber(), oofem::EngngModel::assembleVectorFromDofManagers(), oofem::NodeErrorCheckingRule::check(), oofem::HangingNode::checkConsistency(), oofem::qcNode::checkConsistency(), oofem::Subdivision::createMesh(), oofem::ParmetisLoadBalancer::determineDofManState(), oofem::OutputManager::doDofManOutput(), oofem::Node::drawYourself(), oofem::DirectErrorIndicatorRC::exchangeDofManDensities(), oofem::DirectErrorIndicatorRC::exchangeDofManIndicatorVals(), oofem::Subdivision::exchangeRemoteElements(), oofem::MasterDof::giveUnknown(), oofem::MasterDof::hasBc(), oofem::Natural2GlobalOrdering::init(), oofem::ParallelOrdering::isLocal(), oofem::ParallelOrdering::isShared(), oofem::ParmetisLoadBalancer::labelDofManagers(), oofem::EIPrimaryUnknownMapper::mapAndUpdate(), oofem::LoadBalancer::migrateLoad(), oofem::NonlocalMaterialWTP::packRemoteElements(), oofem::Subdivision::packRemoteElements(), oofem::ParmetisLoadBalancer::packSharedDmanPartitions(), oofem::LoadBalancer::printStatistics(), oofem::SPRNodalRecoveryModel::recoverValues(), oofem::FETICommunicator::setUpCommunicationMaps(), oofem::FETISolver::setUpCommunicationMaps(), oofem::NodeCommunicator::setUpCommunicationMaps(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem2D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem3D(), oofem::HuertaErrorEstimator::solveRefinedPatchProblem(), oofem::NlDEIDynamic::solveYourselfAt(), and oofem::EngngModel::unpackDofManagers().

int oofem::DofManager::givePartitionsConnectivitySize ( )
void oofem::DofManager::givePrescribedUnknownVector ( FloatArray answer,
const IntArray dofMask,
ValueModeType  mode,
TimeStep tStep 
)

Assembles the vector of prescribed unknowns in nodal c.s for given dofs of receiver.

This vector may have size different from number of dofs requested, because some dofs may depend on other dofs. Default implementation uses Dof::giveBcValue and Dof::hasBc service.

Parameters
answerResult (in nodal cs.)
dofMaskArray containing dof mask. This mask containing DofIDItem values (they describe physical meaning of dofs, see cltypes.h) is used to extract only required values. If dof with requested physical meaning does not exist in receiver, an error is generated and execution exits.
modeMode of unknown (e.g, total value, velocity or acceleration of unknown).
tStepTime step when unknown requested. See documentation of particular EngngModel class for valid tStep values (most implementation can return only values for current and possibly for previous time step).
See also
Dof::giveBcValue
Dof::hasBc
Todo:
Remove all usage of this. Just ask for the unknown vector instead, they are the same.

Definition at line 748 of file dofmanager.C.

References oofem::FloatArray::at(), computeL2GTransformation(), oofem::Dof::giveBcValue(), giveDofWithID(), oofem::IntArray::giveSize(), oofem::FloatArray::resize(), and oofem::FloatArray::rotatedWith().

Referenced by oofem::Element::computeVectorOfPrescribed().

void oofem::DofManager::giveUnknownVector ( FloatArray answer,
const IntArray dofMask,
ValueModeType  mode,
TimeStep tStep,
bool  padding = false 
)

Assembles the vector of unknowns in global c.s for given dofs of receiver.

Parameters
answerResult (in global c.s.)
dofMaskArray containing DOF mask. This mask containing DofIDItem values (they describe physical meaning of dofs, see cltypes.h) is used to extract only required values. If dof with requested physical meaning does not exist in receiver, an error is generated and execution exits.
modeMode of unknown (e.g, total value, velocity or acceleration of unknown).
tStepTime step when unknown requested. See documentation of particular EngngModel class for valid tStep values (most implementation can return only values for current and possibly for previous time step).
paddingDetermines if zero value should be inserted when a dof isn't found (otherwise it is just skipped).
See also
Dof::giveUnknown

Definition at line 685 of file dofmanager.C.

References oofem::FloatArray::at(), computeL2GTransformation(), end(), findDofWithDofId(), oofem::IntArray::giveSize(), oofem::FloatArray::resize(), oofem::FloatArray::resizeWithValues(), and oofem::FloatArray::rotatedWith().

Referenced by oofem::PrimaryField::__evaluateAt(), oofem::PrescribedGradientBCWeak::assembleVector(), oofem::Element::computeBoundaryVectorOf(), oofem::Node2NodeContact::computeGap(), oofem::PrescribedGradientBCWeak::computeIntForceGPContrib(), oofem::Element::computeVectorOf(), oofem::LEPlic::doLagrangianPhase(), oofem::Quad1_ht::drawScalar(), oofem::PrescribedMean::giveInternalForcesVector(), oofem::SolutionbasedShapeFunction::giveUnknown(), oofem::HuertaErrorEstimator::solveRefinedElementProblem(), oofem::HuertaErrorEstimator::solveRefinedWholeProblem(), and oofem::LevelSetPCS::updatePosition().

void oofem::DofManager::giveUnknownVector ( FloatArray answer,
const IntArray dofMask,
PrimaryField field,
ValueModeType  mode,
TimeStep tStep,
bool  padding = false 
)

Assembles the vector of unknowns of given filed in global c.s for given dofs of receiver.

Parameters
answerResult (in nodal cs.)
dofMaskArray containing DOF mask. This mask containing DofIDItem values (they describe physical meaning of dofs, see cltypes.h) is used to extract only required values. If dof with requested physical meaning does not exist in receiver, an error is generated and execution exits.
fieldPrimary field.
modeMode of unknown (e.g, total value, velocity or acceleration of unknown).
tStepTime step when unknown is requested. See documentation of particular EngngModel class for valid tStep values (most implementation can return only values for current and possibly for previous time step).
paddingDetermines if zero value should be inserted when a dof isn't found (otherwise it is just skipped).
See also
Dof::giveUnknown

Definition at line 711 of file dofmanager.C.

References oofem::FloatArray::at(), computeL2GTransformation(), end(), findDofWithDofId(), oofem::IntArray::giveSize(), oofem::FloatArray::resize(), oofem::FloatArray::resizeWithValues(), and oofem::FloatArray::rotatedWith().

void oofem::DofManager::giveUnknownVectorOfType ( FloatArray answer,
UnknownType  ut,
ValueModeType  mode,
TimeStep tStep 
)

Constructs the requested vector by assembling e.g.

[D_u, D_v, D_w] or [V_u, V_v, V_w]. If for example D_v or V_w doesn't exist, then zero value is inserted.

Note
Error is produced if the unknown type doesn't represent a vector.
Parameters
answerThe requested vector.
utThe unknown type to assemble.
modeValue mode (total, incremental, etc.)
tStepTime step to evaluate at.

Definition at line 766 of file dofmanager.C.

References oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatArray::beProductOf(), computeL2GTransformation(), oofem::Dof::giveDofID(), oofem::Dof::giveUnknown(), OOFEM_ERROR, oofem::FloatArray::resize(), and oofem::FloatArray::zero().

Referenced by oofem::Node::giveUpdatedCoordinates().

bool oofem::DofManager::hasAnySlaveDofs ( )
virtual
Returns
True if receiver contains slave dofs

Definition at line 830 of file dofmanager.C.

References oofem::Dof::isPrimaryDof().

Referenced by computeM2LTransformation(), oofem::ParmetisLoadBalancer::handleMasterSlaveDofManLinks(), and requiresTransformation().

virtual bool oofem::DofManager::hasCoordinates ( )
inlinevirtual
IRResultType oofem::DofManager::initializeFrom ( InputRecord ir)
virtual

Initializes receiver according to object description stored in input record.

This function is called immediately after creating object using constructor. Input record can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record.

See also
IR_GIVE_FIELD
IR_GIVE_OPTIONAL_FIELD
Parameters
irInput record to initialize from.
Returns
IRResultType
Todo:
This is unnecessary, we should just check if user has supplied a dofidmask field or not and just drop "numberOfDofs"). It is left for now because it will give lots of warnings otherwise, but it is effectively ignored.
Todo:
This should eventually be removed, still here to preserve backwards compatibility:

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::Node, oofem::RigidArmNode, oofem::qcNode, oofem::HangingNode, oofem::PFEMParticle, oofem::ElementSide, oofem::SlaveNode, oofem::InteractionPFEMParticle, oofem::ElementDofManager, and oofem::Particle.

Definition at line 326 of file dofmanager.C.

References _IFT_DofManager_bc, _IFT_DofManager_boundaryflag, _IFT_DofManager_dofidmask, _IFT_DofManager_doftypemask, _IFT_DofManager_ic, _IFT_DofManager_load, _IFT_DofManager_mastermask, _IFT_DofManager_nullflag, _IFT_DofManager_partitions, _IFT_DofManager_remoteflag, _IFT_DofManager_sharedflag, oofem::IntArray::at(), oofem::IntArray::clear(), oofem::IntArray::contains(), dofBCmap, dofICmap, dofidmask, oofem::DofManager_local, oofem::DofManager_null, oofem::DofManager_remote, oofem::DofManager_shared, dofMastermap, dofTypemap, oofem::FEMComponent::domain, oofem::Domain::giveDefaultNodeDofIDArry(), oofem::IntArray::giveSize(), oofem::InputRecord::hasField(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_OK, isBoundaryFlag, loadArray, mBC, OOFEM_ERROR, parallel_mode, and partitions.

Referenced by oofem::ElementDofManager::initializeFrom(), oofem::ElementSide::initializeFrom(), and oofem::Node::initializeFrom().

bool oofem::DofManager::isBoundary ( )
inline
Returns
True if dofmanager is on boundary.

Definition at line 426 of file dofmanager.h.

Referenced by oofem::Subdivision::createMesh(), and oofem::SPRNodalRecoveryModel::determinePatchAssemblyPoints().

virtual bool oofem::DofManager::isDofTypeCompatible ( dofType  type) const
inlinevirtual

Returns true if dof of given type is allowed to be associated to receiver.

Reimplemented in oofem::Node, oofem::RigidArmNode, oofem::ElementSide, oofem::qcNode, oofem::HangingNode, oofem::SlaveNode, and oofem::ElementDofManager.

Definition at line 451 of file dofmanager.h.

Referenced by oofem::Domain::createDofs().

bool oofem::DofManager::isNull ( )
inline
bool oofem::DofManager::isShared ( )
inline
void oofem::DofManager::mergePartitionList ( IntArray _p)

Merges receiver partition list with given lists.

Definition at line 965 of file dofmanager.C.

References oofem::IntArray::at(), oofem::IntArray::giveSize(), oofem::IntArray::insertOnce(), and partitions.

void oofem::DofManager::postInitialize ( )
virtual

Performs post-initialization such like checking if there are any slave dofs etc.

Reimplemented in oofem::RigidArmNode, oofem::qcNode, oofem::HangingNode, and oofem::SlaveNode.

Definition at line 862 of file dofmanager.C.

References hasSlaveDofs, and oofem::Dof::isPrimaryDof().

Referenced by oofem::SlaveNode::postInitialize(), oofem::HangingNode::postInitialize(), oofem::qcNode::postInitialize(), and oofem::RigidArmNode::postInitialize().

void oofem::DofManager::printOutputAt ( FILE *  file,
TimeStep tStep 
)
virtual

Prints output of receiver to stream, for given time step.

This is used for output into the standard output file.

Parameters
fileFile pointer to print to.
tStepTime step to write for.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::qcNode, oofem::PFEMParticle, and oofem::InteractionPFEMParticle.

Definition at line 510 of file dofmanager.C.

References oofem::FEMComponent::giveClassName(), oofem::FEMComponent::giveDomain(), oofem::Domain::giveEngngModel(), giveLabel(), oofem::FEMComponent::giveNumber(), and oofem::EngngModel::printDofOutputAt().

Referenced by oofem::OutputManager::doDofManOutput(), oofem::EngngModel::outputNodes(), oofem::PFEMParticle::printOutputAt(), and oofem::Element::printOutputAt().

void oofem::DofManager::printYourself ( )
virtual

Prints receiver state on stdout. Useful for debugging.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::Node, oofem::ElementSide, and oofem::ElementDofManager.

Definition at line 521 of file dofmanager.C.

References loadArray, oofem::FEMComponent::number, oofem::IntArray::printYourself(), and oofem::Dof::printYourself().

void oofem::DofManager::removeDof ( DofIDItem  id)

Removes Dof with given id from dofArray.

Parameters
id

Definition at line 151 of file dofmanager.C.

References begin(), dofArray, oofem::Dof::giveDofID(), and OOFEM_WARNING.

void oofem::DofManager::removePartitionFromList ( int  _part)
inline

Removes given partition from receiver list.

Definition at line 523 of file dofmanager.h.

References oofem::IntArray::erase(), and oofem::IntArray::findFirstIndexOf().

bool oofem::DofManager::requiresTransformation ( )
virtual

Indicates, whether dofManager requires the transformation.

Returns
Nonzero if transformation is necessary, even for single dof.

Reimplemented in oofem::Node, and oofem::ElementSide.

Definition at line 944 of file dofmanager.C.

References hasAnySlaveDofs().

Referenced by oofem::Element::computeDofTransformationMatrix().

contextIOResultType oofem::DofManager::restoreContext ( DataStream stream,
ContextMode  mode,
void *  obj = NULL 
)
virtual

Restores the receiver state previously written in stream.

See also
saveContext
Parameters
streamInput stream.
modeDetermines amount of info available in stream (state, definition, ...).
objSpecial parameter for sending extra information.
Returns
contextIOResultType.
Exceptions
throwsan ContextIOERR exception if error encountered.
Todo:
Smart pointers would be nicer here

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::Node.

Definition at line 611 of file dofmanager.C.

References appendDof(), oofem::CIO_IOERR, oofem::CIO_OK, oofem::classFactory, CM_Definition, oofem::ClassFactory::createDof(), dofArray, globalNumber, hasSlaveDofs, isBoundaryFlag, loadArray, parallel_mode, partitions, oofem::DataStream::read(), oofem::FEMComponent::restoreContext(), oofem::Dof::restoreContext(), oofem::IntArray::restoreYourself(), and THROW_CIOERR.

Referenced by oofem::Node::restoreContext(), oofem::LoadBalancer::unpackMigratingData(), oofem::NonlocalMaterialWTP::unpackRemoteElements(), and oofem::Subdivision::unpackRemoteElements().

contextIOResultType oofem::DofManager::saveContext ( DataStream stream,
ContextMode  mode,
void *  obj = NULL 
)
virtual

Stores receiver state to output stream.

Parameters
streamOutput stream.
modeDetermines amount of info required in stream (state, definition, ...).
objSpecial parameter, used only to send particular integration point to material class version of this method.
Returns
contextIOResultType.
Exceptions
throwsan ContextIOERR exception if error encountered.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::Node.

Definition at line 543 of file dofmanager.C.

References oofem::CIO_IOERR, oofem::CIO_OK, CM_Definition, globalNumber, hasSlaveDofs, isBoundaryFlag, loadArray, parallel_mode, partitions, oofem::FEMComponent::saveContext(), oofem::IntArray::storeYourself(), THROW_CIOERR, and oofem::DataStream::write().

Referenced by oofem::LoadBalancer::packMigratingData(), oofem::NonlocalMaterialWTP::packRemoteElements(), oofem::Subdivision::packRemoteElements(), and oofem::Node::saveContext().

void oofem::DofManager::setBoundaryFlag ( bool  isBoundary)
inline

Sets the boundary flag.

Parameters
isBoundaryDetermines if receiver is on the boundary.

Definition at line 431 of file dofmanager.h.

Referenced by oofem::Subdivision::createMesh().

void oofem::DofManager::setGlobalNumber ( int  newNumber)
inline

Sets receiver global number.

Parameters
numberNew global number for receiver.

Definition at line 507 of file dofmanager.h.

Referenced by oofem::Subdivision::createMesh(), oofem::LoadBalancer::unpackMigratingData(), and oofem::Subdivision::unpackRemoteElements().

void oofem::DofManager::setLoadArray ( IntArray load)

Sets the array of applied loadings of the receiver.

Parameters
loadArray with indices to the applied loads.

Definition at line 91 of file dofmanager.C.

References loadArray.

Referenced by oofem::Subdivision::createMesh().

void oofem::DofManager::setNumberOfDofs ( int  _ndofs)

Sets number of dofs of the receiver; Deallocates existing DOFs; Resizes the dofArray accordingly.

Definition at line 286 of file dofmanager.C.

References dofArray.

Referenced by oofem::Domain::createDofs(), oofem::Subdivision::createMesh(), and oofem::T3DInterface::t3d_2_OOFEM().

void oofem::DofManager::setPartitionList ( const IntArray _p)
inline
void oofem::DofManager::updateLocalNumbering ( EntityRenumberingFunctor f)
virtual

Local renumbering support.

For some tasks (parallel load balancing, for example) it is necessary to renumber the entities. The various FEM components (such as nodes or elements) typically contain links to other entities in terms of their local numbers, etc. This service allows to update these relations to reflect updated numbering. The renumbering function is passed, which is supposed to return an updated number of specified entity type based on old number.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::RigidArmNode, and oofem::SlaveNode.

Definition at line 950 of file dofmanager.C.

References dofMastermap, oofem::ERS_DofManager, and oofem::Dof::updateLocalNumbering().

Referenced by oofem::SlaveNode::updateLocalNumbering(), and oofem::RigidArmNode::updateLocalNumbering().

void oofem::DofManager::updateYourself ( TimeStep tStep)
virtual

Updates receiver after equilibrium in time step has been reached.

Parameters
tStepActive time step.

Reimplemented in oofem::Node, oofem::PFEMParticle, and oofem::InteractionPFEMParticle.

Definition at line 534 of file dofmanager.C.

References oofem::Dof::updateYourself().

Referenced by oofem::Node::updateYourself().

Member Data Documentation

std::vector< Dof * > oofem::DofManager::dofArray
protected
std :: map< int, int >* oofem::DofManager::dofBCmap
protected

Map from DofIDItem to bc (to be removed).

Definition at line 149 of file dofmanager.h.

Referenced by DofManager(), initializeFrom(), and ~DofManager().

std :: map< int, int >* oofem::DofManager::dofICmap
protected

Map from DofIDItem to ic (to be removed).

Definition at line 151 of file dofmanager.h.

Referenced by DofManager(), initializeFrom(), and ~DofManager().

IntArray* oofem::DofManager::dofidmask
protected
std :: map< int, int >* oofem::DofManager::dofMastermap
protected

Map from DofIDItem to master node.

Definition at line 147 of file dofmanager.h.

Referenced by DofManager(), giveInputRecord(), initializeFrom(), updateLocalNumbering(), and ~DofManager().

std :: map< int, int >* oofem::DofManager::dofTypemap
protected

Map from DofIDItem to dofType.

Definition at line 145 of file dofmanager.h.

Referenced by DofManager(), giveInputRecord(), initializeFrom(), oofem::qcNode::setAsHanging(), and ~DofManager().

int oofem::DofManager::globalNumber
protected

In parallel mode, globalNumber contains globally unique DoFManager number.

The component number, inherited from FEMComponent class contains local domain number.

Definition at line 132 of file dofmanager.h.

Referenced by restoreContext(), and saveContext().

bool oofem::DofManager::hasSlaveDofs
protected
bool oofem::DofManager::isBoundaryFlag
protected

Indicates if dofManager is boundary (true boundary or on boundary between regions) or interior.

This information is required by some recovery techniques.

Definition at line 124 of file dofmanager.h.

Referenced by DofManager(), giveInputRecord(), initializeFrom(), restoreContext(), and saveContext().

IntArray oofem::DofManager::mBC
protected

Definition at line 154 of file dofmanager.h.

Referenced by initializeFrom().

IntArray oofem::DofManager::partitions
protected

List of partition sharing the shared dof manager or remote partition containing remote dofmanager counterpart.

Definition at line 140 of file dofmanager.h.

Referenced by giveInputRecord(), givePartitionsConnectivitySize(), initializeFrom(), isLocal(), mergePartitionList(), restoreContext(), and saveContext().


The documentation for this class was generated from the following files:

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:34 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011