OOFEM 3.0
Loading...
Searching...
No Matches
oofem::DofManager Class Reference

#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)
virtual ~DofManager ()
 Destructor.
Dof management methods
DofgiveDofWithID (int dofID) const
int giveNumberOfDofs () const
void askNewEquationNumbers (TimeStep *tStep)
 Renumbers all contained DOFs.
int giveNumberOfPrimaryMasterDofs (const IntArray &dofIDArray) const
void giveLocationArray (const IntArray &dofIDArry, IntArray &locationArray, const UnknownNumberingScheme &s) const
void giveMasterDofIDArray (const IntArray &dofIDArry, IntArray &masterDofIDs) const
void giveCompleteLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s) const
void giveCompleteMasterDofIDArray (IntArray &dofIDArray) const
std::vector< Dof * >::const_iterator findDofWithDofId (DofIDItem dofID) const
void giveUnknownVector (FloatArray &answer, const IntArray &dofMask, ValueModeType mode, TimeStep *tStep, bool padding=false)
void giveUnknownVector (FloatArray &answer, const IntArray &dofMask, PrimaryField &field, ValueModeType mode, TimeStep *tStep, bool padding=false)
void giveCompleteUnknownVector (FloatArray &answer, ValueModeType mode, TimeStep *tStep)
void giveUnknownVectorOfType (FloatArray &answer, UnknownType ut, ValueModeType mode, TimeStep *tStep)
virtual void givePrescribedUnknownVector (FloatArray &answer, const IntArray &dofMask, ValueModeType mode, TimeStep *tStep)
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)
virtual bool computeL2GTransformation (FloatMatrix &answer, const IntArray &dofIDArry)
virtual bool computeM2LTransformation (FloatMatrix &answer, const IntArray &dofIDArry)
virtual bool requiresTransformation ()
Load related functions
virtual void computeLoadVector (FloatArray &answer, Load *load, CharType type, TimeStep *tStep, ValueModeType mode)
IntArraygiveLoadArray ()
void setLoadArray (IntArray &load)
Position query functions
double giveCoordinate (int i) const
const FloatArraygiveCoordinates () const
void setCoordinates (const FloatArray &coords)
 Set coordinates.
Functions necessary for dof creation. All optional.
const IntArraygiveForcedDofIDs ()
std ::map< int, int > * giveDofTypeMap ()
std ::map< int, int > * giveMasterMap ()
std ::map< int, int > * giveBcMap ()
std ::map< int, int > * giveIcMap ()
void printOutputAt (FILE *file, TimeStep *tStep) override
virtual void updateYourself (TimeStep *tStep)
bool isBoundary ()
void setBoundaryFlag (bool isBoundary)
virtual bool hasAnySlaveDofs ()
virtual bool giveMasterDofMans (IntArray &masters)
void initializeFrom (InputRecord &ir) override
void initializeFrom (InputRecord &ir, int priority) override
void initializeFinish () override
void postInitialize () override
 Performs post initialization steps. Called after all components are created and initialized.
void giveInputRecord (DynamicInputRecord &input) override
void printYourself () override
 Prints receiver state on stdout. Useful for debugging.
void saveContext (DataStream &stream, ContextMode mode) override
void restoreContext (DataStream &stream, ContextMode mode) override
virtual bool isDofTypeCompatible (dofType type) const
 Returns true if dof of given type is allowed to be associated to receiver.
void updateLocalNumbering (EntityRenumberingFunctor &f) override
Advanced functions
void setNumberOfDofs (int _ndofs)
void appendDof (Dof *dof)
void removeDof (DofIDItem id)
bool hasDofID (DofIDItem id) const
virtual void drawYourself (oofegGraphicContext &gc, TimeStep *tStep)
int giveGlobalNumber () const
int giveLabel () const
void setGlobalNumber (int newNumber)
dofManagerParallelMode giveParallelMode () const
void setParallelMode (dofManagerParallelMode _mode)
const IntArraygivePartitionList ()
void setPartitionList (const IntArray *_p)
void removePartitionFromList (int _part)
 Removes given partition from receiver list.
void mergePartitionList (IntArray &_p)
 Merges receiver partition list with given lists.
int givePartitionsConnectivitySize ()
bool isLocal ()
 Returns true if receiver is locally maintained.
bool isShared ()
 Returns true if receiver is shared.
bool isNull ()
 Returns true if receiver is shared.
const char * giveClassName () const override
const char * giveInputRecordName () const override
Public Member Functions inherited from oofem::FEMComponent
 FEMComponent (int n, Domain *d)
virtual ~FEMComponent ()=default
 Virtual destructor.
DomaingiveDomain () const
virtual void setDomain (Domain *d)
int giveNumber () const
void setNumber (int num)
virtual int checkConsistency ()
virtual InterfacegiveInterface (InterfaceType t)
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros).

Static Public Attributes

static ParamKey IPK_DofManager_dofidmask
static ParamKey IPK_DofManager_load
static ParamKey IPK_DofManager_bc
static ParamKey IPK_DofManager_ic
static ParamKey IPK_DofManager_mastermask
static ParamKey IPK_DofManager_doftypemask
static ParamKey IPK_DofManager_boundaryflag
static ParamKey IPK_DofManager_globnum
static ParamKey IPK_DofManager_partitions
static ParamKey IPK_DofManager_sharedflag
static ParamKey IPK_DofManager_remoteflag
static ParamKey IPK_DofManager_nullflag

Protected Attributes

FloatArray coordinates
 Array storing nodal coordinates.
std::vector< Dof * > dofArray
 Array of DOFs.
