MaMiCo 1.2
Loading...
Searching...
No Matches
coupling::interface::MDSolverInterface< LinkedCell, dim > Class Template Referenceabstract

interface to the MD simulation More...

#include <MDSolverInterface.h>

Public Types

using CellIndex_T = I11
 

Public Member Functions

virtual ~MDSolverInterface ()
 
virtual LinkedCell & getLinkedCell (const CellIndex_T &couplingCellIndex, const tarch::la::Vector< dim, unsigned int > &linkedCellInCouplingCell, const tarch::la::Vector< dim, unsigned int > &linkedCellsPerCouplingCell)=0
 
virtual tarch::la::Vector< dim, double > getGlobalMDDomainSize () const =0
 
virtual tarch::la::Vector< dim, double > getGlobalMDDomainOffset () const =0
 
virtual double getMoleculeMass () const =0
 
virtual double getKB () const =0
 
virtual double getMoleculeSigma () const =0
 
virtual double getMoleculeEpsilon () const =0
 
virtual void getInitialVelocity (const tarch::la::Vector< dim, double > &meanVelocity, const double &kB, const double &temperature, tarch::la::Vector< dim, double > &initialVelocity) const =0
 
virtual void deleteMoleculeFromMDSimulation (const coupling::interface::Molecule< dim > &molecule, LinkedCell &cell)=0
 
virtual void addMoleculeToMDSimulation (const coupling::interface::Molecule< dim > &molecule)=0
 
virtual void setupPotentialEnergyLandscape (const tarch::la::Vector< dim, unsigned int > &indexOfFirstCouplingCell, const tarch::la::Vector< dim, unsigned int > &rangeCouplingCells, const tarch::la::Vector< dim, unsigned int > &linkedCellsPerCouplingCell)=0
 
virtual tarch::la::Vector< dim, unsigned int > getLinkedCellIndexForMoleculePosition (const tarch::la::Vector< dim, double > &position)=0
 
virtual void calculateForceAndEnergy (coupling::interface::Molecule< dim > &molecule)=0
 
virtual void synchronizeMoleculesAfterMassModification ()=0
 
virtual void synchronizeMoleculesAfterMomentumModification ()=0
 
virtual double getDt ()=0
 
virtual coupling::interface::MoleculeIterator< LinkedCell, dim > * getMoleculeIterator (LinkedCell &cell)=0
 

Detailed Description

template<class LinkedCell, unsigned int dim>
class coupling::interface::MDSolverInterface< LinkedCell, dim >

interface to the MD simulation

This class provides

Template Parameters
dimNumber of dimensions; it can be 1, 2 or 3
Author
Philipp Neumann

Constructor & Destructor Documentation

◆ ~MDSolverInterface()

template<class LinkedCell, unsigned int dim>
virtual coupling::interface::MDSolverInterface< LinkedCell, dim >::~MDSolverInterface ( )
inlinevirtual

Destructor

Member Function Documentation

◆ addMoleculeToMDSimulation()

template<class LinkedCell, unsigned int dim>
virtual void coupling::interface::MDSolverInterface< LinkedCell, dim >::addMoleculeToMDSimulation ( const coupling::interface::Molecule< dim > & molecule)
pure virtual

This function adds the molecule to the MD simulation.

Parameters
molecule

◆ calculateForceAndEnergy()

template<class LinkedCell, unsigned int dim>
virtual void coupling::interface::MDSolverInterface< LinkedCell, dim >::calculateForceAndEnergy ( coupling::interface::Molecule< dim > & molecule)
pure virtual

This function assumes that a molecule is placed somewhere inside the linked cell at index 'linkedCellIndex' and computes the force and potential energy contributions from all molecules in the same linked cell and the neighboured linked cells onto this molecule. The molecule "molecule" is considered to NOT be part of the simulation domain and thus the linked cells. Therefore, another molecule inside the linked cells may even coincide with "molecule". The results are stored within the molecule.

Parameters
molecule

◆ deleteMoleculeFromMDSimulation()

template<class LinkedCell, unsigned int dim>
virtual void coupling::interface::MDSolverInterface< LinkedCell, dim >::deleteMoleculeFromMDSimulation ( const coupling::interface::Molecule< dim > & molecule,
LinkedCell & cell )
pure virtual

This function deletes the molecule from the MD simulation

Parameters
molecule
cell

◆ getDt()

template<class LinkedCell, unsigned int dim>
virtual double coupling::interface::MDSolverInterface< LinkedCell, dim >::getDt ( )
pure virtual
Returns
the timestep of the MD simulation

◆ getGlobalMDDomainOffset()

template<class LinkedCell, unsigned int dim>
virtual tarch::la::Vector< dim, double > coupling::interface::MDSolverInterface< LinkedCell, dim >::getGlobalMDDomainOffset ( ) const
pure virtual

This function determines the offset (i.e. lower,left corner) of MD domain

Returns
the offset (i.e. lower,left corner) of MD domain

◆ getGlobalMDDomainSize()

template<class LinkedCell, unsigned int dim>
virtual tarch::la::Vector< dim, double > coupling::interface::MDSolverInterface< LinkedCell, dim >::getGlobalMDDomainSize ( ) const
pure virtual

This function specifies the global size of the box-shaped MD domain

Returns
the global size of the box-shaped MD domain

◆ getInitialVelocity()

template<class LinkedCell, unsigned int dim>
virtual void coupling::interface::MDSolverInterface< LinkedCell, dim >::getInitialVelocity ( const tarch::la::Vector< dim, double > & meanVelocity,
const double & kB,
const double & temperature,
tarch::la::Vector< dim, double > & initialVelocity ) const
pure virtual

