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 "[REF].pdb[ref] and [REF].tng[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].tng[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. When writing [REF].tng[ref]",
270 "output the file will contain one molecule type of the correct count",
271 "if the selection name matches the molecule name and the selected atoms",
272 "match all atoms of that molecule. Otherwise the whole selection will",
273 "be treated as one single molecule containing all the selected atoms.[PAR]",
275 "There are two options for fitting the trajectory to a reference",
276 "either for essential dynamics analysis, etc.",
277 "The first option is just plain fitting to a reference structure",
278 "in the structure file. The second option is a progressive fit",
279 "in which the first timeframe is fitted to the reference structure ",
280 "in the structure file to obtain and each subsequent timeframe is ",
281 "fitted to the previously fitted structure. This way a continuous",
282 "trajectory is generated, which might not be the case when using the",
283 "regular fit method, e.g. when your protein undergoes large",
284 "conformational transitions.[PAR]",
286 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
289 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
290 " and requires a run input file to be supplied with [TT]-s[tt].",
291 " * [TT]res[tt] puts the center of mass of residues in the box.",
292 " * [TT]atom[tt] puts all the atoms in the box.",
293 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
294 " them back. This has the effect that all molecules",
295 " will remain whole (provided they were whole in the initial",
296 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
297 " molecules may diffuse out of the box. The starting configuration",
298 " for this procedure is taken from the structure file, if one is",
299 " supplied, otherwise it is the first frame.",
300 " * [TT]cluster[tt] clusters all the atoms in the selected index",
301 " such that they are all closest to the center of mass of the cluster,",
302 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
303 " results if you in fact have a cluster. Luckily that can be checked",
304 " afterwards using a trajectory viewer. Note also that if your molecules",
305 " are broken this will not work either.",
306 " * [TT]whole[tt] only makes broken molecules whole.",
309 "Option [TT]-ur[tt] sets the unit cell representation for options",
310 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
311 "All three options give different results for triclinic boxes and",
312 "identical results for rectangular boxes.",
313 "[TT]rect[tt] is the ordinary brick shape.",
314 "[TT]tric[tt] is the triclinic unit cell.",
315 "[TT]compact[tt] puts all atoms at the closest distance from the center",
316 "of the box. This can be useful for visualizing e.g. truncated octahedra",
317 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
318 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
319 "is set differently.[PAR]",
321 "Option [TT]-center[tt] centers the system in the box. The user can",
322 "select the group which is used to determine the geometrical center.",
323 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
324 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
325 "[TT]tric[tt]: half of the sum of the box vectors,",
326 "[TT]rect[tt]: half of the box diagonal,",
327 "[TT]zero[tt]: zero.",
328 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
329 "want all molecules in the box after the centering.[PAR]",
331 "Option [TT]-box[tt] sets the size of the new box. This option only works",
332 "for leading dimensions and is thus generally only useful for rectangular boxes.",
333 "If you want to modify only some of the dimensions, e.g. when reading from",
334 "a trajectory, you can use -1 for those dimensions that should stay the same",
336 "It is not always possible to use combinations of [TT]-pbc[tt],",
337 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
338 "you want in one call to [THISMODULE]. Consider using multiple",
339 "calls, and check out the GROMACS website for suggestions.[PAR]",
341 "With [TT]-dt[tt], it is possible to reduce the number of ",
342 "frames in the output. This option relies on the accuracy of the times",
343 "in your input trajectory, so if these are inaccurate use the",
344 "[TT]-timestep[tt] option to modify the time (this can be done",
345 "simultaneously). For making smooth movies, the program [gmx-filter]",
346 "can reduce the number of frames while using low-pass frequency",
347 "filtering, this reduces aliasing of high frequency motions.[PAR]",
349 "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
350 "without copying the file. This is useful when a run has crashed",
351 "during disk I/O (i.e. full disk), or when two contiguous",
352 "trajectories must be concatenated without having double frames.[PAR]",
354 "Option [TT]-dump[tt] can be used to extract a frame at or near",
355 "one specific time from your trajectory, but only works reliably",
356 "if the time interval between frames is uniform.[PAR]",
358 "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
359 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
360 "frames with a value below and above the value of the respective options",
361 "will not be written."
377 const char* pbc_opt[epNR + 1] = { nullptr, "none", "mol", "res", "atom",
378 "nojump", "cluster", "whole", nullptr };
381 const char* unitcell_opt[euNR + 1] = { nullptr, "rect", "tric", "compact", nullptr };
391 const char* center_opt[ecNR + 1] = { nullptr, "tric", "rect", "zero", nullptr };
406 const char* fit[efNR + 1] = { nullptr, "none", "rot+trans", "rotxy+transxy",
407 "translation", "transxy", "progressive", nullptr };
409 gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
410 gmx_bool bCenter = FALSE;
411 int skip_nr = 1, ndec = 3, nzero = 0;
412 real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
413 rvec newbox = { 0, 0, 0 }, shift = { 0, 0, 0 }, trans = { 0, 0, 0 };
414 char* exec_command = nullptr;
415 real dropunder = 0, dropover = 0;
416 gmx_bool bRound = FALSE;
419 { "-skip", FALSE, etINT, { &skip_nr }, "Only write every nr-th frame" },
420 { "-dt", FALSE, etTIME, { &delta_t }, "Only write frame when t MOD dt = first time (%t)" },
421 { "-round", FALSE, etBOOL, { &bRound }, "Round measurements to nearest picosecond" },
422 { "-dump", FALSE, etTIME, { &tdump }, "Dump frame nearest specified time (%t)" },
423 { "-t0", FALSE, etTIME, { &tzero }, "Starting time (%t) (default: don't change)" },
424 { "-timestep", FALSE, etTIME, { ×tep }, "Change time step between input frames (%t)" },
425 { "-pbc", FALSE, etENUM, { pbc_opt }, "PBC treatment (see help text for full description)" },
426 { "-ur", FALSE, etENUM, { unitcell_opt }, "Unit-cell representation" },
427 { "-center", FALSE, etBOOL, { &bCenter }, "Center atoms in box" },
428 { "-boxcenter", FALSE, etENUM, { center_opt }, "Center for -pbc and -center" },
429 { "-box", FALSE, etRVEC, { newbox }, "Size for new cubic box (default: read from input)" },
434 "All coordinates will be translated by trans. This "
435 "can advantageously be combined with -pbc mol -ur "
437 { "-shift", FALSE, etRVEC, { shift }, "All coordinates will be shifted by framenr*shift" },
438 { "-fit", FALSE, etENUM, { fit }, "Fit molecule to ref structure in the structure file" },
439 { "-ndec", FALSE, etINT, { &ndec }, "Number of decimal places to write to .xtc output" },
440 { "-vel", FALSE, etBOOL, { &bVels }, "Read and write velocities if possible" },
441 { "-force", FALSE, etBOOL, { &bForce }, "Read and write forces if possible" },
446 "Truncate input trajectory file after this time (%t)" },
451 "Execute command for every output frame with the "
452 "frame number as argument" },
457 "Start writing new file when t MOD split = first "
463 "Write each frame to a separate .gro, .g96 or .pdb "
469 "If the -sep flag is set, use these many digits "
470 "for the file numbers and prepend zeros as needed" },
471 { "-dropunder", FALSE, etREAL, { &dropunder }, "Drop all frames below this value" },
472 { "-dropover", FALSE, etREAL, { &dropover }, "Drop all frames above this value" },
477 "Add CONECT PDB records when writing [REF].pdb[ref] files. Useful "
478 "for visualization of non-standard molecules, e.g. "
479 "coarse grained ones" }
481 #define NPA asize(pa)
484 t_trxstatus* trxout = nullptr;
487 t_trxframe fr, frout;
489 rvec * xmem = nullptr, *vmem = nullptr, *fmem = nullptr;
490 rvec * xp = nullptr, x_shift, hbox;
491 real* w_rls = nullptr;
492 int m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
494 t_topology* top = nullptr;
495 gmx_conect gc = nullptr;
496 PbcType pbcType = PbcType::Unset;
497 t_atoms * atoms = nullptr, useatoms;
499 int * index = nullptr, *cindex = nullptr;
500 char* grpnm = nullptr;
501 int * frindex, nrfri;
506 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
508 real tshift = 0, dt = -1, prec;
509 gmx_bool bFit, bPFit, bReset;
511 gmx_rmpbc_t gpbc = nullptr;
512 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
513 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
514 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetXtcPrec, bNeedPrec;
515 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
516 gmx_bool bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
517 gmx_bool bWriteFrame, bSplitHere;
518 const char *top_file, *in_file, *out_file = nullptr;
519 char out_file2[256], *charpt;
520 char* outf_base = nullptr;
521 const char* outf_ext = nullptr;
522 char top_title[256], timestr[32], stepstr[32], filemode[5];
523 gmx_output_env_t* oenv;
525 t_filenm fnm[] = { { efTRX, "-f", nullptr, ffREAD }, { efTRO, "-o", nullptr, ffWRITE },
526 { efTPS, nullptr, nullptr, ffOPTRD }, { efNDX, nullptr, nullptr, ffOPTRD },
527 { efNDX, "-fr", "frames", ffOPTRD }, { efNDX, "-sub", "cluster", ffOPTRD },
528 { efXVG, "-drop", "drop", ffOPTRD } };
529 #define NFILE asize(fnm)
531 if (!parse_common_args(&argc,
533 PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW | PCA_TIME_UNIT,
547 "Note that major changes are planned in future for "
548 "trjconv, to improve usability and utility.\n");
550 top_file = ftp2fn(efTPS, NFILE, fnm);
552 /* Check command line */
553 in_file = opt2fn("-f", NFILE, fnm);
557 do_trunc(in_file, ttrunc);
561 /* mark active cmdline options */
562 bSetBox = opt2parg_bSet("-box", NPA, pa);
563 bSetTime = opt2parg_bSet("-t0", NPA, pa);
564 bSetXtcPrec = opt2parg_bSet("-ndec", NPA, pa);
565 bSetUR = opt2parg_bSet("-ur", NPA, pa);
566 bExec = opt2parg_bSet("-exec", NPA, pa);
567 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
568 bTDump = opt2parg_bSet("-dump", NPA, pa);
569 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
570 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
571 bTrans = opt2parg_bSet("-trans", NPA, pa);
572 bSplit = (split_t != 0);
574 /* parse enum options */
575 fit_enum = nenum(fit);
576 bFit = (fit_enum == efFit || fit_enum == efFitXY);
577 bReset = (fit_enum == efReset || fit_enum == efResetXY);
578 bPFit = fit_enum == efPFit;
579 pbc_enum = nenum(pbc_opt);
580 bPBCWhole = pbc_enum == epWhole;
581 bPBCcomRes = pbc_enum == epComRes;
582 bPBCcomMol = pbc_enum == epComMol;
583 bPBCcomAtom = pbc_enum == epComAtom;
584 bNoJump = pbc_enum == epNojump;
585 bCluster = pbc_enum == epCluster;
586 bPBC = pbc_enum != epNone;
587 unitcell_enum = nenum(unitcell_opt);
588 ecenter = nenum(center_opt) - ecTric;
590 /* set and check option dependencies */
593 bFit = TRUE; /* for pfit, fit *must* be set */
597 bReset = TRUE; /* for fit, reset *must* be set */
602 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
604 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
608 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
611 "WARNING: Option for unitcell representation (-ur %s)\n"
612 " only has effect in combination with -pbc %s, %s or %s.\n"
613 " Ingoring unitcell representation.\n\n",
623 "PBC condition treatment does not work together with rotational fit.\n"
624 "Please do the PBC condition treatment first and then run trjconv in a "
626 "for the rotational fit.\n"
627 "First doing the rotational fit and then doing the PBC treatment gives "
632 /* ndec for XTC writing is in nr of decimal places, prec is a multiplication factor: */
634 for (i = 0; i < ndec; i++)
639 bIndex = ftp2bSet(efNDX, NFILE, fnm);
642 /* Determine output type */
643 out_file = opt2fn("-o", NFILE, fnm);
644 int ftp = fn2ftp(out_file);
645 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
646 bNeedPrec = (ftp == efXTC);
647 int ftpin = fn2ftp(in_file);
650 /* check if velocities are possible in input and output files */
651 bVels = (ftp == efTRR || ftp == efGRO || ftp == efG96 || ftp == efTNG)
652 && (ftpin == efTRR || ftpin == efGRO || ftpin == efG96 || ftpin == efTNG
655 if (bSeparate || bSplit)
657 outf_ext = std::strrchr(out_file, '.');
658 if (outf_ext == nullptr)
660 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
662 outf_base = gmx_strdup(out_file);
663 outf_base[outf_ext - out_file] = '\0';
666 bool bSubTraj = opt2bSet("-sub", NFILE, fnm);
670 "The -sub option has been removed from gmx trjconv and is now part\n"
671 "of gmx extract-cluster and does nothing here\n");
677 gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
680 std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
682 /* Determine whether to read a topology */
683 bTPS = (ftp2bSet(efTPS, NFILE, fnm) || bRmPBC || bReset || bPBCcomMol || bCluster
684 || (ftp == efGRO) || (ftp == efPDB) || bCONECT);
686 /* Determine if when can read index groups */
687 bIndex = (bIndex || bTPS);
692 read_tps_conf(top_file, top, &pbcType, &xp, nullptr, top_box, bReset || bPBCcomRes);
693 std::strncpy(top_title, *top->name, 255);
694 top_title[255] = '\0';
697 if (0 == top->mols.nr && (bCluster || bPBCcomMol))
699 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
702 /* top_title is only used for gro and pdb,
703 * the header in such a file is top_title, followed by
704 * t= ... and/or step= ...
705 * to prevent double t= or step=, remove it from top_title.
706 * From GROMACS-2018 we only write t/step when the frame actually
707 * has a valid time/step, so we need to check for both separately.
709 if ((charpt = std::strstr(top_title, " t= ")))
713 if ((charpt = std::strstr(top_title, " step= ")))
720 gc = gmx_conect_generate(top);
724 gpbc = gmx_rmpbc_init(&top->idef, pbcType, top->atoms.nr);
728 /* get frame number index */
730 if (opt2bSet("-fr", NFILE, fnm))
732 printf("Select groups of frame number indices:\n");
733 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
736 for (i = 0; i < nrfri; i++)
738 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
743 /* get index groups etc. */
746 printf("Select group for %s fit\n", bFit ? "least squares" : "translational");
747 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit, &ind_fit, &gn_fit);
753 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
757 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
763 printf("Select group for clustering\n");
764 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit, &ind_fit, &gn_fit);
771 printf("Select group for centering\n");
772 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ncent, &cindex, &grpnm);
774 printf("Select group for output\n");
775 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nout, &index, &grpnm);
779 /* no index file, so read natoms from TRX */
780 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
782 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
788 for (i = 0; i < natoms; i++)
802 snew(w_rls, atoms->nr);
803 for (i = 0; (i < ifit); i++)
805 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
808 /* Restore reference structure and set to origin,
809 store original location (to put structure back) */
812 gmx_rmpbc(gpbc, top->atoms.nr, top_box, xp);
814 copy_rvec(xp[index[0]], x_shift);
815 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
816 rvec_dec(x_shift, xp[index[0]]);
823 if (bDropUnder || bDropOver)
825 /* Read the .xvg file with the drop values */
826 fprintf(stderr, "\nReading drop file ...");
827 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
828 fprintf(stderr, " %d time points\n", ndrop);
829 if (ndrop == 0 || ncol < 2)
831 gmx_fatal(FARGS, "Found no data points in %s", opt2fn("-drop", NFILE, fnm));
837 /* Make atoms struct for output in GRO or PDB files */
838 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
840 /* get memory for stuff to go in .pdb file, and initialize
841 * the pdbinfo structure part if the input has it.
843 init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
844 sfree(useatoms.resinfo);
845 useatoms.resinfo = atoms->resinfo;
846 for (i = 0; (i < nout); i++)
848 useatoms.atomname[i] = atoms->atomname[index[i]];
849 useatoms.atom[i] = atoms->atom[index[i]];
850 if (atoms->havePdbInfo)
852 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
854 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind + 1);
858 /* select what to read */
869 flags = flags | TRX_READ_V;
873 flags = flags | TRX_READ_F;
876 /* open trx file for reading */
877 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
880 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1 / fr.prec);
884 if (bSetXtcPrec || !fr.bPrec)
886 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1 / prec);
890 fprintf(stderr, "Using output precision of %g (nm)\n", 1 / prec);
898 // Determine timestep (assuming constant spacing for now) if we
899 // need to dump frames based on time. This is required so we do not
900 // skip the first frame if that was the one that should have been dumped
901 double firstFrameTime = fr.time;
902 if (read_next_frame(oenv, trxin, &fr))
904 dt = fr.time - firstFrameTime;
909 "Warning: Frame times are not incrementing - will dump first "
913 // Now close and reopen so we are at first frame again
916 // Reopen at first frame (We already know it exists if we got here)
917 read_first_frame(oenv, &trxin, in_file, &fr, flags);
920 setTrxFramePbcType(&fr, pbcType);
925 tshift = tzero - fr.time;
935 /* check if index is meaningful */
936 for (i = 0; i < nout; i++)
938 if (index[i] >= natoms)
941 "Index[%d] %d is larger than the number of atoms in the\n"
942 "trajectory file (%d). There is a mismatch in the contents\n"
943 "of your -f, -s and/or -n files.",
948 bCopy = bCopy || (i != index[i]);
952 /* open output for writing */
953 std::strcpy(filemode, "w");
957 trxout = trjtools_gmx_prepare_tng_writing(out_file,
963 gmx::arrayRefFromArray(index, nout),
971 trxout = open_trx(out_file, filemode);
977 if (!bSeparate && !bSplit)
979 out = gmx_ffopen(out_file, filemode);
982 default: gmx_incons("Illegal output file format");
998 /* Start the big loop over frames */
1004 /* Main loop over frames */
1016 /* generate new box */
1021 for (m = 0; m < DIM; m++)
1025 fr.box[m][m] = newbox[m];
1031 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1039 for (i = 0; i < natoms; i++)
1041 rvec_inc(fr.x[i], trans);
1047 // If we could not read two frames or times are not incrementing
1048 // we have almost no idea what to do,
1049 // but dump the first frame so output is not broken.
1050 if (dt <= 0 || !bDTset)
1056 // Dump the frame if we are less than half a frame time
1057 // below it. This will also ensure we at least dump a
1058 // somewhat reasonable frame if the spacing is unequal
1059 // and we have overrun the frame time. Once we dump one
1060 // frame based on time we quit, so it does not matter
1061 // that this might be true for all subsequent frames too.
1062 bDumpFrame = (fr.time > tdump - 0.5 * dt);
1070 /* determine if an atom jumped across the box and reset it if so */
1071 if (bNoJump && (bTPS || frame != 0))
1073 for (d = 0; d < DIM; d++)
1075 hbox[d] = 0.5 * fr.box[d][d];
1077 for (i = 0; i < natoms; i++)
1081 rvec_dec(fr.x[i], x_shift);
1083 for (m = DIM - 1; m >= 0; m--)
1087 while (fr.x[i][m] - xp[i][m] <= -hbox[m])
1089 for (d = 0; d <= m; d++)
1091 fr.x[i][d] += fr.box[m][d];
1094 while (fr.x[i][m] - xp[i][m] > hbox[m])
1096 for (d = 0; d <= m; d++)
1098 fr.x[i][d] -= fr.box[m][d];
1107 calc_pbc_cluster(ecenter, ifit, top, pbcType, fr.x, ind_fit, fr.box);
1112 /* Now modify the coords according to the flags,
1113 for normal fit, this is only done for output frames */
1116 gmx_rmpbc_trxfr(gpbc, &fr);
1119 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1120 do_fit(natoms, w_rls, xp, fr.x);
1123 /* store this set of coordinates for future use */
1124 if (bPFit || bNoJump)
1130 for (i = 0; (i < natoms); i++)
1132 copy_rvec(fr.x[i], xp[i]);
1133 rvec_inc(fr.x[i], x_shift);
1139 /* see if we have a frame from the frame index group */
1140 for (i = 0; i < nrfri && !bDumpFrame; i++)
1142 bDumpFrame = frame == frindex[i];
1145 if (debug && bDumpFrame)
1147 fprintf(debug, "dumping %d\n", frame);
1150 bWriteFrame = ((!bTDump && (frindex == nullptr) && frame % skip_nr == 0) || bDumpFrame);
1152 if (bWriteFrame && (bDropUnder || bDropOver))
1154 while (dropval[0][drop1] < fr.time && drop1 + 1 < ndrop)
1159 if (std::abs(dropval[0][drop0] - fr.time) < std::abs(dropval[0][drop1] - fr.time))
1167 if ((bDropUnder && dropval[1][dropuse] < dropunder)
1168 || (bDropOver && dropval[1][dropuse] > dropover))
1170 bWriteFrame = FALSE;
1176 /* We should avoid modifying the input frame,
1177 * but since here we don't have the output frame yet,
1178 * we introduce a temporary output frame time variable.
1182 frout_time = fr.time;
1187 frout_time = tzero + frame * timestep;
1191 frout_time += tshift;
1197 "\nDumping frame at t= %g %s\n",
1198 output_env_conv_time(oenv, frout_time),
1199 output_env_get_time_unit(oenv).c_str());
1202 /* check for writing at each delta_t */
1203 bDoIt = (delta_t == 0);
1208 bDoIt = bRmod(frout_time, tzero, delta_t);
1212 /* round() is not C89 compatible, so we do this: */
1213 bDoIt = bRmod(std::floor(frout_time + 0.5),
1214 std::floor(tzero + 0.5),
1215 std::floor(delta_t + 0.5));
1219 if (bDoIt || bTDump)
1221 /* print sometimes */
1222 if (((outframe % SKIP) == 0) || (outframe < SKIP))
1225 " -> frame %6d time %8.3f \r",
1227 output_env_conv_time(oenv, frout_time));
1233 /* Now modify the coords according to the flags,
1234 for PFit we did this already! */
1238 gmx_rmpbc_trxfr(gpbc, &fr);
1243 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1246 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1250 for (i = 0; i < natoms; i++)
1252 rvec_inc(fr.x[i], x_shift);
1259 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1263 auto positionsArrayRef =
1264 gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec*>(fr.x), natoms);
1267 switch (unitcell_enum)
1270 put_atoms_in_box(pbcType, fr.box, positionsArrayRef);
1273 put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1276 put_atoms_in_compact_unitcell(
1277 pbcType, ecenter, fr.box, positionsArrayRef);
1283 put_residue_com_in_box(
1284 unitcell_enum, ecenter, natoms, atoms->atom, pbcType, fr.box, fr.x);
1288 put_molecule_com_in_box(
1289 unitcell_enum, ecenter, &top->mols, natoms, atoms->atom, pbcType, fr.box, fr.x);
1291 /* Copy the input trxframe struct to the output trxframe struct */
1293 frout.time = frout_time;
1294 frout.bV = (frout.bV && bVels);
1295 frout.bF = (frout.bF && bForce);
1296 frout.natoms = nout;
1297 if (bNeedPrec && (bSetXtcPrec || !fr.bPrec))
1313 for (i = 0; i < nout; i++)
1315 copy_rvec(fr.x[index[i]], frout.x[i]);
1318 copy_rvec(fr.v[index[i]], frout.v[i]);
1322 copy_rvec(fr.f[index[i]], frout.f[i]);
1327 if (opt2parg_bSet("-shift", NPA, pa))
1329 for (i = 0; i < nout; i++)
1331 for (d = 0; d < DIM; d++)
1333 frout.x[i][d] += outframe * shift[d];
1340 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1344 /* round() is not C89 compatible, so we do this: */
1346 && bRmod(std::floor(frout.time + 0.5),
1347 std::floor(tzero + 0.5),
1348 std::floor(split_t + 0.5));
1350 if (bSeparate || bSplitHere)
1352 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1359 write_tng_frame(trxout, &frout);
1360 // TODO when trjconv behaves better: work how to read and write lambda
1370 trxout = open_trx(out_file2, filemode);
1372 write_trxframe(trxout, &frout, gc);
1377 // Only add a generator statement if title is empty,
1378 // to avoid multiple generated-by statements from various programs
1379 if (std::strlen(top_title) == 0)
1381 sprintf(top_title, "Generated by trjconv");
1385 sprintf(timestr, " t= %9.5f", frout.time);
1389 std::strcpy(timestr, "");
1393 sprintf(stepstr, " step= %" PRId64, frout.step);
1397 std::strcpy(stepstr, "");
1399 title = gmx::formatString("%s%s%s", top_title, timestr, stepstr);
1400 if (bSeparate || bSplitHere)
1402 out = gmx_ffopen(out_file2, "w");
1411 frout.bV ? frout.v : nullptr,
1415 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1416 /* if reading from pdb, we want to keep the original
1417 model numbering else we write the output frame
1418 number plus one, because model 0 is not allowed in pdb */
1419 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1438 const char* outputTitle = "";
1439 if (bSeparate || bTDump)
1441 outputTitle = title.c_str();
1444 frout.bAtoms = TRUE;
1446 frout.atoms = &useatoms;
1447 frout.bStep = FALSE;
1448 frout.bTime = FALSE;
1454 outputTitle = title.c_str();
1456 frout.bAtoms = FALSE;
1460 write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1462 if (bSeparate || bSplitHere)
1468 default: gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1470 if (bSeparate || bSplitHere)
1475 /* execute command */
1479 sprintf(c, "%s %d", exec_command, file_nr - 1);
1480 /*fprintf(stderr,"Executing '%s'\n",c);*/
1483 gmx_fatal(FARGS, "Error executing command: %s", c);
1490 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1491 } while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1494 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1497 "\nWARNING no output, "
1498 "last frame read at t=%g\n",
1501 fprintf(stderr, "\n");
1508 gmx_rmpbc_done(gpbc);
1515 else if (out != nullptr)
1535 do_view(oenv, out_file, nullptr);
1537 output_env_done(oenv);