#include "gromacs/fileio/tpxio.h"
#include "gromacs/fileio/trxio.h"
#include "gromacs/fileio/trnio.h"
+#include "gromacs/fileio/tngio_for_tools.h"
#include "gromacs/commandline/pargs.h"
#include "gromacs/fileio/futil.h"
#include "gromacs/fileio/pdbio.h"
}
#endif
+/*! \brief Read a full molecular topology if useful and available.
+ *
+ * If the input trajectory file is not in TNG format, and the output
+ * file is in TNG format, then we want to try to read a full topology
+ * (if available), so that we can write molecule information to the
+ * output file. The full topology provides better molecule information
+ * than is available from the normal t_topology data used by GROMACS
+ * tools.
+ *
+ * Also, the t_topology is only read under (different) particular
+ * conditions. If both apply, then a .tpr file might be read
+ * twice. Trying to fix this redundancy while trjconv is still an
+ * all-purpose tool does not seem worthwhile.
+ *
+ * Because of the way gmx_prepare_tng_writing is implemented, the case
+ * where the input TNG file has no molecule information will never
+ * lead to an output TNG file having molecule information. Since
+ * molecule information will generally be present if the input TNG
+ * file was written by a GROMACS tool, this seems like reasonable
+ * behaviour. */
+static gmx_mtop_t *read_mtop_for_tng(const char *tps_file,
+ const char *input_file,
+ const char *output_file)
+{
+ gmx_mtop_t *mtop = NULL;
+
+ if (fn2bTPX(tps_file) &&
+ efTNG != fn2ftp(input_file) &&
+ efTNG == fn2ftp(output_file))
+ {
+ int temp_natoms = -1;
+ snew(mtop, 1);
+ read_tpx(tps_file, NULL, NULL, &temp_natoms,
+ NULL, NULL, NULL, mtop);
+ }
+
+ return mtop;
+}
+
int gmx_trjconv(int argc, char *argv[])
{
const char *desc[] = {
FILE *out = NULL;
t_trxstatus *trxout = NULL;
- t_trxstatus *status;
+ t_trxstatus *trxin;
int ftp, ftpin = 0, file_nr;
t_trxframe fr, frout;
int flags;
int m, i, d, frame, outframe, natoms, nout, ncent, nre, newstep = 0, model_nr;
#define SKIP 10
t_topology top;
+ gmx_mtop_t *mtop = NULL;
gmx_conect gc = NULL;
int ePBC = -1;
t_atoms *atoms = NULL, useatoms;
{
}
+ mtop = read_mtop_for_tng(top_file, in_file, out_file);
+
/* Determine whether to read a topology */
bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
bRmPBC || bReset || bPBCcomMol || bCluster ||
else
{
/* no index file, so read natoms from TRX */
- if (!read_first_frame(oenv, &status, in_file, &fr, TRX_DONT_SKIP))
+ if (!read_first_frame(oenv, &trxin, in_file, &fr, TRX_DONT_SKIP))
{
gmx_fatal(FARGS, "Could not read a frame from %s", in_file);
}
natoms = fr.natoms;
- close_trj(status);
+ close_trj(trxin);
sfree(fr.x);
snew(index, natoms);
for (i = 0; i < natoms; i++)
}
/* open trx file for reading */
- bHaveFirstFrame = read_first_frame(oenv, &status, in_file, &fr, flags);
+ bHaveFirstFrame = read_first_frame(oenv, &trxin, in_file, &fr, flags);
if (fr.bPrec)
{
fprintf(stderr, "\nPrecision of %s is %g (nm)\n", in_file, 1/fr.prec);
tzero = fr.time;
}
+ bCopy = FALSE;
+ if (bIndex)
+ {
+ /* check if index is meaningful */
+ for (i = 0; i < nout; i++)
+ {
+ if (index[i] >= natoms)
+ {
+ gmx_fatal(FARGS,
+ "Index[%d] %d is larger than the number of atoms in the\n"
+ "trajectory file (%d). There is a mismatch in the contents\n"
+ "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
+ }
+ bCopy = bCopy || (i != index[i]);
+ }
+ }
+
/* open output for writing */
if ((bAppend) && (gmx_fexist(out_file)))
{
}
switch (ftp)
{
+ case efTNG:
+ trjtools_gmx_prepare_tng_writing(out_file,
+ filemode[0],
+ trxin,
+ &trxout,
+ NULL,
+ nout,
+ mtop,
+ index,
+ grpnm);
+ break;
case efXTC:
case efG87:
case efTRR:
out = ffopen(out_file, filemode);
}
break;
+ default:
+ gmx_incons("Illegal output file format");
}
- bCopy = FALSE;
- if (bIndex)
- {
- /* check if index is meaningful */
- for (i = 0; i < nout; i++)
- {
- if (index[i] >= natoms)
- {
- gmx_fatal(FARGS,
- "Index[%d] %d is larger than the number of atoms in the\n"
- "trajectory file (%d). There is a mismatch in the contents\n"
- "of your -f, -s and/or -n files.", i, index[i]+1, natoms);
- }
- bCopy = bCopy || (i != index[i]);
- }
- }
if (bCopy)
{
snew(xmem, nout);
switch (ftp)
{
+ case efTNG:
+ write_tng_frame(trxout, &frout);
+ // TODO when trjconv behaves better: work how to read and write lambda
+ break;
case efTRJ:
case efTRR:
case efG87:
}
}
frame++;
- bHaveNextFrame = read_next_frame(oenv, status, &fr);
+ bHaveNextFrame = read_next_frame(oenv, trxin, &fr);
}
while (!(bTDump && bDumpFrame) && bHaveNextFrame);
}
}
fprintf(stderr, "\n");
- close_trj(status);
+ close_trj(trxin);
sfree(outf_base);
if (bRmPBC)
}
}
+ sfree(mtop);
+
do_view(oenv, out_file, NULL);
return 0;