Split simulationWork.useGpuBufferOps into separate x and f flags
[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 /*! \brief Return a string from a list of strings
50  *
51  * If index if within 0 .. max_index-1 returns the corresponding string
52  * or "no name defined" otherwise, in other words this is a range-check that does
53  * not crash.
54  * \param[in] index     The index in the array
55  * \param[in] max_index The length of the array
56  * \param[in] names     The array
57  * \return the correct string or "no name defined"
58  */
59 const char* enum_name(int index, int max_index, const char* const names[]);
60
61 /*! \brief Enum for setting answer to yes or no
62  */
63 enum class Boolean : int
64 {
65     No,
66     Yes,
67     Count,
68     Default = No
69 };
70
71 //! Return name of boolean selection.
72 const char* enumValueToString(Boolean enumValue);
73 //! Return name of boolean selection for actual bool.
74 const char* booleanValueToString(bool value);
75
76 //! \brief The two compartments for CompEL setups.
77 enum class Compartment : int
78 {
79     A,
80     B,
81     Count
82 };
83
84 /*! \brief The channels that define with their COM the compartment boundaries in CompEL setups.
85  *
86  * In principle one could also use modified setups with more than two channels.
87  */
88 enum class Channel : int
89 {
90     Zero,
91     One,
92     Count
93 };
94
95 /*! \brief Temperature coupling type
96  *
97  * yes is an alias for berendsen
98  *
99  * Note: Keep `Count` as the second-to-last entry, and `Default` as the last entry -
100  *       this is needed to keep EnumerationWrapper, EnumerationArray and (de)serialization
101  *       working.
102  */
103 enum class TemperatureCoupling : int
104 {
105     No,
106     Berendsen,
107     NoseHoover,
108     Yes,
109     Andersen,
110     AndersenMassive,
111     VRescale,
112     Count,
113     Default = No
114 };
115 //! Return names of temperature coupling schemes
116 const char* enumValueToString(TemperatureCoupling enumValue);
117 //! Return whether this is andersen coupling
118 #define ETC_ANDERSEN(e) \
119     (((e) == TemperatureCoupling::AndersenMassive) || ((e) == TemperatureCoupling::Andersen))
120
121 /*! \brief Pressure coupling types
122  *
123  * isotropic is an alias for berendsen
124  *
125  * Note: Keep `Count` as the second-to-last entry, and `Default` as the last entry -
126  *       this is needed to keep EnumerationWrapper, EnumerationArray and (de)serialization
127  *       working.
128  */
129 enum class PressureCoupling : int
130 {
131     No,
132     Berendsen,
133     ParrinelloRahman,
134     Isotropic,
135     Mttk,
136     CRescale,
137     Count,
138     Default = No
139 };
140 //! Return names of pressure coupling schemes
141 const char* enumValueToString(PressureCoupling enumValue);
142
143 //! Flat-bottom posres geometries
144 enum
145 {
146     efbposresZERO,
147     efbposresSPHERE,
148     efbposresCYLINDER,
149     efbposresX,
150     efbposresY,
151     efbposresZ,
152     efbposresCYLINDERX,
153     efbposresCYLINDERY,
154     efbposresCYLINDERZ,
155     efbposresNR
156 };
157
158 //! Relative coordinate scaling type for position restraints.
159 enum class RefCoordScaling : int
160 {
161     No,
162     All,
163     Com,
164     Count,
165     Default = No
166 };
167
168 //! String corresponding to relative coordinate scaling.
169 const char* enumValueToString(RefCoordScaling enumValue);
170
171 //! Trotter decomposition extended variable parts.
172 enum
173 {
174     etrtNONE,
175     etrtNHC,
176     etrtBAROV,
177     etrtBARONHC,
178     etrtNHC2,
179     etrtBAROV2,
180     etrtBARONHC2,
181     etrtVELOCITY1,
182     etrtVELOCITY2,
183     etrtPOSITION,
184     etrtSKIPALL,
185     etrtNR
186 };
187
188 //! Sequenced parts of the trotter decomposition.
189 enum class TrotterSequence : int
190 {
191     Zero,
192     One,
193     Two,
194     Three,
195     Four,
196     Count
197 };
198
199 //! Pressure coupling type
200 enum class PressureCouplingType : int
201 {
202     Isotropic,
203     SemiIsotropic,
204     Anisotropic,
205     SurfaceTension,
206     Count,
207     Default = Isotropic
208 };
209 //! String corresponding to pressure coupling type
210 const char* enumValueToString(PressureCouplingType enumValue);
211
212 //! \\brief Cutoff scheme
213 enum class CutoffScheme : int
214 {
215     Verlet,
216     Group,
217     Count,
218     Default = Verlet
219 };
220 //! String corresponding to cutoff scheme
221 const char* enumValueToString(CutoffScheme enumValue);
222
223 /*! \brief Coulomb / VdW interaction modifiers.
224  *
225  * grompp replaces eintmodPOTSHIFT_VERLET_UNSUPPORTED by eintmodPOTSHIFT.
226  * Exactcutoff is only used by Reaction-field-zero, and is not user-selectable.
227  */
228 enum class InteractionModifiers : int
229 {
230     PotShiftVerletUnsupported,
231     PotShift,
232     None,
233     PotSwitch,
234     ExactCutoff,
235     ForceSwitch,
236     Count,
237     Default = PotShiftVerletUnsupported
238 };
239 //! String corresponding to interaction modifiers
240 const char* enumValueToString(InteractionModifiers enumValue);
241
242 /*! \brief Cut-off treatment for Coulomb */
243 enum class CoulombInteractionType : int
244 {
245     Cut,
246     RF,
247     GRFNotused,
248     Pme,
249     Ewald,
250     P3mAD,
251     Poisson,
252     Switch,
253     Shift,
254     User,
255     GBNotused,
256     RFNecUnsupported,
257     EncadShiftNotused,
258     PmeUser,
259     PmeSwitch,
260     PmeUserSwitch,
261     RFZero,
262     Count,
263     Default = Cut
264 };
265 //! String corresponding to Coulomb treatment
266 const char* enumValueToString(CoulombInteractionType enumValue);
267
268 //! Ewald geometry.
269 enum class EwaldGeometry : int
270 {
271     ThreeD,
272     ThreeDC,
273     Count,
274     Default = ThreeD
275 };
276 //! String corresponding to Ewald geometry
277 const char* enumValueToString(EwaldGeometry enumValue);
278
279 //! Macro telling us whether we use reaction field
280 #define EEL_RF(e)                                                                   \
281     ((e) == CoulombInteractionType::RF || (e) == CoulombInteractionType::GRFNotused \
282      || (e) == CoulombInteractionType::RFNecUnsupported || (e) == CoulombInteractionType::RFZero)
283
284 //! Macro telling us whether we use PME
285 #define EEL_PME(e)                                                                             \
286     ((e) == CoulombInteractionType::Pme || (e) == CoulombInteractionType::PmeSwitch            \
287      || (e) == CoulombInteractionType::PmeUser || (e) == CoulombInteractionType::PmeUserSwitch \
288      || (e) == CoulombInteractionType::P3mAD)
289 //! Macro telling us whether we use PME or full Ewald
290 #define EEL_PME_EWALD(e) (EEL_PME(e) || (e) == CoulombInteractionType::Ewald)
291 //! Macro telling us whether we use full electrostatics of any sort
292 #define EEL_FULL(e) (EEL_PME_EWALD(e) || (e) == CoulombInteractionType::Poisson)
293 //! Macro telling us whether we use user defined electrostatics
294 #define EEL_USER(e)                                                                \
295     ((e) == CoulombInteractionType::User || (e) == CoulombInteractionType::PmeUser \
296      || (e) == (CoulombInteractionType::PmeUserSwitch))
297
298 //! Van der Waals interaction treatment
299 enum class VanDerWaalsType : int
300 {
301     Cut,
302     Switch,
303     Shift,
304     User,
305     EncadShiftUnused,
306     Pme,
307     Count,
308     Default = Cut
309 };
310 //! String corresponding to Van der Waals treatment
311 const char* enumValueToString(VanDerWaalsType enumValue);
312
313 //! Type of long-range VdW treatment of combination rules
314 enum class LongRangeVdW : int
315 {
316     Geom,
317     LB,
318     Count,
319     Default = Geom
320 };
321 //! String for LJPME combination rule treatment
322 const char* enumValueToString(LongRangeVdW enumValue);
323
324 //! Macro to tell us whether we use LJPME
325 #define EVDW_PME(e) ((e) == VanDerWaalsType::Pme)
326
327 /*! \brief Integrator algorithm
328  *
329  * eiSD2 has been removed, but we keep a renamed enum entry,
330  * so we can refuse to do MD with such .tpr files.
331  * eiVV is normal velocity verlet
332  * eiVVAK uses 1/2*(KE(t-dt/2)+KE(t+dt/2)) as the kinetic energy,
333  * and the half step kinetic energy for temperature control
334  */
335 enum class IntegrationAlgorithm : int
336 {
337     MD,
338     Steep,
339     CG,
340     BD,
341     SD2Removed,
342     NM,
343     LBFGS,
344     TPI,
345     TPIC,
346     SD1,
347     VV,
348     VVAK,
349     Mimic,
350     Count,
351     Default = MD
352 };
353 //! Name of the integrator algorithm
354 const char* enumValueToString(IntegrationAlgorithm enumValue);
355 //! Do we use MiMiC QM/MM?
356 #define EI_MIMIC(e) ((e) == IntegrationAlgorithm::Mimic)
357 //! Do we use velocity Verlet
358 #define EI_VV(e) ((e) == IntegrationAlgorithm::VV || (e) == IntegrationAlgorithm::VVAK)
359 //! Do we use molecular dynamics
360 #define EI_MD(e) ((e) == IntegrationAlgorithm::MD || EI_VV(e) || EI_MIMIC(e))
361 //! Do we use stochastic dynamics
362 #define EI_SD(e) ((e) == IntegrationAlgorithm::SD1)
363 //! Do we use any stochastic integrator
364 #define EI_RANDOM(e) (EI_SD(e) || (e) == IntegrationAlgorithm::BD)
365 /*above integrators may not conserve momenta*/
366 //! Do we use any type of dynamics
367 #define EI_DYNAMICS(e) (EI_MD(e) || EI_RANDOM(e))
368 //! Or do we use minimization
369 #define EI_ENERGY_MINIMIZATION(e)                                          \
370     ((e) == IntegrationAlgorithm::Steep || (e) == IntegrationAlgorithm::CG \
371      || (e) == IntegrationAlgorithm::LBFGS)
372 //! Do we apply test particle insertion
373 #define EI_TPI(e) ((e) == IntegrationAlgorithm::TPI || (e) == IntegrationAlgorithm::TPIC)
374 //! Do we deal with particle velocities
375 #define EI_STATE_VELOCITY(e) (EI_MD(e) || EI_SD(e))
376
377 //! Constraint algorithm
378 enum class ConstraintAlgorithm : int
379 {
380     Lincs,
381     Shake,
382     Count,
383     Default = Lincs
384 };
385 //! String corresponding to constraint algorithm
386 const char* enumValueToString(ConstraintAlgorithm enumValue);
387
388 //! Distance restraint refinement algorithm
389 enum class DistanceRestraintRefinement : int
390 {
391     None,
392     Simple,
393     Ensemble,
394     Count,
395     Default = None
396 };
397 //! String corresponding to distance restraint algorithm
398 const char* enumValueToString(DistanceRestraintRefinement enumValue);
399
400 //! Distance restraints weighting type
401 enum class DistanceRestraintWeighting : int
402 {
403     Conservative,
404     Equal,
405     Count,
406     Default = Conservative
407 };
408 //! String corresponding to distance restraint weighting
409 const char* enumValueToString(DistanceRestraintWeighting enumValue);
410
411 //! Combination rule algorithm.
412 enum class CombinationRule : int
413 {
414     None,
415     Geometric,
416     Arithmetic,
417     GeomSigEps,
418     Count,
419     Default = Geometric
420 };
421 //! String for combination rule algorithm
422 const char* enumValueToString(CombinationRule enumValue);
423
424 //! Van der Waals potential.
425 enum class VanDerWaalsPotential : int
426 {
427     None,
428     LJ,
429     Buckingham,
430     Count,
431     Default = LJ
432 };
433 //! String corresponding to Van der Waals potential
434 const char* enumValueToString(VanDerWaalsPotential enumValue);
435
436 //! Simulated tempering methods.
437 enum class SimulatedTempering : int
438 {
439     Geometric,
440     Exponential,
441     Linear,
442     Count,
443     Default = Geometric
444 };
445 //! String corresponding to simulated tempering
446 const char* enumValueToString(SimulatedTempering enumValue);
447
448 /*! \brief Free energy perturbation type
449  */
450 enum class FreeEnergyPerturbationType : int
451 {
452     //! there are no evaluations at other states
453     No,
454     //! treated equivalently to Static
455     Yes,
456     //! then lambdas do not change during the simulation
457     Static,
458     //! then the states change monotonically throughout the simulation
459     SlowGrowth,
460     //! then expanded ensemble simulations are occuring
461     Expanded,
462     Count,
463     Default = No
464 };
465 //! String corresponding to FEP type.
466 const char* enumValueToString(FreeEnergyPerturbationType enumValue);
467
468 //! Free energy pertubation coupling types.
469 enum class FreeEnergyPerturbationCouplingType : int
470 {
471     Fep,
472     Mass,
473     Coul,
474     Vdw,
475     Bonded,
476     Restraint,
477     Temperature,
478     Count,
479     Default = Fep
480 };
481 //! String for FEP coupling type
482 const char* enumValueToString(FreeEnergyPerturbationCouplingType enumValue);
483 //! String for FEP coupling type, singular mention.
484 const char* enumValueToStringSingular(FreeEnergyPerturbationCouplingType enumValue);
485
486 /*! \brief What to print for free energy calculations
487  *
488  * Printing the energy to the free energy dhdl file.
489  * Yes is an alias to Total, and
490  * will be converted in readir, so we never have to account for it in code.
491  */
492 enum class FreeEnergyPrintEnergy : int
493 {
494     No,
495     Total,
496     Potential,
497     Yes,
498     Count,
499     Default = No
500 };
501 //! String corresponding to printing of free energy
502 const char* enumValueToString(FreeEnergyPrintEnergy enumValue);
503
504 /*! \brief How the lambda weights are calculated
505  */
506 enum class LambdaWeightCalculation : int
507 {
508     //! don't calculate
509     No,
510     //! using the metropolis criteria
511     Metropolis,
512     //! using the Barker critera for transition weights, also called unoptimized Bennett
513     Barker,
514     //! using Barker + minimum variance for weights
515     Minvar,
516     //! Wang-Landu (using visitation counts)
517     WL,
518     //! Weighted Wang-Landau (using optimized Gibbs weighted visitation counts)
519     WWL,
520     Count,
521     Default = No
522 };
523 //! String corresponding to lambda weights
524 const char* enumValueToString(LambdaWeightCalculation enumValue);
525 //! Macro telling us whether we use expanded ensemble
526 #define ELAMSTATS_EXPANDED(e) ((e) > LambdaWeightCalculation::No)
527 //! Macro telling us whether we use some kind of Wang-Landau
528 #define EWL(e) ((e) == LambdaWeightCalculation::WL || (e) == LambdaWeightCalculation::WWL)
529
530 /*! \brief How moves in lambda are calculated
531  */
532 enum class LambdaMoveCalculation : int
533 {
534     //! don't calculate move
535     No,
536     //! using the Metropolis criteria, and 50% up and down
537     Metropolis,
538     //! using the Barker criteria, and 50% up and down
539     Barker,
540     //! computing the transition using the marginalized probabilities of the lambdas
541     Gibbs,
542     /*! \brief
543      * using the metropolized version of Gibbs
544      *
545      * Monte Carlo Strategies in Scientific computing, Liu, p. 134
546      */
547     MetropolisGibbs,
548     Count,
549     Default = No
550 };
551 //! String corresponding to lambda moves
552 const char* enumValueToString(LambdaMoveCalculation enumValue);
553
554 /*! \brief How we decide whether weights have reached equilibrium
555  */
556 enum class LambdaWeightWillReachEquilibrium : int
557 {
558     //! never stop, weights keep going
559     No,
560     //! fix the weights from the beginning; no movement
561     Yes,
562     //! stop when the WL-delta falls below a certain level
563     WLDelta,
564     //! stop when we have a certain number of samples at every step
565     NumAtLambda,
566     //! stop when we've run a certain total number of steps
567     Steps,
568     //! stop when we've run a certain total number of samples
569     Samples,
570     //! stop when the ratio of samples (lowest to highest) is sufficiently large
571     Ratio,
572     Count,
573     Default = No
574 };
575 //! String corresponding to equilibrium algorithm
576 const char* enumValueToString(LambdaWeightWillReachEquilibrium enumValue);
577
578 /*! \brief separate_dhdl_file selection
579  *
580  * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
581  * Why was this done this way, just .........
582  */
583 enum class SeparateDhdlFile : int
584 {
585     Yes,
586     No,
587     Count,
588     Default = Yes
589 };
590 //! String corresponding to separate DHDL file selection
591 const char* enumValueToString(SeparateDhdlFile enumValue);
592
593 /*! \brief dhdl_derivatives selection \
594  *
595  * NOTE: YES is the first one. Do NOT interpret this one as a gmx_bool
596  * Why was this done this way, just .........
597  */
598 enum class DhDlDerivativeCalculation : int
599 {
600     Yes,
601     No,
602     Count,
603     Default = Yes
604 };
605 //! String for DHDL derivatives
606 const char* enumValueToString(DhDlDerivativeCalculation enumValue);
607
608 /*! \brief soft-core function \
609  *
610  * Distinguishes between soft-core functions in the input.
611  */
612 enum class SoftcoreType : int
613 {
614     Beutler,
615     Gapsys,
616     Count,
617     Default = Beutler
618 };
619 //! Strings for softcore function names
620 const char* enumValueToString(SoftcoreType enumValue);
621
622 /*! \brief soft-core function as parameter to the nb-fep kernel/14-interaction.\
623  *
624  * Distinguishes between soft-core functions internally. This is different
625  * from SoftcoreType in that it offers 'None' which is not exposed to the user.
626  */
627 enum class KernelSoftcoreType : int
628 {
629     Beutler,
630     Gapsys,
631     None,
632     Count,
633     Default = Beutler
634 };
635 //! Strings for softcore function names
636 const char* enumValueToString(KernelSoftcoreType enumValue);
637
638 /*! \brief Solvent model
639  *
640  * Distinguishes classical water types with 3 or 4 particles
641  */
642 enum class SolventModel : int
643 {
644     No,
645     Spc,
646     Tip4p,
647     Count,
648     Default = Spc
649 };
650 //! String corresponding to solvent type
651 const char* enumValueToString(SolventModel enumValue);
652
653 //! Dispersion correction.
654 enum class DispersionCorrectionType : int
655 {
656     No,
657     EnerPres,
658     Ener,
659     AllEnerPres,
660     AllEner,
661     Count,
662     Default = No
663 };
664 //! String corresponding to dispersion corrections
665 const char* enumValueToString(DispersionCorrectionType enumValue);
666
667 //! Algorithm for simulated annealing.
668 enum class SimulatedAnnealing : int
669 {
670     No,
671     Single,
672     Periodic,
673     Count,
674     Default = No
675 };
676 //! String for simulated annealing
677 const char* enumValueToString(SimulatedAnnealing enumValue);
678
679 //! Wall types.
680 enum class WallType : int
681 {
682     NineThree,
683     TenFour,
684     Table,
685     TwelveSix,
686     Count,
687     Default = NineThree
688 };
689 //! String corresponding to wall type
690 const char* enumValueToString(WallType enumValue);
691
692 //! Pulling algorithm.
693 enum class PullingAlgorithm : int
694 {
695     Umbrella,
696     Constraint,
697     ConstantForce,
698     FlatBottom,
699     FlatBottomHigh,
700     External,
701     Count,
702     Default = Umbrella
703 };
704 //! String for pulling algorithm
705 const char* enumValueToString(PullingAlgorithm enumValue);
706
707 //! Control of pull groups
708 enum class PullGroupGeometry : int
709 {
710     Distance,
711     Direction,
712     Cylinder,
713     DirectionPBC,
714     DirectionRelative,
715     Angle,
716     Dihedral,
717     AngleAxis,
718     Transformation,
719     Count,
720     Default = Distance
721 };
722 //! String for pull groups
723 const char* enumValueToString(PullGroupGeometry enumValue);
724
725 //! Enforced rotation groups.
726 enum class EnforcedRotationGroupType : int
727 {
728     Iso,
729     Isopf,
730     Pm,
731     Pmpf,
732     Rm,
733     Rmpf,
734     Rm2,
735     Rm2pf,
736     Flex,
737     Flext,
738     Flex2,
739     Flex2t,
740     Count,
741     Default = Iso
742 };
743 //! Rotation group names
744 const char* enumValueToString(EnforcedRotationGroupType enumValue);
745 //! String for rotation group origin names
746 const char* enumValueToLongString(EnforcedRotationGroupType enumValue);
747
748 //! Rotation group fitting type
749 enum class RotationGroupFitting : int
750 {
751     Rmsd,
752     Norm,
753     Pot,
754     Count,
755     Default = Rmsd
756 };
757 //! String corresponding to rotation group fitting
758 const char* enumValueToString(RotationGroupFitting enumValue);
759
760 /*! \brief Direction along which ion/water swaps happen
761  *
762  * Part of "Computational Electrophysiology" (CompEL) setups
763  */
764 enum class SwapType : int
765 {
766     No,
767     X,
768     Y,
769     Z,
770     Count,
771     Default = No
772 };
773 //! Names for swapping
774 const char* enumValueToString(SwapType enumValue);
775
776 /*! \brief Swap group splitting type
777  *
778  * These are just the fixed groups we need for any setup. In t_swap's grp
779  * entry after that follows the variable number of swap groups.
780  */
781 enum class SwapGroupSplittingType : int
782 {
783     Split0,
784     Split1,
785     Solvent,
786     Count,
787     Default = Solvent
788 };
789 //! String for swap group splitting
790 const char* enumValueToString(SwapGroupSplittingType enumValue);
791
792 /*! \brief Types of electrostatics calculations
793  *
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.
797  */
798 enum class NbkernelElecType : int
799 {
800     None,
801     Coulomb,
802     ReactionField,
803     CubicSplineTable,
804     Ewald,
805     Count,
806     Default = None
807 };
808 //! String corresponding to electrostatics kernels
809 const char* enumValueToString(NbkernelElecType enumValue);
810
811 /*! \brief Types of vdw calculations available
812  *
813  * Types of vdw calculations available inside nonbonded kernels.
814  * Note that these do NOT necessarily correspond to the user selections
815  * in the MDP file; many interactions for instance map to tabulated kernels.
816  */
817 enum class NbkernelVdwType : int
818 {
819     None,
820     LennardJones,
821     Buckingham,
822     CubicSplineTable,
823     LJEwald,
824     Count,
825     Default = None
826 };
827 //! String corresponding to VdW kernels
828 const char* enumValueToString(NbkernelVdwType enumValue);
829
830 //! Center of mass motion removal algorithm.
831 enum class ComRemovalAlgorithm : int
832 {
833     Linear,
834     Angular,
835     No,
836     LinearAccelerationCorrection,
837     Count,
838     Default = Linear
839 };
840 //! String corresponding to COM removal
841 const char* enumValueToString(ComRemovalAlgorithm enumValue);
842
843 //! Enumeration that contains all supported periodic boundary setups.
844 enum class PbcType : int
845 {
846     Xyz     = 0, //!< Periodic boundaries in all dimensions.
847     No      = 1, //!< No periodic boundaries.
848     XY      = 2, //!< Only two dimensions are periodic.
849     Screw   = 3, //!< Screw.
850     Unset   = 4, //!< The type of PBC is not set or invalid.
851     Count   = 5,
852     Default = Xyz
853 };
854
855 #endif /* GMX_MDTYPES_MD_ENUMS_H */