IntArray loadArray
 List of applied loads.
bool isBoundaryFlag
bool hasSlaveDofs
 Flag indicating whether receiver has slave DOFs.
int globalNumber
dofManagerParallelMode parallel_mode
IntArray partitions
IntArray dofidmask
 List of additional dof ids to include.
std ::map< int, int > dofTypemap
 Map from DofIDItem to dofType.
std ::map< int, int > dofMastermap
 Map from DofIDItem to master node.
std ::map< int, int > dofBCmap
 Map from DofIDItem to bc (to be removed).
std ::map< int, int > dofICmap
 Map from DofIDItem to ic (to be removed).
IntArray mBC
Protected Attributes inherited from oofem::FEMComponent
int number
 Component number.
Domaindomain
 Link to domain object, useful for communicating with other FEM components.

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 99 of file dofmanager.h.

Constructor & Destructor Documentation

◆ DofManager()

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

◆ ~DofManager()

oofem::DofManager::~DofManager ( )
virtual

Destructor.

Definition at line 82 of file dofmanager.C.

References dofArray.

Member Function Documentation

◆ appendDof()

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 142 of file dofmanager.C.

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

Referenced by oofem::Domain::createDofs(), oofem::MPMSymbolicTerm::initializeCell(), and restoreContext().

◆ askNewEquationNumbers()

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

Renumbers all contained DOFs.

Definition at line 303 of file dofmanager.C.

References oofem::Dof::askNewEquationNumber().

◆ begin() [1/2]

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

Definition at line 161 of file dofmanager.h.

References dofArray.

Referenced by oofem::Shell7BaseXFEM::computeOrderingArray(), findDofWithDofId(), and removeDof().

◆ begin() [2/2]

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

Definition at line 163 of file dofmanager.h.

References dofArray.

◆ computeL2GTransformation()

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 858 of file dofmanager.C.

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

◆ computeLoadVector()

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 105 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().

◆ computeM2GTransformation()

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 833 of file dofmanager.C.

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

◆ 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 864 of file dofmanager.C.

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

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

◆ drawYourself()

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

Reimplemented in oofem::Node.

Definition at line 511 of file dofmanager.h.

References gc.

◆ end() [1/2]

◆ end() [2/2]

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

Definition at line 164 of file dofmanager.h.

References dofArray.

◆ findDofWithDofId()

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

◆ giveBcMap()

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 421 of file dofmanager.h.

References dofBCmap.

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

◆ giveClassName()

const char * oofem::DofManager::giveClassName ( ) const
inlineoverridevirtual

◆ giveCompleteLocationArray()

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

Returns full location array of receiver containing equation numbers of all dofs of receiver. Their order is specific to every DofManager. Mainly used at EngngModel level to assemble DofManager contribution (typically load vector).

Parameters
locationArrayComplete location array of receiver.
sDetermines the equation numbering scheme.

Definition at line 237 of file dofmanager.C.

References oofem::IntArray::followedBy(), oofem::UnknownNumberingScheme::giveDofEquationNumber(), oofem::Dof::giveEquationNumbers(), giveNumberOfDofs(), hasSlaveDofs, oofem::IntArray::resize(), and oofem::IntArray::resizeWithValues().

Referenced by oofem::EngngModel::assembleVectorFromDofManagers(), and oofem::FETICommunicator::setUpCommunicationMaps().

◆ giveCompleteMasterDofIDArray()

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 257 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::VTKXMLLatticeExportModule::exportPrimaryVars(), and oofem::VTKXMLPeriodicExportModule::exportPrimaryVars().

◆ giveCompleteUnknownVector()

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 709 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().

◆ giveCoordinate()

double oofem::DofManager::giveCoordinate ( int i) const
inline
Returns
The i-th coordinate.

Definition at line 383 of file dofmanager.h.

References oofem::FloatArray::at(), and oofem::FloatArray::giveSize().

