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::PairwiseInteraction< PairwiseKernel >, mirheo::PairwiseInteractionWithStress< PairwiseKernel >, mirheo::PairwiseInteraction< mirheo::PairwiseKernel >, mirheo::PairwiseInteraction< mirheo::PairwiseStressWrapper< mirheo::PairwiseKernel > >

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.

Implementation

The factory instantiates one of this templated class. See below for the requirements on the kernels.

template <class PairwiseKernel>
class PairwiseInteraction : public mirheo::BasePairwiseInteraction

Short-range symmetric pairwise interactions.

See the pairwise interaction entry of the developer documentation for the interface requirements of the kernel.

Template Parameters
  • PairwiseKernel: The functor that describes the interaction between two particles (interaction kernel).

Public Functions

PairwiseInteraction(const MirState *state, const std::string &name, real rc, KernelParams pairParams, long seed = 42424242)

Construct a PairwiseInteraction object.

Parameters
  • state: The global state of the system
  • name: The name of the interaction
  • rc: The cut-off radius of the interaction
  • pairParams: The parameters used to construct the interaction kernel
  • seed: used to initialize random number generator (needed to construct some interaction kernels).

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

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 ParticleVector
  • pv2: 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 pv1
  • cl2: cell-list built for the appropriate cut-off radius for pv2
  • stream: 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 ParticleVector
  • pv2: 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 pv1
  • cl2: cell-list built for the appropriate cut-off radius for pv2
  • stream: Execution stream

Stage getStage() const

returns the Stage corresponding to this interaction.

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.

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.

Public Static Functions

static std::string getTypeName()

Return
A string that describes the type of this object

A specific class can be used to compute addtionally the stresses of a given interaction.

template <class PairwiseKernel>
class PairwiseInteractionWithStress : public mirheo::BasePairwiseInteraction

Short-range symmetric pairwise interactions with stress output.

This object manages two interaction: one with stress, which is used every stressPeriod time, and one with no stress wrapper, that is used the rest of the time. This is motivated by the fact that stresses are not needed for the simulation but rather for post processing; thus the stresses may not need to be computed at every time step.

Template Parameters
  • PairwiseKernel: The functor that describes the interaction between two particles (interaction kernel).

Public Functions

PairwiseInteractionWithStress(const MirState *state, const std::string &name, real rc, real stressPeriod, KernelParams pairParams, long seed = 42424242)

Construct a PairwiseInteractionWithStress object.

Parameters
  • state: The global state of the system
  • name: The name of the interaction
  • rc: The cut-off radius of the interaction
  • stressPeriod: The simulation time between two stress computation
  • pairParams: The parameters used to construct the interaction kernel
  • seed: used to initialize random number generator (needed to construct some interaction kernels).

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

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 ParticleVector
  • pv2: 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 pv1
  • cl2: cell-list built for the appropriate cut-off radius for pv2
  • stream: 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 ParticleVector
  • pv2: 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 pv1
  • cl2: cell-list built for the appropriate cut-off radius for pv2
  • stream: 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.

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.

Public Static Functions

static std::string getTypeName()

Return
A string that describes the type of this object

Kernels

Interface

The mirheo::PairwiseInteraction takes a functor that describes a pairwise interaction. This functor may be 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:

  1. 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>
    
  2. A generic constructor from the ParamsType parameters:

    PairwiseKernelType(real rc, const ParamsType& p, real dt, long seed=42424242);
    
  3. Setup function (on Host, for manager only)

    void setup(LocalParticleVector* lpv1, LocalParticleVector* lpv2, CellList* cl1, CellList* cl2, const MirState *state);
    
  4. Handler function (on Host, for manager only)

    const HandlerType& handler() const;
    
  5. 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;
    
  6. Accumulator initializer (on GPU)

    __D__ <Accumulator> getZeroedAccumulator() const;
    
  7. 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;
    
  8. Interacting checker to discard pairs not within cutoff:

    __D__ bool withinCutoff(const ParticleType& src, const ParticleType& dst) const;
    
  9. 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 - add it to the variant 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::PairwiseStressWrapper< mirheo::PairwiseKernel >

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

The rest is directly implemented in the kernels, as no virtual functions are allowed on the device.

Implemented kernels

template <typename DensityKernel>
class PairwiseDensity : 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

using ViewType = PVview

compatible view type

using ParticleType = Particle

compatible particle type

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

Public Static Functions

static std::string getTypeName()

Return
type name string

class PairwiseLJ : public mirheo::PairwiseKernel, public mirheo::ParticleFetcher

Compute Lennard-Jones forces on the device.

Public Types

using ViewType = PVview

Compatible view type.

using ParticleType = Particle

Compatible particle type.

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

Public Static Functions

static std::string getTypeName()

Return
type name string

