3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
70 euSel, euRect, euTric, euCompact, euNR
74 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
75 rvec x[], atom_id index[],
76 rvec clust_com, matrix box, rvec clustercenter)
78 int m, i, j, j0, j1, jj, ai, aj;
81 rvec dx, xtest, box_center;
82 int nmol, imol_center;
84 gmx_bool *bMol, *bTmp;
85 rvec *m_com, *m_shift;
93 calc_box_center(ecenter, box, box_center);
95 /* Initiate the pbc structure */
96 memset(&pbc, 0, sizeof(pbc));
97 set_pbc(&pbc, ePBC, box);
99 /* Convert atom index to molecular */
101 molind = top->mols.index;
107 snew(bTmp, top->atoms.nr);
109 for (i = 0; (i < nrefat); i++)
111 /* Mark all molecules in the index */
114 /* Binary search assuming the molecules are sorted */
119 if (ai < molind[j0+1])
123 else if (ai >= molind[j1])
130 if (ai < molind[jj+1])
142 /* Double check whether all atoms in all molecules that are marked are part
143 * of the cluster. Simultaneously compute the center of geometry.
145 min_dist2 = 10*sqr(trace(box));
148 for (i = 0; i < nmol; i++)
150 for (j = molind[i]; j < molind[i+1]; j++)
152 if (bMol[i] && !bTmp[j])
154 gmx_fatal(FARGS, "Molecule %d marked for clustering but not atom %d in it - check your index!", i+1, j+1);
156 else if (!bMol[i] && bTmp[j])
158 gmx_fatal(FARGS, "Atom %d marked for clustering but not molecule %d - this is an internal error...", j+1, i+1);
162 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
165 pbc_dx(&pbc, x[j], x[j-1], dx);
166 rvec_add(x[j-1], dx, x[j]);
168 /* Compute center of geometry of molecule - m_com[i] was zeroed when we did snew() on it! */
169 rvec_inc(m_com[i], x[j]);
174 /* Normalize center of geometry */
175 fac = 1.0/(molind[i+1]-molind[i]);
176 for (m = 0; (m < DIM); m++)
180 /* Determine which molecule is closest to the center of the box */
181 pbc_dx(&pbc, box_center, m_com[i], dx);
182 tmp_r2 = iprod(dx, dx);
184 if (tmp_r2 < min_dist2)
189 cluster[ncluster++] = i;
196 fprintf(stderr, "No molecules selected in the cluster\n");
199 else if (imol_center == -1)
201 fprintf(stderr, "No central molecules could be found\n");
206 added[nadded++] = imol_center;
207 bMol[imol_center] = FALSE;
209 while (nadded < ncluster)
211 /* Find min distance between cluster molecules and those remaining to be added */
212 min_dist2 = 10*sqr(trace(box));
215 /* Loop over added mols */
216 for (i = 0; i < nadded; i++)
219 /* Loop over all mols */
220 for (j = 0; j < ncluster; j++)
223 /* check those remaining to be added */
226 pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
227 tmp_r2 = iprod(dx, dx);
228 if (tmp_r2 < min_dist2)
238 /* Add the best molecule */
239 added[nadded++] = jmin;
241 /* Calculate the shift from the ai molecule */
242 pbc_dx(&pbc, m_com[jmin], m_com[imin], dx);
243 rvec_add(m_com[imin], dx, xtest);
244 rvec_sub(xtest, m_com[jmin], m_shift[jmin]);
245 rvec_inc(m_com[jmin], m_shift[jmin]);
247 for (j = molind[jmin]; j < molind[jmin+1]; j++)
249 rvec_inc(x[j], m_shift[jmin]);
251 fprintf(stdout, "\rClustering iteration %d of %d...", nadded, ncluster);
261 fprintf(stdout, "\n");
264 static void put_molecule_com_in_box(int unitcell_enum, int ecenter,
266 int natoms, t_atom atom[],
267 int ePBC, matrix box, rvec x[])
271 rvec com, new_com, shift, dx, box_center;
276 calc_box_center(ecenter, box, box_center);
277 set_pbc(&pbc, ePBC, box);
280 gmx_fatal(FARGS, "There are no molecule descriptions. I need a .tpr file for this pbc option.");
282 for (i = 0; (i < mols->nr); i++)
287 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
290 for (d = 0; d < DIM; d++)
296 /* calculate final COM */
297 svmul(1.0/mtot, com, com);
299 /* check if COM is outside box */
300 copy_rvec(com, new_com);
301 switch (unitcell_enum)
304 put_atoms_in_box(ePBC, box, 1, &new_com);
307 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
310 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
313 rvec_sub(new_com, com, shift);
314 if (norm2(shift) > 0)
318 fprintf (debug, "\nShifting position of molecule %d "
319 "by %8.3f %8.3f %8.3f\n", i+1, PR_VEC(shift));
321 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
323 rvec_inc(x[j], shift);
329 static void put_residue_com_in_box(int unitcell_enum, int ecenter,
330 int natoms, t_atom atom[],
331 int ePBC, matrix box, rvec x[])
333 atom_id i, j, res_start, res_end, res_nat;
337 rvec box_center, com, new_com, shift;
339 calc_box_center(ecenter, box, box_center);
345 for (i = 0; i < natoms+1; i++)
347 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
349 /* calculate final COM */
351 res_nat = res_end - res_start;
352 svmul(1.0/mtot, com, com);
354 /* check if COM is outside box */
355 copy_rvec(com, new_com);
356 switch (unitcell_enum)
359 put_atoms_in_box(ePBC, box, 1, &new_com);
362 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
365 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
368 rvec_sub(new_com, com, shift);
373 fprintf (debug, "\nShifting position of residue %d (atoms %u-%u) "
374 "by %g,%g,%g\n", atom[res_start].resind+1,
375 res_start+1, res_end+1, PR_VEC(shift));
377 for (j = res_start; j < res_end; j++)
379 rvec_inc(x[j], shift);
385 /* remember start of new residue */
392 for (d = 0; d < DIM; d++)
398 presnr = atom[i].resind;
403 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, atom_id ci[])
406 rvec cmin, cmax, box_center, dx;
410 copy_rvec(x[ci[0]], cmin);
411 copy_rvec(x[ci[0]], cmax);
412 for (i = 0; i < nc; i++)
415 for (m = 0; m < DIM; m++)
417 if (x[ai][m] < cmin[m])
421 else if (x[ai][m] > cmax[m])
427 calc_box_center(ecenter, box, box_center);
428 for (m = 0; m < DIM; m++)
430 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
433 for (i = 0; i < n; i++)
440 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
446 strcpy(out_file, base);
457 strncat(out_file, "00000000000", ndigit-nd);
459 sprintf(nbuf, "%d.", file_nr);
460 strcat(out_file, nbuf);
461 strcat(out_file, ext);
464 void check_trn(const char *fn)
466 if ((fn2ftp(fn) != efTRJ) && (fn2ftp(fn) != efTRR))
468 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
472 #ifndef GMX_NATIVE_WINDOWS
473 void do_trunc(const char *fn, real t0)
486 gmx_fatal(FARGS, "You forgot to set the truncation time");
489 /* Check whether this is a .trj file */
492 in = open_trn(fn, "r");
493 fp = gmx_fio_getfp(in);
496 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
502 fpos = gmx_fio_ftell(in);
504 while (!bStop && fread_trnheader(in, &sh, &bOK))
506 fread_htrn(in, &sh, NULL, NULL, NULL, NULL);
507 fpos = gmx_ftell(fp);
511 gmx_fseek(fp, fpos, SEEK_SET);
517 fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
518 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
519 fn, j, t, (long int)fpos);
520 if (1 != scanf("%s", yesno))
522 gmx_fatal(FARGS, "Error reading user input");
524 if (strcmp(yesno, "YES") == 0)
526 fprintf(stderr, "Once again, I'm gonna DO this...\n");
528 if (0 != truncate(fn, fpos))
530 gmx_fatal(FARGS, "Error truncating file %s", fn);
535 fprintf(stderr, "Ok, I'll forget about it\n");
540 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
547 int gmx_trjconv(int argc, char *argv[])
549 const char *desc[] = {
550 "[TT]trjconv[tt] can convert trajectory files in many ways:[BR]",
551 "[BB]1.[bb] from one format to another[BR]",
552 "[BB]2.[bb] select a subset of atoms[BR]",
553 "[BB]3.[bb] change the periodicity representation[BR]",
554 "[BB]4.[bb] keep multimeric molecules together[BR]",
555 "[BB]5.[bb] center atoms in the box[BR]",
556 "[BB]6.[bb] fit atoms to reference structure[BR]",
557 "[BB]7.[bb] reduce the number of frames[BR]",
558 "[BB]8.[bb] change the timestamps of the frames ",
559 "([TT]-t0[tt] and [TT]-timestep[tt])[BR]",
560 "[BB]9.[bb] cut the trajectory in small subtrajectories according",
561 "to information in an index file. This allows subsequent analysis of",
562 "the subtrajectories that could, for example, be the result of a",
563 "cluster analysis. Use option [TT]-sub[tt].",
564 "This assumes that the entries in the index file are frame numbers and",
565 "dumps each group in the index file to a separate trajectory file.[BR]",
566 "[BB]10.[bb] select frames within a certain range of a quantity given",
567 "in an [TT].xvg[tt] file.[PAR]",
569 "The program [TT]trjcat[tt] is better suited for concatenating multiple trajectory files.",
572 "Currently seven formats are supported for input and output:",
573 "[TT].xtc[tt], [TT].trr[tt], [TT].trj[tt], [TT].gro[tt], [TT].g96[tt],",
574 "[TT].pdb[tt] and [TT].g87[tt].",
575 "The file formats are detected from the file extension.",
576 "The precision of [TT].xtc[tt] and [TT].gro[tt] output is taken from the",
577 "input file for [TT].xtc[tt], [TT].gro[tt] and [TT].pdb[tt],",
578 "and from the [TT]-ndec[tt] option for other input formats. The precision",
579 "is always taken from [TT]-ndec[tt], when this option is set.",
580 "All other formats have fixed precision. [TT].trr[tt] and [TT].trj[tt]",
581 "output can be single or double precision, depending on the precision",
582 "of the [TT]trjconv[tt] binary.",
583 "Note that velocities are only supported in",
584 "[TT].trr[tt], [TT].trj[tt], [TT].gro[tt] and [TT].g96[tt] files.[PAR]",
586 "Option [TT]-app[tt] can be used to",
587 "append output to an existing trajectory file.",
588 "No checks are performed to ensure integrity",
589 "of the resulting combined trajectory file.[PAR]",
591 "Option [TT]-sep[tt] can be used to write every frame to a separate",
592 "[TT].gro, .g96[tt] or [TT].pdb[tt] file. By default, all frames all written to one file.",
593 "[TT].pdb[tt] files with all frames concatenated can be viewed with",
594 "[TT]rasmol -nmrpdb[tt].[PAR]",
596 "It is possible to select part of your trajectory and write it out",
597 "to a new trajectory file in order to save disk space, e.g. for leaving",
598 "out the water from a trajectory of a protein in water.",
599 "[BB]ALWAYS[bb] put the original trajectory on tape!",
600 "We recommend to use the portable [TT].xtc[tt] format for your analysis",
601 "to save disk space and to have portable files.[PAR]",
603 "There are two options for fitting the trajectory to a reference",
604 "either for essential dynamics analysis, etc.",
605 "The first option is just plain fitting to a reference structure",
606 "in the structure file. The second option is a progressive fit",
607 "in which the first timeframe is fitted to the reference structure ",
608 "in the structure file to obtain and each subsequent timeframe is ",
609 "fitted to the previously fitted structure. This way a continuous",
610 "trajectory is generated, which might not be the case when using the",
611 "regular fit method, e.g. when your protein undergoes large",
612 "conformational transitions.[PAR]",
614 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
616 "[TT]* mol[tt] puts the center of mass of molecules in the box,",
617 "and requires a run input file to be supplied with [TT]-s[tt].[BR]",
618 "[TT]* res[tt] puts the center of mass of residues in the box.[BR]",
619 "[TT]* atom[tt] puts all the atoms in the box.[BR]",
620 "[TT]* nojump[tt] checks if atoms jump across the box and then puts",
621 "them back. This has the effect that all molecules",
622 "will remain whole (provided they were whole in the initial",
623 "conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
624 "molecules may diffuse out of the box. The starting configuration",
625 "for this procedure is taken from the structure file, if one is",
626 "supplied, otherwise it is the first frame.[BR]",
627 "[TT]* cluster[tt] clusters all the atoms in the selected index",
628 "such that they are all closest to the center of mass of the cluster,",
629 "which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
630 "results if you in fact have a cluster. Luckily that can be checked",
631 "afterwards using a trajectory viewer. Note also that if your molecules",
632 "are broken this will not work either.[BR]",
633 "The separate option [TT]-clustercenter[tt] can be used to specify an",
634 "approximate center for the cluster. This is useful e.g. if you have",
635 "two big vesicles, and you want to maintain their relative positions.[BR]",
636 "[TT]* whole[tt] only makes broken molecules whole.[PAR]",
638 "Option [TT]-ur[tt] sets the unit cell representation for options",
639 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
640 "All three options give different results for triclinic boxes and",
641 "identical results for rectangular boxes.",
642 "[TT]rect[tt] is the ordinary brick shape.",
643 "[TT]tric[tt] is the triclinic unit cell.",
644 "[TT]compact[tt] puts all atoms at the closest distance from the center",
645 "of the box. This can be useful for visualizing e.g. truncated octahedra",
646 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
647 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
648 "is set differently.[PAR]",
650 "Option [TT]-center[tt] centers the system in the box. The user can",
651 "select the group which is used to determine the geometrical center.",
652 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
653 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
654 "[TT]tric[tt]: half of the sum of the box vectors,",
655 "[TT]rect[tt]: half of the box diagonal,",
656 "[TT]zero[tt]: zero.",
657 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
658 "want all molecules in the box after the centering.[PAR]",
660 "It is not always possible to use combinations of [TT]-pbc[tt],",
661 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
662 "you want in one call to [TT]trjconv[tt]. Consider using multiple",
663 "calls, and check out the GROMACS website for suggestions.[PAR]",
665 "With [TT]-dt[tt], it is possible to reduce the number of ",
666 "frames in the output. This option relies on the accuracy of the times",
667 "in your input trajectory, so if these are inaccurate use the",
668 "[TT]-timestep[tt] option to modify the time (this can be done",
669 "simultaneously). For making smooth movies, the program [TT]g_filter[tt]",
670 "can reduce the number of frames while using low-pass frequency",
671 "filtering, this reduces aliasing of high frequency motions.[PAR]",
673 "Using [TT]-trunc[tt] [TT]trjconv[tt] can truncate [TT].trj[tt] in place, i.e.",
674 "without copying the file. This is useful when a run has crashed",
675 "during disk I/O (i.e. full disk), or when two contiguous",
676 "trajectories must be concatenated without having double frames.[PAR]",
678 "Option [TT]-dump[tt] can be used to extract a frame at or near",
679 "one specific time from your trajectory.[PAR]",
681 "Option [TT]-drop[tt] reads an [TT].xvg[tt] file with times and values.",
682 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
683 "frames with a value below and above the value of the respective options",
684 "will not be written."
700 const char *pbc_opt[epNR + 1] =
702 NULL, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
707 const char *unitcell_opt[euNR+1] =
708 { NULL, "rect", "tric", "compact", NULL };
712 ecSel, ecTric, ecRect, ecZero, ecNR
714 const char *center_opt[ecNR+1] =
715 { NULL, "tric", "rect", "zero", NULL };
721 efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
723 const char *fit[efNR + 1] =
725 NULL, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
729 static gmx_bool bAppend = FALSE, bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
730 static gmx_bool bCenter = FALSE;
731 static int skip_nr = 1, ndec = 3, nzero = 0;
732 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
733 static rvec newbox = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
734 static char *exec_command = NULL;
735 static real dropunder = 0, dropover = 0;
736 static gmx_bool bRound = FALSE;
737 static rvec clustercenter = {0, 0, 0};
742 { "-skip", FALSE, etINT,
743 { &skip_nr }, "Only write every nr-th frame" },
744 { "-dt", FALSE, etTIME,
746 "Only write frame when t MOD dt = first time (%t)" },
747 { "-round", FALSE, etBOOL,
748 { &bRound }, "Round measurements to nearest picosecond"},
749 { "-dump", FALSE, etTIME,
750 { &tdump }, "Dump frame nearest specified time (%t)" },
751 { "-t0", FALSE, etTIME,
753 "Starting time (%t) (default: don't change)" },
754 { "-timestep", FALSE, etTIME,
756 "Change time step between input frames (%t)" },
757 { "-pbc", FALSE, etENUM,
759 "PBC treatment (see help text for full description)" },
760 { "-ur", FALSE, etENUM,
761 { unitcell_opt }, "Unit-cell representation" },
762 { "-center", FALSE, etBOOL,
763 { &bCenter }, "Center atoms in box" },
764 { "-boxcenter", FALSE, etENUM,
765 { center_opt }, "Center for -pbc and -center" },
766 { "-box", FALSE, etRVEC,
768 "Size for new cubic box (default: read from input)" },
769 { "-clustercenter", FALSE, etRVEC,
771 "Optional starting point for pbc cluster option" },
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 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,
893 top_file = ftp2fn(efTPS, NFILE, fnm);
896 /* Check command line */
897 in_file = opt2fn("-f", NFILE, fnm);
901 #ifndef GMX_NATIVE_WINDOWS
902 do_trunc(in_file, ttrunc);
907 /* mark active cmdline options */
908 bSetBox = opt2parg_bSet("-box", NPA, pa);
909 bSetTime = opt2parg_bSet("-t0", NPA, pa);
910 bSetPrec = opt2parg_bSet("-ndec", NPA, pa);
911 bSetUR = opt2parg_bSet("-ur", NPA, pa);
912 bExec = opt2parg_bSet("-exec", NPA, pa);
913 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
914 bTDump = opt2parg_bSet("-dump", NPA, pa);
915 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
916 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
917 bTrans = opt2parg_bSet("-trans", NPA, pa);
918 bSplit = (split_t != 0);
920 /* parse enum options */
921 fit_enum = nenum(fit);
922 bFit = (fit_enum == efFit || fit_enum == efFitXY);
923 bFitXY = fit_enum == efFitXY;
924 bReset = (fit_enum == efReset || fit_enum == efResetXY);
925 bPFit = fit_enum == efPFit;
926 pbc_enum = nenum(pbc_opt);
927 bPBCWhole = pbc_enum == epWhole;
928 bPBCcomRes = pbc_enum == epComRes;
929 bPBCcomMol = pbc_enum == epComMol;
930 bPBCcomAtom = pbc_enum == epComAtom;
931 bNoJump = pbc_enum == epNojump;
932 bCluster = pbc_enum == epCluster;
933 bPBC = pbc_enum != epNone;
934 unitcell_enum = nenum(unitcell_opt);
935 ecenter = nenum(center_opt) - ecTric;
937 /* set and check option dependencies */
940 bFit = TRUE; /* for pfit, fit *must* be set */
944 bReset = TRUE; /* for fit, reset *must* be set */
949 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
951 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
955 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
958 "WARNING: Option for unitcell representation (-ur %s)\n"
959 " only has effect in combination with -pbc %s, %s or %s.\n"
960 " Ingoring unitcell representation.\n\n",
961 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
967 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
968 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
969 "for the rotational fit.\n"
970 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
974 /* ndec is in nr of decimal places, prec is a multiplication factor: */
976 for (i = 0; i < ndec; i++)
981 bIndex = ftp2bSet(efNDX, NFILE, fnm);
984 /* Determine output type */
985 out_file = opt2fn("-o", NFILE, fnm);
986 ftp = fn2ftp(out_file);
987 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
988 bNeedPrec = (ftp == efXTC || ftp == efGRO);
991 /* check if velocities are possible in input and output files */
992 ftpin = fn2ftp(in_file);
993 bVels = (ftp == efTRR || ftp == efTRJ || ftp == efGRO || ftp == efG96)
994 && (ftpin == efTRR || ftpin == efTRJ || ftpin == efGRO || ftpin == efG96 ||
997 if (bSeparate || bSplit)
999 outf_ext = strrchr(out_file, '.');
1000 if (outf_ext == NULL)
1002 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1004 outf_base = strdup(out_file);
1005 outf_base[outf_ext - out_file] = '\0';
1008 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1011 if ((ftp != efXTC) && (ftp != efTRR))
1013 /* It seems likely that other trajectory file types
1014 * could work here. */
1015 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1018 clust = cluster_index(NULL, opt2fn("-sub", NFILE, fnm));
1020 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1021 * number to check for. In my linux box it is only 16.
1023 if (0 && (clust->clust->nr > FOPEN_MAX-4))
1025 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1026 " trajectories.\ntry splitting the index file in %d parts.\n"
1028 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1030 gmx_warning("The -sub option could require as many open output files as there are\n"
1031 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1032 "try reducing the number of index groups in the file, and perhaps\n"
1033 "using trjconv -sub several times on different chunks of your index file.\n",
1036 snew(clust_status, clust->clust->nr);
1037 snew(clust_status_id, clust->clust->nr);
1038 snew(nfwritten, clust->clust->nr);
1039 for (i = 0; (i < clust->clust->nr); i++)
1041 clust_status[i] = NULL;
1042 clust_status_id[i] = -1;
1044 bSeparate = bSplit = FALSE;
1051 /* Determine whether to read a topology */
1052 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1053 bRmPBC || bReset || bPBCcomMol || bCluster ||
1054 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1056 /* Determine if when can read index groups */
1057 bIndex = (bIndex || bTPS);
1061 read_tps_conf(top_file, top_title, &top, &ePBC, &xp, NULL, top_box,
1062 bReset || bPBCcomRes);
1065 if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1067 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1070 /* top_title is only used for gro and pdb,
1071 * the header in such a file is top_title t= ...
1072 * to prevent a double t=, remove it from top_title
1074 if ((charpt = strstr(top_title, " t= ")))
1081 gc = gmx_conect_generate(&top);
1085 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr, top_box);
1089 /* get frame number index */
1091 if (opt2bSet("-fr", NFILE, fnm))
1093 printf("Select groups of frame number indices:\n");
1094 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (atom_id **)&frindex, &frname);
1097 for (i = 0; i < nrfri; i++)
1099 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1104 /* get index groups etc. */
1107 printf("Select group for %s fit\n",
1108 bFit ? "least squares" : "translational");
1109 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1110 1, &ifit, &ind_fit, &gn_fit);
1116 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1120 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1126 printf("Select group for clustering\n");
1127 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1128 1, &ifit, &ind_fit, &gn_fit);
1135 printf("Select group for centering\n");
1136 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1137 1, &ncent, &cindex, &grpnm);
1139 printf("Select group for output\n");
1140 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1141 1, &nout, &index, &grpnm);
1145 /* no index file, so read natoms from TRX */
1146 if (!read_first_frame(oenv, &status, in_file, &fr, TRX_DONT_SKIP))
1148 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1153 snew(index, natoms);
1154 for (i = 0; i < natoms; i++)
1168 snew(w_rls, atoms->nr);
1169 for (i = 0; (i < ifit); i++)
1171 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1174 /* Restore reference structure and set to origin,
1175 store original location (to put structure back) */
1178 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1180 copy_rvec(xp[index[0]], x_shift);
1181 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, NULL, xp, w_rls);
1182 rvec_dec(x_shift, xp[index[0]]);
1186 clear_rvec(x_shift);
1189 if (bDropUnder || bDropOver)
1191 /* Read the .xvg file with the drop values */
1192 fprintf(stderr, "\nReading drop file ...");
1193 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1194 fprintf(stderr, " %d time points\n", ndrop);
1195 if (ndrop == 0 || ncol < 2)
1197 gmx_fatal(FARGS, "Found no data points in %s",
1198 opt2fn("-drop", NFILE, fnm));
1204 /* Make atoms struct for output in GRO or PDB files */
1205 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1207 /* get memory for stuff to go in .pdb file */
1208 init_t_atoms(&useatoms, atoms->nr, FALSE);
1209 sfree(useatoms.resinfo);
1210 useatoms.resinfo = atoms->resinfo;
1211 for (i = 0; (i < nout); i++)
1213 useatoms.atomname[i] = atoms->atomname[index[i]];
1214 useatoms.atom[i] = atoms->atom[index[i]];
1215 useatoms.nres = max(useatoms.nres, useatoms.atom[i].resind+1);
1219 /* select what to read */
1220 if (ftp == efTRR || ftp == efTRJ)
1230 flags = flags | TRX_READ_V;
1234 flags = flags | TRX_READ_F;
1237 /* open trx file for reading */
1238 bHaveFirstFrame = read_first_frame(oenv, &status, in_file, &fr, flags);
1241 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1245 if (bSetPrec || !fr.bPrec)
1247 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1251 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1255 if (bHaveFirstFrame)
1257 set_trxframe_ePBC(&fr, ePBC);
1263 tshift = tzero-fr.time;
1270 /* open output for writing */
1271 if ((bAppend) && (gmx_fexist(out_file)))
1273 strcpy(filemode, "a");
1274 fprintf(stderr, "APPENDING to existing file %s\n", out_file);
1278 strcpy(filemode, "w");
1287 if (!bSplit && !bSubTraj)
1289 trxout = open_trx(out_file, filemode);
1295 if (( !bSeparate && !bSplit ) && !bSubTraj)
1297 out = ffopen(out_file, filemode);
1305 /* check if index is meaningful */
1306 for (i = 0; i < nout; i++)
1308 if (index[i] >= natoms)
1311 "Index[%d] %d is larger than the number of atoms in the\n"
1312 "trajectory file (%d). There is a mismatch in the contents\n"
1313 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1315 bCopy = bCopy || (i != index[i]);
1333 fprintf(gmx_fio_getfp(trx_get_fileio(trxout)),
1334 "Generated by %s. #atoms=%d, a BOX is"
1335 " stored in this file.\n", Program(), nout);
1338 /* Start the big loop over frames */
1345 /* Main loop over frames */
1356 /*if (frame >= clust->clust->nra)
1357 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1358 if (frame > clust->maxframe)
1364 my_clust = clust->inv_clust[frame];
1366 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1367 (my_clust == NO_ATID))
1375 /* generate new box */
1377 for (m = 0; m < DIM; m++)
1379 fr.box[m][m] = newbox[m];
1385 for (i = 0; i < natoms; i++)
1387 rvec_inc(fr.x[i], trans);
1393 /* determine timestep */
1406 /* This is not very elegant, as one can not dump a frame after
1407 * a timestep with is more than twice as small as the first one. */
1408 bDumpFrame = (fr.time > tdump-0.5*dt) && (fr.time <= tdump+0.5*dt);
1415 /* determine if an atom jumped across the box and reset it if so */
1416 if (bNoJump && (bTPS || frame != 0))
1418 for (d = 0; d < DIM; d++)
1420 hbox[d] = 0.5*fr.box[d][d];
1422 for (i = 0; i < natoms; i++)
1426 rvec_dec(fr.x[i], x_shift);
1428 for (m = DIM-1; m >= 0; m--)
1432 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1434 for (d = 0; d <= m; d++)
1436 fr.x[i][d] += fr.box[m][d];
1439 while (fr.x[i][m]-xp[i][m] > hbox[m])
1441 for (d = 0; d <= m; d++)
1443 fr.x[i][d] -= fr.box[m][d];
1454 calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, com, fr.box, clustercenter);
1459 /* Now modify the coords according to the flags,
1460 for normal fit, this is only done for output frames */
1463 gmx_rmpbc_trxfr(gpbc, &fr);
1466 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1467 do_fit(natoms, w_rls, xp, fr.x);
1470 /* store this set of coordinates for future use */
1471 if (bPFit || bNoJump)
1477 for (i = 0; (i < natoms); i++)
1479 copy_rvec(fr.x[i], xp[i]);
1480 rvec_inc(fr.x[i], x_shift);
1486 /* see if we have a frame from the frame index group */
1487 for (i = 0; i < nrfri && !bDumpFrame; i++)
1489 bDumpFrame = frame == frindex[i];
1492 if (debug && bDumpFrame)
1494 fprintf(debug, "dumping %d\n", frame);
1498 ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1500 if (bWriteFrame && (bDropUnder || bDropOver))
1502 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1507 if (fabs(dropval[0][drop0] - fr.time)
1508 < fabs(dropval[0][drop1] - fr.time))
1516 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1517 (bDropOver && dropval[1][dropuse] > dropover))
1519 bWriteFrame = FALSE;
1529 fr.time = tzero+frame*timestep;
1539 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1540 output_env_conv_time(oenv, fr.time), output_env_get_time_unit(oenv));
1543 /* check for writing at each delta_t */
1544 bDoIt = (delta_t == 0);
1549 bDoIt = bRmod(fr.time, tzero, delta_t);
1553 /* round() is not C89 compatible, so we do this: */
1554 bDoIt = bRmod(floor(fr.time+0.5), floor(tzero+0.5),
1555 floor(delta_t+0.5));
1559 if (bDoIt || bTDump)
1561 /* print sometimes */
1562 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1564 fprintf(stderr, " -> frame %6d time %8.3f \r",
1565 outframe, output_env_conv_time(oenv, fr.time));
1570 /* Now modify the coords according to the flags,
1571 for PFit we did this already! */
1575 gmx_rmpbc_trxfr(gpbc, &fr);
1580 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1583 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1587 for (i = 0; i < natoms; i++)
1589 rvec_inc(fr.x[i], x_shift);
1596 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1602 switch (unitcell_enum)
1605 put_atoms_in_box(ePBC, fr.box, natoms, fr.x);
1608 put_atoms_in_triclinic_unitcell(ecenter, fr.box, natoms, fr.x);
1611 warn = put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1613 if (warn && !bWarnCompact)
1615 fprintf(stderr, "\n%s\n", warn);
1616 bWarnCompact = TRUE;
1623 put_residue_com_in_box(unitcell_enum, ecenter,
1624 natoms, atoms->atom, ePBC, fr.box, fr.x);
1628 put_molecule_com_in_box(unitcell_enum, ecenter,
1630 natoms, atoms->atom, ePBC, fr.box, fr.x);
1632 /* Copy the input trxframe struct to the output trxframe struct */
1634 frout.bV = (frout.bV && bVels);
1635 frout.bF = (frout.bF && bForce);
1636 frout.natoms = nout;
1637 if (bNeedPrec && (bSetPrec || !fr.bPrec))
1653 for (i = 0; i < nout; i++)
1655 copy_rvec(fr.x[index[i]], frout.x[i]);
1658 copy_rvec(fr.v[index[i]], frout.v[i]);
1662 copy_rvec(fr.f[index[i]], frout.f[i]);
1667 if (opt2parg_bSet("-shift", NPA, pa))
1669 for (i = 0; i < nout; i++)
1671 for (d = 0; d < DIM; d++)
1673 frout.x[i][d] += outframe*shift[d];
1680 bSplitHere = bSplit && bRmod(fr.time, tzero, split_t);
1684 /* round() is not C89 compatible, so we do this: */
1685 bSplitHere = bSplit && bRmod(floor(fr.time+0.5),
1687 floor(split_t+0.5));
1689 if (bSeparate || bSplitHere)
1691 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1706 trxout = open_trx(out_file2, filemode);
1713 if (clust_status_id[my_clust] == -1)
1715 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1716 clust_status[my_clust] = open_trx(buf, "w");
1717 clust_status_id[my_clust] = 1;
1720 else if (clust_status_id[my_clust] == -2)
1722 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1723 clust->grpname[my_clust], ntrxopen, frame,
1726 write_trxframe(clust_status[my_clust], &frout, gc);
1727 nfwritten[my_clust]++;
1728 if (nfwritten[my_clust] ==
1729 (clust->clust->index[my_clust+1]-
1730 clust->clust->index[my_clust]))
1732 close_trx(clust_status[my_clust]);
1733 clust_status[my_clust] = NULL;
1734 clust_status_id[my_clust] = -2;
1738 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1745 write_trxframe(trxout, &frout, gc);
1751 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1752 top_title, fr.time);
1753 if (bSeparate || bSplitHere)
1755 out = ffopen(out_file2, "w");
1760 write_hconf_p(out, title, &useatoms, prec2ndec(frout.prec),
1761 frout.x, frout.bV ? frout.v : NULL, frout.box);
1764 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1765 sprintf(title, "%s t= %9.5f", top_title, frout.time);
1766 /* if reading from pdb, we want to keep the original
1767 model numbering else we write the output frame
1768 number plus one, because model 0 is not allowed in pdb */
1769 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1777 write_pdbfile(out, title, &useatoms, frout.x,
1778 frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1781 frout.title = title;
1782 if (bSeparate || bTDump)
1784 frout.bTitle = TRUE;
1787 frout.bAtoms = TRUE;
1789 frout.atoms = &useatoms;
1790 frout.bStep = FALSE;
1791 frout.bTime = FALSE;
1795 frout.bTitle = (outframe == 0);
1796 frout.bAtoms = FALSE;
1800 write_g96_conf(out, &frout, -1, NULL);
1809 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1811 if (bSeparate || bSplitHere)
1816 /* execute command */
1820 sprintf(c, "%s %d", exec_command, file_nr-1);
1821 /*fprintf(stderr,"Executing '%s'\n",c);*/
1822 #ifdef GMX_NO_SYSTEM
1823 printf("Warning-- No calls to system(3) supported on this platform.");
1824 printf("Warning-- Skipping execution of 'system(\"%s\")'.", c);
1828 gmx_fatal(FARGS, "Error executing command: %s", c);
1836 bHaveNextFrame = read_next_frame(oenv, status, &fr);
1838 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1841 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1843 fprintf(stderr, "\nWARNING no output, "
1844 "last frame read at t=%g\n", fr.time);
1846 fprintf(stderr, "\n");
1853 gmx_rmpbc_done(gpbc);
1860 else if (out != NULL)
1866 for (i = 0; (i < clust->clust->nr); i++)
1868 if (clust_status_id[i] >= 0)
1870 close_trx(clust_status[i]);
1876 do_view(oenv, out_file, NULL);