Referenced by oofem::TR1_2D_PFEM::checkConsistency(), oofem::MatlabExportModule::computeArea(), oofem::AxisymElement::computeBHmatrixAt(), oofem::AxisymElement::computeBmatrixAt(), oofem::CohesiveSurface3d::computeBmatrixAt(), oofem::Truss2d::computeBmatrixAt(), oofem::HTSelement::computeCenterOfGravity(), oofem::DelaunayTriangle::computeCircumcircle(), oofem::TR1_2D_PFEM::computeCriticalTimeStep(), oofem::FRCFCMNL::computeElementCentroid(), oofem::BondLink3d::computeGeometryProperties(), oofem::BondLink3dBoundary::computeGeometryProperties(), oofem::Lattice3d::computeGeometryProperties(), oofem::Lattice3d_mt::computeGeometryProperties(), oofem::Lattice3dBoundary::computeGeometryProperties(), oofem::Lattice3dBoundaryTruss::computeGeometryProperties(), oofem::LatticeBeam3d::computeGeometryProperties(), oofem::LatticeBeam3dBoundary::computeGeometryProperties(), oofem::LatticeLink3d::computeGeometryProperties(), oofem::LatticeLink3dBoundary::computeGeometryProperties(), oofem::CCTPlate3d::computeGlobalCoordinates(), oofem::DKTPlate3d::computeGlobalCoordinates(), oofem::LIBeam2dNL::computeGlobalCoordinates(), oofem::LIBeam3d2::computeGlobalCoordinates(), oofem::LIBeam3dNL2::computeGlobalCoordinates(), oofem::LIBeam3dNL::computeGlobalCoordinates(), oofem::MITC4Shell::computeGlobalCoordinates(), oofem::Truss2d::computeGlobalCoordinates(), oofem::Lattice2d_mt::computeInternalSourceRhsVectorAt(), oofem::Line::computeIntersectionPoints(), oofem::Beam2d::computeLength(), oofem::Beam3d::computeLength(), oofem::LIBeam2d::computeLength(), oofem::LIBeam2dNL::computeLength(), oofem::LIBeam3d2::computeLength(), oofem::LIBeam3d::computeLength(), oofem::LIBeam3dNL2::computeLength(), oofem::LIBeam3dNL::computeLength(), oofem::Truss2d::computeLength(), oofem::SurfaceTensionBoundaryCondition::computeLoadVectorFromElement(), oofem::Quasicontinuum::computeStiffnessTensorOf1Link(), oofem::SurfaceTensionBoundaryCondition::computeTangentFromElement(), oofem::LIBeam3dNL2::computeXdVector(), oofem::LIBeam3dNL::computeXdVector(), oofem::Quasicontinuum::createGlobalStiffnesMatrix(), oofem::FreemInterface::createInput(), oofem::T3DInterface::createInput(), oofem::T3DInterface::createInput(), oofem::DofManExportModule::doOutput(), oofem::Element::drawAnnotation(), oofem::CCTPlate::drawRawGeometry(), oofem::CohesiveSurface3d::drawRawGeometry(), oofem::DKTPlate::drawRawGeometry(), oofem::LIBeam2dNL::drawRawGeometry(), oofem::LIBeam3d2::drawRawGeometry(), oofem::LIBeam3dNL2::drawRawGeometry(), oofem::LIBeam3dNL::drawRawGeometry(), oofem::QDKTPlate::drawRawGeometry(), oofem::Truss1d::drawRawGeometry(), oofem::Truss2d::drawRawGeometry(), oofem::Truss3d::drawRawGeometry(), oofem::CCTPlate::drawScalar(), oofem::CohesiveSurface3d::drawScalar(), oofem::DKTPlate::drawScalar(), oofem::LIBeam3d2::drawScalar(), oofem::PFEMParticle::drawScalar(), oofem::QDKTPlate::drawScalar(), oofem::Truss1d::drawScalar(), oofem::Node::drawYourself(), oofem::InsertNode::evaluate(), oofem::CohesiveSurface3d::evaluateCenter(), oofem::CohesiveSurface3d::evaluateLocalCoordinateSystem(), oofem::VTKExportModule::exportIntVarAs(), oofem::VTKXMLPeriodicExportModule::exportPrimaryVars(), oofem::LEPlic::findCellLineConstant(), oofem::VTKBaseExportModule::getNodalVariableFromIS(), oofem::UserDefDirichletBC::give(), oofem::Beam3d::giveCompositeExportData(), oofem::LIBeam3d2::giveCurrentLength(), oofem::Beam3d::giveInternalForcesVectorAtPoint(), oofem::StructuralElement::giveInternalStateAtNode(), oofem::CohesiveSurface3d::giveLength(), oofem::Lattice2d::giveLength(), oofem::Lattice2d_mt::giveLength(), oofem::Lattice2dBoundary::giveLength(), oofem::LIBeam3d2::giveLocalCoordinateSystem(), oofem::LIBeam3d::giveLocalCoordinateSystem(), oofem::LIBeam3dBoundary::giveLocalCoordinateSystem(), oofem::LIBeam3dNL2::giveLocalCoordinateSystem(), oofem::LIBeam3dNL::giveLocalCoordinateSystem(), oofem::Beam2d::givePitch(), oofem::Lattice2d::givePitch(), oofem::Lattice2dBoundary::givePitch(), oofem::LIBeam2d::givePitch(), oofem::LIBeam2dNL::givePitch(), oofem::Truss2d::givePitch(), oofem::HTSelement::giveSideLength(), oofem::Node::giveUpdatedCoordinate(), oofem::LatticeBeam3dBoundary::giveVTKCoordinates(), oofem::TR1_2D_SUPG::initGeometry(), oofem::PolygonLine::intersects(), oofem::PolylineNonlocalBarrier::isActivated(), oofem::FRCFCMNL::isInElementProjection(), oofem::AbaqusUserElement::postInitialize(), oofem::CohesiveSurface3d::postInitialize(), oofem::PFEM::printDofOutputAt(), oofem::LatticeDirichletCouplingNode::printYourself(), oofem::LatticeNeumannCouplingNode::printYourself(), oofem::Node::printYourself(), oofem::PLHoopStressCirc::propagateInterface(), oofem::PLPrincipalStrain::propagateInterface(), oofem::Lattice2dBoundary::recalculateCoordinates(), oofem::Lattice3dBoundary::recalculateCoordinates(), oofem::Lattice3dBoundaryTruss::recalculateCoordinates(), oofem::LatticeLink3dBoundary::recalculateCoordinates(), oofem::QClinearStatic::transformMeshToParticles(), oofem::GeometryBasedEI::updateNodeEnrMarker(), and oofem::XfemElementInterface::XfemElementInterface_createEnrNmatrixAt().

◆ giveCoordinates()

const FloatArray & oofem::DofManager::giveCoordinates ( ) const
inline
Returns
Pointer to node coordinate array.

Definition at line 390 of file dofmanager.h.

