OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Abstract class Dof represents Degree Of Freedom in finite element mesh. More...
#include <dof.h>
Public Member Functions | |
Dof (DofManager *aNode, DofIDItem id=Undef) | |
Constructor. More... | |
virtual | ~Dof () |
Destructor. More... | |
virtual dofType | giveDofType ()=0 |
Returns the type of the receiver. More... | |
virtual const char * | giveClassName () const |
Returns class name of the receiver. More... | |
int | giveDofManNumber () const |
DofManager * | giveDofManager () const |
int | giveDofManGlobalNumber () const |
virtual double | giveBcValue (ValueModeType mode, TimeStep *tStep) |
Returns value of boundary condition of dof if it is prescribed. More... | |
int | giveEquationNumber (const UnknownNumberingScheme &s) |
Returns equation number of receiver for given equation numbering scheme. More... | |
virtual int | __giveEquationNumber () const =0 |
Returns equation number of receiver, usually assigned by emodel. More... | |
virtual void | giveEquationNumbers (IntArray &masterEqNumbers, const UnknownNumberingScheme &s) |
Returns equation number of receiver. More... | |
virtual void | giveDofIDs (IntArray &masterDofIDs) |
As giveEquationNumbers but for dof IDs. More... | |
virtual int | __givePrescribedEquationNumber ()=0 |
Returns prescribed equation number of receiver. More... | |
virtual int | askNewEquationNumber (TimeStep *tStep)=0 |
Asks EngngModel for new equation number. More... | |
virtual double | giveUnknown (ValueModeType mode, TimeStep *tStep)=0 |
The key method of class Dof. More... | |
virtual double | giveUnknown (PrimaryField &field, ValueModeType mode, TimeStep *tStep)=0 |
The key method of class Dof. More... | |
virtual void | giveUnknowns (FloatArray &masterUnknowns, ValueModeType mode, TimeStep *tStep) |
The key method of class Dof. More... | |
virtual void | giveUnknowns (FloatArray &masterUnknowns, PrimaryField &field, ValueModeType mode, TimeStep *tStep) |
The key method of class Dof. More... | |
virtual void | computeDofTransformation (FloatArray &masterContribs) |
Computes dof transformation array, which describes the dependence of receiver value on values of master dofs. More... | |
virtual int | giveNumberOfPrimaryMasterDofs () |
virtual bool | hasBc (TimeStep *tStep)=0 |
Test if Dof has active boundary condition. More... | |
virtual bool | hasIc ()=0 |
Test if Dof has initial condition. More... | |
virtual bool | hasIcOn (ValueModeType u)=0 |
Test if Dof has initial condition of required ValueModeType. More... | |
DofIDItem | giveDofID () const |
Returns DofID value of receiver, which determines type of of unknown connected to receiver (e.g., u-displacement, v-displacement, ...). More... | |
void | setDofID (DofIDItem id) |
Sets the ID of receiver. More... | |
virtual bool | isPrimaryDof () |
Tests if receiver is primary DOF. More... | |
virtual int | giveBcId ()=0 |
Returns the id of associated boundary condition, if there is any. More... | |
virtual int | giveIcId ()=0 |
Returns the id of associated initial condition, if there is any. More... | |
virtual void | giveMasterDofManArray (IntArray &answer) |
virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
Local renumbering support. More... | |
virtual void | printSingleOutputAt (FILE *file, TimeStep *tStep, char ch, ValueModeType mode, double scale=1.0) |
Prints Dof output (it prints value of unknown related to dof at given timeStep). More... | |
virtual void | printMultipleOutputAt (FILE *File, TimeStep *tStep, char *ch, ValueModeType *mode, int nite) |
Prints Dof output (it prints value of unknown related to dof at given timeStep). More... | |
void | printSingleOutputWithAdditionAt (FILE *File, TimeStep *tStep, char ch, ValueModeType mode, double addend) |
virtual void | printYourself () |
Prints the receiver state on stdout. More... | |
virtual void | updateYourself (TimeStep *tStep) |
Updates receiver after finishing time step. More... | |
virtual void | updateUnknownsDictionary (TimeStep *tStep, ValueModeType mode, double dofValue) |
Abstract function, allowing Dof to store its unknowns in its own private dictionary. More... | |
virtual double | giveUnknownsDictionaryValue (TimeStep *tStep, ValueModeType mode) |
Access dictionary value, if not present zero is returned. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). 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 void | setBcId (int bcId) |
Overwrites the boundary condition id (0-inactive BC), intended for specific purposes such as coupling of bc's in multiscale simulations. More... | |
virtual void | setIcId (int icId) |
Overwrites the initial condition id (0-inactive IC) More... | |
virtual void | setEquationNumber (int equationNumber) |
Sets a specific equation number to receiver. More... | |
virtual void | setUnknowns (Dictionary *unknowns) |
Sets the dictionary of unknowns for receiver. More... | |
virtual Dictionary * | giveUnknowns () |
Receives the dictionary of unknowns in receiver. More... | |
virtual int | giveEqn () |
Gives number for equation, negative for prescribed equations. More... | |
virtual int | packUnknowns (DataStream &buff, ValueModeType mode, TimeStep *tStep) |
Packs specific DOF Manager's dofs unknowns into communication buffer. More... | |
virtual int | unpackAndUpdateUnknown (DataStream &buff, ValueModeType mode, TimeStep *tStep) |
Unpacks DOF unknown from communication buffer and updates unknown if necessary. More... | |
Protected Member Functions | |
virtual BoundaryCondition * | giveBc () |
Returns boundary condition of dof if it is prescribed. More... | |
virtual InitialCondition * | giveIc () |
Returns initial condition of dof if it is prescribed. More... | |
Protected Attributes | |
DofManager * | dofManager |
Link to related DofManager. More... | |
DofIDItem | dofID |
Physical meaning of DOF. More... | |
Friends | |
class | SimpleSlaveDof |
Abstract class Dof represents Degree Of Freedom in finite element mesh.
DOFs are possessed by DofManagers (i.e, nodes, sides or whatever) and one DOF belongs to only one DofManager. Dof maintain its related equation or prescribed equation number. This equation number is usually assigned by EngngModel, however, several numbering schemes can exists (see giveEquationNumber and similar services).
It maintains also its physical meaning and reference to related DofManager (reference to DofManager which possess particular DOF). To describe physical meaning of particular Dof, special enum type DofId has been introduced (see cltypes.h). This type is more descriptive than UnknownType, which determines physical meaning for unknowns generally (displacement or temperature). DofId type must distinguish between dofs having displacement unknown, but in different directions, because only some of these may be required by elements.
Dof can be subjected to boundary (BC) or initial (IC) condition. Method for obtaining corresponding DOF unknown value is provided. If no IC condition has been given, zero value IC is assumed otherwise when needed.
Dof class generally supports changes of static system during computation. This feature generally leads to equation renumbering. Then because equation number associated to dof may change, it may become extremely complicated to ask EngngModel for unknown from previous time step (because equation number may have been changed). To overcome this problem, derived class will implement so called unknown dictionary, which is updated after finishing each time step and where unknowns for particular dof are stored. Dof then uses this dictionary for requests for unknowns instead of asking EngngModel for unknowns. Unknowns in dof dictionary are updated by EngngModel automatically (if EngngModel supports changes of static system) after finishing time step.
oofem::Dof::Dof | ( | DofManager * | aNode, |
DofIDItem | id = Undef |
||
) |
Constructor.
Creates DOF with number i, belonging to DofManager aNode and with physical meaning described by id.
aNode | DofManager which possess DOF. |
id | Physical meaning type. |
Definition at line 50 of file dof.C.
References dofID, and dofManager.
|
pure virtual |
Returns equation number of receiver, usually assigned by emodel.
If Dof has active BC, returned equation number is zero. After initializing Dof by calling constructor, Dof has no equation number assigned. When firstly invoked, this function asks EngngModel object for next equation prescribed equation number (this will increase also total number of equation at EngngModel level). Note: By asking nodal code numbers or element code numbers when initializing code numbers in EngngMode, designer should alter equation numbering strategy.
Implemented in oofem::MasterDof, oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Referenced by oofem::SimpleSlaveDof::__giveEquationNumber(), oofem::EModelDefaultEquationNumbering::giveDofEquationNumber(), oofem::VelocityNumberingScheme::giveDofEquationNumber(), oofem::DofIDEquationNumbering::giveDofEquationNumber(), oofem::VelocityEquationNumbering::giveDofEquationNumber(), oofem::PressureEquationNumbering::giveDofEquationNumber(), oofem::NLTransientTransportProblem::giveUnknownComponent(), oofem::StationaryTransportProblem::giveUnknownComponent(), oofem::LinearStatic::giveUnknownComponent(), oofem::FreeWarping::giveUnknownComponent(), oofem::DEIDynamic::giveUnknownComponent(), oofem::EigenValueDynamic::giveUnknownComponent(), oofem::DIIDynamic::giveUnknownComponent(), oofem::AdaptiveNonLinearStatic::giveUnknownComponent(), oofem::LinearStability::giveUnknownComponent(), oofem::NonLinearDynamic::giveUnknownComponent(), oofem::NlDEIDynamic::giveUnknownComponent(), oofem::NonStationaryTransportProblem::giveUnknownComponent(), oofem::NonLinearStatic::giveUnknownComponent(), oofem::SUPG::giveUnknownComponent(), oofem::CBS::giveUnknownComponent(), oofem::HuertaErrorEstimator::solveRefinedElementProblem(), and oofem::HuertaErrorEstimator::solveRefinedWholeProblem().
|
pure virtual |
Returns prescribed equation number of receiver.
If Dof has inactive BC, returned prescribed equation number is zero. If Dof has active BC, then the corresponding prescribed equation number is returned. is zero. After initializing Dof by calling constructor, Dof has no prescribed equation number assigned. When firstly invoked, this function asks EngngModel object for next equation or prescribed equation number (this will increase also total number of equation at EngngModel level). Note: By asking nodal code numbers or element code numbers when initializing code numbers in EngngMode, designer should alter equation numbering strategy.
Implemented in oofem::MasterDof, oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Referenced by oofem::SimpleSlaveDof::__givePrescribedEquationNumber(), oofem::TransportGradientDirichlet::computeCoefficientMatrix(), oofem::PrescribedGenStrainShell7::evaluateHigherOrderContribution(), oofem::EModelDefaultPrescribedEquationNumbering::giveDofEquationNumber(), oofem::VelocityNumberingScheme::giveDofEquationNumber(), oofem::DofIDEquationNumbering::giveDofEquationNumber(), oofem::VelocityEquationNumbering::giveDofEquationNumber(), oofem::PressureEquationNumbering::giveDofEquationNumber(), oofem::CBS::giveTractionPressure(), oofem::PrescribedGradient::updateCoefficientMatrix(), and oofem::PrescribedGenStrainShell7::updateCoefficientMatrix().
|
pure virtual |
Asks EngngModel for new equation number.
Necessary for EngngModels supporting changes of static system during solution. Then it is necessary to force equation renumbering after finishing each time step.
tStep | Time step determining the time. |
Implemented in oofem::MasterDof, oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Referenced by oofem::DofManager::askNewEquationNumbers().
|
virtual |
Computes dof transformation array, which describes the dependence of receiver value on values of master dofs.
For primary dof, this transformation is unity, however, for slave DOFs, this array contains weights, which are multiplied by corresponding master DOF values to obtain slave value.
masterContribs | Contributions of master dofs for receiver. |
Reimplemented in oofem::SlaveDof, and oofem::ActiveDof.
Definition at line 176 of file dof.C.
References oofem::FloatArray::at(), and oofem::FloatArray::resize().
Referenced by oofem::ActiveDof::computeDofTransformation(), oofem::SlaveDof::computeDofTransformation(), and oofem::DofManager::computeM2LTransformation().
std::string oofem::Dof::errorInfo | ( | const char * | func | ) | const |
Returns string for prepending output (used by error reporting macros).
Definition at line 113 of file dof.C.
References dofID, dofManager, giveClassName(), and oofem::FEMComponent::giveNumber().
|
inlineprotectedvirtual |
Returns boundary condition of dof if it is prescribed.
Reimplemented in oofem::MasterDof, and oofem::SimpleSlaveDof.
Definition at line 448 of file dof.h.
Referenced by oofem::SimpleSlaveDof::giveBc(), and giveBcValue().
|
pure virtual |
Returns the id of associated boundary condition, if there is any.
Used only for printing purposes. In general, id could not be used to decide whether BC is active. Use appropriate services instead.
Implemented in oofem::MasterDof, oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Referenced by oofem::SimpleSlaveDof::giveBcId(), oofem::RefinedElement::giveCompatibleBcDofArray(), oofem::DofManager::giveInputRecord(), oofem::GnuplotExportModule::outputReactionForces(), oofem::SolutionbasedShapeFunction::setBoundaryConditionOnDof(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem1D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem2D(), and oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem3D().
|
virtual |
Returns value of boundary condition of dof if it is prescribed.
Use hasBc service to determine, if boundary condition is active. The physical meaning of BC is determined by corresponding DOF.
mode | Unknown char type (if total or incremental value is returned). |
tStep | Time step. |
Reimplemented in oofem::SimpleSlaveDof, and oofem::ActiveDof.
Definition at line 120 of file dof.C.
References oofem::InitialCondition::give(), oofem::BoundaryCondition::give(), giveBc(), giveIc(), hasBc(), hasIcOn(), and oofem::TimeStep::isTheFirstStep().
Referenced by oofem::SimpleSlaveDof::giveBcValue(), oofem::RefinedElement::giveCompatibleBcDofArray(), oofem::InteractionPFEMParticle::givePrescribedUnknownVector(), oofem::DofManager::givePrescribedUnknownVector(), oofem::MasterDof::giveUnknown(), and oofem::PFEM::updateDofUnknownsDictionaryPressure().
|
inlinevirtual |
Returns class name of the receiver.
Reimplemented in oofem::SlaveDof, oofem::MasterDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Definition at line 117 of file dof.h.
Referenced by errorInfo().
|
inline |
Returns DofID value of receiver, which determines type of of unknown connected to receiver (e.g., u-displacement, v-displacement, ...).
Definition at line 276 of file dof.h.
Referenced by oofem::DofManager::appendDof(), oofem::MixedGradientPressureDirichlet::assembleVector(), oofem::SolutionbasedShapeFunction::computeDofTransformation(), oofem::MixedGradientPressureDirichlet::computeDofTransformation(), oofem::Node::computeL2GTransformation(), oofem::RigidArmNode::computeMasterContribution(), oofem::Node::drawYourself(), oofem::DofManager::findDofWithDofId(), oofem::InteractionBoundaryCondition::give(), oofem::PrescribedGradient::give(), oofem::RotatingBoundary::give(), oofem::PrescribedGenStrainShell7::give(), oofem::UserDefDirichletBC::give(), oofem::BoundaryCondition::give(), oofem::PrescribedGradientBCPeriodic::giveBcValue(), oofem::TransportGradientPeriodic::giveBcValue(), oofem::MixedGradientPressureDirichlet::giveBcValue(), oofem::RefinedElement::giveCompatibleBcDofArray(), oofem::DofManager::giveCompleteMasterDofIDArray(), oofem::PressureNumberingScheme::giveDofEquationNumber(), oofem::QuasicontinuumNumberingscheme::giveDofEquationNumber(), oofem::VelocityNumberingScheme::giveDofEquationNumber(), oofem::MicroMaterial::giveDofEquationNumber(), oofem::DofIDEquationNumbering::giveDofEquationNumber(), oofem::VelocityEquationNumbering::giveDofEquationNumber(), oofem::AuxVelocityNumberingScheme::giveDofEquationNumber(), oofem::PressureEquationNumbering::giveDofEquationNumber(), oofem::ActiveDof::giveDofIDs(), giveDofIDs(), oofem::DofManager::giveInputRecord(), oofem::PrescribedGradientBCPeriodic::giveMasterDof(), oofem::TransportGradientPeriodic::giveMasterDof(), oofem::PrescribedGradientBCPeriodic::giveUnknown(), oofem::TransportGradientPeriodic::giveUnknown(), oofem::MixedGradientPressureDirichlet::giveUnknown(), oofem::NLTransientTransportProblem::giveUnknownComponent(), oofem::NonStationaryTransportProblem::giveUnknownComponent(), oofem::CBS::giveUnknownComponent(), oofem::DofManager::giveUnknownVectorOfType(), oofem::Node::giveUpdatedCoordinate(), oofem::DofManager::hasDofID(), oofem::HangingNode::postInitialize(), oofem::RigidArmNode::postInitialize(), oofem::qcNode::postInitializeAsHangingNode(), oofem::SUPG::printDofOutputAt(), oofem::CBS::printDofOutputAt(), oofem::PFEM::printDofOutputAt(), oofem::DofManager::removeDof(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem3D(), and oofem::Node::updateYourself().
|
virtual |
As giveEquationNumbers but for dof IDs.
[out] | masterDofIDs | Dof IDs of master DOFs for receiver. |
Reimplemented in oofem::SlaveDof, and oofem::ActiveDof.
Definition at line 67 of file dof.C.
References giveDofID().
Referenced by oofem::DofManager::giveCompleteMasterDofIDArray(), oofem::ActiveDof::giveDofIDs(), and oofem::SlaveDof::giveDofIDs().
|
inline |
Definition at line 123 of file dof.h.
Referenced by oofem::PrescribedGradientBCPeriodic::computeDofTransformation(), oofem::TransportGradientPeriodic::computeDofTransformation(), oofem::SolutionbasedShapeFunction::computeDofTransformation(), oofem::MixedGradientPressureDirichlet::computeDofTransformation(), oofem::InteractionBoundaryCondition::give(), oofem::PrescribedGradient::give(), oofem::RotatingBoundary::give(), oofem::PrescribedGenStrainShell7::give(), oofem::UserDefDirichletBC::give(), oofem::TransportGradientDirichlet::give(), oofem::AuxVelocityNumberingScheme::giveDofEquationNumber(), oofem::PrescribedGradientBCPeriodic::giveMasterDof(), oofem::TransportGradientPeriodic::giveMasterDof(), oofem::PrescribedGradientBCPeriodic::giveUnknown(), oofem::TransportGradientPeriodic::giveUnknown(), oofem::MixedGradientPressureDirichlet::giveUnknown(), oofem::MixedGradientPressureDirichlet::isDevDof(), oofem::TransportGradientPeriodic::isGradDof(), oofem::PrescribedGradientBCPeriodic::isStrainDof(), oofem::PFEM::printDofOutputAt(), and oofem::SolutionbasedShapeFunction::setBoundaryConditionOnDof().
int oofem::Dof::giveDofManGlobalNumber | ( | ) | const |
Definition at line 74 of file dof.C.
References dofManager, and oofem::DofManager::giveGlobalNumber().
int oofem::Dof::giveDofManNumber | ( | ) | const |
Definition at line 72 of file dof.C.
References dofManager, and oofem::FEMComponent::giveNumber().
Referenced by oofem::PressureNumberingScheme::giveDofEquationNumber(), oofem::QuasicontinuumNumberingscheme::giveDofEquationNumber(), oofem::MicroMaterial::giveDofEquationNumber(), oofem::ActiveDof::giveMasterDofManArray(), and giveMasterDofManArray().
|
pure virtual |
Returns the type of the receiver.
Implemented in oofem::MasterDof, oofem::SimpleSlaveDof, oofem::ActiveDof, and oofem::SlaveDof.
|
inlinevirtual |
Gives number for equation, negative for prescribed equations.
Reimplemented in oofem::MasterDof.
Definition at line 407 of file dof.h.
Referenced by oofem::PrimaryField::applyBoundaryCondition(), oofem::MixedGradientPressureDirichlet::computeTangents(), and oofem::PrimaryField::giveUnknownValue().
int oofem::Dof::giveEquationNumber | ( | const UnknownNumberingScheme & | s | ) |
Returns equation number of receiver for given equation numbering scheme.
s | Numbering scheme used to obtain equation number. |
Definition at line 56 of file dof.C.
References oofem::UnknownNumberingScheme::giveDofEquationNumber().
Referenced by oofem::MixedGradientPressureDirichlet::assembleVector(), oofem::SparseNonLinearSystemNM::convertPertMap(), oofem::VTKXMLExportModule::exportExternalForces(), oofem::ActiveDof::giveEquationNumbers(), oofem::NRSolver::initPrescribedEqs(), and oofem::EIPrimaryUnknownMapper::mapAndUpdate().
|
virtual |
Returns equation number of receiver.
If Dof has active BC, returned equation number is zero. After initializing Dof by calling constructor, Dof has no equation number assigned. When firstly invoked, this function asks EngngModel object for next equation prescribed equation number (this will increase also total number of equation at EngngModel level). Note: By asking nodal code numbers or element code numbers when initializing code numbers in EngngMode, designer should alter equation numbering strategy.
For slave dofs (dependent on other primary dofs) the array of master equation numbers is returned.
[out] | masterEqNumbers | Equation numbers of master DOFs for receiver. |
[in] | s | Numbering scheme used enumeration of equations. |
Reimplemented in oofem::SlaveDof, and oofem::ActiveDof.
Definition at line 61 of file dof.C.
References oofem::IntArray::at(), oofem::UnknownNumberingScheme::giveDofEquationNumber(), and oofem::IntArray::resize().
Referenced by oofem::DofManager::giveCompleteLocationArray(), oofem::ActiveDof::giveEquationNumbers(), and oofem::SlaveDof::giveEquationNumbers().
|
inlineprotectedvirtual |
Returns initial condition of dof if it is prescribed.
Reimplemented in oofem::MasterDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Definition at line 453 of file dof.h.
Referenced by giveBcValue(), and oofem::SimpleSlaveDof::giveIc().
|
pure virtual |
Returns the id of associated initial condition, if there is any.
Used only for printing purposes. In general, id could not be used to decide whether IC is active. Use appropriate services instead.
Implemented in oofem::MasterDof, oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Referenced by oofem::SimpleSlaveDof::giveIcId().
|
virtual |
Reimplemented in oofem::SimpleSlaveDof, oofem::SlaveDof, and oofem::ActiveDof.
Definition at line 183 of file dof.C.
References oofem::IntArray::at(), giveDofManNumber(), and oofem::IntArray::resize().
Referenced by oofem::ActiveDof::giveMasterDofManArray(), oofem::SlaveDof::giveMasterDofManArray(), and oofem::DofManager::giveMasterDofMans().
|
inlinevirtual |
Reimplemented in oofem::SlaveDof, and oofem::ActiveDof.
Definition at line 252 of file dof.h.
Referenced by oofem::DofManager::computeM2LTransformation(), oofem::ActiveDof::giveNumberOfPrimaryMasterDofs(), and oofem::SlaveDof::giveNumberOfPrimaryMasterDofs().
|
pure virtual |
The key method of class Dof.
Returns the value of the unknown of the receiver at given time step. Unknown is characterized by its physical meaning (i.g., displacement) an by its mode (e.g., value of displacement, velocity of displacement or acceleration of displacement). UnknownType of requested unknown must be same as UnknownType of Dof.
mode | Mode of unknown (e.g, total value, velocity or acceleration of unknown). |
tStep | Time 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). |
Implemented in oofem::MasterDof, oofem::SimpleSlaveDof, oofem::SlaveDof, and oofem::ActiveDof.
Referenced by oofem::LinearConstraintBC::assembleVector(), oofem::MacroLSpace::changeMicroBoundaryConditions(), oofem::NodeErrorCheckingRule::check(), oofem::Node2NodeContactL::computeContactTractionAt(), oofem::Tr_Warp::computeEdgeLoadVectorAt(), oofem::MixedGradientPressureDirichlet::computeFields(), oofem::CoupledFieldsElement::computeVectorOfDofIDs(), oofem::ContactElement::ContactElement(), oofem::CohesiveSurface3d::drawDeformedGeometry(), oofem::CohesiveSurface3d::drawScalar(), oofem::Node::drawYourself(), oofem::tet21ghostsolid::EIPrimaryUnknownMI_computePrimaryUnknownVectorAtLocal(), oofem::Tet21Stokes::EIPrimaryUnknownMI_computePrimaryUnknownVectorAtLocal(), oofem::VTKXMLExportModule::getNodalVariableFromPrimaryField(), oofem::DofManager::giveCompleteUnknownVector(), oofem::PrescribedGradientBCPeriodic::giveUnknown(), oofem::SimpleSlaveDof::giveUnknown(), giveUnknowns(), oofem::DofManager::giveUnknownVectorOfType(), oofem::Node::giveUpdatedCoordinate(), oofem::Hexa21Stokes::NodalAveragingRecoveryMI_computeNodalValue(), oofem::tet21ghostsolid::NodalAveragingRecoveryMI_computeNodalValue(), oofem::Tet21Stokes::NodalAveragingRecoveryMI_computeNodalValue(), oofem::Tr21Stokes::NodalAveragingRecoveryMI_computeNodalValue(), oofem::Shell7Base::NodalAveragingRecoveryMI_computeNodalValue(), oofem::GnuplotExportModule::outputBoundaryCondition(), printMultipleOutputAt(), printSingleOutputAt(), printSingleOutputWithAdditionAt(), and oofem::Node::updateYourself().
|
pure virtual |
The key method of class Dof.
Returns the value of the unknown of the receiver at given time step associated to given field.
field | Field used to provide values. |
mode | Mode of unknown. |
tStep | Time 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). |
Implemented in oofem::MasterDof, oofem::SimpleSlaveDof, oofem::SlaveDof, and oofem::ActiveDof.
|
virtual |
The key method of class Dof.
Returns the value of the unknown of the receiver at given time step associated to given field. For primary dof it returns is associated unknown value, for slave dofs it returns an array of master values (in recursive way).
masterUnknowns | Values of master unknowns for receiver. |
mode | Value mode for unknowns. |
tStep | Time step for when unknowns are requested. |
Reimplemented in oofem::SlaveDof, and oofem::ActiveDof.
Definition at line 162 of file dof.C.
References oofem::FloatArray::at(), giveUnknown(), and oofem::FloatArray::resize().
Referenced by oofem::NLTransientTransportProblem::giveUnknownComponent(), oofem::IncrementalLinearStatic::giveUnknownComponent(), oofem::NonStationaryTransportProblem::giveUnknownComponent(), oofem::SUPG::giveUnknownComponent(), oofem::PFEM::giveUnknownComponent(), oofem::ActiveDof::giveUnknowns(), and oofem::SlaveDof::giveUnknowns().
|
virtual |
The key method of class Dof.
Returns the value of the unknown of the receiver at given time step associated to given field. For primary dof it returns is associated unknown value, for slave dofs it returns an array of master values (in recursive way).
masterUnknowns | |
field | The field to pick unknowns from. |
mode | Value mode for unknowns. |
tStep | Time step for when unknowns are requested. |
Reimplemented in oofem::SlaveDof, and oofem::ActiveDof.
Definition at line 169 of file dof.C.
References oofem::FloatArray::at(), giveUnknown(), and oofem::FloatArray::resize().
|
inlinevirtual |
Receives the dictionary of unknowns in receiver.
Reimplemented in oofem::MasterDof.
Definition at line 401 of file dof.h.
Referenced by oofem::SlaveDof::giveUnknown().
|
inlinevirtual |
Access dictionary value, if not present zero is returned.
tStep | Time step. |
mode | Mode of value. |
Reimplemented in oofem::MasterDof.
Definition at line 373 of file dof.h.
References oofem::errorInfo().
Referenced by oofem::StaticStructural::giveUnknownComponent(), and oofem::DofDistributedPrimaryField::giveUnknownValue().
|
pure virtual |
Test if Dof has active boundary condition.
tStep | Time when test is evaluated. |
Implemented in oofem::MasterDof, oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Referenced by oofem::Node::drawYourself(), giveBcValue(), oofem::RefinedElement::giveCompatibleBcDofArray(), oofem::SimpleSlaveDof::hasBc(), oofem::HuertaErrorEstimator::solveRefinedElementProblem(), oofem::HuertaErrorEstimator::solveRefinedWholeProblem(), and oofem::PFEM::updateDofUnknownsDictionaryPressure().
|
pure virtual |
Test if Dof has initial condition.
Implemented in oofem::MasterDof, oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Referenced by oofem::SimpleSlaveDof::hasIc().
|
pure virtual |
Test if Dof has initial condition of required ValueModeType.
u | Type of required IC. |
Implemented in oofem::SlaveDof, oofem::MasterDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Referenced by giveBcValue(), and oofem::SimpleSlaveDof::hasIcOn().
|
inlinevirtual |
Tests if receiver is primary DOF.
Dof is primary if it possess or directly represent certain DOF. If it is linked somehow (rigid arm, doubled node) to other DOF(s) then it is not primary DOF.
Reimplemented in oofem::MasterDof, and oofem::ActiveDof.
Definition at line 287 of file dof.h.
Referenced by oofem::DofDistributedPrimaryField::applyBoundaryCondition(), oofem::Node::drawYourself(), oofem::DofManager::giveMasterDofMans(), oofem::DofManager::hasAnySlaveDofs(), oofem::EIPrimaryUnknownMapper::mapAndUpdate(), and oofem::DofManager::postInitialize().
|
inlinevirtual |
Packs specific DOF Manager's dofs unknowns into communication buffer.
If dof is slave, then no packing is done, this is maintained by master. This requires master be available at same partition as slave.
buff | Communication buffer to pack data. |
mode | Mode of unknown (e.g, total value, velocity or acceleration of unknown). |
tStep | Time step when unknown requested. See documentation of particular EngngModel class for valid tStep values (most implementations can return only values for current and possibly for previous time step). |
|
virtual |
Prints Dof output (it prints value of unknown related to dof at given timeStep).
The format of output depends on analysis type. Called from corresponding e-model.
Definition at line 92 of file dof.C.
References dofID, and giveUnknown().
Referenced by oofem::DEIDynamic::printDofOutputAt(), oofem::DIIDynamic::printDofOutputAt(), oofem::NonLinearDynamic::printDofOutputAt(), and oofem::NlDEIDynamic::printDofOutputAt().
|
virtual |
Prints Dof output (it prints value of unknown related to dof at given timeStep).
The format of output depends on analysis type. Called from corresponding e-model.
Definition at line 76 of file dof.C.
References dofID, and giveUnknown().
Referenced by oofem::SUPG::printDofOutputAt(), oofem::CBS::printDofOutputAt(), oofem::PFEM::printDofOutputAt(), and oofem::EngngModel::printDofOutputAt().
void oofem::Dof::printSingleOutputWithAdditionAt | ( | FILE * | File, |
TimeStep * | tStep, | ||
char | ch, | ||
ValueModeType | mode, | ||
double | addend | ||
) |
Definition at line 84 of file dof.C.
References dofID, and giveUnknown().
|
virtual |
Prints the receiver state on stdout.
Reimplemented in oofem::MasterDof.
Definition at line 107 of file dof.C.
References dofID, dofManager, oofem::FEMComponent::giveClassName(), and oofem::FEMComponent::giveNumber().
Referenced by oofem::ElementDofManager::printYourself(), oofem::ElementSide::printYourself(), oofem::Node::printYourself(), and oofem::DofManager::printYourself().
|
virtual |
Restores the receiver state previously written in stream.
Reimplemented in oofem::MasterDof, oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Definition at line 148 of file dof.C.
References oofem::CIO_IOERR, oofem::CIO_OK, dofID, oofem::DataStream::read(), and THROW_CIOERR.
Referenced by oofem::SimpleSlaveDof::restoreContext(), oofem::SlaveDof::restoreContext(), oofem::MasterDof::restoreContext(), and oofem::DofManager::restoreContext().
|
virtual |
Stores receiver state to output stream.
Reimplemented in oofem::MasterDof, oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Definition at line 136 of file dof.C.
References oofem::CIO_IOERR, oofem::CIO_OK, dofID, THROW_CIOERR, and oofem::DataStream::write().
Referenced by oofem::SimpleSlaveDof::saveContext(), oofem::SlaveDof::saveContext(), and oofem::MasterDof::saveContext().
|
inlinevirtual |
Overwrites the boundary condition id (0-inactive BC), intended for specific purposes such as coupling of bc's in multiscale simulations.
Reimplemented in oofem::MasterDof, and oofem::ActiveDof.
Definition at line 382 of file dof.h.
Referenced by oofem::Domain::createDofs(), oofem::SolutionbasedShapeFunction::setBoundaryConditionOnDof(), oofem::IncrementalLinearStatic::solveYourselfAt(), and oofem::DelaunayTriangulator::writeMesh().
|
inline |
|
inlinevirtual |
Sets a specific equation number to receiver.
equationNumber | New equation number. |
Reimplemented in oofem::MasterDof.
|
inlinevirtual |
Overwrites the initial condition id (0-inactive IC)
Reimplemented in oofem::MasterDof.
Definition at line 384 of file dof.h.
Referenced by oofem::Domain::createDofs().
|
inlinevirtual |
|
inlinevirtual |
Unpacks DOF unknown from communication buffer and updates unknown if necessary.
Unknown is always updated using EngngModel::updateUnknownComponent, if DOFManager to which receiver belongs has DofManager_shared dofManagerParallelMode type. Unknown is unpacked and stored in unknowns dictionary, if DOFManager to which receiver belongs has DofManager_remote dofManagerParallelMode type. There is no reason for invoking this service if DOFManager has DofManager_local mode. If do is slave, then no unpacking and updating is done. This is left on master, which must be available on same partition.
buff | Buffer containing packed message. |
mode | Mode of unknown (e.g, total value, velocity or acceleration of unknown). |
tStep | Time step when unknown requested. See documentation of particular EngngModel class for valid tStep values (most implementations can return only values for current and possibly for previous time step). |
|
inlinevirtual |
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.
f | Function that converts old to new equation number. |
Reimplemented in oofem::SlaveDof, oofem::SimpleSlaveDof, and oofem::ActiveDof.
Definition at line 315 of file dof.h.
Referenced by oofem::DofManager::updateLocalNumbering().
|
inlinevirtual |
Abstract function, allowing Dof to store its unknowns in its own private dictionary.
Dof then uses this dictionary instead of forwarding the requests to EngngModel (with equationNUmber as parameter). If EngngModel does not support changes of static system (see EngngModel::requiresUnknownsDictionaryUpdate method), the dof forwards the requests for its unknowns to EngngModel, where unknowns are naturally kept. This is possible, because dof equation number is same during whole solution. But when changes of static system are allowed, several problem arise. For example by solving simple incremental static with allowed static changes, the incremental displacement vector of structure can not be added to total displacement vector of structure, because equation numbers may have changed, and one can not simply add these vector to obtain new total displacement vector, because incompatible displacement will be added. To solve this problem, unknown dictionary at dof level has been assumed. Dof then keeps its unknowns in its own private dictionary. After computing increment of solution, engngModel updates for each dof its unknowns in its dictionary (using updateUnknownsDictionary function). For aforementioned example engngModel updates incremental values but also total value by asking dof for previous total value (dof will use its dictionary, does not asks back EngngModel) adds corresponding increment and updates total value in dictionary. In fact on EngngModel level only incremental solution is stored, but total values are always stored in dofs dictionaries. Implementation is not provided, only interface declared. Children must implement this method.
tStep | Time step when unknowns are updated. In current version it is unused parameter. It is EngngModel responsibility to update values, and values stored in dictionary are always related to timeStep when they were lastly updated. |
mode | Mode of stored unknown. |
dofValue | Value of unknown. Old value will generally be lost. |
Reimplemented in oofem::MasterDof.
Definition at line 366 of file dof.h.
Referenced by oofem::DofDistributedPrimaryField::applyBoundaryCondition(), oofem::DofDistributedPrimaryField::applyInitialCondition(), and oofem::PFEM::updateDofUnknownsDictionaryPressure().
|
inlinevirtual |
Updates receiver after finishing time step.
tStep | Finished time step. |
Reimplemented in oofem::MasterDof.
Definition at line 336 of file dof.h.
Referenced by oofem::MasterDof::updateYourself(), and oofem::DofManager::updateYourself().
|
friend |
|
protected |
Physical meaning of DOF.
Definition at line 99 of file dof.h.
Referenced by oofem::ActiveDof::askNewEquationNumber(), oofem::MasterDof::askNewEquationNumber(), Dof(), errorInfo(), oofem::SimpleSlaveDof::giveMasterDof(), oofem::SlaveDof::initialize(), printMultipleOutputAt(), printSingleOutputAt(), printSingleOutputWithAdditionAt(), oofem::MasterDof::printYourself(), printYourself(), restoreContext(), and saveContext().
|
protected |
Link to related DofManager.
Definition at line 97 of file dof.h.
Referenced by oofem::ActiveDof::askNewEquationNumber(), oofem::MasterDof::askNewEquationNumber(), Dof(), errorInfo(), oofem::ActiveDof::giveActiveBoundaryCondition(), oofem::MasterDof::giveBc(), giveDofManGlobalNumber(), giveDofManNumber(), oofem::MasterDof::giveIc(), oofem::SimpleSlaveDof::giveMasterDof(), oofem::SlaveDof::giveMasterDof(), oofem::MasterDof::giveUnknown(), oofem::MasterDof::giveUnknownsDictionaryValue(), oofem::MasterDof::hasBc(), oofem::MasterDof::printYourself(), printYourself(), oofem::MasterDof::restoreContext(), oofem::SimpleSlaveDof::saveContext(), oofem::SlaveDof::saveContext(), oofem::MasterDof::saveContext(), and oofem::MasterDof::updateUnknownsDictionary().