35 #ifndef nodalspringelement_h 36 #define nodalspringelement_h 38 #include "../sm/Elements/structuralelement.h" 42 #define _IFT_NodalSpringElement_Name "nodalspring" 43 #define _IFT_NodalSpringElement_dofmask "dofmask" 44 #define _IFT_NodalSpringElement_springConstants "k" 45 #define _IFT_NodalSpringElement_masses "m" 120 #endif // nodalspringelement_h virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
IntArray dofMask
Dof mask.
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
void clear()
Clears receiver (zero size).
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
FloatArray springConstants
Spring constants.
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Class implementing an array of integers.
MatResponseMode
Describes the character of characteristic material matrix.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Abstract base class for all "structural" finite elements.
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
virtual void updateInternalState(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual const char * giveInputRecordName() const
virtual ~NodalSpringElement()
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes initial stress matrix for linear stability problem.
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Class representing vector of real numbers.
virtual int checkConsistency()
Performs consistency check.
Implementation of matrix containing floating point numbers.
IRResultType
Type defining the return values of InputRecord reading operations.
NodalSpringElement(int n, Domain *d)
virtual bool isCast(TimeStep *tStep)
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
This class implements a simple linear spring element connecting the given node and the ground...
FloatArray masses
total mass of the spring; to be distributed to nodes
virtual void computeMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes mass matrix of receiver.
virtual const char * giveClassName() const
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
the oofem namespace is to define a context or scope in which all oofem names are defined.
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
#define _IFT_NodalSpringElement_Name
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Returns equivalent nodal forces vectors.
Class representing integration point in finite element program.
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Class representing solution step.
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
Computes the stress vector of receiver at given integration point, at time step tStep.
virtual int computeNumberOfGlobalDofs()
Computes the total number of element's global dofs.