OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
This class implements a model for concrete creep and shrinkage according to EuroCode 2 The creep is assumed to be linear (formula from section 3.7 is not considered here) with aging. More...
#include <eurocode2creep.h>
Public Member Functions | |
Eurocode2CreepMaterial (int n, Domain *d) | |
virtual | ~Eurocode2CreepMaterial () |
virtual void | giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) |
Computes the real stress vector for given total strain and integration point. More... | |
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. More... | |
virtual const char * | giveClassName () const |
virtual const char * | giveInputRecordName () const |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description stored in input record. More... | |
virtual MaterialStatus * | CreateStatus (GaussPoint *gp) const |
Creates new copy of associated status and inserts it into given integration point. More... | |
virtual double | computeConcreteStrengthAtAge (double age) |
evaluates concrete strength at given age More... | |
virtual double | computeMeanElasticModulusAtAge (double age) |
evaluates concrete mean elastic modulus at given age More... | |
virtual double | computeCreepFunction (double t, double t_prime, GaussPoint *gp, TimeStep *tStep) |
Evaluation of the compliance function. More... | |
virtual double | computeCreepCoefficient (double t, double t_prime, GaussPoint *gp, TimeStep *tStep) |
Evaluation of the compliance function (according to appendix B from the EC) More... | |
Public Member Functions inherited from oofem::KelvinChainMaterial | |
KelvinChainMaterial (int n, Domain *d) | |
virtual | ~KelvinChainMaterial () |
virtual int | hasNonLinearBehaviour () |
Returns nonzero if receiver is non linear. More... | |
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 creep strain). More... | |
void | computeHiddenVars (GaussPoint *gp, TimeStep *tStep) |
Public Member Functions inherited from oofem::RheoChainMaterial | |
RheoChainMaterial (int n, Domain *d) | |
virtual | ~RheoChainMaterial () |
virtual void | giveRealStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More... | |
virtual void | giveRealStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_3d. More... | |
virtual void | giveRealStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_2dBeamLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_PlateLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveThermalDilatationVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis. More... | |
virtual int | hasMaterialModeCapability (MaterialMode mode) |
Tests if material supports material mode. More... | |
virtual int | hasCastingTimeSupport () |
Tests if material supports casting time. More... | |
virtual int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
Returns the integration point corresponding value in Reduced form. More... | |
virtual contextIOResultType | saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
Stores integration point state to output stream. More... | |
virtual contextIOResultType | restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
Reads integration point state to output stream. More... | |
virtual void | give3dMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point. More... | |
virtual void | givePlaneStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing plane stress stiffness matrix of receiver. More... | |
virtual void | givePlaneStrainStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing plane strain stiffness matrix of receiver. More... | |
virtual void | give1dStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 1d stiffness matrix of receiver. More... | |
virtual void | give2dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 2d lattice stiffness matrix of receiver. More... | |
virtual void | give3dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 3d lattice stiffness matrix of receiver. More... | |
virtual void | computeStressIndependentStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode) |
Computes reduced strain vector in given integration point, generated by internal processes in material, which are independent on loading in particular integration point. More... | |
double | giveAlphaOne () const |
double | giveAlphaTwo () const |
virtual bool | isActivated (TimeStep *tStep) |
Public Member Functions inherited from oofem::StructuralMaterial | |
StructuralMaterial (int n, Domain *d) | |
Constructor. More... | |
virtual | ~StructuralMaterial () |
Destructor. More... | |
virtual void | giveInputRecord (DynamicInputRecord &input) |
Setups the input record string of receiver. More... | |
virtual void | giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
Computes the stiffness matrix for giveRealStressVector of receiver in given integration point, respecting its history. More... | |
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· More... | |
virtual void | giveRealStressVector_ShellStressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep) |
virtual void | giveRealStressVector_Warping (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_Fiber (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
virtual void | giveRealStressVector_Lattice2d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
virtual void | giveRealStressVector_Lattice3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
virtual void | giveRealStressVector_2dPlateSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation is not provided. More... | |
virtual void | giveRealStressVector_3dBeamSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
virtual void | giveEshelbyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Prototype for computation of Eshelby stress. More... | |
void | give_dPdF_from (const FloatMatrix &dSdE, FloatMatrix &answer, GaussPoint *gp) |
void | convert_dSdE_2_dPdF (FloatMatrix &answer, const FloatMatrix &dSdE, const FloatArray &S, const FloatArray &F, MaterialMode matMode) |
double | giveReferenceTemperature () |
Returns the reference temperature of receiver. More... | |
virtual void | computeStressIndependentStrainVector_3d (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode) |
virtual void | give3dMaterialStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
virtual void | give3dMaterialStiffnessMatrix_dCde (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
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 temperature, shrinkage, creep and possibly by other phenomena). More... | |
void | giveStressDependentPartOfStrainVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) |
virtual int | setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type) |
Sets the value of a certain variable at a given integration point to the given value. More... | |
virtual void | give2dBeamLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 2d beam layer stiffness matrix of receiver. More... | |
virtual void | givePlateLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 2d plate layer stiffness matrix of receiver. More... | |
virtual void | giveFiberStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 1d fiber stiffness matrix of receiver. More... | |
virtual void | give2dPlateSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing stiffness matrix of plate subsoil model. More... | |
virtual void | give3dBeamSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing stiffness matrix of beam3d subsoil model. More... | |
virtual void | giveFirstPKStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More... | |
virtual void | giveFirstPKStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveFirstPKStressVector_3d. More... | |
virtual void | giveFirstPKStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveFirstPKStressVector_3d. More... | |
virtual void | giveFirstPKStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveFirstPKStressVector_3d. More... | |
virtual void | giveCauchyStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | giveCauchyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | giveCauchyStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | giveCauchyStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | givePlaneStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | givePlaneStrainStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | give1dStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Public Member Functions inherited from oofem::Material | |
Material (int n, Domain *d) | |
Constructor. More... | |
virtual | ~Material () |
Destructor. More... | |
virtual bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true. More... | |
virtual double | give (int aProperty, GaussPoint *gp) |
Returns the value of material property 'aProperty'. More... | |
virtual bool | hasProperty (int aProperty, GaussPoint *gp) |
Returns true if 'aProperty' exists on material. More... | |
virtual void | modifyProperty (int aProperty, double value, GaussPoint *gp) |
Modify 'aProperty', which already exists on material. More... | |
double | giveCastingTime () |
virtual void | printYourself () |
Prints receiver state on stdout. Useful for debugging. More... | |
virtual int | checkConsistency () |
Allows programmer to test some internal data, before computation begins. More... | |
virtual int | initMaterial (Element *element) |
Optional function to call specific procedures when initializing a material. More... | |
virtual MaterialStatus * | giveStatus (GaussPoint *gp) const |
Returns material status of receiver in given integration point. More... | |
virtual int | packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
Pack all necessary data of integration point (according to element parallel_mode) into given communication buffer. More... | |
virtual int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
Unpack and updates all necessary data of given integration point (according to element parallel_mode) into given communication buffer. More... | |
virtual int | estimatePackSize (DataStream &buff, GaussPoint *ip) |
Estimates the necessary pack size to hold all packed data of receiver. More... | |
virtual double | predictRelativeComputationalCost (GaussPoint *gp) |
Returns the weight representing relative computational cost of receiver The reference material model is linear isotropic material - its weight is set to 1.0 The other material models should compare to this reference model. More... | |
virtual double | predictRelativeRedistributionCost (GaussPoint *gp) |
Returns the relative redistribution cost of the receiver. More... | |
virtual void | initTempStatus (GaussPoint *gp) |
Initializes temporary variables stored in integration point status at the beginning of new time step. More... | |
Public Member Functions inherited from oofem::FEMComponent | |
FEMComponent (int n, Domain *d) | |
Regular constructor, creates component with given number and belonging to given domain. More... | |
virtual | ~FEMComponent () |
Virtual destructor. More... | |
Domain * | giveDomain () const |
virtual void | setDomain (Domain *d) |
Sets associated Domain. More... | |
int | giveNumber () const |
void | setNumber (int num) |
Sets number of receiver. More... | |
virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
Local renumbering support. More... | |
virtual contextIOResultType | saveContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Stores receiver state to output stream. More... | |
virtual contextIOResultType | restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Restores the receiver state previously written in stream. More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
Prints output of receiver to stream, for given time step. More... | |
virtual Interface * | giveInterface (InterfaceType t) |
Interface requesting service. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Protected Types | |
enum | ec2ShrinkageType { EC2_NoShrinkage, EC2_TotalShrinkage, EC2_DryingShrinkage, EC2_AutogenousShrinkage } |
shrinkage option More... | |
Protected Member Functions | |
virtual int | hasIncrementalShrinkageFormulation () |
If only incremental shrinkage strain formulation is provided, then total shrinkage strain must be tracked in status in order to be able to compute total value. More... | |
virtual double | giveEModulus (GaussPoint *gp, TimeStep *tStep) |
Evaluation of the incremental modulus. More... | |
virtual double | computeEquivalentAge (GaussPoint *gp, TimeStep *tStep) |
implements B.9 More... | |
virtual double | computeEquivalentMaturity (GaussPoint *gp, TimeStep *tStep) |
implements B.10 More... | |
void | computeElasticityStrengthParams (int) |
sets parameters for elasticity and strength according to formulas from EC More... | |
void | computeShrinkageParams (int, double) |
sets parameters for shrinkage according to formulas from EC More... | |
void | computeCreepParams (int, double) |
sets parameters for creep according to formulas from EC More... | |
virtual void | computeCharTimes () |
computes retardation times of the aging Kelvin chain More... | |
double | computeRetardationTimeCorrection (int mu) |
computes correction factor which multiplies the retardation times More... | |
double | evaluateSpectrumAt (double tau) |
evaluates retardation spectrum at given time (t-t') More... | |
virtual void | computeCharCoefficients (FloatArray &answer, double tPrime, GaussPoint *gp, TimeStep *tStep) |
Evaluation of characteristic moduli of the Kelvin chain. More... | |
void | computeIncrementOfDryingShrinkageVector (FloatArray &answer, GaussPoint *gp, double tNow, double tThen) |
computes increment of drying shrinkage - the shrinkage strain is isotropic More... | |
void | computeIncrementOfAutogenousShrinkageVector (FloatArray &answer, GaussPoint *gp, double tNow, double tThen) |
computes increment of autogenous shrinkage - the shrinkage strain is isotropic More... | |
Protected Member Functions inherited from oofem::KelvinChainMaterial | |
LinearElasticMaterial * | giveLinearElasticMaterial () |
Protected Member Functions inherited from oofem::RheoChainMaterial | |
const FloatArray & | giveDiscreteTimes () |
void | computeDiscreteRelaxationFunction (FloatArray &answer, const FloatArray &tSteps, double t0, double tr, GaussPoint *gp, TimeStep *tSte) |
Evaluation of the relaxation function at given times. More... | |
void | giveUnitComplianceMatrix (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep) |
Evaluation of elastic compliance matrix for unit Young's modulus. More... | |
void | giveUnitStiffnessMatrix (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep) |
Evaluation of elastic stiffness matrix for unit Young's modulus. More... | |
virtual void | updateEparModuli (double tPrime, GaussPoint *gp, TimeStep *tStep) |
Update of partial moduli of individual chain units. More... | |
double | giveEparModulus (int iChain) |
Access to partial modulus of a given unit. More... | |
double | giveCharTime (int) const |
Access to the characteristic time of a given unit. More... | |
virtual double | giveCharTimeExponent (int i) |
Exponent to be used with the char time of a given unit, usually = 1.0. More... | |
LinearElasticMaterial * | giveLinearElasticMaterial () |
Access to the underlying linear elastic material with unit Young's modulus. More... | |
double | giveEndOfTimeOfInterest () |
Access to the time up to which the response should be accurate. More... | |
void | computeTrueStressIndependentStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode) |
Computes, for the given integration point, the strain vector induced by stress-independent internal processes in the material. More... | |
Protected Attributes | |
double | tau1 |
fixed retardation time of the first unit More... | |
double | EspringVal |
stiffness of the zeroth Kelvin unit More... | |
double | fcm28 |
mean compressive strength at 28 days default - to be specified in units of the analysis (e.g. 30.e6 + stiffnessFacotr 1. or 30. + stiffnessFactor 1.e6) More... | |
double | Ecm28 |
Young's modulus at 28 days default [MPa]. More... | |
double | stiffnessFactor |
factor unifying stiffnesses (Ecm is predicted from fcm...) More... | |
double | s |
parameter determined by cement type More... | |
double | phi_RH |
drying creep coefficient More... | |
double | beta_fcm |
drying creep coefficient More... | |
double | beta_H |
drying creep coefficient More... | |
double | h0 |
effective thickness [mm] More... | |
double | alpha_T_cement |
influence of cement type on concrete equivalent age, B.9 in EC2 More... | |
double | t0 |
duration of curing [day, sec, ...] More... | |
double | begOfDrying |
age (absolute) when concrete started drying, must be in days! More... | |
double | kh |
drying shrinkage coefficient More... | |
double | eps_cd_0 |
asymptotic value of drying shrinkage at zero relative humidity, B.11 in EC2 More... | |
double | eps_ca_infty |
asymptotic value of autogenous shrinakge, 3.12 in EC2 More... | |
enum oofem::Eurocode2CreepMaterial::ec2ShrinkageType | shType |
bool | retardationSpectrumApproximation |
If true, analysis of retardation spectrum is used for evaluation of Kelvin units moduli If false, least-squares method is used for evaluation of Kelvin units moduli (default) More... | |
bool | temperatureDependent |
switch for temperature dependence of concrete maturity (default option is off) More... | |
Protected Attributes inherited from oofem::RheoChainMaterial | |
double | talpha |
thermal dilatation coeff. More... | |
int | nUnits |
Number of (Maxwell or Kelvin) units in the rheologic chain. More... | |
double | relMatAge |
Physical age of the material at simulation time = 0. More... | |
bool | lattice |
double | nu |
Poisson's ratio (assumed to be constant, unaffected by creep). More... | |
double | alphaOne |
Parameters for the lattice model. More... | |
double | alphaTwo |
double | EparValTime |
Time for which the partial moduli of individual units have been evaluated. More... | |
double | begOfTimeOfInterest |
Time from which the model should give a good approximation. Optional field. Default value is 0.1 [day]. More... | |
double | endOfTimeOfInterest |
Time (age???) up to which the model should give a good approximation. More... | |
LinearElasticMaterial * | linearElasticMaterial |
Associated linearElasticMaterial, with E = 1. More... | |
FloatArray | EparVal |
Partial moduli of individual units. More... | |
FloatArray | charTimes |
Characteristic times of individual units (relaxation or retardation times). More... | |
FloatArray | discreteTimeScale |
Times at which the errors are evaluated if the least-square method is used. More... | |
double | timeFactor |
Scaling factor transforming the simulation time units into days (gives the number of simulation time units in one day, e.g. More... | |
int | preCastingTimeMat |
Stiffness at time less than casting time - optional parameter, negative by default. More... | |
Protected Attributes inherited from oofem::StructuralMaterial | |
double | referenceTemperature |
Reference temperature (temperature, when material has been built into structure). More... | |
Protected Attributes inherited from oofem::Material | |
Dictionary | propertyDictionary |
Property dictionary. More... | |
double | castingTime |
Casting time. More... | |
Protected Attributes inherited from oofem::FEMComponent | |
int | number |
Component number. More... | |
Domain * | domain |
Link to domain object, useful for communicating with other FEM components. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from oofem::StructuralMaterial | |
static int | giveSymVI (int ind1, int ind2) |
static int | giveVI (int ind1, int ind2) |
static int | giveVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
Returns a mask of the vector indicies corresponding to components in a general (non-symmetric) second order tensor of some stress/strain/deformation measure that performes work. More... | |
static int | giveVoigtSymVectorMask (IntArray &answer, MaterialMode mmode) |
The same as giveVoigtVectorMask but returns a mask corresponding to a symmetric second order tensor. More... | |
static void | giveInvertedVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
Gives the inverted version of giveVoigtVectorMask. More... | |
static int | giveSizeOfVoigtVector (MaterialMode mmode) |
Returns the size of reduced stress/strain vector according to given mode. More... | |
static int | giveSizeOfVoigtSymVector (MaterialMode mmode) |
Returns the size of symmetric part of a reduced stress/strain vector according to given mode. More... | |
static void | giveFullVectorForm (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
Converts the reduced symmetric Voigt vector (2nd order tensor) to full form. More... | |
static void | giveFullVectorFormF (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
Converts the reduced deformation gradient Voigt vector (2nd order tensor). More... | |
static void | giveFullSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form. More... | |
static void | giveReducedVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
Converts the full symmetric Voigt vector (2nd order tensor) to reduced form. More... | |
static void | giveReducedSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form. More... | |
static void | giveFullSymMatrixForm (FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode) |
Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More... | |
static void | giveReducedMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
Converts the full symmetric Voigt matrix (4th order tensor) to reduced form. More... | |
static void | giveReducedSymMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More... | |
static void | transformStrainVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &strainVector, bool transpose=false) |
Transforms 3d strain vector into another coordinate system. More... | |
static void | transformStressVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &stressVector, bool transpose=false) |
Transforms 3d stress vector into another coordinate system. More... | |
static double | computeVonMisesStress (const FloatArray *currentStress) |
Computes equivalent of von Mises stress. More... | |
static void | giveStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
Computes 3d strain vector transformation matrix from standard vector transformation matrix. More... | |
static void | give2DStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
Computes 2d strain vector transformation matrix from standard vector transformation matrix. More... | |
static void | giveStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
Computes 3d stress vector transformation matrix from standard vector transformation matrix. More... | |
static void | givePlaneStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
Computes 2d stress vector transformation matrix from standard vector transformation matrix. More... | |
static void | sortPrincDirAndValCloseTo (FloatArray *pVal, FloatMatrix *pDir, FloatMatrix *toPDir) |
Method for sorting newly computed principal values (pVal) and corresponding principal directions (pDir) to be closed to some (often previous) principal directions (toPDir). More... | |
static void | computePrincipalValues (FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode) |
Common functions for convenience. More... | |
static void | computePrincipalValDir (FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode) |
Computes principal values and directions of stress or strain vector. More... | |
static double | computeDeviatoricVolumetricSplit (FloatArray &dev, const FloatArray &s) |
Computes split of receiver into deviatoric and volumetric part. More... | |
static void | computeDeviatoricVolumetricSum (FloatArray &s, const FloatArray &dev, double mean) |
static void | applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu) |
static void | applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double GModulus) |
static void | applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu) |
static void | applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double GModulus) |
static void | applyElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu) |
static void | applyElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu) |
static double | computeStressNorm (const FloatArray &stress) |
static double | computeFirstInvariant (const FloatArray &s) |
static double | computeSecondStressInvariant (const FloatArray &s) |
static double | computeThirdStressInvariant (const FloatArray &s) |
static double | computeFirstCoordinate (const FloatArray &s) |
static double | computeSecondCoordinate (const FloatArray &s) |
static double | computeThirdCoordinate (const FloatArray &s) |
Static Public Attributes inherited from oofem::StructuralMaterial | |
static std::vector< std::vector< int > > | vIindex |
Voigt index map. More... | |
static std::vector< std::vector< int > > | svIndex |
Symmetric Voigt index map. More... | |
Static Protected Member Functions inherited from oofem::RheoChainMaterial | |
static void | generateLogTimeScale (FloatArray &answer, double from, double to, int nsteps) |
Generates discrete times starting from time "from" to time "to" uniformly distributed in log time scale. More... | |
This class implements a model for concrete creep and shrinkage according to EuroCode 2 The creep is assumed to be linear (formula from section 3.7 is not considered here) with aging.
The linearity assumption is considered OK for stresses below 0.45 of the compressive strength. The shrinkage deformation is additively split into two parts: drying shrinkage ${sh,d}$ and autogenous shrinkage ${sh,a}$. A detailed description of the material model is given in the material manual. See the tests "/sm/EC2creep.in" and "/sm/EC2shrinkage.in" for examples.
Definition at line 100 of file eurocode2creep.h.
|
protected |
shrinkage option
Enumerator | |
---|---|
EC2_NoShrinkage | |
EC2_TotalShrinkage | |
EC2_DryingShrinkage | |
EC2_AutogenousShrinkage |
Definition at line 161 of file eurocode2creep.h.
|
inline |
Definition at line 174 of file eurocode2creep.h.
|
inlinevirtual |
Definition at line 177 of file eurocode2creep.h.
References oofem::RheoChainMaterialStatus::giveShrinkageStrainVector(), and oofem::IntegrationPointStatus::gp.
|
protectedvirtual |
Evaluation of characteristic moduli of the Kelvin chain.
Reimplemented from oofem::KelvinChainMaterial.
Definition at line 512 of file eurocode2creep.C.
References oofem::FloatArray::at(), beta_fcm, oofem::KelvinChainMaterial::computeCharCoefficients(), computeEquivalentAge(), Ecm28, EspringVal, evaluateSpectrumAt(), oofem::RheoChainMaterial::nUnits, phi_RH, oofem::FloatArray::resize(), retardationSpectrumApproximation, tau1, temperatureDependent, oofem::RheoChainMaterial::timeFactor, oofem::FloatArray::times(), and oofem::FloatArray::zero().
|
protectedvirtual |
computes retardation times of the aging Kelvin chain
Reimplemented from oofem::RheoChainMaterial.
Definition at line 408 of file eurocode2creep.C.
References oofem::FloatArray::at(), oofem::RheoChainMaterial::begOfTimeOfInterest, oofem::RheoChainMaterial::charTimes, computeRetardationTimeCorrection(), oofem::RheoChainMaterial::endOfTimeOfInterest, oofem::RheoChainMaterial::nUnits, OOFEM_WARNING, oofem::FloatArray::resize(), retardationSpectrumApproximation, tau1, oofem::RheoChainMaterial::timeFactor, and oofem::FloatArray::zero().
|
virtual |
evaluates concrete strength at given age
Definition at line 234 of file eurocode2creep.C.
References fcm28, s, and oofem::RheoChainMaterial::timeFactor.
Referenced by computeMeanElasticModulusAtAge().
|
virtual |
Evaluation of the compliance function (according to appendix B from the EC)
Definition at line 273 of file eurocode2creep.C.
References beta_fcm, beta_H, computeEquivalentAge(), phi_RH, temperatureDependent, and oofem::RheoChainMaterial::timeFactor.
Referenced by computeCreepFunction().
|
virtual |
Evaluation of the compliance function.
Implements oofem::RheoChainMaterial.
Definition at line 254 of file eurocode2creep.C.
References computeCreepCoefficient(), computeMeanElasticModulusAtAge(), and Ecm28.
|
protected |
sets parameters for creep according to formulas from EC
Definition at line 176 of file eurocode2creep.C.
References alpha_T_cement, beta_fcm, beta_H, fcm28, h0, oofem::min(), OOFEM_ERROR, phi_RH, and stiffnessFactor.
Referenced by initializeFrom().
|
protected |
sets parameters for elasticity and strength according to formulas from EC
Definition at line 105 of file eurocode2creep.C.
References Ecm28, fcm28, OOFEM_ERROR, s, and stiffnessFactor.
Referenced by initializeFrom().
|
protectedvirtual |
implements B.9
Definition at line 345 of file eurocode2creep.C.
References alpha_T_cement, computeEquivalentMaturity(), and oofem::max().
Referenced by computeCharCoefficients(), and computeCreepCoefficient().
|
protectedvirtual |
implements B.10
Definition at line 309 of file eurocode2creep.C.
References oofem::FloatArray::at(), oofem::Material::castingTime, oofem::StructuralElement::computeResultingIPTemperatureAt(), oofem::Eurocode2CreepMaterialStatus::giveConcreteMaturity(), oofem::GaussPoint::giveElement(), oofem::TimeStep::giveIntrinsicTime(), oofem::Material::giveStatus(), oofem::Eurocode2CreepMaterialStatus::giveTemperature(), oofem::TimeStep::giveTimeIncrement(), oofem::TimeStep::isTheFirstStep(), oofem::max(), oofem::min(), oofem::RheoChainMaterial::relMatAge, oofem::Eurocode2CreepMaterialStatus::setTempConcreteMaturity(), and oofem::Eurocode2CreepMaterialStatus::setTempTemperature().
Referenced by computeEquivalentAge().
|
protected |
computes increment of autogenous shrinkage - the shrinkage strain is isotropic
Definition at line 668 of file eurocode2creep.C.
References oofem::FloatArray::at(), eps_ca_infty, oofem::GaussPoint::giveMaterialMode(), oofem::StructuralMaterial::giveReducedSymVectorForm(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by giveShrinkageStrainVector().
|
protected |
computes increment of drying shrinkage - the shrinkage strain is isotropic
Definition at line 634 of file eurocode2creep.C.
References oofem::FloatArray::at(), eps_cd_0, oofem::GaussPoint::giveMaterialMode(), oofem::StructuralMaterial::giveReducedSymVectorForm(), h0, kh, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by giveShrinkageStrainVector().
|
virtual |
evaluates concrete mean elastic modulus at given age
Definition at line 242 of file eurocode2creep.C.
References computeConcreteStrengthAtAge(), Ecm28, and fcm28.
Referenced by computeCreepFunction(), and giveEModulus().
|
protected |
computes correction factor which multiplies the retardation times
Definition at line 476 of file eurocode2creep.C.
References beta_H, tau1, and oofem::RheoChainMaterial::timeFactor.
Referenced by computeCharTimes().
|
protected |
sets parameters for shrinkage according to formulas from EC
Definition at line 125 of file eurocode2creep.C.
References EC2_AutogenousShrinkage, EC2_DryingShrinkage, EC2_TotalShrinkage, eps_ca_infty, eps_cd_0, fcm28, h0, kh, oofem::max(), OOFEM_ERROR, shType, and stiffnessFactor.
Referenced by initializeFrom().
|
virtual |
Creates new copy of associated status and inserts it into given integration point.
gp | Integration point where newly created status will be stored. |
Reimplemented from oofem::KelvinChainMaterial.
Definition at line 703 of file eurocode2creep.C.
References oofem::FEMComponent::giveDomain(), and oofem::RheoChainMaterial::nUnits.
|
protected |
evaluates retardation spectrum at given time (t-t')
Definition at line 486 of file eurocode2creep.C.
References beta_H, and oofem::RheoChainMaterial::timeFactor.
Referenced by computeCharCoefficients().
|
inlinevirtual |
Reimplemented from oofem::KelvinChainMaterial.
Definition at line 184 of file eurocode2creep.h.
|
protectedvirtual |
Evaluation of the incremental modulus.
Reimplemented from oofem::KelvinChainMaterial.
Definition at line 361 of file eurocode2creep.C.
References computeMeanElasticModulusAtAge(), EspringVal, oofem::KelvinChainMaterial::giveEModulus(), oofem::TimeStep::giveTargetTime(), oofem::TimeStep::giveTimeIncrement(), oofem::RheoChainMaterial::isActivated(), OOFEM_ERROR, oofem::RheoChainMaterial::relMatAge, and retardationSpectrumApproximation.
|
inlinevirtual |
Implements oofem::FEMComponent.
Definition at line 185 of file eurocode2creep.h.
References _IFT_Eurocode2CreepMaterial_Name, and oofem::MaterialStatus::initializeFrom().
|
virtual |
Computes the real stress vector for given total strain and integration point.
The total strain is defined as strain computed directly from displacement field at given time. The stress independent parts (temperature, eigenstrains) are subtracted in constitutive driver. The service should use previously reached equilibrium history variables. Also it should update temporary history variables in status according to newly reached state. The temporary history variables are moved into equilibrium ones after global structure equilibrium has been reached by iteration process.
answer | Stress vector in reduced form. For large deformations it is treated as the second Piola-Kirchoff stress. |
gp | Integration point. |
reducedStrain | Strain vector in reduced form. For large deformations it is treated as the Green-Lagrange strain. |
tStep | Current time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::KelvinChainMaterial.
Definition at line 710 of file eurocode2creep.C.
References oofem::KelvinChainMaterial::giveRealStressVector(), and oofem::RheoChainMaterial::isActivated().
|
virtual |
Computes, for the given integration point, the strain vector induced by stress-independent shrinkage.
answer | Returned strain vector. |
gp | Integration point. |
tStep | Time step (most models are able to respond only when tStep is current time step). |
mode | Determines response mode (Total or incremental). |
Reimplemented from oofem::KelvinChainMaterial.
Definition at line 570 of file eurocode2creep.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), oofem::Material::castingTime, computeIncrementOfAutogenousShrinkageVector(), computeIncrementOfDryingShrinkageVector(), EC2_AutogenousShrinkage, EC2_DryingShrinkage, EC2_NoShrinkage, EC2_TotalShrinkage, oofem::GaussPoint::giveMaterialMode(), oofem::StructuralMaterial::giveSizeOfVoigtSymVector(), oofem::TimeStep::giveTargetTime(), oofem::TimeStep::giveTimeIncrement(), oofem::RheoChainMaterial::isActivated(), oofem::max(), OOFEM_ERROR, oofem::RheoChainMaterial::relMatAge, oofem::FloatArray::resize(), shType, t0, oofem::RheoChainMaterial::timeFactor, and oofem::FloatArray::zero().
|
inlineprotectedvirtual |
If only incremental shrinkage strain formulation is provided, then total shrinkage strain must be tracked in status in order to be able to compute total value.
Reimplemented from oofem::KelvinChainMaterial.
Definition at line 205 of file eurocode2creep.h.
|
virtual |
Initializes receiver according to object description stored in input record.
This function is called immediately after creating object using constructor. Input record can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record.
ir | Input record to initialize from. |
Reimplemented from oofem::KelvinChainMaterial.
Definition at line 48 of file eurocode2creep.C.
References _IFT_Eurocode2CreepMaterial_cemType, _IFT_Eurocode2CreepMaterial_fcm28, _IFT_Eurocode2CreepMaterial_h0, _IFT_Eurocode2CreepMaterial_henv, _IFT_Eurocode2CreepMaterial_shType, _IFT_Eurocode2CreepMaterial_spectrum, _IFT_Eurocode2CreepMaterial_stiffnessFactor, _IFT_Eurocode2CreepMaterial_t0, _IFT_Eurocode2CreepMaterial_temperatureDependent, computeCreepParams(), computeElasticityStrengthParams(), computeShrinkageParams(), fcm28, h0, oofem::InputRecord::hasField(), oofem::KelvinChainMaterial::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_OK, OOFEM_ERROR, retardationSpectrumApproximation, shType, stiffnessFactor, t0, and temperatureDependent.
|
protected |
influence of cement type on concrete equivalent age, B.9 in EC2
Definition at line 139 of file eurocode2creep.h.
Referenced by computeCreepParams(), and computeEquivalentAge().
|
protected |
age (absolute) when concrete started drying, must be in days!
Definition at line 147 of file eurocode2creep.h.
|
protected |
drying creep coefficient
Definition at line 130 of file eurocode2creep.h.
Referenced by computeCharCoefficients(), computeCreepCoefficient(), and computeCreepParams().
|
protected |
drying creep coefficient
Definition at line 133 of file eurocode2creep.h.
Referenced by computeCreepCoefficient(), computeCreepParams(), computeRetardationTimeCorrection(), and evaluateSpectrumAt().
|
protected |
Young's modulus at 28 days default [MPa].
Definition at line 117 of file eurocode2creep.h.
Referenced by computeCharCoefficients(), computeCreepFunction(), computeElasticityStrengthParams(), and computeMeanElasticModulusAtAge().
|
protected |
asymptotic value of autogenous shrinakge, 3.12 in EC2
Definition at line 158 of file eurocode2creep.h.
Referenced by computeIncrementOfAutogenousShrinkageVector(), and computeShrinkageParams().
|
protected |
asymptotic value of drying shrinkage at zero relative humidity, B.11 in EC2
Definition at line 153 of file eurocode2creep.h.
Referenced by computeIncrementOfDryingShrinkageVector(), and computeShrinkageParams().
|
protected |
stiffness of the zeroth Kelvin unit
Definition at line 110 of file eurocode2creep.h.
Referenced by computeCharCoefficients(), and giveEModulus().
|
protected |
mean compressive strength at 28 days default - to be specified in units of the analysis (e.g. 30.e6 + stiffnessFacotr 1. or 30. + stiffnessFactor 1.e6)
Definition at line 114 of file eurocode2creep.h.
Referenced by computeConcreteStrengthAtAge(), computeCreepParams(), computeElasticityStrengthParams(), computeMeanElasticModulusAtAge(), computeShrinkageParams(), and initializeFrom().
|
protected |
effective thickness [mm]
Definition at line 136 of file eurocode2creep.h.
Referenced by computeCreepParams(), computeIncrementOfDryingShrinkageVector(), computeShrinkageParams(), and initializeFrom().
|
protected |
drying shrinkage coefficient
Definition at line 150 of file eurocode2creep.h.
Referenced by computeIncrementOfDryingShrinkageVector(), and computeShrinkageParams().
|
protected |
drying creep coefficient
Definition at line 127 of file eurocode2creep.h.
Referenced by computeCharCoefficients(), computeCreepCoefficient(), and computeCreepParams().
|
protected |
If true, analysis of retardation spectrum is used for evaluation of Kelvin units moduli If false, least-squares method is used for evaluation of Kelvin units moduli (default)
Definition at line 167 of file eurocode2creep.h.
Referenced by computeCharCoefficients(), computeCharTimes(), giveEModulus(), and initializeFrom().
|
protected |
parameter determined by cement type
Definition at line 123 of file eurocode2creep.h.
Referenced by computeConcreteStrengthAtAge(), and computeElasticityStrengthParams().
|
protected |
Referenced by computeShrinkageParams(), giveShrinkageStrainVector(), and initializeFrom().
|
protected |
factor unifying stiffnesses (Ecm is predicted from fcm...)
Definition at line 120 of file eurocode2creep.h.
Referenced by computeCreepParams(), computeElasticityStrengthParams(), computeShrinkageParams(), and initializeFrom().
|
protected |
duration of curing [day, sec, ...]
Definition at line 144 of file eurocode2creep.h.
Referenced by giveShrinkageStrainVector(), and initializeFrom().
|
protected |
fixed retardation time of the first unit
Definition at line 107 of file eurocode2creep.h.
Referenced by computeCharCoefficients(), computeCharTimes(), and computeRetardationTimeCorrection().
|
protected |
switch for temperature dependence of concrete maturity (default option is off)
Definition at line 170 of file eurocode2creep.h.
Referenced by computeCharCoefficients(), computeCreepCoefficient(), and initializeFrom().