spam.mesh package#

Submodules#

spam.mesh.objects module#

The objects module offers functions that enables to generate and manipulates geometrical objects in order to represent various phases of materials.

>>> import spam.mesh
>>> spam.mesh.packSpheres()
>>> spam.mesh.packSpheresFromList()

Note

Objects array conventions:

  • Sphere: [radius, centerX, centerY, centerZ, phase]

  • Cylinder: [radius, centerX, centerY, centerZ, directionX, directionY, directionZ, phase]

Warning

This submodule will move to a different module in the near future.

spam.mesh.objects.packSpheres(totalVolumeFraction, rejectionLength, phases, origin=[0.0, 0.0, 0.0], lengths=[1.0, 1.0, 1.0], inside=True, domain='cube', vtk=None)[source]#

This function packs one or several sets (phase) of spheres of different radii and create the corresponding distance fields (one per set).

The packing algorithm is an iterative process based collective rearrangement.

Parameters:
  • totalVolumeFraction (float) – The total volume fraction of all the phases

  • rejectionLength (float) – The minimal distance between two sphere surfaces

  • phases ((nPhase x 3) array) –

    A 2D array containing the phases parameteres. A line corresponds to a phase and a column to a parameter:

    column 0: the minimal ray of the spheres of the phase
    column 1: the maximal ray of the spheres of the phase
    column 2: the relative volume fraction of the phase
    

  • inside (bool, default=True) – Defines whether or not the spheres have to be completly inside the domain or if they can intersect it (centres always remain inside).

  • lengths (array, default=[1, 1, 1]) – The size of the domain the spheres are packed into.

  • origin (array, default=[0, 0, 0]) – The origin of the domain.

  • domain (string, default='cube') –

    The domain type the spheres are packed into. Options are:

    • cube`: which corresponds to a cuboid. lengths is then

    the length of the cuboids.

    • cylinder: which corresponds to a cylinder of diameter lengths[0] and height lengths[2].

  • vtk (string, default=None) – Save vtk files of the spheres for each iterations of the packing algorithm under base name vtk.

Returns:

For each sphere: [radius, ox, oy, oz, phase]

Return type:

(nSpheres x 4) array

>>> import spam.mesh
>>> volumeFraction = 0.1
>>> rejectionLength = 1.0
>>> # phase 1: rmin = 5, rmax = 6.5, volume fraction = 0.6 (of total volume fraction)
>>> # phase 2: rmin = 6.5, rmax = 8, volume fraction = 0.4 (of total volume fraction)
>>> phases = [[5.0, 6.5, 0.6], [6.5, 8.0, 0.4]]
>>> # cylinder going from -5 to 135 in z
>>> # with a base radius 50 and center [50, 60
>>> domain = "cylinder"
>>> lengths = [100, 0, 140]
>>> origin = [0, 10, -5]
>>> # generate and pack spheres
>>> spheres = spam.mesh.packSpheres(volumeFraction, rejectionLength, phases, domain=domain, origin=origin, lengths=lengths, vtk="packing-1")

Note

The output of this function can directly be used by the function spam.mesh.projectObjects.

spam.mesh.objects.packObjectsFromList(objects, rejectionLength, origin=[0.0, 0.0, 0.0], lengths=[1.0, 1.0, 1.0], inside=True, domain='cube', vtk=None)[source]#

This function packs a set of predefine spheres.

The packing algorithm is an iterative process based collective rearrangement.

Parameters:
  • objects ((nSpheres x nPram) array) –

    The list of objects. Each line corresponds to:

    • for spheres: [radius, ox, oy, oz, phase]

  • rejectionLength (float) – The minimal distance between two spheres surfaces

  • inside (bool, default=True) – Defines whether or not the spheres have to be completly inside the domain or if they can intersect it (centres always remain inside).

  • lengths (array, default=[1.0, 1.0, 1.0]) – The size of the domain the spheres are packed into.

  • origin (array, default=[0.0, 0.0, 0.0]) – The origin of the domain.

  • domain (string, default='cube') –

    The domain type the spheres are packed into. Options are:

    • cube`: which corresponds to a cuboid. lengths is then

    the length of the cuboids.

    • cylinder: which corresponds to a cylinder of diameter lengths[0] and height lengths[2].

  • vtk (string, default=None) – Save vtk files of the spheres for each iterations of the packing algorithm under base name vtk.

