underworld.swarm module

This module contains routines relating to swarm type objects.

Module Summary

classes:

underworld.swarm.GaussBorderIntegrationSwarm Integration swarm which creates particles within the boundary faces of an element, at the Gauss points.
underworld.swarm.GaussIntegrationSwarm Integration swarm which creates particles within an element at the Gauss points.
underworld.swarm.IntegrationSwarm Abstract class definition for IntegrationSwarms.
underworld.swarm.PopulationControl This class implements swarm population control mechanism.
underworld.swarm.Swarm The Swarm class supports particle like data structures.
underworld.swarm.SwarmAbstract The SwarmAbstract class supports particle like data structures.
underworld.swarm.SwarmVariable The SwarmVariable class allows users to add data to swarm particles. The data
underworld.swarm.VoronoiIntegrationSwarm Class for an IntegrationSwarm that maps to another Swarm

Module Details

classes:

class underworld.swarm.GaussBorderIntegrationSwarm(mesh, particleCount=None, **kwargs)[source]

Bases: underworld.swarm._integration_swarm.GaussIntegrationSwarm

Integration swarm which creates particles within the boundary faces of an element, at the Gauss points.

See parent class for parameters.

class underworld.swarm.GaussIntegrationSwarm(mesh, particleCount=None, **kwargs)[source]

Bases: underworld.swarm._integration_swarm.IntegrationSwarm

Integration swarm which creates particles within an element at the Gauss points.

Parameters:
  • mesh (underworld.mesh.FeMesh) – The FeMesh the swarm is supported by. See Swarm.mesh property docstring for further information.
  • particleCount (unsigned. Default is 3, unless Q2 mesh which takes default 5.) – Number of gauss particles in each direction. Must take value in [1,5].
class underworld.swarm.IntegrationSwarm(mesh, **kwargs)[source]

Bases: underworld.swarm._swarmabstract.SwarmAbstract

Abstract class definition for IntegrationSwarms.

All IntegrationSwarms have the following SwarmVariables from this class:

  1. localCoordVariable : double (number of particle, dimensions)
    For local element coordinates of the particle
  2. weightVariable : double (number of particles)
    For the integration weight of each particle
particleWeights
Returns:Swarm variable recording the weight of the swarm particles.
Return type:underworld.swarm.SwarmVariable
class underworld.swarm.PopulationControl(swarm, deleteThreshold=0.006, splitThreshold=0.25, maxDeletions=0, maxSplits=3, aggressive=False, aggressiveThreshold=0.8, particlesPerCell=None, **kwargs)[source]

Bases: underworld._stgermain.LeftOverParamsChecker

This class implements swarm population control mechanism. Population control acts on a per element basic, with a discrete voronoi algorithm is used to determine where particles should be added or removed.

Parameters:
  • swarm (underworld.swarm.Swarm) – The swarm for which population control should occur.
  • deleteThreshold (float) – Particle volume fraction threshold below which particle is deleted. i.e if (particleVolume/elementVolume)<deleteThreshold, then the particle is deleted.
  • splitThreshold (float) – Particle volume fraction threshold above which particle is split. i.e if (particleVolume/elementVolume)>splitThreshold, then the particle is split.
  • maxDeletions (int) – maximum number of particle deletions per cell
  • maxSplits (int) – maximum number of particles splits per cell
  • aggressive (bool) – When enabled, this option will invoke aggressive population control in elements where particle counts drop below some threshold.
  • aggressiveThreshold (float) – lower cell particle population threshold beyond which aggressive population control occurs. i.e if (cellParticleCount/particlesPerCell)<aggressiveThreshold, then aggressive pop control will occur. Note that this option is only valid if ‘aggressive’ is enabled.
  • particlesPerCell (int) – This is the desired number of particles each element should contain. Note that this option is only valid if ‘aggressive’ is enabled.

Example

This simple example generates a swarm, then applies population control to split particles.

