2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
41 /* This file is completely threadsafe - please keep it that way! */
53 int pr_indent(FILE *fp, int n)
57 for (i = 0; i < n; i++)
59 (void) fprintf(fp, " ");
64 int available(FILE *fp, void *p, int indent, const char *title)
70 pr_indent(fp, indent);
72 (void) fprintf(fp, "%s: not available\n", title);
77 int pr_title(FILE *fp, int indent, const char *title)
79 (void) pr_indent(fp, indent);
80 (void) fprintf(fp, "%s:\n", title);
81 return (indent+INDENT);
84 int pr_title_n(FILE *fp, int indent, const char *title, int n)
86 (void) pr_indent(fp, indent);
87 (void) fprintf(fp, "%s (%d):\n", title, n);
88 return (indent+INDENT);
91 int pr_title_nxn(FILE *fp, int indent, const char *title, int n1, int n2)
93 (void) pr_indent(fp, indent);
94 (void) fprintf(fp, "%s (%dx%d):\n", title, n1, n2);
95 return (indent+INDENT);
98 void pr_ivec(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
102 if (available(fp, vec, indent, title))
104 indent = pr_title_n(fp, indent, title, n);
105 for (i = 0; i < n; i++)
107 (void) pr_indent(fp, indent);
108 (void) fprintf(fp, "%s[%d]=%d\n", title, bShowNumbers ? i : -1, vec[i]);
113 void pr_ivec_block(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
117 if (available(fp, vec, indent, title))
119 indent = pr_title_n(fp, indent, title, n);
124 while (j < n && vec[j] == vec[j-1]+1)
128 /* Print consecutive groups of 3 or more as blocks */
133 (void) pr_indent(fp, indent);
134 (void) fprintf(fp, "%s[%d]=%d\n",
135 title, bShowNumbers ? i : -1, vec[i]);
141 (void) pr_indent(fp, indent);
142 (void) fprintf(fp, "%s[%d,...,%d] = {%d,...,%d}\n",
144 bShowNumbers ? i : -1,
145 bShowNumbers ? j-1 : -1,
153 void pr_bvec(FILE *fp, int indent, const char *title, gmx_bool vec[], int n, gmx_bool bShowNumbers)
157 if (available(fp, vec, indent, title))
159 indent = pr_title_n(fp, indent, title, n);
160 for (i = 0; i < n; i++)
162 (void) pr_indent(fp, indent);
163 (void) fprintf(fp, "%s[%d]=%s\n", title, bShowNumbers ? i : -1,
169 void pr_ivecs(FILE *fp, int indent, const char *title, ivec vec[], int n, gmx_bool bShowNumbers)
173 if (available(fp, vec, indent, title))
175 indent = pr_title_nxn(fp, indent, title, n, DIM);
176 for (i = 0; i < n; i++)
178 (void) pr_indent(fp, indent);
179 (void) fprintf(fp, "%s[%d]={", title, bShowNumbers ? i : -1);
180 for (j = 0; j < DIM; j++)
184 (void) fprintf(fp, ", ");
186 fprintf(fp, "%d", vec[i][j]);
188 (void) fprintf(fp, "}\n");
193 void pr_rvec(FILE *fp, int indent, const char *title, real vec[], int n, gmx_bool bShowNumbers)
197 if (available(fp, vec, indent, title))
199 indent = pr_title_n(fp, indent, title, n);
200 for (i = 0; i < n; i++)
202 pr_indent(fp, indent);
203 fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
208 void pr_dvec(FILE *fp, int indent, const char *title, double vec[], int n, gmx_bool bShowNumbers)
212 if (available(fp, vec, indent, title))
214 indent = pr_title_n(fp, indent, title, n);
215 for (i = 0; i < n; i++)
217 pr_indent(fp, indent);
218 fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
225 void pr_mat(FILE *fp,int indent,char *title,matrix m)
229 if (available(fp,m,indent,title)) {
230 indent=pr_title_n(fp,indent,title,n);
232 pr_indent(fp,indent);
233 fprintf(fp,"%s[%d]=%12.5e %12.5e %12.5e\n",
234 title,bShowNumbers?i:-1,m[i][XX],m[i][YY],m[i][ZZ]);
240 void pr_rvecs_len(FILE *fp, int indent, const char *title, rvec vec[], int n)
244 if (available(fp, vec, indent, title))
246 indent = pr_title_nxn(fp, indent, title, n, DIM);
247 for (i = 0; i < n; i++)
249 (void) pr_indent(fp, indent);
250 (void) fprintf(fp, "%s[%5d]={", title, i);
251 for (j = 0; j < DIM; j++)
255 (void) fprintf(fp, ", ");
257 (void) fprintf(fp, "%12.5e", vec[i][j]);
259 (void) fprintf(fp, "} len=%12.5e\n", norm(vec[i]));
264 void pr_rvecs(FILE *fp, int indent, const char *title, rvec vec[], int n)
266 const char *fshort = "%12.5e";
267 const char *flong = "%15.8e";
271 if (getenv("LONGFORMAT") != NULL)
280 if (available(fp, vec, indent, title))
282 indent = pr_title_nxn(fp, indent, title, n, DIM);
283 for (i = 0; i < n; i++)
285 (void) pr_indent(fp, indent);
286 (void) fprintf(fp, "%s[%5d]={", title, i);
287 for (j = 0; j < DIM; j++)
291 (void) fprintf(fp, ", ");
293 (void) fprintf(fp, format, vec[i][j]);
295 (void) fprintf(fp, "}\n");
301 void pr_reals(FILE *fp, int indent, const char *title, real *vec, int n)
305 if (available(fp, vec, indent, title))
307 (void) pr_indent(fp, indent);
308 (void) fprintf(fp, "%s:\t", title);
309 for (i = 0; i < n; i++)
311 fprintf(fp, " %10g", vec[i]);
313 (void) fprintf(fp, "\n");
317 void pr_doubles(FILE *fp, int indent, const char *title, double *vec, int n)
321 if (available(fp, vec, indent, title))
323 (void) pr_indent(fp, indent);
324 (void) fprintf(fp, "%s:\t", title);
325 for (i = 0; i < n; i++)
327 fprintf(fp, " %10g", vec[i]);
329 (void) fprintf(fp, "\n");
333 static void pr_int(FILE *fp, int indent, const char *title, int i)
335 pr_indent(fp, indent);
336 fprintf(fp, "%-20s = %d\n", title, i);
339 static void pr_gmx_large_int(FILE *fp, int indent, const char *title, gmx_int64_t i)
341 char buf[STEPSTRSIZE];
343 pr_indent(fp, indent);
344 fprintf(fp, "%-20s = %s\n", title, gmx_step_str(i, buf));
347 static void pr_real(FILE *fp, int indent, const char *title, real r)
349 pr_indent(fp, indent);
350 fprintf(fp, "%-20s = %g\n", title, r);
353 static void pr_double(FILE *fp, int indent, const char *title, double d)
355 pr_indent(fp, indent);
356 fprintf(fp, "%-20s = %g\n", title, d);
359 static void pr_str(FILE *fp, int indent, const char *title, const char *s)
361 pr_indent(fp, indent);
362 fprintf(fp, "%-20s = %s\n", title, s);
365 void pr_qm_opts(FILE *fp, int indent, const char *title, t_grpopts *opts)
369 fprintf(fp, "%s:\n", title);
371 pr_int(fp, indent, "ngQM", opts->ngQM);
374 pr_ivec(fp, indent, "QMmethod", opts->QMmethod, opts->ngQM, FALSE);
375 pr_ivec(fp, indent, "QMbasis", opts->QMbasis, opts->ngQM, FALSE);
376 pr_ivec(fp, indent, "QMcharge", opts->QMcharge, opts->ngQM, FALSE);
377 pr_ivec(fp, indent, "QMmult", opts->QMmult, opts->ngQM, FALSE);
378 pr_bvec(fp, indent, "bSH", opts->bSH, opts->ngQM, FALSE);
379 pr_ivec(fp, indent, "CASorbitals", opts->CASorbitals, opts->ngQM, FALSE);
380 pr_ivec(fp, indent, "CASelectrons", opts->CASelectrons, opts->ngQM, FALSE);
381 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
382 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
383 pr_ivec(fp, indent, "SAsteps", opts->SAsteps, opts->ngQM, FALSE);
384 pr_bvec(fp, indent, "bOPT", opts->bOPT, opts->ngQM, FALSE);
385 pr_bvec(fp, indent, "bTS", opts->bTS, opts->ngQM, FALSE);
389 static void pr_grp_opts(FILE *out, int indent, const char *title, t_grpopts *opts,
396 fprintf(out, "%s:\n", title);
399 pr_indent(out, indent);
400 fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
401 for (i = 0; (i < opts->ngtc); i++)
403 fprintf(out, " %10g", opts->nrdf[i]);
407 pr_indent(out, indent);
408 fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
409 for (i = 0; (i < opts->ngtc); i++)
411 fprintf(out, " %10g", opts->ref_t[i]);
415 pr_indent(out, indent);
416 fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
417 for (i = 0; (i < opts->ngtc); i++)
419 fprintf(out, " %10g", opts->tau_t[i]);
423 /* Pretty-print the simulated annealing info */
424 fprintf(out, "anneal%s", bMDPformat ? " = " : ":");
425 for (i = 0; (i < opts->ngtc); i++)
427 fprintf(out, " %10s", EANNEAL(opts->annealing[i]));
431 fprintf(out, "ann-npoints%s", bMDPformat ? " = " : ":");
432 for (i = 0; (i < opts->ngtc); i++)
434 fprintf(out, " %10d", opts->anneal_npoints[i]);
438 for (i = 0; (i < opts->ngtc); i++)
440 if (opts->anneal_npoints[i] > 0)
442 fprintf(out, "ann. times [%d]:\t", i);
443 for (j = 0; (j < opts->anneal_npoints[i]); j++)
445 fprintf(out, " %10.1f", opts->anneal_time[i][j]);
448 fprintf(out, "ann. temps [%d]:\t", i);
449 for (j = 0; (j < opts->anneal_npoints[i]); j++)
451 fprintf(out, " %10.1f", opts->anneal_temp[i][j]);
457 pr_indent(out, indent);
458 fprintf(out, "acc:\t");
459 for (i = 0; (i < opts->ngacc); i++)
461 for (m = 0; (m < DIM); m++)
463 fprintf(out, " %10g", opts->acc[i][m]);
468 pr_indent(out, indent);
469 fprintf(out, "nfreeze:");
470 for (i = 0; (i < opts->ngfrz); i++)
472 for (m = 0; (m < DIM); m++)
474 fprintf(out, " %10s", opts->nFreeze[i][m] ? "Y" : "N");
480 for (i = 0; (i < opts->ngener); i++)
482 pr_indent(out, indent);
483 fprintf(out, "energygrp-flags[%3d]:", i);
484 for (m = 0; (m < opts->ngener); m++)
486 fprintf(out, " %d", opts->egp_flags[opts->ngener*i+m]);
494 static void pr_matrix(FILE *fp, int indent, const char *title, rvec *m,
499 fprintf(fp, "%-10s = %g %g %g %g %g %g\n", title,
500 m[XX][XX], m[YY][YY], m[ZZ][ZZ], m[XX][YY], m[XX][ZZ], m[YY][ZZ]);
504 pr_rvecs(fp, indent, title, m, DIM);
508 static void pr_cosine(FILE *fp, int indent, const char *title, t_cosines *cos,
515 fprintf(fp, "%s = %d\n", title, cos->n);
519 indent = pr_title(fp, indent, title);
520 (void) pr_indent(fp, indent);
521 fprintf(fp, "n = %d\n", cos->n);
524 (void) pr_indent(fp, indent+2);
526 for (j = 0; (j < cos->n); j++)
528 fprintf(fp, " %e", cos->a[j]);
531 (void) pr_indent(fp, indent+2);
532 fprintf(fp, "phi =");
533 for (j = 0; (j < cos->n); j++)
535 fprintf(fp, " %e", cos->phi[j]);
542 #define PS(t, s) pr_str(fp, indent, t, s)
543 #define PI(t, s) pr_int(fp, indent, t, s)
544 #define PSTEP(t, s) pr_gmx_large_int(fp, indent, t, s)
545 #define PR(t, s) pr_real(fp, indent, t, s)
546 #define PD(t, s) pr_double(fp, indent, t, s)
548 static void pr_pull_group(FILE *fp, int indent, int g, t_pull_group *pgrp)
550 pr_indent(fp, indent);
551 fprintf(fp, "pull-group %d:\n", g);
553 pr_ivec_block(fp, indent, "atom", pgrp->ind, pgrp->nat, TRUE);
554 pr_rvec(fp, indent, "weight", pgrp->weight, pgrp->nweight, TRUE);
555 PI("pbcatom", pgrp->pbcatom);
558 static void pr_pull_coord(FILE *fp, int indent, int c, t_pull_coord *pcrd)
560 pr_indent(fp, indent);
561 fprintf(fp, "pull-coord %d:\n", c);
562 PI("group[0]", pcrd->group[0]);
563 PI("group[1]", pcrd->group[1]);
564 pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
565 pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
566 PR("init", pcrd->init);
567 PR("rate", pcrd->rate);
572 static void pr_simtempvals(FILE *fp, int indent, t_simtemp *simtemp, int n_lambda)
574 PR("simtemp_low", simtemp->simtemp_low);
575 PR("simtemp_high", simtemp->simtemp_high);
576 PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
577 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
580 static void pr_expandedvals(FILE *fp, int indent, t_expanded *expand, int n_lambda)
583 PI("nstexpanded", expand->nstexpanded);
584 PS("lambda-stats", elamstats_names[expand->elamstats]);
585 PS("lambda-mc-move", elmcmove_names[expand->elmcmove]);
586 PI("lmc-repeats", expand->lmc_repeats);
587 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
588 PI("lmc-nstart", expand->lmc_forced_nstart);
589 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
590 PI("nst-transition-matrix", expand->nstTij);
591 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
592 PI("weight-c-range", expand->c_range); /* default is just C=0 */
593 PR("wl-scale", expand->wl_scale);
594 PR("init-wl-delta", expand->init_wl_delta);
595 PR("wl-ratio", expand->wl_ratio);
596 PS("bWLoneovert", EBOOL(expand->bWLoneovert));
597 PI("lmc-seed", expand->lmc_seed);
598 PR("mc-temperature", expand->mc_temp);
599 PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
600 if (expand->elmceq == elmceqNUMATLAM)
602 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
604 if (expand->elmceq == elmceqSAMPLES)
606 PI("weight-equil-number-samples", expand->equil_samples);
608 if (expand->elmceq == elmceqSTEPS)
610 PI("weight-equil-number-steps", expand->equil_steps);
612 if (expand->elmceq == elmceqWLDELTA)
614 PR("weight-equil-wl-delta", expand->equil_wl_delta);
616 if (expand->elmceq == elmceqRATIO)
618 PR("weight-equil-count-ratio", expand->equil_ratio);
621 pr_indent(fp, indent);
622 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
623 PS("init-weights", EBOOL(expand->bInit_weights));
626 static void pr_fepvals(FILE *fp, int indent, t_lambda *fep, gmx_bool bMDPformat)
630 PI("nstdhdl", fep->nstdhdl);
631 PI("init-lambda-state", fep->init_fep_state);
632 PR("init-lambda", fep->init_lambda);
633 PR("delta-lambda", fep->delta_lambda);
636 PI("n-lambdas", fep->n_lambda);
638 if (fep->n_lambda > 0)
640 pr_indent(fp, indent);
641 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
642 for (i = 0; i < efptNR; i++)
644 fprintf(fp, "%18s = ", efpt_names[i]);
645 if (fep->separate_dvdl[i])
647 fprintf(fp, " TRUE");
651 fprintf(fp, " FALSE");
655 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
656 for (i = 0; i < efptNR; i++)
658 fprintf(fp, "%18s = ", efpt_names[i]);
659 for (j = 0; j < fep->n_lambda; j++)
661 fprintf(fp, " %10g", fep->all_lambda[i][j]);
666 PI("calc-lambda-neighbors", fep->lambda_neighbors);
668 PR("sc-alpha", fep->sc_alpha);
669 PS("bScCoul", EBOOL(fep->bScCoul));
670 PS("bScPrintEnergy", EBOOL(fep->bPrintEnergy));
671 PI("sc-power", fep->sc_power);
672 PR("sc-r-power", fep->sc_r_power);
673 PR("sc-sigma", fep->sc_sigma);
674 PR("sc-sigma-min", fep->sc_sigma_min);
675 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
676 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
677 PI("dh-hist-size", fep->dh_hist_size);
678 PD("dh-hist-spacing", fep->dh_hist_spacing);
681 static void pr_pull(FILE *fp, int indent, t_pull *pull)
685 PS("pull-geometry", EPULLGEOM(pull->eGeom));
686 pr_ivec(fp, indent, "pull-dim", pull->dim, DIM, TRUE);
687 PR("pull-r1", pull->cyl_r1);
688 PR("pull-r0", pull->cyl_r0);
689 PR("pull-constr-tol", pull->constr_tol);
690 PS("pull-bPrintRef", EBOOL(pull->bPrintRef));
691 PI("pull-nstxout", pull->nstxout);
692 PI("pull-nstfout", pull->nstfout);
693 PI("pull-ngroup", pull->ngroup);
694 for (g = 0; g < pull->ngroup; g++)
696 pr_pull_group(fp, indent, g, &pull->group[g]);
698 PI("pull-ncoord", pull->ncoord);
699 for (g = 0; g < pull->ncoord; g++)
701 pr_pull_coord(fp, indent, g, &pull->coord[g]);
705 static void pr_rotgrp(FILE *fp, int indent, int g, t_rotgrp *rotg)
707 pr_indent(fp, indent);
708 fprintf(fp, "rotation_group %d:\n", g);
710 PS("type", EROTGEOM(rotg->eType));
711 PS("massw", EBOOL(rotg->bMassW));
712 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
713 pr_rvecs(fp, indent, "x_ref", rotg->x_ref, rotg->nat);
714 pr_rvec(fp, indent, "vec", rotg->vec, DIM, TRUE);
715 pr_rvec(fp, indent, "pivot", rotg->pivot, DIM, TRUE);
716 PR("rate", rotg->rate);
718 PR("slab_dist", rotg->slab_dist);
719 PR("min_gaussian", rotg->min_gaussian);
720 PR("epsilon", rotg->eps);
721 PS("fit_method", EROTFIT(rotg->eFittype));
722 PI("potfitangle_nstep", rotg->PotAngle_nstep);
723 PR("potfitangle_step", rotg->PotAngle_step);
726 static void pr_rot(FILE *fp, int indent, t_rot *rot)
730 PI("rot_nstrout", rot->nstrout);
731 PI("rot_nstsout", rot->nstsout);
732 PI("rot_ngrp", rot->ngrp);
733 for (g = 0; g < rot->ngrp; g++)
735 pr_rotgrp(fp, indent, g, &rot->grp[g]);
739 void pr_inputrec(FILE *fp, int indent, const char *title, t_inputrec *ir,
742 const char *infbuf = "inf";
745 if (available(fp, ir, indent, title))
749 indent = pr_title(fp, indent, title);
751 PS("integrator", EI(ir->eI));
752 PSTEP("nsteps", ir->nsteps);
753 PSTEP("init-step", ir->init_step);
754 PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
755 PS("ns_type", ENS(ir->ns_type));
756 PI("nstlist", ir->nstlist);
757 PI("ndelta", ir->ndelta);
758 PI("nstcomm", ir->nstcomm);
759 PS("comm-mode", ECOM(ir->comm_mode));
760 PI("nstlog", ir->nstlog);
761 PI("nstxout", ir->nstxout);
762 PI("nstvout", ir->nstvout);
763 PI("nstfout", ir->nstfout);
764 PI("nstcalcenergy", ir->nstcalcenergy);
765 PI("nstenergy", ir->nstenergy);
766 PI("nstxtcout", ir->nstxtcout);
767 PR("init-t", ir->init_t);
768 PR("delta-t", ir->delta_t);
770 PR("xtcprec", ir->xtcprec);
771 PR("fourierspacing", ir->fourier_spacing);
775 PI("pme-order", ir->pme_order);
776 PR("ewald-rtol", ir->ewald_rtol);
777 PR("ewald-rtol-lj", ir->ewald_rtol_lj);
778 PR("ewald-geometry", ir->ewald_geometry);
779 PR("epsilon-surface", ir->epsilon_surface);
780 PS("optimize-fft", EBOOL(ir->bOptFFT));
781 PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir->ljpme_combination_rule));
782 PS("ePBC", EPBC(ir->ePBC));
783 PS("bPeriodicMols", EBOOL(ir->bPeriodicMols));
784 PS("bContinuation", EBOOL(ir->bContinuation));
785 PS("bShakeSOR", EBOOL(ir->bShakeSOR));
786 PS("etc", ETCOUPLTYPE(ir->etc));
787 PS("bPrintNHChains", EBOOL(ir->bPrintNHChains));
788 PI("nsttcouple", ir->nsttcouple);
789 PS("epc", EPCOUPLTYPE(ir->epc));
790 PS("epctype", EPCOUPLTYPETYPE(ir->epct));
791 PI("nstpcouple", ir->nstpcouple);
792 PR("tau-p", ir->tau_p);
793 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
794 pr_matrix(fp, indent, "compress", ir->compress, bMDPformat);
795 PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
798 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX],
799 ir->posres_com[YY], ir->posres_com[ZZ]);
803 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
807 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX],
808 ir->posres_comB[YY], ir->posres_comB[ZZ]);
812 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
814 PR("verlet-buffer-tolerance", ir->verletbuf_tol);
815 PR("rlist", ir->rlist);
816 PR("rlistlong", ir->rlistlong);
817 PR("nstcalclr", ir->nstcalclr);
818 PR("rtpi", ir->rtpi);
819 PS("coulombtype", EELTYPE(ir->coulombtype));
820 PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
821 PR("rcoulomb-switch", ir->rcoulomb_switch);
822 PR("rcoulomb", ir->rcoulomb);
823 PS("vdwtype", EVDWTYPE(ir->vdwtype));
824 PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
825 PR("rvdw-switch", ir->rvdw_switch);
826 PR("rvdw", ir->rvdw);
827 if (ir->epsilon_r != 0)
829 PR("epsilon-r", ir->epsilon_r);
833 PS("epsilon-r", infbuf);
835 if (ir->epsilon_rf != 0)
837 PR("epsilon-rf", ir->epsilon_rf);
841 PS("epsilon-rf", infbuf);
843 PR("tabext", ir->tabext);
844 PS("implicit-solvent", EIMPLICITSOL(ir->implicit_solvent));
845 PS("gb-algorithm", EGBALGORITHM(ir->gb_algorithm));
846 PR("gb-epsilon-solvent", ir->gb_epsilon_solvent);
847 PI("nstgbradii", ir->nstgbradii);
848 PR("rgbradii", ir->rgbradii);
849 PR("gb-saltconc", ir->gb_saltconc);
850 PR("gb-obc-alpha", ir->gb_obc_alpha);
851 PR("gb-obc-beta", ir->gb_obc_beta);
852 PR("gb-obc-gamma", ir->gb_obc_gamma);
853 PR("gb-dielectric-offset", ir->gb_dielectric_offset);
854 PS("sa-algorithm", ESAALGORITHM(ir->gb_algorithm));
855 PR("sa-surface-tension", ir->sa_surface_tension);
856 PS("DispCorr", EDISPCORR(ir->eDispCorr));
857 PS("bSimTemp", EBOOL(ir->bSimTemp));
860 pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
862 PS("free-energy", EFEPTYPE(ir->efep));
863 if (ir->efep != efepNO || ir->bSimTemp)
865 pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
869 pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
872 PI("nwall", ir->nwall);
873 PS("wall-type", EWALLTYPE(ir->wall_type));
874 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
875 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
876 PR("wall-density[0]", ir->wall_density[0]);
877 PR("wall-density[1]", ir->wall_density[1]);
878 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
880 PS("pull", EPULLTYPE(ir->ePull));
881 if (ir->ePull != epullNO)
883 pr_pull(fp, indent, ir->pull);
886 PS("rotation", EBOOL(ir->bRot));
889 pr_rot(fp, indent, ir->rot);
892 PS("disre", EDISRETYPE(ir->eDisre));
893 PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
894 PS("disre-mixed", EBOOL(ir->bDisreMixed));
895 PR("dr-fc", ir->dr_fc);
896 PR("dr-tau", ir->dr_tau);
897 PR("nstdisreout", ir->nstdisreout);
898 PR("orires-fc", ir->orires_fc);
899 PR("orires-tau", ir->orires_tau);
900 PR("nstorireout", ir->nstorireout);
902 PR("dihre-fc", ir->dihre_fc);
904 PR("em-stepsize", ir->em_stepsize);
905 PR("em-tol", ir->em_tol);
906 PI("niter", ir->niter);
907 PR("fc-stepsize", ir->fc_stepsize);
908 PI("nstcgsteep", ir->nstcgsteep);
909 PI("nbfgscorr", ir->nbfgscorr);
911 PS("ConstAlg", ECONSTRTYPE(ir->eConstrAlg));
912 PR("shake-tol", ir->shake_tol);
913 PI("lincs-order", ir->nProjOrder);
914 PR("lincs-warnangle", ir->LincsWarnAngle);
915 PI("lincs-iter", ir->nLincsIter);
916 PR("bd-fric", ir->bd_fric);
917 PI("ld-seed", ir->ld_seed);
918 PR("cos-accel", ir->cos_accel);
919 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
921 PS("adress", EBOOL(ir->bAdress));
924 PS("adress_type", EADRESSTYPE(ir->adress->type));
925 PR("adress_const_wf", ir->adress->const_wf);
926 PR("adress_ex_width", ir->adress->ex_width);
927 PR("adress_hy_width", ir->adress->hy_width);
928 PS("adress_interface_correction", EADRESSICTYPE(ir->adress->icor));
929 PS("adress_site", EADRESSSITETYPE(ir->adress->site));
930 PR("adress_ex_force_cap", ir->adress->ex_forcecap);
931 PS("adress_do_hybridpairs", EBOOL(ir->adress->do_hybridpairs));
933 pr_rvec(fp, indent, "adress_reference_coords", ir->adress->refs, DIM, TRUE);
935 PI("userint1", ir->userint1);
936 PI("userint2", ir->userint2);
937 PI("userint3", ir->userint3);
938 PI("userint4", ir->userint4);
939 PR("userreal1", ir->userreal1);
940 PR("userreal2", ir->userreal2);
941 PR("userreal3", ir->userreal3);
942 PR("userreal4", ir->userreal4);
943 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
944 pr_cosine(fp, indent, "efield-x", &(ir->ex[XX]), bMDPformat);
945 pr_cosine(fp, indent, "efield-xt", &(ir->et[XX]), bMDPformat);
946 pr_cosine(fp, indent, "efield-y", &(ir->ex[YY]), bMDPformat);
947 pr_cosine(fp, indent, "efield-yt", &(ir->et[YY]), bMDPformat);
948 pr_cosine(fp, indent, "efield-z", &(ir->ex[ZZ]), bMDPformat);
949 pr_cosine(fp, indent, "efield-zt", &(ir->et[ZZ]), bMDPformat);
950 PS("bQMMM", EBOOL(ir->bQMMM));
951 PI("QMconstraints", ir->QMconstraints);
952 PI("QMMMscheme", ir->QMMMscheme);
953 PR("scalefactor", ir->scalefactor);
954 pr_qm_opts(fp, indent, "qm-opts", &(ir->opts));
961 static void pr_harm(FILE *fp, t_iparams *iparams, const char *r, const char *kr)
963 fprintf(fp, "%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
964 r, iparams->harmonic.rA, kr, iparams->harmonic.krA,
965 r, iparams->harmonic.rB, kr, iparams->harmonic.krB);
968 void pr_iparams(FILE *fp, t_functype ftype, t_iparams *iparams)
971 real VA[4], VB[4], *rbcA, *rbcB;
977 pr_harm(fp, iparams, "th", "ct");
979 case F_CROSS_BOND_BONDS:
980 fprintf(fp, "r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
981 iparams->cross_bb.r1e, iparams->cross_bb.r2e,
982 iparams->cross_bb.krr);
984 case F_CROSS_BOND_ANGLES:
985 fprintf(fp, "r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
986 iparams->cross_ba.r1e, iparams->cross_ba.r2e,
987 iparams->cross_ba.r3e, iparams->cross_ba.krt);
989 case F_LINEAR_ANGLES:
990 fprintf(fp, "klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
991 iparams->linangle.klinA, iparams->linangle.aA,
992 iparams->linangle.klinB, iparams->linangle.aB);
995 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);
997 case F_QUARTIC_ANGLES:
998 fprintf(fp, "theta=%15.8e", iparams->qangle.theta);
999 for (i = 0; i < 5; i++)
1001 fprintf(fp, ", c%c=%15.8e", '0'+i, iparams->qangle.c[i]);
1006 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1007 iparams->bham.a, iparams->bham.b, iparams->bham.c);
1012 pr_harm(fp, iparams, "b0", "cb");
1015 pr_harm(fp, iparams, "xi", "cx");
1018 fprintf(fp, "b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
1019 iparams->morse.b0A, iparams->morse.cbA, iparams->morse.betaA,
1020 iparams->morse.b0B, iparams->morse.cbB, iparams->morse.betaB);
1023 fprintf(fp, "b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
1024 iparams->cubic.b0, iparams->cubic.kb, iparams->cubic.kcub);
1030 fprintf(fp, "bm=%15.8e, kb=%15.8e\n", iparams->fene.bm, iparams->fene.kb);
1033 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",
1034 iparams->restraint.lowA, iparams->restraint.up1A,
1035 iparams->restraint.up2A, iparams->restraint.kA,
1036 iparams->restraint.lowB, iparams->restraint.up1B,
1037 iparams->restraint.up2B, iparams->restraint.kB);
1043 fprintf(fp, "tab=%d, kA=%15.8e, kB=%15.8e\n",
1044 iparams->tab.table, iparams->tab.kA, iparams->tab.kB);
1046 case F_POLARIZATION:
1047 fprintf(fp, "alpha=%15.8e\n", iparams->polarize.alpha);
1050 fprintf(fp, "alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
1051 iparams->anharm_polarize.alpha,
1052 iparams->anharm_polarize.drcut,
1053 iparams->anharm_polarize.khyp);
1056 fprintf(fp, "a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
1057 iparams->thole.a, iparams->thole.alpha1, iparams->thole.alpha2,
1058 iparams->thole.rfac);
1061 fprintf(fp, "al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
1062 iparams->wpol.al_x, iparams->wpol.al_y, iparams->wpol.al_z,
1063 iparams->wpol.rOH, iparams->wpol.rHH, iparams->wpol.rOD);
1066 fprintf(fp, "c6=%15.8e, c12=%15.8e\n", iparams->lj.c6, iparams->lj.c12);
1069 fprintf(fp, "c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
1070 iparams->lj14.c6A, iparams->lj14.c12A,
1071 iparams->lj14.c6B, iparams->lj14.c12B);
1074 fprintf(fp, "fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1076 iparams->ljc14.qi, iparams->ljc14.qj,
1077 iparams->ljc14.c6, iparams->ljc14.c12);
1079 case F_LJC_PAIRS_NB:
1080 fprintf(fp, "qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1081 iparams->ljcnb.qi, iparams->ljcnb.qj,
1082 iparams->ljcnb.c6, iparams->ljcnb.c12);
1088 fprintf(fp, "phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
1089 iparams->pdihs.phiA, iparams->pdihs.cpA,
1090 iparams->pdihs.phiB, iparams->pdihs.cpB,
1091 iparams->pdihs.mult);
1094 fprintf(fp, "label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
1095 iparams->disres.label, iparams->disres.type,
1096 iparams->disres.low, iparams->disres.up1,
1097 iparams->disres.up2, iparams->disres.kfac);
1100 fprintf(fp, "ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
1101 iparams->orires.ex, iparams->orires.label, iparams->orires.power,
1102 iparams->orires.c, iparams->orires.obs, iparams->orires.kfac);
1105 fprintf(fp, "phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
1106 iparams->dihres.phiA, iparams->dihres.dphiA, iparams->dihres.kfacA,
1107 iparams->dihres.phiB, iparams->dihres.dphiB, iparams->dihres.kfacB);
1110 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",
1111 iparams->posres.pos0A[XX], iparams->posres.pos0A[YY],
1112 iparams->posres.pos0A[ZZ], iparams->posres.fcA[XX],
1113 iparams->posres.fcA[YY], iparams->posres.fcA[ZZ],
1114 iparams->posres.pos0B[XX], iparams->posres.pos0B[YY],
1115 iparams->posres.pos0B[ZZ], iparams->posres.fcB[XX],
1116 iparams->posres.fcB[YY], iparams->posres.fcB[ZZ]);
1119 fprintf(fp, "pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e\n",
1120 iparams->fbposres.pos0[XX], iparams->fbposres.pos0[YY],
1121 iparams->fbposres.pos0[ZZ], iparams->fbposres.geom,
1122 iparams->fbposres.r, iparams->fbposres.k);
1125 for (i = 0; i < NR_RBDIHS; i++)
1127 fprintf(fp, "%srbcA[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcA[i]);
1130 for (i = 0; i < NR_RBDIHS; i++)
1132 fprintf(fp, "%srbcB[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcB[i]);
1137 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
1138 * OPLS potential constants back.
1140 rbcA = iparams->rbdihs.rbcA;
1141 rbcB = iparams->rbdihs.rbcB;
1143 VA[3] = -0.25*rbcA[4];
1144 VA[2] = -0.5*rbcA[3];
1145 VA[1] = 4.0*VA[3]-rbcA[2];
1146 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1148 VB[3] = -0.25*rbcB[4];
1149 VB[2] = -0.5*rbcB[3];
1150 VB[1] = 4.0*VB[3]-rbcB[2];
1151 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1153 for (i = 0; i < NR_FOURDIHS; i++)
1155 fprintf(fp, "%sFourA[%d]=%15.8e", i == 0 ? "" : ", ", i, VA[i]);
1158 for (i = 0; i < NR_FOURDIHS; i++)
1160 fprintf(fp, "%sFourB[%d]=%15.8e", i == 0 ? "" : ", ", i, VB[i]);
1167 fprintf(fp, "dA=%15.8e, dB=%15.8e\n", iparams->constr.dA, iparams->constr.dB);
1170 fprintf(fp, "doh=%15.8e, dhh=%15.8e\n", iparams->settle.doh,
1171 iparams->settle.dhh);
1174 fprintf(fp, "a=%15.8e\n", iparams->vsite.a);
1179 fprintf(fp, "a=%15.8e, b=%15.8e\n", iparams->vsite.a, iparams->vsite.b);
1184 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1185 iparams->vsite.a, iparams->vsite.b, iparams->vsite.c);
1188 fprintf(fp, "n=%2d, a=%15.8e\n", iparams->vsiten.n, iparams->vsiten.a);
1193 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);
1196 fprintf(fp, "cmapA=%1d, cmapB=%1d\n", iparams->cmap.cmapA, iparams->cmap.cmapB);
1199 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1200 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1204 void pr_ilist(FILE *fp, int indent, const char *title,
1205 t_functype *functype, t_ilist *ilist, gmx_bool bShowNumbers)
1207 int i, j, k, type, ftype;
1210 if (available(fp, ilist, indent, title) && ilist->nr > 0)
1212 indent = pr_title(fp, indent, title);
1213 (void) pr_indent(fp, indent);
1214 fprintf(fp, "nr: %d\n", ilist->nr);
1217 (void) pr_indent(fp, indent);
1218 fprintf(fp, "iatoms:\n");
1219 iatoms = ilist->iatoms;
1220 for (i = j = 0; i < ilist->nr; )
1223 (void) pr_indent(fp, indent+INDENT);
1225 ftype = functype[type];
1226 (void) fprintf(fp, "%d type=%d (%s)",
1227 bShowNumbers ? j : -1, bShowNumbers ? type : -1,
1228 interaction_function[ftype].name);
1230 for (k = 0; k < interaction_function[ftype].nratoms; k++)
1232 (void) fprintf(fp, " %u", *(iatoms++));
1234 (void) fprintf(fp, "\n");
1235 i += 1+interaction_function[ftype].nratoms;
1237 fprintf(fp, "%5d%5d\n", i, iatoms[i]);
1245 static void pr_cmap(FILE *fp, int indent, const char *title,
1246 gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1251 dx = 360.0 / cmap_grid->grid_spacing;
1252 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1254 if (available(fp, cmap_grid, indent, title))
1256 fprintf(fp, "%s\n", title);
1258 for (i = 0; i < cmap_grid->ngrid; i++)
1261 fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
1263 fprintf(fp, "grid[%3d]={\n", bShowNumbers ? i : -1);
1265 for (j = 0; j < nelem; j++)
1267 if ( (j%cmap_grid->grid_spacing) == 0)
1269 fprintf(fp, "%8.1f\n", idx);
1273 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4]);
1274 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+1]);
1275 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+2]);
1276 fprintf(fp, "%8.3f\n", cmap_grid->cmapdata[i].cmap[j*4+3]);
1284 void pr_ffparams(FILE *fp, int indent, const char *title,
1285 gmx_ffparams_t *ffparams,
1286 gmx_bool bShowNumbers)
1290 indent = pr_title(fp, indent, title);
1291 (void) pr_indent(fp, indent);
1292 (void) fprintf(fp, "atnr=%d\n", ffparams->atnr);
1293 (void) pr_indent(fp, indent);
1294 (void) fprintf(fp, "ntypes=%d\n", ffparams->ntypes);
1295 for (i = 0; i < ffparams->ntypes; i++)
1297 (void) pr_indent(fp, indent+INDENT);
1298 (void) fprintf(fp, "functype[%d]=%s, ",
1299 bShowNumbers ? i : -1,
1300 interaction_function[ffparams->functype[i]].name);
1301 pr_iparams(fp, ffparams->functype[i], &ffparams->iparams[i]);
1303 (void) pr_double(fp, indent, "reppow", ffparams->reppow);
1304 (void) pr_real(fp, indent, "fudgeQQ", ffparams->fudgeQQ);
1305 pr_cmap(fp, indent, "cmap", &ffparams->cmap_grid, bShowNumbers);
1308 void pr_idef(FILE *fp, int indent, const char *title, t_idef *idef, gmx_bool bShowNumbers)
1312 if (available(fp, idef, indent, title))
1314 indent = pr_title(fp, indent, title);
1315 (void) pr_indent(fp, indent);
1316 (void) fprintf(fp, "atnr=%d\n", idef->atnr);
1317 (void) pr_indent(fp, indent);
1318 (void) fprintf(fp, "ntypes=%d\n", idef->ntypes);
1319 for (i = 0; i < idef->ntypes; i++)
1321 (void) pr_indent(fp, indent+INDENT);
1322 (void) fprintf(fp, "functype[%d]=%s, ",
1323 bShowNumbers ? i : -1,
1324 interaction_function[idef->functype[i]].name);
1325 pr_iparams(fp, idef->functype[i], &idef->iparams[i]);
1327 (void) pr_real(fp, indent, "fudgeQQ", idef->fudgeQQ);
1329 for (j = 0; (j < F_NRE); j++)
1331 pr_ilist(fp, indent, interaction_function[j].longname,
1332 idef->functype, &idef->il[j], bShowNumbers);
1337 static int pr_block_title(FILE *fp, int indent, const char *title, t_block *block)
1341 if (available(fp, block, indent, title))
1343 indent = pr_title(fp, indent, title);
1344 (void) pr_indent(fp, indent);
1345 (void) fprintf(fp, "nr=%d\n", block->nr);
1350 static int pr_blocka_title(FILE *fp, int indent, const char *title, t_blocka *block)
1354 if (available(fp, block, indent, title))
1356 indent = pr_title(fp, indent, title);
1357 (void) pr_indent(fp, indent);
1358 (void) fprintf(fp, "nr=%d\n", block->nr);
1359 (void) pr_indent(fp, indent);
1360 (void) fprintf(fp, "nra=%d\n", block->nra);
1365 static void low_pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1369 if (available(fp, block, indent, title))
1371 indent = pr_blocka_title(fp, indent, title, block);
1372 for (i = 0; i <= block->nr; i++)
1374 (void) pr_indent(fp, indent+INDENT);
1375 (void) fprintf(fp, "%s->index[%d]=%u\n",
1376 title, bShowNumbers ? i : -1, block->index[i]);
1378 for (i = 0; i < block->nra; i++)
1380 (void) pr_indent(fp, indent+INDENT);
1381 (void) fprintf(fp, "%s->a[%d]=%u\n",
1382 title, bShowNumbers ? i : -1, block->a[i]);
1387 void pr_block(FILE *fp, int indent, const char *title, t_block *block, gmx_bool bShowNumbers)
1389 int i, j, ok, size, start, end;
1391 if (available(fp, block, indent, title))
1393 indent = pr_block_title(fp, indent, title, block);
1396 if ((ok = (block->index[start] == 0)) == 0)
1398 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1402 for (i = 0; i < block->nr; i++)
1404 end = block->index[i+1];
1405 size = pr_indent(fp, indent);
1408 size += fprintf(fp, "%s[%d]={}\n", title, i);
1412 size += fprintf(fp, "%s[%d]={%d..%d}\n",
1413 title, bShowNumbers ? i : -1,
1414 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1422 void pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1424 int i, j, ok, size, start, end;
1426 if (available(fp, block, indent, title))
1428 indent = pr_blocka_title(fp, indent, title, block);
1431 if ((ok = (block->index[start] == 0)) == 0)
1433 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1437 for (i = 0; i < block->nr; i++)
1439 end = block->index[i+1];
1440 size = pr_indent(fp, indent);
1443 size += fprintf(fp, "%s[%d]={", title, i);
1447 size += fprintf(fp, "%s[%d][%d..%d]={",
1448 title, bShowNumbers ? i : -1,
1449 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1451 for (j = start; j < end; j++)
1455 size += fprintf(fp, ", ");
1457 if ((size) > (USE_WIDTH))
1459 (void) fprintf(fp, "\n");
1460 size = pr_indent(fp, indent+INDENT);
1462 size += fprintf(fp, "%u", block->a[j]);
1464 (void) fprintf(fp, "}\n");
1468 if ((end != block->nra) || (!ok))
1470 (void) pr_indent(fp, indent);
1471 (void) fprintf(fp, "tables inconsistent, dumping complete tables:\n");
1472 low_pr_blocka(fp, indent, title, block, bShowNumbers);
1477 static void pr_strings(FILE *fp, int indent, const char *title, char ***nm, int n, gmx_bool bShowNumbers)
1481 if (available(fp, nm, indent, title))
1483 indent = pr_title_n(fp, indent, title, n);
1484 for (i = 0; i < n; i++)
1486 (void) pr_indent(fp, indent);
1487 (void) fprintf(fp, "%s[%d]={name=\"%s\"}\n",
1488 title, bShowNumbers ? i : -1, *(nm[i]));
1493 static void pr_strings2(FILE *fp, int indent, const char *title,
1494 char ***nm, char ***nmB, int n, gmx_bool bShowNumbers)
1498 if (available(fp, nm, indent, title))
1500 indent = pr_title_n(fp, indent, title, n);
1501 for (i = 0; i < n; i++)
1503 (void) pr_indent(fp, indent);
1504 (void) fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1505 title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
1510 static void pr_resinfo(FILE *fp, int indent, const char *title, t_resinfo *resinfo, int n, gmx_bool bShowNumbers)
1514 if (available(fp, resinfo, indent, title))
1516 indent = pr_title_n(fp, indent, title, n);
1517 for (i = 0; i < n; i++)
1519 (void) pr_indent(fp, indent);
1520 (void) fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1521 title, bShowNumbers ? i : -1,
1522 *(resinfo[i].name), resinfo[i].nr,
1523 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1528 static void pr_atom(FILE *fp, int indent, const char *title, t_atom *atom, int n)
1532 if (available(fp, atom, indent, title))
1534 indent = pr_title_n(fp, indent, title, n);
1535 for (i = 0; i < n; i++)
1537 (void) pr_indent(fp, indent);
1538 fprintf(fp, "%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1539 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1540 title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype],
1541 atom[i].m, atom[i].q, atom[i].mB, atom[i].qB,
1542 atom[i].resind, atom[i].atomnumber);
1547 static void pr_grps(FILE *fp, const char *title, t_grps grps[], char **grpname[])
1551 for (i = 0; (i < egcNR); i++)
1553 fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
1554 for (j = 0; (j < grps[i].nr); j++)
1556 fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
1562 static void pr_groups(FILE *fp, int indent,
1563 gmx_groups_t *groups,
1564 gmx_bool bShowNumbers)
1569 pr_grps(fp, "grp", groups->grps, groups->grpname);
1570 pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
1572 (void) pr_indent(fp, indent);
1573 fprintf(fp, "groups ");
1574 for (g = 0; g < egcNR; g++)
1576 printf(" %5.5s", gtypes[g]);
1580 (void) pr_indent(fp, indent);
1581 fprintf(fp, "allocated ");
1583 for (g = 0; g < egcNR; g++)
1585 printf(" %5d", groups->ngrpnr[g]);
1586 nat_max = max(nat_max, groups->ngrpnr[g]);
1592 (void) pr_indent(fp, indent);
1593 fprintf(fp, "groupnr[%5s] =", "*");
1594 for (g = 0; g < egcNR; g++)
1596 fprintf(fp, " %3d ", 0);
1602 for (i = 0; i < nat_max; i++)
1604 (void) pr_indent(fp, indent);
1605 fprintf(fp, "groupnr[%5d] =", i);
1606 for (g = 0; g < egcNR; g++)
1608 fprintf(fp, " %3d ",
1609 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1616 void pr_atoms(FILE *fp, int indent, const char *title, t_atoms *atoms,
1617 gmx_bool bShownumbers)
1619 if (available(fp, atoms, indent, title))
1621 indent = pr_title(fp, indent, title);
1622 pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
1623 pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
1624 pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
1625 pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
1630 void pr_atomtypes(FILE *fp, int indent, const char *title, t_atomtypes *atomtypes,
1631 gmx_bool bShowNumbers)
1634 if (available(fp, atomtypes, indent, title))
1636 indent = pr_title(fp, indent, title);
1637 for (i = 0; i < atomtypes->nr; i++)
1639 pr_indent(fp, indent);
1641 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1642 bShowNumbers ? i : -1, atomtypes->radius[i], atomtypes->vol[i],
1643 atomtypes->gb_radius[i],
1644 atomtypes->surftens[i], atomtypes->atomnumber[i], atomtypes->S_hct[i]);
1649 static void pr_moltype(FILE *fp, int indent, const char *title,
1650 gmx_moltype_t *molt, int n,
1651 gmx_ffparams_t *ffparams,
1652 gmx_bool bShowNumbers)
1656 indent = pr_title_n(fp, indent, title, n);
1657 (void) pr_indent(fp, indent);
1658 (void) fprintf(fp, "name=\"%s\"\n", *(molt->name));
1659 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
1660 pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
1661 pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
1662 for (j = 0; (j < F_NRE); j++)
1664 pr_ilist(fp, indent, interaction_function[j].longname,
1665 ffparams->functype, &molt->ilist[j], bShowNumbers);
1669 static void pr_molblock(FILE *fp, int indent, const char *title,
1670 gmx_molblock_t *molb, int n,
1671 gmx_moltype_t *molt)
1673 indent = pr_title_n(fp, indent, title, n);
1674 (void) pr_indent(fp, indent);
1675 (void) fprintf(fp, "%-20s = %d \"%s\"\n",
1676 "moltype", molb->type, *(molt[molb->type].name));
1677 pr_int(fp, indent, "#molecules", molb->nmol);
1678 pr_int(fp, indent, "#atoms_mol", molb->natoms_mol);
1679 pr_int(fp, indent, "#posres_xA", molb->nposres_xA);
1680 if (molb->nposres_xA > 0)
1682 pr_rvecs(fp, indent, "posres_xA", molb->posres_xA, molb->nposres_xA);
1684 pr_int(fp, indent, "#posres_xB", molb->nposres_xB);
1685 if (molb->nposres_xB > 0)
1687 pr_rvecs(fp, indent, "posres_xB", molb->posres_xB, molb->nposres_xB);
1691 void pr_mtop(FILE *fp, int indent, const char *title, gmx_mtop_t *mtop,
1692 gmx_bool bShowNumbers)
1696 if (available(fp, mtop, indent, title))
1698 indent = pr_title(fp, indent, title);
1699 (void) pr_indent(fp, indent);
1700 (void) fprintf(fp, "name=\"%s\"\n", *(mtop->name));
1701 pr_int(fp, indent, "#atoms", mtop->natoms);
1702 pr_int(fp, indent, "#molblock", mtop->nmolblock);
1703 for (mb = 0; mb < mtop->nmolblock; mb++)
1705 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
1707 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
1708 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
1709 for (mt = 0; mt < mtop->nmoltype; mt++)
1711 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
1712 &mtop->ffparams, bShowNumbers);
1714 pr_groups(fp, indent, &mtop->groups, bShowNumbers);
1718 void pr_top(FILE *fp, int indent, const char *title, t_topology *top, gmx_bool bShowNumbers)
1720 if (available(fp, top, indent, title))
1722 indent = pr_title(fp, indent, title);
1723 (void) pr_indent(fp, indent);
1724 (void) fprintf(fp, "name=\"%s\"\n", *(top->name));
1725 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
1726 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
1727 pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
1728 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
1729 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
1730 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers);
1734 void pr_header(FILE *fp, int indent, const char *title, t_tpxheader *sh)
1738 if (available(fp, sh, indent, title))
1740 indent = pr_title(fp, indent, title);
1741 pr_indent(fp, indent);
1742 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
1743 pr_indent(fp, indent);
1744 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
1745 pr_indent(fp, indent);
1746 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
1747 pr_indent(fp, indent);
1748 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
1749 pr_indent(fp, indent);
1750 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
1751 pr_indent(fp, indent);
1752 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
1754 pr_indent(fp, indent);
1755 fprintf(fp, "natoms = %d\n", sh->natoms);
1756 pr_indent(fp, indent);
1757 fprintf(fp, "lambda = %e\n", sh->lambda);
1761 void pr_commrec(FILE *fp, int indent, t_commrec *cr)
1763 pr_indent(fp, indent);
1764 fprintf(fp, "commrec:\n");
1766 pr_indent(fp, indent);
1767 fprintf(fp, "nodeid = %d\n", cr->nodeid);
1768 pr_indent(fp, indent);
1769 fprintf(fp, "nnodes = %d\n", cr->nnodes);
1770 pr_indent(fp, indent);
1771 fprintf(fp, "npmenodes = %d\n", cr->npmenodes);
1773 pr_indent(fp,indent);
1774 fprintf(fp,"threadid = %d\n",cr->threadid);
1775 pr_indent(fp,indent);
1776 fprintf(fp,"nthreads = %d\n",cr->nthreads);