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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
49 #include "gromacs/commandline/pargs.h"
50 #include "gromacs/commandline/viewit.h"
51 #include "gromacs/fileio/confio.h"
52 #include "gromacs/fileio/g96io.h"
53 #include "gromacs/fileio/gmxfio.h"
54 #include "gromacs/fileio/groio.h"
55 #include "gromacs/fileio/pdbio.h"
56 #include "gromacs/fileio/tngio.h"
57 #include "gromacs/fileio/tpxio.h"
58 #include "gromacs/fileio/trrio.h"
59 #include "gromacs/fileio/trxio.h"
60 #include "gromacs/fileio/xtcio.h"
61 #include "gromacs/fileio/xvgr.h"
62 #include "gromacs/math/do_fit.h"
63 #include "gromacs/math/functions.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/mdtypes/md_enums.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/pbcutil/pbcmethods.h"
68 #include "gromacs/pbcutil/rmpbc.h"
69 #include "gromacs/topology/index.h"
70 #include "gromacs/topology/topology.h"
71 #include "gromacs/trajectory/trajectoryframe.h"
72 #include "gromacs/utility/arrayref.h"
73 #include "gromacs/utility/arraysize.h"
74 #include "gromacs/utility/fatalerror.h"
75 #include "gromacs/utility/futil.h"
76 #include "gromacs/utility/smalloc.h"
78 static void mk_filenm(char* base, const char* ext, int ndigit, int file_nr, char out_file[])
83 std::strcpy(out_file, base);
93 std::strncat(out_file, "00000000000", ndigit - nd);
95 sprintf(nbuf, "%d.", file_nr);
96 std::strcat(out_file, nbuf);
97 std::strcat(out_file, ext);
100 static void check_trr(const char* fn)
102 if (fn2ftp(fn) != efTRR)
104 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
108 static void do_trunc(const char* fn, real t0)
121 gmx_fatal(FARGS, "You forgot to set the truncation time");
124 /* Check whether this is a .trr file */
127 in = gmx_trr_open(fn, "r");
128 fp = gmx_fio_getfp(in);
131 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
137 fpos = gmx_fio_ftell(in);
139 while (!bStop && gmx_trr_read_frame_header(in, &sh, &bOK))
141 gmx_trr_read_frame_data(in, &sh, nullptr, nullptr, nullptr, nullptr);
142 fpos = gmx_ftell(fp);
146 gmx_fseek(fp, fpos, SEEK_SET);
153 "Do you REALLY want to truncate this trajectory (%s) at:\n"
154 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
158 static_cast<long int>(fpos));
159 if (1 != scanf("%s", yesno))
161 gmx_fatal(FARGS, "Error reading user input");
163 if (std::strcmp(yesno, "YES") == 0)
165 fprintf(stderr, "Once again, I'm gonna DO this...\n");
167 if (0 != gmx_truncate(fn, fpos))
169 gmx_fatal(FARGS, "Error truncating file %s", fn);
174 fprintf(stderr, "Ok, I'll forget about it\n");
179 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
185 /*! \brief Read a full molecular topology if useful and available.
187 * If the input trajectory file is not in TNG format, and the output
188 * file is in TNG format, then we want to try to read a full topology
189 * (if available), so that we can write molecule information to the
190 * output file. The full topology provides better molecule information
191 * than is available from the normal t_topology data used by GROMACS
194 * Also, the t_topology is only read under (different) particular
195 * conditions. If both apply, then a .tpr file might be read
196 * twice. Trying to fix this redundancy while trjconv is still an
197 * all-purpose tool does not seem worthwhile.
199 * Because of the way gmx_prepare_tng_writing is implemented, the case
200 * where the input TNG file has no molecule information will never
201 * lead to an output TNG file having molecule information. Since
202 * molecule information will generally be present if the input TNG
203 * file was written by a GROMACS tool, this seems like reasonable
205 static std::unique_ptr<gmx_mtop_t> read_mtop_for_tng(const char* tps_file,
206 const char* input_file,
207 const char* output_file)
209 std::unique_ptr<gmx_mtop_t> mtop;
211 if (fn2bTPX(tps_file) && efTNG != fn2ftp(input_file) && efTNG == fn2ftp(output_file))
213 int temp_natoms = -1;
214 mtop = std::make_unique<gmx_mtop_t>();
215 read_tpx(tps_file, nullptr, nullptr, &temp_natoms, nullptr, nullptr, mtop.get());
221 int gmx_trjconv(int argc, char* argv[])
223 const char* desc[] = {
224 "[THISMODULE] can convert trajectory files in many ways:",
226 "* from one format to another",
227 "* select a subset of atoms",
228 "* change the periodicity representation",
229 "* keep multimeric molecules together",
230 "* center atoms in the box",
231 "* fit atoms to reference structure",
232 "* reduce the number of frames",
233 "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
234 "* select frames within a certain range of a quantity given",
235 " in an [REF].xvg[ref] file.",
237 "The option to write subtrajectories (-sub) based on the information obtained from",
238 "cluster analysis has been removed from [THISMODULE] and is now part of",
239 "[gmx extract-cluster]",
241 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
244 "The following formats are supported for input and output:",
245 "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
246 "and [REF].pdb[ref].",
247 "The file formats are detected from the file extension.",
248 "The precision of the [REF].xtc[ref] output is taken from the",
249 "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
250 "and from the [TT]-ndec[tt] option for other input formats. The precision",
251 "is always taken from [TT]-ndec[tt], when this option is set.",
252 "All other formats have fixed precision. [REF].trr[ref]",
253 "output can be single or double precision, depending on the precision",
254 "of the [THISMODULE] binary.",
255 "Note that velocities are only supported in",
256 "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
258 "Option [TT]-sep[tt] can be used to write every frame to a separate",
259 "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to ",
261 "[REF].pdb[ref] files with all frames concatenated can be viewed with",
262 "[TT]rasmol -nmrpdb[tt].[PAR]",
264 "It is possible to select part of your trajectory and write it out",
265 "to a new trajectory file in order to save disk space, e.g. for leaving",
266 "out the water from a trajectory of a protein in water.",
267 "[BB]ALWAYS[bb] put the original trajectory on tape!",
268 "We recommend to use the portable [REF].xtc[ref] format for your analysis",
269 "to save disk space and to have portable files.[PAR]",
271 "There are two options for fitting the trajectory to a reference",
272 "either for essential dynamics analysis, etc.",
273 "The first option is just plain fitting to a reference structure",
274 "in the structure file. The second option is a progressive fit",
275 "in which the first timeframe is fitted to the reference structure ",
276 "in the structure file to obtain and each subsequent timeframe is ",
277 "fitted to the previously fitted structure. This way a continuous",
278 "trajectory is generated, which might not be the case when using the",
279 "regular fit method, e.g. when your protein undergoes large",
280 "conformational transitions.[PAR]",
282 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
285 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
286 " and requires a run input file to be supplied with [TT]-s[tt].",
287 " * [TT]res[tt] puts the center of mass of residues in the box.",
288 " * [TT]atom[tt] puts all the atoms in the box.",
289 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
290 " them back. This has the effect that all molecules",
291 " will remain whole (provided they were whole in the initial",
292 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
293 " molecules may diffuse out of the box. The starting configuration",
294 " for this procedure is taken from the structure file, if one is",
295 " supplied, otherwise it is the first frame.",
296 " * [TT]cluster[tt] clusters all the atoms in the selected index",
297 " such that they are all closest to the center of mass of the cluster,",
298 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
299 " results if you in fact have a cluster. Luckily that can be checked",
300 " afterwards using a trajectory viewer. Note also that if your molecules",
301 " are broken this will not work either.",
302 " * [TT]whole[tt] only makes broken molecules whole.",
305 "Option [TT]-ur[tt] sets the unit cell representation for options",
306 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
307 "All three options give different results for triclinic boxes and",
308 "identical results for rectangular boxes.",
309 "[TT]rect[tt] is the ordinary brick shape.",
310 "[TT]tric[tt] is the triclinic unit cell.",
311 "[TT]compact[tt] puts all atoms at the closest distance from the center",
312 "of the box. This can be useful for visualizing e.g. truncated octahedra",
313 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
314 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
315 "is set differently.[PAR]",
317 "Option [TT]-center[tt] centers the system in the box. The user can",
318 "select the group which is used to determine the geometrical center.",
319 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
320 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
321 "[TT]tric[tt]: half of the sum of the box vectors,",
322 "[TT]rect[tt]: half of the box diagonal,",
323 "[TT]zero[tt]: zero.",
324 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
325 "want all molecules in the box after the centering.[PAR]",
327 "Option [TT]-box[tt] sets the size of the new box. This option only works",
328 "for leading dimensions and is thus generally only useful for rectangular boxes.",
329 "If you want to modify only some of the dimensions, e.g. when reading from",
330 "a trajectory, you can use -1 for those dimensions that should stay the same",
332 "It is not always possible to use combinations of [TT]-pbc[tt],",
333 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
334 "you want in one call to [THISMODULE]. Consider using multiple",
335 "calls, and check out the GROMACS website for suggestions.[PAR]",
337 "With [TT]-dt[tt], it is possible to reduce the number of ",
338 "frames in the output. This option relies on the accuracy of the times",
339 "in your input trajectory, so if these are inaccurate use the",
340 "[TT]-timestep[tt] option to modify the time (this can be done",
341 "simultaneously). For making smooth movies, the program [gmx-filter]",
342 "can reduce the number of frames while using low-pass frequency",
343 "filtering, this reduces aliasing of high frequency motions.[PAR]",
345 "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
346 "without copying the file. This is useful when a run has crashed",
347 "during disk I/O (i.e. full disk), or when two contiguous",
348 "trajectories must be concatenated without having double frames.[PAR]",
350 "Option [TT]-dump[tt] can be used to extract a frame at or near",
351 "one specific time from your trajectory, but only works reliably",
352 "if the time interval between frames is uniform.[PAR]",
354 "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
355 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
356 "frames with a value below and above the value of the respective options",
357 "will not be written."
373 const char* pbc_opt[epNR + 1] = { nullptr, "none", "mol", "res", "atom",
374 "nojump", "cluster", "whole", nullptr };
377 const char* unitcell_opt[euNR + 1] = { nullptr, "rect", "tric", "compact", nullptr };
387 const char* center_opt[ecNR + 1] = { nullptr, "tric", "rect", "zero", nullptr };
402 const char* fit[efNR + 1] = { nullptr, "none", "rot+trans", "rotxy+transxy",
403 "translation", "transxy", "progressive", nullptr };
405 gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
406 gmx_bool bCenter = FALSE;
407 int skip_nr = 1, ndec = 3, nzero = 0;
408 real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
409 rvec newbox = { 0, 0, 0 }, shift = { 0, 0, 0 }, trans = { 0, 0, 0 };
410 char* exec_command = nullptr;
411 real dropunder = 0, dropover = 0;
412 gmx_bool bRound = FALSE;
415 { "-skip", FALSE, etINT, { &skip_nr }, "Only write every nr-th frame" },
416 { "-dt", FALSE, etTIME, { &delta_t }, "Only write frame when t MOD dt = first time (%t)" },
417 { "-round", FALSE, etBOOL, { &bRound }, "Round measurements to nearest picosecond" },
418 { "-dump", FALSE, etTIME, { &tdump }, "Dump frame nearest specified time (%t)" },
419 { "-t0", FALSE, etTIME, { &tzero }, "Starting time (%t) (default: don't change)" },
420 { "-timestep", FALSE, etTIME, { ×tep }, "Change time step between input frames (%t)" },
421 { "-pbc", FALSE, etENUM, { pbc_opt }, "PBC treatment (see help text for full description)" },
422 { "-ur", FALSE, etENUM, { unitcell_opt }, "Unit-cell representation" },
423 { "-center", FALSE, etBOOL, { &bCenter }, "Center atoms in box" },
424 { "-boxcenter", FALSE, etENUM, { center_opt }, "Center for -pbc and -center" },
425 { "-box", FALSE, etRVEC, { newbox }, "Size for new cubic box (default: read from input)" },
430 "All coordinates will be translated by trans. This "
431 "can advantageously be combined with -pbc mol -ur "
433 { "-shift", FALSE, etRVEC, { shift }, "All coordinates will be shifted by framenr*shift" },
434 { "-fit", FALSE, etENUM, { fit }, "Fit molecule to ref structure in the structure file" },
435 { "-ndec", FALSE, etINT, { &ndec }, "Number of decimal places to write to .xtc output" },
436 { "-vel", FALSE, etBOOL, { &bVels }, "Read and write velocities if possible" },
437 { "-force", FALSE, etBOOL, { &bForce }, "Read and write forces if possible" },
442 "Truncate input trajectory file after this time (%t)" },
447 "Execute command for every output frame with the "
448 "frame number as argument" },
453 "Start writing new file when t MOD split = first "
459 "Write each frame to a separate .gro, .g96 or .pdb "
465 "If the -sep flag is set, use these many digits "
466 "for the file numbers and prepend zeros as needed" },
467 { "-dropunder", FALSE, etREAL, { &dropunder }, "Drop all frames below this value" },
468 { "-dropover", FALSE, etREAL, { &dropover }, "Drop all frames above this value" },
473 "Add CONECT PDB records when writing [REF].pdb[ref] files. Useful "
474 "for visualization of non-standard molecules, e.g. "
475 "coarse grained ones" }
477 #define NPA asize(pa)
480 t_trxstatus* trxout = nullptr;
483 t_trxframe fr, frout;
485 rvec * xmem = nullptr, *vmem = nullptr, *fmem = nullptr;
486 rvec * xp = nullptr, x_shift, hbox;
487 real* w_rls = nullptr;
488 int m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
490 t_topology* top = nullptr;
491 gmx_conect gc = nullptr;
492 PbcType pbcType = PbcType::Unset;
493 t_atoms * atoms = nullptr, useatoms;
495 int * index = nullptr, *cindex = nullptr;
496 char* grpnm = nullptr;
497 int * frindex, nrfri;
502 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
504 real tshift = 0, dt = -1, prec;
505 gmx_bool bFit, bPFit, bReset;
507 gmx_rmpbc_t gpbc = nullptr;
508 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
509 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
510 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetXtcPrec, bNeedPrec;
511 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
512 gmx_bool bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
513 gmx_bool bWriteFrame, bSplitHere;
514 const char *top_file, *in_file, *out_file = nullptr;
515 char out_file2[256], *charpt;
516 char* outf_base = nullptr;
517 const char* outf_ext = nullptr;
518 char top_title[256], timestr[32], stepstr[32], filemode[5];
519 gmx_output_env_t* oenv;
521 t_filenm fnm[] = { { efTRX, "-f", nullptr, ffREAD }, { efTRO, "-o", nullptr, ffWRITE },
522 { efTPS, nullptr, nullptr, ffOPTRD }, { efNDX, nullptr, nullptr, ffOPTRD },
523 { efNDX, "-fr", "frames", ffOPTRD }, { efNDX, "-sub", "cluster", ffOPTRD },
524 { efXVG, "-drop", "drop", ffOPTRD } };
525 #define NFILE asize(fnm)
527 if (!parse_common_args(&argc,
529 PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW | PCA_TIME_UNIT,
543 "Note that major changes are planned in future for "
544 "trjconv, to improve usability and utility.\n");
546 top_file = ftp2fn(efTPS, NFILE, fnm);
548 /* Check command line */
549 in_file = opt2fn("-f", NFILE, fnm);
553 do_trunc(in_file, ttrunc);
557 /* mark active cmdline options */
558 bSetBox = opt2parg_bSet("-box", NPA, pa);
559 bSetTime = opt2parg_bSet("-t0", NPA, pa);
560 bSetXtcPrec = opt2parg_bSet("-ndec", NPA, pa);
561 bSetUR = opt2parg_bSet("-ur", NPA, pa);
562 bExec = opt2parg_bSet("-exec", NPA, pa);
563 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
564 bTDump = opt2parg_bSet("-dump", NPA, pa);
565 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
566 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
567 bTrans = opt2parg_bSet("-trans", NPA, pa);
568 bSplit = (split_t != 0);
570 /* parse enum options */
571 fit_enum = nenum(fit);
572 bFit = (fit_enum == efFit || fit_enum == efFitXY);
573 bReset = (fit_enum == efReset || fit_enum == efResetXY);
574 bPFit = fit_enum == efPFit;
575 pbc_enum = nenum(pbc_opt);
576 bPBCWhole = pbc_enum == epWhole;
577 bPBCcomRes = pbc_enum == epComRes;
578 bPBCcomMol = pbc_enum == epComMol;
579 bPBCcomAtom = pbc_enum == epComAtom;
580 bNoJump = pbc_enum == epNojump;
581 bCluster = pbc_enum == epCluster;
582 bPBC = pbc_enum != epNone;
583 unitcell_enum = nenum(unitcell_opt);
584 ecenter = nenum(center_opt) - ecTric;
586 /* set and check option dependencies */
589 bFit = TRUE; /* for pfit, fit *must* be set */
593 bReset = TRUE; /* for fit, reset *must* be set */
598 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
600 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
604 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
607 "WARNING: Option for unitcell representation (-ur %s)\n"
608 " only has effect in combination with -pbc %s, %s or %s.\n"
609 " Ingoring unitcell representation.\n\n",
619 "PBC condition treatment does not work together with rotational fit.\n"
620 "Please do the PBC condition treatment first and then run trjconv in a "
622 "for the rotational fit.\n"
623 "First doing the rotational fit and then doing the PBC treatment gives "
628 /* ndec for XTC writing is in nr of decimal places, prec is a multiplication factor: */
630 for (i = 0; i < ndec; i++)
635 bIndex = ftp2bSet(efNDX, NFILE, fnm);
638 /* Determine output type */
639 out_file = opt2fn("-o", NFILE, fnm);
640 int ftp = fn2ftp(out_file);
641 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
642 bNeedPrec = (ftp == efXTC);
643 int ftpin = fn2ftp(in_file);
646 /* check if velocities are possible in input and output files */
647 bVels = (ftp == efTRR || ftp == efGRO || ftp == efG96 || ftp == efTNG)
648 && (ftpin == efTRR || ftpin == efGRO || ftpin == efG96 || ftpin == efTNG
651 if (bSeparate || bSplit)
653 outf_ext = std::strrchr(out_file, '.');
654 if (outf_ext == nullptr)
656 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
658 outf_base = gmx_strdup(out_file);
659 outf_base[outf_ext - out_file] = '\0';
662 bool bSubTraj = opt2bSet("-sub", NFILE, fnm);
666 "The -sub option has been removed from gmx trjconv and is now part\n"
667 "of gmx extract-cluster and does nothing here\n");
673 gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
676 std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
678 /* Determine whether to read a topology */
679 bTPS = (ftp2bSet(efTPS, NFILE, fnm) || bRmPBC || bReset || bPBCcomMol || bCluster
680 || (ftp == efGRO) || (ftp == efPDB) || bCONECT);
682 /* Determine if when can read index groups */
683 bIndex = (bIndex || bTPS);
688 read_tps_conf(top_file, top, &pbcType, &xp, nullptr, top_box, bReset || bPBCcomRes);
689 std::strncpy(top_title, *top->name, 255);
690 top_title[255] = '\0';
693 if (0 == top->mols.nr && (bCluster || bPBCcomMol))
695 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
698 /* top_title is only used for gro and pdb,
699 * the header in such a file is top_title, followed by
700 * t= ... and/or step= ...
701 * to prevent double t= or step=, remove it from top_title.
702 * From GROMACS-2018 we only write t/step when the frame actually
703 * has a valid time/step, so we need to check for both separately.
705 if ((charpt = std::strstr(top_title, " t= ")))
709 if ((charpt = std::strstr(top_title, " step= ")))
716 gc = gmx_conect_generate(top);
720 gpbc = gmx_rmpbc_init(&top->idef, pbcType, top->atoms.nr);
724 /* get frame number index */
726 if (opt2bSet("-fr", NFILE, fnm))
728 printf("Select groups of frame number indices:\n");
729 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
732 for (i = 0; i < nrfri; i++)
734 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
739 /* get index groups etc. */
742 printf("Select group for %s fit\n", bFit ? "least squares" : "translational");
743 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit, &ind_fit, &gn_fit);
749 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
753 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
759 printf("Select group for clustering\n");
760 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit, &ind_fit, &gn_fit);
767 printf("Select group for centering\n");
768 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ncent, &cindex, &grpnm);
770 printf("Select group for output\n");
771 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nout, &index, &grpnm);
775 /* no index file, so read natoms from TRX */
776 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
778 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
784 for (i = 0; i < natoms; i++)
798 snew(w_rls, atoms->nr);
799 for (i = 0; (i < ifit); i++)
801 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
804 /* Restore reference structure and set to origin,
805 store original location (to put structure back) */
808 gmx_rmpbc(gpbc, top->atoms.nr, top_box, xp);
810 copy_rvec(xp[index[0]], x_shift);
811 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
812 rvec_dec(x_shift, xp[index[0]]);
819 if (bDropUnder || bDropOver)
821 /* Read the .xvg file with the drop values */
822 fprintf(stderr, "\nReading drop file ...");
823 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
824 fprintf(stderr, " %d time points\n", ndrop);
825 if (ndrop == 0 || ncol < 2)
827 gmx_fatal(FARGS, "Found no data points in %s", opt2fn("-drop", NFILE, fnm));
833 /* Make atoms struct for output in GRO or PDB files */
834 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
836 /* get memory for stuff to go in .pdb file, and initialize
837 * the pdbinfo structure part if the input has it.
839 init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
840 sfree(useatoms.resinfo);
841 useatoms.resinfo = atoms->resinfo;
842 for (i = 0; (i < nout); i++)
844 useatoms.atomname[i] = atoms->atomname[index[i]];
845 useatoms.atom[i] = atoms->atom[index[i]];
846 if (atoms->havePdbInfo)
848 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
850 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind + 1);
854 /* select what to read */
865 flags = flags | TRX_READ_V;
869 flags = flags | TRX_READ_F;
872 /* open trx file for reading */
873 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
876 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1 / fr.prec);
880 if (bSetXtcPrec || !fr.bPrec)
882 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1 / prec);
886 fprintf(stderr, "Using output precision of %g (nm)\n", 1 / prec);
894 // Determine timestep (assuming constant spacing for now) if we
895 // need to dump frames based on time. This is required so we do not
896 // skip the first frame if that was the one that should have been dumped
897 double firstFrameTime = fr.time;
898 if (read_next_frame(oenv, trxin, &fr))
900 dt = fr.time - firstFrameTime;
905 "Warning: Frame times are not incrementing - will dump first "
909 // Now close and reopen so we are at first frame again
912 // Reopen at first frame (We already know it exists if we got here)
913 read_first_frame(oenv, &trxin, in_file, &fr, flags);
916 setTrxFramePbcType(&fr, pbcType);
921 tshift = tzero - fr.time;
931 /* check if index is meaningful */
932 for (i = 0; i < nout; i++)
934 if (index[i] >= natoms)
937 "Index[%d] %d is larger than the number of atoms in the\n"
938 "trajectory file (%d). There is a mismatch in the contents\n"
939 "of your -f, -s and/or -n files.",
944 bCopy = bCopy || (i != index[i]);
948 /* open output for writing */
949 std::strcpy(filemode, "w");
953 trxout = trjtools_gmx_prepare_tng_writing(out_file,
959 gmx::arrayRefFromArray(index, nout),
967 trxout = open_trx(out_file, filemode);
973 if (!bSeparate && !bSplit)
975 out = gmx_ffopen(out_file, filemode);
978 default: gmx_incons("Illegal output file format");
994 /* Start the big loop over frames */
1000 /* Main loop over frames */
1012 /* generate new box */
1017 for (m = 0; m < DIM; m++)
1021 fr.box[m][m] = newbox[m];
1027 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1035 for (i = 0; i < natoms; i++)
1037 rvec_inc(fr.x[i], trans);
1043 // If we could not read two frames or times are not incrementing
1044 // we have almost no idea what to do,
1045 // but dump the first frame so output is not broken.
1046 if (dt <= 0 || !bDTset)
1052 // Dump the frame if we are less than half a frame time
1053 // below it. This will also ensure we at least dump a
1054 // somewhat reasonable frame if the spacing is unequal
1055 // and we have overrun the frame time. Once we dump one
1056 // frame based on time we quit, so it does not matter
1057 // that this might be true for all subsequent frames too.
1058 bDumpFrame = (fr.time > tdump - 0.5 * dt);
1066 /* determine if an atom jumped across the box and reset it if so */
1067 if (bNoJump && (bTPS || frame != 0))
1069 for (d = 0; d < DIM; d++)
1071 hbox[d] = 0.5 * fr.box[d][d];
1073 for (i = 0; i < natoms; i++)
1077 rvec_dec(fr.x[i], x_shift);
1079 for (m = DIM - 1; m >= 0; m--)
1083 while (fr.x[i][m] - xp[i][m] <= -hbox[m])
1085 for (d = 0; d <= m; d++)
1087 fr.x[i][d] += fr.box[m][d];
1090 while (fr.x[i][m] - xp[i][m] > hbox[m])
1092 for (d = 0; d <= m; d++)
1094 fr.x[i][d] -= fr.box[m][d];
1103 calc_pbc_cluster(ecenter, ifit, top, pbcType, fr.x, ind_fit, fr.box);
1108 /* Now modify the coords according to the flags,
1109 for normal fit, this is only done for output frames */
1112 gmx_rmpbc_trxfr(gpbc, &fr);
1115 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1116 do_fit(natoms, w_rls, xp, fr.x);
1119 /* store this set of coordinates for future use */
1120 if (bPFit || bNoJump)
1126 for (i = 0; (i < natoms); i++)
1128 copy_rvec(fr.x[i], xp[i]);
1129 rvec_inc(fr.x[i], x_shift);
1135 /* see if we have a frame from the frame index group */
1136 for (i = 0; i < nrfri && !bDumpFrame; i++)
1138 bDumpFrame = frame == frindex[i];
1141 if (debug && bDumpFrame)
1143 fprintf(debug, "dumping %d\n", frame);
1146 bWriteFrame = ((!bTDump && (frindex == nullptr) && frame % skip_nr == 0) || bDumpFrame);
1148 if (bWriteFrame && (bDropUnder || bDropOver))
1150 while (dropval[0][drop1] < fr.time && drop1 + 1 < ndrop)
1155 if (std::abs(dropval[0][drop0] - fr.time) < std::abs(dropval[0][drop1] - fr.time))
1163 if ((bDropUnder && dropval[1][dropuse] < dropunder)
1164 || (bDropOver && dropval[1][dropuse] > dropover))
1166 bWriteFrame = FALSE;
1172 /* We should avoid modifying the input frame,
1173 * but since here we don't have the output frame yet,
1174 * we introduce a temporary output frame time variable.
1178 frout_time = fr.time;
1183 frout_time = tzero + frame * timestep;
1187 frout_time += tshift;
1193 "\nDumping frame at t= %g %s\n",
1194 output_env_conv_time(oenv, frout_time),
1195 output_env_get_time_unit(oenv).c_str());
1198 /* check for writing at each delta_t */
1199 bDoIt = (delta_t == 0);
1204 bDoIt = bRmod(frout_time, tzero, delta_t);
1208 /* round() is not C89 compatible, so we do this: */
1209 bDoIt = bRmod(std::floor(frout_time + 0.5),
1210 std::floor(tzero + 0.5),
1211 std::floor(delta_t + 0.5));
1215 if (bDoIt || bTDump)
1217 /* print sometimes */
1218 if (((outframe % SKIP) == 0) || (outframe < SKIP))
1221 " -> frame %6d time %8.3f \r",
1223 output_env_conv_time(oenv, frout_time));
1229 /* Now modify the coords according to the flags,
1230 for PFit we did this already! */
1234 gmx_rmpbc_trxfr(gpbc, &fr);
1239 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1242 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1246 for (i = 0; i < natoms; i++)
1248 rvec_inc(fr.x[i], x_shift);
1255 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1259 auto positionsArrayRef =
1260 gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec*>(fr.x), natoms);
1263 switch (unitcell_enum)
1266 put_atoms_in_box(pbcType, fr.box, positionsArrayRef);
1269 put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1272 put_atoms_in_compact_unitcell(
1273 pbcType, ecenter, fr.box, positionsArrayRef);
1279 put_residue_com_in_box(
1280 unitcell_enum, ecenter, natoms, atoms->atom, pbcType, fr.box, fr.x);
1284 put_molecule_com_in_box(
1285 unitcell_enum, ecenter, &top->mols, natoms, atoms->atom, pbcType, fr.box, fr.x);
1287 /* Copy the input trxframe struct to the output trxframe struct */
1289 frout.time = frout_time;
1290 frout.bV = (frout.bV && bVels);
1291 frout.bF = (frout.bF && bForce);
1292 frout.natoms = nout;
1293 if (bNeedPrec && (bSetXtcPrec || !fr.bPrec))
1309 for (i = 0; i < nout; i++)
1311 copy_rvec(fr.x[index[i]], frout.x[i]);
1314 copy_rvec(fr.v[index[i]], frout.v[i]);
1318 copy_rvec(fr.f[index[i]], frout.f[i]);
1323 if (opt2parg_bSet("-shift", NPA, pa))
1325 for (i = 0; i < nout; i++)
1327 for (d = 0; d < DIM; d++)
1329 frout.x[i][d] += outframe * shift[d];
1336 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1340 /* round() is not C89 compatible, so we do this: */
1342 && bRmod(std::floor(frout.time + 0.5),
1343 std::floor(tzero + 0.5),
1344 std::floor(split_t + 0.5));
1346 if (bSeparate || bSplitHere)
1348 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1355 write_tng_frame(trxout, &frout);
1356 // TODO when trjconv behaves better: work how to read and write lambda
1366 trxout = open_trx(out_file2, filemode);
1368 write_trxframe(trxout, &frout, gc);
1373 // Only add a generator statement if title is empty,
1374 // to avoid multiple generated-by statements from various programs
1375 if (std::strlen(top_title) == 0)
1377 sprintf(top_title, "Generated by trjconv");
1381 sprintf(timestr, " t= %9.5f", frout.time);
1385 std::strcpy(timestr, "");
1389 sprintf(stepstr, " step= %" PRId64, frout.step);
1393 std::strcpy(stepstr, "");
1395 title = gmx::formatString("%s%s%s", top_title, timestr, stepstr);
1396 if (bSeparate || bSplitHere)
1398 out = gmx_ffopen(out_file2, "w");
1407 frout.bV ? frout.v : nullptr,
1411 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1412 /* if reading from pdb, we want to keep the original
1413 model numbering else we write the output frame
1414 number plus one, because model 0 is not allowed in pdb */
1415 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1434 const char* outputTitle = "";
1435 if (bSeparate || bTDump)
1437 outputTitle = title.c_str();
1440 frout.bAtoms = TRUE;
1442 frout.atoms = &useatoms;
1443 frout.bStep = FALSE;
1444 frout.bTime = FALSE;
1450 outputTitle = title.c_str();
1452 frout.bAtoms = FALSE;
1456 write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1458 if (bSeparate || bSplitHere)
1464 default: gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1466 if (bSeparate || bSplitHere)
1471 /* execute command */
1475 sprintf(c, "%s %d", exec_command, file_nr - 1);
1476 /*fprintf(stderr,"Executing '%s'\n",c);*/
1479 gmx_fatal(FARGS, "Error executing command: %s", c);
1486 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1487 } while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1490 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1493 "\nWARNING no output, "
1494 "last frame read at t=%g\n",
1497 fprintf(stderr, "\n");
1504 gmx_rmpbc_done(gpbc);
1511 else if (out != nullptr)
1531 do_view(oenv, out_file, nullptr);
1533 output_env_done(oenv);