38 #include "../sm/Elements/structuralelement.h" 41 #define _IFT_HTSelement_Name "htselement" 84 OOFEM_ERROR(
"Function not defined for this element and should never be called. This is a bug.");
140 #endif // htselement_h virtual int giveNumberOfNodes() const
Returns number of nodes of receiver.
double computeVolumeAroundSide(GaussPoint *gp, int elemSideNumber)
void uv5(FloatArray &answer, double x, double y)
void uv12(FloatArray &answer, double x, double y)
void uv3(FloatArray &answer, double x, double y)
#define _IFT_HTSelement_Name
void uv4(FloatArray &answer, double x, double y)
void sv1(FloatArray &answer, double x, double y)
void computePrescribedDisplacementLoadVectorAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
void uv8(FloatArray &answer, double x, double y)
void computeFMatrixAt(FloatMatrix &answer, FloatMatrix N, GaussPoint *gp, int sideNumber)
void uv6(FloatArray &answer, double x, double y)
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
void uv7(FloatArray &answer, double x, double y)
void sv11(FloatArray &answer, double x, double y)
void sv7(FloatArray &answer, double x, double y)
void uv10(FloatArray &answer, double x, double y)
Implements a Hybrid-Trefftz element See http://en.wikipedia.org/wiki/Trefftz_method for description...
void computeCenterOfGravity()
void sv12(FloatArray &answer, double x, double y)
Class implementing an array of integers.
MatResponseMode
Describes the character of characteristic material matrix.
void computeSvMatrixAt(FloatMatrix &answer, GaussPoint *gp, int sideNumber)
void uv11(FloatArray &answer, double x, double y)
virtual StructuralElement * giveStructuralElement()
void uv9(FloatArray &answer, double x, double y)
virtual const char * giveInputRecordName() const
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int, int)
Computes the geometrical matrix of receiver in given integration point.
Abstract base class for all "structural" finite elements.
void computeAMatrixAt(FloatMatrix &answer, FloatMatrix N, GaussPoint *gp, int sideNumber)
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
Returns equivalent nodal forces vectors.
void sv5(FloatArray &answer, double x, double y)
void sv2(FloatArray &answer, double x, double y)
ElementExtension
Type representing element extension.
void sv8(FloatArray &answer, double x, double y)
void sv3(FloatArray &answer, double x, double y)
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Node * giveSideNode(int elementSideNumber, int nodeNumber)
void sv10(FloatArray &answer, double x, double y)
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.
void computeOutwardNormalMatrix(FloatMatrix &answer, int sideNumber)
void uv25_4(FloatArray &answer, double x, double y)
void sv9(FloatArray &answer, double x, double y)
void sv25_4(FloatArray &answer, double x, double y)
void computePsVectorAt(FloatArray &answer, FloatArray t, GaussPoint *gp)
double u_gammaLin(GaussPoint *gp)
void sv4(FloatArray &answer, double x, double y)
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Class representing vector of real numbers.
void computePuVectorAt(FloatArray &answer, FloatMatrix N, FloatArray u, GaussPoint *gp, int sideNumber)
Implementation of matrix containing floating point numbers.
void computeUgammaMatrixAt(FloatMatrix &answer, GaussPoint *gp)
IRResultType
Type defining the return values of InputRecord reading operations.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
HTSelement(int n, Domain *d)
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
void uv2(FloatArray &answer, double x, double y)
virtual const char * giveClassName() const
void computeUvMatrixAt(FloatMatrix &answer, GaussPoint *gp, int sideNubmer)
double u_gammaConst(GaussPoint *gp)
virtual int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
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...
void uv1(FloatArray &answer, double x, double y)
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
the oofem namespace is to define a context or scope in which all oofem names are defined.
void sv6(FloatArray &answer, double x, double y)
Class implementing node in finite element mesh.
double giveSideLength(int sideNumber)
Class representing integration point in finite element program.
Class representing solution step.
Element extension for edge loads.
void resize(int s)
Resizes receiver towards requested size.