OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Class implementing node in finite element mesh. More...
#include <node.h>
Public Member Functions | |
Node (int n, Domain *aDomain) | |
Constructor. More... | |
virtual | ~Node () |
Destructor. More... | |
virtual bool | hasCoordinates () |
virtual double | giveCoordinate (int i) |
virtual FloatArray * | giveCoordinates () |
const FloatArray & | giveNodeCoordinates () const |
As giveCoordinates, but non-virtual and therefore faster (because it can be inlined). More... | |
void | setCoordinates (FloatArray coords) |
Sets node coordinates to given array. More... | |
virtual double | giveUpdatedCoordinate (int ic, TimeStep *tStep, double scale=1.) |
Returns updated ic-th coordinate of receiver. More... | |
virtual void | giveUpdatedCoordinates (FloatArray &answer, TimeStep *tStep, double scale=1.) |
Returns updated coordinate of receiver. More... | |
bool | hasLocalCS () |
Returns nonzero if node has prescribed local coordinate system. More... | |
FloatMatrix * | giveLocalCoordinateTriplet () |
Returns pointer to local coordinate triplet in node. More... | |
bool | hasSameLCS (Node *remote) |
Returns true, if the local coordinate systems of receiver and given node are the same. More... | |
virtual bool | computeL2GTransformation (FloatMatrix &answer, const IntArray &dofIDArry) |
Computes transformation matrix from global c.s. More... | |
virtual bool | requiresTransformation () |
Indicates, whether dofManager requires the transformation. More... | |
virtual void | computeLoadVector (FloatArray &answer, Load *load, CharType type, TimeStep *tStep, ValueModeType mode) |
Computes the load vector for given load. More... | |
virtual void | updateYourself (TimeStep *tStep) |
Updates receiver at end of time step (i.e. More... | |
virtual const char * | giveClassName () const |
virtual const char * | giveInputRecordName () const |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description stored in input record. More... | |
virtual void | giveInputRecord (DynamicInputRecord &input) |
Setups the input record string of receiver. More... | |
virtual void | printYourself () |
Prints receiver state on stdout. Useful for debugging. More... | |
virtual int | checkConsistency () |
Allows programmer to test some internal data, before computation begins. More... | |
virtual bool | isDofTypeCompatible (dofType type) const |
Returns true if dof of given type is allowed to be associated to receiver. More... | |
virtual int | giveQcNodeType () |
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 | drawYourself (oofegGraphicContext &gc, TimeStep *tStep) |
Public Member Functions inherited from oofem::DofManager | |
std::vector< Dof * >::iterator | begin () |
std::vector< Dof * >::iterator | end () |
std::vector< Dof * >::const_iterator | begin () const |
std::vector< Dof * >::const_iterator | end () const |
DofManager (int n, Domain *aDomain) | |
Constructor. More... | |
virtual | ~DofManager () |
Destructor. More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
Prints output of receiver to stream, for given time step. More... | |
bool | isBoundary () |
void | setBoundaryFlag (bool isBoundary) |
Sets the boundary flag. More... | |
virtual bool | hasAnySlaveDofs () |
virtual bool | giveMasterDofMans (IntArray &masters) |
Returns true if the receiver is linked (its slave DOFs depend on master values) to some other dof managers. More... | |
virtual void | postInitialize () |
Performs post-initialization such like checking if there are any slave dofs etc. More... | |
virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
Local renumbering support. More... | |
void | appendDof (Dof *dof) |
Adds the given Dof into the receiver. More... | |
void | removeDof (DofIDItem id) |
Removes Dof with given id from dofArray. More... | |
bool | hasDofID (DofIDItem id) const |
Checks if receiver contains dof with given ID. More... | |
int | giveGlobalNumber () const |
int | giveLabel () const |
void | setGlobalNumber (int newNumber) |
Sets receiver global number. More... | |
dofManagerParallelMode | giveParallelMode () const |
Return dofManagerParallelMode of receiver. More... | |
void | setParallelMode (dofManagerParallelMode _mode) |
Sets parallel mode of receiver. More... | |
const IntArray * | givePartitionList () |
Returns partition list of receiver. More... | |
void | setPartitionList (const IntArray *_p) |
Sets receiver's partition list. More... | |
void | removePartitionFromList (int _part) |
Removes given partition from receiver list. More... | |
void | mergePartitionList (IntArray &_p) |
Merges receiver partition list with given lists. More... | |
int | givePartitionsConnectivitySize () |
Returns number of partitions sharing given receiver (=number of shared partitions + local one). More... | |
bool | isLocal () |
Returns true if receiver is locally maintained. More... | |
bool | isShared () |
Returns true if receiver is shared. More... | |
bool | isNull () |
Returns true if receiver is shared. More... | |
Dof * | giveDofWithID (int dofID) const |
Returns DOF with given dofID; issues error if not present. More... | |
int | giveNumberOfDofs () const |
void | askNewEquationNumbers (TimeStep *tStep) |
Renumbers all contained DOFs. More... | |
int | giveNumberOfPrimaryMasterDofs (const IntArray &dofIDArray) const |
Returns the number of primary dofs on which receiver dofs (given in dofArray) depend on. More... | |
void | giveLocationArray (const IntArray &dofIDArry, IntArray &locationArray, const UnknownNumberingScheme &s) const |
Returns location array (array containing for each requested dof related equation number) for given numbering scheme. More... | |
void | giveMasterDofIDArray (const IntArray &dofIDArry, IntArray &masterDofIDs) const |
Returns master dof ID array of receiver. More... | |
void | giveCompleteLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s) const |
Returns full location array of receiver containing equation numbers of all dofs of receiver. More... | |
void | giveCompleteMasterDofIDArray (IntArray &dofIDArray) const |
Returns the full dof ID array of receiver. More... | |
std::vector< Dof * >::const_iterator | findDofWithDofId (DofIDItem dofID) const |
Finds index of DOF with required physical meaning of receiver. More... | |
void | giveUnknownVector (FloatArray &answer, const IntArray &dofMask, ValueModeType mode, TimeStep *tStep, bool padding=false) |
Assembles the vector of unknowns in global c.s for given dofs of receiver. More... | |
void | giveUnknownVector (FloatArray &answer, const IntArray &dofMask, PrimaryField &field, ValueModeType mode, TimeStep *tStep, bool padding=false) |
Assembles the vector of unknowns of given filed in global c.s for given dofs of receiver. More... | |
void | giveCompleteUnknownVector (FloatArray &answer, ValueModeType mode, TimeStep *tStep) |
Assembles the complete unknown vector in node. More... | |
void | giveUnknownVectorOfType (FloatArray &answer, UnknownType ut, ValueModeType mode, TimeStep *tStep) |
Constructs the requested vector by assembling e.g. More... | |
void | givePrescribedUnknownVector (FloatArray &answer, const IntArray &dofMask, ValueModeType mode, TimeStep *tStep) |
Assembles the vector of prescribed unknowns in nodal c.s for given dofs of receiver. More... | |
bool | computeM2GTransformation (FloatMatrix &answer, const IntArray &dofIDArry) |
Computes receiver transformation matrix from global CS to dofManager specific coordinate system; . More... | |
virtual bool | computeM2LTransformation (FloatMatrix &answer, const IntArray &dofIDArry) |
Computes transformation matrix from local DOFs to master DOFs; . More... | |
IntArray * | giveLoadArray () |
Returns the array containing applied loadings of the receiver. More... | |
void | setLoadArray (IntArray &load) |
Sets the array of applied loadings of the receiver. More... | |
const IntArray * | giveForcedDofIDs () |
Returns list of specific dofs that should be included in node. More... | |
std::map< int, int > * | giveDofTypeMap () |
Returns map from DofIDItem to dofType. More... | |
std::map< int, int > * | giveMasterMap () |
Returns map from DofIDItem to dofType. More... | |
std::map< int, int > * | giveBcMap () |
Returns map from DofIDItem to dofType. More... | |
std::map< int, int > * | giveIcMap () |
Returns map from DofIDItem to initial condition. More... | |
void | setNumberOfDofs (int _ndofs) |
Sets number of dofs of the receiver; Deallocates existing DOFs; Resizes the dofArray accordingly. More... | |
Public Member Functions inherited from oofem::FEMComponent | |
FEMComponent (int n, Domain *d) | |
Regular constructor, creates component with given number and belonging to given domain. More... | |
virtual | ~FEMComponent () |
Virtual destructor. More... | |
Domain * | giveDomain () const |
virtual void | setDomain (Domain *d) |
Sets associated Domain. More... | |
int | giveNumber () const |
void | setNumber (int num) |
Sets number of receiver. More... | |
virtual Interface * | giveInterface (InterfaceType t) |
Interface requesting service. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Protected Attributes | |
FloatArray | coordinates |
Array storing nodal coordinates. More... | |
FloatMatrix * | localCoordinateSystem |
Triplet defining the local coordinate system in node. More... | |
Protected Attributes inherited from oofem::DofManager | |
std::vector< Dof * > | dofArray |
Array of DOFs. More... | |
IntArray | loadArray |
List of applied loads. More... | |
bool | isBoundaryFlag |
Indicates if dofManager is boundary (true boundary or on boundary between regions) or interior. More... | |
bool | hasSlaveDofs |
Flag indicating whether receiver has slave DOFs. More... | |
int | globalNumber |
In parallel mode, globalNumber contains globally unique DoFManager number. More... | |
dofManagerParallelMode | parallel_mode |
IntArray | partitions |
List of partition sharing the shared dof manager or remote partition containing remote dofmanager counterpart. More... | |
IntArray * | dofidmask |
List of additional dof ids to include. More... | |
std::map< int, int > * | dofTypemap |
Map from DofIDItem to dofType. More... | |
std::map< int, int > * | dofMastermap |
Map from DofIDItem to master node. More... | |
std::map< int, int > * | dofBCmap |
Map from DofIDItem to bc (to be removed). More... | |
std::map< int, int > * | dofICmap |
Map from DofIDItem to ic (to be removed). More... | |
IntArray | mBC |
Protected Attributes inherited from oofem::FEMComponent | |
int | number |
Component number. More... | |
Domain * | domain |
Link to domain object, useful for communicating with other FEM components. More... | |
Class implementing node in finite element mesh.
Node possess degrees of freedom (see base class DofManager). Node is attribute of few elements and it is managed by domain. Node manages its position in space, and if specified local coordinate system in node. If local coordinate system is defined, all equilibrium equations are assembled in this system and therefore all DOFs and applied boundary and initial conditions apply in this local coordinate system. By default, global coordinate system is assumed in each node. For description, how to prescribe local coordinate system in node, see input file description section.
Tasks include:
oofem::Node::Node | ( | int | n, |
Domain * | aDomain | ||
) |
Constructor.
Creates a node belonging to domain.
Definition at line 68 of file node.C.
References localCoordinateSystem.
|
virtual |
|
virtual |
Allows programmer to test some internal data, before computation begins.
For example, one may use this function, to ensure that element has material with required capabilities is assigned to element. This must be done after all mesh components are instanciated.
Reimplemented from oofem::FEMComponent.
Reimplemented in oofem::RigidArmNode, oofem::qcNode, oofem::HangingNode, oofem::PFEMParticle, and oofem::InteractionPFEMParticle.
Definition at line 325 of file node.C.
References oofem::IntArray::at(), oofem::FEMComponent::checkConsistency(), oofem::FEMComponent::domain, oofem::Domain::giveDofManager(), oofem::SimpleSlaveDof::giveMasterDofManagerNum(), hasSameLCS(), and OOFEM_WARNING.
Referenced by oofem::PFEMParticle::checkConsistency(), oofem::HangingNode::checkConsistency(), oofem::qcNode::checkConsistency(), and oofem::RigidArmNode::checkConsistency().
|
virtual |
Computes transformation matrix from global c.s.
to DOF-manager specific c.s; .
answer | Computed transformation matrix. |
dofIDArry | Array containing DofIDItem-type values for which transformation matrix is assembled. If dofIDArry is NULL, then all receiver DOFs are assumed. |
Reimplemented from oofem::DofManager.
Definition at line 424 of file node.C.
References oofem::__DofIDItemToString(), oofem::IntArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::clear(), oofem::Dof::giveDofID(), oofem::DofManager::giveNumberOfDofs(), oofem::IntArray::giveSize(), oofem::IntArray::isEmpty(), localCoordinateSystem, OOFEM_ERROR, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by computeLoadVector(), oofem::RigidArmNode::computeMasterContribution(), drawYourself(), and oofem::InteractionPFEMParticle::givePrescribedUnknownVector().
|
virtual |
Computes the load vector for given load.
answer | Load vector for given load. |
load | Given load. |
type | Characteristic type of the vector. |
tStep | Time step when answer is computed. |
mode | Determines response mode. |
Reimplemented from oofem::DofManager.
Definition at line 178 of file node.C.
References oofem::FloatArray::clear(), oofem::Load::computeComponentArrayAt(), computeL2GTransformation(), oofem::Load::CST_Global, oofem::NodalLoad::giveBCGeoType(), oofem::NodalLoad::giveCoordSystMode(), oofem::GeneralBoundaryCondition::giveDofIDs(), oofem::NodalLoadBGT, OOFEM_ERROR, and oofem::FloatArray::rotatedWith().
Referenced by drawYourself().
|
virtual |
Reimplemented from oofem::DofManager.
Definition at line 640 of file node.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatArray::clear(), computeL2GTransformation(), computeLoadVector(), oofem::DofManager_local, oofem::DofManager_shared, oofem::FEMComponent::domain, oofem::oofegGraphicContext::getBcForceColor(), oofem::oofegGraphicContext::getBcIcColor(), oofem::oofegGraphicContext::getCrackPatternColor(), oofem::oofegGraphicContext::getDeformedElementColor(), oofem::oofegGraphicContext::getDefScale(), oofem::oofegGraphicContext::getNodeColor(), giveCoordinate(), oofem::Dof::giveDofID(), oofem::FEMComponent::giveDomain(), oofem::Domain::giveEngngModel(), oofem::XfemManager::giveEnrichmentItem(), oofem::DofManager::giveGlobalNumber(), oofem::oofegGraphicContext::giveIntVarPlotMode(), oofem::oofegGraphicContext::giveIntVarType(), oofem::Domain::giveLoad(), oofem::DofManager::giveLoadArray(), giveLocalCoordinateTriplet(), oofem::FEMComponent::giveNumber(), oofem::XfemManager::giveNumberOfEnrichmentItems(), oofem::DofManager::giveParallelMode(), oofem::Dof::giveUnknown(), oofem::Domain::giveXfemManager(), oofem::Dof::hasBc(), hasLocalCS(), oofem::EnrichmentItem::isDofManEnriched(), oofem::IntArray::isEmpty(), oofem::Dof::isPrimaryDof(), oofem::OGC_essentialBC, oofem::OGC_naturalBC, oofem::OGC_nodeAnnotation, oofem::OGC_nodeGeometry, oofem::OGC_nodeVectorPlot, OOFEG_BCIC_ANNOTATION_LAYER, OOFEG_DEFORMED_GEOMETRY_LAYER, OOFEG_NATURALBC_LAYER, OOFEG_NODE_ANNOTATION_LAYER, and oofem::FloatArray::rotatedWith().
|
inlinevirtual |
Implements oofem::FEMComponent.
Reimplemented in oofem::RigidArmNode, oofem::qcNode, oofem::PFEMParticle, oofem::HangingNode, oofem::InteractionPFEMParticle, oofem::SlaveNode, and oofem::Particle.
|
virtual |
Reimplemented from oofem::DofManager.
Definition at line 82 of file node.C.
References oofem::FloatArray::at(), coordinates, and oofem::FloatArray::giveSize().
Referenced by oofem::TR1_2D_CBS::checkConsistency(), oofem::TR1_2D_PFEM::checkConsistency(), oofem::IntElLine1PhF::computeAreaAround(), oofem::IntElLine1::computeAreaAround(), oofem::TR1_2D_SUPG_AXI::computeBCRhsTerm_MB(), oofem::TR1_2D_SUPG2_AXI::computeBCRhsTerm_MB(), oofem::TR1_2D_SUPG::computeBCRhsTerm_MB(), oofem::TR1_2D_SUPG2::computeBCRhsTerm_MB(), oofem::AxisymElement::computeBHmatrixAt(), oofem::CohesiveSurface3d::computeBmatrixAt(), oofem::Q4Axisymm::computeBmatrixAt(), oofem::L4Axisymm::computeBmatrixAt(), oofem::Truss2d::computeBmatrixAt(), oofem::Lattice2d::computeBmatrixAt(), oofem::AxisymElement::computeBmatrixAt(), oofem::HTSelement::computeCenterOfGravity(), oofem::IntElLine1IntPen::computeCovarBaseVectorAt(), oofem::IntElLine2IntPen::computeCovarBaseVectorAt(), oofem::IntElLine1PF::computeCovarBaseVectorAt(), oofem::IntElLine1PhF::computeCovarBaseVectorAt(), oofem::IntElLine1::computeCovarBaseVectorAt(), oofem::TR1_2D_PFEM::computeCriticalTimeStep(), oofem::TR1_2D_CBS::computeDensityRhsVelocityTerms(), oofem::TR1_2D_CBS::computeDiffusionTermsI(), oofem::FRCFCMNL::computeElementCentroid(), oofem::InterfaceElement3dTrLin::computeGlobalCoordinates(), oofem::LIBeam2dNL::computeGlobalCoordinates(), oofem::LIBeam2d::computeGlobalCoordinates(), oofem::LIBeam3d::computeGlobalCoordinates(), oofem::Truss2d::computeGlobalCoordinates(), oofem::LIBeam3dNL2::computeGlobalCoordinates(), oofem::LIBeam3dNL::computeGlobalCoordinates(), oofem::LIBeam3d2::computeGlobalCoordinates(), oofem::Beam3d::computeGlobalCoordinates(), oofem::CCTPlate3d::computeGlobalCoordinates(), oofem::DKTPlate3d::computeGlobalCoordinates(), oofem::MITC4Shell::computeGlobalCoordinates(), oofem::InterfaceElem2dLin::computeGtoLRotationMatrix(), oofem::InterfaceElem2dQuad::computeGtoLRotationMatrix(), oofem::Lattice2d_mt::computeInternalSourceRhsVectorAt(), oofem::InterfaceElement3dTrLin::computeLCS(), oofem::LIBeam2d::computeLength(), oofem::LIBeam2dNL::computeLength(), oofem::LIBeam3d::computeLength(), oofem::Truss2d::computeLength(), oofem::LIBeam3dNL::computeLength(), oofem::LIBeam3dNL2::computeLength(), oofem::LIBeam3d2::computeLength(), oofem::Beam2d::computeLength(), oofem::Beam3d::computeLength(), oofem::CCTPlate::computeLoadLEToLRotationMatrix(), oofem::Quad1Mindlin::computeLoadLEToLRotationMatrix(), oofem::DKTPlate::computeLoadLEToLRotationMatrix(), oofem::QDKTPlate::computeLoadLEToLRotationMatrix(), oofem::Quad1MindlinShell3D::computeLoadLEToLRotationMatrix(), oofem::SurfaceTensionBoundaryCondition::computeLoadVectorFromElement(), oofem::InterfaceElem1d::computeLocalSlipDir(), oofem::TR1_2D_SUPG_AXI::computeOutFlowBCTerm_MB(), oofem::TR1_2D_SUPG::computeOutFlowBCTerm_MB(), oofem::HTSelement::computeOutwardNormalMatrix(), oofem::TR1_2D_SUPG_AXI::computePenetrationWithResistanceBCTerm_MB(), oofem::TR1_2D_SUPG::computePenetrationWithResistanceBCTerm_MB(), oofem::TR1_2D_CBS::computePrescribedTractionPressure(), oofem::TR1_2D_SUPG_AXI::computeRadiusAt(), oofem::TR1_2D_SUPG2_AXI::computeRadiusAt(), oofem::TR1_2D_SUPG_AXI::computeSlipWithFrictionBCTerm_MB(), oofem::TR1_2D_SUPG::computeSlipWithFrictionBCTerm_MB(), oofem::HTSelement::computeSvMatrixAt(), oofem::SurfaceTensionBoundaryCondition::computeTangentFromElement(), oofem::HTSelement::computeUvMatrixAt(), oofem::InterfaceElem2dLin::computeVolumeAround(), oofem::InterfaceElem2dQuad::computeVolumeAround(), oofem::LIBeam3dNL::computeXdVector(), oofem::LIBeam3dNL2::computeXdVector(), oofem::FreemInterface::createInput(), oofem::T3DInterface::createInput(), oofem::Element::drawAnnotation(), oofem::QTrPlaneStress2d::drawRawGeometry(), oofem::IntElLine2::drawRawGeometry(), oofem::QPlaneStrain::drawRawGeometry(), oofem::QTrPlaneStrain::drawRawGeometry(), oofem::QPlaneStress2d::drawRawGeometry(), oofem::InterfaceElement3dTrLin::drawRawGeometry(), oofem::InterfaceElem2dLin::drawRawGeometry(), oofem::InterfaceElem2dQuad::drawRawGeometry(), oofem::Quad1_ht::drawRawGeometry(), oofem::Axisymm3d::drawRawGeometry(), oofem::Tetrah1_ht::drawRawGeometry(), oofem::Q4Axisymm::drawRawGeometry(), oofem::L4Axisymm::drawRawGeometry(), oofem::Tet1_3D_SUPG::drawRawGeometry(), oofem::TrPlaneStress2d::drawRawGeometry(), oofem::LTRSpace::drawRawGeometry(), oofem::InterfaceElem1d::drawRawGeometry(), oofem::IntElSurfTr1::drawRawGeometry(), oofem::LumpedMassElement::drawRawGeometry(), oofem::Brick1_ht::drawRawGeometry(), oofem::Quad1PlaneStrain::drawRawGeometry(), oofem::LIBeam2dNL::drawRawGeometry(), oofem::PlaneStress2d::drawRawGeometry(), oofem::CohesiveSurface3d::drawRawGeometry(), oofem::Truss3d::drawRawGeometry(), oofem::LIBeam3d::drawRawGeometry(), oofem::Truss1d::drawRawGeometry(), oofem::TrPlaneStrain::drawRawGeometry(), oofem::LIBeam3dNL2::drawRawGeometry(), oofem::IntElLine1::drawRawGeometry(), oofem::Truss2d::drawRawGeometry(), oofem::LIBeam3dNL::drawRawGeometry(), oofem::LSpace::drawRawGeometry(), oofem::TR_SHELL02::drawRawGeometry(), oofem::TR_SHELL01::drawRawGeometry(), oofem::TR1_2D_SUPG2_AXI::drawRawGeometry(), oofem::Lattice2d::drawRawGeometry(), oofem::IntElPoint::drawRawGeometry(), oofem::TR21_2D_SUPG::drawRawGeometry(), oofem::LIBeam3d2::drawRawGeometry(), oofem::Lattice2d_mt::drawRawGeometry(), oofem::Quad10_2D_SUPG::drawRawGeometry(), oofem::TR1_2D_PFEM::drawRawGeometry(), oofem::Beam2d::drawRawGeometry(), oofem::CCTPlate::drawRawGeometry(), oofem::QDKTPlate::drawRawGeometry(), oofem::TR1_2D_CBS::drawRawGeometry(), oofem::TR1_2D_SUPG2::drawRawGeometry(), oofem::DKTPlate::drawRawGeometry(), oofem::TR1_2D_SUPG::drawRawGeometry(), oofem::Beam3d::drawRawGeometry(), oofem::QTrPlaneStress2d::drawScalar(), oofem::QPlaneStrain::drawScalar(), oofem::QTrPlaneStrain::drawScalar(), oofem::QPlaneStress2d::drawScalar(), oofem::Quad1_ht::drawScalar(), oofem::Tetrah1_ht::drawScalar(), oofem::Axisymm3d::drawScalar(), oofem::L4Axisymm::drawScalar(), oofem::TrPlaneStress2d::drawScalar(), oofem::LTRSpace::drawScalar(), oofem::InterfaceElem1d::drawScalar(), oofem::Brick1_ht::drawScalar(), oofem::Quad1PlaneStrain::drawScalar(), oofem::PlaneStress2d::drawScalar(), oofem::CohesiveSurface3d::drawScalar(), oofem::Truss1d::drawScalar(), oofem::TrPlaneStrain::drawScalar(), oofem::LSpace::drawScalar(), oofem::PFEMParticle::drawScalar(), oofem::TR_SHELL02::drawScalar(), oofem::TR_SHELL01::drawScalar(), oofem::TR1_2D_SUPG2_AXI::drawScalar(), oofem::TR21_2D_SUPG::drawScalar(), oofem::IntElPoint::drawScalar(), oofem::LIBeam3d2::drawScalar(), oofem::Quad10_2D_SUPG::drawScalar(), oofem::TR1_2D_PFEM::drawScalar(), oofem::CCTPlate::drawScalar(), oofem::QDKTPlate::drawScalar(), oofem::TR1_2D_CBS::drawScalar(), oofem::TR1_2D_SUPG2::drawScalar(), oofem::DKTPlate::drawScalar(), oofem::TR1_2D_SUPG::drawScalar(), oofem::TrPlaneStress2d::drawSpecial(), oofem::LTRSpace::drawSpecial(), drawYourself(), oofem::InsertNode::evaluate(), oofem::CohesiveSurface3d::evaluateCenter(), oofem::CohesiveSurface3d::evaluateLocalCoordinateSystem(), oofem::LEPlic::findCellLineConstant(), oofem::TR1_2D_CBS::formMaterialVolumePoly(), oofem::TR1_2D_SUPG2_AXI::formMaterialVolumePoly(), oofem::TR1_2D_SUPG2::formMaterialVolumePoly(), oofem::TR1_2D_SUPG::formMaterialVolumePoly(), oofem::TR1_2D_CBS::formMyVolumePoly(), oofem::TR1_2D_SUPG2_AXI::formMyVolumePoly(), oofem::TR1_2D_SUPG2::formMyVolumePoly(), oofem::TR1_2D_SUPG::formMyVolumePoly(), oofem::TR1_2D_CBS::formVolumeInterfacePoly(), oofem::TR1_2D_SUPG2_AXI::formVolumeInterfacePoly(), oofem::TR1_2D_SUPG2::formVolumeInterfacePoly(), oofem::TR1_2D_SUPG::formVolumeInterfacePoly(), oofem::VTKXMLExportModule::getNodalVariableFromIS(), oofem::Element::giveCharacteristicLengthForAxisymmElements(), oofem::PlaneStress2d::giveCharacteristicSize(), oofem::TrPlaneStress2d::giveCharacteristicSize(), oofem::Beam3d::giveCompositeExportData(), oofem::Lattice2d::giveCrossSectionCoordinates(), oofem::Lattice2d_mt::giveCrossSectionCoordinates(), oofem::LIBeam3d2::giveCurrentLength(), oofem::TR1_2D_CBS::giveElementCenter(), oofem::TR1_2D_SUPG2_AXI::giveElementCenter(), oofem::TR1_2D_SUPG2::giveElementCenter(), oofem::TR1_2D_SUPG::giveElementCenter(), oofem::PolygonLine::giveInputRecord(), oofem::Beam3d::giveInternalForcesVectorAtPoint(), oofem::StructuralElement::giveInternalStateAtNode(), oofem::CohesiveSurface3d::giveLength(), oofem::Lattice2d::giveLength(), oofem::Lattice2d_mt::giveLength(), oofem::LIBeam3d::giveLocalCoordinateSystem(), oofem::LIBeam3d2::giveLocalCoordinateSystem(), oofem::LIBeam3dNL::giveLocalCoordinateSystem(), oofem::LIBeam3dNL2::giveLocalCoordinateSystem(), oofem::LIBeam2d::givePitch(), oofem::LIBeam2dNL::givePitch(), oofem::Truss2d::givePitch(), oofem::Lattice2d::givePitch(), oofem::Beam2d::givePitch(), oofem::HTSelement::giveSideLength(), giveUpdatedCoordinate(), oofem::TR1_2D_SUPG_AXI::initGeometry(), oofem::TR1_2D_SUPG::initGeometry(), oofem::LevelSetPCS::initialize(), oofem::CohesiveSurface3d::initializeFrom(), oofem::PolygonLine::intersects(), oofem::PolylineNonlocalBarrier::isActivated(), oofem::FRCFCMNL::isInElementProjection(), oofem::TR1_2D_SUPG::LS_PCS_computeS(), oofem::Tet1_3D_SUPG::LS_PCS_computeVOFFractions(), oofem::TR1_2D_SUPG_AXI::LS_PCS_computeVOFFractions(), oofem::TR21_2D_SUPG::LS_PCS_computeVOFFractions(), oofem::TR1_2D_SUPG::LS_PCS_computeVOFFractions(), oofem::AbaqusUserElement::postInitialize(), printYourself(), oofem::TR1_2D_SUPG2_AXI::updateIntegrationRules(), oofem::TR1_2D_SUPG2::updateIntegrationRules(), oofem::TR1_2D_SUPG2_AXI::updateVolumePolygons(), and oofem::TR1_2D_SUPG2::updateVolumePolygons().
|
inlinevirtual |
Reimplemented from oofem::DofManager.
Definition at line 114 of file node.h.
Referenced by oofem::OctreeSpatialLocalizer::buildOctreeDataStructure(), oofem::IntElSurfTr1::computeCovarBaseVectorsAt(), oofem::Tet1_3D_SUPG::computeCriticalTimeStep(), oofem::IntElSurfTr1::computeGlobalCoordinates(), oofem::InterfaceElem1d::computeGlobalCoordinates(), oofem::IntElPoint::computeGlobalCoordinates(), oofem::Shell7Base::computeGlobalCoordinates(), oofem::GeometryBasedEI::computeIntersectionPoints(), oofem::Quasicontinuum::computeIntersectionsOfLinkWith2DTringleElements(), oofem::Quasicontinuum::computeIntersectionsOfLinkWith3DTetrahedraElements(), oofem::DKTPlate3d::computeLoadLEToLRotationMatrix(), oofem::TrPlanestressRotAllman::computeLocalNodalCoordinates(), oofem::TrPlanestressRotAllman3d::computeLocalNodalCoordinates(), oofem::LinQuad3DPlaneStress::computeLocalNodalCoordinates(), oofem::IntElPoint::computeLocalSlipDir(), oofem::RigidArmNode::computeMasterContribution(), oofem::Line::computeNumberOfIntersectionPoints(), oofem::MeshQualityErrorEstimator::computeTriangleRadiusError(), oofem::TransportGradientDirichlet::computeXi(), oofem::Subdivision::createMesh(), oofem::PFEM::deactivateTooCloseParticles(), oofem::LEPlic::doLagrangianPhase(), oofem::Shell7Base::evalInitialCovarBaseVectorsAt(), oofem::ClosestNode::evaluate(), oofem::PrescribedGenStrainShell7::evaluateHigherOrderContribution(), oofem::PrescribedGradientBCPeriodic::findSlaveToMasterMap(), oofem::TransportGradientPeriodic::findSlaveToMasterMap(), oofem::VTKXMLExportModule::getNodalVariableFromXFEMST(), oofem::DummySpatialLocalizer::giveAllNodesWithinBox(), oofem::MITC4Shell::giveDirectorVectors(), oofem::Element::giveLengthInDir(), oofem::Beam3d::giveLocalCoordinateSystem(), oofem::DummySpatialLocalizer::giveNodeClosestToPoint(), oofem::OctreeSpatialLocalizer::giveNodeClosestToPointWithinOctant(), oofem::TrPlaneStrRot::giveNodeCoordinates(), oofem::CCTPlate::giveNodeCoordinates(), oofem::QDKTPlate::giveNodeCoordinates(), oofem::DKTPlate::giveNodeCoordinates(), oofem::OctreeSpatialLocalizer::giveNodesWithinBox(), oofem::MITC4Shell::giveThickness(), oofem::FEIIGAElementGeometryWrapper::giveVertexCoordinates(), oofem::PlaneStress2d::HuertaErrorEstimatorI_setupRefinedElementProblem(), oofem::Quad1PlaneStrain::HuertaErrorEstimatorI_setupRefinedElementProblem(), oofem::TrPlaneStrain::HuertaErrorEstimatorI_setupRefinedElementProblem(), oofem::LSpace::HuertaErrorEstimatorI_setupRefinedElementProblem(), oofem::LTRSpace::HuertaErrorEstimatorI_setupRefinedElementProblem(), oofem::TrPlaneStress2d::HuertaErrorEstimatorI_setupRefinedElementProblem(), oofem::Truss1d::HuertaErrorEstimatorI_setupRefinedElementProblem(), oofem::SolutionbasedShapeFunction::init(), oofem::OctreeSpatialLocalizer::initElementIPDataStructure(), oofem::OctreeSpatialLocalizer::insertNodeIntoOctree(), oofem::QCFullsolveddomain::isNodeInside(), oofem::QClinearStatic::nodeInFullSolvedDomainTest(), oofem::XfemElementInterface::partitionEdgeSegment(), oofem::MicroMaterial::setMacroProperties(), oofem::Shell7Base::setupInitialNodeDirectors(), oofem::Shell7Base::setupInitialSolutionVector(), oofem::FreemInterface::smoothNodalDensities(), oofem::SpatialLocalizerInterface::SpatialLocalizerI_giveBBox(), oofem::TR_SHELL01::SpatialLocalizerI_giveBBox(), oofem::TR_SHELL02::SpatialLocalizerI_giveBBox(), oofem::Quasicontinuum::stiffnessAssignment(), oofem::QClinearStatic::transformMeshToParticles(), oofem::PrescribedGenStrainShell7::updateCoefficientMatrix(), oofem::GeometryBasedEI::updateLevelSets(), oofem::FastMarchingMethod::updateTrialValue(), oofem::Shell7Base::vtkEvalInitialGlobalCoordinateAt(), and oofem::Shell7Base::vtkEvalInitialGlobalCZCoordinateAt().
|
virtual |
Setups the input record string of receiver.
input | Dynamic input record to be filled by receiver. |
Reimplemented from oofem::DofManager.
Definition at line 165 of file node.C.
References _IFT_Node_coords, _IFT_Node_lcs, coordinates, oofem::DofManager::giveInputRecord(), localCoordinateSystem, and oofem::DynamicInputRecord::setField().
|
inlinevirtual |
Implements oofem::FEMComponent.
Reimplemented in oofem::RigidArmNode, oofem::qcNode, oofem::PFEMParticle, oofem::HangingNode, oofem::InteractionPFEMParticle, oofem::SlaveNode, and oofem::Particle.
Definition at line 177 of file node.h.
References _IFT_Node_Name.
|
inline |
Returns pointer to local coordinate triplet in node.
If not defined, returns NULL.
Definition at line 158 of file node.h.
Referenced by drawYourself(), oofem::VTKXMLExportModule::getNodalVariableFromPrimaryField(), oofem::RefinedElement::giveCompatibleBcDofArray(), giveUpdatedCoordinate(), hasSameLCS(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem1D(), oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem2D(), and oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem3D().
|
inline |
As giveCoordinates, but non-virtual and therefore faster (because it can be inlined).
/ES
Definition at line 120 of file node.h.
Referenced by oofem::XfemElementInterface::ComputeBOrBHMatrix(), oofem::GeometryBasedEI::computeIntersectionPoints(), oofem::GeometryBasedEI::evaluateEnrFuncInNode(), oofem::Delamination::findInitiationFronts(), oofem::Shell7BaseXFEM::giveCZExportData(), oofem::Shell7BaseXFEM::giveShellExportData(), oofem::FEIElementGeometryWrapper::giveVertexCoordinates(), oofem::HybridEI::interpGradLevelSet(), oofem::HybridEI::interpLevelSet(), oofem::HybridEI::interpLevelSetTangential(), oofem::PLCZdamageRadius::propagateInterface(), oofem::PLnodeRadius::propagateInterface(), and oofem::GeometryBasedEI::updateNodeEnrMarker().
|
inlinevirtual |
Reimplemented in oofem::qcNode.
Definition at line 183 of file node.h.
References gc.
Referenced by oofem::Quasicontinuum::applyApproach2(), oofem::Quasicontinuum::applyApproach3(), and oofem::QClinearStatic::setActivatedNodeList().
|
virtual |
Returns updated ic-th coordinate of receiver.
Return value is computed as coordinate + scale * displacement, where corresponding displacement is obtained from corresponding nodal DOF. Local coordinate system is taken into account. Useful mainly for postprocessing.
ic | Index of coordinate. |
tStep | Time step for the displacement. |
scale | Scaling of displacement. |
Definition at line 245 of file node.C.
References oofem::FloatMatrix::at(), giveCoordinate(), oofem::Dof::giveDofID(), giveLocalCoordinateTriplet(), oofem::FloatMatrix::giveNumberOfRows(), oofem::TimeStep::giveTimeIncrement(), oofem::Dof::giveUnknown(), hasLocalCS(), oofem::TimeStep::isTheCurrentTimeStep(), OOFEM_ERROR, and oofem::FloatArray::zero().
Referenced by oofem::QTrPlaneStress2d::drawDeformedGeometry(), oofem::IntElLine2::drawDeformedGeometry(), oofem::QPlaneStrain::drawDeformedGeometry(), oofem::QTrPlaneStrain::drawDeformedGeometry(), oofem::QPlaneStress2d::drawDeformedGeometry(), oofem::InterfaceElem2dLin::drawDeformedGeometry(), oofem::InterfaceElem2dQuad::drawDeformedGeometry(), oofem::Axisymm3d::drawDeformedGeometry(), oofem::Q4Axisymm::drawDeformedGeometry(), oofem::L4Axisymm::drawDeformedGeometry(), oofem::TrPlaneStress2d::drawDeformedGeometry(), oofem::InterfaceElem1d::drawDeformedGeometry(), oofem::LTRSpace::drawDeformedGeometry(), oofem::LumpedMassElement::drawDeformedGeometry(), oofem::PlaneStress2d::drawDeformedGeometry(), oofem::LIBeam2dNL::drawDeformedGeometry(), oofem::Quad1PlaneStrain::drawDeformedGeometry(), oofem::Truss3d::drawDeformedGeometry(), oofem::CohesiveSurface3d::drawDeformedGeometry(), oofem::LIBeam3d::drawDeformedGeometry(), oofem::Truss1d::drawDeformedGeometry(), oofem::TrPlaneStrain::drawDeformedGeometry(), oofem::LIBeam3dNL2::drawDeformedGeometry(), oofem::IntElLine1::drawDeformedGeometry(), oofem::Truss2d::drawDeformedGeometry(), oofem::LIBeam3dNL::drawDeformedGeometry(), oofem::LSpace::drawDeformedGeometry(), oofem::TR_SHELL02::drawDeformedGeometry(), oofem::TR_SHELL01::drawDeformedGeometry(), oofem::Lattice2d::drawDeformedGeometry(), oofem::IntElPoint::drawDeformedGeometry(), oofem::LIBeam3d2::drawDeformedGeometry(), oofem::CCTPlate::drawDeformedGeometry(), oofem::Beam2d::drawDeformedGeometry(), oofem::QDKTPlate::drawDeformedGeometry(), oofem::DKTPlate::drawDeformedGeometry(), oofem::Beam3d::drawDeformedGeometry(), oofem::QTrPlaneStress2d::drawScalar(), oofem::QPlaneStrain::drawScalar(), oofem::QTrPlaneStrain::drawScalar(), oofem::QPlaneStress2d::drawScalar(), oofem::Axisymm3d::drawScalar(), oofem::L4Axisymm::drawScalar(), oofem::TrPlaneStress2d::drawScalar(), oofem::InterfaceElem1d::drawScalar(), oofem::LTRSpace::drawScalar(), oofem::PlaneStress2d::drawScalar(), oofem::Quad1PlaneStrain::drawScalar(), oofem::CohesiveSurface3d::drawScalar(), oofem::Truss1d::drawScalar(), oofem::TrPlaneStrain::drawScalar(), oofem::LSpace::drawScalar(), oofem::TR_SHELL02::drawScalar(), oofem::TR_SHELL01::drawScalar(), oofem::IntElPoint::drawScalar(), oofem::LIBeam3d2::drawScalar(), oofem::CCTPlate::drawScalar(), oofem::QDKTPlate::drawScalar(), oofem::DKTPlate::drawScalar(), oofem::TrPlaneStress2d::drawSpecial(), oofem::LTRSpace::drawSpecial(), oofem::PlaneStress2d::drawSpecial(), oofem::Quad1PlaneStrain::drawSpecial(), oofem::VTKXMLExportModule::getNodalVariableFromIS(), and oofem::StructuralElement::giveInternalStateAtNode().
|
virtual |
Returns updated coordinate of receiver.
Return value is computed as coordinate + scale * displacement, where corresponding displacement is obtained from corresponding nodal DOF. Local coordinate system is taken into account and the answer is given in global coordinates.
answer | Updated coordinate. |
tStep | Time step for the displacement. |
scale | Scaling of displacement. |
Definition at line 308 of file node.C.
References oofem::FloatArray::at(), coordinates, oofem::FloatArray::giveSize(), oofem::DofManager::giveUnknownVectorOfType(), and oofem::TimeStep::isTheCurrentTimeStep().
|
inlinevirtual |
Reimplemented from oofem::DofManager.
|
inline |
Returns nonzero if node has prescribed local coordinate system.
Definition at line 150 of file node.h.
Referenced by oofem::RigidArmNode::computeMasterContribution(), drawYourself(), oofem::VTKXMLExportModule::getNodalVariableFromPrimaryField(), giveUpdatedCoordinate(), hasSameLCS(), and saveContext().
bool oofem::Node::hasSameLCS | ( | Node * | remote | ) |
Returns true, if the local coordinate systems of receiver and given node are the same.
Definition at line 387 of file node.C.
References oofem::FloatMatrix::at(), giveLocalCoordinateTriplet(), and hasLocalCS().
Referenced by oofem::HangingNode::checkConsistency(), oofem::qcNode::checkConsistency(), and checkConsistency().
|
virtual |
Initializes receiver according to object description stored in input record.
This function is called immediately after creating object using constructor. Input record can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record.
ir | Input record to initialize from. |
Reimplemented from oofem::DofManager.
Reimplemented in oofem::RigidArmNode, oofem::qcNode, oofem::HangingNode, oofem::PFEMParticle, oofem::SlaveNode, oofem::InteractionPFEMParticle, and oofem::Particle.
Definition at line 93 of file node.C.
References _IFT_Node_coords, _IFT_Node_lcs, oofem::FloatArray::at(), oofem::FloatMatrix::at(), coordinates, oofem::FEMComponent::domain, oofem::Domain::giveEngngModel(), oofem::EngngModel::giveEquationScalingFlag(), oofem::FEMComponent::giveNumber(), oofem::FloatArray::giveSize(), oofem::EngngModel::giveVariableScale(), oofem::InputRecord::hasField(), oofem::DofManager::initializeFrom(), IR_GIVE_FIELD, oofem::IRRT_OK, localCoordinateSystem, OOFEM_ERROR, OOFEM_WARNING, oofem::FloatArray::times(), and oofem::VST_Length.
Referenced by oofem::Particle::initializeFrom(), oofem::SlaveNode::initializeFrom(), oofem::PFEMParticle::initializeFrom(), oofem::HangingNode::initializeFrom(), oofem::qcNode::initializeFrom(), and oofem::RigidArmNode::initializeFrom().
|
inlinevirtual |
Returns true if dof of given type is allowed to be associated to receiver.
Reimplemented from oofem::DofManager.
Reimplemented in oofem::RigidArmNode, oofem::qcNode, oofem::HangingNode, and oofem::SlaveNode.
|
virtual |
Prints receiver state on stdout. Useful for debugging.
Reimplemented from oofem::DofManager.
Definition at line 207 of file node.C.
References giveCoordinate(), oofem::DofManager::loadArray, oofem::FEMComponent::number, oofem::IntArray::printYourself(), and oofem::Dof::printYourself().
|
inlinevirtual |
Indicates, whether dofManager requires the transformation.
Reimplemented from oofem::DofManager.
|
virtual |
Restores the receiver state previously written in stream.
stream | Input stream. |
mode | Determines amount of info available in stream (state, definition, ...). |
obj | Special parameter for sending extra information. |
throws | an ContextIOERR exception if error encountered. |
Reimplemented from oofem::DofManager.
Definition at line 599 of file node.C.
References oofem::CIO_IOERR, oofem::CIO_OK, CM_Definition, coordinates, localCoordinateSystem, oofem::DataStream::read(), oofem::DofManager::restoreContext(), oofem::FloatArray::restoreYourself(), oofem::FloatMatrix::restoreYourself(), and THROW_CIOERR.
|
virtual |
Stores receiver state to output stream.
stream | Output stream. |
mode | Determines amount of info required in stream (state, definition, ...). |
obj | Special parameter, used only to send particular integration point to material class version of this method. |
throws | an ContextIOERR exception if error encountered. |
Reimplemented from oofem::DofManager.
Definition at line 565 of file node.C.
References oofem::CIO_IOERR, oofem::CIO_OK, CM_Definition, coordinates, hasLocalCS(), localCoordinateSystem, oofem::DofManager::saveContext(), oofem::FloatArray::storeYourself(), oofem::FloatMatrix::storeYourself(), THROW_CIOERR, and oofem::DataStream::write().
|
inline |
Sets node coordinates to given array.
coords | New coordinates for node. |
Definition at line 126 of file node.h.
Referenced by oofem::DelaunayTriangulator::buildInitialBBXMesh(), oofem::ParticleTopologyDescription::replaceFEMesh(), and oofem::T3DInterface::t3d_2_OOFEM().
|
virtual |
Updates receiver at end of time step (i.e.
after equilibrium has been reached). If EngngModel formulation ( see giveFormulation() member function) returns actualized Lagrange mode, receiver updates its coordinates according to solution.
tStep | Time step for which to update. |
Reimplemented from oofem::DofManager.
Reimplemented in oofem::PFEMParticle, and oofem::InteractionPFEMParticle.
Definition at line 222 of file node.C.
References oofem::AL, oofem::FloatArray::at(), coordinates, oofem::FEMComponent::domain, oofem::Dof::giveDofID(), oofem::Domain::giveEngngModel(), oofem::EngngModel::giveFormulation(), oofem::TimeStep::giveTimeIncrement(), oofem::Dof::giveUnknown(), and oofem::DofManager::updateYourself().
Referenced by oofem::PFEMParticle::updateYourself().
|
protected |
Array storing nodal coordinates.
Definition at line 91 of file node.h.
Referenced by giveCoordinate(), giveInputRecord(), giveUpdatedCoordinates(), initializeFrom(), oofem::HangingNode::postInitialize(), oofem::qcNode::postInitializeAsHangingNode(), restoreContext(), saveContext(), and updateYourself().
|
protected |
Triplet defining the local coordinate system in node.
Value at position (i,j) represents angle between e'(i) and e(j), where e' is base vector of local coordinate system and e is base vector of global c.s.
Definition at line 98 of file node.h.
Referenced by computeL2GTransformation(), oofem::RigidArmNode::computeMasterContribution(), giveInputRecord(), initializeFrom(), Node(), restoreContext(), saveContext(), and ~Node().