Replaced NOTSET from typedefs.h by local solutions.
[gromacs.git] / src / gromacs / gmxana / gmx_trjconv.cpp
blobef2bf12c8d13339956a91f1d852e7141c16ebec7
1 /*
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,2015, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
41 #include <cstring>
43 #include <algorithm>
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/g96io.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/groio.h"
51 #include "gromacs/fileio/pdbio.h"
52 #include "gromacs/fileio/tngio_for_tools.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trrio.h"
55 #include "gromacs/fileio/trx.h"
56 #include "gromacs/fileio/trxio.h"
57 #include "gromacs/fileio/xtcio.h"
58 #include "gromacs/fileio/xvgr.h"
59 #include "gromacs/gmxana/gmx_ana.h"
60 #include "gromacs/legacyheaders/copyrite.h"
61 #include "gromacs/legacyheaders/names.h"
62 #include "gromacs/legacyheaders/typedefs.h"
63 #include "gromacs/math/do_fit.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/pbcutil/rmpbc.h"
67 #include "gromacs/topology/index.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/utility/arraysize.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/futil.h"
72 #include "gromacs/utility/smalloc.h"
74 enum {
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;
83 int imin, jmin;
84 real fac, min_dist2;
85 rvec dx, xtest, box_center;
86 int nmol, imol_center;
87 atom_id *molind;
88 gmx_bool *bMol, *bTmp;
89 rvec *m_com, *m_shift;
90 t_pbc pbc;
91 int *cluster;
92 int *added;
93 int ncluster, nadded;
94 real tmp_r2;
96 calc_box_center(ecenter, box, box_center);
98 /* Initiate the pbc structure */
99 std::memset(&pbc, 0, sizeof(pbc));
100 set_pbc(&pbc, ePBC, box);
102 /* Convert atom index to molecular */
103 nmol = top->mols.nr;
104 molind = top->mols.index;
105 snew(bMol, nmol);
106 snew(m_com, nmol);
107 snew(m_shift, nmol);
108 snew(cluster, nmol);
109 snew(added, nmol);
110 snew(bTmp, top->atoms.nr);
112 for (i = 0; (i < nrefat); i++)
114 /* Mark all molecules in the index */
115 ai = index[i];
116 bTmp[ai] = TRUE;
117 /* Binary search assuming the molecules are sorted */
118 j0 = 0;
119 j1 = nmol-1;
120 while (j0 < j1)
122 if (ai < molind[j0+1])
124 j1 = j0;
126 else if (ai >= molind[j1])
128 j0 = j1;
130 else
132 jj = (j0+j1)/2;
133 if (ai < molind[jj+1])
135 j1 = jj;
137 else
139 j0 = jj;
143 bMol[j0] = TRUE;
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));
149 imol_center = -1;
150 ncluster = 0;
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);
163 else if (bMol[i])
165 /* Make molecule whole, move 2nd and higher atom to same periodicity as 1st atom in molecule */
166 if (j > molind[i])
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]);
175 if (bMol[i])
177 /* Normalize center of geometry */
178 fac = 1.0/(molind[i+1]-molind[i]);
179 for (m = 0; (m < DIM); m++)
181 m_com[i][m] *= fac;
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)
189 min_dist2 = tmp_r2;
190 imol_center = i;
192 cluster[ncluster++] = i;
195 sfree(bTmp);
197 if (ncluster <= 0)
199 fprintf(stderr, "No molecules selected in the cluster\n");
200 return;
202 else if (imol_center == -1)
204 fprintf(stderr, "No central molecules could be found\n");
205 return;
208 nadded = 0;
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));
216 imin = -1;
217 jmin = -1;
218 /* Loop over added mols */
219 for (i = 0; i < nadded; i++)
221 ai = added[i];
222 /* Loop over all mols */
223 for (j = 0; j < ncluster; j++)
225 aj = cluster[j];
226 /* check those remaining to be added */
227 if (bMol[aj])
229 pbc_dx(&pbc, m_com[aj], m_com[ai], dx);
230 tmp_r2 = iprod(dx, dx);
231 if (tmp_r2 < min_dist2)
233 min_dist2 = tmp_r2;
234 imin = ai;
235 jmin = aj;
241 /* Add the best molecule */
242 added[nadded++] = jmin;
243 bMol[jmin] = FALSE;
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);
255 fflush(stdout);
258 sfree(added);
259 sfree(cluster);
260 sfree(bMol);
261 sfree(m_com);
262 sfree(m_shift);
264 fprintf(stdout, "\n");
267 static void put_molecule_com_in_box(int unitcell_enum, int ecenter,
268 t_block *mols,
269 int natoms, t_atom atom[],
270 int ePBC, matrix box, rvec x[])
272 atom_id i, j;
273 int d;
274 rvec com, new_com, shift, box_center;
275 real m;
276 double mtot;
277 t_pbc pbc;
279 calc_box_center(ecenter, box, box_center);
280 set_pbc(&pbc, ePBC, box);
281 if (mols->nr <= 0)
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++)
287 /* calc COM */
288 clear_rvec(com);
289 mtot = 0;
290 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
292 m = atom[j].m;
293 for (d = 0; d < DIM; d++)
295 com[d] += m*x[j][d];
297 mtot += m;
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)
306 case euRect:
307 put_atoms_in_box(ePBC, box, 1, &new_com);
308 break;
309 case euTric:
310 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
311 break;
312 case euCompact:
313 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
314 break;
316 rvec_sub(new_com, com, shift);
317 if (norm2(shift) > 0)
319 if (debug)
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;
338 int d, presnr;
339 real m;
340 double mtot;
341 rvec box_center, com, new_com, shift;
342 static const int NOTSET = -12347;
343 calc_box_center(ecenter, box, box_center);
345 presnr = NOTSET;
346 res_start = 0;
347 clear_rvec(com);
348 mtot = 0;
349 for (i = 0; i < natoms+1; i++)
351 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
353 /* calculate final COM */
354 res_end = i;
355 svmul(1.0/mtot, com, com);
357 /* check if COM is outside box */
358 copy_rvec(com, new_com);
359 switch (unitcell_enum)
361 case euRect:
362 put_atoms_in_box(ePBC, box, 1, &new_com);
363 break;
364 case euTric:
365 put_atoms_in_triclinic_unitcell(ecenter, box, 1, &new_com);
366 break;
367 case euCompact:
368 put_atoms_in_compact_unitcell(ePBC, ecenter, box, 1, &new_com);
369 break;
371 rvec_sub(new_com, com, shift);
372 if (norm2(shift))
374 if (debug)
376 fprintf(debug, "\nShifting position of residue %d (atoms %d-%d) "
377 "by %g,%g,%g\n", atom[res_start].resind+1,
378 res_start+1, res_end+1, shift[XX], shift[YY], shift[ZZ]);
380 for (j = res_start; j < res_end; j++)
382 rvec_inc(x[j], shift);
385 clear_rvec(com);
386 mtot = 0;
388 /* remember start of new residue */
389 res_start = i;
391 if (i < natoms)
393 /* calc COM */
394 m = atom[i].m;
395 for (d = 0; d < DIM; d++)
397 com[d] += m*x[i][d];
399 mtot += m;
401 presnr = atom[i].resind;
406 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, atom_id ci[])
408 int i, m, ai;
409 rvec cmin, cmax, box_center, dx;
411 if (nc > 0)
413 copy_rvec(x[ci[0]], cmin);
414 copy_rvec(x[ci[0]], cmax);
415 for (i = 0; i < nc; i++)
417 ai = ci[i];
418 for (m = 0; m < DIM; m++)
420 if (x[ai][m] < cmin[m])
422 cmin[m] = x[ai][m];
424 else if (x[ai][m] > cmax[m])
426 cmax[m] = x[ai][m];
430 calc_box_center(ecenter, box, box_center);
431 for (m = 0; m < DIM; m++)
433 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
436 for (i = 0; i < n; i++)
438 rvec_inc(x[i], dx);
443 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
444 char out_file[])
446 char nbuf[128];
447 int nd = 0, fnr;
449 std::strcpy(out_file, base);
450 fnr = file_nr;
453 fnr /= 10;
454 nd++;
456 while (fnr > 0);
458 if (nd < ndigit)
460 std::strncat(out_file, "00000000000", ndigit-nd);
462 sprintf(nbuf, "%d.", file_nr);
463 std::strcat(out_file, nbuf);
464 std::strcat(out_file, ext);
467 void check_trr(const char *fn)
469 if (fn2ftp(fn) != efTRR)
471 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
475 void do_trunc(const char *fn, real t0)
477 t_fileio *in;
478 FILE *fp;
479 gmx_bool bStop, bOK;
480 gmx_trr_header_t sh;
481 gmx_off_t fpos;
482 char yesno[256];
483 int j;
484 real t = 0;
486 if (t0 == -1)
488 gmx_fatal(FARGS, "You forgot to set the truncation time");
491 /* Check whether this is a .trr file */
492 check_trr(fn);
494 in = gmx_trr_open(fn, "r");
495 fp = gmx_fio_getfp(in);
496 if (fp == NULL)
498 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
499 gmx_trr_close(in);
501 else
503 j = 0;
504 fpos = gmx_fio_ftell(in);
505 bStop = FALSE;
506 while (!bStop && gmx_trr_read_frame_header(in, &sh, &bOK))
508 gmx_trr_read_frame_data(in, &sh, NULL, NULL, NULL, NULL);
509 fpos = gmx_ftell(fp);
510 t = sh.t;
511 if (t >= t0)
513 gmx_fseek(fp, fpos, SEEK_SET);
514 bStop = TRUE;
517 if (bStop)
519 fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
520 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
521 fn, j, t, (long int)fpos);
522 if (1 != scanf("%s", yesno))
524 gmx_fatal(FARGS, "Error reading user input");
526 if (std::strcmp(yesno, "YES") == 0)
528 fprintf(stderr, "Once again, I'm gonna DO this...\n");
529 gmx_trr_close(in);
530 if (0 != gmx_truncate(fn, fpos))
532 gmx_fatal(FARGS, "Error truncating file %s", fn);
535 else
537 fprintf(stderr, "Ok, I'll forget about it\n");
540 else
542 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
543 gmx_trr_close(in);
548 /*! \brief Read a full molecular topology if useful and available.
550 * If the input trajectory file is not in TNG format, and the output
551 * file is in TNG format, then we want to try to read a full topology
552 * (if available), so that we can write molecule information to the
553 * output file. The full topology provides better molecule information
554 * than is available from the normal t_topology data used by GROMACS
555 * tools.
557 * Also, the t_topology is only read under (different) particular
558 * conditions. If both apply, then a .tpr file might be read
559 * twice. Trying to fix this redundancy while trjconv is still an
560 * all-purpose tool does not seem worthwhile.
562 * Because of the way gmx_prepare_tng_writing is implemented, the case
563 * where the input TNG file has no molecule information will never
564 * lead to an output TNG file having molecule information. Since
565 * molecule information will generally be present if the input TNG
566 * file was written by a GROMACS tool, this seems like reasonable
567 * behaviour. */
568 static gmx_mtop_t *read_mtop_for_tng(const char *tps_file,
569 const char *input_file,
570 const char *output_file)
572 gmx_mtop_t *mtop = NULL;
574 if (fn2bTPX(tps_file) &&
575 efTNG != fn2ftp(input_file) &&
576 efTNG == fn2ftp(output_file))
578 int temp_natoms = -1;
579 snew(mtop, 1);
580 read_tpx(tps_file, NULL, NULL, &temp_natoms,
581 NULL, NULL, mtop);
584 return mtop;
587 int gmx_trjconv(int argc, char *argv[])
589 const char *desc[] = {
590 "[THISMODULE] can convert trajectory files in many ways:",
592 "* from one format to another",
593 "* select a subset of atoms",
594 "* change the periodicity representation",
595 "* keep multimeric molecules together",
596 "* center atoms in the box",
597 "* fit atoms to reference structure",
598 "* reduce the number of frames",
599 "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
600 "* cut the trajectory in small subtrajectories according",
601 " to information in an index file. This allows subsequent analysis of",
602 " the subtrajectories that could, for example, be the result of a",
603 " cluster analysis. Use option [TT]-sub[tt].",
604 " This assumes that the entries in the index file are frame numbers and",
605 " dumps each group in the index file to a separate trajectory file.",
606 "* select frames within a certain range of a quantity given",
607 " in an [REF].xvg[ref] file.",
609 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
610 "[PAR]",
612 "The following formats are supported for input and output:",
613 "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
614 "and [REF].pdb[ref].",
615 "The file formats are detected from the file extension.",
616 "The precision of [REF].xtc[ref] and [REF].gro[ref] output is taken from the",
617 "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
618 "and from the [TT]-ndec[tt] option for other input formats. The precision",
619 "is always taken from [TT]-ndec[tt], when this option is set.",
620 "All other formats have fixed precision. [REF].trr[ref]",
621 "output can be single or double precision, depending on the precision",
622 "of the [THISMODULE] binary.",
623 "Note that velocities are only supported in",
624 "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
626 "Option [TT]-sep[tt] can be used to write every frame to a separate",
627 "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to one file.",
628 "[REF].pdb[ref] files with all frames concatenated can be viewed with",
629 "[TT]rasmol -nmrpdb[tt].[PAR]",
631 "It is possible to select part of your trajectory and write it out",
632 "to a new trajectory file in order to save disk space, e.g. for leaving",
633 "out the water from a trajectory of a protein in water.",
634 "[BB]ALWAYS[bb] put the original trajectory on tape!",
635 "We recommend to use the portable [REF].xtc[ref] format for your analysis",
636 "to save disk space and to have portable files.[PAR]",
638 "There are two options for fitting the trajectory to a reference",
639 "either for essential dynamics analysis, etc.",
640 "The first option is just plain fitting to a reference structure",
641 "in the structure file. The second option is a progressive fit",
642 "in which the first timeframe is fitted to the reference structure ",
643 "in the structure file to obtain and each subsequent timeframe is ",
644 "fitted to the previously fitted structure. This way a continuous",
645 "trajectory is generated, which might not be the case when using the",
646 "regular fit method, e.g. when your protein undergoes large",
647 "conformational transitions.[PAR]",
649 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
650 "treatment:",
652 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
653 " and requires a run input file to be supplied with [TT]-s[tt].",
654 " * [TT]res[tt] puts the center of mass of residues in the box.",
655 " * [TT]atom[tt] puts all the atoms in the box.",
656 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
657 " them back. This has the effect that all molecules",
658 " will remain whole (provided they were whole in the initial",
659 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
660 " molecules may diffuse out of the box. The starting configuration",
661 " for this procedure is taken from the structure file, if one is",
662 " supplied, otherwise it is the first frame.",
663 " * [TT]cluster[tt] clusters all the atoms in the selected index",
664 " such that they are all closest to the center of mass of the cluster,",
665 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
666 " results if you in fact have a cluster. Luckily that can be checked",
667 " afterwards using a trajectory viewer. Note also that if your molecules",
668 " are broken this will not work either.",
670 " The separate option [TT]-clustercenter[tt] can be used to specify an",
671 " approximate center for the cluster. This is useful e.g. if you have",
672 " two big vesicles, and you want to maintain their relative positions.",
673 " * [TT]whole[tt] only makes broken molecules whole.",
676 "Option [TT]-ur[tt] sets the unit cell representation for options",
677 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
678 "All three options give different results for triclinic boxes and",
679 "identical results for rectangular boxes.",
680 "[TT]rect[tt] is the ordinary brick shape.",
681 "[TT]tric[tt] is the triclinic unit cell.",
682 "[TT]compact[tt] puts all atoms at the closest distance from the center",
683 "of the box. This can be useful for visualizing e.g. truncated octahedra",
684 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
685 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
686 "is set differently.[PAR]",
688 "Option [TT]-center[tt] centers the system in the box. The user can",
689 "select the group which is used to determine the geometrical center.",
690 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
691 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
692 "[TT]tric[tt]: half of the sum of the box vectors,",
693 "[TT]rect[tt]: half of the box diagonal,",
694 "[TT]zero[tt]: zero.",
695 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
696 "want all molecules in the box after the centering.[PAR]",
698 "Option [TT]-box[tt] sets the size of the new box. This option only works",
699 "for leading dimensions and is thus generally only useful for rectangular boxes.",
700 "If you want to modify only some of the dimensions, e.g. when reading from",
701 "a trajectory, you can use -1 for those dimensions that should stay the same",
703 "It is not always possible to use combinations of [TT]-pbc[tt],",
704 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
705 "you want in one call to [THISMODULE]. Consider using multiple",
706 "calls, and check out the GROMACS website for suggestions.[PAR]",
708 "With [TT]-dt[tt], it is possible to reduce the number of ",
709 "frames in the output. This option relies on the accuracy of the times",
710 "in your input trajectory, so if these are inaccurate use the",
711 "[TT]-timestep[tt] option to modify the time (this can be done",
712 "simultaneously). For making smooth movies, the program [gmx-filter]",
713 "can reduce the number of frames while using low-pass frequency",
714 "filtering, this reduces aliasing of high frequency motions.[PAR]",
716 "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
717 "without copying the file. This is useful when a run has crashed",
718 "during disk I/O (i.e. full disk), or when two contiguous",
719 "trajectories must be concatenated without having double frames.[PAR]",
721 "Option [TT]-dump[tt] can be used to extract a frame at or near",
722 "one specific time from your trajectory, but only works reliably",
723 "if the time interval between frames is uniform.[PAR]",
725 "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
726 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
727 "frames with a value below and above the value of the respective options",
728 "will not be written."
731 int pbc_enum;
732 enum
734 epSel,
735 epNone,
736 epComMol,
737 epComRes,
738 epComAtom,
739 epNojump,
740 epCluster,
741 epWhole,
742 epNR
744 const char *pbc_opt[epNR + 1] =
746 NULL, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
747 NULL
750 int unitcell_enum;
751 const char *unitcell_opt[euNR+1] =
752 { NULL, "rect", "tric", "compact", NULL };
754 enum
756 ecSel, ecTric, ecRect, ecZero, ecNR
758 const char *center_opt[ecNR+1] =
759 { NULL, "tric", "rect", "zero", NULL };
760 int ecenter;
762 int fit_enum;
763 enum
765 efSel, efNone, efFit, efFitXY, efReset, efResetXY, efPFit, efNR
767 const char *fit[efNR + 1] =
769 NULL, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
770 "progressive", NULL
773 static gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
774 static gmx_bool bCenter = FALSE;
775 static int skip_nr = 1, ndec = 3, nzero = 0;
776 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
777 static rvec newbox = {0, 0, 0}, shift = {0, 0, 0}, trans = {0, 0, 0};
778 static char *exec_command = NULL;
779 static real dropunder = 0, dropover = 0;
780 static gmx_bool bRound = FALSE;
782 t_pargs
783 pa[] =
785 { "-skip", FALSE, etINT,
786 { &skip_nr }, "Only write every nr-th frame" },
787 { "-dt", FALSE, etTIME,
788 { &delta_t },
789 "Only write frame when t MOD dt = first time (%t)" },
790 { "-round", FALSE, etBOOL,
791 { &bRound }, "Round measurements to nearest picosecond"},
792 { "-dump", FALSE, etTIME,
793 { &tdump }, "Dump frame nearest specified time (%t)" },
794 { "-t0", FALSE, etTIME,
795 { &tzero },
796 "Starting time (%t) (default: don't change)" },
797 { "-timestep", FALSE, etTIME,
798 { &timestep },
799 "Change time step between input frames (%t)" },
800 { "-pbc", FALSE, etENUM,
801 { pbc_opt },
802 "PBC treatment (see help text for full description)" },
803 { "-ur", FALSE, etENUM,
804 { unitcell_opt }, "Unit-cell representation" },
805 { "-center", FALSE, etBOOL,
806 { &bCenter }, "Center atoms in box" },
807 { "-boxcenter", FALSE, etENUM,
808 { center_opt }, "Center for -pbc and -center" },
809 { "-box", FALSE, etRVEC,
810 { newbox },
811 "Size for new cubic box (default: read from input)" },
812 { "-trans", FALSE, etRVEC,
813 { trans },
814 "All coordinates will be translated by trans. This "
815 "can advantageously be combined with -pbc mol -ur "
816 "compact." },
817 { "-shift", FALSE, etRVEC,
818 { shift },
819 "All coordinates will be shifted by framenr*shift" },
820 { "-fit", FALSE, etENUM,
821 { fit },
822 "Fit molecule to ref structure in the structure file" },
823 { "-ndec", FALSE, etINT,
824 { &ndec },
825 "Precision for .xtc and .gro writing in number of "
826 "decimal places" },
827 { "-vel", FALSE, etBOOL,
828 { &bVels }, "Read and write velocities if possible" },
829 { "-force", FALSE, etBOOL,
830 { &bForce }, "Read and write forces if possible" },
831 { "-trunc", FALSE, etTIME,
832 { &ttrunc },
833 "Truncate input trajectory file after this time (%t)" },
834 { "-exec", FALSE, etSTR,
835 { &exec_command },
836 "Execute command for every output frame with the "
837 "frame number as argument" },
838 { "-split", FALSE, etTIME,
839 { &split_t },
840 "Start writing new file when t MOD split = first "
841 "time (%t)" },
842 { "-sep", FALSE, etBOOL,
843 { &bSeparate },
844 "Write each frame to a separate .gro, .g96 or .pdb "
845 "file" },
846 { "-nzero", FALSE, etINT,
847 { &nzero },
848 "If the -sep flag is set, use these many digits "
849 "for the file numbers and prepend zeros as needed" },
850 { "-dropunder", FALSE, etREAL,
851 { &dropunder }, "Drop all frames below this value" },
852 { "-dropover", FALSE, etREAL,
853 { &dropover }, "Drop all frames above this value" },
854 { "-conect", FALSE, etBOOL,
855 { &bCONECT },
856 "Add conect records when writing [REF].pdb[ref] files. Useful "
857 "for visualization of non-standard molecules, e.g. "
858 "coarse grained ones" }
860 #define NPA asize(pa)
862 FILE *out = NULL;
863 t_trxstatus *trxout = NULL;
864 t_trxstatus *trxin;
865 int ftp, ftpin = 0, file_nr;
866 t_trxframe fr, frout;
867 int flags;
868 rvec *xmem = NULL, *vmem = NULL, *fmem = NULL;
869 rvec *xp = NULL, x_shift, hbox;
870 real *w_rls = NULL;
871 int m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
872 #define SKIP 10
873 t_topology top;
874 gmx_mtop_t *mtop = NULL;
875 gmx_conect gc = NULL;
876 int ePBC = -1;
877 t_atoms *atoms = NULL, useatoms;
878 matrix top_box;
879 atom_id *index, *cindex;
880 char *grpnm;
881 int *frindex, nrfri;
882 char *frname;
883 int ifit, my_clust = -1;
884 atom_id *ind_fit;
885 char *gn_fit;
886 t_cluster_ndx *clust = NULL;
887 t_trxstatus **clust_status = NULL;
888 int *clust_status_id = NULL;
889 int ntrxopen = 0;
890 int *nfwritten = NULL;
891 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
892 double **dropval;
893 real tshift = 0, t0 = -1, dt = 0.001, prec;
894 gmx_bool bFit, bPFit, bReset;
895 int nfitdim;
896 gmx_rmpbc_t gpbc = NULL;
897 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
898 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
899 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetPrec, bNeedPrec;
900 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
901 gmx_bool bSubTraj = FALSE, bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
902 gmx_bool bWriteFrame, bSplitHere;
903 const char *top_file, *in_file, *out_file = NULL;
904 char out_file2[256], *charpt;
905 char *outf_base = NULL;
906 const char *outf_ext = NULL;
907 char top_title[256], title[256], filemode[5];
908 gmx_output_env_t *oenv;
910 t_filenm fnm[] = {
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,
924 NFILE, fnm, NPA, pa, asize(desc), desc,
925 0, NULL, &oenv))
927 return 0;
930 top_file = ftp2fn(efTPS, NFILE, fnm);
931 init_top(&top);
933 /* Check command line */
934 in_file = opt2fn("-f", NFILE, fnm);
936 if (ttrunc != -1)
938 do_trunc(in_file, ttrunc);
940 else
942 /* mark active cmdline options */
943 bSetBox = opt2parg_bSet("-box", NPA, pa);
944 bSetTime = opt2parg_bSet("-t0", NPA, pa);
945 bSetPrec = opt2parg_bSet("-ndec", NPA, pa);
946 bSetUR = opt2parg_bSet("-ur", NPA, pa);
947 bExec = opt2parg_bSet("-exec", NPA, pa);
948 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
949 bTDump = opt2parg_bSet("-dump", NPA, pa);
950 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
951 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
952 bTrans = opt2parg_bSet("-trans", NPA, pa);
953 bSplit = (split_t != 0);
955 /* parse enum options */
956 fit_enum = nenum(fit);
957 bFit = (fit_enum == efFit || fit_enum == efFitXY);
958 bReset = (fit_enum == efReset || fit_enum == efResetXY);
959 bPFit = fit_enum == efPFit;
960 pbc_enum = nenum(pbc_opt);
961 bPBCWhole = pbc_enum == epWhole;
962 bPBCcomRes = pbc_enum == epComRes;
963 bPBCcomMol = pbc_enum == epComMol;
964 bPBCcomAtom = pbc_enum == epComAtom;
965 bNoJump = pbc_enum == epNojump;
966 bCluster = pbc_enum == epCluster;
967 bPBC = pbc_enum != epNone;
968 unitcell_enum = nenum(unitcell_opt);
969 ecenter = nenum(center_opt) - ecTric;
971 /* set and check option dependencies */
972 if (bPFit)
974 bFit = TRUE; /* for pfit, fit *must* be set */
976 if (bFit)
978 bReset = TRUE; /* for fit, reset *must* be set */
980 nfitdim = 0;
981 if (bFit || bReset)
983 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
985 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
987 if (bSetUR)
989 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
991 fprintf(stderr,
992 "WARNING: Option for unitcell representation (-ur %s)\n"
993 " only has effect in combination with -pbc %s, %s or %s.\n"
994 " Ingoring unitcell representation.\n\n",
995 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
998 if (bFit && bPBC)
1000 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
1001 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1002 "for the rotational fit.\n"
1003 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1004 "results!");
1007 /* ndec is in nr of decimal places, prec is a multiplication factor: */
1008 prec = 1;
1009 for (i = 0; i < ndec; i++)
1011 prec *= 10;
1014 bIndex = ftp2bSet(efNDX, NFILE, fnm);
1017 /* Determine output type */
1018 out_file = opt2fn("-o", NFILE, fnm);
1019 ftp = fn2ftp(out_file);
1020 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1021 bNeedPrec = (ftp == efXTC || ftp == efGRO);
1022 if (bVels)
1024 /* check if velocities are possible in input and output files */
1025 ftpin = fn2ftp(in_file);
1026 bVels = (ftp == efTRR || ftp == efGRO ||
1027 ftp == efG96 || ftp == efTNG)
1028 && (ftpin == efTRR || ftpin == efGRO ||
1029 ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
1031 if (bSeparate || bSplit)
1033 outf_ext = std::strrchr(out_file, '.');
1034 if (outf_ext == NULL)
1036 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1038 outf_base = gmx_strdup(out_file);
1039 outf_base[outf_ext - out_file] = '\0';
1042 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1043 if (bSubTraj)
1045 if ((ftp != efXTC) && (ftp != efTRR))
1047 /* It seems likely that other trajectory file types
1048 * could work here. */
1049 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1050 "xtc and trr");
1052 clust = cluster_index(NULL, opt2fn("-sub", NFILE, fnm));
1054 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1055 * number to check for. In my linux box it is only 16.
1057 if (0 && (clust->clust->nr > FOPEN_MAX-4))
1059 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1060 " trajectories.\ntry splitting the index file in %d parts.\n"
1061 "FOPEN_MAX = %d",
1062 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1064 gmx_warning("The -sub option could require as many open output files as there are\n"
1065 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1066 "try reducing the number of index groups in the file, and perhaps\n"
1067 "using trjconv -sub several times on different chunks of your index file.\n",
1068 clust->clust->nr);
1070 snew(clust_status, clust->clust->nr);
1071 snew(clust_status_id, clust->clust->nr);
1072 snew(nfwritten, clust->clust->nr);
1073 for (i = 0; (i < clust->clust->nr); i++)
1075 clust_status[i] = NULL;
1076 clust_status_id[i] = -1;
1078 bSeparate = bSplit = FALSE;
1080 /* skipping */
1081 if (skip_nr <= 0)
1085 mtop = read_mtop_for_tng(top_file, in_file, out_file);
1087 /* Determine whether to read a topology */
1088 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1089 bRmPBC || bReset || bPBCcomMol || bCluster ||
1090 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1092 /* Determine if when can read index groups */
1093 bIndex = (bIndex || bTPS);
1095 if (bTPS)
1097 read_tps_conf(top_file, &top, &ePBC, &xp, NULL, top_box,
1098 bReset || bPBCcomRes);
1099 std::strncpy(top_title, *top.name, 255);
1100 top_title[255] = '\0';
1101 atoms = &top.atoms;
1103 if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1105 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1108 /* top_title is only used for gro and pdb,
1109 * the header in such a file is top_title t= ...
1110 * to prevent a double t=, remove it from top_title
1112 if ((charpt = std::strstr(top_title, " t= ")))
1114 charpt[0] = '\0';
1117 if (bCONECT)
1119 gc = gmx_conect_generate(&top);
1121 if (bRmPBC)
1123 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1127 /* get frame number index */
1128 frindex = NULL;
1129 if (opt2bSet("-fr", NFILE, fnm))
1131 printf("Select groups of frame number indices:\n");
1132 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (atom_id **)&frindex, &frname);
1133 if (debug)
1135 for (i = 0; i < nrfri; i++)
1137 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1142 /* get index groups etc. */
1143 if (bReset)
1145 printf("Select group for %s fit\n",
1146 bFit ? "least squares" : "translational");
1147 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1148 1, &ifit, &ind_fit, &gn_fit);
1150 if (bFit)
1152 if (ifit < 2)
1154 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1156 else if (ifit == 3)
1158 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1162 else if (bCluster)
1164 printf("Select group for clustering\n");
1165 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1166 1, &ifit, &ind_fit, &gn_fit);
1169 if (bIndex)
1171 if (bCenter)
1173 printf("Select group for centering\n");
1174 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1175 1, &ncent, &cindex, &grpnm);
1177 printf("Select group for output\n");
1178 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1179 1, &nout, &index, &grpnm);
1181 else
1183 /* no index file, so read natoms from TRX */
1184 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1186 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1188 natoms = fr.natoms;
1189 close_trj(trxin);
1190 sfree(fr.x);
1191 snew(index, natoms);
1192 for (i = 0; i < natoms; i++)
1194 index[i] = i;
1196 nout = natoms;
1197 if (bCenter)
1199 ncent = nout;
1200 cindex = index;
1204 if (bReset)
1206 snew(w_rls, atoms->nr);
1207 for (i = 0; (i < ifit); i++)
1209 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1212 /* Restore reference structure and set to origin,
1213 store original location (to put structure back) */
1214 if (bRmPBC)
1216 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1218 copy_rvec(xp[index[0]], x_shift);
1219 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, NULL, xp, w_rls);
1220 rvec_dec(x_shift, xp[index[0]]);
1222 else
1224 clear_rvec(x_shift);
1227 if (bDropUnder || bDropOver)
1229 /* Read the .xvg file with the drop values */
1230 fprintf(stderr, "\nReading drop file ...");
1231 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1232 fprintf(stderr, " %d time points\n", ndrop);
1233 if (ndrop == 0 || ncol < 2)
1235 gmx_fatal(FARGS, "Found no data points in %s",
1236 opt2fn("-drop", NFILE, fnm));
1238 drop0 = 0;
1239 drop1 = 0;
1242 /* Make atoms struct for output in GRO or PDB files */
1243 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1245 /* get memory for stuff to go in .pdb file, and initialize
1246 * the pdbinfo structure part if the input has it.
1248 init_t_atoms(&useatoms, atoms->nr, (atoms->pdbinfo != NULL));
1249 sfree(useatoms.resinfo);
1250 useatoms.resinfo = atoms->resinfo;
1251 for (i = 0; (i < nout); i++)
1253 useatoms.atomname[i] = atoms->atomname[index[i]];
1254 useatoms.atom[i] = atoms->atom[index[i]];
1255 if (atoms->pdbinfo != NULL)
1257 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
1259 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind+1);
1261 useatoms.nr = nout;
1263 /* select what to read */
1264 if (ftp == efTRR)
1266 flags = TRX_READ_X;
1268 else
1270 flags = TRX_NEED_X;
1272 if (bVels)
1274 flags = flags | TRX_READ_V;
1276 if (bForce)
1278 flags = flags | TRX_READ_F;
1281 /* open trx file for reading */
1282 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1283 if (fr.bPrec)
1285 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1287 if (bNeedPrec)
1289 if (bSetPrec || !fr.bPrec)
1291 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1293 else
1295 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1299 if (bHaveFirstFrame)
1301 set_trxframe_ePBC(&fr, ePBC);
1303 natoms = fr.natoms;
1305 if (bSetTime)
1307 tshift = tzero-fr.time;
1309 else
1311 tzero = fr.time;
1314 bCopy = FALSE;
1315 if (bIndex)
1317 /* check if index is meaningful */
1318 for (i = 0; i < nout; i++)
1320 if (index[i] >= natoms)
1322 gmx_fatal(FARGS,
1323 "Index[%d] %d is larger than the number of atoms in the\n"
1324 "trajectory file (%d). There is a mismatch in the contents\n"
1325 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1327 bCopy = bCopy || (i != index[i]);
1331 /* open output for writing */
1332 std::strcpy(filemode, "w");
1333 switch (ftp)
1335 case efTNG:
1336 trjtools_gmx_prepare_tng_writing(out_file,
1337 filemode[0],
1338 trxin,
1339 &trxout,
1340 NULL,
1341 nout,
1342 mtop,
1343 index,
1344 grpnm);
1345 break;
1346 case efXTC:
1347 case efTRR:
1348 out = NULL;
1349 if (!bSplit && !bSubTraj)
1351 trxout = open_trx(out_file, filemode);
1353 break;
1354 case efGRO:
1355 case efG96:
1356 case efPDB:
1357 if (( !bSeparate && !bSplit ) && !bSubTraj)
1359 out = gmx_ffopen(out_file, filemode);
1361 break;
1362 default:
1363 gmx_incons("Illegal output file format");
1366 if (bCopy)
1368 snew(xmem, nout);
1369 if (bVels)
1371 snew(vmem, nout);
1373 if (bForce)
1375 snew(fmem, nout);
1379 /* Start the big loop over frames */
1380 file_nr = 0;
1381 frame = 0;
1382 outframe = 0;
1383 model_nr = 0;
1384 bDTset = FALSE;
1386 /* Main loop over frames */
1389 if (!fr.bStep)
1391 /* set the step */
1392 fr.step = newstep;
1393 newstep++;
1395 if (bSubTraj)
1397 /*if (frame >= clust->clust->nra)
1398 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1399 if (frame > clust->maxframe)
1401 my_clust = -1;
1403 else
1405 my_clust = clust->inv_clust[frame];
1407 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1408 (my_clust == NO_ATID))
1410 my_clust = -1;
1414 if (bSetBox)
1416 /* generate new box */
1417 if (fr.bBox == FALSE)
1419 clear_mat(fr.box);
1421 for (m = 0; m < DIM; m++)
1423 if (newbox[m] >= 0)
1425 fr.box[m][m] = newbox[m];
1427 else
1429 if (fr.bBox == FALSE)
1431 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1437 if (bTrans)
1439 for (i = 0; i < natoms; i++)
1441 rvec_inc(fr.x[i], trans);
1445 if (bTDump)
1447 /* determine timestep */
1448 if (t0 == -1)
1450 t0 = fr.time;
1452 else
1454 if (!bDTset)
1456 dt = fr.time-t0;
1457 bDTset = TRUE;
1460 /* This is not very elegant, as one can not dump a frame after
1461 * a timestep with is more than twice as small as the first one. */
1462 bDumpFrame = (fr.time > tdump-0.5*dt) && (fr.time <= tdump+0.5*dt);
1464 else
1466 bDumpFrame = FALSE;
1469 /* determine if an atom jumped across the box and reset it if so */
1470 if (bNoJump && (bTPS || frame != 0))
1472 for (d = 0; d < DIM; d++)
1474 hbox[d] = 0.5*fr.box[d][d];
1476 for (i = 0; i < natoms; i++)
1478 if (bReset)
1480 rvec_dec(fr.x[i], x_shift);
1482 for (m = DIM-1; m >= 0; m--)
1484 if (hbox[m] > 0)
1486 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1488 for (d = 0; d <= m; d++)
1490 fr.x[i][d] += fr.box[m][d];
1493 while (fr.x[i][m]-xp[i][m] > hbox[m])
1495 for (d = 0; d <= m; d++)
1497 fr.x[i][d] -= fr.box[m][d];
1504 else if (bCluster)
1506 calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1509 if (bPFit)
1511 /* Now modify the coords according to the flags,
1512 for normal fit, this is only done for output frames */
1513 if (bRmPBC)
1515 gmx_rmpbc_trxfr(gpbc, &fr);
1518 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1519 do_fit(natoms, w_rls, xp, fr.x);
1522 /* store this set of coordinates for future use */
1523 if (bPFit || bNoJump)
1525 if (xp == NULL)
1527 snew(xp, natoms);
1529 for (i = 0; (i < natoms); i++)
1531 copy_rvec(fr.x[i], xp[i]);
1532 rvec_inc(fr.x[i], x_shift);
1536 if (frindex)
1538 /* see if we have a frame from the frame index group */
1539 for (i = 0; i < nrfri && !bDumpFrame; i++)
1541 bDumpFrame = frame == frindex[i];
1544 if (debug && bDumpFrame)
1546 fprintf(debug, "dumping %d\n", frame);
1549 bWriteFrame =
1550 ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1552 if (bWriteFrame && (bDropUnder || bDropOver))
1554 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1556 drop0 = drop1;
1557 drop1++;
1559 if (std::abs(dropval[0][drop0] - fr.time)
1560 < std::abs(dropval[0][drop1] - fr.time))
1562 dropuse = drop0;
1564 else
1566 dropuse = drop1;
1568 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1569 (bDropOver && dropval[1][dropuse] > dropover))
1571 bWriteFrame = FALSE;
1575 if (bWriteFrame)
1577 /* We should avoid modifying the input frame,
1578 * but since here we don't have the output frame yet,
1579 * we introduce a temporary output frame time variable.
1581 real frout_time;
1583 frout_time = fr.time;
1585 /* calc new time */
1586 if (bTimeStep)
1588 frout_time = tzero + frame*timestep;
1590 else
1591 if (bSetTime)
1593 frout_time += tshift;
1596 if (bTDump)
1598 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1599 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv));
1602 /* check for writing at each delta_t */
1603 bDoIt = (delta_t == 0);
1604 if (!bDoIt)
1606 if (!bRound)
1608 bDoIt = bRmod(frout_time, tzero, delta_t);
1610 else
1612 /* round() is not C89 compatible, so we do this: */
1613 bDoIt = bRmod(std::floor(frout_time+0.5), std::floor(tzero+0.5),
1614 std::floor(delta_t+0.5));
1618 if (bDoIt || bTDump)
1620 /* print sometimes */
1621 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1623 fprintf(stderr, " -> frame %6d time %8.3f \r",
1624 outframe, output_env_conv_time(oenv, frout_time));
1627 if (!bPFit)
1629 /* Now modify the coords according to the flags,
1630 for PFit we did this already! */
1632 if (bRmPBC)
1634 gmx_rmpbc_trxfr(gpbc, &fr);
1637 if (bReset)
1639 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, NULL, fr.x, w_rls);
1640 if (bFit)
1642 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1644 if (!bCenter)
1646 for (i = 0; i < natoms; i++)
1648 rvec_inc(fr.x[i], x_shift);
1653 if (bCenter)
1655 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1659 if (bPBCcomAtom)
1661 switch (unitcell_enum)
1663 case euRect:
1664 put_atoms_in_box(ePBC, fr.box, natoms, fr.x);
1665 break;
1666 case euTric:
1667 put_atoms_in_triclinic_unitcell(ecenter, fr.box, natoms, fr.x);
1668 break;
1669 case euCompact:
1670 put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1671 natoms, fr.x);
1672 break;
1675 if (bPBCcomRes)
1677 put_residue_com_in_box(unitcell_enum, ecenter,
1678 natoms, atoms->atom, ePBC, fr.box, fr.x);
1680 if (bPBCcomMol)
1682 put_molecule_com_in_box(unitcell_enum, ecenter,
1683 &top.mols,
1684 natoms, atoms->atom, ePBC, fr.box, fr.x);
1686 /* Copy the input trxframe struct to the output trxframe struct */
1687 frout = fr;
1688 frout.time = frout_time;
1689 frout.bV = (frout.bV && bVels);
1690 frout.bF = (frout.bF && bForce);
1691 frout.natoms = nout;
1692 if (bNeedPrec && (bSetPrec || !fr.bPrec))
1694 frout.bPrec = TRUE;
1695 frout.prec = prec;
1697 if (bCopy)
1699 frout.x = xmem;
1700 if (frout.bV)
1702 frout.v = vmem;
1704 if (frout.bF)
1706 frout.f = fmem;
1708 for (i = 0; i < nout; i++)
1710 copy_rvec(fr.x[index[i]], frout.x[i]);
1711 if (frout.bV)
1713 copy_rvec(fr.v[index[i]], frout.v[i]);
1715 if (frout.bF)
1717 copy_rvec(fr.f[index[i]], frout.f[i]);
1722 if (opt2parg_bSet("-shift", NPA, pa))
1724 for (i = 0; i < nout; i++)
1726 for (d = 0; d < DIM; d++)
1728 frout.x[i][d] += outframe*shift[d];
1733 if (!bRound)
1735 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1737 else
1739 /* round() is not C89 compatible, so we do this: */
1740 bSplitHere = bSplit && bRmod(std::floor(frout.time+0.5),
1741 std::floor(tzero+0.5),
1742 std::floor(split_t+0.5));
1744 if (bSeparate || bSplitHere)
1746 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1749 switch (ftp)
1751 case efTNG:
1752 write_tng_frame(trxout, &frout);
1753 // TODO when trjconv behaves better: work how to read and write lambda
1754 break;
1755 case efTRR:
1756 case efXTC:
1757 if (bSplitHere)
1759 if (trxout)
1761 close_trx(trxout);
1763 trxout = open_trx(out_file2, filemode);
1765 if (bSubTraj)
1767 if (my_clust != -1)
1769 char buf[STRLEN];
1770 if (clust_status_id[my_clust] == -1)
1772 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1773 clust_status[my_clust] = open_trx(buf, "w");
1774 clust_status_id[my_clust] = 1;
1775 ntrxopen++;
1777 else if (clust_status_id[my_clust] == -2)
1779 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1780 clust->grpname[my_clust], ntrxopen, frame,
1781 my_clust);
1783 write_trxframe(clust_status[my_clust], &frout, gc);
1784 nfwritten[my_clust]++;
1785 if (nfwritten[my_clust] ==
1786 (clust->clust->index[my_clust+1]-
1787 clust->clust->index[my_clust]))
1789 close_trx(clust_status[my_clust]);
1790 clust_status[my_clust] = NULL;
1791 clust_status_id[my_clust] = -2;
1792 ntrxopen--;
1793 if (ntrxopen < 0)
1795 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1800 else
1802 write_trxframe(trxout, &frout, gc);
1804 break;
1805 case efGRO:
1806 case efG96:
1807 case efPDB:
1808 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1809 top_title, frout.time);
1810 if (bSeparate || bSplitHere)
1812 out = gmx_ffopen(out_file2, "w");
1814 switch (ftp)
1816 case efGRO:
1817 write_hconf_p(out, title, &useatoms, prec2ndec(frout.prec),
1818 frout.x, frout.bV ? frout.v : NULL, frout.box);
1819 break;
1820 case efPDB:
1821 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1822 sprintf(title, "%s t= %9.5f", top_title, frout.time);
1823 /* if reading from pdb, we want to keep the original
1824 model numbering else we write the output frame
1825 number plus one, because model 0 is not allowed in pdb */
1826 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1828 model_nr = fr.step;
1830 else
1832 model_nr++;
1834 write_pdbfile(out, title, &useatoms, frout.x,
1835 frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1836 break;
1837 case efG96:
1838 frout.title = title;
1839 if (bSeparate || bTDump)
1841 frout.bTitle = TRUE;
1842 if (bTPS)
1844 frout.bAtoms = TRUE;
1846 frout.atoms = &useatoms;
1847 frout.bStep = FALSE;
1848 frout.bTime = FALSE;
1850 else
1852 frout.bTitle = (outframe == 0);
1853 frout.bAtoms = FALSE;
1854 frout.bStep = TRUE;
1855 frout.bTime = TRUE;
1857 write_g96_conf(out, &frout, -1, NULL);
1859 if (bSeparate || bSplitHere)
1861 gmx_ffclose(out);
1862 out = NULL;
1864 break;
1865 default:
1866 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1868 if (bSeparate || bSplitHere)
1870 file_nr++;
1873 /* execute command */
1874 if (bExec)
1876 char c[255];
1877 sprintf(c, "%s %d", exec_command, file_nr-1);
1878 /*fprintf(stderr,"Executing '%s'\n",c);*/
1879 if (0 != system(c))
1881 gmx_fatal(FARGS, "Error executing command: %s", c);
1884 outframe++;
1887 frame++;
1888 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1890 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1893 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1895 fprintf(stderr, "\nWARNING no output, "
1896 "last frame read at t=%g\n", fr.time);
1898 fprintf(stderr, "\n");
1900 close_trj(trxin);
1901 sfree(outf_base);
1903 if (bRmPBC)
1905 gmx_rmpbc_done(gpbc);
1908 if (trxout)
1910 close_trx(trxout);
1912 else if (out != NULL)
1914 gmx_ffclose(out);
1916 if (bSubTraj)
1918 for (i = 0; (i < clust->clust->nr); i++)
1920 if (clust_status_id[i] >= 0)
1922 close_trx(clust_status[i]);
1928 sfree(mtop);
1930 do_view(oenv, out_file, NULL);
1932 return 0;