2748a97093af6f6e8dbaf68bcffc5d69735d2a7b
[alexxy/gromacs.git] / src / gromacs / gmxlib / txtdump.cpp
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 <cstdio>
44 #include <cstdlib>
45
46 #include <algorithm>
47
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/legacyheaders/names.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/legacyheaders/types/commrec.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/smalloc.h"
56
57 int pr_indent(FILE *fp, int n)
58 {
59     int i;
60
61     for (i = 0; i < n; i++)
62     {
63         fprintf(fp, " ");
64     }
65     return n;
66 }
67
68 int available(FILE *fp, void *p, int indent, const char *title)
69 {
70     if (!p)
71     {
72         if (indent > 0)
73         {
74             pr_indent(fp, indent);
75         }
76         fprintf(fp, "%s: not available\n", title);
77     }
78     return (p != NULL);
79 }
80
81 int pr_title(FILE *fp, int indent, const char *title)
82 {
83     pr_indent(fp, indent);
84     fprintf(fp, "%s:\n", title);
85     return (indent+INDENT);
86 }
87
88 int pr_title_n(FILE *fp, int indent, const char *title, int n)
89 {
90     pr_indent(fp, indent);
91     fprintf(fp, "%s (%d):\n", title, n);
92     return (indent+INDENT);
93 }
94
95 int pr_title_nxn(FILE *fp, int indent, const char *title, int n1, int n2)
96 {
97     pr_indent(fp, indent);
98     fprintf(fp, "%s (%dx%d):\n", title, n1, n2);
99     return (indent+INDENT);
100 }
101
102 void pr_ivec(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
103 {
104     int i;
105
106     if (available(fp, vec, indent, title))
107     {
108         indent = pr_title_n(fp, indent, title, n);
109         for (i = 0; i < n; i++)
110         {
111             pr_indent(fp, indent);
112             fprintf(fp, "%s[%d]=%d\n", title, bShowNumbers ? i : -1, vec[i]);
113         }
114     }
115 }
116
117 void pr_ivec_block(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
118 {
119     int i, j;
120
121     if (available(fp, vec, indent, title))
122     {
123         indent = pr_title_n(fp, indent, title, n);
124         i      = 0;
125         while (i < n)
126         {
127             j = i+1;
128             while (j < n && vec[j] == vec[j-1]+1)
129             {
130                 j++;
131             }
132             /* Print consecutive groups of 3 or more as blocks */
133             if (j - i < 3)
134             {
135                 while (i < j)
136                 {
137                     pr_indent(fp, indent);
138                     fprintf(fp, "%s[%d]=%d\n",
139                             title, bShowNumbers ? i : -1, vec[i]);
140                     i++;
141                 }
142             }
143             else
144             {
145                 pr_indent(fp, indent);
146                 fprintf(fp, "%s[%d,...,%d] = {%d,...,%d}\n",
147                         title,
148                         bShowNumbers ? i : -1,
149                         bShowNumbers ? j-1 : -1,
150                         vec[i], vec[j-1]);
151                 i = j;
152             }
153         }
154     }
155 }
156
157 void pr_bvec(FILE *fp, int indent, const char *title, gmx_bool vec[], int n, gmx_bool bShowNumbers)
158 {
159     int i;
160
161     if (available(fp, vec, indent, title))
162     {
163         indent = pr_title_n(fp, indent, title, n);
164         for (i = 0; i < n; i++)
165         {
166             pr_indent(fp, indent);
167             fprintf(fp, "%s[%d]=%s\n", title, bShowNumbers ? i : -1,
168                     EBOOL(vec[i]));
169         }
170     }
171 }
172
173 void pr_ivecs(FILE *fp, int indent, const char *title, ivec vec[], int n, gmx_bool bShowNumbers)
174 {
175     int i, j;
176
177     if (available(fp, vec, indent, title))
178     {
179         indent = pr_title_nxn(fp, indent, title, n, DIM);
180         for (i = 0; i < n; i++)
181         {
182             pr_indent(fp, indent);
183             fprintf(fp, "%s[%d]={", title, bShowNumbers ? i : -1);
184             for (j = 0; j < DIM; j++)
185             {
186                 if (j != 0)
187                 {
188                     fprintf(fp, ", ");
189                 }
190                 fprintf(fp, "%d", vec[i][j]);
191             }
192             fprintf(fp, "}\n");
193         }
194     }
195 }
196
197 void pr_rvec(FILE *fp, int indent, const char *title, real vec[], int n, gmx_bool bShowNumbers)
198 {
199     int i;
200
201     if (available(fp, vec, indent, title))
202     {
203         indent = pr_title_n(fp, indent, title, n);
204         for (i = 0; i < n; i++)
205         {
206             pr_indent(fp, indent);
207             fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
208         }
209     }
210 }
211
212 void pr_dvec(FILE *fp, int indent, const char *title, double vec[], int n, gmx_bool bShowNumbers)
213 {
214     int i;
215
216     if (available(fp, vec, indent, title))
217     {
218         indent = pr_title_n(fp, indent, title, n);
219         for (i = 0; i < n; i++)
220         {
221             pr_indent(fp, indent);
222             fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
223         }
224     }
225 }
226
227
228 /*
229    void pr_mat(FILE *fp,int indent,char *title,matrix m)
230    {
231    int i,j;
232
233    if (available(fp,m,indent,title)) {
234     indent=pr_title_n(fp,indent,title,n);
235     for(i=0; i<n; i++) {
236       pr_indent(fp,indent);
237       fprintf(fp,"%s[%d]=%12.5e %12.5e %12.5e\n",
238           title,bShowNumbers?i:-1,m[i][XX],m[i][YY],m[i][ZZ]);
239     }
240    }
241    }
242  */
243
244 void pr_rvecs_len(FILE *fp, int indent, const char *title, rvec vec[], int n)
245 {
246     int i, j;
247
248     if (available(fp, vec, indent, title))
249     {
250         indent = pr_title_nxn(fp, indent, title, n, DIM);
251         for (i = 0; i < n; i++)
252         {
253             pr_indent(fp, indent);
254             fprintf(fp, "%s[%5d]={", title, i);
255             for (j = 0; j < DIM; j++)
256             {
257                 if (j != 0)
258                 {
259                     fprintf(fp, ", ");
260                 }
261                 fprintf(fp, "%12.5e", vec[i][j]);
262             }
263             fprintf(fp, "} len=%12.5e\n", norm(vec[i]));
264         }
265     }
266 }
267
268 void pr_rvecs(FILE *fp, int indent, const char *title, rvec vec[], int n)
269 {
270     const char *fshort = "%12.5e";
271     const char *flong  = "%15.8e";
272     const char *format;
273     int         i, j;
274
275     if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
276     {
277         format = flong;
278     }
279     else
280     {
281         format = fshort;
282     }
283
284     if (available(fp, vec, indent, title))
285     {
286         indent = pr_title_nxn(fp, indent, title, n, DIM);
287         for (i = 0; i < n; i++)
288         {
289             pr_indent(fp, indent);
290             fprintf(fp, "%s[%5d]={", title, i);
291             for (j = 0; j < DIM; j++)
292             {
293                 if (j != 0)
294                 {
295                     fprintf(fp, ", ");
296                 }
297                 fprintf(fp, format, vec[i][j]);
298             }
299             fprintf(fp, "}\n");
300         }
301     }
302 }
303
304
305 void pr_rvecs_of_dim(FILE *fp, int indent, const char *title, rvec vec[], int n, int dim)
306 {
307     const char *fshort = "%12.5e";
308     const char *flong  = "%15.8e";
309     const char *format;
310     int         i, j;
311
312     if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
313     {
314         format = flong;
315     }
316     else
317     {
318         format = fshort;
319     }
320
321     if (available(fp, vec, indent, title))
322     {
323         indent = pr_title_nxn(fp, indent, title, n, dim);
324         for (i = 0; i < n; i++)
325         {
326             pr_indent(fp, indent);
327             fprintf(fp, "%s[%5d]={", title, i);
328             for (j = 0; j < dim; j++)
329             {
330                 if (j != 0)
331                 {
332                     fprintf(fp, ", ");
333                 }
334                 fprintf(fp, format, vec[i][j]);
335             }
336             fprintf(fp, "}\n");
337         }
338     }
339 }
340
341 void pr_reals(FILE *fp, int indent, const char *title, real *vec, int n)
342 {
343     int i;
344
345     if (available(fp, vec, indent, title))
346     {
347         pr_indent(fp, indent);
348         fprintf(fp, "%s:\t", title);
349         for (i = 0; i < n; i++)
350         {
351             fprintf(fp, "  %10g", vec[i]);
352         }
353         fprintf(fp, "\n");
354     }
355 }
356
357 void pr_doubles(FILE *fp, int indent, const char *title, double *vec, int n)
358 {
359     int i;
360
361     if (available(fp, vec, indent, title))
362     {
363         pr_indent(fp, indent);
364         fprintf(fp, "%s:\t", title);
365         for (i = 0; i < n; i++)
366         {
367             fprintf(fp, "  %10g", vec[i]);
368         }
369         fprintf(fp, "\n");
370     }
371 }
372
373 void pr_reals_of_dim(FILE *fp, int indent, const char *title, real *vec, int n, int dim)
374 {
375     int         i, j;
376     const char *fshort = "%12.5e";
377     const char *flong  = "%15.8e";
378     const char *format;
379
380     if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
381     {
382         format = flong;
383     }
384     else
385     {
386         format = fshort;
387     }
388
389     if (available(fp, vec, indent, title))
390     {
391         indent = pr_title_nxn(fp, indent, title, n, dim);
392         for (i = 0; i < n; i++)
393         {
394             pr_indent(fp, indent);
395             fprintf(fp, "%s[%5d]={", title, i);
396             for (j = 0; j < dim; j++)
397             {
398                 if (j != 0)
399                 {
400                     fprintf(fp, ", ");
401                 }
402                 fprintf(fp, format, vec[i * dim  + j]);
403             }
404             fprintf(fp, "}\n");
405         }
406     }
407 }
408
409 static void pr_int(FILE *fp, int indent, const char *title, int i)
410 {
411     pr_indent(fp, indent);
412     fprintf(fp, "%-30s = %d\n", title, i);
413 }
414
415 static void pr_int64(FILE *fp, int indent, const char *title, gmx_int64_t i)
416 {
417     char buf[STEPSTRSIZE];
418
419     pr_indent(fp, indent);
420     fprintf(fp, "%-30s = %s\n", title, gmx_step_str(i, buf));
421 }
422
423 static void pr_real(FILE *fp, int indent, const char *title, real r)
424 {
425     pr_indent(fp, indent);
426     fprintf(fp, "%-30s = %g\n", title, r);
427 }
428
429 static void pr_double(FILE *fp, int indent, const char *title, double d)
430 {
431     pr_indent(fp, indent);
432     fprintf(fp, "%-30s = %g\n", title, d);
433 }
434
435 static void pr_str(FILE *fp, int indent, const char *title, const char *s)
436 {
437     pr_indent(fp, indent);
438     fprintf(fp, "%-30s = %s\n", title, s);
439 }
440
441 void pr_qm_opts(FILE *fp, int indent, const char *title, t_grpopts *opts)
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         pr_indent(fp, indent);
595         fprintf(fp, "n = %d\n", cos->n);
596         if (cos->n > 0)
597         {
598             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             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  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
869     if (available(fp, ir, indent, title))
870     {
871         if (!bMDPformat)
872         {
873             indent = pr_title(fp, indent, title);
874         }
875         /* Try to make this list appear in the same order as the
876          * options are written in the default mdout.mdp, and with
877          * the same user-exposed names to facilitate debugging.
878          */
879         PS("integrator", EI(ir->eI));
880         PR("tinit", ir->init_t);
881         PR("dt", ir->delta_t);
882         PSTEP("nsteps", ir->nsteps);
883         PSTEP("init-step", ir->init_step);
884         PI("simulation-part", ir->simulation_part);
885         PS("comm-mode", ECOM(ir->comm_mode));
886         PI("nstcomm", ir->nstcomm);
887
888         /* Langevin dynamics */
889         PR("bd-fric", ir->bd_fric);
890         PSTEP("ld-seed", ir->ld_seed);
891
892         /* Energy minimization */
893         PR("emtol", ir->em_tol);
894         PR("emstep", ir->em_stepsize);
895         PI("niter", ir->niter);
896         PR("fcstep", ir->fc_stepsize);
897         PI("nstcgsteep", ir->nstcgsteep);
898         PI("nbfgscorr", ir->nbfgscorr);
899
900         /* Test particle insertion */
901         PR("rtpi", ir->rtpi);
902
903         /* Output control */
904         PI("nstxout", ir->nstxout);
905         PI("nstvout", ir->nstvout);
906         PI("nstfout", ir->nstfout);
907         PI("nstlog", ir->nstlog);
908         PI("nstcalcenergy", ir->nstcalcenergy);
909         PI("nstenergy", ir->nstenergy);
910         PI("nstxout-compressed", ir->nstxout_compressed);
911         PR("compressed-x-precision", ir->x_compression_precision);
912
913         /* Neighborsearching parameters */
914         PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
915         PI("nstlist", ir->nstlist);
916         PS("ns-type", ENS(ir->ns_type));
917         PS("pbc", EPBC(ir->ePBC));
918         PS("periodic-molecules", EBOOL(ir->bPeriodicMols));
919         PR("verlet-buffer-tolerance", ir->verletbuf_tol);
920         PR("rlist", ir->rlist);
921         PR("rlistlong", ir->rlistlong);
922         PR("nstcalclr", ir->nstcalclr);
923
924         /* Options for electrostatics and VdW */
925         PS("coulombtype", EELTYPE(ir->coulombtype));
926         PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
927         PR("rcoulomb-switch", ir->rcoulomb_switch);
928         PR("rcoulomb", ir->rcoulomb);
929         if (ir->epsilon_r != 0)
930         {
931             PR("epsilon-r", ir->epsilon_r);
932         }
933         else
934         {
935             PS("epsilon-r", infbuf);
936         }
937         if (ir->epsilon_rf != 0)
938         {
939             PR("epsilon-rf", ir->epsilon_rf);
940         }
941         else
942         {
943             PS("epsilon-rf", infbuf);
944         }
945         PS("vdw-type", EVDWTYPE(ir->vdwtype));
946         PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
947         PR("rvdw-switch", ir->rvdw_switch);
948         PR("rvdw", ir->rvdw);
949         PS("DispCorr", EDISPCORR(ir->eDispCorr));
950         PR("table-extension", ir->tabext);
951
952         PR("fourierspacing", ir->fourier_spacing);
953         PI("fourier-nx", ir->nkx);
954         PI("fourier-ny", ir->nky);
955         PI("fourier-nz", ir->nkz);
956         PI("pme-order", ir->pme_order);
957         PR("ewald-rtol", ir->ewald_rtol);
958         PR("ewald-rtol-lj", ir->ewald_rtol_lj);
959         PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir->ljpme_combination_rule));
960         PR("ewald-geometry", ir->ewald_geometry);
961         PR("epsilon-surface", ir->epsilon_surface);
962
963         /* Implicit solvent */
964         PS("implicit-solvent", EIMPLICITSOL(ir->implicit_solvent));
965
966         /* Generalized born electrostatics */
967         PS("gb-algorithm", EGBALGORITHM(ir->gb_algorithm));
968         PI("nstgbradii", ir->nstgbradii);
969         PR("rgbradii", ir->rgbradii);
970         PR("gb-epsilon-solvent", ir->gb_epsilon_solvent);
971         PR("gb-saltconc", ir->gb_saltconc);
972         PR("gb-obc-alpha", ir->gb_obc_alpha);
973         PR("gb-obc-beta", ir->gb_obc_beta);
974         PR("gb-obc-gamma", ir->gb_obc_gamma);
975         PR("gb-dielectric-offset", ir->gb_dielectric_offset);
976         PS("sa-algorithm", ESAALGORITHM(ir->sa_algorithm));
977         PR("sa-surface-tension", ir->sa_surface_tension);
978
979         /* Options for weak coupling algorithms */
980         PS("tcoupl", ETCOUPLTYPE(ir->etc));
981         PI("nsttcouple", ir->nsttcouple);
982         PI("nh-chain-length", ir->opts.nhchainlength);
983         PS("print-nose-hoover-chain-variables", EBOOL(ir->bPrintNHChains));
984
985         PS("pcoupl", EPCOUPLTYPE(ir->epc));
986         PS("pcoupltype", EPCOUPLTYPETYPE(ir->epct));
987         PI("nstpcouple", ir->nstpcouple);
988         PR("tau-p", ir->tau_p);
989         pr_matrix(fp, indent, "compressibility", ir->compress, bMDPformat);
990         pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
991         PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
992
993         if (bMDPformat)
994         {
995             fprintf(fp, "posres-com  = %g %g %g\n", ir->posres_com[XX],
996                     ir->posres_com[YY], ir->posres_com[ZZ]);
997             fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX],
998                     ir->posres_comB[YY], ir->posres_comB[ZZ]);
999         }
1000         else
1001         {
1002             pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
1003             pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
1004         }
1005
1006         /* QMMM */
1007         PS("QMMM", EBOOL(ir->bQMMM));
1008         PI("QMconstraints", ir->QMconstraints);
1009         PI("QMMMscheme", ir->QMMMscheme);
1010         PR("MMChargeScaleFactor", ir->scalefactor);
1011         pr_qm_opts(fp, indent, "qm-opts", &(ir->opts));
1012
1013         /* CONSTRAINT OPTIONS */
1014         PS("constraint-algorithm", ECONSTRTYPE(ir->eConstrAlg));
1015         PS("continuation", EBOOL(ir->bContinuation));
1016
1017         PS("Shake-SOR", EBOOL(ir->bShakeSOR));
1018         PR("shake-tol", ir->shake_tol);
1019         PI("lincs-order", ir->nProjOrder);
1020         PI("lincs-iter", ir->nLincsIter);
1021         PR("lincs-warnangle", ir->LincsWarnAngle);
1022
1023         /* Walls */
1024         PI("nwall", ir->nwall);
1025         PS("wall-type", EWALLTYPE(ir->wall_type));
1026         PR("wall-r-linpot", ir->wall_r_linpot);
1027         /* wall-atomtype */
1028         PI("wall-atomtype[0]", ir->wall_atomtype[0]);
1029         PI("wall-atomtype[1]", ir->wall_atomtype[1]);
1030         /* wall-density */
1031         PR("wall-density[0]", ir->wall_density[0]);
1032         PR("wall-density[1]", ir->wall_density[1]);
1033         PR("wall-ewald-zfac", ir->wall_ewald_zfac);
1034
1035         /* COM PULLING */
1036         PS("pull", EBOOL(ir->bPull));
1037         if (ir->bPull)
1038         {
1039             pr_pull(fp, indent, ir->pull);
1040         }
1041
1042         /* ENFORCED ROTATION */
1043         PS("rotation", EBOOL(ir->bRot));
1044         if (ir->bRot)
1045         {
1046             pr_rot(fp, indent, ir->rot);
1047         }
1048
1049         /* INTERACTIVE MD */
1050         PS("interactiveMD", EBOOL(ir->bIMD));
1051         if (ir->bIMD)
1052         {
1053             pr_imd(fp, indent, ir->imd);
1054         }
1055
1056         /* NMR refinement stuff */
1057         PS("disre", EDISRETYPE(ir->eDisre));
1058         PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
1059         PS("disre-mixed", EBOOL(ir->bDisreMixed));
1060         PR("dr-fc", ir->dr_fc);
1061         PR("dr-tau", ir->dr_tau);
1062         PR("nstdisreout", ir->nstdisreout);
1063
1064         PR("orire-fc", ir->orires_fc);
1065         PR("orire-tau", ir->orires_tau);
1066         PR("nstorireout", ir->nstorireout);
1067
1068         /* FREE ENERGY VARIABLES */
1069         PS("free-energy", EFEPTYPE(ir->efep));
1070         if (ir->efep != efepNO || ir->bSimTemp)
1071         {
1072             pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
1073         }
1074         if (ir->bExpanded)
1075         {
1076             pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
1077         }
1078
1079         /* NON-equilibrium MD stuff */
1080         PR("cos-acceleration", ir->cos_accel);
1081         pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
1082
1083         /* SIMULATED TEMPERING */
1084         PS("simulated-tempering", EBOOL(ir->bSimTemp));
1085         if (ir->bSimTemp)
1086         {
1087             pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
1088         }
1089
1090         /* ELECTRIC FIELDS */
1091         pr_cosine(fp, indent, "E-x", &(ir->ex[XX]), bMDPformat);
1092         pr_cosine(fp, indent, "E-xt", &(ir->et[XX]), bMDPformat);
1093         pr_cosine(fp, indent, "E-y", &(ir->ex[YY]), bMDPformat);
1094         pr_cosine(fp, indent, "E-yt", &(ir->et[YY]), bMDPformat);
1095         pr_cosine(fp, indent, "E-z", &(ir->ex[ZZ]), bMDPformat);
1096         pr_cosine(fp, indent, "E-zt", &(ir->et[ZZ]), bMDPformat);
1097
1098         /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
1099         PS("swapcoords", ESWAPTYPE(ir->eSwapCoords));
1100         if (ir->eSwapCoords != eswapNO)
1101         {
1102             pr_swap(fp, indent, ir->swap);
1103         }
1104
1105         /* AdResS PARAMETERS */
1106         PS("adress", EBOOL(ir->bAdress));
1107         if (ir->bAdress)
1108         {
1109             PS("adress-type", EADRESSTYPE(ir->adress->type));
1110             PR("adress-const-wf", ir->adress->const_wf);
1111             PR("adress-ex-width", ir->adress->ex_width);
1112             PR("adress-hy-width", ir->adress->hy_width);
1113             PR("adress-ex-forcecap", ir->adress->ex_forcecap);
1114             PS("adress-interface-correction", EADRESSICTYPE(ir->adress->icor));
1115             PS("adress-site", EADRESSSITETYPE(ir->adress->site));
1116             pr_rvec(fp, indent, "adress-reference-coords", ir->adress->refs, DIM, TRUE);
1117             PS("adress-do-hybridpairs", EBOOL(ir->adress->do_hybridpairs));
1118         }
1119
1120         /* USER-DEFINED THINGIES */
1121         PI("userint1", ir->userint1);
1122         PI("userint2", ir->userint2);
1123         PI("userint3", ir->userint3);
1124         PI("userint4", ir->userint4);
1125         PR("userreal1", ir->userreal1);
1126         PR("userreal2", ir->userreal2);
1127         PR("userreal3", ir->userreal3);
1128         PR("userreal4", ir->userreal4);
1129
1130         pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1131     }
1132 }
1133 #undef PS
1134 #undef PR
1135 #undef PI
1136
1137 static void pr_harm(FILE *fp, t_iparams *iparams, const char *r, const char *kr)
1138 {
1139     fprintf(fp, "%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
1140             r, iparams->harmonic.rA, kr, iparams->harmonic.krA,
1141             r, iparams->harmonic.rB, kr, iparams->harmonic.krB);
1142 }
1143
1144 void pr_iparams(FILE *fp, t_functype ftype, t_iparams *iparams)
1145 {
1146     int  i;
1147     real VA[4], VB[4], *rbcA, *rbcB;
1148
1149     switch (ftype)
1150     {
1151         case F_ANGLES:
1152         case F_G96ANGLES:
1153             pr_harm(fp, iparams, "th", "ct");
1154             break;
1155         case F_CROSS_BOND_BONDS:
1156             fprintf(fp, "r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
1157                     iparams->cross_bb.r1e, iparams->cross_bb.r2e,
1158                     iparams->cross_bb.krr);
1159             break;
1160         case F_CROSS_BOND_ANGLES:
1161             fprintf(fp, "r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
1162                     iparams->cross_ba.r1e, iparams->cross_ba.r2e,
1163                     iparams->cross_ba.r3e, iparams->cross_ba.krt);
1164             break;
1165         case F_LINEAR_ANGLES:
1166             fprintf(fp, "klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
1167                     iparams->linangle.klinA, iparams->linangle.aA,
1168                     iparams->linangle.klinB, iparams->linangle.aB);
1169             break;
1170         case F_UREY_BRADLEY:
1171             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);
1172             break;
1173         case F_QUARTIC_ANGLES:
1174             fprintf(fp, "theta=%15.8e", iparams->qangle.theta);
1175             for (i = 0; i < 5; i++)
1176             {
1177                 fprintf(fp, ", c%c=%15.8e", '0'+i, iparams->qangle.c[i]);
1178             }
1179             fprintf(fp, "\n");
1180             break;
1181         case F_BHAM:
1182             fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1183                     iparams->bham.a, iparams->bham.b, iparams->bham.c);
1184             break;
1185         case F_BONDS:
1186         case F_G96BONDS:
1187         case F_HARMONIC:
1188             pr_harm(fp, iparams, "b0", "cb");
1189             break;
1190         case F_IDIHS:
1191             pr_harm(fp, iparams, "xi", "cx");
1192             break;
1193         case F_MORSE:
1194             fprintf(fp, "b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
1195                     iparams->morse.b0A, iparams->morse.cbA, iparams->morse.betaA,
1196                     iparams->morse.b0B, iparams->morse.cbB, iparams->morse.betaB);
1197             break;
1198         case F_CUBICBONDS:
1199             fprintf(fp, "b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
1200                     iparams->cubic.b0, iparams->cubic.kb, iparams->cubic.kcub);
1201             break;
1202         case F_CONNBONDS:
1203             fprintf(fp, "\n");
1204             break;
1205         case F_FENEBONDS:
1206             fprintf(fp, "bm=%15.8e, kb=%15.8e\n", iparams->fene.bm, iparams->fene.kb);
1207             break;
1208         case F_RESTRBONDS:
1209             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",
1210                     iparams->restraint.lowA, iparams->restraint.up1A,
1211                     iparams->restraint.up2A, iparams->restraint.kA,
1212                     iparams->restraint.lowB, iparams->restraint.up1B,
1213                     iparams->restraint.up2B, iparams->restraint.kB);
1214             break;
1215         case F_TABBONDS:
1216         case F_TABBONDSNC:
1217         case F_TABANGLES:
1218         case F_TABDIHS:
1219             fprintf(fp, "tab=%d, kA=%15.8e, kB=%15.8e\n",
1220                     iparams->tab.table, iparams->tab.kA, iparams->tab.kB);
1221             break;
1222         case F_POLARIZATION:
1223             fprintf(fp, "alpha=%15.8e\n", iparams->polarize.alpha);
1224             break;
1225         case F_ANHARM_POL:
1226             fprintf(fp, "alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
1227                     iparams->anharm_polarize.alpha,
1228                     iparams->anharm_polarize.drcut,
1229                     iparams->anharm_polarize.khyp);
1230             break;
1231         case F_THOLE_POL:
1232             fprintf(fp, "a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
1233                     iparams->thole.a, iparams->thole.alpha1, iparams->thole.alpha2,
1234                     iparams->thole.rfac);
1235             break;
1236         case F_WATER_POL:
1237             fprintf(fp, "al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
1238                     iparams->wpol.al_x, iparams->wpol.al_y, iparams->wpol.al_z,
1239                     iparams->wpol.rOH, iparams->wpol.rHH, iparams->wpol.rOD);
1240             break;
1241         case F_LJ:
1242             fprintf(fp, "c6=%15.8e, c12=%15.8e\n", iparams->lj.c6, iparams->lj.c12);
1243             break;
1244         case F_LJ14:
1245             fprintf(fp, "c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
1246                     iparams->lj14.c6A, iparams->lj14.c12A,
1247                     iparams->lj14.c6B, iparams->lj14.c12B);
1248             break;
1249         case F_LJC14_Q:
1250             fprintf(fp, "fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1251                     iparams->ljc14.fqq,
1252                     iparams->ljc14.qi, iparams->ljc14.qj,
1253                     iparams->ljc14.c6, iparams->ljc14.c12);
1254             break;
1255         case F_LJC_PAIRS_NB:
1256             fprintf(fp, "qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1257                     iparams->ljcnb.qi, iparams->ljcnb.qj,
1258                     iparams->ljcnb.c6, iparams->ljcnb.c12);
1259             break;
1260         case F_PDIHS:
1261         case F_PIDIHS:
1262         case F_ANGRES:
1263         case F_ANGRESZ:
1264             fprintf(fp, "phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
1265                     iparams->pdihs.phiA, iparams->pdihs.cpA,
1266                     iparams->pdihs.phiB, iparams->pdihs.cpB,
1267                     iparams->pdihs.mult);
1268             break;
1269         case F_DISRES:
1270             fprintf(fp, "label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
1271                     iparams->disres.label, iparams->disres.type,
1272                     iparams->disres.low, iparams->disres.up1,
1273                     iparams->disres.up2, iparams->disres.kfac);
1274             break;
1275         case F_ORIRES:
1276             fprintf(fp, "ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
1277                     iparams->orires.ex, iparams->orires.label, iparams->orires.power,
1278                     iparams->orires.c, iparams->orires.obs, iparams->orires.kfac);
1279             break;
1280         case F_DIHRES:
1281             fprintf(fp, "phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
1282                     iparams->dihres.phiA, iparams->dihres.dphiA, iparams->dihres.kfacA,
1283                     iparams->dihres.phiB, iparams->dihres.dphiB, iparams->dihres.kfacB);
1284             break;
1285         case F_POSRES:
1286             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",
1287                     iparams->posres.pos0A[XX], iparams->posres.pos0A[YY],
1288                     iparams->posres.pos0A[ZZ], iparams->posres.fcA[XX],
1289                     iparams->posres.fcA[YY], iparams->posres.fcA[ZZ],
1290                     iparams->posres.pos0B[XX], iparams->posres.pos0B[YY],
1291                     iparams->posres.pos0B[ZZ], iparams->posres.fcB[XX],
1292                     iparams->posres.fcB[YY], iparams->posres.fcB[ZZ]);
1293             break;
1294         case F_FBPOSRES:
1295             fprintf(fp, "pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e\n",
1296                     iparams->fbposres.pos0[XX], iparams->fbposres.pos0[YY],
1297                     iparams->fbposres.pos0[ZZ], iparams->fbposres.geom,
1298                     iparams->fbposres.r,        iparams->fbposres.k);
1299             break;
1300         case F_RBDIHS:
1301             for (i = 0; i < NR_RBDIHS; i++)
1302             {
1303                 fprintf(fp, "%srbcA[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcA[i]);
1304             }
1305             fprintf(fp, "\n");
1306             for (i = 0; i < NR_RBDIHS; i++)
1307             {
1308                 fprintf(fp, "%srbcB[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcB[i]);
1309             }
1310             fprintf(fp, "\n");
1311             break;
1312         case F_FOURDIHS:
1313             /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
1314              * OPLS potential constants back.
1315              */
1316             rbcA = iparams->rbdihs.rbcA;
1317             rbcB = iparams->rbdihs.rbcB;
1318
1319             VA[3] = -0.25*rbcA[4];
1320             VA[2] = -0.5*rbcA[3];
1321             VA[1] = 4.0*VA[3]-rbcA[2];
1322             VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1323
1324             VB[3] = -0.25*rbcB[4];
1325             VB[2] = -0.5*rbcB[3];
1326             VB[1] = 4.0*VB[3]-rbcB[2];
1327             VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1328
1329             for (i = 0; i < NR_FOURDIHS; i++)
1330             {
1331                 fprintf(fp, "%sFourA[%d]=%15.8e", i == 0 ? "" : ", ", i, VA[i]);
1332             }
1333             fprintf(fp, "\n");
1334             for (i = 0; i < NR_FOURDIHS; i++)
1335             {
1336                 fprintf(fp, "%sFourB[%d]=%15.8e", i == 0 ? "" : ", ", i, VB[i]);
1337             }
1338             fprintf(fp, "\n");
1339             break;
1340
1341         case F_CONSTR:
1342         case F_CONSTRNC:
1343             fprintf(fp, "dA=%15.8e, dB=%15.8e\n", iparams->constr.dA, iparams->constr.dB);
1344             break;
1345         case F_SETTLE:
1346             fprintf(fp, "doh=%15.8e, dhh=%15.8e\n", iparams->settle.doh,
1347                     iparams->settle.dhh);
1348             break;
1349         case F_VSITE2:
1350             fprintf(fp, "a=%15.8e\n", iparams->vsite.a);
1351             break;
1352         case F_VSITE3:
1353         case F_VSITE3FD:
1354         case F_VSITE3FAD:
1355             fprintf(fp, "a=%15.8e, b=%15.8e\n", iparams->vsite.a, iparams->vsite.b);
1356             break;
1357         case F_VSITE3OUT:
1358         case F_VSITE4FD:
1359         case F_VSITE4FDN:
1360             fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1361                     iparams->vsite.a, iparams->vsite.b, iparams->vsite.c);
1362             break;
1363         case F_VSITEN:
1364             fprintf(fp, "n=%2d, a=%15.8e\n", iparams->vsiten.n, iparams->vsiten.a);
1365             break;
1366         case F_GB12:
1367         case F_GB13:
1368         case F_GB14:
1369             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);
1370             break;
1371         case F_CMAP:
1372             fprintf(fp, "cmapA=%1d, cmapB=%1d\n", iparams->cmap.cmapA, iparams->cmap.cmapB);
1373             break;
1374         case  F_RESTRANGLES:
1375             pr_harm(fp, iparams, "ktheta", "costheta0");
1376             break;
1377         case  F_RESTRDIHS:
1378             fprintf(fp, "phiA=%15.8e, cpA=%15.8e",
1379                     iparams->pdihs.phiA, iparams->pdihs.cpA);
1380             break;
1381         case  F_CBTDIHS:
1382             fprintf(fp, "kphi=%15.8e", iparams->cbtdihs.cbtcA[0]);
1383             for (i = 1; i < NR_CBTDIHS; i++)
1384             {
1385                 fprintf(fp, ", cbtcA[%d]=%15.8e", i-1, iparams->cbtdihs.cbtcA[i]);
1386             }
1387             fprintf(fp, "\n");
1388             break;
1389         default:
1390             gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1391                       ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1392     }
1393 }
1394
1395 void pr_ilist(FILE *fp, int indent, const char *title,
1396               t_functype *functype, t_ilist *ilist, gmx_bool bShowNumbers)
1397 {
1398     int      i, j, k, type, ftype;
1399     t_iatom *iatoms;
1400
1401     if (available(fp, ilist, indent, title) && ilist->nr > 0)
1402     {
1403         indent = pr_title(fp, indent, title);
1404         pr_indent(fp, indent);
1405         fprintf(fp, "nr: %d\n", ilist->nr);
1406         if (ilist->nr > 0)
1407         {
1408             pr_indent(fp, indent);
1409             fprintf(fp, "iatoms:\n");
1410             iatoms = ilist->iatoms;
1411             for (i = j = 0; i < ilist->nr; )
1412             {
1413 #ifndef DEBUG
1414                 pr_indent(fp, indent+INDENT);
1415                 type  = *(iatoms++);
1416                 ftype = functype[type];
1417                 fprintf(fp, "%d type=%d (%s)",
1418                         bShowNumbers ? j : -1, bShowNumbers ? type : -1,
1419                         interaction_function[ftype].name);
1420                 j++;
1421                 for (k = 0; k < interaction_function[ftype].nratoms; k++)
1422                 {
1423                     fprintf(fp, " %d", *(iatoms++));
1424                 }
1425                 fprintf(fp, "\n");
1426                 i += 1+interaction_function[ftype].nratoms;
1427 #else
1428                 fprintf(fp, "%5d%5d\n", i, iatoms[i]);
1429                 i++;
1430 #endif
1431             }
1432         }
1433     }
1434 }
1435
1436 static void pr_cmap(FILE *fp, int indent, const char *title,
1437                     gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1438 {
1439     int  i, j, nelem;
1440     real dx, idx;
1441
1442     dx    = 360.0 / cmap_grid->grid_spacing;
1443     nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1444
1445     if (available(fp, cmap_grid, indent, title))
1446     {
1447         fprintf(fp, "%s\n", title);
1448
1449         for (i = 0; i < cmap_grid->ngrid; i++)
1450         {
1451             idx = -180.0;
1452             fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
1453
1454             fprintf(fp, "grid[%3d]={\n", bShowNumbers ? i : -1);
1455
1456             for (j = 0; j < nelem; j++)
1457             {
1458                 if ( (j%cmap_grid->grid_spacing) == 0)
1459                 {
1460                     fprintf(fp, "%8.1f\n", idx);
1461                     idx += dx;
1462                 }
1463
1464                 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4]);
1465                 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+1]);
1466                 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+2]);
1467                 fprintf(fp, "%8.3f\n", cmap_grid->cmapdata[i].cmap[j*4+3]);
1468             }
1469             fprintf(fp, "\n");
1470         }
1471     }
1472
1473 }
1474
1475 void pr_ffparams(FILE *fp, int indent, const char *title,
1476                  gmx_ffparams_t *ffparams,
1477                  gmx_bool bShowNumbers)
1478 {
1479     int i;
1480
1481     indent = pr_title(fp, indent, title);
1482     pr_indent(fp, indent);
1483     fprintf(fp, "atnr=%d\n", ffparams->atnr);
1484     pr_indent(fp, indent);
1485     fprintf(fp, "ntypes=%d\n", ffparams->ntypes);
1486     for (i = 0; i < ffparams->ntypes; i++)
1487     {
1488         pr_indent(fp, indent+INDENT);
1489         fprintf(fp, "functype[%d]=%s, ",
1490                 bShowNumbers ? i : -1,
1491                 interaction_function[ffparams->functype[i]].name);
1492         pr_iparams(fp, ffparams->functype[i], &ffparams->iparams[i]);
1493     }
1494     pr_double(fp, indent, "reppow", ffparams->reppow);
1495     pr_real(fp, indent, "fudgeQQ", ffparams->fudgeQQ);
1496     pr_cmap(fp, indent, "cmap", &ffparams->cmap_grid, bShowNumbers);
1497 }
1498
1499 void pr_idef(FILE *fp, int indent, const char *title, t_idef *idef, gmx_bool bShowNumbers)
1500 {
1501     int i, j;
1502
1503     if (available(fp, idef, indent, title))
1504     {
1505         indent = pr_title(fp, indent, title);
1506         pr_indent(fp, indent);
1507         fprintf(fp, "atnr=%d\n", idef->atnr);
1508         pr_indent(fp, indent);
1509         fprintf(fp, "ntypes=%d\n", idef->ntypes);
1510         for (i = 0; i < idef->ntypes; i++)
1511         {
1512             pr_indent(fp, indent+INDENT);
1513             fprintf(fp, "functype[%d]=%s, ",
1514                     bShowNumbers ? i : -1,
1515                     interaction_function[idef->functype[i]].name);
1516             pr_iparams(fp, idef->functype[i], &idef->iparams[i]);
1517         }
1518         pr_real(fp, indent, "fudgeQQ", idef->fudgeQQ);
1519
1520         for (j = 0; (j < F_NRE); j++)
1521         {
1522             pr_ilist(fp, indent, interaction_function[j].longname,
1523                      idef->functype, &idef->il[j], bShowNumbers);
1524         }
1525     }
1526 }
1527
1528 static int pr_block_title(FILE *fp, int indent, const char *title, t_block *block)
1529 {
1530     if (available(fp, block, indent, title))
1531     {
1532         indent = pr_title(fp, indent, title);
1533         pr_indent(fp, indent);
1534         fprintf(fp, "nr=%d\n", block->nr);
1535     }
1536     return indent;
1537 }
1538
1539 static int pr_blocka_title(FILE *fp, int indent, const char *title, t_blocka *block)
1540 {
1541     if (available(fp, block, indent, title))
1542     {
1543         indent = pr_title(fp, indent, title);
1544         pr_indent(fp, indent);
1545         fprintf(fp, "nr=%d\n", block->nr);
1546         pr_indent(fp, indent);
1547         fprintf(fp, "nra=%d\n", block->nra);
1548     }
1549     return indent;
1550 }
1551
1552 static void low_pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1553 {
1554     int i;
1555
1556     if (available(fp, block, indent, title))
1557     {
1558         indent = pr_blocka_title(fp, indent, title, block);
1559         for (i = 0; i <= block->nr; i++)
1560         {
1561             pr_indent(fp, indent+INDENT);
1562             fprintf(fp, "%s->index[%d]=%d\n",
1563                     title, bShowNumbers ? i : -1, block->index[i]);
1564         }
1565         for (i = 0; i < block->nra; i++)
1566         {
1567             pr_indent(fp, indent+INDENT);
1568             fprintf(fp, "%s->a[%d]=%d\n",
1569                     title, bShowNumbers ? i : -1, block->a[i]);
1570         }
1571     }
1572 }
1573
1574 void pr_block(FILE *fp, int indent, const char *title, t_block *block, gmx_bool bShowNumbers)
1575 {
1576     int i, start;
1577
1578     if (available(fp, block, indent, title))
1579     {
1580         indent = pr_block_title(fp, indent, title, block);
1581         start  = 0;
1582         if (block->index[start] != 0)
1583         {
1584             fprintf(fp, "block->index[%d] should be 0\n", start);
1585         }
1586         else
1587         {
1588             for (i = 0; i < block->nr; i++)
1589             {
1590                 int end  = block->index[i+1];
1591                 pr_indent(fp, indent);
1592                 if (end <= start)
1593                 {
1594                     fprintf(fp, "%s[%d]={}\n", title, i);
1595                 }
1596                 else
1597                 {
1598                     fprintf(fp, "%s[%d]={%d..%d}\n",
1599                             title, bShowNumbers ? i : -1,
1600                             bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1601                 }
1602                 start = end;
1603             }
1604         }
1605     }
1606 }
1607
1608 void pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1609 {
1610     int i, j, ok, size, start, end;
1611
1612     if (available(fp, block, indent, title))
1613     {
1614         indent = pr_blocka_title(fp, indent, title, block);
1615         start  = 0;
1616         end    = start;
1617         if ((ok = (block->index[start] == 0)) == 0)
1618         {
1619             fprintf(fp, "block->index[%d] should be 0\n", start);
1620         }
1621         else
1622         {
1623             for (i = 0; i < block->nr; i++)
1624             {
1625                 end  = block->index[i+1];
1626                 size = pr_indent(fp, indent);
1627                 if (end <= start)
1628                 {
1629                     size += fprintf(fp, "%s[%d]={", title, i);
1630                 }
1631                 else
1632                 {
1633                     size += fprintf(fp, "%s[%d][%d..%d]={",
1634                                     title, bShowNumbers ? i : -1,
1635                                     bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1636                 }
1637                 for (j = start; j < end; j++)
1638                 {
1639                     if (j > start)
1640                     {
1641                         size += fprintf(fp, ", ");
1642                     }
1643                     if ((size) > (USE_WIDTH))
1644                     {
1645                         fprintf(fp, "\n");
1646                         size = pr_indent(fp, indent+INDENT);
1647                     }
1648                     size += fprintf(fp, "%d", block->a[j]);
1649                 }
1650                 fprintf(fp, "}\n");
1651                 start = end;
1652             }
1653         }
1654         if ((end != block->nra) || (!ok))
1655         {
1656             pr_indent(fp, indent);
1657             fprintf(fp, "tables inconsistent, dumping complete tables:\n");
1658             low_pr_blocka(fp, indent, title, block, bShowNumbers);
1659         }
1660     }
1661 }
1662
1663 static void pr_strings(FILE *fp, int indent, const char *title, char ***nm, int n, gmx_bool bShowNumbers)
1664 {
1665     int i;
1666
1667     if (available(fp, nm, indent, title))
1668     {
1669         indent = pr_title_n(fp, indent, title, n);
1670         for (i = 0; i < n; i++)
1671         {
1672             pr_indent(fp, indent);
1673             fprintf(fp, "%s[%d]={name=\"%s\"}\n",
1674                     title, bShowNumbers ? i : -1, *(nm[i]));
1675         }
1676     }
1677 }
1678
1679 static void pr_strings2(FILE *fp, int indent, const char *title,
1680                         char ***nm, char ***nmB, int n, gmx_bool bShowNumbers)
1681 {
1682     int i;
1683
1684     if (available(fp, nm, indent, title))
1685     {
1686         indent = pr_title_n(fp, indent, title, n);
1687         for (i = 0; i < n; i++)
1688         {
1689             pr_indent(fp, indent);
1690             fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1691                     title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
1692         }
1693     }
1694 }
1695
1696 static void pr_resinfo(FILE *fp, int indent, const char *title, t_resinfo *resinfo, int n, gmx_bool bShowNumbers)
1697 {
1698     int i;
1699
1700     if (available(fp, resinfo, indent, title))
1701     {
1702         indent = pr_title_n(fp, indent, title, n);
1703         for (i = 0; i < n; i++)
1704         {
1705             pr_indent(fp, indent);
1706             fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1707                     title, bShowNumbers ? i : -1,
1708                     *(resinfo[i].name), resinfo[i].nr,
1709                     (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1710         }
1711     }
1712 }
1713
1714 static void pr_atom(FILE *fp, int indent, const char *title, t_atom *atom, int n)
1715 {
1716     int i;
1717
1718     if (available(fp, atom, indent, title))
1719     {
1720         indent = pr_title_n(fp, indent, title, n);
1721         for (i = 0; i < n; i++)
1722         {
1723             pr_indent(fp, indent);
1724             fprintf(fp, "%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1725                     "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1726                     title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype],
1727                     atom[i].m, atom[i].q, atom[i].mB, atom[i].qB,
1728                     atom[i].resind, atom[i].atomnumber);
1729         }
1730     }
1731 }
1732
1733 static void pr_grps(FILE *fp, const char *title, t_grps grps[], char **grpname[])
1734 {
1735     int i, j;
1736
1737     for (i = 0; (i < egcNR); i++)
1738     {
1739         fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
1740         for (j = 0; (j < grps[i].nr); j++)
1741         {
1742             fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
1743         }
1744         fprintf(fp, "]\n");
1745     }
1746 }
1747
1748 static void pr_groups(FILE *fp, int indent,
1749                       gmx_groups_t *groups,
1750                       gmx_bool bShowNumbers)
1751 {
1752     int nat_max, i, g;
1753
1754     pr_grps(fp, "grp", groups->grps, groups->grpname);
1755     pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
1756
1757     pr_indent(fp, indent);
1758     fprintf(fp, "groups          ");
1759     for (g = 0; g < egcNR; g++)
1760     {
1761         printf(" %5.5s", gtypes[g]);
1762     }
1763     printf("\n");
1764
1765     pr_indent(fp, indent);
1766     fprintf(fp, "allocated       ");
1767     nat_max = 0;
1768     for (g = 0; g < egcNR; g++)
1769     {
1770         printf(" %5d", groups->ngrpnr[g]);
1771         nat_max = std::max(nat_max, groups->ngrpnr[g]);
1772     }
1773     printf("\n");
1774
1775     if (nat_max == 0)
1776     {
1777         pr_indent(fp, indent);
1778         fprintf(fp, "groupnr[%5s] =", "*");
1779         for (g = 0; g < egcNR; g++)
1780         {
1781             fprintf(fp, "  %3d ", 0);
1782         }
1783         fprintf(fp, "\n");
1784     }
1785     else
1786     {
1787         for (i = 0; i < nat_max; i++)
1788         {
1789             pr_indent(fp, indent);
1790             fprintf(fp, "groupnr[%5d] =", i);
1791             for (g = 0; g < egcNR; g++)
1792             {
1793                 fprintf(fp, "  %3d ",
1794                         groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1795             }
1796             fprintf(fp, "\n");
1797         }
1798     }
1799 }
1800
1801 void pr_atoms(FILE *fp, int indent, const char *title, t_atoms *atoms,
1802               gmx_bool bShownumbers)
1803 {
1804     if (available(fp, atoms, indent, title))
1805     {
1806         indent = pr_title(fp, indent, title);
1807         pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
1808         pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
1809         pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
1810         pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
1811     }
1812 }
1813
1814
1815 void pr_atomtypes(FILE *fp, int indent, const char *title, t_atomtypes *atomtypes,
1816                   gmx_bool bShowNumbers)
1817 {
1818     int i;
1819     if (available(fp, atomtypes, indent, title))
1820     {
1821         indent = pr_title(fp, indent, title);
1822         for (i = 0; i < atomtypes->nr; i++)
1823         {
1824             pr_indent(fp, indent);
1825             fprintf(fp,
1826                     "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1827                     bShowNumbers ? i : -1, atomtypes->radius[i], atomtypes->vol[i],
1828                     atomtypes->gb_radius[i],
1829                     atomtypes->surftens[i], atomtypes->atomnumber[i], atomtypes->S_hct[i]);
1830         }
1831     }
1832 }
1833
1834 static void pr_moltype(FILE *fp, int indent, const char *title,
1835                        gmx_moltype_t *molt, int n,
1836                        gmx_ffparams_t *ffparams,
1837                        gmx_bool bShowNumbers)
1838 {
1839     int j;
1840
1841     indent = pr_title_n(fp, indent, title, n);
1842     pr_indent(fp, indent);
1843     fprintf(fp, "name=\"%s\"\n", *(molt->name));
1844     pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
1845     pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
1846     pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
1847     for (j = 0; (j < F_NRE); j++)
1848     {
1849         pr_ilist(fp, indent, interaction_function[j].longname,
1850                  ffparams->functype, &molt->ilist[j], bShowNumbers);
1851     }
1852 }
1853
1854 static void pr_molblock(FILE *fp, int indent, const char *title,
1855                         gmx_molblock_t *molb, int n,
1856                         gmx_moltype_t *molt)
1857 {
1858     indent = pr_title_n(fp, indent, title, n);
1859     pr_indent(fp, indent);
1860     fprintf(fp, "%-20s = %d \"%s\"\n",
1861             "moltype", molb->type, *(molt[molb->type].name));
1862     pr_int(fp, indent, "#molecules", molb->nmol);
1863     pr_int(fp, indent, "#atoms_mol", molb->natoms_mol);
1864     pr_int(fp, indent, "#posres_xA", molb->nposres_xA);
1865     if (molb->nposres_xA > 0)
1866     {
1867         pr_rvecs(fp, indent, "posres_xA", molb->posres_xA, molb->nposres_xA);
1868     }
1869     pr_int(fp, indent, "#posres_xB", molb->nposres_xB);
1870     if (molb->nposres_xB > 0)
1871     {
1872         pr_rvecs(fp, indent, "posres_xB", molb->posres_xB, molb->nposres_xB);
1873     }
1874 }
1875
1876 void pr_mtop(FILE *fp, int indent, const char *title, gmx_mtop_t *mtop,
1877              gmx_bool bShowNumbers)
1878 {
1879     int mt, mb, j;
1880
1881     if (available(fp, mtop, indent, title))
1882     {
1883         indent = pr_title(fp, indent, title);
1884         pr_indent(fp, indent);
1885         fprintf(fp, "name=\"%s\"\n", *(mtop->name));
1886         pr_int(fp, indent, "#atoms", mtop->natoms);
1887         pr_int(fp, indent, "#molblock", mtop->nmolblock);
1888         for (mb = 0; mb < mtop->nmolblock; mb++)
1889         {
1890             pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
1891         }
1892         pr_str(fp, indent, "bIntermolecularInteractions", EBOOL(mtop->bIntermolecularInteractions));
1893         if (mtop->bIntermolecularInteractions)
1894         {
1895             for (j = 0; (j < F_NRE); j++)
1896             {
1897                 pr_ilist(fp, indent, interaction_function[j].longname,
1898                          mtop->ffparams.functype,
1899                          &mtop->intermolecular_ilist[j], bShowNumbers);
1900             }
1901         }
1902         pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
1903         pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
1904         for (mt = 0; mt < mtop->nmoltype; mt++)
1905         {
1906             pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
1907                        &mtop->ffparams, bShowNumbers);
1908         }
1909         pr_groups(fp, indent, &mtop->groups, bShowNumbers);
1910     }
1911 }
1912
1913 void pr_top(FILE *fp, int indent, const char *title, t_topology *top, gmx_bool bShowNumbers)
1914 {
1915     if (available(fp, top, indent, title))
1916     {
1917         indent = pr_title(fp, indent, title);
1918         pr_indent(fp, indent);
1919         fprintf(fp, "name=\"%s\"\n", *(top->name));
1920         pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
1921         pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
1922         pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
1923         pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
1924         pr_str(fp, indent, "bIntermolecularInteractions", EBOOL(top->bIntermolecularInteractions));
1925         pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
1926         pr_idef(fp, indent, "idef", &top->idef, bShowNumbers);
1927     }
1928 }
1929
1930 void pr_header(FILE *fp, int indent, const char *title, t_tpxheader *sh)
1931 {
1932     if (available(fp, sh, indent, title))
1933     {
1934         indent = pr_title(fp, indent, title);
1935         pr_indent(fp, indent);
1936         fprintf(fp, "bIr    = %spresent\n", sh->bIr ? "" : "not ");
1937         pr_indent(fp, indent);
1938         fprintf(fp, "bBox   = %spresent\n", sh->bBox ? "" : "not ");
1939         pr_indent(fp, indent);
1940         fprintf(fp, "bTop   = %spresent\n", sh->bTop ? "" : "not ");
1941         pr_indent(fp, indent);
1942         fprintf(fp, "bX     = %spresent\n", sh->bX ? "" : "not ");
1943         pr_indent(fp, indent);
1944         fprintf(fp, "bV     = %spresent\n", sh->bV ? "" : "not ");
1945         pr_indent(fp, indent);
1946         fprintf(fp, "bF     = %spresent\n", sh->bF ? "" : "not ");
1947
1948         pr_indent(fp, indent);
1949         fprintf(fp, "natoms = %d\n", sh->natoms);
1950         pr_indent(fp, indent);
1951         fprintf(fp, "lambda = %e\n", sh->lambda);
1952     }
1953 }
1954
1955 void pr_commrec(FILE *fp, int indent, t_commrec *cr)
1956 {
1957     pr_indent(fp, indent);
1958     fprintf(fp, "commrec:\n");
1959     indent += 2;
1960     pr_indent(fp, indent);
1961     fprintf(fp, "rank      = %d\n", cr->nodeid);
1962     pr_indent(fp, indent);
1963     fprintf(fp, "number of ranks = %d\n", cr->nnodes);
1964     pr_indent(fp, indent);
1965     fprintf(fp, "PME-only ranks = %d\n", cr->npmenodes);
1966     /*
1967        pr_indent(fp,indent);
1968        fprintf(fp,"threadid  = %d\n",cr->threadid);
1969        pr_indent(fp,indent);
1970        fprintf(fp,"nthreads  = %d\n",cr->nthreads);
1971      */
1972 }