OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
structuralmaterial.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 structuralmaterial_h
36 #define structuralmaterial_h
37 
38 #include "material.h"
39 #include "floatarray.h"
40 #include "floatmatrix.h"
41 #include "matconst.h"
42 #include "matstatus.h"
43 #include "stressstrainprincmode.h"
44 #include "valuemodetype.h"
45 #include <vector>
46 
48 
49 #define _IFT_StructuralMaterial_referencetemperature "referencetemperature"
50 #define _IFT_StructuralMaterial_talpha "talpha"
51 
52 
53 namespace oofem {
54 #define STRAIN_STEPS 10.0
55 
56 class GaussPoint;
58 
114 {
115 protected:
118 
119 public:
121  static std::vector< std::vector<int> > vIindex;
122 
124  static std::vector< std::vector<int> > svIndex;
125 
126  static int giveSymVI(int ind1, int ind2) { return svIndex[ind1-1][ind2-1]; }
127  static int giveVI(int ind1, int ind2) { return vIindex[ind1-1][ind2-1]; }
128 
134  StructuralMaterial(int n, Domain *d);
136  virtual ~StructuralMaterial() { }
137 
138  virtual int hasMaterialModeCapability(MaterialMode mode);
139  virtual const char *giveClassName() const { return "StructuralMaterial"; }
140 
142  virtual void giveInputRecord(DynamicInputRecord &input);
143 
153  virtual void giveStiffnessMatrix(FloatMatrix &answer,
154  MatResponseMode mode,
155  GaussPoint *gp,
156  TimeStep *tStep);
157 
171  virtual void giveRealStressVector(FloatArray &answer, GaussPoint *gp,
172  const FloatArray &reducedStrain, TimeStep *tStep);
174  virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
176  virtual void giveRealStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
178  virtual void giveRealStressVector_StressControl(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep);
179  virtual void giveRealStressVector_ShellStressControl(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep);
181  virtual void giveRealStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
183  virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
185  virtual void giveRealStressVector_Warping(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
187  virtual void giveRealStressVector_2dBeamLayer(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
189  virtual void giveRealStressVector_PlateLayer(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
191  virtual void giveRealStressVector_Fiber(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
192  virtual void giveRealStressVector_Lattice2d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
193  virtual void giveRealStressVector_Lattice3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
195  virtual void giveRealStressVector_2dPlateSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
196  virtual void giveRealStressVector_3dBeamSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
197 
216  virtual void giveFirstPKStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep);
219  virtual void giveFirstPKStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep);
221  virtual void giveFirstPKStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep);
223  virtual void giveFirstPKStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep);
225 
245  virtual void giveCauchyStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
246  { OOFEM_ERROR("not implemented "); }
247  virtual void giveCauchyStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
248  { OOFEM_ERROR("not implemented "); }
249  virtual void giveCauchyStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
250  { OOFEM_ERROR("not implemented "); }
251  virtual void giveCauchyStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
252  { OOFEM_ERROR("not implemented "); }
254 
258  virtual void giveEshelbyStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep);
259 
260  void give_dPdF_from(const FloatMatrix &dSdE, FloatMatrix &answer, GaussPoint *gp);
261  void convert_dSdE_2_dPdF(FloatMatrix &answer, const FloatMatrix &dSdE, const FloatArray &S, const FloatArray &F, MaterialMode matMode);
262 
269  virtual void giveThermalDilatationVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
274 
285  GaussPoint *gp, TimeStep *tStep, ValueModeType mode);
287  GaussPoint *gp, TimeStep *tStep, ValueModeType mode);
289 
290 
296  static void computePrincipalValues(FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode);
304  static void computePrincipalValDir(FloatArray &answer, FloatMatrix &dir, const FloatArray &s,
305  stressStrainPrincMode mode);
306 
313  static double computeDeviatoricVolumetricSplit(FloatArray &dev, const FloatArray &s);
314  static void computeDeviatoricVolumetricSum(FloatArray &s, const FloatArray &dev, double mean);
315 
316  static void applyDeviatoricElasticCompliance(FloatArray &strain, const FloatArray &stress, double EModulus, double nu);
317  static void applyDeviatoricElasticCompliance(FloatArray &strain, const FloatArray &stress, double GModulus);
318 
319  static void applyDeviatoricElasticStiffness(FloatArray &stress, const FloatArray &strain, double EModulus, double nu);
320  static void applyDeviatoricElasticStiffness(FloatArray &stress, const FloatArray &strain, double GModulus);
321 
322  static void applyElasticStiffness(FloatArray &stress, const FloatArray &strain, double EModulus, double nu);
323  static void applyElasticCompliance(FloatArray &strain, const FloatArray &stress, double EModulus, double nu);
324 
325  static double computeStressNorm(const FloatArray &stress);
326 
327  static double computeFirstInvariant(const FloatArray &s);
328  static double computeSecondStressInvariant(const FloatArray &s);
329  static double computeThirdStressInvariant(const FloatArray &s);
330 
331  static double computeFirstCoordinate(const FloatArray &s);
332  static double computeSecondCoordinate(const FloatArray &s);
333  static double computeThirdCoordinate(const FloatArray &s);
335 
345  MatResponseMode mode,
346  GaussPoint *gp,
347  TimeStep *tStep)
348  { OOFEM_ERROR("not implemented "); }
349 
350 
351  virtual void give3dMaterialStiffnessMatrix_dPdF(FloatMatrix &answer,
352  MatResponseMode mode,
353  GaussPoint *gp, TimeStep *tStep);
354 
355  virtual void give3dMaterialStiffnessMatrix_dCde(FloatMatrix &answer,
356  MatResponseMode mode,
357  GaussPoint *gp, TimeStep *tStep);
358 
373  static int giveVoigtVectorMask(IntArray &answer, MaterialMode mmode);
374 
392  static int giveVoigtSymVectorMask(IntArray &answer, MaterialMode mmode);
397  static void giveInvertedVoigtVectorMask(IntArray &answer, MaterialMode mmode);
398 
403  static int giveSizeOfVoigtVector(MaterialMode mmode);
408  static int giveSizeOfVoigtSymVector(MaterialMode mmode);
409 
411  static void giveFullVectorForm(FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode);
413  static void giveFullVectorFormF(FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode);
415  static void giveFullSymVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode);
417  static void giveReducedVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode);
419  static void giveReducedSymVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode);
420 
422  static void giveFullSymMatrixForm(FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode);
424  static void giveReducedMatrixForm(FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode);
426  static void giveReducedSymMatrixForm(FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode);
427 
439  void giveStressDependentPartOfStrainVector(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector,
440  TimeStep *tStep, ValueModeType mode);
441  void giveStressDependentPartOfStrainVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector,
442  TimeStep *tStep, ValueModeType mode);
443 
444  virtual int setIPValue(const FloatArray &value, GaussPoint *gp, InternalStateType type);
445  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
446 
460  virtual void givePlaneStressStiffMtrx(FloatMatrix &answer,
461  MatResponseMode mmode, GaussPoint *gp,
462  TimeStep *tStep);
463 
464  virtual void givePlaneStressStiffMtrx_dPdF(FloatMatrix &answer,
465  MatResponseMode mmode, GaussPoint *gp,
466  TimeStep *tStep);
467 
468  virtual void givePlaneStressStiffMtrx_dCde(FloatMatrix &answer,
469  MatResponseMode mmode, GaussPoint *gp,
470  TimeStep *tStep);
472 
491  virtual void givePlaneStrainStiffMtrx(FloatMatrix &answer,
492  MatResponseMode mmode, GaussPoint *gp,
493  TimeStep *tStep);
494 
495  virtual void givePlaneStrainStiffMtrx_dPdF(FloatMatrix &answer,
496  MatResponseMode mmode, GaussPoint *gp,
497  TimeStep *tStep);
498 
499  virtual void givePlaneStrainStiffMtrx_dCde(FloatMatrix &answer,
500  MatResponseMode mmode, GaussPoint *gp,
501  TimeStep *tStep);
503 
517  virtual void give1dStressStiffMtrx(FloatMatrix &answer,
518  MatResponseMode mmode, GaussPoint *gp,
519  TimeStep *tStep);
520 
521  virtual void give1dStressStiffMtrx_dPdF(FloatMatrix &answer,
522  MatResponseMode mmode, GaussPoint *gp,
523  TimeStep *tStep);
524 
525  virtual void give1dStressStiffMtrx_dCde(FloatMatrix &answer,
526  MatResponseMode mmode, GaussPoint *gp,
527  TimeStep *tStep);
529 
542  virtual void give2dBeamLayerStiffMtrx(FloatMatrix &answer,
543  MatResponseMode mmode, GaussPoint *gp,
544  TimeStep *tStep);
557  virtual void givePlateLayerStiffMtrx(FloatMatrix &answer,
558  MatResponseMode mmode, GaussPoint *gp,
559  TimeStep *tStep);
572  virtual void giveFiberStiffMtrx(FloatMatrix &answer,
573  MatResponseMode mmode, GaussPoint *gp,
574  TimeStep *tStep);
575 
576 
584  virtual void give2dLatticeStiffMtrx(FloatMatrix &answer,
585  MatResponseMode mmode, GaussPoint *gp,
586  TimeStep *tStep);
587 
595  virtual void give3dLatticeStiffMtrx(FloatMatrix &answer,
596  MatResponseMode mmode, GaussPoint *gp,
597  TimeStep *tStep);
598 
599 
608  virtual void give2dPlateSubSoilStiffMtrx(FloatMatrix &answer,
609  MatResponseMode mmode, GaussPoint *gp,
610  TimeStep *tStep);
619  virtual void give3dBeamSubSoilStiffMtrx(FloatMatrix &answer,
620  MatResponseMode mmode, GaussPoint *gp,
621  TimeStep *tStep);
631  static void transformStrainVectorTo(FloatArray &answer, const FloatMatrix &base,
632  const FloatArray &strainVector, bool transpose = false);
642  static void transformStressVectorTo(FloatArray &answer, const FloatMatrix &base,
643  const FloatArray &stressVector, bool transpose = false);
644 
649  static double computeVonMisesStress(const FloatArray *currentStress);
650 
658  static void giveStrainVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base,
659  bool transpose = false);
660 
668  static void give2DStrainVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base,
669  bool transpose = false);
670 
678  static void giveStressVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base,
679  bool transpose = false);
687  static void givePlaneStressVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base,
688  bool transpose = false);
698  static void sortPrincDirAndValCloseTo(FloatArray *pVal, FloatMatrix *pDir, FloatMatrix *toPDir);
699 
700  friend class CrossSection;
702  friend class SimpleCrossSection;
703  friend class LayeredCrossSection;
704 };
705 } // end namespace oofem
706 #endif // structuralmaterial_h
static int giveSizeOfVoigtSymVector(MaterialMode mmode)
Returns the size of symmetric part of a reduced stress/strain vector according to given mode...
virtual void giveFirstPKStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Default implementation relies on giveFirstPKStressVector_3d.
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual void giveEshelbyStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Prototype for computation of Eshelby stress.
static void computePrincipalValues(FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode)
Common functions for convenience.
static int giveVoigtSymVectorMask(IntArray &answer, MaterialMode mmode)
The same as giveVoigtVectorMask but returns a mask corresponding to a symmetric second order tensor...
static void givePlaneStressVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
Computes 2d stress vector transformation matrix from standard vector transformation matrix...
virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 1d stiffness matrix of receiver.
Class and object Domain.
Definition: domain.h:115
static void applyDeviatoricElasticCompliance(FloatArray &strain, const FloatArray &stress, double EModulus, double nu)
virtual void giveFirstPKStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
virtual void giveRealStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_3d.
virtual void give1dStressStiffMtrx_dPdF(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual void giveFirstPKStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Default implementation relies on giveFirstPKStressVector_3d.
static int giveVI(int ind1, int ind2)
virtual void giveFiberStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 1d fiber stiffness matrix of receiver.
static double computeSecondCoordinate(const FloatArray &s)
virtual void giveCauchyStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
static void computeDeviatoricVolumetricSum(FloatArray &s, const FloatArray &dev, double mean)
virtual void giveCauchyStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Defines several material constant (respective their representative number).
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
virtual void givePlaneStrainStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane strain stiffness matrix of receiver.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual void computeStressIndependentStrainVector_3d(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode)
static double computeFirstInvariant(const FloatArray &s)
virtual void giveRealStressVector_Lattice2d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
static double computeVonMisesStress(const FloatArray *currentStress)
Computes equivalent of von Mises stress.
static void giveFullSymMatrixForm(FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode)
Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form.
double giveReferenceTemperature()
Returns the reference temperature of receiver.
virtual void give3dMaterialStiffnessMatrix_dCde(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
static int giveSizeOfVoigtVector(MaterialMode mmode)
Returns the size of reduced stress/strain vector according to given mode.
static void applyElasticCompliance(FloatArray &strain, const FloatArray &stress, double EModulus, double nu)
virtual void givePlateLayerStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d plate layer stiffness matrix of receiver.
static void sortPrincDirAndValCloseTo(FloatArray *pVal, FloatMatrix *pDir, FloatMatrix *toPDir)
Method for sorting newly computed principal values (pVal) and corresponding principal directions (pDi...
static std::vector< std::vector< int > > svIndex
Symmetric Voigt index map.
static void giveStrainVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
Computes 3d strain vector transformation matrix from standard vector transformation matrix...
#define S(p)
Definition: mdm.C:481
void give_dPdF_from(const FloatMatrix &dSdE, FloatMatrix &answer, GaussPoint *gp)
void giveStressDependentPartOfStrainVector(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode)
Method for subtracting from reduced space strain vector its stress-independent parts (caused by tempe...
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual void giveRealStressVector_3dBeamSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
static void computePrincipalValDir(FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode)
Computes principal values and directions of stress or strain vector.
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
This class implements a layered cross section in a finite element problem.
virtual void giveRealStressVector(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
Computes the real stress vector for given total strain and integration point.
virtual ~StructuralMaterial()
Destructor.
static void giveFullSymVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form.
Base abstract class representing cross section in finite element mesh.
Definition: crosssection.h:107
virtual void givePlaneStressStiffMtrx_dPdF(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual void give3dLatticeStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 3d lattice stiffness matrix of receiver.
virtual void givePlaneStrainStiffMtrx_dCde(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual void giveRealStressVector_2dPlateSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation is not provided.
virtual void giveCauchyStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual void giveRealStressVector_2dBeamLayer(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void giveStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix for giveRealStressVector of receiver in given integration point...
static void giveFullVectorForm(FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode)
Converts the reduced symmetric Voigt vector (2nd order tensor) to full form.
static void giveReducedSymVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form.
static void giveReducedVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
Converts the full symmetric Voigt vector (2nd order tensor) to reduced form.
virtual void giveRealStressVector_Lattice3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
static std::vector< std::vector< int > > vIindex
Voigt index map.
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual void givePlaneStressStiffMtrx_dCde(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
static void applyElasticStiffness(FloatArray &stress, const FloatArray &strain, double EModulus, double nu)
virtual void give3dMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes full 3d material stiffness matrix at given integration point, time, respecting load history ...
virtual int setIPValue(const FloatArray &value, GaussPoint *gp, InternalStateType type)
Sets the value of a certain variable at a given integration point to the given value.
static double computeFirstCoordinate(const FloatArray &s)
void giveStressDependentPartOfStrainVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode)
static double computeStressNorm(const FloatArray &stress)
static int giveVoigtVectorMask(IntArray &answer, MaterialMode mmode)
Returns a mask of the vector indicies corresponding to components in a general (non-symmetric) second...
Abstract base class for all material models.
Definition: material.h:95
virtual void giveRealStressVector_ShellStressControl(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep)
virtual void givePlaneStrainStiffMtrx_dPdF(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual void giveRealStressVector_StressControl(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep)
Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· ...
virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
static double computeThirdCoordinate(const FloatArray &s)
virtual void giveFirstPKStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Default implementation relies on giveFirstPKStressVector_3d.
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void giveRealStressVector_Warping(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
void convert_dSdE_2_dPdF(FloatMatrix &answer, const FloatMatrix &dSdE, const FloatArray &S, const FloatArray &F, MaterialMode matMode)
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
stressStrainPrincMode
We have only one algorithm for computing eigenvalues and vectors in order to be able to distinguish b...
static void giveStressVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
Computes 3d stress vector transformation matrix from standard vector transformation matrix...
Class implementing "simple" cross section model in finite element problem.
virtual void give2dLatticeStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d lattice stiffness matrix of receiver.
static void transformStrainVectorTo(FloatArray &answer, const FloatMatrix &base, const FloatArray &strainVector, bool transpose=false)
Transforms 3d strain vector into another coordinate system.
double referenceTemperature
Reference temperature (temperature, when material has been built into structure). ...
virtual void giveRealStressVector_Fiber(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
virtual int hasMaterialModeCapability(MaterialMode mode)
Tests if material supports material mode.
static void applyDeviatoricElasticStiffness(FloatArray &stress, const FloatArray &strain, double EModulus, double nu)
static double computeSecondStressInvariant(const FloatArray &s)
virtual void giveThermalDilatationVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Returns a vector of coefficients of thermal dilatation in direction of each material principal (local...
Class representing the general Input Record.
Definition: inputrecord.h:101
static void giveReducedMatrixForm(FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode)
Converts the full symmetric Voigt matrix (4th order tensor) to reduced form.
virtual void give3dBeamSubSoilStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing stiffness matrix of beam3d subsoil model.
Class representing the a dynamic Input Record.
static double computeDeviatoricVolumetricSplit(FloatArray &dev, const FloatArray &s)
Computes split of receiver into deviatoric and volumetric part.
virtual void give2dPlateSubSoilStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing stiffness matrix of plate subsoil model.
virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Abstract base class for all "structural" constitutive models.
static void giveInvertedVoigtVectorMask(IntArray &answer, MaterialMode mmode)
Gives the inverted version of giveVoigtVectorMask.
virtual void giveRealStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
virtual void give2dBeamLayerStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d beam layer stiffness matrix of receiver.
static void give2DStrainVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
Computes 2d strain vector transformation matrix from standard vector transformation matrix...
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
static int giveSymVI(int ind1, int ind2)
virtual const char * giveClassName() const
virtual void give3dMaterialStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
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 giveRealStressVector_PlateLayer(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
static void giveFullVectorFormF(FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode)
Converts the reduced deformation gradient Voigt vector (2nd order tensor).
virtual void givePlaneStressStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane stress stiffness matrix of receiver.
static void transformStressVectorTo(FloatArray &answer, const FloatMatrix &base, const FloatArray &stressVector, bool transpose=false)
Transforms 3d stress vector into another coordinate system.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
static double computeThirdStressInvariant(const FloatArray &s)
virtual void computeStressIndependentStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode)
Computes reduced strain vector in given integration point, generated by internal processes in materia...
virtual void give1dStressStiffMtrx_dCde(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
static void giveReducedSymMatrixForm(FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode)
Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form.
virtual void giveCauchyStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
StructuralMaterial(int n, Domain *d)
Constructor.

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:31 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011