Code beautification with uncrustify
[alexxy/gromacs.git] / src / gromacs / mdlib / stat.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * GROningen Mixture of Alchemy and Childrens' Stories
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <string.h>
41 #include <stdio.h>
42 #include "typedefs.h"
43 #include "sysstuff.h"
44 #include "gmx_fatal.h"
45 #include "network.h"
46 #include "txtdump.h"
47 #include "names.h"
48 #include "physics.h"
49 #include "vec.h"
50 #include "maths.h"
51 #include "mvdata.h"
52 #include "main.h"
53 #include "force.h"
54 #include "vcm.h"
55 #include "smalloc.h"
56 #include "futil.h"
57 #include "network.h"
58 #include "rbin.h"
59 #include "tgroup.h"
60 #include "xtcio.h"
61 #include "gmxfio.h"
62 #include "trnio.h"
63 #include "statutil.h"
64 #include "domdec.h"
65 #include "partdec.h"
66 #include "constr.h"
67 #include "checkpoint.h"
68 #include "xvgr.h"
69 #include "md_support.h"
70 #include "mdrun.h"
71 #include "sim_util.h"
72
73 typedef struct gmx_global_stat
74 {
75     t_bin *rb;
76     int   *itc0;
77     int   *itc1;
78 } t_gmx_global_stat;
79
80 gmx_global_stat_t global_stat_init(t_inputrec *ir)
81 {
82     gmx_global_stat_t gs;
83
84     snew(gs, 1);
85
86     gs->rb = mk_bin();
87     snew(gs->itc0, ir->opts.ngtc);
88     snew(gs->itc1, ir->opts.ngtc);
89
90     return gs;
91 }
92
93 void global_stat_destroy(gmx_global_stat_t gs)
94 {
95     destroy_bin(gs->rb);
96     sfree(gs->itc0);
97     sfree(gs->itc1);
98     sfree(gs);
99 }
100
101 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
102                             gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner)
103 {
104     int i, to, from;
105
106     from = 0;
107     to   = 0;
108     for (i = 0; i < F_NRE; i++)
109     {
110         if (bToBuffer)
111         {
112             from = i;
113         }
114         else
115         {
116             to = i;
117         }
118         switch (i)
119         {
120             case F_EKIN:
121             case F_TEMP:
122             case F_DKDL:
123                 if (bTemp)
124                 {
125                     ato[to++] = afrom[from++];
126                 }
127                 break;
128             case F_PRES:
129             case F_PDISPCORR:
130                 if (bPres)
131                 {
132                     ato[to++] = afrom[from++];
133                 }
134                 break;
135             default:
136                 if (bEner)
137                 {
138                     ato[to++] = afrom[from++];
139                 }
140                 break;
141         }
142     }
143
144     return to;
145 }
146
147 void global_stat(FILE *fplog, gmx_global_stat_t gs,
148                  t_commrec *cr, gmx_enerdata_t *enerd,
149                  tensor fvir, tensor svir, rvec mu_tot,
150                  t_inputrec *inputrec,
151                  gmx_ekindata_t *ekind, gmx_constr_t constr,
152                  t_vcm *vcm,
153                  int nsig, real *sig,
154                  gmx_mtop_t *top_global, t_state *state_local,
155                  gmx_bool bSumEkinhOld, int flags)
156 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
157 {
158     t_bin     *rb;
159     int       *itc0, *itc1;
160     int        ie    = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
161     int        idedl = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
162     int        isig  = -1;
163     int        icj   = -1, ici = -1, icx = -1;
164     int        inn[egNR];
165     real       copyenerd[F_NRE];
166     int        nener, j;
167     real      *rmsd_data = NULL;
168     double     nb;
169     gmx_bool   bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bFirstIterate, bReadEkin;
170
171     bVV           = EI_VV(inputrec->eI);
172     bTemp         = flags & CGLO_TEMPERATURE;
173     bEner         = flags & CGLO_ENERGY;
174     bPres         = (flags & CGLO_PRESSURE);
175     bConstrVir    = (flags & CGLO_CONSTRAINT);
176     bFirstIterate = (flags & CGLO_FIRSTITERATE);
177     bEkinAveVel   = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
178     bReadEkin     = (flags & CGLO_READEKIN);
179
180     rb   = gs->rb;
181     itc0 = gs->itc0;
182     itc1 = gs->itc1;
183
184
185     reset_bin(rb);
186     /* This routine copies all the data to be summed to one big buffer
187      * using the t_bin struct.
188      */
189
190     /* First, we neeed to identify which enerd->term should be
191        communicated.  Temperature and pressure terms should only be
192        communicated and summed when they need to be, to avoid repeating
193        the sums and overcounting. */
194
195     nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
196
197     /* First, the data that needs to be communicated with velocity verlet every time
198        This is just the constraint virial.*/
199     if (bConstrVir)
200     {
201         isv = add_binr(rb, DIM*DIM, svir[0]);
202         where();
203     }
204
205 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
206     if (bTemp || !bVV)
207     {
208         if (ekind)
209         {
210             for (j = 0; (j < inputrec->opts.ngtc); j++)
211             {
212                 if (bSumEkinhOld)
213                 {
214                     itc0[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
215                 }
216                 if (bEkinAveVel && !bReadEkin)
217                 {
218                     itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinf[0]);
219                 }
220                 else if (!bReadEkin)
221                 {
222                     itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh[0]);
223                 }
224             }
225             /* these probably need to be put into one of these categories */
226             where();
227             idedl = add_binr(rb, 1, &(ekind->dekindl));
228             where();
229             ica   = add_binr(rb, 1, &(ekind->cosacc.mvcos));
230             where();
231         }
232     }
233     where();
234
235     if ((bPres || !bVV) && bFirstIterate)
236     {
237         ifv = add_binr(rb, DIM*DIM, fvir[0]);
238     }
239
240
241     if (bEner)
242     {
243         where();
244         if (bFirstIterate)
245         {
246             ie  = add_binr(rb, nener, copyenerd);
247         }
248         where();
249         if (constr)
250         {
251             rmsd_data = constr_rmsd_data(constr);
252             if (rmsd_data)
253             {
254                 irmsd = add_binr(rb, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
255             }
256         }
257         if (!NEED_MUTOT(*inputrec))
258         {
259             imu = add_binr(rb, DIM, mu_tot);
260             where();
261         }
262
263         if (bFirstIterate)
264         {
265             for (j = 0; (j < egNR); j++)
266             {
267                 inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
268             }
269             where();
270             if (inputrec->efep != efepNO)
271             {
272                 idvdll  = add_bind(rb, efptNR, enerd->dvdl_lin);
273                 idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
274                 if (enerd->n_lambda > 0)
275                 {
276                     iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
277                 }
278             }
279         }
280     }
281
282     if (vcm)
283     {
284         icm   = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
285         where();
286         imass = add_binr(rb, vcm->nr, vcm->group_mass);
287         where();
288         if (vcm->mode == ecmANGULAR)
289         {
290             icj   = add_binr(rb, DIM*vcm->nr, vcm->group_j[0]);
291             where();
292             icx   = add_binr(rb, DIM*vcm->nr, vcm->group_x[0]);
293             where();
294             ici   = add_binr(rb, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
295             where();
296         }
297     }
298
299     if (DOMAINDECOMP(cr))
300     {
301         nb  = cr->dd->nbonded_local;
302         inb = add_bind(rb, 1, &nb);
303     }
304     where();
305     if (nsig > 0)
306     {
307         isig = add_binr(rb, nsig, sig);
308     }
309
310     /* Global sum it all */
311     if (debug)
312     {
313         fprintf(debug, "Summing %d energies\n", rb->maxreal);
314     }
315     sum_bin(rb, cr);
316     where();
317
318     /* Extract all the data locally */
319
320     if (bConstrVir)
321     {
322         extract_binr(rb, isv, DIM*DIM, svir[0]);
323     }
324
325     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
326     if (bTemp || !bVV)
327     {
328         if (ekind)
329         {
330             for (j = 0; (j < inputrec->opts.ngtc); j++)
331             {
332                 if (bSumEkinhOld)
333                 {
334                     extract_binr(rb, itc0[j], DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
335                 }
336                 if (bEkinAveVel && !bReadEkin)
337                 {
338                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinf[0]);
339                 }
340                 else if (!bReadEkin)
341                 {
342                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
343                 }
344             }
345             extract_binr(rb, idedl, 1, &(ekind->dekindl));
346             extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
347             where();
348         }
349     }
350     if ((bPres || !bVV) && bFirstIterate)
351     {
352         extract_binr(rb, ifv, DIM*DIM, fvir[0]);
353     }
354
355     if (bEner)
356     {
357         if (bFirstIterate)
358         {
359             extract_binr(rb, ie, nener, copyenerd);
360             if (rmsd_data)
361             {
362                 extract_binr(rb, irmsd, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
363             }
364             if (!NEED_MUTOT(*inputrec))
365             {
366                 extract_binr(rb, imu, DIM, mu_tot);
367             }
368
369             for (j = 0; (j < egNR); j++)
370             {
371                 extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
372             }
373             if (inputrec->efep != efepNO)
374             {
375                 extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
376                 extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
377                 if (enerd->n_lambda > 0)
378                 {
379                     extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
380                 }
381             }
382             if (DOMAINDECOMP(cr))
383             {
384                 extract_bind(rb, inb, 1, &nb);
385                 if ((int)(nb + 0.5) != cr->dd->nbonded_global)
386                 {
387                     dd_print_missing_interactions(fplog, cr, (int)(nb + 0.5), top_global, state_local);
388                 }
389             }
390             where();
391
392             filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
393         }
394     }
395
396     if (vcm)
397     {
398         extract_binr(rb, icm, DIM*vcm->nr, vcm->group_p[0]);
399         where();
400         extract_binr(rb, imass, vcm->nr, vcm->group_mass);
401         where();
402         if (vcm->mode == ecmANGULAR)
403         {
404             extract_binr(rb, icj, DIM*vcm->nr, vcm->group_j[0]);
405             where();
406             extract_binr(rb, icx, DIM*vcm->nr, vcm->group_x[0]);
407             where();
408             extract_binr(rb, ici, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
409             where();
410         }
411     }
412
413     if (nsig > 0)
414     {
415         extract_binr(rb, isig, nsig, sig);
416     }
417     where();
418 }
419
420 int do_per_step(gmx_large_int_t step, gmx_large_int_t nstep)
421 {
422     if (nstep != 0)
423     {
424         return ((step % nstep) == 0);
425     }
426     else
427     {
428         return 0;
429     }
430 }
431
432 static void moveit(t_commrec *cr,
433                    int left, int right, const char *s, rvec xx[])
434 {
435     if (!xx)
436     {
437         return;
438     }
439
440     move_rvecs(cr, FALSE, FALSE, left, right,
441                xx, NULL, (cr->nnodes-cr->npmenodes)-1, NULL);
442 }
443
444 gmx_mdoutf_t *init_mdoutf(int nfile, const t_filenm fnm[], int mdrun_flags,
445                           const t_commrec *cr, const t_inputrec *ir,
446                           const output_env_t oenv)
447 {
448     gmx_mdoutf_t *of;
449     char          filemode[3];
450     gmx_bool      bAppendFiles;
451
452     snew(of, 1);
453
454     of->fp_trn   = NULL;
455     of->fp_ene   = NULL;
456     of->fp_xtc   = NULL;
457     of->fp_dhdl  = NULL;
458     of->fp_field = NULL;
459
460     of->eIntegrator     = ir->eI;
461     of->bExpanded       = ir->bExpanded;
462     of->elamstats       = ir->expandedvals->elamstats;
463     of->simulation_part = ir->simulation_part;
464
465     if (MASTER(cr))
466     {
467         bAppendFiles = (mdrun_flags & MD_APPENDFILES);
468
469         of->bKeepAndNumCPT = (mdrun_flags & MD_KEEPANDNUMCPT);
470
471         sprintf(filemode, bAppendFiles ? "a+" : "w+");
472
473         if ((EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
474 #ifndef GMX_FAHCORE
475             &&
476             !(EI_DYNAMICS(ir->eI) &&
477               ir->nstxout == 0 &&
478               ir->nstvout == 0 &&
479               ir->nstfout == 0)
480 #endif
481             )
482         {
483             of->fp_trn = open_trn(ftp2fn(efTRN, nfile, fnm), filemode);
484         }
485         if (EI_DYNAMICS(ir->eI) &&
486             ir->nstxtcout > 0)
487         {
488             of->fp_xtc   = open_xtc(ftp2fn(efXTC, nfile, fnm), filemode);
489             of->xtc_prec = ir->xtcprec;
490         }
491         if (EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
492         {
493             of->fp_ene = open_enx(ftp2fn(efEDR, nfile, fnm), filemode);
494         }
495         of->fn_cpt = opt2fn("-cpo", nfile, fnm);
496
497         if ((ir->efep != efepNO || ir->bSimTemp) && ir->fepvals->nstdhdl > 0 &&
498             (ir->fepvals->separate_dhdl_file == esepdhdlfileYES ) &&
499             EI_DYNAMICS(ir->eI))
500         {
501             if (bAppendFiles)
502             {
503                 of->fp_dhdl = gmx_fio_fopen(opt2fn("-dhdl", nfile, fnm), filemode);
504             }
505             else
506             {
507                 of->fp_dhdl = open_dhdl(opt2fn("-dhdl", nfile, fnm), ir, oenv);
508             }
509         }
510
511         if (opt2bSet("-field", nfile, fnm) &&
512             (ir->ex[XX].n || ir->ex[YY].n || ir->ex[ZZ].n))
513         {
514             if (bAppendFiles)
515             {
516                 of->fp_dhdl = gmx_fio_fopen(opt2fn("-field", nfile, fnm),
517                                             filemode);
518             }
519             else
520             {
521                 of->fp_field = xvgropen(opt2fn("-field", nfile, fnm),
522                                         "Applied electric field", "Time (ps)",
523                                         "E (V/nm)", oenv);
524             }
525         }
526     }
527
528     return of;
529 }
530
531 void done_mdoutf(gmx_mdoutf_t *of)
532 {
533     if (of->fp_ene != NULL)
534     {
535         close_enx(of->fp_ene);
536     }
537     if (of->fp_xtc)
538     {
539         close_xtc(of->fp_xtc);
540     }
541     if (of->fp_trn)
542     {
543         close_trn(of->fp_trn);
544     }
545     if (of->fp_dhdl != NULL)
546     {
547         gmx_fio_fclose(of->fp_dhdl);
548     }
549     if (of->fp_field != NULL)
550     {
551         gmx_fio_fclose(of->fp_field);
552     }
553
554     sfree(of);
555 }
556
557 void write_traj(FILE *fplog, t_commrec *cr,
558                 gmx_mdoutf_t *of,
559                 int mdof_flags,
560                 gmx_mtop_t *top_global,
561                 gmx_large_int_t step, double t,
562                 t_state *state_local, t_state *state_global,
563                 rvec *f_local, rvec *f_global,
564                 int *n_xtc, rvec **x_xtc)
565 {
566     int           i, j;
567     gmx_groups_t *groups;
568     rvec         *xxtc;
569     rvec         *local_v;
570     rvec         *global_v;
571
572 #define MX(xvf) moveit(cr, GMX_LEFT, GMX_RIGHT,#xvf, xvf)
573
574     /* MRS -- defining these variables is to manage the difference
575      * between half step and full step velocities, but there must be a better way . . . */
576
577     local_v  = state_local->v;
578     global_v = state_global->v;
579
580     if (DOMAINDECOMP(cr))
581     {
582         if (mdof_flags & MDOF_CPT)
583         {
584             dd_collect_state(cr->dd, state_local, state_global);
585         }
586         else
587         {
588             if (mdof_flags & (MDOF_X | MDOF_XTC))
589             {
590                 dd_collect_vec(cr->dd, state_local, state_local->x,
591                                state_global->x);
592             }
593             if (mdof_flags & MDOF_V)
594             {
595                 dd_collect_vec(cr->dd, state_local, local_v,
596                                global_v);
597             }
598         }
599         if (mdof_flags & MDOF_F)
600         {
601             dd_collect_vec(cr->dd, state_local, f_local, f_global);
602         }
603     }
604     else
605     {
606         if (mdof_flags & MDOF_CPT)
607         {
608             /* All pointers in state_local are equal to state_global,
609              * but we need to copy the non-pointer entries.
610              */
611             state_global->lambda = state_local->lambda;
612             state_global->veta   = state_local->veta;
613             state_global->vol0   = state_local->vol0;
614             copy_mat(state_local->box, state_global->box);
615             copy_mat(state_local->boxv, state_global->boxv);
616             copy_mat(state_local->svir_prev, state_global->svir_prev);
617             copy_mat(state_local->fvir_prev, state_global->fvir_prev);
618             copy_mat(state_local->pres_prev, state_global->pres_prev);
619         }
620         if (cr->nnodes > 1)
621         {
622             /* Particle decomposition, collect the data on the master node */
623             if (mdof_flags & MDOF_CPT)
624             {
625                 if (state_local->flags & (1<<estX))
626                 {
627                     MX(state_global->x);
628                 }
629                 if (state_local->flags & (1<<estV))
630                 {
631                     MX(state_global->v);
632                 }
633                 if (state_local->flags & (1<<estSDX))
634                 {
635                     MX(state_global->sd_X);
636                 }
637                 if (state_global->nrngi > 1)
638                 {
639                     if (state_local->flags & (1<<estLD_RNG))
640                     {
641 #ifdef GMX_MPI
642                         MPI_Gather(state_local->ld_rng,
643                                    state_local->nrng*sizeof(state_local->ld_rng[0]), MPI_BYTE,
644                                    state_global->ld_rng,
645                                    state_local->nrng*sizeof(state_local->ld_rng[0]), MPI_BYTE,
646                                    MASTERRANK(cr), cr->mpi_comm_mygroup);
647 #endif
648                     }
649                     if (state_local->flags & (1<<estLD_RNGI))
650                     {
651 #ifdef GMX_MPI
652                         MPI_Gather(state_local->ld_rngi,
653                                    sizeof(state_local->ld_rngi[0]), MPI_BYTE,
654                                    state_global->ld_rngi,
655                                    sizeof(state_local->ld_rngi[0]), MPI_BYTE,
656                                    MASTERRANK(cr), cr->mpi_comm_mygroup);
657 #endif
658                     }
659                 }
660             }
661             else
662             {
663                 if (mdof_flags & (MDOF_X | MDOF_XTC))
664                 {
665                     MX(state_global->x);
666                 }
667                 if (mdof_flags & MDOF_V)
668                 {
669                     MX(global_v);
670                 }
671             }
672             if (mdof_flags & MDOF_F)
673             {
674                 MX(f_global);
675             }
676         }
677     }
678
679     if (MASTER(cr))
680     {
681         if (mdof_flags & MDOF_CPT)
682         {
683             write_checkpoint(of->fn_cpt, of->bKeepAndNumCPT,
684                              fplog, cr, of->eIntegrator, of->simulation_part,
685                              of->bExpanded, of->elamstats, step, t, state_global);
686         }
687
688         if (mdof_flags & (MDOF_X | MDOF_V | MDOF_F))
689         {
690             fwrite_trn(of->fp_trn, step, t, state_local->lambda[efptFEP],
691                        state_local->box, top_global->natoms,
692                        (mdof_flags & MDOF_X) ? state_global->x : NULL,
693                        (mdof_flags & MDOF_V) ? global_v : NULL,
694                        (mdof_flags & MDOF_F) ? f_global : NULL);
695             if (gmx_fio_flush(of->fp_trn) != 0)
696             {
697                 gmx_file("Cannot write trajectory; maybe you are out of disk space?");
698             }
699             gmx_fio_check_file_position(of->fp_trn);
700         }
701         if (mdof_flags & MDOF_XTC)
702         {
703             groups = &top_global->groups;
704             if (*n_xtc == -1)
705             {
706                 *n_xtc = 0;
707                 for (i = 0; (i < top_global->natoms); i++)
708                 {
709                     if (ggrpnr(groups, egcXTC, i) == 0)
710                     {
711                         (*n_xtc)++;
712                     }
713                 }
714                 if (*n_xtc != top_global->natoms)
715                 {
716                     snew(*x_xtc, *n_xtc);
717                 }
718             }
719             if (*n_xtc == top_global->natoms)
720             {
721                 xxtc = state_global->x;
722             }
723             else
724             {
725                 xxtc = *x_xtc;
726                 j    = 0;
727                 for (i = 0; (i < top_global->natoms); i++)
728                 {
729                     if (ggrpnr(groups, egcXTC, i) == 0)
730                     {
731                         copy_rvec(state_global->x[i], xxtc[j++]);
732                     }
733                 }
734             }
735             if (write_xtc(of->fp_xtc, *n_xtc, step, t,
736                           state_local->box, xxtc, of->xtc_prec) == 0)
737             {
738                 gmx_fatal(FARGS, "XTC error - maybe you are out of disk space?");
739             }
740             gmx_fio_check_file_position(of->fp_xtc);
741         }
742     }
743 }