Code reference

module-level

ggf.fiber_distance_capillary(gel_thickness=2e-06, glass_thickness=4e-05, channel_width=4e-05, gel_index=1.449, glass_index=1.474, medium_index=1.335)[source]

Effective distance between the two optical fibers

When the optical stretcher is combined with a microfluidic channel (closed setup), then the effective distance between the two optical fibers (with the the stretched object at the channel center) is defined by the refractive indices of the optical components: index matching gel between fiber and channel wall, microfluidic glass channel wall, and medium inside the channel.

Parameters:
  • gel_thickness (float) – Thickness of index matching gel (distance between fiber and glass wall) [m]
  • glass_thickness (float) – Thickness of glass wall [m]
  • channel_width (float) – Width of the microfluidic channed [m]
  • gel_index (float) – Refractive index of index matching gel
  • glass_index (float) – Refractive index of channel glass wall
  • medium_index (float) – Refractive index of index medium inside channel
Returns:

eff_dist – Effective distance between the fibers

Return type:

float

Notes

The effective distance is computed relative to the medium, i.e. if gel_index == glass_index == medium_index, then eff_dist = 2*`gel_dist` + 2*`glass_dist` + channel_width.

ggf.get_ggf(model, semi_major, semi_minor, object_index, medium_index, effective_fiber_distance=0.0001, mode_field_diameter=3e-06, power_per_fiber=0.6, wavelength=1.064e-06, poisson_ratio=0.5, n_poly=120, use_lut=None, verbose=False)[source]

Model the global geometric factor

Parameters:
  • model (str) – Model to use, one of: boyde2009
  • semi_major (float) – Semi-major axis of an ellipse fit to the object perimeter [m]
  • semi_minor (float) – Semi-minor axis of an ellipse fit to the object perimeter [m]
  • object_index (float) – Refractive index of the object
  • medium_index (float) – Refractive index of the surrounding medium
  • effective_fiber_distance (float) – Effective distance between the two trapping fibers relative to the medium refractive index [m]. For an open setup, this is the physical distance between the fibers. For a closed setup (capillary), this distance takes into account the refractive indices and thicknesses of the glass capillary and index matching gel. For the closed setup, the convenience function ggf.fiber_distance_capillary() can be used.
  • mode_field_diameter (float) – The mode field diameter MFD of the fiber used [m]. Note that the MFD is dependent on the wavelength used. If the manufacturer did not provide a value for the MFD, the MFD can be approximated as 3*wavelenth for a single-mode fiber.
  • power_per_fiber (float) – The laser power coupled into each of the fibers [W]
  • wavelength (float) – The laser wavelength used for the trap [m]
  • poisson_ratio (float) – The Poisson’s ratio of the stretched material. Set this to 0.5 for volume conservation.
  • n_poly (int) – Number of Legendre polynomials to use for computing the GGF. Note that only even Legendre polynomials are used and thus, this number is effectively halved. To reproduce the GGF as computed with the Boyde2009 Matlab script, set this value to None.
  • use_lut (None, str, pathlib.Path or bool) – Use look-up tables to compute the GGF. If set to None, the internal LUTs will be used or the GGF is computed if it cannot be found in a LUT. If True, the internal LUTs will be used and a NotInLUTError will be raised if the GGF cannot be found there. Alternatively, a path to a LUT hdf5 file can be passed.
  • verbose (int) – Increases verbosity
Returns:

ggf – global geometric factor

Return type:

float

ggf.legendre2ggf(coeff, poisson_ratio)[source]

Compute the global geometric factor from Legendre coefficients

The definition of the Legendre coefficients is given in the theory section.

Parameters:
  • coeff (1d ndarray) – Legendre coefficients as defined in [Lure64]
  • poisson_ratio (float) – Poisson’s ratio of the stretched material. Set this to 0.5 for volume conservation.
Returns:

ggf – Global geometric factor

Return type:

float

Notes

All odd Legendre coefficients are assumed to be zero, because the stress profile is symmetric with respect to the stretcher axis.

ggf.stress2ggf(stress, theta, poisson_ratio, n_poly=120)[source]

Compute the GGf from radial stress using Legendre decomposition

Parameters:
  • stress (1d ndarray) – Radial stress profile (in imaging plane)
  • theta (1d ndarray) – Polar angles corresponding to stress
  • poisson_ratio (float) – Poisson’s ratio of the stretched material. Set this to 0.5 for volume conservation.
  • n_poly (int) – Number of Legendre polynomials to use
Returns:

ggf – Global geometric factor

Return type:

float

Notes

