ls1-MarDyn
ls1-MarDyn molecular dynamics code
Public Member Functions | List of all members
Domain Class Reference

This class is used to read in the phasespace and to handle macroscopic values. More...

#include <Domain.h>

Public Member Functions

 Domain (int rank)
 The constructor sets _localRank to rank and initializes all member variables.
 
void readXML (XMLfileUnits &xmlconfig)
 
void writeCheckpoint (std::string filename, ParticleContainer *particleContainer, DomainDecompBase *domainDecomp, double currentTime, bool useBinaryFormat=false)
 writes a checkpoint file that can be used to continue the simulation More...
 
void writeCheckpointHeader (std::string filename, ParticleContainer *particleContainer, const DomainDecompBase *domainDecomp, double currentTime)
 writes a checkpoint file that can be used to continue the simulation More...
 
void writeCheckpointHeaderXML (std::string filename, ParticleContainer *particleContainer, const DomainDecompBase *domainDecomp, double currentTime)
 
void initParameterStreams (double cutoffRadius, double cutoffRadiusLJ)
 initialize far field correction parameters More...
 
void setLocalUpot (double Upot)
 set the potential of the local process
 
double getLocalUpot () const
 get the potential of the local process
 
void setLocalUpotCompSpecific (double UpotCspec)
 set the fluid and fluid-solid potential of the local process
 
void setNumFluidComponents (unsigned nc)
 set the number of fluid phase components (specified in the config-file)
 
unsigned getNumFluidComponents ()
 get the numbr of fluid molecules as specified in the config file (*_1R.cfg)
 
double getLocalUpotCompSpecific ()
 get the fluid and fluid-solid potential of the local process
 
void setLocalVirial (double Virial)
 set the virial of the local process
 
double getLocalVirial () const
 get the virial of the local process
 
double getGlobalBetaTrans ()
 get thermostat scaling for translations
 
double getGlobalBetaTrans (int thermostat)
 
double getGlobalBetaRot ()
 get thermostat scaling for rotations
 
double getGlobalBetaRot (int thermostat)
 
double getGlobalLength (int d) const
 return the length of the domain More...
 
void setGlobalLength (int index, double length)
 set the length of the domain More...
 
double getGlobalCurrentTemperature ()
 get the global temperature for the whole system (i.e. thermostat ID 0)
 
double getCurrentTemperature (int thermostatID)
 
double getTargetTemperature (int thermostatID)
 
void setGlobalTemperature (double T)
 set the global temperature
 
void setTargetTemperature (int thermostatID, double T)
 
std::vector< double > & getmixcoeff ()
 get the mixcoeff
 
double getepsilonRF () const
 get the epsilonRF
 
void setepsilonRF (double erf)
 set the epsilonRF
 
unsigned long getglobalNumMolecules () const
 get globalNumMolecules
 
void setglobalNumMolecules (unsigned long glnummol)
 set globalNumMolecules
 
void updateglobalNumMolecules (ParticleContainer *particleContainer, DomainDecompBase *domainDecomp)
 update globalNumMolecules
 
CommVar< uint64_t > getMaxMoleculeID () const
 get local/global max. moleculeID
 
void updateMaxMoleculeID (ParticleContainer *particleContainer, DomainDecompBase *domainDecomp)
 update max. moleculeID
 
double getGlobalPressure ()
 get the global pressure
 
double getAverageGlobalUpot () const
 get the global average potential per particle More...
 
double getGlobalUpot () const
 
double getAverageGlobalUpotCSpec ()
 by Stefan Becker: return the average global potential of the fluid-fluid and fluid-solid interaction (but NOT solid-solid interaction)
 
unsigned long getNumFluidMolecules ()
 
double getAverageGlobalVirial () const
 get the global average virial per particle More...
 
void setLocalSummv2 (double summv2, int thermostat)
 sets _localSummv2 to the given value
 
void setLocalSumIw2 (double sumIw2, int thermostat)
 sets _localSumIw2 to the given value
 
void setLocalNrotDOF (int thermostat, unsigned long N, unsigned long rotDOF)
 sets _localThermostatN and _localRotationalDOF for thermostat
 
double getglobalRho ()
 get globalRho
 
void setglobalRho (double grho)
 set globalRho
 
unsigned long getglobalRotDOF ()
 get globalRotDOF
 
void setglobalRotDOF (unsigned long grotdof)
 set globalRotDOF
 
Comp2ParamgetComp2Params ()
 get the parameter streams
 
void calculateGlobalValues (DomainDecompBase *domainDecomp, ParticleContainer *particleContainer, bool collectThermostatVelocities, double Tfactor)
 calculate the global macroscopic values More...
 
void calculateGlobalValues (DomainDecompBase *domainDecomp, ParticleContainer *particleContainer)
 calls this->calculateGlobalValues with Tfactor = 1 and without velocity collection
 
void calculateVelocitySums (ParticleContainer *partCont)
 calculate _localSummv2 and _localSumIw2 More...
 
void calculateThermostatDirectedVelocity (ParticleContainer *partCont)
 
bool thermostatIsUndirected (int thermostat)
 
