OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
misesmat.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 misesmat_h
36 #define misesmat_h
37 
38 #include "../sm/Materials/structuralmaterial.h"
39 #include "../sm/Materials/structuralms.h"
41 #include "dictionary.h"
42 #include "floatarray.h"
43 #include "floatmatrix.h"
44 #include "scalarfunction.h"
45 
47 
48 #define _IFT_MisesMat_Name "misesmat"
49 #define _IFT_MisesMat_sig0 "sig0"
50 #define _IFT_MisesMat_h "h"
51 #define _IFT_MisesMat_omega_crit "omega_crit"
52 #define _IFT_MisesMat_a "a"
53 
54 
55 namespace oofem {
56 class GaussPoint;
57 class Domain;
58 
73 {
74 protected:
77 
79  double G;
80 
82  double K;
83 
85  double H;
86 
89 
90  double omega_crit;
91  double a;
92 
93 public:
94  MisesMat(int n, Domain * d);
95  virtual ~MisesMat();
96 
97  void performPlasticityReturn(GaussPoint *gp, const FloatArray &totalStrain, TimeStep *tStep);
98  double computeDamage(GaussPoint *gp, TimeStep *tStep);
99  double computeDamageParam(double tempKappa);
100  double computeDamageParamPrime(double tempKappa);
101  virtual void computeCumPlastStrain(double &kappa, GaussPoint *gp, TimeStep *tStep);
102 
104 
105  virtual int hasNonLinearBehaviour() { return 1; }
106  virtual bool isCharacteristicMtrxSymmetric(MatResponseMode rMode) { return false; }
107 
108  virtual const char *giveInputRecordName() const { return _IFT_MisesMat_Name; }
109  virtual const char *giveClassName() const { return "MisesMat"; }
110 
113 
114  virtual MaterialStatus *CreateStatus(GaussPoint *gp) const;
115 
116  virtual void give3dMaterialStiffnessMatrix(FloatMatrix &answer,
117  MatResponseMode mode,
118  GaussPoint *gp,
119  TimeStep *tStep);
120 
121  virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
122  virtual void give3dMaterialStiffnessMatrix_dPdF(FloatMatrix &answer,
123  MatResponseMode mode,
124  GaussPoint *gp,
125  TimeStep *tStep);
126 
127  virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
128  virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
129 
130  virtual void giveFirstPKStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &vF, TimeStep *tStep);
131 
132  virtual double give(int aProperty, GaussPoint *gp, TimeStep *tStep);
133  double giveTemperature(GaussPoint *gp, TimeStep *tStep);
134 
135 protected:
136  void computeGLPlasticStrain(const FloatMatrix &F, FloatMatrix &Ep, FloatMatrix b, double J);
137 
138  virtual void give3dLSMaterialStiffnessMatrix(FloatMatrix &answer,
139  MatResponseMode mode,
140  GaussPoint *gp,
141  TimeStep *tStep);
142 
143  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
144 };
145 
146 //=============================================================================
147 
148 
150 {
151 protected:
154 
157 
160 
161  /**************************************************/
162  double trialStressV;
163  /**************************************************/
164 
167 
169  double kappa;
170 
172  double tempKappa;
173 
174  /************************/
175  double tempDamage, damage;
176  /******************************/
177 
182 
183 public:
184  MisesMatStatus(int n, Domain * d, GaussPoint * g);
185  virtual ~MisesMatStatus();
186 
187  const FloatArray &givePlasticStrain() { return plasticStrain; }
188 
189  const FloatArray &giveTrialStressDev() { return trialStressD; }
190 
191  /*******************************************/
192  double giveTrialStressVol() { return trialStressV; }
193  /*******************************************/
194  double giveDamage() { return damage; }
195  double giveTempDamage() { return tempDamage; }
196 
197  double giveCumulativePlasticStrain() { return kappa; }
198  double giveTempCumulativePlasticStrain() { return tempKappa; }
199 
200  const FloatMatrix & giveTempLeftCauchyGreen() { return tempLeftCauchyGreen; }
201  const FloatMatrix & giveLeftCauchyGreen() { return leftCauchyGreen; }
202 
203  const FloatArray & giveTempEffectiveStress() { return tempEffStress; }
204  const FloatArray & giveEffectiveStress() { return effStress; }
205 
206  void letTempPlasticStrainBe(const FloatArray &values) { tempPlasticStrain = values; }
207  const FloatArray &getTempPlasticStrain() const { return tempPlasticStrain; }
208 
209  void letTrialStressDevBe(const FloatArray &values) { trialStressD = values; }
210 
211  void letEffectiveStressBe(const FloatArray &values) { effStress = values; }
212 
213  void letTempEffectiveStressBe(const FloatArray &values) { tempEffStress = values; }
214 
215 
216  void setTrialStressVol(double value) { trialStressV = value; }
217 
218  void setTempCumulativePlasticStrain(double value) { tempKappa = value; }
219  /****************************************/
220  void setTempDamage(double value) { tempDamage = value; }
221  /************************************************/
222 
223  void letTempLeftCauchyGreenBe(const FloatMatrix &values) { tempLeftCauchyGreen = values; }
224  void letLeftCauchyGreenBe(const FloatMatrix &values) { leftCauchyGreen = values; }
225 
226  const FloatArray &givePlasDef() { return plasticStrain; }
227 
228  virtual void printOutputAt(FILE *file, TimeStep *tStep);
229 
230  virtual void initTempStatus();
231 
232  virtual void updateYourself(TimeStep *tStep);
233 
234  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
235  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
236 
237  virtual const char *giveClassName() const { return "MisesMatStatus"; }
238 };
239 } // end namespace oofem
240 #endif // misesmat_h
void letTempLeftCauchyGreenBe(const FloatMatrix &values)
Definition: misesmat.h:223
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
FloatArray plasticStrain
Plastic strain (initial).
Definition: misesmat.h:153
This class implements an isotropic elastoplastic material with Mises yield condition, associated flow rule and linear isotropic hardening.
Definition: misesmat.h:72
virtual const char * giveClassName() const
Definition: misesmat.h:109
Class and object Domain.
Definition: domain.h:115
void letLeftCauchyGreenBe(const FloatMatrix &values)
Definition: misesmat.h:224
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 ...
Definition: misesmat.C:424
FloatArray tempPlasticStrain
Plastic strain (final).
Definition: misesmat.h:156
void letEffectiveStressBe(const FloatArray &values)
Definition: misesmat.h:211
double H
Hardening modulus.
Definition: misesmat.h:85
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
virtual bool isCharacteristicMtrxSymmetric(MatResponseMode rMode)
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true...
Definition: misesmat.h:106
MisesMat(int n, Domain *d)
Definition: misesmat.C:54
virtual void giveFirstPKStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &vF, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Definition: misesmat.C:165
This class implements a structural material status information.
Definition: structuralms.h:65
void letTempEffectiveStressBe(const FloatArray &values)
Definition: misesmat.h:213
double giveCumulativePlasticStrain()
Definition: misesmat.h:197
double computeDamage(GaussPoint *gp, TimeStep *tStep)
Definition: misesmat.C:391
double giveTempCumulativePlasticStrain()
Definition: misesmat.h:198
double giveTempDamage()
Definition: misesmat.h:195
void letTrialStressDevBe(const FloatArray &values)
Definition: misesmat.h:209
const FloatMatrix & giveTempLeftCauchyGreen()
Definition: misesmat.h:200
#define _IFT_MisesMat_Name
Definition: misesmat.h:48
const FloatArray & givePlasDef()
Definition: misesmat.h:226
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: femcmpnn.C:51
virtual const char * giveClassName() const
Definition: misesmat.h:237
FloatArray effStress
Definition: misesmat.h:165
double K
Elastic bulk modulus.
Definition: misesmat.h:82
MatResponseMode
Describes the character of characteristic material matrix.
double computeDamageParam(double tempKappa)
Definition: misesmat.C:369
virtual int hasNonLinearBehaviour()
Returns nonzero if receiver is non linear.
Definition: misesmat.h:105
This class is a abstract base class for all linear elastic material models in a finite element proble...
double giveTrialStressVol()
Definition: misesmat.h:192
FloatArray trialStressD
Deviatoric trial stress - needed for tangent stiffness.
Definition: misesmat.h:159
void setTempCumulativePlasticStrain(double value)
Definition: misesmat.h:218
void performPlasticityReturn(GaussPoint *gp, const FloatArray &totalStrain, TimeStep *tStep)
Definition: misesmat.C:294
FloatMatrix tempLeftCauchyGreen
Left Cauchy-Green deformation gradient (final).
Definition: misesmat.h:181
const FloatArray & giveTrialStressDev()
Definition: misesmat.h:189
virtual void computeCumPlastStrain(double &kappa, GaussPoint *gp, TimeStep *tStep)
Definition: misesmat.C:406
LinearElasticMaterial * linearElasticMaterial
Reference to the basic elastic material.
Definition: misesmat.h:76
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
Definition: misesmat.C:692
double computeDamageParamPrime(double tempKappa)
Definition: misesmat.C:379
virtual void give3dLSMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: misesmat.C:511
double omega_crit
Definition: misesmat.h:90
void setTempDamage(double value)
Definition: misesmat.h:220
Abstract base class representing a material status information.
Definition: matstatus.h:84
Class representing vector of real numbers.
Definition: floatarray.h:82
double kappa
Cumulative plastic strain (initial).
Definition: misesmat.h:169
const FloatArray & getTempPlasticStrain() const
Definition: misesmat.h:207
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual double give(int aProperty, GaussPoint *gp, TimeStep *tStep)
Definition: misesmat.C:828
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
Implementation of Scalar function.
const FloatMatrix & giveLeftCauchyGreen()
Definition: misesmat.h:201
Class representing the general Input Record.
Definition: inputrecord.h:101
void setTrialStressVol(double value)
Definition: misesmat.h:216
const FloatArray & giveTempEffectiveStress()
Definition: misesmat.h:203
virtual void give3dMaterialStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: misesmat.C:414
virtual MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
Definition: misesmat.C:101
virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Definition: misesmat.C:107
double giveTemperature(GaussPoint *gp, TimeStep *tStep)
Definition: misesmat.C:840
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: femcmpnn.h:171
const FloatArray & givePlasticStrain()
Definition: misesmat.h:187
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
virtual ~MisesMat()
Definition: misesmat.C:64
Abstract base class for all "structural" constitutive models.
void letTempPlasticStrainBe(const FloatArray &values)
Definition: misesmat.h:206
const FloatArray & giveEffectiveStress()
Definition: misesmat.h:204
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: femcmpnn.C:64
virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 1d stiffness matrix of receiver.
Definition: misesmat.C:480
the oofem namespace is to define a context or scope in which all oofem names are defined.
double G
Elastic shear modulus.
Definition: misesmat.h:79
double tempKappa
Cumulative plastic strain (final).
Definition: misesmat.h:172
virtual void initTempStatus(GaussPoint *gp)
Initializes temporary variables stored in integration point status at the beginning of new time step...
Definition: material.C:267
virtual const char * giveInputRecordName() const
Definition: misesmat.h:108
ScalarFunction sig0
Initial (uniaxial) yield stress.
Definition: misesmat.h:88
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Definition: misesmat.C:143
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: misesmat.C:71
FloatArray tempEffStress
Definition: misesmat.h:166
void computeGLPlasticStrain(const FloatMatrix &F, FloatMatrix &Ep, FloatMatrix b, double J)
Definition: misesmat.C:279
FloatMatrix leftCauchyGreen
Left Cauchy-Green deformation gradient (initial).
Definition: misesmat.h:179
LinearElasticMaterial * giveLinearElasticMaterial()
Returns a reference to the basic elastic material.
Definition: misesmat.h:112

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