All odd Legendre coefficients are assumed to be zero, because the stress profile is symmetric with respect to the stretcher axis. Therefore, only n_poly/2 polynomials are considered.

ggf.stress2legendre(stress, theta, n_poly)[source]

Decompose stress into even Legendre Polynomials

The definition of the Legendre decomposition is given in the theory section.

Parameters:
  • stress (1d ndarray) – Radial stress profile (in imaging plane)
  • theta (1d ndarray) – Polar angles corresponding to stress
  • n_poly (int) – Number of Legendre polynomials to use
Returns:

coeff – Legendre coefficients as defined in [Lure64]

Return type:

1d ndarray

Notes

All odd Legendre coefficients are assumed to be zero, because the stress profile is symmetric with respect to the stretcher axis. Therefore, only n_poly/2 polynomials are considered.

matlab_funcs

Special functions translated from Matlab to Python

ggf.matlab_funcs.besselh(n, z)[source]

Hankel function with k = 1

Parameters:
  • n (int) – real order
  • z (float) – complex argument

Notes

https://de.mathworks.com/help/matlab/ref/besselh.html

ggf.matlab_funcs.besselj(n, z)[source]

Bessel function of first kind

Parameters:
  • n (int) – real order
  • z (float) – complex argument

Notes

https://de.mathworks.com/help/matlab/ref/besselj.html

ggf.matlab_funcs.gammaln(x)[source]

Logarithm of the absolute value of the Gamma function

Notes

https://de.mathworks.com/help/matlab/ref/gammaln.html

See also

scipy.special.gammaln()

ggf.matlab_funcs.legendre(n, x)[source]

Associated Legendre functions

Parameters:
  • n (int) – degree
  • x (ndarray of floats) – argument

Notes

x is treated always as a row vector

The statement legendre(2,0:0.1:0.2) returns the matrix

/ x = 0 x = 0.1 x = 0.2
m = 0 -0.5000 -0.4850 -0.4400
m = 1 0 -0.2985 -0.5879
m = 2 3.0000 2.9700 2.8800

Notes

https://de.mathworks.com/help/matlab/ref/legendre.html

ggf.matlab_funcs.lscov(A, B, w=None)[source]

Least-squares solution in presence of known covariance

\(A \cdot x = B\), that is, \(x\) minimizes \((B - A \cdot x)^T \cdot \text{diag}(w) \cdot (B - A \cdot x)\). The matrix \(w\) typically contains either counts or inverse variances.

Parameters:
  • A (matrix or 2d ndarray) – input matrix
  • B (vector or 1d ndarray) – input vector

Notes

https://de.mathworks.com/help/matlab/ref/lscov.html

ggf.matlab_funcs.quadl(func, a, b)[source]

Numerically evaluate integral with scipy QUADPACK quadrature

Parameters:
  • func (callable) – function to integrate
  • a (float) – interval start
  • b (float) – interval end

Notes

https://de.mathworks.com/help/matlab/ref/quadl.html

sci_funcs

Other scientific functions

ggf.sci_funcs.legendrePlm(m, l, x)[source]

stress

ggf.stress.get_stress(model, semi_major, semi_minor, object_index, medium_index, effective_fiber_distance=0.0001, mode_field_diameter=3e-06, power_per_fiber=0.6, wavelength=1.064e-06, n_points=100, verbose=False)[source]

Compute the optical stress profile in the optical stretcher

Parameters:
  • model (str) – Model to use, one of: boyde2009
  • semi_major (float) – Semi-major axis of an ellipse fit to the object perimeter [m]
  • semi_minor (float) – Semi-minor axis of an ellipse fit to the object perimeter [m]
  • object_index (float) – Refractive index of the object
  • medium_index (float) – Refractive index of the surrounding medium
  • effective_fiber_distance (float) – Effective distance between the two trapping fibers relative to the medium refractive index [m]. For an open setup, this is the physical distance between the fibers. For a closed setup (capillary), this distance takes into account the refractive indices and thicknesses of the glass capillary and index matching gel. For the closed setup, the convenience function ggf.fiber_distance_capillary() can be used.
  • mode_field_diameter (float) – The mode field diameter MFD of the fiber used [m]. Note that the MFD is dependent on the wavelength used. If the manufacturer did not provide a value for the MFD, the MFD can be approximated as 3*wavelenth for a single-mode fiber.
  • power_per_fiber (float) – The laser power coupled into each of the fibers [W]
  • wavelength (float) – The laser wavelength used for the trap [m]
  • n_points (int) – Number of points to compute.
  • verbose (int) – Increases verbosity
Returns:

  • theta (1d ndarray of length n_points) – Polar angles [rad]
  • sigma (1d ndarray of length n_points) – Radial stress profile along theta [Pa]