Returns:

For each sphere: [radius, ox, oy, oz, phase]

Return type:

(nSpheres x 4) array

Note

The output of this function can directly be used by the function spam.mesh.projectObjects.

spam.mesh.projection module#

The spam.mesh.projection module offers a two functions that enables the construction of three dimensional unstructured meshes (four noded tetrahedra) that embbed heterogeneous data (namely, subvolumes and interface orientation) based on eihter:

  • the distance field of an segmented images (see spaM.filters.distanceField),

  • ideal geometrical objects (see spam.filters.objects).

>>> import spam.mesh
>>> spam.mesh.projectField(mesh, fields)
>>> spam.mesh.projectObjects(mesh, objects)

Note

Objects array conventions:

  • Sphere: [radius, centerX, centerY, centerZ, phase]

  • Cylinder: [radius, centerX, centerY, centerZ, directionX, directionY, directionZ, phase]

Distance field and orientation convention:

  • The distance fields inside connected components have positive values

  • The normal vectors are pointing outside (to negative values)

  • The subvolumes returned correspond to the outside (negative part)

spam.mesh.projection.projectObjects(mesh, objects, analyticalOrientation=True, cutoff=1e-06, writeConnectivity=None, vtkMesh=False)[source]#

This function project a set of objects onto an unstructured mesh. Each objects can be attributed a phase.

Parameters:
  • mesh (dict) –

    Dictionary containing points coordinates and mesh connectivity.

    mesh = {
        # numpy array of size number of nodes x 3
        "points": points,
        # numpy array of size number of elements x 4
        "cells": connectivity,
    }
    

  • objects (2D array) –

    The list of objects. Each line corresponds to an object encoded as follow:

    - spheres: radius, oZ, oY, oX, phase
    - cylinders: radius, oZ, oY, oX, nZ, nY, nX, phase
    

  • analyticalOrientation (bool, default=True) –

    Change how the interface’s normals are computed:

    • True: as the direction from the centroid of the tetrahedron to the closest highest value of the object distance field (ie, the center of a sphere, the axis of a cylinder).

    • False: as the normals of the facets which points lie on the object surface.

  • cutoff (float, default=1e-6) – Volume ratio below wich elements with interfaces are ignored and considered being part of a single phase.

  • writeConnectivity (string, default=None) –

    When not None, it writes a text file called writeConnectivity the list of node and the list of elements

    which format is:

    COORdinates ! number of nodes
    nodeId, 0, x, y, z
    ...
    
    ELEMents ! number of elemens
    elemId, 0, elemType, n1, n2, n3, n4, subVol(-), normalX,
    normalY, normalZ
    ...
    

    where:

    • n1, n2, n3, n4 is the connectivity (node numbers).

    • subVol(-) is the sub volume of the terahedron inside the inclusion.

    • normalX, normalY, normalZ are to components of the interface normal vector.

    • elemType is the type of element. Their meaning depends on the their phase.

    Correspondance can be found in the function output after the key word MATE like:

    <projmorpho::set_materials
    .        field 1
    .        .       MATE,1: background
    .        .       MATE,2: phase 1
    .        .       MATE,3: interface phase 1 with
    background
    .        field 2
    .        .       MATE,1: background
    .        .       MATE,4: phase 2
    .        .       MATE,5: interface phase 2 with
    background
    >
    

    Sub volumes and interface vector are only relevant for interfaces. Default value is 0 and [1, 0, 0].

  • vtkMesh (bool, default=False) – Writes the VTK of interpolated fields and materials. The files take the same name as writeConnectivity.

Returns:

For each element: elemType, subVol(-), normalX, normalY, normalZ (see outputFile for details).

Return type:

(nElem x 5) array

Note

Only four noded tetrahedra meshes are supported.

>>> import spam.mesh
>>> # unstructured mesh
>>> points, cells = spam.mesh.createCuboid([1, 1, 1], 0.1)
>>> mesh = {"points": points, "cells": cells}
>>> # one centered sphere (0.5, 0.5, 0.5) of radius 0.2 (phase 1)
>>> sphere = [[0.2, 0.8, 0.1, 0.1, 1]]
>>> # projection
>>> materials = spam.mesh.projectObjects(mesh, sphere, writeConnectivity="mysphere", vtkMesh=True)
spam.mesh.projection.projectField(mesh, fields, thresholds=[0.0], cutoff=1e-06, writeConnectivity=None, vtkMesh=False)[source]#

