Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / tools / convert_tpr.c
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, 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 <math.h>
40
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/fileio/enxio.h"
43 #include "gromacs/fileio/tpxio.h"
44 #include "gromacs/fileio/trnio.h"
45 #include "gromacs/gmxpreprocess/readir.h"
46 #include "gromacs/legacyheaders/checkpoint.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/legacyheaders/names.h"
49 #include "gromacs/legacyheaders/types/inputrec.h"
50 #include "gromacs/legacyheaders/types/simple.h"
51 #include "gromacs/legacyheaders/types/state.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/random/random.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
60
61 #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
63 static gmx_bool *bKeepIt(int gnx, int natoms, atom_id index[])
64 {
65     gmx_bool *b;
66     int       i;
67
68     snew(b, natoms);
69     for (i = 0; (i < gnx); i++)
70     {
71         RANGECHK(index[i], natoms);
72         b[index[i]] = TRUE;
73     }
74
75     return b;
76 }
77
78 static atom_id *invind(int gnx, int natoms, atom_id index[])
79 {
80     atom_id *inv;
81     int      i;
82
83     snew(inv, natoms);
84     for (i = 0; (i < gnx); i++)
85     {
86         RANGECHK(index[i], natoms);
87         inv[index[i]] = i;
88     }
89
90     return inv;
91 }
92
93 static void reduce_block(gmx_bool bKeep[], t_block *block,
94                          const char *name)
95 {
96     atom_id *index;
97     int      i, j, newi, newj;
98
99     snew(index, block->nr);
100
101     newi = newj = 0;
102     for (i = 0; (i < block->nr); i++)
103     {
104         for (j = block->index[i]; (j < block->index[i+1]); j++)
105         {
106             if (bKeep[j])
107             {
108                 newj++;
109             }
110         }
111         if (newj > index[newi])
112         {
113             newi++;
114             index[newi] = newj;
115         }
116     }
117
118     fprintf(stderr, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
119             name, block->nr, newi, block->index[block->nr], newj);
120     block->index = index;
121     block->nr    = newi;
122 }
123
124 static void reduce_blocka(atom_id invindex[], gmx_bool bKeep[], t_blocka *block,
125                           const char *name)
126 {
127     atom_id *index, *a;
128     int      i, j, k, newi, newj;
129
130     snew(index, block->nr);
131     snew(a, block->nra);
132
133     newi = newj = 0;
134     for (i = 0; (i < block->nr); i++)
135     {
136         for (j = block->index[i]; (j < block->index[i+1]); j++)
137         {
138             k = block->a[j];
139             if (bKeep[k])
140             {
141                 a[newj] = invindex[k];
142                 newj++;
143             }
144         }
145         if (newj > index[newi])
146         {
147             newi++;
148             index[newi] = newj;
149         }
150     }
151
152     fprintf(stderr, "Reduced block %8s from %6d to %6d index-, %6d to %6d a-entries\n",
153             name, block->nr, newi, block->nra, newj);
154     block->index = index;
155     block->a     = a;
156     block->nr    = newi;
157     block->nra   = newj;
158 }
159
160 static void reduce_rvec(int gnx, atom_id index[], rvec vv[])
161 {
162     rvec *ptr;
163     int   i;
164
165     snew(ptr, gnx);
166     for (i = 0; (i < gnx); i++)
167     {
168         copy_rvec(vv[index[i]], ptr[i]);
169     }
170     for (i = 0; (i < gnx); i++)
171     {
172         copy_rvec(ptr[i], vv[i]);
173     }
174     sfree(ptr);
175 }
176
177 static void reduce_atom(int gnx, atom_id index[], t_atom atom[], char ***atomname,
178                         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(atom_id invindex[], gmx_bool bKeep[],
218                          t_ilist *il, int nratoms, const char *name)
219 {
220     t_iatom *ia;
221     int      i, j, newnr;
222     gmx_bool bB;
223
224     if (il->nr)
225     {
226         snew(ia, il->nr);
227         newnr = 0;
228         for (i = 0; (i < il->nr); i += nratoms+1)
229         {
230             bB = TRUE;
231             for (j = 1; (j <= nratoms); j++)
232             {
233                 bB = bB && bKeep[il->iatoms[i+j]];
234             }
235             if (bB)
236             {
237                 ia[newnr++] = il->iatoms[i];
238                 for (j = 1; (j <= nratoms); j++)
239                 {
240                     ia[newnr++] = invindex[il->iatoms[i+j]];
241                 }
242             }
243         }
244         fprintf(stderr, "Reduced ilist %8s from %6d to %6d entries\n",
245                 name, il->nr/(nratoms+1),
246                 newnr/(nratoms+1));
247
248         il->nr = newnr;
249         for (i = 0; (i < newnr); i++)
250         {
251             il->iatoms[i] = ia[i];
252         }
253
254         sfree(ia);
255     }
256 }
257
258 static void reduce_topology_x(int gnx, atom_id index[],
259                               gmx_mtop_t *mtop, rvec x[], rvec v[])
260 {
261     t_topology   top;
262     gmx_bool    *bKeep;
263     atom_id     *invindex;
264     int          i;
265
266     top      = gmx_mtop_t_to_t_topology(mtop);
267     bKeep    = bKeepIt(gnx, top.atoms.nr, index);
268     invindex = invind(gnx, top.atoms.nr, index);
269
270     reduce_block(bKeep, &(top.cgs), "cgs");
271     reduce_block(bKeep, &(top.mols), "mols");
272     reduce_blocka(invindex, bKeep, &(top.excls), "excls");
273     reduce_rvec(gnx, index, x);
274     reduce_rvec(gnx, index, v);
275     reduce_atom(gnx, index, top.atoms.atom, top.atoms.atomname,
276                 &(top.atoms.nres), top.atoms.resinfo);
277
278     for (i = 0; (i < F_NRE); i++)
279     {
280         reduce_ilist(invindex, bKeep, &(top.idef.il[i]),
281                      interaction_function[i].nratoms,
282                      interaction_function[i].name);
283     }
284
285     top.atoms.nr = gnx;
286
287     mtop->nmoltype = 1;
288     snew(mtop->moltype, mtop->nmoltype);
289     mtop->moltype[0].name  = mtop->name;
290     mtop->moltype[0].atoms = top.atoms;
291     for (i = 0; i < F_NRE; i++)
292     {
293         mtop->moltype[0].ilist[i] = top.idef.il[i];
294     }
295     mtop->moltype[0].atoms = top.atoms;
296     mtop->moltype[0].cgs   = top.cgs;
297     mtop->moltype[0].excls = top.excls;
298
299     mtop->nmolblock = 1;
300     snew(mtop->molblock, mtop->nmolblock);
301     mtop->molblock[0].type       = 0;
302     mtop->molblock[0].nmol       = 1;
303     mtop->molblock[0].natoms_mol = top.atoms.nr;
304     mtop->molblock[0].nposres_xA = 0;
305     mtop->molblock[0].nposres_xB = 0;
306
307     mtop->natoms                 = top.atoms.nr;
308 }
309
310 static void zeroq(atom_id index[], gmx_mtop_t *mtop)
311 {
312     int mt, i;
313
314     for (mt = 0; mt < mtop->nmoltype; mt++)
315     {
316         for (i = 0; (i < mtop->moltype[mt].atoms.nr); i++)
317         {
318             mtop->moltype[mt].atoms.atom[index[i]].q  = 0;
319             mtop->moltype[mt].atoms.atom[index[i]].qB = 0;
320         }
321     }
322 }
323
324 int gmx_convert_tpr(int argc, char *argv[])
325 {
326     const char       *desc[] = {
327         "[THISMODULE] can edit run input files in four ways.[PAR]",
328         "[BB]1.[bb] by modifying the number of steps in a run input file",
329         "with options [TT]-extend[tt], [TT]-until[tt] or [TT]-nsteps[tt]",
330         "(nsteps=-1 means unlimited number of steps)[PAR]",
331         "[BB]2.[bb] (OBSOLETE) by creating a run input file",
332         "for a continuation run when your simulation has crashed due to e.g.",
333         "a full disk, or by making a continuation run input file.",
334         "This option is obsolete, since mdrun now writes and reads",
335         "checkpoint files.",
336         "[BB]Note[bb] that a frame with coordinates and velocities is needed.",
337         "When pressure and/or Nose-Hoover temperature coupling is used",
338         "an energy file can be supplied to get an exact continuation",
339         "of the original run.[PAR]",
340         "[BB]3.[bb] by creating a [TT].tpx[tt] file for a subset of your original",
341         "tpx file, which is useful when you want to remove the solvent from",
342         "your [TT].tpx[tt] file, or when you want to make e.g. a pure C[GRK]alpha[grk] [TT].tpx[tt] file.",
343         "Note that you may need to use [TT]-nsteps -1[tt] (or similar) to get",
344         "this to work.",
345         "[BB]WARNING: this [TT].tpx[tt] file is not fully functional[bb].[PAR]",
346         "[BB]4.[bb] by setting the charges of a specified group",
347         "to zero. This is useful when doing free energy estimates",
348         "using the LIE (Linear Interaction Energy) method."
349     };
350
351     const char       *top_fn, *frame_fn;
352     t_fileio         *fp;
353     ener_file_t       fp_ener = NULL;
354     t_trnheader       head;
355     int               i;
356     gmx_int64_t       nsteps_req, run_step, frame;
357     double            run_t, state_t;
358     gmx_bool          bOK, bNsteps, bExtend, bUntil, bTime, bTraj;
359     gmx_bool          bFrame, bUse, bSel, bNeedEner, bReadEner, bScanEner, bFepState;
360     gmx_mtop_t        mtop;
361     t_atoms           atoms;
362     t_inputrec       *ir, *irnew = NULL;
363     t_gromppopts     *gopts;
364     t_state           state;
365     rvec             *newx = NULL, *newv = NULL, *tmpx, *tmpv;
366     matrix            newbox;
367     int               gnx;
368     char             *grpname;
369     atom_id          *index = NULL;
370     int               nre;
371     gmx_enxnm_t      *enm     = NULL;
372     t_enxframe       *fr_ener = NULL;
373     char              buf[200], buf2[200];
374     output_env_t      oenv;
375     t_filenm          fnm[] = {
376         { efTPR, NULL,  NULL,    ffREAD  },
377         { efTRN, "-f",  NULL,    ffOPTRD },
378         { efEDR, "-e",  NULL,    ffOPTRD },
379         { efNDX, NULL,  NULL,    ffOPTRD },
380         { efTPR, "-o",  "tprout", ffWRITE }
381     };
382 #define NFILE asize(fnm)
383
384     /* Command line options */
385     static int      nsteps_req_int = 0;
386     static real     start_t        = -1.0, extend_t = 0.0, until_t = 0.0;
387     static int      init_fep_state = 0;
388     static gmx_bool bContinuation  = TRUE, bZeroQ = FALSE, bVel = TRUE;
389     static t_pargs  pa[]           = {
390         { "-extend",        FALSE, etREAL, {&extend_t},
391           "Extend runtime by this amount (ps)" },
392         { "-until",         FALSE, etREAL, {&until_t},
393           "Extend runtime until this ending time (ps)" },
394         { "-nsteps",        FALSE, etINT,  {&nsteps_req_int},
395           "Change the number of steps" },
396         { "-time",          FALSE, etREAL, {&start_t},
397           "Continue from frame at this time (ps) instead of the last frame" },
398         { "-zeroq",         FALSE, etBOOL, {&bZeroQ},
399           "Set the charges of a group (from the index) to zero" },
400         { "-vel",           FALSE, etBOOL, {&bVel},
401           "Require velocities from trajectory" },
402         { "-cont",          FALSE, etBOOL, {&bContinuation},
403           "For exact continuation, the constraints should not be applied before the first step" },
404         { "-init_fep_state", FALSE, etINT, {&init_fep_state},
405           "fep state to initialize from" },
406     };
407     int             nerror = 0;
408
409     /* Parse the command line */
410     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
411                            asize(desc), desc, 0, NULL, &oenv))
412     {
413         return 0;
414     }
415
416     /* Convert int to gmx_int64_t */
417     nsteps_req = nsteps_req_int;
418     bNsteps    = opt2parg_bSet("-nsteps", asize(pa), pa);
419     bExtend    = opt2parg_bSet("-extend", asize(pa), pa);
420     bUntil     = opt2parg_bSet("-until", asize(pa), pa);
421     bFepState  = opt2parg_bSet("-init_fep_state", asize(pa), pa);
422     bTime      = opt2parg_bSet("-time", asize(pa), pa);
423     bTraj      = (opt2bSet("-f", NFILE, fnm) || bTime);
424
425     top_fn = ftp2fn(efTPR, NFILE, fnm);
426     fprintf(stderr, "Reading toplogy and stuff from %s\n", top_fn);
427
428     snew(ir, 1);
429     read_tpx_state(top_fn, ir, &state, NULL, &mtop);
430     run_step = ir->init_step;
431     run_t    = ir->init_step*ir->delta_t + ir->init_t;
432
433     if (!EI_STATE_VELOCITY(ir->eI))
434     {
435         bVel = FALSE;
436     }
437
438     if (bTraj)
439     {
440         fprintf(stderr, "\n"
441                 "NOTE: Reading the state from trajectory is an obsolete feature of gmx convert-tpr.\n"
442                 "      Continuation should be done by loading a checkpoint file with mdrun -cpi\n"
443                 "      This guarantees that all state variables are transferred.\n"
444                 "      gmx convert-tpr is now only useful for increasing nsteps,\n"
445                 "      but even that can often be avoided by using mdrun -maxh\n"
446                 "\n");
447
448         if (ir->bContinuation != bContinuation)
449         {
450             fprintf(stderr, "Modifying ir->bContinuation to %s\n",
451                     bool_names[bContinuation]);
452         }
453         ir->bContinuation = bContinuation;
454
455
456         bNeedEner = (ir->epc == epcPARRINELLORAHMAN || ir->etc == etcNOSEHOOVER);
457         bReadEner = (bNeedEner && ftp2bSet(efEDR, NFILE, fnm));
458         bScanEner = (bReadEner && !bTime);
459
460         if (ir->epc != epcNO || EI_SD(ir->eI) || ir->eI == eiBD)
461         {
462             fprintf(stderr, "NOTE: The simulation uses pressure coupling and/or stochastic dynamics.\n"
463                     "gmx convert-tpr can not provide binary identical continuation.\n"
464                     "If you want that, supply a checkpoint file to mdrun\n\n");
465         }
466
467         if (EI_SD(ir->eI) || ir->eI == eiBD)
468         {
469             fprintf(stderr, "\nChanging ld-seed from %"GMX_PRId64 " ", ir->ld_seed);
470             ir->ld_seed = (gmx_int64_t)gmx_rng_make_seed();
471             fprintf(stderr, "to %"GMX_PRId64 "\n\n", ir->ld_seed);
472         }
473
474         frame_fn = ftp2fn(efTRN, NFILE, fnm);
475
476         if (fn2ftp(frame_fn) == efCPT)
477         {
478             int sim_part;
479
480             fprintf(stderr,
481                     "\nREADING STATE FROM CHECKPOINT %s...\n\n",
482                     frame_fn);
483
484             read_checkpoint_state(frame_fn, &sim_part,
485                                   &run_step, &run_t, &state);
486         }
487         else
488         {
489             fprintf(stderr,
490                     "\nREADING COORDS, VELS AND BOX FROM TRAJECTORY %s...\n\n",
491                     frame_fn);
492
493             fp = open_trn(frame_fn, "r");
494             if (bScanEner)
495             {
496                 fp_ener = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
497                 do_enxnms(fp_ener, &nre, &enm);
498                 snew(fr_ener, 1);
499                 fr_ener->t = -1e-12;
500             }
501
502             /* Now scan until the last set of x and v (step == 0)
503              * or the ones at step step.
504              */
505             bFrame = TRUE;
506             frame  = 0;
507             while (bFrame)
508             {
509                 bFrame = fread_trnheader(fp, &head, &bOK);
510                 if (bOK && frame == 0)
511                 {
512                     if (mtop.natoms != head.natoms)
513                     {
514                         gmx_fatal(FARGS, "Number of atoms in Topology (%d) "
515                                   "is not the same as in Trajectory (%d)\n",
516                                   mtop.natoms, head.natoms);
517                     }
518                     snew(newx, head.natoms);
519                     snew(newv, head.natoms);
520                 }
521                 bFrame = bFrame && bOK;
522                 if (bFrame)
523                 {
524                     bOK = fread_htrn(fp, &head, newbox, newx, newv, NULL);
525                 }
526                 bFrame = bFrame && bOK;
527                 bUse   = FALSE;
528                 if (bFrame &&
529                     (head.x_size) && (head.v_size || !bVel))
530                 {
531                     bUse = TRUE;
532                     if (bScanEner)
533                     {
534                         /* Read until the energy time is >= the trajectory time */
535                         while (fr_ener->t < head.t && do_enx(fp_ener, fr_ener))
536                         {
537                             ;
538                         }
539                         bUse = (fr_ener->t == head.t);
540                     }
541                     if (bUse)
542                     {
543                         tmpx                  = newx;
544                         newx                  = state.x;
545                         state.x               = tmpx;
546                         tmpv                  = newv;
547                         newv                  = state.v;
548                         state.v               = tmpv;
549                         run_t                 = head.t;
550                         run_step              = head.step;
551                         state.fep_state       = head.fep_state;
552                         state.lambda[efptFEP] = head.lambda;
553                         copy_mat(newbox, state.box);
554                     }
555                 }
556                 if (bFrame || !bOK)
557                 {
558                     sprintf(buf, "\r%s %s frame %s%s: step %s%s time %s",
559                             "%s", "%s", "%6", GMX_PRId64, "%6", GMX_PRId64, " %8.3f");
560                     fprintf(stderr, buf,
561                             bUse ? "Read   " : "Skipped", ftp2ext(fn2ftp(frame_fn)),
562                             frame, head.step, head.t);
563                     frame++;
564                     if (bTime && (head.t >= start_t))
565                     {
566                         bFrame = FALSE;
567                     }
568                 }
569             }
570             if (bScanEner)
571             {
572                 close_enx(fp_ener);
573                 free_enxframe(fr_ener);
574                 free_enxnms(nre, enm);
575             }
576             close_trn(fp);
577             fprintf(stderr, "\n");
578
579             if (!bOK)
580             {
581                 fprintf(stderr, "%s frame %s (step %s, time %g) is incomplete\n",
582                         ftp2ext(fn2ftp(frame_fn)), gmx_step_str(frame-1, buf2),
583                         gmx_step_str(head.step, buf), head.t);
584             }
585             fprintf(stderr, "\nUsing frame of step %s time %g\n",
586                     gmx_step_str(run_step, buf), run_t);
587
588             if (bNeedEner)
589             {
590                 if (bReadEner)
591                 {
592                     get_enx_state(ftp2fn(efEDR, NFILE, fnm), run_t, &mtop.groups, ir, &state);
593                 }
594                 else
595                 {
596                     fprintf(stderr, "\nWARNING: The simulation uses %s temperature and/or %s pressure coupling,\n"
597                             "         the continuation will only be exact when an energy file is supplied\n\n",
598                             ETCOUPLTYPE(etcNOSEHOOVER),
599                             EPCOUPLTYPE(epcPARRINELLORAHMAN));
600                 }
601             }
602             if (bFepState)
603             {
604                 ir->fepvals->init_fep_state = init_fep_state;
605             }
606         }
607     }
608
609     if (bNsteps)
610     {
611         fprintf(stderr, "Setting nsteps to %s\n", gmx_step_str(nsteps_req, buf));
612         ir->nsteps = nsteps_req;
613     }
614     else
615     {
616         /* Determine total number of steps remaining */
617         if (bExtend)
618         {
619             ir->nsteps = ir->nsteps - (run_step - ir->init_step) + (gmx_int64_t)(extend_t/ir->delta_t + 0.5);
620             printf("Extending remaining runtime of by %g ps (now %s steps)\n",
621                    extend_t, gmx_step_str(ir->nsteps, buf));
622         }
623         else if (bUntil)
624         {
625             printf("nsteps = %s, run_step = %s, current_t = %g, until = %g\n",
626                    gmx_step_str(ir->nsteps, buf),
627                    gmx_step_str(run_step, buf2),
628                    run_t, until_t);
629             ir->nsteps = (gmx_int64_t)((until_t - run_t)/ir->delta_t + 0.5);
630             printf("Extending remaining runtime until %g ps (now %s steps)\n",
631                    until_t, gmx_step_str(ir->nsteps, buf));
632         }
633         else
634         {
635             ir->nsteps -= run_step - ir->init_step;
636             /* Print message */
637             printf("%s steps (%g ps) remaining from first run.\n",
638                    gmx_step_str(ir->nsteps, buf), ir->nsteps*ir->delta_t);
639         }
640     }
641
642     if (bNsteps || bZeroQ || (ir->nsteps > 0))
643     {
644         ir->init_step = run_step;
645
646         if (ftp2bSet(efNDX, NFILE, fnm) ||
647             !(bNsteps || bExtend || bUntil || bTraj))
648         {
649             atoms = gmx_mtop_global_atoms(&mtop);
650             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1,
651                       &gnx, &index, &grpname);
652             if (!bZeroQ)
653             {
654                 bSel = (gnx != state.natoms);
655                 for (i = 0; ((i < gnx) && (!bSel)); i++)
656                 {
657                     bSel = (i != index[i]);
658                 }
659             }
660             else
661             {
662                 bSel = FALSE;
663             }
664             if (bSel)
665             {
666                 fprintf(stderr, "Will write subset %s of original tpx containing %d "
667                         "atoms\n", grpname, gnx);
668                 reduce_topology_x(gnx, index, &mtop, state.x, state.v);
669                 state.natoms = gnx;
670             }
671             else if (bZeroQ)
672             {
673                 zeroq(index, &mtop);
674                 fprintf(stderr, "Zero-ing charges for group %s\n", grpname);
675             }
676             else
677             {
678                 fprintf(stderr, "Will write full tpx file (no selection)\n");
679             }
680         }
681
682         state_t = ir->init_t + ir->init_step*ir->delta_t;
683         sprintf(buf,   "Writing statusfile with starting step %s%s and length %s%s steps...\n", "%10", GMX_PRId64, "%10", GMX_PRId64);
684         fprintf(stderr, buf, ir->init_step, ir->nsteps);
685         fprintf(stderr, "                                 time %10.3f and length %10.3f ps\n",
686                 state_t, ir->nsteps*ir->delta_t);
687         write_tpx_state(opt2fn("-o", NFILE, fnm), ir, &state, &mtop);
688     }
689     else
690     {
691         printf("You've simulated long enough. Not writing tpr file\n");
692     }
693
694     return 0;
695 }