woo.pack
¶
Creating packings and filling volumes defined by boundary representation or constructive solid geometry.
Warning
Broken links ahead:
For examples, see
-
woo.pack.
SpherePack_fromSimulation
(self, scene, dem=None)[source]¶ Reset this SpherePack object and initialize it from the current simulation; only spherical particles are taken in account. Periodic boundary conditions are supported, but the hSize matrix must be diagonal. Clumps are supported.
-
woo.pack.
SpherePack_toDem
(self, scene, dem=None, rot=Matrix3(1, 0, 0, 0, 1, 0, 0, 0, 1), **kw)[source]¶ Append spheres directly to the simulation. In addition calling
woo.dem.ParticleContainer.add
, this method also appropriately sets periodic cell information of the simulation.>>> from woo import pack; from math import *; from woo.dem import *; from woo.core 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:
>>> scene=Scene(fields=[DemField()]) >>> scene.periodic False
Add generated packing to the simulation, rotated by 45° along +z
>>> sp.toDem(scene,scene.dem,rot=Quaternion((0,0,1),pi/4),color=0) [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, including rotation:
>>> scene.periodic True >>> scene.cell.size Vector3(5,5,5) >>> scene.cell.hSize Matrix3(3.5355339059327373,-3.5355339059327378,0, 3.5355339059327378,3.5355339059327373,0, 0,0,5)
The current state (even if rotated) is taken as mechanically undeformed, i.e. with identity transformation:
>>> scene.cell.trsf Matrix3(1,0,0, 0,1,0, 0,0,1)
- Parameters
rot (Quaternion/Matrix3) – rotation of the packing, which will be applied on spheres and will be used to set
woo.core.Cell.trsf
as well.kw – passed to
woo.utils.sphere
- Returns
list of body ids added (like
woo.dem.ParticleContainer.add
)
-
woo.pack.
SpherePack_toSimulation
(self, scene, rot=Matrix3(1, 0, 0, 0, 1, 0, 0, 0, 1), **kw)[source]¶ Proxy for old arguments, which merely calls
SpherePack.toDem
.
-
woo.pack.
filterSpherePack
(predicate, spherePack, **kw)[source]¶ 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.
-
woo.pack.
gtsSurface2Facets
(surf, shareNodes=True, flex=False, **kw)[source]¶ Construct facets from given GTS surface.
**kw
is passed towoo.dem.Facet.make
orwoo.fem.Mambrane.make
, depending on flex.
-
woo.pack.
gtsSurfaceBestFitOBB
(surf)[source]¶ Return (Vector3 center, Vector3 halfSize, Quaternion orientation) describing best-fit oriented bounding box (OBB) for the given surface. See cloudBestFitOBB for details.
-
woo.pack.
hexaNet
(radius, cornerCoord=[0, 0, 0], xLength=1.0, yLength=0.5, mos=0.08, a=0.04, b=0.04, startAtCorner=True, isSymmetric=False, **kw)[source]¶ Definition of the particles for a hexagonal wire net in the x-y-plane for the WireMatPM.
- Parameters
radius – radius of the particle
cornerCoord – coordinates of the lower left corner of the net
xLenght – net length in x-direction
yLenght – net length in y-direction
mos – mesh opening size
a – length of double-twist
b – height of single wire section
startAtCorner – if true the generation starts with a double-twist at the lower left corner
isSymmetric – defines if the net is symmetric with respect to the y-axis
- Returns
set of spheres which defines the net (net) and exact dimensions of the net (lx,ly).
note:: This packing works for the WireMatPM only. The particles at the corner are always generated first.
-
class
woo.pack.
inGtsSurface_py
(surf, noPad=False)[source]¶ 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.
-
class
woo.pack.
inSpace
(_center=Vector3(0, 0, 0))[source]¶ Predicate returning True for any points, with infinite bounding box.
-
woo.pack.
makeBandFeedPack
(dim, mat, gravity, psd=[], excessWd=None, damping=0.3, porosity=0.5, goal=0.15, dtSafety=0.9, dontBlock=False, memoizeDir=None, botLine=None, leftLine=None, rightLine=None, clumps=[], returnSpherePack=False, useEnergy=True, gen=None, bias=None)[source]¶ Create dense packing periodic in the +x direction, suitable for use with ConveyorInlet. :param useEnergy: use
woo.utils.unbalancedEnergy
instead ofwoo.utils.unbalancedForce
as stop criterion. :param goal: target unbalanced force/energy; if unbalanced energy is used, this value is multiplied by .2. :param psd: particle size distribution :param mat: material for particles :param gravity: gravity acceleration (as Vector3)
-
woo.pack.
makePeriodicFeedPack
(dim, psd, lenAxis=0, damping=0.3, porosity=0.5, goal=0.15, maxNum=- 1, dontBlock=False, returnSpherePack=False, memoizeDir=None, clumps=None, gen=None)[source]¶
-
woo.pack.
randomDensePack
(predicate, radius, mat=- 1, dim=None, cropLayers=0, rRelFuzz=0.0, spheresInCell=0, memoizeDb=None, useOBB=True, memoDbg=False, color=None)[source]¶ Generator of random dense packing with given geometry properties, using TriaxialTest (aperiodic) or PeriIsoCompressor (periodic). The periodicity depens on whether the spheresInCell parameter is given.
- 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.
-
woo.pack.
randomDensePack2
(predicate, generator, settle=0.3, approxLoosePoro=0.1, maxNum=5000, memoizeDir=None, debug=False, porosity=None)[source]¶ Return dense arrangement of particles clipped by predicate. Particle are losely generated in cuboid volume of sufficient dimensions; the cuboid might be made smaller if maxNum were significantly exceeded, and tiled into the predicate volume afterwards.
The algorithm computes initial (loose) volume from predicate bounding box using approxLoosePoro, and fills periodic cell of that size and aspect ratio (or of a part of it, if maxNum would be significantly exceeded) with generated particles. Isotropic compression (via
woo.dem.PeriIsoCompressor
) is subsequently applied with a predefined material (young 10MPa, no friction) down to 100MPa (and stabilization) and then the sample is unloaded to 1MPa and stabilized (those values are currently hard-coded). The resulting packing is tiled (if smaller than predicate), clipped by the predicate and returned as awoo.dem.ShapePack
object (geometry only, no material data).The packing is not completely overlap-free, thus using options like
woo.dem.Law2_L6Geom_FrictPhys_IdealElPl.iniEqlb
might be useful, depending on the stiffness used in the simulation.- Parameters
predicate (woo.pack.Predicate) – volume definition
generator (woo.dem.ParticleGenerator) – particle definition
settle (float) – relative final volume of dense packing, relative to the loose packing used for generation; usually not necessary to adjust for usual particles, but might have to be made smaller is the initial packing is very loose (significantly elongated particles or clumps and similar)
maxNum (int) – maximum number of initial particle; if the number of particles for the initial volume were higher, the predicate volume will be split into smaller periodic cells and the packing will be simply repeated (tiled) to cover the entire predicate volume afterwards.
approxLoosePoro (float) – estimate of loose (random) packing proosity for given particles, used for estimating number of particles. This number is NOT related to the resulting porosity,, it only has effect on how the predicate volume is partitioned, and thus does not need to be set in most cases.
memoizeDir – path where memoized packings (named as hashed parameters, plus the .randomdense extension) are stored; if None, returned packings are not memoized. If memoized result is found, it is returned immediately, instead of running (often costly) simulation of compaction.
debug (bool) – if True, returns
woo.core.Scene
used for compaction, which can be inspected (instead of the usualwoo.dem.ShapePack
with dense packing)
- Returns
woo.dem.ShapePack
, unless debug isTrue
, in which case aScene
is returned.
-
woo.pack.
randomLoosePsd
(predicate, psd, mass=True, discrete=False, maxAttempts=5000, clumps=[], returnSpherePack=False, **kw)[source]¶ Return loose packing based on given PSD.
-
woo.pack.
randomPeriPack
(radius, initSize, rRelFuzz=0.0, memoizeDb=None)[source]¶ 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.
-
woo.pack.
regularHexa
(predicate, radius, gap, retSpherePack=True, **kw)[source]¶ 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.
-
woo.pack.
regularOrtho
(predicate, radius, gap, retSpherePack=True, **kw)[source]¶ 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.
-
woo.pack.
revolutionSurfaceMeridians
(sects, angles, node=<Node @ 0x21cdbb0, at (0, 0, 0)>, origin=None, orientation=None)[source]¶ 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.
-
woo.pack.
sweptPolylines2gtsSurface
(pts, threshold=0, localCoords=None, capStart=False, capEnd=False)[source]¶ Create swept suface (as GTS triangulation) given same-length sequences of points (as 3-tuples). If localCoords is given, it must be a
woo.core.Node
instance and will be used to convert local coordinates in pts to global coordinates.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.
woo._packPredicates
¶
Note
This module is imported into the woo.pack
module automatically; refer to its objects through woo.pack
.
Spatial predicates for volumes (defined analytically or by triangulation).
-
class
woo._packPredicates.
PredicateBoolean
¶ Boolean operation on 2 predicates (abstract class)
-
property
A
¶
-
property
B
¶
-
property
-
class
woo._packPredicates.
PredicateDifference
¶ 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
.
-
class
woo._packPredicates.
PredicateIntersection
¶ Intersection (conjunction) of 2 predicates. A point has to be inside both predicates. Can be constructed using the
&
operator on predicates:pred1 & pred2
.
-
class
woo._packPredicates.
PredicateSymmetricDifference
¶ 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
.
-
class
woo._packPredicates.
PredicateUnion
¶ 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
.
-
class
woo._packPredicates.
inAlignedBox
¶ Axis-aligned box predicate
Ctor taking minumum and maximum points of the box.
-
class
woo._packPredicates.
inAlignedHalfspace
¶ Half-space given by coordinate at some axis.
Ctor taking axis (0,1,2 for \(x\), \(y\), \(z\) respectively), the coordinate along that axis, and whether the lower half (or the upper half, if lower is false) is considered.
-
class
woo._packPredicates.
inAxisRange
¶ Range of coordinate along some axis, effectively defining two axis-aligned enclosing planes.
Ctor taking axis (0,1,2 for \(x\), \(y\), \(z\) respectively), and coordinate range along that axis.
-
class
woo._packPredicates.
inCylSector
¶ Sector of an arbitrarily oriented cylinder in 3d, limiting cylindrical coordinates \(\rho\), \(\theta\), \(z\); all coordinate limits are optional.
Ctor taking position and orientation of the local system, and aligned box in local coordinates.
-
class
woo._packPredicates.
inCylinder
¶ Cylinder predicate
Ctor taking centers of the lateral walls and radius.
-
class
woo._packPredicates.
inEllipsoid
¶ Ellipsoid predicate
Ctor taking center of the ellipsoid and its 3 radii.
-
class
woo._packPredicates.
inGtsSurface
¶ GTS surface predicate
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.
-
property
surf
¶ The associated gts.Surface object.
-
property
-
class
woo._packPredicates.
inHyperboloid
¶ Hyperboloid predicate
Ctor taking centers of the lateral walls, radius at bases and skirt (middle radius).
-
class
woo._packPredicates.
inOrientedBox
¶ Arbitrarily oriented box specified as local coordinates (pos,ori) and aligned box in those local coordinates.
Ctor taking position and orientation of the local system, and aligned box in local coordinates.
-
class
woo._packPredicates.
inParallelepiped
¶ Parallelepiped predicate
Ctor taking four points:
o
(for origin) and thena
,b
,c
which define endpoints of 3 respective edges fromo
.
-
class
woo._packPredicates.
inSphere
¶ Sphere predicate.
Ctor taking center and radius
-
class
woo._packPredicates.
notInNotch
¶ Outside of infinite, rectangle-shaped notch predicate
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.
woo._packSpheres
¶
Note
This module is imported into the woo.pack
module automatically; refer to its objects through woo.pack
.
Creation, manipulation, IO for generic sphere packings.
-
class
woo._packSpheres.
SpherePack
¶ Set of spheres represented as centers and radii. This class is returned by
woo.pack.randomDensePack
,woo.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
Clumps are supported, by using the clumpId parameter to
add
.Special constructors
Construct from list of
[(c1,r1),(c2,r2),…]
. To convert two same-length lists ofcenters
andradii
, construct withzip(centers,radii)
.Empty constructor, optionally taking list [ ((cx,cy,cz),r), … ] for initial data.
-
aabb
()¶ Get axis-aligned bounding box coordinates.
-
add
()¶ Add single sphere to packing, given center as 3-tuple and radius
-
property
appliedPsdScaling
¶ A factor between 0 and 1, uniformly applied on all sizes of the PSD.
-
canonicalize
()¶ Move all sphere’s centers inside cellSize; only works for periodic packings without clumps.
-
cellFill
()¶ 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
()¶ Repeat the packing given number of times in each dimension. Periodicity is retained, cellSize changes. Raises exception for non-periodic packing.
-
property
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
()¶ Return coordinates of the bounding box center.
-
dim
()¶ Return dimensions of the packing in terms of aabb().
-
filtered
()¶ Return new
SpherePack
object, without any spheres which don’t match predicate. Clumps are handled gracefully, i.e. if any of clump’s spheres does not satisfy the predicate, the whole clump is taken away. If recenter is True (and none of the dimensions of the predicate is infinite), the packing will be translated to have center at the same point as the predicate.
-
fromDem
(scene, dem)¶
-
fromList
()¶ Make packing from given list, same format as for constructor. Discards current data. Make packing from given list, same format as for constructor. Discards current data.
-
fromSimulation
(scene, dem=None)¶ Reset this SpherePack object and initialize it from the current simulation; only spherical particles are taken in account. Periodic boundary conditions are supported, but the hSize matrix must be diagonal. Clumps are supported.
-
getClumps
()¶ 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
()¶ Whether this object contains clumps.
-
load
()¶ Load packing from external text file (current data will be discarded).
-
makeCloud
()¶ Create random loose packing enclosed in a parallelepiped. Sphere radius distribution can be specified using one of the following ways:
rMean, rRelFuzz and num gives uniform radius distribution in rMean (1 ± rRelFuzz ). Less than num spheres can be generated if it is too high.
rRelFuzz, num and (optional) porosity, which estimates mean radius so that porosity is attained at the end. rMean must be less than 0 (default). porosity is only an initial guess for the generation algorithm, which will retry with higher porosity until the prescibed num is obtained.
psdSizes and psdCumm, two arrays specifying points of the particle size distribution function. As many spheres as possible are generated.
psdSizes, psdCumm, num, and (optional) porosity, like above but if num is not obtained, psdSizes will be scaled down uniformly, until num is obtained (see
appliedPsdScaling
).
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 withdistributeMass=True
.If num is defined, then sizes generation is deterministic, giving the best fit of target distribution. It enables spheres placement in descending size order, thus giving lower porosity than the random generation.
- Parameters
minCorner (Vector3) – lower corner of an axis-aligned box
maxCorner (Vector3) – upper corner of an axis-aligned box
hSize (Matrix3) – base vectors of a generalized box (arbitrary parallelepiped, typically
woo.core.Cell.hSize
), superseeds minCorner and maxCorner if defined. For periodic boundaries only.rMean (float) – mean radius or spheres
rRelFuzz (float) – dispersion of radius relative to rMean
num (int) – number of spheres to be generated. If negavite (default), generate as many as possible with stochastic sizes, ending after a fixed number of tries to place the sphere in space, else generate exactly num spheres with deterministic size distribution.
periodic (bool) – whether the packing to be generated should be periodic
porosity (float) – initial guess for the iterative generation procedure (if num>1). The algorithm will be retrying until the number of generated spheres is num. The first iteration tries with the provided porosity, but next iterations increase it if necessary (hence an initialy high porosity can speed-up the algorithm). If psdSizes is not defined, rRelFuzz (\(z\)) and num (\(N\)) are used 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 default is :math:`rho`=0.5. The optimal value depends on rRelFuzz or psdSizes.
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.seed – number used to initialize the random number generator.
- Returns
number of created spheres, which can be lower than num depending on the method used.
-
makeOverlapFree
()¶ Scale by 1+maxRelOverlap(), without changing radii.
-
maxRelOverlap
()¶ Return maximum relative overlap of particles.
-
psd
()¶ 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)
, wherecumm
are cummulative fractions for respective diameters andedges
are those diameter values. Dimension of both arrays is equal tobins+1
.
-
relDensity
()¶ Relative packing density, measured as sum of spheres’ volumes / aabb volume. (Sphere overlaps are ignored.)
-
reset
()¶ Re-inialize this object (clear all spheres and reset periodic cell size.
-
rotate
()¶ Rotate all spheres around packing center (in terms of aabb()), given axis and angle of the rotation.
-
save
()¶ Save packing to external text file (will be overwritten).
-
saveTxt
()¶ Identical to
save:
, with newer name for ShapePack compatibility.
-
scale
()¶ Scale the packing around its center (in terms of aabb()) by given factor (may be negative). If periodic, scale around origin. If keepRadius is True, radii are not scaled, only coordinates.
-
solidVolume
()¶ Alias for
sphereVol
, for forward compatibility withwoo.dem.ShapePack
.
-
sphereVol
()¶ Summary volume of spheres, disregarding any overlaps (\(\frac{4}{3}\pi\sum r_i^3\)).
-
toCcRr
()¶ Return packing data as ([c,c,…],[r,r,…]). Raises exception if there are clumps.
-
toDem
(scene, dem=None, rot=Matrix3(1, 0, 0, 0, 1, 0, 0, 0, 1), **kw)¶ Append spheres directly to the simulation. In addition calling
woo.dem.ParticleContainer.add
, this method also appropriately sets periodic cell information of the simulation.>>> from woo import pack; from math import *; from woo.dem import *; from woo.core 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:
>>> scene=Scene(fields=[DemField()]) >>> scene.periodic False
Add generated packing to the simulation, rotated by 45° along +z
>>> sp.toDem(scene,scene.dem,rot=Quaternion((0,0,1),pi/4),color=0) [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, including rotation:
>>> scene.periodic True >>> scene.cell.size Vector3(5,5,5) >>> scene.cell.hSize Matrix3(3.5355339059327373,-3.5355339059327378,0, 3.5355339059327378,3.5355339059327373,0, 0,0,5)
The current state (even if rotated) is taken as mechanically undeformed, i.e. with identity transformation:
>>> scene.cell.trsf Matrix3(1,0,0, 0,1,0, 0,0,1)
- Parameters
rot (Quaternion/Matrix3) – rotation of the packing, which will be applied on spheres and will be used to set
woo.core.Cell.trsf
as well.kw – passed to
woo.utils.sphere
- Returns
list of body ids added (like
woo.dem.ParticleContainer.add
)
-
toList
()¶ Return packing data as [(c,r),(c,r),…].
-
toSimulation
(scene, rot=Matrix3(1, 0, 0, 0, 1, 0, 0, 0, 1), **kw)¶ Proxy for old arguments, which merely calls
SpherePack.toDem
.
-
toSimulation_fast
()¶ Create spheres from this pack in the simulation; unlike
SpherePack.toSimulation
, this method is implemented in c++ (hence rather fast) but lacks some flexibility
-
translate
()¶ Translate all spheres by given vector.
-
property
userData
¶ Arbitrary string (not cotaining newlines) which will be saved and loaded with this object
-
-
class
woo._packSpheres.
SpherePackIterator
¶
woo._packObb
¶
Note
This module is imported into the woo.pack
module automatically; refer to its objects through woo.pack
.
Computation of oriented bounding box for cloud of points.
-
woo._packObb.
cloudBestFitOBB
()¶ 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).
Tip
Report issues or inclarities to github.