#include <boermat.h>
Public Member Functions | |
boermat (void) | |
void | changeparam (atsel &atm, vector &val) |
void | derpotsigma (matrix &sig, matrix &dgds) |
void | deryieldfsigma (matrix &sig, matrix &dfds) |
double | give_consparam (long ipp, long ido) |
void | giveirrstrains (long ipp, long ido, vector &epsp) |
void | matstiff (matrix &d, long ipp, long ido) |
void | nlstresses (long ipp, long im, long ido) |
void | nonloc_nlstresses (long ipp, long im, long ido) |
void | read (XFILE *in) |
void | updateval (long ipp, long im, long ido) |
double | yieldfunction (matrix &sig) |
~boermat (void) | |
Public Attributes | |
double | a |
double | alpha |
double | alpha1 |
double | alpha2 |
double | beta |
double | beta1 |
double | c |
cohesion | |
double | delta |
double | n |
exponent for smoothing | |
double | phi |
friction angle | |
double | psi |
dilation | |
strretalg | sra |
stress return algorithm |
This class implements Boer plasticity model, which is similar to Mohr-Coulomb but the singularities on the cone edges are removed by the smoothing. Smoothing is controled by the attribute n which has default valu 0.229.
Definition at line 15 of file boermat.h.
boermat | ( | void | ) |
~boermat | ( | void | ) |
This destructor is only for the formal purposes.
Definition at line 28 of file boermat.cpp.
Definition at line 357 of file boermat.cpp.
References a, alpha, alpha1, alpha2, atsel::atrib, beta, beta1, c, delta, n, atsel::num, phi, and psi.
This function computes derivatives of plastic potential function with respect of vector sigma.
sig | - stress tensor | |
dgds | - matrix where the resulting derivatives are stored |
Definition at line 170 of file boermat.cpp.
References alpha1, c, cmulm(), deviator(), and normedtensor().
This function computes derivatives of as-th yield function with respect of vector sigma.
sig | - stress tensor | |
dfds | - matrix where the resulting derivatives are stored |
4.1.2002
Definition at line 104 of file boermat.cpp.
References addm(), alpha1, cmulm(), copym(), delta, deviator(), fillm(), Mp, n, second_invar(), subm(), third_invar(), and probdesc::zero.
Referenced by mechmat::dfdsigma(), and mechmat::dgdsigma().
double give_consparam | ( | long | ipp, | |
long | ido | |||
) |
This function extracts consistency parametr gamma for the reached equilibrium state from the integration point other array.
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 346 of file boermat.cpp.
References mechmat::ip, Mm, intpoints::ncompstr, and intpoints::other.
Referenced by mechmat::give_consparam().
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 329 of file boermat.cpp.
References intpoints::eqother, mechmat::ip, Mm, and vector::n.
Referenced by mechmat::giveirrstrains().
void matstiff | ( | matrix & | d, | |
long | ipp, | |||
long | ido | |||
) |
This function computes material stiffnes matrix.
d | - allocated matrix structure for material stiffness matrix | |
ipp | - integration point number | |
ido | - index of internal variables for given material in the ipp other array |
Definition at line 193 of file boermat.cpp.
References mechmat::elmatstiff(), Mm, Mp, probdesc::nlman, and nonlinman::stmat.
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 218 of file boermat.cpp.
References cp, mechmat::cutting_plane(), intpoints::eqother, strretalg::give_err(), strretalg::give_ni(), strretalg::give_tsra(), mechmat::ip, Mm, intpoints::ncompstr, intpoints::other, sra, and intpoints::strain.
Referenced by mechmat::compnonloc_nlstresses(), and mechmat::computenlstresses().
void nonloc_nlstresses | ( | long | ipp, | |
long | im, | |||
long | ido | |||
) |
This function computes stresses at given integration point ipp, depending on the reached averaged nonlocal 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 266 of file boermat.cpp.
References cp, mechmat::cutting_plane(), intpoints::eqother, strretalg::give_err(), strretalg::give_ni(), strretalg::give_tsra(), mechmat::ip, Mm, intpoints::ncompstr, intpoints::nonloc, intpoints::other, sra, and intpoints::strain.
Referenced by mechmat::compnonloc_nlstresses().
void read | ( | XFILE * | in | ) |
This function reads material parameters from the opened text file given by the parameter in. Then it computes material constants alpha, alpha1, alpha2, beta and beta1
in | - pointer to the opned text file |
Definition at line 42 of file boermat.cpp.
References a, alpha, alpha1, alpha2, beta, beta1, c, delta, n, phi, psi, strretalg::read(), sra, and xfscanf().
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 | - index of material type for given ip | |
ido | - index of internal variables for given material in the ipp other array |
Definition at line 309 of file boermat.cpp.
References intpoints::eqother, mechmat::givencompeqother(), mechmat::ip, Mm, n, and intpoints::other.
Referenced by mechmat::updateipvalmat().
double yieldfunction | ( | matrix & | sig | ) |
This function computes the value of yield functions.
sig | - stress tensor |
The | function returns value of yield function for the given stress tensor |
25.3.2002
Definition at line 68 of file boermat.cpp.
References alpha1, beta1, c, delta, deviator(), f, first_invar(), Mp, n, second_invar(), third_invar(), and probdesc::zero.
Referenced by mechmat::yieldfunction().
double a |
Definition at line 44 of file boermat.h.
Referenced by changeparam(), and read().
double alpha |
Definition at line 44 of file boermat.h.
Referenced by changeparam(), and read().
double alpha1 |
Definition at line 44 of file boermat.h.
Referenced by changeparam(), derpotsigma(), deryieldfsigma(), read(), and yieldfunction().
double alpha2 |
Definition at line 44 of file boermat.h.
Referenced by changeparam(), and read().
double beta |
Definition at line 44 of file boermat.h.
Referenced by changeparam(), and read().
double beta1 |
Definition at line 44 of file boermat.h.
Referenced by changeparam(), read(), and yieldfunction().
double c |
cohesion
Definition at line 35 of file boermat.h.
Referenced by boermat(), changeparam(), derpotsigma(), read(), and yieldfunction().
double delta |
Definition at line 44 of file boermat.h.
Referenced by changeparam(), deryieldfsigma(), read(), and yieldfunction().
double n |
exponent for smoothing
Definition at line 39 of file boermat.h.
Referenced by boermat(), changeparam(), deryieldfsigma(), read(), updateval(), and yieldfunction().
double phi |
friction angle
Definition at line 33 of file boermat.h.
Referenced by boermat(), changeparam(), and read().
double psi |
dilation
Definition at line 37 of file boermat.h.
Referenced by boermat(), changeparam(), and read().
stress return algorithm
Definition at line 47 of file boermat.h.
Referenced by nlstresses(), nonloc_nlstresses(), and read().