OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
perfectlyplasticmaterial.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 
36 #include "../sm/CrossSections/structuralcrosssection.h"
37 #include "gausspoint.h"
38 #include "floatmatrix.h"
39 #include "floatarray.h"
40 #include "mathfem.h"
41 #include "datastream.h"
42 #include "contextioerr.h"
43 
44 namespace oofem {
45 int
47 //
48 // returns whether receiver supports given mode
49 //
50 {
51  return mode == _3dMat ||
52  //<RESTRICTED_SECTION>
53  mode == _PlaneStress ||
54  //</RESTRICTED_SECTION>
55  mode == _PlaneStrain ||
56  mode == _1dMat ||
57  mode == _PlateLayer ||
58  mode == _2dBeamLayer;
59 }
60 
61 
62 void
64  GaussPoint *gp,
65  const FloatArray &totalStrain,
66  TimeStep *tStep)
67 //
68 // returns stress vector (in full or reduced form ) of receiver according to
69 // previous level of stress and current
70 // strain increment, the only way, how to correctly update gp records
71 //
72 //
73 // may be good idea for nonlinear materials overload this function in order to
74 // capture strain - stress history for proper modelling receiver behaviour
75 // and in order to capture possible failure of material
76 //
77 // Note: formulated in full stress strain space
78 {
79  FloatArray elasticStressIncrement;
80  FloatArray workStress, *yieldStressGrad, workStress2, stressVector3d;
81  FloatArray mStrainIncrement3d, mStressElasticIncrement3d, PlasticStrainVector3d;
82  FloatArray dSigmaIncrement3d, *stressCorrection, *loadingStressGrad;
83  FloatArray helpArray;
84  FloatArray plasticStrainIncrement3d, strainIncrement, reducedStrain, reducedStrainIncrement;
85  FloatArray statusFullStressVector, statusFullPlasticVector;
86  FloatArray plasticStrainVector;
87  PerfectlyPlasticMaterialStatus *status = static_cast< PerfectlyPlasticMaterialStatus * >( this->giveStatus(gp) );
88 
89  FloatMatrix dp;
90  StructuralCrossSection *crossSection = static_cast< StructuralCrossSection * >
91  ( gp->giveElement()->giveCrossSection() );
92  double f0, f1, f2, help, dLambda, r, r1, m;
93 
94  // init temp variables (of YC,LC,Material) at the beginning of step
95  this->initTempStatus(gp);
96  // subtract stress independent part
97  // note: eigenStrains (temperature) is not contained in mechanical strain stored in gp
98  // therefore it is necessary to subtract always the total eigen strain value
99  this->giveStressDependentPartOfStrainVector(reducedStrain, gp, totalStrain,
100  tStep, VM_Total);
101 
102  reducedStrainIncrement.beDifferenceOf( reducedStrain, status->giveStrainVector() );
103 
104  StructuralMaterial :: giveFullSymVectorForm( strainIncrement, reducedStrainIncrement, gp->giveMaterialMode() );
105  plasticStrainVector = status->givePlasticStrainVector();
106  StructuralMaterial :: giveFullSymVectorForm( statusFullPlasticVector, plasticStrainVector, gp->giveMaterialMode() );
107  StructuralMaterial :: giveFullSymVectorForm( statusFullStressVector, status->giveStressVector(), gp->giveMaterialMode() );
108 
109  //
110  // Note : formulated in full stress strain space
111  //
112  this->computeTrialStressIncrement(elasticStressIncrement, gp, strainIncrement, tStep);
113  //
114  // calculate deltaSigmaPlastic
115  //
116  workStress = statusFullStressVector;
117  workStress.add(elasticStressIncrement);
118 
119  f0 = this->computeYCValueAt(gp, & statusFullStressVector,
120  & statusFullPlasticVector);
121 
122  f1 = this->computeYCValueAt(gp, & workStress,
123  & statusFullPlasticVector);
124 
125  PlasticStrainVector3d = statusFullPlasticVector;
126 
127 
128  if ( f1 >= 0. ) { // loading surface violated by the elastic trial set
129  if ( f0 < 0. ) { // previously in elastic area
130  // element not yielding - set the print status
131  status->setTempYieldFlag(0);
132  // compute scaling factor
133  r1 = -f0 / ( f1 - f0 ); // linear interpolation in f (Zienkiewicz, 1969)
134  yieldStressGrad = this->GiveYCStressGradient(gp, & statusFullStressVector,
135  & statusFullPlasticVector);
136  crossSection->imposeStressConstrainsOnGradient(gp, yieldStressGrad);
137 
138  help = 0.;
139  for ( int i = 1; i <= 6; i++ ) {
140  help += yieldStressGrad->at(i) * elasticStressIncrement.at(i);
141  }
142 
143  delete yieldStressGrad;
144 
145  workStress2 = elasticStressIncrement;
146  workStress2.times(r1);
147  workStress2.add(statusFullStressVector);
148 
149  f2 = this->computeYCValueAt(gp, & workStress2,
150  & statusFullPlasticVector);
151 
152  if ( fabs(help) > 1.e-6 ) {
153  r = r1 - f2 / help; // improved value of r (Nayak & Zienkiewicz, 1972)
154  } else {
155  r = r1;
156  }
157  } else { // f0 should be zero
158  r = 0.;
159  }
160 
161  stressVector3d = elasticStressIncrement;
162  stressVector3d.times(r);
163  stressVector3d.add(statusFullStressVector);
164 
165  m = STRAIN_STEPS;
166  mStrainIncrement3d = strainIncrement;
167  mStrainIncrement3d.times( ( 1.0 - r ) / m );
168  mStressElasticIncrement3d = elasticStressIncrement;
169  mStressElasticIncrement3d.times( ( 1.0 - r ) / m );
170 
171  // test if fracture or failure occurs
172  this->updateIfFailure(gp, & stressVector3d, & PlasticStrainVector3d);
173  // element yielding - set the print status
174  status->setTempYieldFlag(1);
175  // loop over m-steps
176  for ( int i = 1; i <= m; i++ ) {
177  // yieldStressGrad = yieldCriteria->
178  // GiveStressGradient (gp, stressVector3d,
179  // PlasticStrainVector3d,
180  // gp-> giveHardeningParam());
181 
182  // compute dLambda and dp
183  this->computePlasticStiffnessAt(dp, gp,
184  & stressVector3d,
185  & PlasticStrainVector3d,
186  & mStrainIncrement3d,
187  tStep,
188  dLambda);
189  if ( dLambda < 0. ) {
190  dLambda = 0.;
191  }
192 
193  dSigmaIncrement3d.beProductOf(dp, mStrainIncrement3d);
194  dSigmaIncrement3d.negated();
195  dSigmaIncrement3d.add(mStressElasticIncrement3d);
196  // compute normal to loading surface
197  loadingStressGrad = this->GiveLCStressGradient(gp, & stressVector3d,
198  & PlasticStrainVector3d);
199 
200  // update the stress state
201  stressVector3d.add(dSigmaIncrement3d);
202  // compute stress corrections
203 
204  stressCorrection = this->
205  GiveStressCorrectionBackToYieldSurface(gp, & stressVector3d,
206  & PlasticStrainVector3d);
207  // apply the stress correction -> correction is in the direction of
208  // the normal to the yield surface
209  stressVector3d.add(* stressCorrection);
210  // test = yieldCriteria-> computeValueAt(stressVector3d,
211  // PlasticStrainVector3d,
212  // gp-> giveHardeningParam());
213 
214  // calculate plastic strain increment
215 
216  crossSection->imposeStrainConstrainsOnGradient(gp, loadingStressGrad);
217  loadingStressGrad->times(dLambda);
218  PlasticStrainVector3d.add(* loadingStressGrad);
219  // strainVector3d -> add (mStrainIncrement3d);
220 
221  // update yieldCriteria and loading criteria records to newly reached state
222  // in loadingStressGrad is stored current plastic vector increment
223  this->updateTempYC(gp, & stressVector3d, & PlasticStrainVector3d);
224  this->updateTempLC(gp, & stressVector3d, & PlasticStrainVector3d);
225 
226  // test if fracture or failure occurs
227  this->updateIfFailure(gp, & stressVector3d, & PlasticStrainVector3d);
228 
229  delete loadingStressGrad;
230  delete stressCorrection;
231  }
232 
233  // compute total stress increment during strainIncrement
234  // totalStressIncrement = statusFullStressVector;
235  // totalStressIncrement.times (-1.0);
236  // totalStressIncrement.add (stressVector3d);
237  // update gp
238  FloatArray helpArry;
239  StructuralMaterial :: giveReducedSymVectorForm( helpArry, stressVector3d, gp->giveMaterialMode() );
240 
241  status->letTempStressVectorBe(helpArry);
242  status->letTempStrainVectorBe(totalStrain);
243  } else {
244  // element not yielding - set the print status
245  status->setTempYieldFlag(0);
246  stressVector3d = statusFullStressVector;
247  stressVector3d.add(elasticStressIncrement);
248  // test if fracture or failure occurs
249  this->updateIfFailure(gp, & stressVector3d, & PlasticStrainVector3d);
250  // update gp
251  status->letTempStrainVectorBe(totalStrain);
252  StructuralMaterial :: giveReducedSymVectorForm( helpArray, stressVector3d, gp->giveMaterialMode() );
253  status->letTempStressVectorBe(helpArray);
254  }
255 
256  // update gp plastic strain
257  plasticStrainIncrement3d.beDifferenceOf(PlasticStrainVector3d, statusFullPlasticVector);
258  StructuralMaterial :: giveReducedSymVectorForm( helpArray, plasticStrainIncrement3d, gp->giveMaterialMode() );
259  status->letPlasticStrainIncrementVectorBe(helpArray);
260 
261  StructuralMaterial :: giveReducedSymVectorForm( answer, stressVector3d, gp->giveMaterialMode() );
262 }
263 
264 
265 void
267  MatResponseMode mode,
268  GaussPoint *gp,
269  TimeStep *tStep)
270 //
271 //
272 // for case of perfectly plastic material
273 // computes full elastic constitutive matrix for case of gp stress-strain state.
274 // if strainIncrement == NULL a loading is assumed
275 //
276 // we follow terminology based on paper from R. de Borst:
277 // "Smeared Cracking, plasticity, creep - Unified Aproach"
278 //
279 // if derived material would like to implement failure behaviour
280 // it must redefine basic Give3dMaterialStiffnessMatrix function
281 // in order to take possible failure (tension cracking) into account
282 //
283 {
284  // FloatMatrix *de; // elastic matrix respecting fracture or failure
285  StructuralMaterial *lMat = static_cast< StructuralMaterial * >( this->giveLinearElasticMaterial() );
286 
287  if ( lMat->hasMaterialModeCapability( gp->giveMaterialMode() ) ) {
288  FloatMatrix stiff;
289  lMat->giveStiffnessMatrix(stiff, mode, gp, tStep);
290  this->giveFullSymMatrixForm( answer, stiff, gp->giveMaterialMode() );
291  } else {
292  OOFEM_ERROR("giveEffectiveMaterialStiffnessMatrix - unsupported material mode");
293  }
294 }
295 
296 void
298  GaussPoint *gp,
299  TimeStep *tStep)
300 //
301 //
302 //
303 // computes full constitutive matrix for case of gp stress-strain state.
304 // it returns elasto-plastic stiffness material matrix.
305 // if strainIncrement == NULL a loading is assumed
306 // for detailed description see (W.F.Chen: Plasticity in Reinforced Concrete, McGraw-Hill, 1982,
307 // chapter 6.)
308 //
309 // if derived material would like to implement failure behaviour
310 // it must redefine basic Give3dMaterialStiffnessMatrix function
311 // in order to take possible failure (tension cracking) into account
312 //
313 // in this function answer is Material stiffness matrix respecting
314 // current stress-strain mode in gp. This is reached by using
315 // impose constraints functions
316 {
317  FloatArray statusFullStressVector, statusFullPlasticVector, plasticStrainVector;
318  double lambda;
319  PerfectlyPlasticMaterialStatus *status = static_cast< PerfectlyPlasticMaterialStatus * >( this->giveStatus(gp) );
320 
321  // double f = domain->giveYieldCriteria(yieldCriteria)->
322  // computeValueAt(gp->giveStressVector(), gp->givePlasticStrainVector(),
323  // gp-> giveHardeningParam());
324 
325  this->giveEffectiveMaterialStiffnessMatrix(answer, mode, gp, tStep);
326  // if (f < YIELD_BOUNDARY) // linear elastic behaviour
327  // return de;
328 
329  if ( status->giveYieldFlag() == 0 ) {
330  return;
331  }
332 
333  // if secant stiffness requested assume initial elastic matrix
334  if ( mode == SecantStiffness ) {
335  return;
336  }
337 
338  plasticStrainVector = status->givePlasticStrainVector();
339  StructuralMaterial :: giveFullSymVectorForm( statusFullPlasticVector, plasticStrainVector, gp->giveMaterialMode() );
340  StructuralMaterial :: giveFullSymVectorForm( statusFullStressVector, status->giveStressVector(), gp->giveMaterialMode() );
341 
342  // yield condition satisfied
343  FloatMatrix dp;
344  this->computePlasticStiffnessAt(dp, gp, & statusFullStressVector,
345  & statusFullPlasticVector,
346  NULL,
347  tStep,
348  lambda);
349  answer.add(dp);
350 }
351 
352 
353 void
355  MatResponseMode mode,
356  GaussPoint *gp,
357  TimeStep *tStep)
358 //
359 //
360 //
361 // computes full 3d constitutive matrix for case of 3d stress-strain state.
362 // it returns elasto-plastic stiffness material matrix.
363 // if strainIncrement == NULL a loading is assumed
364 // for detailed description see (W.F.Chen: Plasticity in Reinforced Concrete, McGraw-Hill, 1982,
365 // chapter 6.)
366 //
367 // if derived material would like to implement failure behaviour
368 // it must redefine basic Give3dMaterialStiffnessMatrix function
369 // in order to take possible failure (tension cracking) into account
370 //
371 //
372 {
373  MaterialMode originalMode = gp->giveMaterialMode();
374  if ( originalMode != _3dMat ) {
375  OOFEM_ERROR("Different stressStrain mode encountered");
376  }
377 
378  // we can force 3d response, and we obtain correct 3d tangent matrix,
379  // but in fact, stress integration algorithm will not work
380  // because in stress integration algorithm we are unable to recognize
381  // which reduction from 3d case should be performed to obtain correct result.
382  // so for new stressStrain state, instead of programming 3d reduction,
383  // you should enhance imposeConstraints functions for ne state, and
384  // then programming simple inteface function for you stressstrain state
385  // calling GiveMaterailStiffenssMatrix, which imposes constrains correctly.
386 
387  if ( mode == ElasticStiffness ) {
388  this->giveLinearElasticMaterial()->give3dMaterialStiffnessMatrix(answer, mode, gp, tStep);
389  } else {
390  this->giveMaterialStiffnessMatrix(answer, mode, gp, tStep);
391  }
392 }
393 
394 
395 void
397  MatResponseMode mode,
398  GaussPoint *gp,
399  TimeStep *tStep)
400 
401 //
402 // returns receiver's 2dPlaneStressMtrx
403 // (2dPlaneStres ==> sigma_z = tau_xz = tau_yz = 0.)
404 //
405 // standard method from Material Class overloaded, because no inversion is needed.
406 // the reduction from 3d case will not work
407 // this implementation should be faster.
408 {
409  FloatMatrix fullAnswer;
410  if ( mode == ElasticStiffness ) {
411  this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, mode, gp, tStep);
412  } else {
413  this->giveMaterialStiffnessMatrix(fullAnswer, mode, gp, tStep);
415  }
416 }
417 
418 
419 void
421  MatResponseMode mode,
422  GaussPoint *gp,
423  TimeStep *tStep)
424 
425 //
426 // return receiver's 2dPlaneStrainMtrx constructed from
427 // general 3dMatrialStiffnessMatrix
428 // (2dPlaneStrain ==> eps_z = gamma_xz = gamma_yz = 0.)
429 //
430 {
431  FloatMatrix fullAnswer;
432  if ( mode == ElasticStiffness ) {
433  this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, mode, gp, tStep);
434  } else {
435  this->giveMaterialStiffnessMatrix(fullAnswer, mode, gp, tStep);
437  }
438 }
439 
440 
441 void
443  MatResponseMode mode,
444  GaussPoint *gp,
445  TimeStep *tStep)
446 
447 //
448 // returns receiver's 1dMaterialStiffnessMAtrix
449 // (1d case ==> sigma_y = sigma_z = tau_yz = tau_zx = tau_xy = 0.)
450 {
451  FloatMatrix fullAnswer;
452  if ( mode == ElasticStiffness ) {
453  this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, mode, gp, tStep);
454  } else {
455  this->giveMaterialStiffnessMatrix(fullAnswer, mode, gp, tStep);
457  }
458 }
459 
460 
461 
462 void
464  MatResponseMode mode,
465  GaussPoint *gp,
466  TimeStep *tStep)
467 //
468 // returns receiver's 2dBeamLayerStiffMtrx.
469 // (2dPlaneStres ==> sigma_z = tau_xz = tau_yz = 0.)
470 //
471 // standard method from Material Class overloaded, because no inversion is needed.
472 // the reduction from 3d case will not work
473 // this implementation should be faster.
474 {
475  FloatMatrix fullAnswer;
476  if ( mode == ElasticStiffness ) {
477  this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, mode, gp, tStep);
478  } else {
479  this->giveMaterialStiffnessMatrix(fullAnswer, mode, gp, tStep);
481  }
482 }
483 
484 
485 void
487  MatResponseMode mode,
488  GaussPoint *gp,
489  TimeStep *tStep)
490 //
491 // returns receiver's 2dPlateLayerMtrx
492 // (2dPlaneStres ==> sigma_z = tau_xz = tau_yz = 0.)
493 //
494 // standard method from Material Class overloaded, because no inversion is needed.
495 // the reduction from 3d case will not work
496 // this implementation should be faster.
497 {
498  FloatMatrix fullAnswer;
499  if ( mode == ElasticStiffness ) {
500  this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, mode, gp, tStep);
501  } else {
502  this->giveMaterialStiffnessMatrix(fullAnswer, mode, gp, tStep);
504  }
505 }
506 
507 
508 void
510  const FloatArray &strainIncrement,
511  TimeStep *tStep)
512 //
513 // computest the elastic stress increment
514 // from stressIncrement in full stress strain space
515 //
516 {
517  FloatMatrix materialMatrix;
518 
519  if ( strainIncrement.giveSize() == 0 ) {
520  answer.clear();
521  return;
522  }
523 
524  this->giveEffectiveMaterialStiffnessMatrix(materialMatrix, TangentStiffness, gp, tStep);
525  answer.beProductOf(materialMatrix, strainIncrement);
526 }
527 
528 
529 void
531  GaussPoint *gp,
532  FloatArray *currentStressVector,
533  FloatArray *currentPlasticStrainVector,
534  FloatArray *strainIncrement3d,
535  TimeStep *tStep,
536  double &lambda)
537 //
538 // Computes full form of Plastic stiffness Matrix at given state.
539 // gp is used only and only for setting proper MaterialMode ()
540 // returns proportionality factor lambda also if strainIncrement3d != NULL
541 {
542  StructuralCrossSection *crossSection = static_cast< StructuralCrossSection * >
543  ( gp->giveElement()->giveCrossSection() );
544  FloatMatrix de, *yeldStressGradMat, *loadingStressGradMat;
545  FloatMatrix fsde, gsfsde;
546  FloatArray *yeldStressGrad, *loadingStressGrad;
547  FloatArray help;
548  double denominator, nominator;
549  //
550  // force de to be elastic even if gp in plastic state
551  // to do this, a flag in this class exist -> ForceElasticResponce
552  // if this flag is set to nonzero, function this::Give3dMaterialStiffnessMatrix
553  // will be forced to return only elasticPart of 3dMaterialStiffnessMatrix
554  // This function also set the ForceElasticFlagResponce to zero.
555  //
556  this->giveEffectiveMaterialStiffnessMatrix(de, TangentStiffness, gp, tStep);
557 
558  yeldStressGrad = this->GiveYCStressGradient(gp, currentStressVector,
559  currentPlasticStrainVector);
560  crossSection->imposeStressConstrainsOnGradient(gp, yeldStressGrad);
561  yeldStressGradMat = new FloatMatrix(*yeldStressGrad, 1); // transpose
562 
563  loadingStressGrad = this->GiveLCStressGradient(gp, currentStressVector,
564  currentPlasticStrainVector);
565 
566  crossSection->imposeStrainConstrainsOnGradient(gp, loadingStressGrad);
567  loadingStressGradMat = new FloatMatrix(*yeldStressGrad);
568 
569  help.beProductOf(de, * loadingStressGrad);
570  delete loadingStressGrad;
571 
572  fsde.beProductOf(* yeldStressGradMat, de);
573  delete yeldStressGradMat;
574  gsfsde.beProductOf(* loadingStressGradMat, fsde);
575  delete loadingStressGradMat;
576 
577  denominator = 0.;
578  for ( int i = 1; i <= 6; i++ ) {
579  denominator += yeldStressGrad->at(i) * help.at(i);
580  }
581 
582  answer.beProductOf(de, gsfsde);
583 
584  answer.times( -( 1. / denominator ) );
585 
586  if ( strainIncrement3d != NULL ) { // compute proportional factor lambda
587  nominator = 0.;
588  help.beProductOf(de, * strainIncrement3d);
589  for ( int i = 1; i <= 6; i++ ) {
590  nominator += yeldStressGrad->at(i) * help.at(i);
591  }
592 
593  lambda = nominator / denominator;
594  }
595 
596  delete yeldStressGrad;
597 }
598 
599 
600 FloatArray *
602  FloatArray *stressVector3d,
603  FloatArray *plasticVector3d)
604 //
605 // returns the stress correction -> correction is in the direction of
606 // the normal to the yield surface
607 // in full stress strain space
608 {
609  FloatArray *yeldStressGrad, *stressCorrection;
610  StructuralCrossSection *crossSection = static_cast< StructuralCrossSection * >
611  ( gp->giveElement()->giveCrossSection() );
612  double f3, help;
613 
614  yeldStressGrad = this->GiveYCStressGradient(gp,
615  stressVector3d,
616  plasticVector3d);
617  crossSection->imposeStressConstrainsOnGradient(gp, yeldStressGrad);
618 
619  f3 = this->computeYCValueAt(gp, stressVector3d, plasticVector3d);
620 
621  help = 0.;
622  for ( int j = 1; j <= 6; j++ ) {
623  help += yeldStressGrad->at(j) * yeldStressGrad->at(j);
624  }
625 
626  stressCorrection = new FloatArray(*yeldStressGrad);
627  stressCorrection->times(-f3 / help);
628 
629 
630  delete yeldStressGrad;
631  return stressCorrection;
632 }
633 
634 
637 {
639  if ( result != IRRT_OK ) return result;
640  return this->giveLinearElasticMaterial()->initializeFrom(ir);
641 }
642 
643 
644 double
646 // Returns the value of the property aProperty (e.g. the Young's modulus
647 // 'E') of the receiver.
648 {
649  double value = 0.0;
650 
651  if ( propertyDictionary.includes(aProperty) ) {
652  value = propertyDictionary.at(aProperty);
653  } else {
654  if ( linearElasticMaterial ) {
655  value = this->linearElasticMaterial->give(aProperty, gp);
656  } else {
657  OOFEM_ERROR("property not defined");
658  }
659  }
660 
661  return value;
662 }
663 
664 
667 /*
668  * creates new material status corresponding to this class
669  */
670 {
671  return new PerfectlyPlasticMaterialStatus(1, this->giveDomain(), gp);
672 }
673 
674 
675 int
677 {
678  PerfectlyPlasticMaterialStatus *status = static_cast< PerfectlyPlasticMaterialStatus * >( this->giveStatus(gp) );
679  if ( type == IST_PlasticStrainTensor ) {
680  const FloatArray &ep = status->givePlasticStrainVector();
683  return 1;
684  } else if ( type == IST_PrincipalPlasticStrainTensor ) {
685  FloatArray st;
686 
687  const FloatArray &s = status->givePlasticStrainVector();
690 
691  this->computePrincipalValues(answer, st, principal_strain);
692  return 1;
693  } else {
694  return StructuralMaterial :: giveIPValue(answer, gp, type, tStep);
695  }
696 }
697 
698 
699 //##################################################################################################
700 
701 
703  StructuralMaterialStatus(n, d, g), plasticStrainVector(), plasticStrainIncrementVector()
704 {
705  temp_yield_flag = yield_flag = 0; // elastic case at beginning
706 }
707 
708 
710 { }
711 
712 
715 //
716 // saves full information stored in this Status
717 //
718 {
719  contextIOResultType iores;
720 
721  if ( ( iores = StructuralMaterialStatus :: saveContext(stream, mode, obj) ) != CIO_OK ) {
722  THROW_CIOERR(iores);
723  }
724 
725  // write a raw data
726  if ( !stream.write(yield_flag) ) {
728  }
729 
730  if ( ( iores = plasticStrainVector.storeYourself(stream) ) != CIO_OK ) {
731  THROW_CIOERR(iores);
732  }
733 
734  // return result back
735  return CIO_OK;
736 }
737 
738 
741 //
742 // restore state variables from stream
743 //
744 {
745  contextIOResultType iores;
746 
747  if ( ( iores = StructuralMaterialStatus :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
748  THROW_CIOERR(iores);
749  }
750 
751  // read raw data
752  if ( !stream.read(yield_flag) ) {
754  }
755 
756  if ( ( iores = plasticStrainVector.restoreYourself(stream) ) != CIO_OK ) {
757  THROW_CIOERR(iores);
758  }
759 
760  // return result back
761  return CIO_OK;
762 }
763 
764 
765 void
767 {
769  fprintf(file, "status { yield_flag %d}\n", yield_flag);
770 }
771 
772 
773 void
775 //
776 // initialize record at the begining of new load step
777 //
778 {
780 
781  if ( plasticStrainVector.giveSize() == 0 ) {
784  }
785 
786  if ( plasticStrainIncrementVector.giveSize() == 0 ) {
789  } else {
791  }
792 
794 }
795 
796 
797 void
799 {
801 
804 
806 }
807 } // end namespace oofem
static int giveSizeOfVoigtSymVector(MaterialMode mmode)
Returns the size of symmetric part of a reduced stress/strain vector according to given mode...
CrossSection * giveCrossSection()
Definition: element.C:495
virtual void givePlaneStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane stress stiffness matrix of receiver.
virtual FloatArray * imposeStressConstrainsOnGradient(GaussPoint *gp, FloatArray *gradientStressVector3d)
Returns modified gradient of stress vector, which is used to bring stresses back to yield surface...
virtual FloatArray * GiveYCStressGradient(GaussPoint *gp, FloatArray *, FloatArray *)
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Definition: structuralms.C:96
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
static void computePrincipalValues(FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode)
Common functions for convenience.
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition: gausspoint.h:191
const FloatArray & givePlasticStrainVector() const
void letTempStrainVectorBe(const FloatArray &v)
Assigns tempStrainVector to given vector v.
Definition: structuralms.h:137
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Returns material status of receiver in given integration point.
Definition: material.C:244
virtual void updateTempLC(GaussPoint *gp, FloatArray *, FloatArray *)
GaussPoint * gp
Associated integration point.
Class and object Domain.
Definition: domain.h:115
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
For computing principal strains from engineering strains.
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
double & at(int aKey)
Returns the value of the pair which key is aKey.
Definition: dictionary.C:106
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
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
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
virtual FloatArray * imposeStrainConstrainsOnGradient(GaussPoint *gp, FloatArray *gradientStressVector3d)
Returns modified gradient of strain vector, which is used to compute plastic strain increment...
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
virtual void updateTempYC(GaussPoint *gp, FloatArray *, FloatArray *)
void computeTrialStressIncrement(FloatArray &answer, GaussPoint *gp, const FloatArray &strainIncrement, TimeStep *tStep)
virtual void updateIfFailure(GaussPoint *gp, FloatArray *stressVector3d, FloatArray *PlasticStrainVector3d)
This class implements a structural material status information.
Definition: structuralms.h:65
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 ...
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
This class implements associated Material Status to PerfectlyPlasticMaterial.
bool includes(int aKey)
Checks if dictionary includes given key.
Definition: dictionary.C:126
static void giveFullSymMatrixForm(FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode)
Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Dictionary propertyDictionary
Property dictionary.
Definition: material.h:105
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Element * giveElement()
Returns corresponding element to receiver.
Definition: gausspoint.h:188
General IO error.
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Print receiver&#39;s output to given stream.
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 contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
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.
FloatArray * GiveStressCorrectionBackToYieldSurface(GaussPoint *gp, FloatArray *stressVector3d, FloatArray *plasticVector3d)
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 1d stiffness matrix of receiver.
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
Definition: floatarray.C:341
static void giveFullSymVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form.
virtual int hasMaterialModeCapability(MaterialMode mode)
Tests if material supports material mode.
LinearElasticMaterial * linearElasticMaterial
virtual int write(const int *data, int count)=0
Writes count integer values from array pointed by data.
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of material property &#39;aProperty&#39;.
Definition: material.C:52
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 giveReducedSymVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form.
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of material property &#39;aProperty&#39;.
virtual MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual void give2dBeamLayerStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d beam layer stiffness matrix of receiver.
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 ...
void times(double f)
Multiplies receiver by factor f.
Definition: floatmatrix.C:1594
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
contextIOResultType restoreYourself(DataStream &stream)
Definition: floatarray.C:895
LinearElasticMaterial * giveLinearElasticMaterial()
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
#define STRAIN_STEPS
Definition: material.h:56
virtual FloatArray * GiveLCStressGradient(GaussPoint *gp, FloatArray *, FloatArray *)
virtual void givePlaneStrainStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane strain stiffness matrix of receiver.
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: structuralms.C:157
Abstract base class representing a material status information.
Definition: matstatus.h:84
PerfectlyPlasticMaterialStatus(int n, Domain *d, GaussPoint *g)
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void giveRealStressVector(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
Computes the real stress vector for given total strain and integration point.
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
void letTempStressVectorBe(const FloatArray &v)
Assigns tempStressVector to given vector v.
Definition: structuralms.h:135
const FloatArray & giveStressVector() const
Returns the const pointer to receiver&#39;s stress vector.
Definition: structuralms.h:107
virtual void giveEffectiveMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
virtual int hasMaterialModeCapability(MaterialMode mode)
Tests if material supports material mode.
Class representing the general Input Record.
Definition: inputrecord.h:101
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Definition: floatmatrix.C:1023
virtual void giveMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
void computePlasticStiffnessAt(FloatMatrix &answer, GaussPoint *gp, FloatArray *currentStressVector, FloatArray *currentPlasticStrainVector, FloatArray *strainIncrement3d, TimeStep *tStep, double &lambda)
void times(double s)
Multiplies receiver with scalar.
Definition: floatarray.C:818
virtual void givePlateLayerStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d plate layer stiffness matrix of receiver.
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
Abstract base class for all "structural" constitutive models.
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
Definition: structuralms.C:108
Domain * giveDomain() const
Definition: femcmpnn.h:100
virtual double computeYCValueAt(GaussPoint *gp, FloatArray *, FloatArray *)
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:337
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
Abstract base class for all structural cross section models.
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.
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 negated()
Switches the sign of every coefficient of receiver.
Definition: floatarray.C:739
const FloatArray & giveStrainVector() const
Returns the const pointer to receiver&#39;s strain vector.
Definition: structuralms.h:105
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: material.C:89
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
static void giveReducedSymMatrixForm(FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode)
Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form.
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:30 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011