double getThermostatDirectedVelocity (int thermostat, int d)
 returns the directed velocity associated with a thermostat More...
 
bool severalThermostats ()
 returns whether there are several distinct thermostats in the system
 
int getThermostat (int cid)
 thermostat to be applied to component cid More...
 
void disableComponentwiseThermostat ()
 disables the componentwise thermostat (so that a single thermostat is applied to all DOF)
 
void enableComponentwiseThermostat ()
 enables the componentwise thermostat
 
unsigned maxThermostat ()
 returns the ID of the "last" thermostat in the system More...
 
void setComponentThermostat (int cid, int thermostat)
 associates a component with a thermostat More...
 
void enableUndirectedThermostat (int thermostat)
 enables the "undirected" flag for the specified thermostat More...
 
unsigned long N ()
 
void Nadd (unsigned cid, int N, int localN)
 
double getGlobalVolume () const
 
void thermostatOff ()
 
void thermostatOn ()
 
bool NVE ()
 
bool thermostatWarning ()
 
void evaluateRho (unsigned long localN, DomainDecompBase *comm)
 
void submitDU (unsigned cid, double DU, double *r)
 
void setLambda (double lambda)
 
void setDensityCoefficient (float coeff)
 
void setProfiledComponentMass (double m)
 
void init_cv (unsigned N, double U, double UU)
 
void record_cv ()
 
double cv ()
 
double getSigma (unsigned cid, unsigned nthSigma)
 methods implemented by Stefan Becker stefa.nosp@m.n.be.nosp@m.cker@.nosp@m.mv.u.nosp@m.ni-kl.nosp@m..de
 
unsigned getNumberOfComponents ()
 
void setUpotCorr (double upotcorr)
 
void setVirialCorr (double virialcorr)
 
void SetExplosionHeuristics (bool bVal)
 

Detailed Description

This class is used to read in the phasespace and to handle macroscopic values.

Author
Martin Bernreuther bernr.nosp@m.euth.nosp@m.er@hl.nosp@m.rs.d.nosp@m.e et al. (2011)

This class is responsible for all macroscopic values. It is important to differentiate between local and global version of those values Macroscopic values are values that aggregate some property of a set of molecules. As this program is designed to run on a parallel computer, there are typically several processes. Each process has an instance of this class, but with a different subset of molecules. Therefore, also the macroscopic values are only representative for the "local" domain of that process.

member variables that represent "local" macroscopic values begin with _local

properties of the global system begin with _global if only the rank 0 process obtains the correct value (e.g., as a sum over corresponding local properties) and with _universal if the global value must be communicated to all processes.

At some points of the simulation, macroscopic values for the whole set of molecules have to be calculated. Those values are stored in member variables beginning with _global.

Member Function Documentation

◆ calculateGlobalValues()

void Domain::calculateGlobalValues ( DomainDecompBase domainDecomp,
ParticleContainer particleContainer,
bool  collectThermostatVelocities,
double  Tfactor 
)

calculate the global macroscopic values

Parameters
domainDecompdomain decomposition
particleContainerparticle Container
collectThermostatVelocitiesflag stating whether the directed velocity should be collected for the corresponding thermostats
Tfactortemporary factor applied to the temperature during equilibration

Essentially, this method calculates all thermophysical values that require communication between the subdomains of the system.

In particular, it determines:

  • the potential energy and the virial
  • if (collectThermostatVelocities == true), the directed velocity associated with the appropriately marked thermostats
  • the kinetic energy and number of DOF associated with each thermostat
  • on that basis, the correction factors for translation and rotation, also applying some heuristics in case of an explosion

The heuristics applied in case of an explosion were already commented on, in the appropriate place (Domain.cpp), but apparently it should be repeated here (according to the requirements communicated by TUM): Sometimes an "explosion" occurs, for instance, by inserting a particle at an unfavorable position, or due to imprecise integration of the equations of motion. The thermostat would then remove the kinetic energy from the rest of the system, effectively destroying a simulation run that could be saved by a more intelligent version, which is provided by the present heuristics. It is essential that such an intervention does not occur regularly, therefore it is limited to three times every 3000 time steps. That limit is implemented using the properties _universalSelectiveThermostatCounter, _universalSelectiveThermostatWarning, and _universalSelectiveThermostatError.

◆ calculateThermostatDirectedVelocity()

void Domain::calculateThermostatDirectedVelocity ( ParticleContainer partCont)

Calculates the LOCAL directed kinetic energy associated with the appropriately marked thermostats

◆ calculateVelocitySums()

void Domain::calculateVelocitySums ( ParticleContainer partCont)

calculate _localSummv2 and _localSumIw2

The present method calculates the translational and rotational kinetic energies associated with the thermostats of the system, together with the corresponding number of DOF. If a thermostat is marked as applying to the undirected motion only, this is explicitly considered by subtracting the average directed motion corresponding to the thermostat from the motion of each individual particle, for the purpose of aggregating the kinetic energy of the undirected motion.

◆ enableUndirectedThermostat()

void Domain::enableUndirectedThermostat ( int  thermostat)

enables the "undirected" flag for the specified thermostat

Parameters
tstinternal ID of the respective thermostat