Referenced by oofem::PrimaryField::__evaluateAt(), oofem::PrimaryField::applyInitialCondition(), oofem::OctreeSpatialLocalizer::buildOctreeDataStructure(), oofem::ContactElement::computeAABB(), oofem::XfemElementInterface::ComputeBOrBHMatrix(), oofem::Tet1_3D_SUPG::computeCriticalTimeStep(), oofem::MixedGradientPressureDirichlet::computeDofTransformation(), oofem::PrescribedGradientBCPeriodic::computeDofTransformation(), oofem::SolutionbasedShapeFunction::computeDofTransformation(), oofem::TransportGradientPeriodic::computeDofTransformation(), oofem::PrescribedGradientBCWeak::computeDomainBoundingBox(), oofem::Circle::computeIntersectionPoints(), oofem::GeometryBasedEI::computeIntersectionPoints(), oofem::GeometryBasedEI::computeIntersectionPoints(), oofem::Line::computeIntersectionPoints(), oofem::PolygonLine::computeIntersectionPoints(), oofem::Quasicontinuum::computeIntersectionsOfLinkWith2DTringleElements(), oofem::Quasicontinuum::computeIntersectionsOfLinkWith3DTetrahedraElements(), oofem::LatticeNeumannCouplingNode::computeLoadCouplingContribution(), oofem::MITC4Shell::computeLoadLEToLRotationMatrix(), oofem::RigidArmNode::computeMasterContribution(), oofem::MaterialForceEvaluator::computeMaterialForce(), oofem::Quad1Mindlin::computeMidPlaneNormal(), oofem::Quad1MindlinShell3D::computeMidPlaneNormal(), oofem::Line::computeNumberOfIntersectionPoints(), oofem::PrescribedDispSlipBCDirichletRC::computeReinfStress(), oofem::MeshQualityErrorEstimator::computeTriangleRadiusError(), oofem::TransportGradientDirichlet::computeXi(), oofem::DGProblem::constructBoundaryEntities(), oofem::SolutionbasedShapeFunction::copyDofManagersToSurfaceData(), oofem::PFEM::deactivateTooCloseParticles(), oofem::LEPlic::doLagrangianPhase(), oofem::VTKExportModule::doOutput(), oofem::VTKPFEMXMLExportModule::doOutput(), oofem::ClosestNode::evaluate(), oofem::InternalVariableField::evaluateAt(), oofem::UniformGridField::evaluateAt(), oofem::UnstructuredGridField::evaluateAt(), oofem::VoxelVOFField::evaluateAt(), oofem::Shell7BaseXFEM::EvaluateEnrFuncInDofMan(), oofem::GeometryBasedEI::evaluateEnrFuncInNode(), oofem::VTKXMLLatticeExportModule::exportPrimaryVars(), oofem::VTKXMLPeriodicExportModule::exportPrimaryVars(), oofem::PrescribedGradientBCWeak::findCrackBndIntersecCoord(), oofem::PrescribedGradientBCWeak::findHoleCoord(), oofem::QClinearStatic::findNearestParticle(), oofem::DelaunayTriangulator::findNonDelaunayTriangles(), oofem::PrescribedGradientBCPeriodic::findSlaveToMasterMap(), oofem::TransportGradientPeriodic::findSlaveToMasterMap(), oofem::VTKXMLXFemExportModule::getNodalVariableFromXFEMST(), oofem::PrescribedDispSlipBCDirichletRC::give(), oofem::PrescribedGenStrainShell7::give(), oofem::PrescribedGradient::give(), oofem::RotatingBoundary::give(), oofem::TransportGradientDirichlet::give(), oofem::UserDefDirichletBC::give(), oofem::DummySpatialLocalizer::giveAllNodesWithinBox(), oofem::PlaneStress2dXfem::giveCompositeExportData(), oofem::QTrPlaneStress2dXFEM::giveCompositeExportData(), oofem::TrPlaneStress2dXFEM::giveCompositeExportData(), oofem::MITC4Shell::giveDirectorVectors(), oofem::DelaunayTriangle::giveEdgeLength(), oofem::MacroLSpace::giveInternalForcesVector(), oofem::Beam3d::giveLocalCoordinateSystem(), oofem::MicroMaterial::giveMacroStiffnessMatrix(), oofem::PrescribedGradientBCPeriodic::giveMasterDof(), oofem::DummySpatialLocalizer::giveNodeClosestToPoint(), oofem::OctreeSpatialLocalizer::giveNodeClosestToPointWithinOctant(), oofem::CCTPlate::giveNodeCoordinates(), oofem::DKTPlate::giveNodeCoordinates(), oofem::QDKTPlate::giveNodeCoordinates(), oofem::Shell7BaseXFEM::giveShellExportData(), oofem::XfemStructuralElementInterface::giveSubtriangulationCompositeExportData(), oofem::MITC4Shell::giveThickness(), oofem::MixedGradientPressureDirichlet::giveUnknown(), oofem::PrescribedGradientBCPeriodic::giveUnknown(), oofem::TransportGradientPeriodic::giveUnknown(), oofem::Truss1d::HuertaErrorEstimatorI_setupRefinedElementProblem(), oofem::SolutionbasedShapeFunction::init(), oofem::OctreeSpatialLocalizer::initElementIPDataStructure(), oofem::SolutionbasedShapeFunction::initializeSurfaceData(), oofem::OctreeSpatialLocalizer::insertNodeIntoOctree(), oofem::Circle::intersects(), oofem::Circle::isInside(), oofem::QCFullsolveddomain::isNodeInside(), oofem::SolutionbasedShapeFunction::loadProblem(), oofem::EIPrimaryUnknownMapper::mapAndUpdate(), oofem::EnrFrontLinearBranchFuncRadius::MarkNodesAsFront(), oofem::QClinearStatic::nodeInFullSolvedDomainTest(), oofem::GnuplotExportModule::outputBoundaryCondition(), oofem::IntElLine1::postInitialize(), oofem::PLCZdamageRadius::propagateInterface(), oofem::PLHoopStressCirc::propagateInterface(), oofem::PLnodeRadius::propagateInterface(), oofem::PLPrincipalStrain::propagateInterface(), oofem::MicroMaterial::setMacroProperties(), oofem::QuasicontinuumVTKXMLExportModule::setupVTKPiece(), oofem::VTKBaseExportModule::setupVTKPiece(), oofem::VTKXMLLatticeExportModule::setupVTKPiece(), oofem::VTKXMLPeriodicExportModule::setupVTKPiece(), oofem::FreemInterface::smoothNodalDensities(), oofem::SolutionbasedShapeFunction::splitBoundaryNodeIDs(), oofem::Quasicontinuum::stiffnessAssignment(), oofem::QClinearStatic::transformMeshToParticles(), oofem::PrescribedGenStrainShell7::updateCoefficientMatrix(), oofem::GeometryBasedEI::updateLevelSets(), oofem::GeometryBasedEI::updateNodeEnrMarker(), and oofem::XfemElementInterface::XfemElementInterface_createEnrNmatrixAt().

