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,2018, 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 "convert_tpr.h"
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/fileio/checkpoint.h"
45 #include "gromacs/fileio/enxio.h"
46 #include "gromacs/fileio/tpxio.h"
47 #include "gromacs/fileio/trrio.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/mdtypes/state.h"
52 #include "gromacs/random/seed.h"
53 #include "gromacs/topology/ifunc.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/arraysize.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/real.h"
62 #include "gromacs/utility/smalloc.h"
63 #include "gromacs/utility/stringutil.h"
65 #define RANGECHK(i, n) if ((i) >= (n)) gmx_fatal(FARGS, "Your index file contains atomnumbers (e.g. %d)\nthat are larger than the number of atoms in the tpr file (%d)", (i), (n))
67 static gmx_bool
*bKeepIt(int gnx
, int natoms
, int index
[])
73 for (i
= 0; (i
< gnx
); i
++)
75 RANGECHK(index
[i
], natoms
);
82 static int *invind(int gnx
, int natoms
, int index
[])
88 for (i
= 0; (i
< gnx
); i
++)
90 RANGECHK(index
[i
], natoms
);
97 static void reduce_block(const gmx_bool bKeep
[], t_block
*block
,
101 int i
, j
, newi
, newj
;
103 snew(index
, block
->nr
);
106 for (i
= 0; (i
< block
->nr
); i
++)
108 for (j
= block
->index
[i
]; (j
< block
->index
[i
+1]); j
++)
115 if (newj
> index
[newi
])
122 fprintf(stderr
, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
123 name
, block
->nr
, newi
, block
->index
[block
->nr
], newj
);
124 block
->index
= index
;
128 static void reduce_blocka(const int invindex
[], const gmx_bool bKeep
[], t_blocka
*block
,
132 int i
, j
, k
, newi
, newj
;
134 snew(index
, block
->nr
);
138 for (i
= 0; (i
< block
->nr
); i
++)
140 for (j
= block
->index
[i
]; (j
< block
->index
[i
+1]); j
++)
145 a
[newj
] = invindex
[k
];
149 if (newj
> index
[newi
])
156 fprintf(stderr
, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
157 name
, block
->nr
, newi
, block
->nra
, newj
);
158 block
->index
= index
;
164 static void reduce_rvec(int gnx
, const int index
[], rvec vv
[])
170 for (i
= 0; (i
< gnx
); i
++)
172 copy_rvec(vv
[index
[i
]], ptr
[i
]);
174 for (i
= 0; (i
< gnx
); i
++)
176 copy_rvec(ptr
[i
], vv
[i
]);
181 static void reduce_atom(int gnx
, const int index
[], t_atom atom
[], char ***atomname
,
182 int *nres
, t_resinfo
*resinfo
)
191 snew(rinfo
, atom
[index
[gnx
-1]].resind
+1);
192 for (i
= 0; (i
< gnx
); i
++)
194 ptr
[i
] = atom
[index
[i
]];
195 aname
[i
] = atomname
[index
[i
]];
198 for (i
= 0; (i
< gnx
); i
++)
201 atomname
[i
] = aname
[i
];
202 if ((i
== 0) || (atom
[i
].resind
!= atom
[i
-1].resind
))
205 rinfo
[nr
] = resinfo
[atom
[i
].resind
];
210 for (i
= 0; (i
< nr
); i
++)
212 resinfo
[i
] = rinfo
[i
];
221 static void reduce_ilist(const int invindex
[], const gmx_bool bKeep
[],
222 t_ilist
*il
, int nratoms
, const char *name
)
232 for (i
= 0; (i
< il
->nr
); i
+= nratoms
+1)
235 for (j
= 1; (j
<= nratoms
); j
++)
237 bB
= bB
&& bKeep
[il
->iatoms
[i
+j
]];
241 ia
[newnr
++] = il
->iatoms
[i
];
242 for (j
= 1; (j
<= nratoms
); j
++)
244 ia
[newnr
++] = invindex
[il
->iatoms
[i
+j
]];
248 fprintf(stderr
, "Reduced ilist %8s from %6d to %6d entries\n",
249 name
, il
->nr
/(nratoms
+1),
253 for (i
= 0; (i
< newnr
); i
++)
255 il
->iatoms
[i
] = ia
[i
];
262 static void reduce_topology_x(int gnx
, int index
[],
263 gmx_mtop_t
*mtop
, rvec x
[], rvec v
[])
270 top
= gmx_mtop_t_to_t_topology(mtop
, false);
271 bKeep
= bKeepIt(gnx
, top
.atoms
.nr
, index
);
272 invindex
= invind(gnx
, top
.atoms
.nr
, index
);
274 reduce_block(bKeep
, &(top
.cgs
), "cgs");
275 reduce_block(bKeep
, &(top
.mols
), "mols");
276 reduce_blocka(invindex
, bKeep
, &(top
.excls
), "excls");
277 reduce_rvec(gnx
, index
, x
);
278 reduce_rvec(gnx
, index
, v
);
279 reduce_atom(gnx
, index
, top
.atoms
.atom
, top
.atoms
.atomname
,
280 &(top
.atoms
.nres
), top
.atoms
.resinfo
);
282 for (i
= 0; (i
< F_NRE
); i
++)
284 reduce_ilist(invindex
, bKeep
, &(top
.idef
.il
[i
]),
285 interaction_function
[i
].nratoms
,
286 interaction_function
[i
].name
);
291 mtop
->moltype
.resize(1);
292 mtop
->moltype
[0].name
= mtop
->name
;
293 mtop
->moltype
[0].atoms
= top
.atoms
;
294 for (i
= 0; i
< F_NRE
; i
++)
296 InteractionList
&ilist
= mtop
->moltype
[0].ilist
[i
];
297 ilist
.iatoms
.resize(top
.idef
.il
[i
].nr
);
298 for (int j
= 0; j
< top
.idef
.il
[i
].nr
; j
++)
300 ilist
.iatoms
[j
] = top
.idef
.il
[i
].iatoms
[j
];
303 mtop
->moltype
[0].atoms
= top
.atoms
;
304 mtop
->moltype
[0].cgs
= top
.cgs
;
305 mtop
->moltype
[0].excls
= top
.excls
;
307 mtop
->molblock
.resize(1);
308 mtop
->molblock
[0].type
= 0;
309 mtop
->molblock
[0].nmol
= 1;
311 mtop
->natoms
= top
.atoms
.nr
;
314 static void zeroq(const int index
[], gmx_mtop_t
*mtop
)
316 for (gmx_moltype_t
&moltype
: mtop
->moltype
)
318 for (int i
= 0; i
< moltype
.atoms
.nr
; i
++)
320 moltype
.atoms
.atom
[index
[i
]].q
= 0;
321 moltype
.atoms
.atom
[index
[i
]].qB
= 0;
326 int gmx_convert_tpr(int argc
, char *argv
[])
328 const char *desc
[] = {
329 "[THISMODULE] can edit run input files in three ways.[PAR]",
330 "[BB]1.[bb] by modifying the number of steps in a run input file",
331 "with options [TT]-extend[tt], [TT]-until[tt] or [TT]-nsteps[tt]",
332 "(nsteps=-1 means unlimited number of steps)[PAR]",
333 "[BB]2.[bb] by creating a [REF].tpx[ref] file for a subset of your original",
334 "tpx file, which is useful when you want to remove the solvent from",
335 "your [REF].tpx[ref] file, or when you want to make e.g. a pure C[GRK]alpha[grk] [REF].tpx[ref] file.",
336 "Note that you may need to use [TT]-nsteps -1[tt] (or similar) to get",
338 "[BB]WARNING: this [REF].tpx[ref] file is not fully functional[bb].[PAR]",
339 "[BB]3.[bb] by setting the charges of a specified group",
340 "to zero. This is useful when doing free energy estimates",
341 "using the LIE (Linear Interaction Energy) method."
346 int64_t nsteps_req
, run_step
;
347 double run_t
, state_t
;
349 gmx_bool bNsteps
, bExtend
, bUntil
;
355 int *index
= nullptr;
356 char buf
[200], buf2
[200];
357 gmx_output_env_t
*oenv
;
359 { efTPR
, nullptr, nullptr, ffREAD
},
360 { efNDX
, nullptr, nullptr, ffOPTRD
},
361 { efTPR
, "-o", "tprout", ffWRITE
}
363 #define NFILE asize(fnm)
365 /* Command line options */
366 static int nsteps_req_int
= 0;
367 static real extend_t
= 0.0, until_t
= 0.0;
368 static gmx_bool bZeroQ
= FALSE
;
369 static t_pargs pa
[] = {
370 { "-extend", FALSE
, etREAL
, {&extend_t
},
371 "Extend runtime by this amount (ps)" },
372 { "-until", FALSE
, etREAL
, {&until_t
},
373 "Extend runtime until this ending time (ps)" },
374 { "-nsteps", FALSE
, etINT
, {&nsteps_req_int
},
375 "Change the number of steps" },
376 { "-zeroq", FALSE
, etBOOL
, {&bZeroQ
},
377 "Set the charges of a group (from the index) to zero" }
380 /* Parse the command line */
381 if (!parse_common_args(&argc
, argv
, 0, NFILE
, fnm
, asize(pa
), pa
,
382 asize(desc
), desc
, 0, nullptr, &oenv
))
387 /* Convert int to int64_t */
388 nsteps_req
= nsteps_req_int
;
389 bNsteps
= opt2parg_bSet("-nsteps", asize(pa
), pa
);
390 bExtend
= opt2parg_bSet("-extend", asize(pa
), pa
);
391 bUntil
= opt2parg_bSet("-until", asize(pa
), pa
);
393 top_fn
= ftp2fn(efTPR
, NFILE
, fnm
);
394 fprintf(stderr
, "Reading toplogy and stuff from %s\n", top_fn
);
396 t_inputrec irInstance
;
397 t_inputrec
*ir
= &irInstance
;
398 read_tpx_state(top_fn
, ir
, &state
, &mtop
);
399 run_step
= ir
->init_step
;
400 run_t
= ir
->init_step
*ir
->delta_t
+ ir
->init_t
;
404 fprintf(stderr
, "Setting nsteps to %s\n", gmx_step_str(nsteps_req
, buf
));
405 ir
->nsteps
= nsteps_req
;
409 /* Determine total number of steps remaining */
412 ir
->nsteps
= ir
->nsteps
- (run_step
- ir
->init_step
) + gmx::roundToInt64(extend_t
/ir
->delta_t
);
413 printf("Extending remaining runtime of by %g ps (now %s steps)\n",
414 extend_t
, gmx_step_str(ir
->nsteps
, buf
));
418 printf("nsteps = %s, run_step = %s, current_t = %g, until = %g\n",
419 gmx_step_str(ir
->nsteps
, buf
),
420 gmx_step_str(run_step
, buf2
),
422 ir
->nsteps
= gmx::roundToInt64((until_t
- run_t
)/ir
->delta_t
);
423 printf("Extending remaining runtime until %g ps (now %s steps)\n",
424 until_t
, gmx_step_str(ir
->nsteps
, buf
));
428 ir
->nsteps
-= run_step
- ir
->init_step
;
430 printf("%s steps (%g ps) remaining from first run.\n",
431 gmx_step_str(ir
->nsteps
, buf
), ir
->nsteps
*ir
->delta_t
);
435 if (bNsteps
|| bZeroQ
|| (ir
->nsteps
> 0))
437 ir
->init_step
= run_step
;
439 if (ftp2bSet(efNDX
, NFILE
, fnm
) ||
440 !(bNsteps
|| bExtend
|| bUntil
))
442 atoms
= gmx_mtop_global_atoms(&mtop
);
443 get_index(&atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
), 1,
444 &gnx
, &index
, &grpname
);
447 bSel
= (gnx
!= state
.natoms
);
448 for (i
= 0; ((i
< gnx
) && (!bSel
)); i
++)
450 bSel
= (i
!= index
[i
]);
459 fprintf(stderr
, "Will write subset %s of original tpx containing %d "
460 "atoms\n", grpname
, gnx
);
461 reduce_topology_x(gnx
, index
, &mtop
, as_rvec_array(state
.x
.data()), as_rvec_array(state
.v
.data()));
467 fprintf(stderr
, "Zero-ing charges for group %s\n", grpname
);
471 fprintf(stderr
, "Will write full tpx file (no selection)\n");
475 state_t
= ir
->init_t
+ ir
->init_step
*ir
->delta_t
;
476 sprintf(buf
, "Writing statusfile with starting step %s%s and length %s%s steps...\n", "%10", PRId64
, "%10", PRId64
);
477 fprintf(stderr
, buf
, ir
->init_step
, ir
->nsteps
);
478 fprintf(stderr
, " time %10.3f and length %10.3f ps\n",
479 state_t
, ir
->nsteps
*ir
->delta_t
);
480 write_tpx_state(opt2fn("-o", NFILE
, fnm
), ir
, &state
, &mtop
);
484 printf("You've simulated long enough. Not writing tpr file\n");