OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
This class represent a new concept on how to define elements. More...
#include <structuralelementevaluator.h>
Protected Member Functions | |
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 | computeNMatrixAt (FloatMatrix &answer, GaussPoint *gp)=0 |
Computes the matrix for which the unknown field is obtained, typically [N1, 0, N2, 0, ...; 0, N1, 0, N2, ...]. More... | |
virtual void | computeBMatrixAt (FloatMatrix &answer, GaussPoint *gp)=0 |
virtual void | computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) |
virtual double | computeVolumeAround (GaussPoint *gp) |
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) |
virtual void | computeStressVector (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)=0 |
Computes the stress vector. More... | |
virtual void | computeConstitutiveMatrixAt (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)=0 |
Computes constitutive matrix of receiver. More... | |
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... | |
Protected Attributes | |
FloatMatrix | rotationMatrix |
int | rotationMatrixDefined |
Flag indicating if transformation matrix has been already computed. More... | |
Friends | |
void | drawIGAPatchDeformedGeometry (Element *elem, StructuralElementEvaluator *se, oofegGraphicContext &gc, TimeStep *tStep, UnknownType) |
This class represent a new concept on how to define elements.
Traditionally, Elements are derived from problem specific class (structural element, for example) define their interpolation and implement methods to evaluate shape function matrix, geometrical matrix, etc. This new concept, here represented by StructuralElementEvaluator and derived classes allows to define all problem specific methods (including shape function matrix, geometrical matrix evaluation) to be defined once for all elements of the same type (plane stress elements, space3d elements, etc) just relying on services of finite element interpolation classes. Definition of particular element is then done simply by deriving if from evaluator and providing interpolation.
StructuralElementEvaluator - base class of all structural elements Individual elements supposed to be derived from StructuralElementEvaluator and IGAElement
Definition at line 56 of file structuralelementevaluator.h.
|
protected |
Definition at line 51 of file structuralelementevaluator.C.
References oofem::FloatMatrix::clear(), and rotationMatrix.
|
inlineprotectedvirtual |
Definition at line 64 of file structuralelementevaluator.h.
References giveCharacteristicMatrix(), and giveCharacteristicVector().
|
protectedpure virtual |
Implemented in oofem::Space3dStructuralElementEvaluator, and oofem::PlaneStressStructuralElementEvaluator.
Referenced by computeStiffnessMatrix(), computeStrainVector(), giveInternalForcesVector(), and giveMassMtrxIntegrationMask().
|
protectedvirtual |
Computes consistent mass matrix of receiver using numerical integration over element volume.
Mass matrix is computed as , where is displacement approximation matrix. The number of necessary integration points is determined using this->giveNumberOfIPForMassMtrxIntegration service. Only selected degrees of freedom participate in integration of mass matrix. This is described using dof mass integration mask. This mask is obtained from this->giveMassMtrxIntegrationgMask service. The nonzero mask value at i-th position indicates that i-th element DOF participates in mass matrix computation. The result is in element local coordinate system.
answer | Mass matrix. |
tStep | Time step. |
mass | Total mass of receiver. |
Definition at line 206 of file structuralelementevaluator.C.
References oofem::IntArray::at(), oofem::FloatMatrix::at(), computeNMatrixAt(), oofem::Element::computeNumberOfDofs(), computeVolumeAround(), oofem::CrossSection::give(), giveElement(), giveMassMtrxIntegrationMask(), giveMassMtrxIntegrationRule(), oofem::FloatMatrix::giveNumberOfRows(), oofem::StructuralElement::giveStructuralCrossSection(), isActivated(), oofem::IntArray::isEmpty(), OOFEM_ERROR, oofem::FloatMatrix::plusProductSymmUpper(), oofem::FloatMatrix::resize(), oofem::FloatMatrix::symmetrized(), and oofem::FloatMatrix::zero().
Referenced by computeLumpedMassMatrix(), giveCharacteristicMatrix(), and giveMassMtrxIntegrationMask().
|
protectedpure virtual |
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. |
Implemented in oofem::Space3dStructuralElementEvaluator, and oofem::PlaneStressStructuralElementEvaluator.
Referenced by computeStiffnessMatrix(), and isActivated().
|
protectedvirtual |
Computes lumped mass matrix of receiver.
Default implementation returns lumped consistent mass matrix. Then returns lumped mass transformed into nodal coordinate system. The lumping procedure zeros all off-diagonal members and zeros also all diagonal members corresponding to non-displacement DOFs. Such diagonal matrix is then rescaled, to preserve the element mass. Requires the computeNmatrixAt and giveMassMtrxIntegrationgMask services to be implemented.
answer | Lumped mass matrix. |
tStep | Time step. |
Definition at line 147 of file structuralelementevaluator.C.
References oofem::IntArray::at(), oofem::FloatMatrix::at(), computeConsistentMassMatrix(), oofem::Element::computeNumberOfDofs(), oofem::Element::giveDofManDofIDMask(), giveElement(), oofem::Element::giveNumberOfDofManagers(), oofem::FloatMatrix::giveNumberOfRows(), oofem::IntArray::giveSize(), isActivated(), OOFEM_ERROR, oofem::FloatMatrix::resize(), oofem::FloatMatrix::times(), and oofem::FloatMatrix::zero().
Referenced by giveCharacteristicMatrix(), and giveMassMtrxIntegrationMask().
|
protectedpure virtual |
Computes the matrix for which the unknown field is obtained, typically [N1, 0, N2, 0, ...; 0, N1, 0, N2, ...].
Implemented in oofem::Space3dStructuralElementEvaluator, and oofem::PlaneStressStructuralElementEvaluator.
Referenced by computeConsistentMassMatrix(), oofem::drawIGAPatchDeformedGeometry(), and giveMassMtrxIntegrationMask().
|
protectedvirtual |
Definition at line 372 of file structuralelementevaluator.C.
References oofem::FloatMatrix::assemble(), oofem::FloatMatrix::beProductOf(), oofem::FloatMatrix::clear(), computeBMatrixAt(), computeConstitutiveMatrixAt(), oofem::Element::computeNumberOfDofs(), computeVolumeAround(), oofem::Element_remote, oofem::Element::giveCrossSection(), giveElement(), giveIntegrationElementLocalCodeNumbers(), oofem::Element::giveIntegrationRule(), oofem::Element::giveInterpolation(), oofem::Element::giveNumberOfIntegrationRules(), oofem::FEInterpolation::hasSubPatchFormulation(), oofem::StructuralCrossSection::isCharacteristicMtrxSymmetric(), oofem::FloatMatrix::plusProductSymmUpper(), oofem::FloatMatrix::plusProductUnsym(), oofem::FloatMatrix::resize(), oofem::FloatMatrix::symmetrized(), and oofem::FloatMatrix::zero().
Referenced by giveCharacteristicMatrix(), and giveMassMtrxIntegrationMask().
|
protected |
Optimized version, allowing to pass element displacements as parameter.
Standard version has a huge performance leak; in typical IGA element the element vector is VERY large and its querying for each point take more time than strain evaluation. And this has to be done for each integration point. This optimized version allows to assemble displacement vector only once (for all IP) and pass this vector as parameter
Definition at line 312 of file structuralelementevaluator.C.
References oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatArray::beProductOf(), computeBMatrixAt(), giveElement(), giveIntegrationElementLocalCodeNumbers(), oofem::GaussPoint::giveIntegrationRule(), oofem::GaussPoint::giveMaterialMode(), oofem::FloatMatrix::giveNumberOfColumns(), oofem::IntArray::giveSize(), oofem::StructuralMaterial::giveSizeOfVoigtSymVector(), isActivated(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by oofem::BsplinePlaneStressElement::drawScalar(), oofem::NURBSPlaneStressElement::drawScalar(), oofem::TSplinePlaneStressElement::drawScalar(), oofem::NURBSSpace3dElement::drawScalar(), giveInternalForcesVector(), isActivated(), and updateInternalState().
|
protectedpure virtual |
Computes the stress vector.
answer | Stress vector. |
strain | Strain vector. |
gp | Integration point for which stress is computed. |
tStep | Time step. |
Implemented in oofem::Space3dStructuralElementEvaluator, and oofem::PlaneStressStructuralElementEvaluator.
Referenced by giveInternalForcesVector(), isActivated(), and updateInternalState().
|
inlineprotected |
Definition at line 118 of file structuralelementevaluator.h.
References oofem::Element::computeVectorOf(), and giveElement().
Referenced by oofem::drawIGAPatchDeformedGeometry(), oofem::BsplinePlaneStressElement::drawScalar(), oofem::NURBSPlaneStressElement::drawScalar(), oofem::TSplinePlaneStressElement::drawScalar(), and oofem::NURBSSpace3dElement::drawScalar().
|
inlineprotected |
Definition at line 121 of file structuralelementevaluator.h.
References oofem::Element::computeVectorOf(), and giveElement().
|
inlineprotectedvirtual |
Reimplemented in oofem::Space3dStructuralElementEvaluator, and oofem::PlaneStressStructuralElementEvaluator.
Definition at line 116 of file structuralelementevaluator.h.
References giveInternalForcesVector().
Referenced by computeConsistentMassMatrix(), computeStiffnessMatrix(), and giveInternalForcesVector().
|
protectedvirtual |
Reimplemented in oofem::NURBSSpace3dElement, oofem::TSplinePlaneStressElement, oofem::NURBSPlaneStressElement, and oofem::BsplinePlaneStressElement.
Definition at line 128 of file structuralelementevaluator.C.
References oofem::__CharTypeToString(), computeConsistentMassMatrix(), computeLumpedMassMatrix(), computeStiffnessMatrix(), and OOFEM_ERROR.
Referenced by oofem::BsplinePlaneStressElement::giveCharacteristicMatrix(), oofem::NURBSPlaneStressElement::giveCharacteristicMatrix(), oofem::TSplinePlaneStressElement::giveCharacteristicMatrix(), oofem::NURBSSpace3dElement::giveCharacteristicMatrix(), and ~StructuralElementEvaluator().
|
protectedvirtual |
Reimplemented in oofem::NURBSSpace3dElement, oofem::TSplinePlaneStressElement, oofem::NURBSPlaneStressElement, and oofem::BsplinePlaneStressElement.
Definition at line 116 of file structuralelementevaluator.C.
References oofem::FloatArray::clear(), and giveInternalForcesVector().
Referenced by oofem::BsplinePlaneStressElement::giveCharacteristicVector(), oofem::NURBSPlaneStressElement::giveCharacteristicVector(), oofem::TSplinePlaneStressElement::giveCharacteristicVector(), oofem::NURBSSpace3dElement::giveCharacteristicVector(), and ~StructuralElementEvaluator().
|
protectedpure virtual |
Implemented in oofem::NURBSSpace3dElement, oofem::TSplinePlaneStressElement, oofem::NURBSPlaneStressElement, and oofem::BsplinePlaneStressElement.
Referenced by oofem::Space3dStructuralElementEvaluator::computeBMatrixAt(), computeConsistentMassMatrix(), oofem::PlaneStressStructuralElementEvaluator::computeConstitutiveMatrixAt(), oofem::Space3dStructuralElementEvaluator::computeConstitutiveMatrixAt(), computeLumpedMassMatrix(), oofem::Space3dStructuralElementEvaluator::computeNMatrixAt(), computeStiffnessMatrix(), computeStrainVector(), oofem::PlaneStressStructuralElementEvaluator::computeStressVector(), oofem::Space3dStructuralElementEvaluator::computeStressVector(), computeVectorOf(), oofem::PlaneStressStructuralElementEvaluator::computeVolumeAround(), oofem::Space3dStructuralElementEvaluator::computeVolumeAround(), giveInternalForcesVector(), giveMassMtrxIntegrationMask(), and updateInternalState().
|
protectedvirtual |
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.
Definition at line 81 of file structuralelementevaluator.C.
References oofem::IntArray::at(), oofem::IntArray::clear(), oofem::IntArray::followedBy(), oofem::Element::giveDofManDofIDMask(), oofem::Element::giveInterpolation(), oofem::IGAIntegrationElement::giveKnotSpan(), oofem::FEInterpolation::giveKnotSpanBasisFuncMask(), oofem::FEInterpolation::giveNsd(), oofem::IntArray::giveSize(), oofem::FEInterpolation::hasSubPatchFormulation(), and oofem::IntArray::resize().
Referenced by computeStiffnessMatrix(), computeStrainVector(), oofem::drawIGAPatchDeformedGeometry(), giveInternalForcesVector(), and isActivated().
|
protectedvirtual |
Definition at line 260 of file structuralelementevaluator.C.
References oofem::FloatArray::assemble(), oofem::FloatArray::clear(), computeBMatrixAt(), oofem::Element::computeNumberOfDofs(), computeStrainVector(), computeStressVector(), oofem::Element::computeVectorOf(), computeVolumeAround(), giveElement(), giveIntegrationElementLocalCodeNumbers(), oofem::Element::giveIntegrationRule(), oofem::Element::giveInterpolation(), oofem::Element::giveNumberOfIntegrationRules(), oofem::FloatArray::giveSize(), oofem::FEInterpolation::hasSubPatchFormulation(), isActivated(), oofem::FloatArray::plusProduct(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by computeVolumeAround(), and giveCharacteristicVector().
|
inlineprotectedvirtual |
Returns mask indicating, which unknowns (their type and ordering is the same as element unknown vector) participate in mass matrix integration.
Nonzero value at i-th position indicates that corresponding row in interpolation matrix N will participate in mass matrix integration (typically only displacements are taken into account).
answer | Integration mask, if zero sized, all unknowns participate. This is default. |
Definition at line 81 of file structuralelementevaluator.h.
References oofem::IntArray::clear(), computeBMatrixAt(), computeConsistentMassMatrix(), computeLumpedMassMatrix(), computeNMatrixAt(), computeStiffnessMatrix(), and giveElement().
Referenced by computeConsistentMassMatrix().
|
inlineprotectedvirtual |
Returns the integration rule for mass matrices, if relevant.
Definition at line 72 of file structuralelementevaluator.h.
Referenced by computeConsistentMassMatrix().
|
inlineprotected |
Definition at line 124 of file structuralelementevaluator.h.
References computeConstitutiveMatrixAt(), computeStrainVector(), computeStressVector(), drawIGAPatchDeformedGeometry, gc, giveIntegrationElementLocalCodeNumbers(), and updateInternalState().
Referenced by computeConsistentMassMatrix(), computeLumpedMassMatrix(), computeStrainVector(), and giveInternalForcesVector().
|
protected |
Definition at line 340 of file structuralelementevaluator.C.
References computeStrainVector(), computeStressVector(), oofem::Element::computeVectorOf(), oofem::Element_remote, giveElement(), oofem::Element::giveIntegrationRule(), oofem::Element::giveKnotSpanParallelMode(), oofem::Element::giveNumberOfIntegrationRules(), and oofem::FloatArray::subtract().
Referenced by isActivated(), oofem::BsplinePlaneStressElement::updateInternalState(), oofem::NURBSPlaneStressElement::updateInternalState(), oofem::TSplinePlaneStressElement::updateInternalState(), and oofem::NURBSSpace3dElement::updateInternalState().
|
friend |
Definition at line 1006 of file iga.C.
Referenced by oofem::BsplinePlaneStressElement::drawDeformedGeometry(), oofem::NURBSPlaneStressElement::drawDeformedGeometry(), oofem::NURBSSpace3dElement::drawDeformedGeometry(), and isActivated().
|
protected |
Definition at line 59 of file structuralelementevaluator.h.
Referenced by StructuralElementEvaluator().
|
protected |
Flag indicating if transformation matrix has been already computed.
Definition at line 61 of file structuralelementevaluator.h.