Gromacs  2024.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Typedefs | Enumerations | Functions
state.h File Reference
#include <array>
#include <memory>
#include <vector>
#include "gromacs/gpu_utils/hostallocator.h"
#include "gromacs/math/paddedvector.h"
#include "gromacs/math/vectypes.h"
#include "gromacs/mdtypes/md_enums.h"
#include "gromacs/utility/enumerationhelpers.h"
#include "gromacs/utility/real.h"
+ Include dependency graph for state.h:

Description

This file contains the definition of the microstate of the simulated system.

History of observables that needs to be checkpointed should be stored in ObservablesHistory. The state of the mdrun machinery that needs to be checkpointed is also stored elsewhere.

Author
Berk Hess

Classes

class  gmx::CheckpointData< operation >
 } More...
 
class  history_t
 History information for NMR distance and orientation restraints. More...
 
class  ekinstate_t
 Struct used for checkpointing only. More...
 
struct  df_history_t
 Free-energy sampling history struct. More...
 
class  t_state
 The microstate of the system. More...
 

Typedefs

template<class T >
using PaddedHostVector = gmx::PaddedHostVector< T >
 Convenience alias for until all is moved in the gmx namespace.
 

Enumerations

enum  StateEntry : int {
  Lambda, Box, BoxRel, BoxV,
  PressurePrevious, Nhxi, ThermInt, X,
  V, SDxNotSupported, Cgp, LDRngNotSupported,
  LDRngINotSupported, DisreInitF, DisreRm3Tav, OrireInitF,
  OrireDtav, SVirPrev, Nhvxi, Veta,
  Vol0, Nhpresxi, Nhpresvxi, FVirPrev,
  FepState, MCRngNotSupported, MCRngINotSupported, BarosInt,
  PullComPrevStep, Count
}
 Enum for all entries in t_state. More...
 

Functions

const char * enumValueToString (StateEntry enumValue)
 The names of the state entries, defined in src/gromacs/fileio/checkpoint.cpp.
 
template<typename Enum >
int enumValueToBitMask (Enum enumValue)
 Convert enum to bitmask value. More...
 
void init_gtc_state (t_state *state, int ngtc, int nnhpres, int nhchainlength)
 Resizes the T- and P-coupling state variables.
 
void init_dfhist_state (t_state *state, int dfhistNumLambda)
 Allocates memory for free-energy history.
 
void comp_state (const t_state *st1, const t_state *st2, bool bRMSD, real ftol, real abstol)
 Compares two states, write the differences to stdout.
 
rvecmakeRvecArray (gmx::ArrayRef< const gmx::RVec > v, gmx::Index n)
 Allocates an rvec pointer and copy the contents of v to it.
 
void set_box_rel (const t_inputrec *ir, t_state *state)
 Determine the relative box components. More...
 
void preserveBoxShape (const PressureCouplingOptions &pressureCoupling, const tensor deform, matrix box_rel, matrix box)
 Make sure the relative box shape remains the same. More...
 
static gmx::ArrayRef< const
gmx::RVec
positionsFromStatePointer (const t_state *state)
 Returns an arrayRef to the positions in state when state!=null. More...
 
void printLambdaStateToLog (FILE *fplog, gmx::ArrayRef< const real > lambda, bool isInitialOutput)
 Prints the current lambda state to the log file. More...
 
void initialize_lambdas (FILE *fplog, FreeEnergyPerturbationType freeEnergyPerturbationType, bool haveSimulatedTempering, const t_lambda &fep, gmx::ArrayRef< const real > simulatedTemperingTemps, gmx_ekindata_t *ekind, bool isMain, int *fep_state, gmx::ArrayRef< real > lambda)
 Fills fep_state and lambda if needed. More...
 

Enumeration Type Documentation

enum StateEntry : int
strong

Enum for all entries in t_state.

These enums are used in flags as (1<<est...). The order of these enums should not be changed, since that affects the checkpoint (.cpt) file format.

Function Documentation

template<typename Enum >
int enumValueToBitMask ( Enum  enumValue)
inline

Convert enum to bitmask value.

Used for setting flags in checkpoint header and verifying which flags are set.

void initialize_lambdas ( FILE *  fplog,
FreeEnergyPerturbationType  freeEnergyPerturbationType,
bool  haveSimulatedTempering,
const t_lambda &  fep,
gmx::ArrayRef< const real simulatedTemperingTemps,
gmx_ekindata_t *  ekind,
bool  isMain,
int *  fep_state,
gmx::ArrayRef< real lambda 
)

Fills fep_state and lambda if needed.

If FEP or simulated tempering is in use, fills fep_state and lambda on the main rank and sets the reference temperatures in ekind on all ranks.

Reports the initial lambda state to the log file.

static gmx::ArrayRef<const gmx::RVec> positionsFromStatePointer ( const t_state state)
inlinestatic

Returns an arrayRef to the positions in state when state!=null.

When state=nullptr, returns an empty arrayRef.

Note
The size returned is the number of atoms, without padding.
Parameters
[in]stateThe state, can be nullptr
void preserveBoxShape ( const PressureCouplingOptions &  pressureCoupling,
const tensor  deform,
matrix  box_rel,
matrix  box 
)

Make sure the relative box shape remains the same.

This function ensures that the relative box dimensions are preserved, which otherwise might diffuse away due to rounding errors in pressure coupling or the deform option.

Parameters
[in]pressureCouplingThe pressure-coupling options
[in]deformThe box-deformation tensor
[in]box_relRelative box dimensions
[in,out]boxThe corrected actual box dimensions
void printLambdaStateToLog ( FILE *  fplog,
gmx::ArrayRef< const real lambda,
bool  isInitialOutput 
)

Prints the current lambda state to the log file.

Parameters
[in]fplogThe log file. If fplog == nullptr there will be no output.
[in]lambdaThe array of lambda values.
[in]isInitialOutputWhether this output is the initial lambda state or not.
void set_box_rel ( const t_inputrec *  ir,
t_state state 
)

Determine the relative box components.

Set box_rel e.g. used in mdrun state, used to preserve the box shape

Parameters
[in]irInput record
[in,out]stateState