b515578ef3edbe09a73ec781906a5061b859b41d
[alexxy/gromacs.git] / src / gromacs / fileio / trxio.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "trxio.h"
40
41 #include "config.h"
42
43 #include <assert.h>
44 #include <math.h>
45
46 #ifdef GMX_USE_PLUGINS
47 #include "vmdio.h"
48 #endif
49 #include "gmxfio.h"
50 #include "trxio.h"
51 #include "tpxio.h"
52 #include "trnio.h"
53 #include "tngio.h"
54 #include "tngio_for_tools.h"
55 #include "gromacs/legacyheaders/names.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/utility/futil.h"
58 #include "xtcio.h"
59 #include "pdbio.h"
60 #include "confio.h"
61 #include "gromacs/legacyheaders/checkpoint.h"
62 #include "xdrf.h"
63
64 #include "gromacs/fileio/timecontrol.h"
65 #include "gromacs/fileio/trx.h"
66 #include "gromacs/topology/atoms.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/smalloc.h"
69
70 /* defines for frame counter output */
71 #define SKIP1   10
72 #define SKIP2  100
73 #define SKIP3 1000
74
75 struct t_trxstatus
76 {
77     int                     __frame;
78     t_trxframe             *xframe;
79     int                     nxframe;
80     t_fileio               *fio;
81     tng_trajectory_t        tng;
82     int                     NATOMS;
83     double                  DT, BOX[3];
84     gmx_bool                bReadBox;
85     char                   *persistent_line; /* Persistent line for reading g96 trajectories */
86 };
87
88 /* utility functions */
89
90 gmx_bool bRmod_fd(double a, double b, double c, gmx_bool bDouble)
91 {
92     int    iq;
93     double tol;
94
95     tol = 2*(bDouble ? GMX_DOUBLE_EPS : GMX_FLOAT_EPS);
96
97     iq = (int)((a - b + tol*a)/c);
98
99     if (fabs(a - b - c*iq) <= tol*fabs(a))
100     {
101         return TRUE;
102     }
103     else
104     {
105         return FALSE;
106     }
107 }
108
109
110
111 int check_times2(real t, real t0, gmx_bool bDouble)
112 {
113     int  r;
114
115 #ifndef GMX_DOUBLE
116     /* since t is float, we can not use double precision for bRmod */
117     bDouble = FALSE;
118 #endif
119
120     r = -1;
121     if ((!bTimeSet(TBEGIN) || (t >= rTimeValue(TBEGIN)))  &&
122         (!bTimeSet(TEND)   || (t <= rTimeValue(TEND))))
123     {
124         if (bTimeSet(TDELTA) && !bRmod_fd(t, t0, rTimeValue(TDELTA), bDouble))
125         {
126             r = -1;
127         }
128         else
129         {
130             r = 0;
131         }
132     }
133     else if (bTimeSet(TEND) && (t >= rTimeValue(TEND)))
134     {
135         r = 1;
136     }
137     if (debug)
138     {
139         fprintf(debug, "t=%g, t0=%g, b=%g, e=%g, dt=%g: r=%d\n",
140                 t, t0, rTimeValue(TBEGIN), rTimeValue(TEND), rTimeValue(TDELTA), r);
141     }
142     return r;
143 }
144
145 int check_times(real t)
146 {
147     return check_times2(t, t, FALSE);
148 }
149
150 static void initcount(t_trxstatus *status)
151 {
152     status->__frame = -1;
153 }
154
155 static void status_init(t_trxstatus *status)
156 {
157     status->nxframe         = 0;
158     status->xframe          = NULL;
159     status->fio             = NULL;
160     status->__frame         = -1;
161     status->persistent_line = NULL;
162     status->tng             = NULL;
163 }
164
165
166 int nframes_read(t_trxstatus *status)
167 {
168     return status->__frame;
169 }
170
171 static void printcount_(t_trxstatus *status, const output_env_t oenv,
172                         const char *l, real t)
173 {
174     if ((status->__frame < 2*SKIP1 || status->__frame % SKIP1 == 0) &&
175         (status->__frame < 2*SKIP2 || status->__frame % SKIP2 == 0) &&
176         (status->__frame < 2*SKIP3 || status->__frame % SKIP3 == 0))
177     {
178         fprintf(stderr, "\r%-14s %6d time %8.3f   ", l, status->__frame,
179                 output_env_conv_time(oenv, t));
180     }
181 }
182
183 static void printcount(t_trxstatus *status, const output_env_t oenv, real t,
184                        gmx_bool bSkip)
185 {
186     status->__frame++;
187     printcount_(status, oenv, bSkip ? "Skipping frame" : "Reading frame", t);
188 }
189
190 static void printlast(t_trxstatus *status, const output_env_t oenv, real t)
191 {
192     printcount_(status, oenv, "Last frame", t);
193     fprintf(stderr, "\n");
194 }
195
196 static void printincomp(t_trxstatus *status, t_trxframe *fr)
197 {
198     if (fr->not_ok & HEADER_NOT_OK)
199     {
200         fprintf(stderr, "WARNING: Incomplete header: nr %d time %g\n",
201                 status->__frame+1, fr->time);
202     }
203     else if (fr->not_ok)
204     {
205         fprintf(stderr, "WARNING: Incomplete frame: nr %d time %g\n",
206                 status->__frame+1, fr->time);
207     }
208 }
209
210 int prec2ndec(real prec)
211 {
212     if (prec <= 0)
213     {
214         gmx_fatal(FARGS, "DEATH HORROR prec (%g) <= 0 in prec2ndec", prec);
215     }
216
217     return (int)(log(prec)/log(10.0)+0.5);
218 }
219
220 real ndec2prec(int ndec)
221 {
222     return pow(10.0, ndec);
223 }
224
225 t_fileio *trx_get_fileio(t_trxstatus *status)
226 {
227     return status->fio;
228 }
229
230 float trx_get_time_of_final_frame(t_trxstatus *status)
231 {
232     t_fileio *stfio    = trx_get_fileio(status);
233     int       filetype = gmx_fio_getftp(stfio);
234     int       bOK;
235     float     lasttime = -1;
236
237     if (filetype == efXTC)
238     {
239         lasttime =
240             xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio),
241                                         gmx_fio_getxdr(stfio),
242                                         status->xframe->natoms, &bOK);
243         if (!bOK)
244         {
245             gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported." );
246         }
247     }
248     else if (filetype == efTNG)
249     {
250         tng_trajectory_t tng = status->tng;
251         if (!tng)
252         {
253             gmx_fatal(FARGS, "Error opening TNG file.");
254         }
255         lasttime = gmx_tng_get_time_of_final_frame(tng);
256     }
257     else
258     {
259         gmx_incons("Only supported for TNG and XTC");
260     }
261     return lasttime;
262 }
263
264 void clear_trxframe(t_trxframe *fr, gmx_bool bFirst)
265 {
266     fr->not_ok    = 0;
267     fr->bTitle    = FALSE;
268     fr->bStep     = FALSE;
269     fr->bTime     = FALSE;
270     fr->bLambda   = FALSE;
271     fr->bFepState = FALSE;
272     fr->bAtoms    = FALSE;
273     fr->bPrec     = FALSE;
274     fr->bX        = FALSE;
275     fr->bV        = FALSE;
276     fr->bF        = FALSE;
277     fr->bBox      = FALSE;
278     if (bFirst)
279     {
280         fr->flags     = 0;
281         fr->bDouble   = FALSE;
282         fr->natoms    = -1;
283         fr->t0        = 0;
284         fr->tf        = 0;
285         fr->tpf       = 0;
286         fr->tppf      = 0;
287         fr->title     = NULL;
288         fr->step      = 0;
289         fr->time      = 0;
290         fr->lambda    = 0;
291         fr->fep_state = 0;
292         fr->atoms     = NULL;
293         fr->prec      = 0;
294         fr->x         = NULL;
295         fr->v         = NULL;
296         fr->f         = NULL;
297         clear_mat(fr->box);
298         fr->bPBC   = FALSE;
299         fr->ePBC   = -1;
300     }
301 }
302
303 void set_trxframe_ePBC(t_trxframe *fr, int ePBC)
304 {
305     fr->bPBC = (ePBC == -1);
306     fr->ePBC = ePBC;
307 }
308
309 int write_trxframe_indexed(t_trxstatus *status, t_trxframe *fr, int nind,
310                            const atom_id *ind, gmx_conect gc)
311 {
312     char  title[STRLEN];
313     rvec *xout = NULL, *vout = NULL, *fout = NULL;
314     int   i, ftp = -1;
315     real  prec;
316
317     if (fr->bPrec)
318     {
319         prec = fr->prec;
320     }
321     else
322     {
323         prec = 1000.0;
324     }
325
326     if (status->tng)
327     {
328         ftp = efTNG;
329     }
330     else if (status->fio)
331     {
332         ftp = gmx_fio_getftp(status->fio);
333     }
334     else
335     {
336         gmx_incons("No input file available");
337     }
338
339     switch (ftp)
340     {
341         case efTRJ:
342         case efTRR:
343         case efTNG:
344             break;
345         default:
346             if (!fr->bX)
347             {
348                 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory",
349                           ftp2ext(ftp));
350             }
351             break;
352     }
353
354     switch (ftp)
355     {
356         case efTRJ:
357         case efTRR:
358         case efTNG:
359             if (fr->bV)
360             {
361                 snew(vout, nind);
362                 for (i = 0; i < nind; i++)
363                 {
364                     copy_rvec(fr->v[ind[i]], vout[i]);
365                 }
366             }
367             if (fr->bF)
368             {
369                 snew(fout, nind);
370                 for (i = 0; i < nind; i++)
371                 {
372                     copy_rvec(fr->f[ind[i]], fout[i]);
373                 }
374             }
375         /* no break */
376         case efXTC:
377             if (fr->bX)
378             {
379                 snew(xout, nind);
380                 for (i = 0; i < nind; i++)
381                 {
382                     copy_rvec(fr->x[ind[i]], xout[i]);
383                 }
384             }
385             break;
386         default:
387             break;
388     }
389
390     switch (ftp)
391     {
392         case efTNG:
393             gmx_write_tng_from_trxframe(status->tng, fr, nind);
394             break;
395         case efXTC:
396             write_xtc(status->fio, nind, fr->step, fr->time, fr->box, xout, prec);
397             break;
398         case efTRJ:
399         case efTRR:
400             fwrite_trn(status->fio, nframes_read(status),
401                        fr->time, fr->step, fr->box, nind, xout, vout, fout);
402             break;
403         case efGRO:
404         case efPDB:
405         case efBRK:
406         case efENT:
407             if (!fr->bAtoms)
408             {
409                 gmx_fatal(FARGS, "Can not write a %s file without atom names",
410                           ftp2ext(ftp));
411             }
412             sprintf(title, "frame t= %.3f", fr->time);
413             if (ftp == efGRO)
414             {
415                 write_hconf_indexed_p(gmx_fio_getfp(status->fio), title, fr->atoms, nind, ind,
416                                       prec2ndec(prec),
417                                       fr->x, fr->bV ? fr->v : NULL, fr->box);
418             }
419             else
420             {
421                 write_pdbfile_indexed(gmx_fio_getfp(status->fio), title, fr->atoms,
422                                       fr->x, -1, fr->box, ' ', fr->step, nind, ind, gc, TRUE);
423             }
424             break;
425         case efG96:
426             write_g96_conf(gmx_fio_getfp(status->fio), fr, nind, ind);
427             break;
428         default:
429             gmx_fatal(FARGS, "Sorry, write_trxframe_indexed can not write %s",
430                       ftp2ext(ftp));
431             break;
432     }
433
434     switch (ftp)
435     {
436         case efTRN:
437         case efTRJ:
438         case efTRR:
439         case efTNG:
440             if (vout)
441             {
442                 sfree(vout);
443             }
444             if (fout)
445             {
446                 sfree(fout);
447             }
448         /* no break */
449         case efXTC:
450             sfree(xout);
451             break;
452         default:
453             break;
454     }
455
456     return 0;
457 }
458
459 void trjtools_gmx_prepare_tng_writing(const char       *filename,
460                                       char              filemode,
461                                       t_trxstatus      *in,
462                                       t_trxstatus     **out,
463                                       const char       *infile,
464                                       const int         natoms,
465                                       const gmx_mtop_t *mtop,
466                                       const atom_id    *index,
467                                       const char       *index_group_name)
468 {
469     if (filemode != 'w' && filemode != 'a')
470     {
471         gmx_incons("Sorry, can only prepare for TNG output.");
472     }
473
474     if (*out == NULL)
475     {
476         snew((*out), 1);
477     }
478     status_init(*out);
479
480     if (in != NULL)
481     {
482         gmx_prepare_tng_writing(filename,
483                                 filemode,
484                                 &in->tng,
485                                 &(*out)->tng,
486                                 natoms,
487                                 mtop,
488                                 index,
489                                 index_group_name);
490     }
491     else if (efTNG == fn2ftp(infile))
492     {
493         tng_trajectory_t tng_in;
494         gmx_tng_open(infile, 'r', &tng_in);
495
496         gmx_prepare_tng_writing(filename,
497                                 filemode,
498                                 &tng_in,
499                                 &(*out)->tng,
500                                 natoms,
501                                 mtop,
502                                 index,
503                                 index_group_name);
504     }
505 }
506
507 void write_tng_frame(t_trxstatus *status,
508                      t_trxframe  *frame)
509 {
510     gmx_write_tng_from_trxframe(status->tng, frame, -1);
511 }
512
513 int write_trxframe(t_trxstatus *status, t_trxframe *fr, gmx_conect gc)
514 {
515     char title[STRLEN];
516     real prec;
517
518     if (fr->bPrec)
519     {
520         prec = fr->prec;
521     }
522     else
523     {
524         prec = 1000.0;
525     }
526
527     if (status->tng)
528     {
529         gmx_tng_set_compression_precision(status->tng, prec);
530         write_tng_frame(status, fr);
531
532         return 0;
533     }
534
535     switch (gmx_fio_getftp(status->fio))
536     {
537         case efTRJ:
538         case efTRR:
539             break;
540         default:
541             if (!fr->bX)
542             {
543                 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory",
544                           ftp2ext(gmx_fio_getftp(status->fio)));
545             }
546             break;
547     }
548
549     switch (gmx_fio_getftp(status->fio))
550     {
551         case efXTC:
552             write_xtc(status->fio, fr->natoms, fr->step, fr->time, fr->box, fr->x, prec);
553             break;
554         case efTRJ:
555         case efTRR:
556             fwrite_trn(status->fio, fr->step, fr->time, fr->lambda, fr->box, fr->natoms,
557                        fr->bX ? fr->x : NULL, fr->bV ? fr->v : NULL, fr->bF ? fr->f : NULL);
558             break;
559         case efGRO:
560         case efPDB:
561         case efBRK:
562         case efENT:
563             if (!fr->bAtoms)
564             {
565                 gmx_fatal(FARGS, "Can not write a %s file without atom names",
566                           ftp2ext(gmx_fio_getftp(status->fio)));
567             }
568             sprintf(title, "frame t= %.3f", fr->time);
569             if (gmx_fio_getftp(status->fio) == efGRO)
570             {
571                 write_hconf_p(gmx_fio_getfp(status->fio), title, fr->atoms,
572                               prec2ndec(prec), fr->x, fr->bV ? fr->v : NULL, fr->box);
573             }
574             else
575             {
576                 write_pdbfile(gmx_fio_getfp(status->fio), title,
577                               fr->atoms, fr->x, fr->bPBC ? fr->ePBC : -1, fr->box,
578                               ' ', fr->step, gc, TRUE);
579             }
580             break;
581         case efG96:
582             write_g96_conf(gmx_fio_getfp(status->fio), fr, -1, NULL);
583             break;
584         default:
585             gmx_fatal(FARGS, "Sorry, write_trxframe can not write %s",
586                       ftp2ext(gmx_fio_getftp(status->fio)));
587             break;
588     }
589
590     return 0;
591 }
592
593 int write_trx(t_trxstatus *status, int nind, const atom_id *ind, t_atoms *atoms,
594               int step, real time, matrix box, rvec x[], rvec *v,
595               gmx_conect gc)
596 {
597     t_trxframe fr;
598
599     clear_trxframe(&fr, TRUE);
600     fr.bStep  = TRUE;
601     fr.step   = step;
602     fr.bTime  = TRUE;
603     fr.time   = time;
604     fr.bAtoms = atoms != NULL;
605     fr.atoms  = atoms;
606     fr.bX     = TRUE;
607     fr.x      = x;
608     fr.bV     = v != NULL;
609     fr.v      = v;
610     fr.bBox   = TRUE;
611     copy_mat(box, fr.box);
612
613     return write_trxframe_indexed(status, &fr, nind, ind, gc);
614 }
615
616 void close_trx(t_trxstatus *status)
617 {
618     gmx_tng_close(&status->tng);
619     if (status->fio)
620     {
621         gmx_fio_close(status->fio);
622     }
623     sfree(status);
624 }
625
626 t_trxstatus *open_trx(const char *outfile, const char *filemode)
627 {
628     t_trxstatus *stat;
629     if (filemode[0] != 'w' && filemode[0] != 'a' && filemode[1] != '+')
630     {
631         gmx_fatal(FARGS, "Sorry, write_trx can only write");
632     }
633
634     snew(stat, 1);
635     status_init(stat);
636
637     stat->fio = gmx_fio_open(outfile, filemode);
638     return stat;
639 }
640
641 static gmx_bool gmx_next_frame(t_trxstatus *status, t_trxframe *fr)
642 {
643     t_trnheader sh;
644     gmx_bool    bOK, bRet;
645
646     bRet = FALSE;
647
648     if (fread_trnheader(status->fio, &sh, &bOK))
649     {
650         fr->bDouble   = sh.bDouble;
651         fr->natoms    = sh.natoms;
652         fr->bStep     = TRUE;
653         fr->step      = sh.step;
654         fr->bTime     = TRUE;
655         fr->time      = sh.t;
656         fr->bLambda   = TRUE;
657         fr->bFepState = TRUE;
658         fr->lambda    = sh.lambda;
659         fr->bBox      = sh.box_size > 0;
660         if (fr->flags & (TRX_READ_X | TRX_NEED_X))
661         {
662             if (fr->x == NULL)
663             {
664                 snew(fr->x, sh.natoms);
665             }
666             fr->bX = sh.x_size > 0;
667         }
668         if (fr->flags & (TRX_READ_V | TRX_NEED_V))
669         {
670             if (fr->v == NULL)
671             {
672                 snew(fr->v, sh.natoms);
673             }
674             fr->bV = sh.v_size > 0;
675         }
676         if (fr->flags & (TRX_READ_F | TRX_NEED_F))
677         {
678             if (fr->f == NULL)
679             {
680                 snew(fr->f, sh.natoms);
681             }
682             fr->bF = sh.f_size > 0;
683         }
684         if (fread_htrn(status->fio, &sh, fr->box, fr->x, fr->v, fr->f))
685         {
686             bRet = TRUE;
687         }
688         else
689         {
690             fr->not_ok = DATA_NOT_OK;
691         }
692     }
693     else
694     if (!bOK)
695     {
696         fr->not_ok = HEADER_NOT_OK;
697     }
698
699     return bRet;
700 }
701
702 static gmx_bool pdb_next_x(t_trxstatus *status, FILE *fp, t_trxframe *fr)
703 {
704     t_atoms   atoms;
705     matrix    boxpdb;
706     int       ePBC, model_nr, na;
707     char      title[STRLEN], *time;
708     double    dbl;
709
710     atoms.nr      = fr->natoms;
711     atoms.atom    = NULL;
712     atoms.pdbinfo = NULL;
713     /* the other pointers in atoms should not be accessed if these are NULL */
714     model_nr = NOTSET;
715     na       = read_pdbfile(fp, title, &model_nr, &atoms, fr->x, &ePBC, boxpdb, TRUE, NULL);
716     set_trxframe_ePBC(fr, ePBC);
717     if (nframes_read(status) == 0)
718     {
719         fprintf(stderr, " '%s', %d atoms\n", title, fr->natoms);
720     }
721     fr->bPrec = TRUE;
722     fr->prec  = 10000;
723     fr->bX    = TRUE;
724     fr->bBox  = (boxpdb[XX][XX] != 0.0);
725     if (fr->bBox)
726     {
727         copy_mat(boxpdb, fr->box);
728     }
729
730     if (model_nr != NOTSET)
731     {
732         fr->bStep = TRUE;
733         fr->step  = model_nr;
734     }
735     time = strstr(title, " t= ");
736     if (time)
737     {
738         fr->bTime = TRUE;
739         sscanf(time+4, "%lf", &dbl);
740         fr->time = (real)dbl;
741     }
742     else
743     {
744         fr->bTime = FALSE;
745         /* this is a bit dirty, but it will work: if no time is read from
746            comment line in pdb file, set time to current frame number */
747         if (fr->bStep)
748         {
749             fr->time = (real)fr->step;
750         }
751         else
752         {
753             fr->time = (real)nframes_read(status);
754         }
755     }
756     if (na == 0)
757     {
758         return FALSE;
759     }
760     else
761     {
762         if (na != fr->natoms)
763         {
764             gmx_fatal(FARGS, "Number of atoms in pdb frame %d is %d instead of %d",
765                       nframes_read(status), na, fr->natoms);
766         }
767         return TRUE;
768     }
769 }
770
771 static int pdb_first_x(t_trxstatus *status, FILE *fp, t_trxframe *fr)
772 {
773     initcount(status);
774
775     fprintf(stderr, "Reading frames from pdb file");
776     frewind(fp);
777     get_pdb_coordnum(fp, &fr->natoms);
778     if (fr->natoms == 0)
779     {
780         gmx_fatal(FARGS, "\nNo coordinates in pdb file\n");
781     }
782     frewind(fp);
783     snew(fr->x, fr->natoms);
784     pdb_next_x(status, fp, fr);
785
786     return fr->natoms;
787 }
788
789 gmx_bool read_next_frame(const output_env_t oenv, t_trxstatus *status, t_trxframe *fr)
790 {
791     real     pt;
792     int      ct;
793     gmx_bool bOK, bRet, bMissingData = FALSE, bSkip = FALSE;
794     int      dummy = 0;
795     int      ftp;
796
797     bRet = FALSE;
798     pt   = fr->tf;
799
800     do
801     {
802         clear_trxframe(fr, FALSE);
803         fr->tppf = fr->tpf;
804         fr->tpf  = fr->tf;
805
806         if (status->tng)
807         {
808             /* Special treatment for TNG files */
809             ftp = efTNG;
810         }
811         else
812         {
813             ftp = gmx_fio_getftp(status->fio);
814         }
815         switch (ftp)
816         {
817             case efTRJ:
818             case efTRR:
819                 bRet = gmx_next_frame(status, fr);
820                 break;
821             case efCPT:
822                 /* Checkpoint files can not contain mulitple frames */
823                 break;
824             case efG96:
825                 read_g96_conf(gmx_fio_getfp(status->fio), NULL, fr,
826                               status->persistent_line);
827                 bRet = (fr->natoms > 0);
828                 break;
829             case efXTC:
830                 /* B. Hess 2005-4-20
831                  * Sometimes is off by one frame
832                  * and sometimes reports frame not present/file not seekable
833                  */
834                 /* DvdS 2005-05-31: this has been fixed along with the increased
835                  * accuracy of the control over -b and -e options.
836                  */
837                 if (bTimeSet(TBEGIN) && (fr->tf < rTimeValue(TBEGIN)))
838                 {
839                     if (xtc_seek_time(status->fio, rTimeValue(TBEGIN), fr->natoms, TRUE))
840                     {
841                         gmx_fatal(FARGS, "Specified frame (time %f) doesn't exist or file corrupt/inconsistent.",
842                                   rTimeValue(TBEGIN));
843                     }
844                     initcount(status);
845                 }
846                 bRet = read_next_xtc(status->fio, fr->natoms, &fr->step, &fr->time, fr->box,
847                                      fr->x, &fr->prec, &bOK);
848                 fr->bPrec = (bRet && fr->prec > 0);
849                 fr->bStep = bRet;
850                 fr->bTime = bRet;
851                 fr->bX    = bRet;
852                 fr->bBox  = bRet;
853                 if (!bOK)
854                 {
855                     /* Actually the header could also be not ok,
856                        but from bOK from read_next_xtc this can't be distinguished */
857                     fr->not_ok = DATA_NOT_OK;
858                 }
859                 break;
860             case efTNG:
861                 bRet = gmx_read_next_tng_frame(status->tng, fr, NULL, 0);
862                 break;
863             case efPDB:
864                 bRet = pdb_next_x(status, gmx_fio_getfp(status->fio), fr);
865                 break;
866             case efGRO:
867                 bRet = gro_next_x_or_v(gmx_fio_getfp(status->fio), fr);
868                 break;
869             default:
870 #ifdef GMX_USE_PLUGINS
871                 bRet = read_next_vmd_frame(fr);
872 #else
873                 gmx_fatal(FARGS, "DEATH HORROR in read_next_frame ftp=%s,status=%s",
874                           ftp2ext(gmx_fio_getftp(status->fio)),
875                           gmx_fio_getname(status->fio));
876 #endif
877         }
878         fr->tf = fr->time;
879
880         if (bRet)
881         {
882             bMissingData = (((fr->flags & TRX_NEED_X) && !fr->bX) ||
883                             ((fr->flags & TRX_NEED_V) && !fr->bV) ||
884                             ((fr->flags & TRX_NEED_F) && !fr->bF));
885             bSkip = FALSE;
886             if (!bMissingData)
887             {
888                 ct = check_times2(fr->time, fr->t0, fr->bDouble);
889                 if (ct == 0 || ((fr->flags & TRX_DONT_SKIP) && ct < 0))
890                 {
891                     printcount(status, oenv, fr->time, FALSE);
892                 }
893                 else if (ct > 0)
894                 {
895                     bRet = FALSE;
896                 }
897                 else
898                 {
899                     printcount(status, oenv, fr->time, TRUE);
900                     bSkip = TRUE;
901                 }
902             }
903         }
904
905     }
906     while (bRet && (bMissingData || bSkip));
907
908     if (!bRet)
909     {
910         printlast(status, oenv, pt);
911         if (fr->not_ok)
912         {
913             printincomp(status, fr);
914         }
915     }
916
917     return bRet;
918 }
919
920 int read_first_frame(const output_env_t oenv, t_trxstatus **status,
921                      const char *fn, t_trxframe *fr, int flags)
922 {
923     t_fileio      *fio;
924     gmx_bool       bFirst, bOK;
925     int            dummy = 0;
926     int            ftp   = fn2ftp(fn);
927     gmx_int64_t   *tng_ids;
928
929     clear_trxframe(fr, TRUE);
930     fr->flags = flags;
931
932     bFirst = TRUE;
933
934     snew((*status), 1);
935
936     status_init( *status );
937     (*status)->nxframe = 1;
938     initcount(*status);
939
940     if (efTNG == ftp)
941     {
942         /* Special treatment for TNG files */
943         gmx_tng_open(fn, 'r', &(*status)->tng);
944     }
945     else
946     {
947         fio = (*status)->fio = gmx_fio_open(fn, "r");
948     }
949     switch (ftp)
950     {
951         case efTRJ:
952         case efTRR:
953             break;
954         case efCPT:
955             read_checkpoint_trxframe(fio, fr);
956             bFirst = FALSE;
957             break;
958         case efG96:
959             /* Can not rewind a compressed file, so open it twice */
960             if (!(*status)->persistent_line)
961             {
962                 /* allocate the persistent line */
963                 snew((*status)->persistent_line, STRLEN+1);
964             }
965             read_g96_conf(gmx_fio_getfp(fio), fn, fr, (*status)->persistent_line);
966             gmx_fio_close(fio);
967             clear_trxframe(fr, FALSE);
968             if (flags & (TRX_READ_X | TRX_NEED_X))
969             {
970                 snew(fr->x, fr->natoms);
971             }
972             if (flags & (TRX_READ_V | TRX_NEED_V))
973             {
974                 snew(fr->v, fr->natoms);
975             }
976             fio = (*status)->fio = gmx_fio_open(fn, "r");
977             break;
978         case efXTC:
979             if (read_first_xtc(fio, &fr->natoms, &fr->step, &fr->time, fr->box, &fr->x,
980                                &fr->prec, &bOK) == 0)
981             {
982                 assert(!bOK);
983                 fr->not_ok = DATA_NOT_OK;
984             }
985             if (fr->not_ok)
986             {
987                 fr->natoms = 0;
988                 printincomp(*status, fr);
989             }
990             else
991             {
992                 fr->bPrec = (fr->prec > 0);
993                 fr->bStep = TRUE;
994                 fr->bTime = TRUE;
995                 fr->bX    = TRUE;
996                 fr->bBox  = TRUE;
997                 printcount(*status, oenv, fr->time, FALSE);
998             }
999             bFirst = FALSE;
1000             break;
1001         case efTNG:
1002             fr->step = -1;
1003             if (!gmx_read_next_tng_frame((*status)->tng, fr, NULL, 0))
1004             {
1005                 fr->not_ok = DATA_NOT_OK;
1006                 fr->natoms = 0;
1007                 printincomp(*status, fr);
1008             }
1009             else
1010             {
1011                 printcount(*status, oenv, fr->time, FALSE);
1012             }
1013             bFirst = FALSE;
1014             break;
1015         case efPDB:
1016             pdb_first_x(*status, gmx_fio_getfp(fio), fr);
1017             if (fr->natoms)
1018             {
1019                 printcount(*status, oenv, fr->time, FALSE);
1020             }
1021             bFirst = FALSE;
1022             break;
1023         case efGRO:
1024             if (gro_first_x_or_v(gmx_fio_getfp(fio), fr))
1025             {
1026                 printcount(*status, oenv, fr->time, FALSE);
1027             }
1028             bFirst = FALSE;
1029             break;
1030         default:
1031 #ifdef GMX_USE_PLUGINS
1032             fprintf(stderr, "The file format of %s is not a known trajectory format to GROMACS.\n"
1033                     "Please make sure that the file is a trajectory!\n"
1034                     "GROMACS will now assume it to be a trajectory and will try to open it using the VMD plug-ins.\n"
1035                     "This will only work in case the VMD plugins are found and it is a trajectory format supported by VMD.\n", fn);
1036             gmx_fio_fp_close(fio); /*only close the file without removing FIO entry*/
1037             if (!read_first_vmd_frame(fn, fr))
1038             {
1039                 gmx_fatal(FARGS, "Not supported in read_first_frame: %s", fn);
1040             }
1041 #else
1042             gmx_fatal(FARGS, "Not supported in read_first_frame: %s. Please make sure that the file is a trajectory.\n"
1043                       "GROMACS is not compiled with plug-in support. Thus it cannot read non-GROMACS trajectory formats using the VMD plug-ins.\n"
1044                       "Please compile with plug-in support if you want to read non-GROMACS trajectory formats.\n", fn);
1045 #endif
1046             break;
1047     }
1048     fr->tf = fr->time;
1049
1050     /* Return FALSE if we read a frame that's past the set ending time. */
1051     if (!bFirst && (!(fr->flags & TRX_DONT_SKIP) && check_times(fr->time) > 0))
1052     {
1053         fr->t0 = fr->time;
1054         return FALSE;
1055     }
1056
1057     if (bFirst ||
1058         (!(fr->flags & TRX_DONT_SKIP) && check_times(fr->time) < 0))
1059     {
1060         /* Read a frame when no frame was read or the first was skipped */
1061         if (!read_next_frame(oenv, *status, fr))
1062         {
1063             return FALSE;
1064         }
1065     }
1066     fr->t0 = fr->time;
1067
1068     return (fr->natoms > 0);
1069 }
1070
1071 /***** C O O R D I N A T E   S T U F F *****/
1072
1073 int read_first_x(const output_env_t oenv, t_trxstatus **status, const char *fn,
1074                  real *t, rvec **x, matrix box)
1075 {
1076     t_trxframe fr;
1077
1078     read_first_frame(oenv, status, fn, &fr, TRX_NEED_X);
1079
1080     snew((*status)->xframe, 1);
1081     (*status)->nxframe   = 1;
1082     (*(*status)->xframe) = fr;
1083     *t                   = (*status)->xframe->time;
1084     *x                   = (*status)->xframe->x;
1085     copy_mat((*status)->xframe->box, box);
1086
1087     return (*status)->xframe->natoms;
1088 }
1089
1090 gmx_bool read_next_x(const output_env_t oenv, t_trxstatus *status, real *t,
1091                      rvec x[], matrix box)
1092 {
1093     gmx_bool bRet;
1094
1095     status->xframe->x = x;
1096     /*xframe[status].x = x;*/
1097     bRet = read_next_frame(oenv, status, status->xframe);
1098     *t   = status->xframe->time;
1099     copy_mat(status->xframe->box, box);
1100
1101     return bRet;
1102 }
1103
1104 void close_trj(t_trxstatus *status)
1105 {
1106     gmx_tng_close(&status->tng);
1107     if (status->fio)
1108     {
1109         gmx_fio_close(status->fio);
1110     }
1111
1112     /* The memory in status->xframe is lost here,
1113      * but the read_first_x/read_next_x functions are deprecated anyhow.
1114      * read_first_frame/read_next_frame and close_trx should be used.
1115      */
1116     sfree(status);
1117 }
1118
1119 void rewind_trj(t_trxstatus *status)
1120 {
1121     initcount(status);
1122
1123     gmx_fio_rewind(status->fio);
1124 }
1125
1126 /***** T O P O L O G Y   S T U F F ******/
1127
1128 t_topology *read_top(const char *fn, int *ePBC)
1129 {
1130     int         epbc, natoms;
1131     t_topology *top;
1132
1133     snew(top, 1);
1134     epbc = read_tpx_top(fn, NULL, NULL, &natoms, NULL, NULL, NULL, top);
1135     if (ePBC)
1136     {
1137         *ePBC = epbc;
1138     }
1139
1140     return top;
1141 }