Mirheo coordinator

The coordinator class stitches together data containers, Particle Vectors, and all the handlers, and provides functions to manipulate the system components.

One and only one instance of this class should be created in the beginning of any simulation setup.

Note

Creating the coordinator will internally call MPI_Init() function, and its destruction will call MPI_Finalize(). Therefore if using a mpi4py Python module, it should be imported in the following way:

import  mpi4py
mpi4py.rc(initialize=False, finalize=False)
from mpi4py import MPI
class Mirheo

Main coordination class, should only be one instance at a time

Methods

__init__(nranks: int3, domain: float3, dt: float, log_filename: str = 'log', debug_level: int = 3, checkpoint_every: int = 0, checkpoint_folder: str = 'restart/', checkpoint_mode: str = 'PingPong', cuda_aware_mpi: bool = False, no_splash: bool = False, comm_ptr: int = 0) → None

Create the Mirheo coordinator.

Warning

Debug level determines the amount of output produced by each of the simulation processes:

  1. silent: no log output
  2. only report fatal errors
  3. report serious errors
  4. report important information steps of simulation and warnings (this is the default level)
  5. report not critical information
  6. report some debug information
  7. report more debug
  8. report all the debug
  9. force flushing to the file after each message

Debug levels above 4 or 5 may significanlty increase the runtime, they are only recommended to debug errors. Flushing increases the runtime yet more, but it is required in order not to lose any messages in case of abnormal program abort.

Parameters:
  • nranks – number of MPI simulation tasks per axis: x,y,z. If postprocess is enabled, the same number of the postprocess tasks will be running
  • domain – size of the simulation domain in x,y,z. Periodic boundary conditions are applied at the domain boundaries. The domain will be split in equal chunks between the MPI ranks. The largest chunk size that a single MPI rank can have depends on the total number of particles, handlers and hardware, and is typically about \(120^3 - 200^3\).
  • dt – timestep of the simulation
  • log_filename – prefix of the log files that will be created. Logging is implemented in the form of one file per MPI rank, so in the simulation folder NP files with names log_00000.log, log_00001.log, … will be created, where NP is the total number of MPI ranks. Each MPI task (including postprocess) writes messages about itself into his own log file, and the combined log may be created by merging all the individual ones and sorting with respect to time. If this parameter is set to ‘stdout’ or ‘stderr’ standard output or standard error streams will be used instead of the file, however, there is no guarantee that messages from different ranks are synchronized.
  • debug_level – Debug level from 0 to 8, see above.
  • checkpoint_every – save state of the simulation components (particle vectors and handlers like integrators, plugins, etc.)
  • checkpoint_folder – folder where the checkpoint files will reside
  • checkpoint_mode – set to “PingPong” to keep only the last 2 checkpoint states; set to “Incremental” to keep all checkpoint states.
  • cuda_aware_mpi – enable CUDA Aware MPI. The MPI library must support that feature, otherwise it may fail.
  • no_splash – don’t display the splash screen when at the start-up.
  • comm_ptr – pointer to communicator. By default MPI_COMM_WORLD will be used
applyObjectBelongingChecker(checker: ObjectBelongingChecker, pv: ParticleVector, correct_every: int = 0, inside: str = '', outside: str = '') → ParticleVector

Apply the checker to the given particle vector. One and only one of the options inside or outside has to be specified.

Parameters:
  • checker – instance of BelongingChecker
  • pvParticleVector that will be split (source PV)
  • inside – if specified and not “none”, a new ParticleVector with name inside will be returned, that will keep the inner particles of the source PV. If set to “none”, None object will be returned. In any case, the source PV will only contain the outer particles
  • outside – if specified and not “none”, a new ParticleVector with name outside will be returned, that will keep the outer particles of the source PV. If set to “none”, None object will be returned. In any case, the source PV will only contain the inner particles
  • correct_every – If greater than zero, perform correction every this many time-steps. Correction will move e.g. inner particles of outer PV to the :inner PV and viceversa. If one of the PVs was defined as “none”, the ‘wrong’ particles will be just removed.
Returns:

New ParticleVector or None depending on inside and outside options

computeVolumeInsideWalls(walls: List[Wall], nSamplesPerRank: int = 100000) → float

Compute the volume inside the given walls in the whole domain (negative values are the ‘inside’ of the simulation). The computation is made via simple Monte-Carlo.

Parameters:
  • walls – sdf based walls
  • nSamplesPerRank – number of Monte-Carlo samples used per rank
dumpWalls2XDMF(walls: List[Wall], h: float3, filename: str = 'xdmf/wall') → None

Write Signed Distance Function for the intersection of the provided walls (negative values are the ‘inside’ of the simulation)

Parameters:
  • walls – list of walls to dump; the output sdf will be the union of all walls inside
  • h – cell-size of the resulting grid
  • filename – base filename output, will create to files filename.xmf and filename.h5
