OOFEM 3.0
Loading...
Searching...
No Matches
rcsdnl.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 "rcsdnl.h"
36#include "gausspoint.h"
37#include "floatmatrix.h"
38#include "floatarray.h"
39#include "mathfem.h"
40#include "nonlocalmaterialext.h"
41#include "contextioerr.h"
42#include "classfactory.h"
43#include "datastream.h"
44
45namespace oofem {
47
48RCSDNLMaterial :: RCSDNLMaterial(int n, Domain *d) : RCSDEMaterial(n, d), StructuralNonlocalMaterialExtensionInterface(d)
49{}
50
51
53RCSDNLMaterial :: giveInterface(InterfaceType type)
54{
56 return static_cast< StructuralNonlocalMaterialExtensionInterface * >(this);
57 } else {
58 return nullptr;
59 }
60}
61
62
63void
64RCSDNLMaterial :: updateBeforeNonlocAverage(const FloatArray &strainVector, GaussPoint *gp, TimeStep *tStep) const
65{
66 /* Implements the service updating local variables in given integration points,
67 * which take part in nonlocal average process. Actually, no update is necessary,
68 * because the value used for nonlocal averaging is strain vector used for nonlocal secant stiffness
69 * computation. It is therefore necessary only to store local strain in corresponding status.
70 * This service is declared at StructuralNonlocalMaterial level.
71 */
72 auto status = static_cast< RCSDNLMaterialStatus * >( this->giveStatus(gp) );
73
74 this->initTempStatus(gp);
75
76 status->setLocalStrainVectorForAverage(strainVector);
77}
78
79
80
81void
82RCSDNLMaterial :: giveRealStressVector(FloatArray &answer, GaussPoint *gp,
83 const FloatArray &totalStrain,
84 TimeStep *tStep) const
85//
86// returns real stress vector in 3d stress space of receiver according to
87// previous level of stress and current
88// strain increment, the only way, how to correctly update gp records
89//
90{
91 FloatMatrix Ds0;
92 double equivStrain;
93 FloatArray princStress, nonlocalStrain, reducedSpaceStressVector;
94 FloatArray reducedNonlocStrainVector, fullNonlocStrainVector, principalStrain;
95 FloatMatrix tempCrackDirs;
96 RCSDNLMaterialStatus *nonlocStatus, *status = static_cast< RCSDNLMaterialStatus * >( this->giveStatus(gp) );
97
98 FloatArray nonlocalContribution;
99 FloatArray reducedLocalStrainVector, localStrain;
100
101 this->initTempStatus(gp);
102 this->buildNonlocalPointTable(gp);
104
105 // compute nonlocal strain increment first
106 auto list = this->giveIPIntegrationList(gp); // !
107
108 for ( auto &lir: *list ) {
109 nonlocStatus = static_cast< RCSDNLMaterialStatus * >( this->giveStatus(lir.nearGp) );
110 nonlocalContribution = nonlocStatus->giveLocalStrainVectorForAverage();
111 nonlocalContribution.times(lir.weight);
112
113 reducedNonlocStrainVector.add(nonlocalContribution);
114 }
115
116 reducedNonlocStrainVector.times( 1. / status->giveIntegrationScale() );
117 this->endIPNonlocalAverage(gp); // !
118
119 // subtract stress independent part
121 //
122
123 reducedLocalStrainVector = totalStrain;
124
125 // subtract stress independent part
126 // note: eigenStrains (temperature) is not contained in mechanical strain stored in gp
127 // therefore it is necessary to subtract always the total eigen strain value
128 this->giveStressDependentPartOfStrainVector(nonlocalStrain, gp, reducedNonlocStrainVector,
129 tStep, VM_Total);
130 this->giveStressDependentPartOfStrainVector(localStrain, gp, reducedLocalStrainVector,
131 tStep, VM_Total);
132
133 StructuralMaterial :: giveFullSymVectorForm( fullNonlocStrainVector, nonlocalStrain, gp->giveMaterialMode() );
134
135 tempCrackDirs = status->giveTempCrackDirs();
136 this->computePrincipalValDir(principalStrain, tempCrackDirs,
137 fullNonlocStrainVector,
139
140
141 if ( status->giveTempMode() == RCSDEMaterialStatus :: rcMode ) {
142 // rotating crack mode
143
144 this->giveRealPrincipalStressVector3d(princStress, gp, principalStrain, tempCrackDirs, tStep);
145
151
155
156 this->updateCrackStatus(gp, status->giveCrackStrainVector());
157
159 this->giveMaterialStiffnessMatrix(Ds0, SecantStiffness, gp, tStep);
160
165
166 // check for any currently opening crack
167 int anyOpeningCrack = 0;
168 for ( int i = 1; i <= 3; i++ ) {
169 if ( ( status->giveTempCrackStatus(i) == pscm_SOFTENING ) ||
170 ( status->giveTempCrackStatus(i) == pscm_OPEN ) ) {
171 anyOpeningCrack++;
172 }
173 }
174
175 if ( anyOpeningCrack ) {
176 // test if transition to scalar damage mode take place
177 double minSofteningPrincStress = this->Ft, E, Le, CurrFt, Gf, Gf0, Gf1, e0, ef, ef2, damage;
178 // double minSofteningPrincStress = this->Ft, dCoeff, CurrFt, E, ep, ef, damage;
179 int ipos = 0;
180 for ( int i = 1; i <= 3; i++ ) {
181 if ( ( status->giveTempCrackStatus(i) == pscm_SOFTENING ) ||
182 ( status->giveTempCrackStatus(i) == pscm_OPEN ) ) {
183 if ( princStress.at(i) < minSofteningPrincStress ) {
184 minSofteningPrincStress = princStress.at(i);
185 ipos = i;
186 }
187 }
188 }
189
190 CurrFt = this->computeStrength( gp, status->giveCharLength(ipos) );
191
193 // next pasted from rcm2:giveEffectiveMaterialStiffnessMatrix
194 double G, minG, currG, princStressDis, princStrainDis;
195 int ii, jj;
196
197 minG = G = this->give(pscm_G, gp);
198
200 IntArray indx;
201 StructuralMaterial :: giveVoigtSymVectorMask( indx, gp->giveMaterialMode() );
202 for ( int i = 4; i <= 6; i++ ) {
203 if ( indx.contains(i) ) {
204 if ( i == 4 ) {
205 ii = 2;
206 jj = 3;
207 } else if ( i == 5 ) {
208 ii = 1;
209 jj = 3;
210 } else if ( i == 6 ) {
211 ii = 1;
212 jj = 2;
213 } else {
214 continue;
215 }
216
217 princStressDis = princStress.at(ii) -
218 princStress.at(jj);
219 princStrainDis = principalStrain.at(ii) -
220 principalStrain.at(jj);
221
222 if ( fabs(princStrainDis) < rcm_SMALL_STRAIN ) {
223 currG = G;
224 } else {
225 currG = princStressDis / ( 2.0 * princStrainDis );
226 }
227
228 //currG = Ds0.at(indi, indi);
229 minG = min(minG, currG);
230 }
231 }
232
234
235 if ( ( minSofteningPrincStress <= this->SDTransitionCoeff * CurrFt ) ||
236 ( minG <= this->SDTransitionCoeff2 * G ) ) {
237 // printf ("minSofteningPrincStress=%lf, CurrFt=%lf, SDTransitionCoeff=%lf",minSofteningPrincStress, CurrFt, this->SDTransitionCoeff);
238 // printf ("\nminG=%lf, G=%lf, SDTransitionCoeff2=%lf\n",minG, G, this->SDTransitionCoeff2);
239
240 // sd transition takes place
241 if ( ipos == 0 ) {
242 for ( int i = 1; i <= 3; i++ ) {
243 if ( ( status->giveTempCrackStatus(i) == pscm_SOFTENING ) ||
244 ( status->giveTempCrackStatus(i) == pscm_OPEN ) ) {
245 if ( ipos == 0 ) {
246 ipos = i;
247 minSofteningPrincStress = princStress.at(i);
248 }
249
250 if ( princStress.at(i) < minSofteningPrincStress ) {
251 minSofteningPrincStress = princStress.at(i);
252 ipos = i;
253 }
254 }
255 }
256 }
257
258 // test for internal consistency error
259 // we should switch to scalar damage, but no softening take place
260 if ( ipos == 0 ) {
261 //RCSDEMaterial :: OOFEM_ERROR("can not switch to sd mode, while no cracking");
262 OOFEM_ERROR("can not switch to sd mode, while no cracking");
263 }
264
265 //if (minSofteningPrincStress <= this->SDTransitionCoeff * CurrFt) printf (".");
266 //else printf (":");
267 //
268 Le = status->giveCharLength(ipos);
269 E = linearElasticMaterial->give(Ex, gp);
270 Gf = this->give(pscm_Gf, gp) / Le;
271 ef = this->ef;
272 e0 = principalStrain.at(ipos);
273 Gf0 = -CurrFt * ef * ( exp(-status->giveCrackStrain(ipos) / ef) - 1.0 ); // already disipated + 0.5*sigma0*epsilon0
274 Gf1 = Gf - Gf0;
275
276 ef2 = Gf1 / princStress.at(ipos);
277
278 //this->giveMaterialStiffnessMatrix (Ds0, SecantStiffness, gp, tStep);
279 // compute reached equivalent strain
280 equivStrain = this->computeCurrEquivStrain(gp, nonlocalStrain, E, tStep);
281 damage = this->computeDamageCoeff(equivStrain, e0, ef2);
282
283 // printf ("Gf=%lf, Gf0=%lf, damage=%lf, e0=%lf, ef2=%lf, es=%lf\n",Gf, Gf0, damage,e0,ef2,equivStrain);
284
285 status->setTransitionEpsCoeff(e0);
286 status->setEpsF2Coeff(ef2);
287 status->setDs0Matrix(Ds0);
288 status->setTempMaxEquivStrain(equivStrain);
289 status->setTempDamageCoeff(damage);
290 status->setTempMode(RCSDEMaterialStatus :: sdMode);
291 }
292 }
293 } else if ( status->giveTempMode() == RCSDEMaterialStatus :: sdMode ) {
294 // scalar damage mode
295 double E, e0, ef2;
296 // double ep, ef, E, dCoeff;
297 //FloatArray reducedSpaceStressVector;
298 double damage;
299
300 E = linearElasticMaterial->give(Ex, gp);
301 equivStrain = this->computeCurrEquivStrain(gp, nonlocalStrain, E, tStep);
302 equivStrain = max( equivStrain, status->giveTempMaxEquivStrain() );
304 ef2 = status->giveEpsF2Coeff();
305 e0 = status->giveTransitionEpsCoeff();
306 damage = this->computeDamageCoeff(equivStrain, e0, ef2);
307
308 // dCoeff = status->giveDamageStiffCoeff ();
309 // ef = status->giveDamageEpsfCoeff();
310 // ep = status->giveDamageEpspCoeff();
311 // damage = this->computeDamageCoeff (equivStrain, dCoeff, ep, ef);
313 Ds0 = * status->giveDs0Matrix();
314 Ds0.times(1.0 - damage);
317
319
323
324 status->setTempMaxEquivStrain(equivStrain);
325 status->setTempDamageCoeff(damage);
326 }
327
329 reducedSpaceStressVector.beProductOf(Ds0, localStrain);
330
331 answer = reducedSpaceStressVector;
332
333 status->letTempStressVectorBe(reducedSpaceStressVector);
334
335 status->letTempStrainVectorBe(totalStrain);
336 status->setTempNonlocalStrainVector(reducedNonlocStrainVector);
337}
338
339
340void
341RCSDNLMaterial :: initializeFrom(InputRecord &ir)
342{
343 //RCSDEMaterial::instanciateFrom (ir);
344 this->giveLinearElasticMaterial()->initializeFrom(ir);
345 StructuralNonlocalMaterialExtensionInterface :: initializeFrom(ir);
346
350 if ( SDTransitionCoeff2 > 1.0 ) {
351 SDTransitionCoeff2 = 1.0;
352 }
353
354 if ( SDTransitionCoeff2 < 0.0 ) {
355 SDTransitionCoeff2 = 0.0;
356 }
357
359 if ( R < 0.0 ) {
360 R = 0.0;
361 }
362
363 if ( ir.hasField(_IFT_RCSDNLMaterial_ef) ) { // if ef is specified, Gf is computed acordingly
365 this->Gf = this->Ft * this->ef;
366 } else if ( ir.hasField(_IFT_RCSDNLMaterial_gf) ) { // otherwise if Gf is specified, ef is computed acordingly
368 this->ef = this->Gf / this->Ft;
369 } else {
370 throw ValueInputException(ir, "none", "cannot determine Gf and ef from input data");
371 }
372}
373
374
375double
376RCSDNLMaterial :: giveMinCrackStrainsForFullyOpenCrack(GaussPoint *gp, int i) const
377{
378 return 1.e6; //this->ef;
379}
380
381
382double
383RCSDNLMaterial :: computeWeightFunction(const double R, const FloatArray &src, const FloatArray &coord) const
384{
385 // Bell shaped function decaying with the distance.
386
387 double dist = distance(src, coord);
388
389 if ( ( dist >= 0. ) && ( dist <= this->R ) ) {
390 double help = ( 1. - dist * dist / ( R * R ) );
391 return help * help;
392 }
393
394 return 0.0;
395}
396/*
397 * void
398 * RCSDNLMaterial :: updateStatusForNewCrack (GaussPoint* gp, int i, double Le)
399 * //
400 * // updates gp status when new crack-plane i is formed with charLength Le
401 * // updates Le and computes and sets minEffStrainForFullyOpenCrack
402 * //
403 * {
404 * RCSDNLMaterialStatus *status = (RCSDNLMaterialStatus*) this -> giveStatus (gp);
405 *
406 * if (Le <= 0) {
407 * char errMsg [80];
408 * sprintf (errMsg,"Element %d returned zero char length",
409 * gp->giveElement()->giveNumber());
410 * RCSDMaterial::_error(errMsg);
411 * }
412 *
413 * //status -> setCharLength(i, Le);
414 * status -> setCharLength(i, 1.0);
415 * status ->
416 * setMinCrackStrainsForFullyOpenCrack (i, this->giveMinCrackStrainsForFullyOpenCrack(gp,i));
417 * }
418 */
419
420
421
422RCSDNLMaterialStatus :: RCSDNLMaterialStatus(GaussPoint *g) :
425{
426 nonlocalStrainVector.resize( StructuralMaterial :: giveSizeOfVoigtSymVector( g->giveMaterialMode() ) );
427
428 localStrainVectorForAverage.resize( StructuralMaterial :: giveSizeOfVoigtSymVector( g->giveMaterialMode() ) );
429}
430
431
432void
433RCSDNLMaterialStatus :: printOutputAt(FILE *file, TimeStep *tStep) const
434{
435 FloatArray helpVec;
436
437 RCSDEMaterialStatus :: printOutputAt(file, tStep);
438
439 fprintf(file, "nonlocstatus { ");
440 fprintf(file, " nonloc strains ");
441 StructuralMaterial :: giveFullSymVectorForm( helpVec, nonlocalStrainVector, gp->giveMaterialMode() );
442 for ( auto &val : helpVec ) {
443 fprintf( file, " %.4e", val );
444 }
445
446 fprintf(file, "}\n");
447}
448
449
450void
451RCSDNLMaterialStatus :: initTempStatus()
452//
453// initializes temp variables according to variables form previous equlibrium state.
454// builds new crackMap
455//
456{
457 RCSDEMaterialStatus :: initTempStatus();
458
460}
461
462
463void
464RCSDNLMaterialStatus :: updateYourself(TimeStep *tStep)
465//
466// updates variables (nonTemp variables describing situation at previous equilibrium state)
467// after a new equilibrium state has been reached
468// temporary variables are having values corresponding to newly reched equilibrium.
469//
470{
471 RCSDEMaterialStatus :: updateYourself(tStep);
473}
474
475
476void
477RCSDNLMaterialStatus :: saveContext(DataStream &stream, ContextMode mode)
478{
479 RCSDEMaterialStatus :: saveContext(stream, mode);
480
482 if ( ( iores = nonlocalStrainVector.storeYourself(stream) ) != CIO_OK ) {
483 THROW_CIOERR(iores);
484 }
485}
486
487
488void
489RCSDNLMaterialStatus :: restoreContext(DataStream &stream, ContextMode mode)
490{
491 RCSDEMaterialStatus :: restoreContext(stream, mode);
492
494 if ( ( iores = nonlocalStrainVector.restoreYourself(stream) ) != CIO_OK ) {
495 THROW_CIOERR(iores);
496 }
497}
498
499
500Interface *
501RCSDNLMaterialStatus :: giveInterface(InterfaceType type)
502{
504 return this;
505 } else {
506 return nullptr;
507 }
508}
509
510
511int
512RCSDNLMaterial :: packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
513{
514 RCSDNLMaterialStatus *status = static_cast< RCSDNLMaterialStatus * >( this->giveStatus(ip) );
515
516 this->buildNonlocalPointTable(ip);
518
519 return (status->giveLocalStrainVectorForAverage().storeYourself(buff) == CIO_OK);
520}
521
522int
523RCSDNLMaterial :: unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
524{
525 int result;
526 RCSDNLMaterialStatus *status = static_cast< RCSDNLMaterialStatus * >( this->giveStatus(ip) );
527 FloatArray localStrainVectorForAverage;
528
529 result = (localStrainVectorForAverage.restoreYourself(buff) == CIO_OK);
530 status->setLocalStrainVectorForAverage(localStrainVectorForAverage);
531 return result;
532}
533
534int
535RCSDNLMaterial :: estimatePackSize(DataStream &buff, GaussPoint *ip)
536{
537 //
538 // Note: status localStrainVectorForAverage memeber must be properly sized!
539 //
540 RCSDNLMaterialStatus *status = static_cast< RCSDNLMaterialStatus * >( this->giveStatus(ip) );
541
542 return status->giveLocalStrainVectorForAverage().givePackSize(buff);
543}
544
545} // end namespace oofem
#define E(a, b)
#define REGISTER_Material(class)
double & at(Index i)
Definition floatarray.h:202
int givePackSize(DataStream &buff) const
Definition floatarray.C:943
contextIOResultType storeYourself(DataStream &stream) const
Definition floatarray.C:894
contextIOResultType restoreYourself(DataStream &stream)
Definition floatarray.C:917
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void add(const FloatArray &src)
Definition floatarray.C:218
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 bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
bool contains(int value) const
Definition intarray.h:292
GaussPoint * gp
Associated integration point.
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Definition material.C:206
virtual void initTempStatus(GaussPoint *gp) const
Definition material.C:221
void buildNonlocalPointTable(GaussPoint *gp) const
void updateDomainBeforeNonlocAverage(TimeStep *tStep) const
std ::vector< localIntegrationRecord > * giveIPIntegrationList(GaussPoint *gp) const
double giveIntegrationScale()
Returns associated integration scale.
const IntArray & giveTempCrackStatus()
Definition rcm2.h:130
double giveCrackStrain(int icrack) const
Definition rcm2.h:135
const FloatArray & giveCrackStrainVector() const
Definition rcm2.h:134
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
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
RCSDEMaterialStatus(GaussPoint *g)
Definition rcsde.C:448
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
void setTempDamageCoeff(double val)
Definition rcsde.h:80
__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 give(int aProperty, GaussPoint *gp) const override
Definition rcsde.C:263
RCSDEMaterial(int n, Domain *d)
Definition rcsde.C:49
double computeDamageCoeff(double, double, double) const
Definition rcsde.C:223
double computeCurrEquivStrain(GaussPoint *, const FloatArray &, double, TimeStep *tStep) const
Definition rcsde.C:235
double SDTransitionCoeff
Definition rcsde.h:118
const FloatArray & giveLocalStrainVectorForAverage()
Definition rcsdnl.h:74
FloatArray localStrainVectorForAverage
Definition rcsdnl.h:63
void setTempNonlocalStrainVector(FloatArray ls)
Definition rcsdnl.h:72
FloatArray nonlocalStrainVector
Definition rcsdnl.h:63
void setLocalStrainVectorForAverage(FloatArray ls)
Definition rcsdnl.h:75
FloatArray tempNonlocalStrainVector
Definition rcsdnl.h:63
double computeStrength(GaussPoint *, double) const override
Definition rcsdnl.h:153
double SDTransitionCoeff2
Definition rcsdnl.h:102
void letTempStressVectorBe(const FloatArray &v)
Assigns tempStressVector to given vector v.
void letTempStrainVectorBe(const FloatArray &v)
Assigns tempStrainVector to given vector v.
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
#define THROW_CIOERR(e)
#define OOFEM_ERROR(...)
Definition error.h:79
#define pscm_SOFTENING
Definition fcm.h:58
#define IR_GIVE_FIELD(__ir, __value, __id)
Definition inputrecord.h:67
#define Ex
Definition matconst.h:59
long ContextMode
Definition contextmode.h:43
FloatArrayF< N > min(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
FloatArrayF< N > max(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
@ principal_strain
For computing principal strains from engineering strains.
@ NonlocalMaterialStatusExtensionInterfaceType
@ NonlocalMaterialExtensionInterfaceType
double distance(const FloatArray &x, const FloatArray &y)
#define pscm_Gf
Definition rcm2.h:54
#define rcm_SMALL_STRAIN
Definition rcm2.h:71
#define pscm_G
Definition rcm2.h:56
#define pscm_OPEN
Definition rcm2.h:61
#define _IFT_RCSDNLMaterial_sdtransitioncoeff
Definition rcsdnl.h:49
#define _IFT_RCSDNLMaterial_gf
Definition rcsdnl.h:53
#define _IFT_RCSDNLMaterial_sdtransitioncoeff2
Definition rcsdnl.h:50
#define _IFT_RCSDNLMaterial_r
Definition rcsdnl.h:51
#define _IFT_RCSDNLMaterial_ft
Definition rcsdnl.h:48
#define _IFT_RCSDNLMaterial_ef
Definition rcsdnl.h:52

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