Mirheo

The main coordinator class

API

class Mirheo

Coordinator class for a full simulation.

Manages and splits work between Simulation and Postprocess ranks.

Public Functions

Mirheo(int3 nranks3D, real3 globalDomainSize, LogInfo logInfo, CheckpointInfo checkpointInfo, real maxObjHalfLength, bool gpuAwareMPI = false)

Construct a Mirheo object using MPI_COMM_WORLD.

The product of

nranks3D must be equal to the number of available ranks (or hals if postprocess is used)
Note
MPI will be initialized internally. If this constructor is used, the destructor will also finalize MPI.
Parameters
  • nranks3D: Number of ranks along each cartesian direction.
  • globalDomainSize: The full domain dimensions in length units. Must be positive.
  • logInfo: Information about logging
  • checkpointInfo: Information about checkpoint
  • maxObjHalfLength: Half of the maximum length of all objects.
  • gpuAwareMPI: true to use RDMA (must be compile with a MPI version that supports it)

Mirheo(MPI_Comm comm, int3 nranks3D, real3 globalDomainSize, LogInfo logInfo, CheckpointInfo checkpointInfo, real maxObjHalfLength, bool gpuAwareMPI = false)

Construct a Mirheo object using a given communicator.

Note
MPI will be NOT be initialized. If this constructor is used, the destructor will NOT finalize MPI.

void restart(std::string folder = "restart/")

reset the internal state from a checkpoint folder

MPI_Comm getWorldComm() const

Return
the world communicator

bool isComputeTask() const

Return
true if the current rank is a Simulation rank

bool isMasterTask() const

Return
true if the current rank is the root (i.e. rank = 0)

bool isSimulationMasterTask() const

Return
true if the current rank is the root within the simulation communicator

bool isPostprocessMasterTask() const

Return
true if the current rank is the root within the postprocess communicator

void startProfiler()

start profiling for nvvp

void stopProfiler()

stop profiling for nvvp

void dumpDependencyGraphToGraphML(const std::string &fname, bool current) const

dump the task dependency of the simulation in graphML format.

Parameters
  • fname: The file name to dump the graph to (without extension).
  • current: if true, will only dump the current tasks; otherwise, will dump all possible ones.

void run(MirState::StepType niters, real dt)

advance the system for a given number of time steps

Parameters
  • niters: number of interations
  • dt: time step duration

void registerParticleVector(std::shared_ptr<ParticleVector> pv, std::shared_ptr<InitialConditions> ic)

register a ParticleVector in the simulation and initialize it with the gien InitialConditions.

Parameters

void registerInteraction(std::shared_ptr<Interaction> interaction)

register an Interaction

See
setInteraction().
Parameters

void registerIntegrator(std::shared_ptr<Integrator> integrator)

register an Integrator

See
setIntegrator().
Parameters

void registerWall(std::shared_ptr<Wall> wall, int checkEvery = 0)

register a Wall

Parameters
  • wall: The Wall to register
  • checkEvery: The particles that will bounce against this wall will be checked (inside/outside log info) every this number of time steps. 0 means no check.

void registerBouncer(std::shared_ptr<Bouncer> bouncer)

register a Bouncer

See
setBouncer().
Parameters
  • bouncer: the Bouncer to register.

void registerPlugins(std::shared_ptr<SimulationPlugin> simPlugin, std::shared_ptr<PostprocessPlugin> postPlugin)

register a SimulationPlugin

Parameters
  • simPlugin: the SimulationPlugin to register (only relevant if the current rank is a compute task).
  • postPlugin: the PostprocessPlugin to register (only relevant if the current rank is a postprocess task).

void registerPlugins(const PairPlugin &plugins)

More generic version of registerPlugins()

void registerObjectBelongingChecker(std::shared_ptr<ObjectBelongingChecker> checker, ObjectVector *ov)

register a ObjectBelongingChecker

See
applyObjectBelongingChecker()
Parameters

void deregisterIntegrator(Integrator *integrator)

deregister an Integrator

See
registerIntegrator().
Parameters

void deregisterPlugins(SimulationPlugin *simPlugin, PostprocessPlugin *postPlugin)

deregister a Plugin

Parameters
  • simPlugin: the SimulationPlugin to deregister (only relevant if the current rank is a compute task).
  • postPlugin: the PostprocessPlugin to deregister (only relevant if the current rank is a postprocess task).

void setIntegrator(Integrator *integrator, ParticleVector *pv)

Assign a registered Integrator to a registered ParticleVector.

Parameters
  • integrator: The registered integrator (will die if it was not registered)
  • pv: The registered ParticleVector (will die if it was not registered)

void setInteraction(Interaction *interaction, ParticleVector *pv1, ParticleVector *pv2)

Assign two registered Interaction to two registered ParticleVector objects.

This was designed to handle PairwiseInteraction, which needs up to two

ParticleVector. For self interaction cases (such as MembraneInteraction), pv1 and pv2 must be the same.
Parameters
  • interaction: The registered interaction (will die if it is not registered)
  • pv1: The first registered ParticleVector (will die if it is not registered)
  • pv2: The second registered ParticleVector (will die if it is not registered)

