Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdtypes / group.h
blobd39e7c8cb5cefe5e0e547a439cc408a5f975f7d0
1 /*
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,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.
37 #ifndef GMX_MDTYPES_GROUP_H
38 #define GMX_MDTYPES_GROUP_H
40 #include <vector>
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/utility/basedefinitions.h"
44 #include "gromacs/utility/real.h"
45 #include "gromacs/utility/smalloc.h"
47 struct t_grp_tcstat
49 //! Temperature at half step
50 real Th = 0;
51 //! Temperature at full step
52 real T = 0;
53 //! Kinetic energy at half step
54 tensor ekinh = { { 0 } };
55 //! Kinetic energy at old half step
56 tensor ekinh_old = { { 0 } };
57 //! Kinetic energy at full step
58 tensor ekinf = { { 0 } };
59 //! Berendsen coupling lambda
60 real lambda = 0;
61 //! Scaling factor for NHC- full step
62 double ekinscalef_nhc = 0;
63 //! Scaling factor for NHC- half step
64 double ekinscaleh_nhc = 0;
65 //! Scaling factor for NHC- velocity
66 double vscale_nhc = 0;
69 struct t_grp_acc
71 //! Number of atoms in this group
72 int nat = 0;
73 //! Mean velocities of home particles
74 rvec u = { 0 };
75 //! Previous mean velocities of home particles
76 rvec uold = { 0 };
77 //! Mass for topology A
78 double mA = 0;
79 //! Mass for topology B
80 double mB = 0;
83 struct t_cos_acc
85 //! The acceleration for the cosine profile
86 real cos_accel = 0;
87 //! The cos momenta of home particles
88 real mvcos = 0;
89 //! The velocity of the cosine profile
90 real vcos = 0;
93 struct gmx_ekindata_t
95 //! Whether non-equilibrium MD is active (ie. constant or cosine acceleration)
96 gmx_bool bNEMD;
97 //! The number of T-coupling groups
98 int ngtc = 0;
99 //! For size of ekin_work
100 int nthreads = 0;
101 //! T-coupling data
102 std::vector<t_grp_tcstat> tcstat;
103 //! Allocated locations for *_work members
104 tensor** ekin_work_alloc = nullptr;
105 //! Work arrays for tcstat per thread
106 tensor** ekin_work = nullptr;
107 //! Work location for dekindl per thread
108 real** dekindl_work = nullptr;
109 //! The number of acceleration groups
110 int ngacc = 0;
111 //! Acceleration data
112 std::vector<t_grp_acc> grpstat;
113 //! overall kinetic energy
114 tensor ekin = { { 0 } };
115 //! overall 1/2 step kinetic energy
116 tensor ekinh = { { 0 } };
117 //! dEkin/dlambda at half step
118 real dekindl = 0;
119 //! dEkin/dlambda at old half step
120 real dekindl_old = 0;
121 //! Cosine acceleration data
122 t_cos_acc cosacc;
124 ~gmx_ekindata_t();
127 #define GID(igid, jgid, gnr) \
128 (((igid) < (jgid)) ? ((igid) * (gnr) + (jgid)) : ((jgid) * (gnr) + (igid)))
130 #endif