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