OOFEM 3.0
Loading...
Searching...
No Matches
rcsde.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 - 2025 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 "rcsde.h"
36#include "gausspoint.h"
37#include "floatmatrix.h"
38#include "floatarray.h"
39#include "mathfem.h"
41#include "datastream.h"
42#include "contextioerr.h"
43#include "classfactory.h"
44
45
46namespace oofem {
48
49RCSDEMaterial :: RCSDEMaterial(int n, Domain *d) : RCM2Material(n, d)
50{
52}
53
54
55RCSDEMaterial :: ~RCSDEMaterial()
56{
58}
59
60
61void
62RCSDEMaterial :: giveRealStressVector(FloatArray &answer, GaussPoint *gp,
63 const FloatArray &totalStrain,
64 TimeStep *tStep) const
65//
66// returns real stress vector in 3d stress space of receiver according to
67// previous level of stress and current
68// strain increment, the only way, how to correctly update gp records
69//
70{
71 FloatMatrix Ds0;
72 double equivStrain;
73 FloatArray princStress, reducedAnswer;
74 FloatArray reducedStrainVector, strainVector, principalStrain;
75 FloatArray reducedSpaceStressVector;
76 FloatMatrix tempCrackDirs;
77 RCSDEMaterialStatus *status = static_cast< RCSDEMaterialStatus * >( this->giveStatus(gp) );
78
79 this->initTempStatus(gp);
80
81 // subtract stress independent part
82 // note: eigenStrains (temperature) is not contained in mechanical strain stored in gp
83 // therefore it is necessary to subtract always the total eigen strain value
84 this->giveStressDependentPartOfStrainVector(reducedStrainVector, gp, totalStrain,
85 tStep, VM_Total);
86
87 StructuralMaterial :: giveFullSymVectorForm( strainVector, reducedStrainVector, gp->giveMaterialMode() );
88
89 tempCrackDirs = status->giveTempCrackDirs();
90 this->computePrincipalValDir(principalStrain, tempCrackDirs,
91 strainVector,
93
94
95 if ( status->giveTempMode() == RCSDEMaterialStatus :: rcMode ) {
96 // rotating crack mode
97
98 this->giveRealPrincipalStressVector3d(princStress, gp, principalStrain, tempCrackDirs, tStep);
99 princStress.resize(6);
100 tempCrackDirs = status->giveTempCrackDirs();
101 answer = this->transformStressVectorTo(tempCrackDirs, princStress, 1);
102
103 StructuralMaterial :: giveReducedSymVectorForm( reducedSpaceStressVector, answer, gp->giveMaterialMode() );
104 status->letTempStressVectorBe(reducedSpaceStressVector);
105
106 this->updateCrackStatus(gp, status->giveCrackStrainVector());
107
108 StructuralMaterial :: giveReducedSymVectorForm( reducedAnswer, answer, gp->giveMaterialMode() );
109 answer = reducedAnswer;
110
111 // test if transition to scalar damage mode take place
112 double minSofteningPrincStress = this->Ft, E, Le, CurrFt, Ft, Gf, Gf0, Gf1, e0, ef, ef2, damage;
113 int ipos = 0;
114 for ( int i = 1; i <= 3; i++ ) {
115 if ( status->giveTempCrackStatus(i) == pscm_SOFTENING ) {
116 if ( princStress.at(i) < minSofteningPrincStress ) {
117 minSofteningPrincStress = princStress.at(i);
118 ipos = i;
119 }
120 }
121 }
122
123 CurrFt = this->computeStrength( gp, status->giveCharLength(ipos) );
124
125 if ( minSofteningPrincStress <= this->SDTransitionCoeff * CurrFt ) {
126 // sd transition takes place
127
128 Le = status->giveCharLength(ipos);
129 E = linearElasticMaterial->give(Ex, gp);
130 Gf = this->give(pscm_Gf, gp) / Le;
131 Ft = this->computeStrength(gp, Le);
132 ef = Gf / Ft;
133 e0 = principalStrain.at(ipos);
134 Gf0 = -CurrFt * ef * ( exp(-status->giveCrackStrain(ipos) / ef) - 1.0 ); // already disipated + 0.5*sigma0*epsilon0
135 Gf1 = Gf - Gf0;
136
137 ef2 = Gf1 / princStress.at(ipos);
138
139 this->giveMaterialStiffnessMatrix(Ds0, SecantStiffness, gp, tStep);
140 // compute reached equivalent strain
141 equivStrain = this->computeCurrEquivStrain(gp, reducedStrainVector, E, tStep);
142 damage = this->computeDamageCoeff(equivStrain, e0, ef2);
143
144
145 status->setTransitionEpsCoeff(e0);
146 status->setEpsF2Coeff(ef2);
147 status->setDs0Matrix(Ds0);
148 status->setTempMaxEquivStrain(equivStrain);
149 status->setTempDamageCoeff(damage);
150 status->setTempMode(RCSDEMaterialStatus :: sdMode);
151 }
152 } else {
153 // scalar damage mode
154 double E, e0, ef2;
155 double damage;
156 //int ipos;
157
158 E = linearElasticMaterial->give(Ex, gp);
159 equivStrain = this->computeCurrEquivStrain(gp, reducedStrainVector, E, tStep);
160 equivStrain = max( equivStrain, status->giveTempMaxEquivStrain() );
161 reducedSpaceStressVector.beProductOf(* status->giveDs0Matrix(), reducedStrainVector);
162 //dCoeff = status->giveDamageStiffCoeff();
163 ef2 = status->giveEpsF2Coeff();
164 e0 = status->giveTransitionEpsCoeff();
165 damage = this->computeDamageCoeff(equivStrain, e0, ef2);
166 reducedSpaceStressVector.times(1.0 - damage);
167
168 answer = reducedSpaceStressVector;
169
170 status->letTempStressVectorBe(reducedSpaceStressVector);
171
172 status->setTempMaxEquivStrain(equivStrain);
173 status->setTempDamageCoeff(damage);
174 }
175
176 status->letTempStrainVectorBe(totalStrain);
177}
178
179
180void
181RCSDEMaterial :: giveEffectiveMaterialStiffnessMatrix(FloatMatrix &answer,
182 MatResponseMode rMode, GaussPoint *gp,
183 TimeStep *tStep) const
184//
185// returns effective material stiffness matrix in full form
186// for gp stress strain mode
187//
188{
189 RCSDEMaterialStatus *status = static_cast< RCSDEMaterialStatus * >( this->giveStatus(gp) );
190
191 if ( status->giveTempMode() == RCSDEMaterialStatus :: rcMode ) {
192 // rotating crack mode
193
194 RCM2Material :: giveEffectiveMaterialStiffnessMatrix(answer, rMode, gp, tStep);
195 return;
196 } else {
197 // rcsd mode
198
199 if ( ( rMode == TangentStiffness ) || ( rMode == SecantStiffness ) ) {
200 FloatMatrix reducedAnswer;
201 double dCoeff;
202
203 reducedAnswer = * status->giveDs0Matrix();
204 dCoeff = 1.0 - status->giveDamageCoeff();
205 if ( dCoeff < RCSDE_DAMAGE_EPS ) {
206 dCoeff = RCSDE_DAMAGE_EPS;
207 }
208
209 reducedAnswer.times(dCoeff);
210
211 answer = reducedAnswer;
212 } else if ( rMode == ElasticStiffness ) {
213 this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, rMode, gp, tStep);
214 return;
215 } else {
216 OOFEM_ERROR("usupported mode");
217 }
218 }
219}
220
221
222double
223RCSDEMaterial :: computeDamageCoeff(double equivStrain, double e0, double ef2) const
224{
225 double damage = 1. - e0 / equivStrain *exp(-( equivStrain - e0 ) / ef2);
226 if ( damage > 1.0 ) {
227 return 1.0;
228 }
229
230 return damage;
231}
232
233
234double
235RCSDEMaterial :: computeCurrEquivStrain(GaussPoint *gp, const FloatArray &reducedTotalStrainVector, double e, TimeStep *tStep) const
236{
237 FloatArray effStress, princEffStress, fullEffStress;
238 FloatMatrix De;
239 double answer = 0.0;
240
241 linearElasticMaterial->giveStiffnessMatrix(De, TangentStiffness, gp, tStep);
242 effStress.beProductOf(De, reducedTotalStrainVector);
243 StructuralMaterial :: giveFullSymVectorForm( fullEffStress, effStress, gp->giveMaterialMode() );
244
245 this->computePrincipalValues(princEffStress, fullEffStress, principal_stress);
246 for ( int i = 1; i <= 3; i++ ) {
247 answer = max( answer, macbra( princEffStress.at(i) ) );
248 }
249
250 return answer / e;
251}
252
253
254void
255RCSDEMaterial :: initializeFrom(InputRecord &ir)
256{
257 RCM2Material :: initializeFrom(ir);
259}
260
261
262double
263RCSDEMaterial :: give(int aProperty, GaussPoint *gp) const
264// Returns the value of the property aProperty (e.g. the Young's modulus
265// 'E') of the receiver.
266{
267 if ( aProperty == pscm_SDTransitionCoeff ) {
268 return this->SDTransitionCoeff;
269 }
270
271 return RCM2Material :: give(aProperty, gp);
272}
273
274
275int
276RCSDEMaterial :: checkSizeLimit(GaussPoint *gp, double charLength) const
277//
278// checks if element size (charLength) is too big
279// so that tension strength must be reduced followed
280// by sudden stress drop
281//
282{
283 double Ee, Gf, Ft, LeCrit;
284
285 Ee = this->give(pscm_Ee, gp);
286 Gf = this->give(pscm_Gf, gp);
287 Ft = this->give(pscm_Ft, gp);
288
289 LeCrit = 2.0 * Gf * Ee / ( Ft * Ft );
290 return ( charLength < LeCrit );
291}
292
293
294double
295RCSDEMaterial :: computeStrength(GaussPoint *gp, double charLength) const
296//
297// computes strength for given gp,
298// which may be reduced according to length of "fracture process zone"
299// to be energetically correct
300//
301{
302 double Ee, Gf, Ft;
303
304 Ee = this->give(pscm_Ee, gp);
305 Gf = this->give(pscm_Gf, gp);
306 Ft = this->give(pscm_Ft, gp);
307
308 if ( this->checkSizeLimit(gp, charLength) ) { } else {
309 // we reduce Ft and there is no softening but sudden drop
310 Ft = sqrt(2. * Ee * Gf / charLength);
311 //
312 OOFEM_LOG_RELEVANT("Reducing Ft to %f in element %d, gp %d, Le %f\n",
313 Ft, gp->giveElement()->giveNumber(), gp->giveNumber(), charLength);
314 }
315
316 return Ft;
317}
318
319
320double
321RCSDEMaterial :: giveMinCrackStrainsForFullyOpenCrack(GaussPoint *gp, int i) const
322//
323// computes MinCrackStrainsForFullyOpenCrack for given gp and i-th crack
324//
325{
326 /*
327 * RCM2MaterialStatus *status = (RCM2MaterialStatus*) this -> giveStatus (gp);
328 * double Le, Gf, Ft;
329 *
330 * Le = status -> giveCharLength (i);
331 * Gf = this->give(pscm_Gf);
332 * Ft = this->computeStrength (gp, Le);
333 *
334 * return Gf/(Le*Ft); // Exponential softening
335 */
336 return 1.e6;
337}
338
339
340/*
341 * void
342 * RCSDEMaterial :: updateStatusForNewCrack (GaussPoint* gp, int i, double Le)
343 * //
344 * // updates gp status when new crack-plane i is formed with charLength Le
345 * // updates Le and computes and sets minEffStrainForFullyOpenCrack
346 * //
347 * {
348 * RCM2MaterialStatus *status = (RCM2MaterialStatus*) this -> giveStatus (gp);
349 *
350 * if (Le <= 0) {
351 * char errMsg [80];
352 * sprintf (errMsg,"Element %d returned zero char length",
353 * gp->giveElement()->giveNumber());
354 * OOFEM_ERROR(errMsg);
355 * }
356 *
357 * status -> setCharLength(i, Le);
358 * status ->
359 * setMinCrackStrainsForFullyOpenCrack (i, this->giveMinCrackStrainsForFullyOpenCrack(gp,i));
360 * }
361 */
362
363double
364RCSDEMaterial :: giveCrackingModulus(MatResponseMode rMode, GaussPoint *gp,
365 double crackStrain, int i) const
366//
367// returns current cracking modulus according to crackStrain for i-th
368// crackplane
369// now linear softening is implemented
370// see also CreateStatus () function.
371// softening modulus represents a relation between the normal crack strain
372// rate and the normal stress rate.
373//
374{
375 // double Ee, Gf;
376 double Gf, Cf, Ft, Le, ef;
377 RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
378
379 //
380 // now we have to set proper reduced strength and softening modulus Et
381 // in order to obtain energetically correct solution using the concept
382 // of fracture energy Gf, which is a material constant.
383 //
384 //Ee = this->give(pscm_Ee);
385 Gf = this->give(pscm_Gf, gp);
386 Le = status->giveCharLength(i);
387 Ft = this->computeStrength(gp, Le);
388 //minEffStrainForFullyOpenCrack = this->giveMinCrackStrainsForFullyOpenCrack(gp,i);
389 ef = Gf / ( Le * Ft );
390
391 if ( rMode == TangentStiffness ) {
392 if ( this->checkSizeLimit(gp, Le) ) {
393 if ( crackStrain >= status->giveTempMaxCrackStrain(i) ) {
394 // further softening
395 Cf = -Ft / ef *exp(-crackStrain / ef);
396 } else {
397 // unloading or reloading regime
398 Cf = Ft / status->giveTempMaxCrackStrain(i) * exp(-status->giveTempMaxCrackStrain(i) / ef);
399 }
400 } else {
401 Cf = 0.;
402 }
403 } else {
404 if ( this->checkSizeLimit(gp, Le) ) {
405 Cf = Ft / status->giveTempMaxCrackStrain(i) * exp(-status->giveTempMaxCrackStrain(i) / ef);
406 } else {
407 Cf = 0.;
408 }
409 }
410
411 return Cf;
412}
413
414
415double
416RCSDEMaterial :: giveNormalCrackingStress(GaussPoint *gp, double crackStrain, int i) const
417//
418// returns receivers Normal Stress in crack i for given cracking strain
419//
420{
421 double Gf, Ft, Le, answer, ef;
422 RCM2MaterialStatus *status = static_cast< RCM2MaterialStatus * >( this->giveStatus(gp) );
423
424 Le = status->giveCharLength(i);
425 Ft = this->computeStrength(gp, Le);
426 Gf = this->give(pscm_Gf, gp);
427 ef = Gf / ( Le * Ft );
428
429 if ( this->checkSizeLimit(gp, Le) ) {
430 if ( crackStrain >= status->giveTempMaxCrackStrain(i) ) {
431 // further softening
432 answer = Ft * exp(-crackStrain / ef);
433 } else {
434 // crack closing
435 // or unloading or reloading regime
436 answer = Ft * crackStrain / status->giveTempMaxCrackStrain(i) *
437 exp(-status->giveTempMaxCrackStrain(i) / ef);
438 }
439 } else {
440 answer = 0.;
441 }
442
443 return answer;
444}
445
446
447
448RCSDEMaterialStatus :: RCSDEMaterialStatus(GaussPoint *g) :
450{}
451
452
453void
454RCSDEMaterialStatus :: printOutputAt(FILE *file, TimeStep *tStep) const
455{
456 char s [ 11 ];
457
458 StructuralMaterialStatus :: printOutputAt(file, tStep);
459 fprintf(file, "status { ");
460 if ( this->giveTempMode() == rcMode ) {
461 fprintf(file, "mode :rc ");
462
463 if ( this->giveTempAlreadyCrack() ) {
464 for ( int i = 1; i <= 3; i++ ) {
465 switch ( crackStatuses.at(i) ) {
466 case pscm_NONE:
467 strcpy(s, "NONE");
468 break;
469 case pscm_OPEN:
470 strcpy(s, "OPEN");
471 break;
472 case pscm_SOFTENING:
473 strcpy(s, "SOFTENING");
474 break;
475 case pscm_RELOADING:
476 strcpy(s, "RELOADING");
477 break;
478 case pscm_UNLOADING:
479 strcpy(s, "UNLOADING");
480 break;
481 default:
482 strcpy(s, "UNKNOWN");
483 break;
484 }
485
486 fprintf( file, "crack %d {status %s, normal to crackplane { %f %f %f }} ",
487 i, s, crackDirs.at(1, i), crackDirs.at(2, i), crackDirs.at(3, i) );
488 }
489 }
490 } else {
491 // sd mode
492 fprintf( file, "mode :sd, damageCoeff = %f ", this->giveTempDamageCoeff() );
493 }
494
495 fprintf(file, "}\n");
496}
497
498
499void
500RCSDEMaterialStatus :: initTempStatus()
501//
502// initializes temp variables according to variables form previous equlibrium state.
503// builds new crackMap
504//
505{
506 RCM2MaterialStatus :: initTempStatus();
507
511}
512
513
514void
515RCSDEMaterialStatus :: updateYourself(TimeStep *tStep)
516//
517// updates variables (nonTemp variables describing situation at previous equilibrium state)
518// after a new equilibrium state has been reached
519// temporary variables are having values corresponding to newly reched equilibrium.
520//
521{
522 RCM2MaterialStatus :: updateYourself(tStep);
523
527}
528
529
530void
531RCSDEMaterialStatus :: saveContext(DataStream &stream, ContextMode mode)
532{
533 RCM2MaterialStatus :: saveContext(stream, mode);
534
535 if ( !stream.write(maxEquivStrain) ) {
537 }
538
539 if ( !stream.write(damageCoeff) ) {
541 }
542
543 int _mode = rcsdMode;
544 if ( !stream.write(_mode) ) {
546 }
547
548 if ( !stream.write(transitionEps) ) {
550 }
551
552 if ( !stream.write(epsF2) ) {
554 }
555
557 if ( ( iores = Ds0.storeYourself(stream) ) != CIO_OK ) {
558 THROW_CIOERR(iores);
559 }
560}
561
562
563void
564RCSDEMaterialStatus :: restoreContext(DataStream &stream, ContextMode mode)
565{
566 RCM2MaterialStatus :: restoreContext(stream, mode);
567
568 if ( !stream.read(maxEquivStrain) ) {
570 }
571
572 if ( !stream.read(damageCoeff) ) {
574 }
575
576 int _mode;
577 if ( !stream.read(_mode) ) {
579 }
580
581 rcsdMode = ( __rcsdModeType ) _mode;
582 if ( !stream.read(transitionEps) ) {
584 }
585
586 if ( !stream.read(epsF2) ) {
588 }
589
591 if ( ( iores = Ds0.restoreYourself(stream) ) != CIO_OK ) {
592 THROW_CIOERR(iores);
593 }
594}
595} // end namespace oofem
#define E(a, b)
#define REGISTER_Material(class)
virtual int read(int *data, std::size_t count)=0
Reads count integer values into array pointed by data.
virtual int write(const int *data, std::size_t count)=0
Writes count integer values from array pointed by data.
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void times(double s)
Definition floatarray.C:834
void times(double f)
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition gausspoint.h:190
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Definition material.C:206
virtual void initTempStatus(GaussPoint *gp) const
Definition material.C:221
const IntArray & giveTempCrackStatus()
Definition rcm2.h:130
IntArray crackStatuses
One value from (pscm_NONE, pscm_OPEN, pscm_SOFTENING, pscm_RELOADING, pscm_UNLOADING,...
Definition rcm2.h:81
double giveCrackStrain(int icrack) const
Definition rcm2.h:135
RCM2MaterialStatus(GaussPoint *g)
Definition rcm2.C:977
int giveTempAlreadyCrack() const
Definition rcm2.h:120
const FloatArray & giveCrackStrainVector() const
Definition rcm2.h:134
double giveTempMaxCrackStrain(int icrack)
Definition rcm2.h:128
FloatMatrix crackDirs
Storing direction of cracks in columwise format.
Definition rcm2.h:87
double giveCharLength(int icrack) const
Definition rcm2.h:140
const FloatMatrix & giveTempCrackDirs()
Definition rcm2.h:125
void giveRealPrincipalStressVector3d(FloatArray &answer, GaussPoint *, FloatArray &, FloatMatrix &, TimeStep *) const
Definition rcm2.C:150
LinearElasticMaterial * linearElasticMaterial
Definition rcm2.h:178
RCM2Material(int n, Domain *d)
Definition rcm2.C:48
LinearElasticMaterial * giveLinearElasticMaterial() const
Definition rcm2.h:193
virtual void updateCrackStatus(GaussPoint *gp, const FloatArray &crackStrain) const
Definition rcm2.C:429
virtual void giveMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) const
Definition rcm2.C:83
double giveTransitionEpsCoeff() const
Definition rcsde.h:84
void setTempMaxEquivStrain(double val)
Definition rcsde.h:76
void setDs0Matrix(FloatMatrix &mtrx)
Definition rcsde.h:82
void setTransitionEpsCoeff(double val)
Definition rcsde.h:85
double giveEpsF2Coeff() const
Definition rcsde.h:86
double giveTempMaxEquivStrain() const
Definition rcsde.h:75
__rcsdModeType rcsdMode
Definition rcsde.h:68
void setTempDamageCoeff(double val)
Definition rcsde.h:80
double giveDamageCoeff() const
Definition rcsde.h:94
__rcsdModeType giveTempMode() const
Definition rcsde.h:89
void setTempMode(__rcsdModeType mode)
Definition rcsde.h:90
const FloatMatrix * giveDs0Matrix()
Definition rcsde.h:81
void setEpsF2Coeff(double val)
Definition rcsde.h:87
double giveTempDamageCoeff() const
Definition rcsde.h:79
__rcsdModeType tempRcsdMode
Definition rcsde.h:68
int checkSizeLimit(GaussPoint *gp, double) const override
Definition rcsde.C:276
double give(int aProperty, GaussPoint *gp) const override
Definition rcsde.C:263
double computeDamageCoeff(double, double, double) const
Definition rcsde.C:223
double computeCurrEquivStrain(GaussPoint *, const FloatArray &, double, TimeStep *tStep) const
Definition rcsde.C:235
double computeStrength(GaussPoint *gp, double) const override
Definition rcsde.C:295
double SDTransitionCoeff
Definition rcsde.h:118
void letTempStressVectorBe(const FloatArray &v)
Assigns tempStressVector to given vector v.
void letTempStrainVectorBe(const FloatArray &v)
Assigns tempStrainVector to given vector v.
static void computePrincipalValues(FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode)
Common functions for convenience.
static void computePrincipalValDir(FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode)
void giveStressDependentPartOfStrainVector(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) const
static FloatArrayF< 6 > transformStressVectorTo(const FloatMatrixF< 3, 3 > &base, const FloatArrayF< 6 > &stress, bool transpose=false)
#define THROW_CIOERR(e)
#define OOFEM_ERROR(...)
Definition error.h:79
#define pscm_NONE
Definition fcm.h:56
#define pscm_SOFTENING
Definition fcm.h:58
#define IR_GIVE_FIELD(__ir, __value, __id)
Definition inputrecord.h:67
#define OOFEM_LOG_RELEVANT(...)
Definition logger.h:142
#define Ex
Definition matconst.h:59
long ContextMode
Definition contextmode.h:43
double macbra(double x)
Returns the positive part of given float.
Definition mathfem.h:128
FloatArrayF< N > max(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
@ principal_strain
For computing principal strains from engineering strains.
@ principal_stress
For computing principal stresses.
@ CIO_IOERR
General IO error.
#define pscm_Gf
Definition rcm2.h:54
#define pscm_Ee
Definition rcm2.h:52
#define pscm_UNLOADING
Definition rcm2.h:64
#define pscm_OPEN
Definition rcm2.h:61
#define pscm_Ft
Definition rcm2.h:57
#define pscm_RELOADING
Definition rcm2.h:63
#define pscm_SDTransitionCoeff
Definition rcsd.h:52
#define RCSDE_DAMAGE_EPS
Definition rcsde.h:53
#define _IFT_RCSDEMaterial_sdtransitioncoeff
Definition rcsde.h:47

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011