OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
General purpose 3d structural element evaluator. More...
#include <space3delementevaluator.h>
Public Member Functions | |
Space3dStructuralElementEvaluator () | |
virtual | ~Space3dStructuralElementEvaluator () |
Protected Member Functions | |
virtual void | computeNMatrixAt (FloatMatrix &answer, GaussPoint *gp) |
Assemble interpolation matrix at given IP In case of IGAElements, N is assumed to contain only nonzero interpolation functions. More... | |
virtual void | computeBMatrixAt (FloatMatrix &answer, GaussPoint *gp) |
Assembles the strain-displacement matrix of the receiver at given integration point In case of IGAElements, B is assumed to contain only contribution from nonzero interpolation functions. More... | |
virtual double | computeVolumeAround (GaussPoint *gp) |
virtual void | computeStressVector (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) |
Computes the stress vector. More... | |
virtual void | computeConstitutiveMatrixAt (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) |
Computes constitutive matrix of receiver. More... | |
void | giveDofManDofIDMask (int inode, IntArray &answer) const |
Protected Member Functions inherited from oofem::StructuralElementEvaluator | |
StructuralElementEvaluator () | |
virtual | ~StructuralElementEvaluator () |
virtual void | giveCharacteristicMatrix (FloatMatrix &answer, CharType mtrx, TimeStep *tStep) |
virtual void | giveCharacteristicVector (FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep) |
virtual IntegrationRule * | giveMassMtrxIntegrationRule () |
Returns the integration rule for mass matrices, if relevant. More... | |
virtual void | giveMassMtrxIntegrationMask (IntArray &answer) |
Returns mask indicating, which unknowns (their type and ordering is the same as element unknown vector) participate in mass matrix integration. More... | |
virtual void | computeLumpedMassMatrix (FloatMatrix &answer, TimeStep *tStep) |
Computes lumped mass matrix of receiver. More... | |
virtual void | computeConsistentMassMatrix (FloatMatrix &answer, TimeStep *tStep, double &mass) |
Computes consistent mass matrix of receiver using numerical integration over element volume. More... | |
virtual Element * | giveElement ()=0 |
virtual void | computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) |
virtual void | giveInternalForcesVector (FloatArray &answer, TimeStep *tStep, bool useUpdatedGpRecord=false) |
void | computeVectorOf (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
void | computeVectorOf (PrimaryField &field, const IntArray &dofIdMask, ValueModeType u, TimeStep *tStep, FloatArray &answer) |
bool | isActivated (TimeStep *tStep) |
void | updateInternalState (TimeStep *tStep) |
void | computeStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, FloatArray &u) |
Optimized version, allowing to pass element displacements as parameter. More... | |
virtual int | giveIntegrationElementLocalCodeNumbers (IntArray &answer, Element *elem, IntegrationRule *ie) |
Assembles the local element code numbers of given integration element (sub-patch) This is done by obtaining list of nonzero shape functions and by collecting the code numbers of nodes corresponding to these shape functions. More... | |
Additional Inherited Members | |
Protected Attributes inherited from oofem::StructuralElementEvaluator | |
FloatMatrix | rotationMatrix |
int | rotationMatrixDefined |
Flag indicating if transformation matrix has been already computed. More... | |
General purpose 3d structural element evaluator.
Definition at line 44 of file space3delementevaluator.h.
|
inline |
Definition at line 47 of file space3delementevaluator.h.
|
inlinevirtual |
Definition at line 48 of file space3delementevaluator.h.
References computeBMatrixAt(), computeConstitutiveMatrixAt(), computeNMatrixAt(), computeStressVector(), and computeVolumeAround().
|
protectedvirtual |
Assembles the strain-displacement matrix of the receiver at given integration point In case of IGAElements, B is assumed to contain only contribution from nonzero interpolation functions.
Implements oofem::StructuralElementEvaluator.
Definition at line 62 of file space3delementevaluator.C.
References oofem::FloatMatrix::at(), oofem::FEInterpolation::evaldNdx(), oofem::StructuralElementEvaluator::giveElement(), oofem::GaussPoint::giveIntegrationRule(), oofem::Element::giveInterpolation(), oofem::IntegrationRule::giveKnotSpan(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::FloatMatrix::giveNumberOfRows(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by ~Space3dStructuralElementEvaluator().
|
protectedvirtual |
Computes constitutive matrix of receiver.
answer | Constitutive matrix. |
rMode | Material response mode of answer. |
gp | Integration point for which constitutive matrix is computed. |
tStep | Time step. |
Implements oofem::StructuralElementEvaluator.
Definition at line 106 of file space3delementevaluator.C.
References oofem::Element::giveCrossSection(), and oofem::StructuralElementEvaluator::giveElement().
Referenced by ~Space3dStructuralElementEvaluator().
|
protectedvirtual |
Assemble interpolation matrix at given IP In case of IGAElements, N is assumed to contain only nonzero interpolation functions.
Implements oofem::StructuralElementEvaluator.
Definition at line 51 of file space3delementevaluator.C.
References oofem::FloatMatrix::beNMatrixOf(), oofem::FEInterpolation::evalN(), oofem::StructuralElementEvaluator::giveElement(), oofem::GaussPoint::giveIntegrationRule(), oofem::Element::giveInterpolation(), oofem::IntegrationRule::giveKnotSpan(), oofem::GaussPoint::giveNaturalCoordinates(), and N.
Referenced by ~Space3dStructuralElementEvaluator().
|
protectedvirtual |
Computes the stress vector.
answer | Stress vector. |
strain | Strain vector. |
gp | Integration point for which stress is computed. |
tStep | Time step. |
Implements oofem::StructuralElementEvaluator.
Definition at line 100 of file space3delementevaluator.C.
References oofem::Element::giveCrossSection(), and oofem::StructuralElementEvaluator::giveElement().
Referenced by oofem::BsplinePlaneStressElement::drawScalar(), oofem::NURBSPlaneStressElement::drawScalar(), oofem::TSplinePlaneStressElement::drawScalar(), oofem::NURBSSpace3dElement::drawScalar(), and ~Space3dStructuralElementEvaluator().
|
protectedvirtual |
Reimplemented from oofem::StructuralElementEvaluator.
Definition at line 90 of file space3delementevaluator.C.
References oofem::StructuralElementEvaluator::giveElement(), oofem::GaussPoint::giveIntegrationRule(), oofem::IntegrationRule::giveKnotSpan(), oofem::GaussPoint::giveNaturalCoordinates(), and oofem::GaussPoint::giveWeight().
Referenced by ~Space3dStructuralElementEvaluator().
|
inlineprotected |
Definition at line 64 of file space3delementevaluator.h.
Referenced by oofem::NURBSSpace3dElement::giveDofManDofIDMask().