35 #ifndef layeredcrosssection_h 36 #define layeredcrosssection_h 38 #include "../sm/CrossSections/structuralcrosssection.h" 50 #define _IFT_LayeredCrossSection_Name "layeredcs" 51 #define _IFT_LayeredCrossSection_nlayers "nlayers" 52 #define _IFT_LayeredCrossSection_layermaterials "layermaterials" 53 #define _IFT_LayeredCrossSection_interfacematerials "interfacematerials" 54 #define _IFT_LayeredCrossSection_layerRotations "rotations" 55 #define _IFT_LayeredCrossSection_thicks "thicks" 56 #define _IFT_LayeredCrossSection_widths "widths" 57 #define _IFT_LayeredCrossSection_midsurf "midsurf" 58 #define _IFT_LayeredCrossSection_nintegrationpoints "nintegrationpoints" 59 #define _IFT_LayeredCrossSection_initiationlimits "initiationlimits" 63 class StructuralMaterial;
180 return this->layerMaterials.
at(layer);
186 return this->interfacerMaterials.
at(interface);
202 return this->layerMidZ.
at(layer);
205 return this->layerThicks.
at(layer);
297 virtual const char *
giveClassName()
const {
return "LayeredIntegrationRule"; }
308 virtual int SetUpPointsOnWedge(
int nPointsTri,
int nPointsDepth,
MaterialMode mode);
311 #endif // layeredcrosssection_h
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual void giveStiffnessMatrix_PlaneStress(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
virtual Material * giveMaterial(IntegrationPoint *ip)
Returns the material associated with the GP.
double giveMidSurfaceZcoordFromBottom()
int numberOfIntegrationPoints
num integration points per layer
int giveLayerMaterial(int layer)
virtual int checkConsistency()
Allows programmer to test some internal data, before computation begins.
virtual void give2dPlateStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d plate stiffness matrix.
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
double & at(int i)
Coefficient access function.
virtual void giveRealStress_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
virtual void giveStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes the material stiffness matrix dPdF of receiver in a given integration point, respecting its history.
virtual void give3dShellStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 3d shell stiffness matrix.
virtual void giveGeneralizedStress_Plate(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
FloatArray layerThicks
Thickness for each layer.
Abstract base class for all finite elements.
void mapLayerGpCoordsToShellCoords(std::vector< std::unique_ptr< IntegrationRule > > &layerIntegrationRulesArray)
MaterialMode giveCorrespondingSlaveMaterialMode(MaterialMode mode)
virtual const char * giveClassName() const
CrossSectionProperty
List of properties possibly stored in a cross section.
#define _IFT_LayeredCrossSection_Name
MaterialMode
Type representing material mode of integration point.
virtual FloatArray * imposeStressConstrainsOnGradient(GaussPoint *gp, FloatArray *)
Returns modified gradient of stress vector, which is used to bring stresses back to yield surface...
Class implementing an array of integers.
int & at(int i)
Coefficient access function.
MatResponseMode
Describes the character of characteristic material matrix.
void giveInterfaceXiCoords(FloatArray &answer)
virtual int giveIPValue(FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
virtual void giveFirstPKStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep)
Computes the First Piola-Kirchoff stress vector for a given deformation gradient and integration poin...
virtual void giveStiffnessMatrix_3d(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing the stiffness matrix.
Abstract base class representing integration rule.
This class implements a layered cross section in a finite element problem.
double giveMidSurfaceXiCoordFromBottom()
double midSurfaceXiCoordFromBottom
virtual void give2dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix for 2d beams.
IntArray layerMaterials
Material of each layer.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
double computeIntegralThick()
Returns the total thickness of all layers.
LayeredCrossSection(int n, Domain *d)
virtual void printYourself()
Prints receiver state on stdout. Useful for debugging.
virtual void giveGeneralizedStress_PlateSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
virtual const char * giveInputRecordName() const
FloatArray layerMidZ
z-coord of the mid plane for each layer
Material * giveMaterial(int n)
Service for accessing particular domain material model.
int giveLayer(GaussPoint *gp)
virtual FloatArray * imposeStrainConstrainsOnGradient(GaussPoint *gp, FloatArray *)
Returns modified gradient of strain vector, which is used to compute plastic strain increment...
double midSurfaceZcoordFromBottom
double giveLayerThickness(int layer)
int giveInterfaceMaterialNum(int interface)
virtual void giveGeneralizedStress_Beam2d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
Computes the generalized stress vector for given strain and integration point.
void setupLayeredIntegrationRule(std::vector< std::unique_ptr< IntegrationRule > > &layerIntegrationRulesArray, Element *el, int numInPlanePoints)
FloatArray layerWidths
Width for each layer.
virtual const char * giveClassName() const
double giveLayerMidZ(int layer)
virtual void createMaterialStatus(GaussPoint &iGP)
virtual void giveRealStress_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
virtual void giveStiffnessMatrix_PlaneStrain(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
IntArray upperInterfacePoints
virtual contextIOResultType saveIPContext(DataStream &stream, ContextMode mode, GaussPoint *gp)
Stores integration point state to output stream.
void setupLayerMidPlanes()
Abstract base class for all material models.
IntArray interfacerMaterials
Interface (cohesive zone) material for each interface.
virtual void giveRealStress_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
virtual void giveRealStress_3dDegeneratedShell(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
Class representing vector of real numbers.
virtual void give3dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix for 2d beams.
virtual void giveCharMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix of receiver in given integration point, respecting its history...
Implementation of matrix containing floating point numbers.
IRResultType
Type defining the return values of InputRecord reading operations.
virtual contextIOResultType restoreIPContext(DataStream &stream, ContextMode mode, GaussPoint *gp)
Reads integration point state to output stream.
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
virtual bool isCharacteristicMtrxSymmetric(MatResponseMode mode)
Check for symmetry of stiffness matrix.
FloatArray layerRots
Rotation of the material in each layer.
virtual void giveCauchyStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep)
Computes the Cauchy stress vector for a given increment of deformation gradient and given integration...
virtual void give2dPlateSubSoilStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing subsoil stiffness matrix for plates.
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
virtual void giveGeneralizedStress_Shell(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
LayeredCrossSectionInterface()
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
virtual int packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
Pack all necessary data of integration point (according to element parallel_mode) into given communic...
virtual void giveGeneralizedStress_Beam3d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
virtual int unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
Unpack and updates all necessary data of given integration point (according to element parallel_mode)...
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Domain * giveDomain() const
virtual void giveRealStress_Warping(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
The element interface required by LayeredCrossSection.
int giveNumIntegrationPointsInLayer()
virtual void give3dDegeneratedShellStiffMtrx(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Method for computing 3d shell stiffness matrix on degenerated shell elements.
virtual void giveMembraneRotStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing membrane stiffness matrix with added drilling stiffness.
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 giveStiffnessMatrix_dCde(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes the material stiffness matrix dCde of receiver in a given integration point, respecting its history.
virtual void giveGeneralizedStress_MembraneRot(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
virtual ~LayeredCrossSection()
virtual IRResultType initializeFrom(InputRecord *ir)
GaussPoint * giveSlaveGaussPoint(GaussPoint *gp, int slaveIndex)
virtual int estimatePackSize(DataStream &buff, GaussPoint *ip)
Estimates the necessary pack size to hold all packed data of receiver.
Class representing integration point in finite element program.
Class representing solution step.
Material * giveInterfaceMaterial(int interface)
virtual void giveStiffnessMatrix_1d(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
virtual void giveRealStress_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)