Use ListOfLists in gmx_mtop_t and gmx_localtop_t
[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,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.
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)                                                                           \
66     if ((i) >= (n))                                                                              \
67     gmx_fatal(FARGS,                                                                             \
68               "Your index file contains atomnumbers (e.g. %d)\nthat are larger than the number " \
69               "of atoms in the tpr file (%d)",                                                   \
70               (i), (n))
71
72 static gmx_bool* bKeepIt(int gnx, int natoms, int index[])
73 {
74     gmx_bool* b;
75     int       i;
76
77     snew(b, natoms);
78     for (i = 0; (i < gnx); i++)
79     {
80         RANGECHK(index[i], natoms);
81         b[index[i]] = TRUE;
82     }
83
84     return b;
85 }
86
87 static int* invind(int gnx, int natoms, int index[])
88 {
89     int* inv;
90     int  i;
91
92     snew(inv, natoms);
93     for (i = 0; (i < gnx); i++)
94     {
95         RANGECHK(index[i], natoms);
96         inv[index[i]] = i;
97     }
98
99     return inv;
100 }
101
102 static void reduce_block(const gmx_bool bKeep[], t_block* block, const char* name)
103 {
104     int* index;
105     int  i, j, newi, newj;
106
107     snew(index, block->nr);
108
109     newi = newj = 0;
110     for (i = 0; (i < block->nr); i++)
111     {
112         for (j = block->index[i]; (j < block->index[i + 1]); j++)
113         {
114             if (bKeep[j])
115             {
116                 newj++;
117             }
118         }
119         if (newj > index[newi])
120         {
121             newi++;
122             index[newi] = newj;
123         }
124     }
125
126     fprintf(stderr, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n", name,
127             block->nr, newi, block->index[block->nr], newj);
128     block->index = index;
129     block->nr    = newi;
130 }
131
132 static gmx::ListOfLists<int>
133 reduce_blocka(const int invindex[], const gmx_bool bKeep[], const t_blocka& block, const char* name)
134 {
135     gmx::ListOfLists<int> lists;
136
137     std::vector<int> exclusionsForAtom;
138     for (int i = 0; i < block.nr; i++)
139     {
140         if (bKeep[i])
141         {
142             exclusionsForAtom.clear();
143             for (int j = block.index[i]; j < block.index[i + 1]; j++)
144             {
145                 const int k = block.a[j];
146                 if (bKeep[k])
147                 {
148                     exclusionsForAtom.push_back(invindex[k]);
149                 }
150             }
151             lists.pushBack(exclusionsForAtom);
152         }
153     }
154
155     fprintf(stderr, "Reduced block %8s from %6d to %6zu index-, %6d to %6d a-entries\n", name,
156             block.nr, lists.size(), block.nra, lists.numElements());
157
158     return lists;
159 }
160
161 static void reduce_rvec(int gnx, const int index[], rvec vv[])
162 {
163     rvec* ptr;
164     int   i;
165
166     snew(ptr, gnx);
167     for (i = 0; (i < gnx); i++)
168     {
169         copy_rvec(vv[index[i]], ptr[i]);
170     }
171     for (i = 0; (i < gnx); i++)
172     {
173         copy_rvec(ptr[i], vv[i]);
174     }
175     sfree(ptr);
176 }
177
178 static void reduce_atom(int gnx, const int index[], t_atom atom[], char*** atomname, int* nres, t_resinfo* resinfo)
179 {
180     t_atom*    ptr;
181     char***    aname;
182     t_resinfo* rinfo;
183     int        i, nr;
184
185     snew(ptr, gnx);
186     snew(aname, gnx);
187     snew(rinfo, atom[index[gnx - 1]].resind + 1);
188     for (i = 0; (i < gnx); i++)
189     {
190         ptr[i]   = atom[index[i]];
191         aname[i] = atomname[index[i]];
192     }
193     nr = -1;
194     for (i = 0; (i < gnx); i++)
195     {
196         atom[i]     = ptr[i];
197         atomname[i] = aname[i];
198         if ((i == 0) || (atom[i].resind != atom[i - 1].resind))
199         {
200             nr++;
201             rinfo[nr] = resinfo[atom[i].resind];
202         }
203         atom[i].resind = nr;
204     }
205     nr++;
206     for (i = 0; (i < nr); i++)
207     {
208         resinfo[i] = rinfo[i];
209     }
210     *nres = nr;
211
212     sfree(aname);
213     sfree(ptr);
214     sfree(rinfo);
215 }
216
217 static void reduce_ilist(const int invindex[], const gmx_bool bKeep[], t_ilist* il, int nratoms, const char* name)
218 {
219     t_iatom* ia;
220     int      i, j, newnr;
221     gmx_bool bB;
222
223     if (il->nr)
224     {
225         snew(ia, il->nr);
226         newnr = 0;
227         for (i = 0; (i < il->nr); i += nratoms + 1)
228         {
229             bB = TRUE;
230             for (j = 1; (j <= nratoms); j++)
231             {
232                 bB = bB && bKeep[il->iatoms[i + j]];
233             }
234             if (bB)
235             {
236                 ia[newnr++] = il->iatoms[i];
237                 for (j = 1; (j <= nratoms); j++)
238                 {
239                     ia[newnr++] = invindex[il->iatoms[i + j]];
240                 }
241             }
242         }
243         fprintf(stderr, "Reduced ilist %8s from %6d to %6d entries\n", name, il->nr / (nratoms + 1),
244                 newnr / (nratoms + 1));
245
246         il->nr = newnr;
247         for (i = 0; (i < newnr); i++)
248         {
249             il->iatoms[i] = ia[i];
250         }
251
252         sfree(ia);
253     }
254 }
255
256 static void reduce_topology_x(int gnx, int index[], gmx_mtop_t* mtop, rvec x[], rvec v[])
257 {
258     t_topology top;
259     gmx_bool*  bKeep;
260     int*       invindex;
261     int        i;
262
263     top      = gmx_mtop_t_to_t_topology(mtop, false);
264     bKeep    = bKeepIt(gnx, top.atoms.nr, index);
265     invindex = invind(gnx, top.atoms.nr, index);
266
267     reduce_block(bKeep, &(top.mols), "mols");
268     gmx::ListOfLists<int> excls = reduce_blocka(invindex, bKeep, top.excls, "excls");
269     reduce_rvec(gnx, index, x);
270     reduce_rvec(gnx, index, v);
271     reduce_atom(gnx, index, top.atoms.atom, top.atoms.atomname, &(top.atoms.nres), top.atoms.resinfo);
272
273     for (i = 0; (i < F_NRE); i++)
274     {
275         reduce_ilist(invindex, bKeep, &(top.idef.il[i]), interaction_function[i].nratoms,
276                      interaction_function[i].name);
277     }
278
279     top.atoms.nr = gnx;
280
281     mtop->moltype.resize(1);
282     mtop->moltype[0].name  = mtop->name;
283     mtop->moltype[0].atoms = top.atoms;
284     for (i = 0; i < F_NRE; i++)
285     {
286         InteractionList& ilist = mtop->moltype[0].ilist[i];
287         ilist.iatoms.resize(top.idef.il[i].nr);
288         for (int j = 0; j < top.idef.il[i].nr; j++)
289         {
290             ilist.iatoms[j] = top.idef.il[i].iatoms[j];
291         }
292     }
293     mtop->moltype[0].atoms = top.atoms;
294     mtop->moltype[0].excls = excls;
295
296     mtop->molblock.resize(1);
297     mtop->molblock[0].type = 0;
298     mtop->molblock[0].nmol = 1;
299
300     mtop->natoms = top.atoms.nr;
301 }
302
303 static void zeroq(const int index[], gmx_mtop_t* mtop)
304 {
305     for (gmx_moltype_t& moltype : mtop->moltype)
306     {
307         for (int i = 0; i < moltype.atoms.nr; i++)
308         {
309             moltype.atoms.atom[index[i]].q  = 0;
310             moltype.atoms.atom[index[i]].qB = 0;
311         }
312     }
313 }
314
315 int gmx_convert_tpr(int argc, char* argv[])
316 {
317     const char* desc[] = {
318         "[THISMODULE] can edit run input files in three ways.[PAR]",
319         "[BB]1.[bb] by modifying the number of steps in a run input file",
320         "with options [TT]-extend[tt], [TT]-until[tt] or [TT]-nsteps[tt]",
321         "(nsteps=-1 means unlimited number of steps)[PAR]",
322         "[BB]2.[bb] by creating a [REF].tpx[ref] file for a subset of your original",
323         "tpx file, which is useful when you want to remove the solvent from",
324         "your [REF].tpx[ref] file, or when you want to make e.g. a pure C[GRK]alpha[grk] ",
325         "[REF].tpx[ref] file.",
326         "Note that you may need to use [TT]-nsteps -1[tt] (or similar) to get",
327         "this to work.",
328         "[BB]WARNING: this [REF].tpx[ref] file is not fully functional[bb].[PAR]",
329         "[BB]3.[bb] by setting the charges of a specified group",
330         "to zero. This is useful when doing free energy estimates",
331         "using the LIE (Linear Interaction Energy) method."
332     };
333
334     const char*       top_fn;
335     int               i;
336     int64_t           nsteps_req, run_step;
337     double            run_t, state_t;
338     gmx_bool          bSel;
339     gmx_bool          bNsteps, bExtend, bUntil;
340     gmx_mtop_t        mtop;
341     t_atoms           atoms;
342     t_state           state;
343     int               gnx;
344     char*             grpname;
345     int*              index = nullptr;
346     char              buf[200], buf2[200];
347     gmx_output_env_t* oenv;
348     t_filenm          fnm[] = { { efTPR, nullptr, nullptr, ffREAD },
349                        { efNDX, nullptr, nullptr, ffOPTRD },
350                        { efTPR, "-o", "tprout", ffWRITE } };
351 #define NFILE asize(fnm)
352
353     /* Command line options */
354     static int      nsteps_req_int = 0;
355     static real     extend_t = 0.0, until_t = 0.0;
356     static gmx_bool bZeroQ = FALSE;
357     static t_pargs  pa[]   = {
358         { "-extend", FALSE, etREAL, { &extend_t }, "Extend runtime by this amount (ps)" },
359         { "-until", FALSE, etREAL, { &until_t }, "Extend runtime until this ending time (ps)" },
360         { "-nsteps", FALSE, etINT, { &nsteps_req_int }, "Change the number of steps" },
361         { "-zeroq",
362           FALSE,
363           etBOOL,
364           { &bZeroQ },
365           "Set the charges of a group (from the index) to zero" }
366     };
367
368     /* Parse the command line */
369     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
370     {
371         return 0;
372     }
373
374     /* Convert int to int64_t */
375     nsteps_req = nsteps_req_int;
376     bNsteps    = opt2parg_bSet("-nsteps", asize(pa), pa);
377     bExtend    = opt2parg_bSet("-extend", asize(pa), pa);
378     bUntil     = opt2parg_bSet("-until", asize(pa), pa);
379
380     top_fn = ftp2fn(efTPR, NFILE, fnm);
381     fprintf(stderr, "Reading toplogy and stuff from %s\n", top_fn);
382
383     t_inputrec  irInstance;
384     t_inputrec* ir = &irInstance;
385     read_tpx_state(top_fn, ir, &state, &mtop);
386     run_step = ir->init_step;
387     run_t    = ir->init_step * ir->delta_t + ir->init_t;
388
389     if (bNsteps)
390     {
391         fprintf(stderr, "Setting nsteps to %s\n", gmx_step_str(nsteps_req, buf));
392         ir->nsteps = nsteps_req;
393     }
394     else
395     {
396         /* Determine total number of steps remaining */
397         if (bExtend)
398         {
399             ir->nsteps = ir->nsteps - (run_step - ir->init_step) + gmx::roundToInt64(extend_t / ir->delta_t);
400             printf("Extending remaining runtime of by %g ps (now %s steps)\n", extend_t,
401                    gmx_step_str(ir->nsteps, buf));
402         }
403         else if (bUntil)
404         {
405             printf("nsteps = %s, run_step = %s, current_t = %g, until = %g\n",
406                    gmx_step_str(ir->nsteps, buf), gmx_step_str(run_step, buf2), run_t, until_t);
407             ir->nsteps = gmx::roundToInt64((until_t - run_t) / ir->delta_t);
408             printf("Extending remaining runtime until %g ps (now %s steps)\n", until_t,
409                    gmx_step_str(ir->nsteps, buf));
410         }
411         else
412         {
413             ir->nsteps -= run_step - ir->init_step;
414             /* Print message */
415             printf("%s steps (%g ps) remaining from first run.\n", gmx_step_str(ir->nsteps, buf),
416                    ir->nsteps * ir->delta_t);
417         }
418     }
419
420     if (bNsteps || bZeroQ || (ir->nsteps > 0))
421     {
422         ir->init_step = run_step;
423
424         if (ftp2bSet(efNDX, NFILE, fnm) || !(bNsteps || bExtend || bUntil))
425         {
426             atoms = gmx_mtop_global_atoms(&mtop);
427             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
428             if (!bZeroQ)
429             {
430                 bSel = (gnx != state.natoms);
431                 for (i = 0; ((i < gnx) && (!bSel)); i++)
432                 {
433                     bSel = (i != index[i]);
434                 }
435             }
436             else
437             {
438                 bSel = FALSE;
439             }
440             if (bSel)
441             {
442                 fprintf(stderr,
443                         "Will write subset %s of original tpx containing %d "
444                         "atoms\n",
445                         grpname, gnx);
446                 reduce_topology_x(gnx, index, &mtop, state.x.rvec_array(), state.v.rvec_array());
447                 state.natoms = gnx;
448             }
449             else if (bZeroQ)
450             {
451                 zeroq(index, &mtop);
452                 fprintf(stderr, "Zero-ing charges for group %s\n", grpname);
453             }
454             else
455             {
456                 fprintf(stderr, "Will write full tpx file (no selection)\n");
457             }
458         }
459
460         state_t = ir->init_t + ir->init_step * ir->delta_t;
461         sprintf(buf, "Writing statusfile with starting step %s%s and length %s%s steps...\n", "%10",
462                 PRId64, "%10", PRId64);
463         fprintf(stderr, buf, ir->init_step, ir->nsteps);
464         fprintf(stderr, "                                 time %10.3f and length %10.3f ps\n",
465                 state_t, ir->nsteps * ir->delta_t);
466         write_tpx_state(opt2fn("-o", NFILE, fnm), ir, &state, &mtop);
467     }
468     else
469     {
470         printf("You've simulated long enough. Not writing tpr file\n");
471     }
472
473     return 0;
474 }