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,2015, 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.
39 /* This file is completely threadsafe - please keep it that way! */
41 #include "gromacs/legacyheaders/txtdump.h"
46 #include "gromacs/legacyheaders/macros.h"
47 #include "gromacs/legacyheaders/names.h"
48 #include "gromacs/legacyheaders/typedefs.h"
49 #include "gromacs/legacyheaders/types/commrec.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
54 int pr_indent(FILE *fp, int n)
58 for (i = 0; i < n; i++)
60 (void) fprintf(fp, " ");
65 int available(FILE *fp, void *p, int indent, const char *title)
71 pr_indent(fp, indent);
73 (void) fprintf(fp, "%s: not available\n", title);
78 int pr_title(FILE *fp, int indent, const char *title)
80 (void) pr_indent(fp, indent);
81 (void) fprintf(fp, "%s:\n", title);
82 return (indent+INDENT);
85 int pr_title_n(FILE *fp, int indent, const char *title, int n)
87 (void) pr_indent(fp, indent);
88 (void) fprintf(fp, "%s (%d):\n", title, n);
89 return (indent+INDENT);
92 int pr_title_nxn(FILE *fp, int indent, const char *title, int n1, int n2)
94 (void) pr_indent(fp, indent);
95 (void) fprintf(fp, "%s (%dx%d):\n", title, n1, n2);
96 return (indent+INDENT);
99 void pr_ivec(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
103 if (available(fp, vec, indent, title))
105 indent = pr_title_n(fp, indent, title, n);
106 for (i = 0; i < n; i++)
108 (void) pr_indent(fp, indent);
109 (void) fprintf(fp, "%s[%d]=%d\n", title, bShowNumbers ? i : -1, vec[i]);
114 void pr_ivec_block(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
118 if (available(fp, vec, indent, title))
120 indent = pr_title_n(fp, indent, title, n);
125 while (j < n && vec[j] == vec[j-1]+1)
129 /* Print consecutive groups of 3 or more as blocks */
134 (void) pr_indent(fp, indent);
135 (void) fprintf(fp, "%s[%d]=%d\n",
136 title, bShowNumbers ? i : -1, vec[i]);
142 (void) pr_indent(fp, indent);
143 (void) fprintf(fp, "%s[%d,...,%d] = {%d,...,%d}\n",
145 bShowNumbers ? i : -1,
146 bShowNumbers ? j-1 : -1,
154 void pr_bvec(FILE *fp, int indent, const char *title, gmx_bool vec[], int n, gmx_bool bShowNumbers)
158 if (available(fp, vec, indent, title))
160 indent = pr_title_n(fp, indent, title, n);
161 for (i = 0; i < n; i++)
163 (void) pr_indent(fp, indent);
164 (void) fprintf(fp, "%s[%d]=%s\n", title, bShowNumbers ? i : -1,
170 void pr_ivecs(FILE *fp, int indent, const char *title, ivec vec[], int n, gmx_bool bShowNumbers)
174 if (available(fp, vec, indent, title))
176 indent = pr_title_nxn(fp, indent, title, n, DIM);
177 for (i = 0; i < n; i++)
179 (void) pr_indent(fp, indent);
180 (void) fprintf(fp, "%s[%d]={", title, bShowNumbers ? i : -1);
181 for (j = 0; j < DIM; j++)
185 (void) fprintf(fp, ", ");
187 fprintf(fp, "%d", vec[i][j]);
189 (void) fprintf(fp, "}\n");
194 void pr_rvec(FILE *fp, int indent, const char *title, real vec[], int n, gmx_bool bShowNumbers)
198 if (available(fp, vec, indent, title))
200 indent = pr_title_n(fp, indent, title, n);
201 for (i = 0; i < n; i++)
203 pr_indent(fp, indent);
204 fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
209 void pr_dvec(FILE *fp, int indent, const char *title, double vec[], int n, gmx_bool bShowNumbers)
213 if (available(fp, vec, indent, title))
215 indent = pr_title_n(fp, indent, title, n);
216 for (i = 0; i < n; i++)
218 pr_indent(fp, indent);
219 fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
226 void pr_mat(FILE *fp,int indent,char *title,matrix m)
230 if (available(fp,m,indent,title)) {
231 indent=pr_title_n(fp,indent,title,n);
233 pr_indent(fp,indent);
234 fprintf(fp,"%s[%d]=%12.5e %12.5e %12.5e\n",
235 title,bShowNumbers?i:-1,m[i][XX],m[i][YY],m[i][ZZ]);
241 void pr_rvecs_len(FILE *fp, int indent, const char *title, rvec vec[], int n)
245 if (available(fp, vec, indent, title))
247 indent = pr_title_nxn(fp, indent, title, n, DIM);
248 for (i = 0; i < n; i++)
250 (void) pr_indent(fp, indent);
251 (void) fprintf(fp, "%s[%5d]={", title, i);
252 for (j = 0; j < DIM; j++)
256 (void) fprintf(fp, ", ");
258 (void) fprintf(fp, "%12.5e", vec[i][j]);
260 (void) fprintf(fp, "} len=%12.5e\n", norm(vec[i]));
265 void pr_rvecs(FILE *fp, int indent, const char *title, rvec vec[], int n)
267 const char *fshort = "%12.5e";
268 const char *flong = "%15.8e";
272 if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
281 if (available(fp, vec, indent, title))
283 indent = pr_title_nxn(fp, indent, title, n, DIM);
284 for (i = 0; i < n; i++)
286 (void) pr_indent(fp, indent);
287 (void) fprintf(fp, "%s[%5d]={", title, i);
288 for (j = 0; j < DIM; j++)
292 (void) fprintf(fp, ", ");
294 (void) fprintf(fp, format, vec[i][j]);
296 (void) fprintf(fp, "}\n");
302 void pr_rvecs_of_dim(FILE *fp, int indent, const char *title, rvec vec[], int n, int dim)
304 const char *fshort = "%12.5e";
305 const char *flong = "%15.8e";
309 if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
318 if (available(fp, vec, indent, title))
320 indent = pr_title_nxn(fp, indent, title, n, dim);
321 for (i = 0; i < n; i++)
323 (void) pr_indent(fp, indent);
324 (void) fprintf(fp, "%s[%5d]={", title, i);
325 for (j = 0; j < dim; j++)
329 (void) fprintf(fp, ", ");
331 (void) fprintf(fp, format, vec[i][j]);
333 (void) fprintf(fp, "}\n");
338 void pr_reals(FILE *fp, int indent, const char *title, real *vec, int n)
342 if (available(fp, vec, indent, title))
344 (void) pr_indent(fp, indent);
345 (void) fprintf(fp, "%s:\t", title);
346 for (i = 0; i < n; i++)
348 fprintf(fp, " %10g", vec[i]);
350 (void) fprintf(fp, "\n");
354 void pr_doubles(FILE *fp, int indent, const char *title, double *vec, int n)
358 if (available(fp, vec, indent, title))
360 (void) pr_indent(fp, indent);
361 (void) fprintf(fp, "%s:\t", title);
362 for (i = 0; i < n; i++)
364 fprintf(fp, " %10g", vec[i]);
366 (void) fprintf(fp, "\n");
370 void pr_reals_of_dim(FILE *fp, int indent, const char *title, real *vec, int n, int dim)
373 const char *fshort = "%12.5e";
374 const char *flong = "%15.8e";
377 if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
386 if (available(fp, vec, indent, title))
388 indent = pr_title_nxn(fp, indent, title, n, dim);
389 for (i = 0; i < n; i++)
391 (void) pr_indent(fp, indent);
392 (void) fprintf(fp, "%s[%5d]={", title, i);
393 for (j = 0; j < dim; j++)
397 (void) fprintf(fp, ", ");
399 (void) fprintf(fp, format, vec[i * dim + j]);
401 (void) fprintf(fp, "}\n");
406 static void pr_int(FILE *fp, int indent, const char *title, int i)
408 pr_indent(fp, indent);
409 fprintf(fp, "%-30s = %d\n", title, i);
412 static void pr_int64(FILE *fp, int indent, const char *title, gmx_int64_t i)
414 char buf[STEPSTRSIZE];
416 pr_indent(fp, indent);
417 fprintf(fp, "%-30s = %s\n", title, gmx_step_str(i, buf));
420 static void pr_real(FILE *fp, int indent, const char *title, real r)
422 pr_indent(fp, indent);
423 fprintf(fp, "%-30s = %g\n", title, r);
426 static void pr_double(FILE *fp, int indent, const char *title, double d)
428 pr_indent(fp, indent);
429 fprintf(fp, "%-30s = %g\n", title, d);
432 static void pr_str(FILE *fp, int indent, const char *title, const char *s)
434 pr_indent(fp, indent);
435 fprintf(fp, "%-30s = %s\n", title, s);
438 void pr_qm_opts(FILE *fp, int indent, const char *title, t_grpopts *opts)
442 fprintf(fp, "%s:\n", title);
444 pr_int(fp, indent, "ngQM", opts->ngQM);
447 pr_ivec(fp, indent, "QMmethod", opts->QMmethod, opts->ngQM, FALSE);
448 pr_ivec(fp, indent, "QMbasis", opts->QMbasis, opts->ngQM, FALSE);
449 pr_ivec(fp, indent, "QMcharge", opts->QMcharge, opts->ngQM, FALSE);
450 pr_ivec(fp, indent, "QMmult", opts->QMmult, opts->ngQM, FALSE);
451 pr_bvec(fp, indent, "SH", opts->bSH, opts->ngQM, FALSE);
452 pr_ivec(fp, indent, "CASorbitals", opts->CASorbitals, opts->ngQM, FALSE);
453 pr_ivec(fp, indent, "CASelectrons", opts->CASelectrons, opts->ngQM, FALSE);
454 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
455 pr_rvec(fp, indent, "SAoff", opts->SAoff, opts->ngQM, FALSE);
456 pr_ivec(fp, indent, "SAsteps", opts->SAsteps, opts->ngQM, FALSE);
457 pr_bvec(fp, indent, "bOPT", opts->bOPT, opts->ngQM, FALSE);
458 pr_bvec(fp, indent, "bTS", opts->bTS, opts->ngQM, FALSE);
462 static void pr_grp_opts(FILE *out, int indent, const char *title, t_grpopts *opts,
469 fprintf(out, "%s:\n", title);
472 pr_indent(out, indent);
473 fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
474 for (i = 0; (i < opts->ngtc); i++)
476 fprintf(out, " %10g", opts->nrdf[i]);
480 pr_indent(out, indent);
481 fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
482 for (i = 0; (i < opts->ngtc); i++)
484 fprintf(out, " %10g", opts->ref_t[i]);
488 pr_indent(out, indent);
489 fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
490 for (i = 0; (i < opts->ngtc); i++)
492 fprintf(out, " %10g", opts->tau_t[i]);
496 /* Pretty-print the simulated annealing info */
497 fprintf(out, "annealing%s", bMDPformat ? " = " : ":");
498 for (i = 0; (i < opts->ngtc); i++)
500 fprintf(out, " %10s", EANNEAL(opts->annealing[i]));
504 fprintf(out, "annealing-npoints%s", bMDPformat ? " = " : ":");
505 for (i = 0; (i < opts->ngtc); i++)
507 fprintf(out, " %10d", opts->anneal_npoints[i]);
511 for (i = 0; (i < opts->ngtc); i++)
513 if (opts->anneal_npoints[i] > 0)
515 fprintf(out, "annealing-time [%d]:\t", i);
516 for (j = 0; (j < opts->anneal_npoints[i]); j++)
518 fprintf(out, " %10.1f", opts->anneal_time[i][j]);
521 fprintf(out, "annealing-temp [%d]:\t", i);
522 for (j = 0; (j < opts->anneal_npoints[i]); j++)
524 fprintf(out, " %10.1f", opts->anneal_temp[i][j]);
530 pr_indent(out, indent);
531 fprintf(out, "acc:\t");
532 for (i = 0; (i < opts->ngacc); i++)
534 for (m = 0; (m < DIM); m++)
536 fprintf(out, " %10g", opts->acc[i][m]);
541 pr_indent(out, indent);
542 fprintf(out, "nfreeze:");
543 for (i = 0; (i < opts->ngfrz); i++)
545 for (m = 0; (m < DIM); m++)
547 fprintf(out, " %10s", opts->nFreeze[i][m] ? "Y" : "N");
553 for (i = 0; (i < opts->ngener); i++)
555 pr_indent(out, indent);
556 fprintf(out, "energygrp-flags[%3d]:", i);
557 for (m = 0; (m < opts->ngener); m++)
559 fprintf(out, " %d", opts->egp_flags[opts->ngener*i+m]);
567 static void pr_matrix(FILE *fp, int indent, const char *title, rvec *m,
572 fprintf(fp, "%-10s = %g %g %g %g %g %g\n", title,
573 m[XX][XX], m[YY][YY], m[ZZ][ZZ], m[XX][YY], m[XX][ZZ], m[YY][ZZ]);
577 pr_rvecs(fp, indent, title, m, DIM);
581 static void pr_cosine(FILE *fp, int indent, const char *title, t_cosines *cos,
588 fprintf(fp, "%s = %d\n", title, cos->n);
592 indent = pr_title(fp, indent, title);
593 (void) pr_indent(fp, indent);
594 fprintf(fp, "n = %d\n", cos->n);
597 (void) pr_indent(fp, indent+2);
599 for (j = 0; (j < cos->n); j++)
601 fprintf(fp, " %e", cos->a[j]);
604 (void) pr_indent(fp, indent+2);
605 fprintf(fp, "phi =");
606 for (j = 0; (j < cos->n); j++)
608 fprintf(fp, " %e", cos->phi[j]);
615 #define PS(t, s) pr_str(fp, indent, t, s)
616 #define PI(t, s) pr_int(fp, indent, t, s)
617 #define PSTEP(t, s) pr_int64(fp, indent, t, s)
618 #define PR(t, s) pr_real(fp, indent, t, s)
619 #define PD(t, s) pr_double(fp, indent, t, s)
621 static void pr_pull_group(FILE *fp, int indent, int g, t_pull_group *pgrp)
623 pr_indent(fp, indent);
624 fprintf(fp, "pull-group %d:\n", g);
626 pr_ivec_block(fp, indent, "atom", pgrp->ind, pgrp->nat, TRUE);
627 pr_rvec(fp, indent, "weight", pgrp->weight, pgrp->nweight, TRUE);
628 PI("pbcatom", pgrp->pbcatom);
631 static void pr_pull_coord(FILE *fp, int indent, int c, t_pull_coord *pcrd)
633 pr_indent(fp, indent);
634 fprintf(fp, "pull-coord %d:\n", c);
635 PI("group[0]", pcrd->group[0]);
636 PI("group[1]", pcrd->group[1]);
637 if (pcrd->eGeom == epullgDIRRELATIVE)
639 PI("group[2]", pcrd->group[2]);
640 PI("group[3]", pcrd->group[3]);
642 PS("type", EPULLTYPE(pcrd->eType));
643 PS("geometry", EPULLGEOM(pcrd->eGeom));
644 pr_ivec(fp, indent, "dim", pcrd->dim, DIM, TRUE);
645 pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
646 pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
647 PS("start", EBOOL(pcrd->bStart));
648 PR("init", pcrd->init);
649 PR("rate", pcrd->rate);
654 static void pr_simtempvals(FILE *fp, int indent, t_simtemp *simtemp, int n_lambda)
656 PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
657 PR("sim-temp-low", simtemp->simtemp_low);
658 PR("sim-temp-high", simtemp->simtemp_high);
659 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
662 static void pr_expandedvals(FILE *fp, int indent, t_expanded *expand, int n_lambda)
665 PI("nstexpanded", expand->nstexpanded);
666 PS("lmc-stats", elamstats_names[expand->elamstats]);
667 PS("lmc-move", elmcmove_names[expand->elmcmove]);
668 PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
669 if (expand->elmceq == elmceqNUMATLAM)
671 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
673 if (expand->elmceq == elmceqSAMPLES)
675 PI("weight-equil-number-samples", expand->equil_samples);
677 if (expand->elmceq == elmceqSTEPS)
679 PI("weight-equil-number-steps", expand->equil_steps);
681 if (expand->elmceq == elmceqWLDELTA)
683 PR("weight-equil-wl-delta", expand->equil_wl_delta);
685 if (expand->elmceq == elmceqRATIO)
687 PR("weight-equil-count-ratio", expand->equil_ratio);
689 PI("lmc-seed", expand->lmc_seed);
690 PR("mc-temperature", expand->mc_temp);
691 PI("lmc-repeats", expand->lmc_repeats);
692 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
693 PI("lmc-forced-nstart", expand->lmc_forced_nstart);
694 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
695 PI("nst-transition-matrix", expand->nstTij);
696 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
697 PI("weight-c-range", expand->c_range); /* default is just C=0 */
698 PR("wl-scale", expand->wl_scale);
699 PR("wl-ratio", expand->wl_ratio);
700 PR("init-wl-delta", expand->init_wl_delta);
701 PS("wl-oneovert", EBOOL(expand->bWLoneovert));
703 pr_indent(fp, indent);
704 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
705 PS("init-weights", EBOOL(expand->bInit_weights));
708 static void pr_fepvals(FILE *fp, int indent, t_lambda *fep, gmx_bool bMDPformat)
712 PR("init-lambda", fep->init_lambda);
713 PI("init-lambda-state", fep->init_fep_state);
714 PR("delta-lambda", fep->delta_lambda);
715 PI("nstdhdl", fep->nstdhdl);
719 PI("n-lambdas", fep->n_lambda);
721 if (fep->n_lambda > 0)
723 pr_indent(fp, indent);
724 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
725 for (i = 0; i < efptNR; i++)
727 fprintf(fp, "%18s = ", efpt_names[i]);
728 if (fep->separate_dvdl[i])
730 fprintf(fp, " TRUE");
734 fprintf(fp, " FALSE");
738 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
739 for (i = 0; i < efptNR; i++)
741 fprintf(fp, "%18s = ", efpt_names[i]);
742 for (j = 0; j < fep->n_lambda; j++)
744 fprintf(fp, " %10g", fep->all_lambda[i][j]);
749 PI("calc-lambda-neighbors", fep->lambda_neighbors);
750 PS("dhdl-print-energy", edHdLPrintEnergy_names[fep->edHdLPrintEnergy]);
751 PR("sc-alpha", fep->sc_alpha);
752 PI("sc-power", fep->sc_power);
753 PR("sc-r-power", fep->sc_r_power);
754 PR("sc-sigma", fep->sc_sigma);
755 PR("sc-sigma-min", fep->sc_sigma_min);
756 PS("sc-coul", EBOOL(fep->bScCoul));
757 PI("dh-hist-size", fep->dh_hist_size);
758 PD("dh-hist-spacing", fep->dh_hist_spacing);
759 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
760 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
763 static void pr_pull(FILE *fp, int indent, t_pull *pull)
767 PR("pull-cylinder-r", pull->cylinder_r);
768 PR("pull-constr-tol", pull->constr_tol);
769 PS("pull-print-COM1", EBOOL(pull->bPrintCOM1));
770 PS("pull-print-COM2", EBOOL(pull->bPrintCOM2));
771 PS("pull-print-ref-value", EBOOL(pull->bPrintRefValue));
772 PS("pull-print-components", EBOOL(pull->bPrintComp));
773 PI("pull-nstxout", pull->nstxout);
774 PI("pull-nstfout", pull->nstfout);
775 PI("pull-ngroups", pull->ngroup);
776 for (g = 0; g < pull->ngroup; g++)
778 pr_pull_group(fp, indent, g, &pull->group[g]);
780 PI("pull-ncoords", pull->ncoord);
781 for (g = 0; g < pull->ncoord; g++)
783 pr_pull_coord(fp, indent, g, &pull->coord[g]);
787 static void pr_rotgrp(FILE *fp, int indent, int g, t_rotgrp *rotg)
789 pr_indent(fp, indent);
790 fprintf(fp, "rot-group %d:\n", g);
792 PS("rot-type", EROTGEOM(rotg->eType));
793 PS("rot-massw", EBOOL(rotg->bMassW));
794 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
795 pr_rvecs(fp, indent, "x-ref", rotg->x_ref, rotg->nat);
796 pr_rvec(fp, indent, "rot-vec", rotg->vec, DIM, TRUE);
797 pr_rvec(fp, indent, "rot-pivot", rotg->pivot, DIM, TRUE);
798 PR("rot-rate", rotg->rate);
799 PR("rot-k", rotg->k);
800 PR("rot-slab-dist", rotg->slab_dist);
801 PR("rot-min-gauss", rotg->min_gaussian);
802 PR("rot-eps", rotg->eps);
803 PS("rot-fit-method", EROTFIT(rotg->eFittype));
804 PI("rot-potfit-nstep", rotg->PotAngle_nstep);
805 PR("rot-potfit-step", rotg->PotAngle_step);
808 static void pr_rot(FILE *fp, int indent, t_rot *rot)
812 PI("rot-nstrout", rot->nstrout);
813 PI("rot-nstsout", rot->nstsout);
814 PI("rot-ngroups", rot->ngrp);
815 for (g = 0; g < rot->ngrp; g++)
817 pr_rotgrp(fp, indent, g, &rot->grp[g]);
822 static void pr_swap(FILE *fp, int indent, t_swapcoords *swap)
828 PI("swap-frequency", swap->nstswap);
829 for (j = 0; j < 2; j++)
831 sprintf(str, "massw_split%d", j);
832 PS(str, EBOOL(swap->massw_split[j]));
833 sprintf(str, "split atoms group %d", j);
834 pr_ivec_block(fp, indent, str, swap->ind_split[j], swap->nat_split[j], TRUE);
836 pr_ivec_block(fp, indent, "swap atoms", swap->ind, swap->nat, TRUE);
837 pr_ivec_block(fp, indent, "solvent atoms", swap->ind_sol, swap->nat_sol, TRUE);
838 PR("cyl0-r", swap->cyl0r);
839 PR("cyl0-up", swap->cyl0u);
840 PR("cyl0-down", swap->cyl0l);
841 PR("cyl1-r", swap->cyl1r);
842 PR("cyl1-up", swap->cyl1u);
843 PR("cyl1-down", swap->cyl1l);
844 PI("coupl-steps", swap->nAverage);
845 for (j = 0; j < 2; j++)
847 sprintf(str, "anions%c", j+'A');
848 PI(str, swap->nanions[j]);
849 sprintf(str, "cations%c", j+'A');
850 PI(str, swap->ncations[j]);
852 PR("threshold", swap->threshold);
856 static void pr_imd(FILE *fp, int indent, t_IMD *imd)
858 PI("IMD-atoms", imd->nat);
859 pr_ivec_block(fp, indent, "atom", imd->ind, imd->nat, TRUE);
863 void pr_inputrec(FILE *fp, int indent, const char *title, t_inputrec *ir,
866 const char *infbuf = "inf";
869 if (available(fp, ir, indent, title))
873 indent = pr_title(fp, indent, title);
875 /* Try to make this list appear in the same order as the
876 * options are written in the default mdout.mdp, and with
877 * the same user-exposed names to facilitate debugging.
879 PS("integrator", EI(ir->eI));
880 PR("tinit", ir->init_t);
881 PR("dt", ir->delta_t);
882 PSTEP("nsteps", ir->nsteps);
883 PSTEP("init-step", ir->init_step);
884 PI("simulation-part", ir->simulation_part);
885 PS("comm-mode", ECOM(ir->comm_mode));
886 PI("nstcomm", ir->nstcomm);
888 /* Langevin dynamics */
889 PR("bd-fric", ir->bd_fric);
890 PSTEP("ld-seed", ir->ld_seed);
892 /* Energy minimization */
893 PR("emtol", ir->em_tol);
894 PR("emstep", ir->em_stepsize);
895 PI("niter", ir->niter);
896 PR("fcstep", ir->fc_stepsize);
897 PI("nstcgsteep", ir->nstcgsteep);
898 PI("nbfgscorr", ir->nbfgscorr);
900 /* Test particle insertion */
901 PR("rtpi", ir->rtpi);
904 PI("nstxout", ir->nstxout);
905 PI("nstvout", ir->nstvout);
906 PI("nstfout", ir->nstfout);
907 PI("nstlog", ir->nstlog);
908 PI("nstcalcenergy", ir->nstcalcenergy);
909 PI("nstenergy", ir->nstenergy);
910 PI("nstxout-compressed", ir->nstxout_compressed);
911 PR("compressed-x-precision", ir->x_compression_precision);
913 /* Neighborsearching parameters */
914 PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
915 PI("nstlist", ir->nstlist);
916 PS("ns-type", ENS(ir->ns_type));
917 PS("pbc", EPBC(ir->ePBC));
918 PS("periodic-molecules", EBOOL(ir->bPeriodicMols));
919 PR("verlet-buffer-tolerance", ir->verletbuf_tol);
920 PR("rlist", ir->rlist);
921 PR("rlistlong", ir->rlistlong);
922 PR("nstcalclr", ir->nstcalclr);
924 /* Options for electrostatics and VdW */
925 PS("coulombtype", EELTYPE(ir->coulombtype));
926 PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
927 PR("rcoulomb-switch", ir->rcoulomb_switch);
928 PR("rcoulomb", ir->rcoulomb);
929 if (ir->epsilon_r != 0)
931 PR("epsilon-r", ir->epsilon_r);
935 PS("epsilon-r", infbuf);
937 if (ir->epsilon_rf != 0)
939 PR("epsilon-rf", ir->epsilon_rf);
943 PS("epsilon-rf", infbuf);
945 PS("vdw-type", EVDWTYPE(ir->vdwtype));
946 PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
947 PR("rvdw-switch", ir->rvdw_switch);
948 PR("rvdw", ir->rvdw);
949 PS("DispCorr", EDISPCORR(ir->eDispCorr));
950 PR("table-extension", ir->tabext);
952 PR("fourierspacing", ir->fourier_spacing);
953 PI("fourier-nx", ir->nkx);
954 PI("fourier-ny", ir->nky);
955 PI("fourier-nz", ir->nkz);
956 PI("pme-order", ir->pme_order);
957 PR("ewald-rtol", ir->ewald_rtol);
958 PR("ewald-rtol-lj", ir->ewald_rtol_lj);
959 PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir->ljpme_combination_rule));
960 PR("ewald-geometry", ir->ewald_geometry);
961 PR("epsilon-surface", ir->epsilon_surface);
963 /* Implicit solvent */
964 PS("implicit-solvent", EIMPLICITSOL(ir->implicit_solvent));
966 /* Generalized born electrostatics */
967 PS("gb-algorithm", EGBALGORITHM(ir->gb_algorithm));
968 PI("nstgbradii", ir->nstgbradii);
969 PR("rgbradii", ir->rgbradii);
970 PR("gb-epsilon-solvent", ir->gb_epsilon_solvent);
971 PR("gb-saltconc", ir->gb_saltconc);
972 PR("gb-obc-alpha", ir->gb_obc_alpha);
973 PR("gb-obc-beta", ir->gb_obc_beta);
974 PR("gb-obc-gamma", ir->gb_obc_gamma);
975 PR("gb-dielectric-offset", ir->gb_dielectric_offset);
976 PS("sa-algorithm", ESAALGORITHM(ir->gb_algorithm));
977 PR("sa-surface-tension", ir->sa_surface_tension);
979 /* Options for weak coupling algorithms */
980 PS("tcoupl", ETCOUPLTYPE(ir->etc));
981 PI("nsttcouple", ir->nsttcouple);
982 PI("nh-chain-length", ir->opts.nhchainlength);
983 PS("print-nose-hoover-chain-variables", EBOOL(ir->bPrintNHChains));
985 PS("pcoupl", EPCOUPLTYPE(ir->epc));
986 PS("pcoupltype", EPCOUPLTYPETYPE(ir->epct));
987 PI("nstpcouple", ir->nstpcouple);
988 PR("tau-p", ir->tau_p);
989 pr_matrix(fp, indent, "compressibility", ir->compress, bMDPformat);
990 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
991 PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
995 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX],
996 ir->posres_com[YY], ir->posres_com[ZZ]);
997 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX],
998 ir->posres_comB[YY], ir->posres_comB[ZZ]);
1002 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
1003 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
1007 PS("QMMM", EBOOL(ir->bQMMM));
1008 PI("QMconstraints", ir->QMconstraints);
1009 PI("QMMMscheme", ir->QMMMscheme);
1010 PR("MMChargeScaleFactor", ir->scalefactor);
1011 pr_qm_opts(fp, indent, "qm-opts", &(ir->opts));
1013 /* CONSTRAINT OPTIONS */
1014 PS("constraint-algorithm", ECONSTRTYPE(ir->eConstrAlg));
1015 PS("continuation", EBOOL(ir->bContinuation));
1017 PS("Shake-SOR", EBOOL(ir->bShakeSOR));
1018 PR("shake-tol", ir->shake_tol);
1019 PI("lincs-order", ir->nProjOrder);
1020 PI("lincs-iter", ir->nLincsIter);
1021 PR("lincs-warnangle", ir->LincsWarnAngle);
1024 PI("nwall", ir->nwall);
1025 PS("wall-type", EWALLTYPE(ir->wall_type));
1026 PR("wall-r-linpot", ir->wall_r_linpot);
1028 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
1029 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
1031 PR("wall-density[0]", ir->wall_density[0]);
1032 PR("wall-density[1]", ir->wall_density[1]);
1033 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
1036 PS("pull", EBOOL(ir->bPull));
1039 pr_pull(fp, indent, ir->pull);
1042 /* ENFORCED ROTATION */
1043 PS("rotation", EBOOL(ir->bRot));
1046 pr_rot(fp, indent, ir->rot);
1049 /* INTERACTIVE MD */
1050 PS("interactiveMD", EBOOL(ir->bIMD));
1053 pr_imd(fp, indent, ir->imd);
1056 /* NMR refinement stuff */
1057 PS("disre", EDISRETYPE(ir->eDisre));
1058 PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
1059 PS("disre-mixed", EBOOL(ir->bDisreMixed));
1060 PR("dr-fc", ir->dr_fc);
1061 PR("dr-tau", ir->dr_tau);
1062 PR("nstdisreout", ir->nstdisreout);
1064 PR("orire-fc", ir->orires_fc);
1065 PR("orire-tau", ir->orires_tau);
1066 PR("nstorireout", ir->nstorireout);
1068 /* FREE ENERGY VARIABLES */
1069 PS("free-energy", EFEPTYPE(ir->efep));
1070 if (ir->efep != efepNO || ir->bSimTemp)
1072 pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
1076 pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
1079 /* NON-equilibrium MD stuff */
1080 PR("cos-acceleration", ir->cos_accel);
1081 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
1083 /* SIMULATED TEMPERING */
1084 PS("simulated-tempering", EBOOL(ir->bSimTemp));
1087 pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
1090 /* ELECTRIC FIELDS */
1091 pr_cosine(fp, indent, "E-x", &(ir->ex[XX]), bMDPformat);
1092 pr_cosine(fp, indent, "E-xt", &(ir->et[XX]), bMDPformat);
1093 pr_cosine(fp, indent, "E-y", &(ir->ex[YY]), bMDPformat);
1094 pr_cosine(fp, indent, "E-yt", &(ir->et[YY]), bMDPformat);
1095 pr_cosine(fp, indent, "E-z", &(ir->ex[ZZ]), bMDPformat);
1096 pr_cosine(fp, indent, "E-zt", &(ir->et[ZZ]), bMDPformat);
1098 /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
1099 PS("swapcoords", ESWAPTYPE(ir->eSwapCoords));
1100 if (ir->eSwapCoords != eswapNO)
1102 pr_swap(fp, indent, ir->swap);
1105 /* AdResS PARAMETERS */
1106 PS("adress", EBOOL(ir->bAdress));
1109 PS("adress-type", EADRESSTYPE(ir->adress->type));
1110 PR("adress-const-wf", ir->adress->const_wf);
1111 PR("adress-ex-width", ir->adress->ex_width);
1112 PR("adress-hy-width", ir->adress->hy_width);
1113 PR("adress-ex-forcecap", ir->adress->ex_forcecap);
1114 PS("adress-interface-correction", EADRESSICTYPE(ir->adress->icor));
1115 PS("adress-site", EADRESSSITETYPE(ir->adress->site));
1116 pr_rvec(fp, indent, "adress-reference-coords", ir->adress->refs, DIM, TRUE);
1117 PS("adress-do-hybridpairs", EBOOL(ir->adress->do_hybridpairs));
1120 /* USER-DEFINED THINGIES */
1121 PI("userint1", ir->userint1);
1122 PI("userint2", ir->userint2);
1123 PI("userint3", ir->userint3);
1124 PI("userint4", ir->userint4);
1125 PR("userreal1", ir->userreal1);
1126 PR("userreal2", ir->userreal2);
1127 PR("userreal3", ir->userreal3);
1128 PR("userreal4", ir->userreal4);
1130 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1137 static void pr_harm(FILE *fp, t_iparams *iparams, const char *r, const char *kr)
1139 fprintf(fp, "%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
1140 r, iparams->harmonic.rA, kr, iparams->harmonic.krA,
1141 r, iparams->harmonic.rB, kr, iparams->harmonic.krB);
1144 void pr_iparams(FILE *fp, t_functype ftype, t_iparams *iparams)
1147 real VA[4], VB[4], *rbcA, *rbcB;
1153 pr_harm(fp, iparams, "th", "ct");
1155 case F_CROSS_BOND_BONDS:
1156 fprintf(fp, "r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
1157 iparams->cross_bb.r1e, iparams->cross_bb.r2e,
1158 iparams->cross_bb.krr);
1160 case F_CROSS_BOND_ANGLES:
1161 fprintf(fp, "r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
1162 iparams->cross_ba.r1e, iparams->cross_ba.r2e,
1163 iparams->cross_ba.r3e, iparams->cross_ba.krt);
1165 case F_LINEAR_ANGLES:
1166 fprintf(fp, "klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
1167 iparams->linangle.klinA, iparams->linangle.aA,
1168 iparams->linangle.klinB, iparams->linangle.aB);
1170 case F_UREY_BRADLEY:
1171 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);
1173 case F_QUARTIC_ANGLES:
1174 fprintf(fp, "theta=%15.8e", iparams->qangle.theta);
1175 for (i = 0; i < 5; i++)
1177 fprintf(fp, ", c%c=%15.8e", '0'+i, iparams->qangle.c[i]);
1182 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1183 iparams->bham.a, iparams->bham.b, iparams->bham.c);
1188 pr_harm(fp, iparams, "b0", "cb");
1191 pr_harm(fp, iparams, "xi", "cx");
1194 fprintf(fp, "b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
1195 iparams->morse.b0A, iparams->morse.cbA, iparams->morse.betaA,
1196 iparams->morse.b0B, iparams->morse.cbB, iparams->morse.betaB);
1199 fprintf(fp, "b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
1200 iparams->cubic.b0, iparams->cubic.kb, iparams->cubic.kcub);
1206 fprintf(fp, "bm=%15.8e, kb=%15.8e\n", iparams->fene.bm, iparams->fene.kb);
1209 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",
1210 iparams->restraint.lowA, iparams->restraint.up1A,
1211 iparams->restraint.up2A, iparams->restraint.kA,
1212 iparams->restraint.lowB, iparams->restraint.up1B,
1213 iparams->restraint.up2B, iparams->restraint.kB);
1219 fprintf(fp, "tab=%d, kA=%15.8e, kB=%15.8e\n",
1220 iparams->tab.table, iparams->tab.kA, iparams->tab.kB);
1222 case F_POLARIZATION:
1223 fprintf(fp, "alpha=%15.8e\n", iparams->polarize.alpha);
1226 fprintf(fp, "alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
1227 iparams->anharm_polarize.alpha,
1228 iparams->anharm_polarize.drcut,
1229 iparams->anharm_polarize.khyp);
1232 fprintf(fp, "a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
1233 iparams->thole.a, iparams->thole.alpha1, iparams->thole.alpha2,
1234 iparams->thole.rfac);
1237 fprintf(fp, "al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
1238 iparams->wpol.al_x, iparams->wpol.al_y, iparams->wpol.al_z,
1239 iparams->wpol.rOH, iparams->wpol.rHH, iparams->wpol.rOD);
1242 fprintf(fp, "c6=%15.8e, c12=%15.8e\n", iparams->lj.c6, iparams->lj.c12);
1245 fprintf(fp, "c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
1246 iparams->lj14.c6A, iparams->lj14.c12A,
1247 iparams->lj14.c6B, iparams->lj14.c12B);
1250 fprintf(fp, "fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1252 iparams->ljc14.qi, iparams->ljc14.qj,
1253 iparams->ljc14.c6, iparams->ljc14.c12);
1255 case F_LJC_PAIRS_NB:
1256 fprintf(fp, "qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1257 iparams->ljcnb.qi, iparams->ljcnb.qj,
1258 iparams->ljcnb.c6, iparams->ljcnb.c12);
1264 fprintf(fp, "phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
1265 iparams->pdihs.phiA, iparams->pdihs.cpA,
1266 iparams->pdihs.phiB, iparams->pdihs.cpB,
1267 iparams->pdihs.mult);
1270 fprintf(fp, "label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
1271 iparams->disres.label, iparams->disres.type,
1272 iparams->disres.low, iparams->disres.up1,
1273 iparams->disres.up2, iparams->disres.kfac);
1276 fprintf(fp, "ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
1277 iparams->orires.ex, iparams->orires.label, iparams->orires.power,
1278 iparams->orires.c, iparams->orires.obs, iparams->orires.kfac);
1281 fprintf(fp, "phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
1282 iparams->dihres.phiA, iparams->dihres.dphiA, iparams->dihres.kfacA,
1283 iparams->dihres.phiB, iparams->dihres.dphiB, iparams->dihres.kfacB);
1286 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",
1287 iparams->posres.pos0A[XX], iparams->posres.pos0A[YY],
1288 iparams->posres.pos0A[ZZ], iparams->posres.fcA[XX],
1289 iparams->posres.fcA[YY], iparams->posres.fcA[ZZ],
1290 iparams->posres.pos0B[XX], iparams->posres.pos0B[YY],
1291 iparams->posres.pos0B[ZZ], iparams->posres.fcB[XX],
1292 iparams->posres.fcB[YY], iparams->posres.fcB[ZZ]);
1295 fprintf(fp, "pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e\n",
1296 iparams->fbposres.pos0[XX], iparams->fbposres.pos0[YY],
1297 iparams->fbposres.pos0[ZZ], iparams->fbposres.geom,
1298 iparams->fbposres.r, iparams->fbposres.k);
1301 for (i = 0; i < NR_RBDIHS; i++)
1303 fprintf(fp, "%srbcA[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcA[i]);
1306 for (i = 0; i < NR_RBDIHS; i++)
1308 fprintf(fp, "%srbcB[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcB[i]);
1313 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
1314 * OPLS potential constants back.
1316 rbcA = iparams->rbdihs.rbcA;
1317 rbcB = iparams->rbdihs.rbcB;
1319 VA[3] = -0.25*rbcA[4];
1320 VA[2] = -0.5*rbcA[3];
1321 VA[1] = 4.0*VA[3]-rbcA[2];
1322 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1324 VB[3] = -0.25*rbcB[4];
1325 VB[2] = -0.5*rbcB[3];
1326 VB[1] = 4.0*VB[3]-rbcB[2];
1327 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1329 for (i = 0; i < NR_FOURDIHS; i++)
1331 fprintf(fp, "%sFourA[%d]=%15.8e", i == 0 ? "" : ", ", i, VA[i]);
1334 for (i = 0; i < NR_FOURDIHS; i++)
1336 fprintf(fp, "%sFourB[%d]=%15.8e", i == 0 ? "" : ", ", i, VB[i]);
1343 fprintf(fp, "dA=%15.8e, dB=%15.8e\n", iparams->constr.dA, iparams->constr.dB);
1346 fprintf(fp, "doh=%15.8e, dhh=%15.8e\n", iparams->settle.doh,
1347 iparams->settle.dhh);
1350 fprintf(fp, "a=%15.8e\n", iparams->vsite.a);
1355 fprintf(fp, "a=%15.8e, b=%15.8e\n", iparams->vsite.a, iparams->vsite.b);
1360 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1361 iparams->vsite.a, iparams->vsite.b, iparams->vsite.c);
1364 fprintf(fp, "n=%2d, a=%15.8e\n", iparams->vsiten.n, iparams->vsiten.a);
1369 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);
1372 fprintf(fp, "cmapA=%1d, cmapB=%1d\n", iparams->cmap.cmapA, iparams->cmap.cmapB);
1375 pr_harm(fp, iparams, "ktheta", "costheta0");
1378 fprintf(fp, "phiA=%15.8e, cpA=%15.8e",
1379 iparams->pdihs.phiA, iparams->pdihs.cpA);
1382 fprintf(fp, "kphi=%15.8e", iparams->cbtdihs.cbtcA[0]);
1383 for (i = 1; i < NR_CBTDIHS; i++)
1385 fprintf(fp, ", cbtcA[%d]=%15.8e", i-1, iparams->cbtdihs.cbtcA[i]);
1390 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1391 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1395 void pr_ilist(FILE *fp, int indent, const char *title,
1396 t_functype *functype, t_ilist *ilist, gmx_bool bShowNumbers)
1398 int i, j, k, type, ftype;
1401 if (available(fp, ilist, indent, title) && ilist->nr > 0)
1403 indent = pr_title(fp, indent, title);
1404 (void) pr_indent(fp, indent);
1405 fprintf(fp, "nr: %d\n", ilist->nr);
1408 (void) pr_indent(fp, indent);
1409 fprintf(fp, "iatoms:\n");
1410 iatoms = ilist->iatoms;
1411 for (i = j = 0; i < ilist->nr; )
1414 (void) pr_indent(fp, indent+INDENT);
1416 ftype = functype[type];
1417 (void) fprintf(fp, "%d type=%d (%s)",
1418 bShowNumbers ? j : -1, bShowNumbers ? type : -1,
1419 interaction_function[ftype].name);
1421 for (k = 0; k < interaction_function[ftype].nratoms; k++)
1423 (void) fprintf(fp, " %d", *(iatoms++));
1425 (void) fprintf(fp, "\n");
1426 i += 1+interaction_function[ftype].nratoms;
1428 fprintf(fp, "%5d%5d\n", i, iatoms[i]);
1436 static void pr_cmap(FILE *fp, int indent, const char *title,
1437 gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1442 dx = 360.0 / cmap_grid->grid_spacing;
1443 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1445 if (available(fp, cmap_grid, indent, title))
1447 fprintf(fp, "%s\n", title);
1449 for (i = 0; i < cmap_grid->ngrid; i++)
1452 fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
1454 fprintf(fp, "grid[%3d]={\n", bShowNumbers ? i : -1);
1456 for (j = 0; j < nelem; j++)
1458 if ( (j%cmap_grid->grid_spacing) == 0)
1460 fprintf(fp, "%8.1f\n", idx);
1464 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4]);
1465 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+1]);
1466 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+2]);
1467 fprintf(fp, "%8.3f\n", cmap_grid->cmapdata[i].cmap[j*4+3]);
1475 void pr_ffparams(FILE *fp, int indent, const char *title,
1476 gmx_ffparams_t *ffparams,
1477 gmx_bool bShowNumbers)
1481 indent = pr_title(fp, indent, title);
1482 (void) pr_indent(fp, indent);
1483 (void) fprintf(fp, "atnr=%d\n", ffparams->atnr);
1484 (void) pr_indent(fp, indent);
1485 (void) fprintf(fp, "ntypes=%d\n", ffparams->ntypes);
1486 for (i = 0; i < ffparams->ntypes; i++)
1488 (void) pr_indent(fp, indent+INDENT);
1489 (void) fprintf(fp, "functype[%d]=%s, ",
1490 bShowNumbers ? i : -1,
1491 interaction_function[ffparams->functype[i]].name);
1492 pr_iparams(fp, ffparams->functype[i], &ffparams->iparams[i]);
1494 (void) pr_double(fp, indent, "reppow", ffparams->reppow);
1495 (void) pr_real(fp, indent, "fudgeQQ", ffparams->fudgeQQ);
1496 pr_cmap(fp, indent, "cmap", &ffparams->cmap_grid, bShowNumbers);
1499 void pr_idef(FILE *fp, int indent, const char *title, t_idef *idef, gmx_bool bShowNumbers)
1503 if (available(fp, idef, indent, title))
1505 indent = pr_title(fp, indent, title);
1506 (void) pr_indent(fp, indent);
1507 (void) fprintf(fp, "atnr=%d\n", idef->atnr);
1508 (void) pr_indent(fp, indent);
1509 (void) fprintf(fp, "ntypes=%d\n", idef->ntypes);
1510 for (i = 0; i < idef->ntypes; i++)
1512 (void) pr_indent(fp, indent+INDENT);
1513 (void) fprintf(fp, "functype[%d]=%s, ",
1514 bShowNumbers ? i : -1,
1515 interaction_function[idef->functype[i]].name);
1516 pr_iparams(fp, idef->functype[i], &idef->iparams[i]);
1518 (void) pr_real(fp, indent, "fudgeQQ", idef->fudgeQQ);
1520 for (j = 0; (j < F_NRE); j++)
1522 pr_ilist(fp, indent, interaction_function[j].longname,
1523 idef->functype, &idef->il[j], bShowNumbers);
1528 static int pr_block_title(FILE *fp, int indent, const char *title, t_block *block)
1532 if (available(fp, block, indent, title))
1534 indent = pr_title(fp, indent, title);
1535 (void) pr_indent(fp, indent);
1536 (void) fprintf(fp, "nr=%d\n", block->nr);
1541 static int pr_blocka_title(FILE *fp, int indent, const char *title, t_blocka *block)
1545 if (available(fp, block, indent, title))
1547 indent = pr_title(fp, indent, title);
1548 (void) pr_indent(fp, indent);
1549 (void) fprintf(fp, "nr=%d\n", block->nr);
1550 (void) pr_indent(fp, indent);
1551 (void) fprintf(fp, "nra=%d\n", block->nra);
1556 static void low_pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1560 if (available(fp, block, indent, title))
1562 indent = pr_blocka_title(fp, indent, title, block);
1563 for (i = 0; i <= block->nr; i++)
1565 (void) pr_indent(fp, indent+INDENT);
1566 (void) fprintf(fp, "%s->index[%d]=%d\n",
1567 title, bShowNumbers ? i : -1, block->index[i]);
1569 for (i = 0; i < block->nra; i++)
1571 (void) pr_indent(fp, indent+INDENT);
1572 (void) fprintf(fp, "%s->a[%d]=%d\n",
1573 title, bShowNumbers ? i : -1, block->a[i]);
1578 void pr_block(FILE *fp, int indent, const char *title, t_block *block, gmx_bool bShowNumbers)
1580 int i, j, ok, size, start, end;
1582 if (available(fp, block, indent, title))
1584 indent = pr_block_title(fp, indent, title, block);
1587 if ((ok = (block->index[start] == 0)) == 0)
1589 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1593 for (i = 0; i < block->nr; i++)
1595 end = block->index[i+1];
1596 size = pr_indent(fp, indent);
1599 size += fprintf(fp, "%s[%d]={}\n", title, i);
1603 size += fprintf(fp, "%s[%d]={%d..%d}\n",
1604 title, bShowNumbers ? i : -1,
1605 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1613 void pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1615 int i, j, ok, size, start, end;
1617 if (available(fp, block, indent, title))
1619 indent = pr_blocka_title(fp, indent, title, block);
1622 if ((ok = (block->index[start] == 0)) == 0)
1624 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1628 for (i = 0; i < block->nr; i++)
1630 end = block->index[i+1];
1631 size = pr_indent(fp, indent);
1634 size += fprintf(fp, "%s[%d]={", title, i);
1638 size += fprintf(fp, "%s[%d][%d..%d]={",
1639 title, bShowNumbers ? i : -1,
1640 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1642 for (j = start; j < end; j++)
1646 size += fprintf(fp, ", ");
1648 if ((size) > (USE_WIDTH))
1650 (void) fprintf(fp, "\n");
1651 size = pr_indent(fp, indent+INDENT);
1653 size += fprintf(fp, "%d", block->a[j]);
1655 (void) fprintf(fp, "}\n");
1659 if ((end != block->nra) || (!ok))
1661 (void) pr_indent(fp, indent);
1662 (void) fprintf(fp, "tables inconsistent, dumping complete tables:\n");
1663 low_pr_blocka(fp, indent, title, block, bShowNumbers);
1668 static void pr_strings(FILE *fp, int indent, const char *title, char ***nm, int n, gmx_bool bShowNumbers)
1672 if (available(fp, nm, indent, title))
1674 indent = pr_title_n(fp, indent, title, n);
1675 for (i = 0; i < n; i++)
1677 (void) pr_indent(fp, indent);
1678 (void) fprintf(fp, "%s[%d]={name=\"%s\"}\n",
1679 title, bShowNumbers ? i : -1, *(nm[i]));
1684 static void pr_strings2(FILE *fp, int indent, const char *title,
1685 char ***nm, char ***nmB, int n, gmx_bool bShowNumbers)
1689 if (available(fp, nm, indent, title))
1691 indent = pr_title_n(fp, indent, title, n);
1692 for (i = 0; i < n; i++)
1694 (void) pr_indent(fp, indent);
1695 (void) fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1696 title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
1701 static void pr_resinfo(FILE *fp, int indent, const char *title, t_resinfo *resinfo, int n, gmx_bool bShowNumbers)
1705 if (available(fp, resinfo, indent, title))
1707 indent = pr_title_n(fp, indent, title, n);
1708 for (i = 0; i < n; i++)
1710 (void) pr_indent(fp, indent);
1711 (void) fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1712 title, bShowNumbers ? i : -1,
1713 *(resinfo[i].name), resinfo[i].nr,
1714 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1719 static void pr_atom(FILE *fp, int indent, const char *title, t_atom *atom, int n)
1723 if (available(fp, atom, indent, title))
1725 indent = pr_title_n(fp, indent, title, n);
1726 for (i = 0; i < n; i++)
1728 (void) pr_indent(fp, indent);
1729 fprintf(fp, "%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1730 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1731 title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype],
1732 atom[i].m, atom[i].q, atom[i].mB, atom[i].qB,
1733 atom[i].resind, atom[i].atomnumber);
1738 static void pr_grps(FILE *fp, const char *title, t_grps grps[], char **grpname[])
1742 for (i = 0; (i < egcNR); i++)
1744 fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
1745 for (j = 0; (j < grps[i].nr); j++)
1747 fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
1753 static void pr_groups(FILE *fp, int indent,
1754 gmx_groups_t *groups,
1755 gmx_bool bShowNumbers)
1760 pr_grps(fp, "grp", groups->grps, groups->grpname);
1761 pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
1763 (void) pr_indent(fp, indent);
1764 fprintf(fp, "groups ");
1765 for (g = 0; g < egcNR; g++)
1767 printf(" %5.5s", gtypes[g]);
1771 (void) pr_indent(fp, indent);
1772 fprintf(fp, "allocated ");
1774 for (g = 0; g < egcNR; g++)
1776 printf(" %5d", groups->ngrpnr[g]);
1777 nat_max = max(nat_max, groups->ngrpnr[g]);
1783 (void) pr_indent(fp, indent);
1784 fprintf(fp, "groupnr[%5s] =", "*");
1785 for (g = 0; g < egcNR; g++)
1787 fprintf(fp, " %3d ", 0);
1793 for (i = 0; i < nat_max; i++)
1795 (void) pr_indent(fp, indent);
1796 fprintf(fp, "groupnr[%5d] =", i);
1797 for (g = 0; g < egcNR; g++)
1799 fprintf(fp, " %3d ",
1800 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1807 void pr_atoms(FILE *fp, int indent, const char *title, t_atoms *atoms,
1808 gmx_bool bShownumbers)
1810 if (available(fp, atoms, indent, title))
1812 indent = pr_title(fp, indent, title);
1813 pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
1814 pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
1815 pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
1816 pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
1821 void pr_atomtypes(FILE *fp, int indent, const char *title, t_atomtypes *atomtypes,
1822 gmx_bool bShowNumbers)
1825 if (available(fp, atomtypes, indent, title))
1827 indent = pr_title(fp, indent, title);
1828 for (i = 0; i < atomtypes->nr; i++)
1830 pr_indent(fp, indent);
1832 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1833 bShowNumbers ? i : -1, atomtypes->radius[i], atomtypes->vol[i],
1834 atomtypes->gb_radius[i],
1835 atomtypes->surftens[i], atomtypes->atomnumber[i], atomtypes->S_hct[i]);
1840 static void pr_moltype(FILE *fp, int indent, const char *title,
1841 gmx_moltype_t *molt, int n,
1842 gmx_ffparams_t *ffparams,
1843 gmx_bool bShowNumbers)
1847 indent = pr_title_n(fp, indent, title, n);
1848 (void) pr_indent(fp, indent);
1849 (void) fprintf(fp, "name=\"%s\"\n", *(molt->name));
1850 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
1851 pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
1852 pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
1853 for (j = 0; (j < F_NRE); j++)
1855 pr_ilist(fp, indent, interaction_function[j].longname,
1856 ffparams->functype, &molt->ilist[j], bShowNumbers);
1860 static void pr_molblock(FILE *fp, int indent, const char *title,
1861 gmx_molblock_t *molb, int n,
1862 gmx_moltype_t *molt)
1864 indent = pr_title_n(fp, indent, title, n);
1865 (void) pr_indent(fp, indent);
1866 (void) fprintf(fp, "%-20s = %d \"%s\"\n",
1867 "moltype", molb->type, *(molt[molb->type].name));
1868 pr_int(fp, indent, "#molecules", molb->nmol);
1869 pr_int(fp, indent, "#atoms_mol", molb->natoms_mol);
1870 pr_int(fp, indent, "#posres_xA", molb->nposres_xA);
1871 if (molb->nposres_xA > 0)
1873 pr_rvecs(fp, indent, "posres_xA", molb->posres_xA, molb->nposres_xA);
1875 pr_int(fp, indent, "#posres_xB", molb->nposres_xB);
1876 if (molb->nposres_xB > 0)
1878 pr_rvecs(fp, indent, "posres_xB", molb->posres_xB, molb->nposres_xB);
1882 void pr_mtop(FILE *fp, int indent, const char *title, gmx_mtop_t *mtop,
1883 gmx_bool bShowNumbers)
1887 if (available(fp, mtop, indent, title))
1889 indent = pr_title(fp, indent, title);
1890 (void) pr_indent(fp, indent);
1891 (void) fprintf(fp, "name=\"%s\"\n", *(mtop->name));
1892 pr_int(fp, indent, "#atoms", mtop->natoms);
1893 pr_int(fp, indent, "#molblock", mtop->nmolblock);
1894 for (mb = 0; mb < mtop->nmolblock; mb++)
1896 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
1898 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
1899 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
1900 for (mt = 0; mt < mtop->nmoltype; mt++)
1902 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
1903 &mtop->ffparams, bShowNumbers);
1905 pr_groups(fp, indent, &mtop->groups, bShowNumbers);
1909 void pr_top(FILE *fp, int indent, const char *title, t_topology *top, gmx_bool bShowNumbers)
1911 if (available(fp, top, indent, title))
1913 indent = pr_title(fp, indent, title);
1914 (void) pr_indent(fp, indent);
1915 (void) fprintf(fp, "name=\"%s\"\n", *(top->name));
1916 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
1917 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
1918 pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
1919 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
1920 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
1921 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers);
1925 void pr_header(FILE *fp, int indent, const char *title, t_tpxheader *sh)
1929 if (available(fp, sh, indent, title))
1931 indent = pr_title(fp, indent, title);
1932 pr_indent(fp, indent);
1933 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
1934 pr_indent(fp, indent);
1935 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
1936 pr_indent(fp, indent);
1937 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
1938 pr_indent(fp, indent);
1939 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
1940 pr_indent(fp, indent);
1941 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
1942 pr_indent(fp, indent);
1943 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
1945 pr_indent(fp, indent);
1946 fprintf(fp, "natoms = %d\n", sh->natoms);
1947 pr_indent(fp, indent);
1948 fprintf(fp, "lambda = %e\n", sh->lambda);
1952 void pr_commrec(FILE *fp, int indent, t_commrec *cr)
1954 pr_indent(fp, indent);
1955 fprintf(fp, "commrec:\n");
1957 pr_indent(fp, indent);
1958 fprintf(fp, "rank = %d\n", cr->nodeid);
1959 pr_indent(fp, indent);
1960 fprintf(fp, "number of ranks = %d\n", cr->nnodes);
1961 pr_indent(fp, indent);
1962 fprintf(fp, "PME-only ranks = %d\n", cr->npmenodes);
1964 pr_indent(fp,indent);
1965 fprintf(fp,"threadid = %d\n",cr->threadid);
1966 pr_indent(fp,indent);
1967 fprintf(fp,"nthreads = %d\n",cr->nthreads);