#include <shrinkmat.h>
Public Member Functions | |
double | give_actual_fc (long ipp, long im, long ido) |
double | give_actual_ft (long ipp, long im, long ido) |
void | give_reqnmq (long *anmq) |
void | giveirrstrains (long ipp, long im, long ido, vector &epsirr) |
void | givestressincr (long lcid, long ipp, long im, long ido, long fi, vector &sig) |
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 | print (FILE *out) |
void | read (XFILE *in) |
shrinkmat (void) | |
void | updateval (long ipp, long im, long ido) |
~shrinkmat (void) | |
Public Attributes | |
gfunct | beta |
gfunct | humdef |
gfunct | shmeas |
answertype | thumdef |
tshrlaw | tshr |
This class defines material model which computes irreversible shrinkage strains
21. 11. 2013
Definition at line 17 of file shrinkmat.h.
shrinkmat | ( | void | ) |
This constructor inializes attributes to zero values.
Created 21. 11. 2013 by JK+TKo
Definition at line 23 of file shrinkmat.cpp.
References no, shr_measured, thumdef, and tshr.
~shrinkmat | ( | void | ) |
This destructor is only for the formal purposes.
Created 21. 11. 2013 by JK+TKo
Definition at line 36 of file shrinkmat.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 |
Created 21. 11. 2013 by JK+TKo
Definition at line 448 of file shrinkmat.cpp.
References mechmat::give_actual_fc(), mechmat::ip, Mm, and intpoints::ncompstr.
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 |
Created 21. 11. 2013 by JK+TKo
Definition at line 427 of file shrinkmat.cpp.
References mechmat::give_actual_ft(), mechmat::ip, Mm, and intpoints::ncompstr.
Referenced by mechmat::give_actual_ft().
void give_reqnmq | ( | long * | anmq | ) |
The funtion marks required non-mechanical quantities in the array anmq.
anmq | - array with flags for used material types anmq[i] = 1 => qunatity type nonmechquant(i+1) is required anmq[i] = 0 => qunatity type nonmechquant(i+1) is not required |
Definition at line 470 of file shrinkmat.cpp.
References no, thumdef, and vol_moist_cont.
Referenced by mechmat::give_reqnmq().
void giveirrstrains | ( | long | ipp, | |
long | im, | |||
long | ido, | |||
vector & | epsirr | |||
) |
This function returns the vector of irreversible strains
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 | |
epsirr | - vector containing irreversible strain components (output) |
Created 21. 11. 2013 by JK+TKo
Definition at line 403 of file shrinkmat.cpp.
References intpoints::eqother, mechmat::giveirrstrains(), mechmat::ip, Mm, and intpoints::ncompstr.
Referenced by mechmat::giveirrstrains().
void givestressincr | ( | long | lcid, | |
long | ipp, | |||
long | im, | |||
long | ido, | |||
long | fi, | |||
vector & | sig | |||
) |
This function returns the value of tensile strength
lcid | - load case id | |
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 | |
fi | - first index of the required stress increment component | |
sig | - vector containing stress increment components (output) |
Created 21. 11. 2013 by JK+TKo
Definition at line 378 of file shrinkmat.cpp.
References intpoints::eqother, mechmat::givestressincr(), mechmat::ip, Mm, vector::n, and intpoints::ncompstr.
Referenced by mechmat::givestressincr().
void initvalues | ( | long | lcid, | |
long | ipp, | |||
long | im, | |||
long | ido | |||
) |
Function initializes eqother array with initial values. Actual values of 'vol_moist_cont' from array Mm->nonmechq are taken as 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 |
Created 21. 11. 2013 by JK+TKo
Definition at line 348 of file shrinkmat.cpp.
References intpoints::eqother, gfunct::getval(), mechmat::givenonmechq(), humdef, mechmat::initvalues(), mechmat::ip, Mm, Mp, intpoints::ncompstr, no, thumdef, probdesc::time, and vol_moist_cont.
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 |
Created 21. 11. 2013 by JK+TKo
Definition at line 123 of file shrinkmat.cpp.
References mechmat::ip, mechmat::matstiff(), Mm, and intpoints::ncompstr.
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 |
Created 21. 11. 2013 by JK+TKo
Definition at line 218 of file shrinkmat.cpp.
References mechmat::computenlstresses(), for(), mechmat::give_actual_nu(), mechmat::ip, Mm, intpoints::ncompstr, intpoints::other, planestress, intpoints::ssst, and intpoints::strain.
Referenced by 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.
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 21. 11. 2013 by JK+TKo
Definition at line 141 of file shrinkmat.cpp.
References beta, beta_val, mechmat::computenlstressesincr(), mechmat::elmatstiff(), intpoints::eqother, for(), gfunct::getderiv(), gfunct::getval(), mechmat::givenonmechq(), humdef, mechmat::ip, Mm, Mp, mxv(), intpoints::ncompstr, no, intpoints::other, print_err(), shmeas, shr_measured, intpoints::ssst, strain, tensor_vector(), thumdef, probdesc::time, tshr, and vol_moist_cont.
Referenced by mechmat::computenlstressesincr().
void nonloc_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 |
Created 21. 11. 2013 by JK+TKo
Definition at line 269 of file shrinkmat.cpp.
References mechmat::compnonloc_nlstresses(), for(), mechmat::give_actual_nu(), mechmat::ip, Mm, intpoints::ncompstr, intpoints::other, planestress, intpoints::ssst, and intpoints::strain.
void print | ( | FILE * | out | ) |
Function reads input material parameters
Created 21. 11. 2013 by JK+TKo
Definition at line 83 of file shrinkmat.cpp.
References beta, beta_val, humdef, no, gfunct::print(), print_err(), shmeas, shr_measured, thumdef, tshr, and yes.
Referenced by mechmat::printmatchar().
void read | ( | XFILE * | in | ) |
Function reads input material parameters
Created 21. 11. 2013 by JK+TKo
Definition at line 48 of file shrinkmat.cpp.
References answertype_kwdset(), beta, beta_val, humdef, no, print_err(), gfunct::read(), shmeas, shr_measured, thumdef, tshr, tshrlaw_kwdset(), xfscanf(), and yes.
Referenced by mechmat::readmattype().
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 |
Created 21. 11. 2013 by JK+TKo
Definition at line 320 of file shrinkmat.cpp.
References intpoints::eqother, mechmat::ip, Mm, intpoints::ncompstr, intpoints::other, and mechmat::updateipvalmat().
Referenced by mechmat::updateipvalmat().
Definition at line 40 of file shrinkmat.h.
Referenced by nlstressesincr(), print(), and read().
Definition at line 43 of file shrinkmat.h.
Referenced by initvalues(), nlstressesincr(), print(), and read().
Definition at line 41 of file shrinkmat.h.
Referenced by nlstressesincr(), print(), and read().
Definition at line 42 of file shrinkmat.h.
Referenced by give_reqnmq(), initvalues(), nlstressesincr(), print(), read(), and shrinkmat().
Definition at line 39 of file shrinkmat.h.
Referenced by nlstressesincr(), print(), read(), and shrinkmat().