>>> import underworld as uw
>>> import numpy as np
>>> mesh = uw.mesh.FeMesh_Cartesian()
>>> swarm = uw.swarm.Swarm(mesh)
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellGaussLayout(swarm,4))
>>> population_control = uw.swarm.PopulationControl(swarm,deleteThreshold=0.,splitThreshold=0.,maxDeletions=0,maxSplits=9999)
>>> population_control.repopulate()
>>> swarm.particleGlobalCount
512
repopulate()[source]

This method repopulates the swarm.

class underworld.swarm.Swarm(mesh, particleEscape=False, **kwargs)[source]

Bases: underworld.swarm._swarmabstract.SwarmAbstract, underworld.function._function.FunctionInput, underworld._stgermain.Save

The Swarm class supports particle like data structures. Each instance of this class will store a set of unique particles. In this context, particles are data structures which store a location variable, along with any other variables the user requests.

Parameters:
  • mesh (underworld.mesh.FeMesh) – The FeMesh the swarm is supported by. See Swarm.mesh property docstring for further information.
  • particleEscape (bool) – If set to true, particles are deleted when they leave the domain. This may occur during particle advection, or when the mesh is deformed.

Example

Create a swarm with some variables:

>>> # First we need a mesh:
>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> # Create empty swarm:
>>> swarm = uw.swarm.Swarm(mesh)
>>> # Add a variable:
>>> svar = swarm.add_variable("char",1)
>>> # Add another:
>>> svar = swarm.add_variable("double",3)
>>> # Now use a layout to fill with particles
>>> swarm.particleLocalCount
0
>>> layout = uw.swarm.layouts.PerCellGaussLayout(swarm,2)
>>> swarm.populate_using_layout(layout)
>>> swarm.particleLocalCount
1024
>>> swarm.data[0]
array([ 0.0132078,  0.0132078])
>>> swarm.owningCell.data[0]
array([0], dtype=int32)

With particleEscape enabled, particles which are no longer within the mesh domain are deleted.

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> swarm = uw.swarm.Swarm(mesh, particleEscape=True)
>>> swarm.particleLocalCount
0
>>> layout = uw.swarm.layouts.PerCellGaussLayout(swarm,2)
>>> swarm.populate_using_layout(layout)
>>> swarm.particleGlobalCount
1024
>>> with mesh.deform_mesh():
...     mesh.data[:] += (0.5,0.)
>>> swarm.particleGlobalCount
512

Alternatively, moving the particles:

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> swarm = uw.swarm.Swarm(mesh, particleEscape=True)
>>> swarm.particleLocalCount
0
>>> layout = uw.swarm.layouts.PerCellGaussLayout(swarm,2)
>>> swarm.populate_using_layout(layout)
>>> swarm.particleGlobalCount
1024
>>> with swarm.deform_swarm():
...     swarm.data[:] -= (0.5,0.)
>>> swarm.particleGlobalCount
512
add_particles_with_coordinates(coordinatesArray)[source]

This method adds particles to the swarm using particle coordinates provided using a numpy array.

Note that particles with coordinates NOT local to the current processor will be reject/ignored.

Parameters:coordinatesArray (numpy.ndarray) – The numpy array containing the coordinate of the new particles. Array is expected to take shape n*dim, where n is the number of new particles, and dim is the dimensionality of the swarm’s supporting mesh.
Returns:Array containing the local index of the added particles. Rejected particles are denoted with an index of -1.
Return type:numpy.ndarray

Example

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(4,4), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> swarm = uw.swarm.Swarm(mesh)
>>> import numpy as np
>>> arr = np.zeros((5,2))
>>> arr[0] = [0.1,0.1]
>>> arr[1] = [0.2,0.1]
>>> arr[2] = [0.1,0.2]
>>> arr[3] = [-0.1,-0.1]
>>> arr[4] = [0.8,0.8]
>>> swarm.add_particles_with_coordinates(arr)
array([ 0,  1,  2, -1,  3], dtype=int32)
>>> swarm.particleLocalCount
4
>>> swarm.data
array([[ 0.1,  0.1],
       [ 0.2,  0.1],
       [ 0.1,  0.2],
       [ 0.8,  0.8]])
