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",
155 fn, j, t, static_cast<long int>(fpos));
156 if (1 != scanf("%s", yesno))
158 gmx_fatal(FARGS, "Error reading user input");
160 if (std::strcmp(yesno, "YES") == 0)
162 fprintf(stderr, "Once again, I'm gonna DO this...\n");
164 if (0 != gmx_truncate(fn, fpos))
166 gmx_fatal(FARGS, "Error truncating file %s", fn);
171 fprintf(stderr, "Ok, I'll forget about it\n");
176 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
182 /*! \brief Read a full molecular topology if useful and available.
184 * If the input trajectory file is not in TNG format, and the output
185 * file is in TNG format, then we want to try to read a full topology
186 * (if available), so that we can write molecule information to the
187 * output file. The full topology provides better molecule information
188 * than is available from the normal t_topology data used by GROMACS
191 * Also, the t_topology is only read under (different) particular
192 * conditions. If both apply, then a .tpr file might be read
193 * twice. Trying to fix this redundancy while trjconv is still an
194 * all-purpose tool does not seem worthwhile.
196 * Because of the way gmx_prepare_tng_writing is implemented, the case
197 * where the input TNG file has no molecule information will never
198 * lead to an output TNG file having molecule information. Since
199 * molecule information will generally be present if the input TNG
200 * file was written by a GROMACS tool, this seems like reasonable
202 static std::unique_ptr<gmx_mtop_t> read_mtop_for_tng(const char* tps_file,
203 const char* input_file,
204 const char* output_file)
206 std::unique_ptr<gmx_mtop_t> mtop;
208 if (fn2bTPX(tps_file) && efTNG != fn2ftp(input_file) && efTNG == fn2ftp(output_file))
210 int temp_natoms = -1;
211 mtop = std::make_unique<gmx_mtop_t>();
212 read_tpx(tps_file, nullptr, nullptr, &temp_natoms, nullptr, nullptr, mtop.get());
218 int gmx_trjconv(int argc, char* argv[])
220 const char* desc[] = {
221 "[THISMODULE] can convert trajectory files in many ways:",
223 "* from one format to another",
224 "* select a subset of atoms",
225 "* change the periodicity representation",
226 "* keep multimeric molecules together",
227 "* center atoms in the box",
228 "* fit atoms to reference structure",
229 "* reduce the number of frames",
230 "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
231 "* select frames within a certain range of a quantity given",
232 " in an [REF].xvg[ref] file.",
234 "The option to write subtrajectories (-sub) based on the information obtained from",
235 "cluster analysis has been removed from [THISMODULE] and is now part of",
236 "[gmx extract-cluster]",
238 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
241 "The following formats are supported for input and output:",
242 "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
243 "and [REF].pdb[ref].",
244 "The file formats are detected from the file extension.",
245 "The precision of the [REF].xtc[ref] output is taken from the",
246 "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
247 "and from the [TT]-ndec[tt] option for other input formats. The precision",
248 "is always taken from [TT]-ndec[tt], when this option is set.",
249 "All other formats have fixed precision. [REF].trr[ref]",
250 "output can be single or double precision, depending on the precision",
251 "of the [THISMODULE] binary.",
252 "Note that velocities are only supported in",
253 "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
255 "Option [TT]-sep[tt] can be used to write every frame to a separate",
256 "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to ",
258 "[REF].pdb[ref] files with all frames concatenated can be viewed with",
259 "[TT]rasmol -nmrpdb[tt].[PAR]",
261 "It is possible to select part of your trajectory and write it out",
262 "to a new trajectory file in order to save disk space, e.g. for leaving",
263 "out the water from a trajectory of a protein in water.",
264 "[BB]ALWAYS[bb] put the original trajectory on tape!",
265 "We recommend to use the portable [REF].xtc[ref] format for your analysis",
266 "to save disk space and to have portable files.[PAR]",
268 "There are two options for fitting the trajectory to a reference",
269 "either for essential dynamics analysis, etc.",
270 "The first option is just plain fitting to a reference structure",
271 "in the structure file. The second option is a progressive fit",
272 "in which the first timeframe is fitted to the reference structure ",
273 "in the structure file to obtain and each subsequent timeframe is ",
274 "fitted to the previously fitted structure. This way a continuous",
275 "trajectory is generated, which might not be the case when using the",
276 "regular fit method, e.g. when your protein undergoes large",
277 "conformational transitions.[PAR]",
279 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
282 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
283 " and requires a run input file to be supplied with [TT]-s[tt].",
284 " * [TT]res[tt] puts the center of mass of residues in the box.",
285 " * [TT]atom[tt] puts all the atoms in the box.",
286 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
287 " them back. This has the effect that all molecules",
288 " will remain whole (provided they were whole in the initial",
289 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
290 " molecules may diffuse out of the box. The starting configuration",
291 " for this procedure is taken from the structure file, if one is",
292 " supplied, otherwise it is the first frame.",
293 " * [TT]cluster[tt] clusters all the atoms in the selected index",
294 " such that they are all closest to the center of mass of the cluster,",
295 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
296 " results if you in fact have a cluster. Luckily that can be checked",
297 " afterwards using a trajectory viewer. Note also that if your molecules",
298 " are broken this will not work either.",
299 " * [TT]whole[tt] only makes broken molecules whole.",
302 "Option [TT]-ur[tt] sets the unit cell representation for options",
303 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
304 "All three options give different results for triclinic boxes and",
305 "identical results for rectangular boxes.",
306 "[TT]rect[tt] is the ordinary brick shape.",
307 "[TT]tric[tt] is the triclinic unit cell.",
308 "[TT]compact[tt] puts all atoms at the closest distance from the center",
309 "of the box. This can be useful for visualizing e.g. truncated octahedra",
310 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
311 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
312 "is set differently.[PAR]",
314 "Option [TT]-center[tt] centers the system in the box. The user can",
315 "select the group which is used to determine the geometrical center.",
316 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
317 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
318 "[TT]tric[tt]: half of the sum of the box vectors,",
319 "[TT]rect[tt]: half of the box diagonal,",
320 "[TT]zero[tt]: zero.",
321 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
322 "want all molecules in the box after the centering.[PAR]",
324 "Option [TT]-box[tt] sets the size of the new box. This option only works",
325 "for leading dimensions and is thus generally only useful for rectangular boxes.",
326 "If you want to modify only some of the dimensions, e.g. when reading from",
327 "a trajectory, you can use -1 for those dimensions that should stay the same",
329 "It is not always possible to use combinations of [TT]-pbc[tt],",
330 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
331 "you want in one call to [THISMODULE]. Consider using multiple",
332 "calls, and check out the GROMACS website for suggestions.[PAR]",
334 "With [TT]-dt[tt], it is possible to reduce the number of ",
335 "frames in the output. This option relies on the accuracy of the times",
336 "in your input trajectory, so if these are inaccurate use the",
337 "[TT]-timestep[tt] option to modify the time (this can be done",
338 "simultaneously). For making smooth movies, the program [gmx-filter]",
339 "can reduce the number of frames while using low-pass frequency",
340 "filtering, this reduces aliasing of high frequency motions.[PAR]",
342 "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
343 "without copying the file. This is useful when a run has crashed",
344 "during disk I/O (i.e. full disk), or when two contiguous",
345 "trajectories must be concatenated without having double frames.[PAR]",
347 "Option [TT]-dump[tt] can be used to extract a frame at or near",
348 "one specific time from your trajectory, but only works reliably",
349 "if the time interval between frames is uniform.[PAR]",
351 "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
352 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
353 "frames with a value below and above the value of the respective options",
354 "will not be written."
370 const char* pbc_opt[epNR + 1] = { nullptr, "none", "mol", "res", "atom",
371 "nojump", "cluster", "whole", nullptr };
374 const char* unitcell_opt[euNR + 1] = { nullptr, "rect", "tric", "compact", nullptr };
384 const char* center_opt[ecNR + 1] = { nullptr, "tric", "rect", "zero", nullptr };
399 const char* fit[efNR + 1] = { nullptr, "none", "rot+trans", "rotxy+transxy",
400 "translation", "transxy", "progressive", nullptr };
402 static gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
403 static gmx_bool bCenter = FALSE;
404 static int skip_nr = 1, ndec = 3, nzero = 0;
405 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
406 static rvec newbox = { 0, 0, 0 }, shift = { 0, 0, 0 }, trans = { 0, 0, 0 };
407 static char* exec_command = nullptr;
408 static real dropunder = 0, dropover = 0;
409 static gmx_bool bRound = FALSE;
412 { "-skip", FALSE, etINT, { &skip_nr }, "Only write every nr-th frame" },
413 { "-dt", FALSE, etTIME, { &delta_t }, "Only write frame when t MOD dt = first time (%t)" },
414 { "-round", FALSE, etBOOL, { &bRound }, "Round measurements to nearest picosecond" },
415 { "-dump", FALSE, etTIME, { &tdump }, "Dump frame nearest specified time (%t)" },
416 { "-t0", FALSE, etTIME, { &tzero }, "Starting time (%t) (default: don't change)" },
417 { "-timestep", FALSE, etTIME, { ×tep }, "Change time step between input frames (%t)" },
418 { "-pbc", FALSE, etENUM, { pbc_opt }, "PBC treatment (see help text for full description)" },
419 { "-ur", FALSE, etENUM, { unitcell_opt }, "Unit-cell representation" },
420 { "-center", FALSE, etBOOL, { &bCenter }, "Center atoms in box" },
421 { "-boxcenter", FALSE, etENUM, { center_opt }, "Center for -pbc and -center" },
422 { "-box", FALSE, etRVEC, { newbox }, "Size for new cubic box (default: read from input)" },
427 "All coordinates will be translated by trans. This "
428 "can advantageously be combined with -pbc mol -ur "
430 { "-shift", FALSE, etRVEC, { shift }, "All coordinates will be shifted by framenr*shift" },
431 { "-fit", FALSE, etENUM, { fit }, "Fit molecule to ref structure in the structure file" },
432 { "-ndec", FALSE, etINT, { &ndec }, "Number of decimal places to write to .xtc output" },
433 { "-vel", FALSE, etBOOL, { &bVels }, "Read and write velocities if possible" },
434 { "-force", FALSE, etBOOL, { &bForce }, "Read and write forces if possible" },
439 "Truncate input trajectory file after this time (%t)" },
444 "Execute command for every output frame with the "
445 "frame number as argument" },
450 "Start writing new file when t MOD split = first "
456 "Write each frame to a separate .gro, .g96 or .pdb "
462 "If the -sep flag is set, use these many digits "
463 "for the file numbers and prepend zeros as needed" },
464 { "-dropunder", FALSE, etREAL, { &dropunder }, "Drop all frames below this value" },
465 { "-dropover", FALSE, etREAL, { &dropover }, "Drop all frames above this value" },
470 "Add CONECT PDB records when writing [REF].pdb[ref] files. Useful "
471 "for visualization of non-standard molecules, e.g. "
472 "coarse grained ones" }
474 #define NPA asize(pa)
477 t_trxstatus* trxout = nullptr;
480 t_trxframe fr, frout;
482 rvec * xmem = nullptr, *vmem = nullptr, *fmem = nullptr;
483 rvec * xp = nullptr, x_shift, hbox;
484 real* w_rls = nullptr;
485 int m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
487 t_topology* top = nullptr;
488 gmx_conect gc = nullptr;
489 PbcType pbcType = PbcType::Unset;
490 t_atoms * atoms = nullptr, useatoms;
492 int * index = nullptr, *cindex = nullptr;
493 char* grpnm = nullptr;
494 int * frindex, nrfri;
499 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
501 real tshift = 0, dt = -1, prec;
502 gmx_bool bFit, bPFit, bReset;
504 gmx_rmpbc_t gpbc = nullptr;
505 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
506 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
507 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetXtcPrec, bNeedPrec;
508 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
509 gmx_bool bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
510 gmx_bool bWriteFrame, bSplitHere;
511 const char *top_file, *in_file, *out_file = nullptr;
512 char out_file2[256], *charpt;
513 char* outf_base = nullptr;
514 const char* outf_ext = nullptr;
515 char top_title[256], timestr[32], stepstr[32], filemode[5];
516 gmx_output_env_t* oenv;
518 t_filenm fnm[] = { { efTRX, "-f", nullptr, ffREAD }, { efTRO, "-o", nullptr, ffWRITE },
519 { efTPS, nullptr, nullptr, ffOPTRD }, { efNDX, nullptr, nullptr, ffOPTRD },
520 { efNDX, "-fr", "frames", ffOPTRD }, { efNDX, "-sub", "cluster", ffOPTRD },
521 { efXVG, "-drop", "drop", ffOPTRD } };
522 #define NFILE asize(fnm)
524 if (!parse_common_args(&argc, argv, PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW | PCA_TIME_UNIT,
525 NFILE, fnm, NPA, pa, asize(desc), desc, 0, nullptr, &oenv))
530 "Note that major changes are planned in future for "
531 "trjconv, to improve usability and utility.\n");
533 top_file = ftp2fn(efTPS, NFILE, fnm);
535 /* Check command line */
536 in_file = opt2fn("-f", NFILE, fnm);
540 do_trunc(in_file, ttrunc);
544 /* mark active cmdline options */
545 bSetBox = opt2parg_bSet("-box", NPA, pa);
546 bSetTime = opt2parg_bSet("-t0", NPA, pa);
547 bSetXtcPrec = opt2parg_bSet("-ndec", NPA, pa);
548 bSetUR = opt2parg_bSet("-ur", NPA, pa);
549 bExec = opt2parg_bSet("-exec", NPA, pa);
550 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
551 bTDump = opt2parg_bSet("-dump", NPA, pa);
552 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
553 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
554 bTrans = opt2parg_bSet("-trans", NPA, pa);
555 bSplit = (split_t != 0);
557 /* parse enum options */
558 fit_enum = nenum(fit);
559 bFit = (fit_enum == efFit || fit_enum == efFitXY);
560 bReset = (fit_enum == efReset || fit_enum == efResetXY);
561 bPFit = fit_enum == efPFit;
562 pbc_enum = nenum(pbc_opt);
563 bPBCWhole = pbc_enum == epWhole;
564 bPBCcomRes = pbc_enum == epComRes;
565 bPBCcomMol = pbc_enum == epComMol;
566 bPBCcomAtom = pbc_enum == epComAtom;
567 bNoJump = pbc_enum == epNojump;
568 bCluster = pbc_enum == epCluster;
569 bPBC = pbc_enum != epNone;
570 unitcell_enum = nenum(unitcell_opt);
571 ecenter = nenum(center_opt) - ecTric;
573 /* set and check option dependencies */
576 bFit = TRUE; /* for pfit, fit *must* be set */
580 bReset = TRUE; /* for fit, reset *must* be set */
585 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
587 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
591 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
594 "WARNING: Option for unitcell representation (-ur %s)\n"
595 " only has effect in combination with -pbc %s, %s or %s.\n"
596 " Ignoring unitcell representation.\n\n",
597 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
603 "PBC condition treatment does not work together with rotational fit.\n"
604 "Please do the PBC condition treatment first and then run trjconv in a "
606 "for the rotational fit.\n"
607 "First doing the rotational fit and then doing the PBC treatment gives "
612 /* ndec for XTC writing is in nr of decimal places, prec is a multiplication factor: */
614 for (i = 0; i < ndec; i++)
619 bIndex = ftp2bSet(efNDX, NFILE, fnm);
622 /* Determine output type */
623 out_file = opt2fn("-o", NFILE, fnm);
624 int ftp = fn2ftp(out_file);
625 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
626 bNeedPrec = (ftp == efXTC);
627 int ftpin = fn2ftp(in_file);
630 /* check if velocities are possible in input and output files */
631 bVels = (ftp == efTRR || ftp == efGRO || ftp == efG96 || ftp == efTNG)
632 && (ftpin == efTRR || ftpin == efGRO || ftpin == efG96 || ftpin == efTNG
635 if (bSeparate || bSplit)
637 outf_ext = std::strrchr(out_file, '.');
638 if (outf_ext == nullptr)
640 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
642 outf_base = gmx_strdup(out_file);
643 outf_base[outf_ext - out_file] = '\0';
646 bool bSubTraj = opt2bSet("-sub", NFILE, fnm);
650 "The -sub option has been removed from gmx trjconv and is now part\n"
651 "of gmx extract-cluster and does nothing here\n");
657 gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
660 std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
662 /* Determine whether to read a topology */
663 bTPS = (ftp2bSet(efTPS, NFILE, fnm) || bRmPBC || bReset || bPBCcomMol || bCluster
664 || (ftp == efGRO) || (ftp == efPDB) || bCONECT);
666 /* Determine if when can read index groups */
667 bIndex = (bIndex || bTPS);
672 read_tps_conf(top_file, top, &pbcType, &xp, nullptr, top_box, bReset || bPBCcomRes);
673 std::strncpy(top_title, *top->name, 255);
674 top_title[255] = '\0';
677 if (0 == top->mols.nr && (bCluster || bPBCcomMol))
679 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option",
683 /* top_title is only used for gro and pdb,
684 * the header in such a file is top_title, followed by
685 * t= ... and/or step= ...
686 * to prevent double t= or step=, remove it from top_title.
687 * From GROMACS-2018 we only write t/step when the frame actually
688 * has a valid time/step, so we need to check for both separately.
690 if ((charpt = std::strstr(top_title, " t= ")))
694 if ((charpt = std::strstr(top_title, " step= ")))
701 gc = gmx_conect_generate(top);
705 gpbc = gmx_rmpbc_init(&top->idef, pbcType, top->atoms.nr);
709 /* get frame number index */
711 if (opt2bSet("-fr", NFILE, fnm))
713 printf("Select groups of frame number indices:\n");
714 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
717 for (i = 0; i < nrfri; i++)
719 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
724 /* get index groups etc. */
727 printf("Select group for %s fit\n", bFit ? "least squares" : "translational");
728 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit, &ind_fit, &gn_fit);
734 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
738 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
744 printf("Select group for clustering\n");
745 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit, &ind_fit, &gn_fit);
752 printf("Select group for centering\n");
753 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ncent, &cindex, &grpnm);
755 printf("Select group for output\n");
756 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nout, &index, &grpnm);
760 /* no index file, so read natoms from TRX */
761 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
763 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
769 for (i = 0; i < natoms; i++)
783 snew(w_rls, atoms->nr);
784 for (i = 0; (i < ifit); i++)
786 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
789 /* Restore reference structure and set to origin,
790 store original location (to put structure back) */
793 gmx_rmpbc(gpbc, top->atoms.nr, top_box, xp);
795 copy_rvec(xp[index[0]], x_shift);
796 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
797 rvec_dec(x_shift, xp[index[0]]);
804 if (bDropUnder || bDropOver)
806 /* Read the .xvg file with the drop values */
807 fprintf(stderr, "\nReading drop file ...");
808 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
809 fprintf(stderr, " %d time points\n", ndrop);
810 if (ndrop == 0 || ncol < 2)
812 gmx_fatal(FARGS, "Found no data points in %s", opt2fn("-drop", NFILE, fnm));
818 /* Make atoms struct for output in GRO or PDB files */
819 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
821 /* get memory for stuff to go in .pdb file, and initialize
822 * the pdbinfo structure part if the input has it.
824 init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
825 sfree(useatoms.resinfo);
826 useatoms.resinfo = atoms->resinfo;
827 for (i = 0; (i < nout); i++)
829 useatoms.atomname[i] = atoms->atomname[index[i]];
830 useatoms.atom[i] = atoms->atom[index[i]];
831 if (atoms->havePdbInfo)
833 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
835 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind + 1);
839 /* select what to read */
850 flags = flags | TRX_READ_V;
854 flags = flags | TRX_READ_F;
857 /* open trx file for reading */
858 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
861 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1 / fr.prec);
865 if (bSetXtcPrec || !fr.bPrec)
867 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1 / prec);
871 fprintf(stderr, "Using output precision of %g (nm)\n", 1 / prec);
879 // Determine timestep (assuming constant spacing for now) if we
880 // need to dump frames based on time. This is required so we do not
881 // skip the first frame if that was the one that should have been dumped
882 double firstFrameTime = fr.time;
883 if (read_next_frame(oenv, trxin, &fr))
885 dt = fr.time - firstFrameTime;
890 "Warning: Frame times are not incrementing - will dump first "
894 // Now close and reopen so we are at first frame again
897 // Reopen at first frame (We already know it exists if we got here)
898 read_first_frame(oenv, &trxin, in_file, &fr, flags);
901 setTrxFramePbcType(&fr, pbcType);
906 tshift = tzero - fr.time;
916 /* check if index is meaningful */
917 for (i = 0; i < nout; i++)
919 if (index[i] >= natoms)
922 "Index[%d] %d is larger than the number of atoms in the\n"
923 "trajectory file (%d). There is a mismatch in the contents\n"
924 "of your -f, -s and/or -n files.",
925 i, index[i] + 1, natoms);
927 bCopy = bCopy || (i != index[i]);
931 /* open output for writing */
932 std::strcpy(filemode, "w");
936 trxout = trjtools_gmx_prepare_tng_writing(
937 out_file, filemode[0], trxin, nullptr, nout, mtop.get(),
938 gmx::arrayRefFromArray(index, nout), grpnm);
945 trxout = open_trx(out_file, filemode);
951 if (!bSeparate && !bSplit)
953 out = gmx_ffopen(out_file, filemode);
956 default: gmx_incons("Illegal output file format");
972 /* Start the big loop over frames */
978 /* Main loop over frames */
990 /* generate new box */
995 for (m = 0; m < DIM; m++)
999 fr.box[m][m] = newbox[m];
1005 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1013 for (i = 0; i < natoms; i++)
1015 rvec_inc(fr.x[i], trans);
1021 // If we could not read two frames or times are not incrementing
1022 // we have almost no idea what to do,
1023 // but dump the first frame so output is not broken.
1024 if (dt <= 0 || !bDTset)
1030 // Dump the frame if we are less than half a frame time
1031 // below it. This will also ensure we at least dump a
1032 // somewhat reasonable frame if the spacing is unequal
1033 // and we have overrun the frame time. Once we dump one
1034 // frame based on time we quit, so it does not matter
1035 // that this might be true for all subsequent frames too.
1036 bDumpFrame = (fr.time > tdump - 0.5 * dt);
1044 /* determine if an atom jumped across the box and reset it if so */
1045 if (bNoJump && (bTPS || frame != 0))
1047 for (d = 0; d < DIM; d++)
1049 hbox[d] = 0.5 * fr.box[d][d];
1051 for (i = 0; i < natoms; i++)
1055 rvec_dec(fr.x[i], x_shift);
1057 for (m = DIM - 1; m >= 0; m--)
1061 while (fr.x[i][m] - xp[i][m] <= -hbox[m])
1063 for (d = 0; d <= m; d++)
1065 fr.x[i][d] += fr.box[m][d];
1068 while (fr.x[i][m] - xp[i][m] > hbox[m])
1070 for (d = 0; d <= m; d++)
1072 fr.x[i][d] -= fr.box[m][d];
1081 calc_pbc_cluster(ecenter, ifit, top, pbcType, fr.x, ind_fit, fr.box);
1086 /* Now modify the coords according to the flags,
1087 for normal fit, this is only done for output frames */
1090 gmx_rmpbc_trxfr(gpbc, &fr);
1093 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1094 do_fit(natoms, w_rls, xp, fr.x);
1097 /* store this set of coordinates for future use */
1098 if (bPFit || bNoJump)
1104 for (i = 0; (i < natoms); i++)
1106 copy_rvec(fr.x[i], xp[i]);
1107 rvec_inc(fr.x[i], x_shift);
1113 /* see if we have a frame from the frame index group */
1114 for (i = 0; i < nrfri && !bDumpFrame; i++)
1116 bDumpFrame = frame == frindex[i];
1119 if (debug && bDumpFrame)
1121 fprintf(debug, "dumping %d\n", frame);
1124 bWriteFrame = ((!bTDump && (frindex == nullptr) && frame % skip_nr == 0) || bDumpFrame);
1126 if (bWriteFrame && (bDropUnder || bDropOver))
1128 while (dropval[0][drop1] < fr.time && drop1 + 1 < ndrop)
1133 if (std::abs(dropval[0][drop0] - fr.time) < std::abs(dropval[0][drop1] - fr.time))
1141 if ((bDropUnder && dropval[1][dropuse] < dropunder)
1142 || (bDropOver && dropval[1][dropuse] > dropover))
1144 bWriteFrame = FALSE;
1150 /* We should avoid modifying the input frame,
1151 * but since here we don't have the output frame yet,
1152 * we introduce a temporary output frame time variable.
1156 frout_time = fr.time;
1161 frout_time = tzero + frame * timestep;
1165 frout_time += tshift;
1170 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1171 output_env_conv_time(oenv, frout_time),
1172 output_env_get_time_unit(oenv).c_str());
1175 /* check for writing at each delta_t */
1176 bDoIt = (delta_t == 0);
1181 bDoIt = bRmod(frout_time, tzero, delta_t);
1185 /* round() is not C89 compatible, so we do this: */
1186 bDoIt = bRmod(std::floor(frout_time + 0.5), std::floor(tzero + 0.5),
1187 std::floor(delta_t + 0.5));
1191 if (bDoIt || bTDump)
1193 /* print sometimes */
1194 if (((outframe % SKIP) == 0) || (outframe < SKIP))
1196 fprintf(stderr, " -> frame %6d time %8.3f \r", outframe,
1197 output_env_conv_time(oenv, frout_time));
1203 /* Now modify the coords according to the flags,
1204 for PFit we did this already! */
1208 gmx_rmpbc_trxfr(gpbc, &fr);
1213 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1216 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1220 for (i = 0; i < natoms; i++)
1222 rvec_inc(fr.x[i], x_shift);
1229 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1233 auto positionsArrayRef =
1234 gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec*>(fr.x), natoms);
1237 switch (unitcell_enum)
1240 put_atoms_in_box(pbcType, fr.box, positionsArrayRef);
1243 put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1246 put_atoms_in_compact_unitcell(pbcType, ecenter, fr.box,
1253 put_residue_com_in_box(unitcell_enum, ecenter, natoms, atoms->atom,
1254 pbcType, fr.box, fr.x);
1258 put_molecule_com_in_box(unitcell_enum, ecenter, &top->mols, natoms,
1259 atoms->atom, pbcType, fr.box, fr.x);
1261 /* Copy the input trxframe struct to the output trxframe struct */
1263 frout.time = frout_time;
1264 frout.bV = (frout.bV && bVels);
1265 frout.bF = (frout.bF && bForce);
1266 frout.natoms = nout;
1267 if (bNeedPrec && (bSetXtcPrec || !fr.bPrec))
1283 for (i = 0; i < nout; i++)
1285 copy_rvec(fr.x[index[i]], frout.x[i]);
1288 copy_rvec(fr.v[index[i]], frout.v[i]);
1292 copy_rvec(fr.f[index[i]], frout.f[i]);
1297 if (opt2parg_bSet("-shift", NPA, pa))
1299 for (i = 0; i < nout; i++)
1301 for (d = 0; d < DIM; d++)
1303 frout.x[i][d] += outframe * shift[d];
1310 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1314 /* round() is not C89 compatible, so we do this: */
1316 && bRmod(std::floor(frout.time + 0.5),
1317 std::floor(tzero + 0.5), std::floor(split_t + 0.5));
1319 if (bSeparate || bSplitHere)
1321 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1328 write_tng_frame(trxout, &frout);
1329 // TODO when trjconv behaves better: work how to read and write lambda
1339 trxout = open_trx(out_file2, filemode);
1341 write_trxframe(trxout, &frout, gc);
1346 // Only add a generator statement if title is empty,
1347 // to avoid multiple generated-by statements from various programs
1348 if (std::strlen(top_title) == 0)
1350 sprintf(top_title, "Generated by trjconv");
1354 sprintf(timestr, " t= %9.5f", frout.time);
1358 std::strcpy(timestr, "");
1362 sprintf(stepstr, " step= %" PRId64, frout.step);
1366 std::strcpy(stepstr, "");
1368 title = gmx::formatString("%s%s%s", top_title, timestr, stepstr);
1369 if (bSeparate || bSplitHere)
1371 out = gmx_ffopen(out_file2, "w");
1376 write_hconf_p(out, title.c_str(), &useatoms, frout.x,
1377 frout.bV ? frout.v : nullptr, frout.box);
1380 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1381 /* if reading from pdb, we want to keep the original
1382 model numbering else we write the output frame
1383 number plus one, because model 0 is not allowed in pdb */
1384 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1392 write_pdbfile(out, title.c_str(), &useatoms, frout.x,
1393 frout.pbcType, frout.box, ' ', model_nr, gc);
1396 const char* outputTitle = "";
1397 if (bSeparate || bTDump)
1399 outputTitle = title.c_str();
1402 frout.bAtoms = TRUE;
1404 frout.atoms = &useatoms;
1405 frout.bStep = FALSE;
1406 frout.bTime = FALSE;
1412 outputTitle = title.c_str();
1414 frout.bAtoms = FALSE;
1418 write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1420 if (bSeparate || bSplitHere)
1426 default: gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1428 if (bSeparate || bSplitHere)
1433 /* execute command */
1437 sprintf(c, "%s %d", exec_command, file_nr - 1);
1438 /*fprintf(stderr,"Executing '%s'\n",c);*/
1441 gmx_fatal(FARGS, "Error executing command: %s", c);
1448 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1449 } while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1452 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1455 "\nWARNING no output, "
1456 "last frame read at t=%g\n",
1459 fprintf(stderr, "\n");
1466 gmx_rmpbc_done(gpbc);
1473 else if (out != nullptr)
1493 do_view(oenv, out_file, nullptr);
1495 output_env_done(oenv);