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