ContactMechanics.ReferenceSolutions package
Submodules
ContactMechanics.ReferenceSolutions.Cone module
- ContactMechanics.ReferenceSolutions.Cone.contact_radius_and_area(alpha)
Reference: K.L. Johnson, Contact Mechanics, Cambridge University Press. Chapter 5, page 112, fig(a)”
- Parameters:
alpha (float) – half of Cone angle
- Returns:
Ratio_contact_radius (float) – Contact Radius / penetration –> R/D
Ratio_Area (float) – Contact Area / Penetration**2 –> A/D**2
- ContactMechanics.ReferenceSolutions.Cone.deformation(penetration, alpha)
Reference: K.L. Johnson, Contact Mechanics, Cambridge University Press. Chapter 5, page 112, fig(a)”
- Parameters:
penetration (float) – Radius / Penetration –> R / P
alpha (float) – half of Cone angle
- Returns:
Ratio_Deformation – Ratio Deformation –> Deformation / Penetration
- Return type:
float
- ContactMechanics.ReferenceSolutions.Cone.load_and_mean_pressure(alpha)
Reference: K.L. Johnson, Contact Mechanics, Cambridge University Press. Chapter 5, page 112, fig(a)”
- Parameters:
alpha (float) – half of Cone angle
- Returns:
Ratio_ExternalForece (float) – External Force / (Young’s module * contact area) –> F/(E*A)
Ratio_Mean_Pressure (float) – Mean Pressure / Young’s module –> P/E
- ContactMechanics.ReferenceSolutions.Cone.pressure(mean_pressure, alpha)
Reference: K.L. Johnson, Contact Mechanics, Cambridge University Press. Chapter 5, page 112, fig(a)”
- Parameters:
mean_pressure (float) – Ratio of Pressure and Mean Pressure –> Pressure / Mean Pressure
alpha (float) – half of cone angle
- Returns:
Ratio_Pressure – Ratio Pressure / Mean Pressure
- Return type:
float
ContactMechanics.ReferenceSolutions.GreenwoodTripp module
Greenwood-Tripp model for the contact of rough spheres
- ContactMechanics.ReferenceSolutions.GreenwoodTripp.Fn(h, n, phi=<function <lambda>>)
- ContactMechanics.ReferenceSolutions.GreenwoodTripp.GreenwoodTripp(d, μ, rhomax=5, n=100, eps=1e-06, tol=1e-06, mix=0.1, maxiter=1000)
Greenwood-Tripp solution for the contact of rough spheres. See: Greenwood, Tripp, J. Appl. Mech. 34, 153 (1967) Symbols correspond to notation in this paper.
- Parameters:
d (float) – Dimensionless displacement d* of the sphere, d* = d / σ where σ is the standard deviation of the height distribution (i.e. the rms height).
μ (float) – Material parameter, μ=8/3 η σ sqrt(2 B β) where η is the density of asperities, B is the sphere radius and β the asperity radius.
rhomax (float) – Radius cutoff (in units of σ). Pressure is assumed to be zero beyond this radius.
n (int) – Number of grid points.
eps (float) – Small number for numerical regularization.
tol (float) – Max difference between displacements in consecutive solutions.
mix (float) – Mixing of solution between consecutive steps.
- Returns:
w (array) – Displacement in radial direction in units of σ.
p (array) – Pressure in radial direction in units of E* sqrt(σ / 8B) where E* is the contact modulus.
rho (array) – Radial coordinate in units of sqrt(2 B σ)
- ContactMechanics.ReferenceSolutions.GreenwoodTripp.s(ξ, ξ1=1e-05, ξ2=100000.0)
- ContactMechanics.ReferenceSolutions.GreenwoodTripp.ξ(s, s1=1e-05, s2=100000.0)
ContactMechanics.ReferenceSolutions.Hertz module
Helper tools for PyCo
- ContactMechanics.ReferenceSolutions.Hertz.centerline_stress(z, poisson=0.5)
Given distance from the center of the sphere, contact radius and Poisson number contact, compute the stress at the surface.
- Parameters:
z (array_like) – Array of depths (from the center of the sphere in units of contact radius a).
poisson (float) – Poisson number.
- Returns:
srr (array) – Radial stress (in units of maximum pressure p0).
szz (array) – Contact pressure (in units of maximum pressure p0).
- ContactMechanics.ReferenceSolutions.Hertz.elastic_energy(d, R, Es)
- Parameters:
d (float) – Rigid Body Penetration
R (float) – Sphere Radius
Es (float) – Contact modulus: Es = E/(1-nu**2) with Young’s modulus E and Poisson number nu.
- Return type:
Elastic energy
- ContactMechanics.ReferenceSolutions.Hertz.max_pressure__penetration(delta, R=1, Es=1)
when R and Es is not specified, delta is supposed to be in units of R and p0 in units of Es :param delta: :type delta: penetration :param R: :type R: Radius :param Es: :type Es: contact modulus
- Return type:
max pressure p0
- ContactMechanics.ReferenceSolutions.Hertz.normal_load(d, R, Es)
Given rigid Body Penetration, sphere radius and contact modulus compute normal load.
- Parameters:
d (float) – Rigid Body Penetration
R (float) – Sphere Radius
Es (float) – Contact modulus: Es = E/(1-nu**2) with Young’s modulus E and Poisson number nu.
- Returns:
float
normal_load
- ContactMechanics.ReferenceSolutions.Hertz.penetration(N, R, Es)
Given normal load, sphere radius and contact modulus compute rigid body penetration.
- Parameters:
N (float) – Normal force.
R (float) – Sphere radius.
Es (float) – Contact modulus: Es = E/(1-nu**2) with Young’s modulus E and Poisson number nu.
- Returns:
float
Normal load
- ContactMechanics.ReferenceSolutions.Hertz.radius_and_pressure(N, R, Es)
Given normal load, sphere radius and contact modulus compute contact radius and peak pressure.
- Parameters:
N (float) – Normal force.
R (float) – Sphere radius.
Es (float) – Contact modulus: Es = E/(1-nu**2) with Young’s modulus E and Poisson number nu.
- Returns:
a (float) – Contact radius.
p0 (float) – Maximum pressure inside the contacting area (right under the apex).
- ContactMechanics.ReferenceSolutions.Hertz.stress(r, z, poisson=0.5)
Return components of the stress tensor in the interior of the Hertz solid. This is the solution given by: M.T. Huber, Ann. Phys. 319, 153 (1904)
Note that the stress tensor at any point in the solid has the form below in a cylindrical coordinate system centered at the tip apex. Some off-diagonal components are zero by rotational symmetry. stt is the circumferential component, srr the radial component and szz the normal component of the stress tensor.
/ stt 0 0
- s = | 0 srr srz |
0 srz szz /
- Parameters:
r (array_like) – Radial position (in units of the contact radius a).
z (array_like) – Depth (in units of the contact radius a).
poisson (float) – Poisson number.
- Returns:
stt (array) – Circumferential component of the stress tensor (in units of maximum pressure p0).
srr (array) – Radial component of the stress tensor (in units of maximum pressure p0).
szz (array) – Normal component of the stress tensor (in units of maximum pressure p0).
srz (array) – Shear component of the stress tensor (in units of maximum pressure p0).
- ContactMechanics.ReferenceSolutions.Hertz.stress_Cartesian(x, y, z, poisson=0.5)
Return components of the stress tensor in the interior of solid due to normal Hertz loading. This is the solution given by: G.M. Hamilton, Proc. Instn. Mech. Engrs. 197C, 53-59 (1983)
- Parameters:
x (array_like) – In-plane positions (in units of the contact radius a).
y (array_like) – In-plane positions (in units of the contact radius a).
z (array_like) – Depth (in units of the contact radius a).
poisson (float) – Poisson number.
- Returns:
sxx, syy, szz, syz, sxz, sxy – Individual components of the Cartesian stress tensor.
- Return type:
array
- ContactMechanics.ReferenceSolutions.Hertz.stress_for_tangential_loading(x, y, z, poisson=0.5)
Return components of the stress tensor in the interior of solid due to tangential (Hertz) loading. This is the solution given by: G.M. Hamilton, Proc. Instn. Mech. Engrs. 197C, 53-59 (1983)
- Parameters:
x (array_like) – In-plane positions (in units of the contact radius a).
y (array_like) – In-plane positions (in units of the contact radius a).
z (array_like) – Depth (in units of the contact radius a).
poisson (float) – Poisson number.
- Returns:
sxx, syy, szz, syz, sxz, sxy – Individual components of the Cartesian stress tensor.
- Return type:
array
- ContactMechanics.ReferenceSolutions.Hertz.surface_displacements(r)
Return the displacements at the surface due to an indenting sphere. See: K.L. Johnson, Contact Mechanics, p. 61
- Parameters:
r (array_like) – Radial position normalized by contact radius a.
- Returns:
uz – Normal displacements at the surface of the contact (in units of p0/Es * a where p0 is maximum pressure, Es contact modulus and a contact radius).
- Return type:
array
- ContactMechanics.ReferenceSolutions.Hertz.surface_stress(r, poisson=0.5)
Given distance from the center of the sphere, contact radius and Poisson number contact, compute the stress at the surface.
- Parameters:
r (array_like) – Array of distance (from the center of the sphere in units of contact radius a).
poisson (float) – Poisson number.
- Returns:
pz (array) – Contact pressure (in units of maximum pressure p0).
sr (array) – Radial stress (in units of maximum pressure p0).
stheta (array) – Azimuthal stress (in units of maximum pressure p0).
ContactMechanics.ReferenceSolutions.Westergaard module
Westergaard solution for partial contact of a sinusoidal periodic indenter. See: H.M. Westergaard, Trans. ASME, J. Appl. Mech. 6, 49 (1939)
The indenter geometry is
Nommenclature:
\(E^*\): contact modulus
\(\lambda\): period of the sinusoidal indenter
\(h\): amplitude of the sinusoidal indenter. peak to tale distance is \(2h\)
\(a\): half contact width, “contact radius”
- ContactMechanics.ReferenceSolutions.Westergaard.contact_radius(mean_pressure)
Johnson, K. L. International Journal of Solids and Structures 32, 423–430 (1995)
Equation (4)
- Parameters:
mean_pressure (float, array_like) – mean pressure in units of \(\pi E^* h / \lambda\)
- Return type:
contact radius in units of \(\lambda\)
- ContactMechanics.ReferenceSolutions.Westergaard.displacements(x, a)
Displacements u in units of 2h, the geometry is \(2h sin^2(pi x / lambda)\)
- Parameters:
x (float or np.array) – distance from the indenter tip contact point in units of lambda
a (float) – half contact length in units of lambda
- Return type:
Displacents u in units of 2h
Examples
>>> import matplotlib.pyplot as plt >>> fig, ax = plt.subplots() >>> x = np.linspace(-1.,1., 500) >>> l,= ax.plot(x, displacements(x, 0.1), label="displacements") >>> l2, = ax.plot(x, -np.sin(np.pi * x)**2, "--k", label="geometry") >>> ax.invert_yaxis() >>> leg = ax.legend() >>> plt.show(block=True)
- ContactMechanics.ReferenceSolutions.Westergaard.elastic_energy(mean_pressure)
elastic energy along the equilibrium curve
formula for the JKR contact but with alpha = 0
- Parameters:
mean_pressure (float or np.array) – mean pressure in units of \(\pi Es h/ lambda\)
- Returns:
energy per unit area in units of \(h p_{wfc}\)
with \(p_{wfc} = \pi E^* h/\lambda\)
- ContactMechanics.ReferenceSolutions.Westergaard.elastic_energy_a(a)
elastic energy in dependence of a along the equilibrium curve
- Parameters:
a (float or np.array) – half contact width in units of the wavelength \(\lambda\)
- Returns:
energy per unit area in units of \(h p_{wfc}\)
with \(p_{wfc} = pi E^* h/\lambda\)
- ContactMechanics.ReferenceSolutions.Westergaard.gap(x, a)
gap g in units of 2h, the geometry is \(2h sin^2(pi x / lambda)\)
- Parameters:
x (float or np.array) – distance from the first contact point in units of lambda
a (float) – half contact length in units of lambda
- Return type:
Displacents u in units of 2h
- ContactMechanics.ReferenceSolutions.Westergaard.mean_displacement(a)
mean gap would be 1 / 2 + mean displacement
- Parameters:
a (float or np.array) – half contactwidth in units of the wavelength lambda
- Return type:
mean displacement in units of 2 h
Examples
>>> import matplotlib.pyplot as plt >>> a = np.linspace(0,0.5) >>> l,= plt.plot(a, mean_displacement(a)) >>> plt.show()
- ContactMechanics.ReferenceSolutions.Westergaard.mean_pressure(contact_radius)
Johnson, K. L. International Journal of Solids and Structures 32, 423–430 (1995)
Equation (4)
- Parameters:
contact_radius (float) – in units of \(\lambda\)
- Returns:
mean pressure – mean pressure in units of \(\pi E^* h / \lambda\)
- Return type:
float, array_like
Module contents
Analytical or semi-analytical reference solutions