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 "gmx_fatal.h"
45 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/gmxpreprocess/readir.h"
50 #include "gromacs/commandline/pargs.h"
52 #include "mtop_util.h"
53 #include "checkpoint.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trnio.h"
56 #include "gromacs/fileio/enxio.h"
57 #include "gromacs/fileio/futil.h"
58 #include "gromacs/random/random.h"
60 #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))
62 static gmx_bool *bKeepIt(int gnx, int natoms, atom_id index[])
68 for (i = 0; (i < gnx); i++)
70 RANGECHK(index[i], natoms);
77 static atom_id *invind(int gnx, int natoms, atom_id index[])
83 for (i = 0; (i < gnx); i++)
85 RANGECHK(index[i], natoms);
92 static void reduce_block(gmx_bool bKeep[], t_block *block,
98 snew(index, block->nr);
101 for (i = 0; (i < block->nr); i++)
103 for (j = block->index[i]; (j < block->index[i+1]); j++)
110 if (newj > index[newi])
117 fprintf(stderr, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
118 name, block->nr, newi, block->index[block->nr], newj);
119 block->index = index;
123 static void reduce_blocka(atom_id invindex[], gmx_bool bKeep[], t_blocka *block,
127 int i, j, k, newi, newj;
129 snew(index, block->nr);
133 for (i = 0; (i < block->nr); i++)
135 for (j = block->index[i]; (j < block->index[i+1]); j++)
140 a[newj] = invindex[k];
144 if (newj > index[newi])
151 fprintf(stderr, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
152 name, block->nr, newi, block->nra, newj);
153 block->index = index;
159 static void reduce_rvec(int gnx, atom_id index[], rvec vv[])
165 for (i = 0; (i < gnx); i++)
167 copy_rvec(vv[index[i]], ptr[i]);
169 for (i = 0; (i < gnx); i++)
171 copy_rvec(ptr[i], vv[i]);
176 static void reduce_atom(int gnx, atom_id index[], t_atom atom[], char ***atomname,
177 int *nres, t_resinfo *resinfo)
186 snew(rinfo, atom[index[gnx-1]].resind+1);
187 for (i = 0; (i < gnx); i++)
189 ptr[i] = atom[index[i]];
190 aname[i] = atomname[index[i]];
193 for (i = 0; (i < gnx); i++)
196 atomname[i] = aname[i];
197 if ((i == 0) || (atom[i].resind != atom[i-1].resind))
200 rinfo[nr] = resinfo[atom[i].resind];
205 for (i = 0; (i < nr); i++)
207 resinfo[i] = rinfo[i];
216 static void reduce_ilist(atom_id invindex[], gmx_bool bKeep[],
217 t_ilist *il, int nratoms, const char *name)
227 for (i = 0; (i < il->nr); i += nratoms+1)
230 for (j = 1; (j <= nratoms); j++)
232 bB = bB && bKeep[il->iatoms[i+j]];
236 ia[newnr++] = il->iatoms[i];
237 for (j = 1; (j <= nratoms); j++)
239 ia[newnr++] = invindex[il->iatoms[i+j]];
243 fprintf(stderr, "Reduced ilist %8s from %6d to %6d entries\n",
244 name, il->nr/(nratoms+1),
248 for (i = 0; (i < newnr); i++)
250 il->iatoms[i] = ia[i];
257 static void reduce_topology_x(int gnx, atom_id index[],
258 gmx_mtop_t *mtop, rvec x[], rvec v[])
265 top = gmx_mtop_t_to_t_topology(mtop);
266 bKeep = bKeepIt(gnx, top.atoms.nr, index);
267 invindex = invind(gnx, top.atoms.nr, index);
269 reduce_block(bKeep, &(top.cgs), "cgs");
270 reduce_block(bKeep, &(top.mols), "mols");
271 reduce_blocka(invindex, bKeep, &(top.excls), "excls");
272 reduce_rvec(gnx, index, x);
273 reduce_rvec(gnx, index, v);
274 reduce_atom(gnx, index, top.atoms.atom, top.atoms.atomname,
275 &(top.atoms.nres), top.atoms.resinfo);
277 for (i = 0; (i < F_NRE); i++)
279 reduce_ilist(invindex, bKeep, &(top.idef.il[i]),
280 interaction_function[i].nratoms,
281 interaction_function[i].name);
287 snew(mtop->moltype, mtop->nmoltype);
288 mtop->moltype[0].name = mtop->name;
289 mtop->moltype[0].atoms = top.atoms;
290 for (i = 0; i < F_NRE; i++)
292 mtop->moltype[0].ilist[i] = top.idef.il[i];
294 mtop->moltype[0].atoms = top.atoms;
295 mtop->moltype[0].cgs = top.cgs;
296 mtop->moltype[0].excls = top.excls;
299 snew(mtop->molblock, mtop->nmolblock);
300 mtop->molblock[0].type = 0;
301 mtop->molblock[0].nmol = 1;
302 mtop->molblock[0].natoms_mol = top.atoms.nr;
303 mtop->molblock[0].nposres_xA = 0;
304 mtop->molblock[0].nposres_xB = 0;
306 mtop->natoms = top.atoms.nr;
309 static void zeroq(atom_id index[], gmx_mtop_t *mtop)
313 for (mt = 0; mt < mtop->nmoltype; mt++)
315 for (i = 0; (i < mtop->moltype[mt].atoms.nr); i++)
317 mtop->moltype[mt].atoms.atom[index[i]].q = 0;
318 mtop->moltype[mt].atoms.atom[index[i]].qB = 0;
323 int gmx_convert_tpr(int argc, char *argv[])
325 const char *desc[] = {
326 "[THISMODULE] can edit run input files in four ways.[PAR]",
327 "[BB]1.[bb] by modifying the number of steps in a run input file",
328 "with options [TT]-extend[tt], [TT]-until[tt] or [TT]-nsteps[tt]",
329 "(nsteps=-1 means unlimited number of steps)[PAR]",
330 "[BB]2.[bb] (OBSOLETE) by creating a run input file",
331 "for a continuation run when your simulation has crashed due to e.g.",
332 "a full disk, or by making a continuation run input file.",
333 "This option is obsolete, since mdrun now writes and reads",
335 "[BB]Note[bb] that a frame with coordinates and velocities is needed.",
336 "When pressure and/or Nose-Hoover temperature coupling is used",
337 "an energy file can be supplied to get an exact continuation",
338 "of the original run.[PAR]",
339 "[BB]3.[bb] by creating a [TT].tpx[tt] file for a subset of your original",
340 "tpx file, which is useful when you want to remove the solvent from",
341 "your [TT].tpx[tt] file, or when you want to make e.g. a pure C[GRK]alpha[grk] [TT].tpx[tt] file.",
342 "Note that you may need to use [TT]-nsteps -1[tt] (or similar) to get",
344 "[BB]WARNING: this [TT].tpx[tt] file is not fully functional[bb].[PAR]",
345 "[BB]4.[bb] by setting the charges of a specified group",
346 "to zero. This is useful when doing free energy estimates",
347 "using the LIE (Linear Interaction Energy) method."
350 const char *top_fn, *frame_fn;
352 ener_file_t fp_ener = NULL;
355 gmx_int64_t nsteps_req, run_step, frame;
356 double run_t, state_t;
357 gmx_bool bOK, bNsteps, bExtend, bUntil, bTime, bTraj;
358 gmx_bool bFrame, bUse, bSel, bNeedEner, bReadEner, bScanEner, bFepState;
361 t_inputrec *ir, *irnew = NULL;
364 rvec *newx = NULL, *newv = NULL, *tmpx, *tmpv;
368 atom_id *index = NULL;
370 gmx_enxnm_t *enm = NULL;
371 t_enxframe *fr_ener = NULL;
372 char buf[200], buf2[200];
375 { efTPX, NULL, NULL, ffREAD },
376 { efTRN, "-f", NULL, ffOPTRD },
377 { efEDR, "-e", NULL, ffOPTRD },
378 { efNDX, NULL, NULL, ffOPTRD },
379 { efTPX, "-o", "tpxout", ffWRITE }
381 #define NFILE asize(fnm)
383 /* Command line options */
384 static int nsteps_req_int = 0;
385 static real start_t = -1.0, extend_t = 0.0, until_t = 0.0;
386 static int init_fep_state = 0;
387 static gmx_bool bContinuation = TRUE, bZeroQ = FALSE, bVel = TRUE;
388 static t_pargs pa[] = {
389 { "-extend", FALSE, etREAL, {&extend_t},
390 "Extend runtime by this amount (ps)" },
391 { "-until", FALSE, etREAL, {&until_t},
392 "Extend runtime until this ending time (ps)" },
393 { "-nsteps", FALSE, etINT, {&nsteps_req_int},
394 "Change the number of steps" },
395 { "-time", FALSE, etREAL, {&start_t},
396 "Continue from frame at this time (ps) instead of the last frame" },
397 { "-zeroq", FALSE, etBOOL, {&bZeroQ},
398 "Set the charges of a group (from the index) to zero" },
399 { "-vel", FALSE, etBOOL, {&bVel},
400 "Require velocities from trajectory" },
401 { "-cont", FALSE, etBOOL, {&bContinuation},
402 "For exact continuation, the constraints should not be applied before the first step" },
403 { "-init_fep_state", FALSE, etINT, {&init_fep_state},
404 "fep state to initialize from" },
408 /* Parse the command line */
409 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
410 asize(desc), desc, 0, NULL, &oenv))
415 /* Convert int to gmx_int64_t */
416 nsteps_req = nsteps_req_int;
417 bNsteps = opt2parg_bSet("-nsteps", asize(pa), pa);
418 bExtend = opt2parg_bSet("-extend", asize(pa), pa);
419 bUntil = opt2parg_bSet("-until", asize(pa), pa);
420 bFepState = opt2parg_bSet("-init_fep_state", asize(pa), pa);
421 bTime = opt2parg_bSet("-time", asize(pa), pa);
422 bTraj = (opt2bSet("-f", NFILE, fnm) || bTime);
424 top_fn = ftp2fn(efTPX, NFILE, fnm);
425 fprintf(stderr, "Reading toplogy and stuff from %s\n", top_fn);
428 read_tpx_state(top_fn, ir, &state, NULL, &mtop);
429 run_step = ir->init_step;
430 run_t = ir->init_step*ir->delta_t + ir->init_t;
432 if (!EI_STATE_VELOCITY(ir->eI))
440 "NOTE: Reading the state from trajectory is an obsolete feature of gmx convert-tpr.\n"
441 " Continuation should be done by loading a checkpoint file with mdrun -cpi\n"
442 " This guarantees that all state variables are transferred.\n"
443 " gmx convert-tpr is now only useful for increasing nsteps,\n"
444 " but even that can often be avoided by using mdrun -maxh\n"
447 if (ir->bContinuation != bContinuation)
449 fprintf(stderr, "Modifying ir->bContinuation to %s\n",
450 bool_names[bContinuation]);
452 ir->bContinuation = bContinuation;
455 bNeedEner = (ir->epc == epcPARRINELLORAHMAN || ir->etc == etcNOSEHOOVER);
456 bReadEner = (bNeedEner && ftp2bSet(efEDR, NFILE, fnm));
457 bScanEner = (bReadEner && !bTime);
459 if (ir->epc != epcNO || EI_SD(ir->eI) || ir->eI == eiBD)
461 fprintf(stderr, "NOTE: The simulation uses pressure coupling and/or stochastic dynamics.\n"
462 "gmx convert-tpr can not provide binary identical continuation.\n"
463 "If you want that, supply a checkpoint file to mdrun\n\n");
466 if (EI_SD(ir->eI) || ir->eI == eiBD)
468 fprintf(stderr, "\nChanging ld-seed from %"GMX_PRId64 " ", ir->ld_seed);
469 ir->ld_seed = (gmx_int64_t)gmx_rng_make_seed();
470 fprintf(stderr, "to %"GMX_PRId64 "\n\n", ir->ld_seed);
473 frame_fn = ftp2fn(efTRN, NFILE, fnm);
475 if (fn2ftp(frame_fn) == efCPT)
480 "\nREADING STATE FROM CHECKPOINT %s...\n\n",
483 read_checkpoint_state(frame_fn, &sim_part,
484 &run_step, &run_t, &state);
489 "\nREADING COORDS, VELS AND BOX FROM TRAJECTORY %s...\n\n",
492 fp = open_trn(frame_fn, "r");
495 fp_ener = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
496 do_enxnms(fp_ener, &nre, &enm);
501 /* Now scan until the last set of x and v (step == 0)
502 * or the ones at step step.
508 bFrame = fread_trnheader(fp, &head, &bOK);
509 if (bOK && frame == 0)
511 if (mtop.natoms != head.natoms)
513 gmx_fatal(FARGS, "Number of atoms in Topology (%d) "
514 "is not the same as in Trajectory (%d)\n",
515 mtop.natoms, head.natoms);
517 snew(newx, head.natoms);
518 snew(newv, head.natoms);
520 bFrame = bFrame && bOK;
523 bOK = fread_htrn(fp, &head, newbox, newx, newv, NULL);
525 bFrame = bFrame && bOK;
528 (head.x_size) && (head.v_size || !bVel))
533 /* Read until the energy time is >= the trajectory time */
534 while (fr_ener->t < head.t && do_enx(fp_ener, fr_ener))
538 bUse = (fr_ener->t == head.t);
549 run_step = head.step;
550 state.fep_state = head.fep_state;
551 state.lambda[efptFEP] = head.lambda;
552 copy_mat(newbox, state.box);
557 sprintf(buf, "\r%s %s frame %s%s: step %s%s time %s",
558 "%s", "%s", "%6", GMX_PRId64, "%6", GMX_PRId64, " %8.3f");
560 bUse ? "Read " : "Skipped", ftp2ext(fn2ftp(frame_fn)),
561 frame, head.step, head.t);
563 if (bTime && (head.t >= start_t))
572 free_enxframe(fr_ener);
573 free_enxnms(nre, enm);
576 fprintf(stderr, "\n");
580 fprintf(stderr, "%s frame %s (step %s, time %g) is incomplete\n",
581 ftp2ext(fn2ftp(frame_fn)), gmx_step_str(frame-1, buf2),
582 gmx_step_str(head.step, buf), head.t);
584 fprintf(stderr, "\nUsing frame of step %s time %g\n",
585 gmx_step_str(run_step, buf), run_t);
591 get_enx_state(ftp2fn(efEDR, NFILE, fnm), run_t, &mtop.groups, ir, &state);
595 fprintf(stderr, "\nWARNING: The simulation uses %s temperature and/or %s pressure coupling,\n"
596 " the continuation will only be exact when an energy file is supplied\n\n",
597 ETCOUPLTYPE(etcNOSEHOOVER),
598 EPCOUPLTYPE(epcPARRINELLORAHMAN));
603 ir->fepvals->init_fep_state = init_fep_state;
610 fprintf(stderr, "Setting nsteps to %s\n", gmx_step_str(nsteps_req, buf));
611 ir->nsteps = nsteps_req;
615 /* Determine total number of steps remaining */
618 ir->nsteps = ir->nsteps - (run_step - ir->init_step) + (gmx_int64_t)(extend_t/ir->delta_t + 0.5);
619 printf("Extending remaining runtime of by %g ps (now %s steps)\n",
620 extend_t, gmx_step_str(ir->nsteps, buf));
624 printf("nsteps = %s, run_step = %s, current_t = %g, until = %g\n",
625 gmx_step_str(ir->nsteps, buf),
626 gmx_step_str(run_step, buf2),
628 ir->nsteps = (gmx_int64_t)((until_t - run_t)/ir->delta_t + 0.5);
629 printf("Extending remaining runtime until %g ps (now %s steps)\n",
630 until_t, gmx_step_str(ir->nsteps, buf));
634 ir->nsteps -= run_step - ir->init_step;
636 printf("%s steps (%g ps) remaining from first run.\n",
637 gmx_step_str(ir->nsteps, buf), ir->nsteps*ir->delta_t);
641 if (bNsteps || bZeroQ || (ir->nsteps > 0))
643 ir->init_step = run_step;
645 if (ftp2bSet(efNDX, NFILE, fnm) ||
646 !(bNsteps || bExtend || bUntil || bTraj))
648 atoms = gmx_mtop_global_atoms(&mtop);
649 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1,
650 &gnx, &index, &grpname);
653 bSel = (gnx != state.natoms);
654 for (i = 0; ((i < gnx) && (!bSel)); i++)
656 bSel = (i != index[i]);
665 fprintf(stderr, "Will write subset %s of original tpx containing %d "
666 "atoms\n", grpname, gnx);
667 reduce_topology_x(gnx, index, &mtop, state.x, state.v);
673 fprintf(stderr, "Zero-ing charges for group %s\n", grpname);
677 fprintf(stderr, "Will write full tpx file (no selection)\n");
681 state_t = ir->init_t + ir->init_step*ir->delta_t;
682 sprintf(buf, "Writing statusfile with starting step %s%s and length %s%s steps...\n", "%10", GMX_PRId64, "%10", GMX_PRId64);
683 fprintf(stderr, buf, ir->init_step, ir->nsteps);
684 fprintf(stderr, " time %10.3f and length %10.3f ps\n",
685 state_t, ir->nsteps*ir->delta_t);
686 write_tpx_state(opt2fn("-o", NFILE, fnm), ir, &state, &mtop);
690 printf("You've simulated long enough. Not writing tpr file\n");