allow_parallel_nn

By default, parallel nearest neighbour search is disabled as consistent results are not currently guaranteed. Set this attribute to True to allow parallel NN.

deform_swarm(update_owners=True)[source]

Any particle location modifications must occur within this python context manager. This is necessary as it is critical that certain internal objects are updated when particle locations are modified. Note that this method must be called collectively by all processes, irrespective of whether any given process does or does not need to locate any particles.

Parameters:update_owners (bool, default=True) – If this is set to False, particle ownership (which element owns a particular particle) is not updated at the conclusion of the context manager. This is often necessary when both the mesh and particles are advecting simutaneously.

Notes

This method must be called collectively by all processes.

Example

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> swarm = uw.swarm.Swarm(mesh)
>>> layout = uw.swarm.layouts.PerCellGaussLayout(swarm,2)
>>> swarm.populate_using_layout(layout)
>>> swarm.data[0]
array([ 0.0132078,  0.0132078])

Attempted modification without using deform_swarm() should fail:

>>> swarm.data[0] = [0.2,0.2]
Traceback (most recent call last):
...
ValueError: assignment destination is read-only

Within the deform_swarm() context manager, modification is allowed:

>>> with swarm.deform_swarm():
...     swarm.data[0] = [0.2,0.2]
>>> swarm.data[0]
array([ 0.2,  0.2])
fn_particle_found()[source]

This function returns True where a particle is able to be found using the provided input to function evaluation.

Returns:The function object.
Return type:underworld.function.Function

Example

Setup some things:

>>> import numpy as np
>>> mesh = uw.mesh.FeMesh_Cartesian(elementRes=(32,32))
>>> swarm = uw.swarm.Swarm(mesh, particleEscape=True)
>>> layout = uw.swarm.layouts.PerCellGaussLayout(swarm,2)
>>> swarm.populate_using_layout(layout)

Now, evaluate the fn_particle_found function on the swarm.. all should be true

>>> fn_pf = swarm.fn_particle_found()
>>> fn_pf.evaluate(swarm).all()
True

Evalute at arbitrary coord… should return False

>>> fn_pf.evaluate( (0.3,0.9) )
array([[False]], dtype=bool)

Now, lets get rid of all particles outside of a circle, and look to obtain pi/4. First eject particles:

>>> with swarm.deform_swarm():
...    for ind,coord in enumerate(swarm.data):
...        if np.dot(coord,coord)>1.:
...            swarm.data[ind] = (99999.,99999.)

Now integrate and test

>>> incirc = uw.function.branching.conditional( ( (fn_pf,1.),(True,0.)  ) )
>>> np.isclose(uw.utils.Integral(incirc, mesh).evaluate(),np.pi/4.,rtol=2e-2)
array([ True], dtype=bool)
load(filename, collective=False, try_optimise=True)[source]

Load a swarm from disk. Note that this must be called before any SwarmVariable members are loaded. Note also that this method will not replace the currently existing particles, but will instead add all the particles from the provided h5 file to the current swarm. Usually you will want to call this method when the swarm is still empty.

Parameters:
  • filename (str) – The filename for the saved file. Relative or absolute paths may be used.
  • collective (bool) – If True, swarm is loaded MPI collective. This is usually faster, but currently is problematic for passive swarms which may not have representation on all processes.
  • try_optimise (bool, Default=True) – Will speed up the swarm load time but warning - this algorithm assumes the previously saved swarm data was made on an identical mesh and mesh partitioning (number of processors) with respect to the current mesh. If this isn’t the case then the reloaded particle ordering will be broken, leading to an invalid swarms. One can disable this optimisation and revert to a brute force algorithm, much slower, by setting this option to False.

