2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
40 * Declares enumerated types used throughout the code.
42 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
44 * \ingroup module_mdtypes
46 #ifndef GMX_MDTYPES_MD_ENUMS_H
47 #define GMX_MDTYPES_MD_ENUMS_H
49 #include "gromacs/utility/basedefinitions.h"
51 /*! \brief Return a string from a list of strings
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
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"
61 const char* enum_name(int index
, int max_index
, const char* names
[]);
63 //! Boolean strings no or yes
64 extern const char* yesno_names
[BOOL_NR
+ 1];
66 //! \brief The two compartments for CompEL setups.
74 /*! \brief The channels that define with their COM the compartment boundaries in CompEL setups.
76 * In principle one could also use modified setups with more than two channels.
85 /*! \brief Temperature coupling type
87 * yes is an alias for berendsen
100 //! Strings corresponding to temperatyre coupling types
101 extern const char* etcoupl_names
[etcNR
+ 1];
102 //! Macro for selecting t coupling string
103 #define ETCOUPLTYPE(e) enum_name(e, etcNR, etcoupl_names)
104 //! Return whether this is andersen coupling
105 #define ETC_ANDERSEN(e) (((e) == etcANDERSENMASSIVE) || ((e) == etcANDERSEN))
107 /*! \brief Pressure coupling types
109 * isotropic is an alias for berendsen
120 //! String corresponding to pressure coupling algorithm
121 extern const char* epcoupl_names
[epcNR
+ 1];
122 //! Macro to return the correct pcoupling string
123 #define EPCOUPLTYPE(e) enum_name(e, epcNR, epcoupl_names)
125 //! Flat-bottom posres geometries
140 //! Relative coordinate scaling type for position restraints.
148 //! String corresponding to relativ coordinate scaling.
149 extern const char* erefscaling_names
[erscNR
+ 1];
150 //! Macro to select correct coordinate scaling string.
151 #define EREFSCALINGTYPE(e) enum_name(e, erscNR, erefscaling_names)
153 //! Trotter decomposition extended variable parts.
170 //! Sequenced parts of the trotter decomposition.
181 //! Pressure coupling type
190 //! String corresponding to pressure coupling type
191 extern const char* epcoupltype_names
[epctNR
+ 1];
192 //! Macro to select the right string for pcoupl type
193 #define EPCOUPLTYPETYPE(e) enum_name(e, epctNR, epcoupltype_names)
195 //! \\brief Cutoff scheme
202 //! String corresponding to cutoff scheme
203 extern const char* ecutscheme_names
[ecutsNR
+ 1];
204 //! Macro to select the right string for cutoff scheme
205 #define ECUTSCHEME(e) enum_name(e, ecutsNR, ecutscheme_names)
207 /*! \brief Coulomb / VdW interaction modifiers.
209 * grompp replaces eintmodPOTSHIFT_VERLET_UNSUPPORTED by eintmodPOTSHIFT.
210 * Exactcutoff is only used by Reaction-field-zero, and is not user-selectable.
214 eintmodPOTSHIFT_VERLET_UNSUPPORTED
,
222 //! String corresponding to interaction modifiers
223 extern const char* eintmod_names
[eintmodNR
+ 1];
224 //! Macro to select the correct string for modifiers
225 #define INTMODIFIER(e) enum_name(e, eintmodNR, eintmod_names)
227 /*! \brief Cut-off treatment for Coulomb */
241 eelRF_NEC_UNSUPPORTED
,
242 eelENCADSHIFT_NOTUSED
,
249 //! String corresponding to Coulomb treatment
250 extern const char* eel_names
[eelNR
+ 1];
251 //! Macro for correct string for Coulomb treatment
252 #define EELTYPE(e) enum_name(e, eelNR, eel_names)
261 //! String corresponding to Ewald geometry
262 extern const char* eewg_names
[eewgNR
+ 1];
264 //! Macro telling us whether we use reaction field
266 ((e) == eelRF || (e) == eelGRF_NOTUSED || (e) == eelRF_NEC_UNSUPPORTED || (e) == eelRF_ZERO)
268 //! Macro telling us whether we use PME
270 ((e) == eelPME || (e) == eelPMESWITCH || (e) == eelPMEUSER || (e) == eelPMEUSERSWITCH || (e) == eelP3M_AD)
271 //! Macro telling us whether we use PME or full Ewald
272 #define EEL_PME_EWALD(e) (EEL_PME(e) || (e) == eelEWALD)
273 //! Macro telling us whether we use full electrostatics of any sort
274 #define EEL_FULL(e) (EEL_PME_EWALD(e) || (e) == eelPOISSON)
275 //! Macro telling us whether we use user defined electrostatics
276 #define EEL_USER(e) ((e) == eelUSER || (e) == eelPMEUSER || (e) == (eelPMEUSERSWITCH))
278 //! Van der Waals interaction treatment
285 evdwENCADSHIFT_UNUSED
,
289 //! String corresponding to Van der Waals treatment
290 extern const char* evdw_names
[evdwNR
+ 1];
291 //! Macro for selecting correct string for VdW treatment
292 #define EVDWTYPE(e) enum_name(e, evdwNR, evdw_names)
294 //! Type of long-range VdW treatment of combination rules
301 //! String for LJPME combination rule treatment
302 extern const char* eljpme_names
[eljpmeNR
+ 1];
303 //! Macro for correct LJPME comb rule name
304 #define ELJPMECOMBNAMES(e) enum_name(e, eljpmeNR, eljpme_names)
306 //! Macro to tell us whether we use LJPME
307 #define EVDW_PME(e) ((e) == evdwPME)
309 /*! \brief Integrator algorithm
311 * eiSD2 has been removed, but we keep a renamed enum entry,
312 * so we can refuse to do MD with such .tpr files.
313 * eiVV is normal velocity verlet
314 * eiVVAK uses 1/2*(KE(t-dt/2)+KE(t+dt/2)) as the kinetic energy,
315 * and the half step kinetic energy for temperature control
334 //! Name of the integrator algorithm
335 extern const char* ei_names
[eiNR
+ 1];
336 //! Macro returning integrator string
337 #define EI(e) enum_name(e, eiNR, ei_names)
338 //! Do we use MiMiC QM/MM?
339 #define EI_MIMIC(e) ((e) == eiMimic)
340 //! Do we use velocity Verlet
341 #define EI_VV(e) ((e) == eiVV || (e) == eiVVAK)
342 //! Do we use molecular dynamics
343 #define EI_MD(e) ((e) == eiMD || EI_VV(e) || EI_MIMIC(e))
344 //! Do we use stochastic dynamics
345 #define EI_SD(e) ((e) == eiSD1)
346 //! Do we use any stochastic integrator
347 #define EI_RANDOM(e) (EI_SD(e) || (e) == eiBD)
348 /*above integrators may not conserve momenta*/
349 //! Do we use any type of dynamics
350 #define EI_DYNAMICS(e) (EI_MD(e) || EI_RANDOM(e))
351 //! Or do we use minimization
352 #define EI_ENERGY_MINIMIZATION(e) ((e) == eiSteep || (e) == eiCG || (e) == eiLBFGS)
353 //! Do we apply test particle insertion
354 #define EI_TPI(e) ((e) == eiTPI || (e) == eiTPIC)
355 //! Do we deal with particle velocities
356 #define EI_STATE_VELOCITY(e) (EI_MD(e) || EI_SD(e))
358 //! Constraint algorithm
365 //! String corresponding to constraint algorithm
366 extern const char* econstr_names
[econtNR
+ 1];
367 //! Macro to select the correct string
368 #define ECONSTRTYPE(e) enum_name(e, econtNR, econstr_names)
370 //! Distance restraint refinement algorithm
378 //! String corresponding to distance restraint algorithm
379 extern const char* edisre_names
[edrNR
+ 1];
380 //! Macro to select the right disre algorithm string
381 #define EDISRETYPE(e) enum_name(e, edrNR, edisre_names)
383 //! Distance restraints weighting type
390 //! String corresponding to distance restraint weighting
391 extern const char* edisreweighting_names
[edrwNR
+ 1];
392 //! Macro corresponding to dr weighting
393 #define EDISREWEIGHTING(e) enum_name(e, edrwNR, edisreweighting_names)
395 //! Combination rule algorithm.
404 //! String for combination rule algorithm
405 extern const char* ecomb_names
[eCOMB_NR
+ 1];
406 //! Macro to select the comb rule string
407 #define ECOMBNAME(e) enum_name(e, eCOMB_NR, ecomb_names)
409 //! Van der Waals potential.
417 //! String corresponding to Van der Waals potential
418 extern const char* enbf_names
[eNBF_NR
+ 1];
419 //! Macro for correct VdW potential string
420 #define ENBFNAME(e) enum_name(e, eNBF_NR, enbf_names)
422 //! Simulated tempering methods.
430 //! String corresponding to simulated tempering
431 extern const char* esimtemp_names
[esimtempNR
+ 1];
432 //! Macro for correct tempering string
433 #define ESIMTEMP(e) enum_name(e, esimtempNR, esimtemp_names)
435 /*! \brief Free energy perturbation type
437 * efepNO, there are no evaluations at other states.
438 * efepYES, treated equivalently to efepSTATIC.
439 * efepSTATIC, then lambdas do not change during the simulation.
440 * efepSLOWGROWTH, then the states change monotonically
441 * throughout the simulation.
442 * efepEXPANDED, then expanded ensemble simulations are occuring.
453 //! String corresponding to FEP type.
454 extern const char* efep_names
[efepNR
+ 1];
455 //! Macro corresponding to FEP string.
456 #define EFEPTYPE(e) enum_name(e, efepNR, efep_names)
458 //! Free energy pertubation coupling types.
470 //! String for FEP coupling type
471 extern const char* efpt_names
[efptNR
+ 1];
472 //! Long names for FEP coupling type
473 extern const char* efpt_singular_names
[efptNR
+ 1];
475 /*! \brief What to print for free energy calculations
477 * Printing the energy to the free energy dhdl file.
478 * YES is an alias to TOTAL, and
479 * will be converted in readir, so we never have to account for it in code.
484 edHdLPrintEnergyTOTAL
,
485 edHdLPrintEnergyPOTENTIAL
,
489 //! String corresponding to printing of free energy
490 extern const char* edHdLPrintEnergy_names
[edHdLPrintEnergyNR
+ 1];
492 /*! \brief How the lambda weights are calculated
494 * elamstatsMETROPOLIS - using the metropolis criteria
495 * elamstatsBARKER - using the Barker critera for transition weights,
496 * also called unoptimized Bennett
497 * elamstatsMINVAR - using Barker + minimum variance for weights
498 * elamstatsWL - Wang-Landu (using visitation counts)
499 * elamstatsWWL - Weighted Wang-Landau (using optimized Gibbs
500 * weighted visitation counts)
512 //! String corresponding to lambda weights
513 extern const char* elamstats_names
[elamstatsNR
+ 1];
514 //! Macro telling us whether we use expanded ensemble
515 #define ELAMSTATS_EXPANDED(e) ((e) > elamstatsNO)
516 //! Macro telling us whether we use some kind of Wang-Landau
517 #define EWL(e) ((e) == elamstatsWL || (e) == elamstatsWWL)
519 /*! \brief How moves in lambda are calculated
521 * elmovemcMETROPOLIS - using the Metropolis criteria, and 50% up and down
522 * elmovemcBARKER - using the Barker criteria, and 50% up and down
523 * elmovemcGIBBS - computing the transition using the marginalized
524 * probabilities of the lambdas
525 * elmovemcMETGIBBS - computing the transition using the metropolized
526 * version of Gibbs (Monte Carlo Strategies in
527 * Scientific computing, Liu, p. 134)
538 //! String corresponding to lambda moves
539 extern const char* elmcmove_names
[elmcmoveNR
+ 1];
541 /*! \brief How we decide whether weights have reached equilibrium
543 * elmceqNO - never stop, weights keep going
544 * elmceqYES - fix the weights from the beginning; no movement
545 * elmceqWLDELTA - stop when the WL-delta falls below a certain level
546 * elmceqNUMATLAM - stop when we have a certain number of samples at
548 * elmceqSTEPS - stop when we've run a certain total number of steps
549 * elmceqSAMPLES - stop when we've run a certain total number of samples
550 * elmceqRATIO - stop when the ratio of samples (lowest to highest)
551 * is sufficiently large
564 //! String corresponding to equilibrium algorithm
565 extern const char* elmceq_names
[elmceqNR
+ 1];
567 /*! \brief separate_dhdl_file selection
569 * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
577 //! String corresponding to separate DHDL file selection
578 extern const char* separate_dhdl_file_names
[esepdhdlfileNR
+ 1];
579 //! Monster macro for DHDL file selection
580 #define SEPDHDLFILETYPE(e) enum_name(e, esepdhdlfileNR, separate_dhdl_file_names)
582 /*! \brief dhdl_derivatives selection \
584 * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
592 //! String for DHDL derivatives
593 extern const char* dhdl_derivatives_names
[edhdlderivativesNR
+ 1];
594 //! YAMM (Yet another monster macro)
595 #define DHDLDERIVATIVESTYPE(e) enum_name(e, edhdlderivativesNR, dhdl_derivatives_names)
597 /*! \brief Solvent model
599 * Distinguishes classical water types with 3 or 4 particles
608 //! String corresponding to solvent type
609 extern const char* esol_names
[esolNR
+ 1];
610 //! Macro lest we print the wrong solvent model string
611 #define ESOLTYPE(e) enum_name(e, esolNR, esol_names)
613 //! Dispersion correction.
623 //! String corresponding to dispersion corrections
624 extern const char* edispc_names
[edispcNR
+ 1];
625 //! Macro for dispcorr string
626 #define EDISPCORR(e) enum_name(e, edispcNR, edispc_names)
628 //! Center of mass motion removal algorithm.
634 ecmLINEAR_ACCELERATION_CORRECTION
,
637 //! String corresponding to COM removal
638 extern const char* ecm_names
[ecmNR
+ 1];
639 //! Macro for COM removal string
640 #define ECOM(e) enum_name(e, ecmNR, ecm_names)
642 //! Algorithm for simulated annealing.
650 //! String for simulated annealing
651 extern const char* eann_names
[eannNR
+ 1];
652 //! And macro for simulated annealing string
653 #define EANNEAL(e) enum_name(e, eannNR, eann_names)
664 //! String corresponding to wall type
665 extern const char* ewt_names
[ewtNR
+ 1];
666 //! Macro for wall type string
667 #define EWALLTYPE(e) enum_name(e, ewtNR, ewt_names)
669 //! Pulling algorithm.
680 //! String for pulling algorithm
681 extern const char* epull_names
[epullNR
+ 1];
682 //! Macro for pulling string
683 #define EPULLTYPE(e) enum_name(e, epullNR, epull_names)
685 //! Control of pull groups
698 //! String for pull groups
699 extern const char* epullg_names
[epullgNR
+ 1];
700 //! Macro for pull group string
701 #define EPULLGEOM(e) enum_name(e, epullgNR, epullg_names)
703 //! Enforced rotation groups.
720 //! Rotation group names
721 extern const char* erotg_names
[erotgNR
+ 1];
722 //! Macro for rot group names
723 #define EROTGEOM(e) enum_name(e, erotgNR, erotg_names)
724 //! String for rotation group origin names
725 extern const char* erotg_originnames
[erotgNR
+ 1];
726 //! Macro for rot group origin names
727 #define EROTORIGIN(e) enum_name(e, erotgOriginNR, erotg_originnames)
729 //! Rotation group fitting type
737 //! String corresponding to rotation group fitting
738 extern const char* erotg_fitnames
[erotgFitNR
+ 1];
739 //! Macro for rot group fit names
740 #define EROTFIT(e) enum_name(e, erotgFitNR, erotg_fitnames)
742 /*! \brief Direction along which ion/water swaps happen
744 * Part of "Computational Electrophysiology" (CompEL) setups
754 //! Names for swapping
755 extern const char* eSwapTypes_names
[eSwapTypesNR
+ 1];
756 //! Macro for swapping string
757 #define ESWAPTYPE(e) enum_name(e, eSwapTypesNR, eSwapTypes_names)
759 /*! \brief Swap group splitting type
761 * These are just the fixed groups we need for any setup. In t_swap's grp
762 * entry after that follows the variable number of swap groups.
771 //! String for swap group splitting
772 extern const char* eSwapFixedGrp_names
[eSwapFixedGrpNR
+ 1];
774 /*! \brief Neighborlist geometry type.
776 * Kernels will compute interactions between two particles,
777 * 3-center water, 4-center water or coarse-grained beads.
779 enum gmx_nblist_kernel_geometry
781 GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE
,
782 GMX_NBLIST_GEOMETRY_WATER3_PARTICLE
,
783 GMX_NBLIST_GEOMETRY_WATER3_WATER3
,
784 GMX_NBLIST_GEOMETRY_WATER4_PARTICLE
,
785 GMX_NBLIST_GEOMETRY_WATER4_WATER4
,
786 GMX_NBLIST_GEOMETRY_CG_CG
,
787 GMX_NBLIST_GEOMETRY_NR
789 //! String corresponding to nblist geometry names
790 extern const char* gmx_nblist_geometry_names
[GMX_NBLIST_GEOMETRY_NR
+ 1];
792 /*! \brief Types of electrostatics calculations
794 * Types of electrostatics calculations available inside nonbonded kernels.
795 * Note that these do NOT necessarily correspond to the user selections
796 * in the MDP file; many interactions for instance map to tabulated kernels.
798 enum gmx_nbkernel_elec
800 GMX_NBKERNEL_ELEC_NONE
,
801 GMX_NBKERNEL_ELEC_COULOMB
,
802 GMX_NBKERNEL_ELEC_REACTIONFIELD
,
803 GMX_NBKERNEL_ELEC_CUBICSPLINETABLE
,
804 GMX_NBKERNEL_ELEC_EWALD
,
807 //! String corresponding to electrostatics kernels
808 extern const char* gmx_nbkernel_elec_names
[GMX_NBKERNEL_ELEC_NR
+ 1];
810 /*! \brief Types of vdw calculations available
812 * Types of vdw calculations available inside nonbonded kernels.
813 * Note that these do NOT necessarily correspond to the user selections
814 * in the MDP file; many interactions for instance map to tabulated kernels.
816 enum gmx_nbkernel_vdw
818 GMX_NBKERNEL_VDW_NONE
,
819 GMX_NBKERNEL_VDW_LENNARDJONES
,
820 GMX_NBKERNEL_VDW_BUCKINGHAM
,
821 GMX_NBKERNEL_VDW_CUBICSPLINETABLE
,
822 GMX_NBKERNEL_VDW_LJEWALD
,
825 //! String corresponding to VdW kernels
826 extern const char* gmx_nbkernel_vdw_names
[GMX_NBKERNEL_VDW_NR
+ 1];
828 //! \brief Types of interactions inside the neighborlist
829 enum gmx_nblist_interaction_type
831 GMX_NBLIST_INTERACTION_STANDARD
,
832 GMX_NBLIST_INTERACTION_FREE_ENERGY
,
833 GMX_NBLIST_INTERACTION_NR
835 //! String corresponding to interactions in neighborlist code
836 extern const char* gmx_nblist_interaction_names
[GMX_NBLIST_INTERACTION_NR
+ 1];
838 #endif /* GMX_MDTYPES_MD_ENUMS_H */