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.
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/fileio/trnio.h"
53 #include "gromacs/fileio/tngio_for_tools.h"
54 #include "gromacs/commandline/pargs.h"
55 #include "gromacs/fileio/futil.h"
56 #include "gromacs/fileio/pdbio.h"
57 #include "gromacs/fileio/confio.h"
61 #include "gromacs/fileio/xtcio.h"
64 #include "gromacs/fileio/g87io.h"
75 euSel, euRect, euTric, euCompact, euNR
79 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
80 rvec x[], atom_id index[], matrix box)
82 int m, i, j, j0, j1, jj, ai, aj;
85 rvec dx, xtest, box_center;
86 int nmol, imol_center;
88 gmx_bool *bMol, *bTmp;
89 rvec *m_com, *m_shift;
97 calc_box_center(ecenter, box, box_center);
99 /* Initiate the pbc structure */
100 memset(&pbc, 0, sizeof(pbc));
101 set_pbc(&pbc, ePBC, box);
103 /* Convert atom index to molecular */
105 molind = top->mols.index;
111 snew(bTmp, top->atoms.nr);
113 for (i = 0; (i < nrefat); i++)
115 /* Mark all molecules in the index */
118 /* Binary search assuming the molecules are sorted */
123 if (ai < molind[j0+1])
127 else if (ai >= molind[j1])
134 if (ai < molind[jj+1])
146 /* Double check whether all atoms in all molecules that are marked are part
147 * of the cluster. Simultaneously compute the center of geometry.
149 min_dist2 = 10*sqr(trace(box));
152 for (i = 0; i < nmol; i++)
154 for (j = molind[i]; j < molind[i+1]; j++)
156 if (bMol[i] && !bTmp[j])
158 gmx_fatal(FARGS, "Molecule %d marked for clustering but not atom %d in it - check your index!", i+1, j+1);
160 else if (!bMol[i] && bTmp[j])
162 gmx_fatal(FARGS, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j+1, i+1);
166 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
169 pbc_dx(&pbc, x[j], x[j-1], dx);
170 rvec_add(x[j-1], dx, x[j]);
172 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
173 rvec_inc(m_com[i], x[j]);
178 /* Normalize center of geometry */
179 fac = 1.0/(molind[i+1]-molind[i]);
180 for (m = 0; (m < DIM); m++)
184 /* Determine which molecule is closest to the center of the box */
185 pbc_dx(&pbc, box_center, m_com[i], dx);
186 tmp_r2 = iprod(dx, dx);
188 if (tmp_r2 < min_dist2)
193 cluster[ncluster++] = i;
200 fprintf(stderr, "No molecules selected in the cluster\n");
203 else if (imol_center == -1)
205 fprintf(stderr, "No central molecules could be found\n");
210 added[nadded++] = imol_center;
211 bMol[imol_center] = FALSE;
213 while (nadded < ncluster)
215 /* Find min distance between cluster molecules and those remaining to be added */
216 min_dist2 = 10*sqr(trace(box));
219 /* Loop over added mols */
220 for (i = 0; i < nadded; i++)
223 /* Loop over all mols */
224 for (j = 0; j < ncluster; j++)
227 /* check those remaining to be added */
230 pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
231 tmp_r2 = iprod(dx, dx);
232 if (tmp_r2 < min_dist2)
242 /* Add the best molecule */
243 added[nadded++] = jmin;
245 /* Calculate the shift from the ai molecule */
246 pbc_dx(&pbc, m_com[jmin], m_com[imin], dx);
247 rvec_add(m_com[imin], dx, xtest);
248 rvec_sub(xtest, m_com[jmin], m_shift[jmin]);
249 rvec_inc(m_com[jmin], m_shift[jmin]);
251 for (j = molind[jmin]; j < molind[jmin+1]; j++)
253 rvec_inc(x[j], m_shift[jmin]);
255 fprintf(stdout, "\rClustering iteration %d of %d...", nadded, ncluster);
265 fprintf(stdout, "\n");
268 static void put_molecule_com_in_box(int unitcell_enum, int ecenter,
270 int natoms, t_atom atom[],
271 int ePBC, matrix box, rvec x[])
275 rvec com, new_com, shift, dx, box_center;
280 calc_box_center(ecenter, box, box_center);
281 set_pbc(&pbc, ePBC, box);
284 gmx_fatal(FARGS, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
286 for (i = 0; (i < mols->nr); i++)
291 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
294 for (d = 0; d < DIM; d++)
300 /* calculate final COM */
301 svmul(1.0/mtot, com, com);
303 /* check if COM is outside box */
304 copy_rvec(com, new_com);
305 switch (unitcell_enum)
308 put_atoms_in_box(ePBC, box, 1, &new_com);
311 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
314 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
317 rvec_sub(new_com, com, shift);
318 if (norm2(shift) > 0)
322 fprintf(debug, "\nShifting position of molecule %d "
323 "by %8.3f %8.3f %8.3f\n", i+1,
324 shift[XX], shift[YY], shift[ZZ]);
326 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
328 rvec_inc(x[j], shift);
334 static void put_residue_com_in_box(int unitcell_enum, int ecenter,
335 int natoms, t_atom atom[],
336 int ePBC, matrix box, rvec x[])
338 atom_id i, j, res_start, res_end, res_nat;
342 rvec box_center, com, new_com, shift;
344 calc_box_center(ecenter, box, box_center);
350 for (i = 0; i < natoms+1; i++)
352 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
354 /* calculate final COM */
356 res_nat = res_end - res_start;
357 svmul(1.0/mtot, com, com);
359 /* check if COM is outside box */
360 copy_rvec(com, new_com);
361 switch (unitcell_enum)
364 put_atoms_in_box(ePBC, box, 1, &new_com);
367 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
370 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
373 rvec_sub(new_com, com, shift);
378 fprintf(debug, "\nShifting position of residue %d (atoms %u-%u) "
379 "by %g,%g,%g\n", atom[res_start].resind+1,
380 res_start+1, res_end+1, shift[XX], shift[YY], shift[ZZ]);
382 for (j = res_start; j < res_end; j++)
384 rvec_inc(x[j], shift);
390 /* remember start of new residue */
397 for (d = 0; d < DIM; d++)
403 presnr = atom[i].resind;
408 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, atom_id ci[])
411 rvec cmin, cmax, box_center, dx;
415 copy_rvec(x[ci[0]], cmin);
416 copy_rvec(x[ci[0]], cmax);
417 for (i = 0; i < nc; i++)
420 for (m = 0; m < DIM; m++)
422 if (x[ai][m] < cmin[m])
426 else if (x[ai][m] > cmax[m])
432 calc_box_center(ecenter, box, box_center);
433 for (m = 0; m < DIM; m++)
435 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
438 for (i = 0; i < n; i++)
445 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
451 strcpy(out_file, base);
462 strncat(out_file, "00000000000", ndigit-nd);
464 sprintf(nbuf, "%d.", file_nr);
465 strcat(out_file, nbuf);
466 strcat(out_file, ext);
469 void check_trn(const char *fn)
471 if ((fn2ftp(fn) != efTRJ) && (fn2ftp(fn) != efTRR))
473 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
477 #ifndef GMX_NATIVE_WINDOWS
478 void do_trunc(const char *fn, real t0)
491 gmx_fatal(FARGS, "You forgot to set the truncation time");
494 /* Check whether this is a .trj file */
497 in = open_trn(fn, "r");
498 fp = gmx_fio_getfp(in);
501 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
507 fpos = gmx_fio_ftell(in);
509 while (!bStop && fread_trnheader(in, &sh, &bOK))
511 fread_htrn(in, &sh, NULL, NULL, NULL, NULL);
512 fpos = gmx_ftell(fp);
516 gmx_fseek(fp, fpos, SEEK_SET);
522 fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
523 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
524 fn, j, t, (long int)fpos);
525 if (1 != scanf("%s", yesno))
527 gmx_fatal(FARGS, "Error reading user input");
529 if (strcmp(yesno, "YES") == 0)
531 fprintf(stderr, "Once again, I'm gonna DO this...\n");
533 if (0 != truncate(fn, fpos))
535 gmx_fatal(FARGS, "Error truncating file %s", fn);
540 fprintf(stderr, "Ok, I'll forget about it\n");
545 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
552 /*! \brief Read a full molecular topology if useful and available.
554 * If the input trajectory file is not in TNG format, and the output
555 * file is in TNG format, then we want to try to read a full topology
556 * (if available), so that we can write molecule information to the
557 * output file. The full topology provides better molecule information
558 * than is available from the normal t_topology data used by GROMACS
561 * Also, the t_topology is only read under (different) particular
562 * conditions. If both apply, then a .tpr file might be read
563 * twice. Trying to fix this redundancy while trjconv is still an
564 * all-purpose tool does not seem worthwhile.
566 * Because of the way gmx_prepare_tng_writing is implemented, the case
567 * where the input TNG file has no molecule information will never
568 * lead to an output TNG file having molecule information. Since
569 * molecule information will generally be present if the input TNG
570 * file was written by a GROMACS tool, this seems like reasonable
572 static gmx_mtop_t *read_mtop_for_tng(const char *tps_file,
573 const char *input_file,
574 const char *output_file)
576 gmx_mtop_t *mtop = NULL;
578 if (fn2bTPX(tps_file) &&
579 efTNG != fn2ftp(input_file) &&
580 efTNG == fn2ftp(output_file))
582 int temp_natoms = -1;
584 read_tpx(tps_file, NULL, NULL, &temp_natoms,
585 NULL, NULL, NULL, mtop);
591 int gmx_trjconv(int argc, char *argv[])
593 const char *desc[] = {
594 "[THISMODULE] can convert trajectory files in many ways:[BR]",
595 "* from one format to another[BR]",
596 "* select a subset of atoms[BR]",
597 "* change the periodicity representation[BR]",
598 "* keep multimeric molecules together[BR]",
599 "* center atoms in the box[BR]",
600 "* fit atoms to reference structure[BR]",
601 "* reduce the number of frames[BR]",
602 "* change the timestamps of the frames ",
603 "([TT]-t0[tt] and [TT]-timestep[tt])[BR]",
604 "* cut the trajectory in small subtrajectories according",
605 "to information in an index file. This allows subsequent analysis of",
606 "the subtrajectories that could, for example, be the result of a",
607 "cluster analysis. Use option [TT]-sub[tt].",
608 "This assumes that the entries in the index file are frame numbers and",
609 "dumps each group in the index file to a separate trajectory file.[BR]",
610 "* select frames within a certain range of a quantity given",
611 "in an [TT].xvg[tt] file.[PAR]",
613 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
616 "Currently seven formats are supported for input and output:",
617 "[TT].xtc[tt], [TT].trr[tt], [TT].trj[tt], [TT].gro[tt], [TT].g96[tt],",
618 "[TT].pdb[tt] and [TT].g87[tt].",
619 "The file formats are detected from the file extension.",
620 "The precision of [TT].xtc[tt] and [TT].gro[tt] output is taken from the",
621 "input file for [TT].xtc[tt], [TT].gro[tt] and [TT].pdb[tt],",
622 "and from the [TT]-ndec[tt] option for other input formats. The precision",
623 "is always taken from [TT]-ndec[tt], when this option is set.",
624 "All other formats have fixed precision. [TT].trr[tt] and [TT].trj[tt]",
625 "output can be single or double precision, depending on the precision",
626 "of the [THISMODULE] binary.",
627 "Note that velocities are only supported in",
628 "[TT].trr[tt], [TT].trj[tt], [TT].gro[tt] and [TT].g96[tt] files.[PAR]",
630 "Option [TT]-sep[tt] can be used to write every frame to a separate",
631 "[TT].gro, .g96[tt] or [TT].pdb[tt] file. By default, all frames all written to one file.",
632 "[TT].pdb[tt] files with all frames concatenated can be viewed with",
633 "[TT]rasmol -nmrpdb[tt].[PAR]",
635 "It is possible to select part of your trajectory and write it out",
636 "to a new trajectory file in order to save disk space, e.g. for leaving",
637 "out the water from a trajectory of a protein in water.",
638 "[BB]ALWAYS[bb] put the original trajectory on tape!",
639 "We recommend to use the portable [TT].xtc[tt] format for your analysis",
640 "to save disk space and to have portable files.[PAR]",
642 "There are two options for fitting the trajectory to a reference",
643 "either for essential dynamics analysis, etc.",
644 "The first option is just plain fitting to a reference structure",
645 "in the structure file. The second option is a progressive fit",
646 "in which the first timeframe is fitted to the reference structure ",
647 "in the structure file to obtain and each subsequent timeframe is ",
648 "fitted to the previously fitted structure. This way a continuous",
649 "trajectory is generated, which might not be the case when using the",
650 "regular fit method, e.g. when your protein undergoes large",
651 "conformational transitions.[PAR]",
653 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
655 "[TT]* mol[tt] puts the center of mass of molecules in the box,",
656 "and requires a run input file to be supplied with [TT]-s[tt].[BR]",
657 "[TT]* res[tt] puts the center of mass of residues in the box.[BR]",
658 "[TT]* atom[tt] puts all the atoms in the box.[BR]",
659 "[TT]* nojump[tt] checks if atoms jump across the box and then puts",
660 "them back. This has the effect that all molecules",
661 "will remain whole (provided they were whole in the initial",
662 "conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
663 "molecules may diffuse out of the box. The starting configuration",
664 "for this procedure is taken from the structure file, if one is",
665 "supplied, otherwise it is the first frame.[BR]",
666 "[TT]* cluster[tt] clusters all the atoms in the selected index",
667 "such that they are all closest to the center of mass of the cluster,",
668 "which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
669 "results if you in fact have a cluster. Luckily that can be checked",
670 "afterwards using a trajectory viewer. Note also that if your molecules",
671 "are broken this will not work either.[BR]",
672 "The separate option [TT]-clustercenter[tt] can be used to specify an",
673 "approximate center for the cluster. This is useful e.g. if you have",
674 "two big vesicles, and you want to maintain their relative positions.[BR]",
675 "[TT]* whole[tt] only makes broken molecules whole.[PAR]",
677 "Option [TT]-ur[tt] sets the unit cell representation for options",
678 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
679 "All three options give different results for triclinic boxes and",
680 "identical results for rectangular boxes.",
681 "[TT]rect[tt] is the ordinary brick shape.",
682 "[TT]tric[tt] is the triclinic unit cell.",
683 "[TT]compact[tt] puts all atoms at the closest distance from the center",
684 "of the box. This can be useful for visualizing e.g. truncated octahedra",
685 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
686 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
687 "is set differently.[PAR]",
689 "Option [TT]-center[tt] centers the system in the box. The user can",
690 "select the group which is used to determine the geometrical center.",
691 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
692 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
693 "[TT]tric[tt]: half of the sum of the box vectors,",
694 "[TT]rect[tt]: half of the box diagonal,",
695 "[TT]zero[tt]: zero.",
696 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
697 "want all molecules in the box after the centering.[PAR]",
699 "It is not always possible to use combinations of [TT]-pbc[tt],",
700 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
701 "you want in one call to [THISMODULE]. Consider using multiple",
702 "calls, and check out the GROMACS website for suggestions.[PAR]",
704 "With [TT]-dt[tt], it is possible to reduce the number of ",
705 "frames in the output. This option relies on the accuracy of the times",
706 "in your input trajectory, so if these are inaccurate use the",
707 "[TT]-timestep[tt] option to modify the time (this can be done",
708 "simultaneously). For making smooth movies, the program [gmx-filter]",
709 "can reduce the number of frames while using low-pass frequency",
710 "filtering, this reduces aliasing of high frequency motions.[PAR]",
712 "Using [TT]-trunc[tt] [THISMODULE] can truncate [TT].trj[tt] in place, i.e.",
713 "without copying the file. This is useful when a run has crashed",
714 "during disk I/O (i.e. full disk), or when two contiguous",
715 "trajectories must be concatenated without having double frames.[PAR]",
717 "Option [TT]-dump[tt] can be used to extract a frame at or near",
718 "one specific time from your trajectory.[PAR]",
720 "Option [TT]-drop[tt] reads an [TT].xvg[tt] file with times and values.",
721 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
722 "frames with a value below and above the value of the respective options",
723 "will not be written."
739 const char *pbc_opt[epNR + 1] =
741 NULL, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
746 const char *unitcell_opt[euNR+1] =
747 { NULL, "rect", "tric", "compact", NULL };
751 ecSel, ecTric, ecRect, ecZero, ecNR
753 const char *center_opt[ecNR+1] =
754 { NULL, "tric", "rect", "zero", NULL };
760 efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
762 const char *fit[efNR + 1] =
764 NULL, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
768 static gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
769 static gmx_bool bCenter = FALSE;
770 static int skip_nr = 1, ndec = 3, nzero = 0;
771 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
772 static rvec newbox = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
773 static char *exec_command = NULL;
774 static real dropunder = 0, dropover = 0;
775 static gmx_bool bRound = FALSE;
780 { "-skip", FALSE, etINT,
781 { &skip_nr }, "Only write every nr-th frame" },
782 { "-dt", FALSE, etTIME,
784 "Only write frame when t MOD dt = first time (%t)" },
785 { "-round", FALSE, etBOOL,
786 { &bRound }, "Round measurements to nearest picosecond"},
787 { "-dump", FALSE, etTIME,
788 { &tdump }, "Dump frame nearest specified time (%t)" },
789 { "-t0", FALSE, etTIME,
791 "Starting time (%t) (default: don't change)" },
792 { "-timestep", FALSE, etTIME,
794 "Change time step between input frames (%t)" },
795 { "-pbc", FALSE, etENUM,
797 "PBC treatment (see help text for full description)" },
798 { "-ur", FALSE, etENUM,
799 { unitcell_opt }, "Unit-cell representation" },
800 { "-center", FALSE, etBOOL,
801 { &bCenter }, "Center atoms in box" },
802 { "-boxcenter", FALSE, etENUM,
803 { center_opt }, "Center for -pbc and -center" },
804 { "-box", FALSE, etRVEC,
806 "Size for new cubic box (default: read from input)" },
807 { "-trans", FALSE, etRVEC,
809 "All coordinates will be translated by trans. This "
810 "can advantageously be combined with -pbc mol -ur "
812 { "-shift", FALSE, etRVEC,
814 "All coordinates will be shifted by framenr*shift" },
815 { "-fit", FALSE, etENUM,
817 "Fit molecule to ref structure in the structure file" },
818 { "-ndec", FALSE, etINT,
820 "Precision for .xtc and .gro writing in number of "
822 { "-vel", FALSE, etBOOL,
823 { &bVels }, "Read and write velocities if possible" },
824 { "-force", FALSE, etBOOL,
825 { &bForce }, "Read and write forces if possible" },
826 #ifndef GMX_NATIVE_WINDOWS
827 { "-trunc", FALSE, etTIME,
829 "Truncate input trajectory file after this time (%t)" },
831 { "-exec", FALSE, etSTR,
833 "Execute command for every output frame with the "
834 "frame number as argument" },
835 { "-split", FALSE, etTIME,
837 "Start writing new file when t MOD split = first "
839 { "-sep", FALSE, etBOOL,
841 "Write each frame to a separate .gro, .g96 or .pdb "
843 { "-nzero", FALSE, etINT,
845 "If the -sep flag is set, use these many digits "
846 "for the file numbers and prepend zeros as needed" },
847 { "-dropunder", FALSE, etREAL,
848 { &dropunder }, "Drop all frames below this value" },
849 { "-dropover", FALSE, etREAL,
850 { &dropover }, "Drop all frames above this value" },
851 { "-conect", FALSE, etBOOL,
853 "Add conect records when writing [TT].pdb[tt] files. Useful "
854 "for visualization of non-standard molecules, e.g. "
855 "coarse grained ones" }
857 #define NPA asize(pa)
860 t_trxstatus *trxout = NULL;
862 int ftp, ftpin = 0, file_nr;
863 t_trxframe fr, frout;
865 rvec *xmem = NULL, *vmem = NULL, *fmem = NULL;
866 rvec *xp = NULL, x_shift, hbox, box_center, dx;
867 real xtcpr, lambda, *w_rls = NULL;
868 int m, i, d, frame, outframe, natoms, nout, ncent, nre, newstep = 0, model_nr;
871 gmx_mtop_t *mtop = NULL;
872 gmx_conect gc = NULL;
874 t_atoms *atoms = NULL, useatoms;
876 atom_id *index, *cindex;
880 int ifit, irms, my_clust = -1;
881 atom_id *ind_fit, *ind_rms;
882 char *gn_fit, *gn_rms;
883 t_cluster_ndx *clust = NULL;
884 t_trxstatus **clust_status = NULL;
885 int *clust_status_id = NULL;
887 int *nfwritten = NULL;
888 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
890 real tshift = 0, t0 = -1, dt = 0.001, prec;
891 gmx_bool bFit, bFitXY, bPFit, bReset;
893 gmx_rmpbc_t gpbc = NULL;
894 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
895 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
896 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetPrec, bNeedPrec;
897 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
898 gmx_bool bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
899 gmx_bool bWriteFrame, bSplitHere;
900 const char *top_file, *in_file, *out_file = NULL;
901 char out_file2[256], *charpt;
902 char *outf_base = NULL;
903 const char *outf_ext = NULL;
904 char top_title[256], title[256], command[256], filemode[5];
906 gmx_bool bWarnCompact = FALSE;
911 { efTRX, "-f", NULL, ffREAD },
912 { efTRO, "-o", NULL, ffWRITE },
913 { efTPS, NULL, NULL, ffOPTRD },
914 { efNDX, NULL, NULL, ffOPTRD },
915 { efNDX, "-fr", "frames", ffOPTRD },
916 { efNDX, "-sub", "cluster", ffOPTRD },
917 { efXVG, "-drop", "drop", ffOPTRD }
919 #define NFILE asize(fnm)
921 if (!parse_common_args(&argc, argv,
922 PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
923 PCA_TIME_UNIT | PCA_BE_NICE,
924 NFILE, fnm, NPA, pa, asize(desc), desc,
930 top_file = ftp2fn(efTPS, NFILE, fnm);
933 /* Check command line */
934 in_file = opt2fn("-f", NFILE, fnm);
938 #ifndef GMX_NATIVE_WINDOWS
939 do_trunc(in_file, ttrunc);
944 /* mark active cmdline options */
945 bSetBox = opt2parg_bSet("-box", NPA, pa);
946 bSetTime = opt2parg_bSet("-t0", NPA, pa);
947 bSetPrec = opt2parg_bSet("-ndec", NPA, pa);
948 bSetUR = opt2parg_bSet("-ur", NPA, pa);
949 bExec = opt2parg_bSet("-exec", NPA, pa);
950 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
951 bTDump = opt2parg_bSet("-dump", NPA, pa);
952 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
953 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
954 bTrans = opt2parg_bSet("-trans", NPA, pa);
955 bSplit = (split_t != 0);
957 /* parse enum options */
958 fit_enum = nenum(fit);
959 bFit = (fit_enum == efFit || fit_enum == efFitXY);
960 bFitXY = fit_enum == efFitXY;
961 bReset = (fit_enum == efReset || fit_enum == efResetXY);
962 bPFit = fit_enum == efPFit;
963 pbc_enum = nenum(pbc_opt);
964 bPBCWhole = pbc_enum == epWhole;
965 bPBCcomRes = pbc_enum == epComRes;
966 bPBCcomMol = pbc_enum == epComMol;
967 bPBCcomAtom = pbc_enum == epComAtom;
968 bNoJump = pbc_enum == epNojump;
969 bCluster = pbc_enum == epCluster;
970 bPBC = pbc_enum != epNone;
971 unitcell_enum = nenum(unitcell_opt);
972 ecenter = nenum(center_opt) - ecTric;
974 /* set and check option dependencies */
977 bFit = TRUE; /* for pfit, fit *must* be set */
981 bReset = TRUE; /* for fit, reset *must* be set */
986 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
988 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
992 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
995 "WARNING: Option for unitcell representation (-ur %s)\n"
996 " only has effect in combination with -pbc %s, %s or %s.\n"
997 " Ingoring unitcell representation.\n\n",
998 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
1004 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
1005 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1006 "for the rotational fit.\n"
1007 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1011 /* ndec is in nr of decimal places, prec is a multiplication factor: */
1013 for (i = 0; i < ndec; i++)
1018 bIndex = ftp2bSet(efNDX, NFILE, fnm);
1021 /* Determine output type */
1022 out_file = opt2fn("-o", NFILE, fnm);
1023 ftp = fn2ftp(out_file);
1024 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1025 bNeedPrec = (ftp == efXTC || ftp == efGRO);
1028 /* check if velocities are possible in input and output files */
1029 ftpin = fn2ftp(in_file);
1030 bVels = (ftp == efTRR || ftp == efTRJ || ftp == efGRO || ftp == efG96)
1031 && (ftpin == efTRR || ftpin == efTRJ || ftpin == efGRO || ftpin == efG96 ||
1034 if (bSeparate || bSplit)
1036 outf_ext = strrchr(out_file, '.');
1037 if (outf_ext == NULL)
1039 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1041 outf_base = strdup(out_file);
1042 outf_base[outf_ext - out_file] = '\0';
1045 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1048 if ((ftp != efXTC) && (ftp != efTRR))
1050 /* It seems likely that other trajectory file types
1051 * could work here. */
1052 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1055 clust = cluster_index(NULL, opt2fn("-sub", NFILE, fnm));
1057 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1058 * number to check for. In my linux box it is only 16.
1060 if (0 && (clust->clust->nr > FOPEN_MAX-4))
1062 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1063 " trajectories.\ntry splitting the index file in %d parts.\n"
1065 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1067 gmx_warning("The -sub option could require as many open output files as there are\n"
1068 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1069 "try reducing the number of index groups in the file, and perhaps\n"
1070 "using trjconv -sub several times on different chunks of your index file.\n",
1073 snew(clust_status, clust->clust->nr);
1074 snew(clust_status_id, clust->clust->nr);
1075 snew(nfwritten, clust->clust->nr);
1076 for (i = 0; (i < clust->clust->nr); i++)
1078 clust_status[i] = NULL;
1079 clust_status_id[i] = -1;
1081 bSeparate = bSplit = FALSE;
1088 mtop = read_mtop_for_tng(top_file, in_file, out_file);
1090 /* Determine whether to read a topology */
1091 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1092 bRmPBC || bReset || bPBCcomMol || bCluster ||
1093 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1095 /* Determine if when can read index groups */
1096 bIndex = (bIndex || bTPS);
1100 read_tps_conf(top_file, top_title, &top, &ePBC, &xp, NULL, top_box,
1101 bReset || bPBCcomRes);
1104 if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1106 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1109 /* top_title is only used for gro and pdb,
1110 * the header in such a file is top_title t= ...
1111 * to prevent a double t=, remove it from top_title
1113 if ((charpt = strstr(top_title, " t= ")))
1120 gc = gmx_conect_generate(&top);
1124 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1128 /* get frame number index */
1130 if (opt2bSet("-fr", NFILE, fnm))
1132 printf("Select groups of frame number indices:\n");
1133 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (atom_id **)&frindex, &frname);
1136 for (i = 0; i < nrfri; i++)
1138 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1143 /* get index groups etc. */
1146 printf("Select group for %s fit\n",
1147 bFit ? "least squares" : "translational");
1148 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1149 1, &ifit, &ind_fit, &gn_fit);
1155 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1159 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1165 printf("Select group for clustering\n");
1166 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1167 1, &ifit, &ind_fit, &gn_fit);
1174 printf("Select group for centering\n");
1175 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1176 1, &ncent, &cindex, &grpnm);
1178 printf("Select group for output\n");
1179 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1180 1, &nout, &index, &grpnm);
1184 /* no index file, so read natoms from TRX */
1185 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1187 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1192 snew(index, natoms);
1193 for (i = 0; i < natoms; i++)
1207 snew(w_rls, atoms->nr);
1208 for (i = 0; (i < ifit); i++)
1210 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1213 /* Restore reference structure and set to origin,
1214 store original location (to put structure back) */
1217 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1219 copy_rvec(xp[index[0]], x_shift);
1220 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, NULL, xp, w_rls);
1221 rvec_dec(x_shift, xp[index[0]]);
1225 clear_rvec(x_shift);
1228 if (bDropUnder || bDropOver)
1230 /* Read the .xvg file with the drop values */
1231 fprintf(stderr, "\nReading drop file ...");
1232 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1233 fprintf(stderr, " %d time points\n", ndrop);
1234 if (ndrop == 0 || ncol < 2)
1236 gmx_fatal(FARGS, "Found no data points in %s",
1237 opt2fn("-drop", NFILE, fnm));
1243 /* Make atoms struct for output in GRO or PDB files */
1244 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1246 /* get memory for stuff to go in .pdb file */
1247 init_t_atoms(&useatoms, atoms->nr, FALSE);
1248 sfree(useatoms.resinfo);
1249 useatoms.resinfo = atoms->resinfo;
1250 for (i = 0; (i < nout); i++)
1252 useatoms.atomname[i] = atoms->atomname[index[i]];
1253 useatoms.atom[i] = atoms->atom[index[i]];
1254 useatoms.nres = max(useatoms.nres, useatoms.atom[i].resind+1);
1258 /* select what to read */
1259 if (ftp == efTRR || ftp == efTRJ)
1269 flags = flags | TRX_READ_V;
1273 flags = flags | TRX_READ_F;
1276 /* open trx file for reading */
1277 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1280 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1284 if (bSetPrec || !fr.bPrec)
1286 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1290 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1294 if (bHaveFirstFrame)
1296 set_trxframe_ePBC(&fr, ePBC);
1302 tshift = tzero-fr.time;
1312 /* check if index is meaningful */
1313 for (i = 0; i < nout; i++)
1315 if (index[i] >= natoms)
1318 "Index[%d] %d is larger than the number of atoms in the\n"
1319 "trajectory file (%d). There is a mismatch in the contents\n"
1320 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1322 bCopy = bCopy || (i != index[i]);
1326 /* open output for writing */
1327 strcpy(filemode, "w");
1331 trjtools_gmx_prepare_tng_writing(out_file,
1346 if (!bSplit && !bSubTraj)
1348 trxout = open_trx(out_file, filemode);
1354 if (( !bSeparate && !bSplit ) && !bSubTraj)
1356 out = ffopen(out_file, filemode);
1360 gmx_incons("Illegal output file format");
1378 fprintf(gmx_fio_getfp(trx_get_fileio(trxout)),
1379 "Generated by %s. #atoms=%d, a BOX is"
1380 " stored in this file.\n", ShortProgram(), nout);
1383 /* Start the big loop over frames */
1390 /* Main loop over frames */
1401 /*if (frame >= clust->clust->nra)
1402 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1403 if (frame > clust->maxframe)
1409 my_clust = clust->inv_clust[frame];
1411 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1412 (my_clust == NO_ATID))
1420 /* generate new box */
1422 for (m = 0; m < DIM; m++)
1424 fr.box[m][m] = newbox[m];
1430 for (i = 0; i < natoms; i++)
1432 rvec_inc(fr.x[i], trans);
1438 /* determine timestep */
1451 /* This is not very elegant, as one can not dump a frame after
1452 * a timestep with is more than twice as small as the first one. */
1453 bDumpFrame = (fr.time > tdump-0.5*dt) && (fr.time <= tdump+0.5*dt);
1460 /* determine if an atom jumped across the box and reset it if so */
1461 if (bNoJump && (bTPS || frame != 0))
1463 for (d = 0; d < DIM; d++)
1465 hbox[d] = 0.5*fr.box[d][d];
1467 for (i = 0; i < natoms; i++)
1471 rvec_dec(fr.x[i], x_shift);
1473 for (m = DIM-1; m >= 0; m--)
1477 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1479 for (d = 0; d <= m; d++)
1481 fr.x[i][d] += fr.box[m][d];
1484 while (fr.x[i][m]-xp[i][m] > hbox[m])
1486 for (d = 0; d <= m; d++)
1488 fr.x[i][d] -= fr.box[m][d];
1497 calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1502 /* Now modify the coords according to the flags,
1503 for normal fit, this is only done for output frames */
1506 gmx_rmpbc_trxfr(gpbc, &fr);
1509 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1510 do_fit(natoms, w_rls, xp, fr.x);
1513 /* store this set of coordinates for future use */
1514 if (bPFit || bNoJump)
1520 for (i = 0; (i < natoms); i++)
1522 copy_rvec(fr.x[i], xp[i]);
1523 rvec_inc(fr.x[i], x_shift);
1529 /* see if we have a frame from the frame index group */
1530 for (i = 0; i < nrfri && !bDumpFrame; i++)
1532 bDumpFrame = frame == frindex[i];
1535 if (debug && bDumpFrame)
1537 fprintf(debug, "dumping %d\n", frame);
1541 ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1543 if (bWriteFrame && (bDropUnder || bDropOver))
1545 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1550 if (fabs(dropval[0][drop0] - fr.time)
1551 < fabs(dropval[0][drop1] - fr.time))
1559 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1560 (bDropOver && dropval[1][dropuse] > dropover))
1562 bWriteFrame = FALSE;
1572 fr.time = tzero+frame*timestep;
1582 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1583 output_env_conv_time(oenv, fr.time), output_env_get_time_unit(oenv));
1586 /* check for writing at each delta_t */
1587 bDoIt = (delta_t == 0);
1592 bDoIt = bRmod(fr.time, tzero, delta_t);
1596 /* round() is not C89 compatible, so we do this: */
1597 bDoIt = bRmod(floor(fr.time+0.5), floor(tzero+0.5),
1598 floor(delta_t+0.5));
1602 if (bDoIt || bTDump)
1604 /* print sometimes */
1605 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1607 fprintf(stderr, " -> frame %6d time %8.3f \r",
1608 outframe, output_env_conv_time(oenv, fr.time));
1613 /* Now modify the coords according to the flags,
1614 for PFit we did this already! */
1618 gmx_rmpbc_trxfr(gpbc, &fr);
1623 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1626 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1630 for (i = 0; i < natoms; i++)
1632 rvec_inc(fr.x[i], x_shift);
1639 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1645 switch (unitcell_enum)
1648 put_atoms_in_box(ePBC, fr.box, natoms, fr.x);
1651 put_atoms_in_triclinic_unitcell(ecenter, fr.box, natoms, fr.x);
1654 warn = put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1656 if (warn && !bWarnCompact)
1658 fprintf(stderr, "\n%s\n", warn);
1659 bWarnCompact = TRUE;
1666 put_residue_com_in_box(unitcell_enum, ecenter,
1667 natoms, atoms->atom, ePBC, fr.box, fr.x);
1671 put_molecule_com_in_box(unitcell_enum, ecenter,
1673 natoms, atoms->atom, ePBC, fr.box, fr.x);
1675 /* Copy the input trxframe struct to the output trxframe struct */
1677 frout.bV = (frout.bV && bVels);
1678 frout.bF = (frout.bF && bForce);
1679 frout.natoms = nout;
1680 if (bNeedPrec && (bSetPrec || !fr.bPrec))
1696 for (i = 0; i < nout; i++)
1698 copy_rvec(fr.x[index[i]], frout.x[i]);
1701 copy_rvec(fr.v[index[i]], frout.v[i]);
1705 copy_rvec(fr.f[index[i]], frout.f[i]);
1710 if (opt2parg_bSet("-shift", NPA, pa))
1712 for (i = 0; i < nout; i++)
1714 for (d = 0; d < DIM; d++)
1716 frout.x[i][d] += outframe*shift[d];
1723 bSplitHere = bSplit && bRmod(fr.time, tzero, split_t);
1727 /* round() is not C89 compatible, so we do this: */
1728 bSplitHere = bSplit && bRmod(floor(fr.time+0.5),
1730 floor(split_t+0.5));
1732 if (bSeparate || bSplitHere)
1734 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1740 write_tng_frame(trxout, &frout);
1741 // TODO when trjconv behaves better: work how to read and write lambda
1753 trxout = open_trx(out_file2, filemode);
1760 if (clust_status_id[my_clust] == -1)
1762 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1763 clust_status[my_clust] = open_trx(buf, "w");
1764 clust_status_id[my_clust] = 1;
1767 else if (clust_status_id[my_clust] == -2)
1769 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1770 clust->grpname[my_clust], ntrxopen, frame,
1773 write_trxframe(clust_status[my_clust], &frout, gc);
1774 nfwritten[my_clust]++;
1775 if (nfwritten[my_clust] ==
1776 (clust->clust->index[my_clust+1]-
1777 clust->clust->index[my_clust]))
1779 close_trx(clust_status[my_clust]);
1780 clust_status[my_clust] = NULL;
1781 clust_status_id[my_clust] = -2;
1785 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1792 write_trxframe(trxout, &frout, gc);
1798 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1799 top_title, fr.time);
1800 if (bSeparate || bSplitHere)
1802 out = ffopen(out_file2, "w");
1807 write_hconf_p(out, title, &useatoms, prec2ndec(frout.prec),
1808 frout.x, frout.bV ? frout.v : NULL, frout.box);
1811 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1812 sprintf(title, "%s t= %9.5f", top_title, frout.time);
1813 /* if reading from pdb, we want to keep the original
1814 model numbering else we write the output frame
1815 number plus one, because model 0 is not allowed in pdb */
1816 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1824 write_pdbfile(out, title, &useatoms, frout.x,
1825 frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1828 frout.title = title;
1829 if (bSeparate || bTDump)
1831 frout.bTitle = TRUE;
1834 frout.bAtoms = TRUE;
1836 frout.atoms = &useatoms;
1837 frout.bStep = FALSE;
1838 frout.bTime = FALSE;
1842 frout.bTitle = (outframe == 0);
1843 frout.bAtoms = FALSE;
1847 write_g96_conf(out, &frout, -1, NULL);
1856 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1858 if (bSeparate || bSplitHere)
1863 /* execute command */
1867 sprintf(c, "%s %d", exec_command, file_nr-1);
1868 /*fprintf(stderr,"Executing '%s'\n",c);*/
1869 #ifdef GMX_NO_SYSTEM
1870 printf("Warning-- No calls to system(3) supported on this platform.");
1871 printf("Warning-- Skipping execution of 'system(\"%s\")'.", c);
1875 gmx_fatal(FARGS, "Error executing command: %s", c);
1883 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1885 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1888 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1890 fprintf(stderr, "\nWARNING no output, "
1891 "last frame read at t=%g\n", fr.time);
1893 fprintf(stderr, "\n");
1900 gmx_rmpbc_done(gpbc);
1907 else if (out != NULL)
1913 for (i = 0; (i < clust->clust->nr); i++)
1915 if (clust_status_id[i] >= 0)
1917 close_trx(clust_status[i]);
1925 do_view(oenv, out_file, NULL);