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, 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.
67 eCompA
, eCompB
, eCompNR
70 /*! \brief The channels that define with their COM the compartment boundaries in CompEL setups.
72 * In principle one could also use modified setups with more than two channels.
75 eChan0
, eChan1
, eChanNR
78 /*! \brief Temperature coupling type
80 * yes is an alias for berendsen
83 etcNO
, etcBERENDSEN
, etcNOSEHOOVER
, etcYES
, etcANDERSEN
, etcANDERSENMASSIVE
, etcVRESCALE
, etcNR
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))
92 /*! \brief Pressure coupling types
94 * isotropic is an alias for berendsen
97 epcNO
, epcBERENDSEN
, epcPARRINELLORAHMAN
, epcISOTROPIC
, epcMTTK
, epcNR
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)
104 //! Flat-bottom posres geometries
106 efbposresZERO
, efbposresSPHERE
, efbposresCYLINDER
, efbposresX
, efbposresY
, efbposresZ
,
107 efbposresCYLINDERX
, efbposresCYLINDERY
, efbposresCYLINDERZ
, efbposresNR
110 //! Relative coordinate scaling type for position restraints.
112 erscNO
, erscALL
, erscCOM
, erscNR
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)
119 //! Trotter decomposition extended variable parts.
121 etrtNONE
, etrtNHC
, etrtBAROV
, etrtBARONHC
, etrtNHC2
, etrtBAROV2
, etrtBARONHC2
,
122 etrtVELOCITY1
, etrtVELOCITY2
, etrtPOSITION
, etrtSKIPALL
, etrtNR
125 //! Sequenced parts of the trotter decomposition.
127 ettTSEQ0
, ettTSEQ1
, ettTSEQ2
, ettTSEQ3
, ettTSEQ4
, ettTSEQMAX
130 //! Pressure coupling type
132 epctISOTROPIC
, epctSEMIISOTROPIC
, epctANISOTROPIC
,
133 epctSURFACETENSION
, epctNR
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)
140 //! \\brief Cutoff scheme
142 ecutsVERLET
, ecutsGROUP
, ecutsNR
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)
149 /*! \brief Coulomb / VdW interaction modifiers.
151 * grompp replaces eintmodPOTSHIFT_VERLET by eintmodPOTSHIFT or eintmodNONE.
152 * Exactcutoff is only used by Reaction-field-zero, and is not user-selectable.
155 eintmodPOTSHIFT_VERLET
, eintmodPOTSHIFT
, eintmodNONE
, eintmodPOTSWITCH
, eintmodEXACTCUTOFF
, eintmodFORCESWITCH
, eintmodNR
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)
162 /*! \brief Cut-off treatment for Coulomb
164 * eelNOTUSED1 used to be GB, but to enable generalized born with different
165 * forms of electrostatics (RF, switch, etc.) in the future it is now selected
166 * separately (through the implicit_solvent option).
169 eelCUT
, eelRF
, eelGRF
, eelPME
, eelEWALD
, eelP3M_AD
,
170 eelPOISSON
, eelSWITCH
, eelSHIFT
, eelUSER
, eelGB_NOTUSED
, eelRF_NEC_UNSUPPORTED
, eelENCADSHIFT
,
171 eelPMEUSER
, eelPMESWITCH
, eelPMEUSERSWITCH
, eelRF_ZERO
, eelNR
173 //! String corresponding to Coulomb treatment
174 extern const char *eel_names
[eelNR
+1];
175 //! Macro for correct string for Coulomb treatment
176 #define EELTYPE(e) enum_name(e, eelNR, eel_names)
180 eewg3D
, eewg3DC
, eewgNR
182 //! String corresponding to Ewald geometry
183 extern const char *eewg_names
[eewgNR
+1];
185 //! Macro telling us whether we use reaction field
186 #define EEL_RF(e) ((e) == eelRF || (e) == eelGRF || (e) == eelRF_NEC_UNSUPPORTED || (e) == eelRF_ZERO )
188 //! Macro telling us whether we use PME
189 #define EEL_PME(e) ((e) == eelPME || (e) == eelPMESWITCH || (e) == eelPMEUSER || (e) == eelPMEUSERSWITCH || (e) == eelP3M_AD)
190 //! Macro telling us whether we use PME or full Ewald
191 #define EEL_PME_EWALD(e) (EEL_PME(e) || (e) == eelEWALD)
192 //! Macro telling us whether we use full electrostatics of any sort
193 #define EEL_FULL(e) (EEL_PME_EWALD(e) || (e) == eelPOISSON)
194 //! Macro telling us whether we use user defined electrostatics
195 #define EEL_USER(e) ((e) == eelUSER || (e) == eelPMEUSER || (e) == (eelPMEUSERSWITCH))
197 //! Van der Waals interaction treatment
199 evdwCUT
, evdwSWITCH
, evdwSHIFT
, evdwUSER
, evdwENCADSHIFT
,
202 //! String corresponding to Van der Waals treatment
203 extern const char *evdw_names
[evdwNR
+1];
204 //! Macro for selecting correct string for VdW treatment
205 #define EVDWTYPE(e) enum_name(e, evdwNR, evdw_names)
207 //! Type of long-range VdW treatment of combination rules
209 eljpmeGEOM
, eljpmeLB
, eljpmeNR
211 //! String for LJPME combination rule treatment
212 extern const char *eljpme_names
[eljpmeNR
+1];
213 //! Macro for correct LJPME comb rule name
214 #define ELJPMECOMBNAMES(e) enum_name(e, eljpmeNR, eljpme_names)
216 //! Macro to tell us whether we use LJPME
217 #define EVDW_PME(e) ((e) == evdwPME)
219 //! Neighborsearching algorithm
221 ensGRID
, ensSIMPLE
, ensNR
223 //! String corresponding to neighborsearching
224 extern const char *ens_names
[ensNR
+1];
225 //! Macro for correct NS algorithm
226 #define ENS(e) enum_name(e, ensNR, ens_names)
228 /*! \brief Integrator algorithm
230 * eiSD2 has been removed, but we keep a renamed enum entry,
231 * so we can refuse to do MD with such .tpr files.
232 * eiVV is normal velocity verlet
233 * eiVVAK uses 1/2*(KE(t-dt/2)+KE(t+dt/2)) as the kinetic energy,
234 * and the half step kinetic energy for temperature control
237 eiMD
, eiSteep
, eiCG
, eiBD
, eiSD2_REMOVED
, eiNM
, eiLBFGS
, eiTPI
, eiTPIC
, eiSD1
, eiVV
, eiVVAK
, eiNR
239 //! Name of the integrator algorithm
240 extern const char *ei_names
[eiNR
+1];
241 //! Macro returning integrator string
242 #define EI(e) enum_name(e, eiNR, ei_names)
243 //! Do we use velocity Verlet
244 #define EI_VV(e) ((e) == eiVV || (e) == eiVVAK)
245 //! Do we use molecular dynamics
246 #define EI_MD(e) ((e) == eiMD || EI_VV(e))
247 //! Do we use stochastic dynamics
248 #define EI_SD(e) ((e) == eiSD1)
249 //! Do we use any stochastic integrator
250 #define EI_RANDOM(e) (EI_SD(e) || (e) == eiBD)
251 /*above integrators may not conserve momenta*/
252 //! Do we use any type of dynamics
253 #define EI_DYNAMICS(e) (EI_MD(e) || EI_RANDOM(e))
254 //! Or do we use minimization
255 #define EI_ENERGY_MINIMIZATION(e) ((e) == eiSteep || (e) == eiCG || (e) == eiLBFGS)
256 //! Do we apply test particle insertion
257 #define EI_TPI(e) ((e) == eiTPI || (e) == eiTPIC)
258 //! Do we deal with particle velocities
259 #define EI_STATE_VELOCITY(e) (EI_MD(e) || EI_SD(e))
261 //! Constraint algorithm
263 econtLINCS
, econtSHAKE
, econtNR
265 //! String corresponding to constraint algorithm
266 extern const char *econstr_names
[econtNR
+1];
267 //! Macro to select the correct string
268 #define ECONSTRTYPE(e) enum_name(e, econtNR, econstr_names)
270 //! Distance restraint refinement algorithm
272 edrNone
, edrSimple
, edrEnsemble
, edrNR
274 //! String corresponding to distance restraint algorithm
275 extern const char *edisre_names
[edrNR
+1];
276 //! Macro to select the right disre algorithm string
277 #define EDISRETYPE(e) enum_name(e, edrNR, edisre_names)
279 //! Distance restraints weighting type
281 edrwConservative
, edrwEqual
, edrwNR
283 //! String corresponding to distance restraint weighting
284 extern const char *edisreweighting_names
[edrwNR
+1];
285 //! Macro corresponding to dr weighting
286 #define EDISREWEIGHTING(e) enum_name(e, edrwNR, edisreweighting_names)
288 //! Combination rule algorithm.
290 eCOMB_NONE
, eCOMB_GEOMETRIC
, eCOMB_ARITHMETIC
, eCOMB_GEOM_SIG_EPS
, eCOMB_NR
292 //! String for combination rule algorithm
293 extern const char *ecomb_names
[eCOMB_NR
+1];
294 //! Macro to select the comb rule string
295 #define ECOMBNAME(e) enum_name(e, eCOMB_NR, ecomb_names)
297 //! Van der Waals potential.
299 eNBF_NONE
, eNBF_LJ
, eNBF_BHAM
, eNBF_NR
301 //! String corresponding to Van der Waals potential
302 extern const char *enbf_names
[eNBF_NR
+1];
303 //! Macro for correct VdW potential string
304 #define ENBFNAME(e) enum_name(e, eNBF_NR, enbf_names)
306 //! Simulated tempering methods.
308 esimtempGEOMETRIC
, esimtempEXPONENTIAL
, esimtempLINEAR
, esimtempNR
310 //! String corresponding to simulated tempering
311 extern const char *esimtemp_names
[esimtempNR
+1];
312 //! Macro for correct tempering string
313 #define ESIMTEMP(e) enum_name(e, esimtempNR, esimtemp_names)
315 /*! \brief Free energy perturbation type
317 * efepNO, there are no evaluations at other states.
318 * efepYES, treated equivalently to efepSTATIC.
319 * efepSTATIC, then lambdas do not change during the simulation.
320 * efepSLOWGROWTH, then the states change monotonically
321 * throughout the simulation.
322 * efepEXPANDED, then expanded ensemble simulations are occuring.
325 efepNO
, efepYES
, efepSTATIC
, efepSLOWGROWTH
, efepEXPANDED
, efepNR
327 //! String corresponding to FEP type.
328 extern const char *efep_names
[efepNR
+1];
329 //! Macro corresponding to FEP string.
330 #define EFEPTYPE(e) enum_name(e, efepNR, efep_names)
332 //! Free energy pertubation coupling types.
334 efptFEP
, efptMASS
, efptCOUL
, efptVDW
, efptBONDED
, efptRESTRAINT
, efptTEMPERATURE
, efptNR
336 //! String for FEP coupling type
337 extern const char *efpt_names
[efptNR
+1];
338 //! Long names for FEP coupling type
339 extern const char *efpt_singular_names
[efptNR
+1];
341 /*! \brief What to print for free energy calculations
343 * Printing the energy to the free energy dhdl file.
344 * YES is an alias to TOTAL, and
345 * will be converted in readir, so we never have to account for it in code.
348 edHdLPrintEnergyNO
, edHdLPrintEnergyTOTAL
, edHdLPrintEnergyPOTENTIAL
, edHdLPrintEnergyYES
, edHdLPrintEnergyNR
350 //! String corresponding to printing of free energy
351 extern const char *edHdLPrintEnergy_names
[edHdLPrintEnergyNR
+1];
353 /*! \brief How the lambda weights are calculated
355 * elamstatsMETROPOLIS - using the metropolis criteria
356 * elamstatsBARKER - using the Barker critera for transition weights,
357 * also called unoptimized Bennett
358 * elamstatsMINVAR - using Barker + minimum variance for weights
359 * elamstatsWL - Wang-Landu (using visitation counts)
360 * elamstatsWWL - Weighted Wang-Landau (using optimized Gibbs
361 * weighted visitation counts)
364 elamstatsNO
, elamstatsMETROPOLIS
, elamstatsBARKER
, elamstatsMINVAR
, elamstatsWL
, elamstatsWWL
, elamstatsNR
366 //! String corresponding to lambda weights
367 extern const char *elamstats_names
[elamstatsNR
+1];
368 //! Macro telling us whether we use expanded ensemble
369 #define ELAMSTATS_EXPANDED(e) ((e) > elamstatsNO)
370 //! Macro telling us whether we use some kind of Wang-Landau
371 #define EWL(e) ((e) == elamstatsWL || (e) == elamstatsWWL)
373 /*! \brief How moves in lambda are calculated
375 * elmovemcMETROPOLIS - using the Metropolis criteria, and 50% up and down
376 * elmovemcBARKER - using the Barker criteria, and 50% up and down
377 * elmovemcGIBBS - computing the transition using the marginalized
378 * probabilities of the lambdas
379 * elmovemcMETGIBBS - computing the transition using the metropolized
380 * version of Gibbs (Monte Carlo Strategies in
381 * Scientific computing, Liu, p. 134)
384 elmcmoveNO
, elmcmoveMETROPOLIS
, elmcmoveBARKER
, elmcmoveGIBBS
, elmcmoveMETGIBBS
, elmcmoveNR
386 //! String corresponding to lambda moves
387 extern const char *elmcmove_names
[elmcmoveNR
+1];
389 /*! \brief How we decide whether weights have reached equilibrium
391 * elmceqNO - never stop, weights keep going
392 * elmceqYES - fix the weights from the beginning; no movement
393 * elmceqWLDELTA - stop when the WL-delta falls below a certain level
394 * elmceqNUMATLAM - stop when we have a certain number of samples at
396 * elmceqSTEPS - stop when we've run a certain total number of steps
397 * elmceqSAMPLES - stop when we've run a certain total number of samples
398 * elmceqRATIO - stop when the ratio of samples (lowest to highest)
399 * is sufficiently large
402 elmceqNO
, elmceqYES
, elmceqWLDELTA
, elmceqNUMATLAM
, elmceqSTEPS
, elmceqSAMPLES
, elmceqRATIO
, elmceqNR
404 //! String corresponding to equilibrium algorithm
405 extern const char *elmceq_names
[elmceqNR
+1];
407 /*! \brief separate_dhdl_file selection
409 * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
413 esepdhdlfileYES
, esepdhdlfileNO
, esepdhdlfileNR
415 //! String corresponding to separate DHDL file selection
416 extern const char *separate_dhdl_file_names
[esepdhdlfileNR
+1];
417 //! Monster macro for DHDL file selection
418 #define SEPDHDLFILETYPE(e) enum_name(e, esepdhdlfileNR, separate_dhdl_file_names)
420 /*! \brief dhdl_derivatives selection \
422 * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
426 edhdlderivativesYES
, edhdlderivativesNO
, edhdlderivativesNR
428 //! String for DHDL derivatives
429 extern const char *dhdl_derivatives_names
[edhdlderivativesNR
+1];
430 //! YAMM (Yet another monster macro)
431 #define DHDLDERIVATIVESTYPE(e) enum_name(e, edhdlderivativesNR, dhdl_derivatives_names)
433 /*! \brief Solvent model
435 * Distinguishes classical water types with 3 or 4 particles
438 esolNO
, esolSPC
, esolTIP4P
, esolNR
440 //! String corresponding to solvent type
441 extern const char *esol_names
[esolNR
+1];
442 //! Macro lest we print the wrong solvent model string
443 #define ESOLTYPE(e) enum_name(e, esolNR, esol_names)
445 //! Dispersion correction.
447 edispcNO
, edispcEnerPres
, edispcEner
, edispcAllEnerPres
, edispcAllEner
, edispcNR
449 //! String corresponding to dispersion corrections
450 extern const char *edispc_names
[edispcNR
+1];
451 //! Macro for dispcorr string
452 #define EDISPCORR(e) enum_name(e, edispcNR, edispc_names)
454 //! Center of mass motion removal algorithm.
456 ecmLINEAR
, ecmANGULAR
, ecmNO
, ecmNR
458 //! String corresponding to COM removal
459 extern const char *ecm_names
[ecmNR
+1];
460 //! Macro for COM removal string
461 #define ECOM(e) enum_name(e, ecmNR, ecm_names)
463 //! Algorithm for simulated annealing.
465 eannNO
, eannSINGLE
, eannPERIODIC
, eannNR
467 //! String for simulated annealing
468 extern const char *eann_names
[eannNR
+1];
469 //! And macro for simulated annealing string
470 #define EANNEAL(e) enum_name(e, eannNR, eann_names)
472 //! Implicit solvent algorithms.
474 eisNO
, eisGBSA
, eisNR
476 //! String corresponding to implicit solvent.
477 extern const char *eis_names
[eisNR
+1];
478 //! Macro for implicit solvent string.
479 #define EIMPLICITSOL(e) enum_name(e, eisNR, eis_names)
481 //! Algorithms for calculating GB radii.
483 egbSTILL
, egbHCT
, egbOBC
, egbNR
485 //! String for GB algorithm name.
486 extern const char *egb_names
[egbNR
+1];
487 //! Macro for GB string.
488 #define EGBALGORITHM(e) enum_name(e, egbNR, egb_names)
490 //! Surface area algorithm for implicit solvent.
492 esaAPPROX
, esaNO
, esaSTILL
, esaNR
494 //! String corresponding to surface area algorithm.
495 extern const char *esa_names
[esaNR
+1];
496 //! brief Macro for SA algorithm string.
497 #define ESAALGORITHM(e) enum_name(e, esaNR, esa_names)
501 ewt93
, ewt104
, ewtTABLE
, ewt126
, ewtNR
503 //! String corresponding to wall type
504 extern const char *ewt_names
[ewtNR
+1];
505 //! Macro for wall type string
506 #define EWALLTYPE(e) enum_name(e, ewtNR, ewt_names)
508 //! Pulling algorithm.
510 epullUMBRELLA
, epullCONSTRAINT
, epullCONST_F
, epullFLATBOTTOM
, epullFLATBOTTOMHIGH
, epullEXTERNAL
, epullNR
512 //! String for pulling algorithm
513 extern const char *epull_names
[epullNR
+1];
514 //! Macro for pulling string
515 #define EPULLTYPE(e) enum_name(e, epullNR, epull_names)
517 //! Control of pull groups
519 epullgDIST
, epullgDIR
, epullgCYL
, epullgDIRPBC
, epullgDIRRELATIVE
, epullgANGLE
, epullgDIHEDRAL
, epullgANGLEAXIS
, epullgNR
521 //! String for pull groups
522 extern const char *epullg_names
[epullgNR
+1];
523 //! Macro for pull group string
524 #define EPULLGEOM(e) enum_name(e, epullgNR, epullg_names)
526 //! Enforced rotation groups.
528 erotgISO
, erotgISOPF
,
531 erotgRM2
, erotgRM2PF
,
532 erotgFLEX
, erotgFLEXT
,
533 erotgFLEX2
, erotgFLEX2T
,
536 //! Rotation group names
537 extern const char *erotg_names
[erotgNR
+1];
538 //! Macro for rot group names
539 #define EROTGEOM(e) enum_name(e, erotgNR, erotg_names)
540 //! String for rotation group origin names
541 extern const char *erotg_originnames
[erotgNR
+1];
542 //! Macro for rot group origin names
543 #define EROTORIGIN(e) enum_name(e, erotgOriginNR, erotg_originnames)
545 //! Rotation group fitting type
547 erotgFitRMSD
, erotgFitNORM
, erotgFitPOT
, erotgFitNR
549 //! String corresponding to rotation group fitting
550 extern const char *erotg_fitnames
[erotgFitNR
+1];
551 //! Macro for rot group fit names
552 #define EROTFIT(e) enum_name(e, erotgFitNR, erotg_fitnames)
554 /*! \brief Direction along which ion/water swaps happen
556 * Part of "Computational Electrophysiology" (CompEL) setups
559 eswapNO
, eswapX
, eswapY
, eswapZ
, eSwapTypesNR
561 //! Names for swapping
562 extern const char *eSwapTypes_names
[eSwapTypesNR
+1];
563 //! Macro for swapping string
564 #define ESWAPTYPE(e) enum_name(e, eSwapTypesNR, eSwapTypes_names)
566 /*! \brief Swap group splitting type
568 * These are just the fixed groups we need for any setup. In t_swap's grp
569 * entry after that follows the variable number of swap groups.
572 eGrpSplit0
, eGrpSplit1
, eGrpSolvent
, eSwapFixedGrpNR
574 //! String for swap group splitting
575 extern const char *eSwapFixedGrp_names
[eSwapFixedGrpNR
+1];
579 eQMmethodAM1
, eQMmethodPM3
, eQMmethodRHF
,
580 eQMmethodUHF
, eQMmethodDFT
, eQMmethodB3LYP
, eQMmethodMP2
, eQMmethodCASSCF
, eQMmethodB3LYPLAN
,
581 eQMmethodDIRECT
, eQMmethodNR
583 //! String corresponding to QMMM methods
584 extern const char *eQMmethod_names
[eQMmethodNR
+1];
585 //! Macro to pick QMMM method name
586 #define EQMMETHOD(e) enum_name(e, eQMmethodNR, eQMmethod_names)
588 //! QMMM basis function for QM part
590 eQMbasisSTO3G
, eQMbasisSTO3G2
, eQMbasis321G
,
591 eQMbasis321Gp
, eQMbasis321dGp
, eQMbasis621G
,
592 eQMbasis631G
, eQMbasis631Gp
, eQMbasis631dGp
,
593 eQMbasis6311G
, eQMbasisNR
595 //! Name for QMMM basis function
596 extern const char *eQMbasis_names
[eQMbasisNR
+1];
597 //! Macro to pick right basis function string
598 #define EQMBASIS(e) enum_name(e, eQMbasisNR, eQMbasis_names)
602 eQMMMschemenormal
, eQMMMschemeoniom
, eQMMMschemeNR
604 //! QMMMM scheme names
605 extern const char *eQMMMscheme_names
[eQMMMschemeNR
+1];
606 //! Macro to pick QMMMM scheme name
607 #define EQMMMSCHEME(e) enum_name(e, eQMMMschemeNR, eQMMMscheme_names)
609 /*! \brief Neighborlist geometry type.
611 * Kernels will compute interactions between two particles,
612 * 3-center water, 4-center water or coarse-grained beads.
614 enum gmx_nblist_kernel_geometry
616 GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE
,
617 GMX_NBLIST_GEOMETRY_WATER3_PARTICLE
,
618 GMX_NBLIST_GEOMETRY_WATER3_WATER3
,
619 GMX_NBLIST_GEOMETRY_WATER4_PARTICLE
,
620 GMX_NBLIST_GEOMETRY_WATER4_WATER4
,
621 GMX_NBLIST_GEOMETRY_CG_CG
,
622 GMX_NBLIST_GEOMETRY_NR
624 //! String corresponding to nblist geometry names
625 extern const char *gmx_nblist_geometry_names
[GMX_NBLIST_GEOMETRY_NR
+1];
627 /*! \brief Types of electrostatics calculations
629 * Types of electrostatics calculations available inside nonbonded kernels.
630 * Note that these do NOT necessarily correspond to the user selections
631 * in the MDP file; many interactions for instance map to tabulated kernels.
633 enum gmx_nbkernel_elec
635 GMX_NBKERNEL_ELEC_NONE
,
636 GMX_NBKERNEL_ELEC_COULOMB
,
637 GMX_NBKERNEL_ELEC_REACTIONFIELD
,
638 GMX_NBKERNEL_ELEC_CUBICSPLINETABLE
,
639 GMX_NBKERNEL_ELEC_GENERALIZEDBORN
,
640 GMX_NBKERNEL_ELEC_EWALD
,
643 //! String corresponding to electrostatics kernels
644 extern const char *gmx_nbkernel_elec_names
[GMX_NBKERNEL_ELEC_NR
+1];
646 /*! \brief Types of vdw calculations available
648 * Types of vdw calculations available inside nonbonded kernels.
649 * Note that these do NOT necessarily correspond to the user selections
650 * in the MDP file; many interactions for instance map to tabulated kernels.
652 enum gmx_nbkernel_vdw
654 GMX_NBKERNEL_VDW_NONE
,
655 GMX_NBKERNEL_VDW_LENNARDJONES
,
656 GMX_NBKERNEL_VDW_BUCKINGHAM
,
657 GMX_NBKERNEL_VDW_CUBICSPLINETABLE
,
658 GMX_NBKERNEL_VDW_LJEWALD
,
661 //! String corresponding to VdW kernels
662 extern const char *gmx_nbkernel_vdw_names
[GMX_NBKERNEL_VDW_NR
+1];
664 //! \brief Types of interactions inside the neighborlist
665 enum gmx_nblist_interaction_type
667 GMX_NBLIST_INTERACTION_STANDARD
,
668 GMX_NBLIST_INTERACTION_FREE_ENERGY
,
669 GMX_NBLIST_INTERACTION_NR
671 //! String corresponding to interactions in neighborlist code
672 extern const char *gmx_nblist_interaction_names
[GMX_NBLIST_INTERACTION_NR
+1];
674 #endif /* GMX_MDTYPES_MD_ENUMS_H */