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, 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.
45 #include "gromacs/domdec/domdec.h"
46 #include "gromacs/domdec/domdec_struct.h"
47 #include "gromacs/fileio/checkpoint.h"
48 #include "gromacs/fileio/xtcio.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdlib/constr.h"
53 #include "gromacs/mdlib/md_support.h"
54 #include "gromacs/mdlib/rbin.h"
55 #include "gromacs/mdlib/tgroup.h"
56 #include "gromacs/mdlib/vcm.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/mdtypes/enerdata.h"
59 #include "gromacs/mdtypes/group.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/futil.h"
64 #include "gromacs/utility/smalloc.h"
66 typedef struct gmx_global_stat
73 gmx_global_stat_t
global_stat_init(const t_inputrec
* ir
)
80 snew(gs
->itc0
, ir
->opts
.ngtc
);
81 snew(gs
->itc1
, ir
->opts
.ngtc
);
86 void global_stat_destroy(gmx_global_stat_t gs
)
94 static int filter_enerdterm(const real
* afrom
, gmx_bool bToBuffer
, real
* ato
, gmx_bool bTemp
, gmx_bool bPres
, gmx_bool bEner
)
100 for (i
= 0; i
< F_NRE
; i
++)
117 ato
[to
++] = afrom
[from
++];
124 ato
[to
++] = afrom
[from
++];
130 ato
[to
++] = afrom
[from
++];
139 void global_stat(const gmx_global_stat
* gs
,
141 gmx_enerdata_t
* enerd
,
144 const t_inputrec
* inputrec
,
145 gmx_ekindata_t
* ekind
,
146 const gmx::Constraints
* constr
,
150 int* totalNumberOfBondedInteractions
,
151 gmx_bool bSumEkinhOld
,
153 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
157 int ie
= 0, ifv
= 0, isv
= 0, irmsd
= 0;
158 int idedl
= 0, idedlo
= 0, idvdll
= 0, idvdlnl
= 0, iepl
= 0, icm
= 0, imass
= 0, ica
= 0, inb
= 0;
160 int icj
= -1, ici
= -1, icx
= -1;
162 real copyenerd
[F_NRE
];
165 gmx_bool bVV
, bTemp
, bEner
, bPres
, bConstrVir
, bEkinAveVel
, bReadEkin
;
166 bool checkNumberOfBondedInteractions
= (flags
& CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
) != 0;
168 bVV
= EI_VV(inputrec
->eI
);
169 bTemp
= ((flags
& CGLO_TEMPERATURE
) != 0);
170 bEner
= ((flags
& CGLO_ENERGY
) != 0);
171 bPres
= ((flags
& CGLO_PRESSURE
) != 0);
172 bConstrVir
= ((flags
& CGLO_CONSTRAINT
) != 0);
173 bEkinAveVel
= (inputrec
->eI
== eiVV
|| (inputrec
->eI
== eiVVAK
&& bPres
));
174 bReadEkin
= ((flags
& CGLO_READEKIN
) != 0);
182 /* This routine copies all the data to be summed to one big buffer
183 * using the t_bin struct.
186 /* First, we neeed to identify which enerd->term should be
187 communicated. Temperature and pressure terms should only be
188 communicated and summed when they need to be, to avoid repeating
189 the sums and overcounting. */
191 nener
= filter_enerdterm(enerd
->term
, TRUE
, copyenerd
, bTemp
, bPres
, bEner
);
193 /* First, the data that needs to be communicated with velocity verlet every time
194 This is just the constraint virial.*/
197 isv
= add_binr(rb
, DIM
* DIM
, svir
[0]);
200 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
205 for (j
= 0; (j
< inputrec
->opts
.ngtc
); j
++)
209 itc0
[j
] = add_binr(rb
, DIM
* DIM
, ekind
->tcstat
[j
].ekinh_old
[0]);
211 if (bEkinAveVel
&& !bReadEkin
)
213 itc1
[j
] = add_binr(rb
, DIM
* DIM
, ekind
->tcstat
[j
].ekinf
[0]);
217 itc1
[j
] = add_binr(rb
, DIM
* DIM
, ekind
->tcstat
[j
].ekinh
[0]);
220 /* these probably need to be put into one of these categories */
221 idedl
= add_binr(rb
, 1, &(ekind
->dekindl
));
224 idedlo
= add_binr(rb
, 1, &(ekind
->dekindl_old
));
226 if (ekind
->cosacc
.cos_accel
!= 0)
228 ica
= add_binr(rb
, 1, &(ekind
->cosacc
.mvcos
));
235 ifv
= add_binr(rb
, DIM
* DIM
, fvir
[0]);
238 gmx::ArrayRef
<real
> rmsdData
;
241 ie
= add_binr(rb
, nener
, copyenerd
);
244 rmsdData
= constr
->rmsdData();
245 if (!rmsdData
.empty())
247 irmsd
= add_binr(rb
, 2, rmsdData
.data());
251 for (j
= 0; (j
< egNR
); j
++)
253 inn
[j
] = add_binr(rb
, enerd
->grpp
.nener
, enerd
->grpp
.ener
[j
].data());
255 if (inputrec
->efep
!= efepNO
)
257 idvdll
= add_bind(rb
, efptNR
, enerd
->dvdl_lin
);
258 idvdlnl
= add_bind(rb
, efptNR
, enerd
->dvdl_nonlin
);
259 if (enerd
->foreignLambdaTerms
.numLambdas() > 0)
261 iepl
= add_bind(rb
, enerd
->foreignLambdaTerms
.energies().size(),
262 enerd
->foreignLambdaTerms
.energies().data());
269 icm
= add_binr(rb
, DIM
* vcm
->nr
, vcm
->group_p
[0]);
270 imass
= add_binr(rb
, vcm
->nr
, vcm
->group_mass
.data());
271 if (vcm
->mode
== ecmANGULAR
)
273 icj
= add_binr(rb
, DIM
* vcm
->nr
, vcm
->group_j
[0]);
274 icx
= add_binr(rb
, DIM
* vcm
->nr
, vcm
->group_x
[0]);
275 ici
= add_binr(rb
, DIM
* DIM
* vcm
->nr
, vcm
->group_i
[0][0]);
279 if (checkNumberOfBondedInteractions
)
281 nb
= cr
->dd
->nbonded_local
;
282 inb
= add_bind(rb
, 1, &nb
);
286 isig
= add_binr(rb
, nsig
, sig
);
289 /* Global sum it all */
292 fprintf(debug
, "Summing %d energies\n", rb
->maxreal
);
296 /* Extract all the data locally */
300 extract_binr(rb
, isv
, DIM
* DIM
, svir
[0]);
303 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
308 for (j
= 0; (j
< inputrec
->opts
.ngtc
); j
++)
312 extract_binr(rb
, itc0
[j
], DIM
* DIM
, ekind
->tcstat
[j
].ekinh_old
[0]);
314 if (bEkinAveVel
&& !bReadEkin
)
316 extract_binr(rb
, itc1
[j
], DIM
* DIM
, ekind
->tcstat
[j
].ekinf
[0]);
320 extract_binr(rb
, itc1
[j
], DIM
* DIM
, ekind
->tcstat
[j
].ekinh
[0]);
323 extract_binr(rb
, idedl
, 1, &(ekind
->dekindl
));
326 extract_binr(rb
, idedlo
, 1, &(ekind
->dekindl_old
));
328 if (ekind
->cosacc
.cos_accel
!= 0)
330 extract_binr(rb
, ica
, 1, &(ekind
->cosacc
.mvcos
));
336 extract_binr(rb
, ifv
, DIM
* DIM
, fvir
[0]);
341 extract_binr(rb
, ie
, nener
, copyenerd
);
342 if (!rmsdData
.empty())
344 extract_binr(rb
, irmsd
, rmsdData
);
347 for (j
= 0; (j
< egNR
); j
++)
349 extract_binr(rb
, inn
[j
], enerd
->grpp
.nener
, enerd
->grpp
.ener
[j
].data());
351 if (inputrec
->efep
!= efepNO
)
353 extract_bind(rb
, idvdll
, efptNR
, enerd
->dvdl_lin
);
354 extract_bind(rb
, idvdlnl
, efptNR
, enerd
->dvdl_nonlin
);
355 if (enerd
->foreignLambdaTerms
.numLambdas() > 0)
357 extract_bind(rb
, iepl
, enerd
->foreignLambdaTerms
.energies().size(),
358 enerd
->foreignLambdaTerms
.energies().data());
362 filter_enerdterm(copyenerd
, FALSE
, enerd
->term
, bTemp
, bPres
, bEner
);
367 extract_binr(rb
, icm
, DIM
* vcm
->nr
, vcm
->group_p
[0]);
368 extract_binr(rb
, imass
, vcm
->nr
, vcm
->group_mass
.data());
369 if (vcm
->mode
== ecmANGULAR
)
371 extract_binr(rb
, icj
, DIM
* vcm
->nr
, vcm
->group_j
[0]);
372 extract_binr(rb
, icx
, DIM
* vcm
->nr
, vcm
->group_x
[0]);
373 extract_binr(rb
, ici
, DIM
* DIM
* vcm
->nr
, vcm
->group_i
[0][0]);
377 if (checkNumberOfBondedInteractions
)
379 extract_bind(rb
, inb
, 1, &nb
);
380 *totalNumberOfBondedInteractions
= gmx::roundToInt(nb
);
385 extract_binr(rb
, isig
, nsig
, sig
);