Allow for "computational electrophysiology" simulations (CompEl)
[alexxy/gromacs.git] / src / gromacs / legacyheaders / types / enums.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 #ifndef ENUMS_H_
39 #define ENUMS_H_
40
41 #ifdef __cplusplus
42 extern "C" {
43 #endif
44 #if 0
45 } /* fixes auto-indentation problems */
46 #endif
47
48 /* note: these enums should correspond to the names in gmxlib/names.c */
49
50 enum {
51     epbcXYZ, epbcNONE, epbcXY, epbcSCREW, epbcNR
52 };
53
54 enum {
55     etcNO, etcBERENDSEN, etcNOSEHOOVER, etcYES, etcANDERSEN, etcANDERSENMASSIVE, etcVRESCALE, etcNR
56 }; /* yes is an alias for berendsen */
57
58 #define ETC_ANDERSEN(e) (((e) == etcANDERSENMASSIVE) || ((e) == etcANDERSEN))
59
60 enum {
61     epcNO, epcBERENDSEN, epcPARRINELLORAHMAN, epcISOTROPIC, epcMTTK, epcNR
62 }; /* isotropic is an alias for berendsen */
63
64 /* trotter decomposition extended variable parts */
65 enum {
66     etrtNONE, etrtNHC, etrtBAROV, etrtBARONHC, etrtNHC2, etrtBAROV2, etrtBARONHC2,
67     etrtVELOCITY1, etrtVELOCITY2, etrtPOSITION, etrtSKIPALL, etrtNR
68 };
69
70 /* sequenced parts of the trotter decomposition */
71 enum {
72     ettTSEQ0,  ettTSEQ1,  ettTSEQ2,  ettTSEQ3,  ettTSEQ4, ettTSEQMAX
73 };
74
75 enum {
76     epctISOTROPIC, epctSEMIISOTROPIC, epctANISOTROPIC,
77     epctSURFACETENSION, epctNR
78 };
79
80 enum {
81     erscNO, erscALL, erscCOM, erscNR
82 };
83
84 enum {
85     ecutsVERLET, ecutsGROUP, ecutsNR
86 };
87
88 /* Coulomb / VdW interaction modifiers.
89  * grompp replaces eintmodPOTSHIFT_VERLET by eintmodPOTSHIFT or eintmodNONE.
90  * Exactcutoff is only used by Reaction-field-zero, and is not user-selectable.
91  */
92 enum eintmod {
93     eintmodPOTSHIFT_VERLET, eintmodPOTSHIFT, eintmodNONE, eintmodPOTSWITCH, eintmodEXACTCUTOFF, eintmodNR
94 };
95
96 /*
97  * eelNOTUSED1 used to be GB, but to enable generalized born with different
98  * forms of electrostatics (RF, switch, etc.) in the future it is now selected
99  * separately (through the implicit_solvent option).
100  */
101 enum {
102     eelCUT,     eelRF,     eelGRF,   eelPME,  eelEWALD,  eelP3M_AD,
103     eelPOISSON, eelSWITCH, eelSHIFT, eelUSER, eelGB_NOTUSED, eelRF_NEC, eelENCADSHIFT,
104     eelPMEUSER, eelPMESWITCH, eelPMEUSERSWITCH, eelRF_ZERO, eelNR
105 };
106
107 /* Ewald geometry */
108 enum {
109     eewg3D, eewg3DC, eewgNR
110 };
111
112 #define EEL_RF(e) ((e) == eelRF || (e) == eelGRF || (e) == eelRF_NEC || (e) == eelRF_ZERO )
113
114 #define EEL_PME(e)  ((e) == eelPME || (e) == eelPMESWITCH || (e) == eelPMEUSER || (e) == eelPMEUSERSWITCH || (e) == eelP3M_AD)
115 #define EEL_EWALD(e)  (EEL_PME(e) || (e) == eelEWALD)
116 #define EEL_FULL(e) (EEL_PME(e) || (e) == eelPOISSON || (e) == eelEWALD)
117
118 #define EEL_SWITCHED(e) ((e) == eelSWITCH || (e) == eelSHIFT || (e) == eelENCADSHIFT || (e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH)
119
120 #define EEL_USER(e) ((e) == eelUSER || (e) == eelPMEUSER || (e) == (eelPMEUSERSWITCH))
121
122 #define EEL_IS_ZERO_AT_CUTOFF(e) (EEL_SWITCHED(e) || (e) == eelRF_ZERO)
123
124 #define EEL_MIGHT_BE_ZERO_AT_CUTOFF(e) (EEL_IS_ZERO_AT_CUTOFF(e) || (e) == eelUSER || (e) == eelPMEUSER)
125
126 enum {
127     evdwCUT, evdwSWITCH, evdwSHIFT, evdwUSER, evdwENCADSHIFT,
128     evdwPME, evdwNR
129 };
130
131 enum {
132     eljpmeGEOM, eljpmeLB, eljpmeNR
133 };
134
135 #define EVDW_PME(e) ((e) == evdwPME)
136
137 #define EVDW_SWITCHED(e) ((e) == evdwSWITCH || (e) == evdwSHIFT || (e) == evdwENCADSHIFT)
138
139 #define EVDW_IS_ZERO_AT_CUTOFF(e) EVDW_SWITCHED(e)
140
141 #define EVDW_MIGHT_BE_ZERO_AT_CUTOFF(e) (EVDW_IS_ZERO_AT_CUTOFF(e) || (e) == evdwUSER)
142
143 enum {
144     ensGRID, ensSIMPLE, ensNR
145 };
146
147 /* eiVV is normal velocity verlet -- eiVVAK uses 1/2*(KE(t-dt/2)+KE(t+dt/2)) as the kinetic energy, and the half step kinetic
148    energy for temperature control */
149
150 enum {
151     eiMD, eiSteep, eiCG, eiBD, eiSD2, eiNM, eiLBFGS, eiTPI, eiTPIC, eiSD1, eiVV, eiVVAK, eiNR
152 };
153 #define EI_VV(e) ((e) == eiVV || (e) == eiVVAK)
154 #define EI_MD(e) ((e) == eiMD || EI_VV(e))
155 #define EI_SD(e) ((e) == eiSD1 || (e) == eiSD2)
156 #define EI_RANDOM(e) (EI_SD(e) || (e) == eiBD)
157 /*above integrators may not conserve momenta*/
158 #define EI_DYNAMICS(e) (EI_MD(e) || EI_SD(e) || (e) == eiBD)
159 #define EI_ENERGY_MINIMIZATION(e) ((e) == eiSteep || (e) == eiCG || (e) == eiLBFGS)
160 #define EI_TPI(e) ((e) == eiTPI || (e) == eiTPIC)
161
162 #define EI_STATE_VELOCITY(e) (EI_MD(e) || EI_SD(e))
163
164 enum {
165     econtLINCS, econtSHAKE, econtNR
166 };
167
168 enum {
169     edrNone, edrSimple, edrEnsemble, edrNR
170 };
171
172 enum {
173     edrwConservative, edrwEqual, edrwNR
174 };
175
176 /* Combination rule things */
177 enum {
178     eCOMB_NONE, eCOMB_GEOMETRIC, eCOMB_ARITHMETIC, eCOMB_GEOM_SIG_EPS, eCOMB_NR
179 };
180
181 /* NBF selection */
182 enum {
183     eNBF_NONE, eNBF_LJ, eNBF_BHAM, eNBF_NR
184 };
185
186 /* simulated tempering methods */
187 enum {
188     esimtempGEOMETRIC, esimtempEXPONENTIAL, esimtempLINEAR, esimtempNR
189 };
190 /* FEP selection */
191 enum {
192     efepNO, efepYES, efepSTATIC, efepSLOWGROWTH, efepEXPANDED, efepNR
193 };
194 /* if efepNO, there are no evaluations at other states.
195    if efepYES, treated equivalently to efepSTATIC.
196    if efepSTATIC, then lambdas do not change during the simulation.
197    if efepSLOWGROWTH, then the states change monotonically throughout the simulation.
198    if efepEXPANDED, then expanded ensemble simulations are occuring.
199  */
200
201 /* FEP coupling types */
202 enum {
203     efptFEP, efptMASS, efptCOUL, efptVDW, efptBONDED, efptRESTRAINT, efptTEMPERATURE, efptNR
204 };
205
206 /* How the lambda weights are calculated:
207    elamstatsMETROPOLIS = using the metropolis criteria
208    elamstatsBARKER = using the Barker critera for transition weights - also called unoptimized Bennett
209    elamstatsMINVAR = using Barker + minimum variance for weights
210    elamstatsWL = Wang-Landu (using visitation counts)
211    elamstatsWWL = Weighted Wang-Landau (using optimized gibbs weighted visitation counts)
212  */
213 enum {
214     elamstatsNO, elamstatsMETROPOLIS, elamstatsBARKER, elamstatsMINVAR, elamstatsWL, elamstatsWWL, elamstatsNR
215 };
216
217 #define ELAMSTATS_EXPANDED(e) ((e) > elamstatsNO)
218
219 #define EWL(e) ((e) == elamstatsWL || (e) == elamstatsWWL)
220
221 /* How moves in lambda are calculated:
222    elmovemcMETROPOLIS - using the Metropolis criteria, and 50% up and down
223    elmovemcBARKER - using the Barker criteria, and 50% up and down
224    elmovemcGIBBS - computing the transition using the marginalized probabilities of the lambdas
225    elmovemcMETGIBBS - computing the transition using the metropolized version of Gibbs (Monte Carlo Strategies in Scientific computing, Liu, p. 134)
226  */
227 enum {
228     elmcmoveNO, elmcmoveMETROPOLIS, elmcmoveBARKER, elmcmoveGIBBS, elmcmoveMETGIBBS, elmcmoveNR
229 };
230
231 /* how we decide whether weights have reached equilibrium
232    elmceqNO - never stop, weights keep going
233    elmceqYES - fix the weights from the beginning; no movement
234    elmceqWLDELTA - stop when the WL-delta falls below a certain level
235    elmceqNUMATLAM - stop when we have a certain number of samples at every step
236    elmceqSTEPS - stop when we've run a certain total number of steps
237    elmceqSAMPLES - stop when we've run a certain total number of samples
238    elmceqRATIO - stop when the ratio of samples (lowest to highest) is sufficiently large
239  */
240 enum {
241     elmceqNO, elmceqYES, elmceqWLDELTA, elmceqNUMATLAM, elmceqSTEPS, elmceqSAMPLES, elmceqRATIO, elmceqNR
242 };
243
244 /* separate_dhdl_file selection */
245 enum
246 {
247     /* NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool */
248     esepdhdlfileYES, esepdhdlfileNO, esepdhdlfileNR
249 };
250
251 /* dhdl_derivatives selection */
252 enum
253 {
254     /* NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool */
255     edhdlderivativesYES, edhdlderivativesNO, edhdlderivativesNR
256 };
257
258 /* Solvent model */
259 enum {
260     esolNO, esolSPC, esolTIP4P, esolNR
261 };
262
263 /* Dispersion correction */
264 enum {
265     edispcNO, edispcEnerPres, edispcEner, edispcAllEnerPres, edispcAllEner, edispcNR
266 };
267
268 /* Center of mass motion selection */
269 enum {
270     ecmLINEAR, ecmANGULAR, ecmNO, ecmNR
271 };
272
273 /* New version of simulated annealing */
274 enum {
275     eannNO, eannSINGLE, eannPERIODIC, eannNR
276 };
277
278 /* Implicit solvent algorithms */
279 enum {
280     eisNO, eisGBSA, eisNR
281 };
282
283 /* Algorithms for calculating GB radii */
284 enum {
285     egbSTILL, egbHCT, egbOBC, egbNR
286 };
287
288 enum {
289     esaAPPROX, esaNO, esaSTILL, esaNR
290 };
291
292 /* Wall types */
293 enum {
294     ewt93, ewt104, ewtTABLE, ewt126, ewtNR
295 };
296
297 /* Pull stuff */
298 enum {
299     epullNO, epullUMBRELLA, epullCONSTRAINT, epullCONST_F, epullNR
300 };
301
302 enum {
303     epullgDIST, epullgDIR, epullgCYL, epullgDIRPBC, epullgNR
304 };
305
306 #define PULL_CYL(pull) ((pull)->eGeom == epullgCYL)
307
308 /* Enforced rotation groups */
309 enum {
310     erotgISO, erotgISOPF,
311     erotgPM, erotgPMPF,
312     erotgRM, erotgRMPF,
313     erotgRM2, erotgRM2PF,
314     erotgFLEX, erotgFLEXT,
315     erotgFLEX2, erotgFLEX2T,
316     erotgNR
317 };
318
319 enum {
320     erotgFitRMSD, erotgFitNORM, erotgFitPOT, erotgFitNR
321 };
322
323 /* Direction along which ion/water swaps happen in "Computational
324  * Electrophysiology" (CompEL) setups */
325 enum eSwaptype {
326     eswapNO, eswapX, eswapY, eswapZ, eSwapTypesNR
327 };
328 /* The following three enums are mainly used for indexing arrays and when
329  * looping over the available ions, channels, or compartments. This hopefully
330  * adds to the code's readability because it makes clear which object is dealt
331  * with in a block of code.
332  *
333  * The two compartments for CompEL setups */
334 enum eCompartment {
335     eCompA, eCompB, eCompNR
336 };
337 /* The positive and negative ions CompEL setups, future versions of the
338  * protocol might consider more than two types of ions */
339 enum eIontype {
340     eIonNEG, eIonPOS, eIonNR
341 };
342 /* The chanels that define with their COM the compartment boundaries in
343  * CompEL setups. In principle one could also use modified setups with
344  * more than two channels. */
345 enum eChannel {
346     eChan0, eChan1, eChanNR
347 };
348
349 /* QMMM */
350 enum {
351     eQMmethodAM1, eQMmethodPM3, eQMmethodRHF,
352     eQMmethodUHF, eQMmethodDFT, eQMmethodB3LYP, eQMmethodMP2, eQMmethodCASSCF, eQMmethodB3LYPLAN,
353     eQMmethodDIRECT, eQMmethodNR
354 };
355
356 enum {
357     eQMbasisSTO3G, eQMbasisSTO3G2, eQMbasis321G,
358     eQMbasis321Gp, eQMbasis321dGp, eQMbasis621G,
359     eQMbasis631G, eQMbasis631Gp, eQMbasis631dGp,
360     eQMbasis6311G, eQMbasisNR
361 };
362
363 enum {
364     eQMMMschemenormal, eQMMMschemeoniom, eQMMMschemeNR
365 };
366
367 enum {
368     eMultentOptName, eMultentOptNo, eMultentOptLast, eMultentOptNR
369 };
370
371 /* flat-bottom posres geometries */
372 enum {
373     efbposresZERO, efbposresSPHERE, efbposresCYLINDER, efbposresX, efbposresY, efbposresZ,
374     efbposresNR
375 };
376
377 enum {
378     eAdressOff, eAdressConst, eAdressXSplit, eAdressSphere, eAdressNR
379 };
380
381 enum {
382     eAdressICOff, eAdressICThermoForce, eAdressICNR
383 };
384
385 enum {
386     eAdressSITEcom, eAdressSITEcog, eAdressSITEatom, eAdressSITEatomatom, eAdressSITENR
387 };
388
389
390 /* The interactions contained in a (possibly merged) table
391  * for computing electrostatic, VDW repulsion and/or VDW dispersion
392  * contributions.
393  */
394 enum gmx_table_interaction
395 {
396     GMX_TABLE_INTERACTION_ELEC,
397     GMX_TABLE_INTERACTION_VDWREP_VDWDISP,
398     GMX_TABLE_INTERACTION_VDWEXPREP_VDWDISP,
399     GMX_TABLE_INTERACTION_VDWDISP,
400     GMX_TABLE_INTERACTION_ELEC_VDWREP_VDWDISP,
401     GMX_TABLE_INTERACTION_ELEC_VDWEXPREP_VDWDISP,
402     GMX_TABLE_INTERACTION_ELEC_VDWDISP,
403     GMX_TABLE_INTERACTION_NR
404 };
405
406 /* Different formats for table data. Cubic spline tables are typically stored
407  * with the four Y,F,G,H intermediate values (check tables.c for format), which
408  * makes it easy to load with a single 4-way SIMD instruction too.
409  * Linear tables only need one value per table point, or two if both V and F
410  * are calculated. However, with SIMD instructions this makes the loads unaligned,
411  * and in that case we store the data as F, D=F(i+1)-F(i), V, and then a blank value,
412  * which again makes it possible to load as a single instruction.
413  */
414 enum gmx_table_format
415 {
416     GMX_TABLE_FORMAT_CUBICSPLINE_YFGH,
417     GMX_TABLE_FORMAT_LINEAR_VF,
418     GMX_TABLE_FORMAT_LINEAR_V,
419     GMX_TABLE_FORMAT_LINEAR_F,
420     GMX_TABLE_FORMAT_LINEAR_FDV0,
421     GMX_TABLE_FORMAT_NR
422 };
423
424 /* Neighborlist geometry type.
425  * Kernels will compute interactions between two particles,
426  * 3-center water, 4-center water or coarse-grained beads.
427  */
428 enum gmx_nblist_kernel_geometry
429 {
430     GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE,
431     GMX_NBLIST_GEOMETRY_WATER3_PARTICLE,
432     GMX_NBLIST_GEOMETRY_WATER3_WATER3,
433     GMX_NBLIST_GEOMETRY_WATER4_PARTICLE,
434     GMX_NBLIST_GEOMETRY_WATER4_WATER4,
435     GMX_NBLIST_GEOMETRY_CG_CG,
436     GMX_NBLIST_GEOMETRY_NR
437 };
438
439 /* Types of electrostatics calculations available inside nonbonded kernels.
440  * Note that these do NOT necessarily correspond to the user selections in the MDP file;
441  * many interactions for instance map to tabulated kernels.
442  */
443 enum gmx_nbkernel_elec
444 {
445     GMX_NBKERNEL_ELEC_NONE,
446     GMX_NBKERNEL_ELEC_COULOMB,
447     GMX_NBKERNEL_ELEC_REACTIONFIELD,
448     GMX_NBKERNEL_ELEC_CUBICSPLINETABLE,
449     GMX_NBKERNEL_ELEC_GENERALIZEDBORN,
450     GMX_NBKERNEL_ELEC_EWALD,
451     GMX_NBKERNEL_ELEC_NR
452 };
453
454 /* Types of vdw calculations available inside nonbonded kernels.
455  * Note that these do NOT necessarily correspond to the user selections in the MDP file;
456  * many interactions for instance map to tabulated kernels.
457  */
458 enum gmx_nbkernel_vdw
459 {
460     GMX_NBKERNEL_VDW_NONE,
461     GMX_NBKERNEL_VDW_LENNARDJONES,
462     GMX_NBKERNEL_VDW_BUCKINGHAM,
463     GMX_NBKERNEL_VDW_CUBICSPLINETABLE,
464     GMX_NBKERNEL_VDW_NR
465 };
466 /* Types of interactions inside the neighborlist
467  */
468 enum gmx_nblist_interaction_type
469 {
470     GMX_NBLIST_INTERACTION_STANDARD,
471     GMX_NBLIST_INTERACTION_FREE_ENERGY,
472     GMX_NBLIST_INTERACTION_ADRESS,
473     GMX_NBLIST_INTERACTION_NR
474 };
475
476 #ifdef __cplusplus
477 }
478 #endif
479
480 #endif /* ENUMS_H_ */