Make temperature and pressure coupling enums enum classes
[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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /*! \file
39  * \brief
40  * Declares enumerated types used throughout the code.
41  *
42  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
43  * \inpublicapi
44  * \ingroup module_mdtypes
45  */
46 #ifndef GMX_MDTYPES_MD_ENUMS_H
47 #define GMX_MDTYPES_MD_ENUMS_H
48
49 #include "gromacs/utility/basedefinitions.h"
50
51 /*! \brief Return a string from a list of strings
52  *
53  * If index if within 0 .. max_index-1 returns the corresponding string
54  * or "no name defined" otherwise, in other words this is a range-check that does
55  * not crash.
56  * \param[in] index     The index in the array
57  * \param[in] max_index The length of the array
58  * \param[in] names     The array
59  * \return the correct string or "no name defined"
60  */
61 const char* enum_name(int index, int max_index, const char* const names[]);
62
63 //! Boolean strings no or yes
64 extern const char* yesno_names[BOOL_NR + 1];
65
66 //! \brief The two compartments for CompEL setups.
67 enum eCompartment
68 {
69     eCompA,
70     eCompB,
71     eCompNR
72 };
73
74 /*! \brief The channels that define with their COM the compartment boundaries in CompEL setups.
75  *
76  * In principle one could also use modified setups with more than two channels.
77  */
78 enum eChannel
79 {
80     eChan0,
81     eChan1,
82     eChanNR
83 };
84
85 /*! \brief Temperature coupling type
86  *
87  * yes is an alias for berendsen
88  *
89  * Note: Keep `Count` as the second-to-last entry, and `Default` as the last entry -
90  *       this is needed to keep EnumerationWrapper, EnumerationArray and (de)serialization
91  *       working.
92  */
93 enum class TemperatureCoupling : int
94 {
95     No,
96     Berendsen,
97     NoseHoover,
98     Yes,
99     Andersen,
100     AndersenMassive,
101     VRescale,
102     Count,
103     Default = No
104 };
105 //! Return names of temperature coupling schemes
106 const char* enumValueToString(TemperatureCoupling enumValue);
107 //! Return whether this is andersen coupling
108 #define ETC_ANDERSEN(e) \
109     (((e) == TemperatureCoupling::AndersenMassive) || ((e) == TemperatureCoupling::Andersen))
110
111 /*! \brief Pressure coupling types
112  *
113  * isotropic is an alias for berendsen
114  *
115  * Note: Keep `Count` as the second-to-last entry, and `Default` as the last entry -
116  *       this is needed to keep EnumerationWrapper, EnumerationArray and (de)serialization
117  *       working.
118  */
119 enum class PressureCoupling : int
120 {
121     No,
122     Berendsen,
123     ParrinelloRahman,
124     Isotropic,
125     Mttk,
126     CRescale,
127     Count,
128     Default = No
129 };
130 //! Return names of pressure coupling schemes
131 const char* enumValueToString(PressureCoupling enumValue);
132
133 //! Flat-bottom posres geometries
134 enum
135 {
136     efbposresZERO,
137     efbposresSPHERE,
138     efbposresCYLINDER,
139     efbposresX,
140     efbposresY,
141     efbposresZ,
142     efbposresCYLINDERX,
143     efbposresCYLINDERY,
144     efbposresCYLINDERZ,
145     efbposresNR
146 };
147
148 //! Relative coordinate scaling type for position restraints.
149 enum
150 {
151     erscNO,
152     erscALL,
153     erscCOM,
154     erscNR
155 };
156 //! String corresponding to relativ coordinate scaling.
157 extern const char* erefscaling_names[erscNR + 1];
158 //! Macro to select correct coordinate scaling string.
159 #define EREFSCALINGTYPE(e) enum_name(e, erscNR, erefscaling_names)
160
161 //! Trotter decomposition extended variable parts.
162 enum
163 {
164     etrtNONE,
165     etrtNHC,
166     etrtBAROV,
167     etrtBARONHC,
168     etrtNHC2,
169     etrtBAROV2,
170     etrtBARONHC2,
171     etrtVELOCITY1,
172     etrtVELOCITY2,
173     etrtPOSITION,
174     etrtSKIPALL,
175     etrtNR
176 };
177
178 //! Sequenced parts of the trotter decomposition.
179 enum
180 {
181     ettTSEQ0,
182     ettTSEQ1,
183     ettTSEQ2,
184     ettTSEQ3,
185     ettTSEQ4,
186     ettTSEQMAX
187 };
188
189 //! Pressure coupling type
190 enum
191 {
192     epctISOTROPIC,
193     epctSEMIISOTROPIC,
194     epctANISOTROPIC,
195     epctSURFACETENSION,
196     epctNR
197 };
198 //! String corresponding to pressure coupling type
199 extern const char* epcoupltype_names[epctNR + 1];
200 //! Macro to select the right string for pcoupl type
201 #define EPCOUPLTYPETYPE(e) enum_name(e, epctNR, epcoupltype_names)
202
203 //! \\brief Cutoff scheme
204 enum
205 {
206     ecutsVERLET,
207     ecutsGROUP,
208     ecutsNR
209 };
210 //! String corresponding to cutoff scheme
211 extern const char* ecutscheme_names[ecutsNR + 1];
212 //! Macro to select the right string for cutoff scheme
213 #define ECUTSCHEME(e) enum_name(e, ecutsNR, ecutscheme_names)
214
215 /*! \brief Coulomb / VdW interaction modifiers.
216  *
217  * grompp replaces eintmodPOTSHIFT_VERLET_UNSUPPORTED by eintmodPOTSHIFT.
218  * Exactcutoff is only used by Reaction-field-zero, and is not user-selectable.
219  */
220 enum eintmod
221 {
222     eintmodPOTSHIFT_VERLET_UNSUPPORTED,
223     eintmodPOTSHIFT,
224     eintmodNONE,
225     eintmodPOTSWITCH,
226     eintmodEXACTCUTOFF,
227     eintmodFORCESWITCH,
228     eintmodNR
229 };
230 //! String corresponding to interaction modifiers
231 extern const char* eintmod_names[eintmodNR + 1];
232 //! Macro to select the correct string for modifiers
233 #define INTMODIFIER(e) enum_name(e, eintmodNR, eintmod_names)
234
235 /*! \brief Cut-off treatment for Coulomb */
236 enum
237 {
238     eelCUT,
239     eelRF,
240     eelGRF_NOTUSED,
241     eelPME,
242     eelEWALD,
243     eelP3M_AD,
244     eelPOISSON,
245     eelSWITCH,
246     eelSHIFT,
247     eelUSER,
248     eelGB_NOTUSED,
249     eelRF_NEC_UNSUPPORTED,
250     eelENCADSHIFT_NOTUSED,
251     eelPMEUSER,
252     eelPMESWITCH,
253     eelPMEUSERSWITCH,
254     eelRF_ZERO,
255     eelNR
256 };
257 //! String corresponding to Coulomb treatment
258 extern const char* eel_names[eelNR + 1];
259 //! Macro for correct string for Coulomb treatment
260 #define EELTYPE(e) enum_name(e, eelNR, eel_names)
261
262 //! Ewald geometry.
263 enum
264 {
265     eewg3D,
266     eewg3DC,
267     eewgNR
268 };
269 //! String corresponding to Ewald geometry
270 extern const char* eewg_names[eewgNR + 1];
271
272 //! Macro telling us whether we use reaction field
273 #define EEL_RF(e) \
274     ((e) == eelRF || (e) == eelGRF_NOTUSED || (e) == eelRF_NEC_UNSUPPORTED || (e) == eelRF_ZERO)
275
276 //! Macro telling us whether we use PME
277 #define EEL_PME(e) \
278     ((e) == eelPME || (e) == eelPMESWITCH || (e) == eelPMEUSER || (e) == eelPMEUSERSWITCH || (e) == eelP3M_AD)
279 //! Macro telling us whether we use PME or full Ewald
280 #define EEL_PME_EWALD(e) (EEL_PME(e) || (e) == eelEWALD)
281 //! Macro telling us whether we use full electrostatics of any sort
282 #define EEL_FULL(e) (EEL_PME_EWALD(e) || (e) == eelPOISSON)
283 //! Macro telling us whether we use user defined electrostatics
284 #define EEL_USER(e) ((e) == eelUSER || (e) == eelPMEUSER || (e) == (eelPMEUSERSWITCH))
285
286 //! Van der Waals interaction treatment
287 enum
288 {
289     evdwCUT,
290     evdwSWITCH,
291     evdwSHIFT,
292     evdwUSER,
293     evdwENCADSHIFT_UNUSED,
294     evdwPME,
295     evdwNR
296 };
297 //! String corresponding to Van der Waals treatment
298 extern const char* evdw_names[evdwNR + 1];
299 //! Macro for selecting correct string for VdW treatment
300 #define EVDWTYPE(e) enum_name(e, evdwNR, evdw_names)
301
302 //! Type of long-range VdW treatment of combination rules
303 enum
304 {
305     eljpmeGEOM,
306     eljpmeLB,
307     eljpmeNR
308 };
309 //! String for LJPME combination rule treatment
310 extern const char* eljpme_names[eljpmeNR + 1];
311 //! Macro for correct LJPME comb rule name
312 #define ELJPMECOMBNAMES(e) enum_name(e, eljpmeNR, eljpme_names)
313
314 //! Macro to tell us whether we use LJPME
315 #define EVDW_PME(e) ((e) == evdwPME)
316
317 /*! \brief Integrator algorithm
318  *
319  * eiSD2 has been removed, but we keep a renamed enum entry,
320  * so we can refuse to do MD with such .tpr files.
321  * eiVV is normal velocity verlet
322  * eiVVAK uses 1/2*(KE(t-dt/2)+KE(t+dt/2)) as the kinetic energy,
323  * and the half step kinetic energy for temperature control
324  */
325 enum
326 {
327     eiMD,
328     eiSteep,
329     eiCG,
330     eiBD,
331     eiSD2_REMOVED,
332     eiNM,
333     eiLBFGS,
334     eiTPI,
335     eiTPIC,
336     eiSD1,
337     eiVV,
338     eiVVAK,
339     eiMimic,
340     eiNR
341 };
342 //! Name of the integrator algorithm
343 extern const char* ei_names[eiNR + 1];
344 //! Macro returning integrator string
345 #define EI(e) enum_name(e, eiNR, ei_names)
346 //! Do we use MiMiC QM/MM?
347 #define EI_MIMIC(e) ((e) == eiMimic)
348 //! Do we use velocity Verlet
349 #define EI_VV(e) ((e) == eiVV || (e) == eiVVAK)
350 //! Do we use molecular dynamics
351 #define EI_MD(e) ((e) == eiMD || EI_VV(e) || EI_MIMIC(e))
352 //! Do we use stochastic dynamics
353 #define EI_SD(e) ((e) == eiSD1)
354 //! Do we use any stochastic integrator
355 #define EI_RANDOM(e) (EI_SD(e) || (e) == eiBD)
356 /*above integrators may not conserve momenta*/
357 //! Do we use any type of dynamics
358 #define EI_DYNAMICS(e) (EI_MD(e) || EI_RANDOM(e))
359 //! Or do we use minimization
360 #define EI_ENERGY_MINIMIZATION(e) ((e) == eiSteep || (e) == eiCG || (e) == eiLBFGS)
361 //! Do we apply test particle insertion
362 #define EI_TPI(e) ((e) == eiTPI || (e) == eiTPIC)
363 //! Do we deal with particle velocities
364 #define EI_STATE_VELOCITY(e) (EI_MD(e) || EI_SD(e))
365
366 //! Constraint algorithm
367 enum
368 {
369     econtLINCS,
370     econtSHAKE,
371     econtNR
372 };
373 //! String corresponding to constraint algorithm
374 extern const char* econstr_names[econtNR + 1];
375 //! Macro to select the correct string
376 #define ECONSTRTYPE(e) enum_name(e, econtNR, econstr_names)
377
378 //! Distance restraint refinement algorithm
379 enum
380 {
381     edrNone,
382     edrSimple,
383     edrEnsemble,
384     edrNR
385 };
386 //! String corresponding to distance restraint algorithm
387 extern const char* edisre_names[edrNR + 1];
388 //! Macro to select the right disre algorithm string
389 #define EDISRETYPE(e) enum_name(e, edrNR, edisre_names)
390
391 //! Distance restraints weighting type
392 enum
393 {
394     edrwConservative,
395     edrwEqual,
396     edrwNR
397 };
398 //! String corresponding to distance restraint weighting
399 extern const char* edisreweighting_names[edrwNR + 1];
400 //! Macro corresponding to dr weighting
401 #define EDISREWEIGHTING(e) enum_name(e, edrwNR, edisreweighting_names)
402
403 //! Combination rule algorithm.
404 enum
405 {
406     eCOMB_NONE,
407     eCOMB_GEOMETRIC,
408     eCOMB_ARITHMETIC,
409     eCOMB_GEOM_SIG_EPS,
410     eCOMB_NR
411 };
412 //! String for combination rule algorithm
413 extern const char* ecomb_names[eCOMB_NR + 1];
414 //! Macro to select the comb rule string
415 #define ECOMBNAME(e) enum_name(e, eCOMB_NR, ecomb_names)
416
417 //! Van der Waals potential.
418 enum
419 {
420     eNBF_NONE,
421     eNBF_LJ,
422     eNBF_BHAM,
423     eNBF_NR
424 };
425 //! String corresponding to Van der Waals potential
426 extern const char* enbf_names[eNBF_NR + 1];
427 //! Macro for correct VdW potential string
428 #define ENBFNAME(e) enum_name(e, eNBF_NR, enbf_names)
429
430 //! Simulated tempering methods.
431 enum
432 {
433     esimtempGEOMETRIC,
434     esimtempEXPONENTIAL,
435     esimtempLINEAR,
436     esimtempNR
437 };
438 //! String corresponding to simulated tempering
439 extern const char* esimtemp_names[esimtempNR + 1];
440 //! Macro for correct tempering string
441 #define ESIMTEMP(e) enum_name(e, esimtempNR, esimtemp_names)
442
443 /*! \brief Free energy perturbation type
444  *
445  * efepNO, there are no evaluations at other states.
446  * efepYES, treated equivalently to efepSTATIC.
447  * efepSTATIC, then lambdas do not change during the simulation.
448  * efepSLOWGROWTH, then the states change monotonically
449  * throughout the simulation.
450  * efepEXPANDED, then expanded ensemble simulations are occuring.
451  */
452 enum
453 {
454     efepNO,
455     efepYES,
456     efepSTATIC,
457     efepSLOWGROWTH,
458     efepEXPANDED,
459     efepNR
460 };
461 //! String corresponding to FEP type.
462 extern const char* efep_names[efepNR + 1];
463 //! Macro corresponding to FEP string.
464 #define EFEPTYPE(e) enum_name(e, efepNR, efep_names)
465
466 //! Free energy pertubation coupling types.
467 enum
468 {
469     efptFEP,
470     efptMASS,
471     efptCOUL,
472     efptVDW,
473     efptBONDED,
474     efptRESTRAINT,
475     efptTEMPERATURE,
476     efptNR
477 };
478 //! String for FEP coupling type
479 extern const char* efpt_names[efptNR + 1];
480 //! Long names for FEP coupling type
481 extern const char* efpt_singular_names[efptNR + 1];
482
483 /*! \brief What to print for free energy calculations
484  *
485  * Printing the energy to the free energy dhdl file.
486  * YES is an alias to TOTAL, and
487  * will be converted in readir, so we never have to account for it in code.
488  */
489 enum
490 {
491     edHdLPrintEnergyNO,
492     edHdLPrintEnergyTOTAL,
493     edHdLPrintEnergyPOTENTIAL,
494     edHdLPrintEnergyYES,
495     edHdLPrintEnergyNR
496 };
497 //! String corresponding to printing of free energy
498 extern const char* edHdLPrintEnergy_names[edHdLPrintEnergyNR + 1];
499
500 /*! \brief How the lambda weights are calculated
501  *
502  * elamstatsMETROPOLIS - using the metropolis criteria
503  * elamstatsBARKER     - using the Barker critera for transition weights,
504  *                       also called unoptimized Bennett
505  * elamstatsMINVAR     - using Barker + minimum variance for weights
506  * elamstatsWL         - Wang-Landu (using visitation counts)
507  * elamstatsWWL        - Weighted Wang-Landau (using optimized Gibbs
508  *                       weighted visitation counts)
509  */
510 enum
511 {
512     elamstatsNO,
513     elamstatsMETROPOLIS,
514     elamstatsBARKER,
515     elamstatsMINVAR,
516     elamstatsWL,
517     elamstatsWWL,
518     elamstatsNR
519 };
520 //! String corresponding to lambda weights
521 extern const char* elamstats_names[elamstatsNR + 1];
522 //! Macro telling us whether we use expanded ensemble
523 #define ELAMSTATS_EXPANDED(e) ((e) > elamstatsNO)
524 //! Macro telling us whether we use some kind of Wang-Landau
525 #define EWL(e) ((e) == elamstatsWL || (e) == elamstatsWWL)
526
527 /*! \brief How moves in lambda are calculated
528  *
529  * elmovemcMETROPOLIS - using the Metropolis criteria, and 50% up and down
530  * elmovemcBARKER     - using the Barker criteria, and 50% up and down
531  * elmovemcGIBBS      - computing the transition using the marginalized
532  *                      probabilities of the lambdas
533  * elmovemcMETGIBBS   - computing the transition using the metropolized
534  *                      version of Gibbs (Monte Carlo Strategies in
535  *                      Scientific computing, Liu, p. 134)
536  */
537 enum
538 {
539     elmcmoveNO,
540     elmcmoveMETROPOLIS,
541     elmcmoveBARKER,
542     elmcmoveGIBBS,
543     elmcmoveMETGIBBS,
544     elmcmoveNR
545 };
546 //! String corresponding to lambda moves
547 extern const char* elmcmove_names[elmcmoveNR + 1];
548
549 /*! \brief How we decide whether weights have reached equilibrium
550  *
551  * elmceqNO       - never stop, weights keep going
552  * elmceqYES      - fix the weights from the beginning; no movement
553  * elmceqWLDELTA  - stop when the WL-delta falls below a certain level
554  * elmceqNUMATLAM - stop when we have a certain number of samples at
555  *                  every step
556  * elmceqSTEPS    - stop when we've run a certain total number of steps
557  * elmceqSAMPLES  - stop when we've run a certain total number of samples
558  * elmceqRATIO    - stop when the ratio of samples (lowest to highest)
559  *                  is sufficiently large
560  */
561 enum
562 {
563     elmceqNO,
564     elmceqYES,
565     elmceqWLDELTA,
566     elmceqNUMATLAM,
567     elmceqSTEPS,
568     elmceqSAMPLES,
569     elmceqRATIO,
570     elmceqNR
571 };
572 //! String corresponding to equilibrium algorithm
573 extern const char* elmceq_names[elmceqNR + 1];
574
575 /*! \brief separate_dhdl_file selection
576  *
577  * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
578  */
579 enum
580 {
581     esepdhdlfileYES,
582     esepdhdlfileNO,
583     esepdhdlfileNR
584 };
585 //! String corresponding to separate DHDL file selection
586 extern const char* separate_dhdl_file_names[esepdhdlfileNR + 1];
587 //! Monster macro for DHDL file selection
588 #define SEPDHDLFILETYPE(e) enum_name(e, esepdhdlfileNR, separate_dhdl_file_names)
589
590 /*! \brief dhdl_derivatives selection \
591  *
592  * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
593  */
594 enum
595 {
596     edhdlderivativesYES,
597     edhdlderivativesNO,
598     edhdlderivativesNR
599 };
600 //! String for DHDL derivatives
601 extern const char* dhdl_derivatives_names[edhdlderivativesNR + 1];
602 //! YAMM (Yet another monster macro)
603 #define DHDLDERIVATIVESTYPE(e) enum_name(e, edhdlderivativesNR, dhdl_derivatives_names)
604
605 /*! \brief Solvent model
606  *
607  * Distinguishes classical water types with 3 or 4 particles
608  */
609 enum
610 {
611     esolNO,
612     esolSPC,
613     esolTIP4P,
614     esolNR
615 };
616 //! String corresponding to solvent type
617 extern const char* esol_names[esolNR + 1];
618 //! Macro lest we print the wrong solvent model string
619 #define ESOLTYPE(e) enum_name(e, esolNR, esol_names)
620
621 //! Dispersion correction.
622 enum
623 {
624     edispcNO,
625     edispcEnerPres,
626     edispcEner,
627     edispcAllEnerPres,
628     edispcAllEner,
629     edispcNR
630 };
631 //! String corresponding to dispersion corrections
632 extern const char* edispc_names[edispcNR + 1];
633 //! Macro for dispcorr string
634 #define EDISPCORR(e) enum_name(e, edispcNR, edispc_names)
635
636 //! Center of mass motion removal algorithm.
637 enum
638 {
639     ecmLINEAR,
640     ecmANGULAR,
641     ecmNO,
642     ecmLINEAR_ACCELERATION_CORRECTION,
643     ecmNR
644 };
645 //! String corresponding to COM removal
646 extern const char* ecm_names[ecmNR + 1];
647 //! Macro for COM removal string
648 #define ECOM(e) enum_name(e, ecmNR, ecm_names)
649
650 //! Algorithm for simulated annealing.
651 enum
652 {
653     eannNO,
654     eannSINGLE,
655     eannPERIODIC,
656     eannNR
657 };
658 //! String for simulated annealing
659 extern const char* eann_names[eannNR + 1];
660 //! And macro for simulated annealing string
661 #define EANNEAL(e) enum_name(e, eannNR, eann_names)
662
663 //! Wall types.
664 enum
665 {
666     ewt93,
667     ewt104,
668     ewtTABLE,
669     ewt126,
670     ewtNR
671 };
672 //! String corresponding to wall type
673 extern const char* ewt_names[ewtNR + 1];
674 //! Macro for wall type string
675 #define EWALLTYPE(e) enum_name(e, ewtNR, ewt_names)
676
677 //! Pulling algorithm.
678 enum
679 {
680     epullUMBRELLA,
681     epullCONSTRAINT,
682     epullCONST_F,
683     epullFLATBOTTOM,
684     epullFLATBOTTOMHIGH,
685     epullEXTERNAL,
686     epullNR
687 };
688 //! String for pulling algorithm
689 extern const char* epull_names[epullNR + 1];
690 //! Macro for pulling string
691 #define EPULLTYPE(e) enum_name(e, epullNR, epull_names)
692
693 //! Control of pull groups
694 enum
695 {
696     epullgDIST,
697     epullgDIR,
698     epullgCYL,
699     epullgDIRPBC,
700     epullgDIRRELATIVE,
701     epullgANGLE,
702     epullgDIHEDRAL,
703     epullgANGLEAXIS,
704     epullgNR
705 };
706 //! String for pull groups
707 extern const char* epullg_names[epullgNR + 1];
708 //! Macro for pull group string
709 #define EPULLGEOM(e) enum_name(e, epullgNR, epullg_names)
710
711 //! Enforced rotation groups.
712 enum
713 {
714     erotgISO,
715     erotgISOPF,
716     erotgPM,
717     erotgPMPF,
718     erotgRM,
719     erotgRMPF,
720     erotgRM2,
721     erotgRM2PF,
722     erotgFLEX,
723     erotgFLEXT,
724     erotgFLEX2,
725     erotgFLEX2T,
726     erotgNR
727 };
728 //! Rotation group names
729 extern const char* erotg_names[erotgNR + 1];
730 //! Macro for rot group names
731 #define EROTGEOM(e) enum_name(e, erotgNR, erotg_names)
732 //! String for rotation group origin names
733 extern const char* erotg_originnames[erotgNR + 1];
734 //! Macro for rot group origin names
735 #define EROTORIGIN(e) enum_name(e, erotgOriginNR, erotg_originnames)
736
737 //! Rotation group fitting type
738 enum
739 {
740     erotgFitRMSD,
741     erotgFitNORM,
742     erotgFitPOT,
743     erotgFitNR
744 };
745 //! String corresponding to rotation group fitting
746 extern const char* erotg_fitnames[erotgFitNR + 1];
747 //! Macro for rot group fit names
748 #define EROTFIT(e) enum_name(e, erotgFitNR, erotg_fitnames)
749
750 /*! \brief Direction along which ion/water swaps happen
751  *
752  * Part of "Computational Electrophysiology" (CompEL) setups
753  */
754 enum eSwaptype
755 {
756     eswapNO,
757     eswapX,
758     eswapY,
759     eswapZ,
760     eSwapTypesNR
761 };
762 //! Names for swapping
763 extern const char* eSwapTypes_names[eSwapTypesNR + 1];
764 //! Macro for swapping string
765 #define ESWAPTYPE(e) enum_name(e, eSwapTypesNR, eSwapTypes_names)
766
767 /*! \brief Swap group splitting type
768  *
769  * These are just the fixed groups we need for any setup. In t_swap's grp
770  * entry after that follows the variable number of swap groups.
771  */
772 enum
773 {
774     eGrpSplit0,
775     eGrpSplit1,
776     eGrpSolvent,
777     eSwapFixedGrpNR
778 };
779 //! String for swap group splitting
780 extern const char* eSwapFixedGrp_names[eSwapFixedGrpNR + 1];
781
782 /*! \brief Neighborlist geometry type.
783  *
784  * Kernels will compute interactions between two particles,
785  * 3-center water, 4-center water or coarse-grained beads.
786  */
787 enum gmx_nblist_kernel_geometry
788 {
789     GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE,
790     GMX_NBLIST_GEOMETRY_WATER3_PARTICLE,
791     GMX_NBLIST_GEOMETRY_WATER3_WATER3,
792     GMX_NBLIST_GEOMETRY_WATER4_PARTICLE,
793     GMX_NBLIST_GEOMETRY_WATER4_WATER4,
794     GMX_NBLIST_GEOMETRY_CG_CG,
795     GMX_NBLIST_GEOMETRY_NR
796 };
797 //! String corresponding to nblist geometry names
798 extern const char* gmx_nblist_geometry_names[GMX_NBLIST_GEOMETRY_NR + 1];
799
800 /*! \brief Types of electrostatics calculations
801  *
802  * Types of electrostatics calculations available inside nonbonded kernels.
803  * Note that these do NOT necessarily correspond to the user selections
804  * in the MDP file; many interactions for instance map to tabulated kernels.
805  */
806 enum gmx_nbkernel_elec
807 {
808     GMX_NBKERNEL_ELEC_NONE,
809     GMX_NBKERNEL_ELEC_COULOMB,
810     GMX_NBKERNEL_ELEC_REACTIONFIELD,
811     GMX_NBKERNEL_ELEC_CUBICSPLINETABLE,
812     GMX_NBKERNEL_ELEC_EWALD,
813     GMX_NBKERNEL_ELEC_NR
814 };
815 //! String corresponding to electrostatics kernels
816 extern const char* gmx_nbkernel_elec_names[GMX_NBKERNEL_ELEC_NR + 1];
817
818 /*! \brief Types of vdw calculations available
819  *
820  * Types of vdw calculations available inside nonbonded kernels.
821  * Note that these do NOT necessarily correspond to the user selections
822  * in the MDP file; many interactions for instance map to tabulated kernels.
823  */
824 enum gmx_nbkernel_vdw
825 {
826     GMX_NBKERNEL_VDW_NONE,
827     GMX_NBKERNEL_VDW_LENNARDJONES,
828     GMX_NBKERNEL_VDW_BUCKINGHAM,
829     GMX_NBKERNEL_VDW_CUBICSPLINETABLE,
830     GMX_NBKERNEL_VDW_LJEWALD,
831     GMX_NBKERNEL_VDW_NR
832 };
833 //! String corresponding to VdW kernels
834 extern const char* gmx_nbkernel_vdw_names[GMX_NBKERNEL_VDW_NR + 1];
835
836 //! \brief Types of interactions inside the neighborlist
837 enum gmx_nblist_interaction_type
838 {
839     GMX_NBLIST_INTERACTION_STANDARD,
840     GMX_NBLIST_INTERACTION_FREE_ENERGY,
841     GMX_NBLIST_INTERACTION_NR
842 };
843 //! String corresponding to interactions in neighborlist code
844 extern const char* gmx_nblist_interaction_names[GMX_NBLIST_INTERACTION_NR + 1];
845
846 #endif /* GMX_MDTYPES_MD_ENUMS_H */