Remove support for implicit solvation
[alexxy/gromacs.git] / src / gromacs / mdtypes / md_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,2015,2016,2017,2018, 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 /*! \file
38  * \brief
39  * Declares enumerated types used throughout the code.
40  *
41  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42  * \inpublicapi
43  * \ingroup module_mdtypes
44  */
45 #ifndef GMX_MDTYPES_MD_ENUMS_H
46 #define GMX_MDTYPES_MD_ENUMS_H
47
48 #include "gromacs/utility/basedefinitions.h"
49
50 /*! \brief Return a string from a list of strings
51  *
52  * If index if within 0 .. max_index-1 returns the corresponding string
53  * or "no name defined" otherwise, in other words this is a range-check that does
54  * not crash.
55  * \param[in] index     The index in the array
56  * \param[in] max_index The length of the array
57  * \param[in] names     The array
58  * \return the correct string or "no name defined"
59  */
60 const char *enum_name(int index, int max_index, const char *names[]);
61
62 //! Boolean strings no or yes
63 extern const char *yesno_names[BOOL_NR+1];
64
65 //! \brief The two compartments for CompEL setups.
66 enum eCompartment {
67     eCompA, eCompB, eCompNR
68 };
69
70 /*! \brief The channels that define with their COM the compartment boundaries in CompEL setups.
71  *
72  * In principle one could also use modified setups with more than two channels.
73  */
74 enum eChannel {
75     eChan0, eChan1, eChanNR
76 };
77
78 /*! \brief Temperature coupling type
79  *
80  * yes is an alias for berendsen
81  */
82 enum {
83     etcNO, etcBERENDSEN, etcNOSEHOOVER, etcYES, etcANDERSEN, etcANDERSENMASSIVE, etcVRESCALE, etcNR
84 };
85 //! Strings corresponding to temperatyre coupling types
86 extern const char *etcoupl_names[etcNR+1];
87 //! Macro for selecting t coupling string
88 #define ETCOUPLTYPE(e) enum_name(e, etcNR, etcoupl_names)
89 //! Return whether this is andersen coupling
90 #define ETC_ANDERSEN(e) (((e) == etcANDERSENMASSIVE) || ((e) == etcANDERSEN))
91
92 /*! \brief Pressure coupling types
93  *
94  * isotropic is an alias for berendsen
95  */
96 enum {
97     epcNO, epcBERENDSEN, epcPARRINELLORAHMAN, epcISOTROPIC, epcMTTK, epcNR
98 };
99 //! String corresponding to pressure coupling algorithm
100 extern const char *epcoupl_names[epcNR+1];
101 //! Macro to return the correct pcoupling string
102 #define EPCOUPLTYPE(e) enum_name(e, epcNR, epcoupl_names)
103
104 //! Flat-bottom posres geometries
105 enum {
106     efbposresZERO, efbposresSPHERE, efbposresCYLINDER, efbposresX, efbposresY, efbposresZ,
107     efbposresCYLINDERX, efbposresCYLINDERY, efbposresCYLINDERZ, efbposresNR
108 };
109
110 //! Relative coordinate scaling type for position restraints.
111 enum {
112     erscNO, erscALL, erscCOM, erscNR
113 };
114 //! String corresponding to relativ coordinate scaling.
115 extern const char *erefscaling_names[erscNR+1];
116 //! Macro to select correct coordinate scaling string.
117 #define EREFSCALINGTYPE(e) enum_name(e, erscNR, erefscaling_names)
118
119 //! Trotter decomposition extended variable parts.
120 enum {
121     etrtNONE, etrtNHC, etrtBAROV, etrtBARONHC, etrtNHC2, etrtBAROV2, etrtBARONHC2,
122     etrtVELOCITY1, etrtVELOCITY2, etrtPOSITION, etrtSKIPALL, etrtNR
123 };
124
125 //! Sequenced parts of the trotter decomposition.
126 enum {
127     ettTSEQ0,  ettTSEQ1,  ettTSEQ2,  ettTSEQ3,  ettTSEQ4, ettTSEQMAX
128 };
129
130 //! Pressure coupling type
131 enum {
132     epctISOTROPIC, epctSEMIISOTROPIC, epctANISOTROPIC,
133     epctSURFACETENSION, epctNR
134 };
135 //! String corresponding to pressure coupling type
136 extern const char *epcoupltype_names[epctNR+1];
137 //! Macro to select the right string for pcoupl type
138 #define EPCOUPLTYPETYPE(e) enum_name(e, epctNR, epcoupltype_names)
139
140 //! \\brief Cutoff scheme
141 enum {
142     ecutsVERLET, ecutsGROUP, ecutsNR
143 };
144 //! String corresponding to cutoff scheme
145 extern const char *ecutscheme_names[ecutsNR+1];
146 //! Macro to select the right string for cutoff scheme
147 #define ECUTSCHEME(e)  enum_name(e, ecutsNR, ecutscheme_names)
148
149 /*! \brief Coulomb / VdW interaction modifiers.
150  *
151  * grompp replaces eintmodPOTSHIFT_VERLET by eintmodPOTSHIFT or eintmodNONE.
152  * Exactcutoff is only used by Reaction-field-zero, and is not user-selectable.
153  */
154 enum eintmod {
155     eintmodPOTSHIFT_VERLET, eintmodPOTSHIFT, eintmodNONE, eintmodPOTSWITCH, eintmodEXACTCUTOFF, eintmodFORCESWITCH, eintmodNR
156 };
157 //! String corresponding to interaction modifiers
158 extern const char *eintmod_names[eintmodNR+1];
159 //! Macro to select the correct string for modifiers
160 #define INTMODIFIER(e) enum_name(e, eintmodNR, eintmod_names)
161
162 /*! \brief Cut-off treatment for Coulomb */
163 enum {
164     eelCUT,     eelRF,     eelGRF,   eelPME,  eelEWALD,  eelP3M_AD,
165     eelPOISSON, eelSWITCH, eelSHIFT, eelUSER, eelGB_NOTUSED, eelRF_NEC_UNSUPPORTED, eelENCADSHIFT,
166     eelPMEUSER, eelPMESWITCH, eelPMEUSERSWITCH, eelRF_ZERO, eelNR
167 };
168 //! String corresponding to Coulomb treatment
169 extern const char *eel_names[eelNR+1];
170 //! Macro for correct string for Coulomb treatment
171 #define EELTYPE(e)     enum_name(e, eelNR, eel_names)
172
173 //! Ewald geometry.
174 enum {
175     eewg3D, eewg3DC, eewgNR
176 };
177 //! String corresponding to Ewald geometry
178 extern const char *eewg_names[eewgNR+1];
179
180 //! Macro telling us whether we use reaction field
181 #define EEL_RF(e) ((e) == eelRF || (e) == eelGRF || (e) == eelRF_NEC_UNSUPPORTED || (e) == eelRF_ZERO )
182
183 //! Macro telling us whether we use PME
184 #define EEL_PME(e)  ((e) == eelPME || (e) == eelPMESWITCH || (e) == eelPMEUSER || (e) == eelPMEUSERSWITCH || (e) == eelP3M_AD)
185 //! Macro telling us whether we use PME or full Ewald
186 #define EEL_PME_EWALD(e) (EEL_PME(e) || (e) == eelEWALD)
187 //! Macro telling us whether we use full electrostatics of any sort
188 #define EEL_FULL(e) (EEL_PME_EWALD(e) || (e) == eelPOISSON)
189 //! Macro telling us whether we use user defined electrostatics
190 #define EEL_USER(e) ((e) == eelUSER || (e) == eelPMEUSER || (e) == (eelPMEUSERSWITCH))
191
192 //! Van der Waals interaction treatment
193 enum {
194     evdwCUT, evdwSWITCH, evdwSHIFT, evdwUSER, evdwENCADSHIFT,
195     evdwPME, evdwNR
196 };
197 //! String corresponding to Van der Waals treatment
198 extern const char *evdw_names[evdwNR+1];
199 //! Macro for selecting correct string for VdW treatment
200 #define EVDWTYPE(e)    enum_name(e, evdwNR, evdw_names)
201
202 //! Type of long-range VdW treatment of combination rules
203 enum {
204     eljpmeGEOM, eljpmeLB, eljpmeNR
205 };
206 //! String for LJPME combination rule treatment
207 extern const char *eljpme_names[eljpmeNR+1];
208 //! Macro for correct LJPME comb rule name
209 #define ELJPMECOMBNAMES(e) enum_name(e, eljpmeNR, eljpme_names)
210
211 //! Macro to tell us whether we use LJPME
212 #define EVDW_PME(e) ((e) == evdwPME)
213
214 //! Neighborsearching algorithm
215 enum {
216     ensGRID, ensSIMPLE, ensNR
217 };
218 //! String corresponding to neighborsearching
219 extern const char *ens_names[ensNR+1];
220 //! Macro for correct NS algorithm
221 #define ENS(e)         enum_name(e, ensNR, ens_names)
222
223 /*! \brief Integrator algorithm
224  *
225  * eiSD2 has been removed, but we keep a renamed enum entry,
226  * so we can refuse to do MD with such .tpr files.
227  * eiVV is normal velocity verlet
228  * eiVVAK uses 1/2*(KE(t-dt/2)+KE(t+dt/2)) as the kinetic energy,
229  * and the half step kinetic energy for temperature control
230  */
231 enum {
232     eiMD, eiSteep, eiCG, eiBD, eiSD2_REMOVED, eiNM, eiLBFGS, eiTPI, eiTPIC, eiSD1, eiVV, eiVVAK, eiNR
233 };
234 //! Name of the integrator algorithm
235 extern const char *ei_names[eiNR+1];
236 //! Macro returning integrator string
237 #define EI(e)          enum_name(e, eiNR, ei_names)
238 //! Do we use velocity Verlet
239 #define EI_VV(e) ((e) == eiVV || (e) == eiVVAK)
240 //! Do we use molecular dynamics
241 #define EI_MD(e) ((e) == eiMD || EI_VV(e))
242 //! Do we use stochastic dynamics
243 #define EI_SD(e) ((e) == eiSD1)
244 //! Do we use any stochastic integrator
245 #define EI_RANDOM(e) (EI_SD(e) || (e) == eiBD)
246 /*above integrators may not conserve momenta*/
247 //! Do we use any type of dynamics
248 #define EI_DYNAMICS(e) (EI_MD(e) || EI_RANDOM(e))
249 //! Or do we use minimization
250 #define EI_ENERGY_MINIMIZATION(e) ((e) == eiSteep || (e) == eiCG || (e) == eiLBFGS)
251 //! Do we apply test particle insertion
252 #define EI_TPI(e) ((e) == eiTPI || (e) == eiTPIC)
253 //! Do we deal with particle velocities
254 #define EI_STATE_VELOCITY(e) (EI_MD(e) || EI_SD(e))
255
256 //! Constraint algorithm
257 enum {
258     econtLINCS, econtSHAKE, econtNR
259 };
260 //! String corresponding to constraint algorithm
261 extern const char *econstr_names[econtNR+1];
262 //! Macro to select the correct string
263 #define ECONSTRTYPE(e) enum_name(e, econtNR, econstr_names)
264
265 //! Distance restraint refinement algorithm
266 enum {
267     edrNone, edrSimple, edrEnsemble, edrNR
268 };
269 //! String corresponding to distance restraint algorithm
270 extern const char *edisre_names[edrNR+1];
271 //! Macro to select the right disre algorithm string
272 #define EDISRETYPE(e)  enum_name(e, edrNR, edisre_names)
273
274 //! Distance restraints weighting type
275 enum {
276     edrwConservative, edrwEqual, edrwNR
277 };
278 //! String corresponding to distance restraint weighting
279 extern const char *edisreweighting_names[edrwNR+1];
280 //! Macro corresponding to dr weighting
281 #define EDISREWEIGHTING(e)  enum_name(e, edrwNR, edisreweighting_names)
282
283 //! Combination rule algorithm.
284 enum {
285     eCOMB_NONE, eCOMB_GEOMETRIC, eCOMB_ARITHMETIC, eCOMB_GEOM_SIG_EPS, eCOMB_NR
286 };
287 //! String for combination rule algorithm
288 extern const char *ecomb_names[eCOMB_NR+1];
289 //! Macro to select the comb rule string
290 #define ECOMBNAME(e)   enum_name(e, eCOMB_NR, ecomb_names)
291
292 //! Van der Waals potential.
293 enum {
294     eNBF_NONE, eNBF_LJ, eNBF_BHAM, eNBF_NR
295 };
296 //! String corresponding to Van der Waals potential
297 extern const char *enbf_names[eNBF_NR+1];
298 //! Macro for correct VdW potential string
299 #define ENBFNAME(e)    enum_name(e, eNBF_NR, enbf_names)
300
301 //! Simulated tempering methods.
302 enum {
303     esimtempGEOMETRIC, esimtempEXPONENTIAL, esimtempLINEAR, esimtempNR
304 };
305 //! String corresponding to simulated tempering
306 extern const char *esimtemp_names[esimtempNR+1];
307 //! Macro for correct tempering string
308 #define ESIMTEMP(e)    enum_name(e, esimtempNR, esimtemp_names)
309
310 /*! \brief Free energy perturbation type
311  *
312  * efepNO, there are no evaluations at other states.
313  * efepYES, treated equivalently to efepSTATIC.
314  * efepSTATIC, then lambdas do not change during the simulation.
315  * efepSLOWGROWTH, then the states change monotonically
316  * throughout the simulation.
317  * efepEXPANDED, then expanded ensemble simulations are occuring.
318  */
319 enum {
320     efepNO, efepYES, efepSTATIC, efepSLOWGROWTH, efepEXPANDED, efepNR
321 };
322 //! String corresponding to FEP type.
323 extern const char *efep_names[efepNR+1];
324 //! Macro corresponding to FEP string.
325 #define EFEPTYPE(e)    enum_name(e, efepNR, efep_names)
326
327 //! Free energy pertubation coupling types.
328 enum {
329     efptFEP, efptMASS, efptCOUL, efptVDW, efptBONDED, efptRESTRAINT, efptTEMPERATURE, efptNR
330 };
331 //! String for FEP coupling type
332 extern const char *efpt_names[efptNR+1];
333 //! Long names for FEP coupling type
334 extern const char *efpt_singular_names[efptNR+1];
335
336 /*! \brief What to print for free energy calculations
337  *
338  * Printing the energy to the free energy dhdl file.
339  * YES is an alias to TOTAL, and
340  * will be converted in readir, so we never have to account for it in code.
341  */
342 enum {
343     edHdLPrintEnergyNO, edHdLPrintEnergyTOTAL, edHdLPrintEnergyPOTENTIAL, edHdLPrintEnergyYES, edHdLPrintEnergyNR
344 };
345 //! String corresponding to printing of free energy
346 extern const char *edHdLPrintEnergy_names[edHdLPrintEnergyNR+1];
347
348 /*! \brief How the lambda weights are calculated
349  *
350  * elamstatsMETROPOLIS - using the metropolis criteria
351  * elamstatsBARKER     - using the Barker critera for transition weights,
352  *                       also called unoptimized Bennett
353  * elamstatsMINVAR     - using Barker + minimum variance for weights
354  * elamstatsWL         - Wang-Landu (using visitation counts)
355  * elamstatsWWL        - Weighted Wang-Landau (using optimized Gibbs
356  *                       weighted visitation counts)
357  */
358 enum {
359     elamstatsNO, elamstatsMETROPOLIS, elamstatsBARKER, elamstatsMINVAR, elamstatsWL, elamstatsWWL, elamstatsNR
360 };
361 //! String corresponding to lambda weights
362 extern const char *elamstats_names[elamstatsNR+1];
363 //! Macro telling us whether we use expanded ensemble
364 #define ELAMSTATS_EXPANDED(e) ((e) > elamstatsNO)
365 //! Macro telling us whether we use some kind of Wang-Landau
366 #define EWL(e) ((e) == elamstatsWL || (e) == elamstatsWWL)
367
368 /*! \brief How moves in lambda are calculated
369  *
370  * elmovemcMETROPOLIS - using the Metropolis criteria, and 50% up and down
371  * elmovemcBARKER     - using the Barker criteria, and 50% up and down
372  * elmovemcGIBBS      - computing the transition using the marginalized
373  *                      probabilities of the lambdas
374  * elmovemcMETGIBBS   - computing the transition using the metropolized
375  *                      version of Gibbs (Monte Carlo Strategies in
376  *                      Scientific computing, Liu, p. 134)
377  */
378 enum {
379     elmcmoveNO, elmcmoveMETROPOLIS, elmcmoveBARKER, elmcmoveGIBBS, elmcmoveMETGIBBS, elmcmoveNR
380 };
381 //! String corresponding to lambda moves
382 extern const char *elmcmove_names[elmcmoveNR+1];
383
384 /*! \brief How we decide whether weights have reached equilibrium
385  *
386  * elmceqNO       - never stop, weights keep going
387  * elmceqYES      - fix the weights from the beginning; no movement
388  * elmceqWLDELTA  - stop when the WL-delta falls below a certain level
389  * elmceqNUMATLAM - stop when we have a certain number of samples at
390  *                  every step
391  * elmceqSTEPS    - stop when we've run a certain total number of steps
392  * elmceqSAMPLES  - stop when we've run a certain total number of samples
393  * elmceqRATIO    - stop when the ratio of samples (lowest to highest)
394  *                  is sufficiently large
395  */
396 enum {
397     elmceqNO, elmceqYES, elmceqWLDELTA, elmceqNUMATLAM, elmceqSTEPS, elmceqSAMPLES, elmceqRATIO, elmceqNR
398 };
399 //! String corresponding to equilibrium algorithm
400 extern const char *elmceq_names[elmceqNR+1];
401
402 /*! \brief separate_dhdl_file selection
403  *
404  * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
405  */
406 enum
407 {
408     esepdhdlfileYES, esepdhdlfileNO, esepdhdlfileNR
409 };
410 //! String corresponding to separate DHDL file selection
411 extern const char *separate_dhdl_file_names[esepdhdlfileNR+1];
412 //! Monster macro for DHDL file selection
413 #define SEPDHDLFILETYPE(e) enum_name(e, esepdhdlfileNR, separate_dhdl_file_names)
414
415 /*! \brief dhdl_derivatives selection \
416  *
417  * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
418  */
419 enum
420 {
421     edhdlderivativesYES, edhdlderivativesNO, edhdlderivativesNR
422 };
423 //! String for DHDL derivatives
424 extern const char *dhdl_derivatives_names[edhdlderivativesNR+1];
425 //! YAMM (Yet another monster macro)
426 #define DHDLDERIVATIVESTYPE(e) enum_name(e, edhdlderivativesNR, dhdl_derivatives_names)
427
428 /*! \brief Solvent model
429  *
430  * Distinguishes classical water types with 3 or 4 particles
431  */
432 enum {
433     esolNO, esolSPC, esolTIP4P, esolNR
434 };
435 //! String corresponding to solvent type
436 extern const char *esol_names[esolNR+1];
437 //! Macro lest we print the wrong solvent model string
438 #define ESOLTYPE(e)    enum_name(e, esolNR, esol_names)
439
440 //! Dispersion correction.
441 enum {
442     edispcNO, edispcEnerPres, edispcEner, edispcAllEnerPres, edispcAllEner, edispcNR
443 };
444 //! String corresponding to dispersion corrections
445 extern const char *edispc_names[edispcNR+1];
446 //! Macro for dispcorr string
447 #define EDISPCORR(e)   enum_name(e, edispcNR, edispc_names)
448
449 //! Center of mass motion removal algorithm.
450 enum {
451     ecmLINEAR, ecmANGULAR, ecmNO, ecmLINEAR_ACCELERATION_CORRECTION, ecmNR
452 };
453 //! String corresponding to COM removal
454 extern const char *ecm_names[ecmNR+1];
455 //! Macro for COM removal string
456 #define ECOM(e)        enum_name(e, ecmNR, ecm_names)
457
458 //! Algorithm for simulated annealing.
459 enum {
460     eannNO, eannSINGLE, eannPERIODIC, eannNR
461 };
462 //! String for simulated annealing
463 extern const char *eann_names[eannNR+1];
464 //! And macro for simulated annealing string
465 #define EANNEAL(e)      enum_name(e, eannNR, eann_names)
466
467 //! Wall types.
468 enum {
469     ewt93, ewt104, ewtTABLE, ewt126, ewtNR
470 };
471 //! String corresponding to wall type
472 extern const char *ewt_names[ewtNR+1];
473 //! Macro for wall type string
474 #define EWALLTYPE(e)   enum_name(e, ewtNR, ewt_names)
475
476 //! Pulling algorithm.
477 enum {
478     epullUMBRELLA, epullCONSTRAINT, epullCONST_F, epullFLATBOTTOM, epullFLATBOTTOMHIGH, epullEXTERNAL, epullNR
479 };
480 //! String for pulling algorithm
481 extern const char *epull_names[epullNR+1];
482 //! Macro for pulling string
483 #define EPULLTYPE(e)   enum_name(e, epullNR, epull_names)
484
485 //! Control of pull groups
486 enum {
487     epullgDIST, epullgDIR, epullgCYL, epullgDIRPBC, epullgDIRRELATIVE, epullgANGLE, epullgDIHEDRAL, epullgANGLEAXIS, epullgNR
488 };
489 //! String for pull groups
490 extern const char *epullg_names[epullgNR+1];
491 //! Macro for pull group string
492 #define EPULLGEOM(e)   enum_name(e, epullgNR, epullg_names)
493
494 //! Enforced rotation groups.
495 enum {
496     erotgISO, erotgISOPF,
497     erotgPM, erotgPMPF,
498     erotgRM, erotgRMPF,
499     erotgRM2, erotgRM2PF,
500     erotgFLEX, erotgFLEXT,
501     erotgFLEX2, erotgFLEX2T,
502     erotgNR
503 };
504 //! Rotation group names
505 extern const char *erotg_names[erotgNR+1];
506 //! Macro for rot group names
507 #define EROTGEOM(e)    enum_name(e, erotgNR, erotg_names)
508 //! String for rotation group origin names
509 extern const char *erotg_originnames[erotgNR+1];
510 //! Macro for rot group origin names
511 #define EROTORIGIN(e)  enum_name(e, erotgOriginNR, erotg_originnames)
512
513 //! Rotation group fitting type
514 enum {
515     erotgFitRMSD, erotgFitNORM, erotgFitPOT, erotgFitNR
516 };
517 //! String corresponding to rotation group fitting
518 extern const char *erotg_fitnames[erotgFitNR+1];
519 //! Macro for rot group fit names
520 #define EROTFIT(e)     enum_name(e, erotgFitNR, erotg_fitnames)
521
522 /*! \brief Direction along which ion/water swaps happen
523  *
524  * Part of "Computational Electrophysiology" (CompEL) setups
525  */
526 enum eSwaptype {
527     eswapNO, eswapX, eswapY, eswapZ, eSwapTypesNR
528 };
529 //! Names for swapping
530 extern const char *eSwapTypes_names[eSwapTypesNR+1];
531 //! Macro for swapping string
532 #define ESWAPTYPE(e)   enum_name(e, eSwapTypesNR, eSwapTypes_names)
533
534 /*! \brief Swap group splitting type
535  *
536  * These are just the fixed groups we need for any setup. In t_swap's grp
537  * entry after that follows the variable number of swap groups.
538  */
539 enum {
540     eGrpSplit0, eGrpSplit1, eGrpSolvent, eSwapFixedGrpNR
541 };
542 //! String for swap group splitting
543 extern const char *eSwapFixedGrp_names[eSwapFixedGrpNR+1];
544
545 //! QMMM methods.
546 enum {
547     eQMmethodAM1, eQMmethodPM3, eQMmethodRHF,
548     eQMmethodUHF, eQMmethodDFT, eQMmethodB3LYP, eQMmethodMP2, eQMmethodCASSCF, eQMmethodB3LYPLAN,
549     eQMmethodDIRECT, eQMmethodNR
550 };
551 //! String corresponding to QMMM methods
552 extern const char *eQMmethod_names[eQMmethodNR+1];
553 //! Macro to pick QMMM method name
554 #define EQMMETHOD(e)   enum_name(e, eQMmethodNR, eQMmethod_names)
555
556 //! QMMM basis function for QM part
557 enum {
558     eQMbasisSTO3G, eQMbasisSTO3G2, eQMbasis321G,
559     eQMbasis321Gp, eQMbasis321dGp, eQMbasis621G,
560     eQMbasis631G, eQMbasis631Gp, eQMbasis631dGp,
561     eQMbasis6311G, eQMbasisNR
562 };
563 //! Name for QMMM basis function
564 extern const char *eQMbasis_names[eQMbasisNR+1];
565 //! Macro to pick right basis function string
566 #define EQMBASIS(e)    enum_name(e, eQMbasisNR, eQMbasis_names)
567
568 //! QMMM scheme
569 enum {
570     eQMMMschemenormal, eQMMMschemeoniom, eQMMMschemeNR
571 };
572 //! QMMMM scheme names
573 extern const char *eQMMMscheme_names[eQMMMschemeNR+1];
574 //! Macro to pick QMMMM scheme name
575 #define EQMMMSCHEME(e) enum_name(e, eQMMMschemeNR, eQMMMscheme_names)
576
577 /*! \brief Neighborlist geometry type.
578  *
579  * Kernels will compute interactions between two particles,
580  * 3-center water, 4-center water or coarse-grained beads.
581  */
582 enum gmx_nblist_kernel_geometry
583 {
584     GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE,
585     GMX_NBLIST_GEOMETRY_WATER3_PARTICLE,
586     GMX_NBLIST_GEOMETRY_WATER3_WATER3,
587     GMX_NBLIST_GEOMETRY_WATER4_PARTICLE,
588     GMX_NBLIST_GEOMETRY_WATER4_WATER4,
589     GMX_NBLIST_GEOMETRY_CG_CG,
590     GMX_NBLIST_GEOMETRY_NR
591 };
592 //! String corresponding to nblist geometry names
593 extern const char *gmx_nblist_geometry_names[GMX_NBLIST_GEOMETRY_NR+1];
594
595 /*! \brief Types of electrostatics calculations
596  *
597  * Types of electrostatics calculations available inside nonbonded kernels.
598  * Note that these do NOT necessarily correspond to the user selections
599  * in the MDP file; many interactions for instance map to tabulated kernels.
600  */
601 enum gmx_nbkernel_elec
602 {
603     GMX_NBKERNEL_ELEC_NONE,
604     GMX_NBKERNEL_ELEC_COULOMB,
605     GMX_NBKERNEL_ELEC_REACTIONFIELD,
606     GMX_NBKERNEL_ELEC_CUBICSPLINETABLE,
607     GMX_NBKERNEL_ELEC_EWALD,
608     GMX_NBKERNEL_ELEC_NR
609 };
610 //! String corresponding to electrostatics kernels
611 extern const char *gmx_nbkernel_elec_names[GMX_NBKERNEL_ELEC_NR+1];
612
613 /*! \brief Types of vdw calculations available
614  *
615  * Types of vdw calculations available inside nonbonded kernels.
616  * Note that these do NOT necessarily correspond to the user selections
617  * in the MDP file; many interactions for instance map to tabulated kernels.
618  */
619 enum gmx_nbkernel_vdw
620 {
621     GMX_NBKERNEL_VDW_NONE,
622     GMX_NBKERNEL_VDW_LENNARDJONES,
623     GMX_NBKERNEL_VDW_BUCKINGHAM,
624     GMX_NBKERNEL_VDW_CUBICSPLINETABLE,
625     GMX_NBKERNEL_VDW_LJEWALD,
626     GMX_NBKERNEL_VDW_NR
627 };
628 //! String corresponding to VdW kernels
629 extern const char *gmx_nbkernel_vdw_names[GMX_NBKERNEL_VDW_NR+1];
630
631 //! \brief Types of interactions inside the neighborlist
632 enum gmx_nblist_interaction_type
633 {
634     GMX_NBLIST_INTERACTION_STANDARD,
635     GMX_NBLIST_INTERACTION_FREE_ENERGY,
636     GMX_NBLIST_INTERACTION_NR
637 };
638 //! String corresponding to interactions in neighborlist code
639 extern const char *gmx_nblist_interaction_names[GMX_NBLIST_INTERACTION_NR+1];
640
641 #endif /* GMX_MDTYPES_MD_ENUMS_H */