Particle Vectors

See also the user interface.

Particle Vectors

class ParticleVector : public mirheo::MirSimulationObject

Base particles container.

Holds two LocalParticleVector: local and halo. The local one contains the data present in the current subdomain. The halo one is used to exchange particle data with the neighboring ranks.

By default, contains positions, velocities, forces and global ids.

Subclassed by mirheo::ObjectVector

Unnamed Group

void setCoordinates_vector(const std::vector<real3> &coordinates)

Python getters / setters Use default blocking stream.

Public Functions

ParticleVector(const MirState *state, const std::string &name, real mass, int n = 0)

Construct a ParticleVector.

Parameters
  • state: The simulation state
  • name: Name of the pv
  • mass: Mass of one particle
  • n: Number of particles

LocalParticleVector *local()

get the local LocalParticleVector

LocalParticleVector *halo()

get the halo LocalParticleVector

LocalParticleVector *get(ParticleVectorLocality locality)

get the LocalParticleVector corresponding to a given locality

Parameters
  • locality: local or halo

const LocalParticleVector *local() const

get the local LocalParticleVector

const LocalParticleVector *halo() const

get the halo LocalParticleVector

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.

template <typename T>
void requireDataPerParticle(const std::string &name, DataManager::PersistenceMode persistence, DataManager::ShiftMode shift = DataManager::ShiftMode::None)

Add a new channel to hold additional data per particle.

Template Parameters
  • T: The type of data to add
Parameters
  • name: channel name
  • persistence: If the data should stich to the particles or not when exchanged
  • shift: If the data needs to be shifted when exchanged

real getMassPerParticle() const

get the particle mass

Public Members

bool haloValid = {false}

true if the halo is up to date

bool redistValid = {false}

true if the particles are redistributed

int cellListStamp = {0}

stamp that keep track if the cell list is up to date

class ObjectVector : public mirheo::ParticleVector

Base objects container.

Holds two LocalObjectVector: local and halo.

Subclassed by mirheo::ChainVector, mirheo::MembraneVector, mirheo::RigidObjectVector, mirheo::RodVector

Public Functions

ObjectVector(const MirState *state, const std::string &name, real mass, int objSize, int nObjects = 0)

Construct a ObjectVector.

Parameters
  • state: The simulation state
  • name: Name of the pv
  • mass: Mass of one particle
  • objSize: Number of particles per object
  • nObjects: Number of objects

void findExtentAndCOM(cudaStream_t stream, ParticleVectorLocality locality)

Compute Extents and center of mass of each object in the given LocalObjectVector.

Parameters
  • stream: The stream to execute the kernel on.
  • locality: Specify which LocalObjectVector to compute the data

LocalObjectVector *local()

get local LocalObjectVector

LocalObjectVector *halo()

get halo LocalObjectVector

LocalObjectVector *get(ParticleVectorLocality locality)

get LocalObjectVector from locality

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.

template <typename T>
void requireDataPerObject(const std::string &name, DataManager::PersistenceMode persistence, DataManager::ShiftMode shift = DataManager::ShiftMode::None)

Add a new channel to hold additional data per object.

Template Parameters
  • T: The type of data to add
Parameters
  • name: channel name
  • persistence: If the data should stich to the objects or not when exchanging
  • shift: If the data needs to be shifted when exchanged

int getObjectSize() const

get number of particles per object

Public Members

std::shared_ptr<Mesh> mesh

Triangle mesh that can be used to represent the object surface.

class RigidObjectVector : public mirheo::ObjectVector

Rigid objects container.

Holds two LocalRigidObjectVector: local and halo.

Subclassed by mirheo::RigidShapedObjectVector< Shape >

Public Functions

RigidObjectVector(const MirState *state, const std::string &name, real partMass, real3 J, const int objSize, std::shared_ptr<Mesh> mesh, const int nObjects = 0)

Construct a RigidObjectVector.

Parameters
  • state: The simulation state
  • name: Name of the pv
  • partMass: Mass of one frozen particle
  • J: Diagonal entries of the inertia tensor, which must be diagonal.
  • objSize: Number of particles per object
  • mesh: Mesh representing the surface of the object.
  • nObjects: Number of objects

