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