◆ giveDofTypeMap()

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 409 of file dofmanager.h.

References dofTypemap.

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

◆ giveDofWithID()

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 127 of file dofmanager.C.

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

Referenced by oofem::PrimaryField::applyBoundaryCondition(), oofem::InternalElementDofManErrorCheckingRule::check(), oofem::NodeErrorCheckingRule::check(), oofem::LatticeNeumannCouplingNode::computeLoadCouplingContribution(), oofem::CoupledFieldsElement::computeLocationArrayOfDofIDs(), oofem::IntElLine1PF::computeLocationArrayOfDofIDs(), oofem::PhaseFieldElement::computeLocationArrayOfDofIDs(), computeM2LTransformation(), oofem::RigidArmNode::computeMasterContribution(), oofem::PrescribedDispSlipBCDirichletRC::computeReinfStress(), oofem::PrescribedDispSlipBCDirichletRC::computeTransferStress(), oofem::InteractionLoad::computeValueAt(), oofem::CoupledFieldsElement::computeVectorOfDofIDs(), oofem::Domain::createDofs(), oofem::CohesiveSurface3d::drawDeformedGeometry(), oofem::CohesiveSurface3d::drawScalar(), oofem::tet21ghostsolid::EIPrimaryUnknownMI_computePrimaryUnknownVectorAtLocal(), oofem::VTKBaseExportModule::exportExternalForces(), oofem::VTKXMLLatticeExportModule::exportPrimaryVars(), oofem::VTKXMLPeriodicExportModule::exportPrimaryVars(), oofem::VTKExportModule::getDofManPrimaryVariable(), oofem::VTKBaseExportModule::getNodalVariableFromPrimaryField(), oofem::InternalElementDofManErrorCheckingRule::getValue(), oofem::NodeErrorCheckingRule::getValue(), oofem::RefinedElement::giveCompatibleBcDofArray(), oofem::MicroMaterial::giveMacroStiffnessMatrix(), givePrescribedUnknownVector(), oofem::InteractionPFEMParticle::givePrescribedUnknownVector(), oofem::tet21ghostsolid::giveRowTransformationMatrix(), oofem::PrescribedGradientBCPeriodic::giveUnknown(), oofem::PrescribedGradientBCPeriodic::giveUnknown(), oofem::TransportGradientPeriodic::giveUnknown(), oofem::TransportGradientPeriodic::giveUnknown(), oofem::tet21ghostsolid::NodalAveragingRecoveryMI_computeNodalValue(), oofem::GnuplotExportModule::outputBoundaryCondition(), oofem::GnuplotExportModule::outputReactionForces(), oofem::TransientTransportProblem::printOutputAt(), oofem::StructuralEngngModel::printReactionForces(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem1D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem2D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem3D(), oofem::HuertaErrorEstimator::solveRefinedElementProblem(), oofem::HuertaErrorEstimator::solveRefinedWholeProblem(), oofem::FluidStructureProblem::solveYourselfAt(), oofem::IncrementalLinearStatic::solveYourselfAt(), oofem::PrescribedGenStrainShell7::updateCoefficientMatrix(), oofem::PFEM::updateDofUnknownsDictionaryPressure(), and oofem::SolutionbasedShapeFunction::updateModelWithFactors().

◆ giveForcedDofIDs()

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 404 of file dofmanager.h.

References dofidmask.

◆ giveGlobalNumber()

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

Definition at line 515 of file dofmanager.h.

References globalNumber.

Referenced by oofem::Element::addDofManager(), oofem::Domain::BuildDofManPlaceInArrayMap(), oofem::ParmetisLoadBalancer::calculateLoadTransfer(), oofem::XfemElementInterface::ComputeBOrBHMatrix(), oofem::Shell7BaseXFEM::computeEnrichedBmatrixAt(), oofem::GeometryBasedEI::computeIntersectionPoints(), oofem::XfemElementInterface::computeNCohesive(), oofem::PrescribedGradientBCWeak::createTractionMesh(), oofem::LoadBalancer::deleteRemoteDofManagers(), oofem::Node::drawYourself(), oofem::Shell7BaseXFEM::EvaluateEnrFuncInDofMan(), oofem::Subdivision::exchangeRemoteElements(), oofem::PlaneStress2dXfem::giveCompositeExportData(), oofem::QTrPlaneStress2dXFEM::giveCompositeExportData(), oofem::TrPlaneStress2dXFEM::giveCompositeExportData(), oofem::AuxVelocityNumberingScheme::giveDofEquationNumber(), oofem::PlaneStress2dXfem::giveDofManDofIDMask(), oofem::QTrPlaneStress2dXFEM::giveDofManDofIDMask(), oofem::TrPlaneStress2dXFEM::giveDofManDofIDMask(), oofem::MicroMaterial::giveMacroStiffnessMatrix(), oofem::EnrichmentItem::giveNumDofManEnrichments(), oofem::PrescribedDispSlipBCDirichletRC::giveOnSteel(), oofem::XfemStructuralElementInterface::giveSubtriangulationCompositeExportData(), oofem::Natural2GlobalOrdering::init(), oofem::EnrichmentItem::isDofManEnriched(), oofem::QCFullsolveddomain::isNodeInside(), oofem::LoadBalancer::migrateLoad(), oofem::LoadBalancer::packMigratingData(), oofem::NonlocalMaterialWTP::packRemoteElements(), oofem::Subdivision::packRemoteElements(), oofem::ParmetisLoadBalancer::packSharedDmanPartitions(), oofem::LatticeDirichletCouplingNode::printOutputAt(), oofem::qcNode::printOutputAt(), oofem::Domain::py_setDofManager(), oofem::REGISTER_EnrichmentFront(), oofem::MicroMaterial::setMacroProperties(), oofem::QClinearStatic::setRepNodesInVerticesOfInterpolationMesh(), oofem::FETICommunicator::setUpCommunicationMaps(), oofem::GeometryBasedEI::updateNodeEnrMarker(), oofem::XfemElementInterface::XfemElementInterface_createEnrNmatrixAt(), and oofem::XfemElementInterface::XfemElementInterface_giveNumDofManEnrichments().

◆ giveIcMap()

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 427 of file dofmanager.h.

References dofICmap.

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

◆ giveInputRecord()

◆ giveInputRecordName()

◆ giveLabel()

◆ giveLoadArray()

◆ giveLocationArray()

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 185 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::PrescribedGradientBCWeak::assemble(), oofem::PrescribedGradientBCWeak::assembleExtraDisplock(), oofem::EngngModel::assembleVectorFromBC(), oofem::EngngModel::assembleVectorFromDofManagers(), oofem::Integral::getElementTermCodeNumbers(), and oofem::PrescribedGradientBCWeak::giveTractionLocationArray().

◆ giveMasterDofIDArray()

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 217 of file dofmanager.C.

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

Referenced by oofem::Integral::getElementTermCodeNumbers().

◆ giveMasterDofMans()

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 814 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().

◆ giveMasterMap()

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 415 of file dofmanager.h.

References dofMastermap.

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

◆ giveNumberOfDofs()

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

Definition at line 287 of file dofmanager.C.

References dofArray.

Referenced by oofem::PetscSparseMtrx::buildInternalStructure(), oofem::RigidArmNode::checkConsistency(), oofem::Node::computeL2GTransformation(), computeLoadVector(), oofem::CoupledFieldsElement::computeLocationArrayOfDofIDs(), oofem::IntElLine1PF::computeLocationArrayOfDofIDs(), oofem::PhaseFieldElement::computeLocationArrayOfDofIDs(), oofem::StructuralInterfaceElementPhF::computeLocationArrayOfDofIDs(), computeM2LTransformation(), oofem::ElementSide::computeTransformation(), oofem::HuertaErrorEstimator::estimateError(), oofem::VTKXMLLatticeExportModule::exportPrimaryVars(), oofem::VTKXMLPeriodicExportModule::exportPrimaryVars(), oofem::RefinedElement::giveBcDofArray1D(), oofem::RefinedElement::giveBcDofArray2D(), oofem::RefinedElement::giveBcDofArray3D(), giveCompleteLocationArray(), giveCompleteMasterDofIDArray(), giveCompleteUnknownVector(), oofem::MacroLSpace::giveInternalForcesVector(), oofem::BaseMixedPressureElement::giveLocationArrayOfDofIDs(), oofem::GradientDamageElement::giveLocationArrayOfDofIDs(), oofem::Natural2GlobalOrdering::init(), oofem::LSPrimaryVariableMapper::mapPrimaryVariables(), oofem::LatticeDirichletCouplingNode::printOutputAt(), oofem::LatticeDirichletCouplingNode::printYourself(), oofem::LatticeNeumannCouplingNode::printYourself(), oofem::MicroMaterial::setMacroProperties(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem1D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem2D(), and oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem3D().

◆ giveNumberOfPrimaryMasterDofs()

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 311 of file dofmanager.C.

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

Referenced by computeM2LTransformation().

◆ giveParallelMode()

◆ givePartitionList()

◆ givePartitionsConnectivitySize()

int oofem::DofManager::givePartitionsConnectivitySize ( )

Returns number of partitions sharing given receiver (=number of shared partitions + local one).

Definition at line 934 of file dofmanager.C.

References oofem::FEMComponent::giveDomain(), and partitions.

Referenced by oofem::EngngModel::assembleVectorFromDofManagers().

◆ givePrescribedUnknownVector()

void oofem::DofManager::givePrescribedUnknownVector ( FloatArray & answer,
const IntArray & dofMask,
ValueModeType mode,
TimeStep * tStep )
virtual

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.

Reimplemented in oofem::InteractionPFEMParticle.

Definition at line 720 of file dofmanager.C.

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

◆ giveUnknownVector() [1/2]

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 683 of file dofmanager.C.

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

◆ giveUnknownVector() [2/2]

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 657 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::computeIntForceGPContrib(), oofem::LEPlic::doLagrangianPhase(), oofem::MPElement::getBoundaryUnknownVector(), oofem::MPElement::getUnknownVector(), oofem::HuertaErrorEstimator::solveRefinedElementProblem(), and oofem::HuertaErrorEstimator::solveRefinedWholeProblem().

◆ giveUnknownVectorOfType()

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 738 of file dofmanager.C.

References oofem::FloatArray::at(), oofem::IntArray::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().

◆ hasAnySlaveDofs()

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

Definition at line 802 of file dofmanager.C.

References oofem::Dof::isPrimaryDof().

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

◆ hasDofID()

◆ initializeFinish()

void oofem::DofManager::initializeFinish ( )
overridevirtual

Finishes the initialization. Note that initializeFrom may be called multiple times. The initializeFinish typycally performs the input parameter checking (if compulsory parameters set, etc.) After initializeFinish, DOFs and other components may be created.

Reimplemented from oofem::FEMComponent.

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

Definition at line 367 of file dofmanager.C.

References oofem::IntArray::at(), oofem::ParameterManager::checkIfSet(), oofem::IntArray::contains(), dofBCmap, dofICmap, dofidmask, dofMastermap, dofTypemap, oofem::FEMComponent::domain, oofem::ParameterManager::getTempParam(), oofem::FEMComponent::giveDomain(), oofem::IntArray::giveSize(), IPK_DofManager_bc, IPK_DofManager_dofidmask, IPK_DofManager_doftypemask, IPK_DofManager_ic, IPK_DofManager_mastermask, mBC, oofem::FEMComponent::number, and OOFEM_ERROR.

◆ initializeFrom() [1/2]

void oofem::DofManager::initializeFrom ( InputRecord & ir)
inlineoverridevirtual

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. Note that initializeFrom may be called mutiple times.

See also
IR_GIVE_FIELD
IR_GIVE_OPTIONAL_FIELD
Parameters
irInput record to initialize from.
priorityPriority of the input record. This is used to determine the order of initialization

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::Node.

Definition at line 456 of file dofmanager.h.

References initializeFrom().

Referenced by oofem::Domain::initializeFinish(), and initializeFrom().

◆ initializeFrom() [2/2]

◆ isBoundary()

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

Definition at line 439 of file dofmanager.h.

References isBoundaryFlag.

Referenced by setBoundaryFlag().

◆ isDofTypeCompatible()

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::ClonedDofManager, oofem::ElementDofManager, oofem::ElementSide, oofem::GeneralSlaveNode, oofem::HangingNode, oofem::Node, oofem::qcNode, oofem::RigidArmNode, and oofem::SlaveNode.

Definition at line 468 of file dofmanager.h.

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

◆ isLocal()

bool oofem::DofManager::isLocal ( )

Returns true if receiver is locally maintained.

Definition at line 946 of file dofmanager.C.

References oofem::DofManager_local, oofem::DofManager_shared, oofem::FEMComponent::giveDomain(), oofem::min(), parallel_mode, and partitions.

◆ isNull()

bool oofem::DofManager::isNull ( )
inline

Returns true if receiver is shared.

Definition at line 554 of file dofmanager.h.

References oofem::DofManager_null, and parallel_mode.

◆ isShared()

bool oofem::DofManager::isShared ( )
inline

◆ mergePartitionList()

void oofem::DofManager::mergePartitionList ( IntArray & _p)

Merges receiver partition list with given lists.

Definition at line 924 of file dofmanager.C.

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

◆ postInitialize()

void oofem::DofManager::postInitialize ( )
overridevirtual

◆ printOutputAt()

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

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::LatticeDirichletCouplingNode, and oofem::qcNode.

Definition at line 507 of file dofmanager.C.

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

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

◆ printYourself()

void oofem::DofManager::printYourself ( )
overridevirtual

Prints receiver state on stdout. Useful for debugging.

Reimplemented from oofem::FEMComponent.

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

Definition at line 518 of file dofmanager.C.

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

◆ removeDof()

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

Removes Dof with given id from dofArray.

Parameters
id

Definition at line 159 of file dofmanager.C.

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

Referenced by oofem::EnrichmentItem::createEnrichedDofs().

◆ removePartitionFromList()

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

Removes given partition from receiver list.

Definition at line 537 of file dofmanager.h.

References partitions.

◆ requiresTransformation()

bool oofem::DofManager::requiresTransformation ( )
virtual

Indicates, whether dofManager requires the transformation.

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

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

Definition at line 903 of file dofmanager.C.

References hasAnySlaveDofs().

◆ restoreContext()

void oofem::DofManager::restoreContext ( DataStream & stream,
ContextMode mode )
overridevirtual

Restores the receiver state previously written in stream.

See also
saveContext
Parameters
streamInput stream.
modeDetermines amount of info available in stream (state, definition, ...).
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 595 of file dofmanager.C.

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

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

◆ saveContext()

void oofem::DofManager::saveContext ( DataStream & stream,
ContextMode mode )
overridevirtual

Stores receiver state to output stream.

Parameters
streamOutput stream.
modeDetermines amount of info required in stream (state, definition, ...).
Exceptions
throwsan ContextIOERR exception if error encountered.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::Node.

Definition at line 540 of file dofmanager.C.

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

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

◆ setBoundaryFlag()

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

Sets the boundary flag.

Parameters
isBoundaryDetermines if receiver is on the boundary.

Definition at line 444 of file dofmanager.h.

References isBoundary().

◆ setCoordinates()

void oofem::DofManager::setCoordinates ( const FloatArray & coords)
inline

Set coordinates.

Definition at line 394 of file dofmanager.h.

◆ setGlobalNumber()

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

Sets receiver global number.

Parameters
numberNew global number for receiver.

Definition at line 521 of file dofmanager.h.

References globalNumber.

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

◆ setLoadArray()

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 99 of file dofmanager.C.

References loadArray.

◆ setNumberOfDofs()

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

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

Definition at line 294 of file dofmanager.C.

References dofArray.

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

◆ setParallelMode()

void oofem::DofManager::setParallelMode ( dofManagerParallelMode _mode)
inline

◆ setPartitionList()

void oofem::DofManager::setPartitionList ( const IntArray * _p)
inline

Sets receiver's partition list.

Definition at line 535 of file dofmanager.h.

References partitions.

Referenced by oofem::LoadBalancer::deleteRemoteDofManagers(), and oofem::LoadBalancer::unpackMigratingData().

◆ updateLocalNumbering()

void oofem::DofManager::updateLocalNumbering ( EntityRenumberingFunctor & f)
overridevirtual

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::GeneralSlaveNode, oofem::RigidArmNode, and oofem::SlaveNode.

Definition at line 909 of file dofmanager.C.

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

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

◆ updateYourself()

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

Updates receiver after equilibrium in time step has been reached.

Parameters
tStepActive time step.

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

Definition at line 531 of file dofmanager.C.

References oofem::Dof::updateYourself().

Member Data Documentation

◆ coordinates

◆ dofArray

◆ dofBCmap

std :: map< int, int > oofem::DofManager::dofBCmap
protected

Map from DofIDItem to bc (to be removed).

Definition at line 138 of file dofmanager.h.

Referenced by DofManager(), giveBcMap(), and initializeFinish().

◆ dofICmap

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

Map from DofIDItem to ic (to be removed).

Definition at line 140 of file dofmanager.h.

Referenced by DofManager(), giveIcMap(), and initializeFinish().

◆ dofidmask

IntArray oofem::DofManager::dofidmask
protected

List of additional dof ids to include.

Definition at line 132 of file dofmanager.h.

Referenced by oofem::RigidArmNode::computeMasterContribution(), giveForcedDofIDs(), giveInputRecord(), initializeFinish(), and initializeFrom().

◆ dofMastermap

std :: map< int, int > oofem::DofManager::dofMastermap
protected

Map from DofIDItem to master node.

Definition at line 136 of file dofmanager.h.

Referenced by DofManager(), giveInputRecord(), giveMasterMap(), initializeFinish(), and updateLocalNumbering().

◆ dofTypemap

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

◆ globalNumber

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 121 of file dofmanager.h.

Referenced by DofManager(), giveGlobalNumber(), giveLabel(), restoreContext(), saveContext(), and setGlobalNumber().

◆ hasSlaveDofs

bool oofem::DofManager::hasSlaveDofs
protected

◆ IPK_DofManager_bc

◆ IPK_DofManager_boundaryflag

ParamKey oofem::DofManager::IPK_DofManager_boundaryflag
static

Definition at line 152 of file dofmanager.h.

Referenced by giveInputRecord(), and initializeFrom().

◆ IPK_DofManager_dofidmask

◆ IPK_DofManager_doftypemask

ParamKey oofem::DofManager::IPK_DofManager_doftypemask
static

◆ IPK_DofManager_globnum

ParamKey oofem::DofManager::IPK_DofManager_globnum
static

Definition at line 153 of file dofmanager.h.

◆ IPK_DofManager_ic

ParamKey oofem::DofManager::IPK_DofManager_ic
static

Definition at line 149 of file dofmanager.h.

Referenced by initializeFinish(), and initializeFrom().

◆ IPK_DofManager_load

ParamKey oofem::DofManager::IPK_DofManager_load
static

◆ IPK_DofManager_mastermask

ParamKey oofem::DofManager::IPK_DofManager_mastermask
static

◆ IPK_DofManager_nullflag

ParamKey oofem::DofManager::IPK_DofManager_nullflag
static

Definition at line 157 of file dofmanager.h.

Referenced by giveInputRecord(), and initializeFrom().

◆ IPK_DofManager_partitions

ParamKey oofem::DofManager::IPK_DofManager_partitions
static

Definition at line 154 of file dofmanager.h.

Referenced by giveInputRecord(), and initializeFrom().

◆ IPK_DofManager_remoteflag

ParamKey oofem::DofManager::IPK_DofManager_remoteflag
static

Definition at line 156 of file dofmanager.h.

Referenced by giveInputRecord(), and initializeFrom().

◆ IPK_DofManager_sharedflag

ParamKey oofem::DofManager::IPK_DofManager_sharedflag
static

Definition at line 155 of file dofmanager.h.

Referenced by giveInputRecord(), and initializeFrom().

◆ isBoundaryFlag

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 113 of file dofmanager.h.

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

◆ loadArray

◆ mBC

IntArray oofem::DofManager::mBC
protected

Definition at line 143 of file dofmanager.h.

Referenced by initializeFinish(), and initializeFrom().

◆ parallel_mode

◆ partitions

IntArray oofem::DofManager::partitions
protected

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

Definition at line 129 of file dofmanager.h.

Referenced by DofManager(), giveInputRecord(), givePartitionList(), givePartitionsConnectivitySize(), initializeFrom(), isLocal(), mergePartitionList(), removePartitionFromList(), restoreContext(), saveContext(), and setPartitionList().


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

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011