This function project a set of distance fields onto an unstructured mesh. Each distance field corresponds to a phase and the interface between the two phases is set by the thresholds.

Parameters:
  • mesh (dict) –

    Dictionary containing points coordinates and mesh connectivity.

    mesh = {
        # numpy array of size number of nodes x 3
        "points": points,
        # numpy array of size number of elements x 4
        "cells": connectivity,
    }
    

  • fields (list of dicts) –

    The fields should be continuous (like a distance field) for a better projection. They are discretised over a regular mesh (lexicographical order) and eahc one corresponds to a phase. The dictionary containing the fields data is defined as follow

    fields = {
        # coordinates of the origin of the field (3 x 1)
        "origin": origin,
        # lengths of fields domain of definition  (3 x 1)
        "lengths": lengths,
        #  list of fields values (list of 3D arrays)
        "values": [phase1, phase2],
    }
    

  • thresholds (list of floats, default=[0]) – The list of thresholds.

  • cutoff (float, default=1e-6) – Volume ratio below wich elements with interfaces are ignored and considered being part of a single phase.

  • writeConnectivity (string, default=None) –

    When not None, it writes a text file called writeConnectivity the list of node and the list of elements

    which format is:

    COORdinates ! number of nodes
    nodeId, 0, x, y, z
    ...
    
    ELEMents ! number of elemens
    elemId, 0, elemType, n1, n2, n3, n4, subVol(-), normalX,
    normalY, normalZ
    ...
    

    where:

    • n1, n2, n3, n4 is the connectivity (node numbers).

    • subVol(-) is the sub volume of the terahedron inside the inclusion.

    • normalX, normalY, normalZ are to components of the interface normal vector.

    • elemType is the type of element. Their meaning depends on the their phase.

    Correspondance can be found in the function output after the key word MATE like:

    <projmorpho::set_materials
    .        field 1
    .        .       MATE,1: background
    .        .       MATE,2: phase 1
    .        .       MATE,3: interface phase 1 with
    background
    .        field 2
    .        .       MATE,1: background
    .        .       MATE,4: phase 2
    .        .       MATE,5: interface phase 2 with
    background
    >
    

    Sub volumes and interface vector are only relevant for interfaces. Default value is 0 and [1, 0, 0].

  • vtkMesh (bool, default=False) – Writes the VTK of interpolated fields and materials. The files take the same name as writeConnectivity.

Returns:

For each element: elemType, subVol(-), normalX, normalY, normalZ (see outputFile for details).

Return type:

(nElem x 5) array

Note

Only four noded tetrahedra meshes are supported.

>>> import spam.mesh
>>> import spam.kalisphera
>>> # unstructured mesh
>>> points, cells = spam.mesh.createCuboid([1, 1, 1], 0.1)
>>> mesh = {"points": points, "cells": cells}
>>> # create an image of a sphere and its distance field
>>> sphereBinary = numpy.zeros([101] * 3, dtype=float)
>>> spam.kalisphera.makeSphere(sphereBinary, [50.5] * 3, 20)
>>> sphereField = spam.mesh.distanceField(one_sphere_binary.astype(bool).astype(int))
>>> # format field object
>>> fields = {
>>>     # coordinates of the origin of the field (3 x 1)
>>>     "origin": [0] * 3,
>>>     # lengths of fields domain of definition  (3 x 1)
>>>     "lengths": [1] * 3,
>>>     # list of fields
>>>     "values": [sphereField],
>>> }
>>> # projection
>>> materials = spam.mesh.projectField(mesh, fields, writeConnectivity="mysphere", vtkMesh=True)

spam.mesh.structured module#

This module offers a set of tools in order to manipulate structured meshes.

>>> # import module
>>> import spam.mesh

The strucutred VTK files used to save the data have the form:

# vtk DataFile Version 2.0
VTK file from spam: spam.vtk
ASCII

DATASET STRUCTURED_POINTS
DIMENSIONS nx ny nz
ASPECT_RATIO lx ly lz
ORIGIN ox oy oz