LocalRigidObjectVector *local()

get local LocalRigidObjectVector

LocalRigidObjectVector *halo()

get halo LocalRigidObjectVector

LocalRigidObjectVector *get(ParticleVectorLocality locality)

get LocalRigidObjectVector from locality

real3 getInertialTensor() const

get diagonal entries of the inertia tensor

Public Members

PinnedBuffer<real4> initialPositions

Coordinates of the frozen particles in the frame of reference of the object.

template <class Shape>
class RigidShapedObjectVector : public mirheo::RigidObjectVector

RigidObjectVector with analytic shape instead of triangle mesh.

Template Parameters
  • Shape: The analytic shape that represents the object surface in its frame of reference

Public Functions

RigidShapedObjectVector(const MirState *state, const std::string &name, real mass, int objSize, Shape shape, int nObjects = 0)

Construct a RigidShapedObjectVector.

Parameters
  • state: The simulation state
  • name: Name of the pv
  • mass: Mass of one frozen particle
  • objSize: Number of particles per object
  • shape: The shape that represents the surface of the object
  • nObjects: Number of objects

RigidShapedObjectVector(const MirState *state, const std::string &name, real mass, int objSize, Shape shape, std::shared_ptr<Mesh> mesh, int nObjects = 0)

Construct a RigidShapedObjectVector.

Note

The mesh is used only for visualization purpose

Parameters
  • state: The simulation state
  • name: Name of the pv
  • mass: Mass of one frozen particle
  • objSize: Number of particles per object
  • shape: The shape that represents the surface of the object
  • mesh: The mesh that represents the surface, should not used in the simulation.
  • nObjects: Number of objects

const Shape &getShape() const

get the handler that represent the shape of the objects

class RodVector : public mirheo::ObjectVector

Rod objects container.

Holds two LocalRodVector: local and halo.

Public Functions

RodVector(const MirState *state, const std::string &name, real mass, int nSegments, int nObjects = 0)

Construct a ObjectVector.

Parameters
  • state: The simulation state
  • name: Name of the pv
  • mass: Mass of one particle
  • nSegments: Number of segments per rod
  • nObjects: Number of rods

LocalRodVector *local()

get local LocalRodVector

LocalRodVector *halo()

get halo LocalRodVector

LocalRodVector *get(ParticleVectorLocality locality)

get LocalRodVector from locality

template <typename T>
void requireDataPerBisegment(const std::string &name, DataManager::PersistenceMode persistence, DataManager::ShiftMode shift = DataManager::ShiftMode::None)

Add a new channel to hold additional data per bisegment.

Template Parameters
  • T: The type of data to add
Parameters
  • name: channel name
  • persistence: If the data should stich to the object or not when exchanging
  • shift: If the data needs to be shifted when exchanged

class MembraneVector : public mirheo::ObjectVector

Represent a set of membranes.

Each membrane is composed of the same connectivity (stored in mesh) and number of vertices. The particles data correspond to the vertices of the membranes.

Public Functions

MembraneVector(const MirState *state, const std::string &name, real mass, std::shared_ptr<MembraneMesh> mptr, int nObjects = 0)

Construct a MembraneVector.

Parameters
  • state: The simulation state
  • name: Name of the pv
  • mass: Mass of one particle
  • mptr: Triangle mesh which stores the connectivity of a single membrane
  • nObjects: Number of objects

Local Particle Vectors

class LocalParticleVector

Particles container.

This is used to represent local or halo particles in ParticleVector.

Subclassed by mirheo::LocalObjectVector

Public Functions

LocalParticleVector(ParticleVector *pv, int np = 0)

Construct a LocalParticleVector.

Parameters

int size() const

return the number of particles

virtual void resize(int n, cudaStream_t stream)

resize the container, preserving the data.

Parameters
  • n: new number of particles
  • stream: that is used to copy data

virtual void resize_anew(int n)

resize the container, without preserving the data.

Parameters
  • n: new number of particles

PinnedBuffer<Force> &forces()

get forces container reference

PinnedBuffer<real4> &positions()

get positions container reference

PinnedBuffer<real4> &velocities()

get velocities container reference

virtual void computeGlobalIds(MPI_Comm comm, cudaStream_t stream)

