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, 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.
43 #include "gromacs/legacyheaders/types/inputrec.h"
44 #include "gromacs/legacyheaders/types/simple.h"
45 #include "gromacs/legacyheaders/types/state.h"
49 #include "gromacs/gmxpreprocess/readir.h"
50 #include "gromacs/topology/mtop_util.h"
51 #include "checkpoint.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/fileio/trnio.h"
54 #include "gromacs/fileio/enxio.h"
55 #include "gromacs/utility/futil.h"
57 #include "gromacs/commandline/pargs.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/random/random.h"
60 #include "gromacs/topology/topology.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/smalloc.h"
64 #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))
66 static gmx_bool *bKeepIt(int gnx, int natoms, atom_id index[])
72 for (i = 0; (i < gnx); i++)
74 RANGECHK(index[i], natoms);
81 static atom_id *invind(int gnx, int natoms, atom_id index[])
87 for (i = 0; (i < gnx); i++)
89 RANGECHK(index[i], natoms);
96 static void reduce_block(gmx_bool bKeep[], t_block *block,
100 int i, j, newi, newj;
102 snew(index, block->nr);
105 for (i = 0; (i < block->nr); i++)
107 for (j = block->index[i]; (j < block->index[i+1]); j++)
114 if (newj > index[newi])
121 fprintf(stderr, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
122 name, block->nr, newi, block->index[block->nr], newj);
123 block->index = index;
127 static void reduce_blocka(atom_id invindex[], gmx_bool bKeep[], t_blocka *block,
131 int i, j, k, newi, newj;
133 snew(index, block->nr);
137 for (i = 0; (i < block->nr); i++)
139 for (j = block->index[i]; (j < block->index[i+1]); j++)
144 a[newj] = invindex[k];
148 if (newj > index[newi])
155 fprintf(stderr, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
156 name, block->nr, newi, block->nra, newj);
157 block->index = index;
163 static void reduce_rvec(int gnx, atom_id index[], rvec vv[])
169 for (i = 0; (i < gnx); i++)
171 copy_rvec(vv[index[i]], ptr[i]);
173 for (i = 0; (i < gnx); i++)
175 copy_rvec(ptr[i], vv[i]);
180 static void reduce_atom(int gnx, atom_id index[], t_atom atom[], char ***atomname,
181 int *nres, t_resinfo *resinfo)
190 snew(rinfo, atom[index[gnx-1]].resind+1);
191 for (i = 0; (i < gnx); i++)
193 ptr[i] = atom[index[i]];
194 aname[i] = atomname[index[i]];
197 for (i = 0; (i < gnx); i++)
200 atomname[i] = aname[i];
201 if ((i == 0) || (atom[i].resind != atom[i-1].resind))
204 rinfo[nr] = resinfo[atom[i].resind];
209 for (i = 0; (i < nr); i++)
211 resinfo[i] = rinfo[i];
220 static void reduce_ilist(atom_id invindex[], gmx_bool bKeep[],
221 t_ilist *il, int nratoms, const char *name)
231 for (i = 0; (i < il->nr); i += nratoms+1)
234 for (j = 1; (j <= nratoms); j++)
236 bB = bB && bKeep[il->iatoms[i+j]];
240 ia[newnr++] = il->iatoms[i];
241 for (j = 1; (j <= nratoms); j++)
243 ia[newnr++] = invindex[il->iatoms[i+j]];
247 fprintf(stderr, "Reduced ilist %8s from %6d to %6d entries\n",
248 name, il->nr/(nratoms+1),
252 for (i = 0; (i < newnr); i++)
254 il->iatoms[i] = ia[i];
261 static void reduce_topology_x(int gnx, atom_id index[],
262 gmx_mtop_t *mtop, rvec x[], rvec v[])
269 top = gmx_mtop_t_to_t_topology(mtop);
270 bKeep = bKeepIt(gnx, top.atoms.nr, index);
271 invindex = invind(gnx, top.atoms.nr, index);
273 reduce_block(bKeep, &(top.cgs), "cgs");
274 reduce_block(bKeep, &(top.mols), "mols");
275 reduce_blocka(invindex, bKeep, &(top.excls), "excls");
276 reduce_rvec(gnx, index, x);
277 reduce_rvec(gnx, index, v);
278 reduce_atom(gnx, index, top.atoms.atom, top.atoms.atomname,
279 &(top.atoms.nres), top.atoms.resinfo);
281 for (i = 0; (i < F_NRE); i++)
283 reduce_ilist(invindex, bKeep, &(top.idef.il[i]),
284 interaction_function[i].nratoms,
285 interaction_function[i].name);
291 snew(mtop->moltype, mtop->nmoltype);
292 mtop->moltype[0].name = mtop->name;
293 mtop->moltype[0].atoms = top.atoms;
294 for (i = 0; i < F_NRE; i++)
296 mtop->moltype[0].ilist[i] = top.idef.il[i];
298 mtop->moltype[0].atoms = top.atoms;
299 mtop->moltype[0].cgs = top.cgs;
300 mtop->moltype[0].excls = top.excls;
303 snew(mtop->molblock, mtop->nmolblock);
304 mtop->molblock[0].type = 0;
305 mtop->molblock[0].nmol = 1;
306 mtop->molblock[0].natoms_mol = top.atoms.nr;
307 mtop->molblock[0].nposres_xA = 0;
308 mtop->molblock[0].nposres_xB = 0;
310 mtop->natoms = top.atoms.nr;
313 static void zeroq(atom_id index[], gmx_mtop_t *mtop)
317 for (mt = 0; mt < mtop->nmoltype; mt++)
319 for (i = 0; (i < mtop->moltype[mt].atoms.nr); i++)
321 mtop->moltype[mt].atoms.atom[index[i]].q = 0;
322 mtop->moltype[mt].atoms.atom[index[i]].qB = 0;
327 int gmx_convert_tpr(int argc, char *argv[])
329 const char *desc[] = {
330 "[THISMODULE] can edit run input files in four ways.[PAR]",
331 "[BB]1.[bb] by modifying the number of steps in a run input file",
332 "with options [TT]-extend[tt], [TT]-until[tt] or [TT]-nsteps[tt]",
333 "(nsteps=-1 means unlimited number of steps)[PAR]",
334 "[BB]2.[bb] (OBSOLETE) by creating a run input file",
335 "for a continuation run when your simulation has crashed due to e.g.",
336 "a full disk, or by making a continuation run input file.",
337 "This option is obsolete, since mdrun now writes and reads",
339 "[BB]Note[bb] that a frame with coordinates and velocities is needed.",
340 "When pressure and/or Nose-Hoover temperature coupling is used",
341 "an energy file can be supplied to get an exact continuation",
342 "of the original run.[PAR]",
343 "[BB]3.[bb] by creating a [TT].tpx[tt] file for a subset of your original",
344 "tpx file, which is useful when you want to remove the solvent from",
345 "your [TT].tpx[tt] file, or when you want to make e.g. a pure C[GRK]alpha[grk] [TT].tpx[tt] file.",
346 "Note that you may need to use [TT]-nsteps -1[tt] (or similar) to get",
348 "[BB]WARNING: this [TT].tpx[tt] file is not fully functional[bb].[PAR]",
349 "[BB]4.[bb] by setting the charges of a specified group",
350 "to zero. This is useful when doing free energy estimates",
351 "using the LIE (Linear Interaction Energy) method."
354 const char *top_fn, *frame_fn;
356 ener_file_t fp_ener = NULL;
359 gmx_int64_t nsteps_req, run_step, frame;
360 double run_t, state_t;
361 gmx_bool bOK, bNsteps, bExtend, bUntil, bTime, bTraj;
362 gmx_bool bFrame, bUse, bSel, bNeedEner, bReadEner, bScanEner, bFepState;
365 t_inputrec *ir, *irnew = NULL;
368 rvec *newx = NULL, *newv = NULL, *tmpx, *tmpv;
372 atom_id *index = NULL;
374 gmx_enxnm_t *enm = NULL;
375 t_enxframe *fr_ener = NULL;
376 char buf[200], buf2[200];
379 { efTPX, NULL, NULL, ffREAD },
380 { efTRN, "-f", NULL, ffOPTRD },
381 { efEDR, "-e", NULL, ffOPTRD },
382 { efNDX, NULL, NULL, ffOPTRD },
383 { efTPX, "-o", "tpxout", ffWRITE }
385 #define NFILE asize(fnm)
387 /* Command line options */
388 static int nsteps_req_int = 0;
389 static real start_t = -1.0, extend_t = 0.0, until_t = 0.0;
390 static int init_fep_state = 0;
391 static gmx_bool bContinuation = TRUE, bZeroQ = FALSE, bVel = TRUE;
392 static t_pargs pa[] = {
393 { "-extend", FALSE, etREAL, {&extend_t},
394 "Extend runtime by this amount (ps)" },
395 { "-until", FALSE, etREAL, {&until_t},
396 "Extend runtime until this ending time (ps)" },
397 { "-nsteps", FALSE, etINT, {&nsteps_req_int},
398 "Change the number of steps" },
399 { "-time", FALSE, etREAL, {&start_t},
400 "Continue from frame at this time (ps) instead of the last frame" },
401 { "-zeroq", FALSE, etBOOL, {&bZeroQ},
402 "Set the charges of a group (from the index) to zero" },
403 { "-vel", FALSE, etBOOL, {&bVel},
404 "Require velocities from trajectory" },
405 { "-cont", FALSE, etBOOL, {&bContinuation},
406 "For exact continuation, the constraints should not be applied before the first step" },
407 { "-init_fep_state", FALSE, etINT, {&init_fep_state},
408 "fep state to initialize from" },
412 /* Parse the command line */
413 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
414 asize(desc), desc, 0, NULL, &oenv))
419 /* Convert int to gmx_int64_t */
420 nsteps_req = nsteps_req_int;
421 bNsteps = opt2parg_bSet("-nsteps", asize(pa), pa);
422 bExtend = opt2parg_bSet("-extend", asize(pa), pa);
423 bUntil = opt2parg_bSet("-until", asize(pa), pa);
424 bFepState = opt2parg_bSet("-init_fep_state", asize(pa), pa);
425 bTime = opt2parg_bSet("-time", asize(pa), pa);
426 bTraj = (opt2bSet("-f", NFILE, fnm) || bTime);
428 top_fn = ftp2fn(efTPX, NFILE, fnm);
429 fprintf(stderr, "Reading toplogy and stuff from %s\n", top_fn);
432 read_tpx_state(top_fn, ir, &state, NULL, &mtop);
433 run_step = ir->init_step;
434 run_t = ir->init_step*ir->delta_t + ir->init_t;
436 if (!EI_STATE_VELOCITY(ir->eI))
444 "NOTE: Reading the state from trajectory is an obsolete feature of gmx convert-tpr.\n"
445 " Continuation should be done by loading a checkpoint file with mdrun -cpi\n"
446 " This guarantees that all state variables are transferred.\n"
447 " gmx convert-tpr is now only useful for increasing nsteps,\n"
448 " but even that can often be avoided by using mdrun -maxh\n"
451 if (ir->bContinuation != bContinuation)
453 fprintf(stderr, "Modifying ir->bContinuation to %s\n",
454 bool_names[bContinuation]);
456 ir->bContinuation = bContinuation;
459 bNeedEner = (ir->epc == epcPARRINELLORAHMAN || ir->etc == etcNOSEHOOVER);
460 bReadEner = (bNeedEner && ftp2bSet(efEDR, NFILE, fnm));
461 bScanEner = (bReadEner && !bTime);
463 if (ir->epc != epcNO || EI_SD(ir->eI) || ir->eI == eiBD)
465 fprintf(stderr, "NOTE: The simulation uses pressure coupling and/or stochastic dynamics.\n"
466 "gmx convert-tpr can not provide binary identical continuation.\n"
467 "If you want that, supply a checkpoint file to mdrun\n\n");
470 if (EI_SD(ir->eI) || ir->eI == eiBD)
472 fprintf(stderr, "\nChanging ld-seed from %"GMX_PRId64 " ", ir->ld_seed);
473 ir->ld_seed = (gmx_int64_t)gmx_rng_make_seed();
474 fprintf(stderr, "to %"GMX_PRId64 "\n\n", ir->ld_seed);
477 frame_fn = ftp2fn(efTRN, NFILE, fnm);
479 if (fn2ftp(frame_fn) == efCPT)
484 "\nREADING STATE FROM CHECKPOINT %s...\n\n",
487 read_checkpoint_state(frame_fn, &sim_part,
488 &run_step, &run_t, &state);
493 "\nREADING COORDS, VELS AND BOX FROM TRAJECTORY %s...\n\n",
496 fp = open_trn(frame_fn, "r");
499 fp_ener = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
500 do_enxnms(fp_ener, &nre, &enm);
505 /* Now scan until the last set of x and v (step == 0)
506 * or the ones at step step.
512 bFrame = fread_trnheader(fp, &head, &bOK);
513 if (bOK && frame == 0)
515 if (mtop.natoms != head.natoms)
517 gmx_fatal(FARGS, "Number of atoms in Topology (%d) "
518 "is not the same as in Trajectory (%d)\n",
519 mtop.natoms, head.natoms);
521 snew(newx, head.natoms);
522 snew(newv, head.natoms);
524 bFrame = bFrame && bOK;
527 bOK = fread_htrn(fp, &head, newbox, newx, newv, NULL);
529 bFrame = bFrame && bOK;
532 (head.x_size) && (head.v_size || !bVel))
537 /* Read until the energy time is >= the trajectory time */
538 while (fr_ener->t < head.t && do_enx(fp_ener, fr_ener))
542 bUse = (fr_ener->t == head.t);
553 run_step = head.step;
554 state.fep_state = head.fep_state;
555 state.lambda[efptFEP] = head.lambda;
556 copy_mat(newbox, state.box);
561 sprintf(buf, "\r%s %s frame %s%s: step %s%s time %s",
562 "%s", "%s", "%6", GMX_PRId64, "%6", GMX_PRId64, " %8.3f");
564 bUse ? "Read " : "Skipped", ftp2ext(fn2ftp(frame_fn)),
565 frame, head.step, head.t);
567 if (bTime && (head.t >= start_t))
576 free_enxframe(fr_ener);
577 free_enxnms(nre, enm);
580 fprintf(stderr, "\n");
584 fprintf(stderr, "%s frame %s (step %s, time %g) is incomplete\n",
585 ftp2ext(fn2ftp(frame_fn)), gmx_step_str(frame-1, buf2),
586 gmx_step_str(head.step, buf), head.t);
588 fprintf(stderr, "\nUsing frame of step %s time %g\n",
589 gmx_step_str(run_step, buf), run_t);
595 get_enx_state(ftp2fn(efEDR, NFILE, fnm), run_t, &mtop.groups, ir, &state);
599 fprintf(stderr, "\nWARNING: The simulation uses %s temperature and/or %s pressure coupling,\n"
600 " the continuation will only be exact when an energy file is supplied\n\n",
601 ETCOUPLTYPE(etcNOSEHOOVER),
602 EPCOUPLTYPE(epcPARRINELLORAHMAN));
607 ir->fepvals->init_fep_state = init_fep_state;
614 fprintf(stderr, "Setting nsteps to %s\n", gmx_step_str(nsteps_req, buf));
615 ir->nsteps = nsteps_req;
619 /* Determine total number of steps remaining */
622 ir->nsteps = ir->nsteps - (run_step - ir->init_step) + (gmx_int64_t)(extend_t/ir->delta_t + 0.5);
623 printf("Extending remaining runtime of by %g ps (now %s steps)\n",
624 extend_t, gmx_step_str(ir->nsteps, buf));
628 printf("nsteps = %s, run_step = %s, current_t = %g, until = %g\n",
629 gmx_step_str(ir->nsteps, buf),
630 gmx_step_str(run_step, buf2),
632 ir->nsteps = (gmx_int64_t)((until_t - run_t)/ir->delta_t + 0.5);
633 printf("Extending remaining runtime until %g ps (now %s steps)\n",
634 until_t, gmx_step_str(ir->nsteps, buf));
638 ir->nsteps -= run_step - ir->init_step;
640 printf("%s steps (%g ps) remaining from first run.\n",
641 gmx_step_str(ir->nsteps, buf), ir->nsteps*ir->delta_t);
645 if (bNsteps || bZeroQ || (ir->nsteps > 0))
647 ir->init_step = run_step;
649 if (ftp2bSet(efNDX, NFILE, fnm) ||
650 !(bNsteps || bExtend || bUntil || bTraj))
652 atoms = gmx_mtop_global_atoms(&mtop);
653 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1,
654 &gnx, &index, &grpname);
657 bSel = (gnx != state.natoms);
658 for (i = 0; ((i < gnx) && (!bSel)); i++)
660 bSel = (i != index[i]);
669 fprintf(stderr, "Will write subset %s of original tpx containing %d "
670 "atoms\n", grpname, gnx);
671 reduce_topology_x(gnx, index, &mtop, state.x, state.v);
677 fprintf(stderr, "Zero-ing charges for group %s\n", grpname);
681 fprintf(stderr, "Will write full tpx file (no selection)\n");
685 state_t = ir->init_t + ir->init_step*ir->delta_t;
686 sprintf(buf, "Writing statusfile with starting step %s%s and length %s%s steps...\n", "%10", GMX_PRId64, "%10", GMX_PRId64);
687 fprintf(stderr, buf, ir->init_step, ir->nsteps);
688 fprintf(stderr, " time %10.3f and length %10.3f ps\n",
689 state_t, ir->nsteps*ir->delta_t);
690 write_tpx_state(opt2fn("-o", NFILE, fnm), ir, &state, &mtop);
694 printf("You've simulated long enough. Not writing tpr file\n");