Pairwise Interactions¶
Base class¶
This is the visible class that is output of the factory function.
-
class
BasePairwiseInteraction
: public mirheo::Interaction¶ Base class for short-range symmetric pairwise interactions.
Subclassed by mirheo::PairwiseDensityInteraction< DensityKernel >, mirheo::PairwiseDPDInteraction, mirheo::PairwiseGrowingRepulsiveLJInteraction< Awareness >, mirheo::PairwiseLJInteraction, mirheo::PairwiseMDPDInteraction, mirheo::PairwiseMorseInteraction< Awareness >, mirheo::PairwiseNoRandomDPDInteraction, mirheo::PairwiseRepulsiveLJInteraction< Awareness >, mirheo::PairwiseSDPDInteraction< PressureEOS, DensityKernel >, mirheo::PairwiseViscoElasticDPDInteraction, mirheo::PairwiseViscoElasticSmoothVelDPDInteraction
Public Functions
-
BasePairwiseInteraction
(const MirState *state, const std::string &name, real rc)¶ Construct a base pairwise interaction from parameters.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.
-
std::optional<real>
getCutoffRadius
() const¶ - Return
- the cut-off radius of the pairwise interaction.
-
Interactions¶
Here are listed the current available pairwise interactions.
-
class
PairwiseDPDInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute dissipative particle dynamics (DPD) forces.
This interaction needs the particles positions and velocities.
Public Functions
-
PairwiseDPDInteraction
(const MirState *state, const std::string &name, real rc, DPDParams params, std::optional<real> stressPeriod = std::nullopt)¶ Create a PairwiseDPDInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the DPD forces.stressPeriod
: The simulation time between two stress computations. If set tostd::nullopt
, disables stress computation.
-
void
setPrerequisites
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)¶ Add needed properties to the given ParticleVectors for future interactions.
Must be called before any other method of this class.
- Parameters
pv1
: One ParticleVector of the interactionpv2
: The other ParticleVector of that will interactcl1
: CellList of pv1cl2
: CellList of pv2
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
std::vector<InteractionChannel>
getOutputChannels
() const¶ Returns which channels are those output by the interactions.
-
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
PairwiseViscoElasticDPDInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute extended dissipative particle dynamics (DPD) forces with visco-elastic component.
This interaction requires an extra state variable (polymeric chain vector end to end Q) per particle. This interaction needs the particles positions, velocities and polymeric chain vectors. In addition to forces, the time derivative of the polymeric chain vectors (dQ/dt) is computed. This interaction must be used with a special integrator that also evolves the polymeric chain vectors over time.
Public Functions
-
PairwiseViscoElasticDPDInteraction
(const MirState *state, const std::string &name, real rc, ViscoElasticDPDParams params, std::optional<real> stressPeriod = std::nullopt)¶ Create a PairwiseViscoElasticDPDInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the DPD forces.stressPeriod
: The simulation time between two stress computations. If set tostd::nullopt
, disables stress computation.
-
void
setPrerequisites
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)¶ Add needed properties to the given ParticleVectors for future interactions.
Must be called before any other method of this class.
- Parameters
pv1
: One ParticleVector of the interactionpv2
: The other ParticleVector of that will interactcl1
: CellList of pv1cl2
: CellList of pv2
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
std::vector<InteractionChannel>
getInputChannels
() const¶ Returns which channels are required as input.
Positions and velocities are always required and are not listed here; Only other channels must be specified here.
-
std::vector<InteractionChannel>
getOutputChannels
() const¶ Returns which channels are those output by the interactions.
-
-
class
PairwiseViscoElasticSmoothVelDPDInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute extended dissipative particle dynamics (DPD) forces with visco-elastic component.
This interaction requires an extra state variable (polymeric chain vector end to end Q) per particle. This interaction needs the particles positions, velocities and polymeric chain vectors. In addition to forces, the time derivative of the polymeric chain vectors (dQ/dt) is computed. This interaction must be used with a special integrator that also evolves the polymeric chain vectors over time.
Public Functions
-
PairwiseViscoElasticSmoothVelDPDInteraction
(const MirState *state, const std::string &name, real rc, ViscoElasticDPDParams params)¶ Create a PairwiseViscoElasticDPDInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the DPD forces.
-
void
setPrerequisites
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)¶ Add needed properties to the given ParticleVectors for future interactions.
Must be called before any other method of this class.
- Parameters
pv1
: One ParticleVector of the interactionpv2
: The other ParticleVector of that will interactcl1
: CellList of pv1cl2
: CellList of pv2
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
std::vector<InteractionChannel>
getInputChannels
() const¶ Returns which channels are required as input.
Positions and velocities are always required and are not listed here; Only other channels must be specified here.
-
std::vector<InteractionChannel>
getOutputChannels
() const¶ Returns which channels are those output by the interactions.
-
-
template <class PressureEOS, class DensityKernel>
classPairwiseSDPDInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute smooth dissipative particle dynamics (SDPD) forces.
This interaction needs the particles positions velocities and densities.
- Template Parameters
PressureEOS
: The kernel used to compute the equation of state.DensityKernel
: The kernel used to compute the number densities.
Public Functions
-
PairwiseSDPDInteraction
(const MirState *state, const std::string &name, real rc, SDPDParams params, std::optional<real> stressPeriod = std::nullopt)¶ Create a PairwiseSDPDInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the DPD forces.stressPeriod
: The simulation time between two stress computations. If set tostd::nullopt
, disables stress computation.
-
void
setPrerequisites
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)¶ Add needed properties to the given ParticleVectors for future interactions.
Must be called before any other method of this class.
- Parameters
pv1
: One ParticleVector of the interactionpv2
: The other ParticleVector of that will interactcl1
: CellList of pv1cl2
: CellList of pv2
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
std::vector<InteractionChannel>
getInputChannels
() const¶ Returns which channels are required as input.
Positions and velocities are always required and are not listed here; Only other channels must be specified here.
-
std::vector<InteractionChannel>
getOutputChannels
() const¶ Returns which channels are those output by the interactions.
-
template <class DensityKernel>
classPairwiseDensityInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute the density per particle.
This interaction needs particles positions and adds the particles contributions to the number density channel.
- Template Parameters
DensityKernel
: The blob function used to compute the density.
Public Functions
-
PairwiseDensityInteraction
(const MirState *state, const std::string &name, real rc, DensityParams params)¶ Create a PairwiseDensityInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the density kernel.
-
void
setPrerequisites
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)¶ Add needed properties to the given ParticleVectors for future interactions.
Must be called before any other method of this class.
- Parameters
pv1
: One ParticleVector of the interactionpv2
: The other ParticleVector of that will interactcl1
: CellList of pv1cl2
: CellList of pv2
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
Stage
getStage
() const¶ returns the Stage corresponding to this interaction.
-
std::vector<InteractionChannel>
getOutputChannels
() const¶ Returns which channels are those output by the interactions.
-
class
PairwiseLJInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute Lennard-Jones (LJ) forces.
Public Functions
-
PairwiseLJInteraction
(const MirState *state, const std::string &name, real rc, LJParams params, std::optional<real> stressPeriod = std::nullopt)¶ Create a PairwiseLJInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the LJ forces.stressPeriod
: The simulation time between two stress computations. If set tostd::nullopt
, disables stress computation.
-
void
setPrerequisites
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)¶ Add needed properties to the given ParticleVectors for future interactions.
Must be called before any other method of this class.
- Parameters
pv1
: One ParticleVector of the interactionpv2
: The other ParticleVector of that will interactcl1
: CellList of pv1cl2
: CellList of pv2
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
std::vector<InteractionChannel>
getOutputChannels
() const¶ Returns which channels are those output by the interactions.
-
-
template <class Awareness>
classPairwiseRepulsiveLJInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute the repulsive Lennard-Jones (LJ) forces between particles.
- Template Parameters
Awareness
: To control which particles interacxt with which particle (e.g. avoiding interactions between particles of the same object).
Public Functions
-
PairwiseRepulsiveLJInteraction
(const MirState *state, const std::string &name, real rc, RepulsiveLJParams params, std::optional<real> stressPeriod = std::nullopt)¶ Create a PairwiseRepulsiveLJInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the forces.stressPeriod
: The simulation time between two stress computations. If set tostd::nullopt
, disables stress computation.
-
void
setPrerequisites
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)¶ Add needed properties to the given ParticleVectors for future interactions.
Must be called before any other method of this class.
- Parameters
pv1
: One ParticleVector of the interactionpv2
: The other ParticleVector of that will interactcl1
: CellList of pv1cl2
: CellList of pv2
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
std::vector<InteractionChannel>
getOutputChannels
() const¶ Returns which channels are those output by the interactions.
-
template <class Awareness>
classPairwiseGrowingRepulsiveLJInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute the repulsive Lennard-Jones (LJ) forces between particles.
- Template Parameters
Awareness
: To control which particles interact with which particle (e.g. avoiding interactions between particles of the same object).
Public Functions
-
PairwiseGrowingRepulsiveLJInteraction
(const MirState *state, const std::string &name, real rc, GrowingRepulsiveLJParams params, std::optional<real> stressPeriod = std::nullopt)¶ Create a PairwiseGrowingRepulsiveLJInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the forces.stressPeriod
: The simulation time between two stress computations. If set tostd::nullopt
, disables stress computation.
-
void
setPrerequisites
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)¶ Add needed properties to the given ParticleVectors for future interactions.
Must be called before any other method of this class.
- Parameters
pv1
: One ParticleVector of the interactionpv2
: The other ParticleVector of that will interactcl1
: CellList of pv1cl2
: CellList of pv2
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
std::vector<InteractionChannel>
getOutputChannels
() const¶ Returns which channels are those output by the interactions.
-
template <class Awareness>
classPairwiseMorseInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute forces arising from the Morse potential between particles.
- Template Parameters
Awareness
: To control which particles interact with which particle (e.g. avoiding interactions between particles of the same object).
Public Functions
-
PairwiseMorseInteraction
(const MirState *state, const std::string &name, real rc, MorseParams params, std::optional<real> stressPeriod = std::nullopt)¶ Create a PairwiseMorseInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the forces.stressPeriod
: The simulation time between two stress computations. If set tostd::nullopt
, disables stress computation.
-
void
setPrerequisites
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)¶ Add needed properties to the given ParticleVectors for future interactions.
Must be called before any other method of this class.
- Parameters
pv1
: One ParticleVector of the interactionpv2
: The other ParticleVector of that will interactcl1
: CellList of pv1cl2
: CellList of pv2
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
std::vector<InteractionChannel>
getOutputChannels
() const¶ Returns which channels are those output by the interactions.
-
class
PairwiseNoRandomDPDInteraction
: public mirheo::BasePairwiseInteraction¶ Interaction to compute dissipative particle dynamics (DPD) forces without the random term.
Useful for testing.
This interaction needs the particles positions and velocities.
Public Functions
-
PairwiseNoRandomDPDInteraction
(const MirState *state, const std::string &name, real rc, NoRandomDPDParams params)¶ Create a PairwiseNoRandomDPDInteraction object.
- Parameters
state
: The global state of the system.name
: The name of the interaction.rc
: The cutoff radius of the interaction. Must be positive and smaller than the sub-domain size.params
: The parameters of the DPD forces.
-
void
local
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. The order of pv1 and pv2 may change the performance of the interactions.- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
void
halo
(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute interactions between bulk particles and halo particles.
The result of the interaction is
added to the corresponding channel of the ParticleVector. In general, the following interactions will be computed: pv1->halo() <> pv2->local() and pv2->halo() <> pv1->local().- Parameters
pv1
: first interacting ParticleVectorpv2
: second interacting ParticleVector. If it is the same as the pv1, self interactions will be computed.cl1
: cell-list built for the appropriate cut-off radius for pv1cl2
: cell-list built for the appropriate cut-off radius for pv2stream
: Execution stream
-
Implementation¶
A generic pairwise CUDA kernel computes interactions between two ParticleVector
.
The kernel is templated and requires a functor describing the interaction forces (see Kernels)
Stress computation is optional with compatible interactions. A helper class is provided:
-
class
StressManager
¶ Helper class to control when the stress must be computed.
In theory stresses could be computed at every time step, but with a higher cost (higher bandwidth needed and more computation). In practice stress does not need to be computed at every iteration. Instead we compute stresses every given periods of simulation time.
Public Functions
-
StressManager
(real stressPeriod)¶ Construct a StressManager object.
- Parameters
stressPeriod
: Time (in simulation units) between two stress computations.
-
template <class PairwiseKernel>
voidcomputeLocalInteractions
(const MirState *state, PairwiseKernel &pair, PairwiseStressWrapper<PairwiseKernel> &pairWithStress, ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Compute local interactions between two ParticleVector.
The stresses are computed only when needed.
- Template Parameters
PairwiseKernel
: The symmetric pairwise kernel of the interaction.
- Parameters
state
: The global state of the simulation.pair
: The interaction kernel (without stresses).pairWithStress
: The interaction kernel (with stresses).pv1
: First ParticleVector.pv2
: Second ParticleVector.cl1
: CellList of pv1.cl2
: CellList of pv2.stream
: Stream of execution.
-
template <class PairwiseKernel>
voidcomputeHaloInteractions
(const MirState *state, PairwiseKernel &pair, PairwiseStressWrapper<PairwiseKernel> &pairWithStress, ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2, cudaStream_t stream)¶ Computehal halo interactions between two ParticleVector.
The stresses are computed only when needed.
- Template Parameters
PairwiseKernel
: The symmetric pairwise kernel of the interaction.
- Parameters
state
: The global state of the simulation.pair
: The interaction kernel (without stresses).pairWithStress
: The interaction kernel (with stresses).pv1
: First ParticleVector.pv2
: Second ParticleVector.cl1
: CellList of pv1.cl2
: CellList of pv2.stream
: Stream of execution.
-
Interaction::InteractionChannel
getStressPredicate
(const MirState *state) const¶ - Return
- InteractioChannel for the output of a pairwise Interaction with predicate active when stresses are needed.
- Parameters
state
: The global state of the simulation.
-
Kernels¶
Interface¶
The kernel functor describes a pairwise interaction. It is splitted into two parts:
- a handler, that must be usable on the device.
- a manager, that may store extra information on the host. For simple interactions, this can be the same as the handler class.
The interface of the functor must follow the following requirements:
Define a view type to be passed (e.g.
mirheo::PVview
) as well as a particle type to be fetched and the parameter struct used for initialization:using ViewType = <particle vector view type> using ParticleType = <particle type> using HandlerType = <type passed to GPU> using ParamsType = <struct that contains the parameters of this functor>
A generic constructor from the
ParamsType
parameters:PairwiseKernelType(real rc, const ParamsType& p, real dt, long seed=42424242);
Setup function (on Host, for manager only)
void setup(LocalParticleVector* lpv1, LocalParticleVector* lpv2, CellList* cl1, CellList* cl2, const MirState *state);
Handler function (on Host, for manager only)
const HandlerType& handler() const;
Interaction function (output must match with accumulator, see below) (on GPU)
__D__ <OutputType> operator()(const ParticleType dst, int dstId, const ParticleType src, int srcId) const;
Accumulator initializer (on GPU)
__D__ <Accumulator> getZeroedAccumulator() const;
Fetch functions (see in fetchers.h or see the docs):
__D__ ParticleType read(const ViewType& view, int id) const; __D__ ParticleType readNoCache(const ViewType& view, int id) const; __D__ void readCoordinates(ParticleType& p, const ViewType& view, int id) const; __D__ void readExtraData(ParticleType& p, const ViewType& view, int id) const;
Interacting checker to discard pairs not within cutoff:
__D__ bool withinCutoff(const ParticleType& src, const ParticleType& dst) const;
Position getter from generic particle type:
__D__ real3 getPosition(const ParticleType& p) const;
Note
To implement a new kernel, the following must be done: - satisfy the above interface - add a corresponding parameter in parameters.h - if necessary, add type traits specialization in type_traits.h
This is the interface for the host calls:
-
class
PairwiseKernel
¶ Interface of host methods required for a pairwise kernel.
Subclassed by mirheo::PairwiseDensity< DensityKernel >, mirheo::PairwiseDPD, mirheo::PairwiseLJ, mirheo::PairwiseMDPD, mirheo::PairwiseMorse< Awareness >, mirheo::PairwiseNoRandomDPD, mirheo::PairwiseRepulsiveLJ< Awareness >, mirheo::PairwiseSDPD< PressureEOS, DensityKernel >, mirheo::PairwiseViscoElasticDPD, mirheo::PairwiseViscoElasticSmoothVelDPD
Public Functions
-
virtual void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup the internal state of the functor
-
virtual void
writeState
(std::ofstream &fout)¶ write internal state to a stream
-
virtual bool
readState
(std::ifstream &fin)¶ restore internal state from a stream
-
virtual void
The rest is directly implemented in the kernels, as no virtual functions are allowed on the device.
Implemented kernels¶
-
template <typename DensityKernel>
classPairwiseDensity
: public mirheo::PairwiseKernel, public mirheo::ParticleFetcher¶ Compute number density from pairwise kernel.
- Template Parameters
DensityKernel
: The kernel used to evaluate the number density
Public Functions
-
PairwiseDensity
(real rc, DensityKernel densityKernel)¶ construct from density kernel
-
PairwiseDensity
(real rc, const ParamsType &p, long seed = 42424242)¶ generic constructor
-
real
operator()
(const ParticleType dst, int dstId, const ParticleType src, int srcId) const¶ evaluate the number density contribution of this pair
-
DensityAccumulator
getZeroedAccumulator
() const¶ initialize the accumulator
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
class
PairwiseDPDHandler
: public mirheo::ParticleFetcher¶ a GPU compatible functor that computes DPD interactions
Subclassed by mirheo::PairwiseDPD
Public Types
Public Functions
-
PairwiseDPDHandler
(real rc, real a, real gamma, real power)¶ constructor
-
real3
operator()
(const ParticleType dst, int dstId, const ParticleType src, int srcId) const¶ evaluate the force
-
ForceAccumulator
getZeroedAccumulator
() const¶ initialize accumulator
-
-
class
PairwiseDPD
: public mirheo::PairwiseKernel, public mirheo::PairwiseDPDHandler¶ Helper class that constructs PairwiseDPDHandler.
Public Types
-
using
HandlerType
= PairwiseDPDHandler¶ handler type corresponding to this object
-
using
ParamsType
= DPDParams¶ parameters that are used to create this object
Public Functions
-
PairwiseDPD
(real rc, real a, real gamma, real kBT, real power, long seed = 42424242)¶ Constructor.
-
PairwiseDPD
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup the internal state of the functor
-
void
writeState
(std::ofstream &fout)¶ write internal state to a stream
-
bool
readState
(std::ifstream &fin)¶ restore internal state from a stream
-
using
-
class
PairwiseLJ
: public mirheo::PairwiseKernel, public mirheo::ParticleFetcher¶ Compute Lennard-Jones forces on the device.
Public Types
-
using
HandlerType
= PairwiseLJ¶ Corresponding handler.
-
using
ParamsType
= LJParams¶ Corresponding parameters type.
Public Functions
-
PairwiseLJ
(real rc, real epsilon, real sigma)¶ Constructor.
-
PairwiseLJ
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
real3
operator()
(ParticleType dst, int, ParticleType src, int) const¶ Evaluate the force.
-
ForceAccumulator
getZeroedAccumulator
() const¶ initialize accumulator
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
using
-
template <class Awareness>
classPairwiseMorse
: public mirheo::PairwiseKernel, public mirheo::ParticleFetcher¶ Compute Morse forces on the device.
- Template Parameters
Awareness
: A functor that describes which particles pairs interact.
Public Functions
-
PairwiseMorse
(real rc, real De, real r0, real beta, real maxForce, Awareness awareness)¶ Constructor.
-
PairwiseMorse
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
real3
operator()
(ParticleType dst, int dstId, ParticleType src, int srcId) const¶ Evaluate the force.
-
ForceAccumulator
getZeroedAccumulator
() const¶ initialize accumulator
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup the internal state of the functor
-
class
PairwiseMDPDHandler
: public mirheo::ParticleFetcherWithDensity¶ a GPU compatible functor that computes MDPD interactions
Subclassed by mirheo::PairwiseMDPD
Public Types
-
using
ViewType
= PVviewWithDensities¶ compatible view type
-
using
ParticleType
= ParticleWithDensity¶ compatible particle type
Public Functions
-
PairwiseMDPDHandler
(real rc, real rd, real a, real b, real gamma, real power)¶ constructor
-
real3
operator()
(const ParticleType dst, int dstId, const ParticleType src, int srcId) const¶ evaluate the force
-
ForceAccumulator
getZeroedAccumulator
() const¶ initialize accumulator
-
using
-
class
PairwiseMDPD
: public mirheo::PairwiseKernel, public mirheo::PairwiseMDPDHandler¶ Helper class that constructs PairwiseMDPDHandler.
Public Types
-
using
HandlerType
= PairwiseMDPDHandler¶ handler type corresponding to this object
-
using
ParamsType
= MDPDParams¶ parameters that are used to create this object
Public Functions
-
PairwiseMDPD
(real rc, real rd, real a, real b, real gamma, real kBT, real power, long seed = 42424242)¶ Constructor.
-
PairwiseMDPD
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup the internal state of the functor
-
void
writeState
(std::ofstream &fout)¶ write internal state to a stream
-
bool
readState
(std::ifstream &fin)¶ restore internal state from a stream
-
using
-
class
PairwiseNoRandomDPD
: public mirheo::PairwiseKernel, public mirheo::ParticleFetcher¶ a GPU compatible functor that computes DPD interactions without fluctuations.
Used in unit tests
Public Types
-
using
HandlerType
= PairwiseNoRandomDPD¶ handler type corresponding to this object
-
using
ParamsType
= NoRandomDPDParams¶ parameters that are used to create this object
Public Functions
-
PairwiseNoRandomDPD
(real rc, real a, real gamma, real kBT, real power)¶ constructor
-
PairwiseNoRandomDPD
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
real3
operator()
(const ParticleType dst, int dstId, const ParticleType src, int srcId) const¶ evaluate the force
-
ForceAccumulator
getZeroedAccumulator
() const¶ initialize accumulator
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup the internal state of the functor
-
using
-
template <class Awareness>
classPairwiseRepulsiveLJ
: public mirheo::PairwiseKernel¶ Kernel for repulsive LJ forces.
Subclassed by mirheo::PairwiseGrowingRepulsiveLJ< Awareness >
Public Functions
-
PairwiseRepulsiveLJ
(real rc, real epsilon, real sigma, real maxForce, Awareness awareness)¶ Constructor.
-
PairwiseRepulsiveLJ
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup the internal state of the functor
-
-
template <typename PressureEOS, typename DensityKernel>
classPairwiseSDPDHandler
: public mirheo::ParticleFetcherWithDensityAndMass¶ Compute smooth dissipative particle dynamics forces on the device.
- Template Parameters
PressureEos
: The equation of stateDensityJKernel
: The kernel used to compute the density
Subclassed by mirheo::PairwiseSDPD< PressureEOS, DensityKernel >
Public Functions
-
PairwiseSDPDHandler
(real rc, PressureEOS pressure, DensityKernel densityKernel, real viscosity)¶ Constructor.
-
real3
operator()
(const ParticleType dst, int dstId, const ParticleType src, int srcId) const¶ evaluate the force
-
ForceAccumulator
getZeroedAccumulator
() const¶ initialize the accumulator
-
template <typename PressureEOS, typename DensityKernel>
classPairwiseSDPD
: public mirheo::PairwiseKernel, public mirheo::PairwiseSDPDHandler<PressureEOS, DensityKernel>¶ Helper class to create PairwiseSDPDHandler from host.
- Template Parameters
PressureEos
: The equation of stateDensityJKernel
: The kernel used to compute the number density
Public Functions
-
PairwiseSDPD
(real rc, PressureEOS pressure, DensityKernel densityKernel, real viscosity, real kBT, long seed = 42424242)¶ Constructor.
-
PairwiseSDPD
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup the internal state of the functor
-
void
writeState
(std::ofstream &fout)¶ write internal state to a stream
-
bool
readState
(std::ifstream &fin)¶ restore internal state from a stream
-
class
PairwiseViscoElasticDPDHandler
: public mirheo::ParticleFetcherWithPolChainVectors¶ a GPU compatible functor that computes DPD interactions
Subclassed by mirheo::PairwiseViscoElasticDPD
Public Types
-
using
ViewType
= PVviewWithPolChainVector¶ compatible view type
-
using
ParticleType
= ParticleWithPolChainVector¶ compatible particle type
Public Functions
-
PairwiseViscoElasticDPDHandler
(real rc, real a, real gamma, real power, real H, real n0)¶ constructor
-
ForceDerPolChain
operator()
(const ParticleType dst, int dstId, const ParticleType src, int srcId) const¶ evaluate the force
-
ForceDerPolChainAccumulator
getZeroedAccumulator
() const¶ initialize accumulator
-
using
-
class
PairwiseViscoElasticDPD
: public mirheo::PairwiseKernel, public mirheo::PairwiseViscoElasticDPDHandler¶ Helper class that constructs PairwiseViscoElasticDPDHandler.
Public Types
-
using
HandlerType
= PairwiseViscoElasticDPDHandler¶ handler type corresponding to this object
-
using
ParamsType
= ViscoElasticDPDParams¶ parameters that are used to create this object
Public Functions
-
PairwiseViscoElasticDPD
(real rc, real a, real gamma, real kBT, real power, real H, real n0, long seed = 42424242)¶ Constructor.
-
PairwiseViscoElasticDPD
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup the internal state of the functor
-
void
writeState
(std::ofstream &fout)¶ write internal state to a stream
-
bool
readState
(std::ifstream &fin)¶ restore internal state from a stream
-
using
-
class
PairwiseViscoElasticSmoothVelDPDHandler
: public mirheo::ParticleFetcherWithPolChainVectorsAndSmoothVel¶ a GPU compatible functor that computes DPD + visco-elastic interactions
Subclassed by mirheo::PairwiseViscoElasticSmoothVelDPD
Public Types
-
using
ViewType
= PVviewWithPolChainVectorAndSmoothVelocity¶ compatible view type
-
using
ParticleType
= ParticleWithPolChainVectorAndSmoothVel¶ compatible particle type
Public Functions
-
PairwiseViscoElasticSmoothVelDPDHandler
(real rc, real a, real gamma, real power, real H, real n0)¶ constructor
-
ForceDerPolChain
operator()
(const ParticleType dst, int dstId, const ParticleType src, int srcId) const¶ evaluate the force
-
ForceDerPolChainAccumulator
getZeroedAccumulator
() const¶ initialize accumulator
-
using
-
class
PairwiseViscoElasticSmoothVelDPD
: public mirheo::PairwiseKernel, public mirheo::PairwiseViscoElasticSmoothVelDPDHandler¶ Helper class that constructs PairwiseViscoElasticSmoothVelDPDHandler.
Public Types
-
using
HandlerType
= PairwiseViscoElasticSmoothVelDPDHandler¶ handler type corresponding to this object
-
using
ParamsType
= ViscoElasticDPDParams¶ parameters that are used to create this object
Public Functions
-
PairwiseViscoElasticSmoothVelDPD
(real rc, real a, real gamma, real kBT, real power, real H, real n0, long seed = 42424242)¶ Constructor.
-
PairwiseViscoElasticSmoothVelDPD
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
-
void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup the internal state of the functor
-
void
writeState
(std::ofstream &fout)¶ write internal state to a stream
-
bool
readState
(std::ifstream &fin)¶ restore internal state from a stream
-
using
The above kernels that output a force can be wrapped by the stress wrapper:
-
template <typename BasePairwiseForceHandler>
classPairwiseStressWrapperHandler
: public BasePairwiseForceHandler¶ Compute force and stress from a pairwise force kernel.
- Template Parameters
BasePairwiseForceHandler
: The underlying pairwise interaction handler (must output a force)
Public Functions
-
PairwiseStressWrapperHandler
(BasePairwiseForceHandler basicForceHandler)¶ Constructor.
-
__device__ auto
operator()
(const ParticleType dst, int dstId, const ParticleType src, int srcId) const¶ Evaluate the force and the stress.
-
auto
getZeroedAccumulator
() const¶ Initialize the accumulator.
-
template <typename BasePairwiseForce>
classPairwiseStressWrapper
: public BasePairwiseForce¶ Create PairwiseStressWrapperHandler from host.
- Template Parameters
BasePairwiseForceHandler
: The underlying pairwise interaction (must output a force)
Public Functions
-
PairwiseStressWrapper
(BasePairwiseForce basicForce)¶ Constructor.
-
PairwiseStressWrapper
(real rc, const ParamsType &p, long seed = 42424242)¶ Generic constructor.
-
void
setup
(LocalParticleVector *lpv1, LocalParticleVector *lpv2, CellList *cl1, CellList *cl2, const MirState *state)¶ setup internal state
-
const HandlerType &
handler
() const¶ get the handler that can be used on device
Fetchers¶
Fetchers are used to load the correct data needed by the pairwise kernels (e.g. the mirheo::PairwiseRepulsiveLJ
kernel needs only the positions while the mirheo::PairwiseSDPD
kernel needs also velocities and number densities).
-
class
ParticleFetcher
¶ fetch position, velocity and global id
Subclassed by mirheo::PairwiseDensity< DensityKernel >, mirheo::PairwiseDPDHandler, mirheo::PairwiseLJ, mirheo::PairwiseMorse< Awareness >, mirheo::PairwiseNoRandomDPD, mirheo::PairwiseRepulsiveLJHandler< Awareness >, mirheo::ParticleFetcherWithDensity, mirheo::ParticleFetcherWithPolChainVectors, mirheo::ParticleFetcherWithPolChainVectorsAndSmoothVel
Public Types
Public Functions
-
ParticleFetcher
(real rc)¶ - Parameters
rc
: cut-off radius
-
ParticleType
read
(const ViewType &view, int id) const¶ fetch the particle information
- Return
- Particle information
- Parameters
view
: The view pointing to the dataid
: The particle index
-
ParticleType
readNoCache
(const ViewType &view, int id) const¶ read particle information directly from the global memory (without going through the L1 or L2 cache) This may be beneficial if one want to maximize the cahe usage on a concurrent stream
-
void
readCoordinates
(ParticleType &p, const ViewType &view, int id) const¶ read the coordinates only (used for the first pass on the neighbors, discard cut-off radius)
-
void
readExtraData
(ParticleType &p, const ViewType &view, int id) const¶ read the additional data contained in the particle (other than coordinates)
-
bool
withinCutoff
(const ParticleType &src, const ParticleType &dst) const¶ - Return
true
if the particlessrc
anddst
are within the cut-off radius distance;false
otherwise.
-
real3
getPosition
(const ParticleType &p) const¶ Generic converter from the ParticleType type to the common
real3
coordinates.
-
int64_t
getId
(const ParticleType &p) const¶ - Return
- Global id of the particle
-
-
class
ParticleFetcherWithDensity
: public mirheo::ParticleFetcher¶ fetch positions, velocities and number densities
Subclassed by mirheo::PairwiseMDPDHandler, mirheo::ParticleFetcherWithDensityAndMass
Public Types
-
using
ViewType
= PVviewWithDensities¶ compatible view type
-
using
ParticleType
= ParticleWithDensity¶ compatible particle type
Public Functions
-
ParticleFetcherWithDensity
(real rc)¶ - Parameters
rc
: cut-off radius
-
ParticleType
read
(const ViewType &view, int id) const¶ read full particle information
-
ParticleType
readNoCache
(const ViewType &view, int id) const¶ read full particle information through global memory
-
void
readCoordinates
(ParticleType &p, const ViewType &view, int id) const¶ read particle coordinates only
-
void
readExtraData
(ParticleType &p, const ViewType &view, int id) const¶ read velocity and number density of the particle
-
bool
withinCutoff
(const ParticleType &src, const ParticleType &dst) const¶ - Return
true
ifsrc
anddst
are within a cut-off radius distance;false
otherwise
-
real3
getPosition
(const ParticleType &p) const¶ fetch position from the generic particle structure
-
int64_t
getId
(const ParticleType &p) const¶ - Return
- Global id of the particle
-
struct
ParticleWithDensity
¶ contains position, global index, velocity and number density of a particle
-
using
-
class
ParticleFetcherWithDensityAndMass
: public mirheo::ParticleFetcherWithDensity¶ fetch that reads positions, velocities, number densities and mass
Subclassed by mirheo::PairwiseSDPDHandler< PressureEOS, DensityKernel >
Public Types
-
using
ViewType
= PVviewWithDensities¶ Compatible view type.
-
using
ParticleType
= ParticleWithDensityAndMass¶ Compatible particle type.
Public Functions
-
ParticleFetcherWithDensityAndMass
(real rc)¶ - Parameters
rc
: The cut-off radius
-
ParticleType
read
(const ViewType &view, int id) const¶ read full particle information
-
ParticleType
readNoCache
(const ViewType &view, int id) const¶ read full particle information through global memory
-
void
readCoordinates
(ParticleType &p, const ViewType &view, int id) const¶ read particle coordinates only
-
void
readExtraData
(ParticleType &p, const ViewType &view, int id) const¶ read velocity, number density and mass of the particle
-
bool
withinCutoff
(const ParticleType &src, const ParticleType &dst) const¶ - Return
true
ifsrc
anddst
are within a cut-off radius distance;false
otherwise
-
real3
getPosition
(const ParticleType &p) const¶ fetch position from the generic particle structure
-
int64_t
getId
(const ParticleType &p) const¶ - Return
- Global id of the particle
-
struct
ParticleWithDensityAndMass
¶ contains position, velocity, global id, number density and mass of a particle
-
using
-
class
ParticleFetcherWithPolChainVectors
: public mirheo::ParticleFetcher¶ fetch positions, velocities and polymeric chain vectors
Subclassed by mirheo::PairwiseViscoElasticDPDHandler
Public Types
-
using
ViewType
= PVviewWithPolChainVector¶ compatible view type
-
using
ParticleType
= ParticleWithPolChainVector¶ compatible particle type
Public Functions
-
ParticleFetcherWithPolChainVectors
(real rc)¶ - Parameters
rc
: cut-off radius
-
ParticleType
read
(const ViewType &view, int id) const¶ read full particle information
-
ParticleType
readNoCache
(const ViewType &view, int id) const¶ read full particle information through global memory
-
void
readCoordinates
(ParticleType &p, const ViewType &view, int id) const¶ read particle coordinates only
-
void
readExtraData
(ParticleType &p, const ViewType &view, int id) const¶ read velocity and number density of the particle
-
bool
withinCutoff
(const ParticleType &src, const ParticleType &dst) const¶ - Return
true
ifsrc
anddst
are within a cut-off radius distance;false
otherwise
-
real3
getPosition
(const ParticleType &p) const¶ fetch position from the generic particle structure
-
int64_t
getId
(const ParticleType &p) const¶ - Return
- Global id of the particle
-
struct
ParticleWithPolChainVector
¶ contains position, global index, velocity and polymeric chain vector of a particle
-
using
-
class
ParticleFetcherWithPolChainVectorsAndSmoothVel
: public mirheo::ParticleFetcher¶ fetch positions, velocities, polymeric chain vectors and smooth velocities
Subclassed by mirheo::PairwiseViscoElasticSmoothVelDPDHandler
Public Types
-
using
ViewType
= PVviewWithPolChainVectorAndSmoothVelocity¶ compatible view type
-
using
ParticleType
= ParticleWithPolChainVectorAndSmoothVel¶ compatible particle type
Public Functions
-
ParticleFetcherWithPolChainVectorsAndSmoothVel
(real rc)¶ - Parameters
rc
: cut-off radius
-
ParticleType
read
(const ViewType &view, int id) const¶ read full particle information
-
ParticleType
readNoCache
(const ViewType &view, int id) const¶ read full particle information through global memory
-
void
readCoordinates
(ParticleType &p, const ViewType &view, int id) const¶ read particle coordinates only
-
void
readExtraData
(ParticleType &p, const ViewType &view, int id) const¶ read velocity and number density of the particle
-
bool
withinCutoff
(const ParticleType &src, const ParticleType &dst) const¶ - Return
true
ifsrc
anddst
are within a cut-off radius distance;false
otherwise
-
real3
getPosition
(const ParticleType &p) const¶ fetch position from the generic particle structure
-
int64_t
getId
(const ParticleType &p) const¶ - Return
- Global id of the particle
-
struct
ParticleWithPolChainVectorAndSmoothVel
¶ contains position, global index, velocity and polymeric chain vector of a particle
-
using
Accumulators¶
Every interaction kernel must initialize an accumulator that is used to add its output quantity. Depending on the kernel, that quantity may be of different type, and may behave in a different way (e.g. forces and stresses are different).
It must satisfy the following interface requirements (in the following, we denote the type of the local variable as LType
and the view type as ViewType
):
A default constructor which initializes the internal local variable
Atomic accumulator from local value to destination view:
__D__ void atomicAddToDst(LType, ViewType&, int id) const;
Atomic accumulator from local value to source view:
__D__ inline void atomicAddToSrc(LType, ViewType&, int id) const;
Accessor of accumulated value:
__D__ inline LType get() const;
Function to add a value to the accumulator (from output of pairwise kernel):
__D__ inline void add(LType);
The following accumulators are currently implemented:
-
class
DensityAccumulator
¶ Accumulate densities on device.
Public Functions
-
DensityAccumulator
()¶ Initialize the DensityAccumulator.
-
void
atomicAddToDst
(real d, PVviewWithDensities &view, int id) const¶ Atomically add density
d
to the destinationview
at idid
.- Parameters
d
: The value to addview
: The destination containerid
: destination index inview
-
void
atomicAddToSrc
(real d, PVviewWithDensities &view, int id) const¶ Atomically add density
d
to the sourceview
at idid
.- Parameters
d
: The value to addview
: The destination containerid
: destination index inview
-
real
get
() const¶ - Return
- the internal accumulated density
-
void
add
(real d)¶ add
d
to the internal density
-
-
class
ForceAccumulator
¶ Accumulate forces on device.
Public Functions
-
ForceAccumulator
()¶ Initialize the ForceAccumulator.
-
real3
get
() const¶ - Return
- the internal accumulated force
-
void
add
(real3 f)¶ add
f
to the internal force
-
-
class
ForceDerPolChainAccumulator
¶ Accumulate ForceDerPolChain structure on device.
Public Functions
-
ForceDerPolChainAccumulator
()¶ Initialize the ForceDerPolChainAccumulator.
-
ForceDerPolChain
get
() const¶ - Return
- the internal accumulated force and stress
-
void
add
(const ForceDerPolChain &fq)¶ add
fq
to the internal value (only to the dst particle)
Public Static Functions
-
static void
atomicAddToDst
(const ForceDerPolChain &fq, PVviewWithPolChainVector &view, int id)¶ Atomically add the force and polymeric chain vector time derivative
fq
to the destinationview
at idid
.- Parameters
fq
: The force and dQdtview
: The destination containerid
: destination index inview
-
static void
atomicAddToSrc
(const ForceDerPolChain &fq, PVviewWithPolChainVector &view, int id)¶ Atomically add the force and polymeric chain vector time derivative
fq
to the sourceview
at idid
.- Parameters
fq
: The force and dQdtview
: The destination containerid
: destination index inview
-
-
template <class ForceType>
structForceStress
¶ Holds generalized force and stress together.
- Template Parameters
ForceType
: The type of the generalized force.
-
template <class BaseView, class BaseAccumulator>
classForceStressAccumulator
¶ Accumulate ForceStress structure on device.
- Template Parameters
BaseView
: The view type without stress, to enforce the use of the stress view wrapperBaseAccumulator
: The accumulator used to accumulate the force of the base interaction
Public Functions
-
ForceStressAccumulator
()¶ Initialize the ForceStressAccumulator.
-
void
atomicAddToDst
(const ForceStress<BaseForceType> &fs, PVviewWithStresses<BaseView> &view, int id) const¶ Atomically add the force and stress
fs
to the destinationview
at idid
.- Parameters
fs
: The force, directed from src to dst, and the corresponding stressview
: The destination containerid
: destination index inview
-
void
atomicAddToSrc
(const ForceStress<BaseForceType> &fs, PVviewWithStresses<BaseView> &view, int id) const¶ Atomically add the force and stress
fs
to the sourceview
at idid
.- Parameters
fs
: The force, directed from src to dst, and the corresponding stressview
: The destination containerid
: destination index inview
-
ForceStress<BaseForceType>
get
() const¶ - Return
- the internal accumulated force and stress
-
void
add
(const ForceStress<BaseForceType> &fs)¶ add
fs
to the internal force