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