POINT_DATA nx x ny x nz

SCALARS myNodalField1 float
LOOKUP_TABLE default
    nodalValue_1
    nodalValue_2
    nodalValue_3
    ...

VECTORS myNodalField2 float
LOOKUP_TABLE default
    nodalValue_1_X nodalValue_1_Y nodalValue_1_Z
    nodalValue_2_X nodalValue_2_Y nodalValue_2_Z
    nodalValue_3_X nodalValue_3_Y nodalValue_3_Z
    ...

CELL_DATA (nx-1) x (ny-1) x (nz-1)

SCALARS myCellField1 float
LOOKUP_TABLE default
    cellValue_1
    cellValue_2
    cellValue_3
    ...

where nx, ny and nz are the number of nodes in each axis, lx, ly, lz, the mesh length in each axis and ox, oy, oz the spatial postiion of the origin.

spam.mesh.structured.createCylindricalMask(shape, radius, voxSize=1.0, centre=None)[source]#

Create a image mask of a cylinder in the z direction.

Parameters:
  • shape (array, int) – The shape of the array the where the cylinder is saved

  • radius (float) – The radius of the cylinder

  • voxSize (float (default=1.0)) – The physical size of a voxel

  • centre (array of floats of size 2, (default None)) – The center [y,x] of the axis of rotation of the cylinder. If None it is taken to be the centre of the array.

Returns:

cyl – The cylinder

Return type:

array, bool

spam.mesh.structured.structuringElement(radius=1, order=2, margin=0, dim=3)[source]#

This function construct a structural element.

Parameters:
  • radius (int, default=1) –

    The radius of the structural element

    radius = 1 gives 3x3x3 arrays
    radius = 2 gives 5x5x5 arrays
    ...
    radius = n gives (2n+1)x(2n+1)x(2n+1) arrays
    

  • order (int, default=2) –

    Defines the shape of the structuring element by setting the order of the norm used to compute the distance between the centre and the border.

    A representation for the slices of a 5x5x5 element (size=2) from the center to on corner (1/8 of the cube)

    order=numpy.inf: the cube
    1 1 1    1 1 1    1 1 1
    1 1 1    1 1 1    1 1 1
    1 1 1    1 1 1    1 1 1
    
    order=2: the sphere
    1 0 0    0 0 0    0 0 0
    1 1 0    1 1 0    0 0 0
    1 1 1    1 1 0    1 0 0
    
    order=1: the diamond
    1 0 0    0 0 0    0 0 0
    1 1 0    1 0 0    0 0 0
    1 1 1    1 1 0    1 0 0
    

  • margin (int, default=0) – Gives a 0 valued margin of size margin.

  • dim (int, default=3) – Spatial dimension (2 or 3).

Returns:

The structural element

Return type:

array

spam.mesh.structured.createLexicoCoordinates(lengths, nNodes, origin=(0, 0, 0))[source]#

Create a list of coordinates following the lexicographical order.

Parameters:
  • lengths (array of floats) – The length of the cuboids in every directions.

  • nNodes (array of int) – The number of nodes of the mesh in every directions.

  • origin (array of floats) – The coordinates of the origin of the mesh.

Returns:

The list of coordinates. shape=(nx*ny*nz, 3)

Return type:

array

spam.mesh.tetrahedra module#

spam.mesh.tetrahedra.tetCentroid(points)[source]#

Compute coordinates of the centroid of a tetrahedron

Parameters:

points (4x3 array) – Array of the 3D coordinates of the 4 points

Returns:

the coordinates of the centroid

Return type:

3x1 array

spam.mesh.tetrahedra.tetVolume(points)[source]#

Compute the volume of a tetrahedron

Parameters:

points (4x3 array) – Array of the 3D coordinates of the 4 points

Returns:

the volume of the tetrahedron

Return type:

float

spam.mesh.tetrahedra.shapeFunctions(points)[source]#

This function computes the shape coefficients matrux from the coordinates of the 4 nodes of a linear tetrahedron (see 4.11 in Zienkiewicz).

coefficients_matrix = [
        a1 b1 c1 d1
        a2 b2 c2 d2
        a3 b3 c3 d3
        a4 b4 c4 d4
    ]
