yade.pack module

Creating packings and filling volumes defined by boundary representation or constructive solid geometry.

For examples, see

yade.pack.SpherePack_toSimulation(self, rot=Matrix3(1, 0, 0, 0, 1, 0, 0, 0, 1), **kw)

Append spheres directly to the simulation. In addition calling O.bodies.append, this method also appropriately sets periodic cell information of the simulation.

>>> from yade import pack; from math import *
>>> sp=pack.SpherePack()

Create random periodic packing with 20 spheres:

>>> sp.makeCloud((0,0,0),(5,5,5),rMean=.5,rRelFuzz=.5,periodic=True,num=20)
20

Virgin simulation is aperiodic:

>>> O.reset()
>>> O.periodic
False

Add generated packing to the simulation, rotated by 45° along +z

>>> sp.toSimulation(rot=Quaternion((0,0,1),pi/4),color=(0,0,1))
[0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19]

Periodic properties are transferred to the simulation correctly:

>>> O.periodic
True
>>> O.cell.refSize
Vector3(5,5,5)
>>> O.cell.trsf
Matrix3(0.707107,-0.707107,0, 0.707107,0.707107,0, 0,0,1)
Parameters:
  • rot (Quaternion/Matrix3) – rotation of the packing, which will be applied on spheres and will be used to set Cell.trsf as well.
  • **kw

    passed to utils.sphere

Returns:

list of body ids added (like O.bodies.append)

yade.pack.filterSpherePack(predicate, spherePack, **kw)

Using given SpherePack instance, return spheres the satisfy predicate. The packing will be recentered to match the predicate and warning is given if the predicate is larger than the packing.

yade.pack.gtsSurface2Facets(surf, **kw)

Construct facets from given GTS surface. **kw is passed to utils.facet.

yade.pack.gtsSurfaceBestFitOBB(surf)

Return (Vector3 center, Vector3 halfSize, Quaternion orientation) describing best-fit oriented bounding box (OBB) for the given surface. See cloudBestFitOBB for details.

class yade.pack.inGtsSurface_py(inherits Predicate)

This class was re-implemented in c++, but should stay here to serve as reference for implementing Predicates in pure python code. C++ allows us to play dirty tricks in GTS which are not accessible through pygts itself; the performance penalty of pygts comes from fact that if constructs and destructs bb tree for the surface at every invocation of gts.Point().is_inside(). That is cached in the c++ code, provided that the surface is not manipulated with during lifetime of the object (user’s responsibility).

Predicate for GTS surfaces. Constructed using an already existing surfaces, which must be closed.

import gts surf=gts.read(open(‘horse.gts’)) inGtsSurface(surf)

Note

Padding is optionally supported by testing 6 points along the axes in the pad distance. This must be enabled in the ctor by saying doSlowPad=True. If it is not enabled and pad is not zero, warning is issued.

aabb()
class yade.pack.inSpace(inherits Predicate)

Predicate returning True for any points, with infinite bounding box.

aabb()
center()
dim()
yade.pack.randomDensePack(predicate, radius, material=-1, dim=None, cropLayers=0, rRelFuzz=0.0, spheresInCell=0, memoizeDb=None, useOBB=True, memoDbg=False, color=None)

Generator of random dense packing with given geometry properties, using TriaxialTest (aperiodic) or PeriIsoCompressor (periodic). The priodicity depens on whether the spheresInCell parameter is given.

O.switchScene() magic is used to have clean simulation for TriaxialTest without deleting the original simulation. This function therefore should never run in parallel with some code accessing your simulation.

