added mdp parameter for FE dH tables
[gromacs.git] / src / gmxlib / txtdump.c
blob87ab255173f906445d15534b963a47289a55155e
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 /* This file is completely threadsafe - please keep it that way! */
41 #ifdef GMX_THREADS
42 #include <thread_mpi.h>
43 #endif
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"
55 int pr_indent(FILE *fp,int n)
57 int i;
59 for (i=0; i<n; i++) (void) fprintf(fp," ");
60 return n;
63 int available(FILE *fp,void *p,int indent,const char *title)
65 if (!p) {
66 if (indent > 0)
67 pr_indent(fp,indent);
68 (void) fprintf(fp,"%s: not available\n",title);
70 return (p!=NULL);
73 int pr_title(FILE *fp,int indent,const char *title)
75 (void) pr_indent(fp,indent);
76 (void) fprintf(fp,"%s:\n",title);
77 return (indent+INDENT);
80 int pr_title_n(FILE *fp,int indent,const char *title,int n)
82 (void) pr_indent(fp,indent);
83 (void) fprintf(fp,"%s (%d):\n",title,n);
84 return (indent+INDENT);
87 int pr_title_nxn(FILE *fp,int indent,const char *title,int n1,int n2)
89 (void) pr_indent(fp,indent);
90 (void) fprintf(fp,"%s (%dx%d):\n",title,n1,n2);
91 return (indent+INDENT);
94 void pr_ivec(FILE *fp,int indent,const char *title,int vec[],int n, bool bShowNumbers)
96 int i;
98 if (available(fp,vec,indent,title))
100 indent=pr_title_n(fp,indent,title,n);
101 for (i=0; i<n; i++)
103 (void) pr_indent(fp,indent);
104 (void) fprintf(fp,"%s[%d]=%d\n",title,bShowNumbers?i:-1,vec[i]);
109 void pr_ivec_block(FILE *fp,int indent,const char *title,int vec[],int n, bool bShowNumbers)
111 int i,j;
113 if (available(fp,vec,indent,title))
115 indent=pr_title_n(fp,indent,title,n);
116 i = 0;
117 while (i < n)
119 j = i+1;
120 while (j < n && vec[j] == vec[j-1]+1)
122 j++;
124 /* Print consecutive groups of 3 or more as blocks */
125 if (j - i < 3)
127 while(i < j)
129 (void) pr_indent(fp,indent);
130 (void) fprintf(fp,"%s[%d]=%d\n",
131 title,bShowNumbers?i:-1,vec[i]);
132 i++;
135 else
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;
149 void pr_bvec(FILE *fp,int indent,const char *title,bool vec[],int n, bool bShowNumbers)
151 int i;
153 if (available(fp,vec,indent,title))
155 indent=pr_title_n(fp,indent,title,n);
156 for (i=0; i<n; i++)
158 (void) pr_indent(fp,indent);
159 (void) fprintf(fp,"%s[%d]=%s\n",title,bShowNumbers?i:-1,
160 BOOL(vec[i]));
165 void pr_ivecs(FILE *fp,int indent,const char *title,ivec vec[],int n, bool bShowNumbers)
167 int i,j;
169 if (available(fp,vec,indent,title))
171 indent=pr_title_nxn(fp,indent,title,n,DIM);
172 for (i=0; i<n; i++)
174 (void) pr_indent(fp,indent);
175 (void) fprintf(fp,"%s[%d]={",title,bShowNumbers?i:-1);
176 for (j=0; j<DIM; j++)
178 if (j!=0) (void) fprintf(fp,", ");
179 fprintf(fp,"%d",vec[i][j]);
181 (void) fprintf(fp,"}\n");
186 void pr_rvec(FILE *fp,int indent,const char *title,real vec[],int n, bool bShowNumbers)
188 int i;
190 if (available(fp,vec,indent,title))
192 indent=pr_title_n(fp,indent,title,n);
193 for (i=0; i<n; i++)
195 pr_indent(fp,indent);
196 fprintf(fp,"%s[%d]=%12.5e\n",title,bShowNumbers?i:-1,vec[i]);
201 void pr_dvec(FILE *fp,int indent,const char *title,double vec[],int n, bool bShowNumbers)
203 int i;
205 if (available(fp,vec,indent,title))
207 indent=pr_title_n(fp,indent,title,n);
208 for (i=0; i<n; i++)
210 pr_indent(fp,indent);
211 fprintf(fp,"%s[%d]=%12.5e\n",title,bShowNumbers?i:-1,vec[i]);
218 void pr_mat(FILE *fp,int indent,char *title,matrix m)
220 int i,j;
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]);
233 void pr_rvecs_len(FILE *fp,int indent,const char *title,rvec vec[],int n)
235 int i,j;
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]);
247 (void) fprintf(fp,"} len=%12.5e\n",norm(vec[i]));
252 void pr_rvecs(FILE *fp,int indent,const char *title,rvec vec[],int n)
254 const char *fshort = "%12.5e";
255 const char *flong = "%15.8e";
256 const char *format;
257 int i,j;
259 if (getenv("LONGFORMAT") != NULL)
260 format = flong;
261 else
262 format = fshort;
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]);
274 (void) fprintf(fp,"}\n");
280 void pr_reals(FILE *fp,int indent,const char *title,real *vec,int n)
282 int i;
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");
293 void pr_doubles(FILE *fp,int indent,const char *title,double *vec,int n)
295 int i;
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");
306 static void pr_int(FILE *fp,int indent,const char *title,int i)
308 pr_indent(fp,indent);
309 fprintf(fp,"%-20s = %d\n",title,i);
312 static void pr_gmx_large_int(FILE *fp,int indent,const char *title,gmx_large_int_t i)
314 char buf[STEPSTRSIZE];
316 pr_indent(fp,indent);
317 fprintf(fp,"%-20s = %s\n",title,gmx_step_str(i,buf));
320 static void pr_real(FILE *fp,int indent,const char *title,real r)
322 pr_indent(fp,indent);
323 fprintf(fp,"%-20s = %g\n",title,r);
326 static void pr_double(FILE *fp,int indent,const char *title,double d)
328 pr_indent(fp,indent);
329 fprintf(fp,"%-20s = %g\n",title,d);
332 static void pr_str(FILE *fp,int indent,const char *title,const char *s)
334 pr_indent(fp,indent);
335 fprintf(fp,"%-20s = %s\n",title,s);
338 void pr_qm_opts(FILE *fp,int indent,const char *title,t_grpopts *opts)
340 int i,m,j;
342 fprintf(fp,"%s:\n",title);
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);
361 static void pr_grp_opts(FILE *out,int indent,const char *title,t_grpopts *opts,
362 bool bMDPformat)
364 int i,m,j;
366 if (!bMDPformat)
367 fprintf(out,"%s:\n",title);
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");
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");
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");
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");
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");
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");
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");
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");
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");
434 fflush(out);
437 static void pr_matrix(FILE *fp,int indent,const char *title,rvec *m,
438 bool bMDPformat)
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);
447 static void pr_cosine(FILE *fp,int indent,const char *title,t_cosines *cos,
448 bool bMDPformat)
450 int j;
452 if (bMDPformat) {
453 fprintf(fp,"%s = %d\n",title,cos->n);
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");
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)
480 static void pr_pullgrp(FILE *fp,int indent,int g,t_pullgrp *pg)
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);
495 static void pr_pull(FILE *fp,int indent,t_pull *pull)
497 int g;
499 PS("pull_geometry",EPULLGEOM(pull->eGeom));
500 pr_ivec(fp,indent,"pull_dim",pull->dim,DIM,TRUE);
501 PR("pull_r1",pull->cyl_r1);
502 PR("pull_r0",pull->cyl_r0);
503 PR("pull_constr_tol",pull->constr_tol);
504 PI("pull_nstxout",pull->nstxout);
505 PI("pull_nstfout",pull->nstfout);
506 PI("pull_ngrp",pull->ngrp);
507 for(g=0; g<pull->ngrp+1; g++)
508 pr_pullgrp(fp,indent,g,&pull->grp[g]);
511 void pr_inputrec(FILE *fp,int indent,const char *title,t_inputrec *ir,
512 bool bMDPformat)
514 const char *infbuf="inf";
515 int i;
517 if (available(fp,ir,indent,title)) {
518 if (!bMDPformat)
519 indent=pr_title(fp,indent,title);
520 PS("integrator",EI(ir->eI));
521 PSTEP("nsteps",ir->nsteps);
522 PSTEP("init_step",ir->init_step);
523 PI("nstcalcenergy",ir->nstcalcenergy);
524 PS("ns_type",ENS(ir->ns_type));
525 PI("nstlist",ir->nstlist);
526 PI("ndelta",ir->ndelta);
527 PI("nstcomm",ir->nstcomm);
528 PS("comm_mode",ECOM(ir->comm_mode));
529 PI("nstlog",ir->nstlog);
530 PI("nstxout",ir->nstxout);
531 PI("nstvout",ir->nstvout);
532 PI("nstfout",ir->nstfout);
533 PI("nstenergy",ir->nstenergy);
534 PI("nstxtcout",ir->nstxtcout);
535 PR("init_t",ir->init_t);
536 PR("delta_t",ir->delta_t);
538 PR("xtcprec",ir->xtcprec);
539 PI("nkx",ir->nkx);
540 PI("nky",ir->nky);
541 PI("nkz",ir->nkz);
542 PI("pme_order",ir->pme_order);
543 PR("ewald_rtol",ir->ewald_rtol);
544 PR("ewald_geometry",ir->ewald_geometry);
545 PR("epsilon_surface",ir->epsilon_surface);
546 PS("optimize_fft",BOOL(ir->bOptFFT));
547 PS("ePBC",EPBC(ir->ePBC));
548 PS("bPeriodicMols",BOOL(ir->bPeriodicMols));
549 PS("bContinuation",BOOL(ir->bContinuation));
550 PS("bShakeSOR",BOOL(ir->bShakeSOR));
551 PS("etc",ETCOUPLTYPE(ir->etc));
552 PI("nsttcouple",ir->nsttcouple);
553 PS("epc",EPCOUPLTYPE(ir->epc));
554 PS("epctype",EPCOUPLTYPETYPE(ir->epct));
555 PI("nstpcouple",ir->nstpcouple);
556 PR("tau_p",ir->tau_p);
557 pr_matrix(fp,indent,"ref_p",ir->ref_p,bMDPformat);
558 pr_matrix(fp,indent,"compress",ir->compress,bMDPformat);
559 PS("refcoord_scaling",EREFSCALINGTYPE(ir->refcoord_scaling));
560 if (bMDPformat)
561 fprintf(fp,"posres_com = %g %g %g\n",ir->posres_com[XX],
562 ir->posres_com[YY],ir->posres_com[ZZ]);
563 else
564 pr_rvec(fp,indent,"posres_com",ir->posres_com,DIM,TRUE);
565 if (bMDPformat)
566 fprintf(fp,"posres_comB = %g %g %g\n",ir->posres_comB[XX],
567 ir->posres_comB[YY],ir->posres_comB[ZZ]);
568 else
569 pr_rvec(fp,indent,"posres_comB",ir->posres_comB,DIM,TRUE);
570 PI("andersen_seed",ir->andersen_seed);
571 PR("rlist",ir->rlist);
572 PR("rlistlong",ir->rlistlong);
573 PR("rtpi",ir->rtpi);
574 PS("coulombtype",EELTYPE(ir->coulombtype));
575 PR("rcoulomb_switch",ir->rcoulomb_switch);
576 PR("rcoulomb",ir->rcoulomb);
577 PS("vdwtype",EVDWTYPE(ir->vdwtype));
578 PR("rvdw_switch",ir->rvdw_switch);
579 PR("rvdw",ir->rvdw);
580 if (ir->epsilon_r != 0)
581 PR("epsilon_r",ir->epsilon_r);
582 else
583 PS("epsilon_r",infbuf);
584 if (ir->epsilon_rf != 0)
585 PR("epsilon_rf",ir->epsilon_rf);
586 else
587 PS("epsilon_rf",infbuf);
588 PR("tabext",ir->tabext);
589 PS("implicit_solvent",EIMPLICITSOL(ir->implicit_solvent));
590 PS("gb_algorithm",EGBALGORITHM(ir->gb_algorithm));
591 PR("gb_epsilon_solvent",ir->gb_epsilon_solvent);
592 PI("nstgbradii",ir->nstgbradii);
593 PR("rgbradii",ir->rgbradii);
594 PR("gb_saltconc",ir->gb_saltconc);
595 PR("gb_obc_alpha",ir->gb_obc_alpha);
596 PR("gb_obc_beta",ir->gb_obc_beta);
597 PR("gb_obc_gamma",ir->gb_obc_gamma);
598 PR("gb_dielectric_offset",ir->gb_dielectric_offset);
599 PS("sa_algorithm",ESAALGORITHM(ir->gb_algorithm));
600 PR("sa_surface_tension",ir->sa_surface_tension);
602 PS("DispCorr",EDISPCORR(ir->eDispCorr));
603 PS("free_energy",EFEPTYPE(ir->efep));
604 PR("init_lambda",ir->init_lambda);
605 PR("delta_lambda",ir->delta_lambda);
606 if (!bMDPformat)
608 PI("n_foreign_lambda",ir->n_flambda);
610 if (ir->n_flambda > 0)
612 pr_indent(fp,indent);
613 fprintf(fp,"foreign_lambda%s",bMDPformat ? " = " : ":");
614 for(i=0; i<ir->n_flambda; i++)
616 fprintf(fp," %10g",ir->flambda[i]);
618 fprintf(fp,"\n");
620 PR("sc_alpha",ir->sc_alpha);
621 PI("sc_power",ir->sc_power);
622 PR("sc_sigma",ir->sc_sigma);
623 PI("nstdhdl", ir->nstdhdl);
624 PI("dh_table_size", ir->dh_table_size);
625 PD("dh_table_spacing", ir->dh_table_spacing);
627 PI("nwall",ir->nwall);
628 PS("wall_type",EWALLTYPE(ir->wall_type));
629 PI("wall_atomtype[0]",ir->wall_atomtype[0]);
630 PI("wall_atomtype[1]",ir->wall_atomtype[1]);
631 PR("wall_density[0]",ir->wall_density[0]);
632 PR("wall_density[1]",ir->wall_density[1]);
633 PR("wall_ewald_zfac",ir->wall_ewald_zfac);
635 PS("pull",EPULLTYPE(ir->ePull));
636 if (ir->ePull != epullNO)
637 pr_pull(fp,indent,ir->pull);
639 PS("disre",EDISRETYPE(ir->eDisre));
640 PS("disre_weighting",EDISREWEIGHTING(ir->eDisreWeighting));
641 PS("disre_mixed",BOOL(ir->bDisreMixed));
642 PR("dr_fc",ir->dr_fc);
643 PR("dr_tau",ir->dr_tau);
644 PR("nstdisreout",ir->nstdisreout);
645 PR("orires_fc",ir->orires_fc);
646 PR("orires_tau",ir->orires_tau);
647 PR("nstorireout",ir->nstorireout);
649 PR("dihre-fc",ir->dihre_fc);
651 PR("em_stepsize",ir->em_stepsize);
652 PR("em_tol",ir->em_tol);
653 PI("niter",ir->niter);
654 PR("fc_stepsize",ir->fc_stepsize);
655 PI("nstcgsteep",ir->nstcgsteep);
656 PI("nbfgscorr",ir->nbfgscorr);
658 PS("ConstAlg",ECONSTRTYPE(ir->eConstrAlg));
659 PR("shake_tol",ir->shake_tol);
660 PI("lincs_order",ir->nProjOrder);
661 PR("lincs_warnangle",ir->LincsWarnAngle);
662 PI("lincs_iter",ir->nLincsIter);
663 PR("bd_fric",ir->bd_fric);
664 PI("ld_seed",ir->ld_seed);
665 PR("cos_accel",ir->cos_accel);
666 pr_matrix(fp,indent,"deform",ir->deform,bMDPformat);
667 PI("userint1",ir->userint1);
668 PI("userint2",ir->userint2);
669 PI("userint3",ir->userint3);
670 PI("userint4",ir->userint4);
671 PR("userreal1",ir->userreal1);
672 PR("userreal2",ir->userreal2);
673 PR("userreal3",ir->userreal3);
674 PR("userreal4",ir->userreal4);
675 pr_grp_opts(fp,indent,"grpopts",&(ir->opts),bMDPformat);
676 pr_cosine(fp,indent,"efield-x",&(ir->ex[XX]),bMDPformat);
677 pr_cosine(fp,indent,"efield-xt",&(ir->et[XX]),bMDPformat);
678 pr_cosine(fp,indent,"efield-y",&(ir->ex[YY]),bMDPformat);
679 pr_cosine(fp,indent,"efield-yt",&(ir->et[YY]),bMDPformat);
680 pr_cosine(fp,indent,"efield-z",&(ir->ex[ZZ]),bMDPformat);
681 pr_cosine(fp,indent,"efield-zt",&(ir->et[ZZ]),bMDPformat);
682 PS("bQMMM",BOOL(ir->bQMMM));
683 PI("QMconstraints",ir->QMconstraints);
684 PI("QMMMscheme",ir->QMMMscheme);
685 PR("scalefactor",ir->scalefactor);
686 pr_qm_opts(fp,indent,"qm_opts",&(ir->opts));
689 #undef PS
690 #undef PR
691 #undef PI
693 static void pr_harm(FILE *fp,t_iparams *iparams,const char *r,const char *kr)
695 fprintf(fp,"%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
696 r,iparams->harmonic.rA,kr,iparams->harmonic.krA,
697 r,iparams->harmonic.rB,kr,iparams->harmonic.krB);
700 void pr_iparams(FILE *fp,t_functype ftype,t_iparams *iparams)
702 int i;
703 real VA[4],VB[4],*rbcA,*rbcB;
705 switch (ftype) {
706 case F_ANGLES:
707 case F_G96ANGLES:
708 pr_harm(fp,iparams,"th","ct");
709 break;
710 case F_CROSS_BOND_BONDS:
711 fprintf(fp,"r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
712 iparams->cross_bb.r1e,iparams->cross_bb.r2e,
713 iparams->cross_bb.krr);
714 break;
715 case F_CROSS_BOND_ANGLES:
716 fprintf(fp,"r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
717 iparams->cross_ba.r1e,iparams->cross_ba.r2e,
718 iparams->cross_ba.r3e,iparams->cross_ba.krt);
719 break;
720 case F_UREY_BRADLEY:
721 fprintf(fp,"theta=%15.8e, ktheta=%15.8e, r13=%15.8e, kUB=%15.8e\n",
722 iparams->u_b.theta,iparams->u_b.ktheta,iparams->u_b.r13,iparams->u_b.kUB);
723 break;
724 case F_QUARTIC_ANGLES:
725 fprintf(fp,"theta=%15.8e",iparams->qangle.theta);
726 for(i=0; i<5; i++)
727 fprintf(fp,", c%c=%15.8e",'0'+i,iparams->qangle.c[i]);
728 fprintf(fp,"\n");
729 break;
730 case F_BHAM:
731 fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
732 iparams->bham.a,iparams->bham.b,iparams->bham.c);
733 break;
734 case F_BONDS:
735 case F_G96BONDS:
736 case F_HARMONIC:
737 pr_harm(fp,iparams,"b0","cb");
738 break;
739 case F_IDIHS:
740 pr_harm(fp,iparams,"xi","cx");
741 break;
742 case F_MORSE:
743 fprintf(fp,"b0=%15.8e, cb=%15.8e, beta=%15.8e\n",
744 iparams->morse.b0,iparams->morse.cb,iparams->morse.beta);
745 break;
746 case F_CUBICBONDS:
747 fprintf(fp,"b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
748 iparams->cubic.b0,iparams->cubic.kb,iparams->cubic.kcub);
749 break;
750 case F_CONNBONDS:
751 fprintf(fp,"\n");
752 break;
753 case F_FENEBONDS:
754 fprintf(fp,"bm=%15.8e, kb=%15.8e\n",iparams->fene.bm,iparams->fene.kb);
755 break;
756 case F_RESTRBONDS:
757 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",
758 iparams->restraint.lowA,iparams->restraint.up1A,
759 iparams->restraint.up2A,iparams->restraint.kA,
760 iparams->restraint.lowB,iparams->restraint.up1B,
761 iparams->restraint.up2B,iparams->restraint.kB);
762 break;
763 case F_TABBONDS:
764 case F_TABBONDSNC:
765 case F_TABANGLES:
766 case F_TABDIHS:
767 fprintf(fp,"tab=%d, kA=%15.8e, kB=%15.8e\n",
768 iparams->tab.table,iparams->tab.kA,iparams->tab.kB);
769 break;
770 case F_POLARIZATION:
771 fprintf(fp,"alpha=%15.8e\n",iparams->polarize.alpha);
772 break;
773 case F_THOLE_POL:
774 fprintf(fp,"a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
775 iparams->thole.a,iparams->thole.alpha1,iparams->thole.alpha2,
776 iparams->thole.rfac);
777 break;
778 case F_WATER_POL:
779 fprintf(fp,"al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
780 iparams->wpol.al_x,iparams->wpol.al_y,iparams->wpol.al_z,
781 iparams->wpol.rOH,iparams->wpol.rHH,iparams->wpol.rOD);
782 break;
783 case F_LJ:
784 fprintf(fp,"c6=%15.8e, c12=%15.8e\n",iparams->lj.c6,iparams->lj.c12);
785 break;
786 case F_LJ14:
787 fprintf(fp,"c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
788 iparams->lj14.c6A,iparams->lj14.c12A,
789 iparams->lj14.c6B,iparams->lj14.c12B);
790 break;
791 case F_LJC14_Q:
792 fprintf(fp,"fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
793 iparams->ljc14.fqq,
794 iparams->ljc14.qi,iparams->ljc14.qj,
795 iparams->ljc14.c6,iparams->ljc14.c12);
796 break;
797 case F_LJC_PAIRS_NB:
798 fprintf(fp,"qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
799 iparams->ljcnb.qi,iparams->ljcnb.qj,
800 iparams->ljcnb.c6,iparams->ljcnb.c12);
801 break;
802 case F_PDIHS:
803 case F_ANGRES:
804 case F_ANGRESZ:
805 fprintf(fp,"phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
806 iparams->pdihs.phiA,iparams->pdihs.cpA,
807 iparams->pdihs.phiB,iparams->pdihs.cpB,
808 iparams->pdihs.mult);
809 break;
810 case F_DISRES:
811 fprintf(fp,"label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
812 iparams->disres.label,iparams->disres.type,
813 iparams->disres.low,iparams->disres.up1,
814 iparams->disres.up2,iparams->disres.kfac);
815 break;
816 case F_ORIRES:
817 fprintf(fp,"ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
818 iparams->orires.ex,iparams->orires.label,iparams->orires.power,
819 iparams->orires.c,iparams->orires.obs,iparams->orires.kfac);
820 break;
821 case F_DIHRES:
822 fprintf(fp,"label=%d, power=%4d phi=%15.8e, dphi=%15.8e, kfac=%15.8e)\n",
823 iparams->dihres.label,iparams->dihres.power,
824 iparams->dihres.phi,iparams->dihres.dphi,iparams->dihres.kfac);
825 break;
826 case F_POSRES:
827 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",
828 iparams->posres.pos0A[XX],iparams->posres.pos0A[YY],
829 iparams->posres.pos0A[ZZ],iparams->posres.fcA[XX],
830 iparams->posres.fcA[YY],iparams->posres.fcA[ZZ],
831 iparams->posres.pos0B[XX],iparams->posres.pos0B[YY],
832 iparams->posres.pos0B[ZZ],iparams->posres.fcB[XX],
833 iparams->posres.fcB[YY],iparams->posres.fcB[ZZ]);
834 break;
835 case F_RBDIHS:
836 for (i=0; i<NR_RBDIHS; i++)
837 fprintf(fp,"%srbcA[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcA[i]);
838 fprintf(fp,"\n");
839 for (i=0; i<NR_RBDIHS; i++)
840 fprintf(fp,"%srbcB[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcB[i]);
841 fprintf(fp,"\n");
842 break;
843 case F_FOURDIHS:
844 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
845 * OPLS potential constants back.
847 rbcA = iparams->rbdihs.rbcA;
848 rbcB = iparams->rbdihs.rbcB;
850 VA[3] = -0.25*rbcA[4];
851 VA[2] = -0.5*rbcA[3];
852 VA[1] = 4.0*VA[3]-rbcA[2];
853 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
855 VB[3] = -0.25*rbcB[4];
856 VB[2] = -0.5*rbcB[3];
857 VB[1] = 4.0*VB[3]-rbcB[2];
858 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
860 for (i=0; i<NR_FOURDIHS; i++)
861 fprintf(fp,"%sFourA[%d]=%15.8e",i==0?"":", ",i,VA[i]);
862 fprintf(fp,"\n");
863 for (i=0; i<NR_FOURDIHS; i++)
864 fprintf(fp,"%sFourB[%d]=%15.8e",i==0?"":", ",i,VB[i]);
865 fprintf(fp,"\n");
866 break;
868 case F_CONSTR:
869 case F_CONSTRNC:
870 fprintf(fp,"dA=%15.8e, dB=%15.8e\n",iparams->constr.dA,iparams->constr.dB);
871 break;
872 case F_SETTLE:
873 fprintf(fp,"doh=%15.8e, dhh=%15.8e\n",iparams->settle.doh,
874 iparams->settle.dhh);
875 break;
876 case F_VSITE2:
877 fprintf(fp,"a=%15.8e\n",iparams->vsite.a);
878 break;
879 case F_VSITE3:
880 case F_VSITE3FD:
881 case F_VSITE3FAD:
882 fprintf(fp,"a=%15.8e, b=%15.8e\n",iparams->vsite.a,iparams->vsite.b);
883 break;
884 case F_VSITE3OUT:
885 case F_VSITE4FD:
886 case F_VSITE4FDN:
887 fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
888 iparams->vsite.a,iparams->vsite.b,iparams->vsite.c);
889 break;
890 case F_VSITEN:
891 fprintf(fp,"n=%2d, a=%15.8e\n",iparams->vsiten.n,iparams->vsiten.a);
892 break;
893 case F_GB12:
894 case F_GB13:
895 case F_GB14:
896 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);
897 break;
898 case F_CMAP:
899 fprintf(fp, "cmapA=%1d, cmapB=%1d\n",iparams->cmap.cmapA, iparams->cmap.cmapB);
900 break;
901 default:
902 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
903 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
907 void pr_ilist(FILE *fp,int indent,const char *title,
908 t_functype *functype,t_ilist *ilist, bool bShowNumbers)
910 int i,j,k,type,ftype;
911 t_iatom *iatoms;
913 if (available(fp,ilist,indent,title) && ilist->nr > 0)
915 indent=pr_title(fp,indent,title);
916 (void) pr_indent(fp,indent);
917 fprintf(fp,"nr: %d\n",ilist->nr);
918 if (ilist->nr > 0) {
919 (void) pr_indent(fp,indent);
920 fprintf(fp,"iatoms:\n");
921 iatoms=ilist->iatoms;
922 for (i=j=0; i<ilist->nr;) {
923 #ifndef DEBUG
924 (void) pr_indent(fp,indent+INDENT);
925 type=*(iatoms++);
926 ftype=functype[type];
927 (void) fprintf(fp,"%d type=%d (%s)",
928 bShowNumbers?j:-1,bShowNumbers?type:-1,
929 interaction_function[ftype].name);
930 j++;
931 for (k=0; k<interaction_function[ftype].nratoms; k++)
932 (void) fprintf(fp," %u",*(iatoms++));
933 (void) fprintf(fp,"\n");
934 i+=1+interaction_function[ftype].nratoms;
935 #else
936 fprintf(fp,"%5d%5d\n",i,iatoms[i]);
937 i++;
938 #endif
944 static void pr_cmap(FILE *fp, int indent, const char *title,
945 gmx_cmap_t *cmap_grid, bool bShowNumbers)
947 int i,j,nelem;
948 real dx,idx;
950 dx = 360.0 / cmap_grid->grid_spacing;
951 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
953 if(available(fp,cmap_grid,indent,title))
955 fprintf(fp,"%s\n",title);
957 for(i=0;i<cmap_grid->ngrid;i++)
959 idx = -180.0;
960 fprintf(fp,"%8s %8s %8s %8s\n","V","dVdx","dVdy","d2dV");
962 fprintf(fp,"grid[%3d]={\n",bShowNumbers?i:-1);
964 for(j=0;j<nelem;j++)
966 if( (j%cmap_grid->grid_spacing)==0)
968 fprintf(fp,"%8.1f\n",idx);
969 idx+=dx;
972 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4]);
973 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+1]);
974 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+2]);
975 fprintf(fp,"%8.3f\n",cmap_grid->cmapdata[i].cmap[j*4+3]);
977 fprintf(fp,"\n");
983 void pr_ffparams(FILE *fp,int indent,const char *title,
984 gmx_ffparams_t *ffparams,
985 bool bShowNumbers)
987 int i,j;
989 indent=pr_title(fp,indent,title);
990 (void) pr_indent(fp,indent);
991 (void) fprintf(fp,"atnr=%d\n",ffparams->atnr);
992 (void) pr_indent(fp,indent);
993 (void) fprintf(fp,"ntypes=%d\n",ffparams->ntypes);
994 for (i=0; i<ffparams->ntypes; i++) {
995 (void) pr_indent(fp,indent+INDENT);
996 (void) fprintf(fp,"functype[%d]=%s, ",
997 bShowNumbers?i:-1,
998 interaction_function[ffparams->functype[i]].name);
999 pr_iparams(fp,ffparams->functype[i],&ffparams->iparams[i]);
1001 (void) pr_double(fp,indent,"reppow",ffparams->reppow);
1002 (void) pr_real(fp,indent,"fudgeQQ",ffparams->fudgeQQ);
1003 pr_cmap(fp,indent,"cmap",&ffparams->cmap_grid,bShowNumbers);
1006 void pr_idef(FILE *fp,int indent,const char *title,t_idef *idef, bool bShowNumbers)
1008 int i,j;
1010 if (available(fp,idef,indent,title)) {
1011 indent=pr_title(fp,indent,title);
1012 (void) pr_indent(fp,indent);
1013 (void) fprintf(fp,"atnr=%d\n",idef->atnr);
1014 (void) pr_indent(fp,indent);
1015 (void) fprintf(fp,"ntypes=%d\n",idef->ntypes);
1016 for (i=0; i<idef->ntypes; i++) {
1017 (void) pr_indent(fp,indent+INDENT);
1018 (void) fprintf(fp,"functype[%d]=%s, ",
1019 bShowNumbers?i:-1,
1020 interaction_function[idef->functype[i]].name);
1021 pr_iparams(fp,idef->functype[i],&idef->iparams[i]);
1023 (void) pr_real(fp,indent,"fudgeQQ",idef->fudgeQQ);
1025 for(j=0; (j<F_NRE); j++)
1026 pr_ilist(fp,indent,interaction_function[j].longname,
1027 idef->functype,&idef->il[j],bShowNumbers);
1031 static int pr_block_title(FILE *fp,int indent,const char *title,t_block *block)
1033 int i;
1035 if (available(fp,block,indent,title))
1037 indent=pr_title(fp,indent,title);
1038 (void) pr_indent(fp,indent);
1039 (void) fprintf(fp,"nr=%d\n",block->nr);
1041 return indent;
1044 static int pr_blocka_title(FILE *fp,int indent,const char *title,t_blocka *block)
1046 int i;
1048 if (available(fp,block,indent,title))
1050 indent=pr_title(fp,indent,title);
1051 (void) pr_indent(fp,indent);
1052 (void) fprintf(fp,"nr=%d\n",block->nr);
1053 (void) pr_indent(fp,indent);
1054 (void) fprintf(fp,"nra=%d\n",block->nra);
1056 return indent;
1059 static void low_pr_block(FILE *fp,int indent,const char *title,t_block *block, bool bShowNumbers)
1061 int i;
1063 if (available(fp,block,indent,title))
1065 indent=pr_block_title(fp,indent,title,block);
1066 for (i=0; i<=block->nr; i++)
1068 (void) pr_indent(fp,indent+INDENT);
1069 (void) fprintf(fp,"%s->index[%d]=%u\n",
1070 title,bShowNumbers?i:-1,block->index[i]);
1075 static void low_pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block, bool bShowNumbers)
1077 int i;
1079 if (available(fp,block,indent,title))
1081 indent=pr_blocka_title(fp,indent,title,block);
1082 for (i=0; i<=block->nr; i++)
1084 (void) pr_indent(fp,indent+INDENT);
1085 (void) fprintf(fp,"%s->index[%d]=%u\n",
1086 title,bShowNumbers?i:-1,block->index[i]);
1088 for (i=0; i<block->nra; i++)
1090 (void) pr_indent(fp,indent+INDENT);
1091 (void) fprintf(fp,"%s->a[%d]=%u\n",
1092 title,bShowNumbers?i:-1,block->a[i]);
1097 void pr_block(FILE *fp,int indent,const char *title,t_block *block,bool bShowNumbers)
1099 int i,j,ok,size,start,end;
1101 if (available(fp,block,indent,title))
1103 indent=pr_block_title(fp,indent,title,block);
1104 start=0;
1105 end=start;
1106 if ((ok=(block->index[start]==0))==0)
1107 (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1108 else
1109 for (i=0; i<block->nr; i++)
1111 end=block->index[i+1];
1112 size=pr_indent(fp,indent);
1113 if (end<=start)
1114 size+=fprintf(fp,"%s[%d]={}\n",title,i);
1115 else
1116 size+=fprintf(fp,"%s[%d]={%d..%d}\n",
1117 title,bShowNumbers?i:-1,
1118 bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1119 start=end;
1124 void pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block,bool bShowNumbers)
1126 int i,j,ok,size,start,end;
1128 if (available(fp,block,indent,title))
1130 indent=pr_blocka_title(fp,indent,title,block);
1131 start=0;
1132 end=start;
1133 if ((ok=(block->index[start]==0))==0)
1134 (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1135 else
1136 for (i=0; i<block->nr; i++)
1138 end=block->index[i+1];
1139 size=pr_indent(fp,indent);
1140 if (end<=start)
1141 size+=fprintf(fp,"%s[%d]={",title,i);
1142 else
1143 size+=fprintf(fp,"%s[%d][%d..%d]={",
1144 title,bShowNumbers?i:-1,
1145 bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1146 for (j=start; j<end; j++)
1148 if (j>start) size+=fprintf(fp,", ");
1149 if ((size)>(USE_WIDTH))
1151 (void) fprintf(fp,"\n");
1152 size=pr_indent(fp,indent+INDENT);
1154 size+=fprintf(fp,"%u",block->a[j]);
1156 (void) fprintf(fp,"}\n");
1157 start=end;
1159 if ((end!=block->nra)||(!ok))
1161 (void) pr_indent(fp,indent);
1162 (void) fprintf(fp,"tables inconsistent, dumping complete tables:\n");
1163 low_pr_blocka(fp,indent,title,block,bShowNumbers);
1168 static void pr_strings(FILE *fp,int indent,const char *title,char ***nm,int n, bool bShowNumbers)
1170 int i;
1172 if (available(fp,nm,indent,title))
1174 indent=pr_title_n(fp,indent,title,n);
1175 for (i=0; i<n; i++)
1177 (void) pr_indent(fp,indent);
1178 (void) fprintf(fp,"%s[%d]={name=\"%s\"}\n",
1179 title,bShowNumbers?i:-1,*(nm[i]));
1184 static void pr_strings2(FILE *fp,int indent,const char *title,
1185 char ***nm,char ***nmB,int n, bool bShowNumbers)
1187 int i;
1189 if (available(fp,nm,indent,title))
1191 indent=pr_title_n(fp,indent,title,n);
1192 for (i=0; i<n; i++)
1194 (void) pr_indent(fp,indent);
1195 (void) fprintf(fp,"%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1196 title,bShowNumbers?i:-1,*(nm[i]),*(nmB[i]));
1201 static void pr_resinfo(FILE *fp,int indent,const char *title,t_resinfo *resinfo,int n, bool bShowNumbers)
1203 int i;
1205 if (available(fp,resinfo,indent,title))
1207 indent=pr_title_n(fp,indent,title,n);
1208 for (i=0; i<n; i++)
1210 (void) pr_indent(fp,indent);
1211 (void) fprintf(fp,"%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1212 title,bShowNumbers?i:-1,
1213 *(resinfo[i].name),resinfo[i].nr,
1214 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1219 static void pr_atom(FILE *fp,int indent,const char *title,t_atom *atom,int n)
1221 int i,j;
1223 if (available(fp,atom,indent,title)) {
1224 indent=pr_title_n(fp,indent,title,n);
1225 for (i=0; i<n; i++) {
1226 (void) pr_indent(fp,indent);
1227 fprintf(fp,"%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1228 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1229 title,i,atom[i].type,atom[i].typeB,ptype_str[atom[i].ptype],
1230 atom[i].m,atom[i].q,atom[i].mB,atom[i].qB,
1231 atom[i].resind,atom[i].atomnumber);
1236 static void pr_grps(FILE *fp,int indent,const char *title,t_grps grps[],int ngrp,
1237 char **grpname[], bool bShowNumbers)
1239 int i,j;
1241 for(i=0; (i<ngrp); i++) {
1242 fprintf(fp,"%s[%d] nr=%d, name=[",title,bShowNumbers?i:-1,grps[i].nr);
1243 for(j=0; (j<grps[i].nr); j++)
1244 fprintf(fp," %s",*(grpname[grps[i].nm_ind[j]]));
1245 fprintf(fp,"]\n");
1249 static void pr_groups(FILE *fp,int indent,const char *title,
1250 gmx_groups_t *groups,
1251 bool bShowNumbers)
1253 int grpnr[egcNR];
1254 int nat_max,i,g;
1256 pr_grps(fp,indent,"grp",groups->grps,egcNR,groups->grpname,bShowNumbers);
1257 pr_strings(fp,indent,"grpname",groups->grpname,groups->ngrpname,bShowNumbers);
1259 (void) pr_indent(fp,indent);
1260 fprintf(fp,"groups ");
1261 for(g=0; g<egcNR; g++)
1263 printf(" %5.5s",gtypes[g]);
1265 printf("\n");
1267 (void) pr_indent(fp,indent);
1268 fprintf(fp,"allocated ");
1269 nat_max = 0;
1270 for(g=0; g<egcNR; g++)
1272 printf(" %5d",groups->ngrpnr[g]);
1273 nat_max = max(nat_max,groups->ngrpnr[g]);
1275 printf("\n");
1277 if (nat_max == 0)
1279 (void) pr_indent(fp,indent);
1280 fprintf(fp,"groupnr[%5s] =","*");
1281 for(g=0; g<egcNR; g++)
1283 fprintf(fp," %3d ",0);
1285 fprintf(fp,"\n");
1287 else
1289 for(i=0; i<nat_max; i++)
1291 (void) pr_indent(fp,indent);
1292 fprintf(fp,"groupnr[%5d] =",i);
1293 for(g=0; g<egcNR; g++)
1295 fprintf(fp," %3d ",
1296 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1298 fprintf(fp,"\n");
1303 void pr_atoms(FILE *fp,int indent,const char *title,t_atoms *atoms,
1304 bool bShownumbers)
1306 if (available(fp,atoms,indent,title))
1308 indent=pr_title(fp,indent,title);
1309 pr_atom(fp,indent,"atom",atoms->atom,atoms->nr);
1310 pr_strings(fp,indent,"atom",atoms->atomname,atoms->nr,bShownumbers);
1311 pr_strings2(fp,indent,"type",atoms->atomtype,atoms->atomtypeB,atoms->nr,bShownumbers);
1312 pr_resinfo(fp,indent,"residue",atoms->resinfo,atoms->nres,bShownumbers);
1317 void pr_atomtypes(FILE *fp,int indent,const char *title,t_atomtypes *atomtypes,
1318 bool bShowNumbers)
1320 int i;
1321 if (available(fp,atomtypes,indent,title))
1323 indent=pr_title(fp,indent,title);
1324 for(i=0;i<atomtypes->nr;i++) {
1325 pr_indent(fp,indent);
1326 fprintf(fp,
1327 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1328 bShowNumbers?i:-1,atomtypes->radius[i],atomtypes->vol[i],
1329 atomtypes->gb_radius[i],
1330 atomtypes->surftens[i],atomtypes->atomnumber[i],atomtypes->S_hct[i]);
1335 static void pr_moltype(FILE *fp,int indent,const char *title,
1336 gmx_moltype_t *molt,int n,
1337 gmx_ffparams_t *ffparams,
1338 bool bShowNumbers)
1340 int j;
1342 indent = pr_title_n(fp,indent,title,n);
1343 (void) pr_indent(fp,indent);
1344 (void) fprintf(fp,"name=\"%s\"\n",*(molt->name));
1345 pr_atoms(fp,indent,"atoms",&(molt->atoms),bShowNumbers);
1346 pr_block(fp,indent,"cgs",&molt->cgs, bShowNumbers);
1347 pr_blocka(fp,indent,"excls",&molt->excls, bShowNumbers);
1348 for(j=0; (j<F_NRE); j++) {
1349 pr_ilist(fp,indent,interaction_function[j].longname,
1350 ffparams->functype,&molt->ilist[j],bShowNumbers);
1354 static void pr_molblock(FILE *fp,int indent,const char *title,
1355 gmx_molblock_t *molb,int n,
1356 gmx_moltype_t *molt,
1357 bool bShowNumbers)
1359 indent = pr_title_n(fp,indent,title,n);
1360 (void) pr_indent(fp,indent);
1361 (void) fprintf(fp,"%-20s = %d \"%s\"\n",
1362 "moltype",molb->type,*(molt[molb->type].name));
1363 pr_int(fp,indent,"#molecules",molb->nmol);
1364 pr_int(fp,indent,"#atoms_mol",molb->natoms_mol);
1365 pr_int(fp,indent,"#posres_xA",molb->nposres_xA);
1366 if (molb->nposres_xA > 0) {
1367 pr_rvecs(fp,indent,"posres_xA",molb->posres_xA,molb->nposres_xA);
1369 pr_int(fp,indent,"#posres_xB",molb->nposres_xB);
1370 if (molb->nposres_xB > 0) {
1371 pr_rvecs(fp,indent,"posres_xB",molb->posres_xB,molb->nposres_xB);
1375 void pr_mtop(FILE *fp,int indent,const char *title,gmx_mtop_t *mtop,
1376 bool bShowNumbers)
1378 int mt,mb;
1380 if (available(fp,mtop,indent,title)) {
1381 indent=pr_title(fp,indent,title);
1382 (void) pr_indent(fp,indent);
1383 (void) fprintf(fp,"name=\"%s\"\n",*(mtop->name));
1384 pr_int(fp,indent,"#atoms",mtop->natoms);
1385 for(mb=0; mb<mtop->nmolblock; mb++) {
1386 pr_molblock(fp,indent,"molblock",&mtop->molblock[mb],mb,
1387 mtop->moltype,bShowNumbers);
1389 pr_ffparams(fp,indent,"ffparams",&(mtop->ffparams),bShowNumbers);
1390 pr_atomtypes(fp,indent,"atomtypes",&(mtop->atomtypes),bShowNumbers);
1391 for(mt=0; mt<mtop->nmoltype; mt++) {
1392 pr_moltype(fp,indent,"moltype",&mtop->moltype[mt],mt,
1393 &mtop->ffparams,bShowNumbers);
1395 pr_groups(fp,indent,"groups",&mtop->groups,bShowNumbers);
1399 void pr_top(FILE *fp,int indent,const char *title,t_topology *top, bool bShowNumbers)
1401 if (available(fp,top,indent,title)) {
1402 indent=pr_title(fp,indent,title);
1403 (void) pr_indent(fp,indent);
1404 (void) fprintf(fp,"name=\"%s\"\n",*(top->name));
1405 pr_atoms(fp,indent,"atoms",&(top->atoms),bShowNumbers);
1406 pr_atomtypes(fp,indent,"atomtypes",&(top->atomtypes),bShowNumbers);
1407 pr_block(fp,indent,"cgs",&top->cgs, bShowNumbers);
1408 pr_block(fp,indent,"mols",&top->mols, bShowNumbers);
1409 pr_blocka(fp,indent,"excls",&top->excls, bShowNumbers);
1410 pr_idef(fp,indent,"idef",&top->idef,bShowNumbers);
1414 void pr_header(FILE *fp,int indent,const char *title,t_tpxheader *sh)
1416 char buf[22];
1418 if (available(fp,sh,indent,title))
1420 indent=pr_title(fp,indent,title);
1421 pr_indent(fp,indent);
1422 fprintf(fp,"bIr = %spresent\n",sh->bIr?"":"not ");
1423 pr_indent(fp,indent);
1424 fprintf(fp,"bBox = %spresent\n",sh->bBox?"":"not ");
1425 pr_indent(fp,indent);
1426 fprintf(fp,"bTop = %spresent\n",sh->bTop?"":"not ");
1427 pr_indent(fp,indent);
1428 fprintf(fp,"bX = %spresent\n",sh->bX?"":"not ");
1429 pr_indent(fp,indent);
1430 fprintf(fp,"bV = %spresent\n",sh->bV?"":"not ");
1431 pr_indent(fp,indent);
1432 fprintf(fp,"bF = %spresent\n",sh->bF?"":"not ");
1434 pr_indent(fp,indent);
1435 fprintf(fp,"natoms = %d\n",sh->natoms);
1436 pr_indent(fp,indent);
1437 fprintf(fp,"lambda = %e\n",sh->lambda);
1441 void pr_commrec(FILE *fp,int indent,t_commrec *cr)
1443 pr_indent(fp,indent);
1444 fprintf(fp,"commrec:\n");
1445 indent+=2;
1446 pr_indent(fp,indent);
1447 fprintf(fp,"nodeid = %d\n",cr->nodeid);
1448 pr_indent(fp,indent);
1449 fprintf(fp,"nnodes = %d\n",cr->nnodes);
1450 pr_indent(fp,indent);
1451 fprintf(fp,"npmenodes = %d\n",cr->npmenodes);
1452 pr_indent(fp,indent);
1453 fprintf(fp,"threadid = %d\n",cr->threadid);
1454 pr_indent(fp,indent);
1455 fprintf(fp,"nthreads = %d\n",cr->nthreads);