35 #ifndef structuralelement_h 36 #define structuralelement_h 48 #define ALL_STRAINS -1 58 class StructuralCrossSection;
226 TimeStep *tStep,
int useUpdatedGpRecord = 0);
483 int lowerIndx = 1,
int upperIndx =
ALL_STRAINS) = 0;
528 #endif // structuralelement_h InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual void computeResultingIPTemperatureAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
Computes at given time (tStep) the the resulting temperature component array.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Base class for all matrices stored in sparse format.
virtual int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Returns transformation matrix from local surface c.s to element local coordinate system of load vecto...
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
virtual void createMaterialStatus()
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted...
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)=0
Computes constitutive matrix of receiver.
virtual void showExtendedSparseMtrxStructure(CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep)
Shows extended sparse structure (for example, due to nonlocal interactions for tangent stiffness) ...
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...
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Class implementing element body load, acting over whole element volume (e.g., the dead weight)...
Trabecular bone nonlocal material model.
virtual ~StructuralElement()
Destructor.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const
Assembles surface dof mapping mask, which provides mapping between surface local DOFs and "global" el...
virtual void computeEdgeNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
computes edge interpolation matrix
Abstract base class for all finite elements.
virtual void giveNonlocalLocationArray(IntArray &locationArray, const UnknownNumberingScheme &us)
Returns the "nonlocal" location array of receiver.
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
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 c...
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
virtual int adaptiveUpdate(TimeStep *tStep)
Updates the internal state variables stored in all IPs according to already mapped state...
Class implementing an array of integers.
MatResponseMode
Describes the character of characteristic material matrix.
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
virtual void setupIRForMassMtrxIntegration(IntegrationRule &iRule)
Setup Integration Rule Gauss Points for Mass Matrix integration.
virtual void showSparseMtrxStructure(CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep)
Shows sparse structure.
Abstract base class representing integration rule.
Abstract class for gradient formulation of coupled damage-plasticity model(GradDp).
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 m...
Abstract base class for all "structural" finite elements.
virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
Computes the contribution of the given body load (volumetric).
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. ...
Abstract base class representing a boundary load (force, momentum, ...) that acts directly on a bound...
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...
void computeStiffnessMatrix_withIRulesAsSubcells(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
void clear()
Clears the array (zero size).
virtual void giveInternalForcesVector_withIRulesAsSubcells(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void updateInternalState(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Abstract base class allowing to control the way, how equations are assigned to individual DOFs...
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.
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Compute strain vector of receiver evaluated at given integration point at given time step from elemen...
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Class representing vector of real numbers.
Rankine nonlocal material.
virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
Implementation of matrix containing floating point numbers.
IRResultType
Type defining the return values of InputRecord reading operations.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
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.
virtual void computeSurfaceNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
Computes surface interpolation matrix.
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Returns equivalent nodal forces vectors.
StructuralElement(int n, Domain *d)
Constructor.
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.
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes initial stress matrix for linear stability problem.
virtual void computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer)
Computes the unknown vector interpolated at the specified local coordinates.
virtual int checkConsistency()
Performs consistency check.
This class implements a Nonlocal Isotropic Damage Model for Concrete in Tension Model based on nonloc...
virtual void updateBeforeNonlocalAverage(TimeStep *tStep)
Updates internal element state (in all integration points of receiver) before nonlocal averaging take...
virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
Computes volume related to integration point on local surface.
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
virtual int giveNumberOfIPForMassMtrxIntegration()
Return desired number of integration points for consistent mass matrix computation, if required.
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
Load is base abstract class for all loads.
virtual const char * giveClassName() const
Abstract base class for all structural cross section models.
the oofem namespace is to define a context or scope in which all oofem names are defined.
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)...
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
virtual void computeResultingIPEigenstrainAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
Computes at given time the resulting eigenstrain component array.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void addNonlocalStiffnessContributions(SparseMtrx &dest, const UnknownNumberingScheme &s, TimeStep *tStep)
Adds the "nonlocal" contribution to stiffness matrix, to account for nonlocality of material model...
Class representing integration point in finite element program.
Class representing solution step.
InternalStateMode
Determines the mode of internal variable.
virtual IntegrationRule * GetSurfaceIntegrationRule(int order)
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...
virtual void giveMassMtrxIntegrationgMask(IntArray &answer)
Returns mask indicating, which unknowns (their type and ordering is the same as element unknown vecto...