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/commandline/pargs.h"
54 #include "gromacs/fileio/futil.h"
55 #include "gromacs/fileio/pdbio.h"
56 #include "gromacs/fileio/confio.h"
60 #include "gromacs/fileio/xtcio.h"
63 #include "gromacs/fileio/g87io.h"
74 euSel, euRect, euTric, euCompact, euNR
78 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
79 rvec x[], atom_id index[], matrix box)
81 int m, i, j, j0, j1, jj, ai, aj;
84 rvec dx, xtest, box_center;
85 int nmol, imol_center;
87 gmx_bool *bMol, *bTmp;
88 rvec *m_com, *m_shift;
96 calc_box_center(ecenter, box, box_center);
98 /* Initiate the pbc structure */
99 memset(&pbc, 0, sizeof(pbc));
100 set_pbc(&pbc, ePBC, box);
102 /* Convert atom index to molecular */
104 molind = top->mols.index;
110 snew(bTmp, top->atoms.nr);
112 for (i = 0; (i < nrefat); i++)
114 /* Mark all molecules in the index */
117 /* Binary search assuming the molecules are sorted */
122 if (ai < molind[j0+1])
126 else if (ai >= molind[j1])
133 if (ai < molind[jj+1])
145 /* Double check whether all atoms in all molecules that are marked are part
146 * of the cluster. Simultaneously compute the center of geometry.
148 min_dist2 = 10*sqr(trace(box));
151 for (i = 0; i < nmol; i++)
153 for (j = molind[i]; j < molind[i+1]; j++)
155 if (bMol[i] && !bTmp[j])
157 gmx_fatal(FARGS, "Molecule %d marked for clustering but not atom %d in it - check your index!", i+1, j+1);
159 else if (!bMol[i] && bTmp[j])
161 gmx_fatal(FARGS, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j+1, i+1);
165 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
168 pbc_dx(&pbc, x[j], x[j-1], dx);
169 rvec_add(x[j-1], dx, x[j]);
171 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
172 rvec_inc(m_com[i], x[j]);
177 /* Normalize center of geometry */
178 fac = 1.0/(molind[i+1]-molind[i]);
179 for (m = 0; (m < DIM); m++)
183 /* Determine which molecule is closest to the center of the box */
184 pbc_dx(&pbc, box_center, m_com[i], dx);
185 tmp_r2 = iprod(dx, dx);
187 if (tmp_r2 < min_dist2)
192 cluster[ncluster++] = i;
199 fprintf(stderr, "No molecules selected in the cluster\n");
202 else if (imol_center == -1)
204 fprintf(stderr, "No central molecules could be found\n");
209 added[nadded++] = imol_center;
210 bMol[imol_center] = FALSE;
212 while (nadded < ncluster)
214 /* Find min distance between cluster molecules and those remaining to be added */
215 min_dist2 = 10*sqr(trace(box));
218 /* Loop over added mols */
219 for (i = 0; i < nadded; i++)
222 /* Loop over all mols */
223 for (j = 0; j < ncluster; j++)
226 /* check those remaining to be added */
229 pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
230 tmp_r2 = iprod(dx, dx);
231 if (tmp_r2 < min_dist2)
241 /* Add the best molecule */
242 added[nadded++] = jmin;
244 /* Calculate the shift from the ai molecule */
245 pbc_dx(&pbc, m_com[jmin], m_com[imin], dx);
246 rvec_add(m_com[imin], dx, xtest);
247 rvec_sub(xtest, m_com[jmin], m_shift[jmin]);
248 rvec_inc(m_com[jmin], m_shift[jmin]);
250 for (j = molind[jmin]; j < molind[jmin+1]; j++)
252 rvec_inc(x[j], m_shift[jmin]);
254 fprintf(stdout, "\rClustering iteration %d of %d...", nadded, ncluster);
264 fprintf(stdout, "\n");
267 static void put_molecule_com_in_box(int unitcell_enum, int ecenter,
269 int natoms, t_atom atom[],
270 int ePBC, matrix box, rvec x[])
274 rvec com, new_com, shift, dx, box_center;
279 calc_box_center(ecenter, box, box_center);
280 set_pbc(&pbc, ePBC, box);
283 gmx_fatal(FARGS, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
285 for (i = 0; (i < mols->nr); i++)
290 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
293 for (d = 0; d < DIM; d++)
299 /* calculate final COM */
300 svmul(1.0/mtot, com, com);
302 /* check if COM is outside box */
303 copy_rvec(com, new_com);
304 switch (unitcell_enum)
307 put_atoms_in_box(ePBC, box, 1, &new_com);
310 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
313 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
316 rvec_sub(new_com, com, shift);
317 if (norm2(shift) > 0)
321 fprintf(debug, "\nShifting position of molecule %d "
322 "by %8.3f %8.3f %8.3f\n", i+1,
323 shift[XX], shift[YY], shift[ZZ]);
325 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
327 rvec_inc(x[j], shift);
333 static void put_residue_com_in_box(int unitcell_enum, int ecenter,
334 int natoms, t_atom atom[],
335 int ePBC, matrix box, rvec x[])
337 atom_id i, j, res_start, res_end, res_nat;
341 rvec box_center, com, new_com, shift;
343 calc_box_center(ecenter, box, box_center);
349 for (i = 0; i < natoms+1; i++)
351 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
353 /* calculate final COM */
355 res_nat = res_end - res_start;
356 svmul(1.0/mtot, com, com);
358 /* check if COM is outside box */
359 copy_rvec(com, new_com);
360 switch (unitcell_enum)
363 put_atoms_in_box(ePBC, box, 1, &new_com);
366 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
369 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
372 rvec_sub(new_com, com, shift);
377 fprintf(debug, "\nShifting position of residue %d (atoms %u-%u) "
378 "by %g,%g,%g\n", atom[res_start].resind+1,
379 res_start+1, res_end+1, shift[XX], shift[YY], shift[ZZ]);
381 for (j = res_start; j < res_end; j++)
383 rvec_inc(x[j], shift);
389 /* remember start of new residue */
396 for (d = 0; d < DIM; d++)
402 presnr = atom[i].resind;
407 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, atom_id ci[])
410 rvec cmin, cmax, box_center, dx;
414 copy_rvec(x[ci[0]], cmin);
415 copy_rvec(x[ci[0]], cmax);
416 for (i = 0; i < nc; i++)
419 for (m = 0; m < DIM; m++)
421 if (x[ai][m] < cmin[m])
425 else if (x[ai][m] > cmax[m])
431 calc_box_center(ecenter, box, box_center);
432 for (m = 0; m < DIM; m++)
434 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
437 for (i = 0; i < n; i++)
444 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
450 strcpy(out_file, base);
461 strncat(out_file, "00000000000", ndigit-nd);
463 sprintf(nbuf, "%d.", file_nr);
464 strcat(out_file, nbuf);
465 strcat(out_file, ext);
468 void check_trn(const char *fn)
470 if ((fn2ftp(fn) != efTRJ) && (fn2ftp(fn) != efTRR))
472 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
476 #ifndef GMX_NATIVE_WINDOWS
477 void do_trunc(const char *fn, real t0)
490 gmx_fatal(FARGS, "You forgot to set the truncation time");
493 /* Check whether this is a .trj file */
496 in = open_trn(fn, "r");
497 fp = gmx_fio_getfp(in);
500 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
506 fpos = gmx_fio_ftell(in);
508 while (!bStop && fread_trnheader(in, &sh, &bOK))
510 fread_htrn(in, &sh, NULL, NULL, NULL, NULL);
511 fpos = gmx_ftell(fp);
515 gmx_fseek(fp, fpos, SEEK_SET);
521 fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
522 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
523 fn, j, t, (long int)fpos);
524 if (1 != scanf("%s", yesno))
526 gmx_fatal(FARGS, "Error reading user input");
528 if (strcmp(yesno, "YES") == 0)
530 fprintf(stderr, "Once again, I'm gonna DO this...\n");
532 if (0 != truncate(fn, fpos))
534 gmx_fatal(FARGS, "Error truncating file %s", fn);
539 fprintf(stderr, "Ok, I'll forget about it\n");
544 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
551 int gmx_trjconv(int argc, char *argv[])
553 const char *desc[] = {
554 "[THISMODULE] can convert trajectory files in many ways:[BR]",
555 "[BB]1.[bb] from one format to another[BR]",
556 "[BB]2.[bb] select a subset of atoms[BR]",
557 "[BB]3.[bb] change the periodicity representation[BR]",
558 "[BB]4.[bb] keep multimeric molecules together[BR]",
559 "[BB]5.[bb] center atoms in the box[BR]",
560 "[BB]6.[bb] fit atoms to reference structure[BR]",
561 "[BB]7.[bb] reduce the number of frames[BR]",
562 "[BB]8.[bb] change the timestamps of the frames ",
563 "([TT]-t0[tt] and [TT]-timestep[tt])[BR]",
564 "[BB]9.[bb] cut the trajectory in small subtrajectories according",
565 "to information in an index file. This allows subsequent analysis of",
566 "the subtrajectories that could, for example, be the result of a",
567 "cluster analysis. Use option [TT]-sub[tt].",
568 "This assumes that the entries in the index file are frame numbers and",
569 "dumps each group in the index file to a separate trajectory file.[BR]",
570 "[BB]10.[bb] select frames within a certain range of a quantity given",
571 "in an [TT].xvg[tt] file.[PAR]",
573 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
576 "Currently seven formats are supported for input and output:",
577 "[TT].xtc[tt], [TT].trr[tt], [TT].trj[tt], [TT].gro[tt], [TT].g96[tt],",
578 "[TT].pdb[tt] and [TT].g87[tt].",
579 "The file formats are detected from the file extension.",
580 "The precision of [TT].xtc[tt] and [TT].gro[tt] output is taken from the",
581 "input file for [TT].xtc[tt], [TT].gro[tt] and [TT].pdb[tt],",
582 "and from the [TT]-ndec[tt] option for other input formats. The precision",
583 "is always taken from [TT]-ndec[tt], when this option is set.",
584 "All other formats have fixed precision. [TT].trr[tt] and [TT].trj[tt]",
585 "output can be single or double precision, depending on the precision",
586 "of the [THISMODULE] binary.",
587 "Note that velocities are only supported in",
588 "[TT].trr[tt], [TT].trj[tt], [TT].gro[tt] and [TT].g96[tt] files.[PAR]",
590 "Option [TT]-app[tt] can be used to",
591 "append output to an existing trajectory file.",
592 "No checks are performed to ensure integrity",
593 "of the resulting combined trajectory file.[PAR]",
595 "Option [TT]-sep[tt] can be used to write every frame to a separate",
596 "[TT].gro, .g96[tt] or [TT].pdb[tt] file. By default, all frames all written to one file.",
597 "[TT].pdb[tt] files with all frames concatenated can be viewed with",
598 "[TT]rasmol -nmrpdb[tt].[PAR]",
600 "It is possible to select part of your trajectory and write it out",
601 "to a new trajectory file in order to save disk space, e.g. for leaving",
602 "out the water from a trajectory of a protein in water.",
603 "[BB]ALWAYS[bb] put the original trajectory on tape!",
604 "We recommend to use the portable [TT].xtc[tt] format for your analysis",
605 "to save disk space and to have portable files.[PAR]",
607 "There are two options for fitting the trajectory to a reference",
608 "either for essential dynamics analysis, etc.",
609 "The first option is just plain fitting to a reference structure",
610 "in the structure file. The second option is a progressive fit",
611 "in which the first timeframe is fitted to the reference structure ",
612 "in the structure file to obtain and each subsequent timeframe is ",
613 "fitted to the previously fitted structure. This way a continuous",
614 "trajectory is generated, which might not be the case when using the",
615 "regular fit method, e.g. when your protein undergoes large",
616 "conformational transitions.[PAR]",
618 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
620 "[TT]* mol[tt] puts the center of mass of molecules in the box,",
621 "and requires a run input file to be supplied with [TT]-s[tt].[BR]",
622 "[TT]* res[tt] puts the center of mass of residues in the box.[BR]",
623 "[TT]* atom[tt] puts all the atoms in the box.[BR]",
624 "[TT]* nojump[tt] checks if atoms jump across the box and then puts",
625 "them back. This has the effect that all molecules",
626 "will remain whole (provided they were whole in the initial",
627 "conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
628 "molecules may diffuse out of the box. The starting configuration",
629 "for this procedure is taken from the structure file, if one is",
630 "supplied, otherwise it is the first frame.[BR]",
631 "[TT]* cluster[tt] clusters all the atoms in the selected index",
632 "such that they are all closest to the center of mass of the cluster,",
633 "which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
634 "results if you in fact have a cluster. Luckily that can be checked",
635 "afterwards using a trajectory viewer. Note also that if your molecules",
636 "are broken this will not work either.[BR]",
637 "The separate option [TT]-clustercenter[tt] can be used to specify an",
638 "approximate center for the cluster. This is useful e.g. if you have",
639 "two big vesicles, and you want to maintain their relative positions.[BR]",
640 "[TT]* whole[tt] only makes broken molecules whole.[PAR]",
642 "Option [TT]-ur[tt] sets the unit cell representation for options",
643 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
644 "All three options give different results for triclinic boxes and",
645 "identical results for rectangular boxes.",
646 "[TT]rect[tt] is the ordinary brick shape.",
647 "[TT]tric[tt] is the triclinic unit cell.",
648 "[TT]compact[tt] puts all atoms at the closest distance from the center",
649 "of the box. This can be useful for visualizing e.g. truncated octahedra",
650 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
651 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
652 "is set differently.[PAR]",
654 "Option [TT]-center[tt] centers the system in the box. The user can",
655 "select the group which is used to determine the geometrical center.",
656 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
657 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
658 "[TT]tric[tt]: half of the sum of the box vectors,",
659 "[TT]rect[tt]: half of the box diagonal,",
660 "[TT]zero[tt]: zero.",
661 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
662 "want all molecules in the box after the centering.[PAR]",
664 "It is not always possible to use combinations of [TT]-pbc[tt],",
665 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
666 "you want in one call to [THISMODULE]. Consider using multiple",
667 "calls, and check out the GROMACS website for suggestions.[PAR]",
669 "With [TT]-dt[tt], it is possible to reduce the number of ",
670 "frames in the output. This option relies on the accuracy of the times",
671 "in your input trajectory, so if these are inaccurate use the",
672 "[TT]-timestep[tt] option to modify the time (this can be done",
673 "simultaneously). For making smooth movies, the program [gmx-filter]",
674 "can reduce the number of frames while using low-pass frequency",
675 "filtering, this reduces aliasing of high frequency motions.[PAR]",
677 "Using [TT]-trunc[tt] [THISMODULE] can truncate [TT].trj[tt] in place, i.e.",
678 "without copying the file. This is useful when a run has crashed",
679 "during disk I/O (i.e. full disk), or when two contiguous",
680 "trajectories must be concatenated without having double frames.[PAR]",
682 "Option [TT]-dump[tt] can be used to extract a frame at or near",
683 "one specific time from your trajectory.[PAR]",
685 "Option [TT]-drop[tt] reads an [TT].xvg[tt] file with times and values.",
686 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
687 "frames with a value below and above the value of the respective options",
688 "will not be written."
704 const char *pbc_opt[epNR + 1] =
706 NULL, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
711 const char *unitcell_opt[euNR+1] =
712 { NULL, "rect", "tric", "compact", NULL };
716 ecSel, ecTric, ecRect, ecZero, ecNR
718 const char *center_opt[ecNR+1] =
719 { NULL, "tric", "rect", "zero", NULL };
725 efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
727 const char *fit[efNR + 1] =
729 NULL, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
733 static gmx_bool bAppend = FALSE, bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
734 static gmx_bool bCenter = FALSE;
735 static int skip_nr = 1, ndec = 3, nzero = 0;
736 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
737 static rvec newbox = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
738 static char *exec_command = NULL;
739 static real dropunder = 0, dropover = 0;
740 static gmx_bool bRound = FALSE;
745 { "-skip", FALSE, etINT,
746 { &skip_nr }, "Only write every nr-th frame" },
747 { "-dt", FALSE, etTIME,
749 "Only write frame when t MOD dt = first time (%t)" },
750 { "-round", FALSE, etBOOL,
751 { &bRound }, "Round measurements to nearest picosecond"},
752 { "-dump", FALSE, etTIME,
753 { &tdump }, "Dump frame nearest specified time (%t)" },
754 { "-t0", FALSE, etTIME,
756 "Starting time (%t) (default: don't change)" },
757 { "-timestep", FALSE, etTIME,
759 "Change time step between input frames (%t)" },
760 { "-pbc", FALSE, etENUM,
762 "PBC treatment (see help text for full description)" },
763 { "-ur", FALSE, etENUM,
764 { unitcell_opt }, "Unit-cell representation" },
765 { "-center", FALSE, etBOOL,
766 { &bCenter }, "Center atoms in box" },
767 { "-boxcenter", FALSE, etENUM,
768 { center_opt }, "Center for -pbc and -center" },
769 { "-box", FALSE, etRVEC,
771 "Size for new cubic box (default: read from input)" },
772 { "-trans", FALSE, etRVEC,
774 "All coordinates will be translated by trans. This "
775 "can advantageously be combined with -pbc mol -ur "
777 { "-shift", FALSE, etRVEC,
779 "All coordinates will be shifted by framenr*shift" },
780 { "-fit", FALSE, etENUM,
782 "Fit molecule to ref structure in the structure file" },
783 { "-ndec", FALSE, etINT,
785 "Precision for .xtc and .gro writing in number of "
787 { "-vel", FALSE, etBOOL,
788 { &bVels }, "Read and write velocities if possible" },
789 { "-force", FALSE, etBOOL,
790 { &bForce }, "Read and write forces if possible" },
791 #ifndef GMX_NATIVE_WINDOWS
792 { "-trunc", FALSE, etTIME,
794 "Truncate input trajectory file after this time (%t)" },
796 { "-exec", FALSE, etSTR,
798 "Execute command for every output frame with the "
799 "frame number as argument" },
800 { "-app", FALSE, etBOOL,
801 { &bAppend }, "Append output" },
802 { "-split", FALSE, etTIME,
804 "Start writing new file when t MOD split = first "
806 { "-sep", FALSE, etBOOL,
808 "Write each frame to a separate .gro, .g96 or .pdb "
810 { "-nzero", FALSE, etINT,
812 "If the -sep flag is set, use these many digits "
813 "for the file numbers and prepend zeros as needed" },
814 { "-dropunder", FALSE, etREAL,
815 { &dropunder }, "Drop all frames below this value" },
816 { "-dropover", FALSE, etREAL,
817 { &dropover }, "Drop all frames above this value" },
818 { "-conect", FALSE, etBOOL,
820 "Add conect records when writing [TT].pdb[tt] files. Useful "
821 "for visualization of non-standard molecules, e.g. "
822 "coarse grained ones" }
824 #define NPA asize(pa)
827 t_trxstatus *trxout = NULL;
829 int ftp, ftpin = 0, file_nr;
830 t_trxframe fr, frout;
832 rvec *xmem = NULL, *vmem = NULL, *fmem = NULL;
833 rvec *xp = NULL, x_shift, hbox, box_center, dx;
834 real xtcpr, lambda, *w_rls = NULL;
835 int m, i, d, frame, outframe, natoms, nout, ncent, nre, newstep = 0, model_nr;
838 gmx_conect gc = NULL;
840 t_atoms *atoms = NULL, useatoms;
842 atom_id *index, *cindex;
846 int ifit, irms, my_clust = -1;
847 atom_id *ind_fit, *ind_rms;
848 char *gn_fit, *gn_rms;
849 t_cluster_ndx *clust = NULL;
850 t_trxstatus **clust_status = NULL;
851 int *clust_status_id = NULL;
853 int *nfwritten = NULL;
854 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
856 real tshift = 0, t0 = -1, dt = 0.001, prec;
857 gmx_bool bFit, bFitXY, bPFit, bReset;
859 gmx_rmpbc_t gpbc = NULL;
860 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
861 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
862 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetPrec, bNeedPrec;
863 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
864 gmx_bool bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
865 gmx_bool bWriteFrame, bSplitHere;
866 const char *top_file, *in_file, *out_file = NULL;
867 char out_file2[256], *charpt;
868 char *outf_base = NULL;
869 const char *outf_ext = NULL;
870 char top_title[256], title[256], command[256], filemode[5];
872 gmx_bool bWarnCompact = FALSE;
877 { efTRX, "-f", NULL, ffREAD },
878 { efTRO, "-o", NULL, ffWRITE },
879 { efTPS, NULL, NULL, ffOPTRD },
880 { efNDX, NULL, NULL, ffOPTRD },
881 { efNDX, "-fr", "frames", ffOPTRD },
882 { efNDX, "-sub", "cluster", ffOPTRD },
883 { efXVG, "-drop", "drop", ffOPTRD }
885 #define NFILE asize(fnm)
887 if (!parse_common_args(&argc, argv,
888 PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW |
889 PCA_TIME_UNIT | PCA_BE_NICE,
890 NFILE, fnm, NPA, pa, asize(desc), desc,
896 top_file = ftp2fn(efTPS, NFILE, fnm);
899 /* Check command line */
900 in_file = opt2fn("-f", NFILE, fnm);
904 #ifndef GMX_NATIVE_WINDOWS
905 do_trunc(in_file, ttrunc);
910 /* mark active cmdline options */
911 bSetBox = opt2parg_bSet("-box", NPA, pa);
912 bSetTime = opt2parg_bSet("-t0", NPA, pa);
913 bSetPrec = opt2parg_bSet("-ndec", NPA, pa);
914 bSetUR = opt2parg_bSet("-ur", NPA, pa);
915 bExec = opt2parg_bSet("-exec", NPA, pa);
916 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
917 bTDump = opt2parg_bSet("-dump", NPA, pa);
918 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
919 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
920 bTrans = opt2parg_bSet("-trans", NPA, pa);
921 bSplit = (split_t != 0);
923 /* parse enum options */
924 fit_enum = nenum(fit);
925 bFit = (fit_enum == efFit || fit_enum == efFitXY);
926 bFitXY = fit_enum == efFitXY;
927 bReset = (fit_enum == efReset || fit_enum == efResetXY);
928 bPFit = fit_enum == efPFit;
929 pbc_enum = nenum(pbc_opt);
930 bPBCWhole = pbc_enum == epWhole;
931 bPBCcomRes = pbc_enum == epComRes;
932 bPBCcomMol = pbc_enum == epComMol;
933 bPBCcomAtom = pbc_enum == epComAtom;
934 bNoJump = pbc_enum == epNojump;
935 bCluster = pbc_enum == epCluster;
936 bPBC = pbc_enum != epNone;
937 unitcell_enum = nenum(unitcell_opt);
938 ecenter = nenum(center_opt) - ecTric;
940 /* set and check option dependencies */
943 bFit = TRUE; /* for pfit, fit *must* be set */
947 bReset = TRUE; /* for fit, reset *must* be set */
952 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
954 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
958 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
961 "WARNING: Option for unitcell representation (-ur %s)\n"
962 " only has effect in combination with -pbc %s, %s or %s.\n"
963 " Ingoring unitcell representation.\n\n",
964 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
970 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
971 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
972 "for the rotational fit.\n"
973 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
977 /* ndec is in nr of decimal places, prec is a multiplication factor: */
979 for (i = 0; i < ndec; i++)
984 bIndex = ftp2bSet(efNDX, NFILE, fnm);
987 /* Determine output type */
988 out_file = opt2fn("-o", NFILE, fnm);
989 ftp = fn2ftp(out_file);
990 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
991 bNeedPrec = (ftp == efXTC || ftp == efGRO);
994 /* check if velocities are possible in input and output files */
995 ftpin = fn2ftp(in_file);
996 bVels = (ftp == efTRR || ftp == efTRJ || ftp == efGRO || ftp == efG96)
997 && (ftpin == efTRR || ftpin == efTRJ || ftpin == efGRO || ftpin == efG96 ||
1000 if (bSeparate || bSplit)
1002 outf_ext = strrchr(out_file, '.');
1003 if (outf_ext == NULL)
1005 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1007 outf_base = strdup(out_file);
1008 outf_base[outf_ext - out_file] = '\0';
1011 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1014 if ((ftp != efXTC) && (ftp != efTRR))
1016 /* It seems likely that other trajectory file types
1017 * could work here. */
1018 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1021 clust = cluster_index(NULL, opt2fn("-sub", NFILE, fnm));
1023 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1024 * number to check for. In my linux box it is only 16.
1026 if (0 && (clust->clust->nr > FOPEN_MAX-4))
1028 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1029 " trajectories.\ntry splitting the index file in %d parts.\n"
1031 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1033 gmx_warning("The -sub option could require as many open output files as there are\n"
1034 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1035 "try reducing the number of index groups in the file, and perhaps\n"
1036 "using trjconv -sub several times on different chunks of your index file.\n",
1039 snew(clust_status, clust->clust->nr);
1040 snew(clust_status_id, clust->clust->nr);
1041 snew(nfwritten, clust->clust->nr);
1042 for (i = 0; (i < clust->clust->nr); i++)
1044 clust_status[i] = NULL;
1045 clust_status_id[i] = -1;
1047 bSeparate = bSplit = FALSE;
1054 /* Determine whether to read a topology */
1055 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1056 bRmPBC || bReset || bPBCcomMol || bCluster ||
1057 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1059 /* Determine if when can read index groups */
1060 bIndex = (bIndex || bTPS);
1064 read_tps_conf(top_file, top_title, &top, &ePBC, &xp, NULL, top_box,
1065 bReset || bPBCcomRes);
1068 if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1070 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1073 /* top_title is only used for gro and pdb,
1074 * the header in such a file is top_title t= ...
1075 * to prevent a double t=, remove it from top_title
1077 if ((charpt = strstr(top_title, " t= ")))
1084 gc = gmx_conect_generate(&top);
1088 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1092 /* get frame number index */
1094 if (opt2bSet("-fr", NFILE, fnm))
1096 printf("Select groups of frame number indices:\n");
1097 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (atom_id **)&frindex, &frname);
1100 for (i = 0; i < nrfri; i++)
1102 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1107 /* get index groups etc. */
1110 printf("Select group for %s fit\n",
1111 bFit ? "least squares" : "translational");
1112 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1113 1, &ifit, &ind_fit, &gn_fit);
1119 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1123 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1129 printf("Select group for clustering\n");
1130 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1131 1, &ifit, &ind_fit, &gn_fit);
1138 printf("Select group for centering\n");
1139 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1140 1, &ncent, &cindex, &grpnm);
1142 printf("Select group for output\n");
1143 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1144 1, &nout, &index, &grpnm);
1148 /* no index file, so read natoms from TRX */
1149 if (!read_first_frame(oenv, &status, in_file, &fr, TRX_DONT_SKIP))
1151 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1156 snew(index, natoms);
1157 for (i = 0; i < natoms; i++)
1171 snew(w_rls, atoms->nr);
1172 for (i = 0; (i < ifit); i++)
1174 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1177 /* Restore reference structure and set to origin,
1178 store original location (to put structure back) */
1181 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1183 copy_rvec(xp[index[0]], x_shift);
1184 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, NULL, xp, w_rls);
1185 rvec_dec(x_shift, xp[index[0]]);
1189 clear_rvec(x_shift);
1192 if (bDropUnder || bDropOver)
1194 /* Read the .xvg file with the drop values */
1195 fprintf(stderr, "\nReading drop file ...");
1196 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1197 fprintf(stderr, " %d time points\n", ndrop);
1198 if (ndrop == 0 || ncol < 2)
1200 gmx_fatal(FARGS, "Found no data points in %s",
1201 opt2fn("-drop", NFILE, fnm));
1207 /* Make atoms struct for output in GRO or PDB files */
1208 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1210 /* get memory for stuff to go in .pdb file */
1211 init_t_atoms(&useatoms, atoms->nr, FALSE);
1212 sfree(useatoms.resinfo);
1213 useatoms.resinfo = atoms->resinfo;
1214 for (i = 0; (i < nout); i++)
1216 useatoms.atomname[i] = atoms->atomname[index[i]];
1217 useatoms.atom[i] = atoms->atom[index[i]];
1218 useatoms.nres = max(useatoms.nres, useatoms.atom[i].resind+1);
1222 /* select what to read */
1223 if (ftp == efTRR || ftp == efTRJ)
1233 flags = flags | TRX_READ_V;
1237 flags = flags | TRX_READ_F;
1240 /* open trx file for reading */
1241 bHaveFirstFrame = read_first_frame(oenv, &status, in_file, &fr, flags);
1244 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1248 if (bSetPrec || !fr.bPrec)
1250 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1254 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1258 if (bHaveFirstFrame)
1260 set_trxframe_ePBC(&fr, ePBC);
1266 tshift = tzero-fr.time;
1273 /* open output for writing */
1274 if ((bAppend) && (gmx_fexist(out_file)))
1276 strcpy(filemode, "a");
1277 fprintf(stderr, "APPENDING to existing file %s\n", out_file);
1281 strcpy(filemode, "w");
1290 if (!bSplit && !bSubTraj)
1292 trxout = open_trx(out_file, filemode);
1298 if (( !bSeparate && !bSplit ) && !bSubTraj)
1300 out = ffopen(out_file, filemode);
1308 /* check if index is meaningful */
1309 for (i = 0; i < nout; i++)
1311 if (index[i] >= natoms)
1314 "Index[%d] %d is larger than the number of atoms in the\n"
1315 "trajectory file (%d). There is a mismatch in the contents\n"
1316 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1318 bCopy = bCopy || (i != index[i]);
1336 fprintf(gmx_fio_getfp(trx_get_fileio(trxout)),
1337 "Generated by %s. #atoms=%d, a BOX is"
1338 " stored in this file.\n", ShortProgram(), nout);
1341 /* Start the big loop over frames */
1348 /* Main loop over frames */
1359 /*if (frame >= clust->clust->nra)
1360 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1361 if (frame > clust->maxframe)
1367 my_clust = clust->inv_clust[frame];
1369 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1370 (my_clust == NO_ATID))
1378 /* generate new box */
1380 for (m = 0; m < DIM; m++)
1382 fr.box[m][m] = newbox[m];
1388 for (i = 0; i < natoms; i++)
1390 rvec_inc(fr.x[i], trans);
1396 /* determine timestep */
1409 /* This is not very elegant, as one can not dump a frame after
1410 * a timestep with is more than twice as small as the first one. */
1411 bDumpFrame = (fr.time > tdump-0.5*dt) && (fr.time <= tdump+0.5*dt);
1418 /* determine if an atom jumped across the box and reset it if so */
1419 if (bNoJump && (bTPS || frame != 0))
1421 for (d = 0; d < DIM; d++)
1423 hbox[d] = 0.5*fr.box[d][d];
1425 for (i = 0; i < natoms; i++)
1429 rvec_dec(fr.x[i], x_shift);
1431 for (m = DIM-1; m >= 0; m--)
1435 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1437 for (d = 0; d <= m; d++)
1439 fr.x[i][d] += fr.box[m][d];
1442 while (fr.x[i][m]-xp[i][m] > hbox[m])
1444 for (d = 0; d <= m; d++)
1446 fr.x[i][d] -= fr.box[m][d];
1455 calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1460 /* Now modify the coords according to the flags,
1461 for normal fit, this is only done for output frames */
1464 gmx_rmpbc_trxfr(gpbc, &fr);
1467 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1468 do_fit(natoms, w_rls, xp, fr.x);
1471 /* store this set of coordinates for future use */
1472 if (bPFit || bNoJump)
1478 for (i = 0; (i < natoms); i++)
1480 copy_rvec(fr.x[i], xp[i]);
1481 rvec_inc(fr.x[i], x_shift);
1487 /* see if we have a frame from the frame index group */
1488 for (i = 0; i < nrfri && !bDumpFrame; i++)
1490 bDumpFrame = frame == frindex[i];
1493 if (debug && bDumpFrame)
1495 fprintf(debug, "dumping %d\n", frame);
1499 ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1501 if (bWriteFrame && (bDropUnder || bDropOver))
1503 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1508 if (fabs(dropval[0][drop0] - fr.time)
1509 < fabs(dropval[0][drop1] - fr.time))
1517 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1518 (bDropOver && dropval[1][dropuse] > dropover))
1520 bWriteFrame = FALSE;
1530 fr.time = tzero+frame*timestep;
1540 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1541 output_env_conv_time(oenv, fr.time), output_env_get_time_unit(oenv));
1544 /* check for writing at each delta_t */
1545 bDoIt = (delta_t == 0);
1550 bDoIt = bRmod(fr.time, tzero, delta_t);
1554 /* round() is not C89 compatible, so we do this: */
1555 bDoIt = bRmod(floor(fr.time+0.5), floor(tzero+0.5),
1556 floor(delta_t+0.5));
1560 if (bDoIt || bTDump)
1562 /* print sometimes */
1563 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1565 fprintf(stderr, " -> frame %6d time %8.3f \r",
1566 outframe, output_env_conv_time(oenv, fr.time));
1571 /* Now modify the coords according to the flags,
1572 for PFit we did this already! */
1576 gmx_rmpbc_trxfr(gpbc, &fr);
1581 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1584 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1588 for (i = 0; i < natoms; i++)
1590 rvec_inc(fr.x[i], x_shift);
1597 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1603 switch (unitcell_enum)
1606 put_atoms_in_box(ePBC, fr.box, natoms, fr.x);
1609 put_atoms_in_triclinic_unitcell(ecenter, fr.box, natoms, fr.x);
1612 warn = put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1614 if (warn && !bWarnCompact)
1616 fprintf(stderr, "\n%s\n", warn);
1617 bWarnCompact = TRUE;
1624 put_residue_com_in_box(unitcell_enum, ecenter,
1625 natoms, atoms->atom, ePBC, fr.box, fr.x);
1629 put_molecule_com_in_box(unitcell_enum, ecenter,
1631 natoms, atoms->atom, ePBC, fr.box, fr.x);
1633 /* Copy the input trxframe struct to the output trxframe struct */
1635 frout.bV = (frout.bV && bVels);
1636 frout.bF = (frout.bF && bForce);
1637 frout.natoms = nout;
1638 if (bNeedPrec && (bSetPrec || !fr.bPrec))
1654 for (i = 0; i < nout; i++)
1656 copy_rvec(fr.x[index[i]], frout.x[i]);
1659 copy_rvec(fr.v[index[i]], frout.v[i]);
1663 copy_rvec(fr.f[index[i]], frout.f[i]);
1668 if (opt2parg_bSet("-shift", NPA, pa))
1670 for (i = 0; i < nout; i++)
1672 for (d = 0; d < DIM; d++)
1674 frout.x[i][d] += outframe*shift[d];
1681 bSplitHere = bSplit && bRmod(fr.time, tzero, split_t);
1685 /* round() is not C89 compatible, so we do this: */
1686 bSplitHere = bSplit && bRmod(floor(fr.time+0.5),
1688 floor(split_t+0.5));
1690 if (bSeparate || bSplitHere)
1692 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1707 trxout = open_trx(out_file2, filemode);
1714 if (clust_status_id[my_clust] == -1)
1716 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1717 clust_status[my_clust] = open_trx(buf, "w");
1718 clust_status_id[my_clust] = 1;
1721 else if (clust_status_id[my_clust] == -2)
1723 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1724 clust->grpname[my_clust], ntrxopen, frame,
1727 write_trxframe(clust_status[my_clust], &frout, gc);
1728 nfwritten[my_clust]++;
1729 if (nfwritten[my_clust] ==
1730 (clust->clust->index[my_clust+1]-
1731 clust->clust->index[my_clust]))
1733 close_trx(clust_status[my_clust]);
1734 clust_status[my_clust] = NULL;
1735 clust_status_id[my_clust] = -2;
1739 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1746 write_trxframe(trxout, &frout, gc);
1752 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1753 top_title, fr.time);
1754 if (bSeparate || bSplitHere)
1756 out = ffopen(out_file2, "w");
1761 write_hconf_p(out, title, &useatoms, prec2ndec(frout.prec),
1762 frout.x, frout.bV ? frout.v : NULL, frout.box);
1765 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1766 sprintf(title, "%s t= %9.5f", top_title, frout.time);
1767 /* if reading from pdb, we want to keep the original
1768 model numbering else we write the output frame
1769 number plus one, because model 0 is not allowed in pdb */
1770 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1778 write_pdbfile(out, title, &useatoms, frout.x,
1779 frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1782 frout.title = title;
1783 if (bSeparate || bTDump)
1785 frout.bTitle = TRUE;
1788 frout.bAtoms = TRUE;
1790 frout.atoms = &useatoms;
1791 frout.bStep = FALSE;
1792 frout.bTime = FALSE;
1796 frout.bTitle = (outframe == 0);
1797 frout.bAtoms = FALSE;
1801 write_g96_conf(out, &frout, -1, NULL);
1810 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1812 if (bSeparate || bSplitHere)
1817 /* execute command */
1821 sprintf(c, "%s %d", exec_command, file_nr-1);
1822 /*fprintf(stderr,"Executing '%s'\n",c);*/
1823 #ifdef GMX_NO_SYSTEM
1824 printf("Warning-- No calls to system(3) supported on this platform.");
1825 printf("Warning-- Skipping execution of 'system(\"%s\")'.", c);
1829 gmx_fatal(FARGS, "Error executing command: %s", c);
1837 bHaveNextFrame = read_next_frame(oenv, status, &fr);
1839 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1842 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1844 fprintf(stderr, "\nWARNING no output, "
1845 "last frame read at t=%g\n", fr.time);
1847 fprintf(stderr, "\n");
1854 gmx_rmpbc_done(gpbc);
1861 else if (out != NULL)
1867 for (i = 0; (i < clust->clust->nr); i++)
1869 if (clust_status_id[i] >= 0)
1871 close_trx(clust_status[i]);
1877 do_view(oenv, out_file, NULL);