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) 2012,2014,2015,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 #include "compute_io.h"
44 #include "gromacs/mdtypes/inputrec.h"
45 #include "gromacs/mdtypes/md_enums.h"
46 #include "gromacs/mdtypes/pull_params.h"
47 #include "gromacs/topology/topology.h"
49 static int div_nsteps(int nsteps
, int nst
)
53 return (1 + nsteps
+ nst
- 1)/nst
;
61 double compute_io(const t_inputrec
*ir
, int natoms
, const SimulationGroups
&groups
,
62 int nrener
, int nrepl
)
65 int nsteps
= ir
->nsteps
;
67 int nstx
, nstv
, nstf
, nste
, nstlog
, nstxtc
;
70 nstx
= div_nsteps(nsteps
, ir
->nstxout
);
71 nstv
= div_nsteps(nsteps
, ir
->nstvout
);
72 nstf
= div_nsteps(nsteps
, ir
->nstfout
);
73 nstxtc
= div_nsteps(nsteps
, ir
->nstxout_compressed
);
74 if (ir
->nstxout_compressed
> 0)
76 for (int i
= 0; i
< natoms
; i
++)
78 if (groups
.groupNumbers
[SimulationAtomGroupType::CompressedPositionOutput
].empty() ||
79 groups
.groupNumbers
[SimulationAtomGroupType::CompressedPositionOutput
][i
] == 0)
85 nstlog
= div_nsteps(nsteps
, ir
->nstlog
);
86 /* We add 2 for the header */
87 nste
= div_nsteps(2+nsteps
, ir
->nstenergy
);
90 cio
+= (nstx
+nstf
+nstv
)*sizeof(real
)*(3.0*natoms
);
91 cio
+= nstxtc
*(14*4 + nxtcatoms
*5.0); /* roughly 5 bytes per atom */
92 cio
+= nstlog
*(nrener
*16*2.0); /* 16 bytes per energy term plus header */
93 /* t_energy contains doubles, but real is written to edr */
94 cio
+= (1.0*nste
)*nrener
*3*sizeof(real
);
96 if ((ir
->efep
!= efepNO
|| ir
->bSimTemp
) && (ir
->fepvals
->nstdhdl
> 0))
98 int ndh
= ir
->fepvals
->n_lambda
;
102 for (i
= 0; i
< efptNR
; i
++)
104 if (ir
->fepvals
->separate_dvdl
[i
])
110 if (ir
->fepvals
->separate_dhdl_file
== esepdhdlfileYES
)
112 nchars
= 8 + ndhdl
*8 + ndh
*10; /* time data ~8 chars/entry, dH data ~10 chars/entry */
113 if (ir
->expandedvals
->elmcmove
> elmcmoveNO
)
115 nchars
+= 5; /* alchemical state */
118 if (ir
->fepvals
->edHdLPrintEnergy
!= edHdLPrintEnergyNO
)
120 nchars
+= 12; /* energy for dhdl */
122 cio
+= div_nsteps(nsteps
, ir
->fepvals
->nstdhdl
)*nchars
;
126 /* dH output to ener.edr: */
127 if (ir
->fepvals
->dh_hist_size
<= 0)
129 int ndh_tot
= ndh
+ndhdl
;
130 if (ir
->expandedvals
->elmcmove
> elmcmoveNO
)
134 if (ir
->fepvals
->edHdLPrintEnergy
!= edHdLPrintEnergyNO
)
138 /* as data blocks: 1 real per dH point */
139 cio
+= div_nsteps(nsteps
, ir
->fepvals
->nstdhdl
)*ndh_tot
*sizeof(real
);
143 /* as histograms: dh_hist_size ints per histogram */
144 cio
+= div_nsteps(nsteps
, ir
->nstenergy
)*
145 sizeof(int)*ir
->fepvals
->dh_hist_size
*ndh
;
149 if (ir
->pull
!= nullptr)
151 cio
+= div_nsteps(nsteps
, ir
->pull
->nstxout
)*20; /* roughly 20 chars per line */
152 cio
+= div_nsteps(nsteps
, ir
->pull
->nstfout
)*20; /* roughly 20 chars per line */
155 return cio
*nrepl
/(1024*1024);