types
AnyType
module-attribute
AnyType = Union[
ContactLawBase,
DebugVariable,
FabricTensor,
Integration,
Load,
LoadPhase,
Material,
MicromechanicalBase,
Options,
StateVariable,
]
__all__
module-attribute
__all__ = [
"AnyType",
"ContactLawBase",
"DebugVariable",
"ElastoplasticContactLaw",
"FabricTensor",
"Integration",
"LoadPhase",
"Load",
"Material",
"MicromechanicalBase",
"NonlinearElasticContactLaw",
"Options",
"StateVariable",
]
ContactLawBase
Base class for contact laws. Define the relationship between the force and displacement increment
_plasticCorrection
instance-attribute
_plasticCorrection: numba_int = 0
The contact plastic correction step
CSL
CSL(idx: int, sv0: StateVariable, sv: StateVariable) -> float
Calculate critical void ratio. You must define the critical void ratio in the elastoplastic contact law, but it is not necessary for the nonlinear elastic contact law, though you can still define it to use in the nonlinear integration method.
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ec
|
Critical void ratio
TYPE:
|
__init__
Constructor for the base class of all contact laws.
check
check(idx: int, ddisp: ndarray, sv0: StateVariable, sv: StateVariable)
Check the increment is reasonable
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
The displacement increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
elasticStiffness
elasticStiffness(
idx: int, fn: float, sv0: StateVariable, sv: StateVariable
) -> ndarray
Calculate the elastic stiffness matrix, a 3x3 matrix in the following form:
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
fn
|
Normal force
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
Ke
|
Elastic stiffness matrix
TYPE:
|
forceDisplacement
forceDisplacement(
idx: int, ddisp: ndarray, Ke: ndarray, sv0: StateVariable, sv: StateVariable
) -> ndarray
Integrate the force-displacement contact law, update the following variables:
sv.force
: force
The following steps are performed:
- Go to the last step if the displacement increment is zero
preIntegration
: called before the integrationintegrate
: integrate the force-displacement contact law to update the contact forcepostIntegration
: called after the integration
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
The displacement increment
TYPE:
|
Ke
|
Elastic stiffness matrix
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The stiffness matrix |
frictionCoefficient
frictionCoefficient(idx: int, sv0: StateVariable, sv: StateVariable) -> float
Calculate the friction coefficient combined with different criteria.
increment
increment(
idx: int, ddisp: ndarray, Ke: ndarray, sv0: StateVariable, sv: StateVariable
) -> ndarray
Increment the displacement and calculate the force. The following steps are performed:
preIncrement
: called before the increment- Return the stiffness matrix if the displacement increment is zero
disp += ddisp
: update the displacementK = stiffness
: calculate the stiffness matrixforce += K @ ddisp
: update the forcecheck
: check the convergencepostIncrement
: called after the increment
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
The displacement increment
TYPE:
|
Ke
|
Elastic stiffness matrix
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The stiffness matrix |
initialize
initialize(sv: StateVariable)
Initialize the state variables for the contact law
RETURNS | DESCRIPTION |
---|---|
sv
|
Initial state variables
TYPE:
|
integrate
integrate(
idx: int, ddisp: ndarray, Ke: ndarray, sv0: StateVariable, sv: StateVariable
) -> ndarray
Integrate the force-displacement contact law, with a number of increments.
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
The displacement increment calculated from the elastic trial from the force increment
TYPE:
|
Ke
|
Elastic stiffness matrix
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The stiffness matrix at the end of the increment |
postIncrement
postIncrement(
idx: int, ddisp: ndarray, sv0: StateVariable, sv: StateVariable
) -> None
Post-processing after each increment for the contact law
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
The displacement increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
postIntegration
postIntegration(
idx: int, ddisp: ndarray, sv0: StateVariable, sv: StateVariable
) -> None
Post-processing after the contact integration
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
The displacement increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
preIncrement
preIncrement(
idx: int, ddisp: ndarray, sv0: StateVariable, sv: StateVariable
) -> None
Pre-processing before each increment for the contact law
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
The displacement increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
preIntegration
preIntegration(
idx: int, ddisp: ndarray, sv0: StateVariable, sv: StateVariable
) -> None
Pre-processing before the contact integration
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
The displacement increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
stiffness
stiffness(
idx: int, ddisp: ndarray, sv0: StateVariable, sv: StateVariable
) -> ndarray
Stiffness matrix for the contact law. You are required to define the nonlinear relationship between the force and displacement increment:
DebugVariable
enables
instance-attribute
Whether to enable debug information for specific fields
scalars
instance-attribute
scalars: Dict[str, List[Scalar]] = empty(string, ScalarListType)
Scalar debug information
vectors
instance-attribute
vectors: Dict[str, List[Vector]] = empty(string, ScalarListType)
Vector debug information
appendScalar
appendScalar(field: str, index: dict, value: numba_float)
Append scalar debug information.
appendVector
Append vector debug information.
ElastoplasticContactLaw
Bases: ContactLawBase
Base class for elastoplastic contact laws.
CPPM
CPPM(
idx: int, ddisp: ndarray, sv0: StateVariable, sv: StateVariable
) -> ndarray
Solve plastic problem (plastic strain) with the closest point projection method (CPPM)
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
Displacement increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
Current state variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
Kep
|
Elasto-plastic stiffness matrix
TYPE:
|
ExplicitCPA
ExplicitCPA(
idx: int, ddisp: ndarray, sv0: StateVariable, sv: StateVariable
) -> ndarray
Solve plastic problem (plastic strain) with the explicit or cutting plane algorithm (CPA)
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
Displacement increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
Current state variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
Kep
|
Elasto-plastic stiffness matrix
TYPE:
|
dfdforce
dfdforce(idx: int, sv0: StateVariable, sv: StateVariable) -> ndarray
Calculate derivative of yield function f with respect to the force:
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
dfdsig
|
Derivative of yield function f with respect to the force
TYPE:
|
dgdforce
dgdforce(idx: int, sv0: StateVariable, sv: StateVariable) -> ndarray
Calculate derivative of potential function g with respect to the force:
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
dgdsig
|
Derivative of potential function g with respect to the stress
TYPE:
|
hardening
hardening(
idx: int, dgdsig: ndarray, sv0: StateVariable, sv: StateVariable
) -> float
Hardening items in the equation of plastic multiplier, coefficient of plastic multiplier:
make sure to define the derivative of plastic multiplier with respect to plastic strain dkappa_ddispp
(:math:\frac{\partial\kappa}{\partial\boldsymbol{\delta}^p}
) attribute in this class since it will be used in
the :meth:.updateHardeningVariables
method to update hardening variables.
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
dgdsig
|
Derivative of potential function g with respect to the force
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
hardening
|
Hardening items in the equation of plastic multiplier
TYPE:
|
integrate
integrate(
idx: int, ddisp: ndarray, Ke: ndarray, sv0: StateVariable, sv: StateVariable
) -> ndarray
isConverged
isConverged(
idx: int, f: float, trial: bool, sv0: StateVariable, sv: StateVariable
) -> bool
Determine whether the contact force state is converged to the yield surface or the force state is in the elastic region.
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
f
|
Yield function value
TYPE:
|
trial
|
Whether it is for elastic trial
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
Current state variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
bool
|
Whether the contact force state is converged to the yield surface or the force state is in the elastic region |
maintainYieldSurface
maintainYieldSurface(
idx: int,
ddisp: ndarray,
dforce: ndarray,
sv0: StateVariable,
sv: StateVariable,
)
Revert state variables from elastic state to maintain the force state on the yield surface
Notes
This method should be overridden when the force state is already plastic, but the mixed load iteration converges to the elastic state. This is similar to the yield surface correction method, but the force state is in the elastic region initially caused by the mixed load iteration. However, previous mixed load iteration indicates that the force state is already in the plastic region. This method is used to update the state variables to keep the force state on the yield surface.
If you override this method, make sure to initialize the sv.cscalars["plastic"]
state variable to an array
of zeros with the same size as the number of integration points in the initialize
method of the contact law.
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
ddisp
|
Displacement increment
TYPE:
|
dforce
|
Force increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
Current state variables
TYPE:
|
plasticModulus
plasticModulus(
idx: int,
Ke: ndarray,
dfdforce: ndarray,
dgdforce: ndarray,
hardening: float,
sv0: StateVariable,
sv: StateVariable,
) -> float
Calculate the plastic modulus
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
Ke
|
Elastic stiffness matrix
TYPE:
|
dfdforce
|
Derivative of the yield function with respect to the force
TYPE:
|
dgdforce
|
Derivative of the hardening variable with respect to the force
TYPE:
|
hardening
|
Hardening item
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
Current state variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
float
|
Plastic modulus |
plasticMultiplier
plasticMultiplier(
idx: int,
f: float,
Ke: ndarray,
ddisp: ndarray,
dfdforce: ndarray,
Kp: float,
sv0: StateVariable,
sv: StateVariable,
) -> float
Calculate the plastic multiplier
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
f
|
Yield function value
TYPE:
|
Ke
|
Elastic stiffness matrix
TYPE:
|
ddisp
|
Displacement increment
TYPE:
|
dfdforce
|
Derivative of the yield function with respect to the force
TYPE:
|
Kp
|
Plastic modulus
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
Current state variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
float
|
Plastic multiplier |
revertHardeningVariables
revertHardeningVariables(
idx: int,
dlambda: float,
ddispp: ndarray,
sv0: StateVariable,
sv: StateVariable,
)
Revert hardening variables.
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
dlambda
|
Plastic multiplier
TYPE:
|
ddispp
|
Plastic displacement increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
setPlasticMethod
setPlasticMethod(method: str)
Set the plastic method
PARAMETER | DESCRIPTION |
---|---|
method
|
The plastic method
TYPE:
|
updateHardeningVariables
updateHardeningVariables(
idx: int,
dlambda: float,
ddispp: ndarray,
sv0: StateVariable,
sv: StateVariable,
)
Update hardening variables, the derivative of plastic multiplier with respect to plastic strain is available
in the dkappa_ddispp
(:math:\frac{\partial\kappa}{\partial\boldsymbol{\delta}^p}
) attribute:
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
dlambda
|
Plastic multiplier
TYPE:
|
ddispp
|
Plastic displacement increment
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
yieldSurface
yieldSurface(
idx: int, fn: float, fr: float, sv0: StateVariable, sv: StateVariable
) -> float
Calculate yield function. You are required to define the yield function of the contact force and the state variables:
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
fn
|
Normal force
TYPE:
|
fr
|
Shear force
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
float
|
Value of the yield function |
yieldSurfaceCorrection
yieldSurfaceCorrection(idx: int, sv0: StateVariable, sv: StateVariable)
Update state variables after every plastic increment to keep the force state on the yield surface
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
yieldSurfaces
yieldSurfaces(
idx: int, fn: ndarray, fr: ndarray, sv0: StateVariable, sv: StateVariable
) -> ndarray
Calculate yield function for multiple times
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
fn
|
Normal force
TYPE:
|
fr
|
Shear force
TYPE:
|
sv0
|
Initial state variables
TYPE:
|
sv
|
State variables
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
Values of the yield function |
FabricTensor
Base class for fabric tensors.
changed
instance-attribute
changed: FabricTensorChanged = FabricTensorChanged()
fabric tensor changed flags
coefs
instance-attribute
coefs: Dict[str, numba_float] = empty(string, numba_float)
Fabric tensor coefficients
__init__
__init__(
type: AvailableFabricTypes = "chang1990-ext",
evolution: AvailableFabricEvolutionTypes = "zhao2020",
coefs: Dict[str, float] = None,
)
Initialize the fabric tensor.
evolve
evolve(sv0: StateVariable, sv: StateVariable, use_plastic_strain: bool = True)
Update the fabric tensor.
PARAMETER | DESCRIPTION |
---|---|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables at the end of the increment.
TYPE:
|
use_plastic_strain
|
Whether to use the plastic strain for fabric evolution
TYPE:
|
rotate
rotate(rotate: bool = True)
Rotate the fabric tensor, set rotate to False to revert the rotation.
setParameters
setParameters(
type: AvailableFabricTypes = "chang1990-ext",
evolution: AvailableFabricEvolutionTypes = "zhao2020",
coefs: Dict[str, float] = None,
)
Set the fabric tensor parameters.
Integration
_beta
instance-attribute
_beta: VectorXf = zeros(0, dtype=numpy_float)
azimuthal angles of the integration points in the x-y plane
_gamma
instance-attribute
_gamma: VectorXf = zeros(0, dtype=numpy_float)
azimuthal angles of the integration points z-axis
_l
instance-attribute
_l: MatrixXf = zeros((0, 3), dtype=numpy_float)
contact vectors of the integration points
fabric
instance-attribute
fabric: FabricTensor = FabricTensor('chang1990-ext', 'zhao2020', None)
fabric tensor
gauss_weight
instance-attribute
gauss_weight: VectorXf = zeros(0, dtype=numpy_float)
weights of the integration points
n
instance-attribute
n: MatrixXf = zeros((0, 3), dtype=numpy_float)
azimuthal angles of the integration points
Aij
Aij() -> ndarray
The fabric tensor Aij, an array of shape (3, 3).
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The fabric tensor. |
Aij_
Aij_() -> ndarray
Calculate the fabric tensor Aij from the definition, an array of shape (3, 3).
__init__
__init__(npv: float = 1000000000.0, radius=0.00065)
beta
beta() -> ndarray
The azimuthal angles of the integration points in the y-z plane, an array of shape (n,).
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The azimuthal angles of the integration points in the y-z plane. |
gamma
gamma() -> ndarray
The azimuthal angles of the integration points x-axis, an array of shape (n,).
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The azimuthal angles of the integration points x-axis. |
indexOf
Return the indices of the integration points with the given beta and gamma angles and the corresponding
angles. Specify one (and only one) of beta
or gamma
.
PARAMETER | DESCRIPTION |
---|---|
beta
|
The angle between the projection of the azimuth on the YZ plane and the Y axis
TYPE:
|
gamma
|
The angle between the azimuth and the X axis
TYPE:
|
tol
|
The tolerance for the comparison of the angles, by default 1e-6
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The indices of the integration points |
ndarray
|
The corresponding angles, depending on the given angles, gamma when beta is given and vice versa |
l
l() -> ndarray
The contact vectors of the integration points, an array of shape (n, 3).
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The integration points in the x-y plane. |
p
The transformation matrix, an array of shape (3, 3).
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The transformation matrix. |
setAnisotropicParameters
setAnisotropicParameters(
type: AvailableFabricTypes = "chang1990-ext",
evolution: AvailableFabricEvolutionTypes = "zhao2020",
coefs: Dict[str, float] = None,
)
setGauss37
setGauss37(full: bool = False)
Set the Gauss37 integration points.
PARAMETER | DESCRIPTION |
---|---|
full
|
Full 37*2 integration points.
TYPE:
|
setGauss61
setGauss61(full: bool = False)
Set the Gauss61 integration points.
PARAMETER | DESCRIPTION |
---|---|
full
|
Full 61*2 integration points.
TYPE:
|
setNPVFromNormalizedPackingDensity
setNPVFromNormalizedPackingDensity(rho: float)
Set the number of contacts per unit volume from the normalized packing density with the following equation:
PARAMETER | DESCRIPTION |
---|---|
rho
|
normalized packing density
TYPE:
|
setNPVFromVoidRatio
setNPVFromVoidRatio(e: float)
Set the number of contacts per unit volume from the void ratio with the following equation:
PARAMETER | DESCRIPTION |
---|---|
e
|
void ratio
TYPE:
|
setParameters
size
size() -> int
The number of integration points.
RETURNS | DESCRIPTION |
---|---|
int
|
The number of integration points. |
Load
load phases manager.
addConfinementPhase
addDrainedTriaxialPhase
addIsotropicPhase
addLoadPhase
addTrueTriaxialPhase
addTrueTriaxialWithConstantShearPhase
addUndrainedTriaxialPhase
addUndrainedTriaxialWithoutShearPhase
initialize
initialize(sv: StateVariable)
Initialize the load.
PARAMETER | DESCRIPTION |
---|---|
sv
|
The state variables.
TYPE:
|
LoadPhase
A class for a single load phase.
control
instance-attribute
control: VectorXf = control
Control type of the load, 0 for stress control, 1 for strain control, -1 for ignored
total
instance-attribute
total: VectorXf = total
Control whether the load is incremental or total, 0 for incremental, 1 for total
__init__
__init__(
initial: ndarray = zeros(6, dtype=numpy_float),
final: ndarray = zeros(6, dtype=numpy_float),
control: ndarray = zeros(6, dtype=numpy_float),
total: ndarray = zeros(6, dtype=numpy_float),
steps: int = 100,
)
Constructor for the load phase class.
PARAMETER | DESCRIPTION |
---|---|
initial
|
Initial mixed stress/strain.
TYPE:
|
final
|
Final mixed stress/strain.
TYPE:
|
control
|
The control type of the load increment, 0 for stress control, 1 for strain control, -1 for ignored.
TYPE:
|
total
|
Control whether the load is incremental or total, 0 for incremental, 1 for total.
TYPE:
|
steps
|
Number of steps.
TYPE:
|
finalize
Finalize the load phase to ensure the intermediate stress ratio is satisfied.
initialize
initialize(sv: StateVariable)
Initialize the load phase.
PARAMETER | DESCRIPTION |
---|---|
sv
|
The state variables.
TYPE:
|
isCompatibleWith
isCompatibleWith(previous: 'LoadPhase') -> bool
makeCompatibleWith
Make the load phase compatible with the previous load phase.
PARAMETER | DESCRIPTION |
---|---|
previous
|
Previous load phase.
TYPE:
|
Material
Material parameters.
Note
Not using a numba TypedDict because it is slower.
Properties ending with underscores are anisotropic parameters.
Lambda
instance-attribute
Lambda: numba_float = pop('Lambda', 0.06)
compression -> slope of CSL line (lambda)
R
instance-attribute
R: numba_float = pop('R', 0.0)
:math:R
, A material constant, :math:k_n=k_t\times R
b
instance-attribute
b: numba_float = pop('b', 1e-07)
:math:b
, Control the evolution rate of the grain breakage
c
instance-attribute
c: numba_float = pop('c', 0.0)
Usual ratio of extension to compression quantities
custom
instance-attribute
custom: Dict[str, numba_float] = empty(string, numba_float)
Custom parameters
erefu
instance-attribute
erefu: numba_float = pop('erefu', 0.8)
:math:e_{refu}
, Ultimate reference void ratio
fabrics
instance-attribute
fabrics: Dict[str, FabricTensor] = empty(string, FabricTensorType)
Anisotropic fabric tensors
fb0
instance-attribute
fb0: numba_float = pop('fb0', 0.0)
:math:\sigma_{b0}
, Initial adhesive normal stress
fc0
instance-attribute
fc0: numba_float = pop('fc0', 100.0)
:math:\sigma_{c0}
, Initial pre-consolidated stress
fnr
instance-attribute
fnr: numba_float = pop('fnr', 0.0)
Mean effective stress when log(e) = 1 in the LCC line (log(e) - log(p) space)
hp
instance-attribute
hp: numba_float = pop('hp', 5000.0)
:math:h_p
, Related to plastic modulus of bounding surface effect
integration
instance-attribute
integration: Integration = Integration(1000000000.0, 0.00065)
Integration
kp
instance-attribute
kp: numba_float = pop('kp', 0.0)
:math:k_p
, A parameter in the SimSand yield formula
kpr
instance-attribute
kpr: numba_float = pop('kpr', 1.0)
a ratio of normal stiffness (k_n) to peak stiffness (k_p)
krr
instance-attribute
krr: numba_float = pop('krr', 0.4)
a ratio of normal stiffness (k_n) to tangential stiffness (k_r)
m
instance-attribute
m: numba_float = pop('m', 0.0)
:math:m
, Control the distance between NCL and CSL
nd
instance-attribute
nd: numba_float = pop('nd', 0.5)
the exponential in the phase-transformation angle
ng
instance-attribute
ng: numba_float = pop('ng', 0.0)
:math:n_g
, An exponent in the formula of nonlinear shear modulus
np
instance-attribute
np: numba_float = pop('np', 1.0)
the exponential in the dynamic friction angle
nw
instance-attribute
nw: numba_float = pop('nw', 1.0)
:math:n_w
, Control the effect of plastic work to the grain breakage
rho
instance-attribute
rho: numba_float = pop('rho', 0.0)
:math:\rho
, Control the decreasing rate of CSL due to grain breakage
rhoc
instance-attribute
rhoc: numba_float = pop('rhoc', 0.0)
Slope of the LCC line in the log(e) - log(p) space
theta
instance-attribute
theta: numba_float = pop('theta', 0.0)
A constant exponent in the hardening behavior of p0
xi
instance-attribute
xi: numba_float = pop('xi', 0.5)
:math:\xi
, Control the nonlinearity of CSL line
xib
instance-attribute
xib: numba_float = pop('xib', 0.0)
:math:\xi_b
, A soil constant controls hardening rule of :math:\xi_b
xic
instance-attribute
xic: numba_float = pop('xic', 0.0)
:math:\xi_c
, A soil constant controls hardening rule of :math:\xi_c
Mc
Mc(
ba: float = 0.0,
Ra: float = 1.0,
weight_mises: float = 0.0,
weight_tresca: float = 0.0,
weight_mohr_coulomb: float = 1.0,
weight_smp: float = 0.0,
) -> float
Return the critical state slope for compression.
Mc_
Mc_(
idx: int,
ba: float = 0.0,
Ra: float = 1.0,
weight_mises: float = 0.0,
weight_tresca: float = 0.0,
weight_mohr_coulomb: float = 1.0,
weight_smp: float = 0.0,
) -> float
Return the critical state slope for compression.
Me
Me(
ba: float = 0.0,
Ra: float = 1.0,
weight_mises: float = 0.0,
weight_tresca: float = 0.0,
weight_mohr_coulomb: float = 1.0,
weight_smp: float = 0.0,
) -> float
Return the critical state slope for extension.
Me_
Me_(
idx: int,
ba: float = 0.0,
Ra: float = 1.0,
weight_mises: float = 0.0,
weight_tresca: float = 0.0,
weight_mohr_coulomb: float = 1.0,
weight_smp: float = 0.0,
) -> float
Return the critical state slope for extension.
__getitem__
Get the value of a parameter. This method is might be slower than direct getting the attribute.
__setitem__
Set the value of a parameter. This method is might be slower than direct setting the attribute.
anisotropicFactor
get_
mu
mu(
ba: float = 0.0,
Ra: float = 1.0,
weight_mises: float = 0.0,
weight_tresca: float = 0.0,
weight_mohr_coulomb: float = 1.0,
weight_smp: float = 0.0,
) -> float
Return the friction coefficient.
mu_
mu_(
idx: int,
ba: float = 0.0,
Ra: float = 1.0,
weight_mises: float = 0.0,
weight_tresca: float = 0.0,
weight_mohr_coulomb: float = 1.0,
weight_smp: float = 0.0,
) -> float
phi_
phi_(
idx: int,
ba: float = 0.0,
Ra: float = 1.0,
weight_mises: float = 0.0,
weight_tresca: float = 0.0,
weight_mohr_coulomb: float = 1.0,
weight_smp: float = 0.0,
) -> float
phi_combined
phi_combined(
ba: float = 0.0,
Ra: float = 1.0,
weight_mises: float = 0.0,
weight_tresca: float = 0.0,
weight_mohr_coulomb: float = 1.0,
weight_smp: float = 0.0,
) -> float
Friction angle considering combined failure criteria
setAnisotropicParameters
setAnisotropicParameters(
key: str,
type: AvailableFabricTypes = "chang1990-ext",
evolution: AvailableFabricEvolutionTypes = "zhao2020",
coefs: Dict[str, float] = None,
)
Set anisotropic parameters.
PARAMETER | DESCRIPTION |
---|---|
key
|
The key of the anisotropic parameter.
TYPE:
|
type
|
The type of the anisotropic parameter.
TYPE:
|
evolution
|
The evolution type of the anisotropic parameter.
TYPE:
|
coefs
|
The coefficients of the anisotropic parameter. |
sin_phi_mises
Friction angle sine of Von Mises failure criteria
sin_phi_smp
Friction angle sine of SMP failure criteria
sin_phi_tresca
Friction angle sine of Tresca failure criteria
MicromechanicalBase
Base class for all micro-mechanical models.
integration
instance-attribute
integration: Integration = Integration(npv, radius)
The integration points.
stateVars
instance-attribute
stateVars: List[StateVariable] = List(
[StateVariable(0, False) for _ in range(0)]
)
The state variables.
__init__
__init__(
props: Dict[str, float] | None = None,
npv: float = 1000000000.0,
radius: float = 0.00065,
)
averageStrain
averageStrain(disp: ndarray, sv0: StateVariable, sv: StateVariable) -> ndarray
Average the contact displacement into the strain tensor.
PARAMETER | DESCRIPTION |
---|---|
disp
|
The contact displacement.
TYPE:
|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The averaged strain tensor, shape (6,). |
averageStress
averageStress(force: ndarray, sv0: StateVariable, sv: StateVariable) -> ndarray
Average the contact force into the stress tensor.
PARAMETER | DESCRIPTION |
---|---|
force
|
The contact force.
TYPE:
|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The averaged stress tensor, shape (6,). |
contactIntegrate
contactIntegrate(
idx: int,
dforce: ndarray,
ddisp: ndarray,
sv0: StateVariable,
sv: StateVariable,
)
Integrate the contact law.
PARAMETER | DESCRIPTION |
---|---|
dforce
|
The contact force increment.
TYPE:
|
ddisp
|
The contact displacement increment.
TYPE:
|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables.
TYPE:
|
contactIntegrates
contactIntegrates(
dforce: ndarray, ddisp: ndarray, sv0: StateVariable, sv: StateVariable
)
Integrate all the contact laws with cache.
PARAMETER | DESCRIPTION |
---|---|
dforce
|
The contact force increment.
TYPE:
|
ddisp
|
The contact displacement increment.
TYPE:
|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables.
TYPE:
|
createContactLaw
createContactLaw(props: dict[str, float] = None) -> ContactLawBase
The default contact law.
PARAMETER | DESCRIPTION |
---|---|
props
|
The material properties, by default None |
RETURNS | DESCRIPTION |
---|---|
ContactLawBase
|
The default contact law. |
fabricEvolution
fabricEvolution(sv0: StateVariable, sv: StateVariable)
Evolve the fabric tensor.
PARAMETER | DESCRIPTION |
---|---|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables at the end of the increment.
TYPE:
|
increment
increment(
dx: ndarray,
x: ndarray,
S: ndarray,
E: ndarray,
G: ndarray,
sv0: StateVariable,
sv: StateVariable,
) -> ndarray
Run the model for one increment.
PARAMETER | DESCRIPTION |
---|---|
dx
|
The mixed load increment.
TYPE:
|
x
|
The mixed load.
TYPE:
|
S
|
The stress control matrix.
TYPE:
|
E
|
The strain control matrix.
TYPE:
|
G
|
The ignored control matrix.
TYPE:
|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables at the end of the increment.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The updated stress. |
initialize
initialize() -> StateVariable
Initialize the state variables.
RETURNS | DESCRIPTION |
---|---|
StateVariable
|
The state variables. |
localizeStrain
localizeStrain(eps: ndarray, sv0: StateVariable, sv: StateVariable) -> ndarray
Localize the strain to the contact displacement.
PARAMETER | DESCRIPTION |
---|---|
eps
|
The strain or strain increment.
TYPE:
|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The displacement or displacement increment. |
localizeStress
localizeStress(sig: ndarray, sv0: StateVariable, sv: StateVariable) -> ndarray
Localize the stress to the contact force.
PARAMETER | DESCRIPTION |
---|---|
sig
|
The stress or stress increment.
TYPE:
|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The force or force increment. |
macroMicroIntegrate
macroMicroIntegrate(
K: ndarray,
dsig: ndarray,
deps: ndarray,
sv0: StateVariable,
sv: StateVariable,
) -> None
Macro-Micro Integration. The following state variables need to be updated:
sv.force
: the contact force.sv.disp
: the contact displacement.sv.sig
: the global stress.sv.eps
: the global strain.
PARAMETER | DESCRIPTION |
---|---|
K
|
The stiffness matrix.
TYPE:
|
dsig
|
The stress increment.
TYPE:
|
deps
|
The strain increment.
TYPE:
|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variables.
TYPE:
|
postInitialize
postInitialize(sv: StateVariable) -> None
Post-initialize the state variables.
PARAMETER | DESCRIPTION |
---|---|
sv
|
The state variables.
TYPE:
|
postProcess
postProcess(sv: StateVariable) -> None
Post-process the state variables.
PARAMETER | DESCRIPTION |
---|---|
sv
|
The state variables.
TYPE:
|
resetCustomStateVariables
resetCustomStateVariables(sv0: StateVariable, sv: StateVariable)
Reset the custom state variables to zeros existing in the current state variables but not in the initial state variables.
PARAMETER | DESCRIPTION |
---|---|
sv0
|
The initial state variables.
TYPE:
|
sv
|
The current state variables.
TYPE:
|
run
run() -> List[StateVariable]
Run the model.
RETURNS | DESCRIPTION |
---|---|
List[StateVariable]
|
The state variables. |
setContactLaw
Set the contact law.
PARAMETER | DESCRIPTION |
---|---|
contact
|
The contact law.
TYPE:
|
stiffness
stiffness(
sv0: StateVariable,
sv: StateVariable,
kn: ndarray = None,
ks: ndarray = None,
original: bool = False,
) -> ndarray
Calculate the stiffness matrix.
PARAMETER | DESCRIPTION |
---|---|
sv0
|
The state variables at the beginning of the increment.
TYPE:
|
sv
|
The state variable.
TYPE:
|
kn
|
The contact normal stiffness, by default None
TYPE:
|
ks
|
The contact shear stiffness, by default None
TYPE:
|
original
|
Return the original 9x9 stiffness matrix, by default False
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
ndarray
|
The elastic stiffness matrix. |
NonlinearElasticContactLaw
Bases: ContactLawBase
Nonlinear elastic contact law
integrate
integrate(
idx: int, ddisp: ndarray, Ke: ndarray, sv0: StateVariable, sv: StateVariable
) -> ndarray
postIncrement
postIncrement(
idx: int, ddisp: ndarray, sv0: StateVariable, sv: StateVariable
) -> None
Options
Options for micromechanical models
best_fit_strain_averaging
instance-attribute
best_fit_strain_averaging: bool = True
Use the best-fit hypothesis for strain averaging
cache_precision
instance-attribute
cache_precision: numba_int = 15
Cache contact variables precision
contact_check_convergence
instance-attribute
contact_check_convergence: bool = True
Whether to check the convergence
contact_check_plastic_multiplier
instance-attribute
contact_check_plastic_multiplier: bool = True
Whether to check the sign of plastic multiplier
contact_ensure_force_state_on_yield_surface
instance-attribute
contact_ensure_force_state_on_yield_surface: bool = False
Whether to ensure the force state is on the yield surface
contact_friction_angle_combined
instance-attribute
contact_friction_angle_combined: bool = False
Whether to combine the friction angle from different criteria
contact_friction_angle_weight_mises
instance-attribute
contact_friction_angle_weight_mises: numba_float = 0.0
Weight for Von Mises criterion in the friction angle
contact_friction_angle_weight_mohr_coulomb
instance-attribute
contact_friction_angle_weight_mohr_coulomb: numba_float = 1.0
Weight for Mohr-Coulomb criterion in the friction angle
contact_friction_angle_weight_smp
instance-attribute
contact_friction_angle_weight_smp: numba_float = 0.0
Weight for SMP criterion in the friction angle
contact_friction_angle_weight_tresca
instance-attribute
contact_friction_angle_weight_tresca: numba_float = 0.0
Weight for Tresca criterion in the friction angle
contact_linear_critical_state_line
instance-attribute
contact_linear_critical_state_line: bool = False
Whether to use linear critical state line for OSIMSAND contact law
contact_max_increment
instance-attribute
contact_max_increment: numba_float = 0.1
Maximum force ratio increment
contact_max_iterations
instance-attribute
contact_max_iterations: numba_int = 100
Maximum number of iterations
contact_min_normal_force
instance-attribute
contact_min_normal_force: numba_float = 1e-10
Minimum normal force
contact_plastic_method_multiplier
instance-attribute
contact_plastic_method_multiplier: numba_float = 1.0
Multiplier for the contact plastic method
contact_sanisand_constant_hardening_parameter
instance-attribute
contact_sanisand_constant_hardening_parameter: bool = True
Whether to use constant hardening parameter h for the SANISAND contact law
contact_sanisand_maintain_yield_surface
instance-attribute
contact_sanisand_maintain_yield_surface: bool = False
Whether to main yield surface for the SANISAND contact law
disp_steps
instance-attribute
disp_steps: numba_int = -1
Print debug information every disp_steps steps
ignore_distractions
instance-attribute
ignore_distractions: bool = True
Whether to ignore distractions
ignore_fabric_evolution_for_isotropic_loading
instance-attribute
ignore_fabric_evolution_for_isotropic_loading: bool = False
Whether to ignore fabric evolution for isotropic loading
ignore_fabric_rotations_for_isotropic_loading
instance-attribute
ignore_fabric_rotations_for_isotropic_loading: bool = False
Whether to ignore fabric rotations for isotropic loading
increment_plastic_displacement
instance-attribute
increment_plastic_displacement: bool = True
Increment plastic displacement to the total displacement
integration_absolute_tolerance
instance-attribute
integration_absolute_tolerance: numba_float = 1e-10
Absolute tolerance for the unbalanced force
integration_check_absolute_convergence
instance-attribute
integration_check_absolute_convergence: bool = False
Whether to check absolute contact integration error instead of the relative error
integration_check_convergence
instance-attribute
integration_check_convergence: bool = True
Whether to converge the macro-micro integration algorithm when the number of iterations reaches the maximum
integration_loading_ratio
instance-attribute
integration_loading_ratio: numba_float = 1.0
Macro-micro integration loading ratio for the contact force and displacement
integration_max_iterations
instance-attribute
integration_max_iterations: numba_int = 100
Maximum number of macro-micro integration iterations
integration_max_substepping_ratio
instance-attribute
integration_max_substepping_ratio: numba_float = -1.0
Maximum macro-micro integration substepping ratio
integration_min_substepping_ratio
instance-attribute
integration_min_substepping_ratio: numba_float = -1.0
Minimum macro-micro integration substepping ratio
integration_relative_tolerance
instance-attribute
integration_relative_tolerance: numba_float = 0.001
Relative tolerance for the unbalanced force
integration_substepping_relative_to_increment
instance-attribute
integration_substepping_relative_to_increment: bool = False
Whether to use the macro-micro integration substepping relative to the increment of last step
kinematic_hypothesis
instance-attribute
kinematic_hypothesis: bool = False
Whether to use the kinematic hypothesis
mixed_load_absolute_tolerance
instance-attribute
mixed_load_absolute_tolerance: numba_float = 1e-05
Absolute tolerance for the mixed load
mixed_load_check_absolute_convergence
instance-attribute
mixed_load_check_absolute_convergence: bool = False
Whether to check absolute tolerance for the mixed load
mixed_load_check_convergence
instance-attribute
mixed_load_check_convergence: bool = True
Whether to converge the mixed load control algorithm when the number of iterations reaches the maximum
mixed_load_check_strain_convergence
instance-attribute
mixed_load_check_strain_convergence: bool = True
Check strain convergence of the mixed load control algorithm
mixed_load_max_iterations
instance-attribute
mixed_load_max_iterations: numba_int = 100
Maximum number of iterations for the mixed load control algorithm
mixed_load_max_larger_steps
instance-attribute
mixed_load_max_larger_steps: numba_int = -1
Maximum number of larger steps for the mixed load control algorithm
mixed_load_relative_tolerance
instance-attribute
mixed_load_relative_tolerance: numba_float = 0.001
Relative tolerance for the mixed load
precise_incremental_volumetric_strain
instance-attribute
precise_incremental_volumetric_strain: bool = False
Whether to use the precise incremental volumetric strain to update void ratio
run_contact_integration_steps
instance-attribute
run_contact_integration_steps: numba_int = -1
Number of macro-micro integration iterations to run
run_mixed_load_steps
instance-attribute
run_mixed_load_steps: numba_int = -1
Number of mixed load iterations to run
separate_mean_deviatoric_stress_strain
instance-attribute
separate_mean_deviatoric_stress_strain: bool = False
Whether to apply mean stress and deviatoric stress/strain separately in the mixed load control algorithm
update_fabric_coefs
instance-attribute
update_fabric_coefs: bool = False
Whether to update fabric coefficients
update_npv
instance-attribute
update_npv: bool = False
Whether to update npv based on the coordination number and particle radius
update_void_ratio_from_after_isotropic_loading_to_before
instance-attribute
update_void_ratio_from_after_isotropic_loading_to_before: bool = False
Whether to update the void ratio after the isotropic loading to the void ratio before the isotropic loading, useful if the initial void ratio is given after the isotropic loading
use_isotropic_fabric_tensor_in_strain_averaging
instance-attribute
use_isotropic_fabric_tensor_in_strain_averaging: bool = True
Whether to evolve fabric tensor in strain averaging
use_isotropic_fabric_tensor_in_stress_localization
instance-attribute
use_isotropic_fabric_tensor_in_stress_localization: bool = False
Whether to evolve fabric tensor in stress localization
use_plastic_strain_in_fabric_evolution
instance-attribute
use_plastic_strain_in_fabric_evolution: bool = True
Whether to use plastic strain in fabric evolution
StateVariable
contact_min_normal_force
instance-attribute
contact_min_normal_force: numba_float
minimum normal force
cscalars
instance-attribute
custom scalar fields for contacts, a 1D numpy array of size (size,) is expected for each contact
cvectors
instance-attribute
custom vector fields for contacts, a 2D numpy array of size (size, ...) is expected for each contact
__init__
copyContactStateVariables
copyContactStateVariables(
from_idx: int,
to_idx: int,
symmetric: bool = False,
cscalars: List[str] = None,
cvectors: List[str] = None,
)
Copy the contact state variables from one integration point to another
PARAMETER | DESCRIPTION |
---|---|
from_idx
|
Index of the integration point to copy from
TYPE:
|
to_idx
|
Index of the integration point to copy to
TYPE:
|
symmetric
|
Whether to revert the symmetric contact state variables, by default False
TYPE:
|
cscalars
|
Custom scalar fields to be multiplied by -1, by default None |
cvectors
|
Custom vector fields to be multiplied by -1, by default None |
elasticStiffnessMatrix
elasticStiffnessMatrix() -> ndarray
Elastic stiffness matrix, shape (3, 3)
RETURNS | DESCRIPTION |
---|---|
ndarray
|
Elastic stiffness matrix |
fnci
Normal force with minimum cohesion at the integration point with index idx
fsi
Shear force at the first direction at the integration point with index idx
fti
Shear force at the second direction at the integration point with index idx
hash
hash(
dfi: ndarray,
idx: int,
precision: int = 15,
symmetric: bool = False,
cscalars: List[str] = None,
cvectors: List[str] = None,
) -> int
Hash the contact state variables at the integration point with index idx
PARAMETER | DESCRIPTION |
---|---|
dfi
|
Contact force increment
TYPE:
|
idx
|
Index of the integration point
TYPE:
|
precision
|
Precision of the hash value, by default 15
TYPE:
|
symmetric
|
Whether to revert the symmetric contact state variables, by default False
TYPE:
|
cscalars
|
Custom scalar fields to be multiplied by -1, by default None |
cvectors
|
Custom vector fields to be multiplied by -1, by default None |
RETURNS | DESCRIPTION |
---|---|
int
|
Hash value |
nsti
Direction of the shear force in two dimensions with index idx, shape (2,)
plasticDisplacementStrainTensor
plasticDisplacementStrainTensor() -> Tensor2
Plastic strain tensor from displacements
revertSymmetricContactStateVariables
revertSymmetricContactStateVariables(
idx: int, cscalars: List[str] = None, cvectors: List[str] = None
)
Revert contact force, displacement and specified custom variables by multiplying them with -1
PARAMETER | DESCRIPTION |
---|---|
idx
|
Index of the integration point
TYPE:
|
cscalars
|
Custom scalar fields to be multiplied by -1, by default None |
cvectors
|
Custom vector fields to be multiplied by -1, by default None |
ups
ups() -> ndarray
Shear plastic displacement at the first direction for every integration point, shape (size,)
upsi
Shear plastic displacement at the first direction at the integration point with index idx
upsti
Shear plastic displacement at the integration point with index idx, shape (2,)
upt
upt() -> ndarray
Shear plastic displacement at the second direction for every integration point, shape (size,)
upti
Shear plastic displacement at the second direction at the integration point with index idx
us
us() -> ndarray
Shear displacement at the first direction for every integration point, shape (size,)
usi
Shear displacement at the first direction at the integration point with index idx
usti
Shear displacement at the integration point with index idx, shape (2,)
ut
ut() -> ndarray
Shear displacement at the second direction for every integration point, shape (size,)