1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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
34 * GROningen Mixture of Alchemy and Childrens' Stories
40 /* This file is completely threadsafe - please keep it that way! */
42 #include <thread_mpi.h>
56 int pr_indent(FILE *fp,int n)
60 for (i=0; i<n; i++) (void) fprintf(fp," ");
64 int available(FILE *fp,void *p,int indent,const char *title)
69 (void) fprintf(fp,"%s: not available\n",title);
74 int pr_title(FILE *fp,int indent,const char *title)
76 (void) pr_indent(fp,indent);
77 (void) fprintf(fp,"%s:\n",title);
78 return (indent+INDENT);
81 int pr_title_n(FILE *fp,int indent,const char *title,int n)
83 (void) pr_indent(fp,indent);
84 (void) fprintf(fp,"%s (%d):\n",title,n);
85 return (indent+INDENT);
88 int pr_title_nxn(FILE *fp,int indent,const char *title,int n1,int n2)
90 (void) pr_indent(fp,indent);
91 (void) fprintf(fp,"%s (%dx%d):\n",title,n1,n2);
92 return (indent+INDENT);
95 void pr_ivec(FILE *fp,int indent,const char *title,int vec[],int n, gmx_bool bShowNumbers)
99 if (available(fp,vec,indent,title))
101 indent=pr_title_n(fp,indent,title,n);
104 (void) pr_indent(fp,indent);
105 (void) fprintf(fp,"%s[%d]=%d\n",title,bShowNumbers?i:-1,vec[i]);
110 void pr_ivec_block(FILE *fp,int indent,const char *title,int vec[],int n, gmx_bool bShowNumbers)
114 if (available(fp,vec,indent,title))
116 indent=pr_title_n(fp,indent,title,n);
121 while (j < n && vec[j] == vec[j-1]+1)
125 /* Print consecutive groups of 3 or more as blocks */
130 (void) pr_indent(fp,indent);
131 (void) fprintf(fp,"%s[%d]=%d\n",
132 title,bShowNumbers?i:-1,vec[i]);
138 (void) pr_indent(fp,indent);
139 (void) fprintf(fp,"%s[%d,...,%d] = {%d,...,%d}\n",
150 void pr_bvec(FILE *fp,int indent,const char *title,gmx_bool vec[],int n, gmx_bool bShowNumbers)
154 if (available(fp,vec,indent,title))
156 indent=pr_title_n(fp,indent,title,n);
159 (void) pr_indent(fp,indent);
160 (void) fprintf(fp,"%s[%d]=%s\n",title,bShowNumbers?i:-1,
166 void pr_ivecs(FILE *fp,int indent,const char *title,ivec vec[],int n, gmx_bool bShowNumbers)
170 if (available(fp,vec,indent,title))
172 indent=pr_title_nxn(fp,indent,title,n,DIM);
175 (void) pr_indent(fp,indent);
176 (void) fprintf(fp,"%s[%d]={",title,bShowNumbers?i:-1);
177 for (j=0; j<DIM; j++)
179 if (j!=0) (void) fprintf(fp,", ");
180 fprintf(fp,"%d",vec[i][j]);
182 (void) fprintf(fp,"}\n");
187 void pr_rvec(FILE *fp,int indent,const char *title,real vec[],int n, gmx_bool bShowNumbers)
191 if (available(fp,vec,indent,title))
193 indent=pr_title_n(fp,indent,title,n);
196 pr_indent(fp,indent);
197 fprintf(fp,"%s[%d]=%12.5e\n",title,bShowNumbers?i:-1,vec[i]);
202 void pr_dvec(FILE *fp,int indent,const char *title,double vec[],int n, gmx_bool bShowNumbers)
206 if (available(fp,vec,indent,title))
208 indent=pr_title_n(fp,indent,title,n);
211 pr_indent(fp,indent);
212 fprintf(fp,"%s[%d]=%12.5e\n",title,bShowNumbers?i:-1,vec[i]);
219 void pr_mat(FILE *fp,int indent,char *title,matrix m)
223 if (available(fp,m,indent,title)) {
224 indent=pr_title_n(fp,indent,title,n);
226 pr_indent(fp,indent);
227 fprintf(fp,"%s[%d]=%12.5e %12.5e %12.5e\n",
228 title,bShowNumbers?i:-1,m[i][XX],m[i][YY],m[i][ZZ]);
234 void pr_rvecs_len(FILE *fp,int indent,const char *title,rvec vec[],int n)
238 if (available(fp,vec,indent,title)) {
239 indent=pr_title_nxn(fp,indent,title,n,DIM);
240 for (i=0; i<n; i++) {
241 (void) pr_indent(fp,indent);
242 (void) fprintf(fp,"%s[%5d]={",title,i);
243 for (j=0; j<DIM; j++) {
245 (void) fprintf(fp,", ");
246 (void) fprintf(fp,"%12.5e",vec[i][j]);
248 (void) fprintf(fp,"} len=%12.5e\n",norm(vec[i]));
253 void pr_rvecs(FILE *fp,int indent,const char *title,rvec vec[],int n)
255 const char *fshort = "%12.5e";
256 const char *flong = "%15.8e";
260 if (getenv("LONGFORMAT") != NULL)
265 if (available(fp,vec,indent,title)) {
266 indent=pr_title_nxn(fp,indent,title,n,DIM);
267 for (i=0; i<n; i++) {
268 (void) pr_indent(fp,indent);
269 (void) fprintf(fp,"%s[%5d]={",title,i);
270 for (j=0; j<DIM; j++) {
272 (void) fprintf(fp,", ");
273 (void) fprintf(fp,format,vec[i][j]);
275 (void) fprintf(fp,"}\n");
281 void pr_reals(FILE *fp,int indent,const char *title,real *vec,int n)
285 if (available(fp,vec,indent,title)) {
286 (void) pr_indent(fp,indent);
287 (void) fprintf(fp,"%s:\t",title);
289 fprintf(fp," %10g",vec[i]);
290 (void) fprintf(fp,"\n");
294 void pr_doubles(FILE *fp,int indent,const char *title,double *vec,int n)
298 if (available(fp,vec,indent,title)) {
299 (void) pr_indent(fp,indent);
300 (void) fprintf(fp,"%s:\t",title);
302 fprintf(fp," %10g",vec[i]);
303 (void) fprintf(fp,"\n");
307 static void pr_int(FILE *fp,int indent,const char *title,int i)
309 pr_indent(fp,indent);
310 fprintf(fp,"%-20s = %d\n",title,i);
313 static void pr_gmx_large_int(FILE *fp,int indent,const char *title,gmx_large_int_t i)
315 char buf[STEPSTRSIZE];
317 pr_indent(fp,indent);
318 fprintf(fp,"%-20s = %s\n",title,gmx_step_str(i,buf));
321 static void pr_real(FILE *fp,int indent,const char *title,real r)
323 pr_indent(fp,indent);
324 fprintf(fp,"%-20s = %g\n",title,r);
327 static void pr_double(FILE *fp,int indent,const char *title,double d)
329 pr_indent(fp,indent);
330 fprintf(fp,"%-20s = %g\n",title,d);
333 static void pr_str(FILE *fp,int indent,const char *title,const char *s)
335 pr_indent(fp,indent);
336 fprintf(fp,"%-20s = %s\n",title,s);
339 void pr_qm_opts(FILE *fp,int indent,const char *title,t_grpopts *opts)
343 fprintf(fp,"%s:\n",title);
345 pr_int(fp,indent,"ngQM",opts->ngQM);
346 if (opts->ngQM > 0) {
347 pr_ivec(fp,indent,"QMmethod",opts->QMmethod,opts->ngQM,FALSE);
348 pr_ivec(fp,indent,"QMbasis",opts->QMbasis,opts->ngQM,FALSE);
349 pr_ivec(fp,indent,"QMcharge",opts->QMcharge,opts->ngQM,FALSE);
350 pr_ivec(fp,indent,"QMmult",opts->QMmult,opts->ngQM,FALSE);
351 pr_bvec(fp,indent,"bSH",opts->bSH,opts->ngQM,FALSE);
352 pr_ivec(fp,indent,"CASorbitals",opts->CASorbitals,opts->ngQM,FALSE);
353 pr_ivec(fp,indent,"CASelectrons",opts->CASelectrons,opts->ngQM,FALSE);
354 pr_rvec(fp,indent,"SAon",opts->SAon,opts->ngQM,FALSE);
355 pr_rvec(fp,indent,"SAon",opts->SAon,opts->ngQM,FALSE);
356 pr_ivec(fp,indent,"SAsteps",opts->SAsteps,opts->ngQM,FALSE);
357 pr_bvec(fp,indent,"bOPT",opts->bOPT,opts->ngQM,FALSE);
358 pr_bvec(fp,indent,"bTS",opts->bTS,opts->ngQM,FALSE);
362 static void pr_grp_opts(FILE *out,int indent,const char *title,t_grpopts *opts,
368 fprintf(out,"%s:\n",title);
370 pr_indent(out,indent);
371 fprintf(out,"nrdf%s",bMDPformat ? " = " : ":");
372 for(i=0; (i<opts->ngtc); i++)
373 fprintf(out," %10g",opts->nrdf[i]);
376 pr_indent(out,indent);
377 fprintf(out,"ref-t%s",bMDPformat ? " = " : ":");
378 for(i=0; (i<opts->ngtc); i++)
379 fprintf(out," %10g",opts->ref_t[i]);
382 pr_indent(out,indent);
383 fprintf(out,"tau-t%s",bMDPformat ? " = " : ":");
384 for(i=0; (i<opts->ngtc); i++)
385 fprintf(out," %10g",opts->tau_t[i]);
388 /* Pretty-print the simulated annealing info */
389 fprintf(out,"anneal%s",bMDPformat ? " = " : ":");
390 for(i=0; (i<opts->ngtc); i++)
391 fprintf(out," %10s",EANNEAL(opts->annealing[i]));
394 fprintf(out,"ann-npoints%s",bMDPformat ? " = " : ":");
395 for(i=0; (i<opts->ngtc); i++)
396 fprintf(out," %10d",opts->anneal_npoints[i]);
399 for(i=0; (i<opts->ngtc); i++) {
400 if(opts->anneal_npoints[i]>0) {
401 fprintf(out,"ann. times [%d]:\t",i);
402 for(j=0; (j<opts->anneal_npoints[i]); j++)
403 fprintf(out," %10.1f",opts->anneal_time[i][j]);
405 fprintf(out,"ann. temps [%d]:\t",i);
406 for(j=0; (j<opts->anneal_npoints[i]); j++)
407 fprintf(out," %10.1f",opts->anneal_temp[i][j]);
412 pr_indent(out,indent);
413 fprintf(out,"acc:\t");
414 for(i=0; (i<opts->ngacc); i++)
415 for(m=0; (m<DIM); m++)
416 fprintf(out," %10g",opts->acc[i][m]);
419 pr_indent(out,indent);
420 fprintf(out,"nfreeze:");
421 for(i=0; (i<opts->ngfrz); i++)
422 for(m=0; (m<DIM); m++)
423 fprintf(out," %10s",opts->nFreeze[i][m] ? "Y" : "N");
427 for(i=0; (i<opts->ngener); i++) {
428 pr_indent(out,indent);
429 fprintf(out,"energygrp-flags[%3d]:",i);
430 for(m=0; (m<opts->ngener); m++)
431 fprintf(out," %d",opts->egp_flags[opts->ngener*i+m]);
438 static void pr_matrix(FILE *fp,int indent,const char *title,rvec *m,
442 fprintf(fp,"%-10s = %g %g %g %g %g %g\n",title,
443 m[XX][XX],m[YY][YY],m[ZZ][ZZ],m[XX][YY],m[XX][ZZ],m[YY][ZZ]);
445 pr_rvecs(fp,indent,title,m,DIM);
448 static void pr_cosine(FILE *fp,int indent,const char *title,t_cosines *cos,
454 fprintf(fp,"%s = %d\n",title,cos->n);
457 indent=pr_title(fp,indent,title);
458 (void) pr_indent(fp,indent);
459 fprintf(fp,"n = %d\n",cos->n);
461 (void) pr_indent(fp,indent+2);
463 for(j=0; (j<cos->n); j++)
464 fprintf(fp," %e",cos->a[j]);
466 (void) pr_indent(fp,indent+2);
468 for(j=0; (j<cos->n); j++)
469 fprintf(fp," %e",cos->phi[j]);
475 #define PS(t,s) pr_str(fp,indent,t,s)
476 #define PI(t,s) pr_int(fp,indent,t,s)
477 #define PSTEP(t,s) pr_gmx_large_int(fp,indent,t,s)
478 #define PR(t,s) pr_real(fp,indent,t,s)
479 #define PD(t,s) pr_double(fp,indent,t,s)
481 static void pr_pullgrp(FILE *fp,int indent,int g,t_pullgrp *pg)
483 pr_indent(fp,indent);
484 fprintf(fp,"pull-group %d:\n",g);
486 pr_ivec_block(fp,indent,"atom",pg->ind,pg->nat,TRUE);
487 pr_rvec(fp,indent,"weight",pg->weight,pg->nweight,TRUE);
488 PI("pbcatom",pg->pbcatom);
489 pr_rvec(fp,indent,"vec",pg->vec,DIM,TRUE);
490 pr_rvec(fp,indent,"init",pg->init,DIM,TRUE);
496 static void pr_simtempvals(FILE *fp,int indent,t_simtemp *simtemp, int n_lambda, gmx_bool bMDPformat)
498 PR("simtemp_low",simtemp->simtemp_low);
499 PR("simtemp_high",simtemp->simtemp_high);
500 PS("simulated-tempering-scaling",ESIMTEMP(simtemp->eSimTempScale));
501 pr_rvec(fp,indent,"simulated tempering temperatures",simtemp->temperatures,n_lambda,TRUE);
504 static void pr_expandedvals(FILE *fp,int indent,t_expanded *expand, int n_lambda, gmx_bool bMDPformat)
507 PI("nstexpanded", expand->nstexpanded);
508 PS("lambda-stats", elamstats_names[expand->elamstats]);
509 PS("lambda-mc-move", elmcmove_names[expand->elmcmove]);
510 PI("lmc-repeats",expand->lmc_repeats);
511 PI("lmc-gibbsdelta",expand->gibbsdeltalam);
512 PI("lmc-nstart",expand->lmc_forced_nstart);
513 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
514 PI("nst-transition-matrix",expand->nstTij);
515 PI("mininum-var-min",expand->minvarmin); /*default is reasonable */
516 PI("weight-c-range",expand->c_range); /* default is just C=0 */
517 PR("wl-scale",expand->wl_scale);
518 PR("init-wl-delta",expand->init_wl_delta);
519 PR("wl-ratio",expand->wl_ratio);
520 PS("bWLoneovert",EBOOL(expand->bWLoneovert));
521 PI("lmc-seed",expand->lmc_seed);
522 PR("mc-temperature",expand->mc_temp);
523 PS("lmc-weights-equil",elmceq_names[expand->elmceq]);
524 if (expand->elmceq == elmceqNUMATLAM)
526 PI("weight-equil-number-all-lambda",expand->equil_n_at_lam);
528 if (expand->elmceq == elmceqSAMPLES)
530 PI("weight-equil-number-samples",expand->equil_samples);
532 if (expand->elmceq == elmceqSTEPS)
534 PI("weight-equil-number-steps",expand->equil_steps);
536 if (expand->elmceq == elmceqWLDELTA)
538 PR("weight-equil-wl-delta",expand->equil_wl_delta);
540 if (expand->elmceq == elmceqRATIO)
542 PR("weight-equil-count-ratio",expand->equil_ratio);
545 pr_indent(fp,indent);
546 pr_rvec(fp,indent,"init-lambda-weights",expand->init_lambda_weights,n_lambda,TRUE);
547 PS("init-weights",EBOOL(expand->bInit_weights));
550 static void pr_fepvals(FILE *fp,int indent,t_lambda *fep, gmx_bool bMDPformat)
554 PI("nstdhdl",fep->nstdhdl);
555 PI("init-lambda-state",fep->init_fep_state);
556 PR("init-lambda",fep->init_lambda);
557 PR("delta-lambda",fep->delta_lambda);
560 PI("n-lambdas",fep->n_lambda);
562 if (fep->n_lambda > 0)
564 pr_indent(fp,indent);
565 fprintf(fp,"all-lambdas%s\n",bMDPformat ? " = " : ":");
566 for(i=0; i<efptNR; i++) {
567 fprintf(fp,"%18s = ",efpt_names[i]);
568 for(j=0; j<fep->n_lambda; j++)
570 fprintf(fp," %10g",fep->all_lambda[i][j]);
576 PR("sc-alpha",fep->sc_alpha);
577 PS("bScCoul",EBOOL(fep->bScCoul));
578 PS("bScPrintEnergy",EBOOL(fep->bPrintEnergy));
579 PI("sc-power",fep->sc_power);
580 PR("sc-r-power",fep->sc_r_power);
581 PR("sc-sigma",fep->sc_sigma);
582 PR("sc-sigma-min",fep->sc_sigma_min);
583 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
584 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
585 PI("dh-hist-size", fep->dh_hist_size);
586 PD("dh-hist-spacing", fep->dh_hist_spacing);
589 static void pr_pull(FILE *fp,int indent,t_pull *pull)
593 PS("pull-geometry",EPULLGEOM(pull->eGeom));
594 pr_ivec(fp,indent,"pull-dim",pull->dim,DIM,TRUE);
595 PR("pull-r1",pull->cyl_r1);
596 PR("pull-r0",pull->cyl_r0);
597 PR("pull-constr-tol",pull->constr_tol);
598 PI("pull-nstxout",pull->nstxout);
599 PI("pull-nstfout",pull->nstfout);
600 PI("pull-ngrp",pull->ngrp);
601 for(g=0; g<pull->ngrp+1; g++)
602 pr_pullgrp(fp,indent,g,&pull->grp[g]);
605 static void pr_rotgrp(FILE *fp,int indent,int g,t_rotgrp *rotg)
607 pr_indent(fp,indent);
608 fprintf(fp,"rotation_group %d:\n",g);
610 PS("type",EROTGEOM(rotg->eType));
611 PS("massw",EBOOL(rotg->bMassW));
612 pr_ivec_block(fp,indent,"atom",rotg->ind,rotg->nat,TRUE);
613 pr_rvecs(fp,indent,"x_ref",rotg->x_ref,rotg->nat);
614 pr_rvec(fp,indent,"vec",rotg->vec,DIM,TRUE);
615 pr_rvec(fp,indent,"pivot",rotg->pivot,DIM,TRUE);
616 PR("rate",rotg->rate);
618 PR("slab_dist",rotg->slab_dist);
619 PR("min_gaussian",rotg->min_gaussian);
620 PR("epsilon",rotg->eps);
621 PS("fit_method",EROTFIT(rotg->eFittype));
622 PI("potfitangle_nstep",rotg->PotAngle_nstep);
623 PR("potfitangle_step",rotg->PotAngle_step);
626 static void pr_rot(FILE *fp,int indent,t_rot *rot)
630 PI("rot_nstrout",rot->nstrout);
631 PI("rot_nstsout",rot->nstsout);
632 PI("rot_ngrp",rot->ngrp);
633 for(g=0; g<rot->ngrp; g++)
634 pr_rotgrp(fp,indent,g,&rot->grp[g]);
637 void pr_inputrec(FILE *fp,int indent,const char *title,t_inputrec *ir,
640 const char *infbuf="inf";
643 if (available(fp,ir,indent,title)) {
645 indent=pr_title(fp,indent,title);
646 PS("integrator",EI(ir->eI));
647 PSTEP("nsteps",ir->nsteps);
648 PSTEP("init-step",ir->init_step);
649 PS("cutoff-scheme",ECUTSCHEME(ir->cutoff_scheme));
650 PS("ns_type",ENS(ir->ns_type));
651 PI("nstlist",ir->nstlist);
652 PI("ndelta",ir->ndelta);
653 PI("nstcomm",ir->nstcomm);
654 PS("comm-mode",ECOM(ir->comm_mode));
655 PI("nstlog",ir->nstlog);
656 PI("nstxout",ir->nstxout);
657 PI("nstvout",ir->nstvout);
658 PI("nstfout",ir->nstfout);
659 PI("nstcalcenergy",ir->nstcalcenergy);
660 PI("nstenergy",ir->nstenergy);
661 PI("nstxtcout",ir->nstxtcout);
662 PR("init-t",ir->init_t);
663 PR("delta-t",ir->delta_t);
665 PR("xtcprec",ir->xtcprec);
666 PR("fourierspacing",ir->fourier_spacing);
670 PI("pme-order",ir->pme_order);
671 PR("ewald-rtol",ir->ewald_rtol);
672 PR("ewald-geometry",ir->ewald_geometry);
673 PR("epsilon-surface",ir->epsilon_surface);
674 PS("optimize-fft",EBOOL(ir->bOptFFT));
675 PS("ePBC",EPBC(ir->ePBC));
676 PS("bPeriodicMols",EBOOL(ir->bPeriodicMols));
677 PS("bContinuation",EBOOL(ir->bContinuation));
678 PS("bShakeSOR",EBOOL(ir->bShakeSOR));
679 PS("etc",ETCOUPLTYPE(ir->etc));
680 PS("bPrintNHChains",EBOOL(ir->bPrintNHChains));
681 PI("nsttcouple",ir->nsttcouple);
682 PS("epc",EPCOUPLTYPE(ir->epc));
683 PS("epctype",EPCOUPLTYPETYPE(ir->epct));
684 PI("nstpcouple",ir->nstpcouple);
685 PR("tau-p",ir->tau_p);
686 pr_matrix(fp,indent,"ref-p",ir->ref_p,bMDPformat);
687 pr_matrix(fp,indent,"compress",ir->compress,bMDPformat);
688 PS("refcoord-scaling",EREFSCALINGTYPE(ir->refcoord_scaling));
690 fprintf(fp,"posres-com = %g %g %g\n",ir->posres_com[XX],
691 ir->posres_com[YY],ir->posres_com[ZZ]);
693 pr_rvec(fp,indent,"posres-com",ir->posres_com,DIM,TRUE);
695 fprintf(fp,"posres-comB = %g %g %g\n",ir->posres_comB[XX],
696 ir->posres_comB[YY],ir->posres_comB[ZZ]);
698 pr_rvec(fp,indent,"posres-comB",ir->posres_comB,DIM,TRUE);
699 PR("verlet-buffer-drift",ir->verletbuf_drift);
700 PR("rlist",ir->rlist);
701 PR("rlistlong",ir->rlistlong);
702 PR("nstcalclr",ir->nstcalclr);
704 PS("coulombtype",EELTYPE(ir->coulombtype));
705 PS("coulomb-modifier",INTMODIFIER(ir->coulomb_modifier));
706 PR("rcoulomb-switch",ir->rcoulomb_switch);
707 PR("rcoulomb",ir->rcoulomb);
708 PS("vdwtype",EVDWTYPE(ir->vdwtype));
709 PS("vdw-modifier",INTMODIFIER(ir->vdw_modifier));
710 PR("rvdw-switch",ir->rvdw_switch);
712 if (ir->epsilon_r != 0)
713 PR("epsilon-r",ir->epsilon_r);
715 PS("epsilon-r",infbuf);
716 if (ir->epsilon_rf != 0)
717 PR("epsilon-rf",ir->epsilon_rf);
719 PS("epsilon-rf",infbuf);
720 PR("tabext",ir->tabext);
721 PS("implicit-solvent",EIMPLICITSOL(ir->implicit_solvent));
722 PS("gb-algorithm",EGBALGORITHM(ir->gb_algorithm));
723 PR("gb-epsilon-solvent",ir->gb_epsilon_solvent);
724 PI("nstgbradii",ir->nstgbradii);
725 PR("rgbradii",ir->rgbradii);
726 PR("gb-saltconc",ir->gb_saltconc);
727 PR("gb-obc-alpha",ir->gb_obc_alpha);
728 PR("gb-obc-beta",ir->gb_obc_beta);
729 PR("gb-obc-gamma",ir->gb_obc_gamma);
730 PR("gb-dielectric-offset",ir->gb_dielectric_offset);
731 PS("sa-algorithm",ESAALGORITHM(ir->gb_algorithm));
732 PR("sa-surface-tension",ir->sa_surface_tension);
733 PS("DispCorr",EDISPCORR(ir->eDispCorr));
734 PS("bSimTemp",EBOOL(ir->bSimTemp));
736 pr_simtempvals(fp,indent,ir->simtempvals,ir->fepvals->n_lambda,bMDPformat);
738 PS("free-energy",EFEPTYPE(ir->efep));
739 if (ir->efep != efepNO || ir->bSimTemp) {
740 pr_fepvals(fp,indent,ir->fepvals,bMDPformat);
743 pr_expandedvals(fp,indent,ir->expandedvals,ir->fepvals->n_lambda,bMDPformat);
746 PI("nwall",ir->nwall);
747 PS("wall-type",EWALLTYPE(ir->wall_type));
748 PI("wall-atomtype[0]",ir->wall_atomtype[0]);
749 PI("wall-atomtype[1]",ir->wall_atomtype[1]);
750 PR("wall-density[0]",ir->wall_density[0]);
751 PR("wall-density[1]",ir->wall_density[1]);
752 PR("wall-ewald-zfac",ir->wall_ewald_zfac);
754 PS("pull",EPULLTYPE(ir->ePull));
755 if (ir->ePull != epullNO)
756 pr_pull(fp,indent,ir->pull);
758 PS("rotation",EBOOL(ir->bRot));
760 pr_rot(fp,indent,ir->rot);
762 PS("disre",EDISRETYPE(ir->eDisre));
763 PS("disre-weighting",EDISREWEIGHTING(ir->eDisreWeighting));
764 PS("disre-mixed",EBOOL(ir->bDisreMixed));
765 PR("dr-fc",ir->dr_fc);
766 PR("dr-tau",ir->dr_tau);
767 PR("nstdisreout",ir->nstdisreout);
768 PR("orires-fc",ir->orires_fc);
769 PR("orires-tau",ir->orires_tau);
770 PR("nstorireout",ir->nstorireout);
772 PR("dihre-fc",ir->dihre_fc);
774 PR("em-stepsize",ir->em_stepsize);
775 PR("em-tol",ir->em_tol);
776 PI("niter",ir->niter);
777 PR("fc-stepsize",ir->fc_stepsize);
778 PI("nstcgsteep",ir->nstcgsteep);
779 PI("nbfgscorr",ir->nbfgscorr);
781 PS("ConstAlg",ECONSTRTYPE(ir->eConstrAlg));
782 PR("shake-tol",ir->shake_tol);
783 PI("lincs-order",ir->nProjOrder);
784 PR("lincs-warnangle",ir->LincsWarnAngle);
785 PI("lincs-iter",ir->nLincsIter);
786 PR("bd-fric",ir->bd_fric);
787 PI("ld-seed",ir->ld_seed);
788 PR("cos-accel",ir->cos_accel);
789 pr_matrix(fp,indent,"deform",ir->deform,bMDPformat);
791 PS("adress",EBOOL(ir->bAdress));
793 PS("adress_type",EADRESSTYPE(ir->adress->type));
794 PR("adress_const_wf",ir->adress->const_wf);
795 PR("adress_ex_width",ir->adress->ex_width);
796 PR("adress_hy_width",ir->adress->hy_width);
797 PS("adress_interface_correction",EADRESSICTYPE(ir->adress->icor));
798 PS("adress_site",EADRESSSITETYPE(ir->adress->site));
799 PR("adress_ex_force_cap",ir->adress->ex_forcecap);
800 PS("adress_do_hybridpairs", EBOOL(ir->adress->do_hybridpairs));
802 pr_rvec(fp,indent,"adress_reference_coords",ir->adress->refs,DIM,TRUE);
804 PI("userint1",ir->userint1);
805 PI("userint2",ir->userint2);
806 PI("userint3",ir->userint3);
807 PI("userint4",ir->userint4);
808 PR("userreal1",ir->userreal1);
809 PR("userreal2",ir->userreal2);
810 PR("userreal3",ir->userreal3);
811 PR("userreal4",ir->userreal4);
812 pr_grp_opts(fp,indent,"grpopts",&(ir->opts),bMDPformat);
813 pr_cosine(fp,indent,"efield-x",&(ir->ex[XX]),bMDPformat);
814 pr_cosine(fp,indent,"efield-xt",&(ir->et[XX]),bMDPformat);
815 pr_cosine(fp,indent,"efield-y",&(ir->ex[YY]),bMDPformat);
816 pr_cosine(fp,indent,"efield-yt",&(ir->et[YY]),bMDPformat);
817 pr_cosine(fp,indent,"efield-z",&(ir->ex[ZZ]),bMDPformat);
818 pr_cosine(fp,indent,"efield-zt",&(ir->et[ZZ]),bMDPformat);
819 PS("bQMMM",EBOOL(ir->bQMMM));
820 PI("QMconstraints",ir->QMconstraints);
821 PI("QMMMscheme",ir->QMMMscheme);
822 PR("scalefactor",ir->scalefactor);
823 pr_qm_opts(fp,indent,"qm-opts",&(ir->opts));
830 static void pr_harm(FILE *fp,t_iparams *iparams,const char *r,const char *kr)
832 fprintf(fp,"%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
833 r,iparams->harmonic.rA,kr,iparams->harmonic.krA,
834 r,iparams->harmonic.rB,kr,iparams->harmonic.krB);
837 void pr_iparams(FILE *fp,t_functype ftype,t_iparams *iparams)
840 real VA[4],VB[4],*rbcA,*rbcB;
845 pr_harm(fp,iparams,"th","ct");
847 case F_CROSS_BOND_BONDS:
848 fprintf(fp,"r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
849 iparams->cross_bb.r1e,iparams->cross_bb.r2e,
850 iparams->cross_bb.krr);
852 case F_CROSS_BOND_ANGLES:
853 fprintf(fp,"r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
854 iparams->cross_ba.r1e,iparams->cross_ba.r2e,
855 iparams->cross_ba.r3e,iparams->cross_ba.krt);
857 case F_LINEAR_ANGLES:
858 fprintf(fp,"klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
859 iparams->linangle.klinA,iparams->linangle.aA,
860 iparams->linangle.klinB,iparams->linangle.aB);
863 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);
865 case F_QUARTIC_ANGLES:
866 fprintf(fp,"theta=%15.8e",iparams->qangle.theta);
868 fprintf(fp,", c%c=%15.8e",'0'+i,iparams->qangle.c[i]);
872 fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
873 iparams->bham.a,iparams->bham.b,iparams->bham.c);
878 pr_harm(fp,iparams,"b0","cb");
881 pr_harm(fp,iparams,"xi","cx");
884 fprintf(fp,"b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
885 iparams->morse.b0A,iparams->morse.cbA,iparams->morse.betaA,
886 iparams->morse.b0B,iparams->morse.cbB,iparams->morse.betaB);
889 fprintf(fp,"b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
890 iparams->cubic.b0,iparams->cubic.kb,iparams->cubic.kcub);
896 fprintf(fp,"bm=%15.8e, kb=%15.8e\n",iparams->fene.bm,iparams->fene.kb);
899 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",
900 iparams->restraint.lowA,iparams->restraint.up1A,
901 iparams->restraint.up2A,iparams->restraint.kA,
902 iparams->restraint.lowB,iparams->restraint.up1B,
903 iparams->restraint.up2B,iparams->restraint.kB);
909 fprintf(fp,"tab=%d, kA=%15.8e, kB=%15.8e\n",
910 iparams->tab.table,iparams->tab.kA,iparams->tab.kB);
913 fprintf(fp,"alpha=%15.8e\n",iparams->polarize.alpha);
916 fprintf(fp,"alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
917 iparams->anharm_polarize.alpha,
918 iparams->anharm_polarize.drcut,
919 iparams->anharm_polarize.khyp);
922 fprintf(fp,"a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
923 iparams->thole.a,iparams->thole.alpha1,iparams->thole.alpha2,
924 iparams->thole.rfac);
927 fprintf(fp,"al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
928 iparams->wpol.al_x,iparams->wpol.al_y,iparams->wpol.al_z,
929 iparams->wpol.rOH,iparams->wpol.rHH,iparams->wpol.rOD);
932 fprintf(fp,"c6=%15.8e, c12=%15.8e\n",iparams->lj.c6,iparams->lj.c12);
935 fprintf(fp,"c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
936 iparams->lj14.c6A,iparams->lj14.c12A,
937 iparams->lj14.c6B,iparams->lj14.c12B);
940 fprintf(fp,"fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
942 iparams->ljc14.qi,iparams->ljc14.qj,
943 iparams->ljc14.c6,iparams->ljc14.c12);
946 fprintf(fp,"qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
947 iparams->ljcnb.qi,iparams->ljcnb.qj,
948 iparams->ljcnb.c6,iparams->ljcnb.c12);
954 fprintf(fp,"phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
955 iparams->pdihs.phiA,iparams->pdihs.cpA,
956 iparams->pdihs.phiB,iparams->pdihs.cpB,
957 iparams->pdihs.mult);
960 fprintf(fp,"label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
961 iparams->disres.label,iparams->disres.type,
962 iparams->disres.low,iparams->disres.up1,
963 iparams->disres.up2,iparams->disres.kfac);
966 fprintf(fp,"ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
967 iparams->orires.ex,iparams->orires.label,iparams->orires.power,
968 iparams->orires.c,iparams->orires.obs,iparams->orires.kfac);
971 fprintf(fp,"phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
972 iparams->dihres.phiA,iparams->dihres.dphiA,iparams->dihres.kfacA,
973 iparams->dihres.phiB,iparams->dihres.dphiB,iparams->dihres.kfacB);
976 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",
977 iparams->posres.pos0A[XX],iparams->posres.pos0A[YY],
978 iparams->posres.pos0A[ZZ],iparams->posres.fcA[XX],
979 iparams->posres.fcA[YY],iparams->posres.fcA[ZZ],
980 iparams->posres.pos0B[XX],iparams->posres.pos0B[YY],
981 iparams->posres.pos0B[ZZ],iparams->posres.fcB[XX],
982 iparams->posres.fcB[YY],iparams->posres.fcB[ZZ]);
985 fprintf(fp,"pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e\n",
986 iparams->fbposres.pos0[XX], iparams->fbposres.pos0[YY],
987 iparams->fbposres.pos0[ZZ], iparams->fbposres.geom,
988 iparams->fbposres.r, iparams->fbposres.k);
991 for (i=0; i<NR_RBDIHS; i++)
992 fprintf(fp,"%srbcA[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcA[i]);
994 for (i=0; i<NR_RBDIHS; i++)
995 fprintf(fp,"%srbcB[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcB[i]);
999 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
1000 * OPLS potential constants back.
1002 rbcA = iparams->rbdihs.rbcA;
1003 rbcB = iparams->rbdihs.rbcB;
1005 VA[3] = -0.25*rbcA[4];
1006 VA[2] = -0.5*rbcA[3];
1007 VA[1] = 4.0*VA[3]-rbcA[2];
1008 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1010 VB[3] = -0.25*rbcB[4];
1011 VB[2] = -0.5*rbcB[3];
1012 VB[1] = 4.0*VB[3]-rbcB[2];
1013 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1015 for (i=0; i<NR_FOURDIHS; i++)
1016 fprintf(fp,"%sFourA[%d]=%15.8e",i==0?"":", ",i,VA[i]);
1018 for (i=0; i<NR_FOURDIHS; i++)
1019 fprintf(fp,"%sFourB[%d]=%15.8e",i==0?"":", ",i,VB[i]);
1025 fprintf(fp,"dA=%15.8e, dB=%15.8e\n",iparams->constr.dA,iparams->constr.dB);
1028 fprintf(fp,"doh=%15.8e, dhh=%15.8e\n",iparams->settle.doh,
1029 iparams->settle.dhh);
1032 fprintf(fp,"a=%15.8e\n",iparams->vsite.a);
1037 fprintf(fp,"a=%15.8e, b=%15.8e\n",iparams->vsite.a,iparams->vsite.b);
1042 fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
1043 iparams->vsite.a,iparams->vsite.b,iparams->vsite.c);
1046 fprintf(fp,"n=%2d, a=%15.8e\n",iparams->vsiten.n,iparams->vsiten.a);
1051 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);
1054 fprintf(fp, "cmapA=%1d, cmapB=%1d\n",iparams->cmap.cmapA, iparams->cmap.cmapB);
1057 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1058 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1062 void pr_ilist(FILE *fp,int indent,const char *title,
1063 t_functype *functype,t_ilist *ilist, gmx_bool bShowNumbers)
1065 int i,j,k,type,ftype;
1068 if (available(fp,ilist,indent,title) && ilist->nr > 0)
1070 indent=pr_title(fp,indent,title);
1071 (void) pr_indent(fp,indent);
1072 fprintf(fp,"nr: %d\n",ilist->nr);
1073 if (ilist->nr > 0) {
1074 (void) pr_indent(fp,indent);
1075 fprintf(fp,"iatoms:\n");
1076 iatoms=ilist->iatoms;
1077 for (i=j=0; i<ilist->nr;) {
1079 (void) pr_indent(fp,indent+INDENT);
1081 ftype=functype[type];
1082 (void) fprintf(fp,"%d type=%d (%s)",
1083 bShowNumbers?j:-1,bShowNumbers?type:-1,
1084 interaction_function[ftype].name);
1086 for (k=0; k<interaction_function[ftype].nratoms; k++)
1087 (void) fprintf(fp," %u",*(iatoms++));
1088 (void) fprintf(fp,"\n");
1089 i+=1+interaction_function[ftype].nratoms;
1091 fprintf(fp,"%5d%5d\n",i,iatoms[i]);
1099 static void pr_cmap(FILE *fp, int indent, const char *title,
1100 gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1105 dx = 360.0 / cmap_grid->grid_spacing;
1106 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1108 if(available(fp,cmap_grid,indent,title))
1110 fprintf(fp,"%s\n",title);
1112 for(i=0;i<cmap_grid->ngrid;i++)
1115 fprintf(fp,"%8s %8s %8s %8s\n","V","dVdx","dVdy","d2dV");
1117 fprintf(fp,"grid[%3d]={\n",bShowNumbers?i:-1);
1119 for(j=0;j<nelem;j++)
1121 if( (j%cmap_grid->grid_spacing)==0)
1123 fprintf(fp,"%8.1f\n",idx);
1127 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4]);
1128 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+1]);
1129 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+2]);
1130 fprintf(fp,"%8.3f\n",cmap_grid->cmapdata[i].cmap[j*4+3]);
1138 void pr_ffparams(FILE *fp,int indent,const char *title,
1139 gmx_ffparams_t *ffparams,
1140 gmx_bool bShowNumbers)
1144 indent=pr_title(fp,indent,title);
1145 (void) pr_indent(fp,indent);
1146 (void) fprintf(fp,"atnr=%d\n",ffparams->atnr);
1147 (void) pr_indent(fp,indent);
1148 (void) fprintf(fp,"ntypes=%d\n",ffparams->ntypes);
1149 for (i=0; i<ffparams->ntypes; i++) {
1150 (void) pr_indent(fp,indent+INDENT);
1151 (void) fprintf(fp,"functype[%d]=%s, ",
1153 interaction_function[ffparams->functype[i]].name);
1154 pr_iparams(fp,ffparams->functype[i],&ffparams->iparams[i]);
1156 (void) pr_double(fp,indent,"reppow",ffparams->reppow);
1157 (void) pr_real(fp,indent,"fudgeQQ",ffparams->fudgeQQ);
1158 pr_cmap(fp,indent,"cmap",&ffparams->cmap_grid,bShowNumbers);
1161 void pr_idef(FILE *fp,int indent,const char *title,t_idef *idef, gmx_bool bShowNumbers)
1165 if (available(fp,idef,indent,title)) {
1166 indent=pr_title(fp,indent,title);
1167 (void) pr_indent(fp,indent);
1168 (void) fprintf(fp,"atnr=%d\n",idef->atnr);
1169 (void) pr_indent(fp,indent);
1170 (void) fprintf(fp,"ntypes=%d\n",idef->ntypes);
1171 for (i=0; i<idef->ntypes; i++) {
1172 (void) pr_indent(fp,indent+INDENT);
1173 (void) fprintf(fp,"functype[%d]=%s, ",
1175 interaction_function[idef->functype[i]].name);
1176 pr_iparams(fp,idef->functype[i],&idef->iparams[i]);
1178 (void) pr_real(fp,indent,"fudgeQQ",idef->fudgeQQ);
1180 for(j=0; (j<F_NRE); j++)
1181 pr_ilist(fp,indent,interaction_function[j].longname,
1182 idef->functype,&idef->il[j],bShowNumbers);
1186 static int pr_block_title(FILE *fp,int indent,const char *title,t_block *block)
1190 if (available(fp,block,indent,title))
1192 indent=pr_title(fp,indent,title);
1193 (void) pr_indent(fp,indent);
1194 (void) fprintf(fp,"nr=%d\n",block->nr);
1199 static int pr_blocka_title(FILE *fp,int indent,const char *title,t_blocka *block)
1203 if (available(fp,block,indent,title))
1205 indent=pr_title(fp,indent,title);
1206 (void) pr_indent(fp,indent);
1207 (void) fprintf(fp,"nr=%d\n",block->nr);
1208 (void) pr_indent(fp,indent);
1209 (void) fprintf(fp,"nra=%d\n",block->nra);
1214 static void low_pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block, gmx_bool bShowNumbers)
1218 if (available(fp,block,indent,title))
1220 indent=pr_blocka_title(fp,indent,title,block);
1221 for (i=0; i<=block->nr; i++)
1223 (void) pr_indent(fp,indent+INDENT);
1224 (void) fprintf(fp,"%s->index[%d]=%u\n",
1225 title,bShowNumbers?i:-1,block->index[i]);
1227 for (i=0; i<block->nra; i++)
1229 (void) pr_indent(fp,indent+INDENT);
1230 (void) fprintf(fp,"%s->a[%d]=%u\n",
1231 title,bShowNumbers?i:-1,block->a[i]);
1236 void pr_block(FILE *fp,int indent,const char *title,t_block *block,gmx_bool bShowNumbers)
1238 int i,j,ok,size,start,end;
1240 if (available(fp,block,indent,title))
1242 indent=pr_block_title(fp,indent,title,block);
1245 if ((ok=(block->index[start]==0))==0)
1246 (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1248 for (i=0; i<block->nr; i++)
1250 end=block->index[i+1];
1251 size=pr_indent(fp,indent);
1253 size+=fprintf(fp,"%s[%d]={}\n",title,i);
1255 size+=fprintf(fp,"%s[%d]={%d..%d}\n",
1256 title,bShowNumbers?i:-1,
1257 bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1263 void pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block,gmx_bool bShowNumbers)
1265 int i,j,ok,size,start,end;
1267 if (available(fp,block,indent,title))
1269 indent=pr_blocka_title(fp,indent,title,block);
1272 if ((ok=(block->index[start]==0))==0)
1273 (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1275 for (i=0; i<block->nr; i++)
1277 end=block->index[i+1];
1278 size=pr_indent(fp,indent);
1280 size+=fprintf(fp,"%s[%d]={",title,i);
1282 size+=fprintf(fp,"%s[%d][%d..%d]={",
1283 title,bShowNumbers?i:-1,
1284 bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1285 for (j=start; j<end; j++)
1287 if (j>start) size+=fprintf(fp,", ");
1288 if ((size)>(USE_WIDTH))
1290 (void) fprintf(fp,"\n");
1291 size=pr_indent(fp,indent+INDENT);
1293 size+=fprintf(fp,"%u",block->a[j]);
1295 (void) fprintf(fp,"}\n");
1298 if ((end!=block->nra)||(!ok))
1300 (void) pr_indent(fp,indent);
1301 (void) fprintf(fp,"tables inconsistent, dumping complete tables:\n");
1302 low_pr_blocka(fp,indent,title,block,bShowNumbers);
1307 static void pr_strings(FILE *fp,int indent,const char *title,char ***nm,int n, gmx_bool bShowNumbers)
1311 if (available(fp,nm,indent,title))
1313 indent=pr_title_n(fp,indent,title,n);
1316 (void) pr_indent(fp,indent);
1317 (void) fprintf(fp,"%s[%d]={name=\"%s\"}\n",
1318 title,bShowNumbers?i:-1,*(nm[i]));
1323 static void pr_strings2(FILE *fp,int indent,const char *title,
1324 char ***nm,char ***nmB,int n, gmx_bool bShowNumbers)
1328 if (available(fp,nm,indent,title))
1330 indent=pr_title_n(fp,indent,title,n);
1333 (void) pr_indent(fp,indent);
1334 (void) fprintf(fp,"%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1335 title,bShowNumbers?i:-1,*(nm[i]),*(nmB[i]));
1340 static void pr_resinfo(FILE *fp,int indent,const char *title,t_resinfo *resinfo,int n, gmx_bool bShowNumbers)
1344 if (available(fp,resinfo,indent,title))
1346 indent=pr_title_n(fp,indent,title,n);
1349 (void) pr_indent(fp,indent);
1350 (void) fprintf(fp,"%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1351 title,bShowNumbers?i:-1,
1352 *(resinfo[i].name),resinfo[i].nr,
1353 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1358 static void pr_atom(FILE *fp,int indent,const char *title,t_atom *atom,int n)
1362 if (available(fp,atom,indent,title)) {
1363 indent=pr_title_n(fp,indent,title,n);
1364 for (i=0; i<n; i++) {
1365 (void) pr_indent(fp,indent);
1366 fprintf(fp,"%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1367 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1368 title,i,atom[i].type,atom[i].typeB,ptype_str[atom[i].ptype],
1369 atom[i].m,atom[i].q,atom[i].mB,atom[i].qB,
1370 atom[i].resind,atom[i].atomnumber);
1375 static void pr_grps(FILE *fp,int indent,const char *title,t_grps grps[],
1376 char **grpname[], gmx_bool bShowNumbers)
1380 for(i=0; (i<egcNR); i++)
1382 fprintf(fp,"%s[%-12s] nr=%d, name=[",title,gtypes[i],grps[i].nr);
1383 for(j=0; (j<grps[i].nr); j++)
1385 fprintf(fp," %s",*(grpname[grps[i].nm_ind[j]]));
1391 static void pr_groups(FILE *fp,int indent,const char *title,
1392 gmx_groups_t *groups,
1393 gmx_bool bShowNumbers)
1398 pr_grps(fp,indent,"grp",groups->grps,groups->grpname,bShowNumbers);
1399 pr_strings(fp,indent,"grpname",groups->grpname,groups->ngrpname,bShowNumbers);
1401 (void) pr_indent(fp,indent);
1402 fprintf(fp,"groups ");
1403 for(g=0; g<egcNR; g++)
1405 printf(" %5.5s",gtypes[g]);
1409 (void) pr_indent(fp,indent);
1410 fprintf(fp,"allocated ");
1412 for(g=0; g<egcNR; g++)
1414 printf(" %5d",groups->ngrpnr[g]);
1415 nat_max = max(nat_max,groups->ngrpnr[g]);
1421 (void) pr_indent(fp,indent);
1422 fprintf(fp,"groupnr[%5s] =","*");
1423 for(g=0; g<egcNR; g++)
1425 fprintf(fp," %3d ",0);
1431 for(i=0; i<nat_max; i++)
1433 (void) pr_indent(fp,indent);
1434 fprintf(fp,"groupnr[%5d] =",i);
1435 for(g=0; g<egcNR; g++)
1438 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1445 void pr_atoms(FILE *fp,int indent,const char *title,t_atoms *atoms,
1446 gmx_bool bShownumbers)
1448 if (available(fp,atoms,indent,title))
1450 indent=pr_title(fp,indent,title);
1451 pr_atom(fp,indent,"atom",atoms->atom,atoms->nr);
1452 pr_strings(fp,indent,"atom",atoms->atomname,atoms->nr,bShownumbers);
1453 pr_strings2(fp,indent,"type",atoms->atomtype,atoms->atomtypeB,atoms->nr,bShownumbers);
1454 pr_resinfo(fp,indent,"residue",atoms->resinfo,atoms->nres,bShownumbers);
1459 void pr_atomtypes(FILE *fp,int indent,const char *title,t_atomtypes *atomtypes,
1460 gmx_bool bShowNumbers)
1463 if (available(fp,atomtypes,indent,title))
1465 indent=pr_title(fp,indent,title);
1466 for(i=0;i<atomtypes->nr;i++) {
1467 pr_indent(fp,indent);
1469 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1470 bShowNumbers?i:-1,atomtypes->radius[i],atomtypes->vol[i],
1471 atomtypes->gb_radius[i],
1472 atomtypes->surftens[i],atomtypes->atomnumber[i],atomtypes->S_hct[i]);
1477 static void pr_moltype(FILE *fp,int indent,const char *title,
1478 gmx_moltype_t *molt,int n,
1479 gmx_ffparams_t *ffparams,
1480 gmx_bool bShowNumbers)
1484 indent = pr_title_n(fp,indent,title,n);
1485 (void) pr_indent(fp,indent);
1486 (void) fprintf(fp,"name=\"%s\"\n",*(molt->name));
1487 pr_atoms(fp,indent,"atoms",&(molt->atoms),bShowNumbers);
1488 pr_block(fp,indent,"cgs",&molt->cgs, bShowNumbers);
1489 pr_blocka(fp,indent,"excls",&molt->excls, bShowNumbers);
1490 for(j=0; (j<F_NRE); j++) {
1491 pr_ilist(fp,indent,interaction_function[j].longname,
1492 ffparams->functype,&molt->ilist[j],bShowNumbers);
1496 static void pr_molblock(FILE *fp,int indent,const char *title,
1497 gmx_molblock_t *molb,int n,
1498 gmx_moltype_t *molt,
1499 gmx_bool bShowNumbers)
1501 indent = pr_title_n(fp,indent,title,n);
1502 (void) pr_indent(fp,indent);
1503 (void) fprintf(fp,"%-20s = %d \"%s\"\n",
1504 "moltype",molb->type,*(molt[molb->type].name));
1505 pr_int(fp,indent,"#molecules",molb->nmol);
1506 pr_int(fp,indent,"#atoms_mol",molb->natoms_mol);
1507 pr_int(fp,indent,"#posres_xA",molb->nposres_xA);
1508 if (molb->nposres_xA > 0) {
1509 pr_rvecs(fp,indent,"posres_xA",molb->posres_xA,molb->nposres_xA);
1511 pr_int(fp,indent,"#posres_xB",molb->nposres_xB);
1512 if (molb->nposres_xB > 0) {
1513 pr_rvecs(fp,indent,"posres_xB",molb->posres_xB,molb->nposres_xB);
1517 void pr_mtop(FILE *fp,int indent,const char *title,gmx_mtop_t *mtop,
1518 gmx_bool bShowNumbers)
1522 if (available(fp,mtop,indent,title)) {
1523 indent=pr_title(fp,indent,title);
1524 (void) pr_indent(fp,indent);
1525 (void) fprintf(fp,"name=\"%s\"\n",*(mtop->name));
1526 pr_int(fp,indent,"#atoms",mtop->natoms);
1527 pr_int(fp,indent,"#molblock",mtop->nmolblock);
1528 for(mb=0; mb<mtop->nmolblock; mb++) {
1529 pr_molblock(fp,indent,"molblock",&mtop->molblock[mb],mb,
1530 mtop->moltype,bShowNumbers);
1532 pr_ffparams(fp,indent,"ffparams",&(mtop->ffparams),bShowNumbers);
1533 pr_atomtypes(fp,indent,"atomtypes",&(mtop->atomtypes),bShowNumbers);
1534 for(mt=0; mt<mtop->nmoltype; mt++) {
1535 pr_moltype(fp,indent,"moltype",&mtop->moltype[mt],mt,
1536 &mtop->ffparams,bShowNumbers);
1538 pr_groups(fp,indent,"groups",&mtop->groups,bShowNumbers);
1542 void pr_top(FILE *fp,int indent,const char *title,t_topology *top, gmx_bool bShowNumbers)
1544 if (available(fp,top,indent,title)) {
1545 indent=pr_title(fp,indent,title);
1546 (void) pr_indent(fp,indent);
1547 (void) fprintf(fp,"name=\"%s\"\n",*(top->name));
1548 pr_atoms(fp,indent,"atoms",&(top->atoms),bShowNumbers);
1549 pr_atomtypes(fp,indent,"atomtypes",&(top->atomtypes),bShowNumbers);
1550 pr_block(fp,indent,"cgs",&top->cgs, bShowNumbers);
1551 pr_block(fp,indent,"mols",&top->mols, bShowNumbers);
1552 pr_blocka(fp,indent,"excls",&top->excls, bShowNumbers);
1553 pr_idef(fp,indent,"idef",&top->idef,bShowNumbers);
1557 void pr_header(FILE *fp,int indent,const char *title,t_tpxheader *sh)
1561 if (available(fp,sh,indent,title))
1563 indent=pr_title(fp,indent,title);
1564 pr_indent(fp,indent);
1565 fprintf(fp,"bIr = %spresent\n",sh->bIr?"":"not ");
1566 pr_indent(fp,indent);
1567 fprintf(fp,"bBox = %spresent\n",sh->bBox?"":"not ");
1568 pr_indent(fp,indent);
1569 fprintf(fp,"bTop = %spresent\n",sh->bTop?"":"not ");
1570 pr_indent(fp,indent);
1571 fprintf(fp,"bX = %spresent\n",sh->bX?"":"not ");
1572 pr_indent(fp,indent);
1573 fprintf(fp,"bV = %spresent\n",sh->bV?"":"not ");
1574 pr_indent(fp,indent);
1575 fprintf(fp,"bF = %spresent\n",sh->bF?"":"not ");
1577 pr_indent(fp,indent);
1578 fprintf(fp,"natoms = %d\n",sh->natoms);
1579 pr_indent(fp,indent);
1580 fprintf(fp,"lambda = %e\n",sh->lambda);
1584 void pr_commrec(FILE *fp,int indent,t_commrec *cr)
1586 pr_indent(fp,indent);
1587 fprintf(fp,"commrec:\n");
1589 pr_indent(fp,indent);
1590 fprintf(fp,"nodeid = %d\n",cr->nodeid);
1591 pr_indent(fp,indent);
1592 fprintf(fp,"nnodes = %d\n",cr->nnodes);
1593 pr_indent(fp,indent);
1594 fprintf(fp,"npmenodes = %d\n",cr->npmenodes);
1596 pr_indent(fp,indent);
1597 fprintf(fp,"threadid = %d\n",cr->threadid);
1598 pr_indent(fp,indent);
1599 fprintf(fp,"nthreads = %d\n",cr->nthreads);