with
N1 = a1 + b1x + c1y + d1z
N2 = a2 + b2x + c2y + d2z
N3 = a3 + b3x + c3y + d3z
N4 = a4 + b4x + c4y + d4z
Parameters:

points (4x3 array) – Array of the 3D coordinates of the 4 points

Returns:

  • volume (float) – The volume of the tetrahedron

  • coefficients_matrix (4x4 array) – The coefficients 4 coefficients of the 4 shape functions

Note

Pure python function.

spam.mesh.tetrahedra.elementaryStiffnessMatrix(points, young, poisson)[source]#

This function computes elementary stiffness matrix from the coordinates of the 4 nodes of a linear tetrahedron.

D = [something with E and mu] (6x6)
B = [B1, B2, B3, B4]          (4 times 6x3 -> 6x12)
ke = V.BT.D.B                 (12x12)
Parameters:
  • points (4x3 array) – Array of the 3D coordinates of the 4 points

  • young (float) – Young modulus

  • poisson (float) – Poisson ratio

Returns:

stiffness_matrix – The stiffness matrix of the tetrahedron

Return type:

12x12 array

Note

Pure python function.

spam.mesh.tetrahedra.globalStiffnessMatrix(points, connectivity, young, poisson)[source]#

This function assembles elementary stiffness matrices to computes an elastic (3 dof) global stiffness matrix for a 4 noded tetrahedra mesh

Notations#

neq: number of global equations nel: number of elements npo: number of points neq = 3 x npo

param points:

Array of the 3D coordinates of the all nodes points

type points:

npo x 3 array

param connectivity:

Connectivity matrix

type connectivity:

nel x 4 array

param young:

Young modulus

type young:

float

param poisson:

Poisson ratio

type poisson:

float

returns:

stiffness_matrix – The stiffness matrix of the tetrahedron

rtype:

neq x neq array

Note

Pure python function.

spam.mesh.unstructured module#

This module offers a set of tools for unstructured 3D meshes made of tetrahedra.

>>> # import module
>>> import spam.mesh
>>> spam.mesh.createCuboid()
spam.mesh.unstructured.gmshToSpam(elementsByType, nodesByType)[source]#

Helper function Converts gmsh mesh data to SPAM format by: 1. reordering by z y x 2. keeping only tetrahedra

Parameters:
  • elementsByType (array) – Should be the output of gmsh.model.mesh.getElementsByType(4)

  • nodesByType (array) – Should be the output of gmsh.model.mesh.getNodesByType(4)

Returns:

  • points (2D numpy array) – The coordinates of the mesh nodes (zyx) Each line is [zPos, yPos, xPos]

  • connectivity (2D numpy array) – The connectivity matrix of the tetrahedra elements Each line is [node1, node2, node3, node4]

spam.mesh.unstructured.getMeshCharacteristicLength(points, connectivity)[source]#

Computes the average distance between two nodes of the edges of each elements.

Parameters:
  • points (Nx3 array) – List of coordinates of the mesh nodes.

  • connectivity (Mx4 array) – Connectivity matrix of the mesh.

Returns:

float

Return type:

the characteristic length

spam.mesh.unstructured.createCuboid(lengths, lc, origin=[0.0, 0.0, 0.0], periodicity=False, verbosity=1, gmshFile=None, vtkFile=None, binary=False, skipOutput=False)[source]#

Creates an unstructured mesh of tetrahedra inside a cuboid.

Parameters:
  • lengths (1D array) – The lengths of the cuboid in the three directions The axis order is zyx

  • origin (1D array) – The origin of the cuboid (zyx)

  • lc (float) – characteristic length of the elements of the mesh (i.e., the average distance between two nodes)

  • periodicity (bool, optional) – if periodicity is True, the generated cube will have a periodicity of mesh on surfaces Default = False

  • gmshFile (string, optional) – If not None, save the gmsh file with name gmshFile and suffix .msh Default = None

  • vtkFile (string, optional) – If not None, save the vtk file with name vtkFile and suffix .vtk Defaut = None

  • binary (bool, optional) – Save files in binary when possible Default = False

  • skipOutput (bool, optional) – Returns None to save time (only write the files) Default = False

  • verbosity (int, optional) – Level of information printed on the terminal and the message console. 0: silent except for fatal errors 1: +errors 2: +warnings 3: +direct 4: +information 5: +status 99: +debug Default = 1