Notes

This method must be called collectively by all processes.

Example

Refer to example provided for ‘save’ method.

particleGlobalCount
Returns:The global number (across all processes) of particles in the swarm.
Return type:int
save(filename, collective=False)[source]

Save the swarm to disk.

Parameters:
  • filename (str) – The filename for the saved file. Relative or absolute paths may be used, but all directories must exist.
  • collective (bool) – If True, swarm is saved MPI collective. This is usually faster, but currently is problematic for passive swarms which may not have representation on all processes.
Returns:

Data object relating to saved file. This only needs to be retained if you wish to create XDMF files and can be ignored otherwise.

Return type:

underworld.utils.SavedFileData

Notes

This method must be called collectively by all processes.

Example

First create the swarm, and populate with layout:

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> swarm = uw.swarm.Swarm(mesh)
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellGaussLayout(swarm,2))

Save to a file:

>>> ignoreMe = swarm.save("saved_swarm.h5")

Now let’s try and reload. First create an empty swarm, and then load:

>>> clone_swarm = uw.swarm.Swarm(mesh)
>>> clone_swarm.load( "saved_swarm.h5" )

Now check for equality:

>>> import numpy as np
>>> np.allclose(swarm.data,clone_swarm.data)
True
>>> # clean up:
>>> if uw.mpi.rank == 0:
...     import os;
...     os.remove( "saved_swarm.h5" )
shadow_particles_fetch()[source]

When called, neighbouring processor particles which have coordinates within the current processor’s shadow zone will be communicated to the current processor. Ie, the processor shadow zone is populated using particles that are owned by neighbouring processors. After this method has been called, particle shadow data is available via the data_shadow handles of SwarmVariable objects. This data is read only.

Note that you will need to call this whenever neighbouring information has potentially changed, for example after swarm advection, or after you have modified a SwarmVariable object.

Any existing shadow information will be discarded when this is called.

Notes

This method must be called collectively by all processes.

update_particle_owners()[source]

This routine will update particles owners after particles have been moved. This is both in terms of the cell/element that the particle resides within, and also in terms of the parallel processor decomposition (particles belonging on other processors will be sent across).

Users should not generally need to call this as it will be called automatically at the conclusion of a deform_swarm() block.

Notes

This method must be called collectively by all processes.

Example

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> swarm = uw.swarm.Swarm(mesh)
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellGaussLayout(swarm,2))
>>> swarm.data[0]
array([ 0.0132078,  0.0132078])
>>> swarm.owningCell.data[0]
array([0], dtype=int32)
>>> with swarm.deform_swarm():
...     swarm.data[0] = [0.1,0.1]
>>> swarm.owningCell.data[0]
array([17], dtype=int32)
class underworld.swarm.SwarmAbstract(mesh, **kwargs)[source]

Bases: underworld._stgermain.StgCompoundComponent

The SwarmAbstract class supports particle like data structures. Each instance of this class will store a set of unique particles. In this context, particles are data structures which store a location variable, along with any other variables the user requests.

Parameters:mesh (underworld.mesh.FeMesh) – The FeMesh the swarm is supported by. See Swarm.mesh property docstring for further information.
add_variable(dataType, count)[source]

Add a variable to each particle in this swarm. Variables can be added at any point. Removal of variables is however not currently supported. See help(SwarmVariable) for further information.

Parameters:
  • dataType (str) – The data type for the variable. Available types are “char”, “short”, “int”, “float” or “double”.
  • count (unsigned) – The number of values to be stored for each particle.
Returns:

The newly created swarm variable.

Return type:

underworld.swarm.SwarmVariable

Example

