OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
rankinemat.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 rankinemat_h
36 #define rankinemat_h
37 
38 #include "../sm/Materials/structuralmaterial.h"
39 #include "../sm/Materials/structuralms.h"
40 #include "linearelasticmaterial.h"
41 #include "dictionary.h"
42 #include "floatarray.h"
43 #include "floatmatrix.h"
44 
45 // this turns on or off a bunch of internal variables
46 // that allow tracking the distribution of dissipated energy
47 // (can be turned off if such information is not needed)
48 #define keep_track_of_dissipated_energy
49 //#undefine keep_track_of_dissipated_energy
50 
52 
53 #define _IFT_RankineMat_Name "rankmat"
54 #define _IFT_RankineMat_sig0 "sig0"
55 #define _IFT_RankineMat_h "h"
56 #define _IFT_RankineMat_a "a"
57 #define _IFT_RankineMat_plasthardtype "plasthardtype"
58 #define _IFT_RankineMat_delsigy "delsigy"
59 #define _IFT_RankineMat_yieldtol "yieldtol"
60 #define _IFT_RankineMat_gf "gf"
61 #define _IFT_RankineMat_ep "ep"
62 #define _IFT_RankineMat_damlaw "damlaw"
63 #define _IFT_RankineMat_param1 "param1"
64 #define _IFT_RankineMat_param2 "param2"
65 #define _IFT_RankineMat_param3 "param3"
66 #define _IFT_RankineMat_param4 "param4"
67 #define _IFT_RankineMat_param5 "param5"
68 
69 
70 namespace oofem {
71 class RankineMatStatus;
72 
86 {
87 protected:
90 
92  double E;
93 
95  double nu;
96 
98  double H0;
99 
102 
104  double delSigY;
105 
107  double sig0;
108 
110  double yieldtol;
111 
113  double a;
114 
116  double ep;
117 
119  double md;
120 
122  int damlaw;
123 
125  double param1;
126 
128  double param2;
129 
131  double param3;
132 
134  double param4;
135 
137  double param5;
138 
139 
140 public:
141  RankineMat(int n, Domain * d);
142  virtual ~RankineMat() { }
143 
144  double evalYieldFunction(const FloatArray &sigPrinc, const double kappa);
145  double evalYieldStress(const double kappa);
146  double evalPlasticModulus(const double kappa);
147  void performPlasticityReturn(GaussPoint *gp, const FloatArray &totalStrain);
148  double computeDamage(GaussPoint *gp, TimeStep *tStep);
149  double computeDamageParam(double tempKappa);
150  double computeDamageParamPrime(double tempKappa);
151  virtual void computeCumPlastStrain(double &kappa, GaussPoint *gp, TimeStep *tStep);
152 
153  virtual int hasMaterialModeCapability(MaterialMode mode);
154 
156 
157  // identification and auxiliary functions
158  virtual int hasNonLinearBehaviour() { return 1; }
159  virtual const char *giveInputRecordName() const { return _IFT_RankineMat_Name; }
160  virtual const char *giveClassName() const { return "RankineMat"; }
161 
164 
165  virtual bool isCharacteristicMtrxSymmetric(MatResponseMode rMode) { return ( a == 0. ); }
166 
167  virtual MaterialStatus *CreateStatus(GaussPoint *gp) const;
168 
169  virtual void giveRealStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp,
170  const FloatArray &reducesStrain, TimeStep *tStep);
171  virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &totalStrain, TimeStep *tStep);
172 protected:
173  virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep);
174  virtual void givePlaneStressStiffMtrx(FloatMatrix &answer,
175  MatResponseMode mode,
176  GaussPoint *gp,
177  TimeStep *tStep);
183  MatResponseMode mode,
184  GaussPoint *gp,
185  TimeStep *tStep, double gprime);
186 
188  void computeEta(FloatArray &answer, RankineMatStatus *status);
189 
190  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
191 };
192 
193 //=============================================================================
194 
195 
197 {
198 protected:
201 
204 
207 
210 
212  double kappa;
213 
215  double tempKappa;
216 
222  double dKappa1, dKappa2;
223 
225  double damage;
226 
228  double tempDamage;
229 
231  double tanG;
232 
233 #ifdef keep_track_of_dissipated_energy
234  double stressWork;
239  double dissWork;
241  double tempDissWork;
242 #endif
243 
244 public:
245  RankineMatStatus(int n, Domain * d, GaussPoint * g);
246  virtual ~RankineMatStatus();
247 
248  const FloatArray & givePlasticStrain() const { return plasticStrain; }
249 
250  double giveDamage() { return damage; }
251  double giveTempDamage() { return tempDamage; }
252 
253  double giveCumulativePlasticStrain() { return kappa; }
254  double giveTempCumulativePlasticStrain() { return tempKappa; }
255 
256  double giveDKappa(int i)
257  {
258  if ( i == 1 ) {
259  return dKappa1;
260  } else {
261  return dKappa2;
262  }
263  }
264 
266  { return tanG; }
267 
268  const FloatArray &giveEffectiveStress() const { return effStress; }
269  const FloatArray &giveTempEffectiveStress() const { return tempEffStress; }
270 
271  void letTempPlasticStrainBe(FloatArray values) { tempPlasticStrain = std :: move(values); }
272 
273  void letEffectiveStressBe(FloatArray values) { effStress = std :: move(values); }
274 
275  void letTempEffectiveStressBe(FloatArray values) { tempEffStress = std :: move(values); }
276 
277  void setTempCumulativePlasticStrain(double value) { tempKappa = value; }
278 
279  void setDKappa(double val1, double val2) {
280  dKappa1 = val1;
281  dKappa2 = val2;
282  }
283 
284  void setTempDamage(double value) { tempDamage = value; }
285 
286  void setTangentShearStiffness(double value) { tanG = value; }
287 
288  const FloatArray &givePlasDef() { return plasticStrain; }
289 
290  virtual void printOutputAt(FILE *file, TimeStep *tStep);
291 
292  virtual void initTempStatus();
293  virtual void updateYourself(TimeStep *tStep);
294 
295  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
296  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
297 
298 #ifdef keep_track_of_dissipated_energy
299  double giveStressWork() { return stressWork; }
302  double giveTempStressWork() { return tempStressWork; }
304  void setTempStressWork(double w) { tempStressWork = w; }
306  double giveDissWork() { return dissWork; }
308  double giveTempDissWork() { return tempDissWork; }
310  void setTempDissWork(double w) { tempDissWork = w; }
318  void computeWork_PlaneStress(GaussPoint *gp, double gf);
319  void computeWork_1d(GaussPoint *gp, double gf);
320 #endif
321 
322  virtual const char *giveClassName() const { return "RankineMatStatus"; }
323 };
324 } // end namespace oofem
325 #endif // rankinemat_h
FloatArray plasticStrain
Plastic strain (initial).
Definition: rankinemat.h:200
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
This class implements an isotropic elastoplastic material with Rankine yield condition, associated flow rule and linear isotropic softening, and with isotropic damage that leads to softening.
Definition: rankinemat.h:85
virtual void giveRealStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducesStrain, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Definition: rankinemat.C:202
double tempKappa
Cumulative plastic strain (final).
Definition: rankinemat.h:215
Class and object Domain.
Definition: domain.h:115
const FloatArray & giveEffectiveStress() const
Definition: rankinemat.h:268
double tanG
Tangent shear stiffness (needed for tangent matrix).
Definition: rankinemat.h:231
int damlaw
type of damage law (0=exponential, 1=exponential and damage starts after peak stress sig0) ...
Definition: rankinemat.h:122
double giveTempCumulativePlasticStrain()
Definition: rankinemat.h:254
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
double damage
Damage (initial).
Definition: rankinemat.h:225
void evaluatePlaneStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep, double gprime)
Executive method used by local and gradient version.
Definition: rankinemat.C:521
double sig0
Initial (uniaxial) yield stress.
Definition: rankinemat.h:107
const FloatArray & giveTempEffectiveStress() const
Definition: rankinemat.h:269
This class implements a structural material status information.
Definition: structuralms.h:65
virtual bool isCharacteristicMtrxSymmetric(MatResponseMode rMode)
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true...
Definition: rankinemat.h:165
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
Definition: rankinemat.C:650
void setDKappa(double val1, double val2)
Definition: rankinemat.h:279
double yieldtol
Relative tolerance in yield condition.
Definition: rankinemat.h:110
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: femcmpnn.C:51
LinearElasticMaterial * giveLinearElasticMaterial()
Returns a reference to the basic elastic material.
Definition: rankinemat.h:163
double evalYieldStress(const double kappa)
Definition: rankinemat.C:238
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
MatResponseMode
Describes the character of characteristic material matrix.
double giveTangentShearStiffness()
Definition: rankinemat.h:265
double a
Parameter that controls damage evolution (a=0 turns damage off).
Definition: rankinemat.h:113
double ep
Total strain at peak stress sig0–Used only if plasthardtype=2.
Definition: rankinemat.h:116
This class is a abstract base class for all linear elastic material models in a finite element proble...
void setTempDissWork(double w)
Sets the density of dissipated work to given value.
Definition: rankinemat.h:310
double delSigY
Final increment of yield stress (at infinite cumulative plastic strain)
Definition: rankinemat.h:104
FloatArray effStress
Effective stress (initial).
Definition: rankinemat.h:206
double giveTempDissWork()
Returns the density of temp dissipated work.
Definition: rankinemat.h:308
void setTempCumulativePlasticStrain(double value)
Definition: rankinemat.h:277
virtual int hasNonLinearBehaviour()
Returns nonzero if receiver is non linear.
Definition: rankinemat.h:158
void letEffectiveStressBe(FloatArray values)
Definition: rankinemat.h:273
const FloatArray & givePlasticStrain() const
Definition: rankinemat.h:248
const FloatArray & givePlasDef()
Definition: rankinemat.h:288
virtual const char * giveInputRecordName() const
Definition: rankinemat.h:159
double computeDamage(GaussPoint *gp, TimeStep *tStep)
Definition: rankinemat.C:470
virtual MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
Definition: rankinemat.C:166
virtual int hasMaterialModeCapability(MaterialMode mode)
Tests if material supports material mode.
Definition: rankinemat.C:73
double param3
coefficient required when damlaw=2
Definition: rankinemat.h:131
#define _IFT_RankineMat_Name
Definition: rankinemat.h:53
double giveCumulativePlasticStrain()
Definition: rankinemat.h:253
Abstract base class representing a material status information.
Definition: matstatus.h:84
Class representing vector of real numbers.
Definition: floatarray.h:82
void setTempDamage(double value)
Definition: rankinemat.h:284
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
void performPlasticityReturn(GaussPoint *gp, const FloatArray &totalStrain)
Definition: rankinemat.C:283
double param1
coefficient required when damlaw=1 or 2
Definition: rankinemat.h:125
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
int plasthardtype
Type of plastic hardening (0=linear, 1=exponential)
Definition: rankinemat.h:101
virtual void givePlaneStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane stress stiffness matrix of receiver.
Definition: rankinemat.C:492
double tempStressWork
Non-equilibrated density of total work done by stresses on strain increments.
Definition: rankinemat.h:237
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 1d stiffness matrix of receiver.
Definition: rankinemat.C:504
double giveDKappa(int i)
Definition: rankinemat.h:256
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: femcmpnn.h:171
double param4
coefficient required when damlaw=2
Definition: rankinemat.h:134
FloatArray tempEffStress
Effective stress (final).
Definition: rankinemat.h:209
virtual ~RankineMat()
Definition: rankinemat.h:142
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: rankinemat.C:81
LinearElasticMaterial * linearElasticMaterial
Reference to the basic elastic material.
Definition: rankinemat.h:89
double md
Exponent in hardening law–Used only if plasthardtype=2.
Definition: rankinemat.h:119
double computeDamageParamPrime(double tempKappa)
Definition: rankinemat.C:453
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
double E
Young's modulus.
Definition: rankinemat.h:92
double evalYieldFunction(const FloatArray &sigPrinc, const double kappa)
Definition: rankinemat.C:231
Abstract base class for all "structural" constitutive models.
void setTangentShearStiffness(double value)
Definition: rankinemat.h:286
double kappa
Cumulative plastic strain (initial).
Definition: rankinemat.h:212
FloatArray tempPlasticStrain
Plastic strain (final).
Definition: rankinemat.h:203
virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &totalStrain, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Definition: rankinemat.C:176
double computeDamageParam(double tempKappa)
Definition: rankinemat.C:436
double dissWork
Density of dissipated work.
Definition: rankinemat.h:239
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: femcmpnn.C:64
virtual const char * giveClassName() const
Definition: rankinemat.h:322
double param2
coefficient required when damlaw=1 or 2
Definition: rankinemat.h:128
void letTempPlasticStrainBe(FloatArray values)
Definition: rankinemat.h:271
virtual void computeCumPlastStrain(double &kappa, GaussPoint *gp, TimeStep *tStep)
Definition: rankinemat.C:484
RankineMat(int n, Domain *d)
Definition: rankinemat.C:52
the oofem namespace is to define a context or scope in which all oofem names are defined.
double param5
coefficient required when damlaw=2
Definition: rankinemat.h:137
virtual void initTempStatus(GaussPoint *gp)
Initializes temporary variables stored in integration point status at the beginning of new time step...
Definition: material.C:267
void computeEta(FloatArray &answer, RankineMatStatus *status)
Computes derivatives of final kappa with respect to final strain.
Definition: rankinemat.C:612
double evalPlasticModulus(const double kappa)
Definition: rankinemat.C:261
double nu
Poisson's ratio.
Definition: rankinemat.h:95
void setTempStressWork(double w)
Sets the density of total work of stress on strain increments to given value.
Definition: rankinemat.h:304
void letTempEffectiveStressBe(FloatArray values)
Definition: rankinemat.h:275
virtual const char * giveClassName() const
Definition: rankinemat.h:160
double H0
Initial hardening modulus.
Definition: rankinemat.h:98
Class representing integration point in finite element program.
Definition: gausspoint.h:93
double giveTempStressWork()
Returns the temp density of total work of stress on strain increments.
Definition: rankinemat.h:302
Class representing solution step.
Definition: timestep.h:80
double tempDamage
Damage (final).
Definition: rankinemat.h:228
double tempDissWork
Non-equilibrated density of dissipated work.
Definition: rankinemat.h:241
double giveDissWork()
Returns the density of dissipated work.
Definition: rankinemat.h:306

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