Parameters:
  • predicate – solid-defining predicate for which we generate packing
  • spheresInCell – if given, the packing will be periodic, with given number of spheres in the periodic cell.
  • radius – mean radius of spheres
  • rRelFuzz – relative fuzz of the radius – e.g. radius=10, rRelFuzz=.2, then spheres will have radii 10 ± (10*.2)). 0 by default, meaning all spheres will have exactly the same radius.
  • cropLayers – (aperiodic only) how many layers of spheres will be added to the computed dimension of the box so that there no (or not so much, at least) boundary effects at the boundaries of the predicate.
  • dim – dimension of the packing, to override dimensions of the predicate (if it is infinite, for instance)
  • memoizeDb – name of sqlite database (existent or nonexistent) to find an already generated packing or to store the packing that will be generated, if not found (the technique of caching results of expensive computations is known as memoization). Fuzzy matching is used to select suitable candidate – packing will be scaled, rRelFuzz and dimensions compared. Packing that are too small are dictarded. From the remaining candidate, the one with the least number spheres will be loaded and returned.
  • useOBB – effective only if a inGtsSurface predicate is given. If true (default), oriented bounding box will be computed first; it can reduce substantially number of spheres for the triaxial compression (like 10× depending on how much asymmetric the body is), see scripts/test/gts-triax-pack-obb.py.
  • memoDbg – show packigns that are considered and reasons why they are rejected/accepted
Returns:

SpherePack object with spheres, filtered by the predicate.

yade.pack.randomPeriPack(radius, initSize, rRelFuzz=0.0, memoizeDb=None)

Generate periodic dense packing.

A cell of initSize is stuffed with as many spheres as possible, then we run periodic compression with PeriIsoCompressor, just like with randomDensePack.

Parameters:
  • radius – mean sphere radius
  • rRelFuzz – relative fuzz of sphere radius (equal distribution); see the same param for randomDensePack.
  • initSize – initial size of the periodic cell.
Returns:

SpherePack object, which also contains periodicity information.

yade.pack.regularHexa(predicate, radius, gap, **kw)

Return set of spheres in regular hexagonal grid, clipped inside solid given by predicate. Created spheres will have given radius and will be separated by gap space.

yade.pack.regularOrtho(predicate, radius, gap, **kw)

Return set of spheres in regular orthogonal grid, clipped inside solid given by predicate. Created spheres will have given radius and will be separated by gap space.

yade.pack.revolutionSurfaceMeridians(sects, angles, origin=Vector3(0, 0, 0), orientation=Quaternion((1, 0, 0), 0))

Revolution surface given sequences of 2d points and sequence of corresponding angles, returning sequences of 3d points representing meridian sections of the revolution surface. The 2d sections are turned around z-axis, but they can be transformed using the origin and orientation arguments to give arbitrary orientation.

yade.pack.sweptPolylines2gtsSurface(pts, threshold=0, capStart=False, capEnd=False)

Create swept suface (as GTS triangulation) given same-length sequences of points (as 3-tuples).

If threshold is given (>0), then

  • degenerate faces (with edges shorter than threshold) will not be created
  • gts.Surface().cleanup(threshold) will be called before returning, which merges vertices mutually closer than threshold. In case your pts are closed (last point concident with the first one) this will the surface strip of triangles. If you additionally have capStart==True and capEnd==True, the surface will be closed.

Note

capStart and capEnd make the most naive polygon triangulation (diagonals) and will perhaps fail for non-convex sections.

Warning

the algorithm connects points sequentially; if two polylines are mutually rotated or have inverse sense, the algorithm will not detect it and connect them regardless in their given order.

Creation, manipulation, IO for generic sphere packings.

class yade._packSpheres.SpherePack

Set of spheres represented as centers and radii. This class is returned by pack.randomDensePack, pack.randomPeriPack and others. The object supports iteration over spheres, as in

>>> sp=SpherePack()
>>> for center,radius in sp: print center,radius
>>> for sphere in sp: print sphere[0],sphere[1]   ## same, but without unpacking the tuple automatically
>>> for i in range(0,len(sp)): print sp[i][0], sp[i][1]   ## same, but accessing spheres by index

Special constructors

Construct from list of [(c1,r1),(c2,r2),…]. To convert two same-length lists of centers and radii, construct with zip(centers,radii).

__init__([(list)list]) → None

Empty constructor, optionally taking list [ ((cx,cy,cz),r), … ] for initial data.

aabb() → tuple

