#include <camclay.h>
Public Member Functions | |
camclay (void) | |
void | changeparam (atsel &atm, vector &val) |
double | cohesion (vector &qtr) |
void | dderyieldfsigma (matrix &ddfds) |
void | der_q_gamma (long ipp, long ido, matrix &sig, vector &qtr, vector &dqdg) |
void | derpotsigma (matrix &sig, vector &q, matrix &dgds) |
void | deryieldfdqdq (matrix &dfq) |
void | deryieldfdsdq (matrix &dfdsdqt) |
void | deryieldfq (matrix &sig, vector &dq) |
void | deryieldfsigma (matrix &sig, vector &q, matrix &dfds) |
void | dhardfdq (long ipp, long ido, double dgamma, vector &qt, matrix &dqdq) |
void | dqdsigma (matrix &sigt, vector &qt, matrix &dqds) |
double | give_actual_nu (long ipp, long im, long ido) |
double | give_actual_ym (long ipp, long im, long ido) |
double | give_consparam (long ipp, long ido) |
double | give_iniporosity (long ipp, long ido) |
double | give_preconspress (long ipp, long ido) |
double | give_virgporosity (long ipp, long ido) |
void | giveirrstrains (long ipp, long ido, vector &epsp) |
void | initval (long ipp, long ido) |
void | matstiff (matrix &d, long im, long ipp, long ido) |
void | nlstresses (long ipp, long im, long ido) |
double | plasmodscalar (long ipp, long ido, matrix &sig, vector &qtr) |
void | read (XFILE *in) |
void | updateq (long ipp, long ido, vector &eps, vector &epsp, matrix &sig, vector &q) |
void | updateval (long ipp, long ido) |
double | yieldfunction (matrix &sig, vector &q) |
~camclay (void) | |
Public Attributes | |
double | kappa |
slope of swelling line | |
double | lambda |
slope of normal consolidation line | |
double | m |
frictional constant | |
strretalg | sra |
stress return algorithm |
This class defines modified cam-clay material model the material is defined by three material constants : m - is frictional constant and it depends on the frictional angle phi lambda - slope of the normal consolidationa line kappa - slope of the swelling line
In addition to these material constants, several initial values have to be specified at each integration point in the following order : v_kappa1 - initial specific volume at the reference pressure p1 after unloading from initial consolidation pressure p_c0 (v_kappa0) p1 - reference pressure (p_1), it must be negative p_c0 - initial consolidation pressure, it must be negative iepsp_1 - . \ . \ components of initial plastic strains . / iepsp_ncomp -
Order of the eqother array components id | description ------------------+------------------------------------------------ <0; ncompstr-1> | plastic strains ncompstr | gamma - consistency parameter ncompstr+1 | hardening parameter p_c - consolidation pressure ncompstr+2 | v_lambda1 - initial value of specific volume on the NCL at reference pressure p1 . | v_ini - initial specific volume . | i1s - mean stress . | gammap - consistency parameter from the previous step //j2s - the second invariant of stress deviator . | epsv - total volume strain ncompstr+7 | epsvp - plastic volume strain
18.1.2005 TKo
Definition at line 48 of file camclay.h.
camclay | ( | void | ) |
This constructor inializes attributes to zero values.
Definition at line 19 of file camclay.cpp.
~camclay | ( | void | ) |
This destructor is only for the formal purposes.
Definition at line 31 of file camclay.cpp.
The function changes material parameters for stochastic analysis.
atm | - selected material parameters (parameters which are changed) | |
val | - array containing new values of parameters |
Created by Tomas Koudelka,
Definition at line 877 of file camclay.cpp.
References atsel::atrib, kappa, lambda, m, atsel::num, and print_err().
double cohesion | ( | vector & | qtr | ) |
void dderyieldfsigma | ( | matrix & | ddfds | ) |
The function computes the second derivatives of yield function with respect of stress tensor sigma.
ddfds | - tensor of the 4-th order where the are derivatives stored |
19.12.2002
Definition at line 175 of file camclay.cpp.
Referenced by mechmat::dfdsigmadsigma(), mechmat::dgdsigmadsigma(), and matstiff().
This function computes derivatives of hardening paramters with respect of consistency parameter gamma.
ipp | - integration point number | |
ido | - index of internal variables for given material in the ipp other array | |
sig | - stress tensor | |
qtr | - vector of hardening variables | |
dqdg | - matrix where the resulting derivatives are stored |
4.1.2002
Definition at line 351 of file camclay.cpp.
References intpoints::eqother, first_invar(), mechmat::ip, kappa, lambda, Mm, and intpoints::ncompstr.
Referenced by plasmodscalar().
The function computes derivatives of plastic potential function with respect of vector sigma.
sig | - stress tensor | |
q | - vector of the hardening parameter | |
dgds | - matrix where the resulting derivatives are stored |
Definition at line 221 of file camclay.cpp.
References deryieldfsigma().
Referenced by matstiff().
void deryieldfdqdq | ( | matrix & | dfq | ) |
This function computes the second derivatives of yield function with respect to hradening parameter.
dfds | - matrix, where the resulting derivatives are stored |
12.4.2005
Definition at line 259 of file camclay.cpp.
void deryieldfdsdq | ( | matrix & | dfdsdqt | ) |
This function computes the second derivatives of yield function with respect to stresses.
dfds | - tensor, where the resulting derivatives are stored. size of dfds = (6,number_of_hardening_param) |
12.4.2005
Definition at line 277 of file camclay.cpp.
Referenced by mechmat::dfdsigmadq(), and mechmat::dgdsigmadq().
This function computes derivatives of as-th yield function with respect of vector of hradening parameters.
sig | - stress tensor | |
dfds | - vector where the resulting derivatives are stored |
4.1.2002
Definition at line 238 of file camclay.cpp.
References first_invar().
Referenced by mechmat::dfdqpar(), mechmat::dgdqpar(), and plasmodscalar().
This function computes derivatives of as-th yield function with respect of vector sigma.
sig | - stress tensor | |
q | - vector of the hardening parameter | |
dfds | - matrix where the resulting derivatives are stored |
4.1.2002
Definition at line 142 of file camclay.cpp.
References cmulm(), deviator(), first_invar(), m, and volume.
Referenced by derpotsigma(), mechmat::dfdsigma(), mechmat::dgdsigma(), and matstiff().
This function computes the derivatives of hardening function with respect to hardening parameters.
ipp | - integration point number | |
ido | - index of internal variables for given material in the ipp other array | |
dgamma | - increment of consistency parameter gamma | |
qt | - vector of hardening variables | |
dqdq | - tensor, where the resulting derivatives are stored. size of dqdq = (number_of_hardening_param, number_of_hardening_param) |
12.4.2005
Definition at line 325 of file camclay.cpp.
References intpoints::eqother, fillm(), mechmat::ip, kappa, lambda, Mm, and intpoints::ncompstr.
Referenced by mechmat::dhdqpar().
This function computes the derivatives of hardening parameters with respect to stresses.
sigt | - stress tensor | |
qt | - vector of hardening variables | |
dqds | - tensor, where the resulting derivatives are stored. size of dqds = (3,3) |
12.4.2005
Definition at line 297 of file camclay.cpp.
References fillm(), first_invar(), and kappa.
double give_actual_nu | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
The function returns actual Poisson's ratio.
ipp | - integration point number | |
im | - index of material type for given ip | |
ido | - index of internal variables for given material in the ipp other array |
Definition at line 519 of file camclay.cpp.
References intpoints::gemid(), mechmat::give_initial_nu(), mechmat::givencompeqother(), mechmat::ip, and Mm.
Referenced by mechmat::give_actual_nu(), and give_actual_ym().
double give_actual_ym | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
The function returns actual Young modulus value.
ipp | - integration point number | |
im | - index of material type for given ip | |
ido | - index of internal variables for given material in the ipp other array |
Definition at line 455 of file camclay.cpp.
References intpoints::eqother, first_invar(), intpoints::gemid(), give_actual_nu(), mechmat::give_initial_ym(), mechmat::givencompeqother(), mechmat::givestrain(), mechmat::givestress(), mechmat::ic, mechmat::ip, kappa, Mm, intpoints::ncompstr, p, intpoints::ssst, strain, stress, and vector_tensor().
Referenced by mechmat::give_actual_ym(), and nlstresses().
double give_consparam | ( | long | ipp, | |
long | ido | |||
) |
This function extracts consistency parametr gamma for the attained equilibrium state from the integration point eqother array.
ipp | - integration point number in the mechmat ip array. | |
ido | - index of internal variables for given material in the ipp other array |
The | function returns value of consistency parameter. |
Definition at line 781 of file camclay.cpp.
References intpoints::eqother, mechmat::ip, Mm, and intpoints::ncompstr.
Referenced by mechmat::give_consparam().
double give_iniporosity | ( | long | ipp, | |
long | ido | |||
) |
The function extracts initial porosity e_ini for the attained equilibrium state from the integration point eqother array.
ipp | - integration point number in the mechmat ip array. | |
ido | - index of internal variables for given material in the ipp other array |
The | function returns value of initial porosity e_ini. |
Created by Tomas Koudelka, 8.10.2013
Definition at line 853 of file camclay.cpp.
References intpoints::eqother, mechmat::ip, Mm, and intpoints::ncompstr.
Referenced by mechmat::give_iniporosity().
double give_preconspress | ( | long | ipp, | |
long | ido | |||
) |
The function extracts preconsolidation pressure p_c for the attained equilibrium state from the integration point eqother array.
ipp | - integration point number in the mechmat ip array. | |
ido | - index of internal variables for given material in the ipp other array |
The | function returns value of preconsolidation pressure. |
Created by Tomas Koudelka, 8.10.2013
Definition at line 805 of file camclay.cpp.
References intpoints::eqother, mechmat::ip, Mm, and intpoints::ncompstr.
Referenced by mechmat::give_preconspress().
double give_virgporosity | ( | long | ipp, | |
long | ido | |||
) |
The function extracts virgin porosity e_lambda1 for the attained equilibrium state from the integration point eqother array.
ipp | - integration point number in the mechmat ip array. | |
ido | - index of internal variables for given material in the ipp other array |
The | function returns value of virgin porosity e_lambda1. |
Created by Tomas Koudelka, 8.10.2013
Definition at line 829 of file camclay.cpp.
References intpoints::eqother, mechmat::ip, Mm, and intpoints::ncompstr.
Referenced by mechmat::give_virgporosity().
void giveirrstrains | ( | long | ipp, | |
long | ido, | |||
vector & | epsp | |||
) |
Function returns irreversible plastic strains.
ipp | - integration point number in the mechmat ip array. | |
ido | - index of the first internal variable for given material in the ipp other array | |
epsp | - vector of irreversible strains |
Returns vector of irreversible strains via parameter epsp
Definition at line 762 of file camclay.cpp.
References intpoints::eqother, mechmat::ip, Mm, and vector::n.
Referenced by mechmat::giveirrstrains().
void initval | ( | long | ipp, | |
long | ido | |||
) |
This function initializes material model data with respect of consistency parameter gamma.
ipp | - integration point number | |
ido | - index of internal variables for given material in the ipp other array |
12/06/2012 TKo
Definition at line 60 of file camclay.cpp.
References mechmat::elmatstiff(), intpoints::eqother, first_invar(), mechmat::ic, mechmat::ip, kappa, lambda, Mm, mxv(), intpoints::ncompstr, intpoints::ssst, stress, and vector_tensor().
Referenced by mechmat::initvalues().
void matstiff | ( | matrix & | d, | |
long | ipp, | |||
long | im, | |||
long | ido | |||
) |
This function computes material stiffness matrix.
d | - allocated matrix structure for material stiffness matrix | |
ipp | - integration point number | |
im | - index of material type for given ip | |
ido | - index of internal variables for given material in the ipp other array |
Definition at line 539 of file camclay.cpp.
References cmulm(), dderyieldfsigma(), derpotsigma(), deryieldfsigma(), mechmat::elmatstiff(), intpoints::eqother, mechmat::givestress(), invm(), mechmat::ip, matrix::m, Mm, Mp, mxm(), mxv(), matrix::n, intpoints::ncompstr, probdesc::nlman, nlstresses(), intpoints::other, plasmodscalar(), scprd(), intpoints::ssst, nonlinman::stmat, stress, probdesc::stressstate, tensor4_matrix(), tensor_vector(), vector_tensor(), vxm(), vxv(), and yieldfunction().
Referenced by mechmat::matstiff().
void nlstresses | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
This function computes stresses at given integration point ipp, depending on the reached strains. The cutting plane algorithm is used. The stress and the other attribute of given integration point is actualized.
ipp | - integration point number in the mechmat ip array. | |
im | - index of material type for given ip | |
ido | - index of internal variables for given material in the ipp other array |
Definition at line 649 of file camclay.cpp.
References cp, mechmat::cutting_plane(), deviator(), intpoints::eqother, first_invar(), give_actual_ym(), strretalg::give_err(), strretalg::give_ni(), strretalg::give_tsra(), mechmat::givestress(), gsra, mechmat::ic, mechmat::ip, Mm, intpoints::ncompstr, mechmat::newton_stress_return(), intpoints::other, print_err(), second_invar(), sra, intpoints::ssst, strain, intpoints::strain, and vector_tensor().
Referenced by mechmat::computenlstresses(), and matstiff().
This function computes plastic modulus.
ipp | - integration point number | |
ido | - index of internal variables for given material in the ipp other array | |
sig | - stress tensor | |
qtr | - vector of hardening parameters |
The | function returns value of the plastic modulus. |
Definition at line 379 of file camclay.cpp.
References der_q_gamma(), deryieldfq(), and scprd().
Referenced by matstiff(), and mechmat::plasmodscalar().
void read | ( | XFILE * | in | ) |
This function reads material parameters from the opened text file given by the parameter in.
in | - pointer to the opned text file |
Definition at line 45 of file camclay.cpp.
References kappa, lambda, m, strretalg::read(), sra, and xfscanf().
Referenced by mechmat::readmattype().
This function computes new value of the hardening parameter q.
ipp | - integration point pointer | |
ido | - index of internal variables for given material in the ipp other array | |
eps | - vector of the attained strains | |
epsp | - vector of the attained plastic strains | |
sig | - attained stress tensor | |
ssst | - stress/strain state parameter (for the used vector_tensor function) | |
q | - vector of the hardening parameters |
Definition at line 406 of file camclay.cpp.
References intpoints::eqother, first_invar(), mechmat::ic, mechmat::ip, kappa, lambda, Mm, vector::n, intpoints::ncompstr, intpoints::ssst, strain, and vector_tensor().
Referenced by mechmat::updateq().
void updateval | ( | long | ipp, | |
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. | |
ido | - index of internal variables for given material in the ipp other array |
Definition at line 732 of file camclay.cpp.
References intpoints::eqother, mechmat::ip, Mm, intpoints::ncompstr, and intpoints::other.
Referenced by mechmat::updateipvalmat().
This function computes the value of yield functions.
sig | - stress tensor | |
q | - vector of hardening parameter |
The | function returns value of yield function for the given stress tensor |
25.3.2002
Definition at line 112 of file camclay.cpp.
References deviator(), f, first_invar(), m, and second_invar().
Referenced by matstiff(), and mechmat::yieldfunction().
double kappa |
slope of swelling line
Definition at line 89 of file camclay.h.
Referenced by camclay(), changeparam(), der_q_gamma(), dhardfdq(), dqdsigma(), give_actual_ym(), initval(), read(), and updateq().
double lambda |
slope of normal consolidation line
Definition at line 86 of file camclay.h.
Referenced by camclay(), changeparam(), der_q_gamma(), dhardfdq(), initval(), read(), and updateq().
double m |
frictional constant
Definition at line 83 of file camclay.h.
Referenced by camclay(), changeparam(), dderyieldfsigma(), deryieldfsigma(), read(), and yieldfunction().
stress return algorithm
Definition at line 92 of file camclay.h.
Referenced by nlstresses(), and read().