OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
anisodamagemodel.h
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 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 #ifndef anisodamagemodel_h
36 #define anisodamagemodel_h
37 
38 // this turns on or off a bunch of internal variables
39 // that allow tracing the distribution of dissipated energy
40 // (can be turned off if such information is not needed)
41 #define keep_track_of_dissipated_energy
42 
43 #include "material.h"
44 #include "../sm/Materials/structuralmaterial.h"
46 #include "../sm/Materials/structuralms.h"
48 
50 
51 // @todo add input parametres which are read from input file here
52 #define _IFT_AnisotropicDamageMaterial_Name "adm"
53 #define _IFT_AnisotropicDamageMaterial_equivStrainType "equivstraintype"
54 #define _IFT_AnisotropicDamageMaterial_damageLawType "damlaw"
55 #define _IFT_AnisotropicDamageMaterial_kappa0 "kappa0"
56 #define _IFT_AnisotropicDamageMaterial_kappaf "kappaf"
57 #define _IFT_AnisotropicDamageMaterial_aA "aa"
58 
59 
60 namespace oofem {
61 class GaussPoint;
62 
68 {
69 protected:
71  double kappa;
73  double tempKappa;
79  double tempStrainZ;
81  double strainZ;
83  int flag;
84  int tempFlag;
85  double storedFactor;
87 
88 #ifdef keep_track_of_dissipated_energy
89  double stressWork;
94  double dissWork;
96  double tempDissWork;
97 #endif
98 
99 public:
104 
105  virtual void printOutputAt(FILE *file, TimeStep *tStep);
106 
108  double giveKappa() { return kappa; }
110  double giveTempKappa() { return tempKappa; }
112  void setTempKappa(double newKappa) { tempKappa = newKappa; }
114  const FloatMatrix &giveDamage() { return damage; }
116  const FloatMatrix &giveTempDamage() { return tempDamage; }
118  void setTempDamage(const FloatMatrix &d) { tempDamage = d; }
120  double giveStrainZ() { return strainZ; }
122  double giveTempStrainZ() { return tempStrainZ; }
124  void setTempStrainZ(double newStrainZ) { tempStrainZ = newStrainZ; }
126  int giveFlag() { return flag; }
128  void setTempFlag(int newflag) { tempFlag = newflag; }
130  int giveTempFlag() { return tempFlag; }
132  double giveStoredFactor() { return storedFactor; }
134  void setStoredFactor(double newStoredFactor) { storedFactor = newStoredFactor; }
138  void setTempStoredFactor(double newTempStoredFactor) { tempStoredFactor = newTempStoredFactor; }
139 
140 
141 
142 #ifdef keep_track_of_dissipated_energy
143  double giveStressWork() { return stressWork; }
146  double giveTempStressWork() { return tempStressWork; }
148  void setTempStressWork(double w) { tempStressWork = w; }
150  double giveDissWork() { return dissWork; }
152  double giveTempDissWork() { return tempDissWork; }
154  void setTempDissWork(double w) { tempDissWork = w; }
156  void computeWork(GaussPoint *gp);
157 #endif
158 
159  // definition
160  virtual const char *giveClassName() const { return "AnisotropicDamageMaterialModelStatus"; }
161 
162  virtual void initTempStatus();
163  virtual void updateYourself(TimeStep *tStep);
164 
165  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
166  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
167 };
168 
169 
180 {
181 protected:
182 
186  double E;
188  double nu;
190  double kappa0;
192  double kappaf;
194  double aA;
196  enum EquivStrainType { EST_Unknown, EST_Mazars, EST_Rankine_Smooth, EST_Rankine_Standard, EST_ElasticEnergy, EST_ElasticEnergyPositiveStress, EST_ElasticEnergyPositiveStrain, EST_Mises, EST_Griffith };
200  enum DamLawType { DLT_Unknown, DLT_Desmorat1, DLT_Desmorat2, DLT_Linear, DLT_Exponential };
203 
204 public:
208  virtual ~AnisotropicDamageMaterial();
209 
210  virtual int hasNonLinearBehaviour() { return 1; }
211  virtual int hasMaterialModeCapability(MaterialMode mode);
212 
213  // identification and auxiliary functions
214  virtual const char *giveClassName() const { return "AnisotropicDamageMaterial"; }
215  virtual const char *giveInputRecordName() const { return _IFT_AnisotropicDamageMaterial_Name; }
216 
218  IsotropicLinearElasticMaterial *giveLinearElasticMaterial() { return linearElasticMaterial; }
219 
221  virtual void giveRealStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep);
222  void computePrincValDir2D(double &D1, double &D2, double &c, double &s, double Dx, double Dy, double Dxy);
223  bool checkPrincVal2D(double Dx, double Dy, double Dxy);
224  void computeDamage(FloatMatrix &tempDamage, const FloatMatrix &damage, double kappa, double eps1, double eps2, double ceps, double seps, double epsZ);
225  double computeTraceD(double equivStrain);
226  double computeOutOfPlaneStrain(const FloatArray &inplaneStrain, const FloatMatrix &dam, bool tens_flag);
227  double computeDimensionlessOutOfPlaneStress(const FloatArray &inplaneStrain, double epsZ, const FloatMatrix &dam);
228  void computeInplaneStress(FloatArray &inplaneStress, const FloatArray &inplaneStrain, double epsZ, const FloatMatrix &dam);
229 
230 
232  double obtainAlpha1(FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, double damageThreshold);
234  double obtainAlpha2(FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, FloatMatrix ProjMatrix, double damageThreshold);
236  double obtainAlpha3(FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, FloatArray vec3, double damageThreshold);
237 
238  double checkSymmetry(FloatMatrix matrix);
239 
240  void correctBigValues(FloatMatrix &matrix);
241 
242  double computeTraceD(FloatMatrix tempDamageTensor, FloatMatrix strainTensor, GaussPoint *gp);
243 
244  double computeCorrectionFactor(FloatMatrix tempDamageTensor, FloatMatrix strainTensor, GaussPoint *gp);
245 
246  virtual void give3dMaterialStiffnessMatrix(FloatMatrix &answer,
247  MatResponseMode mode,
248  GaussPoint *gp,
249  TimeStep *tStep);
250 
251  virtual void giveRealStressVector(FloatArray &answer, GaussPoint *gp,
252  const FloatArray &reducedStrain, TimeStep *tStep);
253 
254  virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
255  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
256  virtual void giveRealStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
257  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
258  virtual void giveRealStressVector_StressControl(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep)
259  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
260  virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
261  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
262 
263  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *atTime);
264 
265  // virtual InternalStateValueType giveIPValueType(InternalStateType type);
266  //virtual void giveThermalDilatationVector(FloatArray &answer, GaussPoint *, TimeStep *);
267 
275  virtual void computeEquivalentStrain(double &kappa, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep);
276 
278  virtual void giveInputRecord(DynamicInputRecord &input);
279  virtual void computeDamageTensor(FloatMatrix &answer, GaussPoint *gp,
280  const FloatArray &totalStrain, double equivStrain,
281  TimeStep *atTime);
282 
284 
285 protected:
286 
287 
288  virtual void givePlaneStressStiffMtrx(FloatMatrix &answer, MatResponseMode mmode,
289  GaussPoint *gp,
290  TimeStep *tStep);
291 
292  virtual void givePlaneStrainStiffMtrx(FloatMatrix &answer, MatResponseMode mmode,
293  GaussPoint *gp,
294  TimeStep *tStep);
295 
296  virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mmode,
297  GaussPoint *gp,
298  TimeStep *tStep);
299  virtual void computePlaneStressStrain(FloatMatrix &answer, FloatMatrix damageTensor, FloatArray totalStrain, GaussPoint *gp,
300  TimeStep *atTime);
301  virtual void computePlaneStressSigmaZ(double &answer, FloatMatrix damageTensor, FloatArray reducedTotalStrainVector,
302  double epsZ, GaussPoint *gp, TimeStep *atTime);
303 
304  /* virtual void computeDamageTensor(FloatMatrix &answer, GaussPoint *gp,
305  * const FloatArray &totalStrain, double equivStrain,
306  * TimeStep *atTime);
307  */
308 
309  virtual void computeSecantOperator(FloatMatrix &answer, FloatMatrix strainTensor,
310  FloatMatrix damageTensor, GaussPoint *gp);
311 
312  double computeK(GaussPoint *gp);
313 
314  double computeKappa(FloatMatrix damageTensor);
315 };
316 } // end namespace oofem
317 #endif // anisodamagemodel_h
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
double aA
Damage parameter a*A, needed to obtain Kappa(trD), according to eq. 33 in the paper mentioned above...
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
This class implements associated Material Status to AnisotropicDamageMaterial.
int giveTempFlag()
Returns the value of the temporary value of flag.
GaussPoint * gp
Associated integration point.
Class and object Domain.
Definition: domain.h:115
IsotropicLinearElasticMaterial * linearElasticMaterial
Reference to bulk (undamaged) material.
virtual const char * giveClassName() const
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
double tempDissWork
Non-equilibrated density of dissipated work.
double tempStressWork
Non-equilibrated density of total work done by stresses on strain increments.
AnisotropicDamageMaterialStatus(int n, Domain *d, GaussPoint *g)
Constructor.
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
FloatMatrix tempDamage
Non-equilibrated second order damage tensor.
Class representing anisotropic damage model.
This class implements a structural material status information.
Definition: structuralms.h:65
double nu
Poisson's ratio.
virtual ~AnisotropicDamageMaterialStatus()
Destructor.
double giveStressWork()
Returns the density of total work of stress on strain increments.
double giveTempDissWork()
Returns the density of temp dissipated work.
double giveTempStoredFactor()
Returns the last Temp Stored Factor.
double giveStoredFactor()
Returns the last Stored Factor.
DamLawType damageLawType
Parameter specifying the damage law.
void setTempFlag(int newflag)
Sets the value of the temporary value of flag.
#define _IFT_AnisotropicDamageMaterial_Name
double stressWork
Density of total work done by stresses on strain increments.
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Definition: femcmpnn.C:77
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
FloatMatrix damage
Second order damage tensor.
const FloatMatrix & giveDamage()
Returns the last equilibrated second order damage tensor.
double giveTempStressWork()
Returns the temp density of total work of stress on strain increments.
virtual void giveRealStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_3d.
double dissWork
Density of dissipated work.
virtual void giveRealStressVector_StressControl(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep)
Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· ...
double tempKappa
Non-equilibrated scalar measure of the largest strain level.
void setTempStrainZ(double newStrainZ)
Sets the temp scalar measure of the out-of-plane strain to given value (for 2dPlaneStress mode)...
void setTempDamage(const FloatMatrix &d)
Assigns temp. damage tensor to given tensor d.
double giveKappa()
Returns the last equilibrated scalar measure of the largest strain level.
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: matstatus.h:140
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Print receiver's output to given stream.
This class implements an isotropic linear elastic material in a finite element problem.
virtual const char * giveInputRecordName() const
void setTempDissWork(double w)
Sets the density of dissipated work to given value.
void setStoredFactor(double newStoredFactor)
Sets the Stored Factor to given value .
double giveTempStrainZ()
Returns the temp scalar measure of the out-of-plane strain to given value (for 2dPlaneStress mode)...
double giveTempKappa()
Returns the temp. scalar measure of the largest strain level.
double strainZ
Non-equilibrated out-of-plane value for 2dPlaneStress mode.
EquivStrainType
Type characterizing the algorithm used to compute equivalent strain measure.
void computeWork(GaussPoint *gp)
Computes the increment of total stress work and of dissipated work.
IsotropicLinearElasticMaterial * giveLinearElasticMaterial()
Returns reference to undamaged (bulk) material.
double giveDissWork()
Returns the density of dissipated work.
int giveFlag()
Returns the value of the flag.
Abstract base class representing a material status information.
Definition: matstatus.h:84
Class representing vector of real numbers.
Definition: floatarray.h:82
void setTempKappa(double newKappa)
Sets the temp scalar measure of the largest strain level to given value.
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
double kappaf
Damage parameter kappa_f (in the paper denoted as "a")
virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
double kappa0
Damage threshold kappa0, as defined in the paper mentioned above.
MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual const char * giveClassName() const
void setTempStoredFactor(double newTempStoredFactor)
Sets the Temp Stored Factor to given value .
double kappa
Scalar measure of the largest strain level ever reached in material.
Class representing the a dynamic Input Record.
DamLawType
Type characterizing the damage law.
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
Abstract base class for all "structural" constitutive models.
double tempStrainZ
Out-of-plane value for 2dPlaneStress mode.
const FloatMatrix & giveTempDamage()
Returns the temp. second order damage tensor.
virtual int hasNonLinearBehaviour()
Returns nonzero if receiver is non linear.
double E
Young's modulus.
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
int flag
This flag turns into 1 and remains 1 when the trace of the damage tensor is >1 in compression (tr(str...
the oofem namespace is to define a context or scope in which all oofem names are defined.
void setTempStressWork(double w)
Sets the density of total work of stress on strain increments to given value.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
EquivStrainType equivStrainType
Parameter specifying the definition of equivalent strain.
double giveStrainZ()
Returns the last equilibrated scalar measure of the out-of-plane strain to given value (for 2dPlaneSt...

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:27 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011