2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/commandline/viewit.h"
50 #include "gromacs/fileio/confio.h"
51 #include "gromacs/fileio/g96io.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/groio.h"
54 #include "gromacs/fileio/pdbio.h"
55 #include "gromacs/fileio/tngio.h"
56 #include "gromacs/fileio/tpxio.h"
57 #include "gromacs/fileio/trrio.h"
58 #include "gromacs/fileio/trxio.h"
59 #include "gromacs/fileio/xtcio.h"
60 #include "gromacs/fileio/xvgr.h"
61 #include "gromacs/math/do_fit.h"
62 #include "gromacs/math/functions.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/pbcutil/pbcmethods.h"
67 #include "gromacs/pbcutil/rmpbc.h"
68 #include "gromacs/topology/index.h"
69 #include "gromacs/topology/topology.h"
70 #include "gromacs/trajectory/trajectoryframe.h"
71 #include "gromacs/utility/arrayref.h"
72 #include "gromacs/utility/arraysize.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/futil.h"
75 #include "gromacs/utility/smalloc.h"
77 static void mk_filenm(char* base, const char* ext, int ndigit, int file_nr, char out_file[])
82 std::strcpy(out_file, base);
92 std::strncat(out_file, "00000000000", ndigit - nd);
94 sprintf(nbuf, "%d.", file_nr);
95 std::strcat(out_file, nbuf);
96 std::strcat(out_file, ext);
99 static void check_trr(const char* fn)
101 if (fn2ftp(fn) != efTRR)
103 gmx_fatal(FARGS, "%s is not a trajectory file, exiting\n", fn);
107 static void do_trunc(const char* fn, real t0)
120 gmx_fatal(FARGS, "You forgot to set the truncation time");
123 /* Check whether this is a .trr file */
126 in = gmx_trr_open(fn, "r");
127 fp = gmx_fio_getfp(in);
130 fprintf(stderr, "Sorry, can not trunc %s, truncation of this filetype is not supported\n", fn);
136 fpos = gmx_fio_ftell(in);
138 while (!bStop && gmx_trr_read_frame_header(in, &sh, &bOK))
140 gmx_trr_read_frame_data(in, &sh, nullptr, nullptr, nullptr, nullptr);
141 fpos = gmx_ftell(fp);
145 gmx_fseek(fp, fpos, SEEK_SET);
152 "Do you REALLY want to truncate this trajectory (%s) at:\n"
153 "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
154 fn, j, t, static_cast<long int>(fpos));
155 if (1 != scanf("%s", yesno))
157 gmx_fatal(FARGS, "Error reading user input");
159 if (std::strcmp(yesno, "YES") == 0)
161 fprintf(stderr, "Once again, I'm gonna DO this...\n");
163 if (0 != gmx_truncate(fn, fpos))
165 gmx_fatal(FARGS, "Error truncating file %s", fn);
170 fprintf(stderr, "Ok, I'll forget about it\n");
175 fprintf(stderr, "Already at end of file (t=%g)...\n", t);
181 /*! \brief Read a full molecular topology if useful and available.
183 * If the input trajectory file is not in TNG format, and the output
184 * file is in TNG format, then we want to try to read a full topology
185 * (if available), so that we can write molecule information to the
186 * output file. The full topology provides better molecule information
187 * than is available from the normal t_topology data used by GROMACS
190 * Also, the t_topology is only read under (different) particular
191 * conditions. If both apply, then a .tpr file might be read
192 * twice. Trying to fix this redundancy while trjconv is still an
193 * all-purpose tool does not seem worthwhile.
195 * Because of the way gmx_prepare_tng_writing is implemented, the case
196 * where the input TNG file has no molecule information will never
197 * lead to an output TNG file having molecule information. Since
198 * molecule information will generally be present if the input TNG
199 * file was written by a GROMACS tool, this seems like reasonable
201 static std::unique_ptr<gmx_mtop_t> read_mtop_for_tng(const char* tps_file,
202 const char* input_file,
203 const char* output_file)
205 std::unique_ptr<gmx_mtop_t> mtop;
207 if (fn2bTPX(tps_file) && efTNG != fn2ftp(input_file) && efTNG == fn2ftp(output_file))
209 int temp_natoms = -1;
210 mtop = std::make_unique<gmx_mtop_t>();
211 read_tpx(tps_file, nullptr, nullptr, &temp_natoms, nullptr, nullptr, mtop.get());
217 int gmx_trjconv(int argc, char* argv[])
219 const char* desc[] = {
220 "[THISMODULE] can convert trajectory files in many ways:",
222 "* from one format to another",
223 "* select a subset of atoms",
224 "* change the periodicity representation",
225 "* keep multimeric molecules together",
226 "* center atoms in the box",
227 "* fit atoms to reference structure",
228 "* reduce the number of frames",
229 "* change the timestamps of the frames ([TT]-t0[tt] and [TT]-timestep[tt])",
230 "* select frames within a certain range of a quantity given",
231 " in an [REF].xvg[ref] file.",
233 "The option to write subtrajectories (-sub) based on the information obtained from",
234 "cluster analysis has been removed from [THISMODULE] and is now part of",
235 "[gmx extract-cluster]",
237 "[gmx-trjcat] is better suited for concatenating multiple trajectory files.",
240 "The following formats are supported for input and output:",
241 "[REF].xtc[ref], [REF].trr[ref], [REF].gro[ref], [TT].g96[tt]",
242 "and [REF].pdb[ref].",
243 "The file formats are detected from the file extension.",
244 "The precision of the [REF].xtc[ref] output is taken from the",
245 "input file for [REF].xtc[ref], [REF].gro[ref] and [REF].pdb[ref],",
246 "and from the [TT]-ndec[tt] option for other input formats. The precision",
247 "is always taken from [TT]-ndec[tt], when this option is set.",
248 "All other formats have fixed precision. [REF].trr[ref]",
249 "output can be single or double precision, depending on the precision",
250 "of the [THISMODULE] binary.",
251 "Note that velocities are only supported in",
252 "[REF].trr[ref], [REF].gro[ref] and [TT].g96[tt] files.[PAR]",
254 "Option [TT]-sep[tt] can be used to write every frame to a separate",
255 "[TT].gro, .g96[tt] or [REF].pdb[ref] file. By default, all frames all written to ",
257 "[REF].pdb[ref] files with all frames concatenated can be viewed with",
258 "[TT]rasmol -nmrpdb[tt].[PAR]",
260 "It is possible to select part of your trajectory and write it out",
261 "to a new trajectory file in order to save disk space, e.g. for leaving",
262 "out the water from a trajectory of a protein in water.",
263 "[BB]ALWAYS[bb] put the original trajectory on tape!",
264 "We recommend to use the portable [REF].xtc[ref] format for your analysis",
265 "to save disk space and to have portable files.[PAR]",
267 "There are two options for fitting the trajectory to a reference",
268 "either for essential dynamics analysis, etc.",
269 "The first option is just plain fitting to a reference structure",
270 "in the structure file. The second option is a progressive fit",
271 "in which the first timeframe is fitted to the reference structure ",
272 "in the structure file to obtain and each subsequent timeframe is ",
273 "fitted to the previously fitted structure. This way a continuous",
274 "trajectory is generated, which might not be the case when using the",
275 "regular fit method, e.g. when your protein undergoes large",
276 "conformational transitions.[PAR]",
278 "Option [TT]-pbc[tt] sets the type of periodic boundary condition",
281 " * [TT]mol[tt] puts the center of mass of molecules in the box,",
282 " and requires a run input file to be supplied with [TT]-s[tt].",
283 " * [TT]res[tt] puts the center of mass of residues in the box.",
284 " * [TT]atom[tt] puts all the atoms in the box.",
285 " * [TT]nojump[tt] checks if atoms jump across the box and then puts",
286 " them back. This has the effect that all molecules",
287 " will remain whole (provided they were whole in the initial",
288 " conformation). [BB]Note[bb] that this ensures a continuous trajectory but",
289 " molecules may diffuse out of the box. The starting configuration",
290 " for this procedure is taken from the structure file, if one is",
291 " supplied, otherwise it is the first frame.",
292 " * [TT]cluster[tt] clusters all the atoms in the selected index",
293 " such that they are all closest to the center of mass of the cluster,",
294 " which is iteratively updated. [BB]Note[bb] that this will only give meaningful",
295 " results if you in fact have a cluster. Luckily that can be checked",
296 " afterwards using a trajectory viewer. Note also that if your molecules",
297 " are broken this will not work either.",
298 " * [TT]whole[tt] only makes broken molecules whole.",
301 "Option [TT]-ur[tt] sets the unit cell representation for options",
302 "[TT]mol[tt], [TT]res[tt] and [TT]atom[tt] of [TT]-pbc[tt].",
303 "All three options give different results for triclinic boxes and",
304 "identical results for rectangular boxes.",
305 "[TT]rect[tt] is the ordinary brick shape.",
306 "[TT]tric[tt] is the triclinic unit cell.",
307 "[TT]compact[tt] puts all atoms at the closest distance from the center",
308 "of the box. This can be useful for visualizing e.g. truncated octahedra",
309 "or rhombic dodecahedra. The center for options [TT]tric[tt] and [TT]compact[tt]",
310 "is [TT]tric[tt] (see below), unless the option [TT]-boxcenter[tt]",
311 "is set differently.[PAR]",
313 "Option [TT]-center[tt] centers the system in the box. The user can",
314 "select the group which is used to determine the geometrical center.",
315 "Option [TT]-boxcenter[tt] sets the location of the center of the box",
316 "for options [TT]-pbc[tt] and [TT]-center[tt]. The center options are:",
317 "[TT]tric[tt]: half of the sum of the box vectors,",
318 "[TT]rect[tt]: half of the box diagonal,",
319 "[TT]zero[tt]: zero.",
320 "Use option [TT]-pbc mol[tt] in addition to [TT]-center[tt] when you",
321 "want all molecules in the box after the centering.[PAR]",
323 "Option [TT]-box[tt] sets the size of the new box. This option only works",
324 "for leading dimensions and is thus generally only useful for rectangular boxes.",
325 "If you want to modify only some of the dimensions, e.g. when reading from",
326 "a trajectory, you can use -1 for those dimensions that should stay the same",
328 "It is not always possible to use combinations of [TT]-pbc[tt],",
329 "[TT]-fit[tt], [TT]-ur[tt] and [TT]-center[tt] to do exactly what",
330 "you want in one call to [THISMODULE]. Consider using multiple",
331 "calls, and check out the GROMACS website for suggestions.[PAR]",
333 "With [TT]-dt[tt], it is possible to reduce the number of ",
334 "frames in the output. This option relies on the accuracy of the times",
335 "in your input trajectory, so if these are inaccurate use the",
336 "[TT]-timestep[tt] option to modify the time (this can be done",
337 "simultaneously). For making smooth movies, the program [gmx-filter]",
338 "can reduce the number of frames while using low-pass frequency",
339 "filtering, this reduces aliasing of high frequency motions.[PAR]",
341 "Using [TT]-trunc[tt] [THISMODULE] can truncate [REF].trr[ref] in place, i.e.",
342 "without copying the file. This is useful when a run has crashed",
343 "during disk I/O (i.e. full disk), or when two contiguous",
344 "trajectories must be concatenated without having double frames.[PAR]",
346 "Option [TT]-dump[tt] can be used to extract a frame at or near",
347 "one specific time from your trajectory, but only works reliably",
348 "if the time interval between frames is uniform.[PAR]",
350 "Option [TT]-drop[tt] reads an [REF].xvg[ref] file with times and values.",
351 "When options [TT]-dropunder[tt] and/or [TT]-dropover[tt] are set,",
352 "frames with a value below and above the value of the respective options",
353 "will not be written."
369 const char* pbc_opt[epNR + 1] = { nullptr, "none", "mol", "res", "atom",
370 "nojump", "cluster", "whole", nullptr };
373 const char* unitcell_opt[euNR + 1] = { nullptr, "rect", "tric", "compact", nullptr };
383 const char* center_opt[ecNR + 1] = { nullptr, "tric", "rect", "zero", nullptr };
398 const char* fit[efNR + 1] = { nullptr, "none", "rot+trans", "rotxy+transxy",
399 "translation", "transxy", "progressive", nullptr };
401 static gmx_bool bSeparate = FALSE, bVels = TRUE, bForce = FALSE, bCONECT = FALSE;
402 static gmx_bool bCenter = FALSE;
403 static int skip_nr = 1, ndec = 3, nzero = 0;
404 static real tzero = 0, delta_t = 0, timestep = 0, ttrunc = -1, tdump = -1, split_t = 0;
405 static rvec newbox = { 0, 0, 0 }, shift = { 0, 0, 0 }, trans = { 0, 0, 0 };
406 static char* exec_command = nullptr;
407 static real dropunder = 0, dropover = 0;
408 static gmx_bool bRound = FALSE;
411 { "-skip", FALSE, etINT, { &skip_nr }, "Only write every nr-th frame" },
412 { "-dt", FALSE, etTIME, { &delta_t }, "Only write frame when t MOD dt = first time (%t)" },
413 { "-round", FALSE, etBOOL, { &bRound }, "Round measurements to nearest picosecond" },
414 { "-dump", FALSE, etTIME, { &tdump }, "Dump frame nearest specified time (%t)" },
415 { "-t0", FALSE, etTIME, { &tzero }, "Starting time (%t) (default: don't change)" },
416 { "-timestep", FALSE, etTIME, { ×tep }, "Change time step between input frames (%t)" },
417 { "-pbc", FALSE, etENUM, { pbc_opt }, "PBC treatment (see help text for full description)" },
418 { "-ur", FALSE, etENUM, { unitcell_opt }, "Unit-cell representation" },
419 { "-center", FALSE, etBOOL, { &bCenter }, "Center atoms in box" },
420 { "-boxcenter", FALSE, etENUM, { center_opt }, "Center for -pbc and -center" },
421 { "-box", FALSE, etRVEC, { newbox }, "Size for new cubic box (default: read from input)" },
426 "All coordinates will be translated by trans. This "
427 "can advantageously be combined with -pbc mol -ur "
429 { "-shift", FALSE, etRVEC, { shift }, "All coordinates will be shifted by framenr*shift" },
430 { "-fit", FALSE, etENUM, { fit }, "Fit molecule to ref structure in the structure file" },
431 { "-ndec", FALSE, etINT, { &ndec }, "Number of decimal places to write to .xtc output" },
432 { "-vel", FALSE, etBOOL, { &bVels }, "Read and write velocities if possible" },
433 { "-force", FALSE, etBOOL, { &bForce }, "Read and write forces if possible" },
438 "Truncate input trajectory file after this time (%t)" },
443 "Execute command for every output frame with the "
444 "frame number as argument" },
449 "Start writing new file when t MOD split = first "
455 "Write each frame to a separate .gro, .g96 or .pdb "
461 "If the -sep flag is set, use these many digits "
462 "for the file numbers and prepend zeros as needed" },
463 { "-dropunder", FALSE, etREAL, { &dropunder }, "Drop all frames below this value" },
464 { "-dropover", FALSE, etREAL, { &dropover }, "Drop all frames above this value" },
469 "Add conect records when writing [REF].pdb[ref] files. Useful "
470 "for visualization of non-standard molecules, e.g. "
471 "coarse grained ones" }
473 #define NPA asize(pa)
476 t_trxstatus* trxout = nullptr;
479 t_trxframe fr, frout;
481 rvec * xmem = nullptr, *vmem = nullptr, *fmem = nullptr;
482 rvec * xp = nullptr, x_shift, hbox;
483 real* w_rls = nullptr;
484 int m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
486 t_topology* top = nullptr;
487 gmx_conect gc = nullptr;
489 t_atoms * atoms = nullptr, useatoms;
491 int * index = nullptr, *cindex = nullptr;
492 char* grpnm = nullptr;
493 int * frindex, nrfri;
498 int ndrop = 0, ncol, drop0 = 0, drop1 = 0, dropuse = 0;
500 real tshift = 0, dt = -1, prec;
501 gmx_bool bFit, bPFit, bReset;
503 gmx_rmpbc_t gpbc = nullptr;
504 gmx_bool bRmPBC, bPBCWhole, bPBCcomRes, bPBCcomMol, bPBCcomAtom, bPBC, bNoJump, bCluster;
505 gmx_bool bCopy, bDoIt, bIndex, bTDump, bSetTime, bTPS = FALSE, bDTset = FALSE;
506 gmx_bool bExec, bTimeStep = FALSE, bDumpFrame = FALSE, bSetXtcPrec, bNeedPrec;
507 gmx_bool bHaveFirstFrame, bHaveNextFrame, bSetBox, bSetUR, bSplit = FALSE;
508 gmx_bool bDropUnder = FALSE, bDropOver = FALSE, bTrans = FALSE;
509 gmx_bool bWriteFrame, bSplitHere;
510 const char *top_file, *in_file, *out_file = nullptr;
511 char out_file2[256], *charpt;
512 char* outf_base = nullptr;
513 const char* outf_ext = nullptr;
514 char top_title[256], timestr[32], stepstr[32], filemode[5];
515 gmx_output_env_t* oenv;
517 t_filenm fnm[] = { { efTRX, "-f", nullptr, ffREAD }, { efTRO, "-o", nullptr, ffWRITE },
518 { efTPS, nullptr, nullptr, ffOPTRD }, { efNDX, nullptr, nullptr, ffOPTRD },
519 { efNDX, "-fr", "frames", ffOPTRD }, { efNDX, "-sub", "cluster", ffOPTRD },
520 { efXVG, "-drop", "drop", ffOPTRD } };
521 #define NFILE asize(fnm)
523 if (!parse_common_args(&argc, argv, PCA_CAN_BEGIN | PCA_CAN_END | PCA_CAN_VIEW | PCA_TIME_UNIT,
524 NFILE, fnm, NPA, pa, asize(desc), desc, 0, nullptr, &oenv))
529 "Note that major changes are planned in future for "
530 "trjconv, to improve usability and utility.\n");
532 top_file = ftp2fn(efTPS, NFILE, fnm);
534 /* Check command line */
535 in_file = opt2fn("-f", NFILE, fnm);
539 do_trunc(in_file, ttrunc);
543 /* mark active cmdline options */
544 bSetBox = opt2parg_bSet("-box", NPA, pa);
545 bSetTime = opt2parg_bSet("-t0", NPA, pa);
546 bSetXtcPrec = opt2parg_bSet("-ndec", NPA, pa);
547 bSetUR = opt2parg_bSet("-ur", NPA, pa);
548 bExec = opt2parg_bSet("-exec", NPA, pa);
549 bTimeStep = opt2parg_bSet("-timestep", NPA, pa);
550 bTDump = opt2parg_bSet("-dump", NPA, pa);
551 bDropUnder = opt2parg_bSet("-dropunder", NPA, pa);
552 bDropOver = opt2parg_bSet("-dropover", NPA, pa);
553 bTrans = opt2parg_bSet("-trans", NPA, pa);
554 bSplit = (split_t != 0);
556 /* parse enum options */
557 fit_enum = nenum(fit);
558 bFit = (fit_enum == efFit || fit_enum == efFitXY);
559 bReset = (fit_enum == efReset || fit_enum == efResetXY);
560 bPFit = fit_enum == efPFit;
561 pbc_enum = nenum(pbc_opt);
562 bPBCWhole = pbc_enum == epWhole;
563 bPBCcomRes = pbc_enum == epComRes;
564 bPBCcomMol = pbc_enum == epComMol;
565 bPBCcomAtom = pbc_enum == epComAtom;
566 bNoJump = pbc_enum == epNojump;
567 bCluster = pbc_enum == epCluster;
568 bPBC = pbc_enum != epNone;
569 unitcell_enum = nenum(unitcell_opt);
570 ecenter = nenum(center_opt) - ecTric;
572 /* set and check option dependencies */
575 bFit = TRUE; /* for pfit, fit *must* be set */
579 bReset = TRUE; /* for fit, reset *must* be set */
584 nfitdim = (fit_enum == efFitXY || fit_enum == efResetXY) ? 2 : 3;
586 bRmPBC = bFit || bPBCWhole || bPBCcomRes || bPBCcomMol;
590 if (!(bPBCcomRes || bPBCcomMol || bPBCcomAtom))
593 "WARNING: Option for unitcell representation (-ur %s)\n"
594 " only has effect in combination with -pbc %s, %s or %s.\n"
595 " Ingoring unitcell representation.\n\n",
596 unitcell_opt[0], pbc_opt[2], pbc_opt[3], pbc_opt[4]);
602 "PBC condition treatment does not work together with rotational fit.\n"
603 "Please do the PBC condition treatment first and then run trjconv in a "
605 "for the rotational fit.\n"
606 "First doing the rotational fit and then doing the PBC treatment gives "
611 /* ndec for XTC writing is in nr of decimal places, prec is a multiplication factor: */
613 for (i = 0; i < ndec; i++)
618 bIndex = ftp2bSet(efNDX, NFILE, fnm);
621 /* Determine output type */
622 out_file = opt2fn("-o", NFILE, fnm);
623 int ftp = fn2ftp(out_file);
624 fprintf(stderr, "Will write %s: %s\n", ftp2ext(ftp), ftp2desc(ftp));
625 bNeedPrec = (ftp == efXTC);
626 int ftpin = fn2ftp(in_file);
629 /* check if velocities are possible in input and output files */
630 bVels = (ftp == efTRR || ftp == efGRO || ftp == efG96 || ftp == efTNG)
631 && (ftpin == efTRR || ftpin == efGRO || ftpin == efG96 || ftpin == efTNG
634 if (bSeparate || bSplit)
636 outf_ext = std::strrchr(out_file, '.');
637 if (outf_ext == nullptr)
639 gmx_fatal(FARGS, "Output file name '%s' does not contain a '.'", out_file);
641 outf_base = gmx_strdup(out_file);
642 outf_base[outf_ext - out_file] = '\0';
645 bool bSubTraj = opt2bSet("-sub", NFILE, fnm);
649 "The -sub option has been removed from gmx trjconv and is now part\n"
650 "of gmx extract-cluster and does nothing here\n");
656 gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
659 std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
661 /* Determine whether to read a topology */
662 bTPS = (ftp2bSet(efTPS, NFILE, fnm) || bRmPBC || bReset || bPBCcomMol || bCluster
663 || (ftp == efGRO) || (ftp == efPDB) || bCONECT);
665 /* Determine if when can read index groups */
666 bIndex = (bIndex || bTPS);
671 read_tps_conf(top_file, top, &ePBC, &xp, nullptr, top_box, bReset || bPBCcomRes);
672 std::strncpy(top_title, *top->name, 255);
673 top_title[255] = '\0';
676 if (0 == top->mols.nr && (bCluster || bPBCcomMol))
678 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option",
682 /* top_title is only used for gro and pdb,
683 * the header in such a file is top_title, followed by
684 * t= ... and/or step= ...
685 * to prevent double t= or step=, remove it from top_title.
686 * From GROMACS-2018 we only write t/step when the frame actually
687 * has a valid time/step, so we need to check for both separately.
689 if ((charpt = std::strstr(top_title, " t= ")))
693 if ((charpt = std::strstr(top_title, " step= ")))
700 gc = gmx_conect_generate(top);
704 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
708 /* get frame number index */
710 if (opt2bSet("-fr", NFILE, fnm))
712 printf("Select groups of frame number indices:\n");
713 rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
716 for (i = 0; i < nrfri; i++)
718 fprintf(debug, "frindex[%4d]=%4d\n", i, frindex[i]);
723 /* get index groups etc. */
726 printf("Select group for %s fit\n", bFit ? "least squares" : "translational");
727 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit, &ind_fit, &gn_fit);
733 gmx_fatal(FARGS, "Need at least 2 atoms to fit!\n");
737 fprintf(stderr, "WARNING: fitting with only 2 atoms is not unique\n");
743 printf("Select group for clustering\n");
744 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit, &ind_fit, &gn_fit);
751 printf("Select group for centering\n");
752 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ncent, &cindex, &grpnm);
754 printf("Select group for output\n");
755 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nout, &index, &grpnm);
759 /* no index file, so read natoms from TRX */
760 if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
762 gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
768 for (i = 0; i < natoms; i++)
782 snew(w_rls, atoms->nr);
783 for (i = 0; (i < ifit); i++)
785 w_rls[ind_fit[i]] = atoms->atom[ind_fit[i]].m;
788 /* Restore reference structure and set to origin,
789 store original location (to put structure back) */
792 gmx_rmpbc(gpbc, top->atoms.nr, top_box, xp);
794 copy_rvec(xp[index[0]], x_shift);
795 reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
796 rvec_dec(x_shift, xp[index[0]]);
803 if (bDropUnder || bDropOver)
805 /* Read the .xvg file with the drop values */
806 fprintf(stderr, "\nReading drop file ...");
807 ndrop = read_xvg(opt2fn("-drop", NFILE, fnm), &dropval, &ncol);
808 fprintf(stderr, " %d time points\n", ndrop);
809 if (ndrop == 0 || ncol < 2)
811 gmx_fatal(FARGS, "Found no data points in %s", opt2fn("-drop", NFILE, fnm));
817 /* Make atoms struct for output in GRO or PDB files */
818 if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
820 /* get memory for stuff to go in .pdb file, and initialize
821 * the pdbinfo structure part if the input has it.
823 init_t_atoms(&useatoms, atoms->nr, atoms->havePdbInfo);
824 sfree(useatoms.resinfo);
825 useatoms.resinfo = atoms->resinfo;
826 for (i = 0; (i < nout); i++)
828 useatoms.atomname[i] = atoms->atomname[index[i]];
829 useatoms.atom[i] = atoms->atom[index[i]];
830 if (atoms->havePdbInfo)
832 useatoms.pdbinfo[i] = atoms->pdbinfo[index[i]];
834 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind + 1);
838 /* select what to read */
849 flags = flags | TRX_READ_V;
853 flags = flags | TRX_READ_F;
856 /* open trx file for reading */
857 bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
860 fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1 / fr.prec);
864 if (bSetXtcPrec || !fr.bPrec)
866 fprintf(stderr, "\nSetting output precision to %g (nm)\n", 1 / prec);
870 fprintf(stderr, "Using output precision of %g (nm)\n", 1 / prec);
878 // Determine timestep (assuming constant spacing for now) if we
879 // need to dump frames based on time. This is required so we do not
880 // skip the first frame if that was the one that should have been dumped
881 double firstFrameTime = fr.time;
882 if (read_next_frame(oenv, trxin, &fr))
884 dt = fr.time - firstFrameTime;
889 "Warning: Frame times are not incrementing - will dump first "
893 // Now close and reopen so we are at first frame again
896 // Reopen at first frame (We already know it exists if we got here)
897 read_first_frame(oenv, &trxin, in_file, &fr, flags);
900 set_trxframe_ePBC(&fr, ePBC);
905 tshift = tzero - fr.time;
915 /* check if index is meaningful */
916 for (i = 0; i < nout; i++)
918 if (index[i] >= natoms)
921 "Index[%d] %d is larger than the number of atoms in the\n"
922 "trajectory file (%d). There is a mismatch in the contents\n"
923 "of your -f, -s and/or -n files.",
924 i, index[i] + 1, natoms);
926 bCopy = bCopy || (i != index[i]);
930 /* open output for writing */
931 std::strcpy(filemode, "w");
935 trxout = trjtools_gmx_prepare_tng_writing(
936 out_file, filemode[0], trxin, nullptr, nout, mtop.get(),
937 gmx::arrayRefFromArray(index, nout), grpnm);
944 trxout = open_trx(out_file, filemode);
950 if (!bSeparate && !bSplit)
952 out = gmx_ffopen(out_file, filemode);
955 default: gmx_incons("Illegal output file format");
971 /* Start the big loop over frames */
977 /* Main loop over frames */
989 /* generate new box */
994 for (m = 0; m < DIM; m++)
998 fr.box[m][m] = newbox[m];
1004 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
1012 for (i = 0; i < natoms; i++)
1014 rvec_inc(fr.x[i], trans);
1020 // If we could not read two frames or times are not incrementing
1021 // we have almost no idea what to do,
1022 // but dump the first frame so output is not broken.
1023 if (dt <= 0 || !bDTset)
1029 // Dump the frame if we are less than half a frame time
1030 // below it. This will also ensure we at least dump a
1031 // somewhat reasonable frame if the spacing is unequal
1032 // and we have overrun the frame time. Once we dump one
1033 // frame based on time we quit, so it does not matter
1034 // that this might be true for all subsequent frames too.
1035 bDumpFrame = (fr.time > tdump - 0.5 * dt);
1043 /* determine if an atom jumped across the box and reset it if so */
1044 if (bNoJump && (bTPS || frame != 0))
1046 for (d = 0; d < DIM; d++)
1048 hbox[d] = 0.5 * fr.box[d][d];
1050 for (i = 0; i < natoms; i++)
1054 rvec_dec(fr.x[i], x_shift);
1056 for (m = DIM - 1; m >= 0; m--)
1060 while (fr.x[i][m] - xp[i][m] <= -hbox[m])
1062 for (d = 0; d <= m; d++)
1064 fr.x[i][d] += fr.box[m][d];
1067 while (fr.x[i][m] - xp[i][m] > hbox[m])
1069 for (d = 0; d <= m; d++)
1071 fr.x[i][d] -= fr.box[m][d];
1080 calc_pbc_cluster(ecenter, ifit, top, ePBC, fr.x, ind_fit, fr.box);
1085 /* Now modify the coords according to the flags,
1086 for normal fit, this is only done for output frames */
1089 gmx_rmpbc_trxfr(gpbc, &fr);
1092 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1093 do_fit(natoms, w_rls, xp, fr.x);
1096 /* store this set of coordinates for future use */
1097 if (bPFit || bNoJump)
1103 for (i = 0; (i < natoms); i++)
1105 copy_rvec(fr.x[i], xp[i]);
1106 rvec_inc(fr.x[i], x_shift);
1112 /* see if we have a frame from the frame index group */
1113 for (i = 0; i < nrfri && !bDumpFrame; i++)
1115 bDumpFrame = frame == frindex[i];
1118 if (debug && bDumpFrame)
1120 fprintf(debug, "dumping %d\n", frame);
1123 bWriteFrame = ((!bTDump && (frindex == nullptr) && frame % skip_nr == 0) || bDumpFrame);
1125 if (bWriteFrame && (bDropUnder || bDropOver))
1127 while (dropval[0][drop1] < fr.time && drop1 + 1 < ndrop)
1132 if (std::abs(dropval[0][drop0] - fr.time) < std::abs(dropval[0][drop1] - fr.time))
1140 if ((bDropUnder && dropval[1][dropuse] < dropunder)
1141 || (bDropOver && dropval[1][dropuse] > dropover))
1143 bWriteFrame = FALSE;
1149 /* We should avoid modifying the input frame,
1150 * but since here we don't have the output frame yet,
1151 * we introduce a temporary output frame time variable.
1155 frout_time = fr.time;
1160 frout_time = tzero + frame * timestep;
1164 frout_time += tshift;
1169 fprintf(stderr, "\nDumping frame at t= %g %s\n",
1170 output_env_conv_time(oenv, frout_time),
1171 output_env_get_time_unit(oenv).c_str());
1174 /* check for writing at each delta_t */
1175 bDoIt = (delta_t == 0);
1180 bDoIt = bRmod(frout_time, tzero, delta_t);
1184 /* round() is not C89 compatible, so we do this: */
1185 bDoIt = bRmod(std::floor(frout_time + 0.5), std::floor(tzero + 0.5),
1186 std::floor(delta_t + 0.5));
1190 if (bDoIt || bTDump)
1192 /* print sometimes */
1193 if (((outframe % SKIP) == 0) || (outframe < SKIP))
1195 fprintf(stderr, " -> frame %6d time %8.3f \r", outframe,
1196 output_env_conv_time(oenv, frout_time));
1202 /* Now modify the coords according to the flags,
1203 for PFit we did this already! */
1207 gmx_rmpbc_trxfr(gpbc, &fr);
1212 reset_x_ndim(nfitdim, ifit, ind_fit, natoms, nullptr, fr.x, w_rls);
1215 do_fit_ndim(nfitdim, natoms, w_rls, xp, fr.x);
1219 for (i = 0; i < natoms; i++)
1221 rvec_inc(fr.x[i], x_shift);
1228 center_x(ecenter, fr.x, fr.box, natoms, ncent, cindex);
1232 auto positionsArrayRef =
1233 gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec*>(fr.x), natoms);
1236 switch (unitcell_enum)
1239 put_atoms_in_box(ePBC, fr.box, positionsArrayRef);
1242 put_atoms_in_triclinic_unitcell(ecenter, fr.box, positionsArrayRef);
1245 put_atoms_in_compact_unitcell(ePBC, ecenter, fr.box, positionsArrayRef);
1251 put_residue_com_in_box(unitcell_enum, ecenter, natoms, atoms->atom,
1252 ePBC, fr.box, fr.x);
1256 put_molecule_com_in_box(unitcell_enum, ecenter, &top->mols, natoms,
1257 atoms->atom, ePBC, fr.box, fr.x);
1259 /* Copy the input trxframe struct to the output trxframe struct */
1261 frout.time = frout_time;
1262 frout.bV = (frout.bV && bVels);
1263 frout.bF = (frout.bF && bForce);
1264 frout.natoms = nout;
1265 if (bNeedPrec && (bSetXtcPrec || !fr.bPrec))
1281 for (i = 0; i < nout; i++)
1283 copy_rvec(fr.x[index[i]], frout.x[i]);
1286 copy_rvec(fr.v[index[i]], frout.v[i]);
1290 copy_rvec(fr.f[index[i]], frout.f[i]);
1295 if (opt2parg_bSet("-shift", NPA, pa))
1297 for (i = 0; i < nout; i++)
1299 for (d = 0; d < DIM; d++)
1301 frout.x[i][d] += outframe * shift[d];
1308 bSplitHere = bSplit && bRmod(frout.time, tzero, split_t);
1312 /* round() is not C89 compatible, so we do this: */
1314 && bRmod(std::floor(frout.time + 0.5),
1315 std::floor(tzero + 0.5), std::floor(split_t + 0.5));
1317 if (bSeparate || bSplitHere)
1319 mk_filenm(outf_base, ftp2ext(ftp), nzero, file_nr, out_file2);
1326 write_tng_frame(trxout, &frout);
1327 // TODO when trjconv behaves better: work how to read and write lambda
1337 trxout = open_trx(out_file2, filemode);
1339 write_trxframe(trxout, &frout, gc);
1344 // Only add a generator statement if title is empty,
1345 // to avoid multiple generated-by statements from various programs
1346 if (std::strlen(top_title) == 0)
1348 sprintf(top_title, "Generated by trjconv");
1352 sprintf(timestr, " t= %9.5f", frout.time);
1356 std::strcpy(timestr, "");
1360 sprintf(stepstr, " step= %" PRId64, frout.step);
1364 std::strcpy(stepstr, "");
1366 title = gmx::formatString("%s%s%s", top_title, timestr, stepstr);
1367 if (bSeparate || bSplitHere)
1369 out = gmx_ffopen(out_file2, "w");
1374 write_hconf_p(out, title.c_str(), &useatoms, frout.x,
1375 frout.bV ? frout.v : nullptr, frout.box);
1378 fprintf(out, "REMARK GENERATED BY TRJCONV\n");
1379 /* if reading from pdb, we want to keep the original
1380 model numbering else we write the output frame
1381 number plus one, because model 0 is not allowed in pdb */
1382 if (ftpin == efPDB && fr.bStep && fr.step > model_nr)
1390 write_pdbfile(out, title.c_str(), &useatoms, frout.x,
1391 frout.ePBC, frout.box, ' ', model_nr, gc);
1394 const char* outputTitle = "";
1395 if (bSeparate || bTDump)
1397 outputTitle = title.c_str();
1400 frout.bAtoms = TRUE;
1402 frout.atoms = &useatoms;
1403 frout.bStep = FALSE;
1404 frout.bTime = FALSE;
1410 outputTitle = title.c_str();
1412 frout.bAtoms = FALSE;
1416 write_g96_conf(out, outputTitle, &frout, -1, nullptr);
1418 if (bSeparate || bSplitHere)
1424 default: gmx_fatal(FARGS, "DHE, ftp=%d\n", ftp);
1426 if (bSeparate || bSplitHere)
1431 /* execute command */
1435 sprintf(c, "%s %d", exec_command, file_nr - 1);
1436 /*fprintf(stderr,"Executing '%s'\n",c);*/
1439 gmx_fatal(FARGS, "Error executing command: %s", c);
1446 bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
1447 } while (!(bTDump && bDumpFrame) && bHaveNextFrame);
1450 if (!bHaveFirstFrame || (bTDump && !bDumpFrame))
1453 "\nWARNING no output, "
1454 "last frame read at t=%g\n",
1457 fprintf(stderr, "\n");
1464 gmx_rmpbc_done(gpbc);
1471 else if (out != nullptr)
1491 do_view(oenv, out_file, nullptr);
1493 output_env_done(oenv);