Merge branch release-2016
[alexxy/gromacs.git] / src / gromacs / tools / dump.cpp
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-2013, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017, 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 "config.h"
40
41 #include <cassert>
42 #include <cmath>
43 #include <cstdio>
44 #include <cstring>
45
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/checkpoint.h"
48 #include "gromacs/fileio/enxio.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/mtxio.h"
51 #include "gromacs/fileio/tngio.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/fileio/trrio.h"
54 #include "gromacs/fileio/xtcio.h"
55 #include "gromacs/gmxpreprocess/gmxcpp.h"
56 #include "gromacs/linearalgebra/sparsematrix.h"
57 #include "gromacs/math/vecdump.h"
58 #include "gromacs/mdrunutility/mdmodules.h"
59 #include "gromacs/mdtypes/forcerec.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/mdtypes/state.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/topology/topology.h"
65 #include "gromacs/trajectory/trajectoryframe.h"
66 #include "gromacs/utility/arraysize.h"
67 #include "gromacs/utility/basedefinitions.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/futil.h"
70 #include "gromacs/utility/smalloc.h"
71 #include "gromacs/utility/txtdump.h"
72
73 static void list_tpx(const char *fn,
74                      gmx_bool    bShowNumbers,
75                      gmx_bool    bShowParameters,
76                      const char *mdpfn,
77                      gmx_bool    bSysTop,
78                      gmx_bool    bOriginalInputrec)
79 {
80     FILE         *gp;
81     int           indent, i, j, **gcount, atot;
82     t_state       state;
83     t_tpxheader   tpx;
84     gmx_mtop_t    mtop;
85     gmx_groups_t *groups;
86     t_topology    top;
87
88     read_tpxheader(fn, &tpx, TRUE);
89     t_inputrec     ir;
90     read_tpx_state(fn,
91                    tpx.bIr ? &ir : nullptr,
92                    &state,
93                    tpx.bTop ? &mtop : nullptr);
94     if (tpx.bIr && !bOriginalInputrec)
95     {
96         gmx::MDModules mdModules;
97         mdModules.adjustInputrecBasedOnModules(&ir);
98     }
99
100     if (mdpfn && tpx.bIr)
101     {
102         gp = gmx_fio_fopen(mdpfn, "w");
103         pr_inputrec(gp, 0, nullptr, &ir, TRUE);
104         gmx_fio_fclose(gp);
105     }
106
107     if (!mdpfn)
108     {
109         if (bSysTop)
110         {
111             top = gmx_mtop_t_to_t_topology(&mtop, false);
112         }
113
114         if (available(stdout, &tpx, 0, fn))
115         {
116             indent = 0;
117             pr_title(stdout, indent, fn);
118             pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &ir : nullptr, FALSE);
119
120             pr_tpxheader(stdout, indent, "header", &(tpx));
121
122             if (!bSysTop)
123             {
124                 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers, bShowParameters);
125             }
126             else
127             {
128                 pr_top(stdout, indent, "topology", &(top), bShowNumbers, bShowParameters);
129             }
130
131             pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : nullptr, DIM);
132             pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : nullptr, DIM);
133             pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : nullptr, DIM);
134             pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : nullptr, DIM);
135             pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : nullptr, DIM);
136             pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : nullptr, DIM);
137             /* leave nosehoover_xi in for now to match the tpr version */
138             pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi.data(), state.ngtc);
139             /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
140             /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
141             pr_rvecs(stdout, indent, "x", tpx.bX ? as_rvec_array(state.x.data()) : nullptr, state.natoms);
142             pr_rvecs(stdout, indent, "v", tpx.bV ? as_rvec_array(state.v.data()) : nullptr, state.natoms);
143         }
144
145         groups = &mtop.groups;
146
147         snew(gcount, egcNR);
148         for (i = 0; (i < egcNR); i++)
149         {
150             snew(gcount[i], groups->grps[i].nr);
151         }
152
153         for (i = 0; (i < mtop.natoms); i++)
154         {
155             for (j = 0; (j < egcNR); j++)
156             {
157                 gcount[j][ggrpnr(groups, j, i)]++;
158             }
159         }
160         printf("Group statistics\n");
161         for (i = 0; (i < egcNR); i++)
162         {
163             atot = 0;
164             printf("%-12s: ", gtypes[i]);
165             for (j = 0; (j < groups->grps[i].nr); j++)
166             {
167                 printf("  %5d", gcount[i][j]);
168                 atot += gcount[i][j];
169             }
170             printf("  (total %d atoms)\n", atot);
171             sfree(gcount[i]);
172         }
173         sfree(gcount);
174     }
175 }
176
177 static void list_top(const char *fn)
178 {
179     int       status, done;
180 #define BUFLEN 256
181     char      buf[BUFLEN];
182     gmx_cpp_t handle;
183     char     *cppopts[] = { nullptr };
184
185     status = cpp_open_file(fn, &handle, cppopts);
186     if (status != 0)
187     {
188         gmx_fatal(FARGS, cpp_error(&handle, status));
189     }
190     do
191     {
192         status = cpp_read_line(&handle, BUFLEN, buf);
193         done   = (status == eCPP_EOF);
194         if (!done)
195         {
196             if (status != eCPP_OK)
197             {
198                 gmx_fatal(FARGS, cpp_error(&handle, status));
199             }
200             else
201             {
202                 printf("%s\n", buf);
203             }
204         }
205     }
206     while (!done);
207     status = cpp_close_file(&handle);
208     if (status != eCPP_OK)
209     {
210         gmx_fatal(FARGS, cpp_error(&handle, status));
211     }
212 }
213
214 static void list_trr(const char *fn)
215 {
216     t_fileio         *fpread;
217     int               nframe, indent;
218     char              buf[256];
219     rvec             *x, *v, *f;
220     matrix            box;
221     gmx_trr_header_t  trrheader;
222     gmx_bool          bOK;
223
224     fpread  = gmx_trr_open(fn, "r");
225
226     nframe = 0;
227     while (gmx_trr_read_frame_header(fpread, &trrheader, &bOK))
228     {
229         snew(x, trrheader.natoms);
230         snew(v, trrheader.natoms);
231         snew(f, trrheader.natoms);
232         if (gmx_trr_read_frame_data(fpread, &trrheader,
233                                     trrheader.box_size ? box : nullptr,
234                                     trrheader.x_size   ? x : nullptr,
235                                     trrheader.v_size   ? v : nullptr,
236                                     trrheader.f_size   ? f : nullptr))
237         {
238             sprintf(buf, "%s frame %d", fn, nframe);
239             indent = 0;
240             indent = pr_title(stdout, indent, buf);
241             pr_indent(stdout, indent);
242             fprintf(stdout, "natoms=%10d  step=%10" GMX_PRId64 "  time=%12.7e  lambda=%10g\n",
243                     trrheader.natoms, trrheader.step, trrheader.t, trrheader.lambda);
244             if (trrheader.box_size)
245             {
246                 pr_rvecs(stdout, indent, "box", box, DIM);
247             }
248             if (trrheader.x_size)
249             {
250                 pr_rvecs(stdout, indent, "x", x, trrheader.natoms);
251             }
252             if (trrheader.v_size)
253             {
254                 pr_rvecs(stdout, indent, "v", v, trrheader.natoms);
255             }
256             if (trrheader.f_size)
257             {
258                 pr_rvecs(stdout, indent, "f", f, trrheader.natoms);
259             }
260         }
261         else
262         {
263             fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
264                     nframe, trrheader.t);
265         }
266
267         sfree(x);
268         sfree(v);
269         sfree(f);
270         nframe++;
271     }
272     if (!bOK)
273     {
274         fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
275                 nframe, trrheader.t);
276     }
277     gmx_trr_close(fpread);
278 }
279
280 void list_xtc(const char *fn)
281 {
282     t_fileio   *xd;
283     int         indent;
284     char        buf[256];
285     rvec       *x;
286     matrix      box;
287     int         nframe, natoms;
288     gmx_int64_t step;
289     real        prec, time;
290     gmx_bool    bOK;
291
292     xd = open_xtc(fn, "r");
293     read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
294
295     nframe = 0;
296     do
297     {
298         sprintf(buf, "%s frame %d", fn, nframe);
299         indent = 0;
300         indent = pr_title(stdout, indent, buf);
301         pr_indent(stdout, indent);
302         fprintf(stdout, "natoms=%10d  step=%10" GMX_PRId64 "  time=%12.7e  prec=%10g\n",
303                 natoms, step, time, prec);
304         pr_rvecs(stdout, indent, "box", box, DIM);
305         pr_rvecs(stdout, indent, "x", x, natoms);
306         nframe++;
307     }
308     while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK));
309     if (!bOK)
310     {
311         fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
312     }
313     sfree(x);
314     close_xtc(xd);
315 }
316
317 /*! \brief Callback used by list_tng_for_gmx_dump. */
318 static void list_tng_inner(const char *fn,
319                            gmx_bool    bFirstFrame,
320                            real       *values,
321                            gmx_int64_t step,
322                            double      frame_time,
323                            gmx_int64_t n_values_per_frame,
324                            gmx_int64_t n_atoms,
325                            real        prec,
326                            gmx_int64_t nframe,
327                            char       *block_name)
328 {
329     char                 buf[256];
330     int                  indent = 0;
331
332     if (bFirstFrame)
333     {
334         sprintf(buf, "%s frame %" GMX_PRId64, fn, nframe);
335         indent = 0;
336         indent = pr_title(stdout, indent, buf);
337         pr_indent(stdout, indent);
338         fprintf(stdout, "natoms=%10" GMX_PRId64 "  step=%10" GMX_PRId64 "  time=%12.7e",
339                 n_atoms, step, frame_time);
340         if (prec > 0)
341         {
342             fprintf(stdout, "  prec=%10g", prec);
343         }
344         fprintf(stdout, "\n");
345     }
346     pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
347 }
348
349 static void list_tng(const char gmx_unused *fn)
350 {
351 #ifdef GMX_USE_TNG
352     tng_trajectory_t     tng;
353     gmx_int64_t          nframe = 0;
354     gmx_int64_t          i, *block_ids = nullptr, step, ndatablocks;
355     gmx_bool             bOK;
356     real                *values = nullptr;
357
358     gmx_tng_open(fn, 'r', &tng);
359     gmx_print_tng_molecule_system(tng, stdout);
360
361     bOK    = gmx_get_tng_data_block_types_of_next_frame(tng, -1,
362                                                         0,
363                                                         nullptr,
364                                                         &step, &ndatablocks,
365                                                         &block_ids);
366     do
367     {
368         for (i = 0; i < ndatablocks; i++)
369         {
370             double               frame_time;
371             real                 prec;
372             gmx_int64_t          n_values_per_frame, n_atoms;
373             char                 block_name[STRLEN];
374
375             gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values,
376                                                       &step, &frame_time,
377                                                       &n_values_per_frame, &n_atoms,
378                                                       &prec,
379                                                       block_name, STRLEN, &bOK);
380             if (!bOK)
381             {
382                 /* Can't write any output because we don't know what
383                    arrays are valid. */
384                 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
385             }
386             else
387             {
388                 list_tng_inner(fn, (0 == i), values, step, frame_time,
389                                n_values_per_frame, n_atoms, prec, nframe, block_name);
390             }
391         }
392         nframe++;
393     }
394     while (gmx_get_tng_data_block_types_of_next_frame(tng, step,
395                                                       0,
396                                                       nullptr,
397                                                       &step,
398                                                       &ndatablocks,
399                                                       &block_ids));
400
401     if (block_ids)
402     {
403         sfree(block_ids);
404     }
405     sfree(values);
406     gmx_tng_close(&tng);
407 #endif
408 }
409
410 void list_trx(const char *fn)
411 {
412     switch (fn2ftp(fn))
413     {
414         case efXTC:
415             list_xtc(fn);
416             break;
417         case efTRR:
418             list_trr(fn);
419             break;
420         case efTNG:
421             list_tng(fn);
422             break;
423         default:
424             fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
425                     fn, fn);
426     }
427 }
428
429 void list_ene(const char *fn)
430 {
431     ener_file_t    in;
432     gmx_bool       bCont;
433     gmx_enxnm_t   *enm = nullptr;
434     t_enxframe    *fr;
435     int            i, j, nre, b;
436     char           buf[22];
437
438     printf("gmx dump: %s\n", fn);
439     in = open_enx(fn, "r");
440     do_enxnms(in, &nre, &enm);
441     assert(enm);
442
443     printf("energy components:\n");
444     for (i = 0; (i < nre); i++)
445     {
446         printf("%5d  %-24s (%s)\n", i, enm[i].name, enm[i].unit);
447     }
448
449     snew(fr, 1);
450     do
451     {
452         bCont = do_enx(in, fr);
453
454         if (bCont)
455         {
456             printf("\n%24s  %12.5e  %12s  %12s\n", "time:",
457                    fr->t, "step:", gmx_step_str(fr->step, buf));
458             printf("%24s  %12s  %12s  %12s\n",
459                    "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
460             printf("%24s  %12.5e  %12s  %12s\n",
461                    "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
462             if (fr->nre == nre)
463             {
464                 printf("%24s  %12s  %12s  %12s\n",
465                        "Component", "Energy", "Av. Energy", "Sum Energy");
466                 if (fr->nsum > 0)
467                 {
468                     for (i = 0; (i < nre); i++)
469                     {
470                         printf("%24s  %12.5e  %12.5e  %12.5e\n",
471                                enm[i].name, fr->ener[i].e, fr->ener[i].eav,
472                                fr->ener[i].esum);
473                     }
474                 }
475                 else
476                 {
477                     for (i = 0; (i < nre); i++)
478                     {
479                         printf("%24s  %12.5e\n",
480                                enm[i].name, fr->ener[i].e);
481                     }
482                 }
483             }
484             for (b = 0; b < fr->nblock; b++)
485             {
486                 const char *typestr = "";
487
488                 t_enxblock *eb = &(fr->block[b]);
489                 printf("Block data %2d (%3d subblocks, id=%d)\n",
490                        b, eb->nsub, eb->id);
491
492                 if (eb->id < enxNR)
493                 {
494                     typestr = enx_block_id_name[eb->id];
495                 }
496                 printf("  id='%s'\n", typestr);
497                 for (i = 0; i < eb->nsub; i++)
498                 {
499                     t_enxsubblock *sb = &(eb->sub[i]);
500                     printf("  Sub block %3d (%5d elems, type=%s) values:\n",
501                            i, sb->nr, xdr_datatype_names[sb->type]);
502
503                     switch (sb->type)
504                     {
505                         case xdr_datatype_float:
506                             for (j = 0; j < sb->nr; j++)
507                             {
508                                 printf("%14d   %8.4f\n", j, sb->fval[j]);
509                             }
510                             break;
511                         case xdr_datatype_double:
512                             for (j = 0; j < sb->nr; j++)
513                             {
514                                 printf("%14d   %10.6f\n", j, sb->dval[j]);
515                             }
516                             break;
517                         case xdr_datatype_int:
518                             for (j = 0; j < sb->nr; j++)
519                             {
520                                 printf("%14d %10d\n", j, sb->ival[j]);
521                             }
522                             break;
523                         case xdr_datatype_int64:
524                             for (j = 0; j < sb->nr; j++)
525                             {
526                                 printf("%14d %s\n",
527                                        j, gmx_step_str(sb->lval[j], buf));
528                             }
529                             break;
530                         case xdr_datatype_char:
531                             for (j = 0; j < sb->nr; j++)
532                             {
533                                 printf("%14d %1c\n", j, sb->cval[j]);
534                             }
535                             break;
536                         case xdr_datatype_string:
537                             for (j = 0; j < sb->nr; j++)
538                             {
539                                 printf("%14d %80s\n", j, sb->sval[j]);
540                             }
541                             break;
542                         default:
543                             gmx_incons("Unknown subblock type");
544                     }
545                 }
546             }
547         }
548     }
549     while (bCont);
550
551     close_enx(in);
552
553     free_enxframe(fr);
554     sfree(fr);
555     sfree(enm);
556 }
557
558 static void list_mtx(const char *fn)
559 {
560     int                  nrow, ncol, i, j, k;
561     real                *full   = nullptr, value;
562     gmx_sparsematrix_t * sparse = nullptr;
563
564     gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
565
566     if (full == nullptr)
567     {
568         snew(full, nrow*ncol);
569         for (i = 0; i < nrow*ncol; i++)
570         {
571             full[i] = 0;
572         }
573
574         for (i = 0; i < sparse->nrow; i++)
575         {
576             for (j = 0; j < sparse->ndata[i]; j++)
577             {
578                 k              = sparse->data[i][j].col;
579                 value          = sparse->data[i][j].value;
580                 full[i*ncol+k] = value;
581                 full[k*ncol+i] = value;
582             }
583         }
584         gmx_sparsematrix_destroy(sparse);
585     }
586
587     printf("%d %d\n", nrow, ncol);
588     for (i = 0; i < nrow; i++)
589     {
590         for (j = 0; j < ncol; j++)
591         {
592             printf(" %g", full[i*ncol+j]);
593         }
594         printf("\n");
595     }
596
597     sfree(full);
598 }
599
600 int gmx_dump(int argc, char *argv[])
601 {
602     const char *desc[] = {
603         "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
604         "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
605         "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
606         "and prints that to standard output in a readable format.",
607         "This program is essential for checking your run input file in case of",
608         "problems.[PAR]",
609         "The program can also preprocess a topology to help finding problems.",
610         "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
611         "directories used for searching include files.",
612     };
613     const char *bugs[] = {
614         "Position restraint output from -sys -s is broken"
615     };
616     t_filenm    fnm[] = {
617         { efTPR, "-s", nullptr, ffOPTRD },
618         { efTRX, "-f", nullptr, ffOPTRD },
619         { efEDR, "-e", nullptr, ffOPTRD },
620         { efCPT, nullptr, nullptr, ffOPTRD },
621         { efTOP, "-p", nullptr, ffOPTRD },
622         { efMTX, "-mtx", "hessian", ffOPTRD },
623         { efMDP, "-om", nullptr, ffOPTWR }
624     };
625 #define NFILE asize(fnm)
626
627     gmx_output_env_t *oenv;
628     /* Command line options */
629     gmx_bool          bShowNumbers      = TRUE;
630     gmx_bool          bShowParams       = FALSE;
631     gmx_bool          bSysTop           = FALSE;
632     gmx_bool          bOriginalInputrec = FALSE;
633     t_pargs           pa[]              = {
634         { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
635         { "-param", FALSE, etBOOL, {&bShowParams}, "Show parameters for each bonded interaction (for comparing dumps, it is useful to combine this with -nonr)" },
636         { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" },
637         { "-orgir", FALSE, etBOOL, {&bOriginalInputrec}, "Show input parameters from tpr as they were written by the version that produced the file, instead of how the current version reads them" }
638     };
639
640     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
641                            asize(desc), desc, asize(bugs), bugs, &oenv))
642     {
643         return 0;
644     }
645
646
647     if (ftp2bSet(efTPR, NFILE, fnm))
648     {
649         list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers, bShowParams,
650                  ftp2fn_null(efMDP, NFILE, fnm), bSysTop, bOriginalInputrec);
651     }
652     else if (ftp2bSet(efTRX, NFILE, fnm))
653     {
654         list_trx(ftp2fn(efTRX, NFILE, fnm));
655     }
656     else if (ftp2bSet(efEDR, NFILE, fnm))
657     {
658         list_ene(ftp2fn(efEDR, NFILE, fnm));
659     }
660     else if (ftp2bSet(efCPT, NFILE, fnm))
661     {
662         list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
663     }
664     else if (ftp2bSet(efTOP, NFILE, fnm))
665     {
666         list_top(ftp2fn(efTOP, NFILE, fnm));
667     }
668     else if (ftp2bSet(efMTX, NFILE, fnm))
669     {
670         list_mtx(ftp2fn(efMTX, NFILE, fnm));
671     }
672
673     return 0;
674 }