Get axis-aligned bounding box coordinates, as 2 3-tuples.

add((Vector3)arg2, (float)arg3) → None

Add single sphere to packing, given center as 3-tuple and radius

cellFill((Vector3)arg2) → None

Repeat the packing (if periodic) so that the results has dim() >= given size. The packing retains periodicity, but changes cellSize. Raises exception for non-periodic packing.

cellRepeat((Vector3i)arg2) → None

Repeat the packing given number of times in each dimension. Periodicity is retained, cellSize changes. Raises exception for non-periodic packing.

cellSize

Size of periodic cell; is Vector3(0,0,0) if not periodic. (Change this property only if you know what you’re doing).

center() → Vector3

Return coordinates of the bounding box center.

dim() → Vector3

Return dimensions of the packing in terms of aabb(), as a 3-tuple.

fromList((list)arg2) → None
Make packing from given list, same format as for constructor. Discards current data.
fromList( (SpherePack)arg1, (object)centers, (object)radii) → None :
Make packing from given list, same format as for constructor. Discards current data.
fromSimulation() → None

Make packing corresponding to the current simulation. Discards current data.

getClumps() → tuple

Return lists of sphere ids sorted by clumps they belong to. The return value is (standalones,[clump1,clump2,…]), where each item is list of id’s of spheres.

hasClumps() → bool

Whether this object contains clumps.

load((str)fileName) → None

Load packing from external text file (current data will be discarded).

makeCloud((Vector3)minCorner, (Vector3)maxCorner[, (float)rMean=-1[, (float)rRelFuzz=0[, (int)num=-1[, (bool)periodic=False[, (float)porosity=-1[, (object)psdSizes=[][, (object)psdCumm=[][, (bool)distributeMass=False]]]]]]]]) → int

Create random loose packing enclosed in box. Sphere radius distribution can be specified using one of the following ways (they are mutually exclusive):

  1. rMean and rRelFuzz gives uniform radius distribution between rMean*(1 ± *rRelFuzz).
  2. porosity, num and rRelFuzz which estimates mean radius so that porosity is attained at the end
  3. psdSizes and psdCumm, two arrays specifying points of the particle size distribution function.

By default (with distributeMass==False), the distribution is applied to particle radii. The usual sense of “particle size distribution” is the distribution of mass fraction (rather than particle count); this can be achieved with distributeMass=True.

Parameters:
  • minCorner (Vector3) – lower corner of the box
  • maxCorner (Vector3) – upper corner of the box
  • rMean (float) – mean radius or spheres
  • rRelFuzz (float) – dispersion of radius relative to rMean
  • num (int) – number of spheres to be generated (if negavite, generate as many as possible, ending after a fixed number of tries to place the sphere in space)
  • periodic (bool) – whether the packing to be generated should be periodic
  • porosity (float) – if specified, estimate mean radius r_m (rMean must not be given) using rRelFuzz (z) and num (N) so that the porosity given (\rho) is approximately achieved at the end of generation, r_m=\sqrt[3]{\frac{V(1-\rho)}{\frac{4}{3}\pi(1+z^2)N}}. The value of \rho=0.65 is recommended.
  • psdSizes – sieve sizes (particle diameters) when particle size distribution (PSD) is specified
  • psdCumm – cummulative fractions of particle sizes given by psdSizes; must be the same length as psdSizes and should be non-decreasing
  • distributeMass (bool) – if True, given distribution will be used to distribute sphere’s mass rather than radius of them.
Returns:

number of created spheres, which can be lower than num is the packing is too tight.

makeClumpCloud((Vector3)minCorner, (Vector3)maxCorner, (object)clumps[, (bool)periodic=False[, (int)num=-1]]) → int

Create random loose packing of clumps within box given by minCorner and maxCorner. Clumps are selected with equal probability. At most num clumps will be positioned if num is positive; otherwise, as many clumps as possible will be put in space, until maximum number of attempts to place a new clump randomly is attained.

