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/utility/fatalerror.h"
44 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/gmxpreprocess/readir.h"
49 #include "gromacs/commandline/pargs.h"
50 #include "gromacs/math/vec.h"
51 #include "mtop_util.h"
52 #include "checkpoint.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trnio.h"
55 #include "gromacs/fileio/enxio.h"
56 #include "gromacs/utility/futil.h"
57 #include "gromacs/random/random.h"
59 #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))
61 static gmx_bool *bKeepIt(int gnx, int natoms, atom_id index[])
67 for (i = 0; (i < gnx); i++)
69 RANGECHK(index[i], natoms);
76 static atom_id *invind(int gnx, int natoms, atom_id index[])
82 for (i = 0; (i < gnx); i++)
84 RANGECHK(index[i], natoms);
91 static void reduce_block(gmx_bool bKeep[], t_block *block,
97 snew(index, block->nr);
100 for (i = 0; (i < block->nr); i++)
102 for (j = block->index[i]; (j < block->index[i+1]); j++)
109 if (newj > index[newi])
116 fprintf(stderr, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
117 name, block->nr, newi, block->index[block->nr], newj);
118 block->index = index;
122 static void reduce_blocka(atom_id invindex[], gmx_bool bKeep[], t_blocka *block,
126 int i, j, k, newi, newj;
128 snew(index, block->nr);
132 for (i = 0; (i < block->nr); i++)
134 for (j = block->index[i]; (j < block->index[i+1]); j++)
139 a[newj] = invindex[k];
143 if (newj > index[newi])
150 fprintf(stderr, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
151 name, block->nr, newi, block->nra, newj);
152 block->index = index;
158 static void reduce_rvec(int gnx, atom_id index[], rvec vv[])
164 for (i = 0; (i < gnx); i++)
166 copy_rvec(vv[index[i]], ptr[i]);
168 for (i = 0; (i < gnx); i++)
170 copy_rvec(ptr[i], vv[i]);
175 static void reduce_atom(int gnx, atom_id index[], t_atom atom[], char ***atomname,
176 int *nres, t_resinfo *resinfo)
185 snew(rinfo, atom[index[gnx-1]].resind+1);
186 for (i = 0; (i < gnx); i++)
188 ptr[i] = atom[index[i]];
189 aname[i] = atomname[index[i]];
192 for (i = 0; (i < gnx); i++)
195 atomname[i] = aname[i];
196 if ((i == 0) || (atom[i].resind != atom[i-1].resind))
199 rinfo[nr] = resinfo[atom[i].resind];
204 for (i = 0; (i < nr); i++)
206 resinfo[i] = rinfo[i];
215 static void reduce_ilist(atom_id invindex[], gmx_bool bKeep[],
216 t_ilist *il, int nratoms, const char *name)
226 for (i = 0; (i < il->nr); i += nratoms+1)
229 for (j = 1; (j <= nratoms); j++)
231 bB = bB && bKeep[il->iatoms[i+j]];
235 ia[newnr++] = il->iatoms[i];
236 for (j = 1; (j <= nratoms); j++)
238 ia[newnr++] = invindex[il->iatoms[i+j]];
242 fprintf(stderr, "Reduced ilist %8s from %6d to %6d entries\n",
243 name, il->nr/(nratoms+1),
247 for (i = 0; (i < newnr); i++)
249 il->iatoms[i] = ia[i];
256 static void reduce_topology_x(int gnx, atom_id index[],
257 gmx_mtop_t *mtop, rvec x[], rvec v[])
264 top = gmx_mtop_t_to_t_topology(mtop);
265 bKeep = bKeepIt(gnx, top.atoms.nr, index);
266 invindex = invind(gnx, top.atoms.nr, index);
268 reduce_block(bKeep, &(top.cgs), "cgs");
269 reduce_block(bKeep, &(top.mols), "mols");
270 reduce_blocka(invindex, bKeep, &(top.excls), "excls");
271 reduce_rvec(gnx, index, x);
272 reduce_rvec(gnx, index, v);
273 reduce_atom(gnx, index, top.atoms.atom, top.atoms.atomname,
274 &(top.atoms.nres), top.atoms.resinfo);
276 for (i = 0; (i < F_NRE); i++)
278 reduce_ilist(invindex, bKeep, &(top.idef.il[i]),
279 interaction_function[i].nratoms,
280 interaction_function[i].name);
286 snew(mtop->moltype, mtop->nmoltype);
287 mtop->moltype[0].name = mtop->name;
288 mtop->moltype[0].atoms = top.atoms;
289 for (i = 0; i < F_NRE; i++)
291 mtop->moltype[0].ilist[i] = top.idef.il[i];
293 mtop->moltype[0].atoms = top.atoms;
294 mtop->moltype[0].cgs = top.cgs;
295 mtop->moltype[0].excls = top.excls;
298 snew(mtop->molblock, mtop->nmolblock);
299 mtop->molblock[0].type = 0;
300 mtop->molblock[0].nmol = 1;
301 mtop->molblock[0].natoms_mol = top.atoms.nr;
302 mtop->molblock[0].nposres_xA = 0;
303 mtop->molblock[0].nposres_xB = 0;
305 mtop->natoms = top.atoms.nr;
308 static void zeroq(atom_id index[], gmx_mtop_t *mtop)
312 for (mt = 0; mt < mtop->nmoltype; mt++)
314 for (i = 0; (i < mtop->moltype[mt].atoms.nr); i++)
316 mtop->moltype[mt].atoms.atom[index[i]].q = 0;
317 mtop->moltype[mt].atoms.atom[index[i]].qB = 0;
322 int gmx_convert_tpr(int argc, char *argv[])
324 const char *desc[] = {
325 "[THISMODULE] can edit run input files in four ways.[PAR]",
326 "[BB]1.[bb] by modifying the number of steps in a run input file",
327 "with options [TT]-extend[tt], [TT]-until[tt] or [TT]-nsteps[tt]",
328 "(nsteps=-1 means unlimited number of steps)[PAR]",
329 "[BB]2.[bb] (OBSOLETE) by creating a run input file",
330 "for a continuation run when your simulation has crashed due to e.g.",
331 "a full disk, or by making a continuation run input file.",
332 "This option is obsolete, since mdrun now writes and reads",
334 "[BB]Note[bb] that a frame with coordinates and velocities is needed.",
335 "When pressure and/or Nose-Hoover temperature coupling is used",
336 "an energy file can be supplied to get an exact continuation",
337 "of the original run.[PAR]",
338 "[BB]3.[bb] by creating a [TT].tpx[tt] file for a subset of your original",
339 "tpx file, which is useful when you want to remove the solvent from",
340 "your [TT].tpx[tt] file, or when you want to make e.g. a pure C[GRK]alpha[grk] [TT].tpx[tt] file.",
341 "Note that you may need to use [TT]-nsteps -1[tt] (or similar) to get",
343 "[BB]WARNING: this [TT].tpx[tt] file is not fully functional[bb].[PAR]",
344 "[BB]4.[bb] by setting the charges of a specified group",
345 "to zero. This is useful when doing free energy estimates",
346 "using the LIE (Linear Interaction Energy) method."
349 const char *top_fn, *frame_fn;
351 ener_file_t fp_ener = NULL;
354 gmx_int64_t nsteps_req, run_step, frame;
355 double run_t, state_t;
356 gmx_bool bOK, bNsteps, bExtend, bUntil, bTime, bTraj;
357 gmx_bool bFrame, bUse, bSel, bNeedEner, bReadEner, bScanEner, bFepState;
360 t_inputrec *ir, *irnew = NULL;
363 rvec *newx = NULL, *newv = NULL, *tmpx, *tmpv;
367 atom_id *index = NULL;
369 gmx_enxnm_t *enm = NULL;
370 t_enxframe *fr_ener = NULL;
371 char buf[200], buf2[200];
374 { efTPX, NULL, NULL, ffREAD },
375 { efTRN, "-f", NULL, ffOPTRD },
376 { efEDR, "-e", NULL, ffOPTRD },
377 { efNDX, NULL, NULL, ffOPTRD },
378 { efTPX, "-o", "tpxout", ffWRITE }
380 #define NFILE asize(fnm)
382 /* Command line options */
383 static int nsteps_req_int = 0;
384 static real start_t = -1.0, extend_t = 0.0, until_t = 0.0;
385 static int init_fep_state = 0;
386 static gmx_bool bContinuation = TRUE, bZeroQ = FALSE, bVel = TRUE;
387 static t_pargs pa[] = {
388 { "-extend", FALSE, etREAL, {&extend_t},
389 "Extend runtime by this amount (ps)" },
390 { "-until", FALSE, etREAL, {&until_t},
391 "Extend runtime until this ending time (ps)" },
392 { "-nsteps", FALSE, etINT, {&nsteps_req_int},
393 "Change the number of steps" },
394 { "-time", FALSE, etREAL, {&start_t},
395 "Continue from frame at this time (ps) instead of the last frame" },
396 { "-zeroq", FALSE, etBOOL, {&bZeroQ},
397 "Set the charges of a group (from the index) to zero" },
398 { "-vel", FALSE, etBOOL, {&bVel},
399 "Require velocities from trajectory" },
400 { "-cont", FALSE, etBOOL, {&bContinuation},
401 "For exact continuation, the constraints should not be applied before the first step" },
402 { "-init_fep_state", FALSE, etINT, {&init_fep_state},
403 "fep state to initialize from" },
407 /* Parse the command line */
408 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
409 asize(desc), desc, 0, NULL, &oenv))
414 /* Convert int to gmx_int64_t */
415 nsteps_req = nsteps_req_int;
416 bNsteps = opt2parg_bSet("-nsteps", asize(pa), pa);
417 bExtend = opt2parg_bSet("-extend", asize(pa), pa);
418 bUntil = opt2parg_bSet("-until", asize(pa), pa);
419 bFepState = opt2parg_bSet("-init_fep_state", asize(pa), pa);
420 bTime = opt2parg_bSet("-time", asize(pa), pa);
421 bTraj = (opt2bSet("-f", NFILE, fnm) || bTime);
423 top_fn = ftp2fn(efTPX, NFILE, fnm);
424 fprintf(stderr, "Reading toplogy and stuff from %s\n", top_fn);
427 read_tpx_state(top_fn, ir, &state, NULL, &mtop);
428 run_step = ir->init_step;
429 run_t = ir->init_step*ir->delta_t + ir->init_t;
431 if (!EI_STATE_VELOCITY(ir->eI))
439 "NOTE: Reading the state from trajectory is an obsolete feature of gmx convert-tpr.\n"
440 " Continuation should be done by loading a checkpoint file with mdrun -cpi\n"
441 " This guarantees that all state variables are transferred.\n"
442 " gmx convert-tpr is now only useful for increasing nsteps,\n"
443 " but even that can often be avoided by using mdrun -maxh\n"
446 if (ir->bContinuation != bContinuation)
448 fprintf(stderr, "Modifying ir->bContinuation to %s\n",
449 bool_names[bContinuation]);
451 ir->bContinuation = bContinuation;
454 bNeedEner = (ir->epc == epcPARRINELLORAHMAN || ir->etc == etcNOSEHOOVER);
455 bReadEner = (bNeedEner && ftp2bSet(efEDR, NFILE, fnm));
456 bScanEner = (bReadEner && !bTime);
458 if (ir->epc != epcNO || EI_SD(ir->eI) || ir->eI == eiBD)
460 fprintf(stderr, "NOTE: The simulation uses pressure coupling and/or stochastic dynamics.\n"
461 "gmx convert-tpr can not provide binary identical continuation.\n"
462 "If you want that, supply a checkpoint file to mdrun\n\n");
465 if (EI_SD(ir->eI) || ir->eI == eiBD)
467 fprintf(stderr, "\nChanging ld-seed from %"GMX_PRId64 " ", ir->ld_seed);
468 ir->ld_seed = (gmx_int64_t)gmx_rng_make_seed();
469 fprintf(stderr, "to %"GMX_PRId64 "\n\n", ir->ld_seed);
472 frame_fn = ftp2fn(efTRN, NFILE, fnm);
474 if (fn2ftp(frame_fn) == efCPT)
479 "\nREADING STATE FROM CHECKPOINT %s...\n\n",
482 read_checkpoint_state(frame_fn, &sim_part,
483 &run_step, &run_t, &state);
488 "\nREADING COORDS, VELS AND BOX FROM TRAJECTORY %s...\n\n",
491 fp = open_trn(frame_fn, "r");
494 fp_ener = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
495 do_enxnms(fp_ener, &nre, &enm);
500 /* Now scan until the last set of x and v (step == 0)
501 * or the ones at step step.
507 bFrame = fread_trnheader(fp, &head, &bOK);
508 if (bOK && frame == 0)
510 if (mtop.natoms != head.natoms)
512 gmx_fatal(FARGS, "Number of atoms in Topology (%d) "
513 "is not the same as in Trajectory (%d)\n",
514 mtop.natoms, head.natoms);
516 snew(newx, head.natoms);
517 snew(newv, head.natoms);
519 bFrame = bFrame && bOK;
522 bOK = fread_htrn(fp, &head, newbox, newx, newv, NULL);
524 bFrame = bFrame && bOK;
527 (head.x_size) && (head.v_size || !bVel))
532 /* Read until the energy time is >= the trajectory time */
533 while (fr_ener->t < head.t && do_enx(fp_ener, fr_ener))
537 bUse = (fr_ener->t == head.t);
548 run_step = head.step;
549 state.fep_state = head.fep_state;
550 state.lambda[efptFEP] = head.lambda;
551 copy_mat(newbox, state.box);
556 sprintf(buf, "\r%s %s frame %s%s: step %s%s time %s",
557 "%s", "%s", "%6", GMX_PRId64, "%6", GMX_PRId64, " %8.3f");
559 bUse ? "Read " : "Skipped", ftp2ext(fn2ftp(frame_fn)),
560 frame, head.step, head.t);
562 if (bTime && (head.t >= start_t))
571 free_enxframe(fr_ener);
572 free_enxnms(nre, enm);
575 fprintf(stderr, "\n");
579 fprintf(stderr, "%s frame %s (step %s, time %g) is incomplete\n",
580 ftp2ext(fn2ftp(frame_fn)), gmx_step_str(frame-1, buf2),
581 gmx_step_str(head.step, buf), head.t);
583 fprintf(stderr, "\nUsing frame of step %s time %g\n",
584 gmx_step_str(run_step, buf), run_t);
590 get_enx_state(ftp2fn(efEDR, NFILE, fnm), run_t, &mtop.groups, ir, &state);
594 fprintf(stderr, "\nWARNING: The simulation uses %s temperature and/or %s pressure coupling,\n"
595 " the continuation will only be exact when an energy file is supplied\n\n",
596 ETCOUPLTYPE(etcNOSEHOOVER),
597 EPCOUPLTYPE(epcPARRINELLORAHMAN));
602 ir->fepvals->init_fep_state = init_fep_state;
609 fprintf(stderr, "Setting nsteps to %s\n", gmx_step_str(nsteps_req, buf));
610 ir->nsteps = nsteps_req;
614 /* Determine total number of steps remaining */
617 ir->nsteps = ir->nsteps - (run_step - ir->init_step) + (gmx_int64_t)(extend_t/ir->delta_t + 0.5);
618 printf("Extending remaining runtime of by %g ps (now %s steps)\n",
619 extend_t, gmx_step_str(ir->nsteps, buf));
623 printf("nsteps = %s, run_step = %s, current_t = %g, until = %g\n",
624 gmx_step_str(ir->nsteps, buf),
625 gmx_step_str(run_step, buf2),
627 ir->nsteps = (gmx_int64_t)((until_t - run_t)/ir->delta_t + 0.5);
628 printf("Extending remaining runtime until %g ps (now %s steps)\n",
629 until_t, gmx_step_str(ir->nsteps, buf));
633 ir->nsteps -= run_step - ir->init_step;
635 printf("%s steps (%g ps) remaining from first run.\n",
636 gmx_step_str(ir->nsteps, buf), ir->nsteps*ir->delta_t);
640 if (bNsteps || bZeroQ || (ir->nsteps > 0))
642 ir->init_step = run_step;
644 if (ftp2bSet(efNDX, NFILE, fnm) ||
645 !(bNsteps || bExtend || bUntil || bTraj))
647 atoms = gmx_mtop_global_atoms(&mtop);
648 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1,
649 &gnx, &index, &grpname);
652 bSel = (gnx != state.natoms);
653 for (i = 0; ((i < gnx) && (!bSel)); i++)
655 bSel = (i != index[i]);
664 fprintf(stderr, "Will write subset %s of original tpx containing %d "
665 "atoms\n", grpname, gnx);
666 reduce_topology_x(gnx, index, &mtop, state.x, state.v);
672 fprintf(stderr, "Zero-ing charges for group %s\n", grpname);
676 fprintf(stderr, "Will write full tpx file (no selection)\n");
680 state_t = ir->init_t + ir->init_step*ir->delta_t;
681 sprintf(buf, "Writing statusfile with starting step %s%s and length %s%s steps...\n", "%10", GMX_PRId64, "%10", GMX_PRId64);
682 fprintf(stderr, buf, ir->init_step, ir->nsteps);
683 fprintf(stderr, " time %10.3f and length %10.3f ps\n",
684 state_t, ir->nsteps*ir->delta_t);
685 write_tpx_state(opt2fn("-o", NFILE, fnm), ir, &state, &mtop);
689 printf("You've simulated long enough. Not writing tpr file\n");