getState(self: Mirheo) → MirState

Return mirheo state

isComputeTask(self: Mirheo) → bool

Returns True if the current rank is a simulation task and False if it is a postrprocess task

isMasterTask(self: Mirheo) → bool

Returns True if the current rank is the root

log_compile_options(self: Mirheo) → None

output compile times options in the log

makeFrozenRigidParticles(checker: ObjectBelongingChecker, shape: ObjectVector, icShape: InitialConditions, interactions: List[Interaction], integrator: Integrator, number_density: float, nsteps: int = 1000) → ParticleVector

Create particles frozen inside object.

Note

A separate simulation will be run for every call to this function, which may take certain amount of time. If you want to save time, consider using restarting mechanism instead

Parameters:
  • checker – object belonging checker
  • shape – object vector describing the shape of the rigid object
  • icShape – initial conditions for shape
  • interactions – list of Interaction that will be used to construct the equilibrium particles distribution
  • integrator – this Integrator will be used to construct the equilibrium particles distribution
  • number_density – target particle number density
  • nsteps – run this many steps to achieve equilibrium
Returns:

New ParticleVector that will contain particles that are close to the wall boundary, but still inside the wall.

makeFrozenWallParticles(pvName: str, walls: List[Wall], interactions: List[Interaction], integrator: Integrator, number_density: float, nsteps: int = 1000) → ParticleVector

Create particles frozen inside the walls.

Note

A separate simulation will be run for every call to this function, which may take certain amount of time. If you want to save time, consider using restarting mechanism instead

Parameters:
  • pvName – name of the created particle vector
  • walls – array of instances of Wall for which the frozen particles will be generated
  • interactions – list of Interaction that will be used to construct the equilibrium particles distribution
  • integrator – this Integrator will be used to construct the equilibrium particles distribution
  • number_density – target particle number density
  • nsteps – run this many steps to achieve equilibrium
Returns:

New ParticleVector that will contain particles that are close to the wall boundary, but still inside the wall.

registerBouncer(bouncer: Bouncer) → None

Register Object Bouncer

Parameters:bouncer – the Bouncer to register
registerIntegrator(integrator: Integrator) → None

Register an Integrator to the coordinator

Parameters:integrator – the Integrator to register
registerInteraction(interaction: Interaction) → None

Register an Interaction to the coordinator

Parameters:interaction – the Interaction to register
registerObjectBelongingChecker(checker: ObjectBelongingChecker, ov: ObjectVector) → None

Register Object Belonging Checker

Parameters:
registerParticleVector(pv: ParticleVector, ic: InitialConditions = None) → None

Register particle vector

Parameters:
registerPlugins(arg0: SimulationPlugin, arg1: PostprocessPlugin) → None

Register Plugins

registerWall(wall: Wall, check_every: int = 0) → None

Register a Wall.

Parameters:
  • wall – the Wall to register
  • check_every – if positive, check every this many time steps if particles penetrate the walls
restart(folder: str = 'restart/') → None

Restart the simulation. This function should typically be called just before running the simulation. It will read the state of all previously registered instances of ParticleVector, Interaction, etc. If the folder contains no checkpoint file required for one of those, an error occur.

Warning

Certain Plugins may not implement restarting mechanism and will not restart correctly. Please check the documentation for the plugins.

Parameters:folder – folder with the checkpoint files
run(niters: int) → None

Advance the system for a given amount of time steps.

Parameters:niters – number of time steps to advance
save_dependency_graph_graphml(fname: str, current: bool = True) → None

Exports GraphML file with task graph for the current simulation time-step

Parameters:
  • fname – the output filename (without extension)
  • current – if True, save the current non empty tasks; else, save all tasks that can exist in a simulation

Warning

if current is set to True, this must be called after _mirheo.Mirheo.run().

setBouncer(bouncer: Bouncer, ov: ObjectVector, pv: ParticleVector) → None

Assign a Bouncer between an ObjectVector and a ParticleVector.

Parameters:
setIntegrator(integrator: Integrator, pv: ParticleVector) → None

Set a specific Integrator to a given ParticleVector

Parameters:
setInteraction(interaction: Interaction, pv1: ParticleVector, pv2: ParticleVector) → None

Forces between two instances of ParticleVector (they can be the same) will be computed according to the defined interaction.

Parameters:
setWall(wall: Wall, pv: ParticleVector, maximum_part_travel: float = 0.25) → None

Assign a Wall bouncer to a given ParticleVector. The current implementation does not support ObjectVector.

Parameters:
  • wall – the Wall surface which will bounce the particles
  • pv – the ParticleVector to be bounced
  • maximum_part_travel – maximum distance that one particle travels in one time step. this should be as small as possible for performance reasons but large enough for correctness
start_profiler(self: Mirheo) → None

Tells nvprof to start recording timeline

stop_profiler(self: Mirheo) → None

Tells nvprof to stop recording timeline