Returns:

  • points (2D numpy array) – The coordinates of the mesh nodes (zyx) Each line is [zPos, yPos, xPos]

  • connectivity (2D numpy array) – The connectivity matrix of the tetrahedra elements Each line is [node1, node2, node3, node4]

Example

>>> points, connectivity = spam.mesh.createCuboid((1.0,1.5,2.0), 0.5)
create a mesh in a cuboid of size 1,1.5,2 with a characteristic length of 0.5
spam.mesh.unstructured.createCylinder(centre, radius, height, lc, zOrigin=0.0, verbosity=0, gmshFile=None, vtkFile=None, binary=False, skipOutput=False)[source]#

Creates an unstructured mesh of tetrahedra inside a cylinder. The height of the cylinder is along the z axis.

Parameters:
  • centre (1D array) – The two coordinates of the centre of the base disk (yx).

  • radius (float) – The radius of the base disk.

  • height (float) – The height of the cylinder.

  • lc (float) – characteristic length of the elements of the mesh (i.e., the average distance between two nodes)

  • zOrigin (float, default = 0.0) – Translate the points coordinates by zOrigin in the z direction.

  • gmshFile (string, optional) – If not None, save the gmsh file with name gmshFile and suffix .msh Default = None

  • vtkFile (string, optional) – If not None, save the vtk file with name vtkFile and suffix .vtk Defaut = None

  • binary (bool, optional) – Save files in binary when possible Default = False

  • skipOutput (bool, optional) – Returns None to save time (only write the files) Default = False

  • verbosity (int, optional) – Level of information printed on the terminal and the message console. 0: silent except for fatal errors 1: +errors 2: +warnings 3: +direct 4: +information 5: +status 99: +debug Default = 1

Returns:

  • points (2D numpy array) – The coordinates of the mesh nodes (zyx) Each line is [zPos, yPos, xPos]

  • connectivity (2D numpy array) – The connectivity matrix of the tetrahedra elements Each line is [node1, node2, node3, node4]

Example

>>> points, connectivity =  spam.mesh.createCylinder((0.0,0.0), 0.5, 2.0, 0.5)
create a mesh in a cylinder of centre 0,0,0 radius, 0.5 and height 2.0 with a characteristic length of 0.5
spam.mesh.unstructured.triangulate(points, alpha=None, weights=None)[source]#

Takes a series of points and optionally their weights and returns a tetrahedral connectivity matrix.

If a completely regular grid is passed, there will be trouble, add some tiny noise.

This function uses CGAL’s Regular Triangulation in the background (with all weights=1 if they are not passed). Weights are passed to CGAL’s Regular Triangulation directly and so should be squared if they are radii of particles.

Users can optionally pass an alpha value to the function with the goal of carving flat boundary tetrahedra away from the final mesh. Typical use of the alpha shape could be the removal of flat (almost 0 volume) boundary tetrahedra from concave/convex boundaries. Another example might be the removal of tetrahedra from an internal void that exists within the domain. In all cases, the user can imagine the alpha tool as a spherical “scoop” that will remove any tetrahedra that it is capable of entering. It follows that flat tetrahedra have a very large open face which are easily entered by large alpha “scoops”. Thus, the user should imagine using a very large alpha value (try starting with 5*domain size) to avoid letting the alpha “scoop” enter well formed tetrahedra.

Consider the following CGAL analogy: The full mesh is an ice cream sundae and the vertices are the chocolate chips. The value of alpha is the squared radius of the icecream scoop (following the mesh coordinate units) that will go in and remove any tetrahedra that it can pass into. Positive alpha is user defined, negative alpha allows CGAL to automatically select a continuous solid (not recommended). For more information visit: https://doc.cgal.org/latest/Alpha_shapes_3/index.html

Parameters:
  • points (Nx3 2D numpy array of floats) – N points to triangulate (Z, Y, X)

  • weights (numpy vector of N floats) – list of weights associated with points. Default = None (which means all weights = 1).

  • alpha (float) – size of the alpha shape used for carving the mesh. Zero is no carving. Negative is CGAL-autocarve which gives an overcut mesh. positive is a user-selected size, try 5*domain size. Default = 0 (no carving).

