d49f581b7e21a6bd36e0ea0dfe8fe70df80d9760
[alexxy/gromacs.git] / src / gromacs / tools / convert_tpr.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "convert_tpr.h"
40
41 #include <cmath>
42
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"
64
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))
66
67 static gmx_bool *bKeepIt(int gnx, int natoms, int index[])
68 {
69     gmx_bool *b;
70     int       i;
71
72     snew(b, natoms);
73     for (i = 0; (i < gnx); i++)
74     {
75         RANGECHK(index[i], natoms);
76         b[index[i]] = TRUE;
77     }
78
79     return b;
80 }
81
82 static int *invind(int gnx, int natoms, int index[])
83 {
84     int     *inv;
85     int      i;
86
87     snew(inv, natoms);
88     for (i = 0; (i < gnx); i++)
89     {
90         RANGECHK(index[i], natoms);
91         inv[index[i]] = i;
92     }
93
94     return inv;
95 }
96
97 static void reduce_block(const gmx_bool bKeep[], t_block *block,
98                          const char *name)
99 {
100     int     *index;
101     int      i, j, newi, newj;
102
103     snew(index, block->nr);
104
105     newi = newj = 0;
106     for (i = 0; (i < block->nr); i++)
107     {
108         for (j = block->index[i]; (j < block->index[i+1]); j++)
109         {
110             if (bKeep[j])
111             {
112                 newj++;
113             }
114         }
115         if (newj > index[newi])
116         {
117             newi++;
118             index[newi] = newj;
119         }
120     }
121
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;
125     block->nr    = newi;
126 }
127
128 static void reduce_blocka(const int invindex[], const gmx_bool bKeep[], t_blocka *block,
129                           const char *name)
130 {
131     int     *index, *a;
132     int      i, j, k, newi, newj;
133
134     snew(index, block->nr);
135     snew(a, block->nra);
136
137     newi = newj = 0;
138     for (i = 0; (i < block->nr); i++)
139     {
140         for (j = block->index[i]; (j < block->index[i+1]); j++)
141         {
142             k = block->a[j];
143             if (bKeep[k])
144             {
145                 a[newj] = invindex[k];
146                 newj++;
147             }
148         }
149         if (newj > index[newi])
150         {
151             newi++;
152             index[newi] = newj;
153         }
154     }
155
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;
159     block->a     = a;
160     block->nr    = newi;
161     block->nra   = newj;
162 }
163
164 static void reduce_rvec(int gnx, const int index[], rvec vv[])
165 {
166     rvec *ptr;
167     int   i;
168
169     snew(ptr, gnx);
170     for (i = 0; (i < gnx); i++)
171     {
172         copy_rvec(vv[index[i]], ptr[i]);
173     }
174     for (i = 0; (i < gnx); i++)
175     {
176         copy_rvec(ptr[i], vv[i]);
177     }
178     sfree(ptr);
179 }
180
181 static void reduce_atom(int gnx, const int index[], t_atom atom[], char ***atomname,
182                         int *nres, t_resinfo *resinfo)
183 {
184     t_atom    *ptr;
185     char    ***aname;
186     t_resinfo *rinfo;
187     int        i, nr;
188
189     snew(ptr, gnx);
190     snew(aname, gnx);
191     snew(rinfo, atom[index[gnx-1]].resind+1);
192     for (i = 0; (i < gnx); i++)
193     {
194         ptr[i]   = atom[index[i]];
195         aname[i] = atomname[index[i]];
196     }
197     nr = -1;
198     for (i = 0; (i < gnx); i++)
199     {
200         atom[i]     = ptr[i];
201         atomname[i] = aname[i];
202         if ((i == 0) || (atom[i].resind != atom[i-1].resind))
203         {
204             nr++;
205             rinfo[nr] = resinfo[atom[i].resind];
206         }
207         atom[i].resind = nr;
208     }
209     nr++;
210     for (i = 0; (i < nr); i++)
211     {
212         resinfo[i] = rinfo[i];
213     }
214     *nres = nr;
215
216     sfree(aname);
217     sfree(ptr);
218     sfree(rinfo);
219 }
220
221 static void reduce_ilist(const int invindex[], const gmx_bool bKeep[],
222                          t_ilist *il, int nratoms, const char *name)
223 {
224     t_iatom *ia;
225     int      i, j, newnr;
226     gmx_bool bB;
227
228     if (il->nr)
229     {
230         snew(ia, il->nr);
231         newnr = 0;
232         for (i = 0; (i < il->nr); i += nratoms+1)
233         {
234             bB = TRUE;
235             for (j = 1; (j <= nratoms); j++)
236             {
237                 bB = bB && bKeep[il->iatoms[i+j]];
238             }
239             if (bB)
240             {
241                 ia[newnr++] = il->iatoms[i];
242                 for (j = 1; (j <= nratoms); j++)
243                 {
244                     ia[newnr++] = invindex[il->iatoms[i+j]];
245                 }
246             }
247         }
248         fprintf(stderr, "Reduced ilist %8s from %6d to %6d entries\n",
249                 name, il->nr/(nratoms+1),
250                 newnr/(nratoms+1));
251
252         il->nr = newnr;
253         for (i = 0; (i < newnr); i++)
254         {
255             il->iatoms[i] = ia[i];
256         }
257
258         sfree(ia);
259     }
260 }
261
262 static void reduce_topology_x(int gnx, int index[],
263                               gmx_mtop_t *mtop, rvec x[], rvec v[])
264 {
265     t_topology   top;
266     gmx_bool    *bKeep;
267     int         *invindex;
268     int          i;
269
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);
273
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);
281
282     for (i = 0; (i < F_NRE); i++)
283     {
284         reduce_ilist(invindex, bKeep, &(top.idef.il[i]),
285                      interaction_function[i].nratoms,
286                      interaction_function[i].name);
287     }
288
289     top.atoms.nr = gnx;
290
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++)
295     {
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++)
299         {
300             ilist.iatoms[j] = top.idef.il[i].iatoms[j];
301         }
302     }
303     mtop->moltype[0].atoms = top.atoms;
304     mtop->moltype[0].cgs   = top.cgs;
305     mtop->moltype[0].excls = top.excls;
306
307     mtop->molblock.resize(1);
308     mtop->molblock[0].type = 0;
309     mtop->molblock[0].nmol = 1;
310
311     mtop->natoms           = top.atoms.nr;
312 }
313
314 static void zeroq(const int index[], gmx_mtop_t *mtop)
315 {
316     for (gmx_moltype_t &moltype : mtop->moltype)
317     {
318         for (int i = 0; i < moltype.atoms.nr; i++)
319         {
320             moltype.atoms.atom[index[i]].q  = 0;
321             moltype.atoms.atom[index[i]].qB = 0;
322         }
323     }
324 }
325
326 int gmx_convert_tpr(int argc, char *argv[])
327 {
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",
337         "this to work.",
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."
342     };
343
344     const char       *top_fn;
345     int               i;
346     int64_t           nsteps_req, run_step;
347     double            run_t, state_t;
348     gmx_bool          bSel;
349     gmx_bool          bNsteps, bExtend, bUntil;
350     gmx_mtop_t        mtop;
351     t_atoms           atoms;
352     t_state           state;
353     int               gnx;
354     char             *grpname;
355     int              *index = nullptr;
356     char              buf[200], buf2[200];
357     gmx_output_env_t *oenv;
358     t_filenm          fnm[] = {
359         { efTPR, nullptr,  nullptr,    ffREAD  },
360         { efNDX, nullptr,  nullptr,    ffOPTRD },
361         { efTPR, "-o",  "tprout", ffWRITE }
362     };
363 #define NFILE asize(fnm)
364
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" }
378     };
379
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))
383     {
384         return 0;
385     }
386
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);
392
393     top_fn = ftp2fn(efTPR, NFILE, fnm);
394     fprintf(stderr, "Reading toplogy and stuff from %s\n", top_fn);
395
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;
401
402     if (bNsteps)
403     {
404         fprintf(stderr, "Setting nsteps to %s\n", gmx_step_str(nsteps_req, buf));
405         ir->nsteps = nsteps_req;
406     }
407     else
408     {
409         /* Determine total number of steps remaining */
410         if (bExtend)
411         {
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));
415         }
416         else if (bUntil)
417         {
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),
421                    run_t, until_t);
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));
425         }
426         else
427         {
428             ir->nsteps -= run_step - ir->init_step;
429             /* Print message */
430             printf("%s steps (%g ps) remaining from first run.\n",
431                    gmx_step_str(ir->nsteps, buf), ir->nsteps*ir->delta_t);
432         }
433     }
434
435     if (bNsteps || bZeroQ || (ir->nsteps > 0))
436     {
437         ir->init_step = run_step;
438
439         if (ftp2bSet(efNDX, NFILE, fnm) ||
440             !(bNsteps || bExtend || bUntil))
441         {
442             atoms = gmx_mtop_global_atoms(&mtop);
443             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1,
444                       &gnx, &index, &grpname);
445             if (!bZeroQ)
446             {
447                 bSel = (gnx != state.natoms);
448                 for (i = 0; ((i < gnx) && (!bSel)); i++)
449                 {
450                     bSel = (i != index[i]);
451                 }
452             }
453             else
454             {
455                 bSel = FALSE;
456             }
457             if (bSel)
458             {
459                 fprintf(stderr, "Will write subset %s of original tpx containing %d "
460                         "atoms\n", grpname, gnx);
461                 reduce_topology_x(gnx, index, &mtop, state.x.rvec_array(), state.v.rvec_array());
462                 state.natoms = gnx;
463             }
464             else if (bZeroQ)
465             {
466                 zeroq(index, &mtop);
467                 fprintf(stderr, "Zero-ing charges for group %s\n", grpname);
468             }
469             else
470             {
471                 fprintf(stderr, "Will write full tpx file (no selection)\n");
472             }
473         }
474
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);
481     }
482     else
483     {
484         printf("You've simulated long enough. Not writing tpr file\n");
485     }
486
487     return 0;
488 }