Removed unused parameters and cleanup of inputrec.
[alexxy/gromacs.git] / src / gromacs / gmxlib / txtdump.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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 /* This file is completely threadsafe - please keep it that way! */
42
43 #include <stdio.h>
44 #include "typedefs.h"
45 #include "types/commrec.h"
46 #include "names.h"
47 #include "txtdump.h"
48 #include "vec.h"
49 #include "macros.h"
50
51 #include "gmx_fatal.h"
52 #include "gromacs/utility/smalloc.h"
53
54 int pr_indent(FILE *fp, int n)
55 {
56     int i;
57
58     for (i = 0; i < n; i++)
59     {
60         (void) fprintf(fp, " ");
61     }
62     return n;
63 }
64
65 int available(FILE *fp, void *p, int indent, const char *title)
66 {
67     if (!p)
68     {
69         if (indent > 0)
70         {
71             pr_indent(fp, indent);
72         }
73         (void) fprintf(fp, "%s: not available\n", title);
74     }
75     return (p != NULL);
76 }
77
78 int pr_title(FILE *fp, int indent, const char *title)
79 {
80     (void) pr_indent(fp, indent);
81     (void) fprintf(fp, "%s:\n", title);
82     return (indent+INDENT);
83 }
84
85 int pr_title_n(FILE *fp, int indent, const char *title, int n)
86 {
87     (void) pr_indent(fp, indent);
88     (void) fprintf(fp, "%s (%d):\n", title, n);
89     return (indent+INDENT);
90 }
91
92 int pr_title_nxn(FILE *fp, int indent, const char *title, int n1, int n2)
93 {
94     (void) pr_indent(fp, indent);
95     (void) fprintf(fp, "%s (%dx%d):\n", title, n1, n2);
96     return (indent+INDENT);
97 }
98
99 void pr_ivec(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
100 {
101     int i;
102
103     if (available(fp, vec, indent, title))
104     {
105         indent = pr_title_n(fp, indent, title, n);
106         for (i = 0; i < n; i++)
107         {
108             (void) pr_indent(fp, indent);
109             (void) fprintf(fp, "%s[%d]=%d\n", title, bShowNumbers ? i : -1, vec[i]);
110         }
111     }
112 }
113
114 void pr_ivec_block(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
115 {
116     int i, j;
117
118     if (available(fp, vec, indent, title))
119     {
120         indent = pr_title_n(fp, indent, title, n);
121         i      = 0;
122         while (i < n)
123         {
124             j = i+1;
125             while (j < n && vec[j] == vec[j-1]+1)
126             {
127                 j++;
128             }
129             /* Print consecutive groups of 3 or more as blocks */
130             if (j - i < 3)
131             {
132                 while (i < j)
133                 {
134                     (void) pr_indent(fp, indent);
135                     (void) fprintf(fp, "%s[%d]=%d\n",
136                                    title, bShowNumbers ? i : -1, vec[i]);
137                     i++;
138                 }
139             }
140             else
141             {
142                 (void) pr_indent(fp, indent);
143                 (void) fprintf(fp, "%s[%d,...,%d] = {%d,...,%d}\n",
144                                title,
145                                bShowNumbers ? i : -1,
146                                bShowNumbers ? j-1 : -1,
147                                vec[i], vec[j-1]);
148                 i = j;
149             }
150         }
151     }
152 }
153
154 void pr_bvec(FILE *fp, int indent, const char *title, gmx_bool vec[], int n, gmx_bool bShowNumbers)
155 {
156     int i;
157
158     if (available(fp, vec, indent, title))
159     {
160         indent = pr_title_n(fp, indent, title, n);
161         for (i = 0; i < n; i++)
162         {
163             (void) pr_indent(fp, indent);
164             (void) fprintf(fp, "%s[%d]=%s\n", title, bShowNumbers ? i : -1,
165                            EBOOL(vec[i]));
166         }
167     }
168 }
169
170 void pr_ivecs(FILE *fp, int indent, const char *title, ivec vec[], int n, gmx_bool bShowNumbers)
171 {
172     int i, j;
173
174     if (available(fp, vec, indent, title))
175     {
176         indent = pr_title_nxn(fp, indent, title, n, DIM);
177         for (i = 0; i < n; i++)
178         {
179             (void) pr_indent(fp, indent);
180             (void) fprintf(fp, "%s[%d]={", title, bShowNumbers ? i : -1);
181             for (j = 0; j < DIM; j++)
182             {
183                 if (j != 0)
184                 {
185                     (void) fprintf(fp, ", ");
186                 }
187                 fprintf(fp, "%d", vec[i][j]);
188             }
189             (void) fprintf(fp, "}\n");
190         }
191     }
192 }
193
194 void pr_rvec(FILE *fp, int indent, const char *title, real vec[], int n, gmx_bool bShowNumbers)
195 {
196     int i;
197
198     if (available(fp, vec, indent, title))
199     {
200         indent = pr_title_n(fp, indent, title, n);
201         for (i = 0; i < n; i++)
202         {
203             pr_indent(fp, indent);
204             fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
205         }
206     }
207 }
208
209 void pr_dvec(FILE *fp, int indent, const char *title, double vec[], int n, gmx_bool bShowNumbers)
210 {
211     int i;
212
213     if (available(fp, vec, indent, title))
214     {
215         indent = pr_title_n(fp, indent, title, n);
216         for (i = 0; i < n; i++)
217         {
218             pr_indent(fp, indent);
219             fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
220         }
221     }
222 }
223
224
225 /*
226    void pr_mat(FILE *fp,int indent,char *title,matrix m)
227    {
228    int i,j;
229
230    if (available(fp,m,indent,title)) {
231     indent=pr_title_n(fp,indent,title,n);
232     for(i=0; i<n; i++) {
233       pr_indent(fp,indent);
234       fprintf(fp,"%s[%d]=%12.5e %12.5e %12.5e\n",
235           title,bShowNumbers?i:-1,m[i][XX],m[i][YY],m[i][ZZ]);
236     }
237    }
238    }
239  */
240
241 void pr_rvecs_len(FILE *fp, int indent, const char *title, rvec vec[], int n)
242 {
243     int i, j;
244
245     if (available(fp, vec, indent, title))
246     {
247         indent = pr_title_nxn(fp, indent, title, n, DIM);
248         for (i = 0; i < n; i++)
249         {
250             (void) pr_indent(fp, indent);
251             (void) fprintf(fp, "%s[%5d]={", title, i);
252             for (j = 0; j < DIM; j++)
253             {
254                 if (j != 0)
255                 {
256                     (void) fprintf(fp, ", ");
257                 }
258                 (void) fprintf(fp, "%12.5e", vec[i][j]);
259             }
260             (void) fprintf(fp, "} len=%12.5e\n", norm(vec[i]));
261         }
262     }
263 }
264
265 void pr_rvecs(FILE *fp, int indent, const char *title, rvec vec[], int n)
266 {
267     const char *fshort = "%12.5e";
268     const char *flong  = "%15.8e";
269     const char *format;
270     int         i, j;
271
272     if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
273     {
274         format = flong;
275     }
276     else
277     {
278         format = fshort;
279     }
280
281     if (available(fp, vec, indent, title))
282     {
283         indent = pr_title_nxn(fp, indent, title, n, DIM);
284         for (i = 0; i < n; i++)
285         {
286             (void) pr_indent(fp, indent);
287             (void) fprintf(fp, "%s[%5d]={", title, i);
288             for (j = 0; j < DIM; j++)
289             {
290                 if (j != 0)
291                 {
292                     (void) fprintf(fp, ", ");
293                 }
294                 (void) fprintf(fp, format, vec[i][j]);
295             }
296             (void) fprintf(fp, "}\n");
297         }
298     }
299 }
300
301
302 void pr_rvecs_of_dim(FILE *fp, int indent, const char *title, rvec vec[], int n, int dim)
303 {
304     const char *fshort = "%12.5e";
305     const char *flong  = "%15.8e";
306     const char *format;
307     int         i, j;
308
309     if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
310     {
311         format = flong;
312     }
313     else
314     {
315         format = fshort;
316     }
317
318     if (available(fp, vec, indent, title))
319     {
320         indent = pr_title_nxn(fp, indent, title, n, dim);
321         for (i = 0; i < n; i++)
322         {
323             (void) pr_indent(fp, indent);
324             (void) fprintf(fp, "%s[%5d]={", title, i);
325             for (j = 0; j < dim; j++)
326             {
327                 if (j != 0)
328                 {
329                     (void) fprintf(fp, ", ");
330                 }
331                 (void) fprintf(fp, format, vec[i][j]);
332             }
333             (void) fprintf(fp, "}\n");
334         }
335     }
336 }
337
338 void pr_reals(FILE *fp, int indent, const char *title, real *vec, int n)
339 {
340     int i;
341
342     if (available(fp, vec, indent, title))
343     {
344         (void) pr_indent(fp, indent);
345         (void) fprintf(fp, "%s:\t", title);
346         for (i = 0; i < n; i++)
347         {
348             fprintf(fp, "  %10g", vec[i]);
349         }
350         (void) fprintf(fp, "\n");
351     }
352 }
353
354 void pr_doubles(FILE *fp, int indent, const char *title, double *vec, int n)
355 {
356     int i;
357
358     if (available(fp, vec, indent, title))
359     {
360         (void) pr_indent(fp, indent);
361         (void) fprintf(fp, "%s:\t", title);
362         for (i = 0; i < n; i++)
363         {
364             fprintf(fp, "  %10g", vec[i]);
365         }
366         (void) fprintf(fp, "\n");
367     }
368 }
369
370 void pr_reals_of_dim(FILE *fp, int indent, const char *title, real *vec, int n, int dim)
371 {
372     int         i, j;
373     const char *fshort = "%12.5e";
374     const char *flong  = "%15.8e";
375     const char *format;
376
377     if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
378     {
379         format = flong;
380     }
381     else
382     {
383         format = fshort;
384     }
385
386     if (available(fp, vec, indent, title))
387     {
388         indent = pr_title_nxn(fp, indent, title, n, dim);
389         for (i = 0; i < n; i++)
390         {
391             (void) pr_indent(fp, indent);
392             (void) fprintf(fp, "%s[%5d]={", title, i);
393             for (j = 0; j < dim; j++)
394             {
395                 if (j != 0)
396                 {
397                     (void) fprintf(fp, ", ");
398                 }
399                 (void) fprintf(fp, format, vec[i * dim  + j]);
400             }
401             (void) fprintf(fp, "}\n");
402         }
403     }
404 }
405
406 static void pr_int(FILE *fp, int indent, const char *title, int i)
407 {
408     pr_indent(fp, indent);
409     fprintf(fp, "%-30s = %d\n", title, i);
410 }
411
412 static void pr_int64(FILE *fp, int indent, const char *title, gmx_int64_t i)
413 {
414     char buf[STEPSTRSIZE];
415
416     pr_indent(fp, indent);
417     fprintf(fp, "%-30s = %s\n", title, gmx_step_str(i, buf));
418 }
419
420 static void pr_real(FILE *fp, int indent, const char *title, real r)
421 {
422     pr_indent(fp, indent);
423     fprintf(fp, "%-30s = %g\n", title, r);
424 }
425
426 static void pr_double(FILE *fp, int indent, const char *title, double d)
427 {
428     pr_indent(fp, indent);
429     fprintf(fp, "%-30s = %g\n", title, d);
430 }
431
432 static void pr_str(FILE *fp, int indent, const char *title, const char *s)
433 {
434     pr_indent(fp, indent);
435     fprintf(fp, "%-30s = %s\n", title, s);
436 }
437
438 void pr_qm_opts(FILE *fp, int indent, const char *title, t_grpopts *opts)
439 {
440     int i, m, j;
441
442     fprintf(fp, "%s:\n", title);
443
444     pr_int(fp, indent, "ngQM", opts->ngQM);
445     if (opts->ngQM > 0)
446     {
447         pr_ivec(fp, indent, "QMmethod", opts->QMmethod, opts->ngQM, FALSE);
448         pr_ivec(fp, indent, "QMbasis", opts->QMbasis, opts->ngQM, FALSE);
449         pr_ivec(fp, indent, "QMcharge", opts->QMcharge, opts->ngQM, FALSE);
450         pr_ivec(fp, indent, "QMmult", opts->QMmult, opts->ngQM, FALSE);
451         pr_bvec(fp, indent, "bSH", opts->bSH, opts->ngQM, FALSE);
452         pr_ivec(fp, indent, "CASorbitals", opts->CASorbitals, opts->ngQM, FALSE);
453         pr_ivec(fp, indent, "CASelectrons", opts->CASelectrons, opts->ngQM, FALSE);
454         pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
455         pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
456         pr_ivec(fp, indent, "SAsteps", opts->SAsteps, opts->ngQM, FALSE);
457         pr_bvec(fp, indent, "bOPT", opts->bOPT, opts->ngQM, FALSE);
458         pr_bvec(fp, indent, "bTS", opts->bTS, opts->ngQM, FALSE);
459     }
460 }
461
462 static void pr_grp_opts(FILE *out, int indent, const char *title, t_grpopts *opts,
463                         gmx_bool bMDPformat)
464 {
465     int i, m, j;
466
467     if (!bMDPformat)
468     {
469         fprintf(out, "%s:\n", title);
470     }
471
472     pr_indent(out, indent);
473     fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
474     for (i = 0; (i < opts->ngtc); i++)
475     {
476         fprintf(out, "  %10g", opts->nrdf[i]);
477     }
478     fprintf(out, "\n");
479
480     pr_indent(out, indent);
481     fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
482     for (i = 0; (i < opts->ngtc); i++)
483     {
484         fprintf(out, "  %10g", opts->ref_t[i]);
485     }
486     fprintf(out, "\n");
487
488     pr_indent(out, indent);
489     fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
490     for (i = 0; (i < opts->ngtc); i++)
491     {
492         fprintf(out, "  %10g", opts->tau_t[i]);
493     }
494     fprintf(out, "\n");
495
496     /* Pretty-print the simulated annealing info */
497     fprintf(out, "anneal%s", bMDPformat ? " = " : ":");
498     for (i = 0; (i < opts->ngtc); i++)
499     {
500         fprintf(out, "  %10s", EANNEAL(opts->annealing[i]));
501     }
502     fprintf(out, "\n");
503
504     fprintf(out, "ann-npoints%s", bMDPformat ? " = " : ":");
505     for (i = 0; (i < opts->ngtc); i++)
506     {
507         fprintf(out, "  %10d", opts->anneal_npoints[i]);
508     }
509     fprintf(out, "\n");
510
511     for (i = 0; (i < opts->ngtc); i++)
512     {
513         if (opts->anneal_npoints[i] > 0)
514         {
515             fprintf(out, "ann. times [%d]:\t", i);
516             for (j = 0; (j < opts->anneal_npoints[i]); j++)
517             {
518                 fprintf(out, "  %10.1f", opts->anneal_time[i][j]);
519             }
520             fprintf(out, "\n");
521             fprintf(out, "ann. temps [%d]:\t", i);
522             for (j = 0; (j < opts->anneal_npoints[i]); j++)
523             {
524                 fprintf(out, "  %10.1f", opts->anneal_temp[i][j]);
525             }
526             fprintf(out, "\n");
527         }
528     }
529
530     pr_indent(out, indent);
531     fprintf(out, "acc:\t");
532     for (i = 0; (i < opts->ngacc); i++)
533     {
534         for (m = 0; (m < DIM); m++)
535         {
536             fprintf(out, "  %10g", opts->acc[i][m]);
537         }
538     }
539     fprintf(out, "\n");
540
541     pr_indent(out, indent);
542     fprintf(out, "nfreeze:");
543     for (i = 0; (i < opts->ngfrz); i++)
544     {
545         for (m = 0; (m < DIM); m++)
546         {
547             fprintf(out, "  %10s", opts->nFreeze[i][m] ? "Y" : "N");
548         }
549     }
550     fprintf(out, "\n");
551
552
553     for (i = 0; (i < opts->ngener); i++)
554     {
555         pr_indent(out, indent);
556         fprintf(out, "energygrp-flags[%3d]:", i);
557         for (m = 0; (m < opts->ngener); m++)
558         {
559             fprintf(out, " %d", opts->egp_flags[opts->ngener*i+m]);
560         }
561         fprintf(out, "\n");
562     }
563
564     fflush(out);
565 }
566
567 static void pr_matrix(FILE *fp, int indent, const char *title, rvec *m,
568                       gmx_bool bMDPformat)
569 {
570     if (bMDPformat)
571     {
572         fprintf(fp, "%-10s    = %g %g %g %g %g %g\n", title,
573                 m[XX][XX], m[YY][YY], m[ZZ][ZZ], m[XX][YY], m[XX][ZZ], m[YY][ZZ]);
574     }
575     else
576     {
577         pr_rvecs(fp, indent, title, m, DIM);
578     }
579 }
580
581 static void pr_cosine(FILE *fp, int indent, const char *title, t_cosines *cos,
582                       gmx_bool bMDPformat)
583 {
584     int j;
585
586     if (bMDPformat)
587     {
588         fprintf(fp, "%s = %d\n", title, cos->n);
589     }
590     else
591     {
592         indent = pr_title(fp, indent, title);
593         (void) pr_indent(fp, indent);
594         fprintf(fp, "n = %d\n", cos->n);
595         if (cos->n > 0)
596         {
597             (void) pr_indent(fp, indent+2);
598             fprintf(fp, "a =");
599             for (j = 0; (j < cos->n); j++)
600             {
601                 fprintf(fp, " %e", cos->a[j]);
602             }
603             fprintf(fp, "\n");
604             (void) pr_indent(fp, indent+2);
605             fprintf(fp, "phi =");
606             for (j = 0; (j < cos->n); j++)
607             {
608                 fprintf(fp, " %e", cos->phi[j]);
609             }
610             fprintf(fp, "\n");
611         }
612     }
613 }
614
615 #define PS(t, s) pr_str(fp, indent, t, s)
616 #define PI(t, s) pr_int(fp, indent, t, s)
617 #define PSTEP(t, s) pr_int64(fp, indent, t, s)
618 #define PR(t, s) pr_real(fp, indent, t, s)
619 #define PD(t, s) pr_double(fp, indent, t, s)
620
621 static void pr_pull_group(FILE *fp, int indent, int g, t_pull_group *pgrp)
622 {
623     pr_indent(fp, indent);
624     fprintf(fp, "pull-group %d:\n", g);
625     indent += 2;
626     pr_ivec_block(fp, indent, "atom", pgrp->ind, pgrp->nat, TRUE);
627     pr_rvec(fp, indent, "weight", pgrp->weight, pgrp->nweight, TRUE);
628     PI("pbcatom", pgrp->pbcatom);
629 }
630
631 static void pr_pull_coord(FILE *fp, int indent, int c, t_pull_coord *pcrd)
632 {
633     pr_indent(fp, indent);
634     fprintf(fp, "pull-coord %d:\n", c);
635     PI("group[0]", pcrd->group[0]);
636     PI("group[1]", pcrd->group[1]);
637     pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
638     pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
639     PR("init", pcrd->init);
640     PR("rate", pcrd->rate);
641     PR("k", pcrd->k);
642     PR("kB", pcrd->kB);
643 }
644
645 static void pr_simtempvals(FILE *fp, int indent, t_simtemp *simtemp, int n_lambda)
646 {
647     PR("simtemp_low", simtemp->simtemp_low);
648     PR("simtemp_high", simtemp->simtemp_high);
649     PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
650     pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
651 }
652
653 static void pr_expandedvals(FILE *fp, int indent, t_expanded *expand, int n_lambda)
654 {
655
656     PI("nstexpanded", expand->nstexpanded);
657     PS("lambda-stats", elamstats_names[expand->elamstats]);
658     PS("lambda-mc-move", elmcmove_names[expand->elmcmove]);
659     PI("lmc-repeats", expand->lmc_repeats);
660     PI("lmc-gibbsdelta", expand->gibbsdeltalam);
661     PI("lmc-nstart", expand->lmc_forced_nstart);
662     PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
663     PI("nst-transition-matrix", expand->nstTij);
664     PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
665     PI("weight-c-range", expand->c_range);    /* default is just C=0 */
666     PR("wl-scale", expand->wl_scale);
667     PR("init-wl-delta", expand->init_wl_delta);
668     PR("wl-ratio", expand->wl_ratio);
669     PS("bWLoneovert", EBOOL(expand->bWLoneovert));
670     PI("lmc-seed", expand->lmc_seed);
671     PR("mc-temperature", expand->mc_temp);
672     PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
673     if (expand->elmceq == elmceqNUMATLAM)
674     {
675         PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
676     }
677     if (expand->elmceq == elmceqSAMPLES)
678     {
679         PI("weight-equil-number-samples", expand->equil_samples);
680     }
681     if (expand->elmceq == elmceqSTEPS)
682     {
683         PI("weight-equil-number-steps", expand->equil_steps);
684     }
685     if (expand->elmceq == elmceqWLDELTA)
686     {
687         PR("weight-equil-wl-delta", expand->equil_wl_delta);
688     }
689     if (expand->elmceq == elmceqRATIO)
690     {
691         PR("weight-equil-count-ratio", expand->equil_ratio);
692     }
693
694     pr_indent(fp, indent);
695     pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
696     PS("init-weights", EBOOL(expand->bInit_weights));
697 }
698
699 static void pr_fepvals(FILE *fp, int indent, t_lambda *fep, gmx_bool bMDPformat)
700 {
701     int i, j;
702
703     PI("nstdhdl", fep->nstdhdl);
704     PI("init-lambda-state", fep->init_fep_state);
705     PR("init-lambda", fep->init_lambda);
706     PR("delta-lambda", fep->delta_lambda);
707     if (!bMDPformat)
708     {
709         PI("n-lambdas", fep->n_lambda);
710     }
711     if (fep->n_lambda > 0)
712     {
713         pr_indent(fp, indent);
714         fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
715         for (i = 0; i < efptNR; i++)
716         {
717             fprintf(fp, "%18s = ", efpt_names[i]);
718             if (fep->separate_dvdl[i])
719             {
720                 fprintf(fp, "  TRUE");
721             }
722             else
723             {
724                 fprintf(fp, "  FALSE");
725             }
726             fprintf(fp, "\n");
727         }
728         fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
729         for (i = 0; i < efptNR; i++)
730         {
731             fprintf(fp, "%18s = ", efpt_names[i]);
732             for (j = 0; j < fep->n_lambda; j++)
733             {
734                 fprintf(fp, "  %10g", fep->all_lambda[i][j]);
735             }
736             fprintf(fp, "\n");
737         }
738     }
739     PI("calc-lambda-neighbors", fep->lambda_neighbors);
740
741     PR("sc-alpha", fep->sc_alpha);
742     PS("bScCoul", EBOOL(fep->bScCoul));
743     PS("bScPrintEnergy", EBOOL(fep->bPrintEnergy));
744     PI("sc-power", fep->sc_power);
745     PR("sc-r-power", fep->sc_r_power);
746     PR("sc-sigma", fep->sc_sigma);
747     PR("sc-sigma-min", fep->sc_sigma_min);
748     PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
749     PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
750     PI("dh-hist-size", fep->dh_hist_size);
751     PD("dh-hist-spacing", fep->dh_hist_spacing);
752 };
753
754 static void pr_pull(FILE *fp, int indent, t_pull *pull)
755 {
756     int g;
757
758     PS("pull-geometry", EPULLGEOM(pull->eGeom));
759     pr_ivec(fp, indent, "pull-dim", pull->dim, DIM, TRUE);
760     PR("pull-r1", pull->cyl_r1);
761     PR("pull-r0", pull->cyl_r0);
762     PR("pull-constr-tol", pull->constr_tol);
763     PS("pull-bPrintRef", EBOOL(pull->bPrintRef));
764     PI("pull-nstxout", pull->nstxout);
765     PI("pull-nstfout", pull->nstfout);
766     PI("pull-ngroup", pull->ngroup);
767     for (g = 0; g < pull->ngroup; g++)
768     {
769         pr_pull_group(fp, indent, g, &pull->group[g]);
770     }
771     PI("pull-ncoord", pull->ncoord);
772     for (g = 0; g < pull->ncoord; g++)
773     {
774         pr_pull_coord(fp, indent, g, &pull->coord[g]);
775     }
776 }
777
778 static void pr_rotgrp(FILE *fp, int indent, int g, t_rotgrp *rotg)
779 {
780     pr_indent(fp, indent);
781     fprintf(fp, "rotation_group %d:\n", g);
782     indent += 2;
783     PS("type", EROTGEOM(rotg->eType));
784     PS("massw", EBOOL(rotg->bMassW));
785     pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
786     pr_rvecs(fp, indent, "x_ref", rotg->x_ref, rotg->nat);
787     pr_rvec(fp, indent, "vec", rotg->vec, DIM, TRUE);
788     pr_rvec(fp, indent, "pivot", rotg->pivot, DIM, TRUE);
789     PR("rate", rotg->rate);
790     PR("k", rotg->k);
791     PR("slab_dist", rotg->slab_dist);
792     PR("min_gaussian", rotg->min_gaussian);
793     PR("epsilon", rotg->eps);
794     PS("fit_method", EROTFIT(rotg->eFittype));
795     PI("potfitangle_nstep", rotg->PotAngle_nstep);
796     PR("potfitangle_step", rotg->PotAngle_step);
797 }
798
799 static void pr_rot(FILE *fp, int indent, t_rot *rot)
800 {
801     int g;
802
803     PI("rot_nstrout", rot->nstrout);
804     PI("rot_nstsout", rot->nstsout);
805     PI("rot_ngrp", rot->ngrp);
806     for (g = 0; g < rot->ngrp; g++)
807     {
808         pr_rotgrp(fp, indent, g, &rot->grp[g]);
809     }
810 }
811
812
813 static void pr_swap(FILE *fp, int indent, t_swapcoords *swap)
814 {
815     int  i, j;
816     char str[STRLEN];
817
818
819     PI("frequency", swap->nstswap);
820     for (j = 0; j < 2; j++)
821     {
822         sprintf(str, "nanions%c", j+'A');
823         PI(str, swap->nanions[j]);
824         sprintf(str, "ncations%c", j+'A');
825         PI(str, swap->ncations[j]);
826     }
827     PI("coupling_steps", swap->nAverage);
828     PR("threshold", swap->threshold);
829     for (j = 0; j < 2; j++)
830     {
831         sprintf(str, "splitgroup%d_massw", j);
832         PS(str, EBOOL(swap->massw_split[j]));
833         sprintf(str, "split atoms group %d", j);
834         pr_ivec_block(fp, indent, str, swap->ind_split[j], swap->nat_split[j], TRUE);
835     }
836     pr_ivec_block(fp, indent, "swap atoms", swap->ind, swap->nat, TRUE);
837     pr_ivec_block(fp, indent, "solvent atoms", swap->ind_sol, swap->nat_sol, TRUE);
838     PR("cyl0_radius", swap->cyl0r);
839     PR("cyl0_upper", swap->cyl0u);
840     PR("cyl0_lower", swap->cyl0l);
841     PR("cyl1_radius", swap->cyl1r);
842     PR("cyl1_upper", swap->cyl1u);
843     PR("cyl1_lower", swap->cyl1l);
844 }
845
846
847 static void pr_imd(FILE *fp, int indent, t_IMD *imd)
848 {
849     PI("IMD_atoms", imd->nat);
850     pr_ivec_block(fp, indent, "atom", imd->ind, imd->nat, TRUE);
851 }
852
853
854 void pr_inputrec(FILE *fp, int indent, const char *title, t_inputrec *ir,
855                  gmx_bool bMDPformat)
856 {
857     const char *infbuf = "inf";
858     int         i;
859
860     if (available(fp, ir, indent, title))
861     {
862         if (!bMDPformat)
863         {
864             indent = pr_title(fp, indent, title);
865         }
866         /* This strings do not all have a direct correspondence to
867            .mdp entries, but we should follow the same convention of
868            using hyphens in the names users read and write. */
869         PS("integrator", EI(ir->eI));
870         PSTEP("nsteps", ir->nsteps);
871         PSTEP("init-step", ir->init_step);
872         PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
873         PS("ns-type", ENS(ir->ns_type));
874         PI("nstlist", ir->nstlist);
875         PI("ndelta", ir->ndelta);
876         PI("nstcomm", ir->nstcomm);
877         PS("comm-mode", ECOM(ir->comm_mode));
878         PI("nstlog", ir->nstlog);
879         PI("nstxout", ir->nstxout);
880         PI("nstvout", ir->nstvout);
881         PI("nstfout", ir->nstfout);
882         PI("nstcalcenergy", ir->nstcalcenergy);
883         PI("nstenergy", ir->nstenergy);
884         PI("nstxout-compressed", ir->nstxout_compressed);
885         PR("init-t", ir->init_t);
886         PR("delta-t", ir->delta_t);
887
888         PR("x-compression-precision", ir->x_compression_precision);
889         PR("fourierspacing", ir->fourier_spacing);
890         PI("nkx", ir->nkx);
891         PI("nky", ir->nky);
892         PI("nkz", ir->nkz);
893         PI("pme-order", ir->pme_order);
894         PR("ewald-rtol", ir->ewald_rtol);
895         PR("ewald-rtol-lj", ir->ewald_rtol_lj);
896         PR("ewald-geometry", ir->ewald_geometry);
897         PR("epsilon-surface", ir->epsilon_surface);
898         PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir->ljpme_combination_rule));
899         PS("ePBC", EPBC(ir->ePBC));
900         PS("bPeriodicMols", EBOOL(ir->bPeriodicMols));
901         PS("bContinuation", EBOOL(ir->bContinuation));
902         PS("bShakeSOR", EBOOL(ir->bShakeSOR));
903         PS("etc", ETCOUPLTYPE(ir->etc));
904         PS("bPrintNHChains", EBOOL(ir->bPrintNHChains));
905         PI("nsttcouple", ir->nsttcouple);
906         PS("epc", EPCOUPLTYPE(ir->epc));
907         PS("epctype", EPCOUPLTYPETYPE(ir->epct));
908         PI("nstpcouple", ir->nstpcouple);
909         PR("tau-p", ir->tau_p);
910         pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
911         pr_matrix(fp, indent, "compress", ir->compress, bMDPformat);
912         PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
913         if (bMDPformat)
914         {
915             fprintf(fp, "posres-com  = %g %g %g\n", ir->posres_com[XX],
916                     ir->posres_com[YY], ir->posres_com[ZZ]);
917         }
918         else
919         {
920             pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
921         }
922         if (bMDPformat)
923         {
924             fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX],
925                     ir->posres_comB[YY], ir->posres_comB[ZZ]);
926         }
927         else
928         {
929             pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
930         }
931         PR("verlet-buffer-tolerance", ir->verletbuf_tol);
932         PR("rlist", ir->rlist);
933         PR("rlistlong", ir->rlistlong);
934         PR("nstcalclr", ir->nstcalclr);
935         PR("rtpi", ir->rtpi);
936         PS("coulombtype", EELTYPE(ir->coulombtype));
937         PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
938         PR("rcoulomb-switch", ir->rcoulomb_switch);
939         PR("rcoulomb", ir->rcoulomb);
940         PS("vdwtype", EVDWTYPE(ir->vdwtype));
941         PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
942         PR("rvdw-switch", ir->rvdw_switch);
943         PR("rvdw", ir->rvdw);
944         if (ir->epsilon_r != 0)
945         {
946             PR("epsilon-r", ir->epsilon_r);
947         }
948         else
949         {
950             PS("epsilon-r", infbuf);
951         }
952         if (ir->epsilon_rf != 0)
953         {
954             PR("epsilon-rf", ir->epsilon_rf);
955         }
956         else
957         {
958             PS("epsilon-rf", infbuf);
959         }
960         PR("tabext", ir->tabext);
961         PS("implicit-solvent", EIMPLICITSOL(ir->implicit_solvent));
962         PS("gb-algorithm", EGBALGORITHM(ir->gb_algorithm));
963         PR("gb-epsilon-solvent", ir->gb_epsilon_solvent);
964         PI("nstgbradii", ir->nstgbradii);
965         PR("rgbradii", ir->rgbradii);
966         PR("gb-saltconc", ir->gb_saltconc);
967         PR("gb-obc-alpha", ir->gb_obc_alpha);
968         PR("gb-obc-beta", ir->gb_obc_beta);
969         PR("gb-obc-gamma", ir->gb_obc_gamma);
970         PR("gb-dielectric-offset", ir->gb_dielectric_offset);
971         PS("sa-algorithm", ESAALGORITHM(ir->gb_algorithm));
972         PR("sa-surface-tension", ir->sa_surface_tension);
973         PS("DispCorr", EDISPCORR(ir->eDispCorr));
974         PS("bSimTemp", EBOOL(ir->bSimTemp));
975         if (ir->bSimTemp)
976         {
977             pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
978         }
979         PS("free-energy", EFEPTYPE(ir->efep));
980         if (ir->efep != efepNO || ir->bSimTemp)
981         {
982             pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
983         }
984         if (ir->bExpanded)
985         {
986             pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
987         }
988
989         PI("nwall", ir->nwall);
990         PS("wall-type", EWALLTYPE(ir->wall_type));
991         PI("wall-atomtype[0]", ir->wall_atomtype[0]);
992         PI("wall-atomtype[1]", ir->wall_atomtype[1]);
993         PR("wall-density[0]", ir->wall_density[0]);
994         PR("wall-density[1]", ir->wall_density[1]);
995         PR("wall-ewald-zfac", ir->wall_ewald_zfac);
996
997         PS("pull", EPULLTYPE(ir->ePull));
998         if (ir->ePull != epullNO)
999         {
1000             pr_pull(fp, indent, ir->pull);
1001         }
1002
1003         PS("rotation", EBOOL(ir->bRot));
1004         if (ir->bRot)
1005         {
1006             pr_rot(fp, indent, ir->rot);
1007         }
1008
1009         PS("interactiveMD", EBOOL(ir->bIMD));
1010         if (ir->bIMD)
1011         {
1012             pr_imd(fp, indent, ir->imd);
1013         }
1014
1015         PS("disre", EDISRETYPE(ir->eDisre));
1016         PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
1017         PS("disre-mixed", EBOOL(ir->bDisreMixed));
1018         PR("dr-fc", ir->dr_fc);
1019         PR("dr-tau", ir->dr_tau);
1020         PR("nstdisreout", ir->nstdisreout);
1021         PR("orires-fc", ir->orires_fc);
1022         PR("orires-tau", ir->orires_tau);
1023         PR("nstorireout", ir->nstorireout);
1024
1025         PR("em-stepsize", ir->em_stepsize);
1026         PR("em-tol", ir->em_tol);
1027         PI("niter", ir->niter);
1028         PR("fc-stepsize", ir->fc_stepsize);
1029         PI("nstcgsteep", ir->nstcgsteep);
1030         PI("nbfgscorr", ir->nbfgscorr);
1031
1032         PS("ConstAlg", ECONSTRTYPE(ir->eConstrAlg));
1033         PR("shake-tol", ir->shake_tol);
1034         PI("lincs-order", ir->nProjOrder);
1035         PR("lincs-warnangle", ir->LincsWarnAngle);
1036         PI("lincs-iter", ir->nLincsIter);
1037         PR("bd-fric", ir->bd_fric);
1038         PSTEP("ld-seed", ir->ld_seed);
1039         PR("cos-accel", ir->cos_accel);
1040         pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
1041
1042         PS("adress", EBOOL(ir->bAdress));
1043         if (ir->bAdress)
1044         {
1045             PS("adress-type", EADRESSTYPE(ir->adress->type));
1046             PR("adress-const-wf", ir->adress->const_wf);
1047             PR("adress-ex-width", ir->adress->ex_width);
1048             PR("adress-hy-width", ir->adress->hy_width);
1049             PS("adress-interface-correction", EADRESSICTYPE(ir->adress->icor));
1050             PS("adress-site", EADRESSSITETYPE(ir->adress->site));
1051             PR("adress-ex-force-cap", ir->adress->ex_forcecap);
1052             PS("adress-do-hybridpairs", EBOOL(ir->adress->do_hybridpairs));
1053
1054             pr_rvec(fp, indent, "adress-reference-coords", ir->adress->refs, DIM, TRUE);
1055         }
1056         PI("userint1", ir->userint1);
1057         PI("userint2", ir->userint2);
1058         PI("userint3", ir->userint3);
1059         PI("userint4", ir->userint4);
1060         PR("userreal1", ir->userreal1);
1061         PR("userreal2", ir->userreal2);
1062         PR("userreal3", ir->userreal3);
1063         PR("userreal4", ir->userreal4);
1064         pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1065         pr_cosine(fp, indent, "efield-x", &(ir->ex[XX]), bMDPformat);
1066         pr_cosine(fp, indent, "efield-xt", &(ir->et[XX]), bMDPformat);
1067         pr_cosine(fp, indent, "efield-y", &(ir->ex[YY]), bMDPformat);
1068         pr_cosine(fp, indent, "efield-yt", &(ir->et[YY]), bMDPformat);
1069         pr_cosine(fp, indent, "efield-z", &(ir->ex[ZZ]), bMDPformat);
1070         pr_cosine(fp, indent, "efield-zt", &(ir->et[ZZ]), bMDPformat);
1071         PS("eSwapCoords", ESWAPTYPE(ir->eSwapCoords));
1072         if (ir->eSwapCoords != eswapNO)
1073         {
1074             pr_swap(fp, indent, ir->swap);
1075         }
1076         PS("bQMMM", EBOOL(ir->bQMMM));
1077         PI("QMconstraints", ir->QMconstraints);
1078         PI("QMMMscheme", ir->QMMMscheme);
1079         PR("scalefactor", ir->scalefactor);
1080         pr_qm_opts(fp, indent, "qm-opts", &(ir->opts));
1081     }
1082 }
1083 #undef PS
1084 #undef PR
1085 #undef PI
1086
1087 static void pr_harm(FILE *fp, t_iparams *iparams, const char *r, const char *kr)
1088 {
1089     fprintf(fp, "%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
1090             r, iparams->harmonic.rA, kr, iparams->harmonic.krA,
1091             r, iparams->harmonic.rB, kr, iparams->harmonic.krB);
1092 }
1093
1094 void pr_iparams(FILE *fp, t_functype ftype, t_iparams *iparams)
1095 {
1096     int  i;
1097     real VA[4], VB[4], *rbcA, *rbcB;
1098
1099     switch (ftype)
1100     {
1101         case F_ANGLES:
1102         case F_G96ANGLES:
1103             pr_harm(fp, iparams, "th", "ct");
1104             break;
1105         case F_CROSS_BOND_BONDS:
1106             fprintf(fp, "r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
1107                     iparams->cross_bb.r1e, iparams->cross_bb.r2e,
1108                     iparams->cross_bb.krr);
1109             break;
1110         case F_CROSS_BOND_ANGLES:
1111             fprintf(fp, "r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
1112                     iparams->cross_ba.r1e, iparams->cross_ba.r2e,
1113                     iparams->cross_ba.r3e, iparams->cross_ba.krt);
1114             break;
1115         case F_LINEAR_ANGLES:
1116             fprintf(fp, "klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
1117                     iparams->linangle.klinA, iparams->linangle.aA,
1118                     iparams->linangle.klinB, iparams->linangle.aB);
1119             break;
1120         case F_UREY_BRADLEY:
1121             fprintf(fp, "thetaA=%15.8e, kthetaA=%15.8e, r13A=%15.8e, kUBA=%15.8e, thetaB=%15.8e, kthetaB=%15.8e, r13B=%15.8e, kUBB=%15.8e\n", iparams->u_b.thetaA, iparams->u_b.kthetaA, iparams->u_b.r13A, iparams->u_b.kUBA, iparams->u_b.thetaB, iparams->u_b.kthetaB, iparams->u_b.r13B, iparams->u_b.kUBB);
1122             break;
1123         case F_QUARTIC_ANGLES:
1124             fprintf(fp, "theta=%15.8e", iparams->qangle.theta);
1125             for (i = 0; i < 5; i++)
1126             {
1127                 fprintf(fp, ", c%c=%15.8e", '0'+i, iparams->qangle.c[i]);
1128             }
1129             fprintf(fp, "\n");
1130             break;
1131         case F_BHAM:
1132             fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1133                     iparams->bham.a, iparams->bham.b, iparams->bham.c);
1134             break;
1135         case F_BONDS:
1136         case F_G96BONDS:
1137         case F_HARMONIC:
1138             pr_harm(fp, iparams, "b0", "cb");
1139             break;
1140         case F_IDIHS:
1141             pr_harm(fp, iparams, "xi", "cx");
1142             break;
1143         case F_MORSE:
1144             fprintf(fp, "b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
1145                     iparams->morse.b0A, iparams->morse.cbA, iparams->morse.betaA,
1146                     iparams->morse.b0B, iparams->morse.cbB, iparams->morse.betaB);
1147             break;
1148         case F_CUBICBONDS:
1149             fprintf(fp, "b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
1150                     iparams->cubic.b0, iparams->cubic.kb, iparams->cubic.kcub);
1151             break;
1152         case F_CONNBONDS:
1153             fprintf(fp, "\n");
1154             break;
1155         case F_FENEBONDS:
1156             fprintf(fp, "bm=%15.8e, kb=%15.8e\n", iparams->fene.bm, iparams->fene.kb);
1157             break;
1158         case F_RESTRBONDS:
1159             fprintf(fp, "lowA=%15.8e, up1A=%15.8e, up2A=%15.8e, kA=%15.8e, lowB=%15.8e, up1B=%15.8e, up2B=%15.8e, kB=%15.8e,\n",
1160                     iparams->restraint.lowA, iparams->restraint.up1A,
1161                     iparams->restraint.up2A, iparams->restraint.kA,
1162                     iparams->restraint.lowB, iparams->restraint.up1B,
1163                     iparams->restraint.up2B, iparams->restraint.kB);
1164             break;
1165         case F_TABBONDS:
1166         case F_TABBONDSNC:
1167         case F_TABANGLES:
1168         case F_TABDIHS:
1169             fprintf(fp, "tab=%d, kA=%15.8e, kB=%15.8e\n",
1170                     iparams->tab.table, iparams->tab.kA, iparams->tab.kB);
1171             break;
1172         case F_POLARIZATION:
1173             fprintf(fp, "alpha=%15.8e\n", iparams->polarize.alpha);
1174             break;
1175         case F_ANHARM_POL:
1176             fprintf(fp, "alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
1177                     iparams->anharm_polarize.alpha,
1178                     iparams->anharm_polarize.drcut,
1179                     iparams->anharm_polarize.khyp);
1180             break;
1181         case F_THOLE_POL:
1182             fprintf(fp, "a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
1183                     iparams->thole.a, iparams->thole.alpha1, iparams->thole.alpha2,
1184                     iparams->thole.rfac);
1185             break;
1186         case F_WATER_POL:
1187             fprintf(fp, "al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
1188                     iparams->wpol.al_x, iparams->wpol.al_y, iparams->wpol.al_z,
1189                     iparams->wpol.rOH, iparams->wpol.rHH, iparams->wpol.rOD);
1190             break;
1191         case F_LJ:
1192             fprintf(fp, "c6=%15.8e, c12=%15.8e\n", iparams->lj.c6, iparams->lj.c12);
1193             break;
1194         case F_LJ14:
1195             fprintf(fp, "c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
1196                     iparams->lj14.c6A, iparams->lj14.c12A,
1197                     iparams->lj14.c6B, iparams->lj14.c12B);
1198             break;
1199         case F_LJC14_Q:
1200             fprintf(fp, "fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1201                     iparams->ljc14.fqq,
1202                     iparams->ljc14.qi, iparams->ljc14.qj,
1203                     iparams->ljc14.c6, iparams->ljc14.c12);
1204             break;
1205         case F_LJC_PAIRS_NB:
1206             fprintf(fp, "qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1207                     iparams->ljcnb.qi, iparams->ljcnb.qj,
1208                     iparams->ljcnb.c6, iparams->ljcnb.c12);
1209             break;
1210         case F_PDIHS:
1211         case F_PIDIHS:
1212         case F_ANGRES:
1213         case F_ANGRESZ:
1214             fprintf(fp, "phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
1215                     iparams->pdihs.phiA, iparams->pdihs.cpA,
1216                     iparams->pdihs.phiB, iparams->pdihs.cpB,
1217                     iparams->pdihs.mult);
1218             break;
1219         case F_DISRES:
1220             fprintf(fp, "label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
1221                     iparams->disres.label, iparams->disres.type,
1222                     iparams->disres.low, iparams->disres.up1,
1223                     iparams->disres.up2, iparams->disres.kfac);
1224             break;
1225         case F_ORIRES:
1226             fprintf(fp, "ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
1227                     iparams->orires.ex, iparams->orires.label, iparams->orires.power,
1228                     iparams->orires.c, iparams->orires.obs, iparams->orires.kfac);
1229             break;
1230         case F_DIHRES:
1231             fprintf(fp, "phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
1232                     iparams->dihres.phiA, iparams->dihres.dphiA, iparams->dihres.kfacA,
1233                     iparams->dihres.phiB, iparams->dihres.dphiB, iparams->dihres.kfacB);
1234             break;
1235         case F_POSRES:
1236             fprintf(fp, "pos0A=(%15.8e,%15.8e,%15.8e), fcA=(%15.8e,%15.8e,%15.8e), pos0B=(%15.8e,%15.8e,%15.8e), fcB=(%15.8e,%15.8e,%15.8e)\n",
1237                     iparams->posres.pos0A[XX], iparams->posres.pos0A[YY],
1238                     iparams->posres.pos0A[ZZ], iparams->posres.fcA[XX],
1239                     iparams->posres.fcA[YY], iparams->posres.fcA[ZZ],
1240                     iparams->posres.pos0B[XX], iparams->posres.pos0B[YY],
1241                     iparams->posres.pos0B[ZZ], iparams->posres.fcB[XX],
1242                     iparams->posres.fcB[YY], iparams->posres.fcB[ZZ]);
1243             break;
1244         case F_FBPOSRES:
1245             fprintf(fp, "pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e\n",
1246                     iparams->fbposres.pos0[XX], iparams->fbposres.pos0[YY],
1247                     iparams->fbposres.pos0[ZZ], iparams->fbposres.geom,
1248                     iparams->fbposres.r,        iparams->fbposres.k);
1249             break;
1250         case F_RBDIHS:
1251             for (i = 0; i < NR_RBDIHS; i++)
1252             {
1253                 fprintf(fp, "%srbcA[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcA[i]);
1254             }
1255             fprintf(fp, "\n");
1256             for (i = 0; i < NR_RBDIHS; i++)
1257             {
1258                 fprintf(fp, "%srbcB[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcB[i]);
1259             }
1260             fprintf(fp, "\n");
1261             break;
1262         case F_FOURDIHS:
1263             /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
1264              * OPLS potential constants back.
1265              */
1266             rbcA = iparams->rbdihs.rbcA;
1267             rbcB = iparams->rbdihs.rbcB;
1268
1269             VA[3] = -0.25*rbcA[4];
1270             VA[2] = -0.5*rbcA[3];
1271             VA[1] = 4.0*VA[3]-rbcA[2];
1272             VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1273
1274             VB[3] = -0.25*rbcB[4];
1275             VB[2] = -0.5*rbcB[3];
1276             VB[1] = 4.0*VB[3]-rbcB[2];
1277             VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1278
1279             for (i = 0; i < NR_FOURDIHS; i++)
1280             {
1281                 fprintf(fp, "%sFourA[%d]=%15.8e", i == 0 ? "" : ", ", i, VA[i]);
1282             }
1283             fprintf(fp, "\n");
1284             for (i = 0; i < NR_FOURDIHS; i++)
1285             {
1286                 fprintf(fp, "%sFourB[%d]=%15.8e", i == 0 ? "" : ", ", i, VB[i]);
1287             }
1288             fprintf(fp, "\n");
1289             break;
1290
1291         case F_CONSTR:
1292         case F_CONSTRNC:
1293             fprintf(fp, "dA=%15.8e, dB=%15.8e\n", iparams->constr.dA, iparams->constr.dB);
1294             break;
1295         case F_SETTLE:
1296             fprintf(fp, "doh=%15.8e, dhh=%15.8e\n", iparams->settle.doh,
1297                     iparams->settle.dhh);
1298             break;
1299         case F_VSITE2:
1300             fprintf(fp, "a=%15.8e\n", iparams->vsite.a);
1301             break;
1302         case F_VSITE3:
1303         case F_VSITE3FD:
1304         case F_VSITE3FAD:
1305             fprintf(fp, "a=%15.8e, b=%15.8e\n", iparams->vsite.a, iparams->vsite.b);
1306             break;
1307         case F_VSITE3OUT:
1308         case F_VSITE4FD:
1309         case F_VSITE4FDN:
1310             fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1311                     iparams->vsite.a, iparams->vsite.b, iparams->vsite.c);
1312             break;
1313         case F_VSITEN:
1314             fprintf(fp, "n=%2d, a=%15.8e\n", iparams->vsiten.n, iparams->vsiten.a);
1315             break;
1316         case F_GB12:
1317         case F_GB13:
1318         case F_GB14:
1319             fprintf(fp, "sar=%15.8e, st=%15.8e, pi=%15.8e, gbr=%15.8e, bmlt=%15.8e\n", iparams->gb.sar, iparams->gb.st, iparams->gb.pi, iparams->gb.gbr, iparams->gb.bmlt);
1320             break;
1321         case F_CMAP:
1322             fprintf(fp, "cmapA=%1d, cmapB=%1d\n", iparams->cmap.cmapA, iparams->cmap.cmapB);
1323             break;
1324         case  F_RESTRANGLES:
1325             pr_harm(fp, iparams, "ktheta", "costheta0");
1326             break;
1327         case  F_RESTRDIHS:
1328             fprintf(fp, "phiA=%15.8e, cpA=%15.8e",
1329                     iparams->pdihs.phiA, iparams->pdihs.cpA);
1330             break;
1331         case  F_CBTDIHS:
1332             fprintf(fp, "kphi=%15.8e", iparams->cbtdihs.cbtcA[0]);
1333             for (i = 1; i < NR_CBTDIHS; i++)
1334             {
1335                 fprintf(fp, ", cbtcA[%d]=%15.8e", i-1, iparams->cbtdihs.cbtcA[i]);
1336             }
1337             fprintf(fp, "\n");
1338             break;
1339         default:
1340             gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1341                       ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1342     }
1343 }
1344
1345 void pr_ilist(FILE *fp, int indent, const char *title,
1346               t_functype *functype, t_ilist *ilist, gmx_bool bShowNumbers)
1347 {
1348     int      i, j, k, type, ftype;
1349     t_iatom *iatoms;
1350
1351     if (available(fp, ilist, indent, title) && ilist->nr > 0)
1352     {
1353         indent = pr_title(fp, indent, title);
1354         (void) pr_indent(fp, indent);
1355         fprintf(fp, "nr: %d\n", ilist->nr);
1356         if (ilist->nr > 0)
1357         {
1358             (void) pr_indent(fp, indent);
1359             fprintf(fp, "iatoms:\n");
1360             iatoms = ilist->iatoms;
1361             for (i = j = 0; i < ilist->nr; )
1362             {
1363 #ifndef DEBUG
1364                 (void) pr_indent(fp, indent+INDENT);
1365                 type  = *(iatoms++);
1366                 ftype = functype[type];
1367                 (void) fprintf(fp, "%d type=%d (%s)",
1368                                bShowNumbers ? j : -1, bShowNumbers ? type : -1,
1369                                interaction_function[ftype].name);
1370                 j++;
1371                 for (k = 0; k < interaction_function[ftype].nratoms; k++)
1372                 {
1373                     (void) fprintf(fp, " %u", *(iatoms++));
1374                 }
1375                 (void) fprintf(fp, "\n");
1376                 i += 1+interaction_function[ftype].nratoms;
1377 #else
1378                 fprintf(fp, "%5d%5d\n", i, iatoms[i]);
1379                 i++;
1380 #endif
1381             }
1382         }
1383     }
1384 }
1385
1386 static void pr_cmap(FILE *fp, int indent, const char *title,
1387                     gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1388 {
1389     int  i, j, nelem;
1390     real dx, idx;
1391
1392     dx    = 360.0 / cmap_grid->grid_spacing;
1393     nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1394
1395     if (available(fp, cmap_grid, indent, title))
1396     {
1397         fprintf(fp, "%s\n", title);
1398
1399         for (i = 0; i < cmap_grid->ngrid; i++)
1400         {
1401             idx = -180.0;
1402             fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
1403
1404             fprintf(fp, "grid[%3d]={\n", bShowNumbers ? i : -1);
1405
1406             for (j = 0; j < nelem; j++)
1407             {
1408                 if ( (j%cmap_grid->grid_spacing) == 0)
1409                 {
1410                     fprintf(fp, "%8.1f\n", idx);
1411                     idx += dx;
1412                 }
1413
1414                 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4]);
1415                 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+1]);
1416                 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+2]);
1417                 fprintf(fp, "%8.3f\n", cmap_grid->cmapdata[i].cmap[j*4+3]);
1418             }
1419             fprintf(fp, "\n");
1420         }
1421     }
1422
1423 }
1424
1425 void pr_ffparams(FILE *fp, int indent, const char *title,
1426                  gmx_ffparams_t *ffparams,
1427                  gmx_bool bShowNumbers)
1428 {
1429     int i, j;
1430
1431     indent = pr_title(fp, indent, title);
1432     (void) pr_indent(fp, indent);
1433     (void) fprintf(fp, "atnr=%d\n", ffparams->atnr);
1434     (void) pr_indent(fp, indent);
1435     (void) fprintf(fp, "ntypes=%d\n", ffparams->ntypes);
1436     for (i = 0; i < ffparams->ntypes; i++)
1437     {
1438         (void) pr_indent(fp, indent+INDENT);
1439         (void) fprintf(fp, "functype[%d]=%s, ",
1440                        bShowNumbers ? i : -1,
1441                        interaction_function[ffparams->functype[i]].name);
1442         pr_iparams(fp, ffparams->functype[i], &ffparams->iparams[i]);
1443     }
1444     (void) pr_double(fp, indent, "reppow", ffparams->reppow);
1445     (void) pr_real(fp, indent, "fudgeQQ", ffparams->fudgeQQ);
1446     pr_cmap(fp, indent, "cmap", &ffparams->cmap_grid, bShowNumbers);
1447 }
1448
1449 void pr_idef(FILE *fp, int indent, const char *title, t_idef *idef, gmx_bool bShowNumbers)
1450 {
1451     int i, j;
1452
1453     if (available(fp, idef, indent, title))
1454     {
1455         indent = pr_title(fp, indent, title);
1456         (void) pr_indent(fp, indent);
1457         (void) fprintf(fp, "atnr=%d\n", idef->atnr);
1458         (void) pr_indent(fp, indent);
1459         (void) fprintf(fp, "ntypes=%d\n", idef->ntypes);
1460         for (i = 0; i < idef->ntypes; i++)
1461         {
1462             (void) pr_indent(fp, indent+INDENT);
1463             (void) fprintf(fp, "functype[%d]=%s, ",
1464                            bShowNumbers ? i : -1,
1465                            interaction_function[idef->functype[i]].name);
1466             pr_iparams(fp, idef->functype[i], &idef->iparams[i]);
1467         }
1468         (void) pr_real(fp, indent, "fudgeQQ", idef->fudgeQQ);
1469
1470         for (j = 0; (j < F_NRE); j++)
1471         {
1472             pr_ilist(fp, indent, interaction_function[j].longname,
1473                      idef->functype, &idef->il[j], bShowNumbers);
1474         }
1475     }
1476 }
1477
1478 static int pr_block_title(FILE *fp, int indent, const char *title, t_block *block)
1479 {
1480     int i;
1481
1482     if (available(fp, block, indent, title))
1483     {
1484         indent = pr_title(fp, indent, title);
1485         (void) pr_indent(fp, indent);
1486         (void) fprintf(fp, "nr=%d\n", block->nr);
1487     }
1488     return indent;
1489 }
1490
1491 static int pr_blocka_title(FILE *fp, int indent, const char *title, t_blocka *block)
1492 {
1493     int i;
1494
1495     if (available(fp, block, indent, title))
1496     {
1497         indent = pr_title(fp, indent, title);
1498         (void) pr_indent(fp, indent);
1499         (void) fprintf(fp, "nr=%d\n", block->nr);
1500         (void) pr_indent(fp, indent);
1501         (void) fprintf(fp, "nra=%d\n", block->nra);
1502     }
1503     return indent;
1504 }
1505
1506 static void low_pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1507 {
1508     int i;
1509
1510     if (available(fp, block, indent, title))
1511     {
1512         indent = pr_blocka_title(fp, indent, title, block);
1513         for (i = 0; i <= block->nr; i++)
1514         {
1515             (void) pr_indent(fp, indent+INDENT);
1516             (void) fprintf(fp, "%s->index[%d]=%u\n",
1517                            title, bShowNumbers ? i : -1, block->index[i]);
1518         }
1519         for (i = 0; i < block->nra; i++)
1520         {
1521             (void) pr_indent(fp, indent+INDENT);
1522             (void) fprintf(fp, "%s->a[%d]=%u\n",
1523                            title, bShowNumbers ? i : -1, block->a[i]);
1524         }
1525     }
1526 }
1527
1528 void pr_block(FILE *fp, int indent, const char *title, t_block *block, gmx_bool bShowNumbers)
1529 {
1530     int i, j, ok, size, start, end;
1531
1532     if (available(fp, block, indent, title))
1533     {
1534         indent = pr_block_title(fp, indent, title, block);
1535         start  = 0;
1536         end    = start;
1537         if ((ok = (block->index[start] == 0)) == 0)
1538         {
1539             (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1540         }
1541         else
1542         {
1543             for (i = 0; i < block->nr; i++)
1544             {
1545                 end  = block->index[i+1];
1546                 size = pr_indent(fp, indent);
1547                 if (end <= start)
1548                 {
1549                     size += fprintf(fp, "%s[%d]={}\n", title, i);
1550                 }
1551                 else
1552                 {
1553                     size += fprintf(fp, "%s[%d]={%d..%d}\n",
1554                                     title, bShowNumbers ? i : -1,
1555                                     bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1556                 }
1557                 start = end;
1558             }
1559         }
1560     }
1561 }
1562
1563 void pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1564 {
1565     int i, j, ok, size, start, end;
1566
1567     if (available(fp, block, indent, title))
1568     {
1569         indent = pr_blocka_title(fp, indent, title, block);
1570         start  = 0;
1571         end    = start;
1572         if ((ok = (block->index[start] == 0)) == 0)
1573         {
1574             (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1575         }
1576         else
1577         {
1578             for (i = 0; i < block->nr; i++)
1579             {
1580                 end  = block->index[i+1];
1581                 size = pr_indent(fp, indent);
1582                 if (end <= start)
1583                 {
1584                     size += fprintf(fp, "%s[%d]={", title, i);
1585                 }
1586                 else
1587                 {
1588                     size += fprintf(fp, "%s[%d][%d..%d]={",
1589                                     title, bShowNumbers ? i : -1,
1590                                     bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1591                 }
1592                 for (j = start; j < end; j++)
1593                 {
1594                     if (j > start)
1595                     {
1596                         size += fprintf(fp, ", ");
1597                     }
1598                     if ((size) > (USE_WIDTH))
1599                     {
1600                         (void) fprintf(fp, "\n");
1601                         size = pr_indent(fp, indent+INDENT);
1602                     }
1603                     size += fprintf(fp, "%u", block->a[j]);
1604                 }
1605                 (void) fprintf(fp, "}\n");
1606                 start = end;
1607             }
1608         }
1609         if ((end != block->nra) || (!ok))
1610         {
1611             (void) pr_indent(fp, indent);
1612             (void) fprintf(fp, "tables inconsistent, dumping complete tables:\n");
1613             low_pr_blocka(fp, indent, title, block, bShowNumbers);
1614         }
1615     }
1616 }
1617
1618 static void pr_strings(FILE *fp, int indent, const char *title, char ***nm, int n, gmx_bool bShowNumbers)
1619 {
1620     int i;
1621
1622     if (available(fp, nm, indent, title))
1623     {
1624         indent = pr_title_n(fp, indent, title, n);
1625         for (i = 0; i < n; i++)
1626         {
1627             (void) pr_indent(fp, indent);
1628             (void) fprintf(fp, "%s[%d]={name=\"%s\"}\n",
1629                            title, bShowNumbers ? i : -1, *(nm[i]));
1630         }
1631     }
1632 }
1633
1634 static void pr_strings2(FILE *fp, int indent, const char *title,
1635                         char ***nm, char ***nmB, int n, gmx_bool bShowNumbers)
1636 {
1637     int i;
1638
1639     if (available(fp, nm, indent, title))
1640     {
1641         indent = pr_title_n(fp, indent, title, n);
1642         for (i = 0; i < n; i++)
1643         {
1644             (void) pr_indent(fp, indent);
1645             (void) fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1646                            title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
1647         }
1648     }
1649 }
1650
1651 static void pr_resinfo(FILE *fp, int indent, const char *title, t_resinfo *resinfo, int n, gmx_bool bShowNumbers)
1652 {
1653     int i;
1654
1655     if (available(fp, resinfo, indent, title))
1656     {
1657         indent = pr_title_n(fp, indent, title, n);
1658         for (i = 0; i < n; i++)
1659         {
1660             (void) pr_indent(fp, indent);
1661             (void) fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1662                            title, bShowNumbers ? i : -1,
1663                            *(resinfo[i].name), resinfo[i].nr,
1664                            (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1665         }
1666     }
1667 }
1668
1669 static void pr_atom(FILE *fp, int indent, const char *title, t_atom *atom, int n)
1670 {
1671     int i, j;
1672
1673     if (available(fp, atom, indent, title))
1674     {
1675         indent = pr_title_n(fp, indent, title, n);
1676         for (i = 0; i < n; i++)
1677         {
1678             (void) pr_indent(fp, indent);
1679             fprintf(fp, "%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1680                     "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1681                     title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype],
1682                     atom[i].m, atom[i].q, atom[i].mB, atom[i].qB,
1683                     atom[i].resind, atom[i].atomnumber);
1684         }
1685     }
1686 }
1687
1688 static void pr_grps(FILE *fp, const char *title, t_grps grps[], char **grpname[])
1689 {
1690     int i, j;
1691
1692     for (i = 0; (i < egcNR); i++)
1693     {
1694         fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
1695         for (j = 0; (j < grps[i].nr); j++)
1696         {
1697             fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
1698         }
1699         fprintf(fp, "]\n");
1700     }
1701 }
1702
1703 static void pr_groups(FILE *fp, int indent,
1704                       gmx_groups_t *groups,
1705                       gmx_bool bShowNumbers)
1706 {
1707     int grpnr[egcNR];
1708     int nat_max, i, g;
1709
1710     pr_grps(fp, "grp", groups->grps, groups->grpname);
1711     pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
1712
1713     (void) pr_indent(fp, indent);
1714     fprintf(fp, "groups          ");
1715     for (g = 0; g < egcNR; g++)
1716     {
1717         printf(" %5.5s", gtypes[g]);
1718     }
1719     printf("\n");
1720
1721     (void) pr_indent(fp, indent);
1722     fprintf(fp, "allocated       ");
1723     nat_max = 0;
1724     for (g = 0; g < egcNR; g++)
1725     {
1726         printf(" %5d", groups->ngrpnr[g]);
1727         nat_max = max(nat_max, groups->ngrpnr[g]);
1728     }
1729     printf("\n");
1730
1731     if (nat_max == 0)
1732     {
1733         (void) pr_indent(fp, indent);
1734         fprintf(fp, "groupnr[%5s] =", "*");
1735         for (g = 0; g < egcNR; g++)
1736         {
1737             fprintf(fp, "  %3d ", 0);
1738         }
1739         fprintf(fp, "\n");
1740     }
1741     else
1742     {
1743         for (i = 0; i < nat_max; i++)
1744         {
1745             (void) pr_indent(fp, indent);
1746             fprintf(fp, "groupnr[%5d] =", i);
1747             for (g = 0; g < egcNR; g++)
1748             {
1749                 fprintf(fp, "  %3d ",
1750                         groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1751             }
1752             fprintf(fp, "\n");
1753         }
1754     }
1755 }
1756
1757 void pr_atoms(FILE *fp, int indent, const char *title, t_atoms *atoms,
1758               gmx_bool bShownumbers)
1759 {
1760     if (available(fp, atoms, indent, title))
1761     {
1762         indent = pr_title(fp, indent, title);
1763         pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
1764         pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
1765         pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
1766         pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
1767     }
1768 }
1769
1770
1771 void pr_atomtypes(FILE *fp, int indent, const char *title, t_atomtypes *atomtypes,
1772                   gmx_bool bShowNumbers)
1773 {
1774     int i;
1775     if (available(fp, atomtypes, indent, title))
1776     {
1777         indent = pr_title(fp, indent, title);
1778         for (i = 0; i < atomtypes->nr; i++)
1779         {
1780             pr_indent(fp, indent);
1781             fprintf(fp,
1782                     "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1783                     bShowNumbers ? i : -1, atomtypes->radius[i], atomtypes->vol[i],
1784                     atomtypes->gb_radius[i],
1785                     atomtypes->surftens[i], atomtypes->atomnumber[i], atomtypes->S_hct[i]);
1786         }
1787     }
1788 }
1789
1790 static void pr_moltype(FILE *fp, int indent, const char *title,
1791                        gmx_moltype_t *molt, int n,
1792                        gmx_ffparams_t *ffparams,
1793                        gmx_bool bShowNumbers)
1794 {
1795     int j;
1796
1797     indent = pr_title_n(fp, indent, title, n);
1798     (void) pr_indent(fp, indent);
1799     (void) fprintf(fp, "name=\"%s\"\n", *(molt->name));
1800     pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
1801     pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
1802     pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
1803     for (j = 0; (j < F_NRE); j++)
1804     {
1805         pr_ilist(fp, indent, interaction_function[j].longname,
1806                  ffparams->functype, &molt->ilist[j], bShowNumbers);
1807     }
1808 }
1809
1810 static void pr_molblock(FILE *fp, int indent, const char *title,
1811                         gmx_molblock_t *molb, int n,
1812                         gmx_moltype_t *molt)
1813 {
1814     indent = pr_title_n(fp, indent, title, n);
1815     (void) pr_indent(fp, indent);
1816     (void) fprintf(fp, "%-20s = %d \"%s\"\n",
1817                    "moltype", molb->type, *(molt[molb->type].name));
1818     pr_int(fp, indent, "#molecules", molb->nmol);
1819     pr_int(fp, indent, "#atoms_mol", molb->natoms_mol);
1820     pr_int(fp, indent, "#posres_xA", molb->nposres_xA);
1821     if (molb->nposres_xA > 0)
1822     {
1823         pr_rvecs(fp, indent, "posres_xA", molb->posres_xA, molb->nposres_xA);
1824     }
1825     pr_int(fp, indent, "#posres_xB", molb->nposres_xB);
1826     if (molb->nposres_xB > 0)
1827     {
1828         pr_rvecs(fp, indent, "posres_xB", molb->posres_xB, molb->nposres_xB);
1829     }
1830 }
1831
1832 void pr_mtop(FILE *fp, int indent, const char *title, gmx_mtop_t *mtop,
1833              gmx_bool bShowNumbers)
1834 {
1835     int mt, mb;
1836
1837     if (available(fp, mtop, indent, title))
1838     {
1839         indent = pr_title(fp, indent, title);
1840         (void) pr_indent(fp, indent);
1841         (void) fprintf(fp, "name=\"%s\"\n", *(mtop->name));
1842         pr_int(fp, indent, "#atoms", mtop->natoms);
1843         pr_int(fp, indent, "#molblock", mtop->nmolblock);
1844         for (mb = 0; mb < mtop->nmolblock; mb++)
1845         {
1846             pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
1847         }
1848         pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
1849         pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
1850         for (mt = 0; mt < mtop->nmoltype; mt++)
1851         {
1852             pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
1853                        &mtop->ffparams, bShowNumbers);
1854         }
1855         pr_groups(fp, indent, &mtop->groups, bShowNumbers);
1856     }
1857 }
1858
1859 void pr_top(FILE *fp, int indent, const char *title, t_topology *top, gmx_bool bShowNumbers)
1860 {
1861     if (available(fp, top, indent, title))
1862     {
1863         indent = pr_title(fp, indent, title);
1864         (void) pr_indent(fp, indent);
1865         (void) fprintf(fp, "name=\"%s\"\n", *(top->name));
1866         pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
1867         pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
1868         pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
1869         pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
1870         pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
1871         pr_idef(fp, indent, "idef", &top->idef, bShowNumbers);
1872     }
1873 }
1874
1875 void pr_header(FILE *fp, int indent, const char *title, t_tpxheader *sh)
1876 {
1877     char buf[22];
1878
1879     if (available(fp, sh, indent, title))
1880     {
1881         indent = pr_title(fp, indent, title);
1882         pr_indent(fp, indent);
1883         fprintf(fp, "bIr    = %spresent\n", sh->bIr ? "" : "not ");
1884         pr_indent(fp, indent);
1885         fprintf(fp, "bBox   = %spresent\n", sh->bBox ? "" : "not ");
1886         pr_indent(fp, indent);
1887         fprintf(fp, "bTop   = %spresent\n", sh->bTop ? "" : "not ");
1888         pr_indent(fp, indent);
1889         fprintf(fp, "bX     = %spresent\n", sh->bX ? "" : "not ");
1890         pr_indent(fp, indent);
1891         fprintf(fp, "bV     = %spresent\n", sh->bV ? "" : "not ");
1892         pr_indent(fp, indent);
1893         fprintf(fp, "bF     = %spresent\n", sh->bF ? "" : "not ");
1894
1895         pr_indent(fp, indent);
1896         fprintf(fp, "natoms = %d\n", sh->natoms);
1897         pr_indent(fp, indent);
1898         fprintf(fp, "lambda = %e\n", sh->lambda);
1899     }
1900 }
1901
1902 void pr_commrec(FILE *fp, int indent, t_commrec *cr)
1903 {
1904     pr_indent(fp, indent);
1905     fprintf(fp, "commrec:\n");
1906     indent += 2;
1907     pr_indent(fp, indent);
1908     fprintf(fp, "nodeid    = %d\n", cr->nodeid);
1909     pr_indent(fp, indent);
1910     fprintf(fp, "nnodes    = %d\n", cr->nnodes);
1911     pr_indent(fp, indent);
1912     fprintf(fp, "npmenodes = %d\n", cr->npmenodes);
1913     /*
1914        pr_indent(fp,indent);
1915        fprintf(fp,"threadid  = %d\n",cr->threadid);
1916        pr_indent(fp,indent);
1917        fprintf(fp,"nthreads  = %d\n",cr->nthreads);
1918      */
1919 }