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