Set a unique Id for each particle in the simulation.

The ids are stored in the channel ChannelNames::globalIds.

Parameters
  • comm: MPI communicator of the simulation
  • stream: Stream used to transfer data between host and device

ParticleVector *parent()

get parent ParticleVector

const ParticleVector *parent() const

get parent ParticleVector

Public Members

DataManager dataPerParticle

Contains all particle channels.

Friends

void swap(LocalParticleVector&, LocalParticleVector&)

swap two LocalParticleVector

class LocalObjectVector : public mirheo::LocalParticleVector

Objects container.

This is used to represent local or halo objects in ObjectVector. An object is a chunk of particles, each chunk with the same number of particles within an ObjectVector. Additionally, data can be attached to each of those chunks.

Subclassed by mirheo::LocalRigidObjectVector, mirheo::LocalRodVector

Public Functions

LocalObjectVector(ParticleVector *pv, int objSize, int nObjects = 0)

Construct a LocalParticleVector.

Parameters
  • pv: Parent ObjectVector
  • objSize: Number of particles per object
  • nObjects: Number of objects

void resize(int n, cudaStream_t stream)

resize the container, preserving the data.

Parameters
  • n: new number of particles
  • stream: that is used to copy data

void resize_anew(int n)

resize the container, without preserving the data.

Parameters
  • n: new number of particles

void computeGlobalIds(MPI_Comm comm, cudaStream_t stream)

Set a unique Id for each particle in the simulation.

The ids are stored in the channel ChannelNames::globalIds.

Parameters
  • comm: MPI communicator of the simulation
  • stream: Stream used to transfer data between host and device

virtual PinnedBuffer<real4> *getMeshVertices(cudaStream_t stream)

get positions of the mesh vertices

virtual PinnedBuffer<real4> *getOldMeshVertices(cudaStream_t stream)

get positions of the old mesh vertices

virtual PinnedBuffer<Force> *getMeshForces(cudaStream_t stream)

get forces on the mesh vertices

int getObjectSize() const

get number of particles per object

int getNumObjects() const

get number of objects

Public Members

DataManager dataPerObject

contains object data

Friends

void swap(LocalObjectVector&, LocalObjectVector&)

swap two LocalObjectVector

class LocalRigidObjectVector : public mirheo::LocalObjectVector

Rigid objects container.

This is used to represent local or halo objects in RigidObjectVector. A rigid object is composed of frozen particles inside a volume that is represented by a triangle mesh. There is then two sets of particles: mesh vertices and frozen particles. The frozen particles are stored in the particle data manager, while the mesh particles are stored in additional buffers.

Additionally, each rigid object has a RigidMotion datum associated that fully describes its state.

Public Functions

LocalRigidObjectVector(ParticleVector *pv, int objSize, int nObjects = 0)

Construct a LocalRigidObjectVector.

Parameters
  • pv: Parent RigidObjectVector
  • objSize: Number of frozen particles per object
  • nObjects: Number of objects

PinnedBuffer<real4> *getMeshVertices(cudaStream_t stream)

get positions of the mesh vertices

PinnedBuffer<real4> *getOldMeshVertices(cudaStream_t stream)

get positions of the old mesh vertices

PinnedBuffer<Force> *getMeshForces(cudaStream_t stream)

get forces on the mesh vertices

void clearRigidForces(cudaStream_t stream)

set forces in rigid motions to zero

class LocalRodVector : public mirheo::LocalObjectVector

Rod container.

This is used to represent local or halo rods in RodVector. A rod is a chunk of particles connected implicitly in segments with additional 4 particles per edge. The number of particles per rod is then 5*n + 1 if n is the number of segments. Each object (called a rod) within a LocalRodVector has the same number of particles. Additionally to particle and object, data can be attached to each bisegment.

Public Functions

LocalRodVector(ParticleVector *pv, int objSize, int nObjects = 0)

Construct a LocalRodVector.

Parameters
  • pv: Parent RodVector
  • objSize: Number of particles per object
  • nObjects: Number of objects

void resize(int n, cudaStream_t stream)

resize the container, preserving the data.

Parameters
  • n: new number of particles
  • stream: that is used to copy data