void setBouncer(Bouncer *bouncer, ObjectVector *ov, ParticleVector *pv)

Assign a registered Bouncer to registered ObjectVector and ParticleVector.

Parameters
  • bouncer: The registered bouncer (will die if it is not registered)
  • ov: The registered ObjectVector that contains the surface to bounce on (will die if it is not registered)
  • pv: The registered ParticleVector to bounce (will die if it is not registered)

void setWallBounce(Wall *wall, ParticleVector *pv, real maximumPartTravel = 0.25f)

Set a registered ParticleVector to bounce on a registered Wall.

Parameters
  • wall: The registered wall (will die if it is not registered)
  • pv: The registered ParticleVector (will die if it is not registered)
  • maximumPartTravel: Performance parameter. See Wall for more information.

MirState *getState()

Return
the global state of the system

const MirState *getState() const

Return
the global state of the system (const version)

Simulation *getSimulation()

Return
the Simulation object; nullptr on postprocess tasks.

const Simulation *getSimulation() const

see getSimulation(); const version

std::shared_ptr<MirState> getMirState()

see getMirState(); shared_ptr version

void dumpWalls2XDMF(std::vector<std::shared_ptr<Wall>> walls, real3 h, const std::string &filename)

Compute the SDF field from the given walls and dump it to a file in xmf+h5 format.

Parameters
  • walls: List of Wall objects. The union of these walls will be dumped.
  • h: The grid spacing
  • filename: The base name of the dumped files (without extension)

double computeVolumeInsideWalls(std::vector<std::shared_ptr<Wall>> walls, long nSamplesPerRank = 100000)

Compute the volume inside the geometry formed by the given walls with simple Monte-Carlo integration.

Return
The Monte-Carlo estimate of the volume
Parameters
  • walls: List of Wall objects. The union of these walls form the geometry.
  • nSamplesPerRank: The number of Monte-Carlo samples per rank

std::shared_ptr<ParticleVector> makeFrozenWallParticles(std::string pvName, std::vector<std::shared_ptr<Wall>> walls, std::vector<std::shared_ptr<Interaction>> interactions, std::shared_ptr<Integrator> integrator, real numDensity, real mass, real dt, int nsteps)

Create a layer of frozen particles inside the given walls.

This will run a simulation of “bulk” particles and select the particles that are inside the effective cut-off radius of the given list of interactions.

Return
The frozen particles
Parameters
  • pvName: The name of the frozen ParticleVector that will be created
  • walls: The list of registered walls that need frozen particles
  • interactions: List of interactions (not necessarily registered) that will be used to equilibrate the particles
  • integrator: Integrator object used to equilibrate the particles
  • numDensity: The number density used to initialize the particles
  • mass: The mass of one particle
  • dt: Equilibration time step
  • nsteps: Number of equilibration steps

std::shared_ptr<ParticleVector> makeFrozenRigidParticles(std::shared_ptr<ObjectBelongingChecker> checker, std::shared_ptr<ObjectVector> shape, std::shared_ptr<InitialConditions> icShape, std::vector<std::shared_ptr<Interaction>> interactions, std::shared_ptr<Integrator> integrator, real numDensity, real mass, real dt, int nsteps)

Create frozen particles inside the given objects.

This will run a simulation of “bulk” particles and select the particles that are inside

shape.
Return
The frozen particles, with name “inside_” + name of shape
Note
For now, the output ParticleVector has mass 1.0.
Parameters
  • checker: The ObjectBelongingChecker to split inside particles
  • shape: The ObjectVector that will be used to define inside particles
  • icShape: The InitialConditions object used to set the objects positions
  • interactions: List of interactions (not necessarily registered) that will be used to equilibrate the particles
  • integrator: Integrator object used to equilibrate the particles
  • numDensity: The number density used to initialize the particles
  • mass: The mass of one particle
  • dt: Equilibration time step
  • nsteps: Number of equilibration steps

std::shared_ptr<ParticleVector> applyObjectBelongingChecker(ObjectBelongingChecker *checker, ParticleVector *pv, int checkEvery, std::string inside = "", std::string outside = "")

Enable a registered ObjectBelongingChecker to split particles of a registered ParticleVector.

inside or outside can take the reserved value “none”, in which case the corresponding particles will be deleted. Furthermore, exactly one of inside and outside must be the same as pv.

Parameters
  • checker: The ObjectBelongingChecker (will die if it is not registered)
  • pv: The registered ParticleVector that must be split (will die if it is not registered)
  • checkEvery: The particle split will be performed every this amount of time steps.
  • inside: Name of the ParticleVector that will contain the particles of pv that are inside the objects. See below for more information.
  • outside: Name of the ParticleVector that will contain the particles of pv that are outside the objects. See below for more information.

If inside or outside has the name of a ParticleVector that is not registered, this call will create an empty ParticleVector with the given name and register it in the Simulation. Otherwise the already registered ParticleVector will be used.

void logCompileOptions() const

print the list of all compile options and their current value in the logs