>>> # first we need a mesh
>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> # create swarm
>>> swarm = uw.swarm.Swarm(mesh)
>>> # add a variable
>>> svar = swarm.add_variable("char",1)
>>> # add another
>>> svar = swarm.add_variable("double",3)
>>> # add some particles
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellGaussLayout(swarm,2))
>>> # add another variable
>>> svar = swarm.add_variable("double",5)
data
Returns:Numpy proxy array to underlying particle coordinate data. Note that this is an alias to swarm.particleCoordinates.data. Check SwarmVariable.data for further info.
Return type:numpy.ndarray
globalId
Returns:Swarm variable recording a particle global identifier. Not yet implemented.
Return type:underworld.swarm.SwarmVariable
mesh
Returns:Supporting FeMesh for this Swarm. All swarms are required to be supported by mesh (or similar) objects, which provide the data structures necessary for efficient particle locating/tracking, as well as the necessary mechanism for the swarm parallel decomposition.
Return type:underworld.mesh.FeMesh
owningCell
Returns:Swarm variable recording the owning cell of the swarm particles. This will usually correspond to the owning element local id.
Return type:underworld.swarm.SwarmVariable
particleCoordinates
Returns:Swarm variable recording the coordinates of the swarm particles.
Return type:underworld.swarm.SwarmVariable
particleLocalCount
Returns:Number of swarm particles within the current processes local space.
Return type:int
populate_using_layout(layout)[source]

This method uses the provided layout to populate the swarm with particles. Usually layouts add particles across the entire domain. Available layouts may be found in the swarm.layouts module. Note that Layouts can only currently be used on empty swarms. Also note that all numpy arrays associated with swarm variables must be deleted before a layout can be applied.

Parameters:layout (underworld.swarm.layouts._ParticleLayoutAbstract) – The layout which determines where particles are created and added.

Example

>>> # first we need a mesh
>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> # create swarm
>>> swarm = uw.swarm.Swarm(mesh)
>>> # add populate
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellGaussLayout(swarm,2))
stateId
Returns:Swarm state identifier. This is incremented whenever the swarm is modified.
Return type:int
variables
Returns:List of swarm variables associated with this swarm.
Return type:list
class underworld.swarm.SwarmVariable(swarm, dataType, count, writeable=True, **kwargs)[source]

Bases: underworld._stgermain.StgClass, underworld.function._function.Function

The SwarmVariable class allows users to add data to swarm particles. The data
can be of type “char”, “short”, “int”, “long, “float” or “double”.

Note that the swarm allocates one block of contiguous memory for all the particles. The per particle variable datums is then interlaced across this memory block.

The recommended practise is to add all swarm variables before populating the swarm to avoid costly reallocations.

Swarm variables should be added via the add_variable swarm method.

Parameters:
  • swarm (underworld.swarm.Swarm) – The swarm of particles for which we wish to add the variable
  • dataType (str) – The data type for the variable. Available types are “char”, “short”, “int”, “long”, “float” or “double”.
  • count (unsigned) – The number of values to be stored for each particle.
  • writeable (bool) – Signifies if the variable should be writeable.
copy(deepcopy=False)[source]

This method returns a copy of the swarmvariable.

Parameters:deepcopy (bool) – If True, the variable’s data is also copied into new variable.
Returns:The swarm variable copy.
Return type:underworld.swarm.SwarmVariable

Example

>>> mesh = uw.mesh.FeMesh_Cartesian()
>>> swarm = uw.swarm.Swarm(mesh)
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellGaussLayout(swarm, 2))
>>> svar = swarm.add_variable("double", 1)
>>> svar.data[:] = 1.23456
>>> svarCopy = svar.copy()
>>> svarCopy.swarm == svar.swarm
True
>>> svarCopy.dataType == svar.dataType
True
>>> import numpy as np
>>> np.allclose(svar.data,svarCopy.data)
False
>>> svarCopy2 = svar.copy(deepcopy=True)
>>> np.allclose(svar.data,svarCopy2.data)
True
count
Returns:Number of data items for this variable stored on each particle.
Return type:int
data
Returns:Numpy proxy array to underlying variable data. Note that the returned array is a proxy for all the local particle data. As numpy arrays are simply proxys to the underlying memory structures, no data copying is required.
Return type:numpy.ndarray

