Membrane Interactions

Base class

This is the visible class that is output of the factory function.

class BaseMembraneInteraction : public mirheo::Interaction

Base class that represents membrane interactions.

This kind of interactions does not require any cell-lists and is always a “self-interaction”, hence the halo interaction does not do anything. This must be used only with MembraneVector objects.

Subclassed by mirheo::MembraneInteraction< TriangleInteraction, DihedralInteraction, Filter >

Public Functions

BaseMembraneInteraction(const MirState *state, const std::string &name)

Construct a BaseMembraneInteraction object.

Parameters
  • state: The global state of the system
  • name: The name of the interaction

void setPrerequisites(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)

Set the required channels to the concerned ParticleVector that will participate in the interactions.

This method will fail if pv1 is not a

MembraneVector or if pv1 is not the same as pv2.
Parameters
  • pv1: The conserned data that will participate in the interactions.
  • pv2: The conserned data that will participate in the interactions.
  • cl1: Unused
  • cl2: Unused

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

bool isSelfObjectInteraction() const

This is useful to know if we need exchange / cell-lists for that interaction. Example: membrane interactions are internal, all particles of a membrane are always on the same rank thus it does not need halo particles.

Return
boolean describing if the interaction is an internal interaction.

Implementation

The factory instantiates one of this templated class. See Triangle Kernels, Dihedral Kernels and Filters for possible template parameters.

template <class TriangleInteraction, class DihedralInteraction, class Filter>
class MembraneInteraction : public mirheo::BaseMembraneInteraction

Generic implementation of membrane forces.

Template Parameters
  • TriangleInteraction: Describes what forces are applied to triangles
  • DihedralInteraction: Describes what forces are applied to dihedrals
  • Filter: This allows to apply the interactions only to a subset of membranes

Public Functions

MembraneInteraction(const MirState *state, std::string name, CommonMembraneParameters parameters, typename TriangleInteraction::ParametersType triangleParams, typename DihedralInteraction::ParametersType dihedralParams, real initLengthFraction, real growUntil, Filter filter, long seed = 42424242)

Construct a MembraneInteraction object.

More information can be found on

growUntil in _scaleFromTime().
Parameters
  • state: The global state of the system
  • name: The name of the interaction
  • parameters: The common parameters from all kernel forces
  • triangleParams: Parameters that contain the parameters of the triangle forces kernel
  • dihedralParams: Parameters that contain the parameters of the dihedral forces kernel
  • initLengthFraction: The membrane will grow from this fraction of its size to its full size in growUntil time
  • growUntil: The membrane will grow from initLengthFraction fraction of its size to its full size in this amount of time
  • filter: Describes which membranes to apply the interactions
  • seed: Random seed for rng

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 setPrerequisites(ParticleVector *pv1, ParticleVector *pv2, CellList *cl1, CellList *cl2)

Set the required channels to the concerned ParticleVector that will participate in the interactions.

This method will fail if pv1 is not a

MembraneVector or if pv1 is not the same as pv2.
Parameters
  • pv1: The conserned data that will participate in the interactions.
  • pv2: The conserned data that will participate in the interactions.
  • cl1: Unused
  • cl2: Unused

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.

Triangle Kernels

Each thread is mapped to one vertex v1 and loops over all adjacent triangles labeled as follows:

graph triangle {
node [shape=plaintext]
{rank = same; v2; v1}
v3 -- v2
v3 -- v1
v2 -- v1
}

The output of the kernel is the forces of a given dihedral on v1. The forces on v2 and v3 from the same dihedral are computed by the thread mapped on v2 and v3, respectively.

template <StressFreeState stressFreeState>
class TriangleLimForce

Compute shear energy on a given triangle with the Lim model.

Template Parameters
  • stressFreeState: States if there is a stress free mesh associated with the interaction

Public Functions

TriangleLimForce(ParametersType p, const Mesh *mesh, mReal lscale)

Construct the functor.

Parameters
  • p: The parameters of the model
  • mesh: Triangle mesh information
  • lscale: Scaling length factor, applied to all parameters

