Developer Reference

This page gathers all docstrings of all classes (including class hierarchy), its members and separate functions in the submodules of assess.

assess.cylindrical module

Analytical solutions in cylindrical domains.

class assess.cylindrical.AnalyticalStokesSolution(Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: object

Base class for solutions in sperical or cylindrical shell domains.

Parameters:
  • Rp – outer radius
  • Rm – inner radius
  • g – magnitude of source term
  • nu – viscosity
class assess.cylindrical.CylindricalStokesSolution(n, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.cylindrical.AnalyticalStokesSolution

Base class for solutions in cylindrical shell domains.

This implements analytical solutions based on a streamfunction of the form.

\[\psi(r,\varphi) = \psi_r(r)\sin(n\varphi)\]

where \(\psi_r\) should be defined in a method psi_r() and its derivative in dpsi_rdr().

Initialize basic parameters of analytical solution in cylindrical domain.

dpsi_rdr(r)

Abstract method to be implemented by subclass.

Parameters:r – radius
pressure_cartesian(X)

Return pressure solution at Cartesian location.

Parameters:X – 2D Cartesian location
psi_r(r)

Abstract method to be implemented by subclass.

Parameters:r – radius
radial_stress(r, phi)

Return radial component of stress.

Parameters:
  • r – radius
  • phi – angle with x-axis.
radial_stress_cartesian(X)

Return radial component of stress at Cartesian location.

Parameters:X – 2D Cartesian location
tau_rphi(r, phi)

Return shear stress \(\tau_{r\varphi}\).

Parameters:
  • r – radius
  • phi – angle with x-axis.
tau_rr(r, phi)

Return radial component of deviatoric stress.

Parameters:
  • r – radius
  • phi – angle with x-axis.
u_phi(r, phi)

Return tangential component of velocity.

Parameters:
  • r – radius
  • phi – angle with x-axis.
u_r(r, phi)

Return radial component of velocity.

Parameters:
  • r – radius
  • phi – angle with x-axis.
velocity_cartesian(X)

Return Cartesian velocity at Cartesian location.

Parameters:X – 2D Cartesian location
class assess.cylindrical.CylindricalStokesSolutionDelta(ABCD, n, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.cylindrical.CylindricalStokesSolution

Base class for solutions in cylindrical shell domains with delta(r-r’) forcing

This implements the analytical solution in one half (above or below r’) of the domain which is based on a biharmonic streamfunction

\[\psi(r,\varphi) = \psi_r(r)\sin(n\varphi)\]

determinded by 4 coefficients, A, B, C, and D.

Parameters:
  • ABCD – list or tuple of 4 floats, coefficients for the 4 biharmonic solutions of the streamfunction
  • n – wave number
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
dpsi_rdr(r)

Radial derivative of radial part of streamfunction

Parameters:r – radius
dpsi_rdr2(r)

Second radial derivative of radial part of streamfunction

Parameters:r – radius
p(r, phi)

Pressure solution

Parameters:
  • r – radius
  • phi – angle with x-axis
psi_r(r)

Radial part of streamfunction

Parameters:r – radius
class assess.cylindrical.CylindricalStokesSolutionDeltaFreeSlip(n, sign, Rp=2.22, Rm=1.22, rp=1.72, nu=1.0, g=1.0)

Bases: assess.cylindrical.CylindricalStokesSolutionDelta

Analytical Solution in cylindrical domain with delta(r-r’) forcing and free-slip boundary conditions

Parameters:
  • n – wave number
  • sign – +1 for upper half solution r’<r<Rp -1 for lower half solution Rm<r<r’
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
class assess.cylindrical.CylindricalStokesSolutionDeltaZeroSlip(n, sign, Rp=2.22, Rm=1.22, rp=1.72, nu=1.0, g=1.0)

Bases: assess.cylindrical.CylindricalStokesSolutionDelta

Analytical Solution in cylindrical domain with delta(r-r’) forcing and zero-slip boundary conditions

Parameters:
  • n – wave number
  • sign – +1 for upper half solution r’<r<Rp -1 for lower half solution Rm<r<r’
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
class assess.cylindrical.CylindricalStokesSolutionSmooth(ABCDE, k, n, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.cylindrical.CylindricalStokesSolution

Base class for solutions in cylindrical shell domains with r^k forcing

The analytical solution is based on a streamfunction

determinded by 4 coefficients, A, B, C, and D corresponding to 4 independent biharmonic (i.e. homogenous) solutions and a fifth coefficient E that corresponds to the inhomogeneous part.

Parameters:
  • ABCDE – list or tuple of 5 floats, coefficients for the 4 biharmonic and one inhomogenous solutions of the streamfunction
  • k – polynomial order of forcing
  • n – wave number
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
delta_rho(r, phi)

Perturbation density \(\rho'\) in forcing term: \(g\rho'\hat r\)

Parameters:
  • r – radius
  • phi – angle with x-axis
delta_rho_cartesian(X)

Perturbation density \(\rho'\) in forcing term: \(g\rho'\hat r\)

Parameters:X – 2D Cartesian coordinate
dpsi_rdr(r)

Radial derivative of radial part of streamfunction

Parameters:r – radius
dpsi_rdr2(r)

Second radial derivative of radial part of streamfunction

Parameters:r – radius
p(r, phi)

Pressure solution

Parameters:
  • r – radius
  • phi – angle with x-axis
psi_r(r)

Radial part of streamfunction

Parameters:r – radius
class assess.cylindrical.CylindricalStokesSolutionSmoothFreeSlip(n, k, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.cylindrical.CylindricalStokesSolutionSmooth

Analytical Solution in cylindrical domain with smooth r^k forcing and free-slip boundary conditions

Parameters:
  • n – wave number
  • k – polynomial order of forcing
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
class assess.cylindrical.CylindricalStokesSolutionSmoothZeroSlip(n, k, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.cylindrical.CylindricalStokesSolutionSmooth

Analytical Solution in cylindrical domain with smooth r^k forcing and zero-slip boundary conditions

Parameters:
  • n – wave number
  • k – polynomial order of forcing
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength

assess.delta module

This module contains functions that compute the coefficients in analytical solution to the Stokes equations in cylindrical and spherical shell domains with a delta function RHS forcing term, which in the 2D cylindrical shell domain takes the form:

\[f = -g \delta(r-r') cos(n\varphi) \hat{r}\]

where r is the radius (distance from origin), \(\varphi\) is the angle with the x-axis, \(r'\) is the radius of the density anomaly, and \(\hat{r}\) is the outward pointing radial unit vector. The magnitude g, degree l and wave number n can be chosen arbitrarily. Similarly in a 3D spherical domain we consider the forcing term

\[f = -g \delta(r-r') Y_{lm}(\theta,\varphi) \hat{r}\]

where \(\theta\) and \(\varphi\) are the co-latitude and longitude respectively, and \(Y_{lm}\) is Laplace’s spherical harmonic function of degree l and order m.

We consider the following cases:

cylinder_delta_fs: cylindrical, with free-slip boundary condition at r=Rm and r=Rp cylinder_delta_ns: cylindrical, with zero-slip boundary conditions at r=Rm and r=Rp sphere_delta_fs: spherical, with free-slip boundary condition at r=Rm and r=Rp sphere_smooh_ns: spherical, with zero-slip boundary condition at r=Rm and r=Rp

The functions below take the following parameters:

param Rp:outer radius
param Rm:inner radius
param rp:radius of density anomaly with Rm<rp<Rp
param n:wave number (2D cylindrical only)
param l:spherical degree (3D only)
param g:scalar magnitude of forcing
param nu:viscosity
param sign:+1 or -1 for the upper (r>rp) and lower (r<rp) half of the solution

and return four coefficients A, B, C, D for the linear combination of biharmonic spherical functions that constitute the analytical solution.

This module has been automatically generated by extracting the solutions from the latex source of the associated paper. If sign=+1 we use the following substitutions

alpha_pm = Rp/rp, alpha_mp = Rm/rp, pm=+1, mp=-1

i.e. we take the top case for the symbols alpha_pm, alpha_mp, pm, and mp in the paper to obtain \(A_+, B_+, C_+,\) and \(D_+\) corresponding to the upper half of the domain. If sign=-1 we have

alpha_pm = Rm/rp, alpha_mp = Rp/rp, pm=-1, mp=+1

to obtain the \(A_-, B_-, C_-,\) and \(D_-\) coefficients of the lower half of the solution.

assess.delta.coefficients_cylinder_delta_fs(Rp, Rm, rp, n, g, nu, sign)
assess.delta.coefficients_cylinder_delta_ns(Rp, Rm, rp, n, g, nu, sign)
assess.delta.coefficients_sphere_delta_fs(Rp, Rm, rp, l, g, nu, sign)
assess.delta.coefficients_sphere_delta_ns(Rp, Rm, rp, l, g, nu, sign)

assess.smooth module

This module contains functions that compute the coefficients in analytical solution to the Stokes equations in cylindrical and spherical shell domains with a smooth RHS forcing term, which in the 2D cylindrical shell domain takes the form:

\[f = -g (r/Rp)^k cos(n\varphi) \hat{r}\]

where r is the radius (distance from origin), \(\varphi\) is the angle with the x-axis, Rp and Rm are the outer and inner radius of the shell domain, and \(\hat{r}\) is the outward pointing radial unit vector. The magnitude g, degree l and wave number n can be chosen arbitrarily. Similarly in a 3D spherical domain we consider the forcing term

\[f = -g (r/Rp)^k Y_{lm}(\theta,\varphi) \hat{r}\]

where \(\theta\) and \(\varphi\) are the co-latitude and longitude respectively, and \(Y_{lm}\) is Laplace’s spherical harmonic function of degree l and order m.

We consider the following cases:

cylinder_smooth_fs: cylindrical, with free-slip boundary condition at r=Rm and r=Rp cylinder_smooth_ns: cylindrical, with zero-slip boundary conditions at r=Rm and r=Rp sphere_smooth_fs: spherical, with free-slip boundary condition at r=Rm and r=Rp sphere_smooth_ns: spherical, with zero-slip boundary condition at r=Rm and r=Rp

The functions below take the following parameters:

param Rp:outer radius
param Rm:inner radius
param k:polynomial degree of forcing term
param n:wave number (2D cylindrical only)
param l:spherical degree (3D only)
param g:scalar magnitude of forcing
param nu:viscosity

and return the five coefficients A, B, C, D, and E used in the analytical solution.

This module has been automatically generated by extracting the solutions from the latex source of the associated paper.

assess.smooth.coefficients_cylinder_smooth_fs(Rp, Rm, k, n, g, nu)
assess.smooth.coefficients_cylinder_smooth_ns(Rp, Rm, k, n, g, nu)
assess.smooth.coefficients_sphere_smooth_fs(Rp, Rm, k, l, g, nu)
assess.smooth.coefficients_sphere_smooth_ns(Rp, Rm, k, l, g, nu)

assess.spherical module

Analytical solutions in spherical shell domains.

class assess.spherical.SphericalStokesSolution(l, m, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.cylindrical.AnalyticalStokesSolution

Base class for solutions in spherical shell domains.

This implements analytical solutions based on a poloidal function of the form.

\[\mathcal{P}(r,\theta,\varphi) = \mathcal{P}_l(r)Y_{lm}(\theta, \varphi)\]

so that the velocity solution takes the form:

\[\mathbf{u} = \nabla\times\left(\mathbf{r}\times\nabla\mathcal{P}\right)\]

The function \(\mathcal{P}_l\) and its derivative should be implement in methods Pl() and dPldr() in the subclass.

Pl(r)

Abstract method to be implemented by subclass.

Parameters:r – radius
dPldr(r)

Abstract method to be implemented by subclass.

Parameters:r – radius
pressure_cartesian(X)

Return pressure solution at Cartesian location.

Parameters:X – 3D Cartesian location
radial_stress(r, theta, phi)

Return radial component of stress.

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
radial_stress_cartesian(X)

Return radial component of stress at Cartesian location.

Parameters:X – 3D Cartesian location
tau_rphi(r, theta, phi)

Return longitudinal component of shear stress.

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
tau_rr(r, theta, phi)

Return radial component of deviatoric stress.

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
tau_rtheta(r, theta, phi)

Return colatitudinal component of shear stress.

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
u_phi(r, theta, phi)

Return longitudinal (eastward) component of velocity.

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
u_r(r, theta, phi)

Return radial component of velocity.

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
u_theta(r, theta, phi)

Return colatitudinal (southward) component of velocity.

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
velocity_cartesian(X)

Return Cartesian velocity at Cartesian location.

Parameters:X – 3D Cartesian location
class assess.spherical.SphericalStokesSolutionDelta(ABCD, l, m, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.spherical.SphericalStokesSolution

Base class for solutions in spherical shell domains with delta(r-r’) forcing

This implements the analytical solution in one half (above or below r’) of the domain which is based on a poloidal function

\[\mathcal{P}(r,\theta,\varphi) = \mathcal{P}_l(r)Y_{lm}(\theta, \varphi)\]

and velocity

\[\mathbf{u} = \nabla\times\left(\mathbf{r}\times\nabla\mathcal{P}\right)\]

where for biharmonic solutions, \(\mathcal{P}_l(r)\) is determinded by four coefficients, A, B, C, and D.

Parameters:
  • ABCD – list or tuple of 4 floats, coefficients for the 4 independent biharmonic solutions.
  • l – degree of the harmonic
  • m – order of the harmonic
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
Pl(r)

Radial part of poloidal function

Parameters:r – radius
dPldr(r)

Radial derivative of radial part of poloidal function

Parameters:r – radius
dPldr2(r)

Second radial derivative of radial part of poloidal function

Parameters:r – radius
p(r, theta, phi)

Pressure solution

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
class assess.spherical.SphericalStokesSolutionDeltaFreeSlip(l, m, sign, Rp=2.22, Rm=1.22, rp=1.72, nu=1.0, g=1.0)

Bases: assess.spherical.SphericalStokesSolutionDelta

Analytical Solution in cylindrical domain with delta(r-r’) forcing and free-slip boundary conditions

Parameters:
  • l – degree of the harmonic
  • m – order of the harmonic
  • sign – +1 for upper half solution r’<r<Rp -1 for lower half solution Rm<r<r’
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
class assess.spherical.SphericalStokesSolutionDeltaZeroSlip(l, m, sign, Rp=2.22, Rm=1.22, rp=1.72, nu=1.0, g=1.0)

Bases: assess.spherical.SphericalStokesSolutionDelta

Analytical Solution in cylindrical domain with delta(r-r’) forcing and zero-slip boundary conditions

Parameters:
  • l – degree of the harmonic
  • m – order of the harmonic
  • sign – +1 for upper half solution r’<r<Rp -1 for lower half solution Rm<r<r’
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
class assess.spherical.SphericalStokesSolutionSmooth(ABCDE, k, l, m, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.spherical.SphericalStokesSolution

Base class for solutions in spherical shell domains with r^k forcing

This implements the analytical solution in one half (above or below r’) of the domain which is based on a poloidal function

\[\mathcal{P}(r,\theta,\varphi) = \mathcal{P}_l(r)Y_{lm}(\theta, \varphi)\]

and velocity

\[\mathbf{u} = \nabla\times\left(\mathbf{r}\times\nabla\mathcal{P}\right)\]

where the solution \(\mathcal{P}_l(r)\) is determinded by four coefficients, A, B, C, and D of four independent biharmonic solutions, and one coefficient E associated with the inhomogenous part.

Parameters:
  • ABCDE – list or tuple of 5 floats, coefficients for the 4 biharmonic and one inhomogenous solution.
  • l – degree of the harmonic
  • m – order of the harmonic
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
Pl(r)

Radial part of poloidal function

Parameters:r – radius
dPldr(r)

Radial derivative of radial part of poloidal function

Parameters:r – radius
dPldr2(r)

Second radial derivative of radial part of poloidal function

Parameters:r – radius
delta_rho(r, theta, phi)

Perturbation density \(\rho'\) in forcing term: \(g\rho'\hat r\)

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
delta_rho_cartesian(X)

Perturbation density \(\rho'\) in forcing term: \(g\rho'\hat r\)

Parameters:X – 3D Cartesian coordinate
p(r, theta, phi)

Pressure solution

Parameters:
  • r – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
class assess.spherical.SphericalStokesSolutionSmoothFreeSlip(l, m, k, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.spherical.SphericalStokesSolutionSmooth

Analytical Solution in cylindrical domain with smooth r^k forcing and free-slip boundary conditions

Parameters:
  • l – degree of the harmonic
  • m – order of the harmonic
  • k – polynomial order of forcing
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
class assess.spherical.SphericalStokesSolutionSmoothZeroSlip(l, m, k, Rp=2.22, Rm=1.22, nu=1.0, g=1.0)

Bases: assess.spherical.SphericalStokesSolutionSmooth

Analytical Solution in cylindrical domain with smooth r^k forcing and zero-slip boundary conditions

Parameters:
  • l – degree of the harmonic
  • m – order of the harmonic
  • k – polynomial order of forcing
  • Rp – outer radius
  • Rm – inner radius
  • nu – viscosity
  • g – forcing strength
assess.spherical.Y(l, m, theta, phi)

Real-valued spherical harmonic function \(Y_{lm}(\theta, \varphi)\)

This is based on the following definition:

\[Y_{lm}(\theta, \varphi) = \sqrt{\frac{(2l+1)}{4\pi}\frac{(l-m)!}{(l+m)!}}P_l^m(\cos(\theta))\cos(m\varphi)\]

which is equal to the real part of scipy.special.sph_harm.

Parameters:
  • l – degree of the harmonic
  • m – order of the harmonic
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
assess.spherical.Y_cartesian(l, m, X)

Real-valued spherical harmonic function that takes Cartesian coordinates

See Y() :param m: order of the harmonic :param l: degree of the harmonic :param X: Cartesian 3D coordinates

assess.spherical.dYdphi(l, m, theta, phi)

Colatitudinal derivative of spherical harmonic function \(Y_{lm}(\theta, \varphi)\)

Parameters:
  • l – degree of the harmonic
  • m – order of the harmonic
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
assess.spherical.dYdtheta(l, m, theta, phi)

Longitudinal derivative of spherical harmonic function \(Y_{lm}(\theta, \varphi)\)

Parameters:
  • l – degree of the harmonic
  • m – order of the harmonic
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
assess.spherical.from_spherical(r, theta, phi)

Convert spherical r, theta, phi to 3D Cartesian coordinates.

Parameters:
  • r, – radius
  • theta – co-latitude in [0, pi]
  • phi – longitude in [0, 2*pi]
Returns X:

Cartesian 3D coordinates

assess.spherical.to_spherical(X)

Convert Cartesian 3D coordinates X to spherical r, theta, phi

Parameters:X – Cartesian 3D coordinates
Returns r, theta, phi:
 radius, colatitude, longitude