945a31a4d86c734251bb4f8669ae62438371c767
[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_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
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                          EBOOL(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_simtempvals(FILE *fp,int indent,t_simtemp *simtemp, int n_lambda, gmx_bool bMDPformat)
496 {
497     PR("simtemp_low",simtemp->simtemp_low);
498     PR("simtemp_high",simtemp->simtemp_high);
499     PS("simulated-tempering-scaling",ESIMTEMP(simtemp->eSimTempScale));
500     pr_rvec(fp,indent,"simulated tempering temperatures",simtemp->temperatures,n_lambda,TRUE);
501 }
502
503 static void pr_expandedvals(FILE *fp,int indent,t_expanded *expand, int n_lambda, gmx_bool bMDPformat)
504 {
505
506     PI("nstexpanded", expand->nstexpanded);
507     PS("lambda-stats", elamstats_names[expand->elamstats]);
508     PS("lambda-mc-move", elmcmove_names[expand->elmcmove]);
509     PI("lmc-repeats",expand->lmc_repeats);
510     PI("lmc-gibbsdelta",expand->gibbsdeltalam);
511     PI("lmc-nstart",expand->lmc_forced_nstart);
512     PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
513     PI("nst-transition-matrix",expand->nstTij);
514     PI("mininum-var-min",expand->minvarmin); /*default is reasonable */
515     PI("weight-c-range",expand->c_range); /* default is just C=0 */
516     PR("wl-scale",expand->wl_scale);
517     PR("init-wl-delta",expand->init_wl_delta);
518     PR("wl-ratio",expand->wl_ratio);
519     PS("bWLoneovert",EBOOL(expand->bWLoneovert));
520     PI("lmc-seed",expand->lmc_seed);
521     PR("mc-temperature",expand->mc_temp);
522     PS("lmc-weights-equil",elmceq_names[expand->elmceq]);
523     if (expand->elmceq == elmceqNUMATLAM)
524     {
525         PI("weight-equil-number-all-lambda",expand->equil_n_at_lam);
526     }
527     if (expand->elmceq == elmceqSAMPLES)
528     {
529         PI("weight-equil-number-samples",expand->equil_samples);
530     }
531     if (expand->elmceq == elmceqSTEPS)
532     {
533         PI("weight-equil-number-steps",expand->equil_steps);
534     }
535     if (expand->elmceq == elmceqWLDELTA)
536     {
537         PR("weight-equil-wl-delta",expand->equil_wl_delta);
538     }
539     if (expand->elmceq == elmceqRATIO)
540     {
541         PR("weight-equil-count-ratio",expand->equil_ratio);
542     }
543
544     pr_indent(fp,indent);
545     pr_rvec(fp,indent,"init-lambda-weights",expand->init_lambda_weights,n_lambda,TRUE);
546     PS("init-weights",EBOOL(expand->bInit_weights));
547 }
548
549 static void pr_fepvals(FILE *fp,int indent,t_lambda *fep, gmx_bool bMDPformat)
550 {
551     int i,j;
552
553     PI("nstdhdl",fep->nstdhdl);
554     PI("init-lambda-state",fep->init_fep_state);
555     PR("init-lambda",fep->init_lambda);
556     PR("delta-lambda",fep->delta_lambda);
557     if (!bMDPformat)
558     {
559         PI("n-lambdas",fep->n_lambda);
560     }
561     if (fep->n_lambda > 0)
562     {
563         pr_indent(fp,indent);
564         fprintf(fp,"all-lambdas%s\n",bMDPformat ? " = " : ":");
565         for(i=0; i<efptNR; i++) {
566             fprintf(fp,"%18s = ",efpt_names[i]);
567             for(j=0; j<fep->n_lambda; j++)
568             {
569                 fprintf(fp,"  %10g",fep->all_lambda[i][j]);
570             }
571             fprintf(fp,"\n");
572         }
573     }
574
575     PR("sc-alpha",fep->sc_alpha);
576     PS("bScCoul",EBOOL(fep->bScCoul));
577     PS("bScPrintEnergy",EBOOL(fep->bPrintEnergy));
578     PI("sc-power",fep->sc_power);
579     PR("sc-r-power",fep->sc_r_power);
580     PR("sc-sigma",fep->sc_sigma);
581     PR("sc-sigma-min",fep->sc_sigma_min);
582     PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
583     PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
584     PI("dh-hist-size", fep->dh_hist_size);
585     PD("dh-hist-spacing", fep->dh_hist_spacing);
586 };
587
588 static void pr_pull(FILE *fp,int indent,t_pull *pull)
589 {
590   int g;
591
592   PS("pull-geometry",EPULLGEOM(pull->eGeom));
593   pr_ivec(fp,indent,"pull-dim",pull->dim,DIM,TRUE);
594   PR("pull-r1",pull->cyl_r1);
595   PR("pull-r0",pull->cyl_r0);
596   PR("pull-constr-tol",pull->constr_tol);
597   PI("pull-nstxout",pull->nstxout);
598   PI("pull-nstfout",pull->nstfout);
599   PI("pull-ngrp",pull->ngrp);
600   for(g=0; g<pull->ngrp+1; g++)
601     pr_pullgrp(fp,indent,g,&pull->grp[g]);
602 }
603
604 static void pr_rotgrp(FILE *fp,int indent,int g,t_rotgrp *rotg)
605 {
606   pr_indent(fp,indent);
607   fprintf(fp,"rotation_group %d:\n",g);
608   indent += 2;
609   PS("type",EROTGEOM(rotg->eType));
610   PS("massw",EBOOL(rotg->bMassW));
611   pr_ivec_block(fp,indent,"atom",rotg->ind,rotg->nat,TRUE);
612   pr_rvecs(fp,indent,"x_ref",rotg->x_ref,rotg->nat);
613   pr_rvec(fp,indent,"vec",rotg->vec,DIM,TRUE);
614   pr_rvec(fp,indent,"pivot",rotg->pivot,DIM,TRUE);
615   PR("rate",rotg->rate);
616   PR("k",rotg->k);
617   PR("slab_dist",rotg->slab_dist);
618   PR("min_gaussian",rotg->min_gaussian);
619   PR("epsilon",rotg->eps);
620   PS("fit_method",EROTFIT(rotg->eFittype));
621   PI("potfitangle_nstep",rotg->PotAngle_nstep);
622   PR("potfitangle_step",rotg->PotAngle_step);
623 }
624
625 static void pr_rot(FILE *fp,int indent,t_rot *rot)
626 {
627   int g;
628
629   PI("rot_nstrout",rot->nstrout);
630   PI("rot_nstsout",rot->nstsout);
631   PI("rot_ngrp",rot->ngrp);
632   for(g=0; g<rot->ngrp; g++)
633     pr_rotgrp(fp,indent,g,&rot->grp[g]);
634 }
635
636 void pr_inputrec(FILE *fp,int indent,const char *title,t_inputrec *ir,
637                  gmx_bool bMDPformat)
638 {
639   const char *infbuf="inf";
640   int  i;
641   
642   if (available(fp,ir,indent,title)) {
643     if (!bMDPformat)
644       indent=pr_title(fp,indent,title);
645     PS("integrator",EI(ir->eI));
646     PSTEP("nsteps",ir->nsteps);
647     PSTEP("init-step",ir->init_step);
648     PS("cutoff-scheme",ECUTSCHEME(ir->cutoff_scheme));
649     PS("ns_type",ENS(ir->ns_type));
650     PI("nstlist",ir->nstlist);
651     PI("ndelta",ir->ndelta);
652     PI("nstcomm",ir->nstcomm);
653     PS("comm-mode",ECOM(ir->comm_mode));
654     PI("nstlog",ir->nstlog);
655     PI("nstxout",ir->nstxout);
656     PI("nstvout",ir->nstvout);
657     PI("nstfout",ir->nstfout);
658     PI("nstcalcenergy",ir->nstcalcenergy);
659     PI("nstenergy",ir->nstenergy);
660     PI("nstxtcout",ir->nstxtcout);
661     PR("init-t",ir->init_t);
662     PR("delta-t",ir->delta_t);
663     
664     PR("xtcprec",ir->xtcprec);
665     PR("fourierspacing",ir->fourier_spacing);
666     PI("nkx",ir->nkx);
667     PI("nky",ir->nky);
668     PI("nkz",ir->nkz);
669     PI("pme-order",ir->pme_order);
670     PR("ewald-rtol",ir->ewald_rtol);
671     PR("ewald-geometry",ir->ewald_geometry);
672     PR("epsilon-surface",ir->epsilon_surface);
673     PS("optimize-fft",EBOOL(ir->bOptFFT));
674     PS("ePBC",EPBC(ir->ePBC));
675     PS("bPeriodicMols",EBOOL(ir->bPeriodicMols));
676     PS("bContinuation",EBOOL(ir->bContinuation));
677     PS("bShakeSOR",EBOOL(ir->bShakeSOR));
678     PS("etc",ETCOUPLTYPE(ir->etc));
679     PS("bPrintNHChains",EBOOL(ir->bPrintNHChains));
680     PI("nsttcouple",ir->nsttcouple);
681     PS("epc",EPCOUPLTYPE(ir->epc));
682     PS("epctype",EPCOUPLTYPETYPE(ir->epct));
683     PI("nstpcouple",ir->nstpcouple);
684     PR("tau-p",ir->tau_p);
685     pr_matrix(fp,indent,"ref-p",ir->ref_p,bMDPformat);
686     pr_matrix(fp,indent,"compress",ir->compress,bMDPformat);
687     PS("refcoord-scaling",EREFSCALINGTYPE(ir->refcoord_scaling));
688     if (bMDPformat)
689       fprintf(fp,"posres-com  = %g %g %g\n",ir->posres_com[XX],
690               ir->posres_com[YY],ir->posres_com[ZZ]);
691     else
692       pr_rvec(fp,indent,"posres-com",ir->posres_com,DIM,TRUE);
693     if (bMDPformat)
694       fprintf(fp,"posres-comB = %g %g %g\n",ir->posres_comB[XX],
695               ir->posres_comB[YY],ir->posres_comB[ZZ]);
696     else
697       pr_rvec(fp,indent,"posres-comB",ir->posres_comB,DIM,TRUE);
698     PR("verlet-buffer-drift",ir->verletbuf_drift);
699     PR("rlist",ir->rlist);
700     PR("rlistlong",ir->rlistlong);
701     PR("nstcalclr",ir->nstcalclr);
702     PR("rtpi",ir->rtpi);
703     PS("coulombtype",EELTYPE(ir->coulombtype));
704     PS("coulomb-modifier",INTMODIFIER(ir->coulomb_modifier));
705     PR("rcoulomb-switch",ir->rcoulomb_switch);
706     PR("rcoulomb",ir->rcoulomb);
707     PS("vdwtype",EVDWTYPE(ir->vdwtype));
708     PS("vdw-modifier",INTMODIFIER(ir->vdw_modifier));
709     PR("rvdw-switch",ir->rvdw_switch);
710     PR("rvdw",ir->rvdw);
711     if (ir->epsilon_r != 0)
712       PR("epsilon-r",ir->epsilon_r);
713     else
714       PS("epsilon-r",infbuf);
715     if (ir->epsilon_rf != 0)
716       PR("epsilon-rf",ir->epsilon_rf);
717     else
718       PS("epsilon-rf",infbuf);
719     PR("tabext",ir->tabext);
720     PS("implicit-solvent",EIMPLICITSOL(ir->implicit_solvent));
721     PS("gb-algorithm",EGBALGORITHM(ir->gb_algorithm));
722     PR("gb-epsilon-solvent",ir->gb_epsilon_solvent);
723     PI("nstgbradii",ir->nstgbradii);
724     PR("rgbradii",ir->rgbradii);
725     PR("gb-saltconc",ir->gb_saltconc);
726     PR("gb-obc-alpha",ir->gb_obc_alpha);
727     PR("gb-obc-beta",ir->gb_obc_beta);
728     PR("gb-obc-gamma",ir->gb_obc_gamma);
729     PR("gb-dielectric-offset",ir->gb_dielectric_offset);
730     PS("sa-algorithm",ESAALGORITHM(ir->gb_algorithm));
731     PR("sa-surface-tension",ir->sa_surface_tension);
732     PS("DispCorr",EDISPCORR(ir->eDispCorr));
733     PS("bSimTemp",EBOOL(ir->bSimTemp));
734     if (ir->bSimTemp) {
735         pr_simtempvals(fp,indent,ir->simtempvals,ir->fepvals->n_lambda,bMDPformat);
736     }
737     PS("free-energy",EFEPTYPE(ir->efep));
738     if (ir->efep != efepNO || ir->bSimTemp) {
739         pr_fepvals(fp,indent,ir->fepvals,bMDPformat);
740     }
741     if (ir->bExpanded) {
742         pr_expandedvals(fp,indent,ir->expandedvals,ir->fepvals->n_lambda,bMDPformat);
743     }
744
745     PI("nwall",ir->nwall);
746     PS("wall-type",EWALLTYPE(ir->wall_type));
747     PI("wall-atomtype[0]",ir->wall_atomtype[0]);
748     PI("wall-atomtype[1]",ir->wall_atomtype[1]);
749     PR("wall-density[0]",ir->wall_density[0]);
750     PR("wall-density[1]",ir->wall_density[1]);
751     PR("wall-ewald-zfac",ir->wall_ewald_zfac);
752
753     PS("pull",EPULLTYPE(ir->ePull));
754     if (ir->ePull != epullNO)
755       pr_pull(fp,indent,ir->pull);
756     
757     PS("rotation",EBOOL(ir->bRot));
758     if (ir->bRot)
759       pr_rot(fp,indent,ir->rot);
760
761     PS("disre",EDISRETYPE(ir->eDisre));
762     PS("disre-weighting",EDISREWEIGHTING(ir->eDisreWeighting));
763     PS("disre-mixed",EBOOL(ir->bDisreMixed));
764     PR("dr-fc",ir->dr_fc);
765     PR("dr-tau",ir->dr_tau);
766     PR("nstdisreout",ir->nstdisreout);
767     PR("orires-fc",ir->orires_fc);
768     PR("orires-tau",ir->orires_tau);
769     PR("nstorireout",ir->nstorireout);
770
771     PR("dihre-fc",ir->dihre_fc);
772     
773     PR("em-stepsize",ir->em_stepsize);
774     PR("em-tol",ir->em_tol);
775     PI("niter",ir->niter);
776     PR("fc-stepsize",ir->fc_stepsize);
777     PI("nstcgsteep",ir->nstcgsteep);
778     PI("nbfgscorr",ir->nbfgscorr);
779
780     PS("ConstAlg",ECONSTRTYPE(ir->eConstrAlg));
781     PR("shake-tol",ir->shake_tol);
782     PI("lincs-order",ir->nProjOrder);
783     PR("lincs-warnangle",ir->LincsWarnAngle);
784     PI("lincs-iter",ir->nLincsIter);
785     PR("bd-fric",ir->bd_fric);
786     PI("ld-seed",ir->ld_seed);
787     PR("cos-accel",ir->cos_accel);
788     pr_matrix(fp,indent,"deform",ir->deform,bMDPformat);
789
790     PS("adress",EBOOL(ir->bAdress));
791     if (ir->bAdress){
792         PS("adress_type",EADRESSTYPE(ir->adress->type));
793         PR("adress_const_wf",ir->adress->const_wf);
794         PR("adress_ex_width",ir->adress->ex_width);
795         PR("adress_hy_width",ir->adress->hy_width);
796         PS("adress_interface_correction",EADRESSICTYPE(ir->adress->icor));
797         PS("adress_site",EADRESSSITETYPE(ir->adress->site));
798         PR("adress_ex_force_cap",ir->adress->ex_forcecap);
799         PS("adress_do_hybridpairs", EBOOL(ir->adress->do_hybridpairs));
800
801         pr_rvec(fp,indent,"adress_reference_coords",ir->adress->refs,DIM,TRUE);
802     }
803     PI("userint1",ir->userint1);
804     PI("userint2",ir->userint2);
805     PI("userint3",ir->userint3);
806     PI("userint4",ir->userint4);
807     PR("userreal1",ir->userreal1);
808     PR("userreal2",ir->userreal2);
809     PR("userreal3",ir->userreal3);
810     PR("userreal4",ir->userreal4);
811     pr_grp_opts(fp,indent,"grpopts",&(ir->opts),bMDPformat);
812     pr_cosine(fp,indent,"efield-x",&(ir->ex[XX]),bMDPformat);
813     pr_cosine(fp,indent,"efield-xt",&(ir->et[XX]),bMDPformat);
814     pr_cosine(fp,indent,"efield-y",&(ir->ex[YY]),bMDPformat);
815     pr_cosine(fp,indent,"efield-yt",&(ir->et[YY]),bMDPformat);
816     pr_cosine(fp,indent,"efield-z",&(ir->ex[ZZ]),bMDPformat);
817     pr_cosine(fp,indent,"efield-zt",&(ir->et[ZZ]),bMDPformat);
818     PS("bQMMM",EBOOL(ir->bQMMM));
819     PI("QMconstraints",ir->QMconstraints);
820     PI("QMMMscheme",ir->QMMMscheme);
821     PR("scalefactor",ir->scalefactor);
822     pr_qm_opts(fp,indent,"qm-opts",&(ir->opts));
823   }
824 }
825 #undef PS
826 #undef PR
827 #undef PI
828
829 static void pr_harm(FILE *fp,t_iparams *iparams,const char *r,const char *kr)
830 {
831   fprintf(fp,"%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
832           r,iparams->harmonic.rA,kr,iparams->harmonic.krA,
833           r,iparams->harmonic.rB,kr,iparams->harmonic.krB);
834 }
835
836 void pr_iparams(FILE *fp,t_functype ftype,t_iparams *iparams)
837 {
838   int i;
839   real VA[4],VB[4],*rbcA,*rbcB;
840
841   switch (ftype) {
842   case F_ANGLES:
843   case F_G96ANGLES:
844     pr_harm(fp,iparams,"th","ct");
845     break;
846   case F_CROSS_BOND_BONDS:
847     fprintf(fp,"r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
848             iparams->cross_bb.r1e,iparams->cross_bb.r2e,
849             iparams->cross_bb.krr);
850     break;
851   case F_CROSS_BOND_ANGLES:
852     fprintf(fp,"r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
853             iparams->cross_ba.r1e,iparams->cross_ba.r2e,
854             iparams->cross_ba.r3e,iparams->cross_ba.krt);
855     break;
856   case F_LINEAR_ANGLES:
857     fprintf(fp,"klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
858             iparams->linangle.klinA,iparams->linangle.aA,
859             iparams->linangle.klinB,iparams->linangle.aB);
860     break;
861   case F_UREY_BRADLEY:
862       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);
863     break;
864   case F_QUARTIC_ANGLES:
865     fprintf(fp,"theta=%15.8e",iparams->qangle.theta);
866     for(i=0; i<5; i++)
867       fprintf(fp,", c%c=%15.8e",'0'+i,iparams->qangle.c[i]);
868     fprintf(fp,"\n");
869     break;
870   case F_BHAM:
871     fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
872             iparams->bham.a,iparams->bham.b,iparams->bham.c);
873     break;
874   case F_BONDS:
875   case F_G96BONDS:
876   case F_HARMONIC:
877     pr_harm(fp,iparams,"b0","cb");
878     break;
879   case F_IDIHS:
880     pr_harm(fp,iparams,"xi","cx");
881     break;
882   case F_MORSE:
883     fprintf(fp,"b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
884             iparams->morse.b0A,iparams->morse.cbA,iparams->morse.betaA,
885             iparams->morse.b0B,iparams->morse.cbB,iparams->morse.betaB);
886     break;
887   case F_CUBICBONDS:
888     fprintf(fp,"b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
889             iparams->cubic.b0,iparams->cubic.kb,iparams->cubic.kcub);
890     break;
891   case F_CONNBONDS:
892     fprintf(fp,"\n");
893     break;
894   case F_FENEBONDS:
895     fprintf(fp,"bm=%15.8e, kb=%15.8e\n",iparams->fene.bm,iparams->fene.kb);
896     break;
897   case F_RESTRBONDS:
898       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",
899               iparams->restraint.lowA,iparams->restraint.up1A,
900               iparams->restraint.up2A,iparams->restraint.kA,
901               iparams->restraint.lowB,iparams->restraint.up1B,
902               iparams->restraint.up2B,iparams->restraint.kB);
903       break;
904   case F_TABBONDS:
905   case F_TABBONDSNC:
906   case F_TABANGLES:
907   case F_TABDIHS:
908     fprintf(fp,"tab=%d, kA=%15.8e, kB=%15.8e\n",
909             iparams->tab.table,iparams->tab.kA,iparams->tab.kB);
910     break;
911   case F_POLARIZATION:
912     fprintf(fp,"alpha=%15.8e\n",iparams->polarize.alpha);
913     break;
914   case F_ANHARM_POL:
915     fprintf(fp,"alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
916             iparams->anharm_polarize.alpha,
917             iparams->anharm_polarize.drcut,
918             iparams->anharm_polarize.khyp);
919     break;
920   case F_THOLE_POL:
921     fprintf(fp,"a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
922             iparams->thole.a,iparams->thole.alpha1,iparams->thole.alpha2,
923             iparams->thole.rfac);
924     break;
925   case F_WATER_POL:
926     fprintf(fp,"al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
927             iparams->wpol.al_x,iparams->wpol.al_y,iparams->wpol.al_z,
928             iparams->wpol.rOH,iparams->wpol.rHH,iparams->wpol.rOD);
929     break;
930   case F_LJ:
931     fprintf(fp,"c6=%15.8e, c12=%15.8e\n",iparams->lj.c6,iparams->lj.c12);
932     break;
933   case F_LJ14:
934     fprintf(fp,"c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
935             iparams->lj14.c6A,iparams->lj14.c12A,
936             iparams->lj14.c6B,iparams->lj14.c12B);
937     break;
938   case F_LJC14_Q:
939     fprintf(fp,"fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
940             iparams->ljc14.fqq,
941             iparams->ljc14.qi,iparams->ljc14.qj,
942             iparams->ljc14.c6,iparams->ljc14.c12);
943     break;
944   case F_LJC_PAIRS_NB:
945     fprintf(fp,"qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
946             iparams->ljcnb.qi,iparams->ljcnb.qj,
947             iparams->ljcnb.c6,iparams->ljcnb.c12);
948     break;
949   case F_PDIHS:
950   case F_PIDIHS:
951   case F_ANGRES:
952   case F_ANGRESZ:
953     fprintf(fp,"phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
954             iparams->pdihs.phiA,iparams->pdihs.cpA,
955             iparams->pdihs.phiB,iparams->pdihs.cpB,
956             iparams->pdihs.mult);
957     break;
958   case F_DISRES:
959     fprintf(fp,"label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
960             iparams->disres.label,iparams->disres.type,
961             iparams->disres.low,iparams->disres.up1,
962             iparams->disres.up2,iparams->disres.kfac);
963     break;
964   case F_ORIRES:
965     fprintf(fp,"ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
966             iparams->orires.ex,iparams->orires.label,iparams->orires.power,
967             iparams->orires.c,iparams->orires.obs,iparams->orires.kfac);
968     break;
969   case F_DIHRES:
970       fprintf(fp,"phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
971               iparams->dihres.phiA,iparams->dihres.dphiA,iparams->dihres.kfacA,
972               iparams->dihres.phiB,iparams->dihres.dphiB,iparams->dihres.kfacB);
973     break;
974   case F_POSRES:
975     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",
976             iparams->posres.pos0A[XX],iparams->posres.pos0A[YY],
977             iparams->posres.pos0A[ZZ],iparams->posres.fcA[XX],
978             iparams->posres.fcA[YY],iparams->posres.fcA[ZZ],
979             iparams->posres.pos0B[XX],iparams->posres.pos0B[YY],
980             iparams->posres.pos0B[ZZ],iparams->posres.fcB[XX],
981             iparams->posres.fcB[YY],iparams->posres.fcB[ZZ]);
982     break;
983   case F_RBDIHS:
984     for (i=0; i<NR_RBDIHS; i++) 
985       fprintf(fp,"%srbcA[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcA[i]);
986     fprintf(fp,"\n");
987     for (i=0; i<NR_RBDIHS; i++) 
988       fprintf(fp,"%srbcB[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcB[i]);
989     fprintf(fp,"\n");
990     break;
991   case F_FOURDIHS:
992     /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
993      * OPLS potential constants back.
994      */
995     rbcA = iparams->rbdihs.rbcA;
996     rbcB = iparams->rbdihs.rbcB;
997
998     VA[3] = -0.25*rbcA[4];
999     VA[2] = -0.5*rbcA[3];
1000     VA[1] = 4.0*VA[3]-rbcA[2];
1001     VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1002
1003     VB[3] = -0.25*rbcB[4];
1004     VB[2] = -0.5*rbcB[3];
1005     VB[1] = 4.0*VB[3]-rbcB[2];
1006     VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1007
1008     for (i=0; i<NR_FOURDIHS; i++) 
1009       fprintf(fp,"%sFourA[%d]=%15.8e",i==0?"":", ",i,VA[i]);
1010     fprintf(fp,"\n");
1011     for (i=0; i<NR_FOURDIHS; i++) 
1012       fprintf(fp,"%sFourB[%d]=%15.8e",i==0?"":", ",i,VB[i]);
1013     fprintf(fp,"\n");
1014     break;
1015    
1016   case F_CONSTR:
1017   case F_CONSTRNC:
1018     fprintf(fp,"dA=%15.8e, dB=%15.8e\n",iparams->constr.dA,iparams->constr.dB);
1019     break;
1020   case F_SETTLE:
1021     fprintf(fp,"doh=%15.8e, dhh=%15.8e\n",iparams->settle.doh,
1022             iparams->settle.dhh);
1023     break;
1024   case F_VSITE2:
1025     fprintf(fp,"a=%15.8e\n",iparams->vsite.a);
1026     break;
1027   case F_VSITE3:
1028   case F_VSITE3FD:
1029   case F_VSITE3FAD:
1030     fprintf(fp,"a=%15.8e, b=%15.8e\n",iparams->vsite.a,iparams->vsite.b);
1031     break;
1032   case F_VSITE3OUT:
1033   case F_VSITE4FD:
1034   case F_VSITE4FDN:
1035     fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
1036             iparams->vsite.a,iparams->vsite.b,iparams->vsite.c);
1037     break;
1038   case F_VSITEN:
1039     fprintf(fp,"n=%2d, a=%15.8e\n",iparams->vsiten.n,iparams->vsiten.a);
1040     break;
1041   case F_GB12:
1042   case F_GB13:
1043   case F_GB14:
1044     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);
1045     break;                
1046   case F_CMAP:
1047     fprintf(fp, "cmapA=%1d, cmapB=%1d\n",iparams->cmap.cmapA, iparams->cmap.cmapB);
1048     break;                
1049   default:
1050     gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1051               ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1052   }
1053 }
1054
1055 void pr_ilist(FILE *fp,int indent,const char *title,
1056               t_functype *functype,t_ilist *ilist, gmx_bool bShowNumbers)
1057 {
1058     int i,j,k,type,ftype;
1059     t_iatom *iatoms;
1060     
1061     if (available(fp,ilist,indent,title) && ilist->nr > 0)
1062     {  
1063         indent=pr_title(fp,indent,title);
1064         (void) pr_indent(fp,indent);
1065         fprintf(fp,"nr: %d\n",ilist->nr);
1066         if (ilist->nr > 0) {
1067             (void) pr_indent(fp,indent);
1068             fprintf(fp,"iatoms:\n");
1069             iatoms=ilist->iatoms;
1070             for (i=j=0; i<ilist->nr;) {
1071 #ifndef DEBUG
1072                 (void) pr_indent(fp,indent+INDENT);
1073                 type=*(iatoms++);
1074                 ftype=functype[type];
1075                 (void) fprintf(fp,"%d type=%d (%s)",
1076                                bShowNumbers?j:-1,bShowNumbers?type:-1,
1077                                interaction_function[ftype].name);
1078                 j++;
1079                 for (k=0; k<interaction_function[ftype].nratoms; k++)
1080                     (void) fprintf(fp," %u",*(iatoms++));
1081                 (void) fprintf(fp,"\n");
1082                 i+=1+interaction_function[ftype].nratoms;
1083 #else
1084                 fprintf(fp,"%5d%5d\n",i,iatoms[i]);
1085                 i++;
1086 #endif
1087             }
1088         }
1089     }
1090 }
1091
1092 static void pr_cmap(FILE *fp, int indent, const char *title,
1093                     gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1094 {
1095     int i,j,nelem;
1096     real dx,idx;
1097         
1098     dx    = 360.0 / cmap_grid->grid_spacing;
1099     nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1100         
1101     if(available(fp,cmap_grid,indent,title))
1102     {
1103         fprintf(fp,"%s\n",title);
1104                 
1105         for(i=0;i<cmap_grid->ngrid;i++)
1106         {
1107             idx = -180.0;
1108             fprintf(fp,"%8s %8s %8s %8s\n","V","dVdx","dVdy","d2dV");
1109                         
1110             fprintf(fp,"grid[%3d]={\n",bShowNumbers?i:-1);
1111                         
1112             for(j=0;j<nelem;j++)
1113             {
1114                 if( (j%cmap_grid->grid_spacing)==0)
1115                 {
1116                     fprintf(fp,"%8.1f\n",idx);
1117                     idx+=dx;
1118                 }
1119                                 
1120                 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4]);
1121                 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+1]);
1122                 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+2]);
1123                 fprintf(fp,"%8.3f\n",cmap_grid->cmapdata[i].cmap[j*4+3]);
1124             }
1125             fprintf(fp,"\n");
1126         }
1127     }
1128         
1129 }
1130
1131 void pr_ffparams(FILE *fp,int indent,const char *title,
1132                  gmx_ffparams_t *ffparams,
1133                  gmx_bool bShowNumbers)
1134 {
1135   int i,j;
1136   
1137   indent=pr_title(fp,indent,title);
1138   (void) pr_indent(fp,indent);
1139   (void) fprintf(fp,"atnr=%d\n",ffparams->atnr);
1140   (void) pr_indent(fp,indent);
1141   (void) fprintf(fp,"ntypes=%d\n",ffparams->ntypes);
1142   for (i=0; i<ffparams->ntypes; i++) {
1143       (void) pr_indent(fp,indent+INDENT);
1144       (void) fprintf(fp,"functype[%d]=%s, ",
1145                      bShowNumbers?i:-1,
1146                      interaction_function[ffparams->functype[i]].name);
1147       pr_iparams(fp,ffparams->functype[i],&ffparams->iparams[i]);
1148   }
1149   (void) pr_double(fp,indent,"reppow",ffparams->reppow);
1150   (void) pr_real(fp,indent,"fudgeQQ",ffparams->fudgeQQ);
1151   pr_cmap(fp,indent,"cmap",&ffparams->cmap_grid,bShowNumbers);
1152 }
1153
1154 void pr_idef(FILE *fp,int indent,const char *title,t_idef *idef, gmx_bool bShowNumbers)
1155 {
1156   int i,j;
1157   
1158   if (available(fp,idef,indent,title)) {  
1159     indent=pr_title(fp,indent,title);
1160     (void) pr_indent(fp,indent);
1161     (void) fprintf(fp,"atnr=%d\n",idef->atnr);
1162     (void) pr_indent(fp,indent);
1163     (void) fprintf(fp,"ntypes=%d\n",idef->ntypes);
1164     for (i=0; i<idef->ntypes; i++) {
1165       (void) pr_indent(fp,indent+INDENT);
1166       (void) fprintf(fp,"functype[%d]=%s, ",
1167                      bShowNumbers?i:-1,
1168                      interaction_function[idef->functype[i]].name);
1169       pr_iparams(fp,idef->functype[i],&idef->iparams[i]);
1170     }
1171     (void) pr_real(fp,indent,"fudgeQQ",idef->fudgeQQ);
1172
1173     for(j=0; (j<F_NRE); j++)
1174       pr_ilist(fp,indent,interaction_function[j].longname,
1175                idef->functype,&idef->il[j],bShowNumbers);
1176   }
1177 }
1178
1179 static int pr_block_title(FILE *fp,int indent,const char *title,t_block *block)
1180 {
1181   int i;
1182
1183   if (available(fp,block,indent,title))
1184     {
1185       indent=pr_title(fp,indent,title);
1186       (void) pr_indent(fp,indent);
1187       (void) fprintf(fp,"nr=%d\n",block->nr);
1188     }
1189   return indent;
1190 }
1191
1192 static int pr_blocka_title(FILE *fp,int indent,const char *title,t_blocka *block)
1193 {
1194   int i;
1195
1196   if (available(fp,block,indent,title))
1197     {
1198       indent=pr_title(fp,indent,title);
1199       (void) pr_indent(fp,indent);
1200       (void) fprintf(fp,"nr=%d\n",block->nr);
1201       (void) pr_indent(fp,indent);
1202       (void) fprintf(fp,"nra=%d\n",block->nra);
1203     }
1204   return indent;
1205 }
1206
1207 static void low_pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block, gmx_bool bShowNumbers)
1208 {
1209   int i;
1210   
1211   if (available(fp,block,indent,title))
1212     {
1213       indent=pr_blocka_title(fp,indent,title,block);
1214       for (i=0; i<=block->nr; i++)
1215         {
1216           (void) pr_indent(fp,indent+INDENT);
1217           (void) fprintf(fp,"%s->index[%d]=%u\n",
1218                          title,bShowNumbers?i:-1,block->index[i]);
1219         }
1220       for (i=0; i<block->nra; i++)
1221         {
1222           (void) pr_indent(fp,indent+INDENT);
1223           (void) fprintf(fp,"%s->a[%d]=%u\n",
1224                          title,bShowNumbers?i:-1,block->a[i]);
1225         }
1226     }
1227 }
1228
1229 void pr_block(FILE *fp,int indent,const char *title,t_block *block,gmx_bool bShowNumbers)
1230 {
1231   int i,j,ok,size,start,end;
1232   
1233   if (available(fp,block,indent,title))
1234     {
1235       indent=pr_block_title(fp,indent,title,block);
1236       start=0;
1237       end=start;
1238       if ((ok=(block->index[start]==0))==0)
1239         (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1240       else
1241         for (i=0; i<block->nr; i++)
1242           {
1243             end=block->index[i+1];
1244             size=pr_indent(fp,indent);
1245             if (end<=start)
1246               size+=fprintf(fp,"%s[%d]={}\n",title,i);
1247             else
1248               size+=fprintf(fp,"%s[%d]={%d..%d}\n",
1249                             title,bShowNumbers?i:-1,
1250                             bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1251             start=end;
1252           }
1253     }
1254 }
1255
1256 void pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block,gmx_bool bShowNumbers)
1257 {
1258   int i,j,ok,size,start,end;
1259   
1260   if (available(fp,block,indent,title))
1261     {
1262       indent=pr_blocka_title(fp,indent,title,block);
1263       start=0;
1264       end=start;
1265       if ((ok=(block->index[start]==0))==0)
1266         (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1267       else
1268         for (i=0; i<block->nr; i++)
1269           {
1270             end=block->index[i+1];
1271             size=pr_indent(fp,indent);
1272             if (end<=start)
1273               size+=fprintf(fp,"%s[%d]={",title,i);
1274             else
1275               size+=fprintf(fp,"%s[%d][%d..%d]={",
1276                             title,bShowNumbers?i:-1,
1277                             bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1278             for (j=start; j<end; j++)
1279               {
1280                 if (j>start) size+=fprintf(fp,", ");
1281                 if ((size)>(USE_WIDTH))
1282                   {
1283                     (void) fprintf(fp,"\n");
1284                     size=pr_indent(fp,indent+INDENT);
1285                   }
1286                 size+=fprintf(fp,"%u",block->a[j]);
1287               }
1288             (void) fprintf(fp,"}\n");
1289             start=end;
1290           }
1291       if ((end!=block->nra)||(!ok)) 
1292         {
1293           (void) pr_indent(fp,indent);
1294           (void) fprintf(fp,"tables inconsistent, dumping complete tables:\n");
1295           low_pr_blocka(fp,indent,title,block,bShowNumbers);
1296         }
1297     }
1298 }
1299
1300 static void pr_strings(FILE *fp,int indent,const char *title,char ***nm,int n, gmx_bool bShowNumbers)
1301 {
1302   int i;
1303
1304   if (available(fp,nm,indent,title))
1305     {  
1306       indent=pr_title_n(fp,indent,title,n);
1307       for (i=0; i<n; i++)
1308         {
1309           (void) pr_indent(fp,indent);
1310           (void) fprintf(fp,"%s[%d]={name=\"%s\"}\n",
1311                          title,bShowNumbers?i:-1,*(nm[i]));
1312         }
1313     }
1314 }
1315
1316 static void pr_strings2(FILE *fp,int indent,const char *title,
1317                         char ***nm,char ***nmB,int n, gmx_bool bShowNumbers)
1318 {
1319   int i;
1320
1321   if (available(fp,nm,indent,title))
1322     {  
1323       indent=pr_title_n(fp,indent,title,n);
1324       for (i=0; i<n; i++)
1325         {
1326           (void) pr_indent(fp,indent);
1327           (void) fprintf(fp,"%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1328                          title,bShowNumbers?i:-1,*(nm[i]),*(nmB[i]));
1329         }
1330     }
1331 }
1332
1333 static void pr_resinfo(FILE *fp,int indent,const char *title,t_resinfo *resinfo,int n, gmx_bool bShowNumbers)
1334 {
1335     int i;
1336     
1337     if (available(fp,resinfo,indent,title))
1338     {  
1339         indent=pr_title_n(fp,indent,title,n);
1340         for (i=0; i<n; i++)
1341         {
1342             (void) pr_indent(fp,indent);
1343             (void) fprintf(fp,"%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1344                            title,bShowNumbers?i:-1,
1345                            *(resinfo[i].name),resinfo[i].nr,
1346                            (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1347         }
1348     }
1349 }
1350
1351 static void pr_atom(FILE *fp,int indent,const char *title,t_atom *atom,int n)
1352 {
1353   int i,j;
1354   
1355   if (available(fp,atom,indent,title)) {  
1356     indent=pr_title_n(fp,indent,title,n);
1357     for (i=0; i<n; i++) {
1358       (void) pr_indent(fp,indent);
1359       fprintf(fp,"%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1360               "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1361               title,i,atom[i].type,atom[i].typeB,ptype_str[atom[i].ptype],
1362               atom[i].m,atom[i].q,atom[i].mB,atom[i].qB,
1363               atom[i].resind,atom[i].atomnumber);
1364     }
1365   }
1366 }
1367
1368 static void pr_grps(FILE *fp,int indent,const char *title,t_grps grps[],
1369                     char **grpname[], gmx_bool bShowNumbers)
1370 {
1371     int i,j;
1372
1373     for(i=0; (i<egcNR); i++)
1374     {
1375         fprintf(fp,"%s[%-12s] nr=%d, name=[",title,gtypes[i],grps[i].nr);
1376         for(j=0; (j<grps[i].nr); j++)
1377         {
1378             fprintf(fp," %s",*(grpname[grps[i].nm_ind[j]]));
1379         }
1380         fprintf(fp,"]\n");
1381     }
1382 }
1383
1384 static void pr_groups(FILE *fp,int indent,const char *title,
1385                       gmx_groups_t *groups,
1386                       gmx_bool bShowNumbers)
1387 {
1388     int grpnr[egcNR];
1389     int nat_max,i,g;
1390
1391     pr_grps(fp,indent,"grp",groups->grps,groups->grpname,bShowNumbers);
1392     pr_strings(fp,indent,"grpname",groups->grpname,groups->ngrpname,bShowNumbers);
1393
1394     (void) pr_indent(fp,indent);
1395     fprintf(fp,"groups          ");
1396     for(g=0; g<egcNR; g++)
1397     {
1398        printf(" %5.5s",gtypes[g]);
1399     }
1400     printf("\n");
1401
1402     (void) pr_indent(fp,indent);
1403     fprintf(fp,"allocated       ");
1404     nat_max = 0;
1405     for(g=0; g<egcNR; g++)
1406     {
1407         printf(" %5d",groups->ngrpnr[g]);
1408         nat_max = max(nat_max,groups->ngrpnr[g]);
1409     }
1410     printf("\n");
1411
1412     if (nat_max == 0)
1413     {
1414         (void) pr_indent(fp,indent);
1415         fprintf(fp,"groupnr[%5s] =","*");
1416         for(g=0; g<egcNR; g++)
1417         {
1418             fprintf(fp,"  %3d ",0);
1419         }
1420         fprintf(fp,"\n");
1421     }
1422     else
1423     {
1424         for(i=0; i<nat_max; i++)
1425         {
1426             (void) pr_indent(fp,indent);
1427             fprintf(fp,"groupnr[%5d] =",i);
1428             for(g=0; g<egcNR; g++)
1429             {
1430                 fprintf(fp,"  %3d ",
1431                         groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1432             }
1433             fprintf(fp,"\n");
1434         }
1435     }
1436 }
1437
1438 void pr_atoms(FILE *fp,int indent,const char *title,t_atoms *atoms, 
1439               gmx_bool bShownumbers)
1440 {
1441   if (available(fp,atoms,indent,title))
1442     {
1443       indent=pr_title(fp,indent,title);
1444       pr_atom(fp,indent,"atom",atoms->atom,atoms->nr);
1445       pr_strings(fp,indent,"atom",atoms->atomname,atoms->nr,bShownumbers);
1446       pr_strings2(fp,indent,"type",atoms->atomtype,atoms->atomtypeB,atoms->nr,bShownumbers);
1447       pr_resinfo(fp,indent,"residue",atoms->resinfo,atoms->nres,bShownumbers);
1448     }
1449 }
1450
1451
1452 void pr_atomtypes(FILE *fp,int indent,const char *title,t_atomtypes *atomtypes, 
1453                   gmx_bool bShowNumbers)
1454 {
1455   int i;
1456   if (available(fp,atomtypes,indent,title)) 
1457   {
1458     indent=pr_title(fp,indent,title);
1459     for(i=0;i<atomtypes->nr;i++) {
1460       pr_indent(fp,indent);
1461                 fprintf(fp,
1462                                 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1463                                 bShowNumbers?i:-1,atomtypes->radius[i],atomtypes->vol[i],
1464                                 atomtypes->gb_radius[i],
1465                                 atomtypes->surftens[i],atomtypes->atomnumber[i],atomtypes->S_hct[i]);
1466     }
1467   }
1468 }
1469
1470 static void pr_moltype(FILE *fp,int indent,const char *title,
1471                        gmx_moltype_t *molt,int n,
1472                        gmx_ffparams_t *ffparams,
1473                        gmx_bool bShowNumbers)
1474 {
1475     int j;
1476
1477     indent = pr_title_n(fp,indent,title,n);
1478     (void) pr_indent(fp,indent);
1479     (void) fprintf(fp,"name=\"%s\"\n",*(molt->name));
1480     pr_atoms(fp,indent,"atoms",&(molt->atoms),bShowNumbers);
1481     pr_block(fp,indent,"cgs",&molt->cgs, bShowNumbers);
1482     pr_blocka(fp,indent,"excls",&molt->excls, bShowNumbers);
1483     for(j=0; (j<F_NRE); j++) {
1484         pr_ilist(fp,indent,interaction_function[j].longname,
1485                  ffparams->functype,&molt->ilist[j],bShowNumbers);
1486     }
1487 }
1488
1489 static void pr_molblock(FILE *fp,int indent,const char *title,
1490                         gmx_molblock_t *molb,int n,
1491                         gmx_moltype_t *molt,
1492                         gmx_bool bShowNumbers)
1493 {
1494     indent = pr_title_n(fp,indent,title,n);
1495     (void) pr_indent(fp,indent);
1496     (void) fprintf(fp,"%-20s = %d \"%s\"\n",
1497                    "moltype",molb->type,*(molt[molb->type].name));
1498     pr_int(fp,indent,"#molecules",molb->nmol);
1499     pr_int(fp,indent,"#atoms_mol",molb->natoms_mol);
1500     pr_int(fp,indent,"#posres_xA",molb->nposres_xA);
1501     if (molb->nposres_xA > 0) {
1502         pr_rvecs(fp,indent,"posres_xA",molb->posres_xA,molb->nposres_xA);
1503     }
1504     pr_int(fp,indent,"#posres_xB",molb->nposres_xB);
1505     if (molb->nposres_xB > 0) {
1506         pr_rvecs(fp,indent,"posres_xB",molb->posres_xB,molb->nposres_xB);
1507     }
1508 }
1509
1510 void pr_mtop(FILE *fp,int indent,const char *title,gmx_mtop_t *mtop,
1511              gmx_bool bShowNumbers)
1512 {
1513     int mt,mb;
1514
1515     if (available(fp,mtop,indent,title)) {
1516         indent=pr_title(fp,indent,title);
1517         (void) pr_indent(fp,indent);
1518         (void) fprintf(fp,"name=\"%s\"\n",*(mtop->name));
1519         pr_int(fp,indent,"#atoms",mtop->natoms);
1520         pr_int(fp,indent,"#molblock",mtop->nmolblock);
1521         for(mb=0; mb<mtop->nmolblock; mb++) {
1522             pr_molblock(fp,indent,"molblock",&mtop->molblock[mb],mb,
1523                         mtop->moltype,bShowNumbers);
1524         }
1525         pr_ffparams(fp,indent,"ffparams",&(mtop->ffparams),bShowNumbers);
1526         pr_atomtypes(fp,indent,"atomtypes",&(mtop->atomtypes),bShowNumbers);
1527         for(mt=0; mt<mtop->nmoltype; mt++) {
1528             pr_moltype(fp,indent,"moltype",&mtop->moltype[mt],mt,
1529                        &mtop->ffparams,bShowNumbers);
1530         }
1531         pr_groups(fp,indent,"groups",&mtop->groups,bShowNumbers);
1532     }
1533 }
1534
1535 void pr_top(FILE *fp,int indent,const char *title,t_topology *top, gmx_bool bShowNumbers)
1536 {
1537   if (available(fp,top,indent,title)) {
1538     indent=pr_title(fp,indent,title);
1539     (void) pr_indent(fp,indent);
1540     (void) fprintf(fp,"name=\"%s\"\n",*(top->name));
1541     pr_atoms(fp,indent,"atoms",&(top->atoms),bShowNumbers);
1542     pr_atomtypes(fp,indent,"atomtypes",&(top->atomtypes),bShowNumbers);
1543     pr_block(fp,indent,"cgs",&top->cgs, bShowNumbers);
1544     pr_block(fp,indent,"mols",&top->mols, bShowNumbers);
1545     pr_blocka(fp,indent,"excls",&top->excls, bShowNumbers);
1546     pr_idef(fp,indent,"idef",&top->idef,bShowNumbers);
1547   }
1548 }
1549
1550 void pr_header(FILE *fp,int indent,const char *title,t_tpxheader *sh)
1551 {
1552   char buf[22];
1553     
1554   if (available(fp,sh,indent,title))
1555     {
1556       indent=pr_title(fp,indent,title);
1557       pr_indent(fp,indent);
1558       fprintf(fp,"bIr    = %spresent\n",sh->bIr?"":"not ");
1559       pr_indent(fp,indent);
1560       fprintf(fp,"bBox   = %spresent\n",sh->bBox?"":"not ");
1561       pr_indent(fp,indent);
1562       fprintf(fp,"bTop   = %spresent\n",sh->bTop?"":"not ");
1563       pr_indent(fp,indent);
1564       fprintf(fp,"bX     = %spresent\n",sh->bX?"":"not ");
1565       pr_indent(fp,indent);
1566       fprintf(fp,"bV     = %spresent\n",sh->bV?"":"not ");
1567       pr_indent(fp,indent);
1568       fprintf(fp,"bF     = %spresent\n",sh->bF?"":"not ");
1569       
1570       pr_indent(fp,indent);
1571       fprintf(fp,"natoms = %d\n",sh->natoms);
1572       pr_indent(fp,indent);
1573       fprintf(fp,"lambda = %e\n",sh->lambda);
1574     }
1575 }
1576
1577 void pr_commrec(FILE *fp,int indent,t_commrec *cr)
1578 {
1579   pr_indent(fp,indent);
1580   fprintf(fp,"commrec:\n");
1581   indent+=2;
1582   pr_indent(fp,indent);
1583   fprintf(fp,"nodeid    = %d\n",cr->nodeid);
1584   pr_indent(fp,indent);
1585   fprintf(fp,"nnodes    = %d\n",cr->nnodes);
1586   pr_indent(fp,indent);
1587   fprintf(fp,"npmenodes = %d\n",cr->npmenodes);
1588   /*
1589   pr_indent(fp,indent);
1590   fprintf(fp,"threadid  = %d\n",cr->threadid);
1591   pr_indent(fp,indent);
1592   fprintf(fp,"nthreads  = %d\n",cr->nthreads);
1593   */
1594 }