hoomd.comm

Overview

hoomd.comm.barrier Perform a MPI barrier synchronization across all ranks in the partition.
hoomd.comm.barrier_all Perform a MPI barrier synchronization across the whole MPI run.
hoomd.comm.decomposition Set the domain decomposition.
hoomd.comm.get_num_ranks Get the number of ranks in this partition.
hoomd.comm.get_partition Get the current partition index.
hoomd.comm.get_rank Get the current rank.

Details

MPI communication interface

Use methods in this module to query the number of MPI ranks, the current rank, etc…

hoomd.comm.barrier()

Perform a MPI barrier synchronization across all ranks in the partition.

Note

Does nothing in in non-MPI builds.

hoomd.comm.barrier_all()

Perform a MPI barrier synchronization across the whole MPI run.

Note

Does nothing in in non-MPI builds.

class hoomd.comm.decomposition(x=None, y=None, z=None, nx=None, ny=None, nz=None)

Set the domain decomposition.

Parameters:
  • x (list) – First nx-1 fractional domain widths (if nx is None)
  • y (list) – First ny-1 fractional domain widths (if ny is None)
  • z (list) – First nz-1 fractional domain widths (if nz is None)
  • nx (int) – Number of processors to uniformly space in x dimension (if x is None)
  • ny (int) – Number of processors to uniformly space in y dimension (if y is None)
  • nz (int) – Number of processors to uniformly space in z dimension (if z is None)

A single domain decomposition is defined for the simulation. A standard domain decomposition divides the simulation box into equal volumes along the Cartesian axes while minimizing the surface area between domains. This works well for systems where particles are uniformly distributed and there is equal computational load for each domain, and is the default behavior in HOOMD-blue. If no decomposition is specified for an MPI run, a uniform decomposition is automatically constructed on initialization.

In simulations with density gradients, such as a vapor-liquid interface, there can be a considerable imbalance of particles between different ranks. The simulation time then becomes limited by the slowest processor. It may then be advantageous in certain systems to create domains of unequal volume, for example, by increasing the volume of less dense regions of the simulation box in order to balance the number of particles.

The decomposition command allows the user to control the geometry and positions of the decomposition. The fractional width of the first \(n_i - 1\) domains is specified along each dimension, where \(n_i\) is the number of ranks desired along dimension \(i\). If no cut planes are specified, then a uniform spacing is assumed. The number of domains with uniform spacing can also be specified. If the desired decomposition is not commensurate with the number of ranks available (for example, a 3x3x3 decomposition when only 8 ranks are available), then a default uniform spacing is chosen. For the best control, the user should specify the number of ranks in each dimension even if uniform spacing is desired.

decomposition can only be called before the system is initialized, at which point the particles are decomposed. An error is raised if the system is already initialized.

The decomposition can be adjusted dynamically if the best static decomposition is not known, or the system composition is changing dynamically. For this associated command, see update.balance().

Priority is always given to specified arguments over the command line arguments. If one of these is not set but a command line option is, then the command line option is used. Otherwise, a default decomposition is chosen.

Examples:

comm.decomposition(x=0.4, ny=2, nz=2)
comm.decomposition(nx=2, y=0.8, z=[0.2,0.3])

Warning

The decomposition command will override specified command line options.

Warning

This command must be invoked before the system is initialized because particles are decomposed at this time.

Note

The domain size cannot be chosen arbitrarily small. There are restrictions placed on the decomposition by the ghost layer width set by the pair potentials. An error will be raised at run time if the ghost layer width exceeds half the shortest domain size.

Warning

Both fractional widths and the number of processors cannot be set simultaneously, and an error will be raised if both are set.

set_params(x=None, y=None, z=None, nx=None, ny=None, nz=None)

Set parameters for the decomposition before initialization.

Parameters:
  • x (list) – First nx-1 fractional domain widths (if nx is None)
  • y (list) – First ny-1 fractional domain widths (if ny is None)
  • z (list) – First nz-1 fractional domain widths (if nz is None)
  • nx (int) – Number of processors to uniformly space in x dimension (if x is None)
  • ny (int) – Number of processors to uniformly space in y dimension (if y is None)
  • nz (int) – Number of processors to uniformly space in z dimension (if z is None)

Examples:

decomposition.set_params(x=[0.2])
decomposition.set_params(nx=1, y=[0.3,0.4], nz=2)
hoomd.comm.get_num_ranks()

Get the number of ranks in this partition.

Returns:The number of MPI ranks in this partition.

Note

Returns 1 in non-mpi builds.

hoomd.comm.get_partition()

Get the current partition index.

Returns:Index of the current partition.

Note

Always returns 0 in non-mpi builds.

hoomd.comm.get_rank()

Get the current rank.

Returns:Index of the current rank in this partition.

Note

Always returns 0 in non-mpi builds.