OOFEM 3.0
Loading...
Searching...
No Matches
idmnl1.C
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
35#include "idmnl1.h"
36#include "gausspoint.h"
37#include "floatmatrix.h"
38#include "floatarray.h"
39#include "mathfem.h"
41#include "sparsemtrx.h"
42#include "error.h"
43#include "nonlocalmaterialext.h"
44#include "contextioerr.h"
45#include "stressvector.h"
46#include "strainvector.h"
47#include "classfactory.h"
48#include "dynamicinputrecord.h"
49#include "datastream.h"
51
52#ifdef __OOFEG
53 #include "oofeggraphiccontext.h"
54 #include "connectivitytable.h"
55#endif
56
57namespace oofem {
59
65
66
67void
68IDNLMaterial :: updateBeforeNonlocAverage(const FloatArray &strainVector, GaussPoint *gp, TimeStep *tStep) const
69{
70 /* Implements the service updating local variables in given integration points,
71 * which take part in nonlocal average process. Actually, no update is necessary,
72 * because the value used for nonlocal averaging is strain vector used for nonlocal secant stiffness
73 * computation. It is therefore necessary only to store local strain in corresponding status.
74 * This service is declared at StructuralNonlocalMaterial level.
75 */
76 FloatArray SDstrainVector;
77 double equivStrain;
78 IDNLMaterialStatus *nlstatus = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
79
80 this->initTempStatus(gp);
81
82 // subtract stress-independent part
83 // note: eigenStrains (temperature) is not contained in mechanical strain stored in gp
84 // therefore it is necessary to subtract always the total eigenstrain value
85 this->giveStressDependentPartOfStrainVector(SDstrainVector, gp, strainVector, tStep, VM_Total);
86
87 // compute and store the local variable to be averaged
88 // (typically the local equivalent strain)
89 nlstatus->letTempStrainVectorBe(SDstrainVector);
90 equivStrain = this->computeLocalEquivalentStrain(SDstrainVector, gp, tStep);
91
92 // nonstandard formulation based on averaging of compliance parameter gamma
93 // note: gamma is stored in a variable named localEquivalentStrainForAverage, which can be misleading
94 // perhaps this variable should later be renamed
95 if ( averagedVar == AVT_Compliance ) {
96 double gamma = complianceFunction(equivStrain, gp);
98 }
99 // nonstandard formulation based on averaging of damage variable omega
100 // note: omega is stored in a variable named localEquivalentStrainForAverage, which can be misleading
101 // perhaps this variable should later be renamed
102 else if ( averagedVar == AVT_Damage ) {
103 double omega = damageFunction(equivStrain, gp);
104 nlstatus->setLocalEquivalentStrainForAverage(omega);
105 }
106 // standard formulation based on averaging of equivalent strain
107 else {
108 nlstatus->setLocalEquivalentStrainForAverage(equivStrain);
109 }
110
111 // influence of damage on weight function
112 if ( averType >= 2 && averType <= 6 ) {
114 }
115}
116
117double
118IDNLMaterial :: giveNonlocalMetricModifierAt(GaussPoint *gp) const
119{
120 auto status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
121 double damage = status->giveTempDamage();
122 if ( damage == 0. ) {
123 damage = status->giveDamage();
124 }
125 return damage;
126}
127
128void
129IDNLMaterial :: computeAngleAndSigmaRatio(double &nx, double &ny, double &ratio, GaussPoint *gp, bool &flag) const
130{
131 auto status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
132 MaterialMode matMode = gp->giveMaterialMode();
133 if ( matMode != _PlaneStress ) { //Check if the stress-based approach can be applied
134 OOFEM_ERROR("Stress-based nonlocal averaging is implemented for plane stress only");
135 }
136
137 //Get the temporary strain vector
138 FloatArray strainFloatArray = status->giveTempStrainVector();
139
140 /* This old implementation could be generalized more easily to arbitrary types of stress,
141 * but for plane stress we use a more efficient direct implementation
142 * (note that the stress-based nonlocal model is very expensive and these operations are repeated many times)
143 *
144 * //Check if strain vector is zero. In this case this function is not going to modify nonlocal radius
145 * if ( strainFloatArray.computeNorm() == 0 ) {
146 * flag = false;
147 * return;
148 * }
149 * //Convert the FloatArray to StrainVector
150 * StrainVector strain(strainFloatArray, matMode);
151 * //Compute effective Stress tensor
152 * StressVector effectiveStress(matMode);
153 * const double E = this->giveLinearElasticMaterial()->give('E', gp);
154 * const double nu = this->giveLinearElasticMaterial()->give('n', gp);
155 * strain.applyElasticStiffness(effectiveStress, E, nu);
156 * //Compute principal values and eigenvectors of effective stress tensor
157 * FloatArray principalStress;
158 * FloatMatrix princDir;
159 * effectiveStress.computePrincipalValDir(principalStress, princDir);
160 * //Calculate components of the first eigenvector
161 * nx = princDir.at(1, 1);
162 * ny = princDir.at(2, 1);
163 * //Calculate ratio of principal stresses
164 * if ( principalStress.at(2) < 0. && principalStress.at(1) < 0. ) { //Both eigenvalues negative
165 * ratio = 1.;
166 * flag = false; // modification of nonlocal weights not done under biaxial compression
167 * } else if ( principalStress.at(2) < 0. ) { //One eigenvalue positive
168 * ratio = 0.;
169 * } else {
170 * ratio = principalStress.at(2) / principalStress.at(1); //Both eigenvalues positive
171 * }
172 */
173
174 /* New implementation, directly finds the eigenvalues and eigenvector using an optimized scheme for plane stress */
175 double epsx = strainFloatArray.at(1);
176 double epsy = strainFloatArray.at(2);
177 double gamxy = strainFloatArray.at(3);
178 //Check if strain vector is zero. In this case this function is not going to modify nonlocal radius
179 if ( epsx == 0. && epsy == 0. && gamxy == 0. ) {
180 flag = false;
181 return;
182 }
183 double aux = sqrt( ( epsx - epsy ) * ( epsx - epsy ) + gamxy * gamxy );
184 double e1 = epsx + epsy + aux; // e1 = 2 times the maximum principal strain
185 double e2 = epsx + epsy - aux; // e2 = 2 times the minimum principal strain
186 const double nu = this->linearElasticMaterial->give('n', gp);
187 double s1 = e1 + nu * e2; // s1 = 2*(1-nu*nu)/E times the maximum principal stress
188 double s2 = e2 + nu * e1; // s2 = 2*(1-nu*nu)/E times the minimum principal stress
189
190 //Calculate ratio of principal stresses
191 if ( s1 <= 0. ) { //No positive eigenvalue
192 ratio = 1.;
193 flag = false; // modification of nonlocal weights not done under biaxial compression
194 } else if ( s2 <= 0. ) { //One positive eigenvalue
195 ratio = 0.;
196 } else {
197 ratio = s2 / s1; //Two positive eigenvalues
198 }
199
200 //Calculate components of the first eigenvector
201 nx = gamxy;
202 ny = e1 - 2. * epsx;
203 aux = nx * nx + ny * ny;
204 if ( aux == 0. ) {
205 nx = e1 - 2. * epsy;
206 ny = gamxy;
207 aux = nx * nx + ny * ny;
208 if ( aux == 0. ) {
209 nx = 1.;
210 ny = 0.;
211 return;
212 }
213 }
214 aux = sqrt(aux);
215 nx /= aux;
216 ny /= aux;
217}
218
219double
220IDNLMaterial :: computeStressBasedWeight(double cl, double &nx, double &ny, double &ratio, GaussPoint *gp, GaussPoint *jGp, double weight) const
221{
222 // Take into account periodicity, if required
223 if ( this->px > 0. ) {
224 return computeStressBasedWeightForPeriodicCell(cl, nx, ny, ratio, gp, jGp);
225 }
226
227 //Check if source and receiver point coincide
228 if ( gp == jGp ) {
229 return weight;
230 }
231 //Compute distance between source and receiver point
232 FloatArray gpCoords, distance;
235 distance.subtract(gpCoords); // Vector connecting the two Gauss points
236
237 //Compute modified distance
238 double x1 = nx * distance.at(1) + ny *distance.at(2);
239 double x2 = -ny *distance.at(1) + nx *distance.at(2);
240 // Compute axis of ellipse and scale/stretch weak axis so that ellipse is converted to circle
241 double gamma = this->beta + ( 1. - beta ) * ratio * ratio;
242 x2 /= gamma;
243 double modDistance = sqrt(x1 * x1 + x2 * x2);
244
245 //Get new weight
246 double updatedWeight = this->computeWeightFunction(cl, modDistance);
247 updatedWeight *= jGp->giveElement()->computeVolumeAround(jGp); //weight * (Volume where the weight is applied)
248 return updatedWeight;
249}
250
251// This method is a slight modification of IDNLMaterial :: computeStressBasedWeight but is implemented separately,
252// to keep the basic method as simple (and efficient) as possible
253double
254IDNLMaterial :: computeStressBasedWeightForPeriodicCell(double cl, double &nx, double &ny, double &ratio, GaussPoint *gp, GaussPoint *jGp) const
255{
256 double updatedWeight = 0.;
257 FloatArray gpCoords, distance;
259 int ix, nper = 1; // could be increased in the future, if needed
260
261 for ( ix = -nper; ix <= nper; ix++ ) { // loop over periodic images shifted in x-direction
263 distance.at(1) += ix * px; // shift the x-coordinate
264 distance.subtract(gpCoords); // Vector connecting the two Gauss points
265
266 //Compute modified distance
267 double x1 = nx * distance.at(1) + ny *distance.at(2);
268 double x2 = -ny *distance.at(1) + nx *distance.at(2);
269 // Compute axis of ellipse and scale/stretch weak axis so that ellipse is converted to circle
270 double gamma = this->beta + ( 1. - beta ) * ratio * ratio;
271 x2 /= gamma;
272 double modDistance = sqrt(x1 * x1 + x2 * x2);
273
274 //Get new weight
275 double updatedWeightContribution = this->computeWeightFunction(cl, modDistance);
276 if ( updatedWeightContribution > 0. ) {
277 updatedWeightContribution *= jGp->giveElement()->computeVolumeAround(jGp); //weight * (Volume where the weight is applied)
278 updatedWeight += updatedWeightContribution;
279 }
280 }
281 return updatedWeight;
282}
283
284double
285IDNLMaterial :: computeEquivalentStrain(const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) const
286{
287 double nonlocalContribution, nonlocalEquivalentStrain = 0.0;
288 IDNLMaterialStatus *nonlocStatus, *status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
289
290 this->buildNonlocalPointTable(gp);
292
293 // compute nonlocal equivalent strain
294 // or nonlocal compliance variable gamma (depending on averagedVar)
295
296 auto list = this->giveIPIntegrationList(gp); // !
297
298 double sigmaRatio = 0.; //ratio sigma2/sigma1 used for stress-based averaging
299 double nx, ny; //components of the first principal stress direction (for stress-based averaging)
300 double updatedIntegrationVolume = 0.; //new integration volume. Sum of all new weights used for stress-based averaging
301 //Flag to deactivate stress-based nonlocal averaging for zero stress states.
302 // When SBAflag is not set, no stress-based averaging takes place.
303 // When SBAflag is set, stress-based averaging takes place.
304 bool SBAflag = ( this->nlvar == NLVT_StressBased );
305 //Check if Stress based averaging is enforced and calculate the angle of the first eigenvector and the sigmaratio
306 if ( SBAflag ) {
307 computeAngleAndSigmaRatio(nx, ny, sigmaRatio, gp, SBAflag);
308 }
309
310 //Loop over all Gauss points which are in gp's integration domain
311 for ( auto &lir : *list ) {
312 GaussPoint *neargp = lir.nearGp;
313 nonlocStatus = static_cast< IDNLMaterialStatus * >( neargp->giveMaterialStatus() );
314 nonlocalContribution = nonlocStatus->giveLocalEquivalentStrainForAverage();
315
316 if ( SBAflag ) { //Check if Stress Based Averaging is requested and calculate nonlocal contribution
317 double stressBasedWeight = computeStressBasedWeight(cl, nx, ny, sigmaRatio, gp, neargp, lir.weight); //Compute new weight
318 updatedIntegrationVolume += stressBasedWeight;
319 nonlocalContribution *= stressBasedWeight;
320 } else {
321 nonlocalContribution *= lir.weight;
322 }
323
324 nonlocalEquivalentStrain += nonlocalContribution;
325 }
326
327 if ( SBAflag ) { // Nonlocal weights are modified in stress-based averaging. Thus the integration volume needs to be modified
328 //status->setIntegrationScale(updatedIntegrationVolume);
329 nonlocalEquivalentStrain /= updatedIntegrationVolume;
330 } else if ( scaling == ST_Standard ) { // standard rescaling
331 nonlocalEquivalentStrain *= 1. / status->giveIntegrationScale();
332 } else if ( scaling == ST_Borino ) { // Borino modification
333 double scale = status->giveIntegrationScale();
334 if ( scale > 1. ) {
335 nonlocalEquivalentStrain *= 1. / scale;
336 } else {
337 nonlocalEquivalentStrain += ( 1. - scale ) * status->giveLocalEquivalentStrainForAverage();
338 }
339 }
340
341
342 // undernonlocal or overnonlocal formulation
343 if ( mm != 1. ) {
344 double localEquivalentStrain = status->giveLocalEquivalentStrainForAverage();
345 if ( mm >= 0. ) { // harmonic averaging
346 if ( localEquivalentStrain > 0. && nonlocalEquivalentStrain > 0. ) {
347 nonlocalEquivalentStrain = 1. / ( mm / nonlocalEquivalentStrain + ( 1. - mm ) / localEquivalentStrain );
348 } else {
349 nonlocalEquivalentStrain = 0.;
350 }
351 } else { // arithmetic averaging, -mm is used instead of mm
352 nonlocalEquivalentStrain = -mm * nonlocalEquivalentStrain + ( 1. + mm ) * localEquivalentStrain;
353 }
354 }
355
356 this->endIPNonlocalAverage(gp); // ???????????????????
357
358 return nonlocalEquivalentStrain;
359}
360
361Interface *
362IDNLMaterial :: giveInterface(InterfaceType type)
363{
365 return static_cast< StructuralNonlocalMaterialExtensionInterface * >(this);
366 } else if ( type == NonlocalMaterialStiffnessInterfaceType ) {
367 return static_cast< NonlocalMaterialStiffnessInterface * >(this);
368 } else if ( type == MaterialModelMapperInterfaceType ) {
369 return static_cast< MaterialModelMapperInterface * >(this);
370 } else {
371 return NULL;
372 }
373}
374
375
376void
377IDNLMaterial :: initializeFrom(InputRecord &ir)
378{
379 IsotropicDamageMaterial1 :: initializeFrom(ir);
380 StructuralNonlocalMaterialExtensionInterface :: initializeFrom(ir);
381}
382
383
384void
385IDNLMaterial :: giveInputRecord(DynamicInputRecord &input)
386{
387 IsotropicDamageMaterial1 :: giveInputRecord(input);
388 StructuralNonlocalMaterialExtensionInterface :: giveInputRecord(input);
389}
390
391double
392IDNLMaterial :: computeDamageParam(double kappa, const FloatArray &strain, GaussPoint *g) const
393{
394 if ( averagedVar == AVT_Compliance ) {
395 // formulation based on nonlocal gamma (here formally called kappa)
396 return kappa / ( 1. + kappa );
397 } else if ( averagedVar == AVT_Damage ) {
398 // formulation based on nonlocal damage (here formally called kappa)
399 return kappa;
400 } else {
401 // formulation based on nonlocal equivalent strain
402 return damageFunction(kappa, g);
403 }
404}
405
406void
407IDNLMaterial :: NonlocalMaterialStiffnessInterface_addIPContribution(SparseMtrx &dest, const UnknownNumberingScheme &s,
408 GaussPoint *gp, TimeStep *tStep)
409{
410 double coeff;
411 IDNLMaterialStatus *status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
412 auto list = status->giveIntegrationDomainList();
413 IDNLMaterial *rmat;
414 FloatArray rcontrib, lcontrib;
415 IntArray loc, rloc;
416
417 FloatMatrix contrib;
418
419 if ( this->giveLocalNonlocalStiffnessContribution(gp, loc, s, lcontrib, tStep) == 0 ) {
420 return;
421 }
422
423 for ( auto &lir : *list ) {
424 rmat = dynamic_cast< IDNLMaterial * >( lir.nearGp->giveMaterial() );
425 if ( rmat ) {
426 rmat->giveRemoteNonlocalStiffnessContribution(lir.nearGp, rloc, s, rcontrib, tStep);
427 coeff = gp->giveElement()->computeVolumeAround(gp) * lir.weight / status->giveIntegrationScale();
428 // printf ("\nelement %d:", gp->giveElement()->giveNumber());
429 // lcontrib.printYourself();
430 // rcontrib.printYourself();
431 // assemble the contribution
432 // dest.checkSizeTowards (loc, rloc);
433 // dest.assemble (lcontrib,loc, rcontrib,rloc);
434
435 /* local effective assembly
436 * int i,j, r, c;
437 * for (i=1; i<= loc.giveSize(); i++)
438 * for (j=1; j<=rloc.giveSize(); j++) {
439 * r = loc.at(i);
440 * c = rloc.at(j);
441 * if ((r != 0) && (c!=0)) dest.at(r,c) -= (double) (lcontrib.at(i)*rcontrib.at(j)*coeff);
442 * }
443 */
444 contrib.clear();
445 contrib.plusDyadUnsym(lcontrib, rcontrib, -1.0 * coeff);
446 dest.assemble(loc, rloc, contrib);
447 }
448 }
449}
450
451std :: vector< localIntegrationRecord > *
452IDNLMaterial :: NonlocalMaterialStiffnessInterface_giveIntegrationDomainList(GaussPoint *gp)
453{
454 IDNLMaterialStatus *status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
455 this->buildNonlocalPointTable(gp);
456 return status->giveIntegrationDomainList();
457}
458
459
460
461int
462IDNLMaterial :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
463{
464 IDNLMaterialStatus *status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
465 if ( type == IST_LocalEquivalentStrain ) {
466 answer.resize(1);
467 answer.zero();
468 answer.at(1) = status->giveLocalEquivalentStrainForAverage();
469 } else {
470 return IsotropicDamageMaterial1 :: giveIPValue(answer, gp, type, tStep);
471 }
472
473 return 1; // to make the compiler happy
474}
475
476
477#ifdef __OOFEG
478void
479IDNLMaterial :: NonlocalMaterialStiffnessInterface_showSparseMtrxStructure(GaussPoint *gp, oofegGraphicContext &gc, TimeStep *tStep)
480{
481 IntArray loc, rloc;
482 FloatArray strain;
483 double f, equivStrain;
484 IDNLMaterialStatus *status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
485 IDNLMaterial *rmat;
486
487 const double e0 = this->give(e0_ID, gp);
488 //const double ef = this->give(ef_ID, gp);
489
490 strain = status->giveTempStrainVector();
491 // compute equivalent strain
492 equivStrain = this->computeEquivalentStrain(strain, gp, tStep);
493 f = equivStrain - status->giveTempKappa();
494
495 if ( ( equivStrain <= e0 ) || ( f < 0.0 ) ) {
496 return;
497 }
498
499 EASValsSetLineWidth(OOFEG_SPARSE_PROFILE_WIDTH);
500 EASValsSetColor( gc.getExtendedSparseProfileColor() );
501 EASValsSetLayer(OOFEG_SPARSE_PROFILE_LAYER);
502 EASValsSetFillStyle(FILL_SOLID);
503
504 WCRec p [ 4 ];
505 GraphicObj *go;
506
508
509 int n, m;
510 auto list = status->giveIntegrationDomainList();
511 for ( auto &lir : *list ) {
512 rmat = dynamic_cast< IDNLMaterial * >( lir.nearGp->giveMaterial() );
513 if ( rmat ) {
514 lir.nearGp->giveElement()->giveLocationArray( rloc, EModelDefaultEquationNumbering() );
515 } else {
516 continue;
517 }
518
519 n = loc.giveSize();
520 m = rloc.giveSize();
521 for ( int i = 1; i <= n; i++ ) {
522 if ( loc.at(i) == 0 ) {
523 continue;
524 }
525
526 for ( int j = 1; j <= m; j++ ) {
527 if ( rloc.at(j) == 0 ) {
528 continue;
529 }
530
531 if ( gc.getSparseProfileMode() == 0 ) {
532 p [ 0 ].x = ( FPNum ) loc.at(i) - 0.5;
533 p [ 0 ].y = ( FPNum ) rloc.at(j) - 0.5;
534 p [ 0 ].z = 0.;
535 p [ 1 ].x = ( FPNum ) loc.at(i) + 0.5;
536 p [ 1 ].y = ( FPNum ) rloc.at(j) - 0.5;
537 p [ 1 ].z = 0.;
538 p [ 2 ].x = ( FPNum ) loc.at(i) + 0.5;
539 p [ 2 ].y = ( FPNum ) rloc.at(j) + 0.5;
540 p [ 2 ].z = 0.;
541 p [ 3 ].x = ( FPNum ) loc.at(i) - 0.5;
542 p [ 3 ].y = ( FPNum ) rloc.at(j) + 0.5;
543 p [ 3 ].z = 0.;
544 go = CreateQuad3D(p);
545 EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | LAYER_MASK, go);
546 EMAddGraphicsToModel(ESIModel(), go);
547 } else {
548 p [ 0 ].x = ( FPNum ) loc.at(i);
549 p [ 0 ].y = ( FPNum ) rloc.at(j);
550 p [ 0 ].z = 0.;
551
552 EASValsSetMType(SQUARE_MARKER);
553 go = CreateMarker3D(p);
554 EGWithMaskChangeAttributes(COLOR_MASK | LAYER_MASK | VECMTYPE_MASK, go);
555 EMAddGraphicsToModel(ESIModel(), go);
556 }
557 }
558 }
559 }
560}
561#endif
562
563
564
565
566int
567IDNLMaterial :: giveLocalNonlocalStiffnessContribution(GaussPoint *gp, IntArray &loc, const UnknownNumberingScheme &s,
568 FloatArray &lcontrib, TimeStep *tStep)
569{
570 int nrows, nsize;
571 double sum, f, equivStrain;
572 auto status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
573 auto elem = static_cast< StructuralElement * >( gp->giveElement() );
574 FloatMatrix b, de;
575 FloatArray stress, strain;
576
577 const double e0 = this->give(e0_ID, gp);
578 const double ef = this->give(ef_ID, gp);
579
580 /*
581 * if (fabs(status->giveTempDamage()) <= 1.e-10) {
582 * // already eleastic regime
583 * loc.clear();
584 * return 0;
585 * }
586 */
587 strain = status->giveTempStrainVector();
588
589 // compute equivalent strain
590 equivStrain = this->computeEquivalentStrain(strain, gp, tStep);
591 f = equivStrain - status->giveTempKappa();
592
593 if ( ( equivStrain <= e0 ) || ( f < 0.0 ) ) {
594 /*
595 * if (status->lst == IDNLMaterialStatus::LST_elastic)
596 * printf (" ");
597 * else printf ("_");
598 * status->lst = IDNLMaterialStatus::LST_elastic;
599 */
600 loc.clear();
601 return 0;
602 } else {
603 if ( status->giveDamage() >= 1.00 ) {
604 // printf ("f");
605 return 0;
606 }
607
608 /*
609 * if (status->lst == IDNLMaterialStatus::LST_loading)
610 * printf ("o");
611 * else printf ("O");
612 * status->lst = IDNLMaterialStatus::LST_loading;
613 */
614
615 // no support for reduced integration now
616 elem->computeBmatrixAt(gp, b);
617
619 lmat->giveStiffnessMatrix(de, SecantStiffness, gp, tStep);
620 stress.beProductOf(de, strain);
621
622 f = ( e0 / ( equivStrain * equivStrain ) ) * exp( -( equivStrain - e0 ) / ( ef - e0 ) )
623 + ( e0 / equivStrain ) * exp( -( equivStrain - e0 ) / ( ef - e0 ) ) * 1.0 / ( ef - e0 );
624
625 nrows = b.giveNumberOfColumns();
626 nsize = stress.giveSize();
627 lcontrib.resize(nrows);
628 for ( int i = 1; i <= nrows; i++ ) {
629 sum = 0.0;
630 for ( int j = 1; j <= nsize; j++ ) {
631 sum += b.at(j, i) * stress.at(j);
632 }
633
634 lcontrib.at(i) = sum * f;
635 }
636 }
637
638 // request element code numbers
639 elem->giveLocationArray(loc, s);
640
641 return 1;
642}
643
644
645void
646IDNLMaterial :: giveRemoteNonlocalStiffnessContribution(GaussPoint *gp, IntArray &rloc, const UnknownNumberingScheme &s,
647 FloatArray &rcontrib, TimeStep *tStep)
648{
649 int ncols, nsize;
650 double coeff = 0.0, sum;
651 IDNLMaterialStatus *status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
652 StructuralElement *elem = static_cast< StructuralElement * >( gp->giveElement() );
653 FloatMatrix b, de, princDir(3, 3), t;
654 FloatArray stress, fullStress, strain, principalStress, help, nu;
655
656 elem->giveLocationArray(rloc, s);
657 // no support for reduced integration now
658 elem->computeBmatrixAt(gp, b);
659
660 if ( this->equivStrainType == EST_Rankine_Standard ) {
661 FloatArray fullHelp;
663
664 lmat->giveStiffnessMatrix(de, SecantStiffness, gp, tStep);
665 strain = status->giveTempStrainVector();
666 stress.beProductOf(de, strain);
667 StructuralMaterial :: giveFullSymVectorForm( fullStress, stress, gp->giveMaterialMode() );
668 if ( gp->giveMaterialMode() == _1dMat ) {
669 principalStress = fullStress;
670 } else {
671 this->computePrincipalValDir(principalStress, princDir, fullStress, principal_stress);
672 if ( ( gp->giveMaterialMode() == _3dMat ) || ( gp->giveMaterialMode() == _PlaneStrain ) ) {
673 ;
674 } else if ( gp->giveMaterialMode() == _PlaneStress ) {
675 // force the out-of-plane princ. dir to be last one
676 int indx = 0, zeroFlag = 1;
677 double swap;
678 for ( int i = 1; i <= 3; i++ ) {
679 if ( fabs( principalStress.at(i) ) > 1.e-10 ) {
680 zeroFlag = 0;
681 }
682
683 if ( princDir.at(3, i) > 0.90 ) {
684 indx = i;
685 }
686 }
687
688 if ( indx ) {
689 for ( int i = 1; i <= 3; i++ ) {
690 swap = princDir.at(i, indx);
691 princDir.at(i, indx) = princDir.at(i, 3);
692 princDir.at(i, 3) = swap;
693 }
694
695 swap = principalStress.at(indx);
696 principalStress.at(indx) = principalStress.at(3);
697 principalStress.at(3) = swap;
698 } else if ( zeroFlag == 0 ) {
699 OOFEM_ERROR("internal error");
700 }
701 } else {
702 OOFEM_ERROR("equivStrainType not supported");
703 }
704 }
705
706 sum = 0.0;
707 for ( int i = 1; i <= 3; i++ ) {
708 if ( principalStress.at(i) > 0.0 ) {
709 sum += principalStress.at(i) * principalStress.at(i);
710 }
711 }
712
713 if ( sum > 1.e-15 ) {
714 coeff = 1. / ( lmat->give('E', gp) * sqrt(sum) );
715 } else {
716 coeff = 0.0;
717 }
718
719 //
720 if ( gp->giveMaterialMode() != _1dMat ) {
721 t = this->giveStressVectorTranformationMtrx(princDir);
722 }
723
724 //
725 // if (gp->giveMaterialMode() != _1dMat) this->giveStressVectorTranformationMtrx (t, princDir,1);
726
727 // extract positive part
728 for ( int i = 1; i <= principalStress.giveSize(); i++ ) {
729 principalStress.at(i) = max(principalStress.at(i), 0.0);
730 }
731
732#if 0
733 this->giveNormalElasticStiffnessMatrix(den, SecantStiffness, gp, tStep);
734 help.beProductOf(den, principalStress);
735 fullHelp.resize(6);
736 for ( i = 1; i <= 3; i++ ) {
737 fullHelp.at(i) = help.at(i);
738 }
739
740 if ( gp->giveMaterialMode() != _1dMat ) {
741 fullNu.beProductOf(t, fullHelp);
742 //fullNu.beTProductOf (t, fullHelp);
743 crossSection->giveReducedCharacteristicVector(nu, gp, fullNu);
744 } else {
745 nu = help;
746 }
747
748#endif
749
750 /* Plane stress optimized version
751 *
752 *
753 * help.resize (3); help.zero();
754 * for (i=1; i<=3; i++) {
755 * help.at(1) += t.at(i,1)*principalStress.at(i);
756 * help.at(2) += t.at(i,2)*principalStress.at(i);
757 * help.at(3) += t.at(i,6)*principalStress.at(i);
758 * }
759 */
760 FloatArray fullPrincStress(6);
761 fullPrincStress.zero();
762 for ( int i = 1; i <= 3; i++ ) {
763 fullPrincStress.at(i) = principalStress.at(i);
764 }
765
766 fullHelp.beTProductOf(t, fullPrincStress);
767 StructuralMaterial :: giveReducedSymVectorForm( help, fullHelp, gp->giveMaterialMode() );
768
769 nu.beProductOf(de, help);
770 } else if ( this->equivStrainType == EST_ElasticEnergy ) {
771 double equivStrain;
772
774 lmat->giveStiffnessMatrix(de, SecantStiffness, gp, tStep);
775 strain = status->giveTempStrainVector();
776 stress.beProductOf(de, strain);
777 equivStrain = this->computeLocalEquivalentStrain(strain, gp, tStep);
778
779 nu = stress;
780 coeff = 1.0 / ( lmat->give('E', gp) * equivStrain );
781 } else {
782 OOFEM_ERROR("equivStrainType not supported");
783 }
784
785
786 ncols = b.giveNumberOfColumns();
787 nsize = nu.giveSize();
788 rcontrib.resize(ncols);
789 for ( int i = 1; i <= ncols; i++ ) {
790 sum = 0.0;
791 for ( int j = 1; j <= nsize; j++ ) {
792 sum += nu.at(j) * b.at(j, i);
793 }
794
795 rcontrib.at(i) = sum * coeff;
796 }
797}
798
799
800void
801IDNLMaterial :: giveNormalElasticStiffnessMatrix(FloatMatrix &answer,
802 MatResponseMode rMode,
803 GaussPoint *gp, TimeStep *tStep)
804{
805 //
806 // return Elastic Stiffness matrix for normal Stresses
807 auto de = linearElasticMaterial->give3dMaterialStiffnessMatrix(rMode, gp, tStep);
808 // This isn't used? Do we need one with zeroed entries (below) or the general 3d stiffness (above)?
809 //lMat->giveCharacteristicMatrix(de, rMode, gp, tStep);
810 //StructuralMaterial :: giveFullSymMatrixForm( de, deRed, gp->giveMaterialMode());
811
812 answer.resize(3, 3);
813 // copy first 3x3 submatrix to answer
814 for ( int i = 1; i <= 3; i++ ) {
815 for ( int j = 1; j <= 3; j++ ) {
816 answer.at(i, j) = de.at(i, j);
817 }
818 }
819}
820
821
827
828
829void
830IDNLMaterialStatus :: printOutputAt(FILE *file, TimeStep *tStep) const
831{
832 StructuralMaterialStatus :: printOutputAt(file, tStep);
833 fprintf(file, "status { ");
834 if ( this->damage > 0.0 ) {
835 fprintf(file, "nonloc-kappa %f, damage %f ", this->kappa, this->damage);
836
837#ifdef keep_track_of_dissipated_energy
838 fprintf(file, ", dissW %f, freeE %f, stressW %f ", this->dissWork, ( this->stressWork ) - ( this->dissWork ), this->stressWork);
839 } else {
840 fprintf(file, "stressW %f ", this->stressWork);
841#endif
842 }
843
844 fprintf(file, "}\n");
845}
846
847void
848IDNLMaterialStatus :: initTempStatus()
849//
850// initializes temp variables according to variables form previous equilibrium state.
851// builds new crackMap
852//
853{
854 IsotropicDamageMaterial1Status :: initTempStatus();
855}
856
857
858
859void
860IDNLMaterialStatus :: updateYourself(TimeStep *tStep)
861//
862// updates variables (nonTemp variables describing situation at previous equilibrium state)
863// after a new equilibrium state has been reached
864// temporary variables are having values corresponding to newly reched equilibrium.
865//
866{
867 IsotropicDamageMaterial1Status :: updateYourself(tStep);
868}
869
870
871
872void
873IDNLMaterialStatus :: saveContext(DataStream &stream, ContextMode mode)
874{
875 IsotropicDamageMaterial1Status :: saveContext(stream, mode);
876 //if (!stream.write(&localEquivalentStrainForAverage,1)) THROW_CIOERR(CIO_IOERR);
877}
878
879void
880IDNLMaterialStatus :: restoreContext(DataStream &stream, ContextMode mode)
881{
882 IsotropicDamageMaterial1Status :: restoreContext(stream, mode);
883 //if (!stream.read (&localEquivalentStrainForAverage,1)) THROW_CIOERR(CIO_IOERR);
884}
885
886Interface *
887IDNLMaterialStatus :: giveInterface(InterfaceType type)
888{
890 return static_cast< StructuralNonlocalMaterialStatusExtensionInterface * >(this);
891 } else {
892 return IsotropicDamageMaterial1Status :: giveInterface(type);
893 }
894}
895
896
897int
898IDNLMaterial :: packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
899{
900 IDNLMaterialStatus *status = static_cast< IDNLMaterialStatus * >( this->giveStatus(ip) );
901
902 this->buildNonlocalPointTable(ip);
904
905 return buff.write( status->giveLocalEquivalentStrainForAverage() );
906}
907
908int
909IDNLMaterial :: unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
910{
911 int result;
912 IDNLMaterialStatus *status = static_cast< IDNLMaterialStatus * >( this->giveStatus(ip) );
913 double localEquivalentStrainForAverage;
914
915 result = buff.read(localEquivalentStrainForAverage);
916 status->setLocalEquivalentStrainForAverage(localEquivalentStrainForAverage);
917 return result;
918}
919
920int
921IDNLMaterial :: estimatePackSize(DataStream &buff, GaussPoint *ip)
922{
923 //
924 // Note: status localStrainVectorForAverage memeber must be properly sized!
925 //
926
927 //IDNLMaterialStatus *status = (IDNLMaterialStatus*) this -> giveStatus (ip);
928
929 return buff.givePackSizeOfDouble(1);
930}
931
932
933double
934IDNLMaterial :: predictRelativeComputationalCost(GaussPoint *gp)
935{
936 //
937 // The values returned come from measurement
938 // do not change them unless you know what are you doing
939 //
940 double cost = 1.2;
941
942
943 if ( gp->giveMaterialMode() == _3dMat ) {
944 cost = 1.5;
945 }
946
947 IDNLMaterialStatus *status = static_cast< IDNLMaterialStatus * >( this->giveStatus(gp) );
948 indexType size = status->giveIntegrationDomainList()->size();
949 // just a guess (size/10) found optimal
950 // cost *= (1.0 + (size/10)*0.5);
951 cost *= ( 1.0 + size / 15.0 );
952
953 return cost;
954}
955} // end namespace oofem
#define REGISTER_Material(class)
virtual int read(int *data, std::size_t count)=0
Reads count integer values into array pointed by data.
virtual int givePackSizeOfDouble(std::size_t count)=0
virtual int write(const int *data, std::size_t count)=0
Writes count integer values from array pointed by data.
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Definition element.C:1225
void giveLocationArray(IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
Definition element.C:429
virtual double computeVolumeAround(GaussPoint *gp)
Definition element.h:530
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
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:721
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
int giveNumberOfColumns() const
Returns number of columns of receiver.
void plusDyadUnsym(const FloatArray &a, const FloatArray &b, double dV)
double at(std::size_t i, std::size_t j) const
double weight
Integration weight.
Definition gausspoint.h:108
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition gausspoint.h:190
IntegrationPointStatus * giveMaterialStatus(IntegrationPointStatusIDType key=IPSID_Default)
Definition gausspoint.h:204
Element * giveElement()
Returns corresponding element to receiver.
Definition gausspoint.h:187
double giveLocalEquivalentStrainForAverage()
Returns the local equivalent strain to be averaged.
Definition idmnl1.h:74
double localEquivalentStrainForAverage
Equivalent strain for averaging.
Definition idmnl1.h:59
void setLocalEquivalentStrainForAverage(double ls)
Sets the localEquivalentStrainForAverage to given value.
Definition idmnl1.h:76
int giveLocalNonlocalStiffnessContribution(GaussPoint *gp, IntArray &loc, const UnknownNumberingScheme &s, FloatArray &lcontrib, TimeStep *tStep)
Definition idmnl1.C:567
double computeLocalEquivalentStrain(const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) const
Definition idmnl1.h:145
double computeStressBasedWeight(double cl, double &nx, double &ny, double &ratio, GaussPoint *gp, GaussPoint *jGp, double weight) const
Definition idmnl1.C:220
void giveRemoteNonlocalStiffnessContribution(GaussPoint *gp, IntArray &rloc, const UnknownNumberingScheme &s, FloatArray &rcontrib, TimeStep *tStep)
Definition idmnl1.C:646
double computeStressBasedWeightForPeriodicCell(double cl, double &nx, double &ny, double &ratio, GaussPoint *gp, GaussPoint *jGp) const
Definition idmnl1.C:254
void computeAngleAndSigmaRatio(double &nx, double &ny, double &ratio, GaussPoint *gp, bool &flag) const
Definition idmnl1.C:129
void giveNormalElasticStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition idmnl1.C:801
IDNLMaterial(int n, Domain *d)
Constructor.
Definition idmnl1.C:60
double computeEquivalentStrain(const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) const override
Definition idmnl1.C:285
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
IsotropicDamageMaterial1Status(GaussPoint *g)
Constructor.
Definition idm1.C:1492
double e0
Equivalent strain at stress peak (or a similar parameter).
Definition idm1.h:145
IsotropicDamageMaterial1(int n, Domain *d)
Constructor.
Definition idm1.C:70
MaterialStatus * giveStatus(GaussPoint *gp) const override
Definition idm1.C:1392
EquivStrainType equivStrainType
Parameter specifying the definition of equivalent strain.
Definition idm1.h:194
double ef
Determines ductility -> corresponds to fracturing strain.
Definition idm1.h:147
double give(int aProperty, GaussPoint *gp) const override
Definition idm1.C:1324
double damageFunction(double kappa, GaussPoint *gp) const
Definition idm1.C:954
double e1
Parameters used if softType = 7 (extended smooth damage law).
Definition idm1.h:219
double complianceFunction(double kappa, GaussPoint *gp) const
Definition idm1.C:1134
double stressWork
Density of total work done by stresses on strain increments.
double damage
Damage level of material.
double kappa
Scalar measure of the largest strain level ever reached in material.
double giveTempKappa() const
Returns the temp. scalar measure of the largest strain level.
double giveTempDamage() const
Returns the temp. damage level.
double dissWork
Density of dissipated work.
LinearElasticMaterial * linearElasticMaterial
Reference to bulk (undamaged) material.
LinearElasticMaterial * giveLinearElasticMaterial() const
Returns reference to undamaged (bulk) material.
virtual double give(int aProperty, GaussPoint *gp) const
Definition material.C:51
virtual void initTempStatus(GaussPoint *gp) const
Definition material.C:221
virtual double computeWeightFunction(const double cl, const double distance) const
int averType
Parameter specifying how the weight function should be adjusted due to damage.
NlVariationType nlvar
Parameter specifying the type of nonlocal variation.
void buildNonlocalPointTable(GaussPoint *gp) const
ScalingType scaling
Parameter specifying the type of scaling of nonlocal weight function.
double mm
For "undernonlocal" or "overnonlocal" formulation.
void updateDomainBeforeNonlocAverage(TimeStep *tStep) const
AveragedVarType averagedVar
Parameter specifying the type of averaged (nonlocal) variable.
void modifyNonlocalWeightFunctionAround(GaussPoint *gp) const
std ::vector< localIntegrationRecord > * giveIPIntegrationList(GaussPoint *gp) const
std ::vector< localIntegrationRecord > * giveIntegrationDomainList()
double giveIntegrationScale()
Returns associated integration scale.
virtual int assemble(const IntArray &loc, const FloatMatrix &mat)=0
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int lowerIndx=1, int upperIndx=ALL_STRAINS)=0
const FloatArray & giveTempStrainVector() const
Returns the const pointer to receiver's temporary strain vector.
void letTempStrainVectorBe(const FloatArray &v)
Assigns tempStrainVector to given vector v.
static FloatMatrixF< 6, 6 > giveStressVectorTranformationMtrx(const FloatMatrixF< 3, 3 > &base, bool transpose=false)
virtual void giveStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const
static void computePrincipalValDir(FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode)
void giveStressDependentPartOfStrainVector(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) const
#define OOFEM_ERROR(...)
Definition error.h:79
#define e0_ID
Definition matconst.h:83
#define ef_ID
Definition matconst.h:84
long ContextMode
Definition contextmode.h:43
FloatArrayF< N > max(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
@ principal_stress
For computing principal stresses.
double sum(const FloatArray &x)
@ NonlocalMaterialStiffnessInterfaceType
@ MaterialModelMapperInterfaceType
@ NonlocalMaterialStatusExtensionInterfaceType
@ NonlocalMaterialExtensionInterfaceType
double distance(const FloatArray &x, const FloatArray &y)
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_SPARSE_PROFILE_WIDTH
#define OOFEG_SPARSE_PROFILE_LAYER

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