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