template <class Awareness>
class PairwiseMorse : 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, Awareness awareness)

Constructor.

PairwiseMorse(real rc, const ParamsType &p, long seed)

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

Public Static Functions

static std::string getTypeName()

Return
type name string

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

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

Public Static Functions

static std::string getTypeName()

Return
type name string

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 ViewType = PVview

compatible view type

using ParticleType = Particle

compatible particle type

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

Public Static Functions

static std::string getTypeName()

Return
type name string

template <class Awareness>
class PairwiseRepulsiveLJ : 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>
class PairwiseSDPDHandler : public mirheo::ParticleFetcherWithDensityAndMass

Compute smooth dissipative particle dynamics forces on the device.

Template Parameters
  • PressureEos: The equation of state
  • DensityJKernel: 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>
class PairwiseSDPD : public mirheo::PairwiseKernel, public mirheo::PairwiseSDPDHandler<PressureEOS, DensityKernel>

Helper class to create PairwiseSDPDHandler from host.

Template Parameters
  • PressureEos: The equation of state
  • DensityJKernel: 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

The above kernels that output a force can be wrapped by the stress wrapper:

template <typename BasicPairwiseForceHandler>
class PairwiseStressWrapperHandler : public BasicPairwiseForceHandler

Compute force and stress from a pairwise force kernel.

Template Parameters
  • BasicPairwiseForceHandler: The underlying pairwise interaction handler (must output a force)

Public Functions

PairwiseStressWrapperHandler(BasicPairwiseForceHandler basicForceHandler)

Constructor.

__device__ ForceStress operator()(const ParticleType dst, int dstId, const ParticleType src, int srcId) const

Evaluate the force and the stress.

ForceStressAccumulator<BasicViewType> getZeroedAccumulator() const

Initialize the accumulator.

template <typename BasicPairwiseForce>
class PairwiseStressWrapper : public BasicPairwiseForce

Create PairwiseStressWrapperHandler from host.

Template Parameters
  • BasicPairwiseForceHandler: The underlying pairwise interaction (must output a force)

Public Functions

PairwiseStressWrapper(BasicPairwiseForce 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

Public Types

using ViewType = PVview

compatible view type

using ParticleType = Particle

compatible particle type

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 data
  • id: 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 particles src and dst 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 if src and dst 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

Public Members

Particle p

positions, global id, velocity

real d

number density

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 if src and dst 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

Public Members

Particle p

position, global id, velocity

real d

number density

real m

mass

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):

  1. A default constructor which initializes the internal local variable

  2. Atomic accumulator from local value to destination view:

    __D__ void atomicAddToDst(LType, ViewType&, int id) const;
    
  3. Atomic accumulator from local value to source view:

    __D__ inline void atomicAddToSrc(LType, ViewType&, int id) const;
    
  4. Accessor of accumulated value:

    __D__ inline LType get() const;
    
  5. 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 destination view at id id.

Parameters
  • d: The value to add
  • view: The destination container
  • id: destination index in view

void atomicAddToSrc(real d, PVviewWithDensities &view, int id) const

Atomically add density d to the source view at id id.

Parameters
  • d: The value to add
  • view: The destination container
  • id: destination index in view

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.

void atomicAddToDst(real3 f, PVview &view, int id) const

Atomically add the force f to the destination view at id id.

Parameters
  • f: The force, directed from src to dst
  • view: The destination container
  • id: destination index in view

void atomicAddToSrc(real3 f, PVview &view, int id) const

Atomically add the force f to the source view at id id.

Parameters
  • f: The force, directed from src to dst
  • view: The destination container
  • id: destination index in view

real3 get() const

Return
the internal accumulated force

void add(real3 f)

add f to the internal force

struct ForceStress

Holds force and stress together.

Public Members

real3 force

force value

Stress stress

stress value

template <typename BasicView>
class ForceStressAccumulator

Accumulate ForceStress structure on device.

Template Parameters
  • BasicView: The view type without stress, to enforce the use of the stress view wrapper

Public Functions

ForceStressAccumulator()

Initialize the ForceStressAccumulator.

void atomicAddToDst(const ForceStress &fs, PVviewWithStresses<BasicView> &view, int id) const

Atomically add the force and stress fs to the destination view at id id.

Parameters
  • fs: The force, directed from src to dst, and the corresponding stress
  • view: The destination container
  • id: destination index in view

void atomicAddToSrc(const ForceStress &fs, PVviewWithStresses<BasicView> &view, int id) const

Atomically add the force and stress fs to the source view at id id.

Parameters
  • fs: The force, directed from src to dst, and the corresponding stress
  • view: The destination container
  • id: destination index in view

ForceStress get() const

Return
the internal accumulated force and stress

void add(const ForceStress &fs)

add fs to the internal force