Move smalloc.h to utility/
[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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include "index.h"
43 #include "gmx_fatal.h"
44 #include "sysstuff.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "macros.h"
47 #include "names.h"
48 #include "typedefs.h"
49 #include "gromacs/gmxpreprocess/readir.h"
50 #include "gromacs/commandline/pargs.h"
51 #include "vec.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"
59
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))
61
62 static gmx_bool *bKeepIt(int gnx, int natoms, atom_id index[])
63 {
64     gmx_bool *b;
65     int       i;
66
67     snew(b, natoms);
68     for (i = 0; (i < gnx); i++)
69     {
70         RANGECHK(index[i], natoms);
71         b[index[i]] = TRUE;
72     }
73
74     return b;
75 }
76
77 static atom_id *invind(int gnx, int natoms, atom_id index[])
78 {
79     atom_id *inv;
80     int      i;
81
82     snew(inv, natoms);
83     for (i = 0; (i < gnx); i++)
84     {
85         RANGECHK(index[i], natoms);
86         inv[index[i]] = i;
87     }
88
89     return inv;
90 }
91
92 static void reduce_block(gmx_bool bKeep[], t_block *block,
93                          const char *name)
94 {
95     atom_id *index;
96     int      i, j, newi, newj;
97
98     snew(index, block->nr);
99
100     newi = newj = 0;
101     for (i = 0; (i < block->nr); i++)
102     {
103         for (j = block->index[i]; (j < block->index[i+1]); j++)
104         {
105             if (bKeep[j])
106             {
107                 newj++;
108             }
109         }
110         if (newj > index[newi])
111         {
112             newi++;
113             index[newi] = newj;
114         }
115     }
116
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;
120     block->nr    = newi;
121 }
122
123 static void reduce_blocka(atom_id invindex[], gmx_bool bKeep[], t_blocka *block,
124                           const char *name)
125 {
126     atom_id *index, *a;
127     int      i, j, k, newi, newj;
128
129     snew(index, block->nr);
130     snew(a, block->nra);
131
132     newi = newj = 0;
133     for (i = 0; (i < block->nr); i++)
134     {
135         for (j = block->index[i]; (j < block->index[i+1]); j++)
136         {
137             k = block->a[j];
138             if (bKeep[k])
139             {
140                 a[newj] = invindex[k];
141                 newj++;
142             }
143         }
144         if (newj > index[newi])
145         {
146             newi++;
147             index[newi] = newj;
148         }
149     }
150
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;
154     block->a     = a;
155     block->nr    = newi;
156     block->nra   = newj;
157 }
158
159 static void reduce_rvec(int gnx, atom_id index[], rvec vv[])
160 {
161     rvec *ptr;
162     int   i;
163
164     snew(ptr, gnx);
165     for (i = 0; (i < gnx); i++)
166     {
167         copy_rvec(vv[index[i]], ptr[i]);
168     }
169     for (i = 0; (i < gnx); i++)
170     {
171         copy_rvec(ptr[i], vv[i]);
172     }
173     sfree(ptr);
174 }
175
176 static void reduce_atom(int gnx, atom_id index[], t_atom atom[], char ***atomname,
177                         int *nres, t_resinfo *resinfo)
178 {
179     t_atom    *ptr;
180     char    ***aname;
181     t_resinfo *rinfo;
182     int        i, nr;
183
184     snew(ptr, gnx);
185     snew(aname, gnx);
186     snew(rinfo, atom[index[gnx-1]].resind+1);
187     for (i = 0; (i < gnx); i++)
188     {
189         ptr[i]   = atom[index[i]];
190         aname[i] = atomname[index[i]];
191     }
192     nr = -1;
193     for (i = 0; (i < gnx); i++)
194     {
195         atom[i]     = ptr[i];
196         atomname[i] = aname[i];
197         if ((i == 0) || (atom[i].resind != atom[i-1].resind))
198         {
199             nr++;
200             rinfo[nr] = resinfo[atom[i].resind];
201         }
202         atom[i].resind = nr;
203     }
204     nr++;
205     for (i = 0; (i < nr); i++)
206     {
207         resinfo[i] = rinfo[i];
208     }
209     *nres = nr;
210
211     sfree(aname);
212     sfree(ptr);
213     sfree(rinfo);
214 }
215
216 static void reduce_ilist(atom_id invindex[], gmx_bool bKeep[],
217                          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",
244                 name, il->nr/(nratoms+1),
245                 newnr/(nratoms+1));
246
247         il->nr = newnr;
248         for (i = 0; (i < newnr); i++)
249         {
250             il->iatoms[i] = ia[i];
251         }
252
253         sfree(ia);
254     }
255 }
256
257 static void reduce_topology_x(int gnx, atom_id index[],
258                               gmx_mtop_t *mtop, rvec x[], rvec v[])
259 {
260     t_topology   top;
261     gmx_bool    *bKeep;
262     atom_id     *invindex;
263     int          i;
264
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);
268
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);
276
277     for (i = 0; (i < F_NRE); i++)
278     {
279         reduce_ilist(invindex, bKeep, &(top.idef.il[i]),
280                      interaction_function[i].nratoms,
281                      interaction_function[i].name);
282     }
283
284     top.atoms.nr = gnx;
285
286     mtop->nmoltype = 1;
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++)
291     {
292         mtop->moltype[0].ilist[i] = top.idef.il[i];
293     }
294     mtop->moltype[0].atoms = top.atoms;
295     mtop->moltype[0].cgs   = top.cgs;
296     mtop->moltype[0].excls = top.excls;
297
298     mtop->nmolblock = 1;
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;
305
306     mtop->natoms                 = top.atoms.nr;
307 }
308
309 static void zeroq(atom_id index[], gmx_mtop_t *mtop)
310 {
311     int mt, i;
312
313     for (mt = 0; mt < mtop->nmoltype; mt++)
314     {
315         for (i = 0; (i < mtop->moltype[mt].atoms.nr); i++)
316         {
317             mtop->moltype[mt].atoms.atom[index[i]].q  = 0;
318             mtop->moltype[mt].atoms.atom[index[i]].qB = 0;
319         }
320     }
321 }
322
323 int gmx_convert_tpr(int argc, char *argv[])
324 {
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",
334         "checkpoint files.",
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",
343         "this to work.",
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."
348     };
349
350     const char       *top_fn, *frame_fn;
351     t_fileio         *fp;
352     ener_file_t       fp_ener = NULL;
353     t_trnheader       head;
354     int               i;
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;
359     gmx_mtop_t        mtop;
360     t_atoms           atoms;
361     t_inputrec       *ir, *irnew = NULL;
362     t_gromppopts     *gopts;
363     t_state           state;
364     rvec             *newx = NULL, *newv = NULL, *tmpx, *tmpv;
365     matrix            newbox;
366     int               gnx;
367     char             *grpname;
368     atom_id          *index = NULL;
369     int               nre;
370     gmx_enxnm_t      *enm     = NULL;
371     t_enxframe       *fr_ener = NULL;
372     char              buf[200], buf2[200];
373     output_env_t      oenv;
374     t_filenm          fnm[] = {
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 }
380     };
381 #define NFILE asize(fnm)
382
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" },
405     };
406     int             nerror = 0;
407
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))
411     {
412         return 0;
413     }
414
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);
423
424     top_fn = ftp2fn(efTPX, NFILE, fnm);
425     fprintf(stderr, "Reading toplogy and stuff from %s\n", top_fn);
426
427     snew(ir, 1);
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;
431
432     if (!EI_STATE_VELOCITY(ir->eI))
433     {
434         bVel = FALSE;
435     }
436
437     if (bTraj)
438     {
439         fprintf(stderr, "\n"
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"
445                 "\n");
446
447         if (ir->bContinuation != bContinuation)
448         {
449             fprintf(stderr, "Modifying ir->bContinuation to %s\n",
450                     bool_names[bContinuation]);
451         }
452         ir->bContinuation = bContinuation;
453
454
455         bNeedEner = (ir->epc == epcPARRINELLORAHMAN || ir->etc == etcNOSEHOOVER);
456         bReadEner = (bNeedEner && ftp2bSet(efEDR, NFILE, fnm));
457         bScanEner = (bReadEner && !bTime);
458
459         if (ir->epc != epcNO || EI_SD(ir->eI) || ir->eI == eiBD)
460         {
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");
464         }
465
466         if (EI_SD(ir->eI) || ir->eI == eiBD)
467         {
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);
471         }
472
473         frame_fn = ftp2fn(efTRN, NFILE, fnm);
474
475         if (fn2ftp(frame_fn) == efCPT)
476         {
477             int sim_part;
478
479             fprintf(stderr,
480                     "\nREADING STATE FROM CHECKPOINT %s...\n\n",
481                     frame_fn);
482
483             read_checkpoint_state(frame_fn, &sim_part,
484                                   &run_step, &run_t, &state);
485         }
486         else
487         {
488             fprintf(stderr,
489                     "\nREADING COORDS, VELS AND BOX FROM TRAJECTORY %s...\n\n",
490                     frame_fn);
491
492             fp = open_trn(frame_fn, "r");
493             if (bScanEner)
494             {
495                 fp_ener = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
496                 do_enxnms(fp_ener, &nre, &enm);
497                 snew(fr_ener, 1);
498                 fr_ener->t = -1e-12;
499             }
500
501             /* Now scan until the last set of x and v (step == 0)
502              * or the ones at step step.
503              */
504             bFrame = TRUE;
505             frame  = 0;
506             while (bFrame)
507             {
508                 bFrame = fread_trnheader(fp, &head, &bOK);
509                 if (bOK && frame == 0)
510                 {
511                     if (mtop.natoms != head.natoms)
512                     {
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);
516                     }
517                     snew(newx, head.natoms);
518                     snew(newv, head.natoms);
519                 }
520                 bFrame = bFrame && bOK;
521                 if (bFrame)
522                 {
523                     bOK = fread_htrn(fp, &head, newbox, newx, newv, NULL);
524                 }
525                 bFrame = bFrame && bOK;
526                 bUse   = FALSE;
527                 if (bFrame &&
528                     (head.x_size) && (head.v_size || !bVel))
529                 {
530                     bUse = TRUE;
531                     if (bScanEner)
532                     {
533                         /* Read until the energy time is >= the trajectory time */
534                         while (fr_ener->t < head.t && do_enx(fp_ener, fr_ener))
535                         {
536                             ;
537                         }
538                         bUse = (fr_ener->t == head.t);
539                     }
540                     if (bUse)
541                     {
542                         tmpx                  = newx;
543                         newx                  = state.x;
544                         state.x               = tmpx;
545                         tmpv                  = newv;
546                         newv                  = state.v;
547                         state.v               = tmpv;
548                         run_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);
553                     }
554                 }
555                 if (bFrame || !bOK)
556                 {
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");
559                     fprintf(stderr, buf,
560                             bUse ? "Read   " : "Skipped", ftp2ext(fn2ftp(frame_fn)),
561                             frame, head.step, head.t);
562                     frame++;
563                     if (bTime && (head.t >= start_t))
564                     {
565                         bFrame = FALSE;
566                     }
567                 }
568             }
569             if (bScanEner)
570             {
571                 close_enx(fp_ener);
572                 free_enxframe(fr_ener);
573                 free_enxnms(nre, enm);
574             }
575             close_trn(fp);
576             fprintf(stderr, "\n");
577
578             if (!bOK)
579             {
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);
583             }
584             fprintf(stderr, "\nUsing frame of step %s time %g\n",
585                     gmx_step_str(run_step, buf), run_t);
586
587             if (bNeedEner)
588             {
589                 if (bReadEner)
590                 {
591                     get_enx_state(ftp2fn(efEDR, NFILE, fnm), run_t, &mtop.groups, ir, &state);
592                 }
593                 else
594                 {
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));
599                 }
600             }
601             if (bFepState)
602             {
603                 ir->fepvals->init_fep_state = init_fep_state;
604             }
605         }
606     }
607
608     if (bNsteps)
609     {
610         fprintf(stderr, "Setting nsteps to %s\n", gmx_step_str(nsteps_req, buf));
611         ir->nsteps = nsteps_req;
612     }
613     else
614     {
615         /* Determine total number of steps remaining */
616         if (bExtend)
617         {
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));
621         }
622         else if (bUntil)
623         {
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),
627                    run_t, until_t);
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));
631         }
632         else
633         {
634             ir->nsteps -= run_step - ir->init_step;
635             /* Print message */
636             printf("%s steps (%g ps) remaining from first run.\n",
637                    gmx_step_str(ir->nsteps, buf), ir->nsteps*ir->delta_t);
638         }
639     }
640
641     if (bNsteps || bZeroQ || (ir->nsteps > 0))
642     {
643         ir->init_step = run_step;
644
645         if (ftp2bSet(efNDX, NFILE, fnm) ||
646             !(bNsteps || bExtend || bUntil || bTraj))
647         {
648             atoms = gmx_mtop_global_atoms(&mtop);
649             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1,
650                       &gnx, &index, &grpname);
651             if (!bZeroQ)
652             {
653                 bSel = (gnx != state.natoms);
654                 for (i = 0; ((i < gnx) && (!bSel)); i++)
655                 {
656                     bSel = (i != index[i]);
657                 }
658             }
659             else
660             {
661                 bSel = FALSE;
662             }
663             if (bSel)
664             {
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);
668                 state.natoms = gnx;
669             }
670             else if (bZeroQ)
671             {
672                 zeroq(index, &mtop);
673                 fprintf(stderr, "Zero-ing charges for group %s\n", grpname);
674             }
675             else
676             {
677                 fprintf(stderr, "Will write full tpx file (no selection)\n");
678             }
679         }
680
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);
687     }
688     else
689     {
690         printf("You've simulated long enough. Not writing tpr file\n");
691     }
692
693     return 0;
694 }