◆ getAverageGlobalUpot()

double Domain::getAverageGlobalUpot ( ) const

get the global average potential per particle

Before this method is called, it has to be sure that the global potential has been calculated (method calculateGlobalValues)

◆ getAverageGlobalVirial()

double Domain::getAverageGlobalVirial ( ) const

get the global average virial per particle

Before this method is called, it has to be sure that the global virial has been calculated (method calculateGlobalValues)

◆ getGlobalLength()

double Domain::getGlobalLength ( int  d) const
inline

return the length of the domain

Parameters
ddimension for which the length should be returned

◆ getNumFluidMolecules()

unsigned long Domain::getNumFluidMolecules ( )

by Stefan Becker: determine and return the totel number of fluid molecules this method assumes all molecules with a component-ID less than _numFluidComponent to be fluid molecules

◆ getThermostat()

int Domain::getThermostat ( int  cid)
inline

thermostat to be applied to component cid

Parameters
cidID of the respective component

◆ getThermostatDirectedVelocity()

double Domain::getThermostatDirectedVelocity ( int  thermostat,
int  d 
)
inline

returns the directed velocity associated with a thermostat

Parameters
thID of the thermostat
dcoordinate 0 (v_x), 1 (v_y), or 2 (v_z).

It should be obvious that this method only returns sensible values for thermostats marked as "undirected", because otherwise the directed velocity is not explicitly computed.

◆ initParameterStreams()

void Domain::initParameterStreams ( double  cutoffRadius,
double  cutoffRadiusLJ 
)

initialize far field correction parameters

By limiting the calculation to pairs of particles which have less distance than the given cutoff radius, an error is made By calculating approximations for the neglected pairs, the error can be reduced. The way the cutoff correction works is probably explained by all books on molecular simulation, e.g. Allen/Tildesley [1].

The idea is that for radii above the cutoff radius, the radial distribution function is assumed to be exactly one, so that the correction for the internal energy as well as the virial can be obtained immediately from an integral over the pair potential and its derivative, respectively. For quadrupoles, the long-range correction turns out to be zero so that it is not calculated here. Molecules with point charges and a zero net charge interact over long distances according to their dipole moments, which is why effective dipoles are calculated for the far field correction.

If ions appear, so that the net charge of a molecule is distinct from zero, the far field terms become much more complicated, and THAT IS NOT IMPLEMENTED AT PRESENT.

[1] Computer Simulation of Liquids (1989), Clarendon, Oxford.

Parameters
cutoffRadiuscutoff radius for electrostatics
cutoffRadiusLJcutoff radius for the LJ potential

initialize parameter streams

This method should only be called, after the component information and all molecule data have been read in

Parameters
cutoffRadiuscutoff radius

◆ maxThermostat()

unsigned Domain::maxThermostat ( )
inline

returns the ID of the "last" thermostat in the system

The idea is that ID -1 refers to DOF without thermostats, while ID 0 refers to the entire system (or to the unique thermostat if the componentwise thermostat is disabled). In case of componentwise thermostats being applied, these have the IDs from 1 to this->maxThermostat().

◆ readXML()

void Domain::readXML ( XMLfileUnits xmlconfig)
Todo:
reading temperature is performed in the Ensemble, so check and remove

◆ setComponentThermostat()

void Domain::setComponentThermostat ( int  cid,
int  thermostat 
)
inline

associates a component with a thermostat

Parameters
cidinternal ID of the component
thinternal ID of the thermostat

◆ setGlobalLength()

void Domain::setGlobalLength ( int  index,
double  length 
)

set the length of the domain

Parameters
indexdimension for which the length should be set
lengthvalue which should be set

◆ thermostatIsUndirected()

bool Domain::thermostatIsUndirected ( int  thermostat)
inline

A thermostat is referred to as undirected here if it explicitly excludes the kinetic energy associated with the directed motion of the respective components.

◆ writeCheckpoint()

void Domain::writeCheckpoint ( std::string  filename,
ParticleContainer particleContainer,
DomainDecompBase domainDecomp,
double  currentTime,
bool  useBinaryFormat = false 
)

writes a checkpoint file that can be used to continue the simulation

The format of the checkpointfile written by this method is the same as the format of the input file.

Parameters
filenameName of the checkpointfile (including path)
particleContainerThe molecules that have to be written to the file are stored here
domainDecompIn the parallel version, the file has to be written by more than one process. Methods to achieve this are available in domainDecomp
currentTimeThe current time to be printed.
useBinaryFormatindicates wheter binary I/O is used or not

◆ writeCheckpointHeader()

void Domain::writeCheckpointHeader ( std::string  filename,
ParticleContainer particleContainer,
const DomainDecompBase domainDecomp,
double  currentTime 
)

writes a checkpoint file that can be used to continue the simulation

The format of the checkpointfile written by this method is the same as the format of the input file.

Parameters
filenameName of the header file (including path). Dependent on the I/O format, this could be the same file as the data file
particleContainerThe molecules that have to be written to the file are stored here
domainDecompIn the parallel version, the file has to be written by more than one process. Methods to achieve this are available in domainDecomp
currentTimeThe current time to be printed.

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