OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
idmgrad1.C
Go to the documentation of this file.
1 /*
2  *
3  * ##### ##### ###### ###### ### ###
4  * ## ## ## ## ## ## ## ### ##
5  * ## ## ## ## #### #### ## # ##
6  * ## ## ## ## ## ## ## ##
7  * ## ## ## ## ## ## ## ##
8  * ##### ##### ## ###### ## ##
9  *
10  *
11  * OOFEM : Object Oriented Finite Element Code
12  *
13  * Copyright (C) 1993 - 2013 Borek Patzak
14  *
15  *
16  *
17  * Czech Technical University, Faculty of Civil Engineering,
18  * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19  *
20  * This library is free software; you can redistribute it and/or
21  * modify it under the terms of the GNU Lesser General Public
22  * License as published by the Free Software Foundation; either
23  * version 2.1 of the License, or (at your option) any later version.
24  *
25  * This program is distributed in the hope that it will be useful,
26  * but WITHOUT ANY WARRANTY; without even the implied warranty of
27  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28  * Lesser General Public License for more details.
29  *
30  * You should have received a copy of the GNU Lesser General Public
31  * License along with this library; if not, write to the Free Software
32  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33  */
34 
35 #include "idmgrad1.h"
36 #include "gausspoint.h"
37 #include "floatmatrix.h"
38 #include "floatarray.h"
39 #include "mathfem.h"
40 #include "sparsemtrx.h"
42 #include "error.h"
43 #include "nonlocalmaterialext.h"
44 #include "datastream.h"
45 #include "contextioerr.h"
46 #include "stressvector.h"
47 #include "strainvector.h"
48 #include "classfactory.h"
49 
50 namespace oofem {
51 REGISTER_Material(IDGMaterial);
52 
54  //
55  // constructor
56  //
57 {
58  internalLength = 0;
59  averType = 0;
60 }
61 
62 
64 //
65 // destructor
66 //
67 { }
68 
71 {
72  IRResultType result;
74  if ( result != IRRT_OK ) {
75  return result;
76  }
78  if ( result != IRRT_OK ) {
79  return result;
80  }
81  return this->mapper.initializeFrom(ir);
82 }
83 
85 
86 
87 
88 int
90 {
91  return mode == _1dMat || mode == _PlaneStress || mode == _PlaneStrain;
92 }
93 
94 void
96  MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
97 //
98 // Returns characteristic material stiffness matrix of the receiver
99 //
100 {
101  OOFEM_ERROR("Shouldn't be called.");
102 }
103 
104 
105 void
107 {
108  IsotropicDamageMaterialStatus *status = static_cast< IsotropicDamageMaterialStatus * >( this->giveStatus(gp) );
110  double om;
111  om = status->giveTempDamage();
112  om = min(om, maxOmega);
113  answer.resize(1, 1);
114  answer.at(1, 1) = lmat->give('E', gp);
115  answer.times(1.0 - om);
116 }
117 
118 void
120 {
121  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
122 
123  double kappa = status->giveKappa();
124  double tempKappa = status->giveTempKappa();
125  FloatArray totalStrain = status->giveTempStrainVector();
126 
127  answer.resize(1, 1);
128  if ( tempKappa > kappa ) {
129  FloatArray eta;
130  this->computeEta(eta, totalStrain, gp, tStep);
131  answer.at(1, 1) = eta.at(1);
132  }
133 }
134 
135 void
137 {
138  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
140 
141  double damage = status->giveDamage();
142  double tempDamage = status->giveTempDamage();
143  double E = lmat->give('E', gp);
144 
145  answer.resize(1, 1);
146  if ( tempDamage > damage ) {
147  double nlKappa = status->giveNonlocalCumulatedStrain();
148  answer.at(1, 1) = E * status->giveTempStrainVector().at(1);
149  double gPrime = damageFunctionPrime(nlKappa, gp);
150  answer.times(gPrime);
151  } else {
152  answer.zero();
153  }
154 }
155 
156 
157 
158 void
160 {
161  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
162  double tempDamage;
163  if ( mode == ElasticStiffness ) {
164  tempDamage = 0.0;
165  } else {
166  tempDamage = status->giveTempDamage();
167  if ( tempDamage > 0.0 ) {
168  tempDamage = min(tempDamage, maxOmega);
169  }
170  }
171 
172  this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, mode, gp, tStep);
173  answer.times(1.0 - tempDamage);
174 }
175 
176 
177 void
179 {
180  answer.resize(1, 3);
181  answer.zero();
182  if ( mode == TangentStiffness ) {
183  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
184  FloatArray eta;
185  FloatArray totalStrain = status->giveTempStrainVector();
186  this->computeEta(eta, totalStrain, gp, tStep);
187  answer.at(1, 1) = eta.at(1);
188  answer.at(1, 2) = eta.at(2);
189  answer.at(1, 3) = eta.at(3);
190  }
191 }
192 
193 void
195 {
196  answer.resize(3, 1);
197  if ( mode == TangentStiffness ) {
198  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
199  double tempDamage = status->giveTempDamage();
200  tempDamage = min(tempDamage, maxOmega);
201  double nlKappa = status->giveNonlocalCumulatedStrain();
202  double gPrime = this->damageFunctionPrime(nlKappa, gp);
203  FloatArray stress = status->giveTempStressVector();
204 
205  answer.at(1, 1) = stress.at(1) / ( 1 - tempDamage );
206  answer.at(2, 1) = stress.at(2) / ( 1 - tempDamage );
207  answer.at(3, 1) = stress.at(3) / ( 1 - tempDamage );
208  answer.times(gPrime);
209  }
210 }
211 
212 
213 void
215 {
216  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
217  double tempDamage;
218  if ( mode == ElasticStiffness ) {
219  tempDamage = 0.0;
220  } else {
221  tempDamage = status->giveTempDamage();
222  if ( tempDamage > 0.0 ) {
223  tempDamage = min(tempDamage, maxOmega);
224  }
225  }
226 
227  this->giveLinearElasticMaterial()->giveStiffnessMatrix(answer, mode, gp, tStep);
228  answer.times(1.0 - tempDamage);
229 }
230 
231 
232 void
234 {
235  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
236  FloatArray totalStrain = status->giveTempStrainVector();
237  FloatArray eta;
238 
239  answer.resize(1, 4);
240  this->computeEta(eta, totalStrain, gp, tStep);
241  answer.at(1, 1) = eta.at(1);
242  answer.at(1, 2) = eta.at(2);
243  answer.at(1, 3) = eta.at(3);
244  answer.at(1, 4) = eta.at(4);
245 }
246 
247 void
249 {
250  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
251  double damage = status->giveDamage();
252  double tempDamage = status->giveTempDamage();
253  tempDamage = min(tempDamage, maxOmega);
254  answer.resize(4, 1);
255  if ( tempDamage > damage ) {
256  double nlEquivStrain = status->giveNonlocalCumulatedStrain();
257  double gPrime = this->damageFunctionPrime(nlEquivStrain, gp);
258  FloatArray stress = status->giveTempStressVector();
259  answer.at(1, 1) = stress.at(1) / ( 1 - tempDamage );
260  answer.at(2, 1) = stress.at(2) / ( 1 - tempDamage );
261  answer.at(3, 1) = stress.at(3) / ( 1 - tempDamage );
262  answer.at(4, 1) = stress.at(4) / ( 1 - tempDamage );
263  answer.times(gPrime);
264  }
265 }
266 
267 
268 void
270 {
271  if ( averType == 0 ) {
272  answer.resize(1, 1);
273  answer.at(1, 1) = cl;
274  } else if ( averType == 1 ) {
275  answer.resize(1, 1);
276  FloatArray gpCoords;
277  if ( gp->giveElement()->computeGlobalCoordinates( gpCoords, gp->giveNaturalCoordinates() ) == 0 ) {
278  OOFEM_ERROR("computeGlobalCoordinates of GP failed");
279  }
280 
282  answer.at(1, 1) = cl;
283  } else if ( averType == 2 ) {
284  MaterialMode mMode = gp->giveMaterialMode();
285  if ( mMode == _PlaneStress ) {
286  answer.resize(2, 2);
287  answer.zero();
288  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
289  FloatArray stress = status->giveTempStressVector();
290  StressVector fullStress(stress, _PlaneStress);
291  FloatArray sigPrinc;
292  FloatMatrix nPrinc;
293  // get principal stresses (ordered) and principal stress directions
294  fullStress.computePrincipalValDir(sigPrinc, nPrinc);
295  if ( nPrinc.giveNumberOfRows() == 0 ) {
296  nPrinc.resize(2, 2);
297  nPrinc.at(1, 1) = 1;
298  nPrinc.at(2, 2) = 1;
299  }
300 
301  // principal internal lengths
302  double l1 = cl0;
303  double l2 = cl0;
304  double gamma;
305  if ( sigPrinc.at(1) > 0 ) {
306  if ( sigPrinc.at(2) > 0 ) {
307  gamma = beta + ( 1 - beta ) * sigPrinc.at(2) * sigPrinc.at(2) / ( sigPrinc.at(1) * sigPrinc.at(1) );
308  } else {
309  gamma = beta;
310  }
311  } else {
312  gamma = 1;
313  }
314 
315  l2 = l2 * gamma;
316 
317  // compose the internal Length matrix in global coordinates
318  // the first subscript refers to coordinate
319  // the second subscript refers to eigenvalue
320  double n11 = nPrinc.at(1, 1);
321  double n12 = nPrinc.at(1, 2);
322  double n21 = nPrinc.at(2, 1);
323  double n22 = nPrinc.at(2, 2);
324 
325  answer.at(1, 1) = l1 * l1 * n11 * n11 + l2 * l2 * n12 * n12;
326  answer.at(1, 2) = l1 * l1 * n11 * n21 + l2 * l2 * n12 * n22;
327  answer.at(2, 1) = l1 * l1 * n11 * n21 + l2 * l2 * n12 * n22;
328  answer.at(2, 2) = l1 * l1 * n21 * n21 + l2 * l2 * n22 * n22;
329  } else {
330  OOFEM_ERROR("Unknown material mode.");
331  }
332  }
333 }
334 
335 
336 void
338 {
339  MaterialMode mMode = gp->giveMaterialMode();
340  if ( mMode == _PlaneStress ) {
341  answer.resize(5, 5);
342  answer.zero();
343  if ( mode == TangentStiffness ) {
344  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
345  FloatArray stress = status->giveTempStressVector();
346  stress.resize(3);
347  StressVector fullStress(stress, _PlaneStress);
348  FloatArray sigPrinc;
349  FloatMatrix nPrinc;
350  // get principal trial stresses (ordered) and principal stress directions
351  fullStress.computePrincipalValDir(sigPrinc, nPrinc);
352  if ( sigPrinc.at(1) > 0 ) {
353  if ( sigPrinc.at(2) > 0 && sigPrinc.at(1) != sigPrinc.at(2) ) {
354  double gamma = beta + ( 1 - beta ) * sigPrinc.at(2) * sigPrinc.at(2) / ( sigPrinc.at(1) * sigPrinc.at(1) );
355  double dL2dS1 = 4. * gamma * cl * cl * ( beta - 1. ) * sigPrinc.at(2) * sigPrinc.at(2) / ( sigPrinc.at(1) * sigPrinc.at(1) * sigPrinc.at(1) );
356  double dL2dS2 = 4. * gamma * cl * cl * ( 1. - beta ) * sigPrinc.at(2) / ( sigPrinc.at(1) * sigPrinc.at(1) );
357  double dLdN = ( gamma * gamma - 1. ) * cl * cl / ( sigPrinc.at(2) - sigPrinc.at(1) );
358  answer.at(1, 1) = nPrinc.at(1, 1);
359  answer.at(1, 2) = nPrinc.at(1, 2);
360  answer.at(2, 1) = nPrinc.at(2, 1);
361  answer.at(2, 2) = nPrinc.at(2, 2);
362  answer.at(3, 3) = dL2dS1;
363  answer.at(4, 4) = dL2dS2;
364  answer.at(5, 5) = dLdN;
365  }
366  }
367  }
368  } else {
369  OOFEM_ERROR("Unknown material mode.");
370  }
371 }
372 
373 
374 
375 
376 
377 
378 void
379 IDGMaterial :: giveRealStressVectorGrad(FloatArray &answer1, double &answer2, GaussPoint *gp, const FloatArray &totalStrain, double nonlocalCumulatedStrain, TimeStep *tStep)
380 //
381 // returns real stress vector in 3d stress space of receiver according to
382 // previous level of stress and current
383 // strain increment, the only way, how to correctly update gp records
384 //
385 {
386  IDGMaterialStatus *status = static_cast< IDGMaterialStatus * >( this->giveStatus(gp) );
388  FloatArray reducedTotalStrainVector, strain;
389 
390  FloatMatrix de;
391  double f, equivStrain, tempKappa = 0.0, omega = 0.0;
392 
393  this->initTempStatus(gp);
394 
395  // subtract stress independent part
396  // note: eigenStrains (temperature) is not contained in mechanical strain stored in gp
397  // therefore it is necessary to subtract always the total eigen strain value
398  this->giveStressDependentPartOfStrainVector(reducedTotalStrainVector, gp, totalStrain, tStep, VM_Total);
399 
400 
401  // compute equivalent strain
402  this->computeEquivalentStrain(equivStrain, reducedTotalStrainVector, gp, tStep);
403 
404  if ( llcriteria == idm_strainLevelCR ) {
405  // compute value of loading function if strainLevel crit apply
406  f = nonlocalCumulatedStrain - status->giveKappa();
407 
408  if ( f <= 0.0 ) {
409  // damage does not grow
410  tempKappa = status->giveKappa();
411  omega = status->giveDamage();
412  } else {
413  // damage grow
414  this->initDamaged(tempKappa, reducedTotalStrainVector, gp);
415  // evaluate damage parameter
416  this->computeDamageParam(omega, nonlocalCumulatedStrain, reducedTotalStrainVector, gp);
417  }
418  } else if ( llcriteria == idm_damageLevelCR ) {
419  // evaluate damage parameter first
420  tempKappa = nonlocalCumulatedStrain;
421  this->initDamaged(tempKappa, strain, gp);
422  this->computeDamageParam(omega, tempKappa, reducedTotalStrainVector, gp);
423  if ( omega < status->giveDamage() ) {
424  // unloading takes place
425  omega = status->giveDamage();
426  //printf (".");
427  }
428  } else {
429  OOFEM_ERROR("unsupported loading/uloading criteria");
430  }
431 
432 
433  lmat->giveStiffnessMatrix(de, SecantStiffness, gp, tStep);
434  de.times(1.0 - omega);
435  answer1.beProductOf(de, strain);
436  answer2 = equivStrain;
437 
438  // update gp
439 
440  status->letTempStrainVectorBe(totalStrain);
441  status->letTempStressVectorBe(answer1);
442  status->setNonlocalCumulatedStrain(nonlocalCumulatedStrain);
443  status->setTempKappa(tempKappa);
444  status->setTempDamage(omega);
445 }
446 
447 
448 
449 
452 {
453  return new IDGMaterialStatus(1, IDGMaterial :: domain, gp);
454 }
455 
456 
458 { }
459 
460 
462 { }
463 
464 
465 
466 
467 void
469 //
470 // initializes temp variables according to variables form previous equlibrium state.
471 // builds new crackMap
472 //
473 {
475 }
476 
477 
478 
479 void
481 //
482 // updates variables (nonTemp variables describing situation at previous equilibrium state)
483 // after a new equilibrium state has been reached
484 // temporary variables are having values corresponding to newly reched equilibrium.
485 //
486 {
488 }
489 
490 
491 
494 //
495 // saves full information stored in this Status
496 // no temp variables stored
497 //
498 {
499  contextIOResultType iores;
500  // save parent class status
501  if ( ( iores = IsotropicDamageMaterial1Status :: saveContext(stream, mode, obj) ) != CIO_OK ) {
502  THROW_CIOERR(iores);
503  }
504 
505  // write a raw data
506  if ( !stream.write(le) ) {
508  }
509 
510  return CIO_OK;
511 }
512 
515 //
516 // restores full information stored in stream to this Status
517 //
518 {
519  contextIOResultType iores;
520  // read parent class status
521  if ( ( iores = IsotropicDamageMaterialStatus :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
522  THROW_CIOERR(iores);
523  }
524 
525  // read raw data
526  if ( !stream.read(le) ) {
528  }
529 
530  return CIO_OK;
531 }
532 
533 void
535 {
536  MaterialMode mMode = gp->giveMaterialMode();
537  switch ( mMode ) {
538  case _1dMat:
539  give1dStressStiffMtrx(answer, mode, gp, tStep);
540  break;
541  case _PlaneStress:
542  givePlaneStressStiffMtrx(answer, mode, gp, tStep);
543  break;
544  case _PlaneStrain:
545  givePlaneStrainStiffMtrx(answer, mode, gp, tStep);
546  break;
547  default:
548  OOFEM_ERROR("mMode = %d not supported\n", mMode);
549  }
550 }
551 
552 void
554 {
555  MaterialMode mMode = gp->giveMaterialMode();
556  switch ( mMode ) {
557  case _1dMat:
558  give1dKappaMatrix(answer, mode, gp, tStep);
559  break;
560  case _PlaneStress:
561  givePlaneStressKappaMatrix(answer, mode, gp, tStep);
562  break;
563  case _PlaneStrain:
564  givePlaneStrainKappaMatrix(answer, mode, gp, tStep);
565  break;
566  default:
567  OOFEM_ERROR("mMode = %d not supported\n", mMode);
568  }
569 }
570 
571 void
573 {
574  MaterialMode mMode = gp->giveMaterialMode();
575  switch ( mMode ) {
576  case _1dMat:
577  give1dGprime(answer, mode, gp, tStep);
578  break;
579  case _PlaneStress:
580  givePlaneStressGprime(answer, mode, gp, tStep);
581  break;
582  case _PlaneStrain:
583  givePlaneStrainGprime(answer, mode, gp, tStep);
584  break;
585  default:
586  OOFEM_ERROR("mMode = %d not supported\n", mMode);
587  }
588 }
589 
590 void
592 {
593  MaterialMode mMode = gp->giveMaterialMode();
594  switch ( mMode ) {
595  case _1dMat:
596  giveInternalLength(answer, mode, gp, tStep);
597  break;
598  case _PlaneStress:
599  giveInternalLength(answer, mode, gp, tStep);
600  break;
601  case _PlaneStrain:
602  giveInternalLength(answer, mode, gp, tStep);
603  break;
604  default:
605  OOFEM_ERROR("mMode = %d not supported\n", mMode);
606  }
607 }
608 
609 void
611 {
612  MaterialMode mMode = gp->giveMaterialMode();
613  switch ( mMode ) {
614  case _PlaneStress:
615  giveInternalLengthDerivative(answer, mode, gp, tStep);
616  break;
617  default:
618  OOFEM_ERROR("mMode = %d not supported\n", mMode);
619  }
620 }
621 } // end namespace oofem
enum oofem::IsotropicDamageMaterial::loaUnloCriterium llcriteria
double giveDamage()
Returns the last equilibrated damage level.
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition: gausspoint.h:191
double le
Characteristic element length, computed when damage initialized from direction of maximum positive pr...
void letTempStrainVectorBe(const FloatArray &v)
Assigns tempStrainVector to given vector v.
Definition: structuralms.h:137
void giveInternalLength(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: idmgrad1.C:269
GaussPoint * gp
Associated integration point.
Class and object Domain.
Definition: domain.h:115
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
double cl0
Initial(user defined) characteristic length of the nonlocal model (its interpretation depends on the ...
void givePlaneStrainGprime(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: idmgrad1.C:248
virtual void givePDGradMatrix_kk(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Right lower block.
Definition: idmgrad1.C:591
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
Specialization of a floating point array for representing a stress state.
Definition: stressvector.h:48
virtual void initDamaged(double kappa, FloatArray &totalStrainVector, GaussPoint *gp)
Performs initialization, when damage first appear.
Definition: idm1.C:1159
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
virtual void givePDGradMatrix_uu(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Left upper block.
Definition: idmgrad1.C:534
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: idmgrad1.C:493
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
void computePrincipalValDir(FloatArray &answer, FloatMatrix &dir) const
Computes principal values and directions of receiver vector.
Definition: stressvector.C:182
IDGMaterialStatus(int n, Domain *d, GaussPoint *g)
Definition: idmgrad1.C:457
IDGMaterial(int n, Domain *d)
Constructor.
Definition: idmgrad1.C:53
double cl
Characteristic length of the nonlocal model (its interpretation depends on the type of weight functio...
static MMAContainingElementProjection mapper
Mapper used to map internal variables in adaptivity.
Definition: idm1.h:236
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
Definition: idmgrad1.C:468
virtual void giveRealStressVectorGrad(FloatArray &answer1, double &answer2, GaussPoint *gp, const FloatArray &totalStrain, double nonlocalCumulatedStrain, TimeStep *tStep)
gradient - based giveRealStressVector
Definition: idmgrad1.C:379
Element * giveElement()
Returns corresponding element to receiver.
Definition: gausspoint.h:188
General IO error.
virtual MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
Definition: idmgrad1.C:451
This class implements associated Material Status to IsotropicDamageMaterial1.
Definition: idm1.h:117
void giveStressDependentPartOfStrainVector(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode)
Method for subtracting from reduced space strain vector its stress-independent parts (caused by tempe...
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual int read(int *data, int count)=0
Reads count integer values into array pointed by data.
MatResponseMode
Describes the character of characteristic material matrix.
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
This class is a abstract base class for all linear elastic material models in a finite element proble...
virtual int hasMaterialModeCapability(MaterialMode mode)
Tests if material supports material mode.
Definition: idmgrad1.C:89
double beta
Parameters specifying how the length scale parameter l is adjusted.
Definition: idmgrad1.h:55
Material status for gradient-enhanced isotropic damage model for concrete in tension.
Definition: idmgrad1.h:104
virtual void givePDGradMatrix_LD(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Stress-based averaging.
Definition: idmgrad1.C:610
virtual int write(const int *data, int count)=0
Writes count integer values from array pointed by data.
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of material property &#39;aProperty&#39;.
Definition: material.C:52
virtual void giveStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix for giveRealStressVector of receiver in given integration point...
void give1dKappaMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: idmgrad1.C:119
void giveStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix for giveRealStressVector of receiver in given integration point...
Definition: idmgrad1.C:95
virtual IRResultType initializeFrom(InputRecord *ir)
#define E(p)
Definition: mdm.C:368
void setTempDamage(double newDamage)
Sets the temp damage level to given value.
double damageFunctionPrime(double kappa, GaussPoint *gp)
Returns the value of compliance parameter corresponding to a given value of the damage-driving variab...
Definition: idm1.C:1039
double maxOmega
Maximum limit on omega. The purpose is elimination of a too compliant material which may cause conver...
virtual void updateYourself(TimeStep *)
Update equilibrium history variables according to temp-variables.
Definition: idmgrad1.C:480
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: idm1.C:112
double internalLength
Length scale parameter.
Definition: idmgrad1.h:51
void times(double f)
Multiplies receiver by factor f.
Definition: floatmatrix.C:1594
virtual void computeDamageParam(double &omega, double kappa, const FloatArray &strain, GaussPoint *gp)
Computes the value of damage parameter omega, based on given value of equivalent strain.
Definition: idm1.C:831
Material interface for gradient material models.
double giveKappa()
Returns the last equilibrated scalar measure of the largest strain level.
void givePlaneStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane stress stiffness matrix of receiver.
Definition: idmgrad1.C:159
virtual void givePDGradMatrix_ku(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Left lower block.
Definition: idmgrad1.C:553
virtual ~IDGMaterialStatus()
Definition: idmgrad1.C:461
void givePlaneStressKappaMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: idmgrad1.C:178
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
void giveDistanceBasedCharacteristicLength(const FloatArray &gpCoords)
Provides the distance based interaction radius This function is called when averType is set to 1...
const FloatArray & giveTempStrainVector() const
Returns the const pointer to receiver&#39;s temporary strain vector.
Definition: structuralms.h:115
virtual void computeEta(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
Computes derivative of the equivalent strain with regards to strain.
Definition: idm1.C:606
int averType
Parameter specifiying averaging type(0 - classical, 1 - distance based, 2 - stress based) ...
Definition: idmgrad1.h:53
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
LinearElasticMaterial * giveLinearElasticMaterial()
Returns reference to undamaged (bulk) material.
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: idmgrad1.C:514
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Returns material status of receiver in given integration point.
Definition: idm1.C:1354
Abstract base class representing a material status information.
Definition: matstatus.h:84
double giveTempKappa()
Returns the temp. scalar measure of the largest strain level.
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Class representing vector of real numbers.
Definition: floatarray.h:82
void givePlaneStressGprime(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: idmgrad1.C:194
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual void setNonlocalCumulatedStrain(double nonlocalCumulatedStrain)
Definition: idmgrad1.h:119
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
void letTempStressVectorBe(const FloatArray &v)
Assigns tempStressVector to given vector v.
Definition: structuralms.h:135
double giveTempDamage()
Returns the temp. damage level.
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
Class representing the general Input Record.
Definition: inputrecord.h:101
void give1dGprime(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: idmgrad1.C:136
void giveInternalLengthDerivative(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: idmgrad1.C:337
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
This class implements a simple local isotropic damage model for concrete in tension.
Definition: idm1.h:137
virtual IRResultType initializeFrom(InputRecord *ir)
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
int min(int i, int j)
Returns smaller value from two given decimals.
Definition: mathfem.h:59
REGISTER_Material(DummyMaterial)
void givePlaneStrainKappaMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Definition: idmgrad1.C:233
const FloatArray & giveTempStressVector() const
Returns the const pointer to receiver&#39;s temporary stress vector.
Definition: structuralms.h:117
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: idmgrad1.C:70
virtual void initTempStatus(GaussPoint *gp)
Initializes temporary variables stored in integration point status at the beginning of new time step...
Definition: material.C:267
int giveNumberOfRows() const
Returns number of rows of receiver.
Definition: floatmatrix.h:156
void setTempKappa(double newKappa)
Sets the temp scalar measure of the largest strain level to given value.
This class implements associated Material Status to IsotropicDamageMaterial.
virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep)
Method for computing 1d stiffness matrix of receiver.
Definition: idmgrad1.C:106
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
virtual void computeEquivalentStrain(double &kappa, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
Computes the equivalent strain measure from given strain vector (full form).
Definition: idm1.C:472
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: element.C:1207
virtual void givePDGradMatrix_uk(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Right upper block.
Definition: idmgrad1.C:572
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
virtual ~IDGMaterial()
Destructor.
Definition: idmgrad1.C:63
void givePlaneStrainStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane strain stiffness matrix of receiver.
Definition: idmgrad1.C:214
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
virtual double giveNonlocalCumulatedStrain()
Definition: idmgrad1.h:118

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:29 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011