OOFEM 3.0
Loading...
Searching...
No Matches
anisodamagemodel.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 program is free software; you can redistribute it and/or modify
21 * it under the terms of the GNU General Public License as published by
22 * the Free Software Foundation; either version 2 of the License, or
23 * (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
28 * GNU General Public License for more details.
29 *
30 * You should have received a copy of the GNU General Public License
31 * along with this program; if not, write to the Free Software
32 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
33 */
34
35// This code is based on the anisotropic damage model proposed by Desmorat, Gatuingt and Ragueneau in
36// their paper "Nonlocal anisotropic damage model and related computational aspects for quasi-brittle material"
37// published in Engineering Fracture Mechanics 74 (2007) 1539-1560.
38
39#include "anisodamagemodel.h"
41#include "floatmatrix.h"
42#include "floatarray.h"
43#include "mathfem.h"
44#include "datastream.h"
45#include "contextioerr.h"
46#include "dynamicinputrecord.h"
47#include "gausspoint.h"
48#include "classfactory.h"
49
50namespace oofem {
52
53AnisotropicDamageMaterial :: AnisotropicDamageMaterial(int n, Domain *d) : StructuralMaterial(n, d),
55{}
56
57
58bool
59AnisotropicDamageMaterial :: hasMaterialModeCapability(MaterialMode mode) const
60{
61 return mode == _3dMat || mode == _PlaneStress;
62 //return mode == _3dMat || mode == _PlaneStress || mode == _PlaneStrain || mode == _1dMat;
63}
64
65
66//********************************************************
67// plane stress implementation by Milan Jirasek
68//********************************************************
69
71AnisotropicDamageMaterial :: giveRealStressVector_PlaneStress(const FloatArrayF<3> &totalStrain, GaussPoint *gp, TimeStep *atTime) const
72//
73// special version of the stress-evaluation algorithm applicable under plane stress
74// based on the report by Jirasek & Suarez, 25 April 2014
75//
76// uses the following special methods:
77// computePrincValDir2D ... evaluation of eigenvalues and eigenvector of a symmetric 2x2 matrix
78// computeDamage ... application of the damage law
79// computeTraceD(double) ... function "g" evaluating trace(D) from kappa
80// checkPrincVal2D ... checking whether principal values are <= 1
81// computeOutOfPlaneStrain ... evaluation of epsZ for given in-plane strain and damage
82// computeDimensionlessOutOfPlaneStress ... evaluation of scaled sigZ for given strain and damage
83// computeInplaneStress ... evaluation of in-plane stress for given strain and damage
84//
85// Note: Only Mazars' equivalent strain is implemented for this case,
86// but the damage law can be varied and there are 4 cases implemented in computeTraceD(double).
87// There is another method computeTraceD(...), which is used by the 3D algorithm.
88{
89#define AD_TOLERANCE 1.e-10 // convergence tolerance for the internal iteration used under plane stress
90#define AD_MAXITER 20 // maximum number of internal iterations used under plane stress
91 //this->initGpForNewStep(gp);
92 this->initTempStatus(gp);
93 // subtract the stress-independent part of strains (e.g. due to temperature)
94 FloatArray eps;
95 this->giveStressDependentPartOfStrainVector(eps, gp, totalStrain, atTime, VM_Total);
96
97 // compute in-plane principal strains: epsilon1 and epsilon2
98 // and the components of the first principal strain direction: ceps and seps
99 double epsilon1, epsilon2, ceps, seps;
100 this->computePrincValDir2D(epsilon1, epsilon2, ceps, seps, eps.at(1), eps.at(2), eps.at(3) / 2.);
101
102 // compute the equivalent strain resulting from the in-plane strains only: equivStrainInPlane
103 // only Mazars' strain implemented so far
104 double equivStrainInPlane = 0.;
105 if ( epsilon1 > 0. ) {
106 equivStrainInPlane += epsilon1 * epsilon1;
107 if ( epsilon2 > 0. ) { // note that we always have epsilon2<=epsilon1
108 equivStrainInPlane += epsilon2 * epsilon2;
109 }
110 }
111 equivStrainInPlane = sqrt(equivStrainInPlane);
112
113 // compute ez0 = maximum value of the out-of-plane strain that still has no influence on damage
114 // formula (85)
115 double ez0 = 0.;
116 AnisotropicDamageMaterialStatus *status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
117 double kappa = status->giveKappa();
118 if ( equivStrainInPlane < kappa ) {
119 ez0 = sqrt(kappa * kappa - equivStrainInPlane * equivStrainInPlane);
120 }
121
122 // compute the trial out-of-plane strain, assuming that the equivalent strain is equal to equivStrainInPlane
123 FloatMatrix Dn = status->giveDamage();
124 FloatMatrix tempDamage = Dn;
125 // update damage, if needed
126 if ( equivStrainInPlane > kappa ) {
127 // formula (76)
128 this->computeDamage(tempDamage, Dn, kappa, epsilon1, epsilon2, ceps, seps, 0.);
129 }
130 // compute the first trial value (with volumetric compression assumed, i.e., bulk damage deactivated)
131 // formula (79)
132 double ez1 = this->computeOutOfPlaneStrain(eps, tempDamage, false);
133
134 double strainTraceInPlane = eps.at(1) + eps.at(2);
135 if ( ( ez1 + strainTraceInPlane > 0. ) || ( ez1 > ez0 ) ) { // first trial value is not admissible
136 bool iteration_needed = true;
137 if ( ez0 + strainTraceInPlane > 0. ) {
138 // compute the second trial value (with volumetric tension assumed, i.e., bulk damage activated)
139 ez1 = this->computeOutOfPlaneStrain(eps, tempDamage, true);
140 iteration_needed = ( ( ez1 > ez0 ) || ( ez1 + strainTraceInPlane < 0. ) );
141 }
142 if ( iteration_needed ) { // the second trial value is not admissible
143 // iterative procedure - out-of-plane strain does have an effect on damage
144 // scaled formula (78)
145 double s0 = this->computeDimensionlessOutOfPlaneStress(eps, ez0, tempDamage);
146 this->computeDamage(tempDamage, Dn, kappa, epsilon1, epsilon2, ceps, seps, ez1);
147 double s1 = this->computeDimensionlessOutOfPlaneStress(eps, ez1, tempDamage);
148 int count = 0;
149 while ( fabs(s1) > AD_TOLERANCE * this->kappa0 ) {
150 if ( ++count > AD_MAXITER ) {
151 OOFEM_ERROR("No convergence in AnisotropicDamageMaterial :: giveRealStressVector for the plane stress case\n");
152 }
153 double ez2 = ( ez1 * s0 - ez0 * s1 ) / ( s0 - s1 );
154 this->computeDamage(tempDamage, Dn, kappa, epsilon1, epsilon2, ceps, seps, ez2);
155 double s2 = this->computeDimensionlessOutOfPlaneStress(eps, ez2, tempDamage);
156 ez0 = ez1;
157 s0 = s1;
158 ez1 = ez2;
159 s1 = s2;
160 }
161 }
162 }
163 status->setTempStrainZ(ez1);
164 status->setTempDamage(tempDamage);
165 double equivStrain = sqrt( equivStrainInPlane * equivStrainInPlane + macbra(ez1) * macbra(ez1) );
166 if ( equivStrain > kappa ) {
167 status->setTempKappa(equivStrain);
168 } else {
169 status->setTempKappa(kappa);
170 }
171 // formulae (93)-(100)
172 FloatArray answer;
173 computeInplaneStress(answer, eps, ez1, tempDamage);
174 //this->correctBigValues(stressTensor); // ???
175 status->setTempDamage(tempDamage);
176 status->letTempStrainVectorBe(totalStrain);
177 status->letTempStressVectorBe(answer);
178#ifdef keep_track_of_dissipated_energy
179 status->computeWork(gp);
180#endif
181 return answer;
182}
183
184void
185AnisotropicDamageMaterial :: computePrincValDir2D(double &D1, double &D2, double &c, double &s, double Dx, double Dy, double Dxy) const
186//
187// computes the principal values and directions of a symmetric second-order tensor in 2D
188// input: Dx, Dy, Dxy ... components of the tensor wrt global coordinates
189// output: D1, D2 ... ordered principal values, D1>=D2
190// output: c, s ... components of the unit principal vector associated with D1
191// (cosine and sine of the angle between the major principal direction and the global x-axis)
192{
193 // based on formulae (88)-(89) from the report by Jirasek & Suarez, 25 April 2014
194 double aux1 = ( Dx + Dy ) / 2.;
195 double aux2 = ( Dx - Dy ) / 2.;
196 double aux3 = sqrt(aux2 * aux2 + Dxy * Dxy);
197 D1 = aux1 + aux3;
198 D2 = aux1 - aux3;
199 // formulae (90)-(92) and the two cases preceding them
200 c = 1.;
201 s = 0.; // cases 1 and 2a
202 if ( Dxy != 0. ) { // case 3
203 double t = ( D1 - Dx ) / Dxy;
204 c = 1. / sqrt(1. + t * t);
205 s = c * t;
206 } else if ( Dx < Dy ) { // case 2b
207 c = 0.;
208 s = 1.;
209 }
210 return;
211}
212
213bool
214AnisotropicDamageMaterial :: checkPrincVal2D(double Dx, double Dy, double Dxy) const
215//
216// checks whether both eigenvalues of a symmetric second-order tensor in 2D are <= 1
217//
218{
219 if ( Dx + Dy > 2. ) {
220 return false;
221 }
222 if ( Dx + Dy > 1. + Dx * Dy - Dxy * Dxy ) {
223 return false;
224 }
225 return true;
226}
227
228void
229AnisotropicDamageMaterial :: computeDamage(FloatMatrix &tempDamage, const FloatMatrix &damage, double kappa, double eps1, double eps2, double ceps, double seps, double epsZ) const
230//
231// evaluates the final damage "tempDamage" from the initial damage "damage", initial history variable "kappa"
232// and the final in-plane strain
233// (given by principal values eps1>=eps2, components of first principal direction ceps and seps)
234// and out-of-plane strain
235//
236{
237 // set final damage to initial damage
238 tempDamage = damage;
239
240 // evaluate final equivalent strain using Mazars' definition
241 double tempEpsEq = 0.;
242 if ( eps1 > 0. ) {
243 tempEpsEq += eps1 * eps1;
244 if ( eps2 > 0. ) {
245 tempEpsEq += eps2 * eps2;
246 }
247 }
248 if ( epsZ > 0. ) {
249 tempEpsEq += epsZ * epsZ;
250 }
251 tempEpsEq = sqrt(tempEpsEq);
252
253 // check for damage growth
254 if ( tempEpsEq <= kappa ) { // no damage growth
255 return;
256 }
257
258 // apply incremental damage evolution law
259 // formula (75)
260 double deltaLambda = ( computeTraceD(tempEpsEq) - computeTraceD(kappa) ) / tempEpsEq / tempEpsEq;
261 double eps1p = macbra(eps1);
262 double eps2p = macbra(eps2);
263 double epsZp = macbra(epsZ);
264 double aux1 = deltaLambda * eps1p * eps1p;
265 double aux2 = deltaLambda * eps2p * eps2p;
266 // formula (76), with the square of positive strain expressed using the spectral decomposition
267 tempDamage.at(1, 1) += aux1 * ceps * ceps + aux2 * seps * seps;
268 tempDamage.at(1, 2) += ( aux1 - aux2 ) * ceps * seps;
269 tempDamage.at(2, 1) = tempDamage.at(1, 2);
270 tempDamage.at(2, 2) += aux1 * seps * seps + aux2 * ceps * ceps;
271 tempDamage.at(3, 3) += deltaLambda * epsZp * epsZp;
272
273 // treat the case when the out-of-plane damage exceeds 1
274 if ( tempDamage.at(3, 3) > 1. ) {
275 tempDamage.at(3, 3) = 1.;
276 }
277 // check whether in-plane principal damages do not exceed 1
278 if ( this->checkPrincVal2D( tempDamage.at(1, 1), tempDamage.at(2, 2), tempDamage.at(1, 2) ) ) {
279 return;
280 }
281 // treat the case when the in-plane damage exceeds 1
282 double D1, D2, cdam, sdam;
283 // compute principal values and direction of the initial damage
284 this->computePrincValDir2D( D1, D2, cdam, sdam, damage.at(1, 1), damage.at(2, 2), damage.at(1, 2) );
285 // compute damage increment components 11 and 22 in the principal damage coordinates
286 double dDx = tempDamage.at(1, 1) - damage.at(1, 1);
287 double dDy = tempDamage.at(2, 2) - damage.at(2, 2);
288 double dDxy = tempDamage.at(1, 2) - damage.at(1, 2);
289 double dD11 = cdam * cdam * dDx + sdam * sdam * dDy + 2. * cdam * sdam * dDxy;
290 double dD22 = sdam * sdam * dDx + cdam * cdam * dDy - 2. * cdam * sdam * dDxy;
291 // compute new principal damages, truncated to 1
292 D1 += dD11;
293 D2 += dD22;
294 if ( D1 > 1. ) {
295 D1 = 1.;
296 }
297 if ( D2 > 1. ) {
298 D2 = 1.;
299 }
300 // evaluate global components of the new damage tensor
301 tempDamage.at(1, 1) = cdam * cdam * D1 + sdam * sdam * D2;
302 tempDamage.at(2, 2) = sdam * sdam * D1 + cdam * cdam * D2;
303 tempDamage.at(1, 2) = cdam * sdam * ( D1 - D2 );
304 tempDamage.at(2, 1) = tempDamage.at(1, 2);
305}
306
307double
308AnisotropicDamageMaterial :: computeTraceD(double equivStrain) const
309{
310 double knu, aux;
311 double answer = 0.;
312 if ( equivStrain > this->kappa0 ) {
313 switch ( this->damageLawType ) {
314 case DLT_Desmorat1:
315 answer = ( equivStrain - this->kappa0 ) / ( this->kappaf - this->kappa0 );
316 break;
317 case DLT_Desmorat2:
318 answer = this->aA * ( atan(equivStrain / this->kappaf) - atan(this->kappa0 / this->kappaf) );
319 break;
320 case DLT_Linear:
321 knu = 2. * ( 1. + this->nu ) / 9.;
322 answer = this->kappaf * ( equivStrain - this->kappa0 ) / ( equivStrain * ( this->kappaf - this->kappa0 ) - knu * this->kappa0 * ( this->kappaf - equivStrain ) );
323 break;
324 case DLT_Exponential:
325 knu = 2. * ( 1. + this->nu ) / 9.;
326 aux = exp( -( equivStrain - this->kappa0 ) / ( this->kappaf - this->kappa0 ) );
327 answer = ( equivStrain - aux * this->kappa0 ) / ( equivStrain - aux * knu * this->kappa0 );
328 break;
329 default:
330 OOFEM_ERROR("Unknown type of damage law.\n");
331 }
332 }
333 return answer;
334}
335
336double
337AnisotropicDamageMaterial :: computeOutOfPlaneStrain(const FloatArray &inplaneStrain, const FloatMatrix &dam, bool tens_flag) const
338//
339// evaluate the out-of-plane strain from the condition of zero out-of-plane stress for the given damage
340// based on formula (79) from the report by Jirasek & Suarez, 25 April 2014
341//
342{
343 // evaluate principal damage values and directions (in-plane)
344 double D1, D2, cdam, sdam;
345 this->computePrincValDir2D( D1, D2, cdam, sdam, dam.at(1, 1), dam.at(2, 2), dam.at(1, 2) );
346
347 // evaluate normal strain components in directions of principal damage (in-plane)
348 // formulae (93)-(94)
349 double eps11 = cdam * cdam * inplaneStrain.at(1) + cdam *sdam *inplaneStrain.at(3) + sdam *sdam *inplaneStrain.at(2);
350 double eps22 = sdam * sdam * inplaneStrain.at(1) - cdam *sdam *inplaneStrain.at(3) + cdam *cdam *inplaneStrain.at(2);
351
352 // out-of-plane damage
353 double Dz = dam.at(3, 3);
354 // formula (4)
355 double Bv = 1.;
356 if ( tens_flag ) {
357 Bv = macbra(1. - D1 - D2 - Dz);
358 }
359 // evaluate auxiliary constants
360 // Q corresponds to K*B*Bv, Q1 corresponds to 2*G*Bz*B1 and Q2 corresponds to 2*G*Bz*B2
361 // but all of them divided by E and multiplied by 3*(1+nu)*(1-2*nu)
362 double Q = ( 3. - D1 - D2 - Dz ) * Bv * ( 1. + nu );
363 double Q1 = 3. * ( 1. - D1 ) * ( 1. - Dz ) * ( 1. - 2. * nu );
364 double Q2 = 3. * ( 1. - D2 ) * ( 1. - Dz ) * ( 1. - 2. * nu );
365
366 // evaluate out-of-plane strain which would give zero out-of-plane stress (at the given damage)
367 // formula (79)
368 double answer = ( ( Q1 - Q ) * eps11 + ( Q2 - Q ) * eps22 ) / ( Q + Q1 + Q2 );
369 return answer;
370}
371
372double
373AnisotropicDamageMaterial :: computeDimensionlessOutOfPlaneStress(const FloatArray &inplaneStrain, double epsZ, const FloatMatrix &dam) const
374//
375// evaluate the dimensionless out-of-plane stress (i.e., stress divided by E and multiplied by 3*B*(1+nu)*(1-2*nu))
376// for the given final in-plane strain, out-of-plane strain and damage (which is already updated for this strain)
377// based on formula (78) from the report by Jirasek & Suarez, 25 April 2014
378//
379{
380 // evaluate principal damage values and directions (in-plane)
381 double D1, D2, cdam, sdam;
382 this->computePrincValDir2D( D1, D2, cdam, sdam, dam.at(1, 1), dam.at(2, 2), dam.at(1, 2) );
383
384 // evaluate normal strain components in directions of principal damage (in-plane)
385 // formulae (93)-(94)
386 double eps11 = cdam * cdam * inplaneStrain.at(1) + cdam *sdam *inplaneStrain.at(3) + sdam *sdam *inplaneStrain.at(2);
387 double eps22 = sdam * sdam * inplaneStrain.at(1) - cdam *sdam *inplaneStrain.at(3) + cdam *cdam *inplaneStrain.at(2);
388
389 // out-of-plane damage
390 double Dz = dam.at(3, 3);
391 // formula (4)
392 double Bv = 1.;
393 double strainTrace = inplaneStrain.at(1) + inplaneStrain.at(2) + epsZ;
394 if ( strainTrace > 0. ) {
395 Bv = macbra(1. - D1 - D2 - Dz);
396 }
397 // evaluate auxiliary constants
398 // Q corresponds to K*B*Bv, Q1 corresponds to 2*G*Bz*B1 and Q2 corresponds to 2*G*Bz*B2
399 // but all of them divided by E and multiplied by 3*(1+nu)*(1-2*nu)
400 double Q = ( 3. - D1 - D2 - Dz ) * Bv * ( 1. + nu );
401 double Q1 = 3. * ( 1. - D1 ) * ( 1. - Dz ) * ( 1. - 2. * nu );
402 double Q2 = 3. * ( 1. - D2 ) * ( 1. - Dz ) * ( 1. - 2. * nu );
403 // formula (78) divided by E and multiplied by 3*B*(1+nu)*(1-2*nu)
404 double answer = ( Q - Q1 ) * eps11 + ( Q - Q2 ) * eps22 + ( Q + Q1 + Q2 ) * epsZ;
405 //printf("%g %g %g %g %g %g %g\n",epsZ,eps11,eps22,Q,Q1,Q2,answer);
406 return answer;
407}
408
409void
410AnisotropicDamageMaterial :: computeInplaneStress(FloatArray &inplaneStress, const FloatArray &inplaneStrain, double epsZ, const FloatMatrix &dam) const
411//
412// evaluate the in-plane stress components
413// for the given final in-plane strain, out-of-plane strain and damage (which is already updated for this strain)
414//
415{
416 // evaluate principal damage values and directions (in-plane)
417 double D1, D2, cdam, sdam;
418 this->computePrincValDir2D( D1, D2, cdam, sdam, dam.at(1, 1), dam.at(2, 2), dam.at(1, 2) );
419
420 // evaluate strain components with respect to the principal damage coordinates (in-plane)
421 // formulae (93)-(95), with both sides of (95) multiplied by factor 2
422 double eps11 = cdam * cdam * inplaneStrain.at(1) + cdam *sdam *inplaneStrain.at(3) + sdam *sdam *inplaneStrain.at(2);
423 double eps22 = sdam * sdam * inplaneStrain.at(1) - cdam *sdam *inplaneStrain.at(3) + cdam *cdam *inplaneStrain.at(2);
424 double gam12 = 2. * cdam * sdam * ( inplaneStrain.at(2) - inplaneStrain.at(1) ) + ( cdam * cdam - sdam * sdam ) * inplaneStrain.at(3);
425
426 /* OLD STYLE
427 * // evaluate in-plane effective stress
428 * double Eaux = E / ( ( 1 + nu ) * ( 1 - 2 * nu ) );
429 * double sig11eff = Eaux * ( ( 1. - nu ) * eps11 + nu * ( eps22 + epsZ ) );
430 * double sig22eff = Eaux * ( ( 1. - nu ) * eps22 + nu * ( eps11 + epsZ ) );
431 * double sigZeff = Eaux * ( ( 1. - nu ) * epsZ + nu * ( eps11 + eps22 ) );
432 * double G = E / ( 2. * ( 1 + nu ) );
433 * double sig12eff = G * gam12;
434 *
435 * // evaluate auxiliary constants
436 * double Dz = dam.at(3, 3);
437 * double Bv = 1.;
438 * double strainTrace = eps11 + eps22 + epsZ;
439 * double damTrace = D1 + D2 + Dz;
440 * if ( strainTrace > 0. ) {
441 * Bv = macbra(1. - damTrace);
442 * }
443 * double K = E / ( 3. * ( 1 - 2 * nu ) );
444 * double sigm = Bv * K * strainTrace;
445 * double Bsig = ( 1. - D1 ) * sig11eff + ( 1. - D2 ) * sig22eff + ( 1. - Dz ) * sigZeff;
446 * Bsig /= ( 3. - damTrace );
447 *
448 * // evaluate nominal stress from effective stress (in principal damage coordinate system)
449 * double sig11 = ( 1. - D1 ) * ( sig11eff - Bsig ) + sigm;
450 * double sig22 = ( 1. - D2 ) * ( sig22eff - Bsig ) + sigm;
451 * double sig12 = sqrt( ( 1. - D1 ) * ( 1. - D2 ) ) * sig12eff;
452 */
453
454 // evaluate auxiliary constants
455 double K = E / ( 3. * ( 1 - 2 * nu ) );
456 double G = E / ( 2. * ( 1 + nu ) );
457 double Bv = 1.;
458 double strainTrace = eps11 + eps22 + epsZ;
459 double damTrace = D1 + D2 + dam.at(3, 3);
460 if ( strainTrace > 0. ) {
461 Bv = macbra(1. - damTrace);
462 }
463 double B = 3. - damTrace;
464 double B1 = 1. - D1;
465 double B2 = 1. - D2;
466 double Bz = 1. - dam.at(3, 3);
467
468 // evaluate components of nominal in-plane stress wrt principal damage coordinates
469 // formulae (96)-(97)
470 double sigm = Bv * K * strainTrace;
471 double sig11 = 2. * G * B1 * ( B2 * ( eps11 - eps22 ) + Bz * ( eps11 - epsZ ) ) / B + sigm;
472 double sig22 = 2. * G * B2 * ( B1 * ( eps22 - eps11 ) + Bz * ( eps22 - epsZ ) ) / B + sigm;
473 double sig12 = G * sqrt(B1 * B2) * gam12;
474
475 // evaluate global components of nominal in-plane stress
476 // formulae (98)-(100)
477 inplaneStress.resize(3);
478 inplaneStress.at(1) = cdam * cdam * sig11 - 2. * cdam * sdam * sig12 + sdam * sdam * sig22;
479 inplaneStress.at(2) = sdam * sdam * sig11 + 2. * cdam * sdam * sig12 + cdam * cdam * sig22;
480 inplaneStress.at(3) = cdam * sdam * ( sig11 - sig22 ) + ( cdam * cdam - sdam * sdam ) * sig12;
481}
482
483//********************************************************
484// end of the plane stress implementation by Milan Jirasek
485//********************************************************
486
487void
488AnisotropicDamageMaterial :: giveRealStressVector(FloatArray &answer, GaussPoint *gp,
489 const FloatArray &totalStrain,
490 TimeStep *atTime) const
491//
492// returns real stress vector in 3d stress space of receiver
493// computed from the state at the beginning of the step and strain at the end of the step
494//
495{
496 //this->initGpForNewStep(gp);
497 this->initTempStatus(gp);
498 MaterialMode mode = gp->giveMaterialMode();
499 // subtract the stress-independent part of strains (e.g. due to temperature)
500 FloatArray reducedTotalStrainVector;
501 this->giveStressDependentPartOfStrainVector(reducedTotalStrainVector, gp, totalStrain, atTime, VM_Total);
502
503 // evaluate stress under general triaxial stress conditions
504 AnisotropicDamageMaterialStatus *status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
506 FloatMatrix de, tempDamage;
507 double equivStrain, kappa = 0.0, tempKappa = 0.0, traceTempD;
508 FloatArray eVals, effectiveStressVector, fullEffectiveStressVector, stressVector;
509 FloatMatrix eVecs, effectiveStressTensor, stressTensor, strainTensor;
510 FloatMatrix Dn = status->giveDamage();
511
512 this->computeEquivalentStrain(equivStrain, reducedTotalStrainVector, gp, atTime);
513 kappa = this->computeKappa(Dn);
514 if ( equivStrain <= kappa ) { // damage does not grow
515 tempDamage = Dn;
516 tempKappa = kappa;
517 } else {
518 this->computeDamageTensor(tempDamage, gp, reducedTotalStrainVector, equivStrain, atTime);
519 tempKappa = computeKappa(tempDamage);
520 }
521
522 if ( ( equivStrain <= kappa ) && ( kappa <= this->kappa0 ) ) { // elastic behavior
523 lmat.giveStiffnessMatrix(de, ElasticStiffness, gp, atTime);
524 effectiveStressVector.beProductOf(de, reducedTotalStrainVector);
525 answer = effectiveStressVector;
526 } else {
527 // StructuralMaterial :: giveFullSymVectorForm(fullEffectiveStressVector, effectiveStressVector, mode);
528 // effectiveStressTensor.beMatrixForm(fullEffectiveStressVector);
529 strainTensor.resize(3, 3);
530 strainTensor.zero();
531 if ( mode == _PlaneStress ) {
532 /*
533 * this->computePlaneStressStrain(strainTensor, tempDamage, reducedTotalStrainVector, gp, atTime);
534 * effectiveStressTensor.resize(3, 3);
535 * effectiveStressTensor.zero();
536 * double aux;
537 * aux = E / ( ( 1 + nu ) * ( 1 - 2 * nu ) );
538 * effectiveStressTensor.at(1, 1) = aux * ( ( 1 - nu ) * strainTensor.at(1, 1) + nu * ( strainTensor.at(2, 2) + strainTensor.at(3, 3) ) );
539 * effectiveStressTensor.at(2, 2) = aux * ( ( 1 - nu ) * strainTensor.at(2, 2) + nu * ( strainTensor.at(1, 1) + strainTensor.at(3, 3) ) );
540 * effectiveStressTensor.at(3, 3) = aux * ( ( 1 - nu ) * strainTensor.at(3, 3) + nu * ( strainTensor.at(1, 1) + strainTensor.at(2, 2) ) );
541 * effectiveStressTensor.at(1, 2) = ( E / ( 1 + nu ) ) * strainTensor.at(1, 2);
542 * effectiveStressTensor.at(2, 1) = effectiveStressTensor.at(1, 2);
543 */
544 } else {
545 strainTensor.at(1, 1) = reducedTotalStrainVector.at(1);
546 strainTensor.at(2, 2) = reducedTotalStrainVector.at(2);
547 strainTensor.at(3, 3) = reducedTotalStrainVector.at(3);
548 strainTensor.at(2, 3) = strainTensor.at(3, 2) = reducedTotalStrainVector.at(4) / 2.0;
549 strainTensor.at(1, 3) = strainTensor.at(3, 1) = reducedTotalStrainVector.at(5) / 2.0;
550 strainTensor.at(1, 2) = strainTensor.at(2, 1) = reducedTotalStrainVector.at(6) / 2.0;
551 lmat.giveStiffnessMatrix(de, ElasticStiffness, gp, atTime);
552 effectiveStressVector.beProductOf(de, reducedTotalStrainVector);
553 StructuralMaterial :: giveFullSymVectorForm(fullEffectiveStressVector, effectiveStressVector, mode);
554 effectiveStressTensor.beMatrixForm(fullEffectiveStressVector);
555 }
556 /* lmat.giveStiffnessMatrix(de, ElasticStiffness, gp, atTime);
557 * effectiveStressVector.beProductOf(de, reducedTotalStrainVector);
558 * StructuralMaterial :: giveFullSymVectorForm(fullEffectiveStressVector, effectiveStressVector, mode);
559 * effectiveStressTensor.beMatrixForm(fullEffectiveStressVector);*/
560 //traceTempD=tempDamage.at(1,1)+tempDamage.at(2,2)+tempDamage.at(3,3);
561 double effectiveStressTrace = effectiveStressTensor.at(1, 1) + effectiveStressTensor.at(2, 2) + effectiveStressTensor.at(3, 3);
562 FloatMatrix Part1, Part2, Part3;
563 // First term of the equation 53 of the reference paper****************************************************************************************
564 FloatMatrix AuxMatrix;
565 // Compute (1-D) (called ImD here)
566 FloatMatrix ImD, sqrtImD;
567 ImD.resize(3, 3);
568 ImD.zero();
569 ImD.at(1, 1) = ImD.at(2, 2) = ImD.at(3, 3) = 1;
570 ImD.subtract(tempDamage);
571 // Compute the square root of (1-D), needed in the equation 53 of the reference paper
572 //int checker1 = this->checkSymmetry(ImD);
573 ImD.jaco_(eVals, eVecs, 40);
574 sqrtImD.resize(3, 3);
575 sqrtImD.zero();
576 for ( int i = 1; i <= 3; i++ ) {
577 for ( int j = 1; j <= 3; j++ ) {
578 for ( int k = 1; k <= 3; k++ ) {
579 if ( eVals.at(k) < 0.0 ) {
580 eVals.at(k) = 0.0;
581 }
582 }
583 sqrtImD.at(i, j) = sqrt( eVals.at(1) ) * eVecs.at(i, 1) * eVecs.at(j, 1) + sqrt( eVals.at(2) ) * eVecs.at(i, 2) * eVecs.at(j, 2) + sqrt( eVals.at(3) ) * eVecs.at(i, 3) * eVecs.at(j, 3);
584 }
585 }
586
587 AuxMatrix.beProductOf(effectiveStressTensor, sqrtImD);
588 //stressTensor.beProductOf(sqrtImD,AuxMatrix);
589 Part1.beProductOf(sqrtImD, AuxMatrix);
590
591 // Second term of the equation 53 of the reference paper*****************************************************************************************
593 // Correct the trace of the damage tensor if necessary (see section 8.1 of the reference paper)
594 traceTempD = computeTraceD(tempDamage, strainTensor, gp);
595 double scalar = 0;
596 AuxMatrix.zero();
597 for ( int i = 1; i <= 3; i++ ) {
598 for ( int j = 1; j <= 3; j++ ) {
599 scalar += ImD.at(i, j) * effectiveStressTensor.at(i, j);
600 }
601 }
602 if ( ( 3. - traceTempD ) < 0.0001 ) {
603 scalar = 0.0;
604 } else {
605 scalar = scalar / ( 3. - traceTempD );
606 }
607 //scalar /= 3.-traceTempD;
608 AuxMatrix = ImD;
609 AuxMatrix.times(scalar);
610 Part2 = ImD;
611 Part2.times(scalar);
612
613
614 //stressTensor.subtract(AuxMatrix);
615 // Third term of the equation 53 of the reference paper********************************************************************************************
616 AuxMatrix.zero();
617 AuxMatrix.at(1, 1) = AuxMatrix.at(2, 2) = AuxMatrix.at(3, 3) = 1. / 3.;
618 if ( effectiveStressTrace > 0 ) {
619 AuxMatrix.times( ( 1 - traceTempD ) * effectiveStressTrace );
620 } else {
621 AuxMatrix.times(effectiveStressTrace);
622 }
623 Part3 = AuxMatrix;
624 //stressTensor.add(AuxMatrix);
625 stressTensor = Part1;
626 stressTensor.subtract(Part2);
627 stressTensor.add(Part3);
628 double factor;
629 factor = computeCorrectionFactor(tempDamage, strainTensor, gp);
630 stressTensor.times(factor);
631 stressVector.beSymVectorForm(stressTensor);
632 StructuralMaterial :: giveReducedSymVectorForm(answer, stressVector, mode);
633 }
634
635 // update gp
636 this->correctBigValues(stressTensor);
637 // int checker20 = this->checkSymmetry(tempDamage);
638 // int checker21 = this->checkSymmetry(stressTensor);
639 // int checker22 = this->checkSymmetry(strainTensor);
640 status->letTempStrainVectorBe(totalStrain);
641 status->letTempStressVectorBe(answer);
642 status->setTempDamage(tempDamage);
643 status->setTempKappa(tempKappa);
644#ifdef keep_track_of_dissipated_energy
645 status->computeWork(gp);
646#endif
647}
648
649void
650AnisotropicDamageMaterial :: computeEquivalentStrain(double &kappa, const FloatArray &strain, GaussPoint *gp, TimeStep *atTime) const
651{
652 auto status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
653
654 if ( strain.isEmpty() ) {
655 kappa = 0.;
656 return;
657 }
658
659
660 if ( this->equivStrainType == EST_Mazars ) {
661 double posNorm = 0.0;
662 FloatArray principalStrains, fullstrain;
663
664 StructuralMaterial :: giveFullSymVectorForm( fullstrain, strain, gp->giveMaterialMode() );
665
666 // if plane stress mode -> compute strain in z-direction from condition of zero stress in corresponding direction
667 if ( gp->giveMaterialMode() == _PlaneStress ) {
668 // fullstrain.at(3) = -nu * ( fullstrain.at(1) + fullstrain.at(2) ) / ( 1. - nu );
669 fullstrain.at(3) = status->giveTempStrainZ();
670 } else if ( gp->giveMaterialMode() == _1dMat ) {
671 fullstrain.at(2) = -nu *fullstrain.at(1);
672 fullstrain.at(3) = -nu *fullstrain.at(1);
673 }
674
675 this->computePrincipalValues(principalStrains, fullstrain, principal_strain);
676
677 for ( int i = 1; i <= 3; i++ ) {
678 if ( principalStrains.at(i) > 0.0 ) {
679 posNorm += principalStrains.at(i) * principalStrains.at(i);
680 }
681 }
682
683 kappa = sqrt(posNorm);
684 } else {
685 OOFEM_ERROR("computeEquivalentStrain: unknown EquivStrainType");
686 }
687 /*
688 * if ( this->equivStrainType == EST_Mazars ) {
689 * double posNorm = 0.0;
690 * FloatArray principalStrains, fullstrain;
691 *
692 * StructuralMaterial :: giveFullSymVectorForm( fullstrain, strain, gp->giveMaterialMode() );
693 *
694 * // if plane stress mode -> compute strain in z-direction from condition of zero stress in corresponding direction
695 * if ( gp->giveMaterialMode() == _PlaneStress ) {
696 * fullstrain.at(3) = -nu * ( fullstrain.at(1) + fullstrain.at(2) ) / ( 1. - nu );
697 * } else if ( gp->giveMaterialMode() == _1dMat ) {
698 * fullstrain.at(2) = -nu *fullstrain.at(1);
699 * fullstrain.at(3) = -nu *fullstrain.at(1);
700 * }
701 *
702 * this->computePrincipalValues(principalStrains, fullstrain, principal_strain);
703 *
704 * for ( int i = 1; i <= 3; i++ ) {
705 * if ( principalStrains.at(i) > 0.0 ) {
706 * posNorm += principalStrains.at(i) * principalStrains.at(i);
707 * }
708 * }
709 *
710 * kappa = sqrt(posNorm);
711 * } else if ( ( this->equivStrainType == EST_Rankine_Smooth ) || ( this->equivStrainType == EST_Rankine_Standard ) ) {
712 * // EST_Rankine equiv strain measure
713 * double sum = 0.;
714 * FloatArray stress, fullStress, principalStress;
715 * FloatMatrix de;
716 *
717 * lmat.giveStiffnessMatrix(de, SecantStiffness, gp, atTime);
718 * stress.beProductOf(de, strain);
719 * StructuralMaterial :: giveFullSymVectorForm( fullStress, stress, gp->giveMaterialMode() );
720 * this->computePrincipalValues(principalStress, fullStress, principal_stress);
721 * for ( int i = 1; i <= 3; i++ ) {
722 * if ( principalStress.at(i) > 0.0 ) {
723 * if ( this->equivStrainType == EST_Rankine_Smooth ) {
724 * sum += principalStress.at(i) * principalStress.at(i);
725 * } else if ( sum < principalStress.at(i) ) {
726 * sum = principalStress.at(i);
727 * }
728 * } else if ( sum < principalStress.at(i) ) {
729 * sum = principalStress.at(i);
730 * }
731 * }
732 *
733 * if ( this->equivStrainType == EST_Rankine_Smooth ) {
734 * sum = sqrt(sum);
735 * }
736 *
737 * kappa = sum / lmat.give('E', gp);
738 * } else if ( ( this->equivStrainType == EST_ElasticEnergy ) || ( this->equivStrainType == EST_ElasticEnergyPositiveStress ) || ( this->equivStrainType == EST_ElasticEnergyPositiveStrain ) ) {
739 * // equivalent strain expressions based on elastic energy
740 * FloatMatrix de;
741 * FloatArray stress;
742 * double sum;
743 *
744 * lmat.giveStiffnessMatrix(de, SecantStiffness, gp, atTime);
745 * if ( this->equivStrainType == EST_ElasticEnergy ) {
746 * // standard elastic energy
747 * stress.beProductOf(de, strain);
748 * sum = strain.dotProduct(stress);
749 * } else if ( this->equivStrainType == EST_ElasticEnergyPositiveStress ) {
750 * // elastic energy corresponding to positive part of stress
751 * FloatArray fullStress, principalStress;
752 * StructuralMaterial :: giveFullSymVectorForm( fullStress, stress, gp->giveMaterialMode() );
753 * this->computePrincipalValues(principalStress, fullStress, principal_stress);
754 * // TO BE FINISHED
755 * sum = 0.;
756 * OOFEM_ERROR("Elastic energy corresponding to positive part of stress not finished");
757 * } else {
758 * // elastic energy corresponding to positive part of strain
759 * // TO BE DONE
760 * sum = 0.;
761 * OOFEM_ERROR("Elastic energy corresponding to positive part of strain not finished");
762 * }
763 *
764 * kappa = sqrt( sum / lmat.give('E', gp) );
765 * } else if ( this->equivStrainType == EST_Griffith ) {
766 * double sum = 0.;
767 * FloatArray stress, fullStress, principalStress;
768 * FloatMatrix de;
769 *
770 * lmat.giveStiffnessMatrix(de, SecantStiffness, gp, atTime);
771 * stress.beProductOf(de, strain);
772 * StructuralMaterial :: giveFullSymVectorForm( fullStress, stress, gp->giveMaterialMode() );
773 * this->computePrincipalValues(principalStress, fullStress, principal_stress);
774 * for ( int i = 1; i <= 3; i++ ) {
775 * if ( principalStress.at(i) > 0.0 && sum < principalStress.at(i) ) {
776 * sum = principalStress.at(i);
777 * }
778 * }
779 *
780 * //Use Griffith criterion if Rankine not applied
781 * if (sum == 0.){
782 * sum = -pow(principalStress.at(1)-principalStress.at(3),2.)/8./(principalStress.at(1)+principalStress.at(3));
783 * }
784 * sum = max(sum,0.);
785 * kappa = sum / lmat.give('E', gp);
786 * } else {
787 * OOFEM_ERROR("computeEquivalentStrain: unknown EquivStrainType");
788 * }
789 */
790}
791
792// Computes Kappa according to the first damage law proposed in reference paper.
793double
794AnisotropicDamageMaterial :: computeKappa(FloatMatrix damageTensor) const
795{
796 double trace = damageTensor.giveTrace();
797 double answer = ( this->kappaf - this->kappa0 ) * trace + this->kappa0;
798 return answer;
799}
800
801
802// Computes the percentage of deltaD that needs to be added so the first eigenvalue of (tempDamageTensor + alpha*deltaD) reaches Dc
803// To do this, the middle point algorithm is used.
804// @TODO: this algorithm is not particularly efficient and another algorithm could be implemented.
805double
806AnisotropicDamageMaterial :: obtainAlpha1(FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, double damageThreshold) const
807{
808 double alpha_a, alpha_b, newAlpha, eps, maxDamage, size;
809 FloatMatrix deltaD, positiveStrainTensorSquared, resultingDamageTensor, eVecs;
810 FloatArray eVals;
811 int cont;
812 cont = 1;
813 alpha_a = 0;
814 alpha_b = 1;
815 newAlpha = ( alpha_a + alpha_b ) / 2;
816 positiveStrainTensorSquared.beProductOf(positiveStrainTensor, positiveStrainTensor);
817 deltaD = positiveStrainTensorSquared;
818 deltaD.times(newAlpha * deltaLambda);
819 this->correctBigValues(deltaD);
820 resultingDamageTensor = tempDamageTensor;
821 resultingDamageTensor.add(deltaD);
822 //int checker2 = this->checkSymmetry(resultingDamageTensor);
823 resultingDamageTensor.jaco_(eVals, eVecs, 20);
824 size = eVals.giveSize();
825 maxDamage = eVals.at(1);
826 for ( int i = 2; i <= size; i++ ) {
827 if ( eVals.at(i) > maxDamage ) {
828 maxDamage = eVals.at(i);
829 }
830 }
831 eps = maxDamage - damageThreshold;
832 do {
833 if ( eps > 0.0 ) {
834 alpha_b = newAlpha;
835 newAlpha = ( alpha_a + alpha_b ) / 2;
836 deltaD = positiveStrainTensorSquared;
837 deltaD.times(newAlpha * deltaLambda);
838 this->correctBigValues(deltaD);
839 resultingDamageTensor = tempDamageTensor;
840 resultingDamageTensor.add(deltaD);
841 //int checker3 = this->checkSymmetry(resultingDamageTensor);
842 resultingDamageTensor.jaco_(eVals, eVecs, 20);
843 size = eVals.giveSize();
844 maxDamage = eVals.at(1);
845 for ( int i = 2; i <= size; i++ ) {
846 if ( eVals.at(i) > maxDamage ) {
847 maxDamage = eVals.at(i);
848 }
849 }
850 eps = maxDamage - damageThreshold;
851 cont = cont + 1;
852 if ( cont == 100 ) {
853 return newAlpha;
854 }
855 } else {
856 alpha_a = newAlpha;
857 newAlpha = ( alpha_a + alpha_b ) / 2;
858 deltaD = positiveStrainTensorSquared;
859 deltaD.times(newAlpha * deltaLambda);
860 this->correctBigValues(deltaD);
861 resultingDamageTensor = tempDamageTensor;
862 resultingDamageTensor.add(deltaD);
863 //int checker4 = this->checkSymmetry(resultingDamageTensor);
864 resultingDamageTensor.jaco_(eVals, eVecs, 20);
865 size = eVals.giveSize();
866 maxDamage = eVals.at(1);
867 for ( int i = 2; i <= size; i++ ) {
868 if ( eVals.at(i) > maxDamage ) {
869 maxDamage = eVals.at(i);
870 }
871 }
872 eps = maxDamage - damageThreshold;
873 cont = cont + 1;
874 if ( cont == 100 ) {
875 return newAlpha;
876 }
877 }
878 } while ( fabs(eps) > 1.0e-15 );
879 return newAlpha;
880}
881
882// Computes the percentage of deltaD that needs to be added so the second eigenvalue of (tempDamageTensor + alpha*deltaD) reaches Dc
883// To do this, the middle point algorithm is used.
884// @TODO: this algorithm is not particularly efficient and another algorithm could be implemented.
885double
886AnisotropicDamageMaterial :: obtainAlpha2(FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, FloatMatrix projPosStrainTensor, double damageThreshold) const
887{
888 double alpha_a, alpha_b, newAlpha, eps, maxDamage, size, minVal;
889 FloatMatrix deltaD, resultingDamageTensor, eVecs;
890 FloatArray eVals;
891 int cont;
892 maxDamage = 0.0;
893 cont = 1;
894 alpha_a = 0;
895 alpha_b = 1;
896 newAlpha = ( alpha_a + alpha_b ) / 2;
897 deltaD = projPosStrainTensor;
898 deltaD.times(newAlpha * deltaLambda);
899 this->correctBigValues(deltaD);
900 resultingDamageTensor = tempDamageTensor;
901 resultingDamageTensor.add(deltaD);
902 //int checker5 = this->checkSymmetry(resultingDamageTensor);
903 resultingDamageTensor.jaco_(eVals, eVecs, 40);
904 size = eVals.giveSize();
905 minVal = eVals.at(1);
906 for ( int i = 2; i <= size; i++ ) {
907 if ( eVals.at(i) < minVal ) {
908 minVal = eVals.at(i);
909 maxDamage = ( ( eVals.at(1) + eVals.at(2) + eVals.at(3) ) - eVals.at(i) ) / 2;
910 }
911 }
912 eps = maxDamage - damageThreshold;
913 do {
914 if ( eps > 0.0 ) {
915 alpha_b = newAlpha;
916 newAlpha = ( alpha_a + alpha_b ) / 2;
917 deltaD = projPosStrainTensor;
918 deltaD.times(newAlpha * deltaLambda);
919 this->correctBigValues(deltaD);
920 resultingDamageTensor = tempDamageTensor;
921 resultingDamageTensor.add(deltaD);
922 //int checker6 = this->checkSymmetry(resultingDamageTensor);
923 resultingDamageTensor.jaco_(eVals, eVecs, 40);
924 size = eVals.giveSize();
925 minVal = eVals.at(1);
926 for ( int i = 2; i <= size; i++ ) {
927 if ( eVals.at(i) < minVal ) {
928 minVal = eVals.at(i);
929 maxDamage = ( ( eVals.at(1) + eVals.at(2) + eVals.at(3) ) - eVals.at(i) ) / 2;
930 }
931 }
932 eps = maxDamage - damageThreshold;
933 cont = cont + 1;
934 if ( cont == 100 ) {
935 return newAlpha;
936 }
937 } else {
938 alpha_a = newAlpha;
939 newAlpha = ( alpha_a + alpha_b ) / 2;
940 deltaD = projPosStrainTensor;
941 deltaD.times(newAlpha * deltaLambda);
942 this->correctBigValues(deltaD);
943 resultingDamageTensor = tempDamageTensor;
944 resultingDamageTensor.add(deltaD);
945 //int checker7 = this->checkSymmetry(resultingDamageTensor);
946 resultingDamageTensor.jaco_(eVals, eVecs, 40);
947 size = eVals.giveSize();
948 minVal = eVals.at(1);
949 for ( int i = 2; i <= size; i++ ) {
950 if ( eVals.at(i) < minVal ) {
951 minVal = eVals.at(i);
952 maxDamage = ( ( eVals.at(1) + eVals.at(2) + eVals.at(3) ) - eVals.at(i) ) / 2;
953 }
954 }
955 eps = maxDamage - damageThreshold;
956 cont = cont + 1;
957 if ( cont == 100 ) {
958 return newAlpha;
959 }
960 }
961 } while ( fabs(eps) > 1.0e-15 );
962 return newAlpha;
963}
964
965// Computes the percentage of deltaD that needs to be added so the third eigenvalue of (tempDamageTensor + alpha*deltaD) reaches Dc
966// To do this, the middle point algorithm is used.
967// @TODO: this algorithm is not particularly efficient and another algorithm could be implemented.
968double
969AnisotropicDamageMaterial :: obtainAlpha3(FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, FloatArray vec3, double damageThreshold) const
970{
971 double alpha_a, alpha_b, newAlpha, eps, aux = 0;
972 FloatMatrix deltaD, positiveStrainTensorSquared, resultingDamageTensor, eVecs;
973 FloatArray eVals;
974 int cont;
975 cont = 1;
976 alpha_a = 0;
977 alpha_b = 1;
978 newAlpha = ( alpha_a + alpha_b ) / 2;
979 positiveStrainTensorSquared.beProductOf(positiveStrainTensor, positiveStrainTensor);
980 for ( int i = 1; i <= 3; i++ ) {
981 aux = aux + vec3.at(i) * ( positiveStrainTensorSquared.at(i, 1) * vec3.at(1) + positiveStrainTensorSquared.at(i, 2) * vec3.at(2) + positiveStrainTensorSquared.at(i, 3) * vec3.at(3) );
982 }
983 deltaD.beDyadicProductOf(vec3, vec3);
984 deltaD.times(newAlpha * deltaLambda * aux);
985 this->correctBigValues(deltaD);
986 resultingDamageTensor = tempDamageTensor;
987 resultingDamageTensor.add(deltaD);
988 //int checker8 = this->checkSymmetry(resultingDamageTensor);
989 resultingDamageTensor.jaco_(eVals, eVecs, 20);
990 eps = eVals.at(1) + eVals.at(2) + eVals.at(3) - 3 * damageThreshold;
991 do {
992 if ( eps > 0.0 ) {
993 alpha_b = newAlpha;
994 newAlpha = ( alpha_a + alpha_b ) / 2;
995 deltaD.beDyadicProductOf(vec3, vec3);
996 deltaD.times(newAlpha * deltaLambda * aux);
997 this->correctBigValues(deltaD);
998 resultingDamageTensor = tempDamageTensor;
999 resultingDamageTensor.add(deltaD);
1000 //int checker9 = this->checkSymmetry(resultingDamageTensor);
1001 resultingDamageTensor.jaco_(eVals, eVecs, 20);
1002 eps = eVals.at(1) + eVals.at(2) + eVals.at(3) - 3 * damageThreshold;
1003 cont = cont + 1;
1004 if ( cont == 100 ) {
1005 return newAlpha;
1006 }
1007 } else {
1008 alpha_a = newAlpha;
1009 newAlpha = ( alpha_a + alpha_b ) / 2;
1010 deltaD.beDyadicProductOf(vec3, vec3);
1011 deltaD.times(newAlpha * deltaLambda * aux);
1012 this->correctBigValues(deltaD);
1013 resultingDamageTensor = tempDamageTensor;
1014 resultingDamageTensor.add(deltaD);
1015 //int checker10 = this->checkSymmetry(resultingDamageTensor);
1016 resultingDamageTensor.jaco_(eVals, eVecs, 20);
1017 eps = eVals.at(1) + eVals.at(2) + eVals.at(3) - 3 * damageThreshold;
1018 cont = cont + 1;
1019 if ( cont == 100 ) {
1020 return newAlpha;
1021 }
1022 }
1023 } while ( fabs(eps) > 1.0e-15 );
1024 return newAlpha;
1025}
1026//To check symmetry: delete this function when everything works fine
1027double
1028AnisotropicDamageMaterial :: checkSymmetry(FloatMatrix matrix) const
1029{
1030 int a = 0;
1031 int nRows = matrix.giveNumberOfRows();
1032 for ( int i = 1; i <= nRows; i++ ) {
1033 for ( int j = 1; j <= nRows; j++ ) {
1034 if ( fabs( matrix.at(i, j) - matrix.at(j, i) ) < 1.e-6 ) {
1035 ;
1036 } else {
1037 a = 1;
1038 }
1039 }
1040 }
1041 if ( a == 1 ) {
1042 a = 1;
1043 }
1044 return a;
1045}
1046
1047void
1048AnisotropicDamageMaterial :: correctBigValues(FloatMatrix &matrix) const
1049{
1050 int nRows = matrix.giveNumberOfRows();
1051 for ( int i = 1; i <= nRows; i++ ) {
1052 for ( int j = 1; j <= nRows; j++ ) {
1053 if ( matrix.at(i, j) != matrix.at(j, i) ) {
1054 double Aux = ( matrix.at(i, j) + matrix.at(j, i) ) / 2.0;
1055 matrix.at(i, j) = Aux;
1056 matrix.at(j, i) = Aux;
1057 }
1058 }
1059 }
1060}
1061
1062double
1063AnisotropicDamageMaterial :: computeTraceD(FloatMatrix tempDamageTensor, FloatMatrix strainTensor, GaussPoint *gp) const
1064{
1065 auto status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
1066 int flag = status->giveFlag();
1067 //int tempFlag=status->giveTempFlag();
1068 double Dc = 1.00, trD = 0;
1069 // If flag = 0, the trace of the damage tensor has never been greater than 1 before
1070 if ( flag == 0 ) {
1071 if ( ( strainTensor.at(1, 1) + strainTensor.at(2, 2) + strainTensor.at(3, 3) ) < 0 ) { // Compression
1072 trD = tempDamageTensor.at(1, 1) + tempDamageTensor.at(2, 2) + tempDamageTensor.at(3, 3);
1073 if ( trD >= 1 ) {
1074 status->setTempFlag(1);
1075 } // The trace of the damage tensor is greater than 1 for the first time, then, flag turns into 1
1076 } else { // Tension
1077 if ( ( tempDamageTensor.at(1, 1) + tempDamageTensor.at(2, 2) + tempDamageTensor.at(3, 3) ) >= 1 ) {
1078 trD = Dc;
1079 } else {
1080 trD = tempDamageTensor.at(1, 1) + tempDamageTensor.at(2, 2) + tempDamageTensor.at(3, 3);
1081 }
1082 }
1083 }
1084 // If flag = 1, the trace of the damage tensor has become greater than 1 before
1085 if ( flag == 1 ) {
1086 if ( ( strainTensor.at(1, 1) + strainTensor.at(2, 2) + strainTensor.at(3, 3) ) < 0 ) { // Compression
1087 trD = tempDamageTensor.at(1, 1) + tempDamageTensor.at(2, 2) + tempDamageTensor.at(3, 3);
1088 } else {
1089 trD = Dc;
1090 } // Tension
1091 }
1092 return trD;
1093}
1094
1095double
1096AnisotropicDamageMaterial :: computeCorrectionFactor(FloatMatrix tempDamageTensor, FloatMatrix strainTensor, GaussPoint *gp) const
1097{
1098 // In the case that the material has experimented some damaged under compression, this must affect the material behaviour when it is
1099 // under tension in the future
1100 auto status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
1101 int tempFlag = status->giveFlag();
1102 double tempStoredFactor = status->giveStoredFactor();
1103 double Dc = 1.00, trD = 0;
1104 double factor = 0.0;
1105 trD = tempDamageTensor.at(1, 1) + tempDamageTensor.at(2, 2) + tempDamageTensor.at(3, 3);
1106 // If flag = 0, the trace of the damage tensor has never been greater than Dc under compression before
1107 if ( tempFlag == 0 ) {
1108 if ( ( strainTensor.at(1, 1) + strainTensor.at(2, 2) + strainTensor.at(3, 3) ) < 0 ) { // Compression
1109 factor = 1.0;
1110 if ( ( 1.0 - trD ) < tempStoredFactor ) {
1111 status->setTempStoredFactor(1.0 - trD);
1112 }
1113 if ( ( 1.0 - trD ) <= ( 1.0 - Dc ) ) {
1114 tempFlag = 1;
1115 status->setTempFlag(1); // The trace of the damage tensor is greater than Dc for the first time, then, flag turns into 1
1116 status->setTempStoredFactor(1. - Dc);
1117 }
1118 } else { // Tension
1119 factor = tempStoredFactor;
1120 }
1121 }
1122
1123 // If flag = 1, the trace of the damage tensor has become greater than 1 before
1124 if ( tempFlag == 1 ) {
1125 if ( ( strainTensor.at(1, 1) + strainTensor.at(2, 2) + strainTensor.at(3, 3) ) < 0 ) { // Compression
1126 factor = 1.0;
1127 } else {
1128 //{factor=status->giveStoredFactor();} // Tension
1129 factor = 1. - Dc;
1130 }
1131 }
1132 return factor;
1133}
1134
1136AnisotropicDamageMaterial :: give3dMaterialStiffnessMatrix(MatResponseMode mode,
1137 GaussPoint *gp,
1138 TimeStep *atTime) const
1139//
1140// Implementation of the 3D stiffness matrix, according to the equations 56 and 57 of the reference paper.
1141{
1142 auto status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
1143 if ( mode == ElasticStiffness ) {
1144 return linearElasticMaterial.give3dMaterialStiffnessMatrix(mode, gp, atTime);
1145 } else {
1146 const FloatArray &totalStrain = status->giveTempStrainVector();
1147 FloatArray reducedTotalStrainVector;
1148 this->giveStressDependentPartOfStrainVector(reducedTotalStrainVector, gp, totalStrain, atTime, VM_Total);
1149 //FloatArray totalStrain;
1150 // The strain vector is turned into a tensor; for that, the elements that are out of the diagonal
1151 // must be divided by 2
1152 FloatMatrix strainTensor(3,3);
1153 strainTensor.at(1, 1) = reducedTotalStrainVector.at(1);
1154 strainTensor.at(2, 2) = reducedTotalStrainVector.at(2);
1155 strainTensor.at(3, 3) = reducedTotalStrainVector.at(3);
1156 strainTensor.at(2, 3) = reducedTotalStrainVector.at(4) / 2.0;
1157 strainTensor.at(3, 2) = reducedTotalStrainVector.at(4) / 2.0;
1158 strainTensor.at(1, 3) = reducedTotalStrainVector.at(5) / 2.0;
1159 strainTensor.at(3, 1) = reducedTotalStrainVector.at(5) / 2.0;
1160 strainTensor.at(1, 2) = reducedTotalStrainVector.at(6) / 2.0;
1161 strainTensor.at(2, 1) = reducedTotalStrainVector.at(6) / 2.0;
1162 // The damage tensor is read
1163 FloatMatrix answer;
1164 FloatMatrix damageTensor = status->giveTempDamage();
1165 AnisotropicDamageMaterial :: computeSecantOperator(answer, strainTensor, damageTensor, gp);
1166 for ( int j = 4; j <= 6; j++ ) {
1167 for ( int i = 1; i <= 6; i++ ) {
1168 answer.at(i, j) = answer.at(i, j) / 2.0;
1169 }
1170 }
1171 return answer;
1172 }
1173}
1174
1175
1177AnisotropicDamageMaterial :: givePlaneStressStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *atTime) const
1178{
1179 auto status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
1180 if ( mode == ElasticStiffness ) {
1181 return linearElasticMaterial.givePlaneStressStiffMtrx(mode, gp, atTime);
1182 } else {
1183 FloatArray totalStrain = status->giveTempStrainVector();
1184 FloatArray reducedTotalStrainVector;
1185 this->giveStressDependentPartOfStrainVector(reducedTotalStrainVector, gp, totalStrain, atTime, VM_Total);
1186 FloatMatrix damageTensor, strainTensor;
1187 // The damage tensor is read
1188 damageTensor.resize(3, 3);
1189 damageTensor.zero();
1190 damageTensor = status->giveTempDamage();
1191 this->computePlaneStressStrain(strainTensor, damageTensor, reducedTotalStrainVector, gp, atTime);
1192 // The strain vector is turned into a tensor; for that, the elements that are out of the diagonal
1193 // must be divided by 2
1194 FloatMatrix secantOperator;
1195 secantOperator.resize(6, 6);
1196 AnisotropicDamageMaterial :: computeSecantOperator(secantOperator, strainTensor, damageTensor, gp);
1197 double C11, C12, C13, C16, C21, C22, C23, C26, C61, C62, C63, C66, q, r, s;
1198 C11 = secantOperator.at(1, 1);
1199 C12 = secantOperator.at(1, 2);
1200 C13 = secantOperator.at(1, 3);
1201 C16 = secantOperator.at(1, 6);
1202 C21 = secantOperator.at(2, 1);
1203 C22 = secantOperator.at(2, 2);
1204 C23 = secantOperator.at(2, 3);
1205 C26 = secantOperator.at(2, 6);
1206 C61 = secantOperator.at(6, 1);
1207 C62 = secantOperator.at(6, 2);
1208 C63 = secantOperator.at(6, 3);
1209 C66 = secantOperator.at(6, 6);
1210 //Change 14-march-2014
1211 FloatArray stressVector;
1212 FloatMatrix Dn;
1213 stressVector = status->giveTempStressVector();
1214 Dn = status->giveTempDamage();
1215 this->computePlaneStressStrain(strainTensor, Dn, reducedTotalStrainVector, gp, atTime);
1216 if ( ( stressVector.at(1) + stressVector.at(2) ) < 1.0e-5 ) {
1217 q = -nu / E;
1218 } else {
1219 q = strainTensor.at(3, 3) / ( stressVector.at(1) + stressVector.at(2) );
1220 }
1221 //q = strainTensor.at(3,3)/(stressVector.at(1)+stressVector.at(2));
1222 //q = -nu/E;
1223 r = 1. / ( 1. - C13 * q );
1224 s = 1. / ( 1. - q * C23 - C23 * q * q * r * C13 );
1225 FloatMatrixF<3,3> answer;
1226 answer.at(2, 1) = s * ( C21 + C11 * C23 * q * r );
1227 answer.at(2, 2) = s * ( C22 + C12 * C23 * q * r );
1228 answer.at(2, 3) = s * ( C26 + C16 * C23 * q * r ) * 1. / 2.;
1229 answer.at(1, 1) = r * ( C11 + C13 * q * answer.at(2, 1) );
1230 answer.at(1, 2) = r * ( C12 + C13 * q * answer.at(2, 2) );
1231 answer.at(1, 3) = r * ( C16 + C13 * q * answer.at(2, 3) ) * 1 / 2;
1232 answer.at(3, 1) = C61 + C63 * q * ( answer.at(1, 1) + answer.at(2, 1) );
1233 answer.at(3, 2) = C62 + C63 * q * ( answer.at(1, 2) + answer.at(2, 2) );
1234 answer.at(3, 3) = ( C66 + C63 * q * ( answer.at(1, 3) + answer.at(2, 3) ) ) * 1. / 2.;
1235 return answer;
1236 }
1237}
1238
1239
1240FloatMatrixF<1,1> AnisotropicDamageMaterial :: give1dStressStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *atTime) const
1241{
1242 // Implementation of the 3D stiffness matrix
1243 if ( mode == ElasticStiffness ) {
1244 return linearElasticMaterial.give1dStressStiffMtrx(mode, gp, atTime);
1245 } else {
1246 auto status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
1247 FloatArray strain = status->giveTempStrainVector();
1248 if ( ( strain.at(1) + strain.at(2) + strain.at(3) ) > 0 ) {
1249 //@todo eq 56
1250 } else {
1251 //@todo eq 57
1252 }
1253 OOFEM_ERROR("not implemented");
1254 }
1255}
1256
1257void
1258AnisotropicDamageMaterial :: computePlaneStressStrain(FloatMatrix &answer, FloatMatrix damageTensor, FloatArray reducedTotalStrainVector, GaussPoint *gp,
1259 TimeStep *atTime) const
1260//
1261{
1262 FloatMatrix inPlaneStrain, B, eVecs;
1263 FloatArray eVals;
1264 double outOfPlaneStrain;
1265 B.resize(2, 2);
1266 B.at(1, 1) = 1. - damageTensor.at(1, 1);
1267 B.at(1, 2) = 0. - damageTensor.at(1, 2);
1268 B.at(2, 1) = 0. - damageTensor.at(2, 1);
1269 B.at(2, 2) = 1. - damageTensor.at(2, 2);
1270 //The eigenVectors for the change of base must be computed using the damageTensor, NOT the B matrix!!!
1271 //B.jaco_(eVals, eVecs, 40);
1272 FloatMatrix Auxiliar;
1273 Auxiliar.resize(2, 2);
1274 Auxiliar.at(1, 1) = damageTensor.at(1, 1);
1275 Auxiliar.at(2, 2) = damageTensor.at(2, 2);
1276 Auxiliar.at(1, 2) = damageTensor.at(1, 2);
1277 Auxiliar.at(2, 1) = damageTensor.at(2, 1);
1278
1279 Auxiliar.jaco_(eVals, eVecs, 40);
1280
1281 // Change of base of reducedTotalStrainVector from the canonical to the base formed by the eigenvectors of the damageTensor
1282 inPlaneStrain.resize(2, 2);
1283 inPlaneStrain.at(1, 1) = reducedTotalStrainVector.at(1);
1284 inPlaneStrain.at(2, 2) = reducedTotalStrainVector.at(2);
1285 inPlaneStrain.at(1, 2) = reducedTotalStrainVector.at(3) / 2.0;
1286 inPlaneStrain.at(2, 1) = reducedTotalStrainVector.at(3) / 2.0;
1287
1288 double term1, term2, term3, B1, B2, Bz, trD, h, epsilon11, epsilon22;
1289 B1 = 1.0 - eVals.at(1);
1290 B2 = 1.0 - eVals.at(2);
1291 Bz = 1. - damageTensor.at(3, 3);
1292 FloatArray vector1, vector2, auxVector;
1293 vector1.resize(2);
1294 vector2.resize(2);
1295 vector1.at(1) = eVecs.at(1, 1);
1296 vector1.at(2) = eVecs.at(2, 1);
1297 vector2.at(1) = eVecs.at(1, 2);
1298 vector2.at(2) = eVecs.at(2, 2);
1299 auxVector.beProductOf(inPlaneStrain, vector1);
1300 epsilon11 = vector1.at(1) * auxVector.at(1) + vector1.at(2) * auxVector.at(2);
1301 auxVector.beProductOf(inPlaneStrain, vector2);
1302 epsilon22 = vector2.at(1) * auxVector.at(1) + vector2.at(2) * auxVector.at(2);
1303 trD = damageTensor.at(1, 1) + damageTensor.at(2, 2) + damageTensor.at(3, 3);
1304 // Assuming a Tension state --> h = 1.0
1305 h = 1.0;
1306 term1 = 3. * Bz * B1 * ( 1. - 2. * nu ) - ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu );
1307 term2 = 3. * Bz * B2 * ( 1. - 2. * nu ) - ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu );
1308 term3 = 3. * Bz * ( 1. - 2. * nu ) * ( B1 + B2 ) + ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu );
1309 if ( Bz < 0.00001 ) {
1310 outOfPlaneStrain = -( epsilon11 + epsilon22 );
1311 } else {
1312 outOfPlaneStrain = ( term1 * epsilon11 + term2 * epsilon22 ) / term3;
1313 }
1314 /* if (outOfPlaneStrain != outOfPlaneStrain){
1315 * outOfPlaneStrain=-(epsilon11+epsilon22);
1316 * }*/
1317 double trStrain;
1318 trStrain = inPlaneStrain.at(1, 1) + inPlaneStrain.at(2, 2) + outOfPlaneStrain;
1319 // Check if actually under Tension, if not, recalculate term1, term2 and term3 with h = 0.0
1320 if ( trStrain < 0.0 ) {
1321 h = 0.0;
1322 term1 = 3. * Bz * B1 * ( 1. - 2. * nu ) - ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu );
1323 term2 = 3. * Bz * B2 * ( 1. - 2. * nu ) - ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu );
1324 term3 = 3. * Bz * ( 1. - 2. * nu ) * ( B1 + B2 ) + ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu );
1325 if ( Bz < 0.00001 ) {
1326 outOfPlaneStrain = -( epsilon11 + epsilon22 );
1327 } else {
1328 outOfPlaneStrain = ( term1 * epsilon11 + term2 * epsilon22 ) / term3;
1329 }
1330 }
1331 answer.resize(3, 3);
1332 answer.zero();
1333 answer.at(1, 1) = inPlaneStrain.at(1, 1);
1334 answer.at(1, 2) = inPlaneStrain.at(1, 2);
1335 answer.at(2, 1) = inPlaneStrain.at(2, 1);
1336 answer.at(2, 2) = inPlaneStrain.at(2, 2);
1337 answer.at(3, 3) = outOfPlaneStrain;
1338}
1339
1340void
1341AnisotropicDamageMaterial :: computePlaneStressSigmaZ(double &answer, FloatMatrix damageTensor, FloatArray reducedTotalStrainVector,
1342 double epsilonZ, GaussPoint *gp, TimeStep *atTime) const
1343//
1344{
1345 FloatMatrix Auxiliar, inPlaneStrain;
1346 FloatArray eVals;
1347 FloatMatrix eVecs;
1348 Auxiliar.resize(2, 2);
1349 Auxiliar.at(1, 1) = damageTensor.at(1, 1);
1350 Auxiliar.at(2, 2) = damageTensor.at(2, 2);
1351 Auxiliar.at(1, 2) = damageTensor.at(1, 2);
1352 Auxiliar.at(2, 1) = damageTensor.at(2, 1);
1353 Auxiliar.jaco_(eVals, eVecs, 40);
1354 inPlaneStrain.resize(2, 2);
1355 inPlaneStrain.at(1, 1) = reducedTotalStrainVector.at(1);
1356 inPlaneStrain.at(2, 2) = reducedTotalStrainVector.at(2);
1357 inPlaneStrain.at(1, 2) = reducedTotalStrainVector.at(3) / 2.0;
1358 inPlaneStrain.at(2, 1) = reducedTotalStrainVector.at(3) / 2.0;
1359 double term1, term2, termZ, B1, B2, Bz, trD, h, epsilon11, epsilon22;
1360 B1 = 1.0 - eVals.at(1);
1361 B2 = 1.0 - eVals.at(2);
1362 Bz = 1. - damageTensor.at(3, 3);
1363 FloatArray vector1, vector2, auxVector;
1364 vector1.resize(2);
1365 vector2.resize(2);
1366 vector1.at(1) = eVecs.at(1, 1);
1367 vector1.at(2) = eVecs.at(2, 1);
1368 vector2.at(1) = eVecs.at(1, 2);
1369 vector2.at(2) = eVecs.at(2, 2);
1370 auxVector.beProductOf(inPlaneStrain, vector1);
1371 epsilon11 = vector1.at(1) * auxVector.at(1) + vector1.at(2) * auxVector.at(2);
1372 auxVector.beProductOf(inPlaneStrain, vector2);
1373 epsilon22 = vector2.at(1) * auxVector.at(1) + vector2.at(2) * auxVector.at(2);
1374 trD = damageTensor.giveTrace();
1375 double Estar;
1376 Estar = E / ( ( 1. + nu ) * ( 1. - 2. * nu ) );
1377 if ( ( epsilon11 + epsilon22 + epsilonZ ) >= 0. ) {
1378 h = 1.;
1379 } else {
1380 h = 0.;
1381 }
1382 term1 = ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu ) - 3. * Bz * B1 * ( 1. - 2. * nu );
1383 term2 = ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu ) - 3. * Bz * B2 * ( 1. - 2. * nu );
1384 termZ = ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu ) + 3. * Bz * ( 1. - 2. * nu ) * ( B1 + B2 );
1385 // Finally the expression of sigmaZ is composed
1386 answer = ( Estar / ( 3. * ( B1 + B2 + Bz ) ) ) * ( epsilon11 * term1 + epsilon22 * term2 + epsilonZ * termZ );
1387
1388#if 0
1389 LinearElasticMaterial *lmat = this->giveLinearElasticMaterial();
1390 double B1, B2, Bz, eps11, eps22, Estar, term1, term2, termZ, trD, h;
1391 Estar=E / ((1. + nu)*(1. - 2. * nu));
1392 // Compute the eigenvalues of the in-plane damage tensor
1393 double eVal1, eVal2, aux1, aux2;
1394 aux1 = (damageTensor.at(1,1) + damageTensor.at(2,2))/2.0;
1395 aux2 = sqrt(pow((damageTensor.at(1,1) - damageTensor.at(2,2)) / 2. , 2.) + damageTensor.at(1,2) * damageTensor.at(2,1));
1396 eVal1 = aux1 + aux2 ;
1397 eVal2 = aux1 - aux2 ;
1398 B1 = 1. - eVal1;
1399 B2 = 1. - eVal2;
1400 Bz = 1. - damageTensor.at(3, 3);
1401 eps11 = reducedTotalStrainVector.at(1);
1402 eps22 = reducedTotalStrainVector.at(2);
1403 if ((eps11 + eps22 + epsZ)>=0.) {
1404 h = 1.;
1405 } else {
1406 h = 0.;
1407 }
1408 trD = damageTensor.giveTrace();
1409 term1 = ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu ) - 3. * Bz * B1 * ( 1. - 2. * nu ) ;
1410 term2 = ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu ) - 3. * Bz * B2 * ( 1. - 2. * nu ) ;
1411 termZ = ( 3. - trD ) * ( 1. - h * trD ) * ( 1. + nu ) + 3. * Bz * ( 1. - 2. * nu ) * ( B1 + B2 ) ;
1412 // Finally the expression of sigmaZ is composed
1413 answer = (Estar / (3. * (B1 + B2 + Bz))) * (eps11 * term1 + eps22 * term2 + epsZ * termZ);
1414 Estar = Estar + 0.0;
1415#endif
1416}
1417
1418void
1419AnisotropicDamageMaterial :: computeDamageTensor(FloatMatrix &answer, GaussPoint *gp,
1420 const FloatArray &reducedTotalStrainVector, double equivStrain,
1421 TimeStep *atTime) const
1422//
1423{
1424 //
1425 // returns real stress vector in 3d stress space of receiver according to
1426 // previous level of stress and current strain increment, the only way,
1427 // how to correctly update gp records
1428 //
1429 AnisotropicDamageMaterialStatus *status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
1430 double Dc = 1.00;
1431 double Kappa;
1432 FloatMatrix tempDamageTensor;
1433 //this->computeEquivalentStrain(equivStrain, reducedTotalStrainVector, gp, atTime);
1434 FloatMatrix Dn = status->giveDamage();
1435 Kappa = this->computeKappa(Dn);
1436 //damageTensor = status->giveDamage();
1437
1438 if ( equivStrain <= Kappa ) { // damage does not grow. Elastic behaviour
1439 answer.resize(3, 3);
1440 answer.zero();
1441 answer = Dn;
1442 } else { // damage grows
1443 Kappa = equivStrain;
1444 double deltaLambda;
1445 FloatArray eVals, fullStrainVector;
1446 FloatMatrix eVecs, strainTensor, positiveStrainTensor, positiveStrainTensorSquared, tempDamageTensor0;
1447 MaterialMode mode = gp->giveMaterialMode();
1448 // Compute square of positive part of strain tensor
1449 //1.- converts strain vector to full form;
1450 StructuralMaterial :: giveFullSymVectorForm(fullStrainVector, reducedTotalStrainVector, mode);
1451 // The strain vector is turned into a tensor; for that, the elements that are out of the diagonal
1452 // must be divided by 2
1453 strainTensor.resize(3, 3);
1454 strainTensor.zero();
1455 if ( mode == _PlaneStress ) {
1456 this->computePlaneStressStrain(strainTensor, Dn, reducedTotalStrainVector, gp, atTime);
1457 strainTensor.at(3, 3) = status->giveTempStrainZ();
1458 } else {
1459 strainTensor.at(1, 1) = reducedTotalStrainVector.at(1);
1460 strainTensor.at(2, 2) = reducedTotalStrainVector.at(2);
1461 strainTensor.at(3, 3) = reducedTotalStrainVector.at(3);
1462 strainTensor.at(2, 3) = reducedTotalStrainVector.at(4) / 2.0;
1463 strainTensor.at(3, 2) = reducedTotalStrainVector.at(4) / 2.0;
1464 strainTensor.at(1, 3) = reducedTotalStrainVector.at(5) / 2.0;
1465 strainTensor.at(3, 1) = reducedTotalStrainVector.at(5) / 2.0;
1466 strainTensor.at(1, 2) = reducedTotalStrainVector.at(6) / 2.0;
1467 strainTensor.at(2, 1) = reducedTotalStrainVector.at(6) / 2.0;
1468 }
1469 // computes polar decomposition and negative eigenvalues set to zero
1470 //int checker14 = this->checkSymmetry(strainTensor);
1471 strainTensor.jaco_(eVals, eVecs, 40);
1472 for ( int i = 1; i <= 3; i++ ) {
1473 if ( eVals.at(i) < 0 ) {
1474 eVals.at(i) = 0;
1475 }
1476 }
1477 // computes the positive part of the strain tensor
1478 positiveStrainTensor.resize(3, 3);
1479 for ( int i = 1; i <= 3; i++ ) {
1480 for ( int j = 1; j <= 3; j++ ) {
1481 positiveStrainTensor.at(i, j) = eVals.at(1) * eVecs.at(i, 1) * eVecs.at(j, 1) + eVals.at(2) * eVecs.at(i, 2) * eVecs.at(j, 2) + eVals.at(3) * eVecs.at(i, 3) * eVecs.at(j, 3);
1482 }
1483 }
1484 // computes the square of positiveStrainTensor
1485 positiveStrainTensorSquared.beProductOf(positiveStrainTensor, positiveStrainTensor);
1486 //compute delta Lambda
1487 //double traceD = damageTensor.at(1,1) + damageTensor.at(2,2) + damageTensor.at(3,3);
1488 double traceD = Dn.at(1, 1) + Dn.at(2, 2) + Dn.at(3, 3);
1489
1490 double traceTempD = this->computeTraceD(equivStrain);
1491 // equation 50 of the reference paper
1492 deltaLambda = ( traceTempD - traceD ) / equivStrain / equivStrain;
1493 //compute delta D: equation 48 of the reference paper
1494 FloatMatrix deltaD;
1495 deltaD = positiveStrainTensorSquared;
1496 deltaD.times(deltaLambda);
1497 this->correctBigValues(deltaD);
1498 // compute new damage tensor
1499 tempDamageTensor0 = Dn;
1500 tempDamageTensor0.add(deltaD); //damage tensor is equal to tempDamageTensor+deltaD
1501 // The following loop implements the algorithm for a case in which the maximum damage threshold is reached,
1502 // in such a case, the remaining damage is projected on the other two directions available and, finally,
1503 // if the damage threshold is reached in two of the three possible directions, the remaining damage is projected
1504 // on the remaining third direction. If the threshold is reached in the three possible directions, the damage
1505 // tensor remains unchanged in the future and with all their eigenvalues equal to the damage threshold Dc.
1506 // This part of the code is based on the section 8.2 of the reference paper.
1507 //int checker15 = this->checkSymmetry(tempDamageTensor0);
1508 tempDamageTensor0.jaco_(eVals, eVecs, 20);
1509 if ( ( eVals.at(1) > ( Dc * 1.001 ) ) || ( eVals.at(2) > ( Dc * 1.001 ) ) || ( eVals.at(3) > ( Dc * 1.001 ) ) ) {
1510 double alpha = 0, deltaLambda1 = 0, Aux1 = 0, Aux2 = 0, Aux3 = 0;
1511 FloatMatrix deltaD1(3, 3), positiveStrainTensorSquared(3, 3), tempDamageTensor1(3, 3), deltaD2(3, 3), N11(3, 3), N12(3, 3), N13(3, 3), N12sym(3, 3), N13sym(3, 3), projPosStrainTensor(3, 3);
1512 FloatArray auxVals(3), auxVec1(3), auxVec2(3), auxVec3(3);
1513
1514 // the percentage alpha of deltaD that needs to be added so the first eigenvalue reaches Dc is obtained
1515 alpha = obtainAlpha1(Dn, deltaLambda, positiveStrainTensor, Dc);
1516
1517 // deltaD1 is obtained --> deltaD1=alpha*deltaD
1518 positiveStrainTensorSquared.beProductOf(positiveStrainTensor, positiveStrainTensor);
1519 deltaD1 = positiveStrainTensorSquared;
1520 deltaD1.times(alpha * deltaLambda);
1521 this->correctBigValues(deltaD1);
1522 tempDamageTensor1 = Dn;
1523 tempDamageTensor1.add(deltaD1);
1524 // The following lines describe the process to apply the equation 64 of the reference paper. First, the
1525 // eigenvalues and eigenvectors of the damage tensor resulting at the moment when the threshold is reached
1526 // are obtained
1527 // (note: the equation 64 is not correctly written in the paper, it should be as implemented here:
1528 // D_dot = lambda_dot * [ <e>^2-(nI·<e>^2 nI)(nI x nI) - 2(nII·<e>^2 nI)(nI x nII)_sym - 2(nIII·<e>^2 nI)(nI x nIII)_sym ]
1529 //int checker16 = this->checkSymmetry(tempDamageTensor1);
1530
1531 tempDamageTensor1.jaco_(eVals, eVecs, 40);
1532 // The eigenvalues and eigenvectors are ordered, with the maximum eigenvalue being I, as its corresponding
1533 // eigenvector, and the other two being II and III. This is necessary so the equation 64 of the reference
1534 // paper can be applied
1535 if ( eVals.at(1) >= eVals.at(2) && eVals.at(1) >= eVals.at(3) ) {
1536 auxVals.at(1) = eVals.at(1);
1537 auxVec1.at(1) = eVecs.at(1, 1);
1538 auxVec1.at(2) = eVecs.at(2, 1);
1539 auxVec1.at(3) = eVecs.at(3, 1);
1540 auxVals.at(2) = eVals.at(2);
1541 auxVec2.at(1) = eVecs.at(1, 2);
1542 auxVec2.at(2) = eVecs.at(2, 2);
1543 auxVec2.at(3) = eVecs.at(3, 2);
1544 auxVals.at(3) = eVals.at(3);
1545 auxVec3.at(1) = eVecs.at(1, 3);
1546 auxVec3.at(2) = eVecs.at(2, 3);
1547 auxVec3.at(3) = eVecs.at(3, 3);
1548 } else if ( eVals.at(2) >= eVals.at(1) && eVals.at(2) >= eVals.at(3) ) {
1549 auxVals.at(1) = eVals.at(2);
1550 auxVec1.at(1) = eVecs.at(1, 2);
1551 auxVec1.at(2) = eVecs.at(2, 2);
1552 auxVec1.at(3) = eVecs.at(3, 2);
1553 auxVals.at(2) = eVals.at(1);
1554 auxVec2.at(1) = eVecs.at(1, 1);
1555 auxVec2.at(2) = eVecs.at(2, 1);
1556 auxVec2.at(3) = eVecs.at(3, 1);
1557 auxVals.at(3) = eVals.at(3);
1558 auxVec3.at(1) = eVecs.at(1, 3);
1559 auxVec3.at(2) = eVecs.at(2, 3);
1560 auxVec3.at(3) = eVecs.at(3, 3);
1561 } else {
1562 auxVals.at(1) = eVals.at(3);
1563 auxVec1.at(1) = eVecs.at(1, 3);
1564 auxVec1.at(2) = eVecs.at(2, 3);
1565 auxVec1.at(3) = eVecs.at(3, 3);
1566 auxVals.at(2) = eVals.at(2);
1567 auxVec2.at(1) = eVecs.at(1, 2);
1568 auxVec2.at(2) = eVecs.at(2, 2);
1569 auxVec2.at(3) = eVecs.at(3, 2);
1570 auxVals.at(3) = eVals.at(1);
1571 auxVec3.at(1) = eVecs.at(1, 1);
1572 auxVec3.at(2) = eVecs.at(2, 1);
1573 auxVec3.at(3) = eVecs.at(3, 1);
1574 }
1575
1576 // The symmetric part of the dyadic product of eigenvectors n1 and n2 is obtained
1577 N11.beDyadicProductOf(auxVec1, auxVec1);
1578 N12.beDyadicProductOf(auxVec1, auxVec2);
1579 for ( int i = 1; i <= 3; i++ ) {
1580 for ( int j = 1; j <= 3; j++ ) {
1581 N12sym.at(i, j) = 0.5 * ( N12.at(i, j) + N12.at(j, i) );
1582 }
1583 }
1584
1585 // The symmetric part of the dyadic product of eigenvectors n1 and n3 is obtained
1586 N13.beDyadicProductOf(auxVec1, auxVec3);
1587 for ( int i = 1; i <= 3; i++ ) {
1588 for ( int j = 1; j <= 3; j++ ) {
1589 N13sym.at(i, j) = 0.5 * ( N13.at(i, j) + N13.at(j, i) );
1590 }
1591 }
1592
1593 //The projected positive strain tensor is obtained
1594 for ( int i = 1; i <= 3; i++ ) {
1595 Aux1 = Aux1 + auxVec1.at(i) * ( positiveStrainTensorSquared.at(i, 1) * auxVec1.at(1) + positiveStrainTensorSquared.at(i, 2) * auxVec1.at(2) + positiveStrainTensorSquared.at(i, 3) * auxVec1.at(3) ); //==(n1*(<e>^2*n1)) (eq. 64 of the reference paper)
1596 Aux2 = Aux2 + auxVec2.at(i) * ( positiveStrainTensorSquared.at(i, 1) * auxVec1.at(1) + positiveStrainTensorSquared.at(i, 2) * auxVec1.at(2) + positiveStrainTensorSquared.at(i, 3) * auxVec1.at(3) ); //==(n2*(<e>^2*n1))_sym (eq. 64 of the reference paper)
1597 Aux3 = Aux3 + auxVec3.at(i) * ( positiveStrainTensorSquared.at(i, 1) * auxVec1.at(1) + positiveStrainTensorSquared.at(i, 2) * auxVec1.at(2) + positiveStrainTensorSquared.at(i, 3) * auxVec1.at(3) ); //==(n3*(<e>^2*n1))_sym (eq. 64 of the reference paper)
1598 }
1599 N11.times(Aux1);
1600 N12sym.times(2 * Aux2);
1601 N13sym.times(2 * Aux3);
1602
1603 // Finally, the expression between brackets in equation 64 is built up and called projPosStrainTensor
1604 projPosStrainTensor = positiveStrainTensorSquared;
1605 projPosStrainTensor.subtract(N11);
1606 projPosStrainTensor.subtract(N12sym);
1607 projPosStrainTensor.subtract(N13sym);
1608
1609 // The following loop avoids numerical problems in the case that the trace of projPosStrainTensor is very small
1610 if ( ( projPosStrainTensor.at(1, 1) + projPosStrainTensor.at(2, 2) + projPosStrainTensor.at(3, 3) ) < traceTempD * 1e-10 ) {
1611 deltaLambda1 = 0.;
1612 } else {
1613 deltaLambda1 = ( traceTempD - ( tempDamageTensor1.at(1, 1) + tempDamageTensor1.at(2, 2) + tempDamageTensor1.at(3, 3) ) ) / ( projPosStrainTensor.at(1, 1) + projPosStrainTensor.at(2, 2) + projPosStrainTensor.at(3, 3) );
1614 }
1615 //projPosStrainTensor.symmetrized();
1616 deltaD2 = projPosStrainTensor;
1617 deltaD2.times(deltaLambda1);
1618 this->correctBigValues(deltaD2);
1619 tempDamageTensor1.add(deltaD2); //damage tensor is equal to tempDamageTensor+deltaD1+deltaD2
1620 // The following loop checks if after the addition of D2, any other eigenvalue of the damage tensor
1621 // has reached the threshold. If it has, it repeats the process, but this time projecting the
1622 // remaining damage on the direction of the remaining eigenvector
1623 //int checker17 = this->checkSymmetry(tempDamageTensor1);
1624 tempDamageTensor1.jaco_(eVals, eVecs, 40);
1625 if ( ( eVals.at(1) > ( Dc * 1.001 ) ) || ( eVals.at(2) > ( Dc * 1.001 ) ) || ( eVals.at(3) > ( Dc * 1.001 ) ) ) {
1626 FloatMatrix deltaD3(3, 3), projPosStrainTensor_new(3, 3), tempDamageTensor2(3, 3), deltaD4(3, 3);
1627 FloatArray vec3(3);
1628 double alpha2 = 0, deltaLambda2 = 0, Aux4 = 0;
1629 // double val3=0;
1630 // Restoring the value of tempDamageTensor1 = tempDamageTensor + deltaD1
1631 tempDamageTensor1 = Dn;
1632 tempDamageTensor1.add(deltaD1);
1633 // the percentage alpha2 of deltaD2 that needs to be added to (tempDamageTensor+deltaD1) so the second eigenvalue
1634 // reaches Dc is obtained
1635 alpha2 = obtainAlpha2(tempDamageTensor1, deltaLambda1, positiveStrainTensor, projPosStrainTensor, Dc);
1636 deltaD3 = deltaD2;
1637 deltaD3.times(alpha2);
1638 tempDamageTensor2 = tempDamageTensor1;
1639 tempDamageTensor2.add(deltaD3);
1640 // The smallest eigenvalue is detected and its eigenvector is used to build the new projPosStrainTensor
1641 //int checker18 = this->checkSymmetry(tempDamageTensor2);
1642 tempDamageTensor2.jaco_(eVals, eVecs, 40);
1643 if ( eVals.at(1) <= eVals.at(2) && eVals.at(1) <= eVals.at(3) ) {
1644 //val3=eVals.at(1);
1645 vec3.at(1) = eVecs.at(1, 1);
1646 vec3.at(2) = eVecs.at(2, 1);
1647 vec3.at(3) = eVecs.at(3, 1);
1648 } else if ( eVals.at(2) <= eVals.at(1) && eVals.at(2) <= eVals.at(3) ) {
1649 //val3=eVals.at(2);
1650 vec3.at(1) = eVecs.at(1, 2);
1651 vec3.at(2) = eVecs.at(2, 2);
1652 vec3.at(3) = eVecs.at(3, 2);
1653 } else {
1654 //val3=eVals.at(3);
1655 vec3.at(1) = eVecs.at(1, 3);
1656 vec3.at(2) = eVecs.at(2, 3);
1657 vec3.at(3) = eVecs.at(3, 3);
1658 }
1659
1660 // The following loop computes nIII·<e>^2 nIII
1661 for ( int i = 1; i <= 3; i++ ) {
1662 Aux4 = Aux4 + vec3.at(i) * ( positiveStrainTensorSquared.at(i, 1) * vec3.at(1) + positiveStrainTensorSquared.at(i, 2) * vec3.at(2) + positiveStrainTensorSquared.at(i, 3) * vec3.at(3) );
1663 }
1664 projPosStrainTensor_new.beDyadicProductOf(vec3, vec3);
1665 //
1666 projPosStrainTensor_new.times(Aux4);
1667
1668 // The following loop avoids numerical problems in the case that the trace of projPosStrainTensor is very small
1669 if ( ( projPosStrainTensor_new.at(1, 1) + projPosStrainTensor_new.at(2, 2) + projPosStrainTensor_new.at(3, 3) ) < traceTempD * 1e-10 ) {
1670 deltaLambda2 = 0;
1671 } else {
1672 deltaLambda2 = ( traceTempD - ( tempDamageTensor2.at(1, 1) + tempDamageTensor2.at(2, 2) + tempDamageTensor2.at(3, 3) ) ) / ( projPosStrainTensor_new.at(1, 1) + projPosStrainTensor_new.at(2, 2) + projPosStrainTensor_new.at(3, 3) );
1673 }
1674 deltaD4 = projPosStrainTensor_new;
1675 deltaD4.times(deltaLambda2);
1676 tempDamageTensor2.add(deltaD4); //damage tensor is equal to tempDamageTensor+deltaD1+deltaD3+deltaD4
1677
1678 // The following loop checks if after the addition of D4, the remaining eigenvalue of the damage tensor
1679 // has reached the threshold. If it has, it computes a damage tensor with all its eigenvalues equal
1680 // to the damage threshold Dc
1681 //int checker19 = this->checkSymmetry(tempDamageTensor2);
1682 tempDamageTensor2.jaco_(eVals, eVecs, 40);
1683 if ( ( eVals.at(1) > ( Dc * 1.001 ) ) || ( eVals.at(2) > ( Dc * 1.001 ) ) || ( eVals.at(3) > ( Dc * 1.001 ) ) ) {
1684 double alpha3 = 0;
1685 FloatMatrix deltaD5, tempDamageTensor3;
1686 tempDamageTensor2 = Dn;
1687 tempDamageTensor2.add(deltaD1);
1688 tempDamageTensor2.add(deltaD3);
1689 alpha3 = obtainAlpha3(tempDamageTensor2, deltaLambda2, positiveStrainTensor, vec3, Dc);
1690 deltaD5 = deltaD4;
1691 deltaD5.times(alpha3);
1692 tempDamageTensor3 = tempDamageTensor2;
1693 tempDamageTensor3.add(deltaD5);
1694 tempDamageTensor = tempDamageTensor3;
1695 tempDamageTensor3.jaco_(eVals, eVecs, 40);
1696 if ( ( eVals.at(1) > ( Dc * 1.001 ) ) || ( eVals.at(2) > ( Dc * 1.001 ) ) || ( eVals.at(3) > ( Dc * 1.001 ) ) ) {
1697 tempDamageTensor3.zero();
1698 for ( int i = 1; i <= 3; i++ ) {
1699 if ( eVals.at(i) > Dc * 1.001 ) {
1700 eVals.at(i) = Dc;
1701 }
1702 }
1703 for ( int i = 1; i <= 3; i++ ) {
1704 for ( int j = 1; j <= 3; j++ ) {
1705 tempDamageTensor3.at(i, j) = eVals.at(1) * eVecs.at(i, 1) * eVecs.at(j, 1) + eVals.at(2) * eVecs.at(i, 2) * eVecs.at(j, 2) + eVals.at(3) * eVecs.at(i, 3) * eVecs.at(j, 3);
1706 }
1707 }
1708
1709 /*
1710 * tempDamageTensor3.zero();
1711 * tempDamageTensor3.at(1,1)=Dc;
1712 * tempDamageTensor3.at(2,2)=Dc;
1713 * tempDamageTensor3.at(3,3)=Dc;*/
1714 }
1715 tempDamageTensor = tempDamageTensor3;
1716 } else {
1717 tempDamageTensor = tempDamageTensor2;
1718 }
1719 } else {
1720 tempDamageTensor = tempDamageTensor1;
1721 }
1722 } else {
1723 tempDamageTensor = tempDamageTensor0;
1724 }
1725 answer = tempDamageTensor;
1726 }
1727}
1728
1729void
1730AnisotropicDamageMaterial :: computeSecantOperator(FloatMatrix &answer, FloatMatrix strainTensor, FloatMatrix damageTensor, GaussPoint *gp) const
1731//
1732// Implementation of the 3D stiffness matrix, according to the equations 56 and 57 of the reference paper.
1733{
1734 // AnisotropicDamageMaterialStatus *status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
1735 double G, K;
1736 double traceD, Aux;
1737 FloatMatrix ImD, sqrtImD, eVecs, Imatrix;
1738 FloatArray eVals;
1739 //MaterialMode mode = gp->giveMaterialMode();
1740 G = E / ( 2.0 * ( 1.0 + nu ) );
1741 K = E / ( 3.0 * ( 1.0 - 2.0 * nu ) );
1742 //Compute the trace of the damage tensor, correcting it if necessary (see section 8.1 of the reference paper)
1743 traceD = computeTraceD(damageTensor, strainTensor, gp);
1744
1745 if ( fabs(3. - traceD) < 0.001 ) {
1746 Aux = 0.0;
1747 } else {
1748 Aux = ( 1. / ( 3. - traceD ) );
1749 }
1750
1751 // compute square root of (1-D)
1752 ImD.resize(3, 3);
1753 ImD.zero();
1754 ImD.at(1, 1) = ImD.at(2, 2) = ImD.at(3, 3) = 1.0;
1755 ImD.subtract(damageTensor);
1756
1757 // computes square of positive part of strain tensor
1758 //int checker11=this->checkSymmetry(ImD);
1759 ImD.jaco_(eVals, eVecs, 40);
1760 sqrtImD.resize(3, 3);
1761 for ( int i = 1; i <= 3; i++ ) {
1762 for ( int j = 1; j <= 3; j++ ) {
1763 for ( int k = 1; k <= 3; k++ ) {
1764 if ( eVals.at(k) < 0.0 ) {
1765 eVals.at(k) = 0.0;
1766 }
1767 }
1768 sqrtImD.at(i, j) = sqrt( eVals.at(1) ) * eVecs.at(i, 1) * eVecs.at(j, 1) + sqrt( eVals.at(2) ) * eVecs.at(i, 2) * eVecs.at(j, 2) + sqrt( eVals.at(3) ) * eVecs.at(i, 3) * eVecs.at(j, 3);
1769 }
1770 }
1771 // To compute the expresions 56 and 57 of the reference paper, we need to work with fourth order tensors. To do this,
1772 // a structured called fourthOrderTensor is defined. This structure is composed by nine 3x3 FloatMatrix objects
1773 struct fourthOrderTensor
1774 {
1775 FloatMatrix Matrix_11kl;
1776 FloatMatrix Matrix_12kl;
1777 FloatMatrix Matrix_13kl;
1778 FloatMatrix Matrix_21kl;
1779 FloatMatrix Matrix_22kl;
1780 FloatMatrix Matrix_23kl;
1781 FloatMatrix Matrix_31kl;
1782 FloatMatrix Matrix_32kl;
1783 FloatMatrix Matrix_33kl;
1784 };
1785 // Four fourthOrderTensor structures are defined
1786 fourthOrderTensor Block1, Block2, Block3, secantOperator;
1787 Imatrix.resize(3, 3);
1788 Imatrix.zero();
1789 Imatrix.at(1, 1) = Imatrix.at(2, 2) = Imatrix.at(3, 3) = 1.0;
1790
1791 // The fourthOrderTensor structures are initialised
1792 Block1.Matrix_11kl.resize(3, 3);
1793 Block1.Matrix_12kl.resize(3, 3);
1794 Block1.Matrix_13kl.resize(3, 3);
1795 Block1.Matrix_21kl.resize(3, 3);
1796 Block1.Matrix_22kl.resize(3, 3);
1797 Block1.Matrix_23kl.resize(3, 3);
1798 Block1.Matrix_31kl.resize(3, 3);
1799 Block1.Matrix_32kl.resize(3, 3);
1800 Block1.Matrix_33kl.resize(3, 3);
1801 Block2.Matrix_11kl.resize(3, 3);
1802 Block2.Matrix_12kl.resize(3, 3);
1803 Block2.Matrix_13kl.resize(3, 3);
1804 Block2.Matrix_21kl.resize(3, 3);
1805 Block2.Matrix_22kl.resize(3, 3);
1806 Block2.Matrix_23kl.resize(3, 3);
1807 Block2.Matrix_31kl.resize(3, 3);
1808 Block2.Matrix_32kl.resize(3, 3);
1809 Block2.Matrix_33kl.resize(3, 3);
1810 Block3.Matrix_11kl.resize(3, 3);
1811 Block3.Matrix_12kl.resize(3, 3);
1812 Block3.Matrix_13kl.resize(3, 3);
1813 Block3.Matrix_21kl.resize(3, 3);
1814 Block3.Matrix_22kl.resize(3, 3);
1815 Block3.Matrix_23kl.resize(3, 3);
1816 Block3.Matrix_31kl.resize(3, 3);
1817 Block3.Matrix_32kl.resize(3, 3);
1818 Block3.Matrix_33kl.resize(3, 3);
1819 secantOperator.Matrix_11kl.resize(3, 3);
1820 secantOperator.Matrix_12kl.resize(3, 3);
1821 secantOperator.Matrix_13kl.resize(3, 3);
1822 secantOperator.Matrix_21kl.resize(3, 3);
1823 secantOperator.Matrix_22kl.resize(3, 3);
1824 secantOperator.Matrix_23kl.resize(3, 3);
1825 secantOperator.Matrix_31kl.resize(3, 3);
1826 secantOperator.Matrix_32kl.resize(3, 3);
1827 secantOperator.Matrix_33kl.resize(3, 3);
1828
1829 for ( int k = 1; k <= 3; k++ ) {
1830 for ( int l = 1; l <= 3; l++ ) {
1831 //The first block inside the brackets is obtained --> (1-D)^(1/2) _x_ (1-D)^(1/2)
1832 Block1.Matrix_11kl.at(k, l) = sqrtImD.at(1, k) * sqrtImD.at(1, l);
1833 Block1.Matrix_12kl.at(k, l) = sqrtImD.at(1, k) * sqrtImD.at(2, l);
1834 Block1.Matrix_13kl.at(k, l) = sqrtImD.at(1, k) * sqrtImD.at(3, l);
1835 Block1.Matrix_21kl.at(k, l) = sqrtImD.at(2, k) * sqrtImD.at(1, l);
1836 Block1.Matrix_22kl.at(k, l) = sqrtImD.at(2, k) * sqrtImD.at(2, l);
1837 Block1.Matrix_23kl.at(k, l) = sqrtImD.at(2, k) * sqrtImD.at(3, l);
1838 Block1.Matrix_31kl.at(k, l) = sqrtImD.at(3, k) * sqrtImD.at(1, l);
1839 Block1.Matrix_32kl.at(k, l) = sqrtImD.at(3, k) * sqrtImD.at(2, l);
1840 Block1.Matrix_33kl.at(k, l) = sqrtImD.at(3, k) * sqrtImD.at(3, l);
1841 //The second block inside the brackets is obtained --> ((1-D) x (1-D))/(3-trD)
1842 Block2.Matrix_11kl.at(k, l) = ImD.at(1, 1) * ImD.at(k, l) * Aux;
1843 Block2.Matrix_12kl.at(k, l) = ImD.at(1, 2) * ImD.at(k, l) * Aux;
1844 Block2.Matrix_13kl.at(k, l) = ImD.at(1, 3) * ImD.at(k, l) * Aux;
1845 Block2.Matrix_21kl.at(k, l) = ImD.at(2, 1) * ImD.at(k, l) * Aux;
1846 Block2.Matrix_22kl.at(k, l) = ImD.at(2, 2) * ImD.at(k, l) * Aux;
1847 Block2.Matrix_23kl.at(k, l) = ImD.at(2, 3) * ImD.at(k, l) * Aux;
1848 Block2.Matrix_31kl.at(k, l) = ImD.at(3, 1) * ImD.at(k, l) * Aux;
1849 Block2.Matrix_32kl.at(k, l) = ImD.at(3, 2) * ImD.at(k, l) * Aux;
1850 Block2.Matrix_33kl.at(k, l) = ImD.at(3, 3) * ImD.at(k, l) * Aux;
1851 //The crossed-product of two identity tensors is obtained --> 1 x 1
1852 Block3.Matrix_11kl.at(k, l) = Imatrix.at(1, 1) * Imatrix.at(k, l);
1853 Block3.Matrix_12kl.at(k, l) = Imatrix.at(1, 2) * Imatrix.at(k, l);
1854 Block3.Matrix_13kl.at(k, l) = Imatrix.at(1, 3) * Imatrix.at(k, l);
1855 Block3.Matrix_21kl.at(k, l) = Imatrix.at(2, 1) * Imatrix.at(k, l);
1856 Block3.Matrix_22kl.at(k, l) = Imatrix.at(2, 2) * Imatrix.at(k, l);
1857 Block3.Matrix_23kl.at(k, l) = Imatrix.at(2, 3) * Imatrix.at(k, l);
1858 Block3.Matrix_31kl.at(k, l) = Imatrix.at(3, 1) * Imatrix.at(k, l);
1859 Block3.Matrix_32kl.at(k, l) = Imatrix.at(3, 2) * Imatrix.at(k, l);
1860 Block3.Matrix_33kl.at(k, l) = Imatrix.at(3, 3) * Imatrix.at(k, l);
1861 }
1862 }
1863 // equation 56 of the reference paper
1864 if ( ( strainTensor.at(1, 1) + strainTensor.at(2, 2) + strainTensor.at(3, 3) ) > 0.0 ) {
1865 for ( int k = 1; k <= 3; k++ ) {
1866 for ( int l = 1; l <= 3; l++ ) {
1867 secantOperator.Matrix_11kl.at(k, l) = 2 * G * ( Block1.Matrix_11kl.at(k, l) - Block2.Matrix_11kl.at(k, l) ) + K * ( 1 - traceD ) * Block3.Matrix_11kl.at(k, l);
1868 secantOperator.Matrix_12kl.at(k, l) = 2 * G * ( Block1.Matrix_12kl.at(k, l) - Block2.Matrix_12kl.at(k, l) ) + K * ( 1 - traceD ) * Block3.Matrix_12kl.at(k, l);
1869 secantOperator.Matrix_13kl.at(k, l) = 2 * G * ( Block1.Matrix_13kl.at(k, l) - Block2.Matrix_13kl.at(k, l) ) + K * ( 1 - traceD ) * Block3.Matrix_13kl.at(k, l);
1870 secantOperator.Matrix_21kl.at(k, l) = 2 * G * ( Block1.Matrix_21kl.at(k, l) - Block2.Matrix_21kl.at(k, l) ) + K * ( 1 - traceD ) * Block3.Matrix_21kl.at(k, l);
1871 secantOperator.Matrix_22kl.at(k, l) = 2 * G * ( Block1.Matrix_22kl.at(k, l) - Block2.Matrix_22kl.at(k, l) ) + K * ( 1 - traceD ) * Block3.Matrix_22kl.at(k, l);
1872 secantOperator.Matrix_23kl.at(k, l) = 2 * G * ( Block1.Matrix_23kl.at(k, l) - Block2.Matrix_23kl.at(k, l) ) + K * ( 1 - traceD ) * Block3.Matrix_23kl.at(k, l);
1873 secantOperator.Matrix_31kl.at(k, l) = 2 * G * ( Block1.Matrix_31kl.at(k, l) - Block2.Matrix_31kl.at(k, l) ) + K * ( 1 - traceD ) * Block3.Matrix_31kl.at(k, l);
1874 secantOperator.Matrix_32kl.at(k, l) = 2 * G * ( Block1.Matrix_32kl.at(k, l) - Block2.Matrix_32kl.at(k, l) ) + K * ( 1 - traceD ) * Block3.Matrix_32kl.at(k, l);
1875 secantOperator.Matrix_33kl.at(k, l) = 2 * G * ( Block1.Matrix_33kl.at(k, l) - Block2.Matrix_33kl.at(k, l) ) + K * ( 1 - traceD ) * Block3.Matrix_33kl.at(k, l);
1876 }
1877 }
1878 // equation 57 of the reference paper
1879 } else {
1880 for ( int k = 1; k <= 3; k++ ) {
1881 for ( int l = 1; l <= 3; l++ ) {
1882 secantOperator.Matrix_11kl.at(k, l) = 2 * G * ( Block1.Matrix_11kl.at(k, l) - Block2.Matrix_11kl.at(k, l) ) + K *Block3.Matrix_11kl.at(k, l);
1883 secantOperator.Matrix_12kl.at(k, l) = 2 * G * ( Block1.Matrix_12kl.at(k, l) - Block2.Matrix_12kl.at(k, l) ) + K *Block3.Matrix_12kl.at(k, l);
1884 secantOperator.Matrix_13kl.at(k, l) = 2 * G * ( Block1.Matrix_13kl.at(k, l) - Block2.Matrix_13kl.at(k, l) ) + K *Block3.Matrix_13kl.at(k, l);
1885 secantOperator.Matrix_21kl.at(k, l) = 2 * G * ( Block1.Matrix_21kl.at(k, l) - Block2.Matrix_21kl.at(k, l) ) + K *Block3.Matrix_21kl.at(k, l);
1886 secantOperator.Matrix_22kl.at(k, l) = 2 * G * ( Block1.Matrix_22kl.at(k, l) - Block2.Matrix_22kl.at(k, l) ) + K *Block3.Matrix_22kl.at(k, l);
1887 secantOperator.Matrix_23kl.at(k, l) = 2 * G * ( Block1.Matrix_23kl.at(k, l) - Block2.Matrix_23kl.at(k, l) ) + K *Block3.Matrix_23kl.at(k, l);
1888 secantOperator.Matrix_31kl.at(k, l) = 2 * G * ( Block1.Matrix_31kl.at(k, l) - Block2.Matrix_31kl.at(k, l) ) + K *Block3.Matrix_31kl.at(k, l);
1889 secantOperator.Matrix_32kl.at(k, l) = 2 * G * ( Block1.Matrix_32kl.at(k, l) - Block2.Matrix_32kl.at(k, l) ) + K *Block3.Matrix_32kl.at(k, l);
1890 secantOperator.Matrix_33kl.at(k, l) = 2 * G * ( Block1.Matrix_33kl.at(k, l) - Block2.Matrix_33kl.at(k, l) ) + K *Block3.Matrix_33kl.at(k, l);
1891 }
1892 }
1893 }
1894 // The resulting material stiffness matrix is built
1895 answer.resize(6, 6);
1896 answer.at(1, 1) = secantOperator.Matrix_11kl.at(1, 1);
1897 answer.at(1, 2) = secantOperator.Matrix_11kl.at(2, 2);
1898 answer.at(1, 3) = secantOperator.Matrix_11kl.at(3, 3);
1899 answer.at(1, 4) = secantOperator.Matrix_11kl.at(2, 3);
1900 answer.at(1, 5) = secantOperator.Matrix_11kl.at(3, 1);
1901 answer.at(1, 6) = secantOperator.Matrix_11kl.at(1, 2);
1902 answer.at(2, 1) = secantOperator.Matrix_22kl.at(1, 1);
1903 answer.at(2, 2) = secantOperator.Matrix_22kl.at(2, 2);
1904 answer.at(2, 3) = secantOperator.Matrix_22kl.at(3, 3);
1905 answer.at(2, 4) = secantOperator.Matrix_22kl.at(2, 3);
1906 answer.at(2, 5) = secantOperator.Matrix_22kl.at(3, 1);
1907 answer.at(2, 6) = secantOperator.Matrix_22kl.at(1, 2);
1908 answer.at(3, 1) = secantOperator.Matrix_33kl.at(1, 1);
1909 answer.at(3, 2) = secantOperator.Matrix_33kl.at(2, 2);
1910 answer.at(3, 3) = secantOperator.Matrix_33kl.at(3, 3);
1911 answer.at(3, 4) = secantOperator.Matrix_33kl.at(2, 3);
1912 answer.at(3, 5) = secantOperator.Matrix_33kl.at(3, 1);
1913 answer.at(3, 6) = secantOperator.Matrix_33kl.at(1, 2);
1914 answer.at(4, 1) = secantOperator.Matrix_23kl.at(1, 1);
1915 answer.at(4, 2) = secantOperator.Matrix_23kl.at(2, 2);
1916 answer.at(4, 3) = secantOperator.Matrix_23kl.at(3, 3);
1917 answer.at(4, 4) = secantOperator.Matrix_23kl.at(2, 3);
1918 answer.at(4, 5) = secantOperator.Matrix_23kl.at(3, 1);
1919 answer.at(4, 6) = secantOperator.Matrix_23kl.at(1, 2);
1920 answer.at(5, 1) = secantOperator.Matrix_31kl.at(1, 1);
1921 answer.at(5, 2) = secantOperator.Matrix_31kl.at(2, 2);
1922 answer.at(5, 3) = secantOperator.Matrix_31kl.at(3, 3);
1923 answer.at(5, 4) = secantOperator.Matrix_31kl.at(2, 3);
1924 answer.at(5, 5) = secantOperator.Matrix_31kl.at(3, 1);
1925 answer.at(5, 6) = secantOperator.Matrix_31kl.at(1, 2);
1926 answer.at(6, 1) = secantOperator.Matrix_12kl.at(1, 1);
1927 answer.at(6, 2) = secantOperator.Matrix_12kl.at(2, 2);
1928 answer.at(6, 3) = secantOperator.Matrix_12kl.at(3, 3);
1929 answer.at(6, 4) = secantOperator.Matrix_12kl.at(2, 3);
1930 answer.at(6, 5) = secantOperator.Matrix_12kl.at(3, 1);
1931 answer.at(6, 6) = secantOperator.Matrix_12kl.at(1, 2);
1932}
1933
1934int
1935AnisotropicDamageMaterial :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *atTime)
1936{
1937 AnisotropicDamageMaterialStatus *status = static_cast< AnisotropicDamageMaterialStatus * >( this->giveStatus(gp) );
1938 if ( type == IST_DamageScalar ) { // returning the trace of the damage tensor
1939 answer.resize(1);
1940 answer.at(1) = status->giveDamage().at(1, 1) + status->giveDamage().at(2, 2) + status->giveDamage().at(3, 3);
1941 return 1;
1942 } else if ( type == IST_DamageTensor ) {
1943 answer.resize(6);
1944 answer.zero();
1945 answer.at(1) = status->giveDamage().at(1, 1);
1946 answer.at(2) = status->giveDamage().at(2, 2);
1947 answer.at(3) = status->giveDamage().at(3, 3);
1948 answer.at(4) = status->giveDamage().at(2, 3);
1949 answer.at(5) = status->giveDamage().at(1, 3);
1950 answer.at(6) = status->giveDamage().at(1, 2);
1951 return 1;
1952 } else if ( type == IST_PrincipalDamageTensor ) {
1953 //int checker12=this->checkSymmetry(status->giveDamage());
1954 FloatMatrix dam = status->giveDamage();
1955 FloatMatrix eVecs;
1956 dam.jaco_(answer, eVecs, 20);
1957 return 1;
1958 } else if ( type == IST_DamageTensorTemp ) {
1959 answer.resize(6);
1960 answer.zero();
1961 answer.at(1) = status->giveTempDamage().at(1, 1);
1962 answer.at(2) = status->giveTempDamage().at(2, 2);
1963 answer.at(3) = status->giveTempDamage().at(3, 3);
1964 answer.at(4) = status->giveTempDamage().at(2, 3);
1965 answer.at(5) = status->giveTempDamage().at(1, 3);
1966 answer.at(6) = status->giveTempDamage().at(1, 2);
1967 return 1;
1968 } else if ( type == IST_PrincipalDamageTempTensor ) {
1969 //int checker13=this->checkSymmetry(status->giveTempDamage());
1970 FloatMatrix dam = status->giveTempDamage();
1971 FloatMatrix eVecs;
1972 dam.jaco_(answer, eVecs, 20);
1973 return 1;
1974 } else if ( type == IST_MaxEquivalentStrainLevel ) {
1975 answer.resize(1);
1976 answer.at(1) = status->giveKappa();
1977 return 1;
1978
1979#ifdef keep_track_of_dissipated_energy
1980 } else if ( type == IST_StressWorkDensity ) {
1981 answer.resize(1);
1982 answer.at(1) = status->giveStressWork();
1983 return 1;
1984 } else if ( type == IST_DissWorkDensity ) {
1985 answer.resize(1);
1986 answer.at(1) = status->giveDissWork();
1987 } else if ( type == IST_FreeEnergyDensity ) {
1988 answer.resize(1);
1989 answer.at(1) = status->giveStressWork() - status->giveDissWork();
1990 return 1;
1991
1992#endif
1993 } else {
1994 return StructuralMaterial :: giveIPValue(answer, gp, type, atTime);
1995 }
1996
1997 return 1; // to make the compiler happy
1998}
1999
2000
2001void
2002AnisotropicDamageMaterial :: initializeFrom(InputRecord &ir)
2003{
2004 StructuralMaterial :: initializeFrom(ir);
2005 linearElasticMaterial.initializeFrom(ir);
2006 E = linearElasticMaterial.giveYoungsModulus();
2007 nu = linearElasticMaterial.givePoissonsRatio();
2008 int eqStrain = 0;
2009 // specify the type of formula for equivalent strain
2010 // currently only the Mazars formula is allowed !!! (this should be generalized later)
2011 //
2012 // IR_GIVE_OPTIONAL_FIELD(ir, eqStrain, _IFT_AnisotropicDamageMaterial_equivStrainType);
2013 //
2014 switch ( eqStrain ) {
2015 case 1: this->equivStrainType = EST_Rankine_Smooth;
2016 break;
2017 case 2: this->equivStrainType = EST_ElasticEnergy;
2018 break;
2019 case 3: this->equivStrainType = EST_Mises;
2020 // IR_GIVE_FIELD(ir, k, _IFT_IsotropicDamageMaterial1_k);
2021 break;
2022 case 4: this->equivStrainType = EST_Rankine_Standard;
2023 break;
2025 break;
2027 break;
2028 case 7: this->equivStrainType = EST_Griffith;
2029 break;
2030 default: this->equivStrainType = EST_Mazars;
2031 }
2032
2033 int damlaw = 0;
2034 // specify the type of damage law (which affects the shape of the stress-strain curve)
2036 switch ( damlaw ) {
2037 case 1: this->damageLawType = DLT_Desmorat1;
2038 break;
2039 case 2: this->damageLawType = DLT_Desmorat2;
2040 break;
2041 case 3: this->damageLawType = DLT_Linear;
2042 break;
2043 case 4: this->damageLawType = DLT_Exponential;
2044 break;
2045 default: this->damageLawType = DLT_Desmorat2;
2046 }
2047
2050 if ( damageLawType == DLT_Desmorat2 ) {
2052 }
2053}
2054
2055void
2056AnisotropicDamageMaterial :: giveInputRecord(DynamicInputRecord &input)
2057{
2058 StructuralMaterial :: giveInputRecord(input);
2060}
2061
2062
2063AnisotropicDamageMaterialStatus :: AnisotropicDamageMaterialStatus(GaussPoint *g) : StructuralMaterialStatus(g),
2064 damage(3, 3), tempDamage(3, 3)
2065{}
2066
2067void
2068AnisotropicDamageMaterialStatus :: printOutputAt(FILE *file, TimeStep *tStep) const
2069{
2070 MaterialMode mode = gp->giveMaterialMode();
2071 if ( mode == _PlaneStress ) { // special treatment of the out-of-plane strain
2072 FloatArray helpVec;
2073 MaterialStatus :: printOutputAt(file, tStep);
2074 fprintf(file, " strains ");
2075 StructuralMaterial :: giveFullSymVectorForm(helpVec, strainVector, mode);
2076 helpVec.at(3) = this->strainZ;
2077 for ( auto &v : helpVec ) {
2078 fprintf( file, " %.4e", v );
2079 }
2080 fprintf(file, "\n stresses");
2081 StructuralMaterial :: giveFullSymVectorForm(helpVec, stressVector, mode);
2082 for ( auto &v : helpVec ) {
2083 fprintf( file, " %.4e", v );
2084 }
2085 fprintf(file, "\n");
2086 } else {
2087 StructuralMaterialStatus :: printOutputAt(file, tStep); // standard treatment of strains and stresses
2088 }
2089
2090 fprintf(file, "status { ");
2091 fprintf(file, "kappa %g", this->kappa);
2092 double damtrace = tempDamage.giveTrace();
2093 if ( damtrace > 0.0 ) {
2094 fprintf(file, ", damage");
2095 int n = tempDamage.giveNumberOfRows();
2096 for ( int i = 1; i <= n; i++ ) {
2097 for ( int j = 1; j <= n; j++ ) {
2098 fprintf( file, " %g", tempDamage.at(i, j) );
2099 }
2100 }
2101 }
2102
2103#ifdef keep_track_of_dissipated_energy
2104 fprintf(file, ", dissW %f, freeE %f, stressW %f ", this->dissWork, ( this->stressWork ) - ( this->dissWork ), this->stressWork);
2105#endif
2106
2107 fprintf(file, "}\n");
2108}
2109
2110
2111
2112void
2113AnisotropicDamageMaterialStatus :: initTempStatus()
2114{
2115 StructuralMaterialStatus :: initTempStatus();
2116 this->tempKappa = this->kappa;
2117 this->tempDamage = this->damage;
2118 this->tempStrainZ = this->strainZ;
2119 this->tempStoredFactor = this->storedFactor;
2120#ifdef keep_track_of_dissipated_energy
2121 this->tempStressWork = this->stressWork;
2122 this->tempDissWork = this->dissWork;
2123#endif
2124}
2125
2126
2127void
2128AnisotropicDamageMaterialStatus :: updateYourself(TimeStep *atTime)
2129{
2130 StructuralMaterialStatus :: updateYourself(atTime);
2131 this->kappa = this->tempKappa;
2132 this->damage = this->tempDamage;
2133 this->strainZ = this->tempStrainZ;
2134 this->flag = this->tempFlag;
2135 this->storedFactor = this->tempStoredFactor;
2136#ifdef keep_track_of_dissipated_energy
2137 this->stressWork = this->tempStressWork;
2138 this->dissWork = this->tempDissWork;
2139#endif
2140}
2141
2142
2143void
2144AnisotropicDamageMaterialStatus :: saveContext(DataStream &stream, ContextMode mode)
2145{
2146 StructuralMaterialStatus :: saveContext(stream, mode);
2147
2148#if 0
2149 if ( !stream.write(kappa) ) {
2151 }
2152
2153 if ( !stream.write(damage) ) {
2155 }
2156
2157#ifdef keep_track_of_dissipated_energy
2158 if ( !stream.write(stressWork) ) {
2160 }
2161
2162 if ( !stream.write(dissWork) ) {
2164 }
2165#endif
2166#endif
2167
2168 if ( !stream.write(kappa) ) {
2170 }
2171
2172 contextIOResultType iores;
2173 if ( ( iores = damage.storeYourself(stream) ) != CIO_OK ) {
2174 THROW_CIOERR(iores);
2175 }
2176
2177#ifdef keep_track_of_dissipated_energy
2178 if ( !stream.write(stressWork) ) {
2180 }
2181
2182 if ( !stream.write(dissWork) ) {
2184 }
2185#endif
2186}
2187
2188
2189void
2190AnisotropicDamageMaterialStatus :: restoreContext(DataStream &stream, ContextMode mode)
2191{
2192 StructuralMaterialStatus :: restoreContext(stream, mode);
2193
2194#if 0
2195 if ( !stream.read(kappa) ) {
2197 }
2198
2199 if ( !stream.read(damage) ) {
2201 }
2202
2203#ifdef keep_track_of_dissipated_energy
2204 if ( !stream.read(stressWork) ) {
2206 }
2207
2208 if ( !stream.read(dissWork) ) {
2210 }
2211
2212#endif
2213#endif
2214
2215 if ( !stream.read(kappa) ) {
2217 }
2218
2219 contextIOResultType iores;
2220 if ( ( iores = damage.restoreYourself(stream) ) != CIO_OK ) {
2221 THROW_CIOERR(iores);
2222 }
2223
2224#ifdef keep_track_of_dissipated_energy
2225 if ( !stream.read(stressWork) ) {
2227 }
2228
2229 if ( !stream.read(dissWork) ) {
2231 }
2232#endif
2233}
2234
2235#ifdef keep_track_of_dissipated_energy
2236void
2237AnisotropicDamageMaterialStatus :: computeWork(GaussPoint *gp)
2238{
2239 /*
2240 * // strain increment
2241 * FloatArray deps;
2242 * deps.beDifferenceOf(tempStrainVector, strainVector);
2243 *
2244 * // increment of stress work density
2245 * double dSW = ( tempStressVector.dotProduct(deps) + stressVector.dotProduct(deps) ) / 2.;
2246 * tempStressWork = stressWork + dSW;
2247 *
2248 * // elastically stored energy density
2249 * double We = tempStressVector.dotProduct(tempStrainVector) / 2.;
2250 *
2251 * // dissipative work density
2252 * tempDissWork = tempStressWork - We;
2253 */
2254}
2255#endif
2256} // end namespace oofem
#define AD_TOLERANCE
#define AD_MAXITER
#define _IFT_AnisotropicDamageMaterial_aA
#define _IFT_AnisotropicDamageMaterial_damageLawType
#define _IFT_AnisotropicDamageMaterial_kappaf
#define _IFT_AnisotropicDamageMaterial_kappa0
#define REGISTER_Material(class)
double tempStrainZ
Non-equilibrated out-of-plane value for 2dPlaneStress mode.
double strainZ
Out-of-plane value for 2dPlaneStress mode.
double giveStressWork()
Returns the density of total work of stress on strain increments.
void setTempDamage(const FloatMatrix &d)
Assigns temp. damage tensor to given tensor d.
double giveDissWork()
Returns the density of dissipated work.
double giveTempStrainZ()
Returns the temp scalar measure of the out-of-plane strain to given value (for 2dPlaneStress mode).
int giveFlag()
Returns the value of the flag.
double kappa
Scalar measure of the largest strain level ever reached in material.
FloatMatrix tempDamage
Non-equilibrated second order damage tensor.
FloatMatrix damage
Second order damage tensor.
double tempStressWork
Non-equilibrated density of total work done by stresses on strain increments.
double dissWork
Density of dissipated work.
void setTempKappa(double newKappa)
Sets the temp scalar measure of the largest strain level to given value.
void computeWork(GaussPoint *gp)
Computes the increment of total stress work and of dissipated work.
int flag
This flag turns into 1 and remains 1 when the trace of the damage tensor is >1 in compression (tr(str...
double tempKappa
Non-equilibrated scalar measure of the largest strain level.
const FloatMatrix & giveTempDamage()
Returns the temp. second order damage tensor.
double tempDissWork
Non-equilibrated density of dissipated work.
double stressWork
Density of total work done by stresses on strain increments.
void setTempStrainZ(double newStrainZ)
Sets the temp scalar measure of the out-of-plane strain to given value (for 2dPlaneStress mode).
double giveKappa()
Returns the last equilibrated scalar measure of the largest strain level.
const FloatMatrix & giveDamage()
Returns the last equilibrated second order damage tensor.
IsotropicLinearElasticMaterial linearElasticMaterial
Reference to bulk (undamaged) material.
void computeDamage(FloatMatrix &tempDamage, const FloatMatrix &damage, double kappa, double eps1, double eps2, double ceps, double seps, double epsZ) const
double computeTraceD(double equivStrain) const
double obtainAlpha2(FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, FloatMatrix ProjMatrix, double damageThreshold) const
Obtains the proportion of the damage tensor that is needed to get the second eigenvalue equal to the ...
double computeKappa(FloatMatrix damageTensor) const
double computeDimensionlessOutOfPlaneStress(const FloatArray &inplaneStrain, double epsZ, const FloatMatrix &dam) const
double kappa0
Damage threshold kappa0, as defined in the paper mentioned above.
bool checkPrincVal2D(double Dx, double Dy, double Dxy) const
void computeEquivalentStrain(double &kappa, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) const
double computeCorrectionFactor(FloatMatrix tempDamageTensor, FloatMatrix strainTensor, GaussPoint *gp) const
double aA
Damage parameter a*A, needed to obtain Kappa(trD), according to eq. 33 in the paper mentioned above.
double obtainAlpha3(FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, FloatArray vec3, double damageThreshold) const
Obtains the proportion of the damage tensor that is needed to get the third eigenvalue equal to the d...
DamLawType damageLawType
Parameter specifying the damage law.
void correctBigValues(FloatMatrix &matrix) const
double kappaf
Damage parameter kappa_f (in the paper denoted as "a").
void computeDamageTensor(FloatMatrix &answer, GaussPoint *gp, const FloatArray &totalStrain, double equivStrain, TimeStep *atTime) const
void computeInplaneStress(FloatArray &inplaneStress, const FloatArray &inplaneStrain, double epsZ, const FloatMatrix &dam) const
double obtainAlpha1(FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, double damageThreshold) const
Obtains the proportion of the damage tensor that is needed to get the first eigenvalue equal to the d...
void computePrincValDir2D(double &D1, double &D2, double &c, double &s, double Dx, double Dy, double Dxy) const
void computePlaneStressStrain(FloatMatrix &answer, FloatMatrix damageTensor, FloatArray totalStrain, GaussPoint *gp, TimeStep *atTime) const
EquivStrainType equivStrainType
Parameter specifying the definition of equivalent strain.
double computeOutOfPlaneStrain(const FloatArray &inplaneStrain, const FloatMatrix &dam, bool tens_flag) const
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 setField(int item, InputFieldType id)
void beSymVectorForm(const FloatMatrix &aMatrix)
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
bool isEmpty() const
Returns true if receiver is empty.
Definition floatarray.h:265
double at(std::size_t i, std::size_t j) const
void times(double f)
double giveTrace() const
void add(const FloatMatrix &a)
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
bool jaco_(FloatArray &eval, FloatMatrix &v, int nf)
void zero()
Zeroes all coefficient of receiver.
void beDyadicProductOf(const FloatArray &vec1, const FloatArray &vec2)
void beMatrixForm(const FloatArray &aArray)
int giveNumberOfRows() const
Returns number of rows of receiver.
double at(std::size_t i, std::size_t j) const
void subtract(const FloatMatrix &a)
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition gausspoint.h:190
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
StructuralMaterialStatus(GaussPoint *g)
Constructor. Creates new StructuralMaterialStatus with IntegrationPoint g.
FloatArray stressVector
Equilibrated stress vector in reduced form.
FloatArray strainVector
Equilibrated strain vector in reduced form.
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.
StructuralMaterial(int n, Domain *d)
virtual void giveStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const
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 IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Definition inputrecord.h:75
#define IR_GIVE_FIELD(__ir, __value, __id)
Definition inputrecord.h:67
long ContextMode
Definition contextmode.h:43
double macbra(double x)
Returns the positive part of given float.
Definition mathfem.h:128
@ principal_strain
For computing principal strains from engineering strains.
@ CIO_IOERR
General IO error.
double trace(const FloatMatrixF< N, N > &mat)

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