OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
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 - 2013 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 
38 #include "../sm/CrossSections/structuralcrosssection.h"
39 #include "floatarray.h"
40 #include "floatmatrix.h"
41 #include "interface.h"
42 #include "gaussintegrationrule.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_initiationlimits "initiationlimits"
60 
61 
62 namespace oofem {
63 class StructuralMaterial;
64 
92 {
93 protected:
104  double totalThick;
105  double area;
106 
107 public:
108  LayeredCrossSection(int n, Domain * d) : StructuralCrossSection(n, d), layerMaterials(), layerThicks(), layerWidths()
109  {
110  numberOfLayers = 0;
111  totalThick = 0.;
112  area = -1.0;
113  }
114 
115  virtual ~LayeredCrossSection() { }
116 
118  virtual void giveInputRecord(DynamicInputRecord &input);
119 
120  virtual void createMaterialStatus(GaussPoint &iGP); // ES
121 
122  //Create slave integration points
123  virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element);
133  virtual int setupIntegrationPoints(IntegrationRule &irule, int npointsXY, int npointsZ, Element *element);
134 
135  virtual void giveRealStress_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep);
136  virtual void giveRealStress_3dDegeneratedShell(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep);
137  virtual void giveRealStress_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep);
138  virtual void giveRealStress_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep);
139  virtual void giveRealStress_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep);
140  virtual void giveRealStress_Warping(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep);
141 
142  virtual void giveStiffnessMatrix_3d(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
143  virtual void giveStiffnessMatrix_PlaneStress(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
144  virtual void giveStiffnessMatrix_PlaneStrain(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
145  virtual void giveStiffnessMatrix_1d(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
146 
147 
148  virtual void giveGeneralizedStress_Beam2d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep);
149  virtual void giveGeneralizedStress_Beam3d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep);
150  virtual void giveGeneralizedStress_Plate(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep);
151  virtual void giveGeneralizedStress_Shell(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep);
152  virtual void giveGeneralizedStress_MembraneRot(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep);
153  virtual void giveGeneralizedStress_PlateSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep);
154 
155  virtual void giveCharMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
156 
158 
159  virtual void give2dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
160  virtual void give3dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
161  virtual void give2dPlateStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
162  virtual void give3dShellStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
163  virtual void give3dDegeneratedShellStiffMtrx(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep);
164  virtual void giveMembraneRotStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
165  virtual void give2dPlateSubSoilStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
166 
169 
170  virtual double give(CrossSectionProperty a, GaussPoint *gp);
171  virtual double give(CrossSectionProperty a, const FloatArray &coords, Element *elem, bool local);
172  int giveNumberOfLayers();
173  int giveLayer(GaussPoint *gp);
174 
176  double computeIntegralThick();
177  void setupLayerMidPlanes();
178 
179  int giveLayerMaterial(int layer) {
180  return this->layerMaterials.at(layer);
181  }
182 
183  virtual Material *giveMaterial(IntegrationPoint *ip);
184 
185  int giveInterfaceMaterialNum(int interface) {
186  return this->interfacerMaterials.at(interface);
187  }
188 
189  Material *giveInterfaceMaterial(int interface) {
190  int matNum = this->giveInterfaceMaterialNum(interface);
191  if ( matNum ) {
192  return this->giveDomain()->giveMaterial( this->interfacerMaterials.at(interface) );
193  } else {
194  return NULL;
195  }
196  }
197 
198  virtual int checkConsistency();
199 
200  double giveLayerMidZ(int layer) {
201  // Gives the z-coord measured from the geometric midplane of the (total) cross section.
202  return this->layerMidZ.at(layer);
203  }
204  double giveLayerThickness(int layer) {
205  return this->layerThicks.at(layer);
206  }
208  return this->numberOfIntegrationPoints;
209  }
211  return this->midSurfaceZcoordFromBottom;
212  }
214  return this->midSurfaceXiCoordFromBottom;
215  }
216  void giveInterfaceXiCoords(FloatArray &answer);
217 
218  // identification and auxiliary functions
219  virtual const char *giveInputRecordName() const { return _IFT_LayeredCrossSection_Name; }
220  virtual const char *giveClassName() const { return "LayeredCrossSection"; }
221  virtual void printYourself();
222 
224  GaussPoint *giveSlaveGaussPoint(GaussPoint *gp, int slaveIndex);
225 
228 
229 
230  void mapLayerGpCoordsToShellCoords(std :: vector< std :: unique_ptr< IntegrationRule > > &layerIntegrationRulesArray);
231 
232  void setupLayeredIntegrationRule(std :: vector< std :: unique_ptr< IntegrationRule > > &layerIntegrationRulesArray, Element *el, int numInPlanePoints);
233 
234  virtual int giveIPValue(FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep);
235  virtual double give(int aProperty, GaussPoint *gp);
236 
237  virtual int packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
238  {
239  OOFEM_ERROR("not implemented");
240  return 0;
241  }
242 
243  virtual int unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
244  {
245  OOFEM_ERROR("not implemented");
246  return 0;
247  }
248 
249  virtual int estimatePackSize(DataStream &buff, GaussPoint *ip)
250  {
251  OOFEM_ERROR("not implemented");
252  return 0;
253  }
254 
255 
256  virtual void giveFirstPKStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep)
257  { OOFEM_ERROR("not implemented"); }
258  virtual void giveCauchyStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep)
259  { OOFEM_ERROR("not implemented"); }
260  virtual void giveStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
261  { OOFEM_ERROR("not implemented"); }
262  virtual void giveStiffnessMatrix_dCde(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
263  { OOFEM_ERROR("not implemented"); }
264 
265 
266 protected:
267  double giveArea();
268 };
269 
274 {
275 public:
277 
287  virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep) = 0;
288 };
289 
291 {
292 public:
293  LayeredIntegrationRule(int n, Element * e, int startIndx, int endIndx, bool dynamic = false);
294  LayeredIntegrationRule(int n, Element * e);
295  virtual ~LayeredIntegrationRule();
296 
297  virtual const char *giveClassName() const { return "LayeredIntegrationRule"; }
299 
300  //virtual int getRequiredNumberOfIntegrationPoints(integrationDomain dType, int approxOrder);
301 
302  // Stores the ip numbers of the points lying on the lower and upper surface of the element.
303  // Thus they will correspond to points lying on the interface between layers.
304  IntArray lowerInterfacePoints, upperInterfacePoints;
305 
306  //virtual int SetUpPointsOnLine(int, MaterialMode); // could be used for beams
307  //virtual int SetUpPointsOnCube(int, MaterialMode mode); // could be used for plates/shells/solids
308  virtual int SetUpPointsOnWedge(int nPointsTri, int nPointsDepth, MaterialMode mode);
309 };
310 } // end namespace oofem
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.
Class and object Domain.
Definition: domain.h:115
int numberOfIntegrationPoints
num integration points per 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...
Definition: datastream.h:54
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
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.
Definition: element.h:145
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.
Definition: crosssection.h:58
#define _IFT_LayeredCrossSection_Name
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
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.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
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.
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.
Definition: domain.C:281
virtual FloatArray * imposeStrainConstrainsOnGradient(GaussPoint *gp, FloatArray *)
Returns modified gradient of strain vector, which is used to compute plastic strain increment...
#define OOFEM_ERROR(...)
Definition: error.h:61
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
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)
virtual contextIOResultType saveIPContext(DataStream &stream, ContextMode mode, GaussPoint *gp)
Stores integration point state to output stream.
Abstract base class for all material models.
Definition: material.h:95
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.
Definition: floatarray.h:82
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.
Definition: floatmatrix.h:94
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
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...
Class representing the general Input Record.
Definition: inputrecord.h:101
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)
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...
Class Interface.
Definition: interface.h:82
virtual void giveGeneralizedStress_Beam3d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
Class representing the a dynamic Input Record.
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.
Definition: contextmode.h:43
Domain * giveDomain() const
Definition: femcmpnn.h:100
virtual void giveRealStress_Warping(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
The element interface required by LayeredCrossSection.
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 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.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
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)

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:29 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011