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