OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
oofem::StructuralElementEvaluator Class Referenceabstract

This class represent a new concept on how to define elements. More...

#include <structuralelementevaluator.h>

+ Inheritance diagram for oofem::StructuralElementEvaluator:
+ Collaboration diagram for oofem::StructuralElementEvaluator:

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 IntegrationRulegiveMassMtrxIntegrationRule ()
 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 ElementgiveElement ()=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)
 

Detailed Description

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

Todo:
{Class is missing much documentation.}

Definition at line 56 of file structuralelementevaluator.h.

Constructor & Destructor Documentation

oofem::StructuralElementEvaluator::StructuralElementEvaluator ( )
protected

Definition at line 51 of file structuralelementevaluator.C.

References oofem::FloatMatrix::clear(), and rotationMatrix.

virtual oofem::StructuralElementEvaluator::~StructuralElementEvaluator ( )
inlineprotectedvirtual

Member Function Documentation

virtual void oofem::StructuralElementEvaluator::computeBMatrixAt ( FloatMatrix answer,
GaussPoint gp 
)
protectedpure virtual
void oofem::StructuralElementEvaluator::computeConsistentMassMatrix ( FloatMatrix answer,
TimeStep tStep,
double &  mass 
)
protectedvirtual

Computes consistent mass matrix of receiver using numerical integration over element volume.

Mass matrix is computed as $ M = \int_V N^{\mathrm{T}} \rho N dV $, where $ N $ 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.

Parameters
answerMass matrix.
tStepTime step.
massTotal 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().

virtual void oofem::StructuralElementEvaluator::computeConstitutiveMatrixAt ( FloatMatrix answer,
MatResponseMode  rMode,
GaussPoint gp,
TimeStep tStep 
)
protectedpure virtual

Computes constitutive matrix of receiver.

Parameters
answerConstitutive matrix.
rModeMaterial response mode of answer.
gpIntegration point for which constitutive matrix is computed.
tStepTime step.

Implemented in oofem::Space3dStructuralElementEvaluator, and oofem::PlaneStressStructuralElementEvaluator.

Referenced by computeStiffnessMatrix(), and isActivated().

void oofem::StructuralElementEvaluator::computeLumpedMassMatrix ( FloatMatrix answer,
TimeStep tStep 
)
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.

Parameters
answerLumped mass matrix.
tStepTime 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().

virtual void oofem::StructuralElementEvaluator::computeNMatrixAt ( FloatMatrix answer,
GaussPoint gp 
)
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().

void oofem::StructuralElementEvaluator::computeStrainVector ( FloatArray answer,
GaussPoint gp,
TimeStep tStep,
FloatArray u 
)
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().

virtual void oofem::StructuralElementEvaluator::computeStressVector ( FloatArray answer,
const FloatArray strain,
GaussPoint gp,
TimeStep tStep 
)
protectedpure virtual

Computes the stress vector.

Parameters
answerStress vector.
strainStrain vector.
gpIntegration point for which stress is computed.
tStepTime step.

Implemented in oofem::Space3dStructuralElementEvaluator, and oofem::PlaneStressStructuralElementEvaluator.

Referenced by giveInternalForcesVector(), isActivated(), and updateInternalState().

void oofem::StructuralElementEvaluator::computeVectorOf ( PrimaryField field,
const IntArray dofIdMask,
ValueModeType  u,
TimeStep tStep,
FloatArray answer 
)
inlineprotected
virtual double oofem::StructuralElementEvaluator::computeVolumeAround ( GaussPoint gp)
inlineprotectedvirtual
int oofem::StructuralElementEvaluator::giveIntegrationElementLocalCodeNumbers ( IntArray answer,
Element elem,
IntegrationRule ie 
)
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.

Returns
Nonzero if integration rule code numbers differ from element code numbers

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

virtual void oofem::StructuralElementEvaluator::giveMassMtrxIntegrationMask ( IntArray answer)
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).

Parameters
answerIntegration 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().

virtual IntegrationRule* oofem::StructuralElementEvaluator::giveMassMtrxIntegrationRule ( )
inlineprotectedvirtual

Returns the integration rule for mass matrices, if relevant.

Returns
Number of integration points for mass matrix.

Definition at line 72 of file structuralelementevaluator.h.

Referenced by computeConsistentMassMatrix().

Friends And Related Function Documentation

Member Data Documentation

FloatMatrix oofem::StructuralElementEvaluator::rotationMatrix
protected

Definition at line 59 of file structuralelementevaluator.h.

Referenced by StructuralElementEvaluator().

int oofem::StructuralElementEvaluator::rotationMatrixDefined
protected

Flag indicating if transformation matrix has been already computed.

Definition at line 61 of file structuralelementevaluator.h.


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

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:41 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011