List of plugins

Dump Plugins

These plugins do not modify the state of the simulation. They can be used to dump selected parts of the state of the simulation to the disk.

class Average3D : public mirheo::SimulationPlugin

Average particles quantities into spacial bins over a cartesian grid and average it over time.

Useful to compute e.g. velocity or density profiles. The number density is always computed inside each bin, as it is used to compute the averages. Other quantities must be specified by giving the channel names.

This plugin should be used with UniformCartesianDumper on the postprocessing side.

Cannot be used with multiple invocations of Mirheo.run.

Subclassed by mirheo::AverageRelative3D

Public Types

enum ChannelType

Specify the form of a channel data.

Values:

Scalar
Vector_real3
Vector_real4
Tensor6
None

Public Functions

Average3D(const MirState *state, std::string name, std::vector<std::string> pvNames, std::vector<std::string> channelNames, int sampleEvery, int dumpEvery, real3 binSize)

Create an Average3D object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: The list of names of the ParticleVector that will be used when averaging.
  • channelNames: The list of particle data channels to average. Will die if the channel does not exist.
  • sampleEvery: Compute spatial averages every this number of time steps.
  • dumpEvery: Compute time averages and send to the postprocess side every this number of time steps.
  • binSize: Size of one spatial bin along the three axes.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

struct HostChannelsInfo

A helper structure that contains grid average info for all required channels.

Public Members

int n

The number of channels (excluding number density).

std::vector<std::string> names

List of channel names.

PinnedBuffer<ChannelType> types

List of channel data forms.

PinnedBuffer<real *> averagePtrs

List of averages of each channel.

PinnedBuffer<real *> dataPtrs

List of data to average, for each channel.

std::vector<DeviceBuffer<real>> average

data container for the averages, for each channel.

class AverageRelative3D : public mirheo::Average3D

Perform the same task as AverageRelative3D on a grid that moves relatively to a given object’s center of mass in a RigidObjectVector.

Cannot be used with multiple invocations of Mirheo.run.

Public Functions

AverageRelative3D(const MirState *state, std::string name, std::vector<std::string> pvNames, std::vector<std::string> channelNames, int sampleEvery, int dumpEvery, real3 binSize, std::string relativeOVname, int relativeID)

Create an AverageRelative3D object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: The list of names of the ParticleVector that will be used when averaging.
  • channelNames: The list of particle data channels to average. Will die if the channel does not exist.
  • sampleEvery: Compute spatial averages every this number of time steps.
  • dumpEvery: Compute time averages and send to the postprocess side every this number of time steps.
  • binSize: Size of one spatial bin along the three axes.
  • relativeOVname: Name of the RigidObjectVector that contains the reference object.
  • relativeID: Index of the reference object within the RigidObjectVector.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class UniformCartesianDumper : public mirheo::PostprocessPlugin

Postprocessing side of Average3D or AverageRelative3D.

Dump uniform grid data to xmf + hdf5 format.

Public Functions

UniformCartesianDumper(std::string name, std::string path)

Create a UniformCartesianDumper.

Parameters
  • name: The name of the plugin.
  • path: The files will be dumped to pathXXXXX.[xmf,h5], where XXXXX is the time stamp.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

XDMF::Channel getChannelOrDie(std::string chname) const

Get the average channel data.

Return
The channel data.
Parameters
  • chname: The name of the channel.

std::vector<int> getLocalResolution() const

Get the grid size in the local domain.

Return
An array with 3 entries, contains the number of grid points along each direction.

class MeshPlugin : public mirheo::SimulationPlugin

Send mesh information of an object for dump to MeshDumper postprocess plugin.

Public Functions

MeshPlugin(const MirState *state, std::string name, std::string ovName, int dumpEvery)

Create a MeshPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • ovName: The name of the ObjectVector that has a mesh to dump.
  • dumpEvery: Will dump the mesh every this number of timesteps.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class MeshDumper : public mirheo::PostprocessPlugin

Postprocess side of MeshPlugin.

Receives mesh info and dump it to ply format.

Public Functions

MeshDumper(std::string name, std::string path)

Create a MeshDumper object.

