Fix trjconv for tdump<frame timestep
[gromacs.git] / src / gromacs / gmxana / gmx_trjconv.cpp
blobe49700f45439b81d39ef3d0884015fc85c9c0d09
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,2016,2017,2018, 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.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trrio.h"
55 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/fileio/xtcio.h"
57 #include "gromacs/fileio/xvgr.h"
58 #include "gromacs/gmxana/gmx_ana.h"
59 #include "gromacs/math/do_fit.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/pbcutil/rmpbc.h"
65 #include "gromacs/topology/index.h"
66 #include "gromacs/topology/topology.h"
67 #include "gromacs/trajectory/trajectoryframe.h"
68 #include "gromacs/utility/arrayref.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[], int 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 int *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*gmx::square(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*gmx::square(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 int i, j;
273 int d;
274 rvec 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 gmx::RVec newCom;
304 copy_rvec(com, newCom);
305 auto newComArrayRef = gmx::arrayRefFromArray(&newCom, 1);
306 switch (unitcell_enum)
308 case euRect:
309 put_atoms_in_box(ePBC, box, newComArrayRef);
310 break;
311 case euTric:
312 put_atoms_in_triclinic_unitcell(ecenter, box, newComArrayRef);
313 break;
314 case euCompact:
315 put_atoms_in_compact_unitcell(ePBC, ecenter, box, newComArrayRef);
316 break;
318 rvec_sub(newCom, com, shift);
319 if (norm2(shift) > 0)
321 if (debug)
323 fprintf(debug, "\nShifting position of molecule %d "
324 "by %8.3f %8.3f %8.3f\n", i+1,
325 shift[XX], shift[YY], shift[ZZ]);
327 for (j = mols->index[i]; (j < mols->index[i+1] && j < natoms); j++)
329 rvec_inc(x[j], shift);
335 static void put_residue_com_in_box(int unitcell_enum, int ecenter,
336 int natoms, t_atom atom[],
337 int ePBC, matrix box, rvec x[])
339 int i, j, res_start, res_end;
340 int d, presnr;
341 real m;
342 double mtot;
343 rvec box_center, com, shift;
344 static const int NOTSET = -12347;
345 calc_box_center(ecenter, box, box_center);
347 presnr = NOTSET;
348 res_start = 0;
349 clear_rvec(com);
350 mtot = 0;
351 for (i = 0; i < natoms+1; i++)
353 if (i == natoms || (presnr != atom[i].resind && presnr != NOTSET))
355 /* calculate final COM */
356 res_end = i;
357 svmul(1.0/mtot, com, com);
359 /* check if COM is outside box */
360 gmx::RVec newCom;
361 copy_rvec(com, newCom);
362 auto newComArrayRef = gmx::arrayRefFromArray(&newCom, 1);
363 switch (unitcell_enum)
365 case euRect:
366 put_atoms_in_box(ePBC, box, newComArrayRef);
367 break;
368 case euTric:
369 put_atoms_in_triclinic_unitcell(ecenter, box, newComArrayRef);
370 break;
371 case euCompact:
372 put_atoms_in_compact_unitcell(ePBC, ecenter, box, newComArrayRef);
373 break;
375 rvec_sub(newCom, com, shift);
376 if (norm2(shift))
378 if (debug)
380 fprintf(debug, "\nShifting position of residue %d (atoms %d-%d) "
381 "by %g,%g,%g\n", atom[res_start].resind+1,
382 res_start+1, res_end+1, shift[XX], shift[YY], shift[ZZ]);
384 for (j = res_start; j < res_end; j++)
386 rvec_inc(x[j], shift);
389 clear_rvec(com);
390 mtot = 0;
392 /* remember start of new residue */
393 res_start = i;
395 if (i < natoms)
397 /* calc COM */
398 m = atom[i].m;
399 for (d = 0; d < DIM; d++)
401 com[d] += m*x[i][d];
403 mtot += m;
405 presnr = atom[i].resind;
410 static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, int ci[])
412 int i, m, ai;
413 rvec cmin, cmax, box_center, dx;
415 if (nc > 0)
417 copy_rvec(x[ci[0]], cmin);
418 copy_rvec(x[ci[0]], cmax);
419 for (i = 0; i < nc; i++)
421 ai = ci[i];
422 for (m = 0; m < DIM; m++)
424 if (x[ai][m] < cmin[m])
426 cmin[m] = x[ai][m];
428 else if (x[ai][m] > cmax[m])
430 cmax[m] = x[ai][m];
434 calc_box_center(ecenter, box, box_center);
435 for (m = 0; m < DIM; m++)
437 dx[m] = box_center[m]-(cmin[m]+cmax[m])*0.5;
440 for (i = 0; i < n; i++)
442 rvec_inc(x[i], dx);
447 static void mk_filenm(char *base, const char *ext, int ndigit, int file_nr,
448 char out_file[])
450 char nbuf[128];
451 int nd = 0, fnr;
453 std::strcpy(out_file, base);
454 fnr = file_nr;
457 fnr /= 10;
458 nd++;
460 while (fnr > 0);
462 if (nd < ndigit)
464 std::strncat(out_file, "00000000000", ndigit-nd);
466 sprintf(nbuf, "%d.", file_nr);
467 std::strcat(out_file, nbuf);
468 std::strcat(out_file, ext);
471 static void check_trr(const char *fn)
473 if (fn2ftp(fn) != efTRR)
475 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
479 static void do_trunc(const char *fn, real t0)
481 t_fileio *in;
482 FILE *fp;
483 gmx_bool bStop, bOK;
484 gmx_trr_header_t sh;
485 gmx_off_t fpos;
486 char yesno[256];
487 int j;
488 real t = 0;
490 if (t0 == -1)
492 gmx_fatal(FARGS, "You forgot to set the truncation time");
495 /* Check whether this is a .trr file */
496 check_trr(fn);
498 in = gmx_trr_open(fn, "r");
499 fp = gmx_fio_getfp(in);
500 if (fp == nullptr)
502 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
503 gmx_trr_close(in);
505 else
507 j = 0;
508 fpos = gmx_fio_ftell(in);
509 bStop = FALSE;
510 while (!bStop && gmx_trr_read_frame_header(in, &sh, &bOK))
512 gmx_trr_read_frame_data(in, &sh, nullptr, nullptr, nullptr, nullptr);
513 fpos = gmx_ftell(fp);
514 t = sh.t;
515 if (t >= t0)
517 gmx_fseek(fp, fpos, SEEK_SET);
518 bStop = TRUE;
521 if (bStop)
523 fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
524 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
525 fn, j, t, (long int)fpos);
526 if (1 != scanf("%s", yesno))
528 gmx_fatal(FARGS, "Error reading user input");
530 if (std::strcmp(yesno, "YES") == 0)
532 fprintf(stderr, "Once again, I'm gonna DO this...\n");
533 gmx_trr_close(in);
534 if (0 != gmx_truncate(fn, fpos))
536 gmx_fatal(FARGS, "Error truncating file %s", fn);
539 else
541 fprintf(stderr, "Ok, I'll forget about it\n");
544 else
546 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
547 gmx_trr_close(in);
552 /*! \brief Read a full molecular topology if useful and available.
554 * If the input trajectory file is not in TNG format, and the output
555 * file is in TNG format, then we want to try to read a full topology
556 * (if available), so that we can write molecule information to the
557 * output file. The full topology provides better molecule information
558 * than is available from the normal t_topology data used by GROMACS
559 * tools.
561 * Also, the t_topology is only read under (different) particular
562 * conditions. If both apply, then a .tpr file might be read
563 * twice. Trying to fix this redundancy while trjconv is still an
564 * all-purpose tool does not seem worthwhile.
566 * Because of the way gmx_prepare_tng_writing is implemented, the case
567 * where the input TNG file has no molecule information will never
568 * lead to an output TNG file having molecule information. Since
569 * molecule information will generally be present if the input TNG
570 * file was written by a GROMACS tool, this seems like reasonable
571 * behaviour. */
572 static gmx_mtop_t *read_mtop_for_tng(const char *tps_file,
573 const char *input_file,
574 const char *output_file)
576 gmx_mtop_t *mtop = nullptr;
578 if (fn2bTPX(tps_file) &&
579 efTNG != fn2ftp(input_file) &&
580 efTNG == fn2ftp(output_file))
582 int temp_natoms = -1;
583 snew(mtop, 1);
584 read_tpx(tps_file, nullptr, nullptr, &temp_natoms,
585 nullptr, nullptr, mtop);
588 return mtop;
591 int gmx_trjconv(int argc, char *argv[])
593 const char *desc[] = {
594 "[THISMODULE] can convert trajectory files in many ways:",
596 "* from one format to another",
597 "* select a subset of atoms",
598 "* change the periodicity representation",
599 "* keep multimeric molecules together",
600 "* center atoms in the box",
601 "* fit atoms to reference structure",
602 "* reduce the number of frames",
603 "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
604 "* cut the trajectory in small subtrajectories according",
605 " to information in an index file. This allows subsequent analysis of",
606 " the subtrajectories that could, for example, be the result of a",
607 " cluster analysis. Use option [TT]-sub[tt].",
608 " This assumes that the entries in the index file are frame numbers and",
609 " dumps each group in the index file to a separate trajectory file.",
610 "* select frames within a certain range of a quantity given",
611 " in an [REF].xvg[ref] file.",
613 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
614 "[PAR]",
616 "The following formats are supported for input and output:",
617 "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
618 "and [REF].pdb[ref].",
619 "The file formats are detected from the file extension.",
620 "The precision of [REF].xtc[ref] and [REF].gro[ref] output is taken from the",
621 "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
622 "and from the [TT]-ndec[tt] option for other input formats. The precision",
623 "is always taken from [TT]-ndec[tt], when this option is set.",
624 "All other formats have fixed precision. [REF].trr[ref]",
625 "output can be single or double precision, depending on the precision",
626 "of the [THISMODULE] binary.",
627 "Note that velocities are only supported in",
628 "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
630 "Option [TT]-sep[tt] can be used to write every frame to a separate",
631 "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to one file.",
632 "[REF].pdb[ref] files with all frames concatenated can be viewed with",
633 "[TT]rasmol -nmrpdb[tt].[PAR]",
635 "It is possible to select part of your trajectory and write it out",
636 "to a new trajectory file in order to save disk space, e.g. for leaving",
637 "out the water from a trajectory of a protein in water.",
638 "[BB]ALWAYS[bb] put the original trajectory on tape!",
639 "We recommend to use the portable [REF].xtc[ref] format for your analysis",
640 "to save disk space and to have portable files.[PAR]",
642 "There are two options for fitting the trajectory to a reference",
643 "either for essential dynamics analysis, etc.",
644 "The first option is just plain fitting to a reference structure",
645 "in the structure file. The second option is a progressive fit",
646 "in which the first timeframe is fitted to the reference structure ",
647 "in the structure file to obtain and each subsequent timeframe is ",
648 "fitted to the previously fitted structure. This way a continuous",
649 "trajectory is generated, which might not be the case when using the",
650 "regular fit method, e.g. when your protein undergoes large",
651 "conformational transitions.[PAR]",
653 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
654 "treatment:",
656 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
657 " and requires a run input file to be supplied with [TT]-s[tt].",
658 " * [TT]res[tt] puts the center of mass of residues in the box.",
659 " * [TT]atom[tt] puts all the atoms in the box.",
660 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
661 " them back. This has the effect that all molecules",
662 " will remain whole (provided they were whole in the initial",
663 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
664 " molecules may diffuse out of the box. The starting configuration",
665 " for this procedure is taken from the structure file, if one is",
666 " supplied, otherwise it is the first frame.",
667 " * [TT]cluster[tt] clusters all the atoms in the selected index",
668 " such that they are all closest to the center of mass of the cluster,",
669 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
670 " results if you in fact have a cluster. Luckily that can be checked",
671 " afterwards using a trajectory viewer. Note also that if your molecules",
672 " are broken this will not work either.",
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 nullptr, "none", "mol", "res", "atom", "nojump", "cluster", "whole",
747 nullptr
750 int unitcell_enum;
751 const char *unitcell_opt[euNR+1] =
752 { nullptr, "rect", "tric", "compact", nullptr };
754 enum
756 ecSel, ecTric, ecRect, ecZero, ecNR
758 const char *center_opt[ecNR+1] =
759 { nullptr, "tric", "rect", "zero", nullptr };
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 nullptr, "none", "rot+trans", "rotxy+transxy", "translation", "transxy",
770 "progressive", nullptr
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 = nullptr;
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 = nullptr;
863 t_trxstatus *trxout = nullptr;
864 t_trxstatus *trxin;
865 int file_nr;
866 t_trxframe fr, frout;
867 int flags;
868 rvec *xmem = nullptr, *vmem = nullptr, *fmem = nullptr;
869 rvec *xp = nullptr, x_shift, hbox;
870 real *w_rls = nullptr;
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 = nullptr;
875 gmx_conect gc = nullptr;
876 int ePBC = -1;
877 t_atoms *atoms = nullptr, useatoms;
878 matrix top_box;
879 int *index = nullptr, *cindex = nullptr;
880 char *grpnm = nullptr;
881 int *frindex, nrfri;
882 char *frname;
883 int ifit, my_clust = -1;
884 int *ind_fit;
885 char *gn_fit;
886 t_cluster_ndx *clust = nullptr;
887 t_trxstatus **clust_status = nullptr;
888 int *clust_status_id = nullptr;
889 int ntrxopen = 0;
890 int *nfwritten = nullptr;
891 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
892 double **dropval;
893 real tshift = 0, dt = -1, prec;
894 gmx_bool bFit, bPFit, bReset;
895 int nfitdim;
896 gmx_rmpbc_t gpbc = nullptr;
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 = nullptr;
904 char out_file2[256], *charpt;
905 char *outf_base = nullptr;
906 const char *outf_ext = nullptr;
907 char top_title[256], title[256], filemode[5];
908 gmx_output_env_t *oenv;
910 t_filenm fnm[] = {
911 { efTRX, "-f", nullptr, ffREAD },
912 { efTRO, "-o", nullptr, ffWRITE },
913 { efTPS, nullptr, nullptr, ffOPTRD },
914 { efNDX, nullptr, nullptr, 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, nullptr, &oenv))
927 return 0;
930 top_file = ftp2fn(efTPS, NFILE, fnm);
932 /* Check command line */
933 in_file = opt2fn("-f", NFILE, fnm);
935 if (ttrunc != -1)
937 do_trunc(in_file, ttrunc);
939 else
941 /* mark active cmdline options */
942 bSetBox = opt2parg_bSet("-box", NPA, pa);
943 bSetTime = opt2parg_bSet("-t0", NPA, pa);
944 bSetPrec = opt2parg_bSet("-ndec", NPA, pa);
945 bSetUR = opt2parg_bSet("-ur", NPA, pa);
946 bExec = opt2parg_bSet("-exec", NPA, pa);
947 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
948 bTDump = opt2parg_bSet("-dump", NPA, pa);
949 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
950 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
951 bTrans = opt2parg_bSet("-trans", NPA, pa);
952 bSplit = (split_t != 0);
954 /* parse enum options */
955 fit_enum = nenum(fit);
956 bFit = (fit_enum == efFit || fit_enum == efFitXY);
957 bReset = (fit_enum == efReset || fit_enum == efResetXY);
958 bPFit = fit_enum == efPFit;
959 pbc_enum = nenum(pbc_opt);
960 bPBCWhole = pbc_enum == epWhole;
961 bPBCcomRes = pbc_enum == epComRes;
962 bPBCcomMol = pbc_enum == epComMol;
963 bPBCcomAtom = pbc_enum == epComAtom;
964 bNoJump = pbc_enum == epNojump;
965 bCluster = pbc_enum == epCluster;
966 bPBC = pbc_enum != epNone;
967 unitcell_enum = nenum(unitcell_opt);
968 ecenter = nenum(center_opt) - ecTric;
970 /* set and check option dependencies */
971 if (bPFit)
973 bFit = TRUE; /* for pfit, fit *must* be set */
975 if (bFit)
977 bReset = TRUE; /* for fit, reset *must* be set */
979 nfitdim = 0;
980 if (bFit || bReset)
982 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
984 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
986 if (bSetUR)
988 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
990 fprintf(stderr,
991 "WARNING: Option for unitcell representation (-ur %s)\n"
992 " only has effect in combination with -pbc %s, %s or %s.\n"
993 " Ingoring unitcell representation.\n\n",
994 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
997 if (bFit && bPBC)
999 gmx_fatal(FARGS, "PBC condition treatment does not work together with rotational fit.\n"
1000 "Please do the PBC condition treatment first and then run trjconv in a second step\n"
1001 "for the rotational fit.\n"
1002 "First doing the rotational fit and then doing the PBC treatment gives incorrect\n"
1003 "results!");
1006 /* ndec is in nr of decimal places, prec is a multiplication factor: */
1007 prec = 1;
1008 for (i = 0; i < ndec; i++)
1010 prec *= 10;
1013 bIndex = ftp2bSet(efNDX, NFILE, fnm);
1016 /* Determine output type */
1017 out_file = opt2fn("-o", NFILE, fnm);
1018 int ftp = fn2ftp(out_file);
1019 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
1020 bNeedPrec = (ftp == efXTC || ftp == efGRO);
1021 int ftpin = fn2ftp(in_file);
1022 if (bVels)
1024 /* check if velocities are possible in input and output files */
1025 bVels = (ftp == efTRR || ftp == efGRO ||
1026 ftp == efG96 || ftp == efTNG)
1027 && (ftpin == efTRR || ftpin == efGRO ||
1028 ftpin == efG96 || ftpin == efTNG || ftpin == efCPT);
1030 if (bSeparate || bSplit)
1032 outf_ext = std::strrchr(out_file, '.');
1033 if (outf_ext == nullptr)
1035 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
1037 outf_base = gmx_strdup(out_file);
1038 outf_base[outf_ext - out_file] = '\0';
1041 bSubTraj = opt2bSet("-sub", NFILE, fnm);
1042 if (bSubTraj)
1044 if ((ftp != efXTC) && (ftp != efTRR))
1046 /* It seems likely that other trajectory file types
1047 * could work here. */
1048 gmx_fatal(FARGS, "Can only use the sub option with output file types "
1049 "xtc and trr");
1051 clust = cluster_index(nullptr, opt2fn("-sub", NFILE, fnm));
1053 /* Check for number of files disabled, as FOPEN_MAX is not the correct
1054 * number to check for. In my linux box it is only 16.
1056 if (0 && (clust->clust->nr > FOPEN_MAX-4))
1058 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
1059 " trajectories.\ntry splitting the index file in %d parts.\n"
1060 "FOPEN_MAX = %d",
1061 clust->clust->nr, 1+clust->clust->nr/FOPEN_MAX, FOPEN_MAX);
1063 gmx_warning("The -sub option could require as many open output files as there are\n"
1064 "index groups in the file (%d). If you get I/O errors opening new files,\n"
1065 "try reducing the number of index groups in the file, and perhaps\n"
1066 "using trjconv -sub several times on different chunks of your index file.\n",
1067 clust->clust->nr);
1069 snew(clust_status, clust->clust->nr);
1070 snew(clust_status_id, clust->clust->nr);
1071 snew(nfwritten, clust->clust->nr);
1072 for (i = 0; (i < clust->clust->nr); i++)
1074 clust_status[i] = nullptr;
1075 clust_status_id[i] = -1;
1077 bSeparate = bSplit = FALSE;
1079 /* skipping */
1080 if (skip_nr <= 0)
1084 mtop = read_mtop_for_tng(top_file, in_file, out_file);
1086 /* Determine whether to read a topology */
1087 bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
1088 bRmPBC || bReset || bPBCcomMol || bCluster ||
1089 (ftp == efGRO) || (ftp == efPDB) || bCONECT);
1091 /* Determine if when can read index groups */
1092 bIndex = (bIndex || bTPS);
1094 if (bTPS)
1096 read_tps_conf(top_file, &top, &ePBC, &xp, nullptr, top_box,
1097 bReset || bPBCcomRes);
1098 std::strncpy(top_title, *top.name, 255);
1099 top_title[255] = '\0';
1100 atoms = &top.atoms;
1102 if (0 == top.mols.nr && (bCluster || bPBCcomMol))
1104 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
1107 /* top_title is only used for gro and pdb,
1108 * the header in such a file is top_title t= ...
1109 * to prevent a double t=, remove it from top_title
1111 if ((charpt = std::strstr(top_title, " t= ")))
1113 charpt[0] = '\0';
1116 if (bCONECT)
1118 gc = gmx_conect_generate(&top);
1120 if (bRmPBC)
1122 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1126 /* get frame number index */
1127 frindex = nullptr;
1128 if (opt2bSet("-fr", NFILE, fnm))
1130 printf("Select groups of frame number indices:\n");
1131 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (int **)&frindex, &frname);
1132 if (debug)
1134 for (i = 0; i < nrfri; i++)
1136 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
1141 /* get index groups etc. */
1142 if (bReset)
1144 printf("Select group for %s fit\n",
1145 bFit ? "least squares" : "translational");
1146 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1147 1, &ifit, &ind_fit, &gn_fit);
1149 if (bFit)
1151 if (ifit < 2)
1153 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
1155 else if (ifit == 3)
1157 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
1161 else if (bCluster)
1163 printf("Select group for clustering\n");
1164 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1165 1, &ifit, &ind_fit, &gn_fit);
1168 if (bIndex)
1170 if (bCenter)
1172 printf("Select group for centering\n");
1173 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1174 1, &ncent, &cindex, &grpnm);
1176 printf("Select group for output\n");
1177 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm),
1178 1, &nout, &index, &grpnm);
1180 else
1182 /* no index file, so read natoms from TRX */
1183 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
1185 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
1187 natoms = fr.natoms;
1188 close_trx(trxin);
1189 sfree(fr.x);
1190 snew(index, natoms);
1191 for (i = 0; i < natoms; i++)
1193 index[i] = i;
1195 nout = natoms;
1196 if (bCenter)
1198 ncent = nout;
1199 cindex = index;
1203 if (bReset)
1205 snew(w_rls, atoms->nr);
1206 for (i = 0; (i < ifit); i++)
1208 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
1211 /* Restore reference structure and set to origin,
1212 store original location (to put structure back) */
1213 if (bRmPBC)
1215 gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
1217 copy_rvec(xp[index[0]], x_shift);
1218 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
1219 rvec_dec(x_shift, xp[index[0]]);
1221 else
1223 clear_rvec(x_shift);
1226 if (bDropUnder || bDropOver)
1228 /* Read the .xvg file with the drop values */
1229 fprintf(stderr, "\nReading drop file ...");
1230 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
1231 fprintf(stderr, " %d time points\n", ndrop);
1232 if (ndrop == 0 || ncol < 2)
1234 gmx_fatal(FARGS, "Found no data points in %s",
1235 opt2fn("-drop", NFILE, fnm));
1237 drop0 = 0;
1238 drop1 = 0;
1241 /* Make atoms struct for output in GRO or PDB files */
1242 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
1244 /* get memory for stuff to go in .pdb file, and initialize
1245 * the pdbinfo structure part if the input has it.
1247 init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
1248 sfree(useatoms.resinfo);
1249 useatoms.resinfo = atoms->resinfo;
1250 for (i = 0; (i < nout); i++)
1252 useatoms.atomname[i] = atoms->atomname[index[i]];
1253 useatoms.atom[i] = atoms->atom[index[i]];
1254 if (atoms->havePdbInfo)
1256 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
1258 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind+1);
1260 useatoms.nr = nout;
1262 /* select what to read */
1263 if (ftp == efTRR)
1265 flags = TRX_READ_X;
1267 else
1269 flags = TRX_NEED_X;
1271 if (bVels)
1273 flags = flags | TRX_READ_V;
1275 if (bForce)
1277 flags = flags | TRX_READ_F;
1280 /* open trx file for reading */
1281 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
1282 if (fr.bPrec)
1284 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
1286 if (bNeedPrec)
1288 if (bSetPrec || !fr.bPrec)
1290 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1/prec);
1292 else
1294 fprintf(stderr, "Using output precision of %g (nm)\n", 1/prec);
1298 if (bHaveFirstFrame)
1300 if (bTDump)
1302 // Determine timestep (assuming constant spacing for now) if we
1303 // need to dump frames based on time. This is required so we do not
1304 // skip the first frame if that was the one that should have been dumped
1305 double firstFrameTime = fr.time;
1306 if (read_next_frame(oenv, trxin, &fr))
1308 dt = fr.time - firstFrameTime;
1309 bDTset = TRUE;
1310 if (dt <= 0)
1312 fprintf(stderr, "Warning: Frame times are not incrementing - will dump first frame.\n");
1315 // Now close and reopen so we are at first frame again
1316 close_trx(trxin);
1317 done_frame(&fr);
1318 // Reopen at first frame (We already know it exists if we got here)
1319 read_first_frame(oenv, &trxin, in_file, &fr, flags);
1322 set_trxframe_ePBC(&fr, ePBC);
1323 natoms = fr.natoms;
1325 if (bSetTime)
1327 tshift = tzero-fr.time;
1329 else
1331 tzero = fr.time;
1334 bCopy = FALSE;
1335 if (bIndex)
1337 /* check if index is meaningful */
1338 for (i = 0; i < nout; i++)
1340 if (index[i] >= natoms)
1342 gmx_fatal(FARGS,
1343 "Index[%d] %d is larger than the number of atoms in the\n"
1344 "trajectory file (%d). There is a mismatch in the contents\n"
1345 "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
1347 bCopy = bCopy || (i != index[i]);
1351 /* open output for writing */
1352 std::strcpy(filemode, "w");
1353 switch (ftp)
1355 case efTNG:
1356 trjtools_gmx_prepare_tng_writing(out_file,
1357 filemode[0],
1358 trxin,
1359 &trxout,
1360 nullptr,
1361 nout,
1362 mtop,
1363 index,
1364 grpnm);
1365 break;
1366 case efXTC:
1367 case efTRR:
1368 out = nullptr;
1369 if (!bSplit && !bSubTraj)
1371 trxout = open_trx(out_file, filemode);
1373 break;
1374 case efGRO:
1375 case efG96:
1376 case efPDB:
1377 if (( !bSeparate && !bSplit ) && !bSubTraj)
1379 out = gmx_ffopen(out_file, filemode);
1381 break;
1382 default:
1383 gmx_incons("Illegal output file format");
1386 if (bCopy)
1388 snew(xmem, nout);
1389 if (bVels)
1391 snew(vmem, nout);
1393 if (bForce)
1395 snew(fmem, nout);
1399 /* Start the big loop over frames */
1400 file_nr = 0;
1401 frame = 0;
1402 outframe = 0;
1403 model_nr = 0;
1405 /* Main loop over frames */
1408 if (!fr.bStep)
1410 /* set the step */
1411 fr.step = newstep;
1412 newstep++;
1414 if (bSubTraj)
1416 /*if (frame >= clust->clust->nra)
1417 gmx_fatal(FARGS,"There are more frames in the trajectory than in the cluster index file\n");*/
1418 if (frame > clust->maxframe)
1420 my_clust = -1;
1422 else
1424 my_clust = clust->inv_clust[frame];
1426 if ((my_clust < 0) || (my_clust >= clust->clust->nr) ||
1427 (my_clust == -1))
1429 my_clust = -1;
1433 if (bSetBox)
1435 /* generate new box */
1436 if (fr.bBox == FALSE)
1438 clear_mat(fr.box);
1440 for (m = 0; m < DIM; m++)
1442 if (newbox[m] >= 0)
1444 fr.box[m][m] = newbox[m];
1446 else
1448 if (fr.bBox == FALSE)
1450 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1456 if (bTrans)
1458 for (i = 0; i < natoms; i++)
1460 rvec_inc(fr.x[i], trans);
1464 if (bTDump)
1466 // If we could not read two frames or times are not incrementing
1467 // we have almost no idea what to do,
1468 // but dump the first frame so output is not broken.
1469 if (dt <= 0 || !bDTset)
1471 bDumpFrame = true;
1473 else
1475 // Dump the frame if we are less than half a frame time
1476 // below it. This will also ensure we at least dump a
1477 // somewhat reasonable frame if the spacing is unequal
1478 // and we have overrun the frame time. Once we dump one
1479 // frame based on time we quit, so it does not matter
1480 // that this might be true for all subsequent frames too.
1481 bDumpFrame = (fr.time > tdump-0.5*dt);
1484 else
1486 bDumpFrame = FALSE;
1489 /* determine if an atom jumped across the box and reset it if so */
1490 if (bNoJump && (bTPS || frame != 0))
1492 for (d = 0; d < DIM; d++)
1494 hbox[d] = 0.5*fr.box[d][d];
1496 for (i = 0; i < natoms; i++)
1498 if (bReset)
1500 rvec_dec(fr.x[i], x_shift);
1502 for (m = DIM-1; m >= 0; m--)
1504 if (hbox[m] > 0)
1506 while (fr.x[i][m]-xp[i][m] <= -hbox[m])
1508 for (d = 0; d <= m; d++)
1510 fr.x[i][d] += fr.box[m][d];
1513 while (fr.x[i][m]-xp[i][m] > hbox[m])
1515 for (d = 0; d <= m; d++)
1517 fr.x[i][d] -= fr.box[m][d];
1524 else if (bCluster)
1526 calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
1529 if (bPFit)
1531 /* Now modify the coords according to the flags,
1532 for normal fit, this is only done for output frames */
1533 if (bRmPBC)
1535 gmx_rmpbc_trxfr(gpbc, &fr);
1538 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1539 do_fit(natoms, w_rls, xp, fr.x);
1542 /* store this set of coordinates for future use */
1543 if (bPFit || bNoJump)
1545 if (xp == nullptr)
1547 snew(xp, natoms);
1549 for (i = 0; (i < natoms); i++)
1551 copy_rvec(fr.x[i], xp[i]);
1552 rvec_inc(fr.x[i], x_shift);
1556 if (frindex)
1558 /* see if we have a frame from the frame index group */
1559 for (i = 0; i < nrfri && !bDumpFrame; i++)
1561 bDumpFrame = frame == frindex[i];
1564 if (debug && bDumpFrame)
1566 fprintf(debug, "dumping %d\n", frame);
1569 bWriteFrame =
1570 ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
1572 if (bWriteFrame && (bDropUnder || bDropOver))
1574 while (dropval[0][drop1] < fr.time && drop1+1 < ndrop)
1576 drop0 = drop1;
1577 drop1++;
1579 if (std::abs(dropval[0][drop0] - fr.time)
1580 < std::abs(dropval[0][drop1] - fr.time))
1582 dropuse = drop0;
1584 else
1586 dropuse = drop1;
1588 if ((bDropUnder && dropval[1][dropuse] < dropunder) ||
1589 (bDropOver && dropval[1][dropuse] > dropover))
1591 bWriteFrame = FALSE;
1595 if (bWriteFrame)
1597 /* We should avoid modifying the input frame,
1598 * but since here we don't have the output frame yet,
1599 * we introduce a temporary output frame time variable.
1601 real frout_time;
1603 frout_time = fr.time;
1605 /* calc new time */
1606 if (bTimeStep)
1608 frout_time = tzero + frame*timestep;
1610 else
1611 if (bSetTime)
1613 frout_time += tshift;
1616 if (bTDump)
1618 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1619 output_env_conv_time(oenv, frout_time), output_env_get_time_unit(oenv).c_str());
1622 /* check for writing at each delta_t */
1623 bDoIt = (delta_t == 0);
1624 if (!bDoIt)
1626 if (!bRound)
1628 bDoIt = bRmod(frout_time, tzero, delta_t);
1630 else
1632 /* round() is not C89 compatible, so we do this: */
1633 bDoIt = bRmod(std::floor(frout_time+0.5), std::floor(tzero+0.5),
1634 std::floor(delta_t+0.5));
1638 if (bDoIt || bTDump)
1640 /* print sometimes */
1641 if ( ((outframe % SKIP) == 0) || (outframe < SKIP) )
1643 fprintf(stderr, " -> frame %6d time %8.3f \r",
1644 outframe, output_env_conv_time(oenv, frout_time));
1645 fflush(stderr);
1648 if (!bPFit)
1650 /* Now modify the coords according to the flags,
1651 for PFit we did this already! */
1653 if (bRmPBC)
1655 gmx_rmpbc_trxfr(gpbc, &fr);
1658 if (bReset)
1660 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1661 if (bFit)
1663 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1665 if (!bCenter)
1667 for (i = 0; i < natoms; i++)
1669 rvec_inc(fr.x[i], x_shift);
1674 if (bCenter)
1676 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1680 auto positionsArrayRef = gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec *>(fr.x), natoms);
1681 if (bPBCcomAtom)
1683 switch (unitcell_enum)
1685 case euRect:
1686 put_atoms_in_box(ePBC, fr.box, positionsArrayRef);
1687 break;
1688 case euTric:
1689 put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1690 break;
1691 case euCompact:
1692 put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box,
1693 positionsArrayRef);
1694 break;
1697 if (bPBCcomRes)
1699 put_residue_com_in_box(unitcell_enum, ecenter,
1700 natoms, atoms->atom, ePBC, fr.box, fr.x);
1702 if (bPBCcomMol)
1704 put_molecule_com_in_box(unitcell_enum, ecenter,
1705 &top.mols,
1706 natoms, atoms->atom, ePBC, fr.box, fr.x);
1708 /* Copy the input trxframe struct to the output trxframe struct */
1709 frout = fr;
1710 frout.time = frout_time;
1711 frout.bV = (frout.bV && bVels);
1712 frout.bF = (frout.bF && bForce);
1713 frout.natoms = nout;
1714 if (bNeedPrec && (bSetPrec || !fr.bPrec))
1716 frout.bPrec = TRUE;
1717 frout.prec = prec;
1719 if (bCopy)
1721 frout.x = xmem;
1722 if (frout.bV)
1724 frout.v = vmem;
1726 if (frout.bF)
1728 frout.f = fmem;
1730 for (i = 0; i < nout; i++)
1732 copy_rvec(fr.x[index[i]], frout.x[i]);
1733 if (frout.bV)
1735 copy_rvec(fr.v[index[i]], frout.v[i]);
1737 if (frout.bF)
1739 copy_rvec(fr.f[index[i]], frout.f[i]);
1744 if (opt2parg_bSet("-shift", NPA, pa))
1746 for (i = 0; i < nout; i++)
1748 for (d = 0; d < DIM; d++)
1750 frout.x[i][d] += outframe*shift[d];
1755 if (!bRound)
1757 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1759 else
1761 /* round() is not C89 compatible, so we do this: */
1762 bSplitHere = bSplit && bRmod(std::floor(frout.time+0.5),
1763 std::floor(tzero+0.5),
1764 std::floor(split_t+0.5));
1766 if (bSeparate || bSplitHere)
1768 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1771 switch (ftp)
1773 case efTNG:
1774 write_tng_frame(trxout, &frout);
1775 // TODO when trjconv behaves better: work how to read and write lambda
1776 break;
1777 case efTRR:
1778 case efXTC:
1779 if (bSplitHere)
1781 if (trxout)
1783 close_trx(trxout);
1785 trxout = open_trx(out_file2, filemode);
1787 if (bSubTraj)
1789 if (my_clust != -1)
1791 char buf[STRLEN];
1792 if (clust_status_id[my_clust] == -1)
1794 sprintf(buf, "%s.%s", clust->grpname[my_clust], ftp2ext(ftp));
1795 clust_status[my_clust] = open_trx(buf, "w");
1796 clust_status_id[my_clust] = 1;
1797 ntrxopen++;
1799 else if (clust_status_id[my_clust] == -2)
1801 gmx_fatal(FARGS, "File %s.xtc should still be open (%d open .xtc files)\n" "in order to write frame %d. my_clust = %d",
1802 clust->grpname[my_clust], ntrxopen, frame,
1803 my_clust);
1805 write_trxframe(clust_status[my_clust], &frout, gc);
1806 nfwritten[my_clust]++;
1807 if (nfwritten[my_clust] ==
1808 (clust->clust->index[my_clust+1]-
1809 clust->clust->index[my_clust]))
1811 close_trx(clust_status[my_clust]);
1812 clust_status[my_clust] = nullptr;
1813 clust_status_id[my_clust] = -2;
1814 ntrxopen--;
1815 if (ntrxopen < 0)
1817 gmx_fatal(FARGS, "Less than zero open .xtc files!");
1822 else
1824 write_trxframe(trxout, &frout, gc);
1826 break;
1827 case efGRO:
1828 case efG96:
1829 case efPDB:
1830 sprintf(title, "Generated by trjconv : %s t= %9.5f",
1831 top_title, frout.time);
1832 if (bSeparate || bSplitHere)
1834 out = gmx_ffopen(out_file2, "w");
1836 switch (ftp)
1838 case efGRO:
1839 write_hconf_p(out, title, &useatoms,
1840 frout.x, frout.bV ? frout.v : nullptr, frout.box);
1841 break;
1842 case efPDB:
1843 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1844 sprintf(title, "%s t= %9.5f", top_title, frout.time);
1845 /* if reading from pdb, we want to keep the original
1846 model numbering else we write the output frame
1847 number plus one, because model 0 is not allowed in pdb */
1848 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1850 model_nr = fr.step;
1852 else
1854 model_nr++;
1856 write_pdbfile(out, title, &useatoms, frout.x,
1857 frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
1858 break;
1859 case efG96:
1860 const char *outputTitle = "";
1861 if (bSeparate || bTDump)
1863 outputTitle = title;
1864 if (bTPS)
1866 frout.bAtoms = TRUE;
1868 frout.atoms = &useatoms;
1869 frout.bStep = FALSE;
1870 frout.bTime = FALSE;
1872 else
1874 if (outframe == 0)
1876 outputTitle = title;
1878 frout.bAtoms = FALSE;
1879 frout.bStep = TRUE;
1880 frout.bTime = TRUE;
1882 write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1884 if (bSeparate || bSplitHere)
1886 gmx_ffclose(out);
1887 out = nullptr;
1889 break;
1890 default:
1891 gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1893 if (bSeparate || bSplitHere)
1895 file_nr++;
1898 /* execute command */
1899 if (bExec)
1901 char c[255];
1902 sprintf(c, "%s %d", exec_command, file_nr-1);
1903 /*fprintf(stderr,"Executing '%s'\n",c);*/
1904 if (0 != system(c))
1906 gmx_fatal(FARGS, "Error executing command: %s", c);
1909 outframe++;
1912 frame++;
1913 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1915 while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1918 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1920 fprintf(stderr, "\nWARNING no output, "
1921 "last frame read at t=%g\n", fr.time);
1923 fprintf(stderr, "\n");
1925 close_trx(trxin);
1926 sfree(outf_base);
1928 if (bRmPBC)
1930 gmx_rmpbc_done(gpbc);
1933 if (trxout)
1935 close_trx(trxout);
1937 else if (out != nullptr)
1939 gmx_ffclose(out);
1941 if (bSubTraj)
1943 for (i = 0; (i < clust->clust->nr); i++)
1945 if (clust_status_id[i] >= 0)
1947 close_trx(clust_status[i]);
1953 sfree(mtop);
1954 done_top(&top);
1955 sfree(xp);
1956 sfree(xmem);
1957 sfree(vmem);
1958 sfree(fmem);
1959 sfree(grpnm);
1960 sfree(index);
1961 sfree(cindex);
1962 done_filenms(NFILE, fnm);
1963 done_frame(&fr);
1965 do_view(oenv, out_file, nullptr);
1967 output_env_done(oenv);
1968 return 0;