Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / fileio / filenm.c
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-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, 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 "filenm.h"
38
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <stdio.h>
44 #include <stdlib.h>
45 #include <string.h>
46
47 #include "macros.h"
48 #include "types/commrec.h"
49
50 #include "gromacs/utility/basedefinitions.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
55
56 /* XDR should be available on all platforms now,
57  * but we keep the possibility of turning it off...
58  */
59 #define USE_XDR
60
61 /* Use bitflag ... */
62 #define IS_SET(fn) ((fn.flag & ffSET) != 0)
63 #define IS_OPT(fn) ((fn.flag & ffOPT) != 0)
64 #define IS_MULT(fn) ((fn.flag & ffMULT) != 0)
65 #define UN_SET(fn) (fn.flag = (fn.flag & ~ffSET))
66 #define DO_SET(fn) (fn.flag = (fn.flag |  ffSET))
67
68 enum
69 {
70     eftASC, eftBIN, eftXDR, eftTNG, eftGEN, eftNR
71 };
72
73 /* To support multiple file types with one general (eg TRX) we have
74  * these arrays.
75  */
76 static const int trxs[] =
77 {
78 #ifdef USE_XDR
79     efXTC, efTRR, efCPT,
80 #endif
81     efTRJ, efGRO, efG96, efPDB, efTNG
82 };
83 #define NTRXS asize(trxs)
84
85 static const int trcompressed[] =
86 {
87 #ifdef USE_XDR
88     efXTC,
89 #endif
90     efTNG
91 };
92 #define NTRCOMPRESSED asize(trcompressed)
93
94 static const int tros[] =
95 {
96 #ifdef USE_XDR
97     efXTC, efTRR,
98 #endif
99     efTRJ, efGRO, efG96, efPDB, efTNG
100 };
101 #define NTROS asize(tros)
102
103 static const int trns[] =
104 {
105 #ifdef USE_XDR
106     efTRR, efCPT,
107 #endif
108     efTRJ, efTNG
109 };
110 #define NTRNS asize(trns)
111
112 static const int stos[] =
113 { efGRO, efG96, efPDB, efBRK, efENT, efESP };
114 #define NSTOS asize(stos)
115
116 static const int stxs[] =
117 {
118     efGRO, efG96, efPDB, efBRK, efENT, efESP,
119 #ifdef USE_XDR
120     efTPR,
121 #endif
122     efTPB, efTPA
123 };
124 #define NSTXS asize(stxs)
125
126 static const int tpxs[] =
127 {
128 #ifdef USE_XDR
129     efTPR,
130 #endif
131     efTPB, efTPA
132 };
133 #define NTPXS asize(tpxs)
134
135 static const int tpss[] =
136 {
137 #ifdef USE_XDR
138     efTPR,
139 #endif
140     efTPB, efTPA, efGRO, efG96, efPDB, efBRK, efENT
141 };
142 #define NTPSS asize(tpss)
143
144 typedef struct
145 {
146     int         ftype;
147     const char *ext;
148     const char *defnm;
149     const char *defopt;
150     const char *descr;
151     int         ntps;
152     const int  *tps;
153 } t_deffile;
154
155 /* this array should correspond to the enum in filenm.h */
156 static const t_deffile
157     deffile[efNR] =
158 {
159     { eftASC, ".mdp", "grompp", "-f", "grompp input file with MD parameters" },
160     { eftGEN, ".???", "traj", "-f", "Trajectory", NTRXS, trxs },
161     { eftGEN, ".???", "trajout", "-f", "Trajectory", NTROS, tros },
162     { eftGEN, ".???", "traj", NULL,
163       "Full precision trajectory", NTRNS, trns },
164     { eftXDR, ".trr", "traj", NULL, "Trajectory in portable xdr format" },
165     { eftBIN, ".trj", "traj", NULL, "Trajectory file (architecture specific)" },
166     { eftGEN, ".???", "traj_comp", NULL,
167       "Compressed trajectory (tng format or portable xdr format)", NTRCOMPRESSED, trcompressed},
168     { eftXDR, ".xtc", "traj", NULL,
169       "Compressed trajectory (portable xdr format): xtc" },
170     { eftTNG, ".tng", "traj", NULL,
171       "Trajectory file (tng format)" },
172     { eftXDR, ".edr", "ener",   NULL, "Energy file"},
173     { eftGEN, ".???", "conf", "-c", "Structure file", NSTXS, stxs },
174     { eftGEN, ".???", "out", "-o", "Structure file", NSTOS, stos },
175     { eftASC, ".gro", "conf", "-c", "Coordinate file in Gromos-87 format" },
176     { eftASC, ".g96", "conf", "-c", "Coordinate file in Gromos-96 format" },
177     { eftASC, ".pdb", "eiwit",  "-f", "Protein data bank file"},
178     { eftASC, ".brk", "eiwit",  "-f", "Brookhaven data bank file"},
179     { eftASC, ".ent", "eiwit", "-f", "Entry in the protein date bank" },
180     { eftASC, ".esp", "conf", "-f", "Coordinate file in Espresso format" },
181     { eftASC, ".pqr", "state",  "-o", "Coordinate file for MEAD"},
182     { eftXDR, ".cpt", "state",  "-cp", "Checkpoint file"},
183     { eftASC, ".log", "run",    "-l", "Log file"},
184     { eftASC, ".xvg", "graph",  "-o", "xvgr/xmgr file"},
185     { eftASC, ".out", "hello",  "-o", "Generic output file"},
186     { eftASC, ".ndx", "index",  "-n", "Index file", },
187     { eftASC, ".top", "topol",  "-p", "Topology file"},
188     { eftASC, ".itp", "topinc", NULL, "Include file for topology"},
189     { eftGEN, ".???", "topol", "-s", "Run input file", NTPXS, tpxs },
190     { eftGEN, ".???", "topol", "-s", "Structure+mass(db)", NTPSS, tpss },
191     { eftXDR, ".tpr", "topol",  "-s", "Portable xdr run input file"},
192     { eftASC, ".tpa", "topol",  "-s", "Ascii run input file"},
193     { eftBIN, ".tpb", "topol",  "-s", "Binary run input file"},
194     { eftASC, ".tex", "doc",    "-o", "LaTeX file"},
195     { eftASC, ".rtp", "residue", NULL, "Residue Type file used by pdb2gmx" },
196     { eftASC, ".atp", "atomtp", NULL, "Atomtype file used by pdb2gmx" },
197     { eftASC, ".hdb", "polar",  NULL, "Hydrogen data base"},
198     { eftASC, ".dat", "nnnice", NULL, "Generic data file"},
199     { eftASC, ".dlg", "user",   NULL, "Dialog Box data for ngmx"},
200     { eftASC, ".map", "ss", NULL, "File that maps matrix data to colors" },
201     { eftASC, ".eps", "plot", NULL, "Encapsulated PostScript (tm) file" },
202     { eftASC, ".mat", "ss",     NULL, "Matrix Data file"},
203     { eftASC, ".m2p", "ps",     NULL, "Input file for mat2ps"},
204     { eftXDR, ".mtx", "hessian", "-m", "Hessian matrix"},
205     { eftASC, ".edi", "sam",    NULL, "ED sampling input"},
206     { eftASC, ".cub", "pot",  NULL, "Gaussian cube file" },
207     { eftASC, ".xpm", "root", NULL, "X PixMap compatible matrix file" },
208     { eftASC, "", "rundir", NULL, "Run directory" }
209 };
210
211 #define NZEXT 2
212 static const char *z_ext[NZEXT] =
213 { ".gz", ".Z" };
214
215 const char *ftp2ext(int ftp)
216 {
217     if ((0 <= ftp) && (ftp < efNR))
218     {
219         return deffile[ftp].ext[0] != '\0' ? deffile[ftp].ext + 1 : "";
220     }
221     else
222     {
223         return "unknown";
224     }
225 }
226
227 const char *ftp2ext_generic(int ftp)
228 {
229     if ((0 <= ftp) && (ftp < efNR))
230     {
231         switch (ftp)
232         {
233             case efTRX:
234                 return "trx";
235             case efTRN:
236                 return "trn";
237             case efSTO:
238                 return "sto";
239             case efSTX:
240                 return "stx";
241             case efTPX:
242                 return "tpx";
243             case efTPS:
244                 return "tps";
245             default:
246                 return ftp2ext(ftp);
247         }
248     }
249     else
250     {
251         return "unknown";
252     }
253 }
254
255 const char *ftp2ext_with_dot(int ftp)
256 {
257     if ((0 <= ftp) && (ftp < efNR))
258     {
259         return deffile[ftp].ext;
260     }
261     else
262     {
263         return "unknown";
264     }
265 }
266
267 int ftp2generic_count(int ftp)
268 {
269     if ((0 <= ftp) && (ftp < efNR))
270     {
271         return deffile[ftp].ntps;
272     }
273     else
274     {
275         return 0;
276     }
277 }
278
279 const int *ftp2generic_list(int ftp)
280 {
281     if ((0 <= ftp) && (ftp < efNR))
282     {
283         return deffile[ftp].tps;
284     }
285     else
286     {
287         return 0;
288     }
289 }
290
291 const char *ftp2desc(int ftp)
292 {
293     if ((0 <= ftp) && (ftp < efNR))
294     {
295         return deffile[ftp].descr;
296     }
297     else
298     {
299         return "unknown filetype";
300     }
301 }
302
303 const char *ftp2ftype(int ftp)
304 {
305     if ((ftp >= 0) && (ftp < efNR))
306     {
307         switch (deffile[ftp].ftype)
308         {
309             case eftASC:
310                 return "ASCII";
311             case eftBIN:
312                 return "Binary";
313             case eftXDR:
314                 return "XDR portable";
315             case eftTNG:
316                 return "TNG";
317             case eftGEN:
318                 return "";
319             default:
320                 gmx_fatal(FARGS, "Unknown filetype %d in ftp2ftype", deffile[ftp].ftype);
321                 break;
322         }
323     }
324     return "unknown";
325 }
326
327 const char *ftp2defnm(int ftp)
328 {
329     if ((0 <= ftp) && (ftp < efNR))
330     {
331         return deffile[ftp].defnm;
332     }
333     else
334     {
335         return NULL;
336     }
337 }
338
339 static void check_opts(int nf, t_filenm fnm[])
340 {
341     int              i;
342     const t_deffile *df;
343
344     for (i = 0; (i < nf); i++)
345     {
346         df = &(deffile[fnm[i].ftp]);
347         if (fnm[i].opt == NULL)
348         {
349             if (df->defopt == NULL)
350             {
351                 gmx_fatal(FARGS, "No default cmd-line option for %s (type %d)\n",
352                           deffile[fnm[i].ftp].ext, fnm[i].ftp);
353             }
354             else
355             {
356                 fnm[i].opt = df->defopt;
357             }
358         }
359     }
360 }
361
362 int fn2ftp(const char *fn)
363 {
364     int         i, len;
365     const char *feptr;
366     const char *eptr;
367
368     if (!fn)
369     {
370         return efNR;
371     }
372
373     len = strlen(fn);
374     if ((len >= 4) && (fn[len - 4] == '.'))
375     {
376         feptr = &(fn[len - 4]);
377     }
378     else
379     {
380         return efNR;
381     }
382
383     for (i = 0; (i < efNR); i++)
384     {
385         if ((eptr = deffile[i].ext) != NULL)
386         {
387             if (gmx_strcasecmp(feptr, eptr) == 0)
388             {
389                 break;
390             }
391         }
392     }
393
394     return i;
395 }
396
397 static void set_extension(char *buf, int ftp)
398 {
399     int              len, extlen;
400     const t_deffile *df;
401
402     /* check if extension is already at end of filename */
403     df     = &(deffile[ftp]);
404     len    = strlen(buf);
405     extlen = strlen(df->ext);
406     if ((len <= extlen) || (gmx_strcasecmp(&(buf[len - extlen]), df->ext) != 0))
407     {
408         strcat(buf, df->ext);
409     }
410 }
411
412 static void add_filenm(t_filenm *fnm, const char *filenm)
413 {
414     srenew(fnm->fns, fnm->nfiles+1);
415     fnm->fns[fnm->nfiles] = strdup(filenm);
416     fnm->nfiles++;
417 }
418
419 static void set_grpfnm(t_filenm *fnm, const char *name, const char *deffnm)
420 {
421     char       buf[256], buf2[256];
422     int        i, type;
423     gmx_bool   bValidExt;
424     int        nopts;
425     const int *ftps;
426
427     nopts = deffile[fnm->ftp].ntps;
428     ftps  = deffile[fnm->ftp].tps;
429     if ((nopts == 0) || (ftps == NULL))
430     {
431         gmx_fatal(FARGS, "nopts == 0 || ftps == NULL");
432     }
433
434     bValidExt = FALSE;
435     if (name && deffnm == NULL)
436     {
437         strcpy(buf, name);
438         /* First check whether we have a valid filename already */
439         type = fn2ftp(name);
440         if ((fnm->flag & ffREAD) && (fnm->ftp == efTRX))
441         {
442             /*if file exist don't add an extension for trajectory reading*/
443             bValidExt = gmx_fexist(name);
444         }
445         for (i = 0; (i < nopts) && !bValidExt; i++)
446         {
447             if (type == ftps[i])
448             {
449                 bValidExt = TRUE;
450             }
451         }
452     }
453     else if (deffnm != NULL)
454     {
455         strcpy(buf, deffnm);
456     }
457     else
458     {
459         /* No name given, set the default name */
460         strcpy(buf, ftp2defnm(fnm->ftp));
461     }
462
463     if (!bValidExt && (fnm->flag & ffREAD))
464     {
465         /* for input-files only: search for filenames in the directory */
466         for (i = 0; (i < nopts) && !bValidExt; i++)
467         {
468             type = ftps[i];
469             strcpy(buf2, buf);
470             set_extension(buf2, type);
471             if (gmx_fexist(buf2))
472             {
473                 bValidExt = TRUE;
474                 strcpy(buf, buf2);
475             }
476         }
477     }
478
479     if (!bValidExt)
480     {
481         /* Use the first extension type */
482         set_extension(buf, ftps[0]);
483     }
484
485     add_filenm(fnm, buf);
486 }
487
488 static void set_filenm(t_filenm *fnm, const char *name, const char *deffnm,
489                        gmx_bool bReadNode)
490 {
491     /* Set the default filename, extension and option for those fields that
492      * are not already set. An extension is added if not present, if fn = NULL
493      * or empty, the default filename is given.
494      */
495     char buf[256];
496     int  i, len, extlen;
497
498     if ((fnm->flag & ffREAD) && !bReadNode)
499     {
500         return;
501     }
502
503     if ((fnm->ftp < 0) || (fnm->ftp >= efNR))
504     {
505         gmx_fatal(FARGS, "file type out of range (%d)", fnm->ftp);
506     }
507
508     if (name)
509     {
510         strcpy(buf, name);
511     }
512     if ((fnm->flag & ffREAD) && name && gmx_fexist(name))
513     {
514         /* check if filename ends in .gz or .Z, if so remove that: */
515         len = strlen(name);
516         for (i = 0; i < NZEXT; i++)
517         {
518             extlen = strlen(z_ext[i]);
519             if (len > extlen)
520             {
521                 if (gmx_strcasecmp(name+len-extlen, z_ext[i]) == 0)
522                 {
523                     buf[len-extlen] = '\0';
524                     break;
525                 }
526             }
527         }
528     }
529
530     if (deffile[fnm->ftp].ntps)
531     {
532         set_grpfnm(fnm, name ? buf : NULL, deffnm);
533     }
534     else
535     {
536         if (name == NULL || deffnm != NULL)
537         {
538             if (deffnm != NULL)
539             {
540                 strcpy(buf, deffnm);
541             }
542             else
543             {
544                 strcpy(buf, ftp2defnm(fnm->ftp));
545             }
546         }
547         set_extension(buf, fnm->ftp);
548
549         add_filenm(fnm, buf);
550     }
551 }
552
553 static void set_filenms(int nf, t_filenm fnm[], const char *deffnm, gmx_bool bReadNode)
554 {
555     int i;
556
557     for (i = 0; (i < nf); i++)
558     {
559         if (!IS_SET(fnm[i]))
560         {
561             set_filenm(&(fnm[i]), fnm[i].fn, deffnm, bReadNode);
562         }
563     }
564 }
565
566 void parse_file_args(int *argc, char *argv[], int nf, t_filenm fnm[],
567                      const char *deffnm, gmx_bool bReadNode)
568 {
569     int       i, j;
570     gmx_bool *bRemove;
571
572     check_opts(nf, fnm);
573
574     for (i = 0; (i < nf); i++)
575     {
576         UN_SET(fnm[i]);
577     }
578
579     if (*argc > 1)
580     {
581         snew(bRemove, (*argc)+1);
582         i = 1;
583         do
584         {
585             for (j = 0; (j < nf); j++)
586             {
587                 if (strcmp(argv[i], fnm[j].opt) == 0)
588                 {
589                     DO_SET(fnm[j]);
590                     bRemove[i] = TRUE;
591                     i++;
592                     /* check if we are out of arguments for this option */
593                     if ((i >= *argc) || (argv[i][0] == '-'))
594                     {
595                         set_filenm(&fnm[j], fnm[j].fn, deffnm, bReadNode);
596                     }
597                     /* sweep up all file arguments for this option */
598                     while ((i < *argc) && (argv[i][0] != '-'))
599                     {
600                         set_filenm(&fnm[j], argv[i], NULL, bReadNode);
601                         bRemove[i] = TRUE;
602                         i++;
603                         /* only repeat for 'multiple' file options: */
604                         if (!IS_MULT(fnm[j]))
605                         {
606                             break;
607                         }
608                     }
609
610                     break; /* jump out of 'j' loop */
611                 }
612             }
613             /* No file found corresponding to option argv[i] */
614             if (j == nf)
615             {
616                 i++;
617             }
618         }
619         while (i < *argc);
620
621         /* Remove used entries */
622         for (i = j = 0; (i <= *argc); i++)
623         {
624             if (!bRemove[i])
625             {
626                 argv[j++] = argv[i];
627             }
628         }
629         (*argc) = j - 1;
630         sfree(bRemove);
631     }
632
633     set_filenms(nf, fnm, deffnm, bReadNode);
634
635 }
636
637 const char *opt2fn(const char *opt, int nfile, const t_filenm fnm[])
638 {
639     int i;
640
641     for (i = 0; (i < nfile); i++)
642     {
643         if (strcmp(opt, fnm[i].opt) == 0)
644         {
645             return fnm[i].fns[0];
646         }
647     }
648
649     fprintf(stderr, "No option %s\n", opt);
650
651     return NULL;
652 }
653
654 const char *opt2fn_master(const char *opt, int nfile, const t_filenm fnm[],
655                           t_commrec *cr)
656 {
657     return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : NULL;
658 }
659
660 int opt2fns(char **fns[], const char *opt, int nfile, const t_filenm fnm[])
661 {
662     int i;
663
664     for (i = 0; (i < nfile); i++)
665     {
666         if (strcmp(opt, fnm[i].opt) == 0)
667         {
668             *fns = fnm[i].fns;
669             return fnm[i].nfiles;
670         }
671     }
672
673     fprintf(stderr, "No option %s\n", opt);
674     return 0;
675 }
676
677 const char *ftp2fn(int ftp, int nfile, const t_filenm fnm[])
678 {
679     int i;
680
681     for (i = 0; (i < nfile); i++)
682     {
683         if (ftp == fnm[i].ftp)
684         {
685             return fnm[i].fns[0];
686         }
687     }
688
689     fprintf(stderr, "ftp2fn: No filetype %s\n", deffile[ftp].ext);
690     return NULL;
691 }
692
693 int ftp2fns(char **fns[], int ftp, int nfile, const t_filenm fnm[])
694 {
695     int i;
696
697     for (i = 0; (i < nfile); i++)
698     {
699         if (ftp == fnm[i].ftp)
700         {
701             *fns = fnm[i].fns;
702             return fnm[i].nfiles;
703         }
704     }
705
706     fprintf(stderr, "ftp2fn: No filetype %s\n", deffile[ftp].ext);
707     return 0;
708 }
709
710 gmx_bool ftp2bSet(int ftp, int nfile, const t_filenm fnm[])
711 {
712     int i;
713
714     for (i = 0; (i < nfile); i++)
715     {
716         if (ftp == fnm[i].ftp)
717         {
718             return (gmx_bool) IS_SET(fnm[i]);
719         }
720     }
721
722     fprintf(stderr, "ftp2bSet: No filetype %s\n", deffile[ftp].ext);
723
724     return FALSE;
725 }
726
727 gmx_bool opt2bSet(const char *opt, int nfile, const t_filenm fnm[])
728 {
729     int i;
730
731     for (i = 0; (i < nfile); i++)
732     {
733         if (strcmp(opt, fnm[i].opt) == 0)
734         {
735             return (gmx_bool) IS_SET(fnm[i]);
736         }
737     }
738
739     fprintf(stderr, "No option %s\n", opt);
740
741     return FALSE;
742 }
743
744 const char *opt2fn_null(const char *opt, int nfile, const t_filenm fnm[])
745 {
746     int i;
747
748     for (i = 0; (i < nfile); i++)
749     {
750         if (strcmp(opt, fnm[i].opt) == 0)
751         {
752             if (IS_OPT(fnm[i]) && !IS_SET(fnm[i]))
753             {
754                 return NULL;
755             }
756             else
757             {
758                 return fnm[i].fns[0];
759             }
760         }
761     }
762     fprintf(stderr, "No option %s\n", opt);
763     return NULL;
764 }
765
766 const char *ftp2fn_null(int ftp, int nfile, const t_filenm fnm[])
767 {
768     int i;
769
770     for (i = 0; (i < nfile); i++)
771     {
772         if (ftp == fnm[i].ftp)
773         {
774             if (IS_OPT(fnm[i]) && !IS_SET(fnm[i]))
775             {
776                 return NULL;
777             }
778             else
779             {
780                 return fnm[i].fns[0];
781             }
782         }
783     }
784     fprintf(stderr, "ftp2fn: No filetype %s\n", deffile[ftp].ext);
785     return NULL;
786 }
787
788 gmx_bool is_optional(const t_filenm *fnm)
789 {
790     return ((fnm->flag & ffOPT) == ffOPT);
791 }
792
793 gmx_bool is_output(const t_filenm *fnm)
794 {
795     return ((fnm->flag & ffWRITE) == ffWRITE);
796 }
797
798 gmx_bool is_set(const t_filenm *fnm)
799 {
800     return ((fnm->flag & ffSET) == ffSET);
801 }
802
803 int add_suffix_to_output_names(t_filenm *fnm, int nfile, const char *suffix)
804 {
805     int   i, j, pos;
806     char  buf[STRLEN], newname[STRLEN];
807     char *extpos;
808
809     for (i = 0; i < nfile; i++)
810     {
811         if (is_output(&fnm[i]) && fnm[i].ftp != efCPT)
812         {
813             /* We never use multiple _outputs_, but we might as well check
814                for it, just in case... */
815             for (j = 0; j < fnm[i].nfiles; j++)
816             {
817                 strncpy(buf, fnm[i].fns[j], STRLEN - 1);
818                 extpos  = strrchr(buf, '.');
819                 *extpos = '\0';
820                 sprintf(newname, "%s%s.%s", buf, suffix, extpos + 1);
821                 free(fnm[i].fns[j]);
822                 fnm[i].fns[j] = strdup(newname);
823             }
824         }
825     }
826     return 0;
827 }
828
829 t_filenm *dup_tfn(int nf, const t_filenm tfn[])
830 {
831     int       i, j;
832     t_filenm *ret;
833
834     snew(ret, nf);
835     for (i = 0; i < nf; i++)
836     {
837         ret[i] = tfn[i]; /* just directly copy all non-string fields */
838         if (tfn[i].opt)
839         {
840             ret[i].opt = strdup(tfn[i].opt);
841         }
842         else
843         {
844             ret[i].opt = NULL;
845         }
846
847         if (tfn[i].fn)
848         {
849             ret[i].fn = strdup(tfn[i].fn);
850         }
851         else
852         {
853             ret[i].fn = NULL;
854         }
855
856         if (tfn[i].nfiles > 0)
857         {
858             snew(ret[i].fns, tfn[i].nfiles);
859             for (j = 0; j < tfn[i].nfiles; j++)
860             {
861                 ret[i].fns[j] = strdup(tfn[i].fns[j]);
862             }
863         }
864     }
865     return ret;
866 }
867
868 void done_filenms(int nf, t_filenm fnm[])
869 {
870     int i, j;
871
872     for (i = 0; i < nf; ++i)
873     {
874         for (j = 0; j < fnm[i].nfiles; ++j)
875         {
876             sfree(fnm[i].fns[j]);
877         }
878         sfree(fnm[i].fns);
879         fnm[i].fns = NULL;
880     }
881 }