35 #ifndef planstrssxfem_h 36 #define planstrssxfem_h 38 #include "../sm/Elements/PlaneStress/planstrss.h" 39 #include "../sm/xfem/xfemstructuralelementinterface.h" 42 #define _IFT_PlaneStress2dXfem_Name "planestress2dxfem" 65 virtual const char *
giveClassName()
const {
return "PlaneStress2dXfem"; }
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *, TimeStep *tStep)
Computes constitutive matrix of receiver.
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
This class implements an isoparametric four-node quadrilateral plane- stress elasticity finite elemen...
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Element_Geometry_Type
Enumerative type used to classify element geometry Possible values are: EGT_point - point in space EG...
virtual void XfemElementInterface_computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
virtual const char * giveInputRecordName() const
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
virtual void postInitialize()
Performs post initialization steps.
Elements with geometry defined as EGT_Composite are exported using individual pieces.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual void computeDeformationGradientVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Computes the deformation gradient in Voigt form at integration point ip and at time step tStep...
MaterialMode
Type representing material mode of integration point.
Class implementing an array of integers.
MatResponseMode
Describes the character of characteristic material matrix.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
virtual const char * giveClassName() const
virtual int computeNumberOfDofs()
Computes or simply returns total number of element's local DOFs.
virtual void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer)
Computes a matrix which, multiplied by the column matrix of nodal displacements, gives the displaceme...
PlaneStress2dXfem(int n, Domain *d)
Constructor.
Provides Xfem interface for a structural element.
Class representing vector of real numbers.
Implementation of matrix containing floating point numbers.
virtual Element_Geometry_Type giveGeometryType() const
Returns the element geometry type.
IRResultType
Type defining the return values of InputRecord reading operations.
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
Evaluates nodal representation of real internal forces.
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
InterfaceType
Enumerative type, used to identify interface type.
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 void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
virtual ~PlaneStress2dXfem()
Destructor.
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...
Class representing integration point in finite element program.
Class representing solution step.
virtual void giveCompositeExportData(std::vector< VTKPiece > &vtkPieces, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep)
VTK Interface.
int numberOfDofMans
Number of dofmanagers.
Temporary class for testing in the usual case instead of PlaneStress2dXfem there will be the standard...