OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
rheoChM.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 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 #ifndef rheochm_h
36 #define rheochm_h
37 
39 #define keep_track_of_strains
40 
41 #include "../sm/Materials/structuralmaterial.h"
43 #include "floatarray.h"
44 #include "floatmatrix.h"
45 
46 #include "matconst.h"
47 #include "../sm/Elements/structuralelement.h"
48 #include "../sm/Materials/structuralms.h"
49 
51 
52 #define _IFT_RheoChainMaterial_n "n"
53 #define _IFT_RheoChainMaterial_alphaOne "a1"
54 #define _IFT_RheoChainMaterial_alphaTwo "a2"
55 #define _IFT_RheoChainMaterial_lattice "lattice"
56 #define _IFT_RheoChainMaterial_relmatage "relmatage"
57 #define _IFT_RheoChainMaterial_begoftimeofinterest "begoftimeofinterest"
58 #define _IFT_RheoChainMaterial_endoftimeofinterest "endoftimeofinterest"
59 #define _IFT_RheoChainMaterial_timefactor "timefactor"
60 #define _IFT_RheoChainMaterial_talpha "talpha"
61 #define _IFT_RheoChainMaterial_preCastingTimeMat "precastingtimemat"
62 
63 
64 namespace oofem {
65 #define MNC_NPOINTS 30
66 #define TIME_DIFF 1.e-10
67 
72 {
73 protected:
75  int nUnits;
77  std :: vector< FloatArray >hiddenVars;
78  std :: vector< FloatArray >tempHiddenVars;
79 
85 
86 #ifdef keep_track_of_strains
87  double thermalStrain;
89 #endif
90 
91 public:
92  RheoChainMaterialStatus(int n, Domain *d, GaussPoint *g, int nunits);
93  virtual ~RheoChainMaterialStatus();
94 
95  virtual void printOutputAt(FILE *file, TimeStep *tStep);
96 
97  virtual const FloatArray &giveViscoelasticStressVector() const { return stressVector; }
98 
99  FloatArray &giveHiddenVarsVector(int i) { return hiddenVars [ i - 1 ]; }
100  FloatArray &giveTempHiddenVarsVector(int i) { return tempHiddenVars [ i - 1 ]; }
102  void letTempHiddenVarsVectorBe(int i, FloatArray &valueArray);
103 
105  void setShrinkageStrainVector(FloatArray src) { shrinkageStrain = std :: move(src); }
106 
107  virtual void initTempStatus();
108  virtual void updateYourself(TimeStep *tStep);
109 
110  virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj = NULL);
111  virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj = NULL);
112 
113 
114 #ifdef keep_track_of_strains
115  void setTempThermalStrain(double src) { tempThermalStrain = src; }
116  double giveTempThermalStrain(void) { return tempThermalStrain; }
117  double giveThermalStrain(void) { return thermalStrain; }
118 #endif
119 
120  // definition
121  virtual const char *giveClassName() const { return "RheoChainMaterialStatus"; }
122 };
123 
124 
131 {
132 protected:
134  double talpha;
136  int nUnits;
138  double relMatAge;
139 
140  bool lattice;
142  double nu;
144  double alphaOne, alphaTwo;
146  double EparValTime;
147 
149  double begOfTimeOfInterest; // local one or taken from e-model
151  double endOfTimeOfInterest; // local one or taken from e-model
156  //FloatArray relaxationTimes;
161 
167  double timeFactor;
168 
170  //double zeroStiffness;
171 
173 
174 public:
175  RheoChainMaterial(int n, Domain *d);
176  virtual ~RheoChainMaterial();
177 
178  virtual void giveRealStressVector(FloatArray &answer, GaussPoint *gp,
179  const FloatArray &reducedStrain, TimeStep *tStep);
180 
181  virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
182  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
183  virtual void giveRealStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
184  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
185  virtual void giveRealStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
186  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
187  virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
188  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
189  virtual void giveRealStressVector_2dBeamLayer(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
190  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
191  virtual void giveRealStressVector_PlateLayer(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
192  { this->giveRealStressVector(answer, gp, reducedE, tStep); }
193 
194  virtual void giveThermalDilatationVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep);
195 
196  /* virtual void giveThermalDilatationVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
197  * { answer.clear(); }*/
198 
200  virtual double giveEModulus(GaussPoint *gp, TimeStep *tStep) = 0;
201 
202  /* virtual double giveIncrementalModulus(GaussPoint *gp, TimeStep *tStep) {
203  * if ( (tStep->giveIntrinsicTime() < this->castingTime) && ( this->zeroStiffness > 0. ) ) {
204  * return this->zeroStiffness;
205  * } else {
206  * return this->giveEModulus(gp, tStep);
207  * }
208  * }*/
209 
211  virtual void computeCharCoefficients(FloatArray &answer, double tPrime, GaussPoint *gp, TimeStep *tStep) = 0;
212 
213  // identification and auxiliary functions
214  virtual int hasNonLinearBehaviour() { return 0; }
215  virtual int hasMaterialModeCapability(MaterialMode mode);
216 
217  virtual int hasCastingTimeSupport() { return 1; }
218 
219  virtual const char *giveClassName() const { return "RheoChainMaterial"; }
221 
222  virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep);
223 
224  // store & restore context functions
225  virtual contextIOResultType saveIPContext(DataStream &stream, ContextMode mode, GaussPoint *gp);
226  virtual contextIOResultType restoreIPContext(DataStream &stream, ContextMode mode, GaussPoint *gp);
227 
228  virtual void give3dMaterialStiffnessMatrix(FloatMatrix &answer,
229  MatResponseMode mode,
230  GaussPoint *gp,
231  TimeStep *tStep);
232 
233  virtual void givePlaneStressStiffMtrx(FloatMatrix &answer,
234  MatResponseMode mmode,
235  GaussPoint *gp,
236  TimeStep *tStep);
237  virtual void givePlaneStrainStiffMtrx(FloatMatrix &answer,
238  MatResponseMode mmode,
239  GaussPoint *gp,
240  TimeStep *tStep);
241  virtual void give1dStressStiffMtrx(FloatMatrix &answer,
242  MatResponseMode mmode,
243  GaussPoint *gp,
244  TimeStep *tStep);
245 
246  virtual void give2dLatticeStiffMtrx(FloatMatrix &answer,
247  MatResponseMode mmode,
248  GaussPoint *gp,
249  TimeStep *tStep);
250 
251  virtual void give3dLatticeStiffMtrx(FloatMatrix &answer,
252  MatResponseMode mmode,
253  GaussPoint *gp,
254  TimeStep *tStep);
255 
256  virtual void computeStressIndependentStrainVector(FloatArray &answer,
257  GaussPoint *gp, TimeStep *tStep, ValueModeType mode);
258 
267  virtual void giveShrinkageStrainVector(FloatArray &answer,
268  GaussPoint *gp,
269  TimeStep *tStep,
270  ValueModeType mode)
271  { answer.clear(); }
272 
273  // Note: must take LoadResponseMode into account
282  virtual void giveEigenStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode) { }
283 
284  virtual MaterialStatus *CreateStatus(GaussPoint *gp) const;
285 
286  double giveAlphaOne() const { return this->alphaOne; }
287  double giveAlphaTwo() const { return this->alphaTwo; }
288 
290  virtual double computeCreepFunction(double t, double t_prime, GaussPoint *gp, TimeStep *tStep) = 0;
291 
292  virtual bool isActivated(TimeStep *tStep) {
293  if ( this->preCastingTimeMat > 0 ) {
294  return true;
295  } else {
296  return Material :: isActivated(tStep);
297  }
298  }
299 
300 
301 protected:
306  virtual int hasIncrementalShrinkageFormulation() { return 0; }
307 
318  static void generateLogTimeScale(FloatArray &answer, double from, double to, int nsteps);
319  const FloatArray &giveDiscreteTimes();
320 
336  void computeDiscreteRelaxationFunction(FloatArray &answer, const FloatArray &tSteps, double t0, double tr, GaussPoint *gp, TimeStep *tSte);
337 
339  void giveUnitComplianceMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep);
341  void giveUnitStiffnessMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep);
342 
344  virtual void updateEparModuli(double tPrime, GaussPoint *gp, TimeStep *tStep);
345 
347  double giveEparModulus(int iChain);
348 
350  virtual void computeCharTimes();
351 
353  double giveCharTime(int) const;
354 
356  virtual double giveCharTimeExponent(int i) { return 1.0; }
357 
359  LinearElasticMaterial *giveLinearElasticMaterial();
360 
362  double giveEndOfTimeOfInterest();
363 
374  void computeTrueStressIndependentStrainVector(FloatArray &answer, GaussPoint *gp,
375  TimeStep *tStep, ValueModeType mode);
376 };
377 } // end namespace oofem
378 #endif // rheochm_h
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
double relMatAge
Physical age of the material at simulation time = 0.
Definition: rheoChM.h:138
virtual bool isActivated(TimeStep *tStep)
Definition: rheoChM.h:292
GaussPoint * gp
Associated integration point.
Class and object Domain.
Definition: domain.h:115
virtual const char * giveClassName() const
Definition: rheoChM.h:219
double nu
Poisson&#39;s ratio (assumed to be constant, unaffected by creep).
Definition: rheoChM.h:142
std::vector< FloatArray > tempHiddenVars
Definition: rheoChM.h:78
virtual ~RheoChainMaterialStatus()
Definition: rheoChM.C:735
virtual void updateYourself(TimeStep *tStep)
Update equilibrium history variables according to temp-variables.
Definition: rheoChM.C:781
virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Definition: rheoChM.h:187
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Print receiver&#39;s output to given stream.
Definition: rheoChM.C:751
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
Defines several material constant (respective their representative number).
double giveThermalStrain(void)
Definition: rheoChM.h:117
virtual int hasIncrementalShrinkageFormulation()
If only incremental shrinkage strain formulation is provided, then total shrinkage strain must be tra...
Definition: rheoChM.h:306
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual void giveShrinkageStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode)
Computes, for the given integration point, the strain vector induced by stress-independent shrinkage...
Definition: rheoChM.h:267
This class implements a structural material status information.
Definition: structuralms.h:65
double begOfTimeOfInterest
Time from which the model should give a good approximation. Optional field. Default value is 0...
Definition: rheoChM.h:149
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: rheoChM.C:802
This class implements associated Material Status to RheoChainMaterial.
Definition: rheoChM.h:71
FloatArray EparVal
Partial moduli of individual units.
Definition: rheoChM.h:155
int preCastingTimeMat
Stiffness at time less than casting time - optional parameter, negative by default.
Definition: rheoChM.h:172
FloatArray stressVector
Equilibrated stress vector in reduced form.
Definition: structuralms.h:71
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: rheoChM.C:829
double giveAlphaTwo() const
Definition: rheoChM.h:287
double giveTempThermalStrain(void)
Definition: rheoChM.h:116
void setTempThermalStrain(double src)
Definition: rheoChM.h:115
This class implements a rheologic chain model describing a viscoelastic material. ...
Definition: rheoChM.h:130
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
MatResponseMode
Describes the character of characteristic material matrix.
virtual int hasNonLinearBehaviour()
Returns nonzero if receiver is non linear.
Definition: rheoChM.h:214
This class is a abstract base class for all linear elastic material models in a finite element proble...
FloatArray charTimes
Characteristic times of individual units (relaxation or retardation times).
Definition: rheoChM.h:158
void letTempHiddenVarsVectorBe(int i, FloatArray &valueArray)
Definition: rheoChM.C:739
double timeFactor
Scaling factor transforming the simulation time units into days (gives the number of simulation time ...
Definition: rheoChM.h:167
RheoChainMaterialStatus(int n, Domain *d, GaussPoint *g, int nunits)
Definition: rheoChM.C:718
double endOfTimeOfInterest
Time (age???) up to which the model should give a good approximation.
Definition: rheoChM.h:151
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: matstatus.h:140
virtual int hasCastingTimeSupport()
Tests if material supports casting time.
Definition: rheoChM.h:217
FloatArray * letHiddenVarsVectorBe(int i, FloatArray *)
FloatArray shrinkageStrain
Total shrinkage strain (needed only when the shrinkage evolution is described in the incremental form...
Definition: rheoChM.h:84
virtual void giveRealStressVector_PlateLayer(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Definition: rheoChM.h:191
FloatArray & giveHiddenVarsVector(int i)
Definition: rheoChM.h:99
std::vector< FloatArray > hiddenVars
Hidden (internal) variables, the meaning of which depends on the type of chain.
Definition: rheoChM.h:77
Abstract base class representing a material status information.
Definition: matstatus.h:84
virtual double giveCharTimeExponent(int i)
Exponent to be used with the char time of a given unit, usually = 1.0.
Definition: rheoChM.h:356
Class representing vector of real numbers.
Definition: floatarray.h:82
FloatArray & giveTempHiddenVarsVector(int i)
Definition: rheoChM.h:100
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
virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Definition: rheoChM.h:181
double giveAlphaOne() const
Definition: rheoChM.h:286
Class representing the general Input Record.
Definition: inputrecord.h:101
FloatArray * giveShrinkageStrainVector()
Definition: rheoChM.h:104
virtual void initTempStatus()
Initializes the temporary internal variables, describing the current state according to previously re...
Definition: rheoChM.C:796
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.
virtual void giveRealStressVector_2dBeamLayer(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Definition: rheoChM.h:189
virtual void giveRealStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Definition: rheoChM.h:185
virtual bool isActivated(TimeStep *tStep)
Definition: material.h:161
FloatArray discreteTimeScale
Times at which the errors are evaluated if the least-square method is used.
Definition: rheoChM.h:160
double talpha
thermal dilatation coeff.
Definition: rheoChM.h:134
virtual const FloatArray & giveViscoelasticStressVector() const
Definition: rheoChM.h:97
the oofem namespace is to define a context or scope in which all oofem names are defined.
int nUnits
Number of (Maxwell or Kelvin) units in the rheologic chain.
Definition: rheoChM.h:136
void setShrinkageStrainVector(FloatArray src)
Definition: rheoChM.h:105
virtual void giveEigenStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode)
Computes, for the given integration point, the strain vector induced by the stress history (typically...
Definition: rheoChM.h:282
LinearElasticMaterial * linearElasticMaterial
Associated linearElasticMaterial, with E = 1.
Definition: rheoChM.h:153
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
double EparValTime
Time for which the partial moduli of individual units have been evaluated.
Definition: rheoChM.h:146
virtual void giveRealStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_3d.
Definition: rheoChM.h:183
int nUnits
Number of units in the chain.
Definition: rheoChM.h:75
virtual const char * giveClassName() const
Definition: rheoChM.h:121

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