Example

>>> # create mesh
>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> # create empty swarm
>>> swarm = uw.swarm.Swarm(mesh)
>>> # add a variable
>>> svar = swarm.add_variable("int",1)
>>> # add particles
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellGaussLayout(swarm,2))
>>> swarm.particleLocalCount
1024
>>> len(svar.data)  # should be the same as particle local count
1024
>>> swarm.owningCell.data  # check particle owning cells/elements.
array([[  0],
       [  0],
       [  0],
       ...,
       [255],
       [255],
       [255]], dtype=int32)
>>> # particle coords
>>> swarm.particleCoordinates.data[0]
array([ 0.0132078,  0.0132078])
>>> # move the particle
>>> with swarm.deform_swarm():
...     swarm.particleCoordinates.data[0] = [0.2,0.2]
>>> swarm.particleCoordinates.data[0]
array([ 0.2,  0.2])
dataType
Returns:Data type for variable. Supported types are ‘char’, ‘short’, ‘int’, ‘long’, ‘float’ and ‘double’.
Return type:str
data_shadow
Returns:Numpy proxy array to underlying variable shadow data.
Return type:numpy.ndarray

Example

Refer to example provided for ‘data’ property(/method).

load(filename, collective=False)[source]

Load the swarm variable from disk. This must be called after the swarm.load().

Parameters:
  • filename (str) – The filename for the saved file. Relative or absolute paths may be used, but all directories must exist.
  • collective (bool) – If True, variable is loaded MPI collective. This is usually faster, but currently is problematic for passive swarms which may not have representation on all processes.

Notes

This method must be called collectively by all processes.

Example

Refer to example provided for ‘save’ method.

save(filename, collective=False, swarmHandle=None)[source]

Save the swarm variable to disk.

Parameters:
  • filename (str) – The filename for the saved file. Relative or absolute paths may be used, but all directories must exist.
  • collective (bool) – If True, variable is saved MPI collective. This is usually faster, but currently is problematic for passive type swarms which may not have representation on all processes.
  • swarmHandle (underworld.utils.SavedFileData) – The saved swarm file handle. If provided, a link is created within the swarmvariable file to this saved swarm file. Optional.
Returns:

Data object relating to saved file. This only needs to be retained if you wish to create XDMF files and can be ignored otherwise.

Return type:

underworld.utils.SavedFileData

Notes

This method must be called collectively by all processes.

Example

First create the swarm, populate, then add a variable:

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> swarm = uw.swarm.Swarm(mesh)
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellGaussLayout(swarm,2))
>>> svar = swarm.add_variable("int",1)

Write something to variable

>>> import numpy as np
>>> svar.data[:,0] = np.arange(swarm.particleLocalCount)

Save to a file:

>>> ignoreMe = swarm.save("saved_swarm.h5")
>>> ignoreMe = svar.save("saved_swarm_variable.h5")

Now let’s try and reload. First create a new swarm and swarm variable, and then load both:

>>> clone_swarm = uw.swarm.Swarm(mesh)
>>> clone_svar = clone_swarm.add_variable("int",1)
>>> clone_swarm.load("saved_swarm.h5")
>>> clone_svar.load("saved_swarm_variable.h5")

Now check for equality:

>>> import numpy as np
>>> np.allclose(svar.data,clone_svar.data)
True
>>> # clean up:
>>> if uw.mpi.rank == 0:
...     import os;
...     os.remove( "saved_swarm.h5" )
...     os.remove( "saved_swarm_variable.h5" )
swarm
Returns:The swarm this variable belongs to.
Return type:underworld.swarm.Swarm
xdmf(filename, varSavedData, varname, swarmSavedData, swarmname, modeltime=0.0)[source]

Creates an xdmf file, filename, associating the varSavedData file on the swarmSavedData file

Notes

