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, 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.
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/fileio/tpxio.h"
50 #include "gromacs/fileio/trxio.h"
51 #include "gromacs/fileio/trnio.h"
52 #include "gromacs/fileio/tngio_for_tools.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/fileio/pdbio.h"
55 #include "gromacs/fileio/confio.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/fileio/xtcio.h"
63 #include "gromacs/commandline/pargs.h"
64 #include "gromacs/fileio/xvgr.h"
65 #include "gromacs/math/do_fit.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/pbcutil/rmpbc.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/smalloc.h"
77 euSel, euRect, euTric, euCompact, euNR
81 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
82 rvec x[], atom_id index[], matrix box)
84 int m, i, j, j0, j1, jj, ai, aj;
87 rvec dx, xtest, box_center;
88 int nmol, imol_center;
90 gmx_bool *bMol, *bTmp;
91 rvec *m_com, *m_shift;
99 calc_box_center(ecenter, box, box_center);
101 /* Initiate the pbc structure */
102 memset(&pbc, 0, sizeof(pbc));
103 set_pbc(&pbc, ePBC, box);
105 /* Convert atom index to molecular */
107 molind = top->mols.index;
113 snew(bTmp, top->atoms.nr);
115 for (i = 0; (i < nrefat); i++)
117 /* Mark all molecules in the index */
120 /* Binary search assuming the molecules are sorted */
125 if (ai < molind[j0+1])
129 else if (ai >= molind[j1])
136 if (ai < molind[jj+1])
148 /* Double check whether all atoms in all molecules that are marked are part
149 * of the cluster. Simultaneously compute the center of geometry.
151 min_dist2 = 10*sqr(trace(box));
154 for (i = 0; i < nmol; i++)
156 for (j = molind[i]; j < molind[i+1]; j++)
158 if (bMol[i] && !bTmp[j])
160 gmx_fatal(FARGS, "Molecule %d marked for clustering but not atom %d in it - check your index!", i+1, j+1);
162 else if (!bMol[i] && bTmp[j])
164 gmx_fatal(FARGS, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j+1, i+1);
168 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
171 pbc_dx(&pbc, x[j], x[j-1], dx);
172 rvec_add(x[j-1], dx, x[j]);
174 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
175 rvec_inc(m_com[i], x[j]);
180 /* Normalize center of geometry */
181 fac = 1.0/(molind[i+1]-molind[i]);
182 for (m = 0; (m < DIM); m++)
186 /* Determine which molecule is closest to the center of the box */
187 pbc_dx(&pbc, box_center, m_com[i], dx);
188 tmp_r2 = iprod(dx, dx);
190 if (tmp_r2 < min_dist2)
195 cluster[ncluster++] = i;
202 fprintf(stderr, "No molecules selected in the cluster\n");
205 else if (imol_center == -1)
207 fprintf(stderr, "No central molecules could be found\n");
212 added[nadded++] = imol_center;
213 bMol[imol_center] = FALSE;
215 while (nadded < ncluster)
217 /* Find min distance between cluster molecules and those remaining to be added */
218 min_dist2 = 10*sqr(trace(box));
221 /* Loop over added mols */
222 for (i = 0; i < nadded; i++)
225 /* Loop over all mols */
226 for (j = 0; j < ncluster; j++)
229 /* check those remaining to be added */
232 pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
233 tmp_r2 = iprod(dx, dx);
234 if (tmp_r2 < min_dist2)
244 /* Add the best molecule */
245 added[nadded++] = jmin;
247 /* Calculate the shift from the ai molecule */
248 pbc_dx(&pbc, m_com[jmin], m_com[imin], dx);
249 rvec_add(m_com[imin], dx, xtest);
250 rvec_sub(xtest, m_com[jmin], m_shift[jmin]);
251 rvec_inc(m_com[jmin], m_shift[jmin]);
253 for (j = molind[jmin]; j < molind[jmin+1]; j++)
255 rvec_inc(x[j], m_shift[jmin]);
257 fprintf(stdout, "\rClustering iteration %d of %d...", nadded, ncluster);
267 fprintf(stdout, "\n");
270 static void put_molecule_com_in_box(int unitcell_enum, int ecenter,
272 int natoms, t_atom atom[],
273 int ePBC, matrix box, rvec x[])
277 rvec com, new_com, shift, dx, box_center;
282 calc_box_center(ecenter, box, box_center);
283 set_pbc(&pbc, ePBC, box);
286 gmx_fatal(FARGS, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
288 for (i = 0; (i < mols->nr); i++)
293 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
296 for (d = 0; d < DIM; d++)
302 /* calculate final COM */
303 svmul(1.0/mtot, com, com);
305 /* check if COM is outside box */
306 copy_rvec(com, new_com);
307 switch (unitcell_enum)
310 put_atoms_in_box(ePBC, box, 1, &new_com);
313 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
316 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
319 rvec_sub(new_com, com, shift);
320 if (norm2(shift) > 0)
324 fprintf(debug, "\nShifting position of molecule %d "
325 "by %8.3f %8.3f %8.3f\n", i+1,
326 shift[XX], shift[YY], shift[ZZ]);
328 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
330 rvec_inc(x[j], shift);
336 static void put_residue_com_in_box(int unitcell_enum, int ecenter,
337 int natoms, t_atom atom[],
338 int ePBC, matrix box, rvec x[])
340 atom_id i, j, res_start, res_end, res_nat;
344 rvec box_center, com, new_com, shift;
346 calc_box_center(ecenter, box, box_center);
352 for (i = 0; i < natoms+1; i++)
354 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
356 /* calculate final COM */
358 res_nat = res_end - res_start;
359 svmul(1.0/mtot, com, com);
361 /* check if COM is outside box */
362 copy_rvec(com, new_com);
363 switch (unitcell_enum)
366 put_atoms_in_box(ePBC, box, 1, &new_com);
369 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
372 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
375 rvec_sub(new_com, com, shift);
380 fprintf(debug, "\nShifting position of residue %d (atoms %u-%u) "
381 "by %g,%g,%g\n", atom[res_start].resind+1,
382 res_start+1, res_end+1, shift[XX], shift[YY], shift[ZZ]);
384 for (j = res_start; j < res_end; j++)
386 rvec_inc(x[j], shift);
392 /* remember start of new residue */
399 for (d = 0; d < DIM; d++)
405 presnr = atom[i].resind;
410 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, atom_id ci[])
413 rvec cmin, cmax, box_center, dx;
417 copy_rvec(x[ci[0]], cmin);
418 copy_rvec(x[ci[0]], cmax);
419 for (i = 0; i < nc; i++)
422 for (m = 0; m < DIM; m++)
424 if (x[ai][m] < cmin[m])
428 else if (x[ai][m] > cmax[m])
434 calc_box_center(ecenter, box, box_center);
435 for (m = 0; m < DIM; m++)
437 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
440 for (i = 0; i < n; i++)
447 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
453 strcpy(out_file, base);
464 strncat(out_file, "00000000000", ndigit-nd);
466 sprintf(nbuf, "%d.", file_nr);
467 strcat(out_file, nbuf);
468 strcat(out_file, ext);
471 void check_trn(const char *fn)
473 if ((fn2ftp(fn) != efTRJ) && (fn2ftp(fn) != efTRR))
475 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
479 #ifndef GMX_NATIVE_WINDOWS
480 void do_trunc(const char *fn, real t0)
493 gmx_fatal(FARGS, "You forgot to set the truncation time");
496 /* Check whether this is a .trj file */
499 in = open_trn(fn, "r");
500 fp = gmx_fio_getfp(in);
503 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
509 fpos = gmx_fio_ftell(in);
511 while (!bStop && fread_trnheader(in, &sh, &bOK))
513 fread_htrn(in, &sh, NULL, NULL, NULL, NULL);
514 fpos = gmx_ftell(fp);
518 gmx_fseek(fp, fpos, SEEK_SET);
524 fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
525 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
526 fn, j, t, (long int)fpos);
527 if (1 != scanf("%s", yesno))
529 gmx_fatal(FARGS, "Error reading user input");
531 if (strcmp(yesno, "YES") == 0)
533 fprintf(stderr, "Once again, I'm gonna DO this...\n");
535 if (0 != truncate(fn, fpos))
537 gmx_fatal(FARGS, "Error truncating file %s", fn);
542 fprintf(stderr, "Ok, I'll forget about it\n");
547 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
554 /*! \brief Read a full molecular topology if useful and available.
556 * If the input trajectory file is not in TNG format, and the output
557 * file is in TNG format, then we want to try to read a full topology
558 * (if available), so that we can write molecule information to the
559 * output file. The full topology provides better molecule information
560 * than is available from the normal t_topology data used by GROMACS
563 * Also, the t_topology is only read under (different) particular
564 * conditions. If both apply, then a .tpr file might be read
565 * twice. Trying to fix this redundancy while trjconv is still an
566 * all-purpose tool does not seem worthwhile.
568 * Because of the way gmx_prepare_tng_writing is implemented, the case
569 * where the input TNG file has no molecule information will never
570 * lead to an output TNG file having molecule information. Since
571 * molecule information will generally be present if the input TNG
572 * file was written by a GROMACS tool, this seems like reasonable
574 static gmx_mtop_t *read_mtop_for_tng(const char *tps_file,
575 const char *input_file,
576 const char *output_file)
578 gmx_mtop_t *mtop = NULL;
580 if (fn2bTPX(tps_file) &&
581 efTNG != fn2ftp(input_file) &&
582 efTNG == fn2ftp(output_file))
584 int temp_natoms = -1;
586 read_tpx(tps_file, NULL, NULL, &temp_natoms,
587 NULL, NULL, NULL, mtop);
593 int gmx_trjconv(int argc, char *argv[])
595 const char *desc[] = {
596 "[THISMODULE] can convert trajectory files in many ways:[BR]",
597 "* from one format to another[BR]",
598 "* select a subset of atoms[BR]",
599 "* change the periodicity representation[BR]",
600 "* keep multimeric molecules together[BR]",
601 "* center atoms in the box[BR]",
602 "* fit atoms to reference structure[BR]",
603 "* reduce the number of frames[BR]",
604 "* change the timestamps of the frames ",
605 "([TT]-t0[tt] and [TT]-timestep[tt])[BR]",
606 "* cut the trajectory in small subtrajectories according",
607 "to information in an index file. This allows subsequent analysis of",
608 "the subtrajectories that could, for example, be the result of a",
609 "cluster analysis. Use option [TT]-sub[tt].",
610 "This assumes that the entries in the index file are frame numbers and",
611 "dumps each group in the index file to a separate trajectory file.[BR]",
612 "* select frames within a certain range of a quantity given",
613 "in an [TT].xvg[tt] file.[PAR]",
615 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
618 "The following formats are supported for input and output:",
619 "[TT].xtc[tt], [TT].trr[tt], [TT].trj[tt], [TT].gro[tt], [TT].g96[tt]",
621 "The file formats are detected from the file extension.",
622 "The precision of [TT].xtc[tt] and [TT].gro[tt] output is taken from the",
623 "input file for [TT].xtc[tt], [TT].gro[tt] and [TT].pdb[tt],",
624 "and from the [TT]-ndec[tt] option for other input formats. The precision",
625 "is always taken from [TT]-ndec[tt], when this option is set.",
626 "All other formats have fixed precision. [TT].trr[tt] and [TT].trj[tt]",
627 "output can be single or double precision, depending on the precision",
628 "of the [THISMODULE] binary.",
629 "Note that velocities are only supported in",
630 "[TT].trr[tt], [TT].trj[tt], [TT].gro[tt] and [TT].g96[tt] files.[PAR]",
632 "Option [TT]-sep[tt] can be used to write every frame to a separate",
633 "[TT].gro, .g96[tt] or [TT].pdb[tt] file. By default, all frames all written to one file.",
634 "[TT].pdb[tt] files with all frames concatenated can be viewed with",
635 "[TT]rasmol -nmrpdb[tt].[PAR]",
637 "It is possible to select part of your trajectory and write it out",
638 "to a new trajectory file in order to save disk space, e.g. for leaving",
639 "out the water from a trajectory of a protein in water.",
640 "[BB]ALWAYS[bb] put the original trajectory on tape!",
641 "We recommend to use the portable [TT].xtc[tt] format for your analysis",
642 "to save disk space and to have portable files.[PAR]",
644 "There are two options for fitting the trajectory to a reference",
645 "either for essential dynamics analysis, etc.",
646 "The first option is just plain fitting to a reference structure",
647 "in the structure file. The second option is a progressive fit",
648 "in which the first timeframe is fitted to the reference structure ",
649 "in the structure file to obtain and each subsequent timeframe is ",
650 "fitted to the previously fitted structure. This way a continuous",
651 "trajectory is generated, which might not be the case when using the",
652 "regular fit method, e.g. when your protein undergoes large",
653 "conformational transitions.[PAR]",
655 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
657 "[TT]* mol[tt] puts the center of mass of molecules in the box,",
658 "and requires a run input file to be supplied with [TT]-s[tt].[BR]",
659 "[TT]* res[tt] puts the center of mass of residues in the box.[BR]",
660 "[TT]* atom[tt] puts all the atoms in the box.[BR]",
661 "[TT]* nojump[tt] checks if atoms jump across the box and then puts",
662 "them back. This has the effect that all molecules",
663 "will remain whole (provided they were whole in the initial",
664 "conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
665 "molecules may diffuse out of the box. The starting configuration",
666 "for this procedure is taken from the structure file, if one is",
667 "supplied, otherwise it is the first frame.[BR]",
668 "[TT]* cluster[tt] clusters all the atoms in the selected index",
669 "such that they are all closest to the center of mass of the cluster,",
670 "which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
671 "results if you in fact have a cluster. Luckily that can be checked",
672 "afterwards using a trajectory viewer. Note also that if your molecules",
673 "are broken this will not work either.[BR]",
674 "The separate option [TT]-clustercenter[tt] can be used to specify an",
675 "approximate center for the cluster. This is useful e.g. if you have",
676 "two big vesicles, and you want to maintain their relative positions.[BR]",
677 "[TT]* whole[tt] only makes broken molecules whole.[PAR]",
679 "Option [TT]-ur[tt] sets the unit cell representation for options",
680 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
681 "All three options give different results for triclinic boxes and",
682 "identical results for rectangular boxes.",
683 "[TT]rect[tt] is the ordinary brick shape.",
684 "[TT]tric[tt] is the triclinic unit cell.",
685 "[TT]compact[tt] puts all atoms at the closest distance from the center",
686 "of the box. This can be useful for visualizing e.g. truncated octahedra",
687 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
688 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
689 "is set differently.[PAR]",
691 "Option [TT]-center[tt] centers the system in the box. The user can",
692 "select the group which is used to determine the geometrical center.",
693 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
694 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
695 "[TT]tric[tt]: half of the sum of the box vectors,",
696 "[TT]rect[tt]: half of the box diagonal,",
697 "[TT]zero[tt]: zero.",
698 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
699 "want all molecules in the box after the centering.[PAR]",
701 "It is not always possible to use combinations of [TT]-pbc[tt],",
702 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
703 "you want in one call to [THISMODULE]. Consider using multiple",
704 "calls, and check out the GROMACS website for suggestions.[PAR]",
706 "With [TT]-dt[tt], it is possible to reduce the number of ",
707 "frames in the output. This option relies on the accuracy of the times",
708 "in your input trajectory, so if these are inaccurate use the",
709 "[TT]-timestep[tt] option to modify the time (this can be done",
710 "simultaneously). For making smooth movies, the program [gmx-filter]",
711 "can reduce the number of frames while using low-pass frequency",
712 "filtering, this reduces aliasing of high frequency motions.[PAR]",
714 "Using [TT]-trunc[tt] [THISMODULE] can truncate [TT].trj[tt] in place, i.e.",
715 "without copying the file. This is useful when a run has crashed",
716 "during disk I/O (i.e. full disk), or when two contiguous",
717 "trajectories must be concatenated without having double frames.[PAR]",
719 "Option [TT]-dump[tt] can be used to extract a frame at or near",
720 "one specific time from your trajectory.[PAR]",
722 "Option [TT]-drop[tt] reads an [TT].xvg[tt] file with times and values.",
723 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
724 "frames with a value below and above the value of the respective options",
725 "will not be written."
741 const char *pbc_opt[epNR + 1] =
743 NULL, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
748 const char *unitcell_opt[euNR+1] =
749 { NULL, "rect", "tric", "compact", NULL };
753 ecSel, ecTric, ecRect, ecZero, ecNR
755 const char *center_opt[ecNR+1] =
756 { NULL, "tric", "rect", "zero", NULL };
762 efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
764 const char *fit[efNR + 1] =
766 NULL, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
770 static gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
771 static gmx_bool bCenter = FALSE;
772 static int skip_nr = 1, ndec = 3, nzero = 0;
773 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
774 static rvec newbox = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
775 static char *exec_command = NULL;
776 static real dropunder = 0, dropover = 0;
777 static gmx_bool bRound = FALSE;
782 { "-skip", FALSE, etINT,
783 { &skip_nr }, "Only write every nr-th frame" },
784 { "-dt", FALSE, etTIME,
786 "Only write frame when t MOD dt = first time (%t)" },
787 { "-round", FALSE, etBOOL,
788 { &bRound }, "Round measurements to nearest picosecond"},
789 { "-dump", FALSE, etTIME,
790 { &tdump }, "Dump frame nearest specified time (%t)" },
791 { "-t0", FALSE, etTIME,
793 "Starting time (%t) (default: don't change)" },
794 { "-timestep", FALSE, etTIME,
796 "Change time step between input frames (%t)" },
797 { "-pbc", FALSE, etENUM,
799 "PBC treatment (see help text for full description)" },
800 { "-ur", FALSE, etENUM,
801 { unitcell_opt }, "Unit-cell representation" },
802 { "-center", FALSE, etBOOL,
803 { &bCenter }, "Center atoms in box" },
804 { "-boxcenter", FALSE, etENUM,
805 { center_opt }, "Center for -pbc and -center" },
806 { "-box", FALSE, etRVEC,
808 "Size for new cubic box (default: read from input)" },
809 { "-trans", FALSE, etRVEC,
811 "All coordinates will be translated by trans. This "
812 "can advantageously be combined with -pbc mol -ur "
814 { "-shift", FALSE, etRVEC,
816 "All coordinates will be shifted by framenr*shift" },
817 { "-fit", FALSE, etENUM,
819 "Fit molecule to ref structure in the structure file" },
820 { "-ndec", FALSE, etINT,
822 "Precision for .xtc and .gro writing in number of "
824 { "-vel", FALSE, etBOOL,
825 { &bVels }, "Read and write velocities if possible" },
826 { "-force", FALSE, etBOOL,
827 { &bForce }, "Read and write forces if possible" },
828 #ifndef GMX_NATIVE_WINDOWS
829 { "-trunc", FALSE, etTIME,
831 "Truncate input trajectory file after this time (%t)" },
833 { "-exec", FALSE, etSTR,
835 "Execute command for every output frame with the "
836 "frame number as argument" },
837 { "-split", FALSE, etTIME,
839 "Start writing new file when t MOD split = first "
841 { "-sep", FALSE, etBOOL,
843 "Write each frame to a separate .gro, .g96 or .pdb "
845 { "-nzero", FALSE, etINT,
847 "If the -sep flag is set, use these many digits "
848 "for the file numbers and prepend zeros as needed" },
849 { "-dropunder", FALSE, etREAL,
850 { &dropunder }, "Drop all frames below this value" },
851 { "-dropover", FALSE, etREAL,
852 { &dropover }, "Drop all frames above this value" },
853 { "-conect", FALSE, etBOOL,
855 "Add conect records when writing [TT].pdb[tt] files. Useful "
856 "for visualization of non-standard molecules, e.g. "
857 "coarse grained ones" }
859 #define NPA asize(pa)
862 t_trxstatus *trxout = NULL;
864 int ftp, ftpin = 0, file_nr;
865 t_trxframe fr, frout;
867 rvec *xmem = NULL, *vmem = NULL, *fmem = NULL;
868 rvec *xp = NULL, x_shift, hbox, box_center, dx;
869 real xtcpr, lambda, *w_rls = NULL;
870 int m, i, d, frame, outframe, natoms, nout, ncent, nre, newstep = 0, model_nr;
873 gmx_mtop_t *mtop = NULL;
874 gmx_conect gc = NULL;
876 t_atoms *atoms = NULL, useatoms;
878 atom_id *index, *cindex;
882 int ifit, irms, my_clust = -1;
883 atom_id *ind_fit, *ind_rms;
884 char *gn_fit, *gn_rms;
885 t_cluster_ndx *clust = NULL;
886 t_trxstatus **clust_status = NULL;
887 int *clust_status_id = NULL;
889 int *nfwritten = NULL;
890 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
892 real tshift = 0, t0 = -1, dt = 0.001, prec;
893 gmx_bool bFit, bFitXY, bPFit, bReset;
895 gmx_rmpbc_t gpbc = NULL;
896 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
897 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
898 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetPrec, bNeedPrec;
899 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
900 gmx_bool bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
901 gmx_bool bWriteFrame, bSplitHere;
902 const char *top_file, *in_file, *out_file = NULL;
903 char out_file2[256], *charpt;
904 char *outf_base = NULL;
905 const char *outf_ext = NULL;
906 char top_title[256], title[256], command[256], filemode[5];
908 gmx_bool bWarnCompact = FALSE;
913 { efTRX, "-f", NULL, ffREAD },
914 { efTRO, "-o", NULL, ffWRITE },
915 { efTPS, NULL, NULL, ffOPTRD },
916 { efNDX, NULL, NULL, ffOPTRD },
917 { efNDX, "-fr", "frames", ffOPTRD },
918 { efNDX, "-sub", "cluster", ffOPTRD },
919 { efXVG, "-drop", "drop", ffOPTRD }
921 #define NFILE asize(fnm)
923 if (!parse_common_args(&argc, argv,
924 PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
925 PCA_TIME_UNIT | PCA_BE_NICE,
926 NFILE, fnm, NPA, pa, asize(desc), desc,
932 top_file = ftp2fn(efTPS, NFILE, fnm);
935 /* Check command line */
936 in_file = opt2fn("-f", NFILE, fnm);
940 #ifndef GMX_NATIVE_WINDOWS
941 do_trunc(in_file, ttrunc);
946 /* mark active cmdline options */
947 bSetBox = opt2parg_bSet("-box", NPA, pa);
948 bSetTime = opt2parg_bSet("-t0", NPA, pa);
949 bSetPrec = opt2parg_bSet("-ndec", NPA, pa);
950 bSetUR = opt2parg_bSet("-ur", NPA, pa);
951 bExec = opt2parg_bSet("-exec", NPA, pa);
952 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
953 bTDump = opt2parg_bSet("-dump", NPA, pa);
954 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
955 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
956 bTrans = opt2parg_bSet("-trans", NPA, pa);
957 bSplit = (split_t != 0);
959 /* parse enum options */
960 fit_enum = nenum(fit);
961 bFit = (fit_enum == efFit || fit_enum == efFitXY);
962 bFitXY = fit_enum == efFitXY;
963 bReset = (fit_enum == efReset || fit_enum == efResetXY);
964 bPFit = fit_enum == efPFit;
965 pbc_enum = nenum(pbc_opt);
966 bPBCWhole = pbc_enum == epWhole;
967 bPBCcomRes = pbc_enum == epComRes;
968 bPBCcomMol = pbc_enum == epComMol;
969 bPBCcomAtom = pbc_enum == epComAtom;
970 bNoJump = pbc_enum == epNojump;
971 bCluster = pbc_enum == epCluster;
972 bPBC = pbc_enum != epNone;
973 unitcell_enum = nenum(unitcell_opt);
974 ecenter = nenum(center_opt) - ecTric;
976 /* set and check option dependencies */
979 bFit = TRUE; /* for pfit, fit *must* be set */
983 bReset = TRUE; /* for fit, reset *must* be set */
988 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
990 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
994 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
997 "WARNING: Option for unitcell representation (-ur %s)\n"
998 " only has effect in combination with -pbc %s, %s or %s.\n"
999 " Ingoring unitcell representation.\n\n",
1000 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
1006 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
1007 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1008 "for the rotational fit.\n"
1009 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1013 /* ndec is in nr of decimal places, prec is a multiplication factor: */
1015 for (i = 0; i < ndec; i++)
1020 bIndex = ftp2bSet(efNDX, NFILE, fnm);
1023 /* Determine output type */
1024 out_file = opt2fn("-o", NFILE, fnm);
1025 ftp = fn2ftp(out_file);
1026 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1027 bNeedPrec = (ftp == efXTC || ftp == efGRO);
1030 /* check if velocities are possible in input and output files */
1031 ftpin = fn2ftp(in_file);
1032 bVels = (ftp == efTRR || ftp == efTRJ || ftp == efGRO || ftp == efG96)
1033 && (ftpin == efTRR || ftpin == efTRJ || ftpin == efGRO || ftpin == efG96 ||
1036 if (bSeparate || bSplit)
1038 outf_ext = strrchr(out_file, '.');
1039 if (outf_ext == NULL)
1041 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1043 outf_base = strdup(out_file);
1044 outf_base[outf_ext - out_file] = '\0';
1047 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1050 if ((ftp != efXTC) && (ftp != efTRR))
1052 /* It seems likely that other trajectory file types
1053 * could work here. */
1054 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1057 clust = cluster_index(NULL, opt2fn("-sub", NFILE, fnm));
1059 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1060 * number to check for. In my linux box it is only 16.
1062 if (0 && (clust->clust->nr > FOPEN_MAX-4))
1064 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1065 " trajectories.\ntry splitting the index file in %d parts.\n"
1067 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1069 gmx_warning("The -sub option could require as many open output files as there are\n"
1070 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1071 "try reducing the number of index groups in the file, and perhaps\n"
1072 "using trjconv -sub several times on different chunks of your index file.\n",
1075 snew(clust_status, clust->clust->nr);
1076 snew(clust_status_id, clust->clust->nr);
1077 snew(nfwritten, clust->clust->nr);
1078 for (i = 0; (i < clust->clust->nr); i++)
1080 clust_status[i] = NULL;
1081 clust_status_id[i] = -1;
1083 bSeparate = bSplit = FALSE;
1090 mtop = read_mtop_for_tng(top_file, in_file, out_file);
1092 /* Determine whether to read a topology */
1093 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1094 bRmPBC || bReset || bPBCcomMol || bCluster ||
1095 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1097 /* Determine if when can read index groups */
1098 bIndex = (bIndex || bTPS);
1102 read_tps_conf(top_file, top_title, &top, &ePBC, &xp, NULL, top_box,
1103 bReset || bPBCcomRes);
1106 if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1108 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1111 /* top_title is only used for gro and pdb,
1112 * the header in such a file is top_title t= ...
1113 * to prevent a double t=, remove it from top_title
1115 if ((charpt = strstr(top_title, " t= ")))
1122 gc = gmx_conect_generate(&top);
1126 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1130 /* get frame number index */
1132 if (opt2bSet("-fr", NFILE, fnm))
1134 printf("Select groups of frame number indices:\n");
1135 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (atom_id **)&frindex, &frname);
1138 for (i = 0; i < nrfri; i++)
1140 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1145 /* get index groups etc. */
1148 printf("Select group for %s fit\n",
1149 bFit ? "least squares" : "translational");
1150 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1151 1, &ifit, &ind_fit, &gn_fit);
1157 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1161 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1167 printf("Select group for clustering\n");
1168 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1169 1, &ifit, &ind_fit, &gn_fit);
1176 printf("Select group for centering\n");
1177 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1178 1, &ncent, &cindex, &grpnm);
1180 printf("Select group for output\n");
1181 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1182 1, &nout, &index, &grpnm);
1186 /* no index file, so read natoms from TRX */
1187 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1189 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1194 snew(index, natoms);
1195 for (i = 0; i < natoms; i++)
1209 snew(w_rls, atoms->nr);
1210 for (i = 0; (i < ifit); i++)
1212 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1215 /* Restore reference structure and set to origin,
1216 store original location (to put structure back) */
1219 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1221 copy_rvec(xp[index[0]], x_shift);
1222 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, NULL, xp, w_rls);
1223 rvec_dec(x_shift, xp[index[0]]);
1227 clear_rvec(x_shift);
1230 if (bDropUnder || bDropOver)
1232 /* Read the .xvg file with the drop values */
1233 fprintf(stderr, "\nReading drop file ...");
1234 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1235 fprintf(stderr, " %d time points\n", ndrop);
1236 if (ndrop == 0 || ncol < 2)
1238 gmx_fatal(FARGS, "Found no data points in %s",
1239 opt2fn("-drop", NFILE, fnm));
1245 /* Make atoms struct for output in GRO or PDB files */
1246 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1248 /* get memory for stuff to go in .pdb file */
1249 init_t_atoms(&useatoms, atoms->nr, FALSE);
1250 sfree(useatoms.resinfo);
1251 useatoms.resinfo = atoms->resinfo;
1252 for (i = 0; (i < nout); i++)
1254 useatoms.atomname[i] = atoms->atomname[index[i]];
1255 useatoms.atom[i] = atoms->atom[index[i]];
1256 useatoms.nres = max(useatoms.nres, useatoms.atom[i].resind+1);
1260 /* select what to read */
1261 if (ftp == efTRR || ftp == efTRJ)
1271 flags = flags | TRX_READ_V;
1275 flags = flags | TRX_READ_F;
1278 /* open trx file for reading */
1279 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1282 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1286 if (bSetPrec || !fr.bPrec)
1288 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1292 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1296 if (bHaveFirstFrame)
1298 set_trxframe_ePBC(&fr, ePBC);
1304 tshift = tzero-fr.time;
1314 /* check if index is meaningful */
1315 for (i = 0; i < nout; i++)
1317 if (index[i] >= natoms)
1320 "Index[%d] %d is larger than the number of atoms in the\n"
1321 "trajectory file (%d). There is a mismatch in the contents\n"
1322 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1324 bCopy = bCopy || (i != index[i]);
1328 /* open output for writing */
1329 strcpy(filemode, "w");
1333 trjtools_gmx_prepare_tng_writing(out_file,
1347 if (!bSplit && !bSubTraj)
1349 trxout = open_trx(out_file, filemode);
1355 if (( !bSeparate && !bSplit ) && !bSubTraj)
1357 out = gmx_ffopen(out_file, filemode);
1361 gmx_incons("Illegal output file format");
1377 /* Start the big loop over frames */
1384 /* Main loop over frames */
1395 /*if (frame >= clust->clust->nra)
1396 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1397 if (frame > clust->maxframe)
1403 my_clust = clust->inv_clust[frame];
1405 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1406 (my_clust == NO_ATID))
1414 /* generate new box */
1416 for (m = 0; m < DIM; m++)
1418 fr.box[m][m] = newbox[m];
1424 for (i = 0; i < natoms; i++)
1426 rvec_inc(fr.x[i], trans);
1432 /* determine timestep */
1445 /* This is not very elegant, as one can not dump a frame after
1446 * a timestep with is more than twice as small as the first one. */
1447 bDumpFrame = (fr.time > tdump-0.5*dt) && (fr.time <= tdump+0.5*dt);
1454 /* determine if an atom jumped across the box and reset it if so */
1455 if (bNoJump && (bTPS || frame != 0))
1457 for (d = 0; d < DIM; d++)
1459 hbox[d] = 0.5*fr.box[d][d];
1461 for (i = 0; i < natoms; i++)
1465 rvec_dec(fr.x[i], x_shift);
1467 for (m = DIM-1; m >= 0; m--)
1471 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1473 for (d = 0; d <= m; d++)
1475 fr.x[i][d] += fr.box[m][d];
1478 while (fr.x[i][m]-xp[i][m] > hbox[m])
1480 for (d = 0; d <= m; d++)
1482 fr.x[i][d] -= fr.box[m][d];
1491 calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1496 /* Now modify the coords according to the flags,
1497 for normal fit, this is only done for output frames */
1500 gmx_rmpbc_trxfr(gpbc, &fr);
1503 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1504 do_fit(natoms, w_rls, xp, fr.x);
1507 /* store this set of coordinates for future use */
1508 if (bPFit || bNoJump)
1514 for (i = 0; (i < natoms); i++)
1516 copy_rvec(fr.x[i], xp[i]);
1517 rvec_inc(fr.x[i], x_shift);
1523 /* see if we have a frame from the frame index group */
1524 for (i = 0; i < nrfri && !bDumpFrame; i++)
1526 bDumpFrame = frame == frindex[i];
1529 if (debug && bDumpFrame)
1531 fprintf(debug, "dumping %d\n", frame);
1535 ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1537 if (bWriteFrame && (bDropUnder || bDropOver))
1539 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1544 if (fabs(dropval[0][drop0] - fr.time)
1545 < fabs(dropval[0][drop1] - fr.time))
1553 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1554 (bDropOver && dropval[1][dropuse] > dropover))
1556 bWriteFrame = FALSE;
1566 fr.time = tzero+frame*timestep;
1576 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1577 output_env_conv_time(oenv, fr.time), output_env_get_time_unit(oenv));
1580 /* check for writing at each delta_t */
1581 bDoIt = (delta_t == 0);
1586 bDoIt = bRmod(fr.time, tzero, delta_t);
1590 /* round() is not C89 compatible, so we do this: */
1591 bDoIt = bRmod(floor(fr.time+0.5), floor(tzero+0.5),
1592 floor(delta_t+0.5));
1596 if (bDoIt || bTDump)
1598 /* print sometimes */
1599 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1601 fprintf(stderr, " -> frame %6d time %8.3f \r",
1602 outframe, output_env_conv_time(oenv, fr.time));
1607 /* Now modify the coords according to the flags,
1608 for PFit we did this already! */
1612 gmx_rmpbc_trxfr(gpbc, &fr);
1617 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1620 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1624 for (i = 0; i < natoms; i++)
1626 rvec_inc(fr.x[i], x_shift);
1633 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1639 switch (unitcell_enum)
1642 put_atoms_in_box(ePBC, fr.box, natoms, fr.x);
1645 put_atoms_in_triclinic_unitcell(ecenter, fr.box, natoms, fr.x);
1648 warn = put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1650 if (warn && !bWarnCompact)
1652 fprintf(stderr, "\n%s\n", warn);
1653 bWarnCompact = TRUE;
1660 put_residue_com_in_box(unitcell_enum, ecenter,
1661 natoms, atoms->atom, ePBC, fr.box, fr.x);
1665 put_molecule_com_in_box(unitcell_enum, ecenter,
1667 natoms, atoms->atom, ePBC, fr.box, fr.x);
1669 /* Copy the input trxframe struct to the output trxframe struct */
1671 frout.bV = (frout.bV && bVels);
1672 frout.bF = (frout.bF && bForce);
1673 frout.natoms = nout;
1674 if (bNeedPrec && (bSetPrec || !fr.bPrec))
1690 for (i = 0; i < nout; i++)
1692 copy_rvec(fr.x[index[i]], frout.x[i]);
1695 copy_rvec(fr.v[index[i]], frout.v[i]);
1699 copy_rvec(fr.f[index[i]], frout.f[i]);
1704 if (opt2parg_bSet("-shift", NPA, pa))
1706 for (i = 0; i < nout; i++)
1708 for (d = 0; d < DIM; d++)
1710 frout.x[i][d] += outframe*shift[d];
1717 bSplitHere = bSplit && bRmod(fr.time, tzero, split_t);
1721 /* round() is not C89 compatible, so we do this: */
1722 bSplitHere = bSplit && bRmod(floor(fr.time+0.5),
1724 floor(split_t+0.5));
1726 if (bSeparate || bSplitHere)
1728 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1734 write_tng_frame(trxout, &frout);
1735 // TODO when trjconv behaves better: work how to read and write lambda
1746 trxout = open_trx(out_file2, filemode);
1753 if (clust_status_id[my_clust] == -1)
1755 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1756 clust_status[my_clust] = open_trx(buf, "w");
1757 clust_status_id[my_clust] = 1;
1760 else if (clust_status_id[my_clust] == -2)
1762 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1763 clust->grpname[my_clust], ntrxopen, frame,
1766 write_trxframe(clust_status[my_clust], &frout, gc);
1767 nfwritten[my_clust]++;
1768 if (nfwritten[my_clust] ==
1769 (clust->clust->index[my_clust+1]-
1770 clust->clust->index[my_clust]))
1772 close_trx(clust_status[my_clust]);
1773 clust_status[my_clust] = NULL;
1774 clust_status_id[my_clust] = -2;
1778 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1785 write_trxframe(trxout, &frout, gc);
1791 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1792 top_title, fr.time);
1793 if (bSeparate || bSplitHere)
1795 out = gmx_ffopen(out_file2, "w");
1800 write_hconf_p(out, title, &useatoms, prec2ndec(frout.prec),
1801 frout.x, frout.bV ? frout.v : NULL, frout.box);
1804 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1805 sprintf(title, "%s t= %9.5f", top_title, frout.time);
1806 /* if reading from pdb, we want to keep the original
1807 model numbering else we write the output frame
1808 number plus one, because model 0 is not allowed in pdb */
1809 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1817 write_pdbfile(out, title, &useatoms, frout.x,
1818 frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1821 frout.title = title;
1822 if (bSeparate || bTDump)
1824 frout.bTitle = TRUE;
1827 frout.bAtoms = TRUE;
1829 frout.atoms = &useatoms;
1830 frout.bStep = FALSE;
1831 frout.bTime = FALSE;
1835 frout.bTitle = (outframe == 0);
1836 frout.bAtoms = FALSE;
1840 write_g96_conf(out, &frout, -1, NULL);
1849 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1851 if (bSeparate || bSplitHere)
1856 /* execute command */
1860 sprintf(c, "%s %d", exec_command, file_nr-1);
1861 /*fprintf(stderr,"Executing '%s'\n",c);*/
1862 #ifdef GMX_NO_SYSTEM
1863 printf("Warning-- No calls to system(3) supported on this platform.");
1864 printf("Warning-- Skipping execution of 'system(\"%s\")'.", c);
1868 gmx_fatal(FARGS, "Error executing command: %s", c);
1876 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1878 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1881 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1883 fprintf(stderr, "\nWARNING no output, "
1884 "last frame read at t=%g\n", fr.time);
1886 fprintf(stderr, "\n");
1893 gmx_rmpbc_done(gpbc);
1900 else if (out != NULL)
1906 for (i = 0; (i < clust->clust->nr); i++)
1908 if (clust_status_id[i] >= 0)
1910 close_trx(clust_status[i]);
1918 do_view(oenv, out_file, NULL);