Replace gmx_large_int_t with gmx_int64_t
[alexxy/gromacs.git] / src / gromacs / fileio / trajectory_writing_low_level.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 #include "trajectory_writing.h"
37 #include "gromacs/legacyheaders/mvdata.h"
38 #include "gromacs/legacyheaders/domdec.h"
39 #include "checkpoint.h"
40 #include "trnio.h"
41 #include "xtcio.h"
42 #include "gromacs/legacyheaders/smalloc.h"
43
44 static void moveit(t_commrec *cr, rvec xx[])
45 {
46     if (!xx)
47     {
48         return;
49     }
50
51     move_rvecs(cr, FALSE, FALSE, xx, NULL, (cr->nnodes-cr->npmenodes)-1, NULL);
52 }
53
54 void write_traj(FILE *fplog, t_commrec *cr,
55                 gmx_mdoutf_t *of,
56                 int mdof_flags,
57                 gmx_mtop_t *top_global,
58                 gmx_int64_t step, double t,
59                 t_state *state_local, t_state *state_global,
60                 rvec *f_local, rvec *f_global,
61                 int *n_xtc, rvec **x_xtc)
62 {
63     int           i, j;
64     gmx_groups_t *groups;
65     rvec         *xxtc;
66     rvec         *local_v;
67     rvec         *global_v;
68
69 #define MX(xvf) moveit(cr, xvf)
70
71     /* MRS -- defining these variables is to manage the difference
72      * between half step and full step velocities, but there must be a better way . . . */
73
74     local_v  = state_local->v;
75     global_v = state_global->v;
76
77     if (DOMAINDECOMP(cr))
78     {
79         if (mdof_flags & MDOF_CPT)
80         {
81             dd_collect_state(cr->dd, state_local, state_global);
82         }
83         else
84         {
85             if (mdof_flags & (MDOF_X | MDOF_XTC))
86             {
87                 dd_collect_vec(cr->dd, state_local, state_local->x,
88                                state_global->x);
89             }
90             if (mdof_flags & MDOF_V)
91             {
92                 dd_collect_vec(cr->dd, state_local, local_v,
93                                global_v);
94             }
95         }
96         if (mdof_flags & MDOF_F)
97         {
98             dd_collect_vec(cr->dd, state_local, f_local, f_global);
99         }
100     }
101     else
102     {
103         if (mdof_flags & MDOF_CPT)
104         {
105             /* All pointers in state_local are equal to state_global,
106              * but we need to copy the non-pointer entries.
107              */
108             state_global->lambda = state_local->lambda;
109             state_global->veta   = state_local->veta;
110             state_global->vol0   = state_local->vol0;
111             copy_mat(state_local->box, state_global->box);
112             copy_mat(state_local->boxv, state_global->boxv);
113             copy_mat(state_local->svir_prev, state_global->svir_prev);
114             copy_mat(state_local->fvir_prev, state_global->fvir_prev);
115             copy_mat(state_local->pres_prev, state_global->pres_prev);
116         }
117         if (cr->nnodes > 1)
118         {
119             /* Particle decomposition, collect the data on the master node */
120             if (mdof_flags & MDOF_CPT)
121             {
122                 if (state_local->flags & (1<<estX))
123                 {
124                     MX(state_global->x);
125                 }
126                 if (state_local->flags & (1<<estV))
127                 {
128                     MX(state_global->v);
129                 }
130                 if (state_local->flags & (1<<estSDX))
131                 {
132                     MX(state_global->sd_X);
133                 }
134                 if (state_global->nrngi > 1)
135                 {
136                     if (state_local->flags & (1<<estLD_RNG))
137                     {
138 #ifdef GMX_MPI
139                         MPI_Gather(state_local->ld_rng,
140                                    state_local->nrng*sizeof(state_local->ld_rng[0]), MPI_BYTE,
141                                    state_global->ld_rng,
142                                    state_local->nrng*sizeof(state_local->ld_rng[0]), MPI_BYTE,
143                                    MASTERRANK(cr), cr->mpi_comm_mygroup);
144 #endif
145                     }
146                     if (state_local->flags & (1<<estLD_RNGI))
147                     {
148 #ifdef GMX_MPI
149                         MPI_Gather(state_local->ld_rngi,
150                                    sizeof(state_local->ld_rngi[0]), MPI_BYTE,
151                                    state_global->ld_rngi,
152                                    sizeof(state_local->ld_rngi[0]), MPI_BYTE,
153                                    MASTERRANK(cr), cr->mpi_comm_mygroup);
154 #endif
155                     }
156                 }
157             }
158             else
159             {
160                 if (mdof_flags & (MDOF_X | MDOF_XTC))
161                 {
162                     MX(state_global->x);
163                 }
164                 if (mdof_flags & MDOF_V)
165                 {
166                     MX(global_v);
167                 }
168             }
169             if (mdof_flags & MDOF_F)
170             {
171                 MX(f_global);
172             }
173         }
174     }
175
176     if (MASTER(cr))
177     {
178         if (mdof_flags & MDOF_CPT)
179         {
180             write_checkpoint(of->fn_cpt, of->bKeepAndNumCPT,
181                              fplog, cr, of->eIntegrator, of->simulation_part,
182                              of->bExpanded, of->elamstats, step, t, state_global);
183         }
184
185         if (mdof_flags & (MDOF_X | MDOF_V | MDOF_F))
186         {
187             fwrite_trn(of->fp_trn, step, t, state_local->lambda[efptFEP],
188                        state_local->box, top_global->natoms,
189                        (mdof_flags & MDOF_X) ? state_global->x : NULL,
190                        (mdof_flags & MDOF_V) ? global_v : NULL,
191                        (mdof_flags & MDOF_F) ? f_global : NULL);
192             if (gmx_fio_flush(of->fp_trn) != 0)
193             {
194                 gmx_file("Cannot write trajectory; maybe you are out of disk space?");
195             }
196             gmx_fio_check_file_position(of->fp_trn);
197         }
198         if (mdof_flags & MDOF_XTC)
199         {
200             groups = &top_global->groups;
201             if (*n_xtc == -1)
202             {
203                 *n_xtc = 0;
204                 for (i = 0; (i < top_global->natoms); i++)
205                 {
206                     if (ggrpnr(groups, egcXTC, i) == 0)
207                     {
208                         (*n_xtc)++;
209                     }
210                 }
211                 if (*n_xtc != top_global->natoms)
212                 {
213                     snew(*x_xtc, *n_xtc);
214                 }
215             }
216             if (*n_xtc == top_global->natoms)
217             {
218                 xxtc = state_global->x;
219             }
220             else
221             {
222                 xxtc = *x_xtc;
223                 j    = 0;
224                 for (i = 0; (i < top_global->natoms); i++)
225                 {
226                     if (ggrpnr(groups, egcXTC, i) == 0)
227                     {
228                         copy_rvec(state_global->x[i], xxtc[j++]);
229                     }
230                 }
231             }
232             if (write_xtc(of->fp_xtc, *n_xtc, step, t,
233                           state_local->box, xxtc, of->xtc_prec) == 0)
234             {
235                 gmx_fatal(FARGS, "XTC error - maybe you are out of disk space?");
236             }
237             gmx_fio_check_file_position(of->fp_xtc);
238         }
239     }
240 }