void applyLengthScalingFactor(mReal lscale)

Scale length-dependent parameters.

EquilibriumTriangleDesc getEquilibriumDesc(const MembraneMeshView &mesh, int i0, int i1) const

Get the reference triangle information.

Return
The reference triangle information.
Parameters
  • mesh: Mesh view that contains the reference mesh. Only used when stressFreeState is Active.
  • i0: Index (in the adjacent vertex ids space, see Mesh) of the first adjacent vertex
  • i1: Index (in the adjacent vertex ids space, see Mesh) of the second adjacent vertex

mReal3 operator()(mReal3 v1, mReal3 v2, mReal3 v3, EquilibriumTriangleDesc eq) const

Compute the triangle force on v1.

See Developer docs for more details.

Return
The triangle force acting on v1
Parameters
  • v1: vertex 1
  • v2: vertex 2
  • v3: vertex 3
  • eq: The reference triangle information

template <StressFreeState stressFreeState>
class TriangleWLCForce

Compute shear energy on a given triangle with the Lim model.

Template Parameters
  • stressFreeState: States if there is a stress free mesh associated with the interaction

Public Functions

TriangleWLCForce(ParametersType p, const Mesh *mesh, mReal lscale)

Construct the functor.

Parameters
  • p: The parameters of the model
  • mesh: Triangle mesh information
  • lscale: Scaling length factor, applied to all parameters

void applyLengthScalingFactor(mReal lscale)

Scale length-dependent parameters.

EquilibriumTriangleDesc getEquilibriumDesc(const MembraneMeshView &mesh, int i0, int i1) const

Get the reference triangle information.

Return
The reference triangle information.
Parameters
  • mesh: Mesh view that contains the reference mesh. Only used when stressFreeState is Active.
  • i0: Index (in the adjacent vertex ids space, see Mesh) of the first adjacent vertex
  • i1: Index (in the adjacent vertex ids space, see Mesh) of the second adjacent vertex

mReal3 operator()(mReal3 v1, mReal3 v2, mReal3 v3, EquilibriumTriangleDesc eq) const

Compute the triangle force on v1.

See Developer docs for more details.

Return
The triangle force acting on v1
Parameters
  • v1: vertex 1
  • v2: vertex 2
  • v3: vertex 3
  • eq: The reference triangle information

Dihedral Kernels

Each thread is mapped to one vertex v0 and loops over all adjacent dihedrals labeled as follows:

graph dihedral {
node [shape=plaintext]
{rank = same; v2 ; v0}
v3 -- v2
v3 -- v0
v2 -- v0
v2 -- v1
v0 -- v1
}

The output of the kernel is the forces of a given dihedral on v0 and v1. The forces on v2 and v3 from the same dihedral are computed by the thread mapped on v3.

class DihedralJuelicher : public mirheo::VertexFetcherWithMeanCurvatures

Compute bending forces from the extended Juelicher model.

Public Types

using ParametersType = JuelicherBendingParameters

Type of parameters that describe the kernel.

Public Functions

DihedralJuelicher(ParametersType p, mReal lscale)

Initialize the functor.

Parameters
  • p: The parameters of the functor
  • lscale: Scaling length factor, applied to all parameters

void applyLengthScalingFactor(mReal lscale)

Scale length-dependent parameters.

void computeInternalCommonQuantities(const ViewType &view, int rbcId)

Precompute internal values that are common to all vertices in the cell.

Parameters
  • view: The view that contains the required object channels
  • rbcId: The index of the membrane in the view

mReal3 operator()(VertexType v0, VertexType v1, VertexType v2, VertexType v3, mReal3 &f1) const

Compute the dihedral forces.

See Developer docs for more details.

Return
The dihedral force acting on v0
Parameters
  • v0: vertex 0
  • v1: vertex 1
  • v2: vertex 2
  • v3: vertex 3
  • f1: force acting on v1; this method will add (not set) the dihedral force to that quantity.

class DihedralKantor : public mirheo::VertexFetcher