void resize_anew(int n)

resize the container, without preserving the data.

Parameters
  • n: new number of particles

int getNumSegmentsPerRod() const

get the number of segment per rod

Public Members

DataManager dataPerBisegment

contains bisegment data

Views

struct PVview

GPU-compatible struct that contains particle basic data.

Contains particle positions, velocities, forces, and mass info.

Subclassed by mirheo::OVview, mirheo::PVviewWithDensities, mirheo::PVviewWithOldParticles, mirheo::PVviewWithPolChainVector

Public Types

using PVType = ParticleVector

Particle Vector compatible type.

using LPVType = LocalParticleVector

Local Particle Vector compatible type.

Public Functions

PVview(ParticleVector *pv, LocalParticleVector *lpv)

Construct a PVview.

Parameters

real4 readPosition(int id) const

fetch position from given particle index

real4 readVelocity(int id) const

fetch velocity from given particle index

void readPosition(Particle &p, int id) const

fetch position from given particle index and store it into p

void readVelocity(Particle &p, int id) const

fetch velocity from given particle index and store it into p

Particle readParticle(int id) const

fetch particle from given particle index

real4 readPositionNoCache(int id) const

fetch position from given particle index without going through the L1/L2 cache This can be useful to reduce the cache pressure on concurrent kernels

Particle readParticleNoCache(int id) const

fetch particle from given particle index without going through the L1/L2 cache This can be useful to reduce the cache pressure on concurrent kernels

void writePosition(int id, const real4 &r)

Store position at the given particle id.

void writeVelocity(int id, const real4 &u)

Store velocity at the given particle id.

void writeParticle(int id, const Particle &p)

Store particle at the given particle id.

Public Members

int size = {0}

number of particles

real4 *positions = {nullptr}

particle positions in local coordinates

real4 *velocities = {nullptr}

particle velocities

real4 *forces = {nullptr}

particle forces

real mass = {0._r}

mass of one particle

real invMass = {0._r}

1 / mass

struct PVviewWithOldParticles : public mirheo::PVview

PVview with additionally positions from previous time steps

Public Functions

PVviewWithOldParticles(ParticleVector *pv, LocalParticleVector *lpv)

Construct a PVviewWithOldParticles.

Note

if pv does not have old positions channel, this will be ignored and oldPositions will be set to nullptr.

Parameters

real3 readOldPosition(int id) const

fetch positions at previous time step

Public Members

real4 *oldPositions = {nullptr}

particle positions from previous time steps

struct PVviewWithDensities : public mirheo::PVview

PVview with additionally densities data

Public Functions

PVviewWithDensities(ParticleVector *pv, LocalParticleVector *lpv)

Construct a PVviewWithOldParticles.

Warning

The pv must hold a density channel.

Parameters

Public Members

real *densities = {nullptr}

particle densities

template <typename BasicView>
struct PVviewWithStresses : public BasicView

A View with additional stress info.

Template Parameters
  • BasicView: The pv view to extend with stresses

Public Functions

PVviewWithStresses(PVType *pv, LPVType *lpv)

Construct a PVviewWithStresses.

Warning

The pv must hold a stress per particle channel.

Parameters

Public Members

Stress *stresses = {nullptr}

stresses per particle

struct OVview : public mirheo::PVview

A PVview with additionally basic object data.

Contains object ids, object extents.

Subclassed by mirheo::OVviewWithAreaVolume, mirheo::OVviewWithNewOldVertices, mirheo::ROVview, mirheo::RVview

Public Types

using PVType = ObjectVector

Particle Vector compatible type.

using LPVType = LocalObjectVector

Local Particle Vector compatible type.

Public Functions

OVview(ObjectVector *ov, LocalObjectVector *lov)

Construct a OVview.

Parameters

Public Members

int nObjects = {0}

number of objects

int objSize = {0}

number of particles per object

real objMass = {0._r}

mass of one object

real invObjMass = {0._r}

1 / objMass

COMandExtent *comAndExtents = {nullptr}

center of mass and extents of the objects

int64_t *ids = {nullptr}

global ids of objects

struct OVviewWithAreaVolume : public mirheo::OVview

A OVview with additionally area and volume information.

