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

Abstract base class for all "structural" finite elements. More...

#include <structuralelement.h>

+ Inheritance diagram for oofem::StructuralElement:
+ Collaboration diagram for oofem::StructuralElement:

Public Member Functions

 StructuralElement (int n, Domain *d)
 Constructor. More...
 
virtual ~StructuralElement ()
 Destructor. More...
 
virtual void giveCharacteristicMatrix (FloatMatrix &answer, CharType, TimeStep *tStep)
 Computes characteristic matrix of receiver of requested type in given time step. More...
 
virtual void giveCharacteristicVector (FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
 Computes characteristic vector of receiver of requested type in given time step. More...
 
virtual void computeMassMatrix (FloatMatrix &answer, TimeStep *tStep)
 Computes mass matrix of receiver. More...
 
virtual void computeLumpedMassMatrix (FloatMatrix &answer, TimeStep *tStep)
 Computes lumped mass matrix of receiver. More...
 
virtual void computeConsistentMassMatrix (FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
 Computes consistent mass matrix of receiver using numerical integration over element volume. More...
 
virtual void giveMassMtrxIntegrationgMask (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 computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
 Computes numerically stiffness matrix of receiver. More...
 
void computeStiffnessMatrix_withIRulesAsSubcells (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
 
virtual void computeInitialStressMatrix (FloatMatrix &answer, TimeStep *tStep)
 Computes initial stress matrix for linear stability problem. More...
 
virtual void computeField (ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer)
 Computes the unknown vector interpolated at the specified local coordinates. More...
 
virtual void giveInternalForcesVector (FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
 Returns equivalent nodal forces vectors. More...
 
virtual void giveInternalForcesVector_withIRulesAsSubcells (FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
 
virtual void computeStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
 Compute strain vector of receiver evaluated at given integration point at given time step from element displacement vector. More...
 
virtual int giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
 Returns the integration point corresponding value in full form. More...
 
virtual void computeResultingIPTemperatureAt (FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
 Computes at given time (tStep) the the resulting temperature component array. More...
 
virtual void computeResultingIPEigenstrainAt (FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
 Computes at given time the resulting eigenstrain component array. More...
 
virtual int adaptiveUpdate (TimeStep *tStep)
 Updates the internal state variables stored in all IPs according to already mapped state. More...
 
virtual void updateInternalState (TimeStep *tStep)
 Updates element state after equilibrium in time step has been reached. More...
 
virtual void updateYourself (TimeStep *tStep)
 Updates element state after equilibrium in time step has been reached. More...
 
virtual int checkConsistency ()
 Performs consistency check. More...
 
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 const char * giveClassName () const
 
virtual int giveInternalStateAtNode (FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep)
 Returns internal state variable (like stress,strain) at node of element in Reduced form, the way how is obtained is dependent on InternalValueType. More...
 
virtual void showSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep)
 Shows sparse structure. More...
 
virtual void showExtendedSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep)
 Shows extended sparse structure (for example, due to nonlocal interactions for tangent stiffness) More...
 
virtual void computeLoadVector (FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
 Computes the contribution of the given body load (volumetric). More...
 
virtual void computeBoundarySurfaceLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true)
 Computes the contribution of the given load at the given boundary surface in global coordinate system. More...
 
virtual void computeBoundaryEdgeLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true)
 Computes the contribution of the given load at the given boundary edge. More...
 
virtual void computeEdgeNMatrix (FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
 computes edge interpolation matrix More...
 
virtual void computeSurfaceNMatrix (FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
 Computes surface interpolation matrix. More...
 
virtual void computeConstitutiveMatrixAt (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)=0
 Computes constitutive matrix of receiver. More...
 
StructuralCrossSectiongiveStructuralCrossSection ()
 Helper function which returns the structural cross-section for the element. More...
 
virtual void createMaterialStatus ()
 
virtual void computeStressVector (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)=0
 Computes the stress vector of receiver at given integration point, at time step tStep. More...
 
virtual void computeBmatrixAt (GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)=0
 Computes the geometrical matrix of receiver in given integration point. More...
 
virtual void computeNmatrixAt (const FloatArray &iLocCoord, FloatMatrix &answer)
 Computes interpolation matrix for element unknowns. More...
 
Methods related to nonlocal models
virtual void updateBeforeNonlocalAverage (TimeStep *tStep)
 Updates internal element state (in all integration points of receiver) before nonlocal averaging takes place. More...
 
virtual void giveNonlocalLocationArray (IntArray &locationArray, const UnknownNumberingScheme &us)
 Returns the "nonlocal" location array of receiver. More...
 
virtual void addNonlocalStiffnessContributions (SparseMtrx &dest, const UnknownNumberingScheme &s, TimeStep *tStep)
 Adds the "nonlocal" contribution to stiffness matrix, to account for nonlocality of material model. More...
 
- Public Member Functions inherited from oofem::Element
 Element (int n, Domain *aDomain)
 Constructor. More...
 
 Element (const Element &src)=delete
 
Elementoperator= (const Element &src)=delete
 
virtual ~Element ()
 Virtual destructor. More...
 
virtual void drawYourself (oofegGraphicContext &gc, TimeStep *tStep)
 
virtual void drawAnnotation (oofegGraphicContext &gc, TimeStep *tStep)
 
virtual void drawRawGeometry (oofegGraphicContext &gc, TimeStep *tStep)
 
virtual void drawDeformedGeometry (oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
 
virtual void drawScalar (oofegGraphicContext &gc, TimeStep *tStep)
 
virtual void drawSpecial (oofegGraphicContext &gc, TimeStep *tStep)
 
virtual void giveLocalIntVarMaxMin (oofegGraphicContext &gc, TimeStep *tStep, double &emin, double &emax)
 
virtual int giveInternalStateAtSide (FloatArray &answer, InternalStateType type, InternalStateMode mode, int side, TimeStep *tStep)
 Returns internal state variable (like stress,strain) at side of element in Reduced form If side is possessing DOFs, otherwise recover techniques will not work due to absence of side-shape functions. More...
 
int giveLabel () const
 
int giveGlobalNumber () const
 
void setGlobalNumber (int num)
 Sets receiver globally unique number. More...
 
elementParallelMode giveParallelMode () const
 Return elementParallelMode of receiver. More...
 
void setParallelMode (elementParallelMode _mode)
 Sets parallel mode of element. More...
 
virtual elementParallelMode giveKnotSpanParallelMode (int) const
 Returns the parallel mode for particular knot span of the receiver. More...
 
int packUnknowns (DataStream &buff, TimeStep *tStep)
 Pack all necessary data of element (according to its parallel_mode) integration points into given communication buffer. More...
 
int unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep)
 Unpack and updates all necessary data of element (according to its parallel_mode) integration points into given communication buffer. More...
 
int estimatePackSize (DataStream &buff)
 Estimates the necessary pack size to hold all packed data of receiver. More...
 
const IntArraygivePartitionList () const
 Returns partition list of receiver. More...
 
void setPartitionList (IntArray &pl)
 Sets partition list of receiver. More...
 
virtual double predictRelativeComputationalCost ()
 Returns the weight representing relative computational cost of receiver The reference element is triangular plane stress element with linear approximation, single integration point and linear isotropic material. More...
 
virtual double giveRelativeSelfComputationalCost ()
 Returns the weight representing relative computational cost of receiver The reference element is triangular plane stress element. More...
 
virtual double predictRelativeRedistributionCost ()
 Returns the relative redistribution cost of the receiver. More...
 
IntArraygiveBodyLoadArray ()
 Returns array containing load numbers of loads acting on element. More...
 
IntArraygiveBoundaryLoadArray ()
 Returns array containing load numbers of boundary loads acting on element. 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 printOutputAt (FILE *file, TimeStep *tStep)
 Prints output of receiver to stream, for given time step. More...
 
void giveLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
 Returns the location array (array of code numbers) of receiver for given numbering scheme. More...
 
void giveLocationArray (IntArray &locationArray, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
 
virtual void giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL)
 Returns the location array for the boundary of the element. More...
 
virtual void giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL)
 
virtual int giveNumberOfDofs ()
 
virtual int giveNumberOfInternalDofManagers () const
 
virtual DofManagergiveInternalDofManager (int i) const
 Returns i-th internal element dof manager of the receiver. More...
 
virtual double giveCharacteristicValue (CharType type, TimeStep *tStep)
 Computes characteristic value of receiver of requested type in given time step. More...
 
virtual void computeTangentFromSurfaceLoad (FloatMatrix &answer, SurfaceLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep)
 Computes the tangent contribution of the given load at the given boundary. More...
 
virtual void computeTangentFromEdgeLoad (FloatMatrix &answer, EdgeLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep)
 Computes the tangent contribution of the given load at the given boundary. More...
 
const IntArraygiveBodyLoadList () const
 Returns receiver list of bodyloads. More...
 
const IntArraygiveBoundaryLoadList () const
 Returns receiver list of boundary loads. More...
 
void computeVectorOf (ValueModeType u, TimeStep *tStep, FloatArray &answer)
 Returns local vector of unknowns. More...
 
void computeVectorOf (const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false)
 
void computeBoundaryVectorOf (const IntArray &bNodes, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false)
 Boundary version of computeVectorOf. More...
 
void computeVectorOf (PrimaryField &field, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false)
 Returns local vector of unknowns. More...
 
void computeVectorOfPrescribed (ValueModeType u, TimeStep *tStep, FloatArray &answer)
 Returns local vector of prescribed unknowns. More...
 
void computeVectorOfPrescribed (const IntArray &dofIDMask, ValueModeType type, TimeStep *tStep, FloatArray &answer)
 Returns local vector of prescribed unknowns. More...
 
virtual int computeNumberOfDofs ()
 Computes or simply returns total number of element's local DOFs. More...
 
virtual int computeNumberOfGlobalDofs ()
 Computes the total number of element's global dofs. More...
 
int computeNumberOfPrimaryMasterDofs ()
 Computes the total number of element's primary master DOFs. More...
 
virtual bool computeGtoLRotationMatrix (FloatMatrix &answer)
 Returns transformation matrix from global c.s. More...
 
virtual bool giveRotationMatrix (FloatMatrix &answer)
 Transformation matrices updates rotation matrix between element-local and primary DOFs, taking into account nodal c.s. More...
 
virtual bool computeDofTransformationMatrix (FloatMatrix &answer, const IntArray &nodes, bool includeInternal)
 Returns transformation matrix for DOFs from global coordinate system to local coordinate system in nodes. More...
 
virtual void giveDofManDofIDMask (int inode, IntArray &answer) const
 Returns dofmanager dof mask for node. More...
 
virtual void giveInternalDofManDofIDMask (int inode, IntArray &answer) const
 Returns internal dofmanager dof mask for node. More...
 
virtual void giveElementDofIDMask (IntArray &answer) const
 Returns element dof mask for node. More...
 
virtual double computeVolumeAround (GaussPoint *gp)
 Returns volume related to given integration point. More...
 
virtual double computeVolumeAreaOrLength ()
 Computes the volume, area or length of the element depending on its spatial dimension. More...
 
double computeMeanSize ()
 Computes the size of the element defined as its length. More...
 
virtual double computeVolume ()
 Computes the volume. More...
 
virtual double computeArea ()
 Computes the area (zero for all but 2d geometries). More...
 
virtual double computeLength ()
 Computes the length (zero for all but 1D geometries) More...
 
virtual void giveBoundaryEdgeNodes (IntArray &bNodes, int boundary)
 Returns list of receiver boundary nodes for given edge. More...
 
virtual void giveBoundarySurfaceNodes (IntArray &bNodes, int boundary)
 Returns list of receiver boundary nodes for given surface. More...
 
virtual IntegrationRulegiveBoundaryEdgeIntegrationRule (int order, int boundary)
 Returns boundary edge integration rule. More...
 
virtual IntegrationRulegiveBoundarySurfaceIntegrationRule (int order, int boundary)
 Returns boundary surface integration rule. More...
 
int giveDofManagerNumber (int i) const
 Translates local to global indices for dof managers. More...
 
const IntArraygiveDofManArray () const
 
void addDofManager (DofManager *dMan)
 
DofManagergiveDofManager (int i) const
 
NodegiveNode (int i) const
 Returns reference to the i-th node of element. More...
 
virtual ElementSidegiveSide (int i) const
 Returns reference to the i-th side of element. More...
 
virtual FEInterpolationgiveInterpolation () const
 
virtual FEInterpolationgiveInterpolation (DofIDItem id) const
 Returns the interpolation for the specific dof id. More...
 
virtual MaterialgiveMaterial ()
 
int giveMaterialNumber () const
 
CrossSectiongiveCrossSection ()
 
void setMaterial (int matIndx)
 Sets the material of receiver. More...
 
virtual void setCrossSection (int csIndx)
 Sets the cross section model of receiver. More...
 
virtual int giveNumberOfDofManagers () const
 
virtual int giveNumberOfNodes () const
 Returns number of nodes of receiver. More...
 
void setDofManagers (const IntArray &dmans)
 Sets receiver dofManagers. More...
 
void setBodyLoads (const IntArray &bodyLoads)
 Sets receiver bodyLoadArray. More...
 
void setIntegrationRules (std::vector< std::unique_ptr< IntegrationRule > > irlist)
 Sets integration rules. More...
 
virtual integrationDomain giveIntegrationDomain () const
 Returns integration domain for receiver, used to initialize integration point over receiver volume. More...
 
virtual MaterialMode giveMaterialMode ()
 Returns material mode for receiver integration points. More...
 
virtual int giveIntegrationRuleLocalCodeNumbers (IntArray &answer, IntegrationRule &ie)
 Assembles the 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...
 
int giveRegionNumber ()
 
virtual void postInitialize ()
 Performs post initialization steps. More...
 
virtual void initializeYourself (TimeStep *timeStepWhenICApply)
 Initialization according to state given by initial conditions. More...
 
virtual bool isActivated (TimeStep *tStep)
 
virtual bool isCast (TimeStep *tStep)
 
virtual void initForNewStep ()
 Initializes receivers state to new time step. More...
 
virtual Element_Geometry_Type giveGeometryType () const
 Returns the element geometry type. More...
 
virtual int giveSpatialDimension ()
 Returns the element spatial dimension (1, 2, or 3). More...
 
virtual int giveNumberOfBoundarySides ()
 
virtual int giveDefaultIntegrationRule () const
 Returns id of default integration rule. More...
 
virtual IntegrationRulegiveDefaultIntegrationRulePtr ()
 Access method for default integration rule. More...
 
int giveNumberOfIntegrationRules ()
 
virtual IntegrationRulegiveIntegrationRule (int i)
 
std::vector< std::unique_ptr< IntegrationRule > > & giveIntegrationRulesArray ()
 
virtual int testElementExtension (ElementExtension ext)
 Tests if the element implements required extension. More...
 
int giveGlobalIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
 
virtual double giveLengthInDir (const FloatArray &normalToCrackPlane)
 Default implementation returns length of element projection into specified direction. More...
 
virtual double giveCharacteristicLength (const FloatArray &normalToCrackPlane)
 Returns the size of element in the given direction, in some cases adjusted (e.g. More...
 
double giveCharacteristicLengthForPlaneElements (const FloatArray &normalToCrackPlane)
 Returns the size of element in the given direction if the direction is in the XY plane, otherwise gives the mean size defined as the square root of the element area. More...
 
double giveCharacteristicLengthForAxisymmElements (const FloatArray &normalToCrackPlane)
 Returns the size of an axisymmetric element in the given direction if the direction is in the XY plane, otherwise gives the mean distance vrom the symmetry axis multiplied by pi. More...
 
virtual double giveCharacteristicSize (GaussPoint *gp, FloatArray &normalToCrackPlane, ElementCharSizeMethod method)
 Returns characteristic element size for a given integration point and given direction. More...
 
virtual double giveParentElSize () const
 Returns the size (length, area or volume depending on element type) of the parent element. More...
 
virtual int computeGlobalCoordinates (FloatArray &answer, const FloatArray &lcoords)
 Computes the global coordinates from given element's local coordinates. More...
 
virtual bool computeLocalCoordinates (FloatArray &answer, const FloatArray &gcoords)
 Computes the element local coordinates from given global coordinates. More...
 
virtual int giveLocalCoordinateSystem (FloatMatrix &answer)
 Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy. More...
 
virtual void computeMidPlaneNormal (FloatArray &answer, const GaussPoint *gp)
 Computes mid-plane normal of receiver at integration point. More...
 
virtual int adaptiveMap (Domain *oldd, TimeStep *tStep)
 Initializes the internal state variables stored in all IPs according to state in given domain. More...
 
virtual int mapStateVariables (Domain &iOldDom, const TimeStep &iTStep)
 Maps the internal state variables stored in all IPs from the old domain to the new domain. More...
 
virtual int adaptiveFinish (TimeStep *tStep)
 Finishes the mapping for given time step. More...
 
virtual void updateLocalNumbering (EntityRenumberingFunctor &f)
 Local renumbering support. More...
 
template<class T >
void ipEvaluator (T *src, void(T::*f)(GaussPoint *gp))
 Integration point evaluator, loops over receiver IP's and calls given function (passed as f parameter) on them. The IP is parameter to function f. More...
 
template<class T , class S >
void ipEvaluator (T *src, void(T::*f)(GaussPoint *, S &), S &_val)
 Integration point evaluator, loops over receiver IP's and calls given function (passed as f parameter) on them. The IP is parameter to function f as well as additional array. 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...
 
virtual const char * giveInputRecordName () const =0
 
DomaingiveDomain () const
 
virtual void setDomain (Domain *d)
 Sets associated Domain. More...
 
int giveNumber () const
 
void setNumber (int num)
 Sets number of receiver. More...
 
virtual void printYourself ()
 Prints receiver state on stdout. Useful for debugging. More...
 
virtual InterfacegiveInterface (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 Member Functions

virtual void computeBodyLoadVectorAt (FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
 Computes the load vector due to body load acting on receiver, at given time step. More...
 
virtual int giveNumberOfIPForMassMtrxIntegration ()
 Return desired number of integration points for consistent mass matrix computation, if required. More...
 
void condense (FloatMatrix *stiff, FloatMatrix *mass, FloatArray *load, IntArray *what)
 General service for condensation of stiffness and optionally load vector and mass or initial stress matrices of receiver. More...
 
virtual void setupIRForMassMtrxIntegration (IntegrationRule &iRule)
 Setup Integration Rule Gauss Points for Mass Matrix integration. More...
 
Edge and surface load support services.
virtual void computePointLoadVectorAt (FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, bool global=true)
 Computes point load vector contribution of receiver for given load (should has BoundaryLoad Base). More...
 
virtual void giveEdgeDofMapping (IntArray &answer, int iEdge) const
 Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element DOFs. More...
 
virtual void giveSurfaceDofMapping (IntArray &answer, int iSurf) const
 Assembles surface dof mapping mask, which provides mapping between surface local DOFs and "global" element DOFs. More...
 
virtual IntegrationRuleGetSurfaceIntegrationRule (int order)
 
virtual double computeEdgeVolumeAround (GaussPoint *gp, int iEdge)
 Computes volume related to integration point on local edge. More...
 
virtual double computeSurfaceVolumeAround (GaussPoint *gp, int iSurf)
 Computes volume related to integration point on local surface. More...
 
virtual int computeLoadGToLRotationMtrx (FloatMatrix &answer)
 Returns transformation matrix from global coordinate system to local element coordinate system for element load vector components. More...
 
virtual int computeLoadLEToLRotationMatrix (FloatMatrix &answer, int iEdge, GaussPoint *gp)
 Returns transformation matrix from local edge c.s to element local coordinate system of load vector components. More...
 
virtual int computeLoadLSToLRotationMatrix (FloatMatrix &answer, int iSurf, GaussPoint *gp)
 Returns transformation matrix from local surface c.s to element local coordinate system of load vector components. More...
 
- Protected Member Functions inherited from oofem::Element
virtual void computeGaussPoints ()
 Initializes the array of integration rules member variable. More...
 

Protected Attributes

std::unique_ptr< FloatArrayinitialDisplacements
 Initial displacement vector, describes the initial nodal displacements when element has been casted. More...
 
- Protected Attributes inherited from oofem::Element
int numberOfDofMans
 Number of dofmanagers. More...
 
IntArray dofManArray
 Array containing dofmanager numbers. More...
 
int material
 Number of associated material. More...
 
int crossSection
 Number of associated cross section. More...
 
IntArray bodyLoadArray
 Array containing indexes of loads (body loads and boundary loads are kept separately), that apply on receiver. More...
 
IntArray boundaryLoadArray
 
std::vector< std::unique_ptr< IntegrationRule > > integrationRulesArray
 List of integration rules of receiver (each integration rule contains associated integration points also). More...
 
FloatMatrix elemLocalCS
 Transformation material matrix, used in orthotropic and anisotropic materials, global->local transformation. More...
 
int activityTimeFunction
 Element activity time function. If defined, nonzero value indicates active receiver, zero value inactive element. More...
 
int globalNumber
 In parallel mode, globalNumber contains globally unique DoFManager number. More...
 
int numberOfGaussPoints
 Number of integration points as specified by nip. More...
 
elementParallelMode parallel_mode
 Determines the parallel mode of the element. More...
 
IntArray partitions
 List of partition sharing the shared element or remote partition containing remote element counterpart. More...
 
- Protected Attributes inherited from oofem::FEMComponent
int number
 Component number. More...
 
Domaindomain
 Link to domain object, useful for communicating with other FEM components. More...
 

Friends

class IDNLMaterial
 
class TrabBoneNL3D
 
class MisesMatNl
 
class RankineMatNl
 
class GradDpElement
 

Detailed Description

Abstract base class for all "structural" finite elements.

It declares common interface provided by all derived elements. The implementation of these services is partly left on derived classes, some general services are implemented here generally (But they can be overload by more efficient element implementation). The general implementation provided here is intended for both linear and nonlinear computations. At this level, only material model nonlinearities are taken into account. If particular element type will participate in geometrically nonlinear computations, it should be derived from derived NLStructuralElement class, which provide support for this cases.

The basic data of an element are the numbers of its 'numberOfNodes' nodes, stored in 'nodeArray', of its 'material', of its body loads (eg, the dead weight) stored in 'loadArray'. These data are obtained from the domain. The element contains logical reference to Volume object ( which stores material,geometrical characteristics, integration points, state of material and so on) The calculated data of an element are its 'massMatrix', its 'stiffnessMatrix' its 'locationArray'. Since the load vector is recalculated at every time step, it is not given the status of attribute.

Tasks:

  • Obtaining its basic data, the element reads in the data file the number of these objects, then obtains the data from the domain (methods 'giveNode', 'giveMaterial',etc).
  • Calculating its contribution to the problem : Calculating its mass matrix M, its stiffness matrix K, its load vector f, its location array. Calculating its contribution to the LHS and RHS of the linear system, using Static,Newmark,etc, formula. These contributions are usually combinations of M,K,f.
  • Performing end-of-step operations; Calculating the strains and stresses at its Gauss points.
  • Printing its output in the data file and updating itself.

Definition at line 95 of file structuralelement.h.

Constructor & Destructor Documentation

oofem::StructuralElement::StructuralElement ( int  n,
Domain d 
)

Constructor.

Creates structural element with given number, belonging to given domain.

Parameters
nElement number.
dDomain to which new material will belong.

Definition at line 67 of file structuralelement.C.

oofem::StructuralElement::~StructuralElement ( )
virtual

Member Function Documentation

int oofem::StructuralElement::adaptiveUpdate ( TimeStep tStep)
virtual

Updates the internal state variables stored in all IPs according to already mapped state.

Parameters
tStepTime step.
Returns
Nonzero if o.k, otherwise zero.

Reimplemented from oofem::Element.

Definition at line 1202 of file structuralelement.C.

References computeStrainVector(), giveStructuralCrossSection(), oofem::Element::integrationRulesArray, and oofem::MaterialModelMapperInterfaceType.

Referenced by computeInitialStressMatrix().

void oofem::StructuralElement::addNonlocalStiffnessContributions ( SparseMtrx dest,
const UnknownNumberingScheme s,
TimeStep tStep 
)
virtual

Adds the "nonlocal" contribution to stiffness matrix, to account for nonlocality of material model.

Typically, this contribution is obtained by summing up mutual IP contributions.

Todo:
Take into account cross section model (slaves)

Definition at line 1180 of file structuralelement.C.

References oofem::Element::giveDefaultIntegrationRulePtr(), giveStructuralCrossSection(), oofem::Element::isActivated(), and oofem::NonlocalMaterialStiffnessInterfaceType.

Referenced by computeInitialStressMatrix().

int oofem::StructuralElement::checkConsistency ( )
virtual

Performs consistency check.

This method is called at startup for all elements in particular domain. This method is intended to check data compatibility. Particular element types should test if compatible material and crossSection both with required capabilities are specified. Derived classes should provide their own analysis specific tests. Some printed input if incompatibility is found should be provided (error or warning member functions). Method can be also used to initialize some variables, since this is invoked after all domain components are instanciated.

Returns
Zero value if check fail, otherwise nonzero.

Reimplemented from oofem::Element.

Reimplemented in oofem::NLStructuralElement, oofem::Shell7BaseXFEM, oofem::SpringElement, oofem::SolidShell, oofem::LineDistributedSpring, oofem::NodalSpringElement, oofem::LumpedMassElement, oofem::Shell7Base, and oofem::TrPlaneStress2dXFEM.

Definition at line 1007 of file structuralelement.C.

References oofem::CS_StructuralCapability, giveClassName(), oofem::Element::giveCrossSection(), and OOFEM_WARNING.

Referenced by oofem::LumpedMassElement::checkConsistency(), and computeInitialStressMatrix().

virtual void oofem::StructuralElement::computeBmatrixAt ( GaussPoint gp,
FloatMatrix answer,
int  lowerIndx = 1,
int  upperIndx = ALL_STRAINS 
)
pure virtual

Computes the geometrical matrix of receiver in given integration point.

The product of this matrix (assembled at given integration point) and element displacement vector is element strain vector. If lowerIndx and upperIndx parameters are specified, answer is formed only for strains within this interval. This will affects the size of answer.

Parameters
gpIntegration point for which answer is computed.
answerGeometric matrix of receiver.
lowerIndxIf specified, answer is formed only for strain with index equal and greater than lowerIndx. This parameter has default value 1 (answer is formed from first strain).
upperIndxIf specified, answer is formed only for strain with index less and equal than upperIndx. This parameter has default value ALL_STRAINS (answer is formed for all strains).

Implemented in oofem::Shell7Base, oofem::AbaqusUserElement, oofem::Beam3d, oofem::Beam2d, oofem::AxisymElement, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::Quad1MindlinShell3D, oofem::LIBeam3d2, oofem::PlaneStrainElement, oofem::Lattice2d, oofem::RerShell, oofem::SpringElement, oofem::LIBeam3dNL, oofem::LIBeam3dNL2, oofem::MITC4Shell, oofem::PlaneStressElement, oofem::Truss1d, oofem::Truss2d, oofem::NodalSpringElement, oofem::tet21ghostsolid, oofem::Truss3d, oofem::Quad1Mindlin, oofem::LIBeam3d, oofem::LIBeam2dNL, oofem::LumpedMassElement, oofem::PlaneStress2d, oofem::Quad1PlaneStrain, oofem::LIBeam2d, oofem::LineDistributedSpring, oofem::InterfaceElem1d, oofem::DKTPlate, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::Structural2DElement, oofem::InterfaceElem2dQuad, oofem::InterfaceElement3dTrLin, oofem::InterfaceElem2dLin, oofem::L4Axisymm, oofem::QDKTPlate, oofem::CCTPlate, oofem::Structural3DElement, oofem::QTruss1d, oofem::TrPlaneStress2dXFEM, oofem::QTrPlaneStress2dXFEM, oofem::Quad2PlateSubSoil, oofem::Q4Axisymm, oofem::Tr_Warp, oofem::CohesiveSurface3d, oofem::PlaneStress2dXfem, oofem::TrPlanestressRotAllman, oofem::HTSelement, oofem::TrPlaneStrRot, oofem::SolidShell, and oofem::LSpaceBB.

Referenced by oofem::PhaseFieldElement::computeBStress_u(), oofem::NLStructuralElement::computeInitialStressMatrix(), computeLoadLSToLRotationMatrix(), oofem::GradDpElement::computeLocalStrainVector(), oofem::NLStructuralElement::computeStiffnessMatrix(), computeStiffnessMatrix(), oofem::GradDpElement::computeStiffnessMatrix_ku(), oofem::PhaseFieldElement::computeStiffnessMatrix_ud(), oofem::GradDpElement::computeStiffnessMatrix_uk(), oofem::GradDpElement::computeStiffnessMatrix_uu(), oofem::PhaseFieldElement::computeStiffnessMatrix_uu(), oofem::NLStructuralElement::computeStiffnessMatrix_withIRulesAsSubcells(), computeStiffnessMatrix_withIRulesAsSubcells(), computeStrainVector(), oofem::GradDpElement::giveInternalForcesVector(), oofem::NLStructuralElement::giveInternalForcesVector(), giveInternalForcesVector(), oofem::PhaseFieldElement::giveInternalForcesVector_u(), oofem::NLStructuralElement::giveInternalForcesVector_withIRulesAsSubcells(), giveInternalForcesVector_withIRulesAsSubcells(), oofem::GradDpElement::giveLocalInternalForcesVector(), oofem::TrabBoneNL3D::giveLocalNonlocalStiffnessContribution(), oofem::MisesMatNl::giveLocalNonlocalStiffnessContribution(), oofem::RankineMatNl::giveLocalNonlocalStiffnessContribution(), oofem::IDNLMaterial::giveLocalNonlocalStiffnessContribution(), oofem::TrabBoneNL3D::giveRemoteNonlocalStiffnessContribution(), oofem::MisesMatNl::giveRemoteNonlocalStiffnessContribution(), oofem::RankineMatNl::giveRemoteNonlocalStiffnessContribution(), oofem::IDNLMaterial::giveRemoteNonlocalStiffnessContribution(), and oofem::LinearizedDilationForceAssembler::vectorFromElement().

void oofem::StructuralElement::computeBodyLoadVectorAt ( FloatArray answer,
Load load,
TimeStep tStep,
ValueModeType  mode 
)
protectedvirtual

Computes the load vector due to body load acting on receiver, at given time step.

Default implementation computes body load vector numerically as $ l=\int_V N^{\mathrm{T}} f \rho\;\mathrm{d}V $ using default integration rule. Result is transformed to global c.s.

Parameters
answerComputed load vector due to body load
loadBody load which contribution is computed.
tStepTime step.
modedetermines the response mode

Reimplemented in oofem::Beam3d, oofem::Shell7Base, oofem::Beam2d, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::Quad1MindlinShell3D, oofem::LIBeam3d2, oofem::RerShell, oofem::LIBeam3dNL, oofem::LIBeam3dNL2, oofem::Quad1Mindlin, oofem::LIBeam3d, oofem::LIBeam2dNL, oofem::LIBeam2d, oofem::LineDistributedSpring, oofem::DKTPlate, oofem::QDKTPlate, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::CCTPlate3d, oofem::DKTPlate3d, oofem::TrPlanestressRotAllman3d, oofem::TrPlaneStrRot3d, oofem::CCTPlate, and oofem::TrPlaneStrRot.

Definition at line 259 of file structuralelement.C.

References oofem::FloatArray::add(), oofem::FloatArray::beTProductOf(), oofem::BodyLoadBGT, oofem::FloatArray::clear(), oofem::Load::computeComponentArrayAt(), computeLoadGToLRotationMtrx(), computeNmatrixAt(), oofem::Element::computeVolumeAround(), oofem::ForceLoadBVT, oofem::CrossSection::give(), oofem::GeneralBoundaryCondition::giveBCGeoType(), oofem::GeneralBoundaryCondition::giveBCValType(), oofem::Element::giveCrossSection(), oofem::Element::giveDefaultIntegrationRulePtr(), oofem::FloatArray::giveSize(), oofem::GaussPoint::giveSubPatchCoordinates(), OOFEM_ERROR, and oofem::FloatArray::rotatedWith().

Referenced by oofem::LIBeam2d::computeBodyLoadVectorAt(), oofem::LIBeam2dNL::computeBodyLoadVectorAt(), oofem::LIBeam3d::computeBodyLoadVectorAt(), oofem::LIBeam3dNL2::computeBodyLoadVectorAt(), oofem::LIBeam3dNL::computeBodyLoadVectorAt(), oofem::LIBeam3d2::computeBodyLoadVectorAt(), oofem::Beam2d::computeBodyLoadVectorAt(), oofem::Beam3d::computeBodyLoadVectorAt(), computeLoadVector(), oofem::BeamBaseElement::computeLocalForceLoadVector(), and giveClassName().

void oofem::StructuralElement::computeBoundaryEdgeLoadVector ( FloatArray answer,
BoundaryLoad load,
int  edge,
CharType  type,
ValueModeType  mode,
TimeStep tStep,
bool  global = true 
)
virtual

Computes the contribution of the given load at the given boundary edge.

In general, the answer should include only relevant DOFs at the edge. The related is giveBoundaryLocationArray method, which should return corresponding code numbers..

Note
Elements which do not have an contribution should resize the vector to be empty.
Parameters
answerRequested contribution of load (in Global c.s.).
loadLoad to compute contribution from.
edgeEdge number.
typeType of the contribution.
modeDetermines mode of answer.
tStepTime step when answer is computed.
globalif true (default) then contribution is in global c.s., when false then contribution is in element local c.s.
Todo:
Make sure this part is correct.
Todo:
Support this...
Todo:
Change the load type to "BoundaryEdgeLoad" maybe?

Reimplemented from oofem::Element.

Reimplemented in oofem::Beam3d, oofem::Beam2d, oofem::Shell7BaseXFEM, oofem::TrPlanestressRotAllman, and oofem::Shell7Base.

Definition at line 183 of file structuralelement.C.

References oofem::FEInterpolation::boundaryEdgeLocal2Global(), oofem::FloatArray::clear(), computeEdgeNMatrix(), computeEdgeVolumeAround(), computeLoadGToLRotationMtrx(), computeLoadLEToLRotationMatrix(), oofem::BoundaryLoad::computeValueAt(), oofem::Load::CST_Global, oofem::Load::FT_Entity, oofem::BoundaryLoad::giveApproxOrder(), oofem::Element::giveBoundaryEdgeIntegrationRule(), oofem::BoundaryLoad::giveCoordSystMode(), oofem::Load::giveFormulationType(), oofem::Element::giveInterpolation(), oofem::GaussPoint::giveNaturalCoordinates(), OOFEM_ERROR, oofem::FloatArray::plusProduct(), and oofem::FloatArray::rotatedWith().

Referenced by oofem::BeamBaseElement::computeLocalForceLoadVector(), and giveClassName().

void oofem::StructuralElement::computeBoundarySurfaceLoadVector ( FloatArray answer,
BoundaryLoad load,
int  boundary,
CharType  type,
ValueModeType  mode,
TimeStep tStep,
bool  global = true 
)
virtual

Computes the contribution of the given load at the given boundary surface in global coordinate system.

In general, the answer should include only relevant DOFs at the edge. The related is giveBoundaryLocationArray method, which should return corresponding code numbers.

Note
Elements which do not have an contribution should resize the vector to be empty.
Parameters
answerRequested contribution of load.
loadLoad to compute contribution from.
boundaryBoundary number.
typeType of the contribution.
modeDetermines mode of answer.
tStepTime step when answer is computed.
globalif true (default) then contribution is in global c.s., when false then contribution is in element local c.s.
Todo:
Make sure this part is correct.
Todo:
Support this...
Todo:
Some way to ask for the thickness at a global coordinate maybe?

Reimplemented from oofem::Element.

Reimplemented in oofem::tet21ghostsolid.

Definition at line 108 of file structuralelement.C.

References oofem::FEInterpolation::boundaryLocal2Global(), oofem::FloatArray::clear(), computeLoadGToLRotationMtrx(), computeSurfaceNMatrix(), computeSurfaceVolumeAround(), oofem::BoundaryLoad::computeValueAt(), oofem::Load::CST_Global, oofem::Load::FT_Entity, oofem::BoundaryLoad::giveApproxOrder(), oofem::Element::giveBoundarySurfaceIntegrationRule(), oofem::BoundaryLoad::giveCoordSystMode(), oofem::Load::giveFormulationType(), oofem::Element::giveInterpolation(), oofem::GaussPoint::giveNaturalCoordinates(), OOFEM_ERROR, oofem::FloatArray::plusProduct(), and oofem::FloatArray::rotatedWith().

Referenced by oofem::BeamBaseElement::computeLocalForceLoadVector(), and giveClassName().

void oofem::StructuralElement::computeConsistentMassMatrix ( FloatMatrix answer,
TimeStep tStep,
double &  mass,
const double *  ipDensity = NULL 
)
virtual

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.

Reimplemented in oofem::AbaqusUserElement, oofem::Beam3d, oofem::TrPlaneStress2dXFEM, oofem::QTrPlaneStress2dXFEM, oofem::PlaneStress2dXfem, oofem::Beam2d, oofem::LIBeam3dNL2, and oofem::LIBeam3dNL.

Definition at line 334 of file structuralelement.C.

References oofem::IntArray::at(), oofem::FloatMatrix::at(), computeNmatrixAt(), oofem::Element::computeNumberOfDofs(), oofem::Element::computeVolumeAround(), oofem::CrossSection::give(), oofem::Element::giveCrossSection(), giveMassMtrxIntegrationgMask(), oofem::FloatMatrix::giveNumberOfRows(), oofem::GaussPoint::giveSubPatchCoordinates(), oofem::Element::isActivated(), oofem::IntArray::isEmpty(), oofem::FloatMatrix::plusProductSymmUpper(), oofem::FloatMatrix::resize(), setupIRForMassMtrxIntegration(), oofem::FloatMatrix::symmetrized(), and oofem::FloatMatrix::zero().

Referenced by computeLumpedMassMatrix(), computeMassMatrix(), and oofem::LSPrimaryVariableMapper::mapPrimaryVariables().

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

Computes constitutive matrix of receiver.

Default implementation uses element cross section giveCharMaterialStiffnessMatrix service.

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

Implemented in oofem::AbaqusUserElement, oofem::Beam3d, oofem::Beam2d, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::AxisymElement, oofem::Quad1MindlinShell3D, oofem::PlaneStrainElement, oofem::Lattice2d, oofem::RerShell, oofem::LIBeam3dNL2, oofem::SpringElement, oofem::MITC4Shell, oofem::PlaneStressElement, oofem::NodalSpringElement, oofem::Quad1Mindlin, oofem::LIBeam3d, oofem::LIBeam2d, oofem::LineDistributedSpring, oofem::DKTPlate, oofem::CohesiveSurface3d, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::Truss3d, oofem::LumpedMassElement, oofem::LIBeam3d2, oofem::InterfaceElem1d, oofem::QDKTPlate, oofem::LIBeam3dNL, oofem::CCTPlate, oofem::InterfaceElem2dQuad, oofem::InterfaceElement3dTrLin, oofem::Tr2Shell7, oofem::Truss1d, oofem::Truss2d, oofem::InterfaceElem2dLin, oofem::HTSelement, oofem::QTrPlaneStress2dXFEM, oofem::TrPlaneStress2dXFEM, oofem::Structural3DElement, oofem::Tr2Shell7XFEM, oofem::PlaneStress2dXfem, oofem::tet21ghostsolid, oofem::Tr_Warp, oofem::TrPlaneStrRot, oofem::LIBeam2dNL, oofem::QTruss1dGrad, and oofem::QTruss1d.

Referenced by oofem::NLStructuralElement::computeStiffnessMatrix(), computeStiffnessMatrix(), oofem::NLStructuralElement::computeStiffnessMatrix_withIRulesAsSubcells(), computeStiffnessMatrix_withIRulesAsSubcells(), giveClassName(), oofem::PlaneStressElement::giveMaterialMode(), oofem::PlaneStrainElement::giveMaterialMode(), oofem::AxisymElement::giveMaterialMode(), oofem::LinearizedDilationForceAssembler::vectorFromElement(), oofem::ZZErrorEstimatorInterface::ZZErrorEstimatorI_computeElementContributions(), and ~StructuralElement().

void oofem::StructuralElement::computeEdgeNMatrix ( FloatMatrix answer,
int  boundaryID,
const FloatArray lcoords 
)
virtual
void oofem::StructuralElement::computeField ( ValueModeType  mode,
TimeStep tStep,
const FloatArray lcoords,
FloatArray answer 
)
virtual

Computes the unknown vector interpolated at the specified local coordinates.

Used for exporting data and mapping fields.

See also
giveElementDofIDMask The unknown vector should match the element field as specified by the element dof IDs.
Parameters
modeMode (total, increment, etc) of the output
tStepTime step to evaluate at
lcoordsLocal coordinates to evaluate at
answerResults

Reimplemented from oofem::Element.

Reimplemented in oofem::AbaqusUserElement, oofem::TrPlaneStress2dXFEM, oofem::QTrPlaneStress2dXFEM, oofem::MacroLSpace, and oofem::QTruss1dGrad.

Definition at line 573 of file structuralelement.C.

References oofem::FloatArray::beProductOf(), computeNmatrixAt(), and oofem::Element::computeVectorOf().

Referenced by computeInitialStressMatrix().

virtual void oofem::StructuralElement::computeInitialStressMatrix ( FloatMatrix answer,
TimeStep tStep 
)
inlinevirtual

Computes initial stress matrix for linear stability problem.

Default implementation is not provided. Please note, that initial stress matrix depends on normal forces of element, corresponding engineering model must take this into account.

Parameters
answerComputed initial stress matrix.
tStepTime step.

Reimplemented in oofem::NLStructuralElement, oofem::Beam3d, oofem::SpringElement, oofem::Beam2d, oofem::NodalSpringElement, oofem::LumpedMassElement, and oofem::LIBeam2dNL.

Definition at line 197 of file structuralelement.h.

References adaptiveUpdate(), addNonlocalStiffnessContributions(), checkConsistency(), computeField(), computeResultingIPEigenstrainAt(), computeResultingIPTemperatureAt(), computeStrainVector(), giveInputRecord(), giveInternalForcesVector(), giveInternalForcesVector_withIRulesAsSubcells(), giveIPValue(), giveNonlocalLocationArray(), initializeFrom(), OOFEM_ERROR, updateBeforeNonlocalAverage(), updateInternalState(), and updateYourself().

Referenced by giveCharacteristicMatrix().

virtual int oofem::StructuralElement::computeLoadGToLRotationMtrx ( FloatMatrix answer)
inlineprotectedvirtual

Returns transformation matrix from global coordinate system to local element coordinate system for element load vector components.

If no transformation is necessary, answer is empty matrix (default);

Parameters
answerTransformation matrix.
Returns
Nonzero if transformation matrix is not empty matrix, zero otherwise.

Reimplemented in oofem::Beam3d, oofem::LIBeam3d2, oofem::LIBeam3dNL, oofem::LIBeam3dNL2, oofem::MITC4Shell, oofem::LIBeam3d, oofem::LIBeam2dNL, oofem::LIBeam2d, oofem::LinQuad3DPlaneStress, oofem::CCTPlate3d, oofem::DKTPlate3d, oofem::TrPlanestressRotAllman3d, and oofem::TrPlaneStrRot3d.

Definition at line 421 of file structuralelement.h.

References oofem::FloatMatrix::clear().

Referenced by oofem::QDKTPlate::computeBodyLoadVectorAt(), computeBodyLoadVectorAt(), oofem::TrPlanestressRotAllman::computeBoundaryEdgeLoadVector(), oofem::Beam2d::computeBoundaryEdgeLoadVector(), computeBoundaryEdgeLoadVector(), oofem::tet21ghostsolid::computeBoundarySurfaceLoadVector(), computeBoundarySurfaceLoadVector(), and computePointLoadVectorAt().

virtual int oofem::StructuralElement::computeLoadLEToLRotationMatrix ( FloatMatrix answer,
int  iEdge,
GaussPoint gp 
)
inlineprotectedvirtual

Returns transformation matrix from local edge c.s to element local coordinate system of load vector components.

Necessary, because integration must be done in local coordinate system of entity (edge or surface). If no transformation is necessary, answer is empty matrix (default);

Parameters
answerComputed rotation matrix.
iEdgeEdge number.
gpIntegration point (point, where transformation is computed, useful for curved edges).
Returns
Nonzero if transformation matrix is not empty matrix, zero otherwise.

Reimplemented in oofem::Quad1MindlinShell3D, oofem::LIBeam3d2, oofem::DKTPlate, oofem::QDKTPlate, oofem::Truss2d, oofem::Quad1Mindlin, oofem::LIBeam3dNL, oofem::LIBeam3dNL2, oofem::MITC4Shell, oofem::Truss3d, oofem::CCTPlate, oofem::LIBeam3d, oofem::Structural2DElement, oofem::LIBeam2dNL, oofem::DKTPlate3d, oofem::LIBeam2d, and oofem::Structural3DElement.

Definition at line 438 of file structuralelement.h.

References oofem::FloatMatrix::clear().

Referenced by computeBoundaryEdgeLoadVector().

virtual int oofem::StructuralElement::computeLoadLSToLRotationMatrix ( FloatMatrix answer,
int  iSurf,
GaussPoint gp 
)
inlineprotectedvirtual

Returns transformation matrix from local surface c.s to element local coordinate system of load vector components.

Necessary, because integration must be done in local coordinate system of entity (edge or surface). If no transformation is necessary, answer is empty matrix (default);

Parameters
answerComputed rotation matrix.
iSurfSurface number.
gpIntegration point (point, where transformation is computed, useful for curved surfaces)
Returns
Nonzero if transformation matrix is not empty matrix, zero otherwise.

Reimplemented in oofem::DKTPlate, oofem::LSpace, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::DKTPlate3d, oofem::QDKTPlate, oofem::CCTPlate3d, oofem::Structural3DElement, oofem::TrPlaneStrRot3d, oofem::QSpace, and oofem::Q27Space.

Definition at line 452 of file structuralelement.h.

References ALL_STRAINS, oofem::FloatMatrix::clear(), computeBmatrixAt(), computeNmatrixAt(), computeStressVector(), condense(), giveNumberOfIPForMassMtrxIntegration(), and setupIRForMassMtrxIntegration().

void oofem::StructuralElement::computeLoadVector ( FloatArray answer,
BodyLoad load,
CharType  type,
ValueModeType  mode,
TimeStep tStep 
)
virtual

Computes the contribution of the given body load (volumetric).

Parameters
answerRequested contribution of load.
loadLoad to compute contribution from.
typeType of the contribution.
modeDetermines mode of answer.
tStepTime step when answer is computed.
Todo:
This assumption of dead-weight loads needs to be lifted. We can have other body loads, such as

Reimplemented from oofem::Element.

Reimplemented in oofem::tet21ghostsolid.

Definition at line 89 of file structuralelement.C.

References oofem::FloatArray::clear(), computeBodyLoadVectorAt(), oofem::Element::computeLocalCoordinates(), computePointLoadVectorAt(), and oofem::PointLoad::giveCoordinates().

Referenced by giveClassName().

void oofem::StructuralElement::computeLumpedMassMatrix ( FloatMatrix answer,
TimeStep tStep 
)
virtual

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 zeroes all off-diagonal members and zeroes 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.

Reimplemented in oofem::Shell7Base, oofem::DKTPlate, oofem::CCTPlate, oofem::Quad1MindlinShell3D, oofem::Quad1Mindlin, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::SpringElement, oofem::LineDistributedSpring, oofem::LIBeam3d2, oofem::LIBeam3dNL2, oofem::NodalSpringElement, oofem::LIBeam3dNL, oofem::Truss2d, oofem::RerShell, oofem::LTRSpace, oofem::Truss1d, oofem::Truss3d, oofem::LumpedMassElement, oofem::InterfaceElem1d, oofem::LIBeam2d, oofem::LIBeam3d, and oofem::LIBeam2dNL.

Definition at line 410 of file structuralelement.C.

References oofem::IntArray::at(), oofem::FloatMatrix::at(), computeConsistentMassMatrix(), oofem::Element::computeNumberOfDofs(), oofem::Element::giveDofManDofIDMask(), oofem::Element::giveNumberOfDofManagers(), oofem::FloatMatrix::giveNumberOfRows(), oofem::IntArray::giveSize(), oofem::Element::isActivated(), OOFEM_ERROR, oofem::FloatMatrix::resize(), oofem::FloatMatrix::times(), and oofem::FloatMatrix::zero().

Referenced by giveCharacteristicMatrix(), and giveCharacteristicVector().

void oofem::StructuralElement::computeMassMatrix ( FloatMatrix answer,
TimeStep tStep 
)
virtual

Computes mass matrix of receiver.

Default implementation returns consistent mass matrix and uses numerical integration. Returns result of this->computeConsistentMassMatrix service, transformed into nodal coordinate system. Requires the computeNmatrixAt and giveMassMtrxIntegrationgMask services to be implemented.

Parameters
answerMass matrix.
tStepTime step.

Reimplemented in oofem::Shell7Base, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::DKTPlate, oofem::CCTPlate, oofem::Quad1MindlinShell3D, oofem::Quad1Mindlin, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::SpringElement, oofem::LineDistributedSpring, oofem::LIBeam3d2, oofem::NodalSpringElement, oofem::Truss2d, oofem::RerShell, oofem::Truss1d, oofem::Truss3d, oofem::LumpedMassElement, oofem::InterfaceElem1d, oofem::LIBeam2d, oofem::LIBeam3d, and oofem::LIBeam2dNL.

Definition at line 401 of file structuralelement.C.

References computeConsistentMassMatrix().

Referenced by giveCharacteristicMatrix().

void oofem::StructuralElement::computeNmatrixAt ( const FloatArray iLocCoord,
FloatMatrix answer 
)
virtual

Computes interpolation matrix for element unknowns.

The order and meaning of unknowns is element dependent.

Parameters
iLocCoordLocal coordinates.
answerInterpolation matrix evaluated at gp.

Reimplemented in oofem::Shell7Base, oofem::Beam3d, oofem::Beam2d, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::LIBeam3d2, oofem::SpringElement, oofem::RerShell, oofem::LIBeam3dNL, oofem::LIBeam3dNL2, oofem::Truss1d, oofem::Truss2d, oofem::MITC4Shell, oofem::NodalSpringElement, oofem::LatticeStructuralElement, oofem::Truss3d, oofem::LIBeam2dNL, oofem::LIBeam3d, oofem::LumpedMassElement, oofem::LIBeam2d, oofem::InterfaceElem1d, oofem::DKTPlate, oofem::InterfaceElem2dQuad, oofem::InterfaceElement3dTrLin, oofem::InterfaceElem2dLin, oofem::QDKTPlate, oofem::CCTPlate, oofem::CohesiveSurface3d, oofem::TrPlaneStress2dXFEM, oofem::QTrPlaneStress2dXFEM, oofem::Tr_Warp, oofem::TrPlanestressRotAllman, oofem::PlaneStress2dXfem, oofem::HTSelement, and oofem::TrPlaneStrRot.

Definition at line 1024 of file structuralelement.C.

References oofem::FloatMatrix::beNMatrixOf(), oofem::FEInterpolation::evalN(), oofem::Element::giveInterpolation(), oofem::Element::giveNumberOfDofManagers(), oofem::Element::giveSpatialDimension(), N, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().

Referenced by computeBodyLoadVectorAt(), computeConsistentMassMatrix(), computeField(), computeLoadLSToLRotationMatrix(), computePointLoadVectorAt(), oofem::FreeWarping::computeResultAtCenterOfGravity(), oofem::Quad1PlateSubSoil::computeSurfaceNMatrixAt(), oofem::Tria1PlateSubSoil::computeSurfaceNMatrixAt(), giveIPValue(), oofem::Quad1PlaneStrain::HuertaErrorEstimatorI_computeNmatrixAt(), oofem::PlaneStress2d::HuertaErrorEstimatorI_computeNmatrixAt(), oofem::TrPlaneStrain::HuertaErrorEstimatorI_computeNmatrixAt(), oofem::LSpace::HuertaErrorEstimatorI_computeNmatrixAt(), oofem::LTRSpace::HuertaErrorEstimatorI_computeNmatrixAt(), oofem::TrPlaneStress2d::HuertaErrorEstimatorI_computeNmatrixAt(), oofem::LSPrimaryVariableMapper::mapPrimaryVariables(), and oofem::XfemStructuralElementInterface::XfemElementInterface_computeConsistentMassMatrix().

void oofem::StructuralElement::computePointLoadVectorAt ( FloatArray answer,
Load load,
TimeStep tStep,
ValueModeType  mode,
bool  global = true 
)
protectedvirtual

Computes point load vector contribution of receiver for given load (should has BoundaryLoad Base).

Parameters
answerComputed load vector.
loadPoint load which contribution is computed.
tStepTime step.
modedetermines response mode

Definition at line 298 of file structuralelement.C.

References oofem::FloatArray::beTProductOf(), computeLoadGToLRotationMtrx(), oofem::Element::computeLocalCoordinates(), computeNmatrixAt(), oofem::PointLoad::computeValueAt(), oofem::Load::CST_Global, oofem::PointLoad::giveCoordinates(), oofem::PointLoad::giveCoordSystMode(), OOFEM_WARNING, and oofem::FloatArray::rotatedWith().

Referenced by computeLoadVector(), oofem::BeamBaseElement::computeLocalForceLoadVector(), and giveClassName().

void oofem::StructuralElement::computeResultingIPTemperatureAt ( FloatArray answer,
TimeStep tStep,
GaussPoint gp,
ValueModeType  mode 
)
virtual
void oofem::StructuralElement::computeStiffnessMatrix ( FloatMatrix answer,
MatResponseMode  rMode,
TimeStep tStep 
)
virtual

Computes numerically stiffness matrix of receiver.

Default implementation computes element stiffness using $ K=\int_v B^{\mathrm{T}} D B \mathrm{d}V $ formulae, where $ B $ is element geometric matrix and $ D $ is material stiffness matrix. No geometrical nonlinearity is taken into account. NUmerical integration procedure uses integrationRulesArray for numerical integration. Support for reduced or selected integration is implemented. The individual integration rules are assumed to correspond to different terms from which the overall matrix is assembled.

If there is one integration rule, the full integration of all coefficients is performed. Otherwise, integration is performed using following rules. Each integration rule can specify start and end strain index of strain vector components for which is valid. It is necessary to ensure that these start and end indexes, dividing geometrical matrix into blocks, are not overlapping and that each strain component is included.

Then stiffness matrix is obtained as summation of integrals $ I_{ij}=\int_v B^{\mathrm{T}}_i D_{ij} B_j \mathrm{d}V $ where $ B_i $ is i-th block of geometrical matrix and $ D_{ij} $ is corresponding constitutive sub-matrix. The geometrical matrix is obtained using computeBmatrixAt service and the constitutive matrix is obtained using computeConstitutiveMatrixAt service. The $ I_{ij} $ integral is evaluated using such integration rule, which is valid for i-th or j-th block and has smaller number of integration points.

For higher numerical performance, only one half of stiffness matrix is computed and answer is then symmetrized. Therefore, if element matrix will be generally nonsymmetric, one must specialize this method. Finally, the result is transformed into global coordinate system (or nodal coordinate system, if it is defined).

Parameters
answerComputed stiffness matrix (symmetric).
rModeResponse mode.
tStepTime step.

Reimplemented in oofem::Shell7Base, oofem::AbaqusUserElement, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::NLStructuralElement, oofem::Lattice2d, oofem::MITC4Shell, oofem::Beam3d, oofem::LIBeam3dNL, oofem::Quad1MindlinShell3D, oofem::LIBeam3dNL2, oofem::TrPlanestressRotAllman, oofem::Shell7BaseXFEM, oofem::SpringElement, oofem::Beam2d, oofem::LIBeam3d2, oofem::QTrPlaneStress2dXFEM, oofem::TrPlaneStress2dXFEM, oofem::NodalSpringElement, oofem::PlaneStress2dXfem, oofem::QTRSpaceGrad, oofem::QWedgeGrad, oofem::tet21ghostsolid, oofem::SolidShell, oofem::QTruss1dGrad, oofem::HTSelement, oofem::MacroLSpace, oofem::LumpedMassElement, oofem::LIBeam3d, oofem::PlaneStressPhF2d, oofem::QPlaneStressPhF2d, oofem::LIBeam2d, oofem::QPlaneStrainGrad, oofem::QTrPlaneStressGrad, oofem::CoupledFieldsElement, oofem::QPlaneStressGrad, and oofem::QTrPlaneStrainGrad.

Definition at line 583 of file structuralelement.C.

References oofem::FloatMatrix::beProductOf(), oofem::FloatMatrix::beSubMatrixOf(), oofem::FloatMatrix::clear(), computeBmatrixAt(), computeConstitutiveMatrixAt(), oofem::Element::computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::Element::giveDefaultIntegrationRulePtr(), oofem::Element::integrationRulesArray, oofem::Element::isActivated(), oofem::CrossSection::isCharacteristicMtrxSymmetric(), oofem::FloatMatrix::plusProductSymmUpper(), oofem::FloatMatrix::plusProductUnsym(), and oofem::FloatMatrix::symmetrized().

Referenced by oofem::LIBeam2d::computeStiffnessMatrix(), oofem::LIBeam3d::computeStiffnessMatrix(), oofem::LIBeam3d2::computeStiffnessMatrix(), giveCharacteristicMatrix(), oofem::LineDistributedSpring::giveInternalForcesVector(), and giveMassMtrxIntegrationgMask().

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

Computes the stress vector of receiver at given integration point, at time step tStep.

The nature of these stresses depends on the element's type.

Parameters
answerStress vector.
strainStrain vector.
gpIntegration point.
tStepTime step.

Implemented in oofem::AbaqusUserElement, oofem::Beam3d, oofem::Beam2d, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::AxisymElement, oofem::Quad1MindlinShell3D, oofem::PlaneStrainElement, oofem::Lattice2d, oofem::RerShell, oofem::LIBeam3dNL2, oofem::MITC4Shell, oofem::SpringElement, oofem::PlaneStressElement, oofem::LIBeam3d, oofem::NodalSpringElement, oofem::Quad1Mindlin, oofem::LIBeam2d, oofem::LineDistributedSpring, oofem::DKTPlate, oofem::CohesiveSurface3d, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::Truss3d, oofem::LIBeam3d2, oofem::LumpedMassElement, oofem::InterfaceElem1d, oofem::QDKTPlate, oofem::LIBeam3dNL, oofem::CCTPlate, oofem::InterfaceElem2dQuad, oofem::InterfaceElement3dTrLin, oofem::Truss1d, oofem::Truss2d, oofem::InterfaceElem2dLin, oofem::QTrPlaneStress2dXFEM, oofem::TrPlaneStress2dXFEM, oofem::Tr2Shell7, oofem::HTSelement, oofem::Structural3DElement, oofem::PlaneStress2dXfem, oofem::Tr2Shell7XFEM, oofem::tet21ghostsolid, oofem::Tr_Warp, oofem::TrPlaneStrRot, oofem::LIBeam2dNL, and oofem::QTruss1d.

Referenced by oofem::PhaseFieldElement::computeBStress_u(), computeLoadLSToLRotationMatrix(), computeStrainVector(), oofem::NLStructuralElement::giveInternalForcesVector(), giveInternalForcesVector(), oofem::NLStructuralElement::giveInternalForcesVector_withIRulesAsSubcells(), giveInternalForcesVector_withIRulesAsSubcells(), oofem::PlaneStressElement::giveMaterialMode(), oofem::PlaneStrainElement::giveMaterialMode(), oofem::AxisymElement::giveMaterialMode(), and updateInternalState().

void oofem::StructuralElement::computeSurfaceNMatrix ( FloatMatrix answer,
int  boundaryID,
const FloatArray lcoords 
)
virtual

Computes surface interpolation matrix.

Interpolation matrix provide way, how to compute local surface unknowns (nonzero element unknowns on surface) at any integration point of surface, based on local unknowns in surface nodes. Local coordinate system of surface edge and element surface numbering is element dependent. The integration point is specified using two-dimensional iso coordinates, or using area coordinates for triangular surface.

Parameters
answerInterpolation matrix of surface.
boundaryIDSurface number.
localcoordinates

Reimplemented in oofem::MITC4Shell, oofem::CCTPlate, and oofem::QDKTPlate.

Definition at line 175 of file structuralelement.C.

References oofem::FloatMatrix::beNMatrixOf(), oofem::FEInterpolation::boundarySurfaceEvalN(), and oofem::Element::giveInterpolation().

Referenced by computeBoundarySurfaceLoadVector(), and giveClassName().

double oofem::StructuralElement::computeSurfaceVolumeAround ( GaussPoint gp,
int  iSurf 
)
protectedvirtual
void oofem::StructuralElement::condense ( FloatMatrix stiff,
FloatMatrix mass,
FloatArray load,
IntArray what 
)
protected

General service for condensation of stiffness and optionally load vector and mass or initial stress matrices of receiver.

Parameters
stiffStiffness matrix to be condensed. Must be specified.
massMass or initial stress matrix. If parameter is NULL, only stiffness and/or load is condensed.
loadLoad vector of receiver. If specified then it is condensed. If no load vector condensation is necessary set parameter to NULL pointer.
whatInteger array. If at i-th position is nonzero, then i-th component is condensed.

Definition at line 1039 of file structuralelement.C.

References oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::giveNumberOfRows(), oofem::IntArray::giveSize(), oofem::FloatArray::giveSize(), oofem::FloatMatrix::isSquare(), OOFEM_ERROR, and oofem::FloatArray::resize().

Referenced by computeLoadLSToLRotationMatrix().

void oofem::StructuralElement::createMaterialStatus ( )
virtual
virtual IntegrationRule* oofem::StructuralElement::GetSurfaceIntegrationRule ( int  order)
inlineprotectedvirtual
void oofem::StructuralElement::giveCharacteristicMatrix ( FloatMatrix answer,
CharType  type,
TimeStep tStep 
)
virtual

Computes characteristic matrix of receiver of requested type in given time step.

Parameters
answerRequested characteristic matrix (stiffness, tangent, ...). If element has no capability to compute requested type of characteristic matrix error function is invoked.
typeId of characteristic component requested.
tStepTime step when answer is computed.

Reimplemented from oofem::Element.

Reimplemented in oofem::TR_SHELL01, and oofem::TR_SHELL02.

Definition at line 884 of file structuralelement.C.

References oofem::__CharTypeToString(), computeInitialStressMatrix(), computeLumpedMassMatrix(), computeMassMatrix(), computeStiffnessMatrix(), and OOFEM_ERROR.

void oofem::StructuralElement::giveCharacteristicVector ( FloatArray answer,
CharType  type,
ValueModeType  mode,
TimeStep tStep 
)
virtual

Computes characteristic vector of receiver of requested type in given time step.

If element has no capability to compute requested type of characteristic vector error function is invoked.

Parameters
answerRequested characteristic vector.
typeId of characteristic component requested.
modeDetermines mode of answer.
tStepTime step when answer is computed.

Reimplemented from oofem::Element.

Reimplemented in oofem::TR_SHELL01, oofem::TR_SHELL02, and oofem::Tr_Warp.

Definition at line 910 of file structuralelement.C.

References oofem::__CharTypeToString(), computeLumpedMassMatrix(), giveInternalForcesVector(), oofem::FloatMatrix::giveNumberOfColumns(), OOFEM_ERROR, and oofem::FloatArray::resize().

Referenced by oofem::Tr_Warp::giveCharacteristicVector().

virtual const char* oofem::StructuralElement::giveClassName ( ) const
inlinevirtual
Returns
Class name of the receiver.

Reimplemented from oofem::Element.

Reimplemented in oofem::NLStructuralElement, oofem::AbaqusUserElement, oofem::Shell7BaseXFEM, oofem::Beam3d, oofem::MITC4Shell, oofem::Beam2d, oofem::DKTPlate, oofem::QDKTPlate, oofem::RerShell, oofem::DKTPlate3d, oofem::SpringElement, oofem::CCTPlate3d, oofem::LIBeam3d2, oofem::CCTPlate, oofem::Lattice2d, oofem::TrPlaneStrRot3d, oofem::Quad1MindlinShell3D, oofem::Tr2Shell7, oofem::NodalSpringElement, oofem::Truss2d, oofem::LinQuad3DPlaneStress, oofem::TrPlaneStrain, oofem::Truss1d, oofem::Tr2Shell7XFEM, oofem::LSpace, oofem::Truss3d, oofem::Quad1PlaneStrain, oofem::PlaneStress2d, oofem::LIBeam3dNL, oofem::LumpedMassElement, oofem::TrPlanestressRotAllman3d, oofem::Axisymm3d, oofem::tet21ghostsolid, oofem::LIBeam3dNL2, oofem::InterfaceElem1d, oofem::TrPlaneStress2d, oofem::QTRSpace, oofem::LIBeam3d, oofem::TrPlanestressRotAllman, oofem::TR_SHELL01, oofem::LIBeam2d, oofem::TR_SHELL02, oofem::Quad1Mindlin, oofem::LIBeam2dNL, oofem::InterfaceElem2dQuad, oofem::Shell7Base, oofem::QSpace, oofem::QWedge, oofem::CohesiveSurface3d, oofem::SolidShell, oofem::InterfaceElem2dLin, oofem::TrPlaneStrRot, oofem::LWedge, oofem::InterfaceElement3dTrLin, oofem::QTrPlaneStrain, oofem::LTRSpace, oofem::Q27Space, oofem::BasicElement, oofem::LineDistributedSpring, oofem::QTrPlaneStress2d, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::QTruss1d, oofem::TrPlaneStress2dXFEM, oofem::Q4Axisymm, oofem::QTrPlaneStress2dXFEM, oofem::BasicElementQuad, oofem::Quad2PlateSubSoil, oofem::L4Axisymm, oofem::MacroLSpace, oofem::HTSelement, oofem::Tr_Warp, oofem::QSpaceGrad, oofem::PlaneStressPhF2d, oofem::PlaneStress2dXfem, oofem::QPlaneStress2d, oofem::QTRSpaceGrad, oofem::QWedgeGrad, oofem::Q9PlaneStress2d, oofem::QPlaneStressPhF2d, oofem::QPlaneStrain, oofem::QTruss1dGrad, oofem::QPlaneStrainGrad, oofem::QPlaneStressGrad, oofem::QTrPlaneStressGrad, oofem::QTrPlaneStrainGrad, and oofem::LSpaceBB.

Definition at line 295 of file structuralelement.h.

References computeBodyLoadVectorAt(), computeBoundaryEdgeLoadVector(), computeBoundarySurfaceLoadVector(), computeConstitutiveMatrixAt(), computeEdgeNMatrix(), computeLoadVector(), computePointLoadVectorAt(), computeSurfaceNMatrix(), createMaterialStatus(), gc, giveInternalStateAtNode(), giveStructuralCrossSection(), showExtendedSparseMtrxStructure(), and showSparseMtrxStructure().

Referenced by checkConsistency().

virtual void oofem::StructuralElement::giveEdgeDofMapping ( IntArray answer,
int  iEdge 
) const
inlineprotectedvirtual

Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element DOFs.

Mask can be imagined as local edge code numbers used to localize local edge DOFs to element DOFs.

Parameters
answerEdge DOF mapping.
iEdgeEdge number.

Reimplemented in oofem::MITC4Shell, oofem::Quad1MindlinShell3D, oofem::Shell7Base, oofem::LIBeam3d2, oofem::DKTPlate, oofem::QDKTPlate, oofem::Truss2d, oofem::Quad1Mindlin, oofem::LIBeam3dNL, oofem::LIBeam3dNL2, oofem::TrPlanestressRotAllman, oofem::Truss3d, oofem::LIBeam3d, oofem::Structural2DElement, oofem::LIBeam2dNL, oofem::CCTPlate, oofem::DKTPlate3d, oofem::LIBeam2d, oofem::Structural3DElement, oofem::Tr2Shell7, oofem::Tr2Shell7XFEM, and oofem::Tr_Warp.

Definition at line 388 of file structuralelement.h.

References oofem::IntArray::clear().

void oofem::StructuralElement::giveInputRecord ( DynamicInputRecord input)
virtual

Setups the input record string of receiver.

Parameters
inputDynamic input record to be filled by receiver.

TODO: Should initialDisplacements be stored? /ES

Reimplemented from oofem::Element.

Reimplemented in oofem::NLStructuralElement, oofem::AbaqusUserElement, oofem::TrPlaneStress2dXFEM, oofem::PlaneStress2dXfem, and oofem::QTrPlaneStress2dXFEM.

Definition at line 1229 of file structuralelement.C.

References oofem::Element::giveInputRecord().

Referenced by computeInitialStressMatrix(), oofem::AbaqusUserElement::giveInputRecord(), and oofem::NLStructuralElement::giveInputRecord().

void oofem::StructuralElement::giveInternalForcesVector ( FloatArray answer,
TimeStep tStep,
int  useUpdatedGpRecord = 0 
)
virtual

Returns equivalent nodal forces vectors.

Useful for nonlinear analysis. Default implementation computes result as $ F=\int_v B^{\mathrm{T}} \sigma \mathrm{d}V $, where $ \sigma $ is the real element stress vector obtained using computeStressVector service (if useUpdatedGpRecord=0) or (if useUpdatedGpRecord=1) from integration point status. The geometric matrix is obtained using computeBmatrixAt service. Integration is performed using default integration rule, which should produce always valid results, assuming that strains used for computation of stresses are valid.

Parameters
answerInternal nodal forces vector.
tStepTime step.
useUpdatedGpRecordIf equal to zero, the stresses in integration points are computed (slow but safe), else if nonzero the stresses are taken directly from integration point status (should be derived from StructuralMaterialStatus) (fast, but engineering model must ensure valid status data in each integration point).

Reimplemented in oofem::Shell7Base, oofem::NLStructuralElement, oofem::AbaqusUserElement, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::Beam3d, oofem::LIBeam3dNL, oofem::MITC4Shell, oofem::Quad1MindlinShell3D, oofem::LIBeam3dNL2, oofem::SpringElement, oofem::QTrPlaneStress2dXFEM, oofem::TrPlaneStress2dXFEM, oofem::Beam2d, oofem::HTSelement, oofem::LineDistributedSpring, oofem::NodalSpringElement, oofem::PlaneStress2dXfem, oofem::MacroLSpace, oofem::tet21ghostsolid, oofem::Shell7BaseXFEM, oofem::QTRSpaceGrad, oofem::QWedgeGrad, oofem::LumpedMassElement, oofem::QTruss1dGrad, oofem::SolidShell, oofem::PlaneStressPhF2d, oofem::QPlaneStressPhF2d, oofem::QPlaneStrainGrad, oofem::QTrPlaneStressGrad, oofem::CoupledFieldsElement, oofem::QPlaneStressGrad, and oofem::QTrPlaneStrainGrad.

Definition at line 732 of file structuralelement.C.

References oofem::FloatArray::beProductOf(), oofem::FloatArray::clear(), computeBmatrixAt(), computeStressVector(), oofem::Element::computeVectorOf(), oofem::Element::computeVolumeAround(), oofem::Element::giveDefaultIntegrationRulePtr(), oofem::GaussPoint::giveMaterialMode(), oofem::GaussPoint::giveMaterialStatus(), oofem::StructuralMaterial::giveReducedSymVectorForm(), oofem::FloatArray::giveSize(), oofem::StructuralMaterial::giveSizeOfVoigtSymVector(), oofem::StructuralMaterialStatus::giveStressVector(), initialDisplacements, oofem::Element::isActivated(), oofem::FloatArray::plusProduct(), oofem::FloatArray::resize(), oofem::FloatArray::subtract(), and oofem::FloatArray::zero().

Referenced by computeInitialStressMatrix(), giveCharacteristicVector(), oofem::Beam2d::giveInternalForcesVector(), and oofem::Beam3d::giveInternalForcesVector().

int oofem::StructuralElement::giveInternalStateAtNode ( FloatArray answer,
InternalStateType  type,
InternalStateMode  mode,
int  node,
TimeStep tStep 
)
virtual

Returns internal state variable (like stress,strain) at node of element in Reduced form, the way how is obtained is dependent on InternalValueType.

The value may be local, or smoothed using some recovery technique / returns zero if element is unable to respond to request.

Parameters
answerContains result, zero sized if not supported.
typeDetermines the internal variable requested (physical meaning).
modeDetermines the mode of variable (recovered, local, ...).
nodeNode number, for which variable is required.
tStepTime step.
Returns
Nonzero if o.k, zero otherwise.

Reimplemented from oofem::Element.

Definition at line 1257 of file structuralelement.C.

References oofem::FloatArray::at(), oofem::Node::giveCoordinate(), oofem::Element::giveInternalStateAtNode(), oofem::Element::giveNode(), oofem::Node::giveUpdatedCoordinate(), and oofem::FloatArray::resize().

Referenced by oofem::QTrPlaneStress2d::drawScalar(), oofem::QTrPlaneStrain::drawScalar(), oofem::QPlaneStress2d::drawScalar(), oofem::Axisymm3d::drawScalar(), oofem::L4Axisymm::drawScalar(), oofem::TrPlaneStress2d::drawScalar(), oofem::LTRSpace::drawScalar(), oofem::Quad1PlaneStrain::drawScalar(), oofem::PlaneStress2d::drawScalar(), oofem::Truss1d::drawScalar(), oofem::TrPlaneStrain::drawScalar(), oofem::LSpace::drawScalar(), oofem::TR_SHELL02::drawScalar(), oofem::TR_SHELL01::drawScalar(), oofem::CCTPlate::drawScalar(), oofem::QDKTPlate::drawScalar(), oofem::DKTPlate::drawScalar(), and giveClassName().

int oofem::StructuralElement::giveIPValue ( FloatArray answer,
GaussPoint gp,
InternalStateType  type,
TimeStep tStep 
)
virtual

Returns the integration point corresponding value in full form.

Parameters
answerContain corresponding integration point value, zero sized if not available.
gpIntegration point to check.
typeDetermines the type of internal variable.
tStepTime step.
Returns
Nonzero if o.k, zero otherwise.
Todo:
Which "error type" should be used? Why are there several? I don't see the point of this enum when there could be different function calls just as well (and different IST values)

Reimplemented from oofem::Element.

Reimplemented in oofem::Shell7Base, oofem::AbaqusUserElement, oofem::Beam3d, oofem::DKTPlate, oofem::Beam2d, oofem::QDKTPlate, oofem::MITC4Shell, oofem::CCTPlate, oofem::DKTPlate3d, oofem::Quad1MindlinShell3D, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::Quad1Mindlin, oofem::RerShell, oofem::TrPlanestressRotAllman3d, oofem::LinQuad3DPlaneStress, oofem::CCTPlate3d, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::LineDistributedSpring, oofem::tet21ghostsolid, oofem::TrPlaneStrRot, oofem::LIBeam3dNL, oofem::TrPlaneStrRot3d, oofem::LIBeam2dNL, oofem::LIBeam3d, and oofem::LIBeam2d.

Definition at line 1134 of file structuralelement.C.

References oofem::FloatArray::beProductOf(), computeNmatrixAt(), oofem::Element::computeVectorOf(), oofem::Element::giveIPValue(), oofem::GaussPoint::giveSubPatchCoordinates(), and N.

Referenced by oofem::NLStructuralElement::computeInitialStressMatrix(), computeInitialStressMatrix(), oofem::QTrPlaneStress2d::drawScalar(), oofem::QPlaneStrain::drawScalar(), oofem::QPlaneStress2d::drawScalar(), oofem::InterfaceElem2dLin::drawScalar(), oofem::InterfaceElem2dQuad::drawScalar(), oofem::Axisymm3d::drawScalar(), oofem::TrPlaneStress2d::drawScalar(), oofem::L4Axisymm::drawScalar(), oofem::LTRSpace::drawScalar(), oofem::InterfaceElem1d::drawScalar(), oofem::Quad1PlaneStrain::drawScalar(), oofem::PlaneStress2d::drawScalar(), oofem::CohesiveSurface3d::drawScalar(), oofem::PlaneStress2dXfem::drawScalar(), oofem::Truss1d::drawScalar(), oofem::TrPlaneStrain::drawScalar(), oofem::TrPlaneStress2dXFEM::drawScalar(), oofem::LIBeam3d2::drawScalar(), oofem::TrPlaneStress2d::drawSpecial(), oofem::LTRSpace::drawSpecial(), oofem::Quad1PlaneStrain::drawSpecial(), oofem::PlaneStress2d::drawSpecial(), oofem::LSpace::drawSpecial(), oofem::Lattice2d::drawSpecial(), oofem::LIBeam2d::giveIPValue(), oofem::LIBeam3d::giveIPValue(), oofem::LIBeam2dNL::giveIPValue(), oofem::LIBeam3dNL::giveIPValue(), oofem::TrPlaneStrRot::giveIPValue(), oofem::tet21ghostsolid::giveIPValue(), oofem::LineDistributedSpring::giveIPValue(), oofem::Quad1PlateSubSoil::giveIPValue(), oofem::Tria1PlateSubSoil::giveIPValue(), oofem::CCTPlate3d::giveIPValue(), oofem::LinQuad3DPlaneStress::giveIPValue(), oofem::TrPlanestressRotAllman3d::giveIPValue(), oofem::Quad1Mindlin::giveIPValue(), oofem::RerShell::giveIPValue(), oofem::TR_SHELL02::giveIPValue(), oofem::TR_SHELL01::giveIPValue(), oofem::Quad1MindlinShell3D::giveIPValue(), oofem::DKTPlate3d::giveIPValue(), oofem::CCTPlate::giveIPValue(), oofem::MITC4Shell::giveIPValue(), oofem::QDKTPlate::giveIPValue(), oofem::Beam2d::giveIPValue(), oofem::DKTPlate::giveIPValue(), oofem::Beam3d::giveIPValue(), oofem::Q9PlaneStress2d::NodalAveragingRecoveryMI_computeNodalValue(), oofem::QPlaneStress2d::NodalAveragingRecoveryMI_computeNodalValue(), oofem::TrPlaneStrain::NodalAveragingRecoveryMI_computeNodalValue(), oofem::LSpace::NodalAveragingRecoveryMI_computeNodalValue(), oofem::Axisymm3d::NodalAveragingRecoveryMI_computeNodalValue(), oofem::LTRSpace::NodalAveragingRecoveryMI_computeNodalValue(), oofem::TrPlaneStress2d::NodalAveragingRecoveryMI_computeNodalValue(), and oofem::Truss1d::NodalAveragingRecoveryMI_computeNodalValue().

virtual void oofem::StructuralElement::giveMassMtrxIntegrationgMask ( IntArray answer)
inlinevirtual

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 155 of file structuralelement.h.

References oofem::IntArray::clear(), computeStiffnessMatrix(), and computeStiffnessMatrix_withIRulesAsSubcells().

Referenced by computeConsistentMassMatrix(), and oofem::XfemStructuralElementInterface::XfemElementInterface_computeConsistentMassMatrix().

void oofem::StructuralElement::giveNonlocalLocationArray ( IntArray locationArray,
const UnknownNumberingScheme us 
)
virtual

Returns the "nonlocal" location array of receiver.

This is necessary, when stiffness matrix of nonlocal model is assembled. Since model is nonlocal, the value at given IP depends on other IP (generally belonging to different elements) and as a consequence leads to increase of stiffness matrix profile, to take into account this "remote" dependency.

Parameters
locationArrayLocation arrays from neighboring elements.
usUnknown numbering scheme.

Definition at line 1149 of file structuralelement.C.

References oofem::IntArray::clear(), oofem::IntArray::followedBy(), oofem::Element::giveDefaultIntegrationRulePtr(), giveStructuralCrossSection(), and oofem::NonlocalMaterialStiffnessInterfaceType.

Referenced by computeInitialStressMatrix().

int oofem::StructuralElement::giveNumberOfIPForMassMtrxIntegration ( )
protectedvirtual
StructuralCrossSection * oofem::StructuralElement::giveStructuralCrossSection ( )

Helper function which returns the structural cross-section for the element.

Definition at line 1237 of file structuralelement.C.

References oofem::Element::giveCrossSection().

Referenced by adaptiveUpdate(), addNonlocalStiffnessContributions(), oofem::TrPlaneStrRot::computeBodyLoadVectorAt(), oofem::CCTPlate::computeBodyLoadVectorAt(), oofem::TrPlaneStrRot3d::computeBodyLoadVectorAt(), oofem::TrPlanestressRotAllman3d::computeBodyLoadVectorAt(), oofem::DKTPlate3d::computeBodyLoadVectorAt(), oofem::CCTPlate3d::computeBodyLoadVectorAt(), oofem::DKTPlate::computeBodyLoadVectorAt(), oofem::Quad1Mindlin::computeBodyLoadVectorAt(), oofem::RerShell::computeBodyLoadVectorAt(), oofem::Quad1MindlinShell3D::computeBodyLoadVectorAt(), oofem::NLStructuralElement::computeCauchyStressVector(), oofem::Beam2d::computeConsistentMassMatrix(), oofem::Beam3d::computeConsistentMassMatrix(), oofem::StructuralElementEvaluator::computeConsistentMassMatrix(), oofem::QTruss1d::computeConstitutiveMatrixAt(), oofem::QTruss1dGrad::computeConstitutiveMatrixAt(), oofem::LIBeam2dNL::computeConstitutiveMatrixAt(), oofem::TrPlaneStrRot::computeConstitutiveMatrixAt(), oofem::Tr_Warp::computeConstitutiveMatrixAt(), oofem::tet21ghostsolid::computeConstitutiveMatrixAt(), oofem::Structural3DElement::computeConstitutiveMatrixAt(), oofem::Truss1d::computeConstitutiveMatrixAt(), oofem::Truss2d::computeConstitutiveMatrixAt(), oofem::CCTPlate::computeConstitutiveMatrixAt(), oofem::LIBeam3dNL::computeConstitutiveMatrixAt(), oofem::QDKTPlate::computeConstitutiveMatrixAt(), oofem::LIBeam3d2::computeConstitutiveMatrixAt(), oofem::Truss3d::computeConstitutiveMatrixAt(), oofem::Quad1PlateSubSoil::computeConstitutiveMatrixAt(), oofem::Tria1PlateSubSoil::computeConstitutiveMatrixAt(), oofem::LIBeam2d::computeConstitutiveMatrixAt(), oofem::DKTPlate::computeConstitutiveMatrixAt(), oofem::LIBeam3d::computeConstitutiveMatrixAt(), oofem::Quad1Mindlin::computeConstitutiveMatrixAt(), oofem::PlaneStressElement::computeConstitutiveMatrixAt(), oofem::LIBeam3dNL2::computeConstitutiveMatrixAt(), oofem::RerShell::computeConstitutiveMatrixAt(), oofem::Lattice2d::computeConstitutiveMatrixAt(), oofem::PlaneStrainElement::computeConstitutiveMatrixAt(), oofem::Quad1MindlinShell3D::computeConstitutiveMatrixAt(), oofem::AxisymElement::computeConstitutiveMatrixAt(), oofem::Beam2d::computeConstitutiveMatrixAt(), oofem::Beam3d::computeConstitutiveMatrixAt(), oofem::NLStructuralElement::computeFirstPKStressVector(), oofem::RerShell::computeLocalCoordinates(), oofem::CCTPlate::computeLocalCoordinates(), oofem::QDKTPlate::computeLocalCoordinates(), oofem::DKTPlate::computeLocalCoordinates(), oofem::LIBeam2dNL::computeLumpedMassMatrix(), oofem::LIBeam3d::computeLumpedMassMatrix(), oofem::LIBeam2d::computeLumpedMassMatrix(), oofem::Truss3d::computeLumpedMassMatrix(), oofem::Truss1d::computeLumpedMassMatrix(), oofem::LTRSpace::computeLumpedMassMatrix(), oofem::Truss2d::computeLumpedMassMatrix(), oofem::RerShell::computeLumpedMassMatrix(), oofem::LIBeam3dNL::computeLumpedMassMatrix(), oofem::LIBeam3dNL2::computeLumpedMassMatrix(), oofem::LIBeam3d2::computeLumpedMassMatrix(), oofem::Quad1Mindlin::computeLumpedMassMatrix(), oofem::Quad1MindlinShell3D::computeLumpedMassMatrix(), oofem::CCTPlate::computeLumpedMassMatrix(), oofem::DKTPlate::computeLumpedMassMatrix(), oofem::SolidShell::computeStiffnessMatrix(), oofem::Quad1MindlinShell3D::computeStiffnessMatrix(), oofem::MITC4Shell::computeStiffnessMatrix(), oofem::NLStructuralElement::computeStiffnessMatrix(), oofem::GradDpElement::computeStiffnessMatrix_kk(), oofem::GradDpElement::computeStiffnessMatrix_ku(), oofem::GradDpElement::computeStiffnessMatrix_uk(), oofem::GradDpElement::computeStiffnessMatrix_uu(), oofem::NLStructuralElement::computeStiffnessMatrix_withIRulesAsSubcells(), oofem::TrPlanestressRotAllman::computeStiffnessMatrixZeroEnergyStabilization(), computeStrainVector(), oofem::QTruss1d::computeStressVector(), oofem::LIBeam2dNL::computeStressVector(), oofem::TrPlaneStrRot::computeStressVector(), oofem::Tr_Warp::computeStressVector(), oofem::tet21ghostsolid::computeStressVector(), oofem::Structural3DElement::computeStressVector(), oofem::Truss1d::computeStressVector(), oofem::Truss2d::computeStressVector(), oofem::CCTPlate::computeStressVector(), oofem::LIBeam3dNL::computeStressVector(), oofem::QDKTPlate::computeStressVector(), oofem::LIBeam3d2::computeStressVector(), oofem::Truss3d::computeStressVector(), oofem::Quad1PlateSubSoil::computeStressVector(), oofem::Tria1PlateSubSoil::computeStressVector(), oofem::DKTPlate::computeStressVector(), oofem::LIBeam2d::computeStressVector(), oofem::Quad1Mindlin::computeStressVector(), oofem::LIBeam3d::computeStressVector(), oofem::PlaneStressElement::computeStressVector(), oofem::MITC4Shell::computeStressVector(), oofem::LIBeam3dNL2::computeStressVector(), oofem::Lattice2d::computeStressVector(), oofem::RerShell::computeStressVector(), oofem::PlaneStrainElement::computeStressVector(), oofem::Quad1MindlinShell3D::computeStressVector(), oofem::AxisymElement::computeStressVector(), oofem::Beam2d::computeStressVector(), oofem::Beam3d::computeStressVector(), oofem::GradDpElement::computeStressVectorAndLocalCumulatedStrain(), oofem::DKTPlate::computeVertexBendingMoments(), oofem::TrPlaneStrRot3d::computeVolumeAround(), oofem::TrPlanestressRotAllman3d::computeVolumeAround(), createMaterialStatus(), oofem::TrPlaneStrRot3d::giveCharacteristicTensor(), oofem::MITC4Shell::giveCharacteristicTensor(), giveClassName(), oofem::SolidShell::giveInternalForcesVector(), oofem::GradDpElement::giveInternalForcesVector(), oofem::Quad1MindlinShell3D::giveInternalForcesVector(), oofem::MITC4Shell::giveInternalForcesVector(), oofem::GradDpElement::giveLocalInternalForcesVector(), oofem::MITC4Shell::giveMidplaneIPValue(), giveNonlocalLocationArray(), showExtendedSparseMtrxStructure(), updateBeforeNonlocalAverage(), oofem::LinearizedDilationForceAssembler::vectorFromElement(), oofem::XfemStructuralElementInterface::XfemElementInterface_computeConsistentMassMatrix(), oofem::XfemStructuralElementInterface::XfemElementInterface_updateIntegrationRule(), and ~StructuralElement().

virtual void oofem::StructuralElement::giveSurfaceDofMapping ( IntArray answer,
int  iSurf 
) const
inlineprotectedvirtual

Assembles surface dof mapping mask, which provides mapping between surface local DOFs and "global" element DOFs.

Mask can be imagined as local surface code numbers used to localize local DOFs to element DOFs.

Parameters
answerSurface DOF mapping.
iSurfSurface number

Reimplemented in oofem::MITC4Shell, oofem::Shell7Base, oofem::DKTPlate, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::DKTPlate3d, oofem::QDKTPlate, oofem::CCTPlate3d, oofem::Structural3DElement, oofem::TrPlaneStrRot3d, oofem::Tr2Shell7, and oofem::Tr2Shell7XFEM.

Definition at line 396 of file structuralelement.h.

References oofem::IntArray::clear().

IRResultType oofem::StructuralElement::initializeFrom ( InputRecord ir)
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.

See also
IR_GIVE_FIELD
IR_GIVE_OPTIONAL_FIELD
Parameters
irInput record to initialize from.
Returns
IRResultType

Reimplemented from oofem::Element.

Reimplemented in oofem::NLStructuralElement, oofem::AbaqusUserElement, oofem::Shell7BaseXFEM, oofem::Beam3d, oofem::MITC4Shell, oofem::Shell7Base, oofem::Beam2d, oofem::DKTPlate, oofem::QDKTPlate, oofem::RerShell, oofem::SpringElement, oofem::CCTPlate, oofem::Lattice2d, oofem::Quad1MindlinShell3D, oofem::NodalSpringElement, oofem::Truss2d, oofem::TrPlaneStress2dXFEM, oofem::TrPlaneStrain, oofem::LSpace, oofem::Truss3d, oofem::Quad1PlaneStrain, oofem::PlaneStress2d, oofem::PlaneStress2dXfem, oofem::Axisymm3d, oofem::LIBeam3dNL, oofem::LumpedMassElement, oofem::QTrPlaneStress2dXFEM, oofem::LIBeam3dNL2, oofem::InterfaceElem1d, oofem::TrPlanestressRotAllman, oofem::TR_SHELL01, oofem::LIBeam2d, oofem::TR_SHELL02, oofem::CohesiveSurface3d, oofem::Quad1Mindlin, oofem::Structural2DElement, oofem::LIBeam2dNL, oofem::InterfaceElem2dQuad, oofem::SolidShell, oofem::InterfaceElem2dLin, oofem::TrPlaneStrRot, oofem::InterfaceElement3dTrLin, oofem::LIBeam3d2, oofem::LineDistributedSpring, oofem::QTrPlaneStress2d, oofem::Quad1PlateSubSoil, oofem::Tria1PlateSubSoil, oofem::QTRSpace, oofem::Q4Axisymm, oofem::L4Axisymm, oofem::Quad2PlateSubSoil, oofem::Structural3DElement, oofem::Tr_Warp, oofem::MacroLSpace, oofem::QSpace, oofem::QWedge, oofem::LWedge, oofem::Q27Space, oofem::QTruss1dGrad, oofem::HTSelement, oofem::LIBeam3d, oofem::QSpaceGrad, oofem::QTRSpaceGrad, oofem::QWedgeGrad, oofem::LatticeStructuralElement, oofem::QTruss1d, oofem::CoupledFieldsElement, oofem::QPlaneStrainGrad, oofem::QPlaneStressGrad, oofem::QTrPlaneStressGrad, and oofem::QTrPlaneStrainGrad.

Definition at line 1224 of file structuralelement.C.

References oofem::Element::initializeFrom().

Referenced by computeInitialStressMatrix(), oofem::QTruss1d::initializeFrom(), oofem::LatticeStructuralElement::initializeFrom(), oofem::LIBeam3d::initializeFrom(), oofem::HTSelement::initializeFrom(), oofem::QTruss1dGrad::initializeFrom(), oofem::Tr_Warp::initializeFrom(), oofem::Quad2PlateSubSoil::initializeFrom(), oofem::Q4Axisymm::initializeFrom(), oofem::Quad1PlateSubSoil::initializeFrom(), oofem::Tria1PlateSubSoil::initializeFrom(), oofem::LineDistributedSpring::initializeFrom(), oofem::InterfaceElement3dTrLin::initializeFrom(), oofem::InterfaceElem2dLin::initializeFrom(), oofem::InterfaceElem2dQuad::initializeFrom(), oofem::CohesiveSurface3d::initializeFrom(), oofem::LIBeam2d::initializeFrom(), oofem::TR_SHELL02::initializeFrom(), oofem::TR_SHELL01::initializeFrom(), oofem::InterfaceElem1d::initializeFrom(), oofem::LumpedMassElement::initializeFrom(), oofem::Axisymm3d::initializeFrom(), oofem::SpringElement::initializeFrom(), oofem::RerShell::initializeFrom(), oofem::Beam2d::initializeFrom(), oofem::Beam3d::initializeFrom(), oofem::AbaqusUserElement::initializeFrom(), and oofem::NLStructuralElement::initializeFrom().

void oofem::StructuralElement::setupIRForMassMtrxIntegration ( IntegrationRule iRule)
protectedvirtual
void oofem::StructuralElement::showExtendedSparseMtrxStructure ( CharType  mtrx,
oofegGraphicContext gc,
TimeStep tStep 
)
virtual

Shows extended sparse structure (for example, due to nonlocal interactions for tangent stiffness)

Reimplemented from oofem::Element.

Definition at line 1358 of file structuralelement.C.

References oofem::Element::giveDefaultIntegrationRulePtr(), giveStructuralCrossSection(), and oofem::NonlocalMaterialStiffnessInterfaceType.

Referenced by giveClassName().

void oofem::StructuralElement::updateBeforeNonlocalAverage ( TimeStep tStep)
virtual

Updates internal element state (in all integration points of receiver) before nonlocal averaging takes place.

Used by so nonlocal materials, because their response in particular point depends not only on state in this point, but depends also on state in point's neighborhood. Nonlocal quantity is computed as nonlocal average of local quantities. Therefore, before updating integration point state depending on nonlocal quantity (or quantities), local quantities in all integration points must be updated in advance. This function updates local quantities of material model using updateBeforeNonlocalAverage member function of structural nonlocal material class.

Parameters
tStepTime step.

Reimplemented from oofem::Element.

Definition at line 970 of file structuralelement.C.

References computeStrainVector(), oofem::Element_remote, oofem::Element::giveParallelMode(), giveStructuralCrossSection(), oofem::Element::integrationRulesArray, oofem::NonlocalMaterialExtensionInterfaceType, and oofem::StructuralNonlocalMaterialExtensionInterface::updateBeforeNonlocAverage().

Referenced by computeInitialStressMatrix().

void oofem::StructuralElement::updateInternalState ( TimeStep tStep)
virtual

Updates element state after equilibrium in time step has been reached.

Default implementation updates all integration rules defined by integrationRulesArray member variable. Doing this, all integration points and their material statuses are updated also. All temporary history variables, which now describe equilibrium state are copied into equilibrium ones. The existing internal state is used for update.

Parameters
tStepTime step for newly reached state.
See also
Material::updateYourself
IntegrationRule::updateYourself
GaussPoint::updateYourself
Element::updateInternalState

Reimplemented from oofem::Element.

Reimplemented in oofem::AbaqusUserElement, oofem::SpringElement, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::NodalSpringElement, and oofem::LumpedMassElement.

Definition at line 955 of file structuralelement.C.

References computeStrainVector(), computeStressVector(), and oofem::Element::integrationRulesArray.

Referenced by computeInitialStressMatrix().

void oofem::StructuralElement::updateYourself ( TimeStep tStep)
virtual

Updates element state after equilibrium in time step has been reached.

Default implementation updates all integration rules defined by integrationRulesArray member variable. Doing this, all integration points and their material statuses are updated also. All temporary history variables, which now describe equilibrium state are copied into equilibrium ones. The existing internal state is used for update.

Parameters
tStepTime step for newly reached state.
See also
Material::updateYourself
IntegrationRule::updateYourself
GaussPoint::updateYourself
Element::updateInternalState

Reimplemented from oofem::Element.

Reimplemented in oofem::AbaqusUserElement, oofem::LIBeam3dNL, oofem::LIBeam3dNL2, oofem::SpringElement, oofem::TR_SHELL01, oofem::TR_SHELL02, oofem::NodalSpringElement, oofem::MacroLSpace, oofem::LumpedMassElement, oofem::LIBeam3d2, oofem::Shell7BaseXFEM, oofem::QTrPlaneStress2dXFEM, oofem::TrPlaneStress2dXFEM, and oofem::PlaneStress2dXfem.

Definition at line 939 of file structuralelement.C.

References oofem::Element::activityTimeFunction, oofem::Element::computeVectorOf(), initialDisplacements, oofem::Element::isActivated(), and oofem::Element::updateYourself().

Referenced by computeInitialStressMatrix(), oofem::REGISTER_Element(), oofem::TrPlaneStress2dXFEM::updateYourself(), oofem::QTrPlaneStress2dXFEM::updateYourself(), oofem::Shell7BaseXFEM::updateYourself(), oofem::LIBeam3d2::updateYourself(), oofem::TR_SHELL02::updateYourself(), oofem::TR_SHELL01::updateYourself(), oofem::LIBeam3dNL2::updateYourself(), oofem::LIBeam3dNL::updateYourself(), and oofem::AbaqusUserElement::updateYourself().

Friends And Related Function Documentation

friend class GradDpElement
friend

Definition at line 525 of file structuralelement.h.

friend class IDNLMaterial
friend

Definition at line 521 of file structuralelement.h.

friend class MisesMatNl
friend

Definition at line 523 of file structuralelement.h.

friend class RankineMatNl
friend

Definition at line 524 of file structuralelement.h.

friend class TrabBoneNL3D
friend

Definition at line 522 of file structuralelement.h.

Member Data Documentation


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