This function sets a random velocity in the vector 'initialVelocity'. This velocity is sampled from a Maxwellian assuming a mean flow velocity 'meanVelocity' and a temperature 'temperature' of the fluid.

Parameters
meanVelocity
kB
temperature
initialVelocity

◆ getKB()

template<class LinkedCell, unsigned int dim>
virtual double coupling::interface::MDSolverInterface< LinkedCell, dim >::getKB ( ) const
pure virtual
Returns
Boltzmann's constant

◆ getLinkedCell()

template<class LinkedCell, unsigned int dim>
virtual LinkedCell & coupling::interface::MDSolverInterface< LinkedCell, dim >::getLinkedCell ( const CellIndex_T & couplingCellIndex,
const tarch::la::Vector< dim, unsigned int > & linkedCellInCouplingCell,
const tarch::la::Vector< dim, unsigned int > & linkedCellsPerCouplingCell )
pure virtual

This function specifies a particular linked cell inside a coupling cell. The coupling cells are currently located on the same process as the respective linked cells. However, several linked cells may be part of a coupling cell. The coupling cells also contain a ghost layer which surrounds each local domain; the very first coupling cell inside the global MD domain (or local MD domain) is thus given by coordinates (1,1,1) (or (1,1) in 2D, respectively). The index linkedCellInCouplingCell corresponds to the coordinates of the linked cell inside the given coupling cell. These coordinates thus lie in a range (0,linkedCellsPerCouplingCell-1).

Parameters
couplingCellIndex
linkedCellInCouplingCell
linkedCellsPerCouplingCell
Returns
a particular linked cell inside a coupling cell.

◆ getLinkedCellIndexForMoleculePosition()

template<class LinkedCell, unsigned int dim>
virtual tarch::la::Vector< dim, unsigned int > coupling::interface::MDSolverInterface< LinkedCell, dim >::getLinkedCellIndexForMoleculePosition ( const tarch::la::Vector< dim, double > & position)
pure virtual

This function specifies the local index vector (w.r.t. lexicographic ordering of the linked cells in the MD simulation, for the linked cell that the position 'position' belongs to. This method is crucial for the USHER particle insertion scheme as we need to loop over all neighboured linked cells and the cell itself to determine potential energy and forces acting on the molecule at position 'position'.

Parameters
position
Returns
the local index vector for the linked cell that the position 'position' belongs to
See also
getLinkedCell()

◆ getMoleculeEpsilon()

template<class LinkedCell, unsigned int dim>
virtual double coupling::interface::MDSolverInterface< LinkedCell, dim >::getMoleculeEpsilon ( ) const
pure virtual
Returns
the epsilon parameter of the Lennard-Jones potential

◆ getMoleculeIterator()

template<class LinkedCell, unsigned int dim>
virtual coupling::interface::MoleculeIterator< LinkedCell, dim > * coupling::interface::MDSolverInterface< LinkedCell, dim >::getMoleculeIterator ( LinkedCell & cell)
pure virtual
Parameters
cell
Returns
a new molecule iterator for a certain linked cell

◆ getMoleculeMass()

template<class LinkedCell, unsigned int dim>
virtual double coupling::interface::MDSolverInterface< LinkedCell, dim >::getMoleculeMass ( ) const
pure virtual

This function specifies the mass of a single fluid molecule

Returns
the mass of a single fluid molecule

◆ getMoleculeSigma()

template<class LinkedCell, unsigned int dim>
virtual double coupling::interface::MDSolverInterface< LinkedCell, dim >::getMoleculeSigma ( ) const
pure virtual
Returns
the sigma parameter of the Lennard-Jones potential

◆ setupPotentialEnergyLandscape()

template<class LinkedCell, unsigned int dim>
virtual void coupling::interface::MDSolverInterface< LinkedCell, dim >::setupPotentialEnergyLandscape ( const tarch::la::Vector< dim, unsigned int > & indexOfFirstCouplingCell,
const tarch::la::Vector< dim, unsigned int > & rangeCouplingCells,
const tarch::la::Vector< dim, unsigned int > & linkedCellsPerCouplingCell )
pure virtual

This function sets up the potential energy landscape over the domain spanned by indexOfFirstCouplingCell and rangeCoordinates. The first vector denotes the position of the lower,left,front corner of the domain, rangeCoordinates the number of coupling cells in each spatial direction in which the potential energy needs to be computed.

Parameters
indexOfFirstCouplingCell
rangeCouplingCells
linkedCellsPerCouplingCell

◆ synchronizeMoleculesAfterMassModification()

template<class LinkedCell, unsigned int dim>
virtual void coupling::interface::MDSolverInterface< LinkedCell, dim >::synchronizeMoleculesAfterMassModification ( )
pure virtual

This function is called each time when MaMiCo tried to insert/ delete molecules from the MD simulation. As a consequence, a synchronization between processes or with local boundary data might be necessary. Example: If in the lower left cell of a 2D MD simulation a molecule was inserted and we use periodic boundary conditions, we need to provide the inserted molecule (amongst others) to the upper right cell by adding this molecule to the ghost boundary layer (when using the builtin MD simulation). For the builtin MD simulation, the implementation of this method clears all molecules from the ghost layers and re-fills the ghost layers again.

◆ synchronizeMoleculesAfterMomentumModification()

template<class LinkedCell, unsigned int dim>
virtual void coupling::interface::MDSolverInterface< LinkedCell, dim >::synchronizeMoleculesAfterMomentumModification ( )
pure virtual

This function is called each time when MaMiCo tried to insert momentum in the MD simulation. For the builtin MD simulation, this method is empty as the simulation does not need to synchronize changing momentum over the processes within a timestep (only the positions and numbers of molecules in possible ghost layers matter!). However, other solvers might need to implement something in here.


The documentation for this class was generated from the following file: