OOFEM 3.0
Loading...
Searching...
No Matches
layeredcrosssection.h
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
35#ifndef layeredcrosssection_h
36#define layeredcrosssection_h
37
39#include "floatarray.h"
40#include "floatmatrix.h"
41#include "interface.h"
43#include "domain.h"
44
45#include <vector>
46#include <memory>
47
49
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_nlayerintegrationpoints "layerintegrationpoints"
60#define _IFT_LayeredCrossSection_initiationlimits "initiationlimits"
61#define _IFT_LayeredCrossSection_shearcoeff_xz "beamshearcoeffxz"
63
64namespace oofem {
65class StructuralMaterial;
66
94{
95protected:
107 double totalThick = 0.;
108 double area = 0.;
109 double beamShearCoeffxz = 1.0;
110public:
113 { }
114
115 void initializeFrom(InputRecord &ir) override;
116 void giveInputRecord(DynamicInputRecord &input) override;
117
118 void createMaterialStatus(GaussPoint &iGP) override;
119
120 //Create slave integration points
121 int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element) override;
131 int setupIntegrationPoints(IntegrationRule &irule, int npointsXY, int npointsZ, Element *element) override;
132
133 FloatArrayF< 6 >giveRealStress_3d(const FloatArrayF< 6 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override;
134 FloatArrayF< 6 >giveRealStress_3dDegeneratedShell(const FloatArrayF< 6 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override;
135 FloatArrayF< 4 >giveRealStress_PlaneStrain(const FloatArrayF< 4 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override;
136 FloatArrayF< 3 >giveRealStress_PlaneStress(const FloatArrayF< 3 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override;
137 FloatArrayF< 1 >giveRealStress_1d(const FloatArrayF< 1 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override;
138 FloatArrayF< 2 >giveRealStress_Warping(const FloatArrayF< 2 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override;
139
140 FloatMatrixF< 6, 6 >giveStiffnessMatrix_3d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
141 FloatMatrixF< 3, 3 >giveStiffnessMatrix_PlaneStress(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
142 FloatMatrixF< 4, 4 >giveStiffnessMatrix_PlaneStrain(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
143 FloatMatrixF< 1, 1 >giveStiffnessMatrix_1d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
144
145 FloatArrayF< 3 >giveGeneralizedStress_Beam2d(const FloatArrayF< 3 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override;
146 FloatArrayF< 6 >giveGeneralizedStress_Beam3d(const FloatArrayF< 6 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override;
147 FloatArrayF< 5 >giveGeneralizedStress_Plate(const FloatArrayF< 5 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override;
148 FloatArrayF< 8 >giveGeneralizedStress_Shell(const FloatArrayF< 8 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override;
149 FloatArrayF< 9 >giveGeneralizedStress_ShellRot(const FloatArrayF< 9 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override;
150 FloatArrayF< 4 >giveGeneralizedStress_MembraneRot(const FloatArrayF< 4 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override;
151 FloatArrayF< 3 >giveGeneralizedStress_PlateSubSoil(const FloatArrayF< 3 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override;
152
153 void giveCharMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override;
154
155 bool isCharacteristicMtrxSymmetric(MatResponseMode mode) const override;
156
157 FloatMatrixF< 3, 3 >give2dBeamStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
158 FloatMatrixF< 6, 6 >give3dBeamStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
159 FloatMatrixF< 5, 5 >give2dPlateStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
160 FloatMatrixF< 8, 8 >give3dShellStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
161 FloatMatrixF< 9, 9 >give3dShellRotStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
162 FloatMatrixF< 6, 6 >give3dDegeneratedShellStiffMtrx(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const override;
163 FloatMatrixF< 4, 4 >giveMembraneRotStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
164 FloatMatrixF< 3, 3 >give2dPlateSubSoilStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
165
168
169 double give(CrossSectionProperty a, GaussPoint *gp) const override;
170 double give(CrossSectionProperty a, const FloatArray &coords, Element *elem, bool local) const override;
171 int giveNumberOfLayers() const;
172 int giveLayer(GaussPoint *gp) const;
173
175 double computeIntegralThick() const;
176 void setupLayerMidPlanes();
177
178 int giveLayerMaterial(int layer) const {
179 return this->layerMaterials.at(layer);
180 }
181
182 Material *giveMaterial(IntegrationPoint *ip) const override;
183
184 int giveInterfaceMaterialNum(int interface) {
185 return this->interfacerMaterials.at(interface);
186 }
187
189 int matNum = this->giveInterfaceMaterialNum(interface);
190 if ( matNum ) {
191 return this->giveDomain()->giveMaterial(this->interfacerMaterials.at(interface) );
192 } else {
193 return nullptr;
194 }
195 }
196
197 int checkConsistency() override;
198
199 double giveLayerMidZ(int layer) const {
200 // Gives the z-coord measured from the geometric midplane of the (total) cross section.
201 return this->layerMidZ.at(layer);
202 }
203 double giveLayerThickness(int layer) const {
204 return this->layerThicks.at(layer);
205 }
207 return this->numberOfIntegrationPoints;
208 }
210 return this->midSurfaceZcoordFromBottom;
211 }
213 return this->midSurfaceXiCoordFromBottom;
214 }
215 void giveInterfaceXiCoords(FloatArray &answer) const;
216
217 // identification and auxiliary functions
218 const char *giveInputRecordName() const override { return _IFT_LayeredCrossSection_Name; }
219 const char *giveClassName() const override { return "LayeredCrossSection"; }
220 void printYourself() override;
221
222 static MaterialMode giveCorrespondingSlaveMaterialMode(MaterialMode mode);
223 GaussPoint *giveSlaveGaussPoint(GaussPoint *gp, int layer, int igp) const;
224
225 void saveIPContext(DataStream &stream, ContextMode mode, GaussPoint *gp) override;
226 void restoreIPContext(DataStream &stream, ContextMode mode, GaussPoint *gp) override;
227
228
229 void mapLayerGpCoordsToShellCoords(std::vector< std::unique_ptr< IntegrationRule > > &layerIntegrationRulesArray);
230
231 void setupLayeredIntegrationRule(std::vector< std::unique_ptr< IntegrationRule > > &layerIntegrationRulesArray, Element *el, int numInPlanePoints);
232
233 int giveIPValue(FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep) override;
234 double give(int aProperty, GaussPoint *gp) const override;
235
236 int packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip) override
237 {
238 OOFEM_ERROR("not implemented");
239 }
240
241 int unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip) override
242 {
243 OOFEM_ERROR("not implemented");
244 }
245
246 int estimatePackSize(DataStream &buff, GaussPoint *ip) override
247 {
248 OOFEM_ERROR("not implemented");
249 }
250
251
252
253 FloatArrayF< 9 >giveFirstPKStress_3d(const FloatArrayF< 9 > &reducedvF, GaussPoint *gp, TimeStep *tStep) const override;
254 FloatArrayF< 5 >giveFirstPKStress_PlaneStrain(const FloatArrayF< 5 > &reducedvF, GaussPoint *gp, TimeStep *tStep) const override;
255 FloatArrayF< 4 >giveFirstPKStress_PlaneStress(const FloatArrayF< 4 > &reducedvF, GaussPoint *gp, TimeStep *tStep) const override;
256 FloatArrayF< 1 >giveFirstPKStress_1d(const FloatArrayF< 1 > &reducedvF, GaussPoint *gp, TimeStep *tStep) const override;
257
258
259 void giveCharMaterialStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override;
260 FloatMatrixF< 9, 9 >giveStiffnessMatrix_dPdF_3d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
261 FloatMatrixF< 5, 5 >giveStiffnessMatrix_dPdF_PlaneStrain(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
262 FloatMatrixF< 4, 4 >giveStiffnessMatrix_dPdF_PlaneStress(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
263 FloatMatrixF< 1, 1 >giveStiffnessMatrix_dPdF_1d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
264
265 void giveCauchyStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep) override
266 { OOFEM_ERROR("not implemented"); }
267 void giveStiffnessMatrix_dCde(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
268 { OOFEM_ERROR("not implemented"); }
269
270protected:
271 double giveArea() const;
272 int giveSlaveGPIndex(int ilayer, int igp) const;
273};
274
279{
280public:
282
292 virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep) = 0;
293};
294
296{
297public:
298 LayeredIntegrationRule(int n, Element *e, int startIndx, int endIndx, bool dynamic = false);
300
301 const char *giveClassName() const override { return "LayeredIntegrationRule"; }
302
303 //virtual int getRequiredNumberOfIntegrationPoints(integrationDomain dType, int approxOrder);
304
305 // Stores the ip numbers of the points lying on the lower and upper surface of the element.
306 // Thus they will correspond to points lying on the interface between layers.
308
309 //int SetUpPointsOnLine(int, MaterialMode) override; // could be used for beams
310 //int SetUpPointsOnCube(int, MaterialMode mode) override; // could be used for plates/shells/solids
311 int SetUpPointsOnWedge(int nPointsTri, int nPointsDepth, MaterialMode mode) override;
312};
313} // end namespace oofem
314#endif // layeredcrosssection_h
Material * giveMaterial(int n)
Definition domain.C:284
Domain * giveDomain() const
Definition femcmpnn.h:97
double & at(Index i)
Definition floatarray.h:202
int & at(std::size_t i)
Definition intarray.h:104
IntegrationRule(int n, Element *e, int startIndx, int endIndx, bool dynamic)
Interface()
Constructor.
Definition interface.h:86
virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)=0
void printYourself() override
Prints receiver state on stdout. Useful for debugging.
int giveInterfaceMaterialNum(int interface)
FloatArray * imposeStressConstrainsOnGradient(GaussPoint *gp, FloatArray *) override
FloatArrayF< 2 > giveRealStress_Warping(const FloatArrayF< 2 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 3, 3 > give2dBeamStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
void giveInputRecord(DynamicInputRecord &input) override
void createMaterialStatus(GaussPoint &iGP) override
double giveMidSurfaceXiCoordFromBottom() const
FloatArrayF< 9 > giveGeneralizedStress_ShellRot(const FloatArrayF< 9 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatArray layerRots
Rotation of the material in each layer.
FloatMatrixF< 5, 5 > giveStiffnessMatrix_dPdF_PlaneStrain(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
Material * giveInterfaceMaterial(int interface)
FloatMatrixF< 6, 6 > give3dDegeneratedShellStiffMtrx(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const override
void setupLayeredIntegrationRule(std::vector< std::unique_ptr< IntegrationRule > > &layerIntegrationRulesArray, Element *el, int numInPlanePoints)
FloatMatrixF< 3, 3 > give2dPlateSubSoilStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatArray layerThicks
Thickness for each layer.
int giveLayer(GaussPoint *gp) const
double giveLayerMidZ(int layer) const
FloatMatrixF< 4, 4 > giveStiffnessMatrix_PlaneStrain(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 5 > giveGeneralizedStress_Plate(const FloatArrayF< 5 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 8 > giveGeneralizedStress_Shell(const FloatArrayF< 8 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
void giveCharMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override
FloatArrayF< 3 > giveGeneralizedStress_Beam2d(const FloatArrayF< 3 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
double giveLayerThickness(int layer) const
void saveIPContext(DataStream &stream, ContextMode mode, GaussPoint *gp) override
FloatMatrixF< 4, 4 > giveStiffnessMatrix_dPdF_PlaneStress(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 3, 3 > giveStiffnessMatrix_PlaneStress(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
double giveMidSurfaceZcoordFromBottom() const
IntArray layerMaterials
Material of each layer.
FloatMatrixF< 1, 1 > giveStiffnessMatrix_1d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatArray * imposeStrainConstrainsOnGradient(GaussPoint *gp, FloatArray *) override
FloatArrayF< 3 > giveRealStress_PlaneStress(const FloatArrayF< 3 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 4 > giveRealStress_PlaneStrain(const FloatArrayF< 4 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element) override
const char * giveClassName() const override
FloatArrayF< 4 > giveGeneralizedStress_MembraneRot(const FloatArrayF< 4 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
int packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip) override
FloatArrayF< 6 > giveGeneralizedStress_Beam3d(const FloatArrayF< 6 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
IntArray interfacerMaterials
Interface (cohesive zone) material for each interface.
int unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip) override
static MaterialMode giveCorrespondingSlaveMaterialMode(MaterialMode mode)
FloatArrayF< 1 > giveFirstPKStress_1d(const FloatArrayF< 1 > &reducedvF, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 4, 4 > giveMembraneRotStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
void mapLayerGpCoordsToShellCoords(std::vector< std::unique_ptr< IntegrationRule > > &layerIntegrationRulesArray)
double computeIntegralThick() const
Returns the total thickness of all layers.
LayeredCrossSection(int n, Domain *d)
FloatArrayF< 3 > giveGeneralizedStress_PlateSubSoil(const FloatArrayF< 3 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 9, 9 > giveStiffnessMatrix_dPdF_3d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 6, 6 > give3dBeamStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 5, 5 > give2dPlateStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 1 > giveRealStress_1d(const FloatArrayF< 1 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
void initializeFrom(InputRecord &ir) override
double give(CrossSectionProperty a, GaussPoint *gp) const override
int giveSlaveGPIndex(int ilayer, int igp) const
void giveStiffnessMatrix_dCde(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
void restoreIPContext(DataStream &stream, ContextMode mode, GaussPoint *gp) override
void giveCharMaterialStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
void giveCauchyStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep) override
FloatArrayF< 9 > giveFirstPKStress_3d(const FloatArrayF< 9 > &reducedvF, GaussPoint *gp, TimeStep *tStep) const override
bool isCharacteristicMtrxSymmetric(MatResponseMode mode) const override
FloatArrayF< 5 > giveFirstPKStress_PlaneStrain(const FloatArrayF< 5 > &reducedvF, GaussPoint *gp, TimeStep *tStep) const override
int numberOfIntegrationPoints
number of integration points per layer (for 3D elements)
int giveLayerMaterial(int layer) const
FloatArrayF< 6 > giveRealStress_3d(const FloatArrayF< 6 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
FloatArray layerWidths
Width for each layer.
void giveInterfaceXiCoords(FloatArray &answer) const
int giveIPValue(FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep) override
FloatMatrixF< 6, 6 > giveStiffnessMatrix_3d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 8, 8 > give3dShellStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 9, 9 > give3dShellRotStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 6 > giveRealStress_3dDegeneratedShell(const FloatArrayF< 6 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
Material * giveMaterial(IntegrationPoint *ip) const override
hidden by virtual oofem::Material* TransportCrossSection::giveMaterial() const
GaussPoint * giveSlaveGaussPoint(GaussPoint *gp, int layer, int igp) const
int estimatePackSize(DataStream &buff, GaussPoint *ip) override
FloatMatrixF< 1, 1 > giveStiffnessMatrix_dPdF_1d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 4 > giveFirstPKStress_PlaneStress(const FloatArrayF< 4 > &reducedvF, GaussPoint *gp, TimeStep *tStep) const override
const char * giveInputRecordName() const override
FloatArray layerMidZ
z-coord of the mid plane for each layer
LayeredIntegrationRule(int n, Element *e, int startIndx, int endIndx, bool dynamic=false)
const char * giveClassName() const override
int SetUpPointsOnWedge(int nPointsTri, int nPointsDepth, MaterialMode mode) override
#define OOFEM_ERROR(...)
Definition error.h:79
#define _IFT_LayeredCrossSection_Name
long ContextMode
Definition contextmode.h:43
CrossSectionProperty
List of properties possibly stored in a cross section.
GaussPoint IntegrationPoint
Definition gausspoint.h:272

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011