OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
concrete2.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 "concrete2.h"
36 #include "domain.h"
37 #include "floatmatrix.h"
38 #include "gausspoint.h"
39 #include "mathfem.h"
40 #include "timestep.h"
41 #include "datastream.h"
42 #include "contextioerr.h"
43 #include "classfactory.h"
44 
45 namespace oofem {
46 REGISTER_Material(Concrete2);
47 
49 {
51 }
52 
53 
55 {
56  delete linearElasticMaterial;
57 }
58 
61 {
62  IRResultType result; // Required by IR_GIVE_FIELD macro
63 
75 
76  // stirrups constants
83 
84  result = this->linearElasticMaterial->initializeFrom(ir);
85  if ( result != IRRT_OK ) {
86  return result;
87  }
88  return Material :: initializeFrom(ir);
89 }
90 
91 
92 double
93 Concrete2 :: give(int aProperty, GaussPoint *gp)
94 // Returns the value of the property aProperty (e.g. the Young's modulus
95 // 'E') of the receiver.
96 {
97  double value;
98 
99  switch ( aProperty ) {
100  case c2_SCCC:
101  return this->SCCC;
102 
103  case c2_SCCT:
104  return this->SCCT;
105 
106  case c2_EPP:
107  return this->EPP;
108 
109  case c2_EPU:
110  return this->EPU;
111 
112  case c2_EOPP:
113  return this->EOPP;
114 
115  case c2_EOPU:
116  return this->EOPU;
117 
118  case c2_SHEARTOL:
119  return this->SHEARTOL;
120 
121  case c2_E:
122  return this->E;
123 
124  case c2_n:
125  return this->n;
126 
127  case stirr_E:
128  return this->stirrE;
129 
130  case stirr_Ft:
131  return this->stirrFt;
132 
133  case stirr_A:
134  return this->stirrA;
135 
136  case stirr_TOL:
137  return this->stirrTOL;
138 
139  case stirr_EREF:
140  return this->stirrEREF;
141 
142  case stirr_LAMBDA:
143  return this->stirrLAMBDA;
144 
145  case c2_IS_PLASTIC_FLOW:
146  return this->IS_PLASTIC_FLOW;
147 
148  case c2_IFAD:
149  return this->IFAD;
150 
151  default:
152  if ( propertyDictionary.includes(aProperty) ) {
153  value = propertyDictionary.at(aProperty);
154  return value;
155  } else {
156  return this->linearElasticMaterial->give(aProperty, gp);
157  // error ("give: property not defined");
158  }
159  }
160 }
161 
162 
163 void
165  GaussPoint *gp,
166  const FloatArray &totalStrain,
167  TimeStep *tStep)
168 //
169 // returns total stress vector of receiver according to
170 // previous level of stress and current
171 // strain increment, the only way, how to correctly update gp records
172 //
173 // Plane stress or uniaxial stress + transverse shear in
174 // concrete layers with transverse stirrups.
175 //
176 //
177 // gp strains - EXX,EYY,2*EXY,2*EYZ,2*EXZ
178 //
179 // material constant of concrete
180 // stored in this c2_property:
181 //
182 // E - Young modulus
183 // ny - Poisson ratio
184 // Ro - specific mass
185 // SCCC<=0 - pressure strength
186 // SCCT>=0 - tension strength
187 // EPP>=0 - treshold eff. plastic strain for softening in compress.
188 // EPU>=epp- ultimate eff. pl. strain
189 // EOPP>=0 - threshold vlumetric plastic strain for soft. in tension
190 // EOPU>=EOPP ultimate vol. pl. strain
191 // SHEARTOL threshold value of the relative shear deformation
192 // (psi**2/eef) at which shear is considered in layers. for
193 // lower r.s.d. the transverse shear remains elastic decoupled
194 // from bending. default value SHEARTOL = 0.01
195 // IS_PLASTIC_FLOW indicates that plastic flow (not deformation theory)
196 // is used in pressure.
197 // IFAD IFAD<=0 STATE VARIABLES WILL NOT BE UPDATED
198 // >0 UPDATE S.V.
199 //
200 //
201 // state varibles of the layer are stored in corresponding Material Status
202 //
203 // SCCM current pressure strenght
204 // EPM max. eff. plastic strain
205 // SCTM current tension strenght
206 // E0PM max. vol. plastic strain
207 // SRF current stress in stirrups
208 // SEZ current strain in transverse (z) direction.
209 // if flow teory of compresion then also available are:
210 // EXX_PLAST, EYY_PLAST, EZZ_PLAST, EXY_PLAST=2.EXY_PLAST, EYZ_PLAST (2*),EXZ_PLAST
211 // - components of plastic strain associated with epp
212 //
213 // Note: formulated in Full stress strain space
214 {
215  double us, SCC, SCT, psi2, comp, tol = 0.0, crit = 0.0, nez = 0.0, effstold = 0.0, dezold = 0.0, effst = 0.0, dez = 0.0;
216  double epsult = 0.0, ep1 = -1.0, ep2 = -1.0, ep3 = -1.0, relus, G, ez;
217  int ifsh, i, ifplas, ifupd;
218  FloatArray plasticStrainVector, plasticStrainIncrementVector;
219 
220  Concrete2MaterialStatus *status = static_cast< Concrete2MaterialStatus * >( this->giveStatus(gp) );
221  FloatArray currentStress, currentStrain, currentEffStrain(6), pVal, ep, pStress, strainIncr, plasticStrain, help, reducedStrain, helpR;
222  FloatMatrix pDir;
223 
224  pDir.resize(3, 3); // principal directions
225  pStress.resize(3); // principal stresses
226  strainIncr.resize(1);
227 
228  this->initTempStatus(gp);
229 
230  // subtract stress independent part
231  // note: eigenStrains (temperature) is not contained in mechanical strain stored in gp
232  // therefore it is necessary to subtract always the total eigen strain value
233  this->giveStressDependentPartOfStrainVector(reducedStrain, gp,
234  totalStrain, tStep, VM_Total);
235 
236  StructuralMaterial :: giveFullSymVectorForm( currentStrain, reducedStrain, gp->giveMaterialMode() );
238 
239  if ( status->giveTempCurrentTensionStrength() < 0. ) {
240  //
241  // called for the firs time - init flags
242  //
243  status->giveTempCurrentPressureStrength() = this->give(c2_SCCC, gp);
244  status->giveTempCurrentTensionStrength() = this->give(c2_SCCT, gp);
245 
246  status->giveCurrentPressureStrength() = this->give(c2_SCCC, gp);
247  status->giveCurrentTensionStrength() = this->give(c2_SCCT, gp);
248  }
249 
250  plasticStrainVector = status->givePlasticStrainVector();
251  StructuralMaterial :: giveFullSymVectorForm( plasticStrain, plasticStrainVector, gp->giveMaterialMode() );
252 
253  ep.resize(3);
254 
255  //for (i =1; i <= 6; i++) currentStrain.at(i) += strainIncrement.at(i);
256 
257  ez = status->giveTempCurrentStrainInZDir();
258  us = 0.;
259  // eep = 0.;
260  //
261  // Concrete strenghts, SCC and SCT are actual strenghts
262  // which may change in the iteration if there is any
263  //
264  SCC = this->give(c2_SCCC, gp);
265  if ( this->give(c2_EPP, gp) != 0. ) {
266  SCC = status->giveTempCurrentPressureStrength();
267  }
268 
269  SCT = this->give(c2_SCCT, gp);
270  if ( this->give(c2_EOPP, gp) != 0. ) {
271  SCT = status->giveTempCurrentTensionStrength();
272  }
273 
274  //
275  // create full 3d strain tensor
276  //
277 
278  if ( this->give(c2_IS_PLASTIC_FLOW, gp) > 0. ) {
279  currentEffStrain.at(1) = currentStrain.at(1) - plasticStrain.at(1);
280  currentEffStrain.at(2) = currentStrain.at(2) - plasticStrain.at(2);
281  currentEffStrain.at(3) = ez - plasticStrain.at(3);
282  currentEffStrain.at(4) = currentStrain.at(4) - plasticStrain.at(4);
283  currentEffStrain.at(5) = currentStrain.at(5) - plasticStrain.at(5);
284  currentEffStrain.at(6) = currentStrain.at(6) - plasticStrain.at(6);
285  } else {
286  currentEffStrain.at(1) = currentStrain.at(1);
287  currentEffStrain.at(2) = currentStrain.at(2);
288  currentEffStrain.at(3) = ez;
289  currentEffStrain.at(4) = currentStrain.at(4);
290  currentEffStrain.at(5) = currentStrain.at(5);
291  currentEffStrain.at(6) = currentStrain.at(6);
292  }
293 
294  //
295  // check for shear to be considered
296  //
297  ifsh = 0;
298  psi2 = currentEffStrain.at(4) * currentEffStrain.at(4) +
299  currentEffStrain.at(5) * currentEffStrain.at(5);
300  comp = sqrt(currentEffStrain.at(1) * currentEffStrain.at(1) +
301  currentEffStrain.at(2) * currentEffStrain.at(2) +
302  currentEffStrain.at(3) * currentEffStrain.at(3) +
303  currentEffStrain.at(6) * currentEffStrain.at(6) +
304  psi2);
305 
306  if ( comp == 0. ) { // strange
307  for ( i = 1; i <= 6; i++ ) { // directly modify gp-record
308  currentStress.at(i) = 0.;
309  }
310 
311  // note in status are reduced space variables
312  status->letTempStrainVectorBe(totalStrain);
313 
315 
316  status->letTempStressVectorBe(helpR);
317 
319 
320  plasticStrainIncrementVector = status->givePlasticStrainIncrementVector();
321  plasticStrainIncrementVector.subtract(plasticStrainVector);
322  plasticStrainIncrementVector.add(help);
323  status->letPlasticStrainIncrementVectorBe(plasticStrainIncrementVector);
324  // plasticStrain->negated()->add (status->givePlasticStrainVector());
325  // status->givePlasticStrainIncrementVector()-> add(plasticStrain);
326 
328  return;
329  }
330 
331  //
332  // Treshold relative shear strain to trigger shear computation
333  // is stored in SHEARTOL material constant. Once triggered (EZ.ne 0)
334  // then shear computation cannot be abandoned even if the
335  // shear strain drops bellow the limit.
336  //
337  if ( ( ez != 0. ) && ( ( sqrt(psi2) / comp ) < this->give(c2_SHEARTOL, gp) ) ) {
338  //
339  // no transverse shear - 2d case
340  //
341 
342  this->computePrincipalValDir(pVal, pDir, currentEffStrain, principal_strain);
343  //
344  // dtp2 procedure expects first two princ values to be in xy plane
345  // so sort results acordingly
346  //
347  int j = 1;
348  if ( pDir.at(3, 2) > pDir.at(3, j) ) {
349  j = 2; // find the eigen direction
350  }
351 
352  if ( pDir.at(3, 3) > pDir.at(3, j) ) {
353  j = 3; // normal to xy plane
354  }
355 
356  if ( j != 3 ) {
357  // swap values
358  int k;
359  double swap;
360  for ( k = 1; k <= 3; k++ ) { // swap eigenvectors
361  swap = pDir.at(k, 3);
362  pDir.at(k, 3) = pDir.at(k, j);
363  pDir.at(k, j) = swap;
364  }
365 
366  swap = pVal.at(3);
367  pVal.at(3) = pVal.at(j);
368  pVal.at(j) = swap;
369  }
370 
371  this->dtp2(gp, pVal, pStress, ep, SCC, SCT, & ifplas);
372  } else {
373  ifsh = 1;
374  //
375  // optional iteration to achieve equilibrium in the z - direction.
376  // // iteration toleration limit tol (=0. -- no iteration)
377  // iteration only possible when stirrups present
378  // Number of loops nez
379  tol = 0.;
380  if ( this->give(stirr_E, gp) > 0 ) {
381  tol = this->give(stirr_TOL, gp);
382  }
383 
384  crit = fabs( this->give(c2_SCCC, gp) ) + this->give(c2_SCCT, gp);
385  nez = 0;
386  effstold = 0.;
387  dezold = 0.;
388 label18:
389  //
390  // find principal strains and corresponding principal direction
391  //
392  this->computePrincipalValDir(pVal, pDir, currentEffStrain, principal_strain);
393  //
394  // plasticity condition
395  //
396  this->dtp3(gp, pVal, pStress, ep, SCC, SCT, & ifplas);
397  //
398  // transverse equilibrium
399  // transv strain ez (flags -> at(SEZ)) its incr. dez ,
400  // stirrup stress flags->at(SRF)
401  //
402  // unbalanced transverse stress us
403  // concrete component
404  //
405  us = ( pDir.at(3, 1) * pStress.at(1) * pDir.at(3, 1) +
406  pDir.at(3, 2) * pStress.at(2) * pDir.at(3, 2) +
407  pDir.at(3, 3) * pStress.at(3) * pDir.at(3, 3) );
408  //
409  // Stirrups component
410  //
411  if ( this->give(stirr_E, gp) > 0. ) {
412  us += status->giveTempCurrentStressInStirrups() * this->give(stirr_A, gp);
413  }
414 
415  //
416  // increment of transverse strain ez
417  // effective modulus of concrete in z direction
418  // f(1)*red**4 where red is sinus of the angle
419  // between z axis and max. principal strain dir.
420  // when us becomes positive (may happen if the
421  // effective modulus was too small and, hence, the
422  // last change dez of the transverse strain ez
423  // too large in the previous increment), it loses sense.
424  // the Jacobi iteration delivers just approximate princ.
425  // vectors which causes spurious pulses in us and fluctuations
426  // in sxz,syz,qx,qz. they decrease
427  // when toli->0. they are successfully managed by z equil.
428  // iteration and do not harm the solution.
429  //
430  if ( ( ifplas ) && ( us < 0. ) ) {
431  // transverse strain increment if concrete is plastic in tension
432  // and further Z-extension is necessary (US<0).
433  // effst is an iteration factor to get transverse equilibrium
434  if ( pVal.at(1) > pVal.at(2) ) {
435  if ( pVal.at(1) > pVal.at(3) ) {
436  effst = pDir.at(3, 1);
437  } else {
438  effst = pDir.at(3, 3);
439  }
440  } else {
441  if ( pVal.at(2) > pVal.at(3) ) {
442  effst = pDir.at(3, 2);
443  } else {
444  effst = pDir.at(3, 3);
445  }
446  }
447 
448  effst = ( 1. - effst * effst );
449  effst *= effst;
450  if ( effst < 1.e-6 ) {
451  effst = 1.e-6;
452  }
453 
454  effst *= this->give(c2_E, gp);
455  if ( this->give(stirr_E, gp) > 0. ) {
456  effst += this->give(stirr_E, gp) * this->give(stirr_A, gp);
457  }
458 
459  effst *= 2.0;
460  } else {
461  // concrete is elastic
462  effst = this->give(c2_E, gp) * 1.1;
463  }
464 
465  dez = -us / effst;
466  // if the increment of ez has changed sign in this z-equil.
467  // iteration loop, raise the iteration factor by adding up
468  // the two last values
469  if ( nez > 1 ) {
470  if ( dez * dezold < 0. ) {
471  effst += effstold;
472  dez = -us / effst;
473  }
474 
475  dezold = dez;
476  effstold = effst;
477  }
478 
479  //
480  ez += dez;
481  if ( this->give(stirr_E, gp) > 0. ) {
482  strainIncr.at(1) = dez;
483  this->updateStirrups(gp, strainIncr, tStep);
484  // updates flags->at(SEZ)
485  // stirr(dez,sv(l3),fst);
486  }
487  }
488 
489  //
490  // strain softening
491  //
492  ifupd = 0;
493  if ( this->give(c2_IFAD, gp) ) {
494  if ( this->give(c2_EOPP, gp) >= 0. ) {
495  if ( this->give(c2_EOPU, gp) > 0. ) {
496  epsult = this->give(c2_EOPU, gp);
497  } else {
498  // position dependent strain softening
499  // NOT IMPLEMENTED
500  epsult = this->give(c2_EOPU, gp);
501  }
502  }
503 
504  this->strsoft(gp, epsult, ep, ep1, ep2, ep3, SCC, SCT, ifupd);
505  }
506 
507  //
508  // End of the Ez iteration loop
509  //
510  if ( ifsh ) {
511  relus = fabs(us) / crit;
512  if ( ( tol > 0. ) && ( relus > tol ) ) {
513  //
514  // add the increment of EZ to the actual strain tenzor
515  //
516  currentEffStrain.at(1) += pDir.at(3, 1) * dez * pDir.at(3, 1);
517  currentEffStrain.at(2) += pDir.at(3, 2) * dez * pDir.at(3, 2);
518  currentEffStrain.at(3) += pDir.at(3, 3) * dez * pDir.at(3, 3);
519  currentEffStrain.at(4) += pDir.at(3, 2) * dez * pDir.at(3, 3);
520  currentEffStrain.at(5) += pDir.at(3, 1) * dez * pDir.at(3, 3);
521  currentEffStrain.at(6) += pDir.at(3, 1) * dez * pDir.at(3, 2);
522  // back to EZ iteration start
523  //
524  if ( this->give(c2_IS_PLASTIC_FLOW, gp) && this->give(c2_IFAD, gp) ) {
525  plasticStrain.at(1) += ( pDir.at(1, 1) * ep1 * pDir.at(1, 1) +
526  pDir.at(1, 2) * ep2 * pDir.at(1, 2) +
527  pDir.at(1, 3) * ep3 * pDir.at(1, 3) );
528  plasticStrain.at(2) += ( pDir.at(2, 1) * ep1 * pDir.at(2, 1) +
529  pDir.at(2, 2) * ep2 * pDir.at(2, 2) +
530  pDir.at(2, 3) * ep3 * pDir.at(2, 3) );
531  plasticStrain.at(3) += ( pDir.at(3, 1) * ep1 * pDir.at(3, 1) +
532  pDir.at(3, 2) * ep2 * pDir.at(3, 2) +
533  pDir.at(3, 3) * ep3 * pDir.at(3, 3) );
534  plasticStrain.at(4) += 2. * ( pDir.at(2, 1) * ep1 * pDir.at(3, 1) +
535  pDir.at(2, 2) * ep2 * pDir.at(3, 2) +
536  pDir.at(2, 3) * ep3 * pDir.at(3, 3) );
537  plasticStrain.at(5) += 2. * ( pDir.at(1, 1) * ep1 * pDir.at(3, 1) +
538  pDir.at(1, 2) * ep2 * pDir.at(3, 2) +
539  pDir.at(1, 3) * ep3 * pDir.at(3, 3) );
540  plasticStrain.at(6) += 2. * ( pDir.at(1, 1) * ep1 * pDir.at(2, 1) +
541  pDir.at(1, 2) * ep2 * pDir.at(2, 2) +
542  pDir.at(1, 3) * ep3 * pDir.at(2, 3) );
543  }
544 
545  pVal.zero();
546  goto label18;
547  }
548  }
549 
550  //
551  // update compressive plastic strains.
552  //
553  if ( ifupd ) {
554  if ( this->give(c2_IFAD, gp) ) {
555  if ( this->give(c2_IS_PLASTIC_FLOW, gp) ) {
556  // -
557  if ( ifsh == 0 ) {
558  //
559  // no transverse shear, 1-2 to x-y
560  //
561  plasticStrain.at(1) += ( pDir.at(1, 1) * ep1 * pDir.at(1, 1) +
562  pDir.at(1, 2) * ep2 * pDir.at(1, 2) );
563  plasticStrain.at(2) += ( pDir.at(2, 1) * ep1 * pDir.at(2, 1) +
564  pDir.at(2, 2) * ep2 * pDir.at(2, 2) );
565  plasticStrain.at(6) += 2. * ( pDir.at(1, 1) * ep1 * pDir.at(2, 1) +
566  pDir.at(1, 2) * ep2 * pDir.at(2, 2) );
567  } else {
568  //
569  // transverse shear present , transform back to x,y,z
570  //
571  plasticStrain.at(1) += ( pDir.at(1, 1) * ep1 * pDir.at(1, 1) +
572  pDir.at(1, 2) * ep2 * pDir.at(1, 2) +
573  pDir.at(1, 3) * ep3 * pDir.at(1, 3) );
574  plasticStrain.at(2) += ( pDir.at(2, 1) * ep1 * pDir.at(2, 1) +
575  pDir.at(2, 2) * ep2 * pDir.at(2, 2) +
576  pDir.at(2, 3) * ep3 * pDir.at(2, 3) );
577  plasticStrain.at(3) += ( pDir.at(3, 1) * ep1 * pDir.at(3, 1) +
578  pDir.at(3, 2) * ep2 * pDir.at(3, 2) +
579  pDir.at(3, 3) * ep3 * pDir.at(3, 3) );
580  plasticStrain.at(4) += 2. * ( pDir.at(2, 1) * ep1 * pDir.at(3, 1) +
581  pDir.at(2, 2) * ep2 * pDir.at(3, 2) +
582  pDir.at(2, 3) * ep3 * pDir.at(3, 3) );
583  plasticStrain.at(5) += 2. * ( pDir.at(1, 1) * ep1 * pDir.at(3, 1) +
584  pDir.at(1, 2) * ep2 * pDir.at(3, 2) +
585  pDir.at(1, 3) * ep3 * pDir.at(3, 3) );
586  plasticStrain.at(6) += 2. * ( pDir.at(1, 1) * ep1 * pDir.at(2, 1) +
587  pDir.at(1, 2) * ep2 * pDir.at(2, 2) +
588  pDir.at(1, 3) * ep3 * pDir.at(2, 3) );
589  }
590  }
591  }
592  }
593 
594  //
595  // stresses back , 1-2 to x-z , sxz= tau=-sigxz
596  //
597  // 36
598  if ( ifsh == 0 ) {
599  G = 0.5 * this->give(c2_E, gp) / ( 1.0 + this->give(c2_n, gp) );
600 
601  currentStress.at(1) = ( pDir.at(1, 1) * pStress.at(1) * pDir.at(1, 1) +
602  pDir.at(1, 2) * pStress.at(2) * pDir.at(1, 2) );
603  currentStress.at(2) = ( pDir.at(2, 1) * pStress.at(1) * pDir.at(2, 1) +
604  pDir.at(2, 2) * pStress.at(2) * pDir.at(2, 2) );
605  currentStress.at(3) = 0.;
606  currentStress.at(4) = currentStrain.at(4) * G;
607  currentStress.at(5) = currentStrain.at(5) * G;
608  currentStress.at(6) = ( pDir.at(1, 1) * pStress.at(1) * pDir.at(2, 1) +
609  pDir.at(1, 2) * pStress.at(2) * pDir.at(2, 2) );
610  } else {
611  currentStress.at(1) = ( pDir.at(1, 1) * pStress.at(1) * pDir.at(1, 1) +
612  pDir.at(1, 2) * pStress.at(2) * pDir.at(1, 2) +
613  pDir.at(1, 3) * pStress.at(3) * pDir.at(1, 3) );
614  currentStress.at(2) = ( pDir.at(2, 1) * pStress.at(1) * pDir.at(2, 1) +
615  pDir.at(2, 2) * pStress.at(2) * pDir.at(2, 2) +
616  pDir.at(2, 3) * pStress.at(3) * pDir.at(2, 3) );
617  currentStress.at(3) = us;
618  currentStress.at(4) = ( pDir.at(2, 1) * pStress.at(1) * pDir.at(3, 1) +
619  pDir.at(2, 2) * pStress.at(2) * pDir.at(3, 2) +
620  pDir.at(2, 3) * pStress.at(3) * pDir.at(3, 3) );
621  currentStress.at(5) = ( pDir.at(1, 1) * pStress.at(1) * pDir.at(3, 1) +
622  pDir.at(1, 2) * pStress.at(2) * pDir.at(3, 2) +
623  pDir.at(1, 3) * pStress.at(3) * pDir.at(3, 3) );
624  currentStress.at(6) = ( pDir.at(1, 1) * pStress.at(1) * pDir.at(2, 1) +
625  pDir.at(1, 2) * pStress.at(2) * pDir.at(2, 2) +
626  pDir.at(1, 3) * pStress.at(3) * pDir.at(2, 3) );
627  }
628 
629  if ( this->give(c2_IFAD, gp) ) {
630  status->giveTempCurrentStrainInZDir() = ez;
631  }
632 
633  //totalStrainVector = status->giveStrainVector();
634  //totalStrainVector.add (totalStrainIncrement);
635  status->letTempStrainVectorBe(totalStrain);
636 
638  status->letTempStressVectorBe(helpR);
639 
641  plasticStrainIncrementVector = status->givePlasticStrainIncrementVector();
642  plasticStrainIncrementVector.subtract(plasticStrainVector);
643  plasticStrainIncrementVector.add(help);
644  status->letPlasticStrainIncrementVectorBe(plasticStrainIncrementVector);
645 
646  //plasticStrain->negated()->add (status->givePlasticStrainVector());
647  //status->givePlasticStrainIncrementVector()-> add(plasticStrain);
648 
649  StructuralMaterial :: giveFullSymVectorForm( answer, currentStress, gp->giveMaterialMode() );
650 }
651 
652 
653 void
655  double SCC, double SCT, int *ifplas)
656 //
657 // DEFORMATION THEORY OF PLASTICITY, PRNCIPAL STRESSES
658 // CONDITION OF PLASTICITY.
659 //
660 // e - principal strains
661 // s - principal stresses (output)
662 // SCT,SCC - tensile and compressive strengths
663 // ifplas - =0 if material is elastic
664 // =1 if plastic in tension
665 //
666 //
667 // state varibles of the layer are stored
668 // in corresponding gp mat Status
669 //
670 // SCCM current pressure strenght
671 // EPM max. eff. plastic strain
672 // SCTM current tension strenght
673 // E0PM max. vol. plastic strain
674 // SRF current stress in stirrups
675 // SEZ current strain in transverse (z) direction.
676 // if flow teory of compresion then also available are:
677 // EXX_PLAST, EYY_PLAST, EZZ_PLAST, EXY_PLAST=2.EXY_PLAST, EYZ_PLAST (2*),EXZ_PLAST
678 // - components of plastic strain associated with epp in gp->plasticStrainVector
679 //
680 {
681  int i, ii = 0, j, k, jj = 0, kk = 0;
682  double yy, ey, e0, yy1, yy2, s0, sci = 0.0, scj = 0.0, sck = 0.0;
683 
684  if ( e.giveSize() != 3 ) {
685  OOFEM_ERROR("principal strains size mismatch");
686  }
687 
688  if ( s.giveSize() != 3 ) {
689  OOFEM_ERROR("principal stress size mismatch");
690  }
691 
692  if ( ep.giveSize() != 3 ) {
693  OOFEM_ERROR("plastic strains size mismatch");
694  }
695 
696  * ifplas = 0;
697 
698  yy = this->give(c2_n, gp);
699  ey = this->give(c2_E, gp);
700  e0 = this->give(c2_E, gp) / ( 1. + yy ) / ( 1. - yy - yy );
701  yy1 = 1. - yy;
702  yy2 = 1. + yy;
703 
704  if ( this->give(c2_E, gp) < 1.e-6 ) {
705  for ( i = 1; i <= 3; i++ ) {
706  s.at(i) = 0.;
707  ep.at(i) = 0.;
708  }
709 
710  return;
711  }
712 
713  if ( this->give(c2_n, gp) >= 0.01 ) {
714  s.at(1) = e0 * ( yy1 * e.at(1) + yy * ( e.at(2) + e.at(3) ) );
715  s.at(2) = e0 * ( yy1 * e.at(2) + yy * ( e.at(1) + e.at(3) ) );
716  s.at(3) = e0 * ( yy1 * e.at(3) + yy * ( e.at(2) + e.at(1) ) );
717  //
718  // find intersection with side i , sig(i)=sci
719  //
720  s0 = 0.;
721  for ( i = 1; i <= 3; i++ ) {
722  if ( ( s.at(i) - SCT - s0 ) > 0. ) {
723  ii = i;
724  sci = SCT;
725  s0 = s.at(i) - SCT;
726  * ifplas = 1;
727  } else {
728  if ( ( SCC - s.at(i) - s0 ) > 0. ) {
729  ii = i;
730  sci = SCC;
731  s0 = SCC - s.at(i);
732  }
733  }
734  }
735 
736  if ( s0 <= 0. ) {
737  //
738  // linear elastic
739  //
740  ep.at(1) = 0.;
741  ep.at(2) = 0.;
742  ep.at(3) = 0.;
743  return;
744  }
745 
746  //
747  // Strain vector outside elastic cube
748  // Plastic state
749  //
750  j = iperm(ii, 3);
751  k = iperm(j, 3);
752  i = ii;
753  s.at(i) = sci;
754  s.at(j) = ( yy * sci + ey / yy2 * ( e.at(j) + yy * e.at(k) ) ) / yy1;
755  s.at(k) = ( yy * sci + ey / yy2 * ( e.at(k) + yy * e.at(j) ) ) / yy1;
756  s0 = 0.;
757  //
758  // find intersection with edges j and k
759  //
760  if ( ( s.at(j) - SCT - s0 ) > 0. ) {
761  jj = j;
762  kk = k;
763  scj = SCT;
764  s0 = s.at(j) - SCT;
765  }
766 
767  if ( ( SCC - s.at(j) - s0 ) > 0. ) {
768  jj = j;
769  kk = k;
770  scj = SCC;
771  s0 = SCC - s.at(j);
772  }
773 
774  if ( ( s.at(k) - SCT - s0 ) > 0. ) {
775  jj = k;
776  kk = j;
777  scj = SCT;
778  s0 = s.at(k) - SCT;
779  }
780 
781  if ( ( SCC - s.at(k) - s0 ) > 0. ) {
782  jj = k;
783  kk = j;
784  scj = SCC;
785  s0 = SCC - s.at(k);
786  } else {
787  //
788  // staying on the side i
789  //
790  if ( s0 <= 0. ) {
791  ep.at(i) = e.at(i) - sci / yy1 / e0 + yy / yy1 * ( e.at(j) + e.at(k) );
792  ep.at(j) = 0.;
793  ep.at(k) = 0.;
794  return;
795  }
796  }
797 
798  //
799  // edge j on the side i
800  //
801  j = jj;
802  k = kk;
803  s.at(j) = scj;
804  s.at(k) = yy * ( sci + scj ) + ey *e.at(k);
805  //
806  // corner ?
807  //
808  if ( ( s.at(k) - SCT ) > 0. ) {
809  sck = SCT;
810  } else {
811  if ( ( SCC - s.at(k) ) > 0. ) {
812  sck = SCC;
813  } else {
814  //
815  // staying on the edge j - i
816  //
817  ep.at(i) = e.at(i) + yy *e.at(k) + yy2 / ey * ( yy * scj - yy1 * sci );
818  ep.at(j) = e.at(j) + yy *e.at(k) + yy2 / ey * ( yy * sci - yy1 * scj );
819  ep.at(k) = 0.;
820  return;
821  }
822  }
823 
824  //
825  // corner sck on the edge j of the side i
826  //
827  s.at(k) = sck;
828  ep.at(i) = e.at(i) - ( sci - yy * ( scj + sck ) ) / ey;
829  ep.at(j) = e.at(j) - ( scj - yy * ( sck + sci ) ) / ey;
830  ep.at(k) = e.at(k) - ( sck - yy * ( sci + scj ) ) / ey;
831  return;
832  } else {
833  //
834  // Poisson ratio = 0.
835  //
836  for ( i = 1; i <= 3; i++ ) {
837  sci = 0.;
838  ep.at(i) = 0.;
839  s.at(i) = e.at(i) * ey;
840  if ( s.at(i) > SCT ) {
841  sci = SCT;
842  * ifplas = 1;
843  s.at(i) = sci;
844  ep.at(i) -= sci / ey;
845  } else {
846  if ( s.at(i) < SCC ) {
847  sci = SCC;
848  s.at(i) = sci;
849  ep.at(i) -= sci / ey;
850  }
851  }
852  }
853  }
854 }
855 
856 void
858  double SCC, double SCT, int *ifplas)
859 //
860 // DEFORMATION THEORY OF PLASTICITY, PRNCIPAL STRESSES
861 // CONDITION OF PLASTICITY. - PLANE STRESS
862 //
863 // e - principal strains
864 // s - principal stresses (output)
865 // SCT,SCC - tensile and compressive strengths
866 // ifplas - =0 if material is elastic
867 // =1 if plastic in tension
868 //
869 //
870 // state varibles of the layer are stored
871 // in gp->matInfo->flagDictionary = flags
872 //
873 // SCCM current pressure strenght
874 // EPM max. eff. plastic strain
875 // SCTM current tension strenght
876 // E0PM max. vol. plastic strain
877 // SRF current stress in stirrups
878 // SEZ current strain in transverse (z) direction.
879 // if flow teory of compresion then also available are:
880 // EXX_PLAST, EYY_PLAST, EZZ_PLAST, EXY_PLAST=2.EXY_PLAST, EYZ_PLAST (2*),EXZ_PLAST
881 // - components of plastic strain associated with epp
882 //
883 {
884  int i, j, k, jj = 0;
885  double yy, ey, e0, so, scj = 0.0, sck = 0.0, sci = 0.0;
886  // Dictionary *flags = gp->matInfo->flagDictionary;
887 
888  if ( e.giveSize() != 3 ) {
889  OOFEM_ERROR("principal strains size mismatch");
890  }
891 
892  if ( s.giveSize() != 3 ) {
893  OOFEM_ERROR("principal stress size mismatch");
894  }
895 
896  if ( ep.giveSize() != 3 ) {
897  OOFEM_ERROR("plastic strains size mismatch");
898  }
899 
900  * ifplas = 0;
901 
902  yy = this->give(c2_n, gp);
903  ey = this->give(c2_E, gp);
904 
905  * ifplas = 0;
906  if ( ey == 0. ) {
907  s.at(1) = s.at(2) = 0.;
908  ep.at(1) = ep.at(2) = 0.;
909  return;
910  }
911 
912  if ( yy > 0.01 ) {
913  e0 = ey / ( 1. - yy * yy );
914 
915  s.at(1) = e0 * ( e.at(1) + yy * e.at(2) );
916  s.at(2) = e0 * ( e.at(2) + yy * e.at(1) );
917  so = 0.;
918  //
919  // find intersection with the edges of the square plast. co
920  // dition - edge j
921  //
922  for ( j = 1; j <= 2; j++ ) {
923  if ( ( s.at(j) - SCT - so ) > 0. ) {
924  jj = j;
925  scj = SCT;
926  so = s.at(j) - SCT;
927  * ifplas = 1;
928  } else {
929  if ( ( SCC - s.at(j) - so ) > 0. ) {
930  jj = j;
931  scj = SCC;
932  so = SCC - s.at(j);
933  }
934  }
935  }
936 
937  if ( so <= 0. ) {
938  //
939  // LINEAR ELASTIC - NO INTERSECTION
940  //
941  ep.at(1) = 0.;
942  ep.at(2) = 0.;
943  return;
944  }
945 
946  //
947  // edge jj,j
948  //
949  k = iperm(jj, 2);
950  j = jj;
951  s.at(j) = scj;
952  s.at(k) = ey * e.at(k) + scj * yy;
953  //
954  // any corner?
955  //
956  if ( ( s.at(k) - SCT ) > 0. ) {
957  sck = SCT;
958  } else {
959  if ( ( SCC - s.at(k) ) > 0. ) {
960  sck = SCC;
961  } else {
962  //
963  // staying on the edge j
964  //
965  ep.at(j) = e.at(j) + yy *e.at(k) - scj / e0;
966  ep.at(k) = 0.;
967  return;
968  }
969  }
970 
971  //
972  // corner sck on the edge j
973  //
974  s.at(k) = sck;
975  ep.at(j) = e.at(j) - ( s.at(j) - yy * s.at(k) ) / ey;
976  ep.at(k) = e.at(k) - ( s.at(k) - yy * s.at(j) ) / ey;
977  return;
978  } else {
979  //
980  // poisson ratio =0.
981  //
982  for ( i = 1; i <= 2; i++ ) {
983  sci = 0.;
984  ep.at(i) = 0.;
985  s.at(i) = e.at(i) * ey;
986  if ( s.at(i) > SCT ) {
987  sci = SCT;
988  * ifplas = 1;
989  } else {
990  sci = SCC;
991  }
992 
993  s.at(i) = sci;
994  ep.at(i) = e.at(i) - sci / ey;
995  }
996 
997  return;
998  }
999 }
1000 
1001 
1002 
1003 
1004 void
1005 Concrete2 :: strsoft(GaussPoint *gp, double epsult, FloatArray &ep, double &ep1, double &ep2, double &ep3,
1006  double SCC, double SCT, int &ifupd)
1007 // material constant of concrete
1008 // stored in this.propertyDictionary
1009 //
1010 // E - Young modulus
1011 // ny - Poisson ratio
1012 // Ro - specific mass
1013 // SCCC<=0 - pressure strength
1014 // SCCT>=0 - tension strength
1015 // EPP>=0 - treshold eff. plastic strain for softening in compress.
1016 // EPU>=epp- ultimate eff. pl. strain
1017 // EOPP>=0 - threshold vlumetric plastic strain for soft. in tension
1018 // EOPU>=EOPP ultimate vol. pl. strain
1019 // SHEARTOL threshold value of the relative shear deformation
1020 // (psi**2/eef) at which shear is considered in layers. for
1021 // lower r.s.d. the transverse shear remains elastic decoupled
1022 // from bending. default value SHEARTOL = 0.01
1023 // IS_PLASTIC_FLOW indicates that plastic flow (not deformation theory)
1024 // is used in pressure.
1025 // IFAD IFAD<=0 STATE VARIABLES WILL NOT BE UPDATED
1026 // >0 UPDATE S.V.
1027 //
1028 //
1029 // state varibles of the layer are stored
1030 // in gp->matInfo->flagDictionary = flags
1031 //
1032 // SCCM current pressure strenght
1033 // EPM max. eff. plastic strain
1034 // SCTM current tension strenght
1035 // E0PM max. vol. plastic strain
1036 // SRF current stress in stirrups
1037 // SEZ current strain in transverse (z) direction.
1038 // if flow teory of compresion then also available are:
1039 // EXX_PLAST, EYY_PLAST, EZZ_PLAST, EXY_PLAST=2.EXY_PLAST, EYZ_PLAST (2*),EXZ_PLAST
1040 // - components of plastic strain associated with epp in gp->plasticStrainVector.
1041 //
1042 {
1043  Concrete2MaterialStatus *status = static_cast< Concrete2MaterialStatus * >( this->giveStatus(gp) );
1044  double eop, eopr, dep, d, eep, eepr;
1045  //
1046  // strain softening
1047  // tension
1048  //
1049  if ( this->give(c2_EOPP, gp) != 0. ) {
1050  eop = max(ep.at(1), 0.) + max(ep.at(2), 0.) +
1051  max(ep.at(3), 0.);
1052  if ( eop > status->giveTempMaxVolPlasticStrain() ) {
1053  //
1054  // Current plastic volum. strain is greater than max.
1055  // reached up to now.
1056  //
1057  // reduced plastic strain withe respect to the elastic
1058  // limit strain of the virgin concrete E0PR
1059  //
1060  eopr = eop - ( this->give(c2_SCCT, gp) - SCT ) / this->give(c2_E, gp);
1061  if ( eopr <= this->give(c2_EOPP, gp) ) {
1062  status->giveTempMaxVolPlasticStrain() = eop;
1063  goto label14;
1064  }
1065 
1066  //
1067  // softening
1068  //
1069  if ( eop > epsult ) {
1070  SCT = 0.;
1071  } else {
1072  SCT = this->give(c2_SCCT, gp) * ( 1. - ( eopr - this->give(c2_EOPP, gp) ) / ( epsult - this->give(c2_EOPP, gp) ) );
1073  }
1074 
1075  // When actual strength is lowered the eff. plastic strain
1076  // increases though no additional strain occurred
1077  status->giveTempMaxVolPlasticStrain() = eop + ( status->giveTempCurrentTensionStrength() - SCT ) /
1078  this->give(c2_E, gp);
1079  status->giveTempCurrentTensionStrength() = SCT;
1080  }
1081  }
1082 
1083  //
1084  // COMPRESSION - EFF. PLASTIC COMPR. STRAIN EEP
1085  // IFUPD=1 INDICATES THAT COMPRESSIVE PLASTIC STRAINS ARE UPDATED
1086  //
1087  //
1088 
1089 label14:
1090 
1091  ifupd = 0;
1092  if ( !( ( this->give(c2_EPP, gp) == 0. ) && ( !this->give(c2_IS_PLASTIC_FLOW, gp) ) ) ) {
1093  ep1 = min(ep.at(1), 0.);
1094  ep2 = min(ep.at(2), 0.);
1095  ep3 = min(ep.at(3), 0.);
1096  dep = ep1 + ep2 + ep3;
1097  //
1098  if ( dep < 0. ) {
1099  // There is a compressive plastic deformation.
1100  // If the plastic
1101  // flow theory is used then DEP is the
1102  // increment of the plastic
1103  // strain.
1104  if ( this->give(c2_EPP, gp) == 0. ) {
1105  // No compression softening.
1106  // For flow theory the plastic strain component
1107  // must be updated
1108  if ( this->give(c2_IS_PLASTIC_FLOW, gp) ) {
1109  ifupd = 1;
1110  return;
1111  }
1112  }
1113 
1114  d = 1.5 * ( ep1 * ep1 + ep2 * ep2 + ep3 * ep3 ) - 0.5 * dep * dep;
1115  dep = sqrt(d);
1116  if ( this->give(c2_IS_PLASTIC_FLOW, gp) ) {
1117  // flow theory
1118  eep = status->giveTempMaxEffPlasticStrain() + dep;
1119  } else {
1120  // deformation theory
1121  eep = dep;
1122  }
1123 
1124  if ( eep > status->giveTempMaxEffPlasticStrain() ) {
1125  //
1126  // current plast.def. EEP is greater than the max. reached
1127  // upto now EPM - strength must be lowered and
1128  // plastic strain components updated if flow rule.
1129  if ( this->give(c2_IS_PLASTIC_FLOW, gp) ) {
1130  ifupd = 1;
1131  }
1132 
1133  //
1134  // SOFTENING?
1135  //
1136  if ( this->give(c2_EPP, gp) != 0. ) {
1137  //
1138  // reduced plastic strain with respect to the elastic
1139  // limit of the virgin concrete EEPR
1140  //
1141  eepr = eep + ( this->give(c2_SCCC, gp) - SCC ) / this->give(c2_E, gp);
1142  if ( eepr <= this->give(c2_EPP, gp) ) {
1143  status->giveTempMaxEffPlasticStrain() = eep;
1144  } else {
1145  // softening
1146  if ( eepr >= this->give(c2_EPU, gp) ) {
1147  SCC = 0.;
1148  } else {
1149  SCC = this->give(c2_SCCC, gp) * ( 1. - ( eepr - this->give(c2_EPP, gp) ) /
1150  ( this->give(c2_EPU, gp) - this->give(c2_EPP, gp) ) );
1151  }
1152 
1153  //
1154  // When actual strength is lowered the eff. plastic strain
1155  // increases though no additional strain occurred
1156  status->giveTempMaxEffPlasticStrain() = eep - ( status->giveTempCurrentPressureStrength() - SCC ) /
1157  this->give(c2_E, gp);
1158  status->giveTempCurrentPressureStrength() = SCC;
1159  }
1160  }
1161  }
1162  }
1163  }
1164 }
1165 
1166 
1167 
1168 
1169 
1170 void
1172 // stirr (double dez, double srf)
1173 //
1174 //
1175 // EVALUATES THE PLASTIC OR VISCOPLSTIC STRESS- STRAIN
1176 // RELATIONS FOR STIRRUPS AND UPDATES DATA IN MAT STATUS AND GP DATA
1177 // ACORDINGLY.
1178 //
1179 // INPUT
1180 // strainincrement->at(1) = DEZ - INCREMENT OF STRAIN
1181 // gp->status->giveCurrentStressInStirrups() - CURENT STRESS, ALSO OUTPUT (INCREMENTED)
1182 //
1183 // this has - MATERIAL CONSTANTS ARRAY:
1184 // STIRR_E - YOUNG MOD.
1185 // STIRR_R - UNIAX. STRENGTH=ELAST. LIMIT
1186 // STIRR_A - STIRRUPS AREA/UNIT LENGHT (BEAM) OR /UNIT AREA (SHELL)
1187 // STIRR_TOL- TOLERANCE OF THE EQUIL. IN THE Z-DIRECTION (=0 NO ITER.)
1188 // STIRR_EREF - EREF - REFRENCE STRAIN RATE FOR PERZYNA'S MATER.
1189 // STIRR_LAMBDA- COEFF. FOR THAT MATER:
1190 // SHTIRR_H- ISOTROPIC HARDENING FACTOR
1191 //
1192 // OUTPUT
1193 // SRF - REAL STIRRUP STRESS IF IFAD=>1, OTHERWISE STRESS
1194 // BEFORE INCR.
1195 // AT PRESENT NO HARDENING AVAILABLE
1196 //
1197 // PLASTIC STRAIN INCREMENT
1198 //
1199 {
1200  Concrete2MaterialStatus *status = static_cast< Concrete2MaterialStatus * >( this->giveStatus(gp) );
1201 
1202  double dep, srf, dez, ovs, s, dt;
1203 
1204  dez = strainIncrement.at(1);
1205  srf = status->giveTempCurrentStressInStirrups();
1206 
1207  dep = 0.;
1208  ovs = fabs(srf) / this->give(stirr_Ft, gp) - 1.;
1209  if ( ( ovs ) > 0. ) {
1210  if ( this->give(stirr_EREF, gp) <= 0. ) {
1211  if ( ( dez * srf ) > 0. ) {
1212  dep = dez;
1213  }
1214  } else {
1215  dt = tStep->giveTimeIncrement();
1216  dep = this->give(stirr_EREF, gp) * exp( ovs / this->give(stirr_LAMBDA, gp) ) * dt;
1217  }
1218  }
1219 
1220  s = srf + this->give(stirr_E, gp) * ( dez - dep );
1221  if ( this->give(c2_IFAD, gp) ) {
1222  srf = s;
1223  }
1224 
1225  status->giveTempCurrentStressInStirrups() = srf;
1226 }
1227 
1228 
1229 
1230 void
1232  MatResponseMode rMode,
1233  GaussPoint *gp,
1234  TimeStep *tStep)
1235 //
1236 // This material is currently unable compute material stiffness
1237 // so it uses slave material (linearElasticMaterial ) to perform this work
1238 // but this material is currently assumed to be linear elastic ->
1239 // plasticity is not taken into account !.
1240 //
1241 {
1242  // error ("givePlateLayerStiffMtrx: unable to compute");
1243  linearElasticMaterial->givePlateLayerStiffMtrx(answer, rMode, gp, tStep);
1244 }
1245 
1246 
1249 /*
1250  * creates new material status corresponding to this class
1251  */
1252 {
1253  Concrete2MaterialStatus *status;
1254 
1255  status = new Concrete2MaterialStatus(1, this->giveDomain(), gp);
1256  return status;
1257 }
1258 
1259 
1260 
1262  StructuralMaterialStatus(n, d, g), plasticStrainVector(), plasticStrainIncrementVector()
1263  //
1264  // constructor
1265  //
1266 {
1267  SCCM = EPM = E0PM = SRF = SEZ = 0.0;
1268  SCTM = -1.0; // init status if SCTM < 0.;
1269 }
1270 
1271 
1272 
1274 //
1275 // DEstructor
1276 //
1277 { }
1278 
1281 //
1282 // saves full information stored in this Status
1283 //
1284 {
1285  contextIOResultType iores;
1286 
1287  if ( ( iores = StructuralMaterialStatus :: saveContext(stream, mode, obj) ) != CIO_OK ) {
1288  THROW_CIOERR(iores);
1289  }
1290 
1291  // write a raw data
1292  if ( !stream.write(SCCM) ) {
1294  }
1295 
1296  if ( !stream.write(EPM) ) {
1298  }
1299 
1300  if ( !stream.write(SCTM) ) {
1302  }
1303 
1304  if ( !stream.write(E0PM) ) {
1306  }
1307 
1308  if ( !stream.write(SRF) ) {
1310  }
1311 
1312  if ( !stream.write(SEZ) ) {
1314  }
1315 
1316  if ( ( iores = plasticStrainVector.storeYourself(stream) ) != CIO_OK ) {
1317  THROW_CIOERR(iores);
1318  }
1319 
1320  // return result back
1321  return CIO_OK;
1322 }
1323 
1324 
1327 //
1328 // restore state variables from stream
1329 //
1330 {
1331  contextIOResultType iores;
1332 
1333  if ( ( iores = StructuralMaterialStatus :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
1334  THROW_CIOERR(iores);
1335  }
1336 
1337  // read raw data
1338  if ( !stream.read(SCCM) ) {
1340  }
1341 
1342  if ( !stream.read(EPM) ) {
1344  }
1345 
1346  if ( !stream.read(SCTM) ) {
1348  }
1349 
1350  if ( !stream.read(E0PM) ) {
1352  }
1353 
1354  if ( !stream.read(SRF) ) {
1356  }
1357 
1358  if ( !stream.read(SEZ) ) {
1360  }
1361 
1362  if ( ( iores = plasticStrainVector.restoreYourself(stream) ) != CIO_OK ) {
1363  THROW_CIOERR(iores);
1364  }
1365 
1366  // return result back
1367  return CIO_OK;
1368 }
1369 
1370 
1371 
1372 void
1374 //
1375 // initializes temp variables according to variables form previous equlibrium state.
1376 //
1377 {
1379 
1380  tempSCCM = SCCM;
1381  tempEPM = EPM;
1382  tempSCTM = SCTM;
1383  tempE0PM = E0PM;
1384  tempSRF = SRF;
1385  tempSEZ = SEZ;
1386 
1387  if ( plasticStrainVector.giveSize() == 0 ) {
1390  }
1391 
1392  if ( plasticStrainIncrementVector.giveSize() == 0 ) {
1395  } else {
1397  }
1398 }
1399 
1400 
1401 
1402 void
1404 //
1405 // updates variables (nonTemp variables describing situation at previous equilibrium state)
1406 // after a new equilibrium state has been reached
1407 // temporary variables are having values corresponding to newly reched equilibrium.
1408 //
1409 {
1411 
1412  SCCM = tempSCCM;
1413  EPM = tempEPM;
1414  SCTM = tempSCTM;
1415  E0PM = tempE0PM;
1416  SRF = tempSRF;
1417  SEZ = tempSEZ;
1418 
1420 }
1421 } // end namespace oofem
double SCTM
Current tension strength.
Definition: concrete2.h:100
static int giveSizeOfVoigtSymVector(MaterialMode mmode)
Returns the size of symmetric part of a reduced stress/strain vector according to given mode...
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Definition: structuralms.C:96
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition: gausspoint.h:191
double SHEARTOL
Threshold value of the relative shear deformation (psi^2/eef) at which shear is considered in layers...
Definition: concrete2.h:167
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
int IFAD
Determines if state variables should be updated or not (>0 updates).
Definition: concrete2.h:175
double stirrLAMBDA
Definition: concrete2.h:171
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
Definition: concrete2.C:1373
void letTempStrainVectorBe(const FloatArray &v)
Assigns tempStrainVector to given vector v.
Definition: structuralms.h:137
Concrete2MaterialStatus(int n, Domain *d, GaussPoint *g)
Definition: concrete2.C:1261
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Returns material status of receiver in given integration point.
Definition: material.C:244
void updateStirrups(GaussPoint *gp, FloatArray &strainIncrement, TimeStep *tStep)
Definition: concrete2.C:1171
GaussPoint * gp
Associated integration point.
Class and object Domain.
Definition: domain.h:115
double stirrEREF
Definition: concrete2.h:171
For computing principal strains from engineering strains.
#define stirr_E
Definition: concrete2.h:78
void dtp2(GaussPoint *gp, FloatArray &e, FloatArray &s, FloatArray &ep, double SCC, double SCT, int *ifplas)
Definition: concrete2.C:857
double & at(int aKey)
Returns the value of the pair which key is aKey.
Definition: dictionary.C:106
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 & at(int i)
Coefficient access function.
Definition: floatarray.h:131
#define stirr_LAMBDA
Definition: concrete2.h:83
int max(int i, int j)
Returns bigger value form two given decimals.
Definition: mathfem.h:71
double & giveTempMaxVolPlasticStrain()
Definition: concrete2.h:124
This class implements a structural material status information.
Definition: structuralms.h:65
double SEZ
Current strain in transverse (z) direction.
Definition: concrete2.h:103
bool includes(int aKey)
Checks if dictionary includes given key.
Definition: dictionary.C:126
virtual ~Concrete2()
Definition: concrete2.C:54
#define c2_EOPU
Definition: concrete2.h:74
const FloatArray & givePlasticStrainIncrementVector() const
Definition: concrete2.h:112
virtual void givePlateLayerStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d plate layer stiffness matrix of receiver.
double SCCM
Current pressure strength.
Definition: concrete2.h:98
Dictionary propertyDictionary
Property dictionary.
Definition: material.h:105
#define _IFT_Concrete2_eopp
Definition: concrete2.h:55
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
void strsoft(GaussPoint *gp, double epsult, FloatArray &ep, double &ep1, double &ep2, double &ep3, double SCC, double SCT, int &ifupd)
Definition: concrete2.C:1005
#define _IFT_Concrete2_stirr_tol
Definition: concrete2.h:63
#define _IFT_Concrete2_stirr_a
Definition: concrete2.h:62
General IO error.
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...
double SCCT
Tension strength.
Definition: concrete2.h:156
double SCCC
Pressure strength.
Definition: concrete2.h:155
static void computePrincipalValDir(FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode)
Computes principal values and directions of stress or strain vector.
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.
#define c2_EOPP
Definition: concrete2.h:73
virtual MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
Definition: concrete2.C:1248
#define c2_SCCC
Definition: concrete2.h:69
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
#define _IFT_Concrete2_stirr_ft
Definition: concrete2.h:61
#define _IFT_Concrete2_stirr_e
Definition: concrete2.h:60
static void giveFullSymVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form.
double giveTimeIncrement()
Returns solution step associated time increment.
Definition: timestep.h:150
#define c2_SHEARTOL
Definition: concrete2.h:75
#define _IFT_Concrete2_is_plastic_flow
Definition: concrete2.h:58
#define _IFT_Concrete2_scct
Definition: concrete2.h:52
virtual int write(const int *data, int count)=0
Writes count integer values from array pointed by data.
double EPP
Threshold eff. plastic strain for softening in compress.
Definition: concrete2.h:157
#define _IFT_Concrete2_epp
Definition: concrete2.h:53
#define c2_n
Definition: concrete2.h:77
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of material property &#39;aProperty&#39;.
Definition: material.C:52
#define stirr_A
Definition: concrete2.h:80
double & giveTempCurrentTensionStrength()
Definition: concrete2.h:121
LinearElasticMaterial * linearElasticMaterial
Definition: concrete2.h:177
static void giveReducedSymVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form.
#define stirr_TOL
Definition: concrete2.h:81
This class implements an abstract base material, which behaves according to deformation theory...
This class implements an isotropic linear elastic material in a finite element problem.
int iperm(int val, int rank)
Returns iperm of val, in specific rank.
Definition: mathfem.C:260
double EPU
Ultimate eff. pl. strain.
Definition: concrete2.h:158
#define OOFEM_ERROR(...)
Definition: error.h:61
#define _IFT_Concrete2_n
Definition: concrete2.h:50
FloatArray plasticStrainIncrementVector
Definition: concrete2.h:94
#define stirr_Ft
Definition: concrete2.h:79
#define _IFT_Concrete2_stirr_lambda
Definition: concrete2.h:65
double & giveTempMaxEffPlasticStrain()
Definition: concrete2.h:120
contextIOResultType restoreYourself(DataStream &stream)
Definition: floatarray.C:895
#define c2_EPU
Definition: concrete2.h:72
double & giveTempCurrentPressureStrength()
Definition: concrete2.h:119
#define c2_E
Definition: concrete2.h:76
Concrete2(int n, Domain *d)
Definition: concrete2.C:48
#define _IFT_Concrete2_ifad
Definition: concrete2.h:59
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Definition: concrete2.C:1403
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
double E0PM
Max. vol. plastic strain.
Definition: concrete2.h:101
#define c2_IFAD
Definition: concrete2.h:85
virtual void giveRealStressVector_PlateLayer(FloatArray &answer, GaussPoint *gp, const FloatArray &, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Definition: concrete2.C:164
#define _IFT_Concrete2_sccc
Definition: concrete2.h:51
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: structuralms.C:157
virtual void givePlateLayerStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d plate layer stiffness matrix of receiver.
Definition: concrete2.C:1231
#define c2_IS_PLASTIC_FLOW
Definition: concrete2.h:84
#define c2_SCCT
Definition: concrete2.h:70
const FloatArray & givePlasticStrainVector() const
Definition: concrete2.h:111
Abstract base class representing a material status information.
Definition: matstatus.h:84
Class representing vector of real numbers.
Definition: floatarray.h:82
#define _IFT_Concrete2_eopu
Definition: concrete2.h:56
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
void letPlasticStrainIncrementVectorBe(FloatArray v)
Definition: concrete2.h:116
double SRF
current stress in stirrups.
Definition: concrete2.h:102
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: concrete2.C:1326
double & giveTempCurrentStrainInZDir()
Definition: concrete2.h:123
double EOPU
Ultimate vol. pl. strain.
Definition: concrete2.h:160
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
#define _IFT_Concrete2_sheartol
Definition: concrete2.h:57
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
#define _IFT_Concrete2_epu
Definition: concrete2.h:54
#define stirr_EREF
Definition: concrete2.h:82
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
double & giveCurrentPressureStrength()
Definition: concrete2.h:127
void dtp3(GaussPoint *gp, FloatArray &e, FloatArray &s, FloatArray &ep, double SCC, double SCT, int *ifplas)
Definition: concrete2.C:654
double & giveTempCurrentStressInStirrups()
Definition: concrete2.h:122
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
#define _IFT_Concrete2_e
Definition: concrete2.h:49
double EOPP
Threshold volumetric plastic strain for soft. in tension.
Definition: concrete2.h:159
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
int min(int i, int j)
Returns smaller value from two given decimals.
Definition: mathfem.h:59
REGISTER_Material(DummyMaterial)
This class implements associated Material Status to Concrete2Material.
Definition: concrete2.h:90
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
double & giveCurrentTensionStrength()
Definition: concrete2.h:129
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.
#define IR_GIVE_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:69
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: concrete2.C:1280
virtual void initTempStatus(GaussPoint *gp)
Initializes temporary variables stored in integration point status at the beginning of new time step...
Definition: material.C:267
#define _IFT_Concrete2_stirr_eref
Definition: concrete2.h:64
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: concrete2.C:60
virtual double give(int, GaussPoint *gp)
Returns the value of material property &#39;aProperty&#39;.
Definition: concrete2.C:93
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: material.C:89
double EPM
Max. eff. plastic strain.
Definition: concrete2.h:99
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
int IS_PLASTIC_FLOW
Indicates that plastic flow (not deformation theory) is used in pressure.
Definition: concrete2.h:173
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
#define c2_EPP
Definition: concrete2.h:71

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