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,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.
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* const 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
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
93 enum class TemperatureCoupling : int
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))
111 /*! \brief Pressure coupling types
113 * isotropic is an alias for berendsen
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
119 enum class PressureCoupling : int
130 //! Return names of pressure coupling schemes
131 const char* enumValueToString(PressureCoupling enumValue);
133 //! Flat-bottom posres geometries
148 //! Relative coordinate scaling type for position restraints.
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)
161 //! Trotter decomposition extended variable parts.
178 //! Sequenced parts of the trotter decomposition.
189 //! Pressure coupling type
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)
203 //! \\brief Cutoff scheme
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)
215 /*! \brief Coulomb / VdW interaction modifiers.
217 * grompp replaces eintmodPOTSHIFT_VERLET_UNSUPPORTED by eintmodPOTSHIFT.
218 * Exactcutoff is only used by Reaction-field-zero, and is not user-selectable.
222 eintmodPOTSHIFT_VERLET_UNSUPPORTED,
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)
235 /*! \brief Cut-off treatment for Coulomb */
249 eelRF_NEC_UNSUPPORTED,
250 eelENCADSHIFT_NOTUSED,
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)
269 //! String corresponding to Ewald geometry
270 extern const char* eewg_names[eewgNR + 1];
272 //! Macro telling us whether we use reaction field
274 ((e) == eelRF || (e) == eelGRF_NOTUSED || (e) == eelRF_NEC_UNSUPPORTED || (e) == eelRF_ZERO)
276 //! Macro telling us whether we use PME
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))
286 //! Van der Waals interaction treatment
293 evdwENCADSHIFT_UNUSED,
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)
302 //! Type of long-range VdW treatment of combination rules
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)
314 //! Macro to tell us whether we use LJPME
315 #define EVDW_PME(e) ((e) == evdwPME)
317 /*! \brief Integrator algorithm
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
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))
366 //! Constraint algorithm
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)
378 //! Distance restraint refinement algorithm
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)
391 //! Distance restraints weighting type
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)
403 //! Combination rule algorithm.
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)
417 //! Van der Waals potential.
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)
430 //! Simulated tempering methods.
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)
443 /*! \brief Free energy perturbation type
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.
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)
466 //! Free energy pertubation coupling types.
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];
483 /*! \brief What to print for free energy calculations
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.
492 edHdLPrintEnergyTOTAL,
493 edHdLPrintEnergyPOTENTIAL,
497 //! String corresponding to printing of free energy
498 extern const char* edHdLPrintEnergy_names[edHdLPrintEnergyNR + 1];
500 /*! \brief How the lambda weights are calculated
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)
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)
527 /*! \brief How moves in lambda are calculated
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)
546 //! String corresponding to lambda moves
547 extern const char* elmcmove_names[elmcmoveNR + 1];
549 /*! \brief How we decide whether weights have reached equilibrium
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
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
572 //! String corresponding to equilibrium algorithm
573 extern const char* elmceq_names[elmceqNR + 1];
575 /*! \brief separate_dhdl_file selection
577 * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
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)
590 /*! \brief dhdl_derivatives selection \
592 * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
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)
605 /*! \brief Solvent model
607 * Distinguishes classical water types with 3 or 4 particles
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)
621 //! Dispersion correction.
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)
636 //! Center of mass motion removal algorithm.
642 ecmLINEAR_ACCELERATION_CORRECTION,
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)
650 //! Algorithm for simulated annealing.
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)
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)
677 //! Pulling algorithm.
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)
693 //! Control of pull groups
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)
711 //! Enforced rotation groups.
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)
737 //! Rotation group fitting type
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)
750 /*! \brief Direction along which ion/water swaps happen
752 * Part of "Computational Electrophysiology" (CompEL) setups
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)
767 /*! \brief Swap group splitting type
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.
779 //! String for swap group splitting
780 extern const char* eSwapFixedGrp_names[eSwapFixedGrpNR + 1];
782 /*! \brief Types of electrostatics calculations
784 * Types of electrostatics calculations available inside nonbonded kernels.
785 * Note that these do NOT necessarily correspond to the user selections
786 * in the MDP file; many interactions for instance map to tabulated kernels.
788 enum gmx_nbkernel_elec
790 GMX_NBKERNEL_ELEC_NONE,
791 GMX_NBKERNEL_ELEC_COULOMB,
792 GMX_NBKERNEL_ELEC_REACTIONFIELD,
793 GMX_NBKERNEL_ELEC_CUBICSPLINETABLE,
794 GMX_NBKERNEL_ELEC_EWALD,
797 //! String corresponding to electrostatics kernels
798 extern const char* gmx_nbkernel_elec_names[GMX_NBKERNEL_ELEC_NR + 1];
800 /*! \brief Types of vdw calculations available
802 * Types of vdw 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.
806 enum gmx_nbkernel_vdw
808 GMX_NBKERNEL_VDW_NONE,
809 GMX_NBKERNEL_VDW_LENNARDJONES,
810 GMX_NBKERNEL_VDW_BUCKINGHAM,
811 GMX_NBKERNEL_VDW_CUBICSPLINETABLE,
812 GMX_NBKERNEL_VDW_LJEWALD,
815 //! String corresponding to VdW kernels
816 extern const char* gmx_nbkernel_vdw_names[GMX_NBKERNEL_VDW_NR + 1];
818 #endif /* GMX_MDTYPES_MD_ENUMS_H */