#include <creepdam.h>
Public Member Functions | |
creepdam (void) | |
double | give_actual_fc (long ipp, long im, long ido) |
double | give_actual_ft (long ipp, long im, long ido) |
void | initvalues (long lcid, long ipp, long im, long ido) |
void | matstiff (matrix &d, long ipp, long im, long ido) |
void | nlstresses (long ipp, long im, long ido) |
void | nlstressesincr (long ipp, long im, long ido) |
void | nonloc_nlstresses (long ipp, long im, long ido) |
void | updateval (long ipp, long im, long ido) |
~creepdam (void) |
This class defines artificial material model which combines a creep model with damage model
Order of internal variables in the other array : ----------------------------------------------- The first components of eqother array belong to creep model, followed by the damage eqother components and optionally followed by the thermisomatttime components
Exact order of the components depends on the used material models.
15.1.2004
Definition at line 24 of file creepdam.h.
creepdam | ( | void | ) |
This constructor inializes attributes to zero values.
Definition at line 18 of file creepdam.cpp.
~creepdam | ( | void | ) |
This destructor is only for the formal purposes.
Definition at line 26 of file creepdam.cpp.
double give_actual_fc | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
This function returns the value of compressive strength
ipp | - integration point number in the mechmat ip array. | |
im | - material index | |
ido | - index of internal variables for given material in the ipp other array |
10.2008 TKo
Definition at line 286 of file creepdam.cpp.
References mechmat::give_actual_fc(), mechmat::givencompeqother(), and Mm.
Referenced by mechmat::give_actual_fc().
double give_actual_ft | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
This function returns the value of tensile strength
ipp | - integration point number in the mechmat ip array. | |
im | - material index | |
ido | - index of internal variables for given material in the ipp other array |
TKo
Definition at line 268 of file creepdam.cpp.
References mechmat::give_actual_ft(), and Mm.
Referenced by mechmat::give_actual_ft().
void initvalues | ( | long | lcid, | |
long | ipp, | |||
long | im, | |||
long | ido | |||
) |
Function initializes eqother array with initial values. Actual values of quantities 'temperature' and 'rel_hum' from array Mm->nonmechq are taken as initial temperature and initial moisture.
lcid | - load case id | |
ipp | - integration point pointer | |
im | - material index | |
ido | - index of internal variables for given material in the ipp other array |
7.6.2005, TKo
Definition at line 250 of file creepdam.cpp.
References mechmat::givencompeqother(), mechmat::initvalues(), and Mm.
Referenced by mechmat::initvalues().
void matstiff | ( | matrix & | d, | |
long | ipp, | |||
long | im, | |||
long | ido | |||
) |
This function computes material stiffnes matrix.
d | - allocated matrix structure for material stiffness matrix | |
ipp | - integration point number | |
im | - material index | |
ido | - index of internal variables for given material in the ipp eqother array |
Definition at line 41 of file creepdam.cpp.
References mechmat::givencompeqother(), mechmat::matstiff(), and Mm.
Referenced by mechmat::matstiff().
void nlstresses | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
This function computes correct stresses in the integration point and stores them into ip stress array.
ipp | - integration point pointer | |
im | - index of the material in the ipp tm and idm arrays. The standard value is zero. | |
ido | - index of internal variables for given material in the ipp other array |
7.10.2001
Definition at line 127 of file creepdam.cpp.
References mechmat::compnonloc_nlstresses(), mechmat::computenlstresses(), mechmat::eliso, intpoints::gemid(), mechmat::giveirrstrains(), mechmat::givencompeqother(), intpoints::hmt, intpoints::idm, mechmat::ip, Mm, Mp, intpoints::ncompstr, probdesc::nonlocphase, elastisomat::nu, planestress, intpoints::ssst, and intpoints::strain.
Referenced by mechmat::compnonloc_nlstresses(), and mechmat::computenlstresses().
void nlstressesincr | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
This function computes correct stresses in the integration point and stores them into ip stress array. It was used for Hinkley computations and supposes creep strains to be read from backup file.
ipp | - integration point pointer 7.10.2001 |
void creepdam::nlstresses (long ipp) { long i; long ncomp = Mm->ip[ipp].ncompstr; vector epsback(ncomp);
backup of the total strains for (i=0;i<ncomp;i++){ epsback[i] = Mm->ip[ipp].strain[i]; }
compute creep - actually it is supposed the creep strains are read from backup file to the eqother array. initial values of elastic strain for damage = total strain - creep strain for (i=0;i<ncomp;i++){ Mm->ip[ipp].strain[i] -= Mm->ip[ipp].eqother[i]; }
compute damage, it has index 1, eqother values for damage starts at ncomp index Mm->computenlstresses(ipp, 1, 2*ncomp);
recovery of the total strains for (i=0;i<ncomp;i++){ Mm->ip[ipp].strain[i] = epsback[i]; }
} This function computes stresses increments fro creep model in the integration point and stores them into ip stress array.
ipp | - integration point pointer | |
im | - index of the material in the ipp tm and idm arrays. The standard value is zero. | |
ido | - index of internal variables for given material in the ipp other array |
Created by Tomas Koudelka 7.2008
Definition at line 99 of file creepdam.cpp.
References mechmat::computenlstressesincr(), mechmat::give_actual_nu(), mechmat::ip, Mm, planestress, intpoints::ssst, and intpoints::strain.
Referenced by mechmat::computenlstressesincr().
void nonloc_nlstresses | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
This function computes correct stresses for the nonlocal damage models in the integration point and stores them into ip stress array.
ipp | - integration point pointer | |
im | - index of the material in the ipp tm and idm arrays. The standard value is zero. | |
ido | - index of internal variables for given material in the ipp other array |
7.10.2001
Definition at line 177 of file creepdam.cpp.
References mechmat::compnonloc_nlstresses(), mechmat::eliso, intpoints::gemid(), mechmat::giveirrstrains(), mechmat::givencompeqother(), intpoints::hmt, intpoints::idm, mechmat::ip, Mm, Mp, intpoints::ncompstr, probdesc::nonlocphase, elastisomat::nu, planestress, intpoints::ssst, and intpoints::strain.
void updateval | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
This function updates values in the other array reached in the previous equlibrium state to values reached in the new actual equilibrium state.
ipp | - integration point number in the mechmat ip array. | |
im | - material index | |
ido | - index of internal variables for given material in the ipp other array |
Definition at line 226 of file creepdam.cpp.
References mechmat::givencompeqother(), Mm, and mechmat::updateipvalmat().
Referenced by mechmat::updateipvalmat().