Compute bending forces from the Kantor model.

Public Types

using ParametersType = KantorBendingParameters

Type of parameters that describe the kernel.

Public Functions

DihedralKantor(ParametersType p, mReal lscale)

Initialize the functor.

Parameters
  • p: The parameters of the functor
  • lscale: Scaling length factor, applied to all parameters

void applyLengthScalingFactor(mReal lscale)

Scale length-dependent parameters.

void computeInternalCommonQuantities(const ViewType &view, int rbcId)

Precompute internal values that are common to all vertices in the cell.

mReal3 operator()(VertexType v0, VertexType v1, VertexType v2, VertexType v3, mReal3 &f1) const

Compute the dihedral forces.

See Developer docs for more details.

Return
The dihedral force acting on v0
Parameters
  • v0: vertex 0
  • v1: vertex 1
  • v2: vertex 2
  • v3: vertex 3
  • f1: force acting on v1; this method will add (not set) the dihedral force to that quantity.

Filters

The membrane interactions can be applied to only a subset of the given mirheo::MembraneVector. This can be convenient to have different interaction parameters for different membranes with the same mesh topology. Furthermore, reducing the numper of mirheo::ParticleVector is beneficial for performance (less interaction kernel launches so overhead for e.g. FSI).

class FilterKeepAll

Filter that keeps all the membranes.

Public Functions

void setPrerequisites(MembraneVector *mv) const

set required properties to mv

void setup(MembraneVector *mv)

Set internal state of the object.

bool inWhiteList(long membraneId) const

States if the given membrane must be kept or not.

Return
true if the membrane should be kept, false otherwise.
Parameters
  • membraneId: The index of the membrane to keep or not

class FilterKeepByTypeId

Keep membranes that have a given typeId.

The typeId of each membrane is stored in the object channel ChannelNames::membraneTypeId.

Public Functions

FilterKeepByTypeId(int whiteListTypeId)

Construct FilterKeepByTypeId that wil keep only the membranes with type id whiteListTypeId.

Parameters
  • whiteListTypeId: The type id of the membranes to keep

void setPrerequisites(MembraneVector *mv) const

set required properties to mv

Parameters

void setup(MembraneVector *mv)

Set internal state of the object.

This must be called after every change of

mv DataManager
Parameters

bool inWhiteList(long membraneId) const

States if the given membrane must be kept or not.

Return
true if the membrane should be kept, false otherwise.
Parameters
  • membraneId: The index of the membrane to keep or not

Fetchers

Fetchers are used to load generic data that is needed for kernel computation. In most cases, only vertex coordinates are sufficient (see mirheo::VertexFetcher). Additional data attached to each vertex may be needed, such as mean curvature in e.g. mirheo::DihedralJuelicher (see mirheo::VertexFetcherWithMeanCurvatures).

class VertexFetcher

Fetch a vertex for a given view.

Subclassed by mirheo::DihedralKantor, mirheo::VertexFetcherWithMeanCurvatures

Public Types

using VertexType = mReal3

info contained in the fetched data

using ViewType = OVview

compatible view

Public Functions

VertexType fetchVertex(const ViewType &view, int i) const

fetch a vertex coordinates from a view

Return
The vertex coordinates
Parameters
  • view: The view from which to fetch the vertex
  • i: The index of the vertex in view

class VertexFetcherWithMeanCurvatures : public mirheo::VertexFetcher

Fetch vertex coordinates and mean curvature for a given view.

Subclassed by mirheo::DihedralJuelicher

Public Types

using VertexType = VertexWithMeanCurvature

info contained in the fetched data

using ViewType = OVviewWithJuelicherQuants

compatible view

Public Functions

VertexType fetchVertex(const ViewType &view, int i) const

fetch a vertex coordinates and its mean curvature from a view

Return
The vertex coordinates
Parameters
  • view: The view from which to fetch the vertex
  • i: The index of the vertex in view

struct VertexWithMeanCurvature

holds vertex coordinates and mean curvature

Public Members

mReal3 r

vertex coordinates

mReal H

mean curvature