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