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,2019, 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) \
68 "Your index file contains atomnumbers (e.g. %d)\nthat are larger than the number " \
69 "of atoms in the tpr file (%d)", \
72 static std::vector<bool> bKeepIt(int gnx, int natoms, int index[])
74 std::vector<bool> b(natoms);
76 for (int i = 0; (i < gnx); i++)
78 RANGECHK(index[i], natoms);
85 static std::vector<int> invind(int gnx, int natoms, int index[])
87 std::vector<int> inv(natoms);
89 for (int i = 0; (i < gnx); i++)
91 RANGECHK(index[i], natoms);
98 static gmx::ListOfLists<int> reduce_listoflists(gmx::ArrayRef<const int> invindex,
99 const std::vector<bool>& bKeep,
100 const gmx::ListOfLists<int>& src,
103 gmx::ListOfLists<int> lists;
105 std::vector<int> exclusionsForAtom;
106 for (gmx::index i = 0; i < src.ssize(); i++)
110 exclusionsForAtom.clear();
111 for (const int j : src[i])
115 exclusionsForAtom.push_back(invindex[j]);
118 lists.pushBack(exclusionsForAtom);
122 fprintf(stderr, "Reduced block %8s from %6zu to %6zu index-, %6d to %6d a-entries\n", name,
123 src.size(), lists.size(), src.numElements(), lists.numElements());
128 static void reduce_rvec(int gnx, const int index[], rvec vv[])
134 for (i = 0; (i < gnx); i++)
136 copy_rvec(vv[index[i]], ptr[i]);
138 for (i = 0; (i < gnx); i++)
140 copy_rvec(ptr[i], vv[i]);
145 static void reduce_atom(int gnx, const int index[], t_atom atom[], char*** atomname, int* nres, t_resinfo* resinfo)
154 snew(rinfo, atom[index[gnx - 1]].resind + 1);
155 for (i = 0; (i < gnx); i++)
157 ptr[i] = atom[index[i]];
158 aname[i] = atomname[index[i]];
161 for (i = 0; (i < gnx); i++)
164 atomname[i] = aname[i];
165 if ((i == 0) || (atom[i].resind != atom[i - 1].resind))
168 rinfo[nr] = resinfo[atom[i].resind];
173 for (i = 0; (i < nr); i++)
175 resinfo[i] = rinfo[i];
184 static void reduce_ilist(gmx::ArrayRef<const int> invindex,
185 const std::vector<bool>& bKeep,
198 for (i = 0; (i < il->nr); i += nratoms + 1)
201 for (j = 1; (j <= nratoms); j++)
203 bB = bB && bKeep[il->iatoms[i + j]];
207 ia[newnr++] = il->iatoms[i];
208 for (j = 1; (j <= nratoms); j++)
210 ia[newnr++] = invindex[il->iatoms[i + j]];
214 fprintf(stderr, "Reduced ilist %8s from %6d to %6d entries\n", name, il->nr / (nratoms + 1),
215 newnr / (nratoms + 1));
218 for (i = 0; (i < newnr); i++)
220 il->iatoms[i] = ia[i];
227 static void reduce_topology_x(int gnx, int index[], gmx_mtop_t* mtop, rvec x[], rvec v[])
230 gmx_mtop_generate_local_top(*mtop, &top, false);
231 t_atoms atoms = gmx_mtop_global_atoms(mtop);
233 const std::vector<bool> bKeep = bKeepIt(gnx, atoms.nr, index);
234 const std::vector<int> invindex = invind(gnx, atoms.nr, index);
236 reduce_rvec(gnx, index, x);
237 reduce_rvec(gnx, index, v);
238 reduce_atom(gnx, index, atoms.atom, atoms.atomname, &(atoms.nres), atoms.resinfo);
240 for (int i = 0; (i < F_NRE); i++)
242 reduce_ilist(invindex, bKeep, &(top.idef.il[i]), interaction_function[i].nratoms,
243 interaction_function[i].name);
248 mtop->moltype.resize(1);
249 mtop->moltype[0].name = mtop->name;
250 mtop->moltype[0].atoms = atoms;
251 mtop->moltype[0].excls = reduce_listoflists(invindex, bKeep, top.excls, "excls");
252 for (int i = 0; i < F_NRE; i++)
254 InteractionList& ilist = mtop->moltype[0].ilist[i];
255 ilist.iatoms.resize(top.idef.il[i].nr);
256 for (int j = 0; j < top.idef.il[i].nr; j++)
258 ilist.iatoms[j] = top.idef.il[i].iatoms[j];
262 mtop->molblock.resize(1);
263 mtop->molblock[0].type = 0;
264 mtop->molblock[0].nmol = 1;
266 mtop->natoms = atoms.nr;
269 static void zeroq(const int index[], gmx_mtop_t* mtop)
271 for (gmx_moltype_t& moltype : mtop->moltype)
273 for (int i = 0; i < moltype.atoms.nr; i++)
275 moltype.atoms.atom[index[i]].q = 0;
276 moltype.atoms.atom[index[i]].qB = 0;
281 int gmx_convert_tpr(int argc, char* argv[])
283 const char* desc[] = {
284 "[THISMODULE] can edit run input files in three ways.[PAR]",
285 "[BB]1.[bb] by modifying the number of steps in a run input file",
286 "with options [TT]-extend[tt], [TT]-until[tt] or [TT]-nsteps[tt]",
287 "(nsteps=-1 means unlimited number of steps)[PAR]",
288 "[BB]2.[bb] by creating a [REF].tpx[ref] file for a subset of your original",
289 "tpx file, which is useful when you want to remove the solvent from",
290 "your [REF].tpx[ref] file, or when you want to make e.g. a pure C[GRK]alpha[grk] ",
291 "[REF].tpx[ref] file.",
292 "Note that you may need to use [TT]-nsteps -1[tt] (or similar) to get",
294 "[BB]WARNING: this [REF].tpx[ref] file is not fully functional[bb].[PAR]",
295 "[BB]3.[bb] by setting the charges of a specified group",
296 "to zero. This is useful when doing free energy estimates",
297 "using the LIE (Linear Interaction Energy) method."
302 int64_t nsteps_req, run_step;
303 double run_t, state_t;
305 gmx_bool bNsteps, bExtend, bUntil;
311 int* index = nullptr;
312 char buf[200], buf2[200];
313 gmx_output_env_t* oenv;
314 t_filenm fnm[] = { { efTPR, nullptr, nullptr, ffREAD },
315 { efNDX, nullptr, nullptr, ffOPTRD },
316 { efTPR, "-o", "tprout", ffWRITE } };
317 #define NFILE asize(fnm)
319 /* Command line options */
320 static int nsteps_req_int = 0;
321 static real extend_t = 0.0, until_t = 0.0;
322 static gmx_bool bZeroQ = FALSE;
323 static t_pargs pa[] = {
324 { "-extend", FALSE, etREAL, { &extend_t }, "Extend runtime by this amount (ps)" },
325 { "-until", FALSE, etREAL, { &until_t }, "Extend runtime until this ending time (ps)" },
326 { "-nsteps", FALSE, etINT, { &nsteps_req_int }, "Change the number of steps" },
331 "Set the charges of a group (from the index) to zero" }
334 /* Parse the command line */
335 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
340 /* Convert int to int64_t */
341 nsteps_req = nsteps_req_int;
342 bNsteps = opt2parg_bSet("-nsteps", asize(pa), pa);
343 bExtend = opt2parg_bSet("-extend", asize(pa), pa);
344 bUntil = opt2parg_bSet("-until", asize(pa), pa);
346 top_fn = ftp2fn(efTPR, NFILE, fnm);
347 fprintf(stderr, "Reading toplogy and stuff from %s\n", top_fn);
349 t_inputrec irInstance;
350 t_inputrec* ir = &irInstance;
351 read_tpx_state(top_fn, ir, &state, &mtop);
352 run_step = ir->init_step;
353 run_t = ir->init_step * ir->delta_t + ir->init_t;
357 fprintf(stderr, "Setting nsteps to %s\n", gmx_step_str(nsteps_req, buf));
358 ir->nsteps = nsteps_req;
362 /* Determine total number of steps remaining */
365 ir->nsteps = ir->nsteps - (run_step - ir->init_step) + gmx::roundToInt64(extend_t / ir->delta_t);
366 printf("Extending remaining runtime of by %g ps (now %s steps)\n", extend_t,
367 gmx_step_str(ir->nsteps, buf));
371 printf("nsteps = %s, run_step = %s, current_t = %g, until = %g\n",
372 gmx_step_str(ir->nsteps, buf), gmx_step_str(run_step, buf2), run_t, until_t);
373 ir->nsteps = gmx::roundToInt64((until_t - run_t) / ir->delta_t);
374 printf("Extending remaining runtime until %g ps (now %s steps)\n", until_t,
375 gmx_step_str(ir->nsteps, buf));
379 ir->nsteps -= run_step - ir->init_step;
381 printf("%s steps (%g ps) remaining from first run.\n", gmx_step_str(ir->nsteps, buf),
382 ir->nsteps * ir->delta_t);
386 if (bNsteps || bZeroQ || (ir->nsteps > 0))
388 ir->init_step = run_step;
390 if (ftp2bSet(efNDX, NFILE, fnm) || !(bNsteps || bExtend || bUntil))
392 atoms = gmx_mtop_global_atoms(&mtop);
393 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
396 bSel = (gnx != state.natoms);
397 for (i = 0; ((i < gnx) && (!bSel)); i++)
399 bSel = (i != index[i]);
409 "Will write subset %s of original tpx containing %d "
412 reduce_topology_x(gnx, index, &mtop, state.x.rvec_array(), state.v.rvec_array());
418 fprintf(stderr, "Zero-ing charges for group %s\n", grpname);
422 fprintf(stderr, "Will write full tpx file (no selection)\n");
426 state_t = ir->init_t + ir->init_step * ir->delta_t;
427 sprintf(buf, "Writing statusfile with starting step %s%s and length %s%s steps...\n", "%10",
428 PRId64, "%10", PRId64);
429 fprintf(stderr, buf, ir->init_step, ir->nsteps);
430 fprintf(stderr, " time %10.3f and length %10.3f ps\n",
431 state_t, ir->nsteps * ir->delta_t);
432 write_tpx_state(opt2fn("-o", NFILE, fnm), ir, &state, &mtop);
436 printf("You've simulated long enough. Not writing tpr file\n");