Removed refs to charmm files in share/top/Makefile.am
[gromacs/qmmm-gamess-us.git] / include / types / state.h
blobf3a80c2b536aaff654e51553853e7729334772c5
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GRoups of Organic Molecules in ACtion for Science
35 #ifndef _state_h_
36 #define _state_h_
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
43 * The t_state struct should contain all the (possibly) non-static
44 * information required to define the state of the system.
45 * Currently the random seeds for SD and BD are missing.
48 /* These enums are used in flags as (1<<est...).
49 * The order of these enums should not be changed,
50 * since that affects the checkpoint (.cpt) file format.
52 enum { estLAMBDA,
53 estBOX, estBOX_REL, estBOXV, estPRES_PREV, estNH_XI, estTC_INT,
54 estX, estV, estSDX, estCGP, estLD_RNG, estLD_RNGI,
55 estDISRE_INITF, estDISRE_RM3TAV,
56 estORIRE_INITF, estORIRE_DTAV,
57 estNR };
59 /* The names of the state entries, defined in src/gmxib/checkpoint.c */
60 extern const char *est_names[estNR];
62 typedef struct
64 real disre_initf; /* The scaling factor for initializing the time av. */
65 int ndisrepairs; /* The number of distance restraints */
66 real *disre_rm3tav; /* The r^-3 time averaged pair distances */
67 real orire_initf; /* The scaling factor for initializing the time av. */
68 int norire_Dtav; /* The number of matrix element in dtav (npair*5) */
69 real *orire_Dtav; /* The time averaged orientation tensors */
70 } history_t;
72 /* Struct used for checkpointing only.
73 * This struct would not be required with unlimited precision.
74 * But because of limited precision, the COM motion removal implementation
75 * can cause the kinetic energy in the MD loop to differ by a few bits from
76 * the kinetic energy one would determine from state.v.
78 typedef struct
80 bool bUpToDate;
81 int ekinh_n;
82 tensor * ekinh;
83 real dekindl;
84 real mvcos;
85 } ekinstate_t;
87 typedef struct
89 gmx_large_int_t nsteps; /* The number of steps in the history */
90 gmx_large_int_t nsum; /* The nr. of steps in the ener_ave and ener_sum */
91 double * ener_ave; /* Energy term history sum to get fluctuations */
92 double * ener_sum; /* Energy term history sum to get fluctuations */
93 int nener; /* Number of energy terms in two previous arrays */
94 gmx_large_int_t nsteps_sim; /* The number of steps in ener_sum_sim */
95 gmx_large_int_t nsum_sim; /* The number of frames in ener_sum_sim */
96 double * ener_sum_sim; /* Energy term history sum of the whole sim */
98 energyhistory_t;
100 typedef struct
102 int natoms;
103 int ngtc;
104 int nrng;
105 int nrngi;
106 int flags; /* Flags telling which entries are present */
107 real lambda; /* the free energy switching parameter */
108 matrix box; /* box vector coordinates */
109 matrix box_rel; /* Relitaive box vectors to preserve shape */
110 matrix boxv; /* box velocitites for Parrinello-Rahman pcoupl */
111 matrix pres_prev; /* Pressure of the previous step for pcoupl */
112 real *nosehoover_xi; /* for Nose-Hoover tcoupl (ngtc) */
113 double *therm_integral; /* for N-H/V-rescale tcoupl (ngtc) */
114 int nalloc; /* Allocation size for x, v and sd_x when !=NULL*/
115 rvec *x; /* the coordinates (natoms) */
116 rvec *v; /* the velocities (natoms) */
117 rvec *sd_X; /* random part of the x update for stoch. dyn. */
118 rvec *cg_p; /* p vector for conjugate gradient minimization */
120 unsigned int *ld_rng; /* RNG random state */
121 int *ld_rngi; /* RNG index */
123 history_t hist; /* Time history for restraints */
125 ekinstate_t ekinstate; /* The state of the kinetic energy data */
127 energyhistory_t enerhist; /* Energy history for statistics */
129 int ddp_count; /* The DD partitioning count for this state */
130 int ddp_count_cg_gl; /* The DD part. count for index_gl */
131 int ncg_gl; /* The number of local charge groups */
132 int *cg_gl; /* The global cg number of the local cgs */
133 int cg_gl_nalloc; /* Allocation size of cg_gl; */
134 } t_state;
136 #endif /* _state_h_ */