Returns:

connectivity – delaunay triangulation with each row containing point numbers

Return type:

Mx4 numpy array of unsigned ints

Note

Function contributed by Robert Caulk (Laboratoire 3SR, Grenoble)

spam.mesh.unstructured.projectTetFieldToGrains(points, connectivity, tetField)[source]#

Projects/coarse-grains any field defined on tetrahedra onto grains by volume-averaging over all tetrahedra for which a given grain is a node. This can be useful for smoothing out a noisy strain field and will not affect the overall agreement between the average of strains and the macroscopically observed strains (R.C. Hurley has verified this in a 2017 paper).

Parameters:
  • points (m x 3 numpy array) – M Particles’ coordinates (in deformed configuration for strain field)

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

  • tetField (n x 3 x 3 numpy array) – Any field defined on tetrahedra (e.g., Bagi strains from bagiStrain).

Returns:

grainField – array containing (3x3) arrays of strain

Return type:

m x 3 x 3

Example

grainStrain = spam.mesh.projectBagiStrainToGrains(connectivity,bagiStrain[0],points)

Returns strains for each grain.

Notes

Function contributed by Ryan Hurley (Johns Hopkins University)

spam.mesh.unstructured.BCFieldFromDVCField(points, dvcField, mask=None, pixelSize=1.0, meshType='cube', centre=[0, 0], radius=1.0, topBottom=False, tol=1e-06, neighbours=4)[source]#

This function imposes boundary conditions coming from a DVC result to the nodes of an unstructured FE mesh.

Parameters:
  • points (2D numpy array) – Array of n node positions of the unstructured mesh Each line is [nodeNumber, z, y, x]

  • dvcField (2D numpy array) – Array of m points of the dvc field Each line is [zPos, yPos, xPos, zDisp, yDisp, xDisp]

  • mask (2D numpy array, optional) – Boolean array of m points of the dvc field Points with 0 will be ignored for the field interpolation Default = None (i.e. interpolate based on all of the dvc points)

  • pixelSize (float, optional) – physical size of a pixel (i.e. 1mm/px) Default = 1.0

  • meshType (string, optional) – For the moment, valid inputs are cube and cylinder The axis of a cylinder is considered to be z Note that if a cylindrical mesh is passed, centre and radius are highly recommended Default = cube

  • centre (float, optional) – The centre of the cylinder [y, x] in physical units (i.e. mm) Default = [0, 0]

  • radius (float, optional) – The radius of the cylinder in physical units (i.e. mm) Default = 1.0

  • topBottom (bool, optional) – If boundary conditions are passed only for the top (i.e. z=zmax) and bottom (i.e. z=zmin) surfaces of the mesh Default = False

  • tol (float, optional) – Numerical tolerance for floats equality Default = 1e-6

  • neighbours (int, , optional) – Neighbours for field interpolation Default = 4

Returns:

bc – Boundary node displacements Each line is [nodeNumber, zPos, yPos, xPos, zDisp, yDisp, xDisp]

Return type:

2D numpy array

Warning

  1. All coordinates and displacement arrays are z, y, x

  2. The axis of a cylinder is considered to be z

spam.mesh.unstructured.tetVolumes(points, connectivity)[source]#

This function computes volumes of the tetrahedra passed with a connectivity matrix. Using algorithm in https://en.wikipedia.org/wiki/Tetrahedron#Volume

Parameters:
  • points (Nx3 array) – Array of N coordinates of the points

  • connectivity (Mx4 array) – Array of M none numbers that are connected as tetrahedra (e.g., the output from triangulate)

Returns:

volumes – Volumes of tetrahedra

Return type:

vector of length M

Note

Pure python function.

spam.mesh.unstructured.rankPoints(points, neighbourRadius=20, verbose=True)[source]#

This function ranks an array of points around the top point

Parameters:
  • points (numpy array N x 3) – Coordinates (zyx) of the points

  • neighbourRadius (float, optional) – Distance from the current point to include neighbours If no neighbour is found, then the closest point is taken Default: 20

Returns:

  • pointsRanked (numpy array N x 3) – Coordinates (zyx) of the ranked points

  • rowNumbers (1D numpy array N of ints) – Reorganised row numbers from input

Module contents#