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,2018,2019, 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.
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.
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.
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.
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.
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.
39 * Declares enumerated types used throughout the code.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
43 * \ingroup module_mdtypes
45 #ifndef GMX_MDTYPES_MD_ENUMS_H
46 #define GMX_MDTYPES_MD_ENUMS_H
48 #include "gromacs/utility/basedefinitions.h"
50 /*! \brief Return a string from a list of strings
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
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"
60 const char* enum_name(int index, int max_index, const char* names[]);
62 //! Boolean strings no or yes
63 extern const char* yesno_names[BOOL_NR + 1];
65 //! \brief The two compartments for CompEL setups.
73 /*! \brief The channels that define with their COM the compartment boundaries in CompEL setups.
75 * In principle one could also use modified setups with more than two channels.
84 /*! \brief Temperature coupling type
86 * yes is an alias for berendsen
99 //! Strings corresponding to temperatyre coupling types
100 extern const char* etcoupl_names[etcNR + 1];
101 //! Macro for selecting t coupling string
102 #define ETCOUPLTYPE(e) enum_name(e, etcNR, etcoupl_names)
103 //! Return whether this is andersen coupling
104 #define ETC_ANDERSEN(e) (((e) == etcANDERSENMASSIVE) || ((e) == etcANDERSEN))
106 /*! \brief Pressure coupling types
108 * isotropic is an alias for berendsen
119 //! String corresponding to pressure coupling algorithm
120 extern const char* epcoupl_names[epcNR + 1];
121 //! Macro to return the correct pcoupling string
122 #define EPCOUPLTYPE(e) enum_name(e, epcNR, epcoupl_names)
124 //! Flat-bottom posres geometries
139 //! Relative coordinate scaling type for position restraints.
147 //! String corresponding to relativ coordinate scaling.
148 extern const char* erefscaling_names[erscNR + 1];
149 //! Macro to select correct coordinate scaling string.
150 #define EREFSCALINGTYPE(e) enum_name(e, erscNR, erefscaling_names)
152 //! Trotter decomposition extended variable parts.
169 //! Sequenced parts of the trotter decomposition.
180 //! Pressure coupling type
189 //! String corresponding to pressure coupling type
190 extern const char* epcoupltype_names[epctNR + 1];
191 //! Macro to select the right string for pcoupl type
192 #define EPCOUPLTYPETYPE(e) enum_name(e, epctNR, epcoupltype_names)
194 //! \\brief Cutoff scheme
201 //! String corresponding to cutoff scheme
202 extern const char* ecutscheme_names[ecutsNR + 1];
203 //! Macro to select the right string for cutoff scheme
204 #define ECUTSCHEME(e) enum_name(e, ecutsNR, ecutscheme_names)
206 /*! \brief Coulomb / VdW interaction modifiers.
208 * grompp replaces eintmodPOTSHIFT_VERLET_UNSUPPORTED by eintmodPOTSHIFT.
209 * Exactcutoff is only used by Reaction-field-zero, and is not user-selectable.
213 eintmodPOTSHIFT_VERLET_UNSUPPORTED,
221 //! String corresponding to interaction modifiers
222 extern const char* eintmod_names[eintmodNR + 1];
223 //! Macro to select the correct string for modifiers
224 #define INTMODIFIER(e) enum_name(e, eintmodNR, eintmod_names)
226 /*! \brief Cut-off treatment for Coulomb */
240 eelRF_NEC_UNSUPPORTED,
248 //! String corresponding to Coulomb treatment
249 extern const char* eel_names[eelNR + 1];
250 //! Macro for correct string for Coulomb treatment
251 #define EELTYPE(e) enum_name(e, eelNR, eel_names)
260 //! String corresponding to Ewald geometry
261 extern const char* eewg_names[eewgNR + 1];
263 //! Macro telling us whether we use reaction field
265 ((e) == eelRF || (e) == eelGRF_NOTUSED || (e) == eelRF_NEC_UNSUPPORTED || (e) == eelRF_ZERO)
267 //! Macro telling us whether we use PME
269 ((e) == eelPME || (e) == eelPMESWITCH || (e) == eelPMEUSER || (e) == eelPMEUSERSWITCH || (e) == eelP3M_AD)
270 //! Macro telling us whether we use PME or full Ewald
271 #define EEL_PME_EWALD(e) (EEL_PME(e) || (e) == eelEWALD)
272 //! Macro telling us whether we use full electrostatics of any sort
273 #define EEL_FULL(e) (EEL_PME_EWALD(e) || (e) == eelPOISSON)
274 //! Macro telling us whether we use user defined electrostatics
275 #define EEL_USER(e) ((e) == eelUSER || (e) == eelPMEUSER || (e) == (eelPMEUSERSWITCH))
277 //! Van der Waals interaction treatment
288 //! String corresponding to Van der Waals treatment
289 extern const char* evdw_names[evdwNR + 1];
290 //! Macro for selecting correct string for VdW treatment
291 #define EVDWTYPE(e) enum_name(e, evdwNR, evdw_names)
293 //! Type of long-range VdW treatment of combination rules
300 //! String for LJPME combination rule treatment
301 extern const char* eljpme_names[eljpmeNR + 1];
302 //! Macro for correct LJPME comb rule name
303 #define ELJPMECOMBNAMES(e) enum_name(e, eljpmeNR, eljpme_names)
305 //! Macro to tell us whether we use LJPME
306 #define EVDW_PME(e) ((e) == evdwPME)
308 /*! \brief Integrator algorithm
310 * eiSD2 has been removed, but we keep a renamed enum entry,
311 * so we can refuse to do MD with such .tpr files.
312 * eiVV is normal velocity verlet
313 * eiVVAK uses 1/2*(KE(t-dt/2)+KE(t+dt/2)) as the kinetic energy,
314 * and the half step kinetic energy for temperature control
333 //! Name of the integrator algorithm
334 extern const char* ei_names[eiNR + 1];
335 //! Macro returning integrator string
336 #define EI(e) enum_name(e, eiNR, ei_names)
337 //! Do we use MiMiC QM/MM?
338 #define EI_MIMIC(e) ((e) == eiMimic)
339 //! Do we use velocity Verlet
340 #define EI_VV(e) ((e) == eiVV || (e) == eiVVAK)
341 //! Do we use molecular dynamics
342 #define EI_MD(e) ((e) == eiMD || EI_VV(e) || EI_MIMIC(e))
343 //! Do we use stochastic dynamics
344 #define EI_SD(e) ((e) == eiSD1)
345 //! Do we use any stochastic integrator
346 #define EI_RANDOM(e) (EI_SD(e) || (e) == eiBD)
347 /*above integrators may not conserve momenta*/
348 //! Do we use any type of dynamics
349 #define EI_DYNAMICS(e) (EI_MD(e) || EI_RANDOM(e))
350 //! Or do we use minimization
351 #define EI_ENERGY_MINIMIZATION(e) ((e) == eiSteep || (e) == eiCG || (e) == eiLBFGS)
352 //! Do we apply test particle insertion
353 #define EI_TPI(e) ((e) == eiTPI || (e) == eiTPIC)
354 //! Do we deal with particle velocities
355 #define EI_STATE_VELOCITY(e) (EI_MD(e) || EI_SD(e))
357 //! Constraint algorithm
364 //! String corresponding to constraint algorithm
365 extern const char* econstr_names[econtNR + 1];
366 //! Macro to select the correct string
367 #define ECONSTRTYPE(e) enum_name(e, econtNR, econstr_names)
369 //! Distance restraint refinement algorithm
377 //! String corresponding to distance restraint algorithm
378 extern const char* edisre_names[edrNR + 1];
379 //! Macro to select the right disre algorithm string
380 #define EDISRETYPE(e) enum_name(e, edrNR, edisre_names)
382 //! Distance restraints weighting type
389 //! String corresponding to distance restraint weighting
390 extern const char* edisreweighting_names[edrwNR + 1];
391 //! Macro corresponding to dr weighting
392 #define EDISREWEIGHTING(e) enum_name(e, edrwNR, edisreweighting_names)
394 //! Combination rule algorithm.
403 //! String for combination rule algorithm
404 extern const char* ecomb_names[eCOMB_NR + 1];
405 //! Macro to select the comb rule string
406 #define ECOMBNAME(e) enum_name(e, eCOMB_NR, ecomb_names)
408 //! Van der Waals potential.
416 //! String corresponding to Van der Waals potential
417 extern const char* enbf_names[eNBF_NR + 1];
418 //! Macro for correct VdW potential string
419 #define ENBFNAME(e) enum_name(e, eNBF_NR, enbf_names)
421 //! Simulated tempering methods.
429 //! String corresponding to simulated tempering
430 extern const char* esimtemp_names[esimtempNR + 1];
431 //! Macro for correct tempering string
432 #define ESIMTEMP(e) enum_name(e, esimtempNR, esimtemp_names)
434 /*! \brief Free energy perturbation type
436 * efepNO, there are no evaluations at other states.
437 * efepYES, treated equivalently to efepSTATIC.
438 * efepSTATIC, then lambdas do not change during the simulation.
439 * efepSLOWGROWTH, then the states change monotonically
440 * throughout the simulation.
441 * efepEXPANDED, then expanded ensemble simulations are occuring.
452 //! String corresponding to FEP type.
453 extern const char* efep_names[efepNR + 1];
454 //! Macro corresponding to FEP string.
455 #define EFEPTYPE(e) enum_name(e, efepNR, efep_names)
457 //! Free energy pertubation coupling types.
469 //! String for FEP coupling type
470 extern const char* efpt_names[efptNR + 1];
471 //! Long names for FEP coupling type
472 extern const char* efpt_singular_names[efptNR + 1];
474 /*! \brief What to print for free energy calculations
476 * Printing the energy to the free energy dhdl file.
477 * YES is an alias to TOTAL, and
478 * will be converted in readir, so we never have to account for it in code.
483 edHdLPrintEnergyTOTAL,
484 edHdLPrintEnergyPOTENTIAL,
488 //! String corresponding to printing of free energy
489 extern const char* edHdLPrintEnergy_names[edHdLPrintEnergyNR + 1];
491 /*! \brief How the lambda weights are calculated
493 * elamstatsMETROPOLIS - using the metropolis criteria
494 * elamstatsBARKER - using the Barker critera for transition weights,
495 * also called unoptimized Bennett
496 * elamstatsMINVAR - using Barker + minimum variance for weights
497 * elamstatsWL - Wang-Landu (using visitation counts)
498 * elamstatsWWL - Weighted Wang-Landau (using optimized Gibbs
499 * weighted visitation counts)
511 //! String corresponding to lambda weights
512 extern const char* elamstats_names[elamstatsNR + 1];
513 //! Macro telling us whether we use expanded ensemble
514 #define ELAMSTATS_EXPANDED(e) ((e) > elamstatsNO)
515 //! Macro telling us whether we use some kind of Wang-Landau
516 #define EWL(e) ((e) == elamstatsWL || (e) == elamstatsWWL)
518 /*! \brief How moves in lambda are calculated
520 * elmovemcMETROPOLIS - using the Metropolis criteria, and 50% up and down
521 * elmovemcBARKER - using the Barker criteria, and 50% up and down
522 * elmovemcGIBBS - computing the transition using the marginalized
523 * probabilities of the lambdas
524 * elmovemcMETGIBBS - computing the transition using the metropolized
525 * version of Gibbs (Monte Carlo Strategies in
526 * Scientific computing, Liu, p. 134)
537 //! String corresponding to lambda moves
538 extern const char* elmcmove_names[elmcmoveNR + 1];
540 /*! \brief How we decide whether weights have reached equilibrium
542 * elmceqNO - never stop, weights keep going
543 * elmceqYES - fix the weights from the beginning; no movement
544 * elmceqWLDELTA - stop when the WL-delta falls below a certain level
545 * elmceqNUMATLAM - stop when we have a certain number of samples at
547 * elmceqSTEPS - stop when we've run a certain total number of steps
548 * elmceqSAMPLES - stop when we've run a certain total number of samples
549 * elmceqRATIO - stop when the ratio of samples (lowest to highest)
550 * is sufficiently large
563 //! String corresponding to equilibrium algorithm
564 extern const char* elmceq_names[elmceqNR + 1];
566 /*! \brief separate_dhdl_file selection
568 * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
576 //! String corresponding to separate DHDL file selection
577 extern const char* separate_dhdl_file_names[esepdhdlfileNR + 1];
578 //! Monster macro for DHDL file selection
579 #define SEPDHDLFILETYPE(e) enum_name(e, esepdhdlfileNR, separate_dhdl_file_names)
581 /*! \brief dhdl_derivatives selection \
583 * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
591 //! String for DHDL derivatives
592 extern const char* dhdl_derivatives_names[edhdlderivativesNR + 1];
593 //! YAMM (Yet another monster macro)
594 #define DHDLDERIVATIVESTYPE(e) enum_name(e, edhdlderivativesNR, dhdl_derivatives_names)
596 /*! \brief Solvent model
598 * Distinguishes classical water types with 3 or 4 particles
607 //! String corresponding to solvent type
608 extern const char* esol_names[esolNR + 1];
609 //! Macro lest we print the wrong solvent model string
610 #define ESOLTYPE(e) enum_name(e, esolNR, esol_names)
612 //! Dispersion correction.
622 //! String corresponding to dispersion corrections
623 extern const char* edispc_names[edispcNR + 1];
624 //! Macro for dispcorr string
625 #define EDISPCORR(e) enum_name(e, edispcNR, edispc_names)
627 //! Center of mass motion removal algorithm.
633 ecmLINEAR_ACCELERATION_CORRECTION,
636 //! String corresponding to COM removal
637 extern const char* ecm_names[ecmNR + 1];
638 //! Macro for COM removal string
639 #define ECOM(e) enum_name(e, ecmNR, ecm_names)
641 //! Algorithm for simulated annealing.
649 //! String for simulated annealing
650 extern const char* eann_names[eannNR + 1];
651 //! And macro for simulated annealing string
652 #define EANNEAL(e) enum_name(e, eannNR, eann_names)
663 //! String corresponding to wall type
664 extern const char* ewt_names[ewtNR + 1];
665 //! Macro for wall type string
666 #define EWALLTYPE(e) enum_name(e, ewtNR, ewt_names)
668 //! Pulling algorithm.
679 //! String for pulling algorithm
680 extern const char* epull_names[epullNR + 1];
681 //! Macro for pulling string
682 #define EPULLTYPE(e) enum_name(e, epullNR, epull_names)
684 //! Control of pull groups
697 //! String for pull groups
698 extern const char* epullg_names[epullgNR + 1];
699 //! Macro for pull group string
700 #define EPULLGEOM(e) enum_name(e, epullgNR, epullg_names)
702 //! Enforced rotation groups.
719 //! Rotation group names
720 extern const char* erotg_names[erotgNR + 1];
721 //! Macro for rot group names
722 #define EROTGEOM(e) enum_name(e, erotgNR, erotg_names)
723 //! String for rotation group origin names
724 extern const char* erotg_originnames[erotgNR + 1];
725 //! Macro for rot group origin names
726 #define EROTORIGIN(e) enum_name(e, erotgOriginNR, erotg_originnames)
728 //! Rotation group fitting type
736 //! String corresponding to rotation group fitting
737 extern const char* erotg_fitnames[erotgFitNR + 1];
738 //! Macro for rot group fit names
739 #define EROTFIT(e) enum_name(e, erotgFitNR, erotg_fitnames)
741 /*! \brief Direction along which ion/water swaps happen
743 * Part of "Computational Electrophysiology" (CompEL) setups
753 //! Names for swapping
754 extern const char* eSwapTypes_names[eSwapTypesNR + 1];
755 //! Macro for swapping string
756 #define ESWAPTYPE(e) enum_name(e, eSwapTypesNR, eSwapTypes_names)
758 /*! \brief Swap group splitting type
760 * These are just the fixed groups we need for any setup. In t_swap's grp
761 * entry after that follows the variable number of swap groups.
770 //! String for swap group splitting
771 extern const char* eSwapFixedGrp_names[eSwapFixedGrpNR + 1];
788 //! String corresponding to QMMM methods
789 extern const char* eQMmethod_names[eQMmethodNR + 1];
790 //! Macro to pick QMMM method name
791 #define EQMMETHOD(e) enum_name(e, eQMmethodNR, eQMmethod_names)
793 //! QMMM basis function for QM part
808 //! Name for QMMM basis function
809 extern const char* eQMbasis_names[eQMbasisNR + 1];
810 //! Macro to pick right basis function string
811 #define EQMBASIS(e) enum_name(e, eQMbasisNR, eQMbasis_names)
820 //! QMMMM scheme names
821 extern const char* eQMMMscheme_names[eQMMMschemeNR + 1];
822 //! Macro to pick QMMMM scheme name
823 #define EQMMMSCHEME(e) enum_name(e, eQMMMschemeNR, eQMMMscheme_names)
825 /*! \brief Neighborlist geometry type.
827 * Kernels will compute interactions between two particles,
828 * 3-center water, 4-center water or coarse-grained beads.
830 enum gmx_nblist_kernel_geometry
832 GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE,
833 GMX_NBLIST_GEOMETRY_WATER3_PARTICLE,
834 GMX_NBLIST_GEOMETRY_WATER3_WATER3,
835 GMX_NBLIST_GEOMETRY_WATER4_PARTICLE,
836 GMX_NBLIST_GEOMETRY_WATER4_WATER4,
837 GMX_NBLIST_GEOMETRY_CG_CG,
838 GMX_NBLIST_GEOMETRY_NR
840 //! String corresponding to nblist geometry names
841 extern const char* gmx_nblist_geometry_names[GMX_NBLIST_GEOMETRY_NR + 1];
843 /*! \brief Types of electrostatics calculations
845 * Types of electrostatics calculations available inside nonbonded kernels.
846 * Note that these do NOT necessarily correspond to the user selections
847 * in the MDP file; many interactions for instance map to tabulated kernels.
849 enum gmx_nbkernel_elec
851 GMX_NBKERNEL_ELEC_NONE,
852 GMX_NBKERNEL_ELEC_COULOMB,
853 GMX_NBKERNEL_ELEC_REACTIONFIELD,
854 GMX_NBKERNEL_ELEC_CUBICSPLINETABLE,
855 GMX_NBKERNEL_ELEC_EWALD,
858 //! String corresponding to electrostatics kernels
859 extern const char* gmx_nbkernel_elec_names[GMX_NBKERNEL_ELEC_NR + 1];
861 /*! \brief Types of vdw calculations available
863 * Types of vdw calculations available inside nonbonded kernels.
864 * Note that these do NOT necessarily correspond to the user selections
865 * in the MDP file; many interactions for instance map to tabulated kernels.
867 enum gmx_nbkernel_vdw
869 GMX_NBKERNEL_VDW_NONE,
870 GMX_NBKERNEL_VDW_LENNARDJONES,
871 GMX_NBKERNEL_VDW_BUCKINGHAM,
872 GMX_NBKERNEL_VDW_CUBICSPLINETABLE,
873 GMX_NBKERNEL_VDW_LJEWALD,
876 //! String corresponding to VdW kernels
877 extern const char* gmx_nbkernel_vdw_names[GMX_NBKERNEL_VDW_NR + 1];
879 //! \brief Types of interactions inside the neighborlist
880 enum gmx_nblist_interaction_type
882 GMX_NBLIST_INTERACTION_STANDARD,
883 GMX_NBLIST_INTERACTION_FREE_ENERGY,
884 GMX_NBLIST_INTERACTION_NR
886 //! String corresponding to interactions in neighborlist code
887 extern const char* gmx_nblist_interaction_names[GMX_NBLIST_INTERACTION_NR + 1];
890 //! \brief QM/MM mode
891 enum struct GmxQmmmMode
896 #endif /* GMX_MDTYPES_MD_ENUMS_H */