Subclassed by mirheo::OVviewWithJuelicherQuants

Public Functions

OVviewWithAreaVolume(ObjectVector *ov, LocalObjectVector *lov)

Construct a OVviewWithAreaVolume.

Warning

The ov must hold a areaVolumes channel.

Parameters

Public Members

real2 *area_volumes = {nullptr}

area and volume per object

struct OVviewWithJuelicherQuants : public mirheo::OVviewWithAreaVolume

A OVviewWithAreaVolume with additional curvature related quantities.

Public Functions

OVviewWithJuelicherQuants(ObjectVector *ov, LocalObjectVector *lov)

Construct a OVviewWithJuelicherQuants.

Warning

The ov must hold areaVolumes and lenThetaTot object channels and vertex areas, meanCurvatures particle channels.

Parameters

Public Members

real *vertexAreas = {nullptr}

area per vertex (defined on a triangle mesh)

real *vertexMeanCurvatures = {nullptr}

mean curvature vertex (defined on a triangle mesh)

real *lenThetaTot = {nullptr}

helper quantity to compute Juelicher bending energy

struct OVviewWithNewOldVertices : public mirheo::OVview

A OVview with additionally vertices information.

Public Functions

OVviewWithNewOldVertices(ObjectVector *ov, LocalObjectVector *lov, cudaStream_t stream)

Construct a OVviewWithNewOldVertices.

Parameters
  • ov: The ObjectVector that the view represents
  • lov: The LocalObjectVector that the view represents
  • stream: Stream used to create mesh vertices if not already present

Public Members

real4 *vertices = {nullptr}

vertex positions

real4 *old_vertices = {nullptr}

vertex positions at previous time step

real4 *vertexForces = {nullptr}

vertex forces

int nvertices = {0}

number of vertices

struct ROVview : public mirheo::OVview

A OVview with additional rigid object infos.

Subclassed by mirheo::ROVviewWithOldMotion, mirheo::RSOVview< Shape >

Public Functions

ROVview(RigidObjectVector *rov, LocalRigidObjectVector *lrov)

Construct a ROVview.

Parameters

Public Members

RigidMotion *motions = {nullptr}

rigid object states

real3 J = {0._r, 0._r, 0._r}

diagonal entries of inertia tensor

real3 J_1 = {0._r ,0._r, 0._r}

diagonal entries of the inverse inertia tensor

struct ROVviewWithOldMotion : public mirheo::ROVview

A OVview with additional rigid object info from previous time step.

Public Functions

ROVviewWithOldMotion(RigidObjectVector *rov, LocalRigidObjectVector *lrov)

Construct a ROVview.

Warning

The rov must hold old motions channel.

Parameters

Public Members

RigidMotion *old_motions = {nullptr}

rigid object states at previous time step

template <class Shape>
struct RSOVview : public mirheo::ROVview

A ROVview with additional analytic shape infos.

Template Parameters
  • Shape: the analytical shape that represents the object shape

Subclassed by mirheo::RSOVviewWithOldMotion< Shape >

Public Functions

RSOVview(RigidShapedObjectVector<Shape> *rsov, LocalRigidObjectVector *lrov)

Construct a RSOVview.

Parameters

Public Members

Shape shape

Represents the object shape.

template <class Shape>
struct RSOVviewWithOldMotion : public mirheo::RSOVview<Shape>

A RSOVview with additional rigid object info from previous time step.

Template Parameters
  • Shape: the analytical shape that represents the object shape

Public Functions

RSOVviewWithOldMotion(RigidShapedObjectVector<Shape> *rsov, LocalRigidObjectVector *lrov)

Construct a RSOVview.

Warning

The rov must hold old motions channel.

Parameters

Public Members

RigidMotion *old_motions = {nullptr}

rigid object states at previous time step

struct RVview : public mirheo::OVview

A OVview with additional rod object infos.

Subclassed by mirheo::RVviewWithOldParticles

Public Functions

RVview(RodVector *rv, LocalRodVector *lrv)

Construct a RVview.

Parameters

Public Members

int nSegments = {0}

number of segments per rod

int *states = {nullptr}

polymorphic states per bisegment

real *energies = {nullptr}

energies per bisegment

