weac.mixins module¶
Mixins for the elastic analysis of layered snow slabs.
- class weac.mixins.AnalysisMixin[source]¶
Bases:
object
Mixin for the analysis of model outputs.
Provides methods for the analysis of layered slabs on compliant elastic foundations.
- Sxx(Z, phi, dz=2, unit='kPa')[source]¶
Compute axial normal stress in slab layers.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x), psi’(x)]^T
phi (float) – Inclination (degrees). Counterclockwise positive.
dz (float, optional) – Element size along z-axis (mm). Default is 2 mm.
unit ({'kPa', 'MPa'}, optional) – Desired output unit. Default is ‘kPa’.
- Returns:
Axial slab normal stress in specified unit.
- Return type:
ndarray, float
- Szz(Z, phi, dz=2, unit='kPa')[source]¶
Compute transverse normal stress in slab layers.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x), psi’(x)]^T
phi (float) – Inclination (degrees). Counterclockwise positive.
dz (float, optional) – Element size along z-axis (mm). Default is 2 mm.
unit ({'kPa', 'MPa'}, optional) – Desired output unit. Default is ‘kPa’.
- Returns:
Transverse normal stress at grid points in the slab in specified unit.
- Return type:
ndarray, float
- Txz(Z, phi, dz=2, unit='kPa')[source]¶
Compute shear stress in slab layers.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x), psi’(x)]^T
phi (float) – Inclination (degrees). Counterclockwise positive.
dz (float, optional) – Element size along z-axis (mm). Default is 2 mm.
unit ({'kPa', 'MPa'}, optional) – Desired output unit. Default is ‘kPa’.
- Returns:
Shear stress at grid points in the slab in specified unit.
- Return type:
ndarray
- gdif(C, phi, li, ki, unit='kJ/m^2', **kwargs)[source]¶
Compute differential energy release rate of all crack tips.
- Parameters:
C (ndarray) – Free constants of the solution.
phi (float) – Inclination (degress).
li (ndarray) – List of segment lengths.
ki (ndarray) – List of booleans indicating whether segment lies on a foundation or not in the cracked configuration.
- Returns:
List of total, mode I, and mode II energy release rates.
- Return type:
ndarray
- get_zmesh(dz=2)[source]¶
Get z-coordinates of grid points and corresponding elastic properties.
- Parameters:
dz (float, optional) – Element size along z-axis (mm). Default is 2 mm.
- Returns:
mesh – Mesh along z-axis. Columns are a list of z-coordinates (mm) of grid points along z-axis with at least two grid points (top, bottom) per layer, Young’s modulus of each grid point, shear modulus of each grid point, and Poisson’s ratio of each grid point.
- Return type:
ndarray
- ginc(C0, C1, phi, li, ki, k0, **kwargs)[source]¶
Compute incremental energy relase rate of of all cracks.
- Parameters:
C0 (ndarray) – Free constants of uncracked solution.
C1 (ndarray) – Free constants of cracked solution.
phi (float) – Inclination (degress).
li (ndarray) – List of segment lengths.
ki (ndarray) – List of booleans indicating whether segment lies on a foundation or not in the cracked configuration.
k0 (ndarray) – List of booleans indicating whether segment lies on a foundation or not in the uncracked configuration.
- Returns:
List of total, mode I, and mode II energy release rates.
- Return type:
ndarray
- principal_stress_slab(Z, phi, dz=2, unit='kPa', val='max', normalize=False)[source]¶
Compute maxium or minimum principal stress in slab layers.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x), psi’(x)]^T
phi (float) – Inclination (degrees). Counterclockwise positive.
dz (float, optional) – Element size along z-axis (mm). Default is 2 mm.
unit ({'kPa', 'MPa'}, optional) – Desired output unit. Default is ‘kPa’.
val (str, optional) – Maximum ‘max’ or minimum ‘min’ principal stress. Default is ‘max’.
normalize (bool) – Toggle layerwise normalization to strength.
- Returns:
Maximum or minimum principal stress in specified unit.
- Return type:
ndarray
- Raises:
ValueError – If specified principal stress component is neither ‘max’ nor ‘min’, or if normalization of compressive principal stress is requested.
- principal_stress_weaklayer(Z, sc=2.6, unit='kPa', val='min', normalize=False)[source]¶
Compute maxium or minimum principal stress in the weak layer.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x), psi’(x)]^T
sc (float) – Weak-layer compressive strength. Default is 2.6 kPa.
unit ({'kPa', 'MPa'}, optional) – Desired output unit. Default is ‘kPa’.
val (str, optional) – Maximum ‘max’ or minimum ‘min’ principal stress. Default is ‘min’.
normalize (bool) – Toggle layerwise normalization to strength.
- Returns:
Maximum or minimum principal stress in specified unit.
- Return type:
ndarray
- Raises:
ValueError – If specified principal stress component is neither ‘max’ nor ‘min’, or if normalization of tensile principal stress is requested.
- rasterize_solution(C: ndarray, phi: float, li: list[float] | bool, ki: list[bool] | bool, num: int = 250, **kwargs)[source]¶
Compute rasterized solution vector.
- Parameters:
C (ndarray) – Vector of free constants.
phi (float) – Inclination (degrees).
li (ndarray) – List of segment lengths (mm).
ki (ndarray) – List of booleans indicating whether segment lies on a foundation or not.
num (int) – Number of grid points.
- Returns:
xq (ndarray) – Grid point x-coordinates at which solution vector is discretized.
zq (ndarray) – Matrix with solution vectors as colums at grid points xq.
xb (ndarray) – Grid point x-coordinates that lie on a foundation.
- class weac.mixins.FieldQuantitiesMixin[source]¶
Bases:
object
Mixin for field quantities.
Provides methods for the computation of displacements, stresses, strains, and energy release rates from the solution vector.
- Gi(Ztip, unit='kJ/m^2')[source]¶
Get mode I differential energy release rate at crack tip.
- Parameters:
Ztip (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T at the crack tip.
unit ({'N/mm', 'kJ/m^2', 'J/m^2'}, optional) – Desired output unit. Default is kJ/m^2.
- Returns:
Mode I differential energy release rate (N/mm = kJ/m^2 or J/m^2) at the crack tip.
- Return type:
float
- Gii(Ztip, unit='kJ/m^2')[source]¶
Get mode II differential energy release rate at crack tip.
- Parameters:
Ztip (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T at the crack tip.
unit ({'N/mm', 'kJ/m^2', 'J/m^2'}, optional) – Desired output unit. Default is kJ/m^2 = N/mm.
- Returns:
Mode II differential energy release rate (N/mm = kJ/m^2 or J/m^2) at the crack tip.
- Return type:
float
- M(Z)[source]¶
Get bending moment M = B11 u’ + D11 psi’ in the slab.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
- Returns:
Bending moment M (Nmm) in the slab.
- Return type:
float
- N(Z)[source]¶
Get the axial normal force N = A11 u’ + B11 psi’ in the slab.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
- Returns:
Axial normal force N (N) in the slab.
- Return type:
float
- V(Z)[source]¶
Get vertical shear force V = kA55(w’ + psi).
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
- Returns:
Vertical shear force V (N) in the slab.
- Return type:
float
- dpsi_dx(Z)[source]¶
Get first derivative psi’ of the midplane rotation.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
- Returns:
First derivative psi’ of the midplane rotation (radians/mm) of the slab.
- Return type:
float
- dpsi_dxdx(z, phi)[source]¶
Get second derivative of the cross-section rotation psi’’(x).
- Parameters:
z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
phi (float) – Inclination (degrees). Counterclockwise positive.
- Returns:
Second derivative of the cross-section rotation psi’’(x) (1/mm^2).
- Return type:
ndarray, float
- dpsi_dxdxdx(z, phi)[source]¶
Get third derivative of the cross-section rotation psi’’’(x).
- Parameters:
z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
phi (float) – Inclination (degrees). Counterclockwise positive.
- Returns:
Third derivative of the cross-section rotation psi’’’(x) (1/mm^3).
- Return type:
ndarray, float
- du0_dxdx(z, phi)[source]¶
Get second derivative of the horiz. centerline displacement u0’’(x).
- Parameters:
z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
phi (float) – Inclination (degrees). Counterclockwise positive.
- Returns:
Second derivative of the horizontal centerline displacement u0’’(x) (1/mm).
- Return type:
ndarray, float
- du0_dxdxdx(z, phi)[source]¶
Get third derivative of the horiz. centerline displacement u0’’’(x).
- Parameters:
z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
phi (float) – Inclination (degrees). Counterclockwise positive.
- Returns:
Third derivative of the horizontal centerline displacement u0’’’(x) (1/mm^2).
- Return type:
ndarray, float
- du_dx(Z, z0)[source]¶
Get first derivative of the horizontal displacement.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
z0 (float) – Z-coordinate (mm) where u is to be evaluated.
- Returns:
First derivative u’ = u0’ + z0 psi’ of the horizontal displacement of the slab.
- Return type:
float
- dw_dx(Z)[source]¶
Get first derivative w’ of the centerline deflection.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
- Returns:
First derivative w’ of the deflection of the slab.
- Return type:
float
- dz_dx(z, phi)[source]¶
Get first derivative z’(x) = K*z(x) + q of the solution vector.
z’(x) = [u’(x) u’’(x) w’(x) w’’(x) psi’(x), psi’’(x)]^T
- Parameters:
z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x), psi’(x)]^T
phi (float) – Inclination (degrees). Counterclockwise positive.
- Returns:
First derivative z’(x) for the solution vector (6x1).
- Return type:
ndarray, float
- dz_dxdx(z, phi)[source]¶
Get second derivative z’’(x) = K*z’(x) of the solution vector.
z’’(x) = [u’’(x) u’’’(x) w’’(x) w’’’(x) psi’’(x), psi’’’(x)]^T
- Parameters:
z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x), psi’(x)]^T
phi (float) – Inclination (degrees). Counterclockwise positive.
- Returns:
Second derivative z’’(x) = (K*z(x) + q)’ = K*z’(x) = K*(K*z(x) + q) of the solution vector (6x1).
- Return type:
ndarray, float
- eps(Z)[source]¶
Get weak-layer normal strain.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
- Returns:
Weak-layer normal strain epsilon.
- Return type:
float
- gamma(Z)[source]¶
Get weak-layer shear strain.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
- Returns:
Weak-layer shear strain gamma.
- Return type:
float
- int1(x, z0, z1)[source]¶
Get mode I crack opening integrand at integration points xi.
- Parameters:
x (float, ndarray) – X-coordinate where integrand is to be evaluated (mm).
z0 (callable) – Function that returns the solution vector of the uncracked configuration.
z1 (callable) – Function that returns the solution vector of the cracked configuration.
- Returns:
Integrant of the mode I crack opening integral.
- Return type:
float or ndarray
- int2(x, z0, z1)[source]¶
Get mode II crack opening integrand at integration points xi.
- Parameters:
x (float, ndarray) – X-coordinate where integrand is to be evaluated (mm).
z0 (callable) – Function that returns the solution vector of the uncracked configuration.
z1 (callable) – Function that returns the solution vector of the cracked configuration.
- Returns:
Integrant of the mode II crack opening integral.
- Return type:
float or ndarray
- psi(Z, unit='rad')[source]¶
Get midplane rotation psi.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
unit ({'deg', 'degrees', 'rad', 'radians'}, optional) – Desired output unit. Default is radians.
- Returns:
psi – Cross-section rotation psi (radians) of the slab.
- Return type:
float
- sig(Z, unit='MPa')[source]¶
Get weak-layer normal stress.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
unit ({'MPa', 'kPa'}, optional) – Desired output unit. Default is MPa.
- Returns:
Weak-layer normal stress sigma (in specified unit).
- Return type:
float
- tau(Z, unit='MPa')[source]¶
Get weak-layer shear stress.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
unit ({'MPa', 'kPa'}, optional) – Desired output unit. Default is MPa.
- Returns:
Weak-layer shear stress tau (in specified unit).
- Return type:
float
- u(Z, z0, unit='mm')[source]¶
Get horizontal displacement u = u0 + z0 psi.
- Parameters:
Z (ndarray) – Solution vector [u(x) u’(x) w(x) w’(x) psi(x) psi’(x)]^T.
z0 (float) – Z-coordinate (mm) where u is to be evaluated.
unit ({'m', 'cm', 'mm', 'um'}, optional) – Desired output unit. Default is mm.
- Returns:
Horizontal displacement u (unit) of the slab.
- Return type:
float
- class weac.mixins.OutputMixin[source]¶
Bases:
object
Mixin for outputs.
Provides convenience methods for the assembly of output lists such as rasterized displacements or rasterized stresses.
- external_potential(C, phi, L, **segments)[source]¶
Compute total external potential (pst only).
- Parameters:
C (ndarray) – Matrix(6xN) of solution constants for a system of N segements. Columns contain the 6 constants of each segement.
phi (float) – Inclination of the slab (°).
L (float, optional) – Total length of model (mm).
segments (dict) – Dictionary with lists of touchdown booleans (tdi), segement lengths (li), skier weights (mi), and foundation booleans in the cracked (ki) and uncracked (k0) configurations.
- Returns:
Pi_ext – Total external potential (Nmm).
- Return type:
float
- get_slab_deflection(x, z, unit='mm')[source]¶
Compute vertical slab displacement.
- Parameters:
x (ndarray) – Discretized x-coordinates (mm) where coordinates of unsupported (no foundation) segments are NaNs.
z (ndarray) – Solution vectors at positions x as columns of matrix z. Default is mid.
unit ({'m', 'cm', 'mm', 'um'}, optional) – Displacement output unit. Default is mm.
- Returns:
x (ndarray) – Horizontal coordinates (cm).
ndarray – Vertical deflections (unit input).
- get_slab_displacement(x, z, loc='mid', unit='mm')[source]¶
Compute horizontal slab displacement.
- Parameters:
x (ndarray) – Discretized x-coordinates (mm) where coordinates of unsupported (no foundation) segments are NaNs.
z (ndarray) – Solution vectors at positions x as columns of matrix z.
loc ({'top', 'mid', 'bot'}) – Get displacements of top, midplane or bottom of slab. Default is mid.
unit ({'m', 'cm', 'mm', 'um'}, optional) – Displacement output unit. Default is mm.
- Returns:
x (ndarray) – Horizontal coordinates (cm).
ndarray – Horizontal displacements (unit input).
- get_slab_rotation(x, z, unit='degrees')[source]¶
Compute slab cross-section rotation angle.
- Parameters:
x (ndarray) – Discretized x-coordinates (mm) where coordinates of unsupported (no foundation) segments are NaNs.
z (ndarray) – Solution vectors at positions x as columns of matrix z. Default is mid.
unit ({'deg', degrees', 'rad', 'radians'}, optional) – Rotation angle output unit. Default is degrees.
- Returns:
x (ndarray) – Horizontal coordinates (cm).
ndarray – Cross section rotations (unit input).
- get_weaklayer_normalstress(x, z, unit='MPa', removeNaNs=False)[source]¶
Compute weak-layer normal stress.
- Parameters:
x (ndarray) – Discretized x-coordinates (mm) where coordinates of unsupported (no foundation) segments are NaNs.
z (ndarray) – Solution vectors at positions x as columns of matrix z.
unit ({'MPa', 'kPa'}, optional) – Stress output unit. Default is MPa.
keepNaNs (bool) – If set, do not remove
- Returns:
x (ndarray) – Horizontal coordinates (cm).
sig (ndarray) – Normal stress (stress unit input).
- get_weaklayer_shearstress(x, z, unit='MPa', removeNaNs=False)[source]¶
Compute weak-layer shear stress.
- Parameters:
x (ndarray) – Discretized x-coordinates (mm) where coordinates of unsupported (no foundation) segments are NaNs.
z (ndarray) – Solution vectors at positions x as columns of matrix z.
unit ({'MPa', 'kPa'}, optional) – Stress output unit. Default is MPa.
keepNaNs (bool) – If set, do not remove
- Returns:
x (ndarray) – Horizontal coordinates (cm).
sig (ndarray) – Normal stress (stress unit input).
- internal_potential(C, phi, L, **segments)[source]¶
Compute total internal potential (pst only).
- Parameters:
C (ndarray) – Matrix(6xN) of solution constants for a system of N segements. Columns contain the 6 constants of each segement.
phi (float) – Inclination of the slab (°).
L (float, optional) – Total length of model (mm).
segments (dict) – Dictionary with lists of touchdown booleans (tdi), segement lengths (li), skier weights (mi), and foundation booleans in the cracked (ki) and uncracked (k0) configurations.
- Returns:
Pi_int – Total internal potential (Nmm).
- Return type:
float
- total_potential(C, phi, L, **segments)[source]¶
Returns total differential potential
- Parameters:
C (ndarray) – Matrix(6xN) of solution constants for a system of N segements. Columns contain the 6 constants of each segement.
phi (float) – Inclination of the slab (°).
L (float, optional) – Total length of model (mm).
segments (dict) – Dictionary with lists of touchdown booleans (tdi), segement lengths (li), skier weights (mi), and foundation booleans in the cracked (ki) and uncracked (k0) configurations.
- Returns:
Pi – Total differential potential (Nmm).
- Return type:
float
- class weac.mixins.SlabContactMixin[source]¶
Bases:
object
Mixin for handling the touchdown situation in a PST.
Provides Methods for the calculation of substitute spring stiffnesses, cracklength-tresholds and element lengths.
- calc_a1()[source]¶
Calc transition lengths a1 (aAB).
- Returns:
a1 – Length of the crack for transition of stage A to stage B (mm).
- Return type:
float
- calc_a2()[source]¶
Calc transition lengths a2 (aBC).
- Returns:
a2 – Length of the crack for transition of stage B to stage C (mm).
- Return type:
float
- calc_qn()[source]¶
Calc total surface normal load.
- Returns:
Total surface normal load (N/mm).
- Return type:
float
- calc_qt()[source]¶
Calc total surface normal load.
- Returns:
Total surface normal load (N/mm).
- Return type:
float
- set_phi(phi)[source]¶
Set inclination of the slab.
- Parameters:
phi (float) – Inclination of the slab (°).
- set_stiffness_ratio(ratio=1000)[source]¶
Set ratio between collapsed and uncollapsed weak-layer stiffness.
- Parameters:
ratio (int, optional) – Stiffness ratio between collapsed and uncollapsed weak layer. Default is 1000.
- set_tc(cf)[source]¶
Set height of the crack.
- Parameters:
cf (float) – Collapse-factor. Ratio of the crack height to the uncollapsed weak-layer height.
- set_touchdown_attributes(L, a, cf, phi, ratio)[source]¶
Set class attributes for touchdown consideration
- substitute_stiffness(L, support='rested', dof='rot')[source]¶
Calc substitute stiffness for beam on elastic foundation.
- Parameters:
L (float) – Total length of the PST-column (mm).
support (string) – Type of segment foundation. Defaults to ‘rested’.
dof (string) – Type of substitute spring, either ‘rot’ or ‘trans’. Defaults to ‘rot’.
- Returns:
k
- Return type:
stiffness of substitute spring.
- class weac.mixins.SolutionMixin[source]¶
Bases:
object
Mixin for the solution of boundary value problems.
Provides methods for the assembly of the system of equations and for the computation of the free constants.
- assemble_and_solve(phi, li, mi, ki)[source]¶
Compute free constants for arbitrary beam assembly.
Assemble LHS from supported and unsupported segments in the form [ ] [ zh1 0 0 … 0 0 0 ][ ] [ ] [ ] left [ ] [ zh1 zh2 0 … 0 0 0 ][ ] [ ] [ ] mid [ ] [ 0 zh2 zh3 … 0 0 0 ][ ] [ ] [ ] mid [z0] = [ … … … … … … … ][ C ] + [ zp ] = [ rhs ] mid [ ] [ 0 0 0 … zhL zhM 0 ][ ] [ ] [ ] mid [ ] [ 0 0 0 … 0 zhM zhN ][ ] [ ] [ ] mid [ ] [ 0 0 0 … 0 0 zhN ][ ] [ ] [ ] right and solve for constants C.
- Parameters:
phi (float) – Inclination (degrees).
li (ndarray) – List of lengths of segements (mm).
mi (ndarray) – List of skier weigths (kg) at segement boundaries.
ki (ndarray) – List of one bool per segement indicating whether segement has foundation (True) or not (False).
- Returns:
C – Matrix(6xN) of solution constants for a system of N segements. Columns contain the 6 constants of each segement.
- Return type:
ndarray
- bc(z, k=False, pos='mid')[source]¶
Provide equations for free (pst) or infinite (skiers) ends.
- Parameters:
z (ndarray) – Solution vector (6x1) at a certain position x.
l (float, optional) – Length of the segment in consideration. Default is zero.
k (boolean) – Indicates whether segment has foundation(True) or not (False). Default is False.
pos ({'left', 'mid', 'right', 'l', 'm', 'r'}, optional) – Determines whether the segement under consideration is a left boundary segement (left, l), one of the center segement (mid, m), or a right boundary segement (right, r). Default is ‘mid’.
- Returns:
bc – Boundary condition vector (lenght 3) at position x.
- Return type:
ndarray
- calc_segments(li: list[float] | list[int] | bool = False, mi: list[float] | list[int] | bool = False, ki: list[bool] | bool = False, k0: list[bool] | bool = False, L: float = 10000.0, a: float = 0, m: float = 0, phi: float = 0, cf: float = 0.5, ratio: float = 1000, **kwargs)[source]¶
Assemble lists defining the segments.
This includes length (li), foundation (ki, k0), and skier weight (mi).
- Parameters:
li (squence, optional) – List of lengths of segements(mm). Used for system ‘skiers’.
mi (squence, optional) – List of skier weigths (kg) at segement boundaries. Used for system ‘skiers’.
ki (squence, optional) – List of one bool per segement indicating whether segement has foundation (True) or not (False) in the cracked state. Used for system ‘skiers’.
k0 (squence, optional) – List of one bool per segement indicating whether segement has foundation(True) or not (False) in the uncracked state. Used for system ‘skiers’.
L (float, optional) – Total length of model (mm). Used for systems ‘pst-’, ‘-pst’, ‘vpst-’, ‘-vpst’, and ‘skier’.
a (float, optional) – Crack length (mm). Used for systems ‘pst-’, ‘-pst’, ‘pst-‘, ‘-pst’, and ‘skier’.
phi (float, optional) – Inclination (degree).
m (float, optional) – Weight of skier (kg) in the axial center of the model. Used for system ‘skier’.
cf (float, optional) – Collapse factor. Ratio of the crack height to the uncollapsed weak-layer height. Used for systems ‘pst-’, ‘-pst’. Default is 0.5.
ratio (float, optional) – Stiffness ratio between collapsed and uncollapsed weak layer. Default is 1000.
- Returns:
segments – Dictionary with lists of touchdown booleans (tdi), segement lengths (li), skier weights (mi), and foundation booleans in the cracked (ki) and uncracked (k0) configurations.
- Return type:
dict
- eqs(zl, zr, k=False, pos='mid')[source]¶
Provide boundary or transmission conditions for beam segments.
- Parameters:
zl (ndarray) – Solution vector (6x1) at left end of beam segement.
zr (ndarray) – Solution vector (6x1) at right end of beam segement.
k (boolean) – Indicates whether segment has foundation(True) or not (False). Default is False.
pos ({'left', 'mid', 'right', 'l', 'm', 'r'}, optional) – Determines whether the segement under consideration is a left boundary segement (left, l), one of the center segement (mid, m), or a right boundary segement (right, r). Default is ‘mid’.
- Returns:
eqs – Vector (of length 9) of boundary conditions (3) and transmission conditions (6) for boundary segements or vector of transmission conditions (of length 6+6) for center segments.
- Return type:
ndarray