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,2014, 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! */
46 #include "types/commrec.h"
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_rvecs_of_dim(FILE *fp, int indent, const char *title, rvec vec[], int n, int dim)
303 const char *fshort = "%12.5e";
304 const char *flong = "%15.8e";
308 if (getenv("LONGFORMAT") != NULL)
317 if (available(fp, vec, indent, title))
319 indent = pr_title_nxn(fp, indent, title, n, dim);
320 for (i = 0; i < n; i++)
322 (void) pr_indent(fp, indent);
323 (void) fprintf(fp, "%s[%5d]={", title, i);
324 for (j = 0; j < dim; j++)
328 (void) fprintf(fp, ", ");
330 (void) fprintf(fp, format, vec[i][j]);
332 (void) fprintf(fp, "}\n");
337 void pr_reals(FILE *fp, int indent, const char *title, real *vec, int n)
341 if (available(fp, vec, indent, title))
343 (void) pr_indent(fp, indent);
344 (void) fprintf(fp, "%s:\t", title);
345 for (i = 0; i < n; i++)
347 fprintf(fp, " %10g", vec[i]);
349 (void) fprintf(fp, "\n");
353 void pr_doubles(FILE *fp, int indent, const char *title, double *vec, int n)
357 if (available(fp, vec, indent, title))
359 (void) pr_indent(fp, indent);
360 (void) fprintf(fp, "%s:\t", title);
361 for (i = 0; i < n; i++)
363 fprintf(fp, " %10g", vec[i]);
365 (void) fprintf(fp, "\n");
369 void pr_reals_of_dim(FILE *fp, int indent, const char *title, real *vec, int n, int dim)
372 const char *fshort = "%12.5e";
373 const char *flong = "%15.8e";
376 if (getenv("LONGFORMAT") != NULL)
385 if (available(fp, vec, indent, title))
387 indent = pr_title_nxn(fp, indent, title, n, dim);
388 for (i = 0; i < n; i++)
390 (void) pr_indent(fp, indent);
391 (void) fprintf(fp, "%s[%5d]={", title, i);
392 for (j = 0; j < dim; j++)
396 (void) fprintf(fp, ", ");
398 (void) fprintf(fp, format, vec[i * dim + j]);
400 (void) fprintf(fp, "}\n");
405 static void pr_int(FILE *fp, int indent, const char *title, int i)
407 pr_indent(fp, indent);
408 fprintf(fp, "%-20s = %d\n", title, i);
411 static void pr_int64(FILE *fp, int indent, const char *title, gmx_int64_t i)
413 char buf[STEPSTRSIZE];
415 pr_indent(fp, indent);
416 fprintf(fp, "%-20s = %s\n", title, gmx_step_str(i, buf));
419 static void pr_real(FILE *fp, int indent, const char *title, real r)
421 pr_indent(fp, indent);
422 fprintf(fp, "%-20s = %g\n", title, r);
425 static void pr_double(FILE *fp, int indent, const char *title, double d)
427 pr_indent(fp, indent);
428 fprintf(fp, "%-20s = %g\n", title, d);
431 static void pr_str(FILE *fp, int indent, const char *title, const char *s)
433 pr_indent(fp, indent);
434 fprintf(fp, "%-30s = %s\n", title, s);
437 void pr_qm_opts(FILE *fp, int indent, const char *title, t_grpopts *opts)
441 fprintf(fp, "%s:\n", title);
443 pr_int(fp, indent, "ngQM", opts->ngQM);
446 pr_ivec(fp, indent, "QMmethod", opts->QMmethod, opts->ngQM, FALSE);
447 pr_ivec(fp, indent, "QMbasis", opts->QMbasis, opts->ngQM, FALSE);
448 pr_ivec(fp, indent, "QMcharge", opts->QMcharge, opts->ngQM, FALSE);
449 pr_ivec(fp, indent, "QMmult", opts->QMmult, opts->ngQM, FALSE);
450 pr_bvec(fp, indent, "bSH", opts->bSH, opts->ngQM, FALSE);
451 pr_ivec(fp, indent, "CASorbitals", opts->CASorbitals, opts->ngQM, FALSE);
452 pr_ivec(fp, indent, "CASelectrons", opts->CASelectrons, opts->ngQM, FALSE);
453 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
454 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
455 pr_ivec(fp, indent, "SAsteps", opts->SAsteps, opts->ngQM, FALSE);
456 pr_bvec(fp, indent, "bOPT", opts->bOPT, opts->ngQM, FALSE);
457 pr_bvec(fp, indent, "bTS", opts->bTS, opts->ngQM, FALSE);
461 static void pr_grp_opts(FILE *out, int indent, const char *title, t_grpopts *opts,
468 fprintf(out, "%s:\n", title);
471 pr_indent(out, indent);
472 fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
473 for (i = 0; (i < opts->ngtc); i++)
475 fprintf(out, " %10g", opts->nrdf[i]);
479 pr_indent(out, indent);
480 fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
481 for (i = 0; (i < opts->ngtc); i++)
483 fprintf(out, " %10g", opts->ref_t[i]);
487 pr_indent(out, indent);
488 fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
489 for (i = 0; (i < opts->ngtc); i++)
491 fprintf(out, " %10g", opts->tau_t[i]);
495 /* Pretty-print the simulated annealing info */
496 fprintf(out, "anneal%s", bMDPformat ? " = " : ":");
497 for (i = 0; (i < opts->ngtc); i++)
499 fprintf(out, " %10s", EANNEAL(opts->annealing[i]));
503 fprintf(out, "ann-npoints%s", bMDPformat ? " = " : ":");
504 for (i = 0; (i < opts->ngtc); i++)
506 fprintf(out, " %10d", opts->anneal_npoints[i]);
510 for (i = 0; (i < opts->ngtc); i++)
512 if (opts->anneal_npoints[i] > 0)
514 fprintf(out, "ann. times [%d]:\t", i);
515 for (j = 0; (j < opts->anneal_npoints[i]); j++)
517 fprintf(out, " %10.1f", opts->anneal_time[i][j]);
520 fprintf(out, "ann. temps [%d]:\t", i);
521 for (j = 0; (j < opts->anneal_npoints[i]); j++)
523 fprintf(out, " %10.1f", opts->anneal_temp[i][j]);
529 pr_indent(out, indent);
530 fprintf(out, "acc:\t");
531 for (i = 0; (i < opts->ngacc); i++)
533 for (m = 0; (m < DIM); m++)
535 fprintf(out, " %10g", opts->acc[i][m]);
540 pr_indent(out, indent);
541 fprintf(out, "nfreeze:");
542 for (i = 0; (i < opts->ngfrz); i++)
544 for (m = 0; (m < DIM); m++)
546 fprintf(out, " %10s", opts->nFreeze[i][m] ? "Y" : "N");
552 for (i = 0; (i < opts->ngener); i++)
554 pr_indent(out, indent);
555 fprintf(out, "energygrp-flags[%3d]:", i);
556 for (m = 0; (m < opts->ngener); m++)
558 fprintf(out, " %d", opts->egp_flags[opts->ngener*i+m]);
566 static void pr_matrix(FILE *fp, int indent, const char *title, rvec *m,
571 fprintf(fp, "%-10s = %g %g %g %g %g %g\n", title,
572 m[XX][XX], m[YY][YY], m[ZZ][ZZ], m[XX][YY], m[XX][ZZ], m[YY][ZZ]);
576 pr_rvecs(fp, indent, title, m, DIM);
580 static void pr_cosine(FILE *fp, int indent, const char *title, t_cosines *cos,
587 fprintf(fp, "%s = %d\n", title, cos->n);
591 indent = pr_title(fp, indent, title);
592 (void) pr_indent(fp, indent);
593 fprintf(fp, "n = %d\n", cos->n);
596 (void) pr_indent(fp, indent+2);
598 for (j = 0; (j < cos->n); j++)
600 fprintf(fp, " %e", cos->a[j]);
603 (void) pr_indent(fp, indent+2);
604 fprintf(fp, "phi =");
605 for (j = 0; (j < cos->n); j++)
607 fprintf(fp, " %e", cos->phi[j]);
614 #define PS(t, s) pr_str(fp, indent, t, s)
615 #define PI(t, s) pr_int(fp, indent, t, s)
616 #define PSTEP(t, s) pr_int64(fp, indent, t, s)
617 #define PR(t, s) pr_real(fp, indent, t, s)
618 #define PD(t, s) pr_double(fp, indent, t, s)
620 static void pr_pull_group(FILE *fp, int indent, int g, t_pull_group *pgrp)
622 pr_indent(fp, indent);
623 fprintf(fp, "pull-group %d:\n", g);
625 pr_ivec_block(fp, indent, "atom", pgrp->ind, pgrp->nat, TRUE);
626 pr_rvec(fp, indent, "weight", pgrp->weight, pgrp->nweight, TRUE);
627 PI("pbcatom", pgrp->pbcatom);
630 static void pr_pull_coord(FILE *fp, int indent, int c, t_pull_coord *pcrd)
632 pr_indent(fp, indent);
633 fprintf(fp, "pull-coord %d:\n", c);
634 PI("group[0]", pcrd->group[0]);
635 PI("group[1]", pcrd->group[1]);
636 pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
637 pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
638 PR("init", pcrd->init);
639 PR("rate", pcrd->rate);
644 static void pr_simtempvals(FILE *fp, int indent, t_simtemp *simtemp, int n_lambda)
646 PR("simtemp_low", simtemp->simtemp_low);
647 PR("simtemp_high", simtemp->simtemp_high);
648 PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
649 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
652 static void pr_expandedvals(FILE *fp, int indent, t_expanded *expand, int n_lambda)
655 PI("nstexpanded", expand->nstexpanded);
656 PS("lambda-stats", elamstats_names[expand->elamstats]);
657 PS("lambda-mc-move", elmcmove_names[expand->elmcmove]);
658 PI("lmc-repeats", expand->lmc_repeats);
659 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
660 PI("lmc-nstart", expand->lmc_forced_nstart);
661 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
662 PI("nst-transition-matrix", expand->nstTij);
663 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
664 PI("weight-c-range", expand->c_range); /* default is just C=0 */
665 PR("wl-scale", expand->wl_scale);
666 PR("init-wl-delta", expand->init_wl_delta);
667 PR("wl-ratio", expand->wl_ratio);
668 PS("bWLoneovert", EBOOL(expand->bWLoneovert));
669 PI("lmc-seed", expand->lmc_seed);
670 PR("mc-temperature", expand->mc_temp);
671 PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
672 if (expand->elmceq == elmceqNUMATLAM)
674 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
676 if (expand->elmceq == elmceqSAMPLES)
678 PI("weight-equil-number-samples", expand->equil_samples);
680 if (expand->elmceq == elmceqSTEPS)
682 PI("weight-equil-number-steps", expand->equil_steps);
684 if (expand->elmceq == elmceqWLDELTA)
686 PR("weight-equil-wl-delta", expand->equil_wl_delta);
688 if (expand->elmceq == elmceqRATIO)
690 PR("weight-equil-count-ratio", expand->equil_ratio);
693 pr_indent(fp, indent);
694 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
695 PS("init-weights", EBOOL(expand->bInit_weights));
698 static void pr_fepvals(FILE *fp, int indent, t_lambda *fep, gmx_bool bMDPformat)
702 PI("nstdhdl", fep->nstdhdl);
703 PI("init-lambda-state", fep->init_fep_state);
704 PR("init-lambda", fep->init_lambda);
705 PR("delta-lambda", fep->delta_lambda);
708 PI("n-lambdas", fep->n_lambda);
710 if (fep->n_lambda > 0)
712 pr_indent(fp, indent);
713 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
714 for (i = 0; i < efptNR; i++)
716 fprintf(fp, "%18s = ", efpt_names[i]);
717 if (fep->separate_dvdl[i])
719 fprintf(fp, " TRUE");
723 fprintf(fp, " FALSE");
727 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
728 for (i = 0; i < efptNR; i++)
730 fprintf(fp, "%18s = ", efpt_names[i]);
731 for (j = 0; j < fep->n_lambda; j++)
733 fprintf(fp, " %10g", fep->all_lambda[i][j]);
738 PI("calc-lambda-neighbors", fep->lambda_neighbors);
740 PR("sc-alpha", fep->sc_alpha);
741 PS("bScCoul", EBOOL(fep->bScCoul));
742 PS("bScPrintEnergy", EBOOL(fep->bPrintEnergy));
743 PI("sc-power", fep->sc_power);
744 PR("sc-r-power", fep->sc_r_power);
745 PR("sc-sigma", fep->sc_sigma);
746 PR("sc-sigma-min", fep->sc_sigma_min);
747 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
748 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
749 PI("dh-hist-size", fep->dh_hist_size);
750 PD("dh-hist-spacing", fep->dh_hist_spacing);
753 static void pr_pull(FILE *fp, int indent, t_pull *pull)
757 PS("pull-geometry", EPULLGEOM(pull->eGeom));
758 pr_ivec(fp, indent, "pull-dim", pull->dim, DIM, TRUE);
759 PR("pull-r1", pull->cyl_r1);
760 PR("pull-r0", pull->cyl_r0);
761 PR("pull-constr-tol", pull->constr_tol);
762 PS("pull-bPrintRef", EBOOL(pull->bPrintRef));
763 PI("pull-nstxout", pull->nstxout);
764 PI("pull-nstfout", pull->nstfout);
765 PI("pull-ngroup", pull->ngroup);
766 for (g = 0; g < pull->ngroup; g++)
768 pr_pull_group(fp, indent, g, &pull->group[g]);
770 PI("pull-ncoord", pull->ncoord);
771 for (g = 0; g < pull->ncoord; g++)
773 pr_pull_coord(fp, indent, g, &pull->coord[g]);
777 static void pr_rotgrp(FILE *fp, int indent, int g, t_rotgrp *rotg)
779 pr_indent(fp, indent);
780 fprintf(fp, "rotation_group %d:\n", g);
782 PS("type", EROTGEOM(rotg->eType));
783 PS("massw", EBOOL(rotg->bMassW));
784 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
785 pr_rvecs(fp, indent, "x_ref", rotg->x_ref, rotg->nat);
786 pr_rvec(fp, indent, "vec", rotg->vec, DIM, TRUE);
787 pr_rvec(fp, indent, "pivot", rotg->pivot, DIM, TRUE);
788 PR("rate", rotg->rate);
790 PR("slab_dist", rotg->slab_dist);
791 PR("min_gaussian", rotg->min_gaussian);
792 PR("epsilon", rotg->eps);
793 PS("fit_method", EROTFIT(rotg->eFittype));
794 PI("potfitangle_nstep", rotg->PotAngle_nstep);
795 PR("potfitangle_step", rotg->PotAngle_step);
798 static void pr_rot(FILE *fp, int indent, t_rot *rot)
802 PI("rot_nstrout", rot->nstrout);
803 PI("rot_nstsout", rot->nstsout);
804 PI("rot_ngrp", rot->ngrp);
805 for (g = 0; g < rot->ngrp; g++)
807 pr_rotgrp(fp, indent, g, &rot->grp[g]);
812 static void pr_swap(FILE *fp, int indent, t_swapcoords *swap)
818 PI("frequency", swap->nstswap);
819 for (j = 0; j < 2; j++)
821 sprintf(str, "nanions%c", j+'A');
822 PI(str, swap->nanions[j]);
823 sprintf(str, "ncations%c", j+'A');
824 PI(str, swap->ncations[j]);
826 PI("coupling_steps", swap->nAverage);
827 PR("threshold", swap->threshold);
828 for (j = 0; j < 2; j++)
830 sprintf(str, "splitgroup%d_massw", j);
831 PS(str, EBOOL(swap->massw_split[j]));
832 sprintf(str, "split atoms group %d", j);
833 pr_ivec_block(fp, indent, str, swap->ind_split[j], swap->nat_split[j], TRUE);
835 pr_ivec_block(fp, indent, "swap atoms", swap->ind, swap->nat, TRUE);
836 pr_ivec_block(fp, indent, "solvent atoms", swap->ind_sol, swap->nat_sol, TRUE);
837 PR("cyl0_radius", swap->cyl0r);
838 PR("cyl0_upper", swap->cyl0u);
839 PR("cyl0_lower", swap->cyl0l);
840 PR("cyl1_radius", swap->cyl1r);
841 PR("cyl1_upper", swap->cyl1u);
842 PR("cyl1_lower", swap->cyl1l);
846 static void pr_imd(FILE *fp, int indent, t_IMD *imd)
848 PI("IMD_atoms", imd->nat);
849 pr_ivec_block(fp, indent, "atom", imd->ind, imd->nat, TRUE);
853 void pr_inputrec(FILE *fp, int indent, const char *title, t_inputrec *ir,
856 const char *infbuf = "inf";
859 if (available(fp, ir, indent, title))
863 indent = pr_title(fp, indent, title);
865 /* This strings do not all have a direct correspondence to
866 .mdp entries, but we should follow the same convention of
867 using hyphens in the names users read and write. */
868 PS("integrator", EI(ir->eI));
869 PSTEP("nsteps", ir->nsteps);
870 PSTEP("init-step", ir->init_step);
871 PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
872 PS("ns-type", ENS(ir->ns_type));
873 PI("nstlist", ir->nstlist);
874 PI("ndelta", ir->ndelta);
875 PI("nstcomm", ir->nstcomm);
876 PS("comm-mode", ECOM(ir->comm_mode));
877 PI("nstlog", ir->nstlog);
878 PI("nstxout", ir->nstxout);
879 PI("nstvout", ir->nstvout);
880 PI("nstfout", ir->nstfout);
881 PI("nstcalcenergy", ir->nstcalcenergy);
882 PI("nstenergy", ir->nstenergy);
883 PI("nstxout-compressed", ir->nstxout_compressed);
884 PR("init-t", ir->init_t);
885 PR("delta-t", ir->delta_t);
887 PR("x-compression-precision", ir->x_compression_precision);
888 PR("fourierspacing", ir->fourier_spacing);
892 PI("pme-order", ir->pme_order);
893 PR("ewald-rtol", ir->ewald_rtol);
894 PR("ewald-rtol-lj", ir->ewald_rtol_lj);
895 PR("ewald-geometry", ir->ewald_geometry);
896 PR("epsilon-surface", ir->epsilon_surface);
897 PS("optimize-fft", EBOOL(ir->bOptFFT));
898 PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir->ljpme_combination_rule));
899 PS("ePBC", EPBC(ir->ePBC));
900 PS("bPeriodicMols", EBOOL(ir->bPeriodicMols));
901 PS("bContinuation", EBOOL(ir->bContinuation));
902 PS("bShakeSOR", EBOOL(ir->bShakeSOR));
903 PS("etc", ETCOUPLTYPE(ir->etc));
904 PS("bPrintNHChains", EBOOL(ir->bPrintNHChains));
905 PI("nsttcouple", ir->nsttcouple);
906 PS("epc", EPCOUPLTYPE(ir->epc));
907 PS("epctype", EPCOUPLTYPETYPE(ir->epct));
908 PI("nstpcouple", ir->nstpcouple);
909 PR("tau-p", ir->tau_p);
910 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
911 pr_matrix(fp, indent, "compress", ir->compress, bMDPformat);
912 PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
915 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX],
916 ir->posres_com[YY], ir->posres_com[ZZ]);
920 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
924 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX],
925 ir->posres_comB[YY], ir->posres_comB[ZZ]);
929 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
931 PR("verlet-buffer-tolerance", ir->verletbuf_tol);
932 PR("rlist", ir->rlist);
933 PR("rlistlong", ir->rlistlong);
934 PR("nstcalclr", ir->nstcalclr);
935 PR("rtpi", ir->rtpi);
936 PS("coulombtype", EELTYPE(ir->coulombtype));
937 PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
938 PR("rcoulomb-switch", ir->rcoulomb_switch);
939 PR("rcoulomb", ir->rcoulomb);
940 PS("vdwtype", EVDWTYPE(ir->vdwtype));
941 PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
942 PR("rvdw-switch", ir->rvdw_switch);
943 PR("rvdw", ir->rvdw);
944 if (ir->epsilon_r != 0)
946 PR("epsilon-r", ir->epsilon_r);
950 PS("epsilon-r", infbuf);
952 if (ir->epsilon_rf != 0)
954 PR("epsilon-rf", ir->epsilon_rf);
958 PS("epsilon-rf", infbuf);
960 PR("tabext", ir->tabext);
961 PS("implicit-solvent", EIMPLICITSOL(ir->implicit_solvent));
962 PS("gb-algorithm", EGBALGORITHM(ir->gb_algorithm));
963 PR("gb-epsilon-solvent", ir->gb_epsilon_solvent);
964 PI("nstgbradii", ir->nstgbradii);
965 PR("rgbradii", ir->rgbradii);
966 PR("gb-saltconc", ir->gb_saltconc);
967 PR("gb-obc-alpha", ir->gb_obc_alpha);
968 PR("gb-obc-beta", ir->gb_obc_beta);
969 PR("gb-obc-gamma", ir->gb_obc_gamma);
970 PR("gb-dielectric-offset", ir->gb_dielectric_offset);
971 PS("sa-algorithm", ESAALGORITHM(ir->gb_algorithm));
972 PR("sa-surface-tension", ir->sa_surface_tension);
973 PS("DispCorr", EDISPCORR(ir->eDispCorr));
974 PS("bSimTemp", EBOOL(ir->bSimTemp));
977 pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
979 PS("free-energy", EFEPTYPE(ir->efep));
980 if (ir->efep != efepNO || ir->bSimTemp)
982 pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
986 pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
989 PI("nwall", ir->nwall);
990 PS("wall-type", EWALLTYPE(ir->wall_type));
991 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
992 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
993 PR("wall-density[0]", ir->wall_density[0]);
994 PR("wall-density[1]", ir->wall_density[1]);
995 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
997 PS("pull", EPULLTYPE(ir->ePull));
998 if (ir->ePull != epullNO)
1000 pr_pull(fp, indent, ir->pull);
1003 PS("rotation", EBOOL(ir->bRot));
1006 pr_rot(fp, indent, ir->rot);
1009 PS("interactiveMD", EBOOL(ir->bIMD));
1012 pr_imd(fp, indent, ir->imd);
1015 PS("disre", EDISRETYPE(ir->eDisre));
1016 PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
1017 PS("disre-mixed", EBOOL(ir->bDisreMixed));
1018 PR("dr-fc", ir->dr_fc);
1019 PR("dr-tau", ir->dr_tau);
1020 PR("nstdisreout", ir->nstdisreout);
1021 PR("orires-fc", ir->orires_fc);
1022 PR("orires-tau", ir->orires_tau);
1023 PR("nstorireout", ir->nstorireout);
1025 PR("dihre-fc", ir->dihre_fc);
1027 PR("em-stepsize", ir->em_stepsize);
1028 PR("em-tol", ir->em_tol);
1029 PI("niter", ir->niter);
1030 PR("fc-stepsize", ir->fc_stepsize);
1031 PI("nstcgsteep", ir->nstcgsteep);
1032 PI("nbfgscorr", ir->nbfgscorr);
1034 PS("ConstAlg", ECONSTRTYPE(ir->eConstrAlg));
1035 PR("shake-tol", ir->shake_tol);
1036 PI("lincs-order", ir->nProjOrder);
1037 PR("lincs-warnangle", ir->LincsWarnAngle);
1038 PI("lincs-iter", ir->nLincsIter);
1039 PR("bd-fric", ir->bd_fric);
1040 PSTEP("ld-seed", ir->ld_seed);
1041 PR("cos-accel", ir->cos_accel);
1042 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
1044 PS("adress", EBOOL(ir->bAdress));
1047 PS("adress-type", EADRESSTYPE(ir->adress->type));
1048 PR("adress-const-wf", ir->adress->const_wf);
1049 PR("adress-ex-width", ir->adress->ex_width);
1050 PR("adress-hy-width", ir->adress->hy_width);
1051 PS("adress-interface-correction", EADRESSICTYPE(ir->adress->icor));
1052 PS("adress-site", EADRESSSITETYPE(ir->adress->site));
1053 PR("adress-ex-force-cap", ir->adress->ex_forcecap);
1054 PS("adress-do-hybridpairs", EBOOL(ir->adress->do_hybridpairs));
1056 pr_rvec(fp, indent, "adress-reference-coords", ir->adress->refs, DIM, TRUE);
1058 PI("userint1", ir->userint1);
1059 PI("userint2", ir->userint2);
1060 PI("userint3", ir->userint3);
1061 PI("userint4", ir->userint4);
1062 PR("userreal1", ir->userreal1);
1063 PR("userreal2", ir->userreal2);
1064 PR("userreal3", ir->userreal3);
1065 PR("userreal4", ir->userreal4);
1066 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1067 pr_cosine(fp, indent, "efield-x", &(ir->ex[XX]), bMDPformat);
1068 pr_cosine(fp, indent, "efield-xt", &(ir->et[XX]), bMDPformat);
1069 pr_cosine(fp, indent, "efield-y", &(ir->ex[YY]), bMDPformat);
1070 pr_cosine(fp, indent, "efield-yt", &(ir->et[YY]), bMDPformat);
1071 pr_cosine(fp, indent, "efield-z", &(ir->ex[ZZ]), bMDPformat);
1072 pr_cosine(fp, indent, "efield-zt", &(ir->et[ZZ]), bMDPformat);
1073 PS("eSwapCoords", ESWAPTYPE(ir->eSwapCoords));
1074 if (ir->eSwapCoords != eswapNO)
1076 pr_swap(fp, indent, ir->swap);
1078 PS("bQMMM", EBOOL(ir->bQMMM));
1079 PI("QMconstraints", ir->QMconstraints);
1080 PI("QMMMscheme", ir->QMMMscheme);
1081 PR("scalefactor", ir->scalefactor);
1082 pr_qm_opts(fp, indent, "qm-opts", &(ir->opts));
1089 static void pr_harm(FILE *fp, t_iparams *iparams, const char *r, const char *kr)
1091 fprintf(fp, "%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
1092 r, iparams->harmonic.rA, kr, iparams->harmonic.krA,
1093 r, iparams->harmonic.rB, kr, iparams->harmonic.krB);
1096 void pr_iparams(FILE *fp, t_functype ftype, t_iparams *iparams)
1099 real VA[4], VB[4], *rbcA, *rbcB;
1105 pr_harm(fp, iparams, "th", "ct");
1107 case F_CROSS_BOND_BONDS:
1108 fprintf(fp, "r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
1109 iparams->cross_bb.r1e, iparams->cross_bb.r2e,
1110 iparams->cross_bb.krr);
1112 case F_CROSS_BOND_ANGLES:
1113 fprintf(fp, "r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
1114 iparams->cross_ba.r1e, iparams->cross_ba.r2e,
1115 iparams->cross_ba.r3e, iparams->cross_ba.krt);
1117 case F_LINEAR_ANGLES:
1118 fprintf(fp, "klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
1119 iparams->linangle.klinA, iparams->linangle.aA,
1120 iparams->linangle.klinB, iparams->linangle.aB);
1122 case F_UREY_BRADLEY:
1123 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);
1125 case F_QUARTIC_ANGLES:
1126 fprintf(fp, "theta=%15.8e", iparams->qangle.theta);
1127 for (i = 0; i < 5; i++)
1129 fprintf(fp, ", c%c=%15.8e", '0'+i, iparams->qangle.c[i]);
1134 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1135 iparams->bham.a, iparams->bham.b, iparams->bham.c);
1140 pr_harm(fp, iparams, "b0", "cb");
1143 pr_harm(fp, iparams, "xi", "cx");
1146 fprintf(fp, "b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
1147 iparams->morse.b0A, iparams->morse.cbA, iparams->morse.betaA,
1148 iparams->morse.b0B, iparams->morse.cbB, iparams->morse.betaB);
1151 fprintf(fp, "b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
1152 iparams->cubic.b0, iparams->cubic.kb, iparams->cubic.kcub);
1158 fprintf(fp, "bm=%15.8e, kb=%15.8e\n", iparams->fene.bm, iparams->fene.kb);
1161 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",
1162 iparams->restraint.lowA, iparams->restraint.up1A,
1163 iparams->restraint.up2A, iparams->restraint.kA,
1164 iparams->restraint.lowB, iparams->restraint.up1B,
1165 iparams->restraint.up2B, iparams->restraint.kB);
1171 fprintf(fp, "tab=%d, kA=%15.8e, kB=%15.8e\n",
1172 iparams->tab.table, iparams->tab.kA, iparams->tab.kB);
1174 case F_POLARIZATION:
1175 fprintf(fp, "alpha=%15.8e\n", iparams->polarize.alpha);
1178 fprintf(fp, "alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
1179 iparams->anharm_polarize.alpha,
1180 iparams->anharm_polarize.drcut,
1181 iparams->anharm_polarize.khyp);
1184 fprintf(fp, "a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
1185 iparams->thole.a, iparams->thole.alpha1, iparams->thole.alpha2,
1186 iparams->thole.rfac);
1189 fprintf(fp, "al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
1190 iparams->wpol.al_x, iparams->wpol.al_y, iparams->wpol.al_z,
1191 iparams->wpol.rOH, iparams->wpol.rHH, iparams->wpol.rOD);
1194 fprintf(fp, "c6=%15.8e, c12=%15.8e\n", iparams->lj.c6, iparams->lj.c12);
1197 fprintf(fp, "c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
1198 iparams->lj14.c6A, iparams->lj14.c12A,
1199 iparams->lj14.c6B, iparams->lj14.c12B);
1202 fprintf(fp, "fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1204 iparams->ljc14.qi, iparams->ljc14.qj,
1205 iparams->ljc14.c6, iparams->ljc14.c12);
1207 case F_LJC_PAIRS_NB:
1208 fprintf(fp, "qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1209 iparams->ljcnb.qi, iparams->ljcnb.qj,
1210 iparams->ljcnb.c6, iparams->ljcnb.c12);
1216 fprintf(fp, "phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
1217 iparams->pdihs.phiA, iparams->pdihs.cpA,
1218 iparams->pdihs.phiB, iparams->pdihs.cpB,
1219 iparams->pdihs.mult);
1222 fprintf(fp, "label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
1223 iparams->disres.label, iparams->disres.type,
1224 iparams->disres.low, iparams->disres.up1,
1225 iparams->disres.up2, iparams->disres.kfac);
1228 fprintf(fp, "ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
1229 iparams->orires.ex, iparams->orires.label, iparams->orires.power,
1230 iparams->orires.c, iparams->orires.obs, iparams->orires.kfac);
1233 fprintf(fp, "phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
1234 iparams->dihres.phiA, iparams->dihres.dphiA, iparams->dihres.kfacA,
1235 iparams->dihres.phiB, iparams->dihres.dphiB, iparams->dihres.kfacB);
1238 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",
1239 iparams->posres.pos0A[XX], iparams->posres.pos0A[YY],
1240 iparams->posres.pos0A[ZZ], iparams->posres.fcA[XX],
1241 iparams->posres.fcA[YY], iparams->posres.fcA[ZZ],
1242 iparams->posres.pos0B[XX], iparams->posres.pos0B[YY],
1243 iparams->posres.pos0B[ZZ], iparams->posres.fcB[XX],
1244 iparams->posres.fcB[YY], iparams->posres.fcB[ZZ]);
1247 fprintf(fp, "pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e\n",
1248 iparams->fbposres.pos0[XX], iparams->fbposres.pos0[YY],
1249 iparams->fbposres.pos0[ZZ], iparams->fbposres.geom,
1250 iparams->fbposres.r, iparams->fbposres.k);
1253 for (i = 0; i < NR_RBDIHS; i++)
1255 fprintf(fp, "%srbcA[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcA[i]);
1258 for (i = 0; i < NR_RBDIHS; i++)
1260 fprintf(fp, "%srbcB[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcB[i]);
1265 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
1266 * OPLS potential constants back.
1268 rbcA = iparams->rbdihs.rbcA;
1269 rbcB = iparams->rbdihs.rbcB;
1271 VA[3] = -0.25*rbcA[4];
1272 VA[2] = -0.5*rbcA[3];
1273 VA[1] = 4.0*VA[3]-rbcA[2];
1274 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1276 VB[3] = -0.25*rbcB[4];
1277 VB[2] = -0.5*rbcB[3];
1278 VB[1] = 4.0*VB[3]-rbcB[2];
1279 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1281 for (i = 0; i < NR_FOURDIHS; i++)
1283 fprintf(fp, "%sFourA[%d]=%15.8e", i == 0 ? "" : ", ", i, VA[i]);
1286 for (i = 0; i < NR_FOURDIHS; i++)
1288 fprintf(fp, "%sFourB[%d]=%15.8e", i == 0 ? "" : ", ", i, VB[i]);
1295 fprintf(fp, "dA=%15.8e, dB=%15.8e\n", iparams->constr.dA, iparams->constr.dB);
1298 fprintf(fp, "doh=%15.8e, dhh=%15.8e\n", iparams->settle.doh,
1299 iparams->settle.dhh);
1302 fprintf(fp, "a=%15.8e\n", iparams->vsite.a);
1307 fprintf(fp, "a=%15.8e, b=%15.8e\n", iparams->vsite.a, iparams->vsite.b);
1312 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1313 iparams->vsite.a, iparams->vsite.b, iparams->vsite.c);
1316 fprintf(fp, "n=%2d, a=%15.8e\n", iparams->vsiten.n, iparams->vsiten.a);
1321 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);
1324 fprintf(fp, "cmapA=%1d, cmapB=%1d\n", iparams->cmap.cmapA, iparams->cmap.cmapB);
1327 pr_harm(fp, iparams, "ktheta", "costheta0");
1330 fprintf(fp, "phiA=%15.8e, cpA=%15.8e",
1331 iparams->pdihs.phiA, iparams->pdihs.cpA);
1334 fprintf(fp, "kphi=%15.8e", iparams->cbtdihs.cbtcA[0]);
1335 for (i = 1; i < NR_CBTDIHS; i++)
1337 fprintf(fp, ", cbtcA[%d]=%15.8e", i-1, iparams->cbtdihs.cbtcA[i]);
1342 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1343 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1347 void pr_ilist(FILE *fp, int indent, const char *title,
1348 t_functype *functype, t_ilist *ilist, gmx_bool bShowNumbers)
1350 int i, j, k, type, ftype;
1353 if (available(fp, ilist, indent, title) && ilist->nr > 0)
1355 indent = pr_title(fp, indent, title);
1356 (void) pr_indent(fp, indent);
1357 fprintf(fp, "nr: %d\n", ilist->nr);
1360 (void) pr_indent(fp, indent);
1361 fprintf(fp, "iatoms:\n");
1362 iatoms = ilist->iatoms;
1363 for (i = j = 0; i < ilist->nr; )
1366 (void) pr_indent(fp, indent+INDENT);
1368 ftype = functype[type];
1369 (void) fprintf(fp, "%d type=%d (%s)",
1370 bShowNumbers ? j : -1, bShowNumbers ? type : -1,
1371 interaction_function[ftype].name);
1373 for (k = 0; k < interaction_function[ftype].nratoms; k++)
1375 (void) fprintf(fp, " %u", *(iatoms++));
1377 (void) fprintf(fp, "\n");
1378 i += 1+interaction_function[ftype].nratoms;
1380 fprintf(fp, "%5d%5d\n", i, iatoms[i]);
1388 static void pr_cmap(FILE *fp, int indent, const char *title,
1389 gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1394 dx = 360.0 / cmap_grid->grid_spacing;
1395 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1397 if (available(fp, cmap_grid, indent, title))
1399 fprintf(fp, "%s\n", title);
1401 for (i = 0; i < cmap_grid->ngrid; i++)
1404 fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
1406 fprintf(fp, "grid[%3d]={\n", bShowNumbers ? i : -1);
1408 for (j = 0; j < nelem; j++)
1410 if ( (j%cmap_grid->grid_spacing) == 0)
1412 fprintf(fp, "%8.1f\n", idx);
1416 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4]);
1417 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+1]);
1418 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+2]);
1419 fprintf(fp, "%8.3f\n", cmap_grid->cmapdata[i].cmap[j*4+3]);
1427 void pr_ffparams(FILE *fp, int indent, const char *title,
1428 gmx_ffparams_t *ffparams,
1429 gmx_bool bShowNumbers)
1433 indent = pr_title(fp, indent, title);
1434 (void) pr_indent(fp, indent);
1435 (void) fprintf(fp, "atnr=%d\n", ffparams->atnr);
1436 (void) pr_indent(fp, indent);
1437 (void) fprintf(fp, "ntypes=%d\n", ffparams->ntypes);
1438 for (i = 0; i < ffparams->ntypes; i++)
1440 (void) pr_indent(fp, indent+INDENT);
1441 (void) fprintf(fp, "functype[%d]=%s, ",
1442 bShowNumbers ? i : -1,
1443 interaction_function[ffparams->functype[i]].name);
1444 pr_iparams(fp, ffparams->functype[i], &ffparams->iparams[i]);
1446 (void) pr_double(fp, indent, "reppow", ffparams->reppow);
1447 (void) pr_real(fp, indent, "fudgeQQ", ffparams->fudgeQQ);
1448 pr_cmap(fp, indent, "cmap", &ffparams->cmap_grid, bShowNumbers);
1451 void pr_idef(FILE *fp, int indent, const char *title, t_idef *idef, gmx_bool bShowNumbers)
1455 if (available(fp, idef, indent, title))
1457 indent = pr_title(fp, indent, title);
1458 (void) pr_indent(fp, indent);
1459 (void) fprintf(fp, "atnr=%d\n", idef->atnr);
1460 (void) pr_indent(fp, indent);
1461 (void) fprintf(fp, "ntypes=%d\n", idef->ntypes);
1462 for (i = 0; i < idef->ntypes; i++)
1464 (void) pr_indent(fp, indent+INDENT);
1465 (void) fprintf(fp, "functype[%d]=%s, ",
1466 bShowNumbers ? i : -1,
1467 interaction_function[idef->functype[i]].name);
1468 pr_iparams(fp, idef->functype[i], &idef->iparams[i]);
1470 (void) pr_real(fp, indent, "fudgeQQ", idef->fudgeQQ);
1472 for (j = 0; (j < F_NRE); j++)
1474 pr_ilist(fp, indent, interaction_function[j].longname,
1475 idef->functype, &idef->il[j], bShowNumbers);
1480 static int pr_block_title(FILE *fp, int indent, const char *title, t_block *block)
1484 if (available(fp, block, indent, title))
1486 indent = pr_title(fp, indent, title);
1487 (void) pr_indent(fp, indent);
1488 (void) fprintf(fp, "nr=%d\n", block->nr);
1493 static int pr_blocka_title(FILE *fp, int indent, const char *title, t_blocka *block)
1497 if (available(fp, block, indent, title))
1499 indent = pr_title(fp, indent, title);
1500 (void) pr_indent(fp, indent);
1501 (void) fprintf(fp, "nr=%d\n", block->nr);
1502 (void) pr_indent(fp, indent);
1503 (void) fprintf(fp, "nra=%d\n", block->nra);
1508 static void low_pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1512 if (available(fp, block, indent, title))
1514 indent = pr_blocka_title(fp, indent, title, block);
1515 for (i = 0; i <= block->nr; i++)
1517 (void) pr_indent(fp, indent+INDENT);
1518 (void) fprintf(fp, "%s->index[%d]=%u\n",
1519 title, bShowNumbers ? i : -1, block->index[i]);
1521 for (i = 0; i < block->nra; i++)
1523 (void) pr_indent(fp, indent+INDENT);
1524 (void) fprintf(fp, "%s->a[%d]=%u\n",
1525 title, bShowNumbers ? i : -1, block->a[i]);
1530 void pr_block(FILE *fp, int indent, const char *title, t_block *block, gmx_bool bShowNumbers)
1532 int i, j, ok, size, start, end;
1534 if (available(fp, block, indent, title))
1536 indent = pr_block_title(fp, indent, title, block);
1539 if ((ok = (block->index[start] == 0)) == 0)
1541 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1545 for (i = 0; i < block->nr; i++)
1547 end = block->index[i+1];
1548 size = pr_indent(fp, indent);
1551 size += fprintf(fp, "%s[%d]={}\n", title, i);
1555 size += fprintf(fp, "%s[%d]={%d..%d}\n",
1556 title, bShowNumbers ? i : -1,
1557 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1565 void pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1567 int i, j, ok, size, start, end;
1569 if (available(fp, block, indent, title))
1571 indent = pr_blocka_title(fp, indent, title, block);
1574 if ((ok = (block->index[start] == 0)) == 0)
1576 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1580 for (i = 0; i < block->nr; i++)
1582 end = block->index[i+1];
1583 size = pr_indent(fp, indent);
1586 size += fprintf(fp, "%s[%d]={", title, i);
1590 size += fprintf(fp, "%s[%d][%d..%d]={",
1591 title, bShowNumbers ? i : -1,
1592 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1594 for (j = start; j < end; j++)
1598 size += fprintf(fp, ", ");
1600 if ((size) > (USE_WIDTH))
1602 (void) fprintf(fp, "\n");
1603 size = pr_indent(fp, indent+INDENT);
1605 size += fprintf(fp, "%u", block->a[j]);
1607 (void) fprintf(fp, "}\n");
1611 if ((end != block->nra) || (!ok))
1613 (void) pr_indent(fp, indent);
1614 (void) fprintf(fp, "tables inconsistent, dumping complete tables:\n");
1615 low_pr_blocka(fp, indent, title, block, bShowNumbers);
1620 static void pr_strings(FILE *fp, int indent, const char *title, char ***nm, int n, gmx_bool bShowNumbers)
1624 if (available(fp, nm, indent, title))
1626 indent = pr_title_n(fp, indent, title, n);
1627 for (i = 0; i < n; i++)
1629 (void) pr_indent(fp, indent);
1630 (void) fprintf(fp, "%s[%d]={name=\"%s\"}\n",
1631 title, bShowNumbers ? i : -1, *(nm[i]));
1636 static void pr_strings2(FILE *fp, int indent, const char *title,
1637 char ***nm, char ***nmB, int n, gmx_bool bShowNumbers)
1641 if (available(fp, nm, indent, title))
1643 indent = pr_title_n(fp, indent, title, n);
1644 for (i = 0; i < n; i++)
1646 (void) pr_indent(fp, indent);
1647 (void) fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1648 title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
1653 static void pr_resinfo(FILE *fp, int indent, const char *title, t_resinfo *resinfo, int n, gmx_bool bShowNumbers)
1657 if (available(fp, resinfo, indent, title))
1659 indent = pr_title_n(fp, indent, title, n);
1660 for (i = 0; i < n; i++)
1662 (void) pr_indent(fp, indent);
1663 (void) fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1664 title, bShowNumbers ? i : -1,
1665 *(resinfo[i].name), resinfo[i].nr,
1666 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1671 static void pr_atom(FILE *fp, int indent, const char *title, t_atom *atom, int n)
1675 if (available(fp, atom, indent, title))
1677 indent = pr_title_n(fp, indent, title, n);
1678 for (i = 0; i < n; i++)
1680 (void) pr_indent(fp, indent);
1681 fprintf(fp, "%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1682 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1683 title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype],
1684 atom[i].m, atom[i].q, atom[i].mB, atom[i].qB,
1685 atom[i].resind, atom[i].atomnumber);
1690 static void pr_grps(FILE *fp, const char *title, t_grps grps[], char **grpname[])
1694 for (i = 0; (i < egcNR); i++)
1696 fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
1697 for (j = 0; (j < grps[i].nr); j++)
1699 fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
1705 static void pr_groups(FILE *fp, int indent,
1706 gmx_groups_t *groups,
1707 gmx_bool bShowNumbers)
1712 pr_grps(fp, "grp", groups->grps, groups->grpname);
1713 pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
1715 (void) pr_indent(fp, indent);
1716 fprintf(fp, "groups ");
1717 for (g = 0; g < egcNR; g++)
1719 printf(" %5.5s", gtypes[g]);
1723 (void) pr_indent(fp, indent);
1724 fprintf(fp, "allocated ");
1726 for (g = 0; g < egcNR; g++)
1728 printf(" %5d", groups->ngrpnr[g]);
1729 nat_max = max(nat_max, groups->ngrpnr[g]);
1735 (void) pr_indent(fp, indent);
1736 fprintf(fp, "groupnr[%5s] =", "*");
1737 for (g = 0; g < egcNR; g++)
1739 fprintf(fp, " %3d ", 0);
1745 for (i = 0; i < nat_max; i++)
1747 (void) pr_indent(fp, indent);
1748 fprintf(fp, "groupnr[%5d] =", i);
1749 for (g = 0; g < egcNR; g++)
1751 fprintf(fp, " %3d ",
1752 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1759 void pr_atoms(FILE *fp, int indent, const char *title, t_atoms *atoms,
1760 gmx_bool bShownumbers)
1762 if (available(fp, atoms, indent, title))
1764 indent = pr_title(fp, indent, title);
1765 pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
1766 pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
1767 pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
1768 pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
1773 void pr_atomtypes(FILE *fp, int indent, const char *title, t_atomtypes *atomtypes,
1774 gmx_bool bShowNumbers)
1777 if (available(fp, atomtypes, indent, title))
1779 indent = pr_title(fp, indent, title);
1780 for (i = 0; i < atomtypes->nr; i++)
1782 pr_indent(fp, indent);
1784 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1785 bShowNumbers ? i : -1, atomtypes->radius[i], atomtypes->vol[i],
1786 atomtypes->gb_radius[i],
1787 atomtypes->surftens[i], atomtypes->atomnumber[i], atomtypes->S_hct[i]);
1792 static void pr_moltype(FILE *fp, int indent, const char *title,
1793 gmx_moltype_t *molt, int n,
1794 gmx_ffparams_t *ffparams,
1795 gmx_bool bShowNumbers)
1799 indent = pr_title_n(fp, indent, title, n);
1800 (void) pr_indent(fp, indent);
1801 (void) fprintf(fp, "name=\"%s\"\n", *(molt->name));
1802 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
1803 pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
1804 pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
1805 for (j = 0; (j < F_NRE); j++)
1807 pr_ilist(fp, indent, interaction_function[j].longname,
1808 ffparams->functype, &molt->ilist[j], bShowNumbers);
1812 static void pr_molblock(FILE *fp, int indent, const char *title,
1813 gmx_molblock_t *molb, int n,
1814 gmx_moltype_t *molt)
1816 indent = pr_title_n(fp, indent, title, n);
1817 (void) pr_indent(fp, indent);
1818 (void) fprintf(fp, "%-20s = %d \"%s\"\n",
1819 "moltype", molb->type, *(molt[molb->type].name));
1820 pr_int(fp, indent, "#molecules", molb->nmol);
1821 pr_int(fp, indent, "#atoms_mol", molb->natoms_mol);
1822 pr_int(fp, indent, "#posres_xA", molb->nposres_xA);
1823 if (molb->nposres_xA > 0)
1825 pr_rvecs(fp, indent, "posres_xA", molb->posres_xA, molb->nposres_xA);
1827 pr_int(fp, indent, "#posres_xB", molb->nposres_xB);
1828 if (molb->nposres_xB > 0)
1830 pr_rvecs(fp, indent, "posres_xB", molb->posres_xB, molb->nposres_xB);
1834 void pr_mtop(FILE *fp, int indent, const char *title, gmx_mtop_t *mtop,
1835 gmx_bool bShowNumbers)
1839 if (available(fp, mtop, indent, title))
1841 indent = pr_title(fp, indent, title);
1842 (void) pr_indent(fp, indent);
1843 (void) fprintf(fp, "name=\"%s\"\n", *(mtop->name));
1844 pr_int(fp, indent, "#atoms", mtop->natoms);
1845 pr_int(fp, indent, "#molblock", mtop->nmolblock);
1846 for (mb = 0; mb < mtop->nmolblock; mb++)
1848 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
1850 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
1851 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
1852 for (mt = 0; mt < mtop->nmoltype; mt++)
1854 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
1855 &mtop->ffparams, bShowNumbers);
1857 pr_groups(fp, indent, &mtop->groups, bShowNumbers);
1861 void pr_top(FILE *fp, int indent, const char *title, t_topology *top, gmx_bool bShowNumbers)
1863 if (available(fp, top, indent, title))
1865 indent = pr_title(fp, indent, title);
1866 (void) pr_indent(fp, indent);
1867 (void) fprintf(fp, "name=\"%s\"\n", *(top->name));
1868 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
1869 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
1870 pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
1871 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
1872 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
1873 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers);
1877 void pr_header(FILE *fp, int indent, const char *title, t_tpxheader *sh)
1881 if (available(fp, sh, indent, title))
1883 indent = pr_title(fp, indent, title);
1884 pr_indent(fp, indent);
1885 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
1886 pr_indent(fp, indent);
1887 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
1888 pr_indent(fp, indent);
1889 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
1890 pr_indent(fp, indent);
1891 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
1892 pr_indent(fp, indent);
1893 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
1894 pr_indent(fp, indent);
1895 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
1897 pr_indent(fp, indent);
1898 fprintf(fp, "natoms = %d\n", sh->natoms);
1899 pr_indent(fp, indent);
1900 fprintf(fp, "lambda = %e\n", sh->lambda);
1904 void pr_commrec(FILE *fp, int indent, t_commrec *cr)
1906 pr_indent(fp, indent);
1907 fprintf(fp, "commrec:\n");
1909 pr_indent(fp, indent);
1910 fprintf(fp, "nodeid = %d\n", cr->nodeid);
1911 pr_indent(fp, indent);
1912 fprintf(fp, "nnodes = %d\n", cr->nnodes);
1913 pr_indent(fp, indent);
1914 fprintf(fp, "npmenodes = %d\n", cr->npmenodes);
1916 pr_indent(fp,indent);
1917 fprintf(fp,"threadid = %d\n",cr->threadid);
1918 pr_indent(fp,indent);
1919 fprintf(fp,"nthreads = %d\n",cr->nthreads);