stress.boyde2009

stress.boyde2009.core

ggf.stress.boyde2009.core.boundary(costheta, a=1, epsilon=0.1, nu=0)[source]

Projected boundary of a prolate spheroid

Compute the boundary according to equation (4) in [BCG09] with the addition of the Poisson’s ratio of the object.

\[B(\theta) = a (1+\epsilon) \left[ (1+\epsilon)^2 - \epsilon (1+\nu) (2+\epsilon (1-\nu)) \cos^2 \theta \right]^{-1/2}\]

This boundary function was derived for a prolate spheroid under the assumption that the semi-major axis \(a\) and the semi-minor axes \(b=c\) are defined as

\[a = b \cdot \frac{1+ \epsilon}{1- \nu \epsilon}\]

The boundary function \(B(\theta)\) can be derived with the above relation using the equation for a prolate spheroid.

Parameters:
  • costheta (float or np.ndarray) – Cosine of polar coordinates \(\theta\) at which to compute the boundary.
  • a (float) – Equatorial radii of prolate spheroid (semi-minor axis).
  • epsilon (float) – Stretch ratio; defines size of semi-major axis: \(a = (1+\epsilon) b\). Note that this is not the eccentricity of the prolate spheroid.
  • nu (float) – Poisson’s ratio \(\nu\) of the material.
Returns:

B – Radial object boundary in dependence of theta \(B(\theta)\).

Return type:

1d ndarray

Notes

For \(\nu=0\), the above equation becomes equation (4) in [BCG09].

ggf.stress.boyde2009.core.get_hgc[source]

Load hypergeometric coefficients from hypergeomdata2.dat.

These coefficients were computed by Lars Boyde using Wolfram Mathematica.

ggf.stress.boyde2009.core.stress(object_index=1.41, medium_index=1.3465, poisson_ratio=0.45, semi_minor=2.8466e-06, stretch_ratio=0.1, wavelength=7.8e-07, beam_waist=3, power_left=0.6, power_right=0.6, dist=0.0001, n_points=100, theta_max=<Mock name='mock.pi' id='140388817295624'>, field_approx='davis', ret_legendre_decomp=False, verbose=False)[source]

Compute the stress acting on a prolate spheroid

The prolate spheroid has semi-major axis \(a\) and semi-minor axis \(b=c\).

Parameters:
  • object_index (float) – Refractive index of the spheroid
  • medium_index (float) – Refractive index of the surrounding medium
  • poisson_ratio (float) – Poisson’s ratio of the spheroid material
  • semi_minor (float) – Semi-minor axis (inner) radius of the stretched object \(b=c\).
  • stretch_ratio (float) – Measure of the deformation, defined as \((a - b) / b\)
  • wavelength (float) – Wavelenth of the gaussian beam [m]
  • beam_waist (float) – Beam waist radius of the gaussian beam [wavelengths]
  • power_left (float) – Laser power of the left beam [W]
  • power_right (float) – Laser power of the right beam [W]
  • dist (float) – Distance between beam waist and object center [m]
  • n_points (int) – Number of points to compute stresses for
  • theta_max (float) – Maximum angle to compute stressed for
  • field_approx (str) – TODO
  • ret_legendre_decomp (bool) – If True, return coefficients of decomposition of stress into Legendre polynomials
  • verbose (int) – Increase verbosity
Returns:

  • theta (1d ndarray) – Angles for which stresses are computed
  • sigma_rr (1d ndarray) – Radial stress corresponding to angles
  • coeff (1d ndarray) – If ret_legendre_decomp is True, return compositions of decomposition of stress into Legendre polynomials.

Notes

  • The angles theta are computed on a grid that does not include zero and theta_max.
  • This implementation was first presented in [BCG09].

stress.boyde2009.globgeomfact

Computation of the global geometric factor

ggf.stress.boyde2009.globgeomfact.coeff2ggf(coeff, poisson_ratio=0.45)[source]

Compute the global geometric factor from stress coefficients

The radial displacements of an elastic sphere can be expressed in terms of Legendre polynomials (see [Lure64] equation 6.2.9) whose coefficients are computed from the Legendre decomposition of the radial stress.

Notes

  • For a \(\sigma_0 \cos^n(\theta)\) stress profile, the GGF already includes the peak stress according to:

    \[\text{GGF} = \sigma_0 F_\text{G}.\]
  • This is a conversion of the Matlab script GGF.m to Python. The code solves a linear system of equations to determine all Legendre coefficients. The new implementation in ggf.legendre2ggf() uses the direct solution and thus should be preferred.