Parameters
  • name: The name of the plugin.
  • path: The files will be dumped to path-XXXXX.ply.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void setup(const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the PostprocessPlugin.

This method must be called before any other function call.

Parameters
  • comm: Contains all postprocess ranks
  • interComm: used to communicate with the simulation ranks

class ParticleSenderPlugin : public mirheo::SimulationPlugin

Send particle data to ParticleDumperPlugin.

Subclassed by mirheo::ParticleWithMeshSenderPlugin, mirheo::ParticleWithPolylinesSenderPlugin

Public Functions

ParticleSenderPlugin(const MirState *state, std::string name, std::string pvName, int dumpEvery, const std::vector<std::string> &channelNames)

Create a ParticleSenderPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to dump.
  • dumpEvery: Send the data to the postprocess side every this number of steps.
  • channelNames: The list of channels to send, additionally to the default positions, velocities and global ids.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

Public Members

std::string pvName_

name of the ParticleVector to dump.

ParticleVector *pv_

pointer to the ParticleVector to dump.

std::vector<char> sendBuffer_

Buffer used to send the data to the postprocess side.

class ParticleDumperPlugin : public mirheo::PostprocessPlugin

Postprocess side of ParticleSenderPlugin.

Dump particles data to xmf + hdf5 format.

Subclassed by mirheo::ParticleWithMeshDumperPlugin, mirheo::ParticleWithPolylinesDumperPlugin

Public Functions

ParticleDumperPlugin(std::string name, std::string path)

Create a ParticleDumperPlugin object.

Parameters
  • name: The name of the plugin.
  • path: Particle data will be dumped to pathXXXXX.[xmf,h5].

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

class ParticleWithMeshSenderPlugin : public mirheo::ParticleSenderPlugin

Send particle data to ParticleWithMeshSenderPlugin.

Does the same as ParticleSenderPlugin with additional Mesh connectivity information. This is compatible only with ObjectVector.

Public Functions

ParticleWithMeshSenderPlugin(const MirState *state, std::string name, std::string pvName, int dumpEvery, const std::vector<std::string> &channelNames)

Create a ParticleWithMeshSenderPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to dump.
  • dumpEvery: Send the data to the postprocess side every this number of steps.
  • channelNames: The list of channels to send, additionally to the default positions, velocities and global ids.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

class ParticleWithMeshDumperPlugin : public mirheo::ParticleDumperPlugin

Postprocess side of ParticleWithMeshSenderPlugin.

Dump particles data with connectivity to xmf + hdf5 format.

Public Functions

ParticleWithMeshDumperPlugin(std::string name, std::string path)

Create a ParticleWithMeshDumperPlugin object.

Parameters
  • name: The name of the plugin.
  • path: Data will be dumped to pathXXXXX.[xmf,h5].

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

class XYZPlugin : public mirheo::SimulationPlugin

Send particle positions to XYZDumper.

Public Functions

XYZPlugin(const MirState *state, std::string name, std::string pvName, int dumpEvery)

Create a XYZPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to dump.
  • dumpEvery: Send the data to the postprocess side every this number of steps.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class XYZDumper : public mirheo::PostprocessPlugin

Postprocess side of XYZPlugin.

Dump the particle positions to simple .xyz format.

Public Functions

XYZDumper(std::string name, std::string path)

Create a XYZDumper object.

Parameters
  • name: The name of the plugin.
  • path: Data will be dumped to pathXXXXX.xyz.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void setup(const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the PostprocessPlugin.

This method must be called before any other function call.

Parameters
  • comm: Contains all postprocess ranks
  • interComm: used to communicate with the simulation ranks

Statistics and In-situ analysis Plugins

These plugins do not modify the state of the simulation. They are used to measure properties of the simulation that can be processed directly at runtime. Their output is generally much lighter than dump plugins. The prefered format is csv, to allow clean postprocessing from e.g. python.

class ExponentialMovingAveragePlugin : public mirheo::SimulationPlugin

Compute the exponential moving average (EMA) of the particles velocity.

Public Functions

ExponentialMovingAveragePlugin(const MirState *state, std::string name, std::string pvName, real alpha, std::string srcChannelName, std::string emaChannelName)

Create a ExponentialMovingAveragePlugin.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector.
  • alpha: Coefficient in [0,1] to perform the EMA.
  • srcChannelName: The name of the channel to average. Will fail if it does not exist.
  • emaChannelName: The name of the channel containing the EMA.

void beforeIntegration(cudaStream_t stream)

hook before integrating the particle vectors but after the forces are computed

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

class MsdPlugin : public mirheo::SimulationPlugin

Compute the mean squared distance (MSD) of a given ParticleVector.

The MSD is computed every dumpEvery steps on the time interval [startTime, endTime].

Each particle stores the total displacement from startTime. To compute this, it also stores its position at each step.

Public Functions

MsdPlugin(const MirState *state, std::string name, std::string pvName, MirState::TimeType startTime, MirState::TimeType endTime, int dumpEvery)

Create a MsdPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector from which to measure the MSD.
  • startTime: MSD will use this time as origin.
  • endTime: The MSD will be reported only on [startTime, endTime].
  • dumpEvery: Will send the MSD to the postprocess side every this number of steps, only during the valid time interval.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class MsdDumper : public mirheo::PostprocessPlugin

Postprocess side of MsdPlugin.

Dumps the VACF in a csv file.

Public Functions

MsdDumper(std::string name, std::string path)

Create a MsdDumper object.

Parameters
  • name: The name of the plugin.
  • path: The folder that will contain the vacf csv file.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void setup(const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the PostprocessPlugin.

This method must be called before any other function call.

Parameters
  • comm: Contains all postprocess ranks
  • interComm: used to communicate with the simulation ranks

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

class ObjStatsPlugin : public mirheo::SimulationPlugin

Send object information to ObjStatsDumper.

Used to track the center of mass, linear and angular velocities, orintation, forces and torques of an ObjectVector.

Public Functions

ObjStatsPlugin(const MirState *state, std::string name, std::string ovName, int dumpEvery)

Create a ObjStatsPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • ovName: The name of the ObjectVector to extract the information from.
  • dumpEvery: Send the information to the postprocess side every this number of steps.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class ObjStatsDumper : public mirheo::PostprocessPlugin

Postprocess side of ObjStatsPlugin.

Dump object information to a csv file.

Public Functions

ObjStatsDumper(std::string name, std::string filename)

Create a ObjStatsDumper object.

Parameters
  • name: The name of the plugin.
  • filename: The name of the csv file to dump to.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void setup(const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the PostprocessPlugin.

This method must be called before any other function call.

Parameters
  • comm: Contains all postprocess ranks
  • interComm: used to communicate with the simulation ranks

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

void checkpoint(MPI_Comm comm, const std::string &path, int checkPointId)

Save the state of the object on disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.
  • checkPointId: The id of the dump.

void restart(MPI_Comm comm, const std::string &path)

Load the state of the object from the disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.

class RdfPlugin : public mirheo::SimulationPlugin

Measure the radial distribution function (RDF) of a ParticleVector.

The RDF is estimated periodically from ensemble averages. See RdfDump for the I/O.

Public Functions

RdfPlugin(const MirState *state, std::string name, std::string pvName, real maxDist, int nbins, int computeEvery)

Create a RdfPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector from which to measure the RDF.
  • maxDist: The RDF will be measured on [0, maxDist].
  • nbins: The number of bins in the interval [0, maxDist].
  • computeEvery: The number of time steps between two RDF evaluations and dump.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class RdfDump : public mirheo::PostprocessPlugin

Postprocess side of RdfPlugin.

Dump the RDF to a csv file.

Public Functions

RdfDump(std::string name, std::string basename)

Create a RdfDump object.

Parameters
  • name: The name of the plugin.
  • basename: The RDF will be dumped to basenameXXXXX.csv.

void setup(const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the PostprocessPlugin.

This method must be called before any other function call.

Parameters
  • comm: Contains all postprocess ranks
  • interComm: used to communicate with the simulation ranks

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

class SimulationStats : public mirheo::SimulationPlugin

Collect global statistics of the simulation and send it to the postprocess ranks.

Compute total linear momentum and estimate of temperature. Furthermore, measures average wall time of time steps.

Public Functions

SimulationStats(const MirState *state, std::string name, int every, std::vector<std::string> pvNames)

Create a SimulationPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • every: Compute the statistics every this number of steps.
  • pvNames: List of names of the pvs to compute statistics from. If empty, will take all the pvs in the simulation.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class PostprocessStats : public mirheo::PostprocessPlugin

Dump the stats sent by SimulationStats to a csv file and to the console output.

Public Functions

PostprocessStats(std::string name, std::string filename = std::string())

Construct a PostprocessStats plugin.

Parameters
  • name: The name of the plugin.
  • filename: The csv file name that will be dumped.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void checkpoint(MPI_Comm comm, const std::string &path, int checkPointId)

Save the state of the object on disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.
  • checkPointId: The id of the dump.

void restart(MPI_Comm comm, const std::string &path)

Load the state of the object from the disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.

class VacfPlugin : public mirheo::SimulationPlugin

Compute the velocity autocorrelation function (VACF) of a given ParticleVector.

The VACF is computed every dumpEvery steps on the time interval [startTime, endTime].

Public Functions

VacfPlugin(const MirState *state, std::string name, std::string pvName, MirState::TimeType startTime, MirState::TimeType endTime, int dumpEvery)

Create a VacfPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector from which to measure the VACF.
  • startTime: VACF will use this time as origin.
  • endTime: The VACF will be reported only on [startTime, endTime].
  • dumpEvery: Will send the VACF to the postprocess side every this number of steps, only during the valid time interval.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class VacfDumper : public mirheo::PostprocessPlugin

Postprocess side of VacfPlugin.

Dumps the VACF in a csv file.

Public Functions

VacfDumper(std::string name, std::string path)

Create a VacfDumper object.

Parameters
  • name: The name of the plugin.
  • path: The folder that will contain the vacf csv file.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void setup(const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the PostprocessPlugin.

This method must be called before any other function call.

Parameters
  • comm: Contains all postprocess ranks
  • interComm: used to communicate with the simulation ranks

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

class VirialPressurePlugin : public mirheo::SimulationPlugin

Compute the pressure in a given region from the virial theorem and send it to the VirialPressureDumper.

Public Functions

VirialPressurePlugin(const MirState *state, std::string name, std::string pvName, ScalarFieldFunction func, real3 h, int dumpEvery)

Create a VirialPressurePlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to add the particles to.
  • func: The scalar field is negative in the region of interest and positive outside.
  • h: The grid size used to discretize the field.
  • dumpEvery: Will compute and send the pressure every this number of steps.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class VirialPressureDumper : public mirheo::PostprocessPlugin

Postprocess side of VirialPressurePlugin.

Recieves and dump the virial pressure.

Public Functions

VirialPressureDumper(std::string name, std::string path)

Create a VirialPressureDumper.

Parameters
  • name: The name of the plugin.
  • path: The csv file to which the data will be dumped.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void setup(const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the PostprocessPlugin.

This method must be called before any other function call.

Parameters
  • comm: Contains all postprocess ranks
  • interComm: used to communicate with the simulation ranks

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

class WallForceCollectorPlugin : public mirheo::SimulationPlugin

Compute the force exerted by particles on the walls.

It has two contributions:

  • Interactio forces with frozen particles
  • bounce-back.

Public Functions

WallForceCollectorPlugin(const MirState *state, std::string name, std::string wallName, std::string frozenPvName, int sampleEvery, int dumpEvery)

Create a WallForceCollectorPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • wallName: The name of the Wall to collect the forces from.
  • frozenPvName: The name of the frozen ParticleVector assigned to the wall.
  • sampleEvery: Compute forces every this number of steps, and average it in time.
  • dumpEvery: Send the average forces to the postprocessing side every this number of steps.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class WallForceDumperPlugin : public mirheo::PostprocessPlugin

Postprocess side of WallForceCollectorPlugin.

Dump the forces to a txt file.

Public Functions

WallForceDumperPlugin(std::string name, std::string filename, bool detailedDump)

Create a WallForceDumperPlugin.

Parameters
  • name: The name of the plugin.
  • filename: The file to dump the stats to.
  • detailedDump: If true, the file will contain the bounce contribution and particle interactions contributions instead of the sum.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

Modifier plugins

These plugins add more functionalities to the simulation.

class AddForcePlugin : public mirheo::SimulationPlugin

Add a constant force to every particle of a given ParticleVector at every time step.

The force is added at the beforeForce() stage.

Public Functions

AddForcePlugin(const MirState *state, const std::string &name, const std::string &pvName, real3 force)

Create a AddForcePlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to which the force should be applied.
  • force: The force to apply.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class AddFourRollMillForcePlugin : public mirheo::SimulationPlugin

Add the “four roll mill” force to every particle of a given ParticleVector at every time step.

The force is added at the beforeForce() stage. The force has the form f = A * (sin(x) cos(y), -cos(x) sin(y), 0) where A is the intensity of the force and x, y are scaled so that they cover one period over the domain.

Public Functions

AddFourRollMillForcePlugin(const MirState *state, const std::string &name, const std::string &pvName, real intensity)

Create a AddFourRollMillForcePlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to which the force should be applied.
  • intensity: The intensity A of the force to apply.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class AddReversePoiseuilleForcePlugin : public mirheo::SimulationPlugin

Add a constant force to every particle of a given ParticleVector at every time step.

The force is added at the beforeForce() stage. The force has a reversed sign in half of the domain.

Public Functions

AddReversePoiseuilleForcePlugin(const MirState *state, const std::string &name, const std::string &pvName, real3 force, char flipDirection)

Create a AddReversePoiseuilleForcePlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to which the force should be applied.
  • force: The force to apply.
  • flipDirection: x, y or z. Indicates along which direction the force flips sign.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class AddSinusoidalForcePlugin : public mirheo::SimulationPlugin

Add a sinusoidal (Kolmogorov) force to every particle of a given ParticleVector at every time step.

The force is added at the beforeForce() stage. The force has a reversed sign in half of the domain.

Public Functions

AddSinusoidalForcePlugin(const MirState *state, const std::string &name, const std::string &pvName, real magnitude, int waveNumber)

Create a AddSinusoidalForcePlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to which the force should be applied.
  • magnitude: Magnitude of the force.
  • waveNumber: How many wavelengths along y.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class AddTorquePlugin : public mirheo::SimulationPlugin

Add a constant torque to every object of a given RigidObjectVector at every time step.

The torque is added at the beforeForce() stage.

Public Functions

AddTorquePlugin(const MirState *state, const std::string &name, const std::string &rovName, real3 torque)

Create a AddTorquePlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • rovName: The name of the RigidObjectVector to which the torque should be applied.
  • torque: The torque to apply.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class AnchorParticlesPlugin : public mirheo::SimulationPlugin

Add constraints on the positions and velocities of given particles of a ParticleVector.

The forces required to keep the particles along the given constrains are recorded and reported via AnchorParticlesStatsPlugin.

Note
This should not be used with RigidObjectVector.
Note
This was designed to work with ObjectVectors containing a single object, on a single rank. Using a plain ParticleVector might not work since particles will be reordered.

Public Functions

AnchorParticlesPlugin(const MirState *state, std::string name, std::string pvName, FuncTime3D positions, FuncTime3D velocities, std::vector<int> pids, int reportEvery)

Create a AnchorParticlesPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector that contains the particles of interest.
  • positions: The constrains on the positions.
  • velocities: The constrains on the velocities.
  • pids: The concerned particle ids (starting from 0). See the restrictions in the class docs.
  • reportEvery: Statistics (forces) will be sent to the AnchorParticlesStatsPlugin every this number of steps.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class AnchorParticlesStatsPlugin : public mirheo::PostprocessPlugin

Postprocessing side of AnchorParticlesPlugin.

Reports the forces required to achieve the constrains in a csv file.

Public Functions

AnchorParticlesStatsPlugin(std::string name, std::string path)

Create a AnchorParticlesStatsPlugin object.

Parameters
  • name: The name of the plugin.
  • path: The directory to which the stats will be dumped. Will create a single file path/<pv_name>.csv.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void setup(const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the PostprocessPlugin.

This method must be called before any other function call.

Parameters
  • comm: Contains all postprocess ranks
  • interComm: used to communicate with the simulation ranks

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

class BerendsenThermostatPlugin : public mirheo::SimulationPlugin

Apply Berendsen thermostat to the given particles.

Public Functions

BerendsenThermostatPlugin(const MirState *state, std::string name, std::vector<std::string> pvNames, real kBT, real tau, bool increaseIfLower)

Create a BerendsenThermostatPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: The list of names of the concerned ParticleVector s.
  • kBT: The target temperature, in energy units.
  • tau: The relaxation time.
  • increaseIfLower: Whether to increase the temperature if it’s lower than the target temperature.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class DensityControlPlugin : public mirheo::SimulationPlugin

Apply forces on particles in order to keep the number density constant within layers in a field.

The layers are determined by the level sets of the field. Forces are perpendicular to these layers; their magnitude is computed from PID controllers.

Cannot be used with multiple invocations of Mirheo.run.

Public Types

using RegionFunc = std::function<real(real3)>

functor that describes the region in terms of level sets.

Public Functions

DensityControlPlugin(const MirState *state, std::string name, std::vector<std::string> pvNames, real targetDensity, RegionFunc region, real3 resolution, real levelLo, real levelHi, real levelSpace, real Kp, real Ki, real Kd, int tuneEvery, int dumpEvery, int sampleEvery)

Create a DensityControlPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: The names of the ParticleVector that have the target density..
  • targetDensity: The target number density.
  • region: The field used to partition the space.
  • resolution: The grid spacing used to discretized region
  • levelLo: The minimum level set of the region to control.
  • levelHi: The maximum level set of the region to control.
  • levelSpace: Determines the difference between 2 consecutive layers in the partition of space.
  • Kp: “Proportional” coefficient of the PID.
  • Ki: “Integral” coefficient of the PID.
  • Kd: “Derivative” coefficient of the PID.
  • tuneEvery: Update th PID controllers every this number of steps.
  • dumpEvery: Dump statistics every this number of steps. See also PostprocessDensityControl.
  • sampleEvery: Sample statistics every this number of steps. Used by PIDs.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

void checkpoint(MPI_Comm comm, const std::string &path, int checkPointId)

Save the state of the object on disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.
  • checkPointId: The id of the dump.

void restart(MPI_Comm comm, const std::string &path)

Load the state of the object from the disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.

struct LevelBounds

Helper structure to partition the space.

Public Members

real lo

Smallest level set.

real hi

Largest level set.

real space

Difference between two level sets.

class ExternalMagneticTorquePlugin : public mirheo::SimulationPlugin

Apply a magnetic torque on given a RigidObjectVector.

Public Types

using UniformMagneticFunc = std::function<real3(real)>

Time varying uniform field.

Public Functions

ExternalMagneticTorquePlugin(const MirState *state, std::string name, std::string rovName, real3 moment, UniformMagneticFunc magneticFunction)

Create a ExternalMagneticTorquePlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • rovName: The name of the RigidObjectVector to apply the torque to.
  • moment: The constant magnetic moment of one object, in its frame of reference.
  • magneticFunction: The external uniform magnetic field which possibly varies in time.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class PostprocessDensityControl : public mirheo::PostprocessPlugin

Postprocessing side of DensityControlPlugin.

Dumps the density and force in each layer of the space partition.

Public Functions

PostprocessDensityControl(std::string name, std::string filename)

Create a PostprocessDensityControl object.

Parameters
  • name: The name of the plugin.
  • filename: The txt file that will contain the density and corresponding force magnitudes in each layer.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

class ParticleDisplacementPlugin : public mirheo::SimulationPlugin

Compute the dispacement of particles between a given number of time steps.

Public Functions

ParticleDisplacementPlugin(const MirState *state, std::string name, std::string pvName, int updateEvery)

Create a ParticleDisplacementPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the concerned ParticleVector.
  • updateEvery: The number of steps between two steps used to compute the displacement.

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class ExchangePVSFluxPlanePlugin : public mirheo::SimulationPlugin

Transfer particles from one ParticleVector to another when they cross a given plane.

Public Functions

ExchangePVSFluxPlanePlugin(const MirState *state, std::string name, std::string pv1Name, std::string pv2Name, real4 plane)

Create a ExchangePVSFluxPlanePlugin object.

The particle has crossed the plane if a *x + b * y + c * z + d goes from negative to positive.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pv1Name: The name of the source ParticleVector. Only particles from this ParticleVector are transfered.
  • pv2Name: The name of the destination ParticleVector.
  • plane: Coefficients of the plane to be crossed, (a, b, c, d).

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeCellLists(cudaStream_t stream)

hook before building the cell lists

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class ForceSaverPlugin : public mirheo::SimulationPlugin

Copies the forces of a given ParticleVector to a new channel at every time step.

This allows to dump the forces since they are reset to zero at every time step.

Public Functions

ForceSaverPlugin(const MirState *state, std::string name, std::string pvName)

Create a ForceSaverPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to save forces from and to.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

void beforeIntegration(cudaStream_t stream)

hook before integrating the particle vectors but after the forces are computed

class ImposeProfilePlugin : public mirheo::SimulationPlugin

Set the velocity to a given one in a box.

The velocity is set to a constant plus a random velocity that has a Maxwell distribution.

Public Functions

ImposeProfilePlugin(const MirState *state, std::string name, std::string pvName, real3 low, real3 high, real3 targetVel, real kBT)

Create a ImposeProfilePlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to modify.
  • low: Lower coordinates of the region of interest.
  • high: Upper coordinates of the region of interest.
  • targetVel: The constant part of the new velocity.
  • kBT: Temperature used to draw the velocity from the maxwellian distribution.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

class ImposeVelocityPlugin : public mirheo::SimulationPlugin

Add a constant to the velocity of particles in a given region such that it matches a given average.

Public Functions

ImposeVelocityPlugin(const MirState *state, std::string name, std::vector<std::string> pvNames, real3 low, real3 high, real3 targetVel, int every)

Create a ImposeVelocityPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: The name of the (list of) ParticleVector to modify.
  • low: Lower coordinates of the region of interest.
  • high: Upper coordinates of the region of interest.
  • targetVel: The target average velocity in the region.
  • every: Correct the velocity every this number of time steps.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

void setTargetVelocity(real3 v)

Change the target velocity to a new value.

Parameters
  • v: The new target velocity.

class MagneticDipoleInteractionsPlugin : public mirheo::SimulationPlugin

Compute the magnetic dipole-dipole forces and torques induced by the interactions between rigid objects that have a magnetic moment.

Public Functions

MagneticDipoleInteractionsPlugin(const MirState *state, std::string name, std::string rovName, real3 moment, real mu0, bool periodic)

Create a MagneticDipoleInteractionsPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • rovName: The name of the RigidObjectVector interacting.
  • moment: The constant magnetic moment of one object, in its frame of reference.
  • mu0: The magnetic permeability of the medium.
  • periodic: If true, computes interactions with the closest periodic image of the objects

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeCellLists(cudaStream_t stream)

hook before building the cell lists

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class OutletPlugin : public mirheo::SimulationPlugin

Base class for outlet Plugins.

Outlet plugins delete particles of given a ParticleVector list in a region.

Subclassed by mirheo::PlaneOutletPlugin, mirheo::RegionOutletPlugin

Public Functions

OutletPlugin(const MirState *state, std::string name, std::vector<std::string> pvNames)

Create a OutletPlugin.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: List of names of the ParticleVector that the outlet will be applied to.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class PlaneOutletPlugin : public mirheo::OutletPlugin

Delete all particles that cross a given plane.

Public Functions

PlaneOutletPlugin(const MirState *state, std::string name, std::vector<std::string> pvNames, real4 plane)

Create a PlaneOutletPlugin.

A particle crosses the plane if a*x + b*y + c*z + d goes from nbegative to postive accross one time step.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: List of names of the ParticleVector that the outlet will be applied to.
  • plane: Coefficients (a, b, c, d) of the plane.

void beforeCellLists(cudaStream_t stream)

hook before building the cell lists

class RegionOutletPlugin : public mirheo::OutletPlugin

Delete all particles in a given region, defined implicitly by a field.

A particle is considered inside the region if the given field is negative at the particle’s position.

Subclassed by mirheo::DensityOutletPlugin, mirheo::RateOutletPlugin

Public Types

using RegionFunc = std::function<real(real3)>

A scalar field to represent inside (negative) / outside (positive) region.

Public Functions

RegionOutletPlugin(const MirState *state, std::string name, std::vector<std::string> pvNames, RegionFunc region, real3 resolution)

Create a RegionOutletPlugin.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: List of names of the ParticleVector that the outlet will be applied to.
  • region: The field that describes the region. This will be sampled on a uniform grid and uploaded to the GPU.
  • resolution: The grid space used to discretize region.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class DensityOutletPlugin : public mirheo::RegionOutletPlugin

Delete particles located in a given region if the number density is higher than a target one.

Public Functions

DensityOutletPlugin(const MirState *state, std::string name, std::vector<std::string> pvNames, real numberDensity, RegionFunc region, real3 resolution)

Create a DensityOutletPlugin.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: List of names of the ParticleVector that the outlet will be applied to.
  • numberDensity: The target number density.
  • region: The field that describes the region. This will be sampled on a uniform grid and uploaded to the GPU.
  • resolution: The grid space used to discretize region.

void beforeCellLists(cudaStream_t stream)

hook before building the cell lists

class RateOutletPlugin : public mirheo::RegionOutletPlugin

Delete particles located in a given region at a given rate.

Public Functions

RateOutletPlugin(const MirState *state, std::string name, std::vector<std::string> pvNames, real rate, RegionFunc region, real3 resolution)

Create a RateOutletPlugin.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: List of names of the ParticleVector that the outlet will be applied to.
  • rate: The rate of deletion of particles.
  • region: The field that describes the region. This will be sampled on a uniform grid and uploaded to the GPU.
  • resolution: The grid space used to discretize region.

void beforeCellLists(cudaStream_t stream)

hook before building the cell lists

class ParticleChannelAveragerPlugin : public mirheo::SimulationPlugin

Average over time a particle vector channel.

Public Functions

ParticleChannelAveragerPlugin(const MirState *state, std::string name, std::string pvName, std::string channelName, std::string averageName, real updateEvery)

Create a ParticleChannelAveragerPlugin.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector.
  • channelName: The name of the channel to average. Will fail if it does not exist.
  • averageName: The name of the new channel, that will contain the time-averaged quantity..
  • updateEvery: Will reset the averaged channel every this number of steps. Must be positive.

void beforeIntegration(cudaStream_t stream)

hook before integrating the particle vectors but after the forces are computed

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

class ParticleChannelSaverPlugin : public mirheo::SimulationPlugin

Copies a given channel to another one that will “stick” to the particle vector.

This is useful to collect statistics on non permanent quantities (e.g. stresses).

Public Functions

ParticleChannelSaverPlugin(const MirState *state, std::string name, std::string pvName, std::string channelName, std::string savedName)

Create a ParticleChannelSaverPlugin.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector.
  • channelName: The name of the channel to save at every time step. Will fail if it does not exist.
  • savedName: The name of the new channel.

void beforeIntegration(cudaStream_t stream)

hook before integrating the particle vectors but after the forces are computed

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

class ParticleDragPlugin : public mirheo::SimulationPlugin

Apply a drag force proportional to the velocity of every particle in a ParticleVector.

Public Functions

ParticleDragPlugin(const MirState *state, std::string name, std::string pvName, real drag)

Create a ParticleDragPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to which the force should be applied.
  • drag: The drag coefficient applied to each particle.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class PinObjectPlugin : public mirheo::SimulationPlugin

Add constraints on objects of an ObjectVector.

This modifies the velocities and forces on the objects in order to satisfy the given constraints on linear and angular velocities.

This plugin also collects statistics on the required forces and torques used to maintain the constraints. This may be useful to e.g. measure the drag around objects. See ReportPinObjectPlugin.

Public Functions

PinObjectPlugin(const MirState *state, std::string name, std::string ovName, real3 translation, real3 rotation, int reportEvery)

Create a PinObjectPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • ovName: The name of the ObjectVector that will be subjected of the constraints.
  • translation: The target linear velocity. Components set to Unrestricted will not be constrained.
  • rotation: The target angular velocity. Components set to Unrestricted will not be constrained.
  • reportEvery: Send forces and torques stats to the postprocess side every this number of time steps.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeIntegration(cudaStream_t stream)

hook before integrating the particle vectors but after the forces are computed

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

void handshake()

Used to communicate initial information between a SimulationPlugin and a PostprocessPlugin.

Does not do anything by default.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

Public Static Attributes

constexpr real Unrestricted = std::numeric_limits<real>::infinity()

Special value reserved to represent unrestricted components.

class PinRodExtremityPlugin : public mirheo::SimulationPlugin

Add alignment force on a rod segment.

Public Functions

PinRodExtremityPlugin(const MirState *state, std::string name, std::string rvName, int segmentId, real fmagn, real3 targetDirection)

Create a PinRodExtremityPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • rvName: The name of the RodVector to which the force should be applied.
  • segmentId: The segment that will be constrained.
  • fmagn: The force coefficient.
  • targetDirection: The target direction.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeIntegration(cudaStream_t stream)

hook before integrating the particle vectors but after the forces are computed

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class TemperaturizePlugin : public mirheo::SimulationPlugin

Add or set maxwellian drawn velocities to the particles of a given ParticleVector.

Public Functions

TemperaturizePlugin(const MirState *state, std::string name, std::string pvName, real kBT, bool keepVelocity)

Create a TemperaturizePlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to modify.
  • kBT: Target temperature.
  • keepVelocity: Wether to add or reset the velocities.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class SimulationVelocityControl : public mirheo::SimulationPlugin

Apply a force in a box region to all particles.

The force is controled by a PID controller that has a target mean velocity in that same region.

Public Functions

SimulationVelocityControl(const MirState *state, std::string name, std::vector<std::string> pvNames, real3 low, real3 high, int sampleEvery, int tuneEvery, int dumpEvery, real3 targetVel, real Kp, real Ki, real Kd)

Create a SimulationVelocityControl object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvNames: The list of names of the ParticleVector to control.
  • low: The lower coordinates of the control region.
  • high: The upper coordinates of the control region.
  • sampleEvery: Sample the velocity average every this number of steps.
  • tuneEvery: Update the PID controller every this number of steps.
  • dumpEvery: Send statistics of the PID to the postprocess plugin every this number of steps.
  • targetVel: The target mean velocity in the region of interest.
  • Kp: “Proportional” coefficient of the PID.
  • Ki: “Integral” coefficient of the PID.
  • Kd: “Derivative” coefficient of the PID.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeForces(cudaStream_t stream)

hook before computing the forces and after the cell lists are created

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

void serializeAndSend(cudaStream_t stream)

Pack and send data to the postprocess rank.

Happens between beforeForces() and beforeIntegration().

Note
This may happens while computing the forces.

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

void checkpoint(MPI_Comm comm, const std::string &path, int checkPointId)

Save the state of the object on disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.
  • checkPointId: The id of the dump.

void restart(MPI_Comm comm, const std::string &path)

Load the state of the object from the disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.

class PostprocessVelocityControl : public mirheo::PostprocessPlugin

Postprocess side of SimulationVelocityControl.

Receives and dump the PID stats to a csv file.

Public Functions

PostprocessVelocityControl(std::string name, std::string filename)

Create a SimulationVelocityControl object.

Parameters
  • name: The name of the plugin.
  • filename: The csv file to which the statistics will be dumped.

void deserialize()

Perform the action implemented by the plugin using the data received from the SimulationPlugin.

void checkpoint(MPI_Comm comm, const std::string &path, int checkPointId)

Save the state of the object on disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.
  • checkPointId: The id of the dump.

void restart(MPI_Comm comm, const std::string &path)

Load the state of the object from the disk.

Parameters
  • comm: MPI communicator to perform the I/O.
  • path: The directory path to store the object state.

class VelocityInletPlugin : public mirheo::SimulationPlugin

Add particles to a given ParticleVector.

The particles are injected on a given surface at a given influx rate.

Public Types

using ImplicitSurfaceFunc = std::function<real(real3)>

Representation of a surface from a scalar field.

using VelocityFieldFunc = std::function<real3(real3)>

Velocity field used to describe the inflow.

Public Functions

VelocityInletPlugin(const MirState *state, std::string name, std::string pvName, ImplicitSurfaceFunc implicitSurface, VelocityFieldFunc velocityField, real3 resolution, real numberDensity, real kBT)

Create a VelocityInletPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector to add the particles to.
  • implicitSurface: The scalar field that has the desired surface as zero level set.
  • velocityField: The inflow velocity. Only relevant on the surface.
  • resolution: Grid size used to sample the fields.
  • numberDensity: The target number density of injection.
  • kBT: The temperature of the injected particles.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeCellLists(cudaStream_t stream)

hook before building the cell lists

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

class WallRepulsionPlugin : public mirheo::SimulationPlugin

Add a force that pushes particles away from the wall surfaces.

The magnitude of the force decreases linearly down to zero at a given distance h. Furthermore, the force can be capped.

Public Functions

WallRepulsionPlugin(const MirState *state, std::string name, std::string pvName, std::string wallName, real C, real h, real maxForce)

Create a WallRepulsionPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • pvName: The name of the ParticleVector that will be subject to the force.
  • wallName: The name of the Wall.
  • C: Force coefficient.
  • h: Force maximum distance.
  • maxForce: Maximum force magnitude.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeIntegration(cudaStream_t stream)

hook before integrating the particle vectors but after the forces are computed

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

Debugging plugins

class ParticleCheckerPlugin : public mirheo::SimulationPlugin

Check the validity of all ParticleVector in the simulation:

  • Check that the positions are within reasonable bounds.
  • Check that the velocities are within reasonable bounds.
  • Check that forces do not contain NaN or Inf values.

If either of the above is not satisfied, the plugin will make the code die with an informative error.

Public Types

enum Info

Encode error type if there is any.

Values:

Ok
Out
Nan

Public Functions

ParticleCheckerPlugin(const MirState *state, std::string name, int checkEvery)

Create a ParticleCheckerPlugin object.

Parameters
  • state: The global state of the simulation.
  • name: The name of the plugin.
  • checkEvery: Will check the states of particles every this number of steps.

void setup(Simulation *simulation, const MPI_Comm &comm, const MPI_Comm &interComm)

setup the internal state of the SimulationPlugin.

This method must be called before any of the hooks of the plugin. This is the place to fetch reference to objects held by the simulation.

Parameters
  • simulation: The simulation to which the plugin is registered.
  • comm: Contains all simulation ranks
  • interComm: used to communicate with the postprocess ranks

void beforeIntegration(cudaStream_t stream)

hook before integrating the particle vectors but after the forces are computed

void afterIntegration(cudaStream_t stream)

hook after the ObjectVector objects are integrated but before redistribution and bounce back

bool needPostproc()

Return
true if this plugin needs a postprocess side; false otherwise.
Note
The plugin can have a postprocess side but not need it.

struct Status

Helper to encode problematic particles.

Public Members

int id

The Index of the potential problematic particle.

Info info

What is problematic.