Gromacs
2022.2
|
#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"
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.
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 | state_change_natoms (t_state *state, int natoms) |
Change the number of atoms represented by this state, allocating memory as needed. | |
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. | |
rvec * | makeRvecArray (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 | preserve_box_shape (const t_inputrec *ir, 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::ArrayRef< real > ref_t, bool isMaster, int *fep_state, gmx::ArrayRef< real > lambda) |
Fills fep_state and lambda if needed. More... | |
|
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.
|
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::ArrayRef< real > | ref_t, | ||
bool | isMaster, | ||
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 master rank.
Reports the initial lambda state to the log file.
|
inlinestatic |
Returns an arrayRef to the positions in state
when state!=null
.
When state=nullptr
, returns an empty arrayRef.
[in] | state | The state, can be nullptr |
void preserve_box_shape | ( | const t_inputrec * | ir, |
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.
[in] | ir | Input record |
[in] | box_rel | Relative box dimensions |
[in,out] | box | The corrected actual box dimensions |
void printLambdaStateToLog | ( | FILE * | fplog, |
gmx::ArrayRef< const real > | lambda, | ||
bool | isInitialOutput | ||
) |
Prints the current lambda state to the log file.
[in] | fplog | The log file. If fplog == nullptr there will be no output. |
[in] | lambda | The array of lambda values. |
[in] | isInitialOutput | Whether 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
[in] | ir | Input record |
[in,out] | state | State |