OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
rankinemat.C
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 #include "rankinemat.h"
37 #include "gausspoint.h"
38 #include "floatmatrix.h"
39 #include "floatarray.h"
40 #include "intarray.h"
41 #include "stressvector.h"
42 #include "strainvector.h"
43 #include "mathfem.h"
44 #include "contextioerr.h"
45 #include "datastream.h"
46 #include "classfactory.h"
47 
48 namespace oofem {
49 REGISTER_Material(RankineMat);
50 
51 // constructor
53 {
55  E = 0.;
56  nu = 0.;
57  H0 = 0.;
58  sig0 = 0.;
59  delSigY = 0.;
60  ep = 0.;
61  md = 0.;
62  damlaw = 1;
63  param1 = 0.;
64  param2 = 0.;
65  param3 = 0.;
66  param4 = 0.;
67  param5 = 0.;
68 }
69 
70 
71 // specifies whether a given material mode is supported by this model
72 int
74 {
75  return ( ( mode == _PlaneStress ) || ( mode == _1dMat ) );
76 }
77 
78 
79 // reads the model parameters from the input file
82 {
83  IRResultType result; // required by IR_GIVE_FIELD macro
84 
86  if ( result != IRRT_OK ) return result;
87 
88  result = linearElasticMaterial->initializeFrom(ir); // takes care of elastic constants
89  if ( result != IRRT_OK ) return result;
90 
91  E = static_cast< IsotropicLinearElasticMaterial * >(linearElasticMaterial)->giveYoungsModulus();
92  nu = static_cast< IsotropicLinearElasticMaterial * >(linearElasticMaterial)->givePoissonsRatio();
93 
94  IR_GIVE_FIELD(ir, sig0, _IFT_RankineMat_sig0); // uniaxial yield stress
95 
96  H0 = 0.;
97  IR_GIVE_OPTIONAL_FIELD(ir, H0, _IFT_RankineMat_h); // hardening modulus
98 
99  plasthardtype = 0;
100  IR_GIVE_OPTIONAL_FIELD(ir, plasthardtype, _IFT_RankineMat_plasthardtype); // type of plastic hardening (0=linear, 1=exponential, 2= prepeak hardening + linear softening )
101 
102  delSigY = 0.;
103  if ( plasthardtype == 0 ) {
104  //no extra required variables
105  } else if ( plasthardtype == 1 ) {
106  IR_GIVE_FIELD(ir, delSigY, _IFT_RankineMat_delsigy); // final increment of yield stress (at infinite cumulative plastic strain)
107  } else if ( plasthardtype == 2 ) {
109  ep = ep - sig0 / E; // user input is strain at peak stress sig0 and is converted to plastic strain at peak stress sig0
110  md = 1. / log(50. * E * ep / sig0); // exponent used on the 1st plasticity branch
111  } else {
112  OOFEM_WARNING("Plasticity hardening type number %d is unknown", plasthardtype);
113  return IRRT_BAD_FORMAT;
114  }
115 
116  yieldtol = 1.e-10;
117  IR_GIVE_OPTIONAL_FIELD(ir, yieldtol, _IFT_RankineMat_yieldtol); // relative tolerance in yield condition
118 
119  damlaw = 0;
120  IR_GIVE_OPTIONAL_FIELD(ir, damlaw, _IFT_RankineMat_damlaw); // type of damage law (0=exponential, 1=exponential and damage starts after peak stress sig0)
121 
122  if ( damlaw == 0 ) {
123  a = 0.;
124  IR_GIVE_OPTIONAL_FIELD(ir, a, _IFT_RankineMat_a); // coefficient in damage law
125  } else if ( damlaw == 1 ) {
126  IR_GIVE_FIELD(ir, param1, _IFT_RankineMat_param1); // coefficient in damage law
127  IR_GIVE_FIELD(ir, param2, _IFT_RankineMat_param2); // coefficient in damage law. If b<1 use only stiffmode=1
128  } else if ( damlaw == 2 ) {
129  IR_GIVE_FIELD(ir, param1, _IFT_RankineMat_param1); // coefficients in damage law
134  } else {
135  OOFEM_WARNING("Damage law number %d is unknown", damlaw);
136  return IRRT_BAD_FORMAT;
137  }
138 
139  double gf = 0.;
140  IR_GIVE_OPTIONAL_FIELD(ir, gf, _IFT_RankineMat_gf); // dissipated energy per unit VOLUME
141 
142  if ( ( a != 0. ) && ( gf != 0 ) ) {
143  OOFEM_WARNING("parameters a and gf cannot be prescribed simultaneously");
144  return IRRT_BAD_FORMAT;
145  }
146 
147  if ( gf > 0. ) {
148  // evaluate parameter "a" from given "gf"
149  double A = H0 * ( 1. + H0 / E );
150  double B = sig0 * ( 1. + H0 / E );
151  double C = 0.5 * sig0 * sig0 / E - gf;
152  if ( C >= 0. ) {
153  OOFEM_ERROR("parameter gf is too low");
154  }
155 
156  double kappaf = ( -B + sqrt(B * B - 4. * A * C) ) / ( 2. * A );
157  a = 1. / kappaf;
158  }
159 
160  return IRRT_OK;
161 }
162 
163 
164 // creates a new material status corresponding to this class
167 {
168  RankineMatStatus *status;
169  status = new RankineMatStatus(1, this->giveDomain(), gp);
170  return status;
171 }
172 
173 
174 // computes the stress vector corresponding to given (final) strain
175 void
177 {
178  RankineMatStatus *status = static_cast< RankineMatStatus * >( this->giveStatus(gp) );
179 
180  // initialization
181  this->initTempStatus(gp);
182 
183  // elastoplasticity
184  this->performPlasticityReturn(gp, totalStrain);
185 
186  // damage
187  double omega = computeDamage(gp, tStep);
188  answer.beScaled(1. - omega, status->giveTempEffectiveStress());
189 
190  // store variables in status
191  status->setTempDamage(omega);
192  status->letTempStrainVectorBe(totalStrain);
193  status->letTempStressVectorBe(answer);
194 #ifdef keep_track_of_dissipated_energy
195  double gf = sig0 * sig0 / E; // only estimated, but OK for this purpose
196  status->computeWork_1d(gp, gf);
197 #endif
198 }
199 
200 
201 void
203  GaussPoint *gp,
204  const FloatArray &totalStrain,
205  TimeStep *tStep)
206 {
207  RankineMatStatus *status = static_cast< RankineMatStatus * >( this->giveStatus(gp) );
208 
209  // initialization
210  this->initTempStatus(gp);
211 
212  // elastoplasticity
213  this->performPlasticityReturn(gp, totalStrain);
214 
215  // damage
216  double omega = computeDamage(gp, tStep);
217  answer.beScaled(1. - omega, status->giveTempEffectiveStress());
218 
219  // store variables in status
220  status->setTempDamage(omega);
221  status->letTempStrainVectorBe(totalStrain);
222  status->letTempStressVectorBe(answer);
223 #ifdef keep_track_of_dissipated_energy
224  double gf = sig0 * sig0 / E; // only estimated, but OK for this purpose
225  status->computeWork_PlaneStress(gp, gf);
226 #endif
227 }
228 
229 
230 double
231 RankineMat :: evalYieldFunction(const FloatArray &sigPrinc, const double kappa)
232 {
233  return sigPrinc.at(1) - evalYieldStress(kappa);
234 }
235 
236 
237 double
238 RankineMat :: evalYieldStress(const double kappa)
239 {
240  double yieldStress = 0.;
241  if ( plasthardtype == 0 ) { // linear hardening
242  yieldStress = sig0 + H0 * kappa;
243  } else if ( plasthardtype == 1 ) { // exponential hardening
244  if ( delSigY == 0. ) {
245  yieldStress = sig0;
246  } else {
247  yieldStress = sig0 + delSigY * ( 1. - exp(-H0 * kappa / delSigY) );
248  }
249  } else if ( plasthardtype == 2 ) { // exponential hardening before the peak stress sig0 and linear after the peak stress sig0
250  if ( kappa <= ep ) {
251  //1st branch in the rankine variation 2 trying to match the 1st branch of the smooth extended damage law reported by Grassl and Jirasek (2010)
252  yieldStress = 50. *E *kappa *exp(-pow ( kappa / ep, md ) / md);
253  } else { //linear hardening branch
254  yieldStress = sig0 + H0 * kappa;
255  }
256  }
257  return yieldStress;
258 }
259 
260 double
262 {
263  double plasticModulus = 0.;
264  if ( plasthardtype == 0 ) { // linear hardening
265  plasticModulus = H0;
266  } else if ( plasthardtype == 1 ) { // exponential hardening
267  plasticModulus = H0 * exp(-H0 * kappa / delSigY);
268  } else if ( plasthardtype == 2 ) { // exponential hardening before the peak stress sig0 and linear after the peak stress sig0
269  if ( kappa <= ep ) { //1st branch of yield stress
270  double aux = pow(kappa / ep, md);
271  plasticModulus = 50. *E *exp(-aux / md) * ( 1. - aux );
272  } else { //2nd branch of yield stress
273  plasticModulus = H0;
274  }
275  }
276  return plasticModulus;
277 }
278 
279 
280 // computes the stress according to elastoplasticity
281 // (return of trial stress to the yield surface)
282 void
284 {
285  double kappa, tempKappa, H;
286  RankineMatStatus *status = static_cast< RankineMatStatus * >( this->giveStatus(gp) );
287  MaterialMode mode = gp->giveMaterialMode();
288  // get the initial plastic strain and initial kappa from the status
289  FloatArray tempPlasticStrain = status->givePlasticStrain();
290  kappa = tempKappa = status->giveCumulativePlasticStrain();
291 
292  // elastic predictor
293  StrainVector elStrain(totalStrain, mode);
294  elStrain.subtract(tempPlasticStrain);
295  StressVector finalStress(mode);
296  elStrain.applyElasticStiffness(finalStress, E, nu);
297  FloatArray sigPrinc;
298  FloatMatrix nPrinc;
299  // get principal trial stresses (ordered) and principal stress directions
300  finalStress.computePrincipalValDir(sigPrinc, nPrinc);
301  double ftrial = evalYieldFunction(sigPrinc, tempKappa);
302  if ( mode == _1dMat ) { //1d case
304  if ( ftrial > 0. ) {
305  double f = ftrial;
306  // calculate the increment of cumulative plastic strain
307  int i = 1;
308  do {
309  i = i + 1;
310  if ( i > 1000 ) {
311  printf("kappa, ftrial: %g %g\n", kappa, ftrial);
312  OOFEM_ERROR("no convergence of regular stress return algorithm");
313  }
314 
315  double ddKappa = f / ( E + evalPlasticModulus(tempKappa) );
316  finalStress.at(1) -= E * ddKappa;
317  tempKappa += ddKappa;
318  f = finalStress.at(1) - evalYieldStress(tempKappa);
319  } while ( fabs(f) > yieldtol * sig0 );
320  }
321 
322  tempPlasticStrain.at(1) = tempKappa;
323  //End of Dimitris change
324  } else { //Plane stress case
325  double difPrincTrialStresses = sigPrinc.at(1) - sigPrinc.at(2);
326  double tanG = E / ( 2. * ( 1. + nu ) );
327 
328  // plastic corrector - regular case
329  bool vertex_case = false;
330  if ( ftrial > 0. ) {
331  double f = ftrial;
332  double Enu = E / ( 1. - nu * nu );
333  // calculate the increment of cumulative plastic strain
334  int i = 1;
335  do {
336  if ( i++ > 50 ) {
337  finalStress.computePrincipalValDir(sigPrinc, nPrinc);
338  sigPrinc.pY();
339  printf("kappa, ftrial: %g %g\n", kappa, ftrial);
340  OOFEM_ERROR("no convergence of regular stress return algorithm");
341  }
342 
343  H = evalPlasticModulus(tempKappa);
344  double ddKappa = f / ( Enu + H );
345  sigPrinc.at(1) -= Enu * ddKappa;
346  sigPrinc.at(2) -= nu * Enu * ddKappa;
347  tempKappa += ddKappa;
348  f = evalYieldFunction(sigPrinc, tempKappa);
349  } while ( fabs(f) > yieldtol * sig0 );
350 
351  if ( sigPrinc.at(2) > sigPrinc.at(1) ) {
352  // plastic corrector - vertex case
353  // ---------------------------------
354  vertex_case = true;
355  // recompute trial principal stresses
356  finalStress.computePrincipalValDir(sigPrinc, nPrinc);
357  // calculate the increment of cumulative plastic strain
358  double sigstar = ( sigPrinc.at(1) - nu * sigPrinc.at(2) ) / ( 1. - nu );
359  double alpha = E / ( 1. - nu );
360  double dkap0 = ( sigPrinc.at(1) - sigPrinc.at(2) ) * ( 1. + nu ) / E;
361  tempKappa = kappa;
362  f = sigstar - evalYieldStress(tempKappa);
363  double dkap1 = 0.;
364  H = evalPlasticModulus(tempKappa);
365  double C = alpha + H * ( 1. + sqrt(2.) ) / 2.;
366  i = 1;
367  do {
368  if ( i++ > 20 ) {
369  finalStress.computePrincipalValDir(sigPrinc, nPrinc);
370  sigPrinc.pY();
371  printf("kappa, ftrial: %g %g\n", kappa, ftrial);
372  OOFEM_ERROR("no convergence of vertex stress return algorithm");
373  }
374 
375  dkap1 += f / C;
376  tempKappa = kappa + sqrt( dkap1 * dkap1 + ( dkap1 - dkap0 ) * ( dkap1 - dkap0 ) );
377  f = sigstar - evalYieldStress(tempKappa) - alpha * dkap1;
378  double aux = dkap1 * dkap1 + ( dkap1 - dkap0 ) * ( dkap1 - dkap0 );
379  H = evalPlasticModulus(tempKappa);
380  if ( aux > 0. ) {
381  C = alpha + H * ( 2. * dkap1 - dkap0 ) / sqrt(aux);
382  } else {
383  C = alpha + H *sqrt(2.);
384  }
385  } while ( fabs(f) > yieldtol * sig0 );
386 
387  sigPrinc.at(1) = sigPrinc.at(2) = sigstar - alpha * dkap1;
388  status->setDKappa(dkap1, dkap1 - dkap0);
389  }
390 
391  // principal stresses
392  double sig1 = sigPrinc.at(1);
393  double sig2 = sigPrinc.at(2);
394  // compose the stress in global coordinates
395  // the first subscript refers to coordinate
396  // the second subscript refers to eigenvalue
397  double n11 = nPrinc.at(1, 1);
398  double n12 = nPrinc.at(1, 2);
399  double n21 = nPrinc.at(2, 1);
400  double n22 = nPrinc.at(2, 2);
401  finalStress.at(1) = sig1 * n11 * n11 + sig2 * n12 * n12;
402  finalStress.at(2) = sig1 * n21 * n21 + sig2 * n22 * n22;
403  finalStress.at(3) = sig1 * n11 * n21 + sig2 * n12 * n22;
404  // add the increment of plastic strain
405  if ( !vertex_case ) {
406  tempPlasticStrain.at(1) += ( tempKappa - kappa ) * n11 * n11;
407  tempPlasticStrain.at(2) += ( tempKappa - kappa ) * n21 * n21;
408  tempPlasticStrain.at(3) += 2. * ( tempKappa - kappa ) * n11 * n21;
409  } else {
410  double dkap1 = status->giveDKappa(1);
411  double dkap2 = status->giveDKappa(2);
412  tempPlasticStrain.at(1) += dkap1 * n11 * n11 + dkap2 * n12 * n12;
413  tempPlasticStrain.at(2) += dkap1 * n21 * n21 + dkap2 * n22 * n22;
414  tempPlasticStrain.at(3) += 2. * ( dkap1 * n11 * n21 + dkap2 * n12 * n22 );
415  }
416 
417  // evaluate the tangent shear stiffness
418  if ( difPrincTrialStresses != 0. ) {
419  double factor = ( sig1 - sig2 ) / difPrincTrialStresses;
420  if ( factor > 0. && factor <= 1. ) {
421  tanG *= factor;
422  }
423  }
424  }
425 
426  status->setTangentShearStiffness(tanG); // store shear stiffness.Used in 2d/3d cases
427  }
428 
429  // store the effective stress, plastic strain and cumulative plastic strain
430  status->letTempEffectiveStressBe(finalStress);
431  status->letTempPlasticStrainBe(tempPlasticStrain);
432  status->setTempCumulativePlasticStrain(tempKappa);
433 }
434 
435 double
437 {
438  double tempDam = 0.;
439  if ( tempKappa > 0. ) {
440  if ( damlaw == 0 ) {
441  tempDam = 1.0 - exp(-a * tempKappa);
442  } else if ( damlaw == 1 && tempKappa > ep ) {
443  tempDam = 1.0 - exp( -param1 * pow( ( tempKappa - ep ) / ep, param2 ) );
444  } else if ( damlaw == 2 && tempKappa > ep ) {
445  tempDam = 1.0 - param5 *exp( -param1 *pow ( ( tempKappa - ep ) / ep, param2 ) ) - ( 1. - param5 ) * exp( -param3 * pow( ( tempKappa - ep ) / ep, param4 ) );
446  }
447  }
448 
449  return tempDam;
450 }
451 
452 double
454 {
455  double tempDam = 0.;
456  if ( tempKappa >= 0. ) {
457  if ( damlaw == 0 ) {
458  tempDam = a * exp(-a * tempKappa);
459  } else if ( damlaw == 1 && tempKappa >= ep ) {
460  tempDam = param1 * param2 * pow( ( tempKappa - ep ) / ep, param2 - 1 ) / ep *exp( -param1 *pow ( ( tempKappa - ep ) / ep, param2 ) );
461  } else if ( damlaw == 2 && tempKappa >= ep ) {
462  tempDam = param5 * param1 * param2 * pow( ( tempKappa - ep ) / ep, param2 - 1 ) / ep *exp( -param1 *pow ( ( tempKappa - ep ) / ep, param2 ) ) + ( 1. - param5 ) * param3 * param4 * pow( ( tempKappa - ep ) / ep, param4 - 1 ) / ep *exp( -param3 *pow ( ( tempKappa - ep ) / ep, param4 ) );
463  }
464  }
465 
466  return tempDam;
467 }
468 
469 double
471 {
472  double tempKappa, dam;
473  RankineMatStatus *status = static_cast< RankineMatStatus * >( this->giveStatus(gp) );
474  dam = status->giveDamage();
475  computeCumPlastStrain(tempKappa, gp, tStep);
476  double tempDam = computeDamageParam(tempKappa);
477  if ( dam > tempDam ) {
478  tempDam = dam;
479  }
480 
481  return tempDam;
482 }
483 
484 void RankineMat :: computeCumPlastStrain(double &tempKappa, GaussPoint *gp, TimeStep *tStep)
485 {
486  RankineMatStatus *status = static_cast< RankineMatStatus * >( this->giveStatus(gp) );
487  tempKappa = status->giveTempCumulativePlasticStrain();
488 }
489 
490 // returns the consistent (algorithmic) stiffness matrix
491 void
493  MatResponseMode mode,
494  GaussPoint *gp,
495  TimeStep *tStep)
496 {
497  RankineMatStatus *status = static_cast< RankineMatStatus * >( this->giveStatus(gp) );
498  double tempKappa = status->giveTempCumulativePlasticStrain();
499  double gprime = computeDamageParamPrime(tempKappa);
500  evaluatePlaneStressStiffMtrx(answer, mode, gp, tStep, gprime);
501 }
502 
503 void
505 {
506  RankineMatStatus *status = static_cast< RankineMatStatus * >( this->giveStatus(gp) );
507  answer.resize(1, 1);
508  answer.at(1, 1) = this->E;
509  if ( mode == ElasticStiffness ) {
510  } else if ( mode == SecantStiffness ) {
511  double om = status->giveTempDamage();
512  answer.times(1.0 - om);
513  } else {
514  OOFEM_ERROR("unknown type of stiffness (secant stiffness not implemented for 1d)");
515  }
516 }
517 
518 // this method is also used by the gradient version,
519 // with gprime replaced by gprime*m and evaluated for kappa hat
520 void
522  MatResponseMode mode,
523  GaussPoint *gp,
524  TimeStep *tStep, double gprime)
525 {
526  RankineMatStatus *status = static_cast< RankineMatStatus * >( this->giveStatus(gp) );
527  if ( mode == ElasticStiffness || mode == SecantStiffness ) {
528  // start from the elastic stiffness
529  this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, mode, gp, tStep);
530  if ( mode == SecantStiffness ) {
531  // transform to secant stiffness
532  double damage = status->giveTempDamage();
533  answer.times(1. - damage);
534  }
535 
536  return;
537  }
538 
539  // check the unloading condition
540  double kappa = status->giveCumulativePlasticStrain();
541  double tempKappa = status->giveTempCumulativePlasticStrain();
542  if ( tempKappa <= kappa ) { // tangent matrix requested, but unloading takes place - use secant
543  this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, mode, gp, tStep);
544  double damage = status->giveTempDamage();
545  answer.times(1. - damage);
546  return;
547  }
548 
549  // tangent stiffness requested, loading
550 
551  // elastoplastic tangent matrix in principal stress coordinates
552  answer.resize(3, 3);
553  answer.zero(); // just to be sure (components 13,23,31,32 must be zero)
554 
555  double eta1, eta2, dkap2;
556  double dkap1 = status->giveDKappa(1);
557  double H = evalPlasticModulus(tempKappa);
558 
559  if ( dkap1 == 0. ) {
560  // regular case
561 
562  double Enu = E / ( 1. - nu * nu );
563  double aux = Enu / ( Enu + H );
564  answer.at(1, 1) = aux * H;
565  answer.at(1, 2) = answer.at(2, 1) = nu * aux * H;
566  answer.at(2, 2) = aux * ( E + H );
567  answer.at(3, 3) = status->giveTangentShearStiffness();
568  eta1 = aux;
569  eta2 = nu * aux;
570  } else {
571  // vertex case
572 
573  dkap2 = status->giveDKappa(2);
574  double denom = E * dkap1 + H * ( 1. - nu ) * ( dkap1 + dkap2 );
575  eta1 = E * dkap1 / denom;
576  eta2 = E * dkap2 / denom;
577  answer.at(1, 1) = answer.at(2, 1) = H * eta1;
578  answer.at(1, 2) = answer.at(2, 2) = H * eta2;
579  answer.at(3, 3) = 0.; //E*1.e-6; // small stiffness, to suppress singularity
580  }
581 
582  // the elastoplastic stiffness has been computed
583  // now add the effect of damage
584 
585  double damage = status->giveTempDamage();
586  answer.times(1. - damage);
587 
588  FloatArray sigPrinc(2);
589  FloatMatrix nPrinc(2, 2);
590  StressVector effStress(status->giveTempEffectiveStress(), _PlaneStress);
591  effStress.computePrincipalValDir(sigPrinc, nPrinc);
592  // sometimes the method is called with gprime=0., then we can save some work
593  if ( gprime != 0. ) {
594  FloatMatrix correction(3, 3);
595  correction.zero();
596  correction.at(1, 1) = sigPrinc.at(1) * eta1;
597  correction.at(1, 2) = sigPrinc.at(1) * eta2;
598  correction.at(2, 1) = sigPrinc.at(2) * eta1;
599  correction.at(2, 2) = sigPrinc.at(2) * eta2;
600  correction.times(gprime); // input parameter gprime used here
601  answer.subtract(correction);
602  }
603 
604  // transform to global coordinates
605  FloatMatrix T;
607  answer.rotatedWith(T, 't');
608 }
609 
610 // derivatives of final kappa with respect to final strain
611 void
613 {
614  FloatArray eta(3);
615  double dkap1 = status->giveDKappa(1);
616  double kap = status->giveTempCumulativePlasticStrain();
617  double H = evalPlasticModulus(kap);
618 
619  // evaluate in principal coordinates
620 
621  if ( dkap1 == 0. ) {
622  // regular case
623  double Estar = E / ( 1. - nu * nu );
624  double aux = Estar / ( H + Estar );
625  eta.at(1) = aux;
626  eta.at(2) = nu * aux;
627  eta.at(3) = 0.;
628  } else {
629  // vertex case
630  double dkap2 = status->giveDKappa(2);
631  double denom = E * dkap1 + H * ( 1. - nu ) * ( dkap1 + dkap2 );
632  eta.at(1) = E * dkap1 / denom;
633  eta.at(2) = E * dkap2 / denom;
634  eta.at(3) = 0.;
635  }
636 
637  // transform to global coordinates
638 
639  FloatArray sigPrinc(2);
640  FloatMatrix nPrinc(2, 2);
641  StressVector effStress(status->giveTempEffectiveStress(), _PlaneStress);
642  effStress.computePrincipalValDir(sigPrinc, nPrinc);
643 
644  FloatMatrix T(3, 3);
646  answer.beProductOf(T, eta);
647 }
648 
649 int
651 {
652  RankineMatStatus *status = static_cast< RankineMatStatus * >( this->giveStatus(gp) );
653  if ( type == IST_PlasticStrainTensor ) {
654  const FloatArray ep = status->givePlasDef();
655  answer.resize(6);
656  answer.at(1) = ep.at(1);
657  answer.at(2) = ep.at(2);
658  answer.at(3) = 0.;
659  answer.at(4) = 0.;
660  answer.at(5) = 0.;
661  answer.at(6) = ep.at(3);
662  return 1;
663  } else if ( type == IST_CumPlasticStrain ) {
664  answer.resize(1);
665  answer.at(1) = status->giveCumulativePlasticStrain();
666  return 1;
667  } else if ( type == IST_DamageScalar ) {
668  answer.resize(1);
669  answer.at(1) = status->giveDamage();
670  return 1;
671  } else if ( type == IST_DamageTensor ) {
672  answer.resize(6);
673  answer.zero();
674  answer.at(1) = answer.at(2) = answer.at(3) = status->giveDamage();
675  return 1;
676 
677 #ifdef keep_track_of_dissipated_energy
678  } else if ( type == IST_StressWorkDensity ) {
679  answer.resize(1);
680  answer.at(1) = status->giveStressWork();
681  return 1;
682  } else if ( type == IST_DissWorkDensity ) {
683  answer.resize(1);
684  answer.at(1) = status->giveDissWork();
685  return 1;
686  } else if ( type == IST_FreeEnergyDensity ) {
687  answer.resize(1);
688  answer.at(1) = status->giveStressWork() - status->giveDissWork();
689  return 1;
690 
691 #endif
692  } else {
693  return StructuralMaterial :: giveIPValue(answer, gp, type, tStep);
694  }
695 }
696 
697 
698 //=============================================================================
699 // RANKINE MATERIAL STATUS
700 //=============================================================================
701 
703  StructuralMaterialStatus(n, d, g), plasticStrain(), tempPlasticStrain()
704 {
705  damage = tempDamage = 0.;
706  kappa = tempKappa = 0.;
707  dKappa1 = dKappa2 = 0.;
708  tanG = 0.;
709  //effStress.resize(3);
710  //tempEffStress.resize(3);
711 #ifdef keep_track_of_dissipated_energy
712  stressWork = tempStressWork = 0.0;
713  dissWork = tempDissWork = 0.0;
714 #endif
715 }
716 
717 
719 { }
720 
721 
722 void
724 {
725  //int i, n;
726 
728 
729  fprintf(file, "status { ");
730  /*
731  * // this would not be correct, since the status is already updated and kappa is now the "final" value
732  * if ( tempKappa > kappa ) {
733  * fprintf(file, " Yielding, ");
734  * } else {
735  * fprintf(file, " Unloading, ");
736  * }
737  */
738  fprintf(file, "kappa %g, damage %g ", kappa, tempDamage);
739 #ifdef keep_track_of_dissipated_energy
740  fprintf(file, ", dissW %g, freeE %g, stressW %g ", this->dissWork, ( this->stressWork ) - ( this->dissWork ), this->stressWork);
741 #endif
742 #if 0
743  // print the plastic strain
744  fprintf(file, " plastic_strains ");
745  for ( auto &val : plasticStrain ) {
746  fprintf( file, " %.4e", val );
747  }
748 #endif
749  // print the cumulative plastic strain
750  fprintf(file, "}\n");
751 }
752 
753 
754 // initializes temporary variables based on their values at the previous equlibrium state
756 {
758 
759  if ( plasticStrain.giveSize() == 0 ) {
762  }
763 
764  tempDamage = damage;
766  tempKappa = kappa;
767 #ifdef keep_track_of_dissipated_energy
770 #endif
771 }
772 
773 
774 // updates internal variables when equilibrium is reached
775 void
777 {
779 
781  kappa = tempKappa;
782  damage = tempDamage;
783 #ifdef keep_track_of_dissipated_energy
786 #endif
787 }
788 
789 
790 // saves full information stored in this status
791 // temporary variables are NOT stored
794 {
795  contextIOResultType iores;
796 
797  // save parent class status
798  if ( ( iores = StructuralMaterialStatus :: saveContext(stream, mode, obj) ) != CIO_OK ) {
799  THROW_CIOERR(iores);
800  }
801 
802  // write raw data
803 
804  // write plastic strain (vector)
805  if ( ( iores = plasticStrain.storeYourself(stream) ) != CIO_OK ) {
806  THROW_CIOERR(iores);
807  }
808 
809  // write cumulative plastic strain (scalar)
810  if ( !stream.write(kappa) ) {
812  }
813 
814  // write damage (scalar)
815  if ( !stream.write(damage) ) {
817  }
818 
819 #ifdef keep_track_of_dissipated_energy
820  if ( !stream.write(stressWork) ) {
822  }
823 
824  if ( !stream.write(dissWork) ) {
826  }
827 
828 #endif
829 
830  return CIO_OK;
831 }
832 
833 
836 //
837 // restores full information stored in stream to this Status
838 //
839 {
840  contextIOResultType iores;
841 
842  // read parent class status
843  if ( ( iores = StructuralMaterialStatus :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
844  THROW_CIOERR(iores);
845  }
846 
847  // read plastic strain (vector)
848  if ( ( iores = plasticStrain.restoreYourself(stream) ) != CIO_OK ) {
849  THROW_CIOERR(iores);
850  }
851 
852  // read cumulative plastic strain (scalar)
853  if ( !stream.read(kappa) ) {
855  }
856 
857  // read damage (scalar)
858  if ( !stream.read(damage) ) {
860  }
861 
862 #ifdef keep_track_of_dissipated_energy
863  if ( !stream.read(stressWork) ) {
865  }
866 
867  if ( !stream.read(dissWork) ) {
869  }
870 
871 #endif
872 
873  return CIO_OK; // return success
874 }
875 
876 
877 #ifdef keep_track_of_dissipated_energy
878 void
880 {
881  // int n = deps.giveSize(); // would not work for gradient version
882  int n = 3;
883 
884  // strain increment
885  FloatArray deps;
887 
888  // increment of stress work density
889  double dSW = ( tempStressVector.dotProduct(deps, n) + stressVector.dotProduct(deps, n) ) / 2.;
890  tempStressWork = stressWork + dSW;
891 
892  // elastically stored energy density
893  FloatArray tempElasticStrainVector;
894  tempElasticStrainVector.beDifferenceOf(tempStrainVector, tempPlasticStrain, n);
895  double We = tempStressVector.dotProduct(tempElasticStrainVector, n) / 2.;
896 
897  // dissipative work density
899  // to avoid extremely small negative dissipation due to round-off error
900  // (note: gf is the dissipation density at complete failure, per unit volume)
901  if ( fabs(tempDissWork) < 1.e-12 * gf ) {
902  tempDissWork = 0.;
903  }
904 }
905 
906 void
908 {
909  // int n = deps.giveSize(); // would not work for gradient version
910  int n = 1;
911 
912 
913  // strain increment
914  FloatArray deps;
916 
917  // increment of stress work density
918  double dSW = ( tempStressVector.dotProduct(deps, n) + stressVector.dotProduct(deps, n) ) / 2.;
919  tempStressWork = stressWork + dSW;
920 
921  // elastically stored energy density
922  FloatArray tempElasticStrainVector;
923  tempElasticStrainVector.beDifferenceOf(tempStrainVector, tempPlasticStrain, n);
924  double We = tempStressVector.dotProduct(tempElasticStrainVector, n) / 2.;
925 
926  // dissipative work density
928  // to avoid extremely small negative dissipation due to round-off error
929  // (note: gf is the dissipation density at complete failure, per unit volume)
930  if ( fabs(tempDissWork) < 1.e-12 * gf ) {
931  tempDissWork = 0.;
932  }
933 }
934 
935 #endif
936 } // end namespace oofem
static int giveSizeOfVoigtSymVector(MaterialMode mmode)
Returns the size of symmetric part of a reduced stress/strain vector according to given mode...
void computeWork_1d(GaussPoint *gp, double gf)
Definition: rankinemat.C:907
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Definition: structuralms.C:96
FloatArray plasticStrain
Plastic strain (initial).
Definition: rankinemat.h:200
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition: gausspoint.h:191
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
void letTempStrainVectorBe(const FloatArray &v)
Assigns tempStrainVector to given vector v.
Definition: structuralms.h:137
#define _IFT_RankineMat_param4
Definition: rankinemat.h:66
static void givePlaneStressVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
Computes 2d stress vector transformation matrix from standard vector transformation matrix...
virtual void giveRealStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducesStrain, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Definition: rankinemat.C:202
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Definition: rankinemat.C:776
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Returns material status of receiver in given integration point.
Definition: material.C:244
#define _IFT_RankineMat_a
Definition: rankinemat.h:56
double tempKappa
Cumulative plastic strain (final).
Definition: rankinemat.h:215
GaussPoint * gp
Associated integration point.
Class and object Domain.
Definition: domain.h:115
double tanG
Tangent shear stiffness (needed for tangent matrix).
Definition: rankinemat.h:231
void applyElasticStiffness(StressVector &stress, const double EModulus, const double nu) const
Applies the elastic stiffness to the strain.
Definition: strainvector.C:350
FloatArray tempStrainVector
Temporary strain vector in reduced form (to find balanced state)
Definition: structuralms.h:75
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
contextIOResultType storeYourself(DataStream &stream) const
Definition: floatarray.C:872
double damage
Damage (initial).
Definition: rankinemat.h:225
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
void evaluatePlaneStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep, double gprime)
Executive method used by local and gradient version.
Definition: rankinemat.C:521
Specialization of a floating point array for representing a stress state.
Definition: stressvector.h:48
double sig0
Initial (uniaxial) yield stress.
Definition: rankinemat.h:107
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
const FloatArray & giveTempEffectiveStress() const
Definition: rankinemat.h:269
This class implements a structural material status information.
Definition: structuralms.h:65
void computePrincipalValDir(FloatArray &answer, FloatMatrix &dir) const
Computes principal values and directions of receiver vector.
Definition: stressvector.C:182
#define _IFT_RankineMat_plasthardtype
Definition: rankinemat.h:57
FloatArray stressVector
Equilibrated stress vector in reduced form.
Definition: structuralms.h:71
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
Definition: rankinemat.C:650
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: rankinemat.C:835
void setDKappa(double val1, double val2)
Definition: rankinemat.h:279
double yieldtol
Relative tolerance in yield condition.
Definition: rankinemat.h:110
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
General IO error.
LinearElasticMaterial * giveLinearElasticMaterial()
Returns a reference to the basic elastic material.
Definition: rankinemat.h:163
double evalYieldStress(const double kappa)
Definition: rankinemat.C:238
Specialization of a floating point array for representing a strain state.
Definition: strainvector.h:45
double stressWork
Density of total work done by stresses on strain increments.
Definition: rankinemat.h:235
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Print receiver&#39;s output to given stream.
Definition: structuralms.C:73
virtual int read(int *data, int count)=0
Reads count integer values into array pointed by data.
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
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
double ep
Total strain at peak stress sig0–Used only if plasthardtype=2.
Definition: rankinemat.h:116
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
Definition: floatarray.C:341
double delSigY
Final increment of yield stress (at infinite cumulative plastic strain)
Definition: rankinemat.h:104
void computeWork_PlaneStress(GaussPoint *gp, double gf)
Computes the increment of total stress work and of dissipated work (gf is the dissipation density per...
Definition: rankinemat.C:879
virtual int write(const int *data, int count)=0
Writes count integer values from array pointed by data.
void beScaled(double s, const FloatArray &b)
Sets receiver to be .
Definition: floatarray.C:146
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Print receiver&#39;s output to given stream.
Definition: rankinemat.C:723
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...
RankineMatStatus(int n, Domain *d, GaussPoint *g)
Definition: rankinemat.C:702
void setTempCumulativePlasticStrain(double value)
Definition: rankinemat.h:277
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
Definition: rankinemat.C:755
#define _IFT_RankineMat_gf
Definition: rankinemat.h:60
double dotProduct(const FloatArray &x) const
Computes the dot product (or inner product) of receiver and argument.
Definition: floatarray.C:463
const FloatArray & givePlasticStrain() const
Definition: rankinemat.h:248
This class implements an isotropic linear elastic material in a finite element problem.
#define OOFEM_ERROR(...)
Definition: error.h:61
void times(double f)
Multiplies receiver by factor f.
Definition: floatmatrix.C:1594
void rotatedWith(const FloatMatrix &r, char mode= 'n')
Returns the receiver &#39;a&#39; transformed using give transformation matrix r.
Definition: floatmatrix.C:1557
contextIOResultType restoreYourself(DataStream &stream)
Definition: floatarray.C:895
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
const FloatArray & givePlasDef()
Definition: rankinemat.h:288
#define _IFT_RankineMat_damlaw
Definition: rankinemat.h:62
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
double giveStressWork()
Returns the density of total work of stress on strain increments.
Definition: rankinemat.h:300
#define _IFT_RankineMat_yieldtol
Definition: rankinemat.h:59
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
virtual int hasMaterialModeCapability(MaterialMode mode)
Tests if material supports material mode.
Definition: rankinemat.C:73
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: structuralms.C:157
double param3
coefficient required when damlaw=2
Definition: rankinemat.h:131
double giveCumulativePlasticStrain()
Definition: rankinemat.h:253
void subtract(const FloatMatrix &a)
Subtracts matrix from the receiver.
Definition: floatmatrix.C:1084
Abstract base class representing a material status information.
Definition: matstatus.h:84
Class representing vector of real numbers.
Definition: floatarray.h:82
#define _IFT_RankineMat_param2
Definition: rankinemat.h:64
void setTempDamage(double value)
Definition: rankinemat.h:284
virtual ~RankineMatStatus()
Definition: rankinemat.C:718
FloatArray strainVector
Equilibrated strain vector in reduced form.
Definition: structuralms.h:69
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 dKappa1
Increments of cumulative plastic strain associated with the first and secomnd principal stress (used ...
Definition: rankinemat.h:222
void letTempStressVectorBe(const FloatArray &v)
Assigns tempStressVector to given vector v.
Definition: structuralms.h:135
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: rankinemat.C:793
double tempStressWork
Non-equilibrated density of total work done by stresses on strain increments.
Definition: rankinemat.h:237
#define _IFT_RankineMat_sig0
Definition: rankinemat.h:54
virtual void pY() const
Print receiver on stdout with high accuracy.
Definition: floatarray.C:787
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
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
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
double param4
coefficient required when damlaw=2
Definition: rankinemat.h:134
#define _IFT_RankineMat_param1
Definition: rankinemat.h:63
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&#39;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
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
Definition: structuralms.C:108
FloatArray tempPlasticStrain
Plastic strain (final).
Definition: rankinemat.h:203
Domain * giveDomain() const
Definition: femcmpnn.h:100
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
#define _IFT_RankineMat_ep
Definition: rankinemat.h:61
REGISTER_Material(DummyMaterial)
double param2
coefficient required when damlaw=1 or 2
Definition: rankinemat.h:128
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:78
void letTempPlasticStrainBe(FloatArray values)
Definition: rankinemat.h:271
#define _IFT_RankineMat_param5
Definition: rankinemat.h:67
virtual void computeCumPlastStrain(double &kappa, GaussPoint *gp, TimeStep *tStep)
Definition: rankinemat.C:484
RankineMat(int n, Domain *d)
Definition: rankinemat.C:52
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: structuralms.C:133
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
#define IR_GIVE_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:69
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&#39;s ratio.
Definition: rankinemat.h:95
FloatArray tempStressVector
Temporary stress vector in reduced form (increments are used mainly in nonlinear analysis) ...
Definition: structuralms.h:73
void letTempEffectiveStressBe(FloatArray values)
Definition: rankinemat.h:275
double H0
Initial hardening modulus.
Definition: rankinemat.h:98
#define _IFT_RankineMat_delsigy
Definition: rankinemat.h:58
Class representing integration point in finite element program.
Definition: gausspoint.h:93
#define OOFEM_WARNING(...)
Definition: error.h:62
Class representing solution step.
Definition: timestep.h:80
#define _IFT_RankineMat_param3
Definition: rankinemat.h:65
double tempDamage
Damage (final).
Definition: rankinemat.h:228
double tempDissWork
Non-equilibrated density of dissipated work.
Definition: rankinemat.h:241
#define _IFT_RankineMat_h
Definition: rankinemat.h:55
double giveDissWork()
Returns the density of dissipated work.
Definition: rankinemat.h:306
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

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