OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
rcm2.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 "rcm2.h"
36 #include "gausspoint.h"
37 #include "floatmatrix.h"
38 #include "floatarray.h"
39 #include "contextioerr.h"
40 #include "mathfem.h"
41 
42 #include <cstring>
43 
44 namespace oofem {
45 #define rcm_STRESSRELERROR 1.e-5
46 #define rcm_RESIDUALSTIFFFACTOR 1.e-3
47 
49  //
50  // constructor
51  //
52 {
53  linearElasticMaterial = NULL;
54  Gf = Ft = 0.;
55 }
56 
57 
59 //
60 // destructor
61 //
62 {
63  //delete linearElasticMaterial;
64 }
65 
66 int
68 //
69 // returns whether receiver supports given mode
70 //
71 {
72  return mode == _3dMat || mode == _PlaneStress ||
73  mode == _PlaneStrain || mode == _1dMat ||
74  mode == _PlateLayer || mode == _2dBeamLayer;
75 }
76 
77 
78 /*
79  * void
80  * RCM2Material :: computeTrialStressIncrement (FloatArray& answer, GaussPoint *gp,
81  * const FloatArray &strainIncrement,
82  * TimeStep* tStep)
83  * //
84  * // returns trial stress according to strainIncrement
85  * //
86  * {
87  * FloatMatrix def;
88  *
89  * this -> giveEffectiveMaterialStiffnessMatrix (def, TangentStiffness,gp,tStep);
90  * answer.beProductOf (def, strainIncrement);
91  * }
92  */
93 
94 
95 void
97  MatResponseMode mode,
98  GaussPoint *gp,
99  TimeStep *tStep)
100 //
101 // computes full constitutive matrix for case of gp stress-strain state.
102 //
103 {
104  /*
105  * FloatMatrix *def;
106  * RCM2MaterialStatus *status =
107  * (RCM2MaterialStatus*) this -> giveStatus (gp);
108  *
109  * CrossSection *crossSection = gp -> giveElement()->giveCrossSection();
110  */
111  this->giveEffectiveMaterialStiffnessMatrix(answer, mode, gp, tStep);
112  // def is full matrix for current gp stress strain mode.
113 }
114 
115 
116 void
118  const FloatArray &totalStrain,
119  TimeStep *tStep)
120 //
121 // returns real stress vector in 3d stress space of receiver according to
122 // previous level of stress and current
123 // strain increment, the only way, how to correctly update gp records
124 //
125 {
126  FloatArray princStress, reducedStressVector, reducedAnswer;
127  FloatArray reducedTotalStrainVector, principalStrain, strainVector;
128  FloatMatrix tempCrackDirs;
129  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
130 
131  this->initTempStatus(gp);
132 
133  // subtract stress independent part
134  // note: eigenStrains (temperature) is not contained in mechanical strain stored in gp
135  // therefore it is necessary to subtract always the total eigen strain value
136  this->giveStressDependentPartOfStrainVector(reducedTotalStrainVector, gp, totalStrain, tStep, VM_Total);
137  //
138 
139  StructuralMaterial :: giveFullSymVectorForm( strainVector, reducedTotalStrainVector, gp->giveMaterialMode() );
140 
141  tempCrackDirs = status->giveTempCrackDirs();
142  this->computePrincipalValDir(principalStrain, tempCrackDirs,
143  strainVector,
145 
146  this->giveRealPrincipalStressVector3d(princStress, gp, principalStrain, tempCrackDirs, tStep);
147  princStress.resizeWithValues(6);
148 
149  tempCrackDirs = status->giveTempCrackDirs();
150  this->transformStressVectorTo(answer, tempCrackDirs, princStress, 1);
151 
152  status->letTempStrainVectorBe(totalStrain);
153  StructuralMaterial :: giveReducedSymVectorForm( reducedStressVector, answer, gp->giveMaterialMode() );
154  status->letTempStressVectorBe(reducedStressVector);
155  this->updateCrackStatus(gp, status->giveCrackStrainVector());
156 
158  answer = reducedAnswer;
159 }
160 
161 
162 void
164  FloatArray &principalStrain,
165  FloatMatrix &tempCrackDirs,
166  TimeStep *tStep)
167 //
168 // returns real principal stress vector in 3d stress space of receiver according to
169 // previous level of stress and current
170 // strain increment, the only way, how to correctly update gp records
171 // updates principal strain and stress of the receiver's status.
172 //
173 {
174  int ind;
175  double maxErr;
176  FloatArray crackStrainVector;
177  FloatArray strainIncrement, crackStrainIterativeIncrement;
178  FloatArray prevPrincipalStrain;
179  FloatArray dSigma;
180  FloatArray elastStrain, sigmaEl, sigmaCr(3);
181  FloatArray fullDSigma;
182  IntArray activatedCracks, crackMapping;
183  FloatMatrix dcr, de, decr, fullDecr, crackDirs;
184  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
185 
186  /*
187  * // totalStressVector = gp -> giveStressVector()->GiveCopy();
188  * reducedTotalStrainVector = status -> giveStrainVector();
189  * reducedTotalStrainVector.add(fullStrainIncrement);
190  * crossSection->giveFullCharacteristicVector(totalStrainVector, gp, reducedTotalStrainVector);
191  * //delete reducedTotalStrainVector;
192  * // plasticStrainVector = status -> givePlasticStrainVector()->GiveCopy();
193  *
194  *
195  * // already cracked - next directions are determined
196  * // according to principal strain directions
197  * status->giveTempCrackDirs(tempCrackDirs);
198  * this->computePrincipalValDir (principalStrain, tempCrackDirs,
199  * totalStrainVector,
200  * principal_strain);
201  * status->letTempCrackDirsBe (tempCrackDirs);
202  */
203  crackStrainVector = status->giveCrackStrainVector(); // local one
204  crackDirs = status->giveCrackDirs();
205  if ( principalStrain.containsOnlyZeroes() ) {
206  // keep old principal values
207  status->letTempCrackDirsBe(crackDirs);
208  } else {
209  this->sortPrincDirAndValCloseTo(& principalStrain,
210  & tempCrackDirs, & crackDirs);
211  status->letTempCrackDirsBe(tempCrackDirs);
212  }
213 
214  // compute de in local system
215  // for iso materials no transformation if stiffness required
216  //
217  // local strain increment
218  prevPrincipalStrain = status->givePrevPrincStrainVector();
219  strainIncrement.beDifferenceOf(principalStrain, prevPrincipalStrain);
220  status->letPrincipalStrainVectorBe(principalStrain);
221 
222  this->giveNormalElasticStiffnessMatrix(de, false, TangentStiffness,
223  gp, tStep, tempCrackDirs);
224  //
225  // construct mapping matrix of active cracks
226  // this mapping will dynamically change as
227  // some crack can unlo or reload
228  //
229  this->updateActiveCrackMap(gp);
230  crackMapping = status->giveCrackMap();
231  // start iteration until stress computed from elastic increment
232  // is equal to stress computed from cracking strain increment
233  // we do this computation in reduced stress strain space
234  dSigma.clear();
235  for ( int iter = 1; iter <= 20; iter++ ) {
236  //
237  // first check if already cracked
238  //
239  if ( status->giveNumberOfTempActiveCracks() ) {
240  // active crack exist
241  this->giveCrackedStiffnessMatrix(dcr, TangentStiffness, gp, tStep);
242  fullDecr = de;
243  fullDecr.add(dcr);
244  decr.resize( crackMapping.maximum(), crackMapping.maximum() );
245  decr.zero();
246  decr.assemble(fullDecr, crackMapping, crackMapping);
247 
248  if ( dSigma.giveSize() == 0 ) {
249  fullDSigma.beProductOf(de, strainIncrement);
250  dSigma.resize( crackMapping.maximum() );
251  dSigma.zero();
252  dSigma.assemble(fullDSigma, crackMapping);
253  //dSigma.beSubArrayOf(fullDSigma, crackMapping); //@todo fix! used old version of method beSubArrayOf
254  }
255 
256  decr.solveForRhs(dSigma, crackStrainIterativeIncrement);
257  for ( int i = 1; i <= 3; i++ ) {
258  if ( ( ind = crackMapping.at(i) ) ) {
259  crackStrainVector.at(i) += crackStrainIterativeIncrement.at(ind);
260  }
261  }
262 
263  // check for crack closing, updates also cracking map
264  this->checkIfClosedCracks(gp, crackStrainVector, crackMapping);
265 
266  // elastic strain component
267  elastStrain.beDifferenceOf(principalStrain, crackStrainVector);
268  sigmaEl.beProductOf(de, elastStrain);
269 
270  // Stress in cracks
271  for ( int i = 1; i <= 3; i++ ) {
272  if ( crackMapping.at(i) ) {
273  sigmaCr.at(i) = giveNormalCrackingStress(gp, crackStrainVector.at(i), i);
274  }
275  }
276 
277  // update status
278  status->letCrackStrainVectorBe(crackStrainVector);
279  } else {
280  //
281  // no active crack exist - elastic behaviour
282  //
283  elastStrain.beDifferenceOf(principalStrain, crackStrainVector);
284  sigmaEl.beProductOf(de, elastStrain);
285  sigmaCr.zero();
286  }
287 
288  // check for new cracks
289  // and update crack map if necessary
290  // when we update map, we need to add new crack at end
291  // because sigmaCr is build
292  this->checkForNewActiveCracks(activatedCracks, gp, crackStrainVector,
293  sigmaEl, sigmaCr, principalStrain);
294  if ( activatedCracks.giveSize() ) {
295  // update crack map also
296  this->updateActiveCrackMap(gp, & activatedCracks);
297  crackMapping = status->giveCrackMap(); // update crackMap
298  }
299 
300  //
301 
302  // compute unbalanced stress
303  // dSigma = sigmaEl - sigmaCr for active cracks
304  fullDSigma = sigmaEl;
305  fullDSigma.subtract(sigmaCr);
306  //fullDSigma.printYourself();
307  //crackMapping.printYourself();
308 
309  dSigma.resize( crackMapping.maximum() );
310  dSigma.zero();
311  dSigma.assemble(fullDSigma, crackMapping);
312  //dSigma.beSubArrayOf(fullDSigma, crackMapping); //@todo fix! used old version of method beSubArrayOf
313  if ( dSigma.giveSize() ) {
314  // fullDSigma.printYourself();
315  //crackMapping.printYourself();
316  //dSigma.printYourself();
317  }
318  //dSigma.printYourself();
319 
320  // find max error in dSigma
321  // if max err < allovedErr -> stop iteration
322  // allowed Err is computed relative to Ft;
323 
324  // check only for active cracks
325  maxErr = 0.;
326  for ( int i = 1; i <= dSigma.giveSize(); i++ ) {
327  if ( fabs( dSigma.at(i) ) > maxErr ) {
328  maxErr = fabs( dSigma.at(i) );
329  }
330  }
331 
332  if ( maxErr < rcm_STRESSRELERROR * this->give(pscm_Ft, gp) ) {
333  status->letPrincipalStressVectorBe(sigmaEl);
334  answer = sigmaEl;
335  return;
336  }
337  } // loop
338 
339  // convergence not reached
340  OOFEM_ERROR("convergence not reached");
341 }
342 
343 
344 void
346  const FloatArray &crackStrain,
347  const FloatArray &princStressVector,
348  FloatArray &crackStressVector,
349  const FloatArray &princStrainVector)
350 //
351 // returns int_array flag showing if some crack
352 // is newly activated or
353 // closed crack is reopened
354 // return 0 if no crack is activated or reactivated.
355 // modifies crackStressVector for newly activated crack.
356 //
357 {
358  double initStress, Le = 0.0;
359  int upd, activationFlag = 0;
360  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
361  FloatArray localStress;
362  FloatMatrix tempCrackDirs;
363  const IntArray &crackMap = status->giveCrackMap();
364 
365  answer.resize(3);
366  answer.zero();
367  localStress = princStressVector;
368  //
369  // local stress is updated according to reached local crack strain
370  //
371  for ( int i = 1; i <= 3; i++ ) { // loop over each possible crack plane
372  // test previous status of each possible crack plane
373  upd = 0;
374  IntArray indx;
376  if ( crackMap.at(i) == 0 && indx.contains(i) ) {
377  if ( status->giveTempMaxCrackStrain(i) > 0. ) {
378  // if (status->giveTempCrackStatus()->at(i) != pscm_NONE) {
379  //
380  // previously cracked direction
381  //
382  initStress = 0.;
383  } else {
384  // newer cracked, so we compute principal stresses
385  // and compare them with reduced tension strength
386  // as a criterion for crack initiation
387 
388  FloatArray crackPlaneNormal(3);
389  tempCrackDirs = status->giveTempCrackDirs();
390  for ( int j = 1; j <= 3; j++ ) {
391  crackPlaneNormal.at(j) = tempCrackDirs.at(j, i);
392  }
393 
394  Le = this->giveCharacteristicElementLength(gp, crackPlaneNormal);
395  initStress = this->computeStrength(gp, Le);
396  upd = 1;
397  }
398 
399  if ( localStress.at(i) > initStress ) {
400  crackStressVector.at(i) = initStress;
401  answer.at(i) = 1;
402  activationFlag = 1;
403  if ( upd ) {
404  this->updateStatusForNewCrack(gp, i, Le);
405  }
406  }
407  } // end of tested crack
408  } // end of loop over are possible directions
409 
410  if ( activationFlag ) {
411  return;
412  }
413 
414  answer.clear();
415 }
416 
417 
418 void
420 //
421 // updates gp status when new crack-plane i is formed with charLength Le
422 // updates Le and computes and sets minEffStrainForFullyOpenCrack
423 //
424 {
425  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
426 
427  if ( Le <= 0 ) {
428  OOFEM_ERROR("Element %d returned zero char length", gp->giveElement()->giveNumber() );
429  }
430 
431  status->setCharLength(i, Le);
432  //status -> setMinCrackStrainsForFullyOpenCrack (i, this->giveMinCrackStrainsForFullyOpenCrack(gp,i));
433 }
434 
435 double
437 {
438  return gp->giveElement()->giveCharacteristicLength(crackPlaneNormal);
439 }
440 
441 void
443 //
444 // updates gp records and MatStatus due to cracking.
445 //
446 // Updates MatStatus (its temporary variables) to current
447 // reached status during integrating incremental constitutive relations
448 // Temporary variables are used, because we may integrate constitutive
449 // realtions many times for different strainIncrement in order to
450 // reach equilibrium state, so we don't want to change other variables
451 // which describing previously reached equlibrium. After a new equilibrium
452 // is reached, this->updateYourself is called which invokes matStatus->
453 // updateYourself(), which copies temporary variables to variables
454 // describing equilibrium.
455 //
456 //
457 {
458  double minCrackStrainsForFullyOpenCrack;
459  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
460  const IntArray &crackMap = status->giveCrackMap();
461 
462  // check if material previously cracked
463  // and compute possible crack planes
464  // or if newer cracked, so we compute principal stresses
465  // (for iso mat. coincide with princ strains)
466  // and compare them with reduced tension strength
467  // as a criterion for crack initiation
468 
469  // principalStrain = status->givePrincipalStrainVector();
470  // transform stresses to principal directions of strains
471  //
472  // local stress is updated according to reached local crack strain
473  //
474  IntArray indx;
476  for ( int i = 1; i <= 3; i++ ) { // loop over each possible crack plane
477  if ( crackMap.at(i) != 0 && indx.contains(i) ) {
478  if ( status->giveTempMaxCrackStrain(i) < crackStrain.at(i) ) {
479  status->setTempMaxCrackStrain( i, crackStrain.at(i) );
480  }
481 
482  minCrackStrainsForFullyOpenCrack = this->giveMinCrackStrainsForFullyOpenCrack(gp, i);
483  //if ((crackStrain.at(i) >= status->giveMinCrackStrainsForFullyOpenCrack(i)) &&
484  if ( ( crackStrain.at(i) >= minCrackStrainsForFullyOpenCrack ) &&
485  ( crackStrain.at(i) >= status->giveTempMaxCrackStrain(i) ) ) {
486  //
487  // fully open crack
488  //
489  status->setTempCrackStatus(i, pscm_OPEN);
490  } else if ( crackStrain.at(i) >= status->giveTempMaxCrackStrain(i) ) {
491  //
492  // further softening of crack
493  //
495  // status->giveTempReachedSofteningStress()->at(i) = localStress->at(i);
496  } else if ( crackStrain.at(i) <= 0. ) {
497  if ( status->giveTempCrackStatus(i) != pscm_NONE ) {
498  // previously active crack becomes closed
499  status->setTempCrackStatus(i, pscm_CLOSED);
500  }
501  } else {
502  //
503  // crack unloading or reloading
504  //
506  }
507  } // end for possible direction
508  } // end loop over prin directions
509 }
510 
511 
512 void
514  IntArray &crackMap)
515 //
516 // Check if crack closing occurs
517 // if yes updates crackStrainVector and gp-status accordingly
518 //
519 {
520  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
521  int isClosing = 0;
522 
523  for ( int i = 1; i <= 3; i++ ) {
524  if ( crackMap.at(i) ) {
525  if ( crackStrainVector.at(i) < 0. ) {
526  // crack closing occur
527  crackStrainVector.at(i) = 0.;
528  crackMap.at(i) = 0;
529  // status->giveTempCrackStatus()->at(i) = pscm_CLOSED;
530  isClosing = 1;
531  }
532  }
533  }
534 
535  status->letCrackMapBe(crackMap);
536  if ( isClosing ) {
537  this->updateActiveCrackMap(gp);
538  }
539 
540  crackMap = status->giveCrackMap(); // update crack Map
541 }
542 
543 
544 void
546  bool reduce, MatResponseMode rMode,
547  GaussPoint *gp, TimeStep *tStep,
548  const FloatMatrix &dir)
549 //
550 // return Elastic Stiffness matrix for normal Stresses
551 // taking into account the directions of normal stresses
552 // (not supported now)
553 //
554 {
556  FloatMatrix de, fullAnswer(3, 3);
557  IntArray mask;
558  int sd;
559 
560  FloatMatrix stiff;
561  lMat->giveStiffnessMatrix(stiff, rMode, gp, tStep);
562  this->giveFullSymMatrixForm( de, stiff, gp->giveMaterialMode() );
563 
564  // copy first 3x3 submatrix to answer
565  for ( int i = 1; i <= 3; i++ ) {
566  for ( int j = 1; j <= 3; j++ ) {
567  fullAnswer.at(i, j) = de.at(i, j);
568  }
569  }
570 
571  if ( !reduce ) { // 3x3 full form required
572  answer = fullAnswer;
573  } else {
574  // reduced form for only
575  // nonzero normal stresses
576  //
577  // first find spatial dimension of problem
578  sd = 0;
579 
580  IntArray indx;
582  for ( int i = 1; i <= 3; i++ ) {
583  if ( indx.contains(i) ) {
584  sd++;
585  }
586  }
587 
588  answer.resize(sd, sd);
589  answer.zero();
590 
591  // copy fullAnswer to reduced one
593  for ( int i = 1; i <= sd; i++ ) {
594  int iidx = mask.findFirstIndexOf(i);
595  for ( int j = 1; j <= sd; j++ ) {
596  int jidx = mask.findFirstIndexOf(j);
597  if ( iidx && jidx ) {
598  answer.at(i, j) = fullAnswer.at(iidx, jidx);
599  }
600  }
601  }
602  }
603 }
604 
605 
606 void
608  MatResponseMode rMode, GaussPoint *gp,
609  TimeStep *tStep)
610 //
611 // returns effective material stiffness matrix in full form
612 // for gp stress strain mode
613 //
614 {
615  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
616  StructuralMaterial *lMat = static_cast< StructuralMaterial * >( this->giveLinearElasticMaterial() );
617  int numberOfActiveCracks = status->giveNumberOfTempActiveCracks();
618  int indi, indj, ii, jj;
619  double G, princStressDis, princStrainDis;
620  FloatMatrix de, invDe, compliance, dcr, d, df, t, tempCrackDirs;
621  FloatArray principalStressVector, principalStrainVector;
622  IntArray mask;
623 
624  if ( ( rMode == ElasticStiffness ) || ( numberOfActiveCracks == 0 ) ) {
625  lMat->giveStiffnessMatrix(answer, rMode, gp, tStep);
626  return;
627  }
628 
629  // this->updateActiveCrackMap(gp) must be done after restart.
630  this->updateActiveCrackMap(gp);
631  tempCrackDirs = status->giveTempCrackDirs();
632  this->giveNormalElasticStiffnessMatrix(de, true, rMode, gp, tStep,
633  tempCrackDirs);
634  invDe.beInverseOf(de);
635  this->giveCrackedStiffnessMatrix(dcr, rMode, gp, tStep);
637  compliance.resize( mask.giveSize(), mask.giveSize() );
638 
639  // we will set
640  // first we set compliances for normal streses in
641  // local coordinate system defined by crackplane
642  for ( int i = 1; i <= 3; i++ ) {
643  IntArray indx;
645  if ( ( indi = indx.findFirstIndexOf(i) ) ) {
646  for ( int j = 1; j <= 3; j++ ) {
647  if ( ( indj = indx.findFirstIndexOf(j) ) ) {
648  compliance.at(indi, indj) += invDe.at(i, j);
649  }
650  }
651 
652  if ( status->isCrackActive(i) ) {
653  if ( dcr.at(i, i) <= 1.e-8 ) {
654  compliance.at(indi, indi) *= rcm2_BIGNUMBER;
655  } else {
656  compliance.at(indi, indi) += 1. / dcr.at(i, i);
657  }
658  }
659  }
660  }
661 
662  principalStressVector = status->getPrincipalStressVector();
663  principalStrainVector = status->getPrincipalStrainVector();
664 
665  // now remain to set shears
666  G = this->give(pscm_G, gp);
667 
668  IntArray indx;
670  for ( int i = 4; i <= 6; i++ ) {
671  if ( ( indi = indx.findFirstIndexOf(i) ) ) {
672  if ( i == 4 ) {
673  ii = 2;
674  jj = 3;
675  } else if ( i == 5 ) {
676  ii = 1;
677  jj = 3;
678  } else {
679  ii = 1;
680  jj = 2;
681  }
682 
683  princStressDis = principalStressVector.at(ii) -
684  principalStressVector.at(jj);
685  princStrainDis = principalStrainVector.at(ii) -
686  principalStrainVector.at(jj);
687  if ( fabs(princStrainDis) < rcm_SMALL_STRAIN ) {
688  compliance.at(indi, indi) = 1. / G;
689  } else if ( fabs(princStressDis) < 1.e-8 ) {
690  compliance.at(indi, indi) = rcm2_BIGNUMBER;
691  } else {
692  compliance.at(indi, indi) = 2 * princStrainDis / princStressDis;
693  }
694  }
695  }
696 
697  // now we invert compliance to get stiffness in reduced space
698  d.beInverseOf(compliance);
699  // delete compliance;
700  //
701  // now let d to grow to Full Format
702  //
704  df.resize(6, 6);
705  df.zero();
706  df.assemble(d, mask);
707  //
708  // final step - transform stiffnes to global c.s
709  //
710 
711  this->giveStressVectorTranformationMtrx(t, tempCrackDirs, 1);
712  df.rotatedWith(t, 't');
713 
715 }
716 
717 
718 void
720  MatResponseMode rMode,
721  GaussPoint *gp,
722  TimeStep *tStep)
723 //
724 //
725 // Returns material incremental stiffness matrix for cracked concrete.
726 // This matrix is composed from submatrix, each of them corresponds to
727 // one active crack in material point.
728 // when constructing submatrix, following assumptions are made:
729 //
730 // A salient characteriastic of crack formation concerns the fact that in most general case
731 // of tree-dimensional solid only 3 out of 6 components of the crack strain rate vector
732 // are possibly non-zero.(the normal and two shear strain rates).
733 // We therefore assume that the stress-strain law for the crack has a structure
734 // such that the other strains rate components vanish. Moreover we assume that
735 // the novanishing strain rate components are only related to corresponding
736 // stress rate components (submatrix has dimensions 3x3).
737 //
738 // if strainIncrement is defined (not null) then we take care about possible unlo&reloading
739 // we don't teke care about possible cracking (or non-linear softening) during strain increment
740 // this correction is made in this -> updateCrackStatus (gp);
741 {
742  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
743  int numberOfActiveCracks = status->giveNumberOfTempActiveCracks();
744  const IntArray &crackMap = status->giveCrackMap();
745  if ( numberOfActiveCracks == 0 ) {
746  answer.clear();
747  return;
748  }
749 
750  answer.resize(3, 3);
751  answer.zero();
752 
753  // loop over each active crack plane
754  for ( int i = 1; i <= 3; i++ ) {
755  if ( crackMap.at(i) ) {
756  // obtain incremental law for one crack
757  answer.at(i, i) = this->giveCrackingModulus(rMode, gp,
758  status->giveCrackStrain(i),
759  i);
760  }
761  }
762 }
763 
764 
765 void
767 //
768 //
769 // updates mapping matrix of active cracks
770 //
771 {
772  int indx = 1;
773  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
774  IntArray crackMap = status->giveCrackMap();
775 
776  //if (crackMap == NULL) OOFEM_ERROR("NULL pointer encountered");
777 
778  for ( int i = 1; i <= 3; i++ ) {
779  if ( status->isCrackActive(i) ) {
780  crackMap.at(i) = indx++;
781  } else if ( activatedCracks ) {
782  if ( ( activatedCracks->at(i) != 0 ) ) {
783  crackMap.at(i) = indx++;
784  }
785  } else {
786  crackMap.at(i) = 0;
787  }
788  }
789 
790  // store modified map into status
791  status->letCrackMapBe(crackMap);
792 }
793 
794 
797 {
798  IRResultType result; // Required by IR_GIVE_FIELD macro
799 
802 
803  return this->giveLinearElasticMaterial()->initializeFrom(ir);
804 }
805 
806 double
807 RCM2Material :: give(int aProperty, GaussPoint *gp)
808 // Returns the value of the property aProperty (e.g. the Young's modulus
809 // 'E') of the receiver.
810 {
811  double value = 0.0;
812 
813  if ( aProperty == pscm_Gf ) {
814  return this->Gf;
815  }
816 
817  // if (aProperty == pscm_Beta) return this->beta;
818  if ( aProperty == pscm_Ft ) {
819  return this->Ft;
820  }
821 
822  if ( aProperty == pscm_Ee ) {
823  aProperty = 'E';
824  }
825 
826  if ( aProperty == pscm_G ) {
827  aProperty = 'G';
828  }
829 
830  if ( propertyDictionary.includes(aProperty) ) {
831  value = propertyDictionary.at(aProperty);
832  } else {
833  if ( linearElasticMaterial ) {
834  value = this->linearElasticMaterial->give(aProperty, gp);
835  } else {
836  OOFEM_ERROR("property not defined");
837  }
838  }
839 
840  return value;
841 }
842 
843 
844 int
846 {
847  RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
848  if ( type == IST_CrackedFlag ) {
849  answer.resize(1);
850  answer.at(1) = status->giveAlreadyCrack();
851  return 1;
852  } else if ( type == IST_CrackDirs ) {
853  const FloatMatrix &help = status->giveCrackDirs();
854  answer.resize(9);
855  for ( int i = 1; i <= 3; i++ ) {
856  answer.at(i) = help.at(1, i);
857  answer.at(i + 3) = help.at(2, i);
858  answer.at(i + 6) = help.at(3, i);
859  }
860 
861  return 1;
862  } else if ( type == IST_CrackStatuses ) {
863  const IntArray &crackStatus = status->giveCrackStatus();
864  answer.resize(3);
865  for ( int i = 1; i <= 3; i++ ) {
866  answer.at(i) = crackStatus.at(i);
867  }
868 
869  return 1;
870  } else {
871  return StructuralMaterial :: giveIPValue(answer, gp, type, tStep);
872  }
873 }
874 
875 void
877  MatResponseMode mode,
878  GaussPoint *gp,
879  TimeStep *tStep)
880 {
881  //
882  // returns receiver 3d material matrix
883  //
884  this->giveMaterialStiffnessMatrix(answer, mode, gp, tStep);
885 }
886 
887 
888 void
890  MatResponseMode mode,
891  GaussPoint *gp,
892  TimeStep *tStep)
893 
894 //
895 // returns receiver's 2dPlaneStressMtrx
896 // (2dPlaneStres ==> sigma_z = tau_xz = tau_yz = 0.)
897 //
898 // standard method from Material Class overloaded, because no inversion is needed.
899 // the reduction from 3d case will not work
900 // this implementation should be faster.
901 {
902  this->giveMaterialStiffnessMatrix(answer, mode, gp, tStep);
903 }
904 
905 
906 void
908  MatResponseMode mode,
909  GaussPoint *gp,
910  TimeStep *tStep)
911 
912 //
913 // return receiver's 2dPlaneStrainMtrx constructed from
914 // general 3dMatrialStiffnessMatrix
915 // (2dPlaneStrain ==> eps_z = gamma_xz = gamma_yz = 0.)
916 //
917 {
918  this->giveMaterialStiffnessMatrix(answer, mode, gp, tStep);
919 }
920 
921 
922 void
924  MatResponseMode mode,
925  GaussPoint *gp,
926  TimeStep *tStep)
927 
928 //
929 // returns receiver's 1dMaterialStiffnessMAtrix
930 // (1d case ==> sigma_y = sigma_z = tau_yz = tau_zx = tau_xy = 0.)
931 {
932  this->giveMaterialStiffnessMatrix(answer, mode, gp, tStep);
933 }
934 
935 
936 void
938  MatResponseMode mode,
939  GaussPoint *gp,
940  TimeStep *tStep)
941 //
942 // returns receiver's 2dBeamLayerStiffMtrx.
943 // (2dPlaneStres ==> sigma_z = tau_xz = tau_yz = 0.)
944 //
945 // standard method from Material Class overloaded, because no inversion is needed.
946 // the reduction from 3d case will not work
947 // this implementation should be faster.
948 {
949  this->giveMaterialStiffnessMatrix(answer, mode, gp, tStep);
950 }
951 
952 
953 void
955  MatResponseMode mode,
956  GaussPoint *gp,
957  TimeStep *tStep)
958 //
959 // returns receiver's 2dPlateLayerMtrx
960 // (2dPlaneStres ==> sigma_z = tau_xz = tau_yz = 0.)
961 //
962 // standard method from Material Class overloaded, because no inversion is needed.
963 // the reduction from 3d case will not work
964 // this implementation should be faster.
965 {
966  this->giveMaterialStiffnessMatrix(answer, mode, gp, tStep);
967 }
968 
969 
970 
971 
973  StructuralMaterialStatus(n, d, g), crackStatuses(3), tempCrackStatuses(3),
974  maxCrackStrains(3), tempMaxCrackStrains(3), crackStrainVector(3),
975  oldCrackStrainVector(3), crackDirs(3, 3), tempCrackDirs(3, 3), charLengths(3),
976  //minEffStrainsForFullyOpenCrack(3),
977  principalStrain(3), oldPrincipalStrain(3),
978  principalStress(3), oldPrincipalStress(3), crackMap(3)
979 {
980  for ( int i = 1; i <= 3; i++ ) {
981  crackDirs.at(i, i) = tempCrackDirs.at(i, i) = 1.0;
982  }
983 }
984 
985 
987 { }
988 
989 
990 
991 int
993 //
994 // return nonzero if crack is active
995 //
996 {
997  if ( crackStrainVector.at(i) > 0. ) {
998  return 1;
999  }
1000 
1001  // handle situation for new crack with crackStrainVector->at(i) = 0.
1002  if ( crackMap.at(i) ) {
1003  return 1;
1004  }
1005 
1006  return 0;
1007 }
1008 
1009 void
1011 {
1012  int i;
1013  char s [ 11 ];
1014 
1016  fprintf(file, "status { ");
1017  if ( this->giveTempAlreadyCrack() ) {
1018  for ( i = 1; i <= 3; i++ ) {
1019  switch ( crackStatuses.at(i) ) {
1020  case pscm_NONE:
1021  strcpy(s, "NONE");
1022  break;
1023  case pscm_OPEN:
1024  strcpy(s, "OPEN");
1025  break;
1026  case pscm_SOFTENING:
1027  strcpy(s, "SOFTENING");
1028  break;
1029  case pscm_RELOADING:
1030  strcpy(s, "RELOADING");
1031  break;
1032  case pscm_UNLOADING:
1033  strcpy(s, "UNLOADING");
1034  break;
1035  default:
1036  strcpy(s, "UNKNOWN");
1037  break;
1038  }
1039 
1040  fprintf( file, "crack %d {status %s, normal to crackplane { %f %f %f }} ",
1041  i, s, crackDirs.at(1, i), crackDirs.at(2, i), crackDirs.at(3, i) );
1042  }
1043  }
1044 
1045  fprintf(file, "}\n");
1046 }
1047 
1048 int
1050 //
1051 // return number of active cracks, computed upon CrackStatus
1052 //
1053 {
1054  int answer = 0;
1055 
1056  for ( int i = 1; i <= 3; i++ ) {
1057  if ( oldCrackStrainVector.at(i) > 0. ) {
1058  answer++;
1059  }
1060  }
1061 
1062  return answer;
1063 }
1064 
1065 
1066 int
1068 //
1069 // return number of active cracks, computed upon tempCrackStatus
1070 //
1071 {
1072  int answer = 0;
1073 
1074  for ( int i = 1; i <= 3; i++ ) {
1075  if ( isCrackActive(i) ) {
1076  answer++;
1077  }
1078  }
1079 
1080  return answer;
1081 }
1082 
1083 
1084 
1085 void
1087 //
1088 // initializes temp variables according to variables form previous equlibrium state.
1089 // builds new crackMap
1090 //
1091 {
1093 
1098 
1101 
1102  for ( int i = 1; i <= 3; i++ ) {
1103  crackMap.at(i) = 0;
1104  }
1105 }
1106 
1107 
1108 
1109 void
1111 //
1112 // updates variables (nonTemp variables describing situation at previous equilibrium state)
1113 // after a new equilibrium state has been reached
1114 // temporary variables are having values corresponding to newly reched equilibrium.
1115 //
1116 {
1118 
1123 
1126 }
1127 
1128 
1129 
1132 //
1133 // saves full information stored in this Status
1134 // no temp variables stored
1135 //
1136 {
1137  contextIOResultType iores;
1138 
1139  // save parent class status
1140  if ( ( iores = StructuralMaterialStatus :: saveContext(stream, mode, obj) ) != CIO_OK ) {
1141  THROW_CIOERR(iores);
1142  }
1143 
1144  // write a raw data
1145 
1146  if ( ( iores = crackStatuses.storeYourself(stream) ) != CIO_OK ) {
1147  THROW_CIOERR(iores);
1148  }
1149 
1150  if ( ( iores = maxCrackStrains.storeYourself(stream) ) != CIO_OK ) {
1151  THROW_CIOERR(iores);
1152  }
1153 
1154  if ( ( iores = crackDirs.storeYourself(stream) ) != CIO_OK ) {
1155  THROW_CIOERR(iores);
1156  }
1157 
1158  //if ((iores = minEffStrainsForFullyOpenCrack.storeYourself(stream)) != CIO_OK) THROW_CIOERR(iores);
1159  if ( ( iores = charLengths.storeYourself(stream) ) != CIO_OK ) {
1160  THROW_CIOERR(iores);
1161  }
1162 
1163  if ( ( iores = crackStrainVector.storeYourself(stream) ) != CIO_OK ) {
1164  THROW_CIOERR(iores);
1165  }
1166 
1167  if ( ( iores = oldCrackStrainVector.storeYourself(stream) ) != CIO_OK ) {
1168  THROW_CIOERR(iores);
1169  }
1170 
1171  if ( ( iores = oldPrincipalStrain.storeYourself(stream) ) != CIO_OK ) {
1172  THROW_CIOERR(iores);
1173  }
1174 
1175  if ( ( iores = principalStrain.storeYourself(stream) ) != CIO_OK ) {
1176  THROW_CIOERR(iores);
1177  }
1178 
1179  if ( ( iores = oldPrincipalStress.storeYourself(stream) ) != CIO_OK ) {
1180  THROW_CIOERR(iores);
1181  }
1182 
1183  if ( ( iores = principalStress.storeYourself(stream) ) != CIO_OK ) {
1184  THROW_CIOERR(iores);
1185  }
1186 
1187  if ( ( iores = crackMap.storeYourself(stream) ) != CIO_OK ) {
1188  THROW_CIOERR(iores);
1189  }
1190 
1191  return CIO_OK;
1192 }
1193 
1196 //
1197 // restores full information stored in stream to this Status
1198 //
1199 {
1200  contextIOResultType iores;
1201 
1202  // read parent class status
1203  if ( ( iores = StructuralMaterialStatus :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
1204  THROW_CIOERR(iores);
1205  }
1206 
1207  // read raw data
1208 
1209  if ( ( iores = crackStatuses.restoreYourself(stream) ) != CIO_OK ) {
1210  THROW_CIOERR(iores);
1211  }
1212 
1213  if ( ( iores = maxCrackStrains.restoreYourself(stream) ) != CIO_OK ) {
1214  THROW_CIOERR(iores);
1215  }
1216 
1217  if ( ( iores = crackDirs.restoreYourself(stream) ) != CIO_OK ) {
1218  THROW_CIOERR(iores);
1219  }
1220 
1221  //if ((iores = minEffStrainsForFullyOpenCrack.restoreYourself(stream)) != CIO_OK) THROW_CIOERR(iores);
1222  if ( ( iores = charLengths.restoreYourself(stream) ) != CIO_OK ) {
1223  THROW_CIOERR(iores);
1224  }
1225 
1226  if ( ( iores = crackStrainVector.restoreYourself(stream) ) != CIO_OK ) {
1227  THROW_CIOERR(iores);
1228  }
1229 
1230  if ( ( iores = oldCrackStrainVector.restoreYourself(stream) ) != CIO_OK ) {
1231  THROW_CIOERR(iores);
1232  }
1233 
1234  if ( ( iores = oldPrincipalStrain.restoreYourself(stream) ) != CIO_OK ) {
1235  THROW_CIOERR(iores);
1236  }
1237 
1238  if ( ( iores = principalStrain.restoreYourself(stream) ) != CIO_OK ) {
1239  THROW_CIOERR(iores);
1240  }
1241 
1242  if ( ( iores = oldPrincipalStress.restoreYourself(stream) ) != CIO_OK ) {
1243  THROW_CIOERR(iores);
1244  }
1245 
1246  if ( ( iores = principalStress.restoreYourself(stream) ) != CIO_OK ) {
1247  THROW_CIOERR(iores);
1248  }
1249 
1250  if ( ( iores = crackMap.restoreYourself(stream) ) != CIO_OK ) {
1251  THROW_CIOERR(iores);
1252  }
1253 
1254  return CIO_OK; // return succes
1255 }
1256 } // end namespace oofem
bool contains(int value) const
Definition: intarray.h:283
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Definition: structuralms.C:96
const FloatArray & getPrincipalStrainVector() const
Definition: rcm2.h:109
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
contextIOResultType storeYourself(DataStream &stream) const
Stores array to output stream.
Definition: intarray.C:289
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
static int giveVoigtSymVectorMask(IntArray &answer, MaterialMode mmode)
The same as giveVoigtVectorMask but returns a mask corresponding to a symmetric second order tensor...
void giveNormalElasticStiffnessMatrix(FloatMatrix &answer, bool reduce, MatResponseMode, GaussPoint *, TimeStep *tStep, const FloatMatrix &)
Definition: rcm2.C:545
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
Class and object Domain.
Definition: domain.h:115
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: rcm2.C:796
For computing principal strains from engineering strains.
virtual double computeStrength(GaussPoint *gp, double)=0
#define _IFT_RCM2Material_gf
Definition: rcm2.h:46
bool solveForRhs(const FloatArray &b, FloatArray &answer, bool transpose=false)
Solves the system of linear equations .
Definition: floatmatrix.C:1112
void letTempCrackDirsBe(FloatMatrix a)
Definition: rcm2.h:127
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: rcm2.C:1195
double & at(int aKey)
Returns the value of the pair which key is aKey.
Definition: dictionary.C:106
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Definition: rcm2.C:1110
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
FloatMatrix tempCrackDirs
Definition: rcm2.h:87
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.
const IntArray & giveCrackStatus()
Definition: rcm2.h:152
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Print receiver&#39;s output to given stream.
Definition: rcm2.C:1010
void zero()
Sets all component to zero.
Definition: intarray.C:52
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
#define pscm_CLOSED
Definition: fcm.h:58
LinearElasticMaterial * linearElasticMaterial
Definition: rcm2.h:181
virtual void giveEffectiveMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Definition: rcm2.C:607
virtual double giveCharacteristicElementLength(GaussPoint *gp, const FloatArray &crackPlaneNormal)
Definition: rcm2.C:436
LinearElasticMaterial * giveLinearElasticMaterial()
Definition: rcm2.h:199
This class implements a structural material status information.
Definition: structuralms.h:65
#define pscm_NONE
Definition: fcm.h:54
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
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.
const IntArray & giveTempCrackStatus()
Definition: rcm2.h:131
#define rcm2_BIGNUMBER
Definition: rcm2.h:72
double giveCrackStrain(int icrack) const
Definition: rcm2.h:136
static void sortPrincDirAndValCloseTo(FloatArray *pVal, FloatMatrix *pDir, FloatMatrix *toPDir)
Method for sorting newly computed principal values (pVal) and corresponding principal directions (pDi...
virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 1d stiffness matrix of receiver.
Definition: rcm2.C:923
#define _IFT_RCM2Material_ft
Definition: rcm2.h:47
FloatArray oldCrackStrainVector
Definition: rcm2.h:85
Dictionary propertyDictionary
Property dictionary.
Definition: material.h:105
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
void setTempCrackStatus(int icrack, int val)
Definition: rcm2.h:133
Element * giveElement()
Returns corresponding element to receiver.
Definition: gausspoint.h:188
void letCrackMapBe(IntArray map)
Definition: rcm2.h:117
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
FloatArray charLengths
Definition: rcm2.h:89
virtual void givePlateLayerStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d plate layer stiffness matrix of receiver.
Definition: rcm2.C:954
const FloatArray & givePrevPrincStrainVector() const
Definition: rcm2.h:111
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Print receiver&#39;s output to given stream.
Definition: structuralms.C:73
virtual void giveRealStressVector(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
Computes the real stress vector for given total strain and integration point.
Definition: rcm2.C:117
static void computePrincipalValDir(FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode)
Computes principal values and directions of stress or strain vector.
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
MatResponseMode
Describes the character of characteristic material matrix.
virtual void givePlaneStressStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane stress stiffness matrix of receiver.
Definition: rcm2.C:889
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
const FloatMatrix & giveCrackDirs()
Definition: rcm2.h:151
contextIOResultType storeYourself(DataStream &stream) const
Definition: floatmatrix.C:1852
int giveAlreadyCrack() const
Definition: rcm2.h:153
void updateActiveCrackMap(GaussPoint *gp, const IntArray *activatedCracks=NULL)
Definition: rcm2.C:766
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.
RCM2MaterialStatus(int n, Domain *d, GaussPoint *g)
Definition: rcm2.C:972
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...
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of material property &#39;aProperty&#39;.
Definition: rcm2.C:807
const FloatArray & getPrincipalStressVector() const
Definition: rcm2.h:110
#define pscm_Ft
Definition: rcm2.h:57
virtual int giveNumberOfActiveCracks() const
Definition: rcm2.C:1049
virtual ~RCM2Material()
Definition: rcm2.C:58
static void giveReducedSymVectorForm(FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form.
void setCharLength(int icrack, double val)
Definition: rcm2.h:148
#define pscm_OPEN
Definition: rcm2.h:61
FloatArray maxCrackStrains
Max crack strain reached.
Definition: rcm2.h:83
#define OOFEM_ERROR(...)
Definition: error.h:61
void giveCrackedStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Definition: rcm2.C:719
void clear()
Clears the array (zero size).
Definition: intarray.h:177
const IntArray & giveCrackMap() const
Definition: rcm2.h:116
FloatArray crackStrainVector
Components of crack strain vector.
Definition: rcm2.h:85
This class implements associated Material Status to SmearedCrackingMaterail.
Definition: rcm2.h:77
virtual int hasMaterialModeCapability(MaterialMode mode)
Tests if material supports material mode.
Definition: rcm2.C:67
virtual void givePlaneStrainStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane strain stiffness matrix of receiver.
Definition: rcm2.C:907
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
int maximum() const
Finds the maximum component in the array.
Definition: intarray.C:195
virtual void give2dBeamLayerStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d beam layer stiffness matrix of receiver.
Definition: rcm2.C:937
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 & giveCrackStrainVector() const
Definition: rcm2.h:135
virtual void giveMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep)
Definition: rcm2.C:96
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: rcm2.C:1131
void resizeWithValues(int s, int allocChunk=0)
Checks size of receiver towards requested bounds.
Definition: floatarray.C:615
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
int giveTempAlreadyCrack() const
Definition: rcm2.h:121
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: structuralms.C:157
#define pscm_Gf
Definition: rcm2.h:54
contextIOResultType restoreYourself(DataStream &stream)
Restores array from image on stream.
Definition: intarray.C:305
virtual ~RCM2MaterialStatus()
Definition: rcm2.C:986
void letCrackStrainVectorBe(FloatArray a)
Definition: rcm2.h:138
contextIOResultType restoreYourself(DataStream &stream)
Definition: floatmatrix.C:1877
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
#define pscm_UNLOADING
Definition: rcm2.h:64
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual void checkIfClosedCracks(GaussPoint *gp, FloatArray &crackStrainVector, IntArray &)
Definition: rcm2.C:513
static void giveStressVectorTranformationMtrx(FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
Computes 3d stress vector transformation matrix from standard vector transformation matrix...
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
Definition: rcm2.C:1086
void letTempStressVectorBe(const FloatArray &v)
Assigns tempStressVector to given vector v.
Definition: structuralms.h:135
virtual void give3dMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes full 3d material stiffness matrix at given integration point, time, respecting load history ...
Definition: rcm2.C:876
void assemble(const FloatArray &fe, const IntArray &loc)
Assembles the array fe (typically, the load vector of a finite element) into the receiver, using loc as location array.
Definition: floatarray.C:551
virtual double giveMinCrackStrainsForFullyOpenCrack(GaussPoint *gp, int i)=0
virtual int giveNumberOfTempActiveCracks() const
Definition: rcm2.C:1067
void giveRealPrincipalStressVector3d(FloatArray &answer, GaussPoint *, FloatArray &, FloatMatrix &, TimeStep *)
Definition: rcm2.C:163
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 checkForNewActiveCracks(IntArray &answer, GaussPoint *gp, const FloatArray &, const FloatArray &, FloatArray &, const FloatArray &)
Definition: rcm2.C:345
void setTempMaxCrackStrain(int icrack, double val)
Definition: rcm2.h:130
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Definition: floatmatrix.C:1023
RCM2Material(int n, Domain *d)
Definition: rcm2.C:48
#define pscm_G
Definition: rcm2.h:56
IntArray tempCrackStatuses
Definition: rcm2.h:81
void letPrincipalStressVectorBe(FloatArray pv)
Definition: rcm2.h:114
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual double giveCharacteristicLength(const FloatArray &normalToCrackPlane)
Returns the size of element in the given direction, in some cases adjusted (e.g.
Definition: element.h:874
#define pscm_Ee
Definition: rcm2.h:52
FloatArray principalStrain
Definition: rcm2.h:99
FloatArray principalStress
Definition: rcm2.h:100
FloatArray tempMaxCrackStrains
Definition: rcm2.h:83
FloatArray oldPrincipalStrain
Definition: rcm2.h:99
FloatArray oldPrincipalStress
Definition: rcm2.h:100
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.
static void giveInvertedVoigtVectorMask(IntArray &answer, MaterialMode mmode)
Gives the inverted version of giveVoigtVectorMask.
FloatMatrix crackDirs
Storing direction of cracks in columwise format.
Definition: rcm2.h:87
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
const FloatMatrix & giveTempCrackDirs()
Definition: rcm2.h:126
int giveSize() const
Definition: intarray.h:203
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
bool containsOnlyZeroes() const
Returns nonzero if all coefficients of the receiver are 0, else returns zero.
Definition: floatarray.C:646
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 int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
Definition: rcm2.C:845
virtual int isCrackActive(int i) const
Definition: rcm2.C:992
void assemble(const FloatMatrix &src, const IntArray &loc)
Assembles the contribution using localization array into receiver.
Definition: floatmatrix.C:251
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Definition: floatmatrix.h:516
#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
double giveTempMaxCrackStrain(int icrack)
Definition: rcm2.h:129
int giveNumber() const
Definition: femcmpnn.h:107
void beInverseOf(const FloatMatrix &src)
Modifies receiver to become inverse of given parameter.
Definition: floatmatrix.C:835
virtual double giveNormalCrackingStress(GaussPoint *gp, double eps_cr, int i)=0
#define pscm_SOFTENING
Definition: fcm.h:56
static void transformStressVectorTo(FloatArray &answer, const FloatMatrix &base, const FloatArray &stressVector, bool transpose=false)
Transforms 3d stress vector into another coordinate system.
IntArray crackStatuses
One value from (pscm_NONE, pscm_OPEN, pscm_SOFTENING, pscm_RELOADING, pscm_UNLOADING, pscm_CLOSED.
Definition: rcm2.h:81
Class representing integration point in finite element program.
Definition: gausspoint.h:93
#define pscm_RELOADING
Definition: rcm2.h:63
Class representing solution step.
Definition: timestep.h:80
virtual double giveCrackingModulus(MatResponseMode rMode, GaussPoint *gp, double effStrain, int i)
Definition: rcm2.h:242
virtual void updateCrackStatus(GaussPoint *gp, const FloatArray &crackStrain)
Definition: rcm2.C:442
int findFirstIndexOf(int value) const
Finds index of first occurrence of given value in array.
Definition: intarray.C:331
void letPrincipalStrainVectorBe(FloatArray pv)
Definition: rcm2.h:113
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
#define rcm_SMALL_STRAIN
Definition: rcm2.h:71
virtual void updateStatusForNewCrack(GaussPoint *, int, double)
Definition: rcm2.C:419

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