particleSD((Vector3)minCorner, (Vector3)maxCorner, (float)rMean, (bool)periodic=False, (str)name, (int)numSph[, (object)radii=[][, (object)passing=[][, (bool)passingIsNotPercentageButCount=False]]]) → int

Create random packing enclosed in box given by minCorner and maxCorner, containing numSph spheres. Returns number of created spheres, which can be < num if the packing is too tight. The computation is done according to the given psd.

psd([(int)bins=10[, (bool)mass=False]]) → tuple

Return particle size distribution of the packing. :param int bins: number of bins between minimum and maximum diameter :param mass: Compute relative mass rather than relative particle count for each bin. Corresponds to distributeMass parameter for makeCloud. :returns: tuple of (cumm,edges), where cumm are cummulative fractions for respective diameters and edges are those diameter values. Dimension of both arrays is equal to bins+1.

psdScaleExponent

Exponent used to scale cummulative distribution function values so that standard uniform distribution is mapped to uniform mass distribution. Theoretically, this value is 3, but that usually overfavors small particles. The default value is 2.5.

relDensity() → float

Relative packing density, measured as sum of spheres’ volumes / aabb volume. (Sphere overlaps are ignored.)

rotate((Vector3)axis, (float)angle) → None

Rotate all spheres around packing center (in terms of aabb()), given axis and angle of the rotation.

save((str)fileName) → None

Save packing to external text file (will be overwritten).

scale((float)arg2) → None

Scale the packing around its center (in terms of aabb()) by given factor (may be negative).

toList() → list

Return packing data as python list.

toSimulation()

Append spheres directly to the simulation. In addition calling O.bodies.append, this method also appropriately sets periodic cell information of the simulation.

>>> from yade import pack; from math import * >>> sp=pack.SpherePack()

Create random periodic packing with 20 spheres:

>>> sp.makeCloud((0,0,0),(5,5,5),rMean=.5,rRelFuzz=.5,periodic=True,num=20) 20

Virgin simulation is aperiodic:

>>> O.reset() >>> O.periodic False

Add generated packing to the simulation, rotated by 45° along +z

>>> sp.toSimulation(rot=Quaternion((0,0,1),pi/4),color=(0,0,1)) [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19]

Periodic properties are transferred to the simulation correctly:

>>> O.periodic True >>> O.cell.refSize Vector3(5,5,5) >>> O.cell.trsf Matrix3(0.707107,-0.707107,0, 0.707107,0.707107,0, 0,0,1)

Parameters:
  • rot (Quaternion/Matrix3) – rotation of the packing, which will be applied on spheres and will be used to set Cell.trsf as well.
  • **kw

    passed to utils.sphere

Returns:

list of body ids added (like O.bodies.append)

translate((Vector3)arg2) → None

Translate all spheres by given vector.

class yade._packSpheres.SpherePackIterator
__init__((SpherePackIterator)arg2) → None
next() → tuple

Spatial predicates for volumes (defined analytically or by triangulation).

class yade._packPredicates.Predicate
aabb() → tuple

aabb( (Predicate)arg1) → None

center() → Vector3
dim() → Vector3
class yade._packPredicates.PredicateBoolean(inherits Predicate)

Boolean operation on 2 predicates (abstract class)

A
B
__init__()

Raises an exception This class cannot be instantiated from Python

class yade._packPredicates.PredicateDifference(inherits PredicateBoolean → Predicate)

Difference (conjunction with negative predicate) of 2 predicates. A point has to be inside the first and outside the second predicate. Can be constructed using the - operator on predicates: pred1 - pred2.

__init__((object)arg2, (object)arg3) → None
class yade._packPredicates.PredicateIntersection(inherits PredicateBoolean → Predicate)

Intersection (conjunction) of 2 predicates. A point has to be inside both predicates. Can be constructed using the & operator on predicates: pred1 & pred2.

__init__((object)arg2, (object)arg3) → None
class yade._packPredicates.PredicateSymmetricDifference(inherits PredicateBoolean → Predicate)