struct RVviewWithOldParticles : public mirheo::RVview

A RVview with additional particles from previous time step.

Public Functions

RVviewWithOldParticles(RodVector *rv, LocalRodVector *lrv)

Construct a RVview.

Parameters

Public Members

real4 *oldPositions = {nullptr}

positions o the particles at previous time step

Data Manager

This is the building block to create particle vectors.

class DataManager

Container for multiple channels on device and host.

Used by ParticleVector and ObjectVector to hold data per particle and per object correspondingly. All channels are stored as PinnedBuffer, which allows to easily transfer the data between host and device. Channels can hold data of types listed in VarPinnedBufferPtr variant. See ChannelDescription for the description of one channel.

Unnamed Group

DataManager(const DataManager &b)

copy and move constructors

Unnamed Group

ChannelDescription *getChannelDesc(const std::string &name)

Get channel from its name or nullptr if it is not found.

Unnamed Group

ChannelDescription &getChannelDescOrDie(const std::string &name)

Get channel from its name or die if it is not found.

Public Types

using NamedChannelDesc = std::pair<std::string, const ChannelDescription *>

The full description of a channel, contains its name and description.

Public Functions

void copyChannelMap(const DataManager&)

Copy channel names and their types from a given DataManager.

Does not copy data or resize buffers. New buffers are empty.

template <typename T>
void createData(const std::string &name, int size = 0)

Allocate a new channel.

This method will die if a channel with different type but same name already exists. If a channel with the same name and same type exists, this method will not allocate a new channel.

Template Parameters
  • T: datatype of the buffer element. sizeof(T) should be compatible with VarPinnedBufferPtr
Parameters
  • name: buffer name
  • size: resize buffer to size elements

void setPersistenceMode(const std::string &name, PersistenceMode persistence)

Set the persistence mode of the data.

This method will die if the required name does not exist.

Warning

This method can only increase the persistence. If the channel is already persistent, this method can not set its persistent mode to None.

Parameters
  • name: The name of the channel to modify
  • persistence: Persistence mode to add to the channel.

void setShiftMode(const std::string &name, ShiftMode shift)

Set the shift mode of the data.

This method will die if the required name does not exist.

Warning

This method can only increase the shift mode. If the channel already needs shift, this method can not set its shift mode to None.

Parameters
  • name: The name of the channel to modify
  • shift: Shift mode to add to the channel.

GPUcontainer *getGenericData(const std::string &name)

Get gpu buffer by name.

This method will die if the required name does not exist.

Return
pointer to GPUcontainer corresponding to the given name
Parameters
  • name: buffer name

template <typename T>
PinnedBuffer<T> *getData(const std::string &name)

Get buffer by name.

This method will die if the required name does not exist or if T is of the wrong type.

Return
pointer to PinnedBuffer<T> corresponding to the given name
Parameters
  • name: buffer name
Template Parameters

void *getGenericPtr(const std::string &name)

Get device buffer pointer regardless of its type.

This method will die if the required name does not exist.

Return
pointer to device data held by the corresponding PinnedBuffer
Parameters
  • name: buffer name

bool checkChannelExists(const std::string &name) const

true if channel with given name exists, false otherwise

const std::vector<NamedChannelDesc> &getSortedChannels() const

Return
vector of channels sorted (descending) by size of their elements (and then name)

bool checkPersistence(const std::string &name) const

Return
true if the channel is persistent

void resize(int n, cudaStream_t stream)

Resize all the channels and preserve their data.

void resize_anew(int n)

Resize all the channels without preserving the data.

Friends

void swap(DataManager &a, DataManager &b)

swap two DataManager

struct ChannelDescription

Holds information and data of a single channel.

A channel has a type, persistence mode and shift mode.

Public Functions

bool needShift() const

returns true if the channel’s data needs to be shifted when exchanged or redistributed.

Public Members

std::unique_ptr<GPUcontainer> container

The data stored in the channel. Internally stored as a PinnedBuffer.

VarPinnedBufferPtr varDataPtr

Pointer to container that holds the correct type.

PersistenceMode persistence = {PersistenceMode::None}

The persistence mode of the channel.

ShiftMode shift = {ShiftMode::None}

The shift mode of the channel.