Remove gmx custom fixed int (e.g. gmx_int64_t) types
[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   = (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);
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));
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 /*! \brief Callback used by list_tng_for_gmx_dump. */
320 static void list_tng_inner(const char *fn,
321                            gmx_bool    bFirstFrame,
322                            real       *values,
323                            int64_t     step,
324                            double      frame_time,
325                            int64_t     n_values_per_frame,
326                            int64_t     n_atoms,
327                            real        prec,
328                            int64_t     nframe,
329                            char       *block_name)
330 {
331     char                 buf[256];
332     int                  indent = 0;
333
334     if (bFirstFrame)
335     {
336         sprintf(buf, "%s frame %" PRId64, fn, nframe);
337         indent = 0;
338         indent = pr_title(stdout, indent, buf);
339         pr_indent(stdout, indent);
340         fprintf(stdout, "natoms=%10" PRId64 "  step=%10" PRId64 "  time=%12.7e",
341                 n_atoms, step, frame_time);
342         if (prec > 0)
343         {
344             fprintf(stdout, "  prec=%10g", prec);
345         }
346         fprintf(stdout, "\n");
347     }
348     pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
349 }
350
351 static void list_tng(const char gmx_unused *fn)
352 {
353 #ifdef GMX_USE_TNG
354     gmx_tng_trajectory_t tng;
355     int64_t              nframe = 0;
356     int64_t              i, *block_ids = nullptr, step, ndatablocks;
357     gmx_bool             bOK;
358     real                *values = nullptr;
359
360     gmx_tng_open(fn, 'r', &tng);
361     gmx_print_tng_molecule_system(tng, stdout);
362
363     bOK    = gmx_get_tng_data_block_types_of_next_frame(tng, -1,
364                                                         0,
365                                                         nullptr,
366                                                         &step, &ndatablocks,
367                                                         &block_ids);
368     do
369     {
370         for (i = 0; i < ndatablocks; i++)
371         {
372             double               frame_time;
373             real                 prec;
374             int64_t              n_values_per_frame, n_atoms;
375             char                 block_name[STRLEN];
376
377             gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values,
378                                                       &step, &frame_time,
379                                                       &n_values_per_frame, &n_atoms,
380                                                       &prec,
381                                                       block_name, STRLEN, &bOK);
382             if (!bOK)
383             {
384                 /* Can't write any output because we don't know what
385                    arrays are valid. */
386                 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
387             }
388             else
389             {
390                 list_tng_inner(fn, (0 == i), values, step, frame_time,
391                                n_values_per_frame, n_atoms, prec, nframe, block_name);
392             }
393         }
394         nframe++;
395     }
396     while (gmx_get_tng_data_block_types_of_next_frame(tng, step,
397                                                       0,
398                                                       nullptr,
399                                                       &step,
400                                                       &ndatablocks,
401                                                       &block_ids));
402
403     if (block_ids)
404     {
405         sfree(block_ids);
406     }
407     sfree(values);
408     gmx_tng_close(&tng);
409 #endif
410 }
411
412 static void list_trx(const char *fn)
413 {
414     switch (fn2ftp(fn))
415     {
416         case efXTC:
417             list_xtc(fn);
418             break;
419         case efTRR:
420             list_trr(fn);
421             break;
422         case efTNG:
423             list_tng(fn);
424             break;
425         default:
426             fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
427                     fn, fn);
428     }
429 }
430
431 static void list_ene(const char *fn)
432 {
433     ener_file_t    in;
434     gmx_bool       bCont;
435     gmx_enxnm_t   *enm = nullptr;
436     t_enxframe    *fr;
437     int            i, j, nre, b;
438     char           buf[22];
439
440     printf("gmx dump: %s\n", fn);
441     in = open_enx(fn, "r");
442     do_enxnms(in, &nre, &enm);
443     assert(enm);
444
445     printf("energy components:\n");
446     for (i = 0; (i < nre); i++)
447     {
448         printf("%5d  %-24s (%s)\n", i, enm[i].name, enm[i].unit);
449     }
450
451     snew(fr, 1);
452     do
453     {
454         bCont = do_enx(in, fr);
455
456         if (bCont)
457         {
458             printf("\n%24s  %12.5e  %12s  %12s\n", "time:",
459                    fr->t, "step:", gmx_step_str(fr->step, buf));
460             printf("%24s  %12s  %12s  %12s\n",
461                    "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
462             printf("%24s  %12.5e  %12s  %12s\n",
463                    "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
464             if (fr->nre == nre)
465             {
466                 printf("%24s  %12s  %12s  %12s\n",
467                        "Component", "Energy", "Av. Energy", "Sum Energy");
468                 if (fr->nsum > 0)
469                 {
470                     for (i = 0; (i < nre); i++)
471                     {
472                         printf("%24s  %12.5e  %12.5e  %12.5e\n",
473                                enm[i].name, fr->ener[i].e, fr->ener[i].eav,
474                                fr->ener[i].esum);
475                     }
476                 }
477                 else
478                 {
479                     for (i = 0; (i < nre); i++)
480                     {
481                         printf("%24s  %12.5e\n",
482                                enm[i].name, fr->ener[i].e);
483                     }
484                 }
485             }
486             for (b = 0; b < fr->nblock; b++)
487             {
488                 const char *typestr = "";
489
490                 t_enxblock *eb = &(fr->block[b]);
491                 printf("Block data %2d (%3d subblocks, id=%d)\n",
492                        b, eb->nsub, eb->id);
493
494                 if (eb->id < enxNR)
495                 {
496                     typestr = enx_block_id_name[eb->id];
497                 }
498                 printf("  id='%s'\n", typestr);
499                 for (i = 0; i < eb->nsub; i++)
500                 {
501                     t_enxsubblock *sb = &(eb->sub[i]);
502                     printf("  Sub block %3d (%5d elems, type=%s) values:\n",
503                            i, sb->nr, xdr_datatype_names[sb->type]);
504
505                     switch (sb->type)
506                     {
507                         case xdr_datatype_float:
508                             for (j = 0; j < sb->nr; j++)
509                             {
510                                 printf("%14d   %8.4f\n", j, sb->fval[j]);
511                             }
512                             break;
513                         case xdr_datatype_double:
514                             for (j = 0; j < sb->nr; j++)
515                             {
516                                 printf("%14d   %10.6f\n", j, sb->dval[j]);
517                             }
518                             break;
519                         case xdr_datatype_int:
520                             for (j = 0; j < sb->nr; j++)
521                             {
522                                 printf("%14d %10d\n", j, sb->ival[j]);
523                             }
524                             break;
525                         case xdr_datatype_int64:
526                             for (j = 0; j < sb->nr; j++)
527                             {
528                                 printf("%14d %s\n",
529                                        j, gmx_step_str(sb->lval[j], buf));
530                             }
531                             break;
532                         case xdr_datatype_char:
533                             for (j = 0; j < sb->nr; j++)
534                             {
535                                 printf("%14d %1c\n", j, sb->cval[j]);
536                             }
537                             break;
538                         case xdr_datatype_string:
539                             for (j = 0; j < sb->nr; j++)
540                             {
541                                 printf("%14d %80s\n", j, sb->sval[j]);
542                             }
543                             break;
544                         default:
545                             gmx_incons("Unknown subblock type");
546                     }
547                 }
548             }
549         }
550     }
551     while (bCont);
552
553     close_enx(in);
554
555     free_enxframe(fr);
556     sfree(fr);
557     sfree(enm);
558 }
559
560 static void list_mtx(const char *fn)
561 {
562     int                  nrow, ncol, i, j, k;
563     real                *full   = nullptr, value;
564     gmx_sparsematrix_t * sparse = nullptr;
565
566     gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
567
568     if (full == nullptr)
569     {
570         snew(full, nrow*ncol);
571         for (i = 0; i < nrow*ncol; i++)
572         {
573             full[i] = 0;
574         }
575
576         for (i = 0; i < sparse->nrow; i++)
577         {
578             for (j = 0; j < sparse->ndata[i]; j++)
579             {
580                 k              = sparse->data[i][j].col;
581                 value          = sparse->data[i][j].value;
582                 full[i*ncol+k] = value;
583                 full[k*ncol+i] = value;
584             }
585         }
586         gmx_sparsematrix_destroy(sparse);
587     }
588
589     printf("%d %d\n", nrow, ncol);
590     for (i = 0; i < nrow; i++)
591     {
592         for (j = 0; j < ncol; j++)
593         {
594             printf(" %g", full[i*ncol+j]);
595         }
596         printf("\n");
597     }
598
599     sfree(full);
600 }
601
602 int gmx_dump(int argc, char *argv[])
603 {
604     const char *desc[] = {
605         "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
606         "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
607         "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
608         "and prints that to standard output in a readable format.",
609         "This program is essential for checking your run input file in case of",
610         "problems.[PAR]",
611         "The program can also preprocess a topology to help finding problems.",
612         "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
613         "directories used for searching include files.",
614     };
615     const char *bugs[] = {
616         "Position restraint output from -sys -s is broken"
617     };
618     t_filenm    fnm[] = {
619         { efTPR, "-s", nullptr, ffOPTRD },
620         { efTRX, "-f", nullptr, ffOPTRD },
621         { efEDR, "-e", nullptr, ffOPTRD },
622         { efCPT, nullptr, nullptr, ffOPTRD },
623         { efTOP, "-p", nullptr, ffOPTRD },
624         { efMTX, "-mtx", "hessian", ffOPTRD },
625         { efMDP, "-om", nullptr, ffOPTWR }
626     };
627 #define NFILE asize(fnm)
628
629     gmx_output_env_t *oenv;
630     /* Command line options */
631     gmx_bool          bShowNumbers      = TRUE;
632     gmx_bool          bShowParams       = FALSE;
633     gmx_bool          bSysTop           = FALSE;
634     gmx_bool          bOriginalInputrec = FALSE;
635     t_pargs           pa[]              = {
636         { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
637         { "-param", FALSE, etBOOL, {&bShowParams}, "Show parameters for each bonded interaction (for comparing dumps, it is useful to combine this with -nonr)" },
638         { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" },
639         { "-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" }
640     };
641
642     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
643                            asize(desc), desc, asize(bugs), bugs, &oenv))
644     {
645         return 0;
646     }
647
648
649     if (ftp2bSet(efTPR, NFILE, fnm))
650     {
651         list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers, bShowParams,
652                  ftp2fn_null(efMDP, NFILE, fnm), bSysTop, bOriginalInputrec);
653     }
654     else if (ftp2bSet(efTRX, NFILE, fnm))
655     {
656         list_trx(ftp2fn(efTRX, NFILE, fnm));
657     }
658     else if (ftp2bSet(efEDR, NFILE, fnm))
659     {
660         list_ene(ftp2fn(efEDR, NFILE, fnm));
661     }
662     else if (ftp2bSet(efCPT, NFILE, fnm))
663     {
664         list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
665     }
666     else if (ftp2bSet(efTOP, NFILE, fnm))
667     {
668         list_top(ftp2fn(efTOP, NFILE, fnm));
669     }
670     else if (ftp2bSet(efMTX, NFILE, fnm))
671     {
672         list_mtx(ftp2fn(efMTX, NFILE, fnm));
673     }
674
675     return 0;
676 }