SymmetricDifference (exclusive disjunction) of 2 predicates. A point has to be in exactly one predicate of the two. Can be constructed using the ^ operator on predicates: pred1 ^ pred2.

__init__((object)arg2, (object)arg3) → None
class yade._packPredicates.PredicateUnion(inherits PredicateBoolean → Predicate)

Union (non-exclusive disjunction) of 2 predicates. A point has to be inside any of the two predicates to be inside. Can be constructed using the | operator on predicates: pred1 | pred2.

__init__((object)arg2, (object)arg3) → None
class yade._packPredicates.inAlignedBox(inherits Predicate)

Axis-aligned box predicate

__init__((Vector3)minAABB, (Vector3)maxAABB) → None

Ctor taking minumum and maximum points of the box (as 3-tuples).

class yade._packPredicates.inCylinder(inherits Predicate)

Cylinder predicate

__init__((Vector3)centerBottom, (Vector3)centerTop, (float)radius) → None

Ctor taking centers of the lateral walls (as 3-tuples) and radius.

class yade._packPredicates.inEllipsoid(inherits Predicate)

Ellipsoid predicate

__init__((Vector3)centerPoint, (Vector3)abc) → None

Ctor taking center of the ellipsoid (3-tuple) and its 3 radii (3-tuple).

class yade._packPredicates.inGtsSurface(inherits Predicate)

GTS surface predicate

__init__((object)surface[, (bool)noPad]) → None

Ctor taking a gts.Surface() instance, which must not be modified during instance lifetime. The optional noPad can disable padding (if set to True), which speeds up calls several times. Note: padding checks inclusion of 6 points along +- cardinal directions in the pad distance from given point, which is not exact.

surf

The associated gts.Surface object.

class yade._packPredicates.inHyperboloid(inherits Predicate)

Hyperboloid predicate

__init__((Vector3)centerBottom, (Vector3)centerTop, (float)radius, (float)skirt) → None

Ctor taking centers of the lateral walls (as 3-tuples), radius at bases and skirt (middle radius).

class yade._packPredicates.inSphere(inherits Predicate)

Sphere predicate.

__init__((Vector3)center, (float)radius) → None

Ctor taking center (as a 3-tuple) and radius

class yade._packPredicates.notInNotch(inherits Predicate)

Outside of infinite, rectangle-shaped notch predicate

__init__((Vector3)centerPoint, (Vector3)edge, (Vector3)normal, (float)aperture) → None

Ctor taking point in the symmetry plane, vector pointing along the edge, plane normal and aperture size. The side inside the notch is edge×normal. Normal is made perpendicular to the edge. All vectors are normalized at construction time.

Computation of oriented bounding box for cloud of points.

yade._packObb.cloudBestFitOBB((tuple)arg1) → tuple

Return (Vector3 center, Vector3 halfSize, Quaternion orientation) of best-fit oriented bounding-box for given tuple of points (uses brute-force velome minimization, do not use for very large clouds).

class yade._packSpherePadder.SpherePadder

Geometrical algorithm for filling tetrahedral mesh with spheres; the algorithm was designed by Jean-François Jerier and is described in [Jerier2009].

__init__((str)fileName[, (str)meshType='']) → None

Initialize using tetrahedral mesh stored in fileName. Type of file is determined by extension: .gmsh implies meshType*=’GMSH’, .inp implies *meshType*=’INP’. If the extension is different, specify *meshType explicitly. Possible values are ‘GMSH’ and ‘INP’.

asSpherePack() → SpherePack
densify() → None
insert_sphere((float)arg2, (float)arg3, (float)arg4, (float)arg5) → None
maxNumberOfSpheres
maxOverlapRate
maxSolidFractioninProbe
meanSolidFraction
numberOfSpheres
pad_5() → None
place_virtual_spheres() → None
radiusRange
radiusRatio
save_mgpost((str)arg2) → None
setRadiusRatio((float)arg2, (float)arg3) → None

Like radiusRatio, but taking 2nd parameter.

virtualRadiusFactor

Previous topic

yade.log module

Next topic

yade.plot module

This Page