xdmf contain 2 files: an .xml and a .h5 file. See http://www.xdmf.org/index.php/Main_Page This method only needs to be called by the master process, all other processes return quietly.

Parameters:
  • filename (str) – The output path to write the xdmf file. Relative or absolute paths may be used, but all directories must exist.
  • varname (str) – The xdmf name to give the swarmVariable
  • swarmname (str) – The xdmf name to give the swarm
  • swarmSavedData (underworld.utils.SaveFileData) – Handler returned for saving a swarm. underworld.swarm.Swarm.save(xxx)
  • varSavedData (underworld.utils.SavedFileData) – Handler returned from saving a SwarmVariable. underworld.swarm.SwarmVariable.save(xxx)
  • modeltime (float (default 0.0)) – The time recorded in the xdmf output file

Example

First create the swarm and add a variable:

>>> mesh = uw.mesh.FeMesh_Cartesian( elementType='Q1/dQ0', elementRes=(16,16), minCoord=(0.,0.), maxCoord=(1.,1.) )
>>> swarm = uw.swarm.Swarm( mesh=mesh )
>>> swarmLayout = uw.swarm.layouts.PerCellGaussLayout(swarm,2)
>>> swarm.populate_using_layout( layout=swarmLayout )
>>> swarmVar = swarm.add_variable( dataType="int", count=1 )

Write something to variable

>>> import numpy as np
>>> swarmVar.data[:,0] = np.arange(swarmVar.data.shape[0])

Save mesh and var to a file:

>>> swarmDat = swarm.save("saved_swarm.h5")
>>> swarmVarDat = swarmVar.save("saved_swarmvariable.h5")

Now let’s create the xdmf file

>>> swarmVar.xdmf("TESTxdmf", swarmVarDat, "var1", swarmDat, "MrSwarm" )

Does file exist?

>>> import os
>>> if uw.mpi.rank == 0: os.path.isfile("TESTxdmf.xdmf")
True
>>> # clean up:
>>> if uw.mpi.rank == 0:
...     import os;
...     os.remove( "saved_swarm.h5" )
...     os.remove( "saved_swarmvariable.h5" )
...     os.remove( "TESTxdmf.xdmf" )
class underworld.swarm.VoronoiIntegrationSwarm(swarm, **kwargs)[source]

Bases: underworld.swarm._integration_swarm.IntegrationSwarm, underworld.function._function.FunctionInput

Class for an IntegrationSwarm that maps to another Swarm

Note that this swarm can act as a function input. In this capacity, the fundamental function input type is the FEMCoordinate (ie, the particle local coordinate, the owning mesh, and the owning element). This input can be reduced to the global coordinate when returned within python. The FEMCoordinate particle representation is useful when deforming a mesh, as it is possible to deform the mesh, and then use the FEMCoordinate to reset the particles within the moved mesh.

Parameters:swarm (underworld.swarm.Swarm) – The PIC integration swarm maps to this user provided swarm.

Example

This simple example checks that the true global coordiante, and that derived from the local coordinate, are close to equal. Note that the VoronoiIntegrationSwarm uses a voronoi centroid algorithm so we do not expect particle to exactly coincide.

>>> import underworld as uw
>>> import numpy as np
>>> mesh = uw.mesh.FeMesh_Cartesian()
>>> swarm = uw.swarm.Swarm(mesh)
>>> swarm.populate_using_layout(uw.swarm.layouts.PerCellGaussLayout(swarm,4))
>>> vswarm = uw.swarm.VoronoiIntegrationSwarm(swarm)
>>> vswarm.repopulate()
>>> np.allclose(swarm.particleCoordinates.data, uw.function.input().evaluate(vswarm),atol=1e-1)
True
repopulate(weights_calculator=None)[source]

This method repopulates the voronoi swarm using the provided global swarm. The weights are also recalculated.

Parameters:weights_calculator (underworld.swarm.Weights) – The weights calculator for the Voronoi swarm. If none is provided, a default DVCWeights calculator is used.