spam.deformation package#

Submodules#

spam.deformation.deformationField module#

spam.deformation.deformationField.FfieldRegularQ8(displacementField, nodeSpacing, nProcesses=32, verbose=False)[source]#

This function computes the transformation gradient field F from a given displacement field. Please note: the transformation gradient tensor: F = I + du/dx.

This function computes du/dx in the centre of an 8-node cell (Q8 in Finite Elements terminology) using order one (linear) shape functions.

Parameters:
  • displacementField (4D array of floats) – The vector field to compute the derivatives. #Its shape is (nz, ny, nx, 3)

  • nodeSpacing (3-component list of floats) – Length between two nodes in every direction (i.e., size of a cell)

  • nProcesses (integer, optional) – Number of processes for multiprocessing Default = number of CPUs in the system

  • verbose (boolean, optional) – Print progress? Default = True

Returns:

F – The field of the transformation gradient tensors

Return type:

(nz-1) x (ny-1) x (nx-1) x 3x3 array of n cells

spam.deformation.deformationField.FfieldRegularGeers(displacementField, nodeSpacing, neighbourRadius=1.5, nProcesses=32, verbose=False)[source]#

This function computes the transformation gradient field F from a given displacement field. Please note: the transformation gradient tensor: F = I + du/dx.

This function computes du/dx as a weighted function of neighbouring points. Here is implemented the linear model proposed in: “Computing strain fields from discrete displacement fields in 2D-solids”, Geers et al., 1996

Parameters:
  • displacementField (4D array of floats) – The vector field to compute the derivatives. Its shape is (nz, ny, nx, 3).

  • nodeSpacing (3-component list of floats) – Length between two nodes in every direction (i.e., size of a cell)

  • neighbourRadius (float, optional) – Distance in nodeSpacings to include neighbours in the strain calcuation. Default = 1.5*nodeSpacing which will give radius = 1.5*min(nodeSpacing)

  • mask (bool, optional) – Avoid non-correlated NaN points in the displacement field? Default = True

  • nProcesses (integer, optional) – Number of processes for multiprocessing Default = number of CPUs in the system

  • verbose (boolean, optional) – Print progress? Default = True

Returns:

Ffield – The field of the transformation gradient tensors

Return type:

nz x ny x nx x 3x3 array of n cells

Note

Taken from the implementation in “TomoWarp2: A local digital volume correlation code”, Tudisco et al., 2017

spam.deformation.deformationField.FfieldBagi(points, connectivity, displacements, nProcesses=32, verbose=False)[source]#

Calculates transformation gradient function using Bagi’s 1996 paper, especially equation 3 on page 174. Required inputs are connectivity matrix for tetrahedra (for example from spam.mesh.triangulate) and nodal positions in reference and deformed configurations.

Parameters:
  • points (m x 3 numpy array) – M Particles’ points in reference configuration

  • connectivity (n x 4 numpy array) – Delaunay triangulation connectivity generated by spam.mesh.triangulate for example

  • displacements (m x 3 numpy array) – M Particles’ displacement

  • nProcesses (integer, optional) – Number of processes for multiprocessing Default = number of CPUs in the system

  • verbose (boolean, optional) – Print progress? Default = True

Returns:

Ffield – The field of the transformation gradient tensors

Return type:

nx3x3 array of n cells

spam.deformation.deformationField.decomposeFfield(Ffield, components, twoD=False, nProcesses=32, verbose=False)[source]#

This function takes in an F field (from either FfieldRegularQ8, FfieldRegularGeers, FfieldBagi) and returns fields of desired transformation components.

Parameters:
  • Ffield (multidimensional x 3 x 3 numpy array of floats) – Spatial field of Fs

  • components (list of strings) – These indicate the desired components consistent with spam.deformation.decomposeF or decomposePhi

  • twoD (bool, optional) – Is the Ffield in 2D? This changes the strain calculation. Default = False

  • nProcesses (integer, optional) – Number of processes for multiprocessing Default = number of CPUs in the system

  • verbose (boolean, optional) – Print progress? Default = True

Returns:

  • Dictionary containing appropriately reshaped fields of the transformation components requested.

  • Keys – vol, dev, volss, devss are scalars r and z are 3-component vectors e and U are 3x3 tensors

spam.deformation.deformationField.decomposePhiField(PhiField, components, twoD=False, nProcesses=32, verbose=False)[source]#

This function takes in a Phi field (from readCorrelationTSV?) and returns fields of desired transformation components.

Parameters:
  • PhiField (multidimensional x 4 x 4 numpy array of floats) – Spatial field of Phis

  • components (list of strings) – These indicate the desired components consistent with spam.deformation.decomposePhi

  • twoD (bool, optional) – Is the PhiField in 2D? This changes the strain calculation. Default = False

  • nProcesses (integer, optional) – Number of processes for multiprocessing Default = number of CPUs in the system

  • verbose (boolean, optional) – Print progress? Default = True

Returns:

  • Dictionary containing appropriately reshaped fields of the transformation components requested.

  • Keys – vol, dev, volss, devss are scalars t, r, and z are 3-component vectors e and U are 3x3 tensors

spam.deformation.deformationFunction module#

spam.deformation.deformationFunction.computePhi(transformation, PhiCentre=[0.0, 0.0, 0.0], PhiPoint=[0.0, 0.0, 0.0])[source]#

Builds “Phi”, a 4x4 deformation function from a dictionary of transformation parameters (translation, rotation, zoom, shear). Phi can be used to deform coordinates as follows: $$ Phi.x = x’$$

Parameters:
  • transformation (dictionary of 3x1 arrays) –

    Input to computeTransformationOperator is a “transformation” dictionary where all items are optional.

    Keys

    t = translation (z,y,x). Note: (0, 0, 0) does nothing r = rotation in “rotation vector” format. Note: (0, 0, 0) does nothing z = zoom. Note: (1, 1, 1) does nothing s = “shear”. Note: (0, 0, 0) does nothing U = Right stretch tensor

  • PhiCentre (3x1 array, optional) – Point where Phi is centered (centre of rotation)

  • PhiPoint (3x1 array, optional) – Point where Phi is going to be applied

Returns:

Phi – Phi, deformation function

Return type:

4x4 array of floats

Note

Useful reference: Chapter 2 – Rigid Body Registration – John Ashburner & Karl J. Friston, although we use a symmetric shear. Here we follow the common choice of applying components to Phi in the following order: U or (z or s), r, t

spam.deformation.deformationFunction.decomposeF(F, twoD=False, verbose=False)[source]#

Get components out of a transformation gradient tensor F

Parameters:
  • F (3x3 array) – The transformation gradient tensor F (I + du/dx)

  • twoD (bool, optional) – is the F defined in 2D? This applies corrections to the strain outputs Default = False

  • verbose (boolean, optional) – Print errors? Default = True

Returns:

transformation

  • r = 3x1 numpy array: Rotation in “rotation vector” format

  • z = 3x1 numpy array: Zoom in “zoom vector” format (z, y, x)

  • U = 3x3 numpy array: Right stretch tensor, U - numpy.eye(3) is the strain tensor in large strains

  • e = 3x3 numpy array: Strain tensor in small strains (symmetric)

  • vol = float: First invariant of the strain tensor (“Volumetric Strain”), det(F)-1

  • dev = float: Second invariant of the strain tensor (“Deviatoric Strain”)

  • volss = float: First invariant of the strain tensor (“Volumetric Strain”) in small strains, trace(e)/3

  • devss = float: Second invariant of the strain tensor (“Deviatoric Strain”) in small strains

Return type:

dictionary of arrays

spam.deformation.deformationFunction.decomposePhi(Phi, PhiCentre=[0.0, 0.0, 0.0], PhiPoint=[0.0, 0.0, 0.0], twoD=False, verbose=False)[source]#

Get components out of a linear deformation function “Phi”

Parameters:
  • Phi (4x4 array) – The deformation function operator “Phi”

  • PhiCentre (3x1 array, optional) – Point where Phi was calculated

  • PhiPoint (3x1 array, optional) – Point where Phi is going to be applied

  • twoD (bool, optional) – is the Phi defined in 2D? This applies corrections to the strain outputs Default = False

  • verbose (boolean, optional) – Print errors? Default = True

Returns:

transformation

  • t = 3x1 numpy array: Translation vector (z, y, x)

  • r = 3x1 numpy array: Rotation in “rotation vector” format

  • z = 3x1 numpy array: Zoom in “zoom vector” format (z, y, x)

  • U = 3x3 numpy array: Right stretch tensor, U - numpy.eye(3) is the strain tensor in large strains

  • e = 3x3 numpy array: Strain tensor in small strains (symmetric)

  • vol = float: First invariant of the strain tensor (“Volumetric Strain”), det(F)-1

  • dev = float: Second invariant of the strain tensor (“Deviatoric Strain”)

  • volss = float: First invariant of the strain tensor (“Volumetric Strain”) in small strains, trace(e)/3

  • devss = float: Second invariant of the strain tensor (“Deviatoric Strain”) in small strains

Return type:

dictionary of arrays

spam.deformation.deformationFunction.computeRigidPhi(Phi)[source]#

Returns only the rigid part of the passed Phi

Parameters:

Phi (4x4 numpy array of floats) – Phi, deformation function

Returns:

PhiRigid – Phi with only the rotation and translation parts on inputted Phi

Return type:

4x4 numpy array of floats

Module contents#