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 PS("type", EPULLTYPE(pcrd->eType));
638 PS("geometry", EPULLGEOM(pcrd->eGeom));
639 pr_ivec(fp, indent, "dim", pcrd->dim, DIM, TRUE);
640 pr_rvec(fp, indent, "origin", pcrd->origin, DIM, TRUE);
641 pr_rvec(fp, indent, "vec", pcrd->vec, DIM, TRUE);
642 PS("start", EBOOL(pcrd->bStart));
643 PR("init", pcrd->init);
644 PR("rate", pcrd->rate);
649 static void pr_simtempvals(FILE *fp, int indent, t_simtemp *simtemp, int n_lambda)
651 PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
652 PR("sim-temp-low", simtemp->simtemp_low);
653 PR("sim-temp-high", simtemp->simtemp_high);
654 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
657 static void pr_expandedvals(FILE *fp, int indent, t_expanded *expand, int n_lambda)
660 PI("nstexpanded", expand->nstexpanded);
661 PS("lmc-stats", elamstats_names[expand->elamstats]);
662 PS("lmc-move", elmcmove_names[expand->elmcmove]);
663 PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
664 if (expand->elmceq == elmceqNUMATLAM)
666 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
668 if (expand->elmceq == elmceqSAMPLES)
670 PI("weight-equil-number-samples", expand->equil_samples);
672 if (expand->elmceq == elmceqSTEPS)
674 PI("weight-equil-number-steps", expand->equil_steps);
676 if (expand->elmceq == elmceqWLDELTA)
678 PR("weight-equil-wl-delta", expand->equil_wl_delta);
680 if (expand->elmceq == elmceqRATIO)
682 PR("weight-equil-count-ratio", expand->equil_ratio);
684 PI("lmc-seed", expand->lmc_seed);
685 PR("mc-temperature", expand->mc_temp);
686 PI("lmc-repeats", expand->lmc_repeats);
687 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
688 PI("lmc-forced-nstart", expand->lmc_forced_nstart);
689 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
690 PI("nst-transition-matrix", expand->nstTij);
691 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
692 PI("weight-c-range", expand->c_range); /* default is just C=0 */
693 PR("wl-scale", expand->wl_scale);
694 PR("wl-ratio", expand->wl_ratio);
695 PR("init-wl-delta", expand->init_wl_delta);
696 PS("wl-oneovert", EBOOL(expand->bWLoneovert));
698 pr_indent(fp, indent);
699 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
700 PS("init-weights", EBOOL(expand->bInit_weights));
703 static void pr_fepvals(FILE *fp, int indent, t_lambda *fep, gmx_bool bMDPformat)
707 PR("init-lambda", fep->init_lambda);
708 PI("init-lambda-state", fep->init_fep_state);
709 PR("delta-lambda", fep->delta_lambda);
710 PI("nstdhdl", fep->nstdhdl);
714 PI("n-lambdas", fep->n_lambda);
716 if (fep->n_lambda > 0)
718 pr_indent(fp, indent);
719 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
720 for (i = 0; i < efptNR; i++)
722 fprintf(fp, "%18s = ", efpt_names[i]);
723 if (fep->separate_dvdl[i])
725 fprintf(fp, " TRUE");
729 fprintf(fp, " FALSE");
733 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
734 for (i = 0; i < efptNR; i++)
736 fprintf(fp, "%18s = ", efpt_names[i]);
737 for (j = 0; j < fep->n_lambda; j++)
739 fprintf(fp, " %10g", fep->all_lambda[i][j]);
744 PI("calc-lambda-neighbors", fep->lambda_neighbors);
745 PS("dhdl-print-energy", edHdLPrintEnergy_names[fep->edHdLPrintEnergy]);
746 PR("sc-alpha", fep->sc_alpha);
747 PI("sc-power", fep->sc_power);
748 PR("sc-r-power", fep->sc_r_power);
749 PR("sc-sigma", fep->sc_sigma);
750 PR("sc-sigma-min", fep->sc_sigma_min);
751 PS("sc-coul", EBOOL(fep->bScCoul));
752 PI("dh-hist-size", fep->dh_hist_size);
753 PD("dh-hist-spacing", fep->dh_hist_spacing);
754 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
755 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
758 static void pr_pull(FILE *fp, int indent, t_pull *pull)
762 PR("pull-cylinder-r", pull->cylinder_r);
763 PR("pull-constr-tol", pull->constr_tol);
764 PS("pull-print-COM1", EBOOL(pull->bPrintCOM1));
765 PS("pull-print-COM2", EBOOL(pull->bPrintCOM2));
766 PS("pull-print-ref-value", EBOOL(pull->bPrintRefValue));
767 PS("pull-print-components", EBOOL(pull->bPrintComp));
768 PI("pull-nstxout", pull->nstxout);
769 PI("pull-nstfout", pull->nstfout);
770 PI("pull-ngroups", pull->ngroup);
771 for (g = 0; g < pull->ngroup; g++)
773 pr_pull_group(fp, indent, g, &pull->group[g]);
775 PI("pull-ncoords", pull->ncoord);
776 for (g = 0; g < pull->ncoord; g++)
778 pr_pull_coord(fp, indent, g, &pull->coord[g]);
782 static void pr_rotgrp(FILE *fp, int indent, int g, t_rotgrp *rotg)
784 pr_indent(fp, indent);
785 fprintf(fp, "rot-group %d:\n", g);
787 PS("rot-type", EROTGEOM(rotg->eType));
788 PS("rot-massw", EBOOL(rotg->bMassW));
789 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
790 pr_rvecs(fp, indent, "x-ref", rotg->x_ref, rotg->nat);
791 pr_rvec(fp, indent, "rot-vec", rotg->vec, DIM, TRUE);
792 pr_rvec(fp, indent, "rot-pivot", rotg->pivot, DIM, TRUE);
793 PR("rot-rate", rotg->rate);
794 PR("rot-k", rotg->k);
795 PR("rot-slab-dist", rotg->slab_dist);
796 PR("rot-min-gauss", rotg->min_gaussian);
797 PR("rot-eps", rotg->eps);
798 PS("rot-fit-method", EROTFIT(rotg->eFittype));
799 PI("rot-potfit-nstep", rotg->PotAngle_nstep);
800 PR("rot-potfit-step", rotg->PotAngle_step);
803 static void pr_rot(FILE *fp, int indent, t_rot *rot)
807 PI("rot-nstrout", rot->nstrout);
808 PI("rot-nstsout", rot->nstsout);
809 PI("rot-ngroups", rot->ngrp);
810 for (g = 0; g < rot->ngrp; g++)
812 pr_rotgrp(fp, indent, g, &rot->grp[g]);
817 static void pr_swap(FILE *fp, int indent, t_swapcoords *swap)
823 PI("swap-frequency", swap->nstswap);
824 for (j = 0; j < 2; j++)
826 sprintf(str, "massw_split%d", j);
827 PS(str, EBOOL(swap->massw_split[j]));
828 sprintf(str, "split atoms group %d", j);
829 pr_ivec_block(fp, indent, str, swap->ind_split[j], swap->nat_split[j], TRUE);
831 pr_ivec_block(fp, indent, "swap atoms", swap->ind, swap->nat, TRUE);
832 pr_ivec_block(fp, indent, "solvent atoms", swap->ind_sol, swap->nat_sol, TRUE);
833 PR("cyl0-r", swap->cyl0r);
834 PR("cyl0-up", swap->cyl0u);
835 PR("cyl0-down", swap->cyl0l);
836 PR("cyl1-r", swap->cyl1r);
837 PR("cyl1-up", swap->cyl1u);
838 PR("cyl1-down", swap->cyl1l);
839 PI("coupl-steps", swap->nAverage);
840 for (j = 0; j < 2; j++)
842 sprintf(str, "anions%c", j+'A');
843 PI(str, swap->nanions[j]);
844 sprintf(str, "cations%c", j+'A');
845 PI(str, swap->ncations[j]);
847 PR("threshold", swap->threshold);
851 static void pr_imd(FILE *fp, int indent, t_IMD *imd)
853 PI("IMD-atoms", imd->nat);
854 pr_ivec_block(fp, indent, "atom", imd->ind, imd->nat, TRUE);
858 void pr_inputrec(FILE *fp, int indent, const char *title, t_inputrec *ir,
861 const char *infbuf = "inf";
864 if (available(fp, ir, indent, title))
868 indent = pr_title(fp, indent, title);
870 /* Try to make this list appear in the same order as the
871 * options are written in the default mdout.mdp, and with
872 * the same user-exposed names to facilitate debugging.
874 PS("integrator", EI(ir->eI));
875 PR("tinit", ir->init_t);
876 PR("dt", ir->delta_t);
877 PSTEP("nsteps", ir->nsteps);
878 PSTEP("init-step", ir->init_step);
879 PI("simulation-part", ir->simulation_part);
880 PS("comm-mode", ECOM(ir->comm_mode));
881 PI("nstcomm", ir->nstcomm);
883 /* Langevin dynamics */
884 PR("bd-fric", ir->bd_fric);
885 PSTEP("ld-seed", ir->ld_seed);
887 /* Energy minimization */
888 PR("emtol", ir->em_tol);
889 PR("emstep", ir->em_stepsize);
890 PI("niter", ir->niter);
891 PR("fcstep", ir->fc_stepsize);
892 PI("nstcgsteep", ir->nstcgsteep);
893 PI("nbfgscorr", ir->nbfgscorr);
895 /* Test particle insertion */
896 PR("rtpi", ir->rtpi);
899 PI("nstxout", ir->nstxout);
900 PI("nstvout", ir->nstvout);
901 PI("nstfout", ir->nstfout);
902 PI("nstlog", ir->nstlog);
903 PI("nstcalcenergy", ir->nstcalcenergy);
904 PI("nstenergy", ir->nstenergy);
905 PI("nstxout-compressed", ir->nstxout_compressed);
906 PR("compressed-x-precision", ir->x_compression_precision);
908 /* Neighborsearching parameters */
909 PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
910 PI("nstlist", ir->nstlist);
911 PS("ns-type", ENS(ir->ns_type));
912 PS("pbc", EPBC(ir->ePBC));
913 PS("periodic-molecules", EBOOL(ir->bPeriodicMols));
914 PR("verlet-buffer-tolerance", ir->verletbuf_tol);
915 PR("rlist", ir->rlist);
916 PR("rlistlong", ir->rlistlong);
917 PR("nstcalclr", ir->nstcalclr);
919 /* Options for electrostatics and VdW */
920 PS("coulombtype", EELTYPE(ir->coulombtype));
921 PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
922 PR("rcoulomb-switch", ir->rcoulomb_switch);
923 PR("rcoulomb", ir->rcoulomb);
924 if (ir->epsilon_r != 0)
926 PR("epsilon-r", ir->epsilon_r);
930 PS("epsilon-r", infbuf);
932 if (ir->epsilon_rf != 0)
934 PR("epsilon-rf", ir->epsilon_rf);
938 PS("epsilon-rf", infbuf);
940 PS("vdw-type", EVDWTYPE(ir->vdwtype));
941 PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
942 PR("rvdw-switch", ir->rvdw_switch);
943 PR("rvdw", ir->rvdw);
944 PS("DispCorr", EDISPCORR(ir->eDispCorr));
945 PR("table-extension", ir->tabext);
947 PR("fourierspacing", ir->fourier_spacing);
948 PI("fourier-nx", ir->nkx);
949 PI("fourier-ny", ir->nky);
950 PI("fourier-nz", ir->nkz);
951 PI("pme-order", ir->pme_order);
952 PR("ewald-rtol", ir->ewald_rtol);
953 PR("ewald-rtol-lj", ir->ewald_rtol_lj);
954 PS("lj-pme-comb-rule", ELJPMECOMBNAMES(ir->ljpme_combination_rule));
955 PR("ewald-geometry", ir->ewald_geometry);
956 PR("epsilon-surface", ir->epsilon_surface);
958 /* Implicit solvent */
959 PS("implicit-solvent", EIMPLICITSOL(ir->implicit_solvent));
961 /* Generalized born electrostatics */
962 PS("gb-algorithm", EGBALGORITHM(ir->gb_algorithm));
963 PI("nstgbradii", ir->nstgbradii);
964 PR("rgbradii", ir->rgbradii);
965 PR("gb-epsilon-solvent", ir->gb_epsilon_solvent);
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);
974 /* Options for weak coupling algorithms */
975 PS("tcoupl", ETCOUPLTYPE(ir->etc));
976 PI("nsttcouple", ir->nsttcouple);
977 PI("nh-chain-length", ir->opts.nhchainlength);
978 PS("print-nose-hoover-chain-variables", EBOOL(ir->bPrintNHChains));
980 PS("pcoupl", EPCOUPLTYPE(ir->epc));
981 PS("pcoupltype", EPCOUPLTYPETYPE(ir->epct));
982 PI("nstpcouple", ir->nstpcouple);
983 PR("tau-p", ir->tau_p);
984 pr_matrix(fp, indent, "compressibility", ir->compress, bMDPformat);
985 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
986 PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
990 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX],
991 ir->posres_com[YY], ir->posres_com[ZZ]);
992 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX],
993 ir->posres_comB[YY], ir->posres_comB[ZZ]);
997 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
998 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
1002 PS("QMMM", EBOOL(ir->bQMMM));
1003 PI("QMconstraints", ir->QMconstraints);
1004 PI("QMMMscheme", ir->QMMMscheme);
1005 PR("MMChargeScaleFactor", ir->scalefactor);
1006 pr_qm_opts(fp, indent, "qm-opts", &(ir->opts));
1008 /* CONSTRAINT OPTIONS */
1009 PS("constraint-algorithm", ECONSTRTYPE(ir->eConstrAlg));
1010 PS("continuation", EBOOL(ir->bContinuation));
1012 PS("Shake-SOR", EBOOL(ir->bShakeSOR));
1013 PR("shake-tol", ir->shake_tol);
1014 PI("lincs-order", ir->nProjOrder);
1015 PI("lincs-iter", ir->nLincsIter);
1016 PR("lincs-warnangle", ir->LincsWarnAngle);
1019 PI("nwall", ir->nwall);
1020 PS("wall-type", EWALLTYPE(ir->wall_type));
1021 PR("wall-r-linpot", ir->wall_r_linpot);
1023 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
1024 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
1026 PR("wall-density[0]", ir->wall_density[0]);
1027 PR("wall-density[1]", ir->wall_density[1]);
1028 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
1031 PS("pull", EBOOL(ir->bPull));
1034 pr_pull(fp, indent, ir->pull);
1037 /* ENFORCED ROTATION */
1038 PS("rotation", EBOOL(ir->bRot));
1041 pr_rot(fp, indent, ir->rot);
1044 /* INTERACTIVE MD */
1045 PS("interactiveMD", EBOOL(ir->bIMD));
1048 pr_imd(fp, indent, ir->imd);
1051 /* NMR refinement stuff */
1052 PS("disre", EDISRETYPE(ir->eDisre));
1053 PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
1054 PS("disre-mixed", EBOOL(ir->bDisreMixed));
1055 PR("dr-fc", ir->dr_fc);
1056 PR("dr-tau", ir->dr_tau);
1057 PR("nstdisreout", ir->nstdisreout);
1059 PR("orire-fc", ir->orires_fc);
1060 PR("orire-tau", ir->orires_tau);
1061 PR("nstorireout", ir->nstorireout);
1063 /* FREE ENERGY VARIABLES */
1064 PS("free-energy", EFEPTYPE(ir->efep));
1065 if (ir->efep != efepNO || ir->bSimTemp)
1067 pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
1071 pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
1074 /* NON-equilibrium MD stuff */
1075 PR("cos-acceleration", ir->cos_accel);
1076 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
1078 /* SIMULATED TEMPERING */
1079 PS("simulated-tempering", EBOOL(ir->bSimTemp));
1082 pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
1085 /* ELECTRIC FIELDS */
1086 pr_cosine(fp, indent, "E-x", &(ir->ex[XX]), bMDPformat);
1087 pr_cosine(fp, indent, "E-xt", &(ir->et[XX]), bMDPformat);
1088 pr_cosine(fp, indent, "E-y", &(ir->ex[YY]), bMDPformat);
1089 pr_cosine(fp, indent, "E-yt", &(ir->et[YY]), bMDPformat);
1090 pr_cosine(fp, indent, "E-z", &(ir->ex[ZZ]), bMDPformat);
1091 pr_cosine(fp, indent, "E-zt", &(ir->et[ZZ]), bMDPformat);
1093 /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
1094 PS("swapcoords", ESWAPTYPE(ir->eSwapCoords));
1095 if (ir->eSwapCoords != eswapNO)
1097 pr_swap(fp, indent, ir->swap);
1100 /* AdResS PARAMETERS */
1101 PS("adress", EBOOL(ir->bAdress));
1104 PS("adress-type", EADRESSTYPE(ir->adress->type));
1105 PR("adress-const-wf", ir->adress->const_wf);
1106 PR("adress-ex-width", ir->adress->ex_width);
1107 PR("adress-hy-width", ir->adress->hy_width);
1108 PR("adress-ex-forcecap", ir->adress->ex_forcecap);
1109 PS("adress-interface-correction", EADRESSICTYPE(ir->adress->icor));
1110 PS("adress-site", EADRESSSITETYPE(ir->adress->site));
1111 pr_rvec(fp, indent, "adress-reference-coords", ir->adress->refs, DIM, TRUE);
1112 PS("adress-do-hybridpairs", EBOOL(ir->adress->do_hybridpairs));
1115 /* USER-DEFINED THINGIES */
1116 PI("userint1", ir->userint1);
1117 PI("userint2", ir->userint2);
1118 PI("userint3", ir->userint3);
1119 PI("userint4", ir->userint4);
1120 PR("userreal1", ir->userreal1);
1121 PR("userreal2", ir->userreal2);
1122 PR("userreal3", ir->userreal3);
1123 PR("userreal4", ir->userreal4);
1125 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
1132 static void pr_harm(FILE *fp, t_iparams *iparams, const char *r, const char *kr)
1134 fprintf(fp, "%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
1135 r, iparams->harmonic.rA, kr, iparams->harmonic.krA,
1136 r, iparams->harmonic.rB, kr, iparams->harmonic.krB);
1139 void pr_iparams(FILE *fp, t_functype ftype, t_iparams *iparams)
1142 real VA[4], VB[4], *rbcA, *rbcB;
1148 pr_harm(fp, iparams, "th", "ct");
1150 case F_CROSS_BOND_BONDS:
1151 fprintf(fp, "r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
1152 iparams->cross_bb.r1e, iparams->cross_bb.r2e,
1153 iparams->cross_bb.krr);
1155 case F_CROSS_BOND_ANGLES:
1156 fprintf(fp, "r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
1157 iparams->cross_ba.r1e, iparams->cross_ba.r2e,
1158 iparams->cross_ba.r3e, iparams->cross_ba.krt);
1160 case F_LINEAR_ANGLES:
1161 fprintf(fp, "klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
1162 iparams->linangle.klinA, iparams->linangle.aA,
1163 iparams->linangle.klinB, iparams->linangle.aB);
1165 case F_UREY_BRADLEY:
1166 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);
1168 case F_QUARTIC_ANGLES:
1169 fprintf(fp, "theta=%15.8e", iparams->qangle.theta);
1170 for (i = 0; i < 5; i++)
1172 fprintf(fp, ", c%c=%15.8e", '0'+i, iparams->qangle.c[i]);
1177 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1178 iparams->bham.a, iparams->bham.b, iparams->bham.c);
1183 pr_harm(fp, iparams, "b0", "cb");
1186 pr_harm(fp, iparams, "xi", "cx");
1189 fprintf(fp, "b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
1190 iparams->morse.b0A, iparams->morse.cbA, iparams->morse.betaA,
1191 iparams->morse.b0B, iparams->morse.cbB, iparams->morse.betaB);
1194 fprintf(fp, "b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
1195 iparams->cubic.b0, iparams->cubic.kb, iparams->cubic.kcub);
1201 fprintf(fp, "bm=%15.8e, kb=%15.8e\n", iparams->fene.bm, iparams->fene.kb);
1204 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",
1205 iparams->restraint.lowA, iparams->restraint.up1A,
1206 iparams->restraint.up2A, iparams->restraint.kA,
1207 iparams->restraint.lowB, iparams->restraint.up1B,
1208 iparams->restraint.up2B, iparams->restraint.kB);
1214 fprintf(fp, "tab=%d, kA=%15.8e, kB=%15.8e\n",
1215 iparams->tab.table, iparams->tab.kA, iparams->tab.kB);
1217 case F_POLARIZATION:
1218 fprintf(fp, "alpha=%15.8e\n", iparams->polarize.alpha);
1221 fprintf(fp, "alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
1222 iparams->anharm_polarize.alpha,
1223 iparams->anharm_polarize.drcut,
1224 iparams->anharm_polarize.khyp);
1227 fprintf(fp, "a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
1228 iparams->thole.a, iparams->thole.alpha1, iparams->thole.alpha2,
1229 iparams->thole.rfac);
1232 fprintf(fp, "al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
1233 iparams->wpol.al_x, iparams->wpol.al_y, iparams->wpol.al_z,
1234 iparams->wpol.rOH, iparams->wpol.rHH, iparams->wpol.rOD);
1237 fprintf(fp, "c6=%15.8e, c12=%15.8e\n", iparams->lj.c6, iparams->lj.c12);
1240 fprintf(fp, "c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
1241 iparams->lj14.c6A, iparams->lj14.c12A,
1242 iparams->lj14.c6B, iparams->lj14.c12B);
1245 fprintf(fp, "fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1247 iparams->ljc14.qi, iparams->ljc14.qj,
1248 iparams->ljc14.c6, iparams->ljc14.c12);
1250 case F_LJC_PAIRS_NB:
1251 fprintf(fp, "qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1252 iparams->ljcnb.qi, iparams->ljcnb.qj,
1253 iparams->ljcnb.c6, iparams->ljcnb.c12);
1259 fprintf(fp, "phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
1260 iparams->pdihs.phiA, iparams->pdihs.cpA,
1261 iparams->pdihs.phiB, iparams->pdihs.cpB,
1262 iparams->pdihs.mult);
1265 fprintf(fp, "label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
1266 iparams->disres.label, iparams->disres.type,
1267 iparams->disres.low, iparams->disres.up1,
1268 iparams->disres.up2, iparams->disres.kfac);
1271 fprintf(fp, "ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
1272 iparams->orires.ex, iparams->orires.label, iparams->orires.power,
1273 iparams->orires.c, iparams->orires.obs, iparams->orires.kfac);
1276 fprintf(fp, "phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
1277 iparams->dihres.phiA, iparams->dihres.dphiA, iparams->dihres.kfacA,
1278 iparams->dihres.phiB, iparams->dihres.dphiB, iparams->dihres.kfacB);
1281 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",
1282 iparams->posres.pos0A[XX], iparams->posres.pos0A[YY],
1283 iparams->posres.pos0A[ZZ], iparams->posres.fcA[XX],
1284 iparams->posres.fcA[YY], iparams->posres.fcA[ZZ],
1285 iparams->posres.pos0B[XX], iparams->posres.pos0B[YY],
1286 iparams->posres.pos0B[ZZ], iparams->posres.fcB[XX],
1287 iparams->posres.fcB[YY], iparams->posres.fcB[ZZ]);
1290 fprintf(fp, "pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e\n",
1291 iparams->fbposres.pos0[XX], iparams->fbposres.pos0[YY],
1292 iparams->fbposres.pos0[ZZ], iparams->fbposres.geom,
1293 iparams->fbposres.r, iparams->fbposres.k);
1296 for (i = 0; i < NR_RBDIHS; i++)
1298 fprintf(fp, "%srbcA[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcA[i]);
1301 for (i = 0; i < NR_RBDIHS; i++)
1303 fprintf(fp, "%srbcB[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcB[i]);
1308 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
1309 * OPLS potential constants back.
1311 rbcA = iparams->rbdihs.rbcA;
1312 rbcB = iparams->rbdihs.rbcB;
1314 VA[3] = -0.25*rbcA[4];
1315 VA[2] = -0.5*rbcA[3];
1316 VA[1] = 4.0*VA[3]-rbcA[2];
1317 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1319 VB[3] = -0.25*rbcB[4];
1320 VB[2] = -0.5*rbcB[3];
1321 VB[1] = 4.0*VB[3]-rbcB[2];
1322 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1324 for (i = 0; i < NR_FOURDIHS; i++)
1326 fprintf(fp, "%sFourA[%d]=%15.8e", i == 0 ? "" : ", ", i, VA[i]);
1329 for (i = 0; i < NR_FOURDIHS; i++)
1331 fprintf(fp, "%sFourB[%d]=%15.8e", i == 0 ? "" : ", ", i, VB[i]);
1338 fprintf(fp, "dA=%15.8e, dB=%15.8e\n", iparams->constr.dA, iparams->constr.dB);
1341 fprintf(fp, "doh=%15.8e, dhh=%15.8e\n", iparams->settle.doh,
1342 iparams->settle.dhh);
1345 fprintf(fp, "a=%15.8e\n", iparams->vsite.a);
1350 fprintf(fp, "a=%15.8e, b=%15.8e\n", iparams->vsite.a, iparams->vsite.b);
1355 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1356 iparams->vsite.a, iparams->vsite.b, iparams->vsite.c);
1359 fprintf(fp, "n=%2d, a=%15.8e\n", iparams->vsiten.n, iparams->vsiten.a);
1364 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);
1367 fprintf(fp, "cmapA=%1d, cmapB=%1d\n", iparams->cmap.cmapA, iparams->cmap.cmapB);
1370 pr_harm(fp, iparams, "ktheta", "costheta0");
1373 fprintf(fp, "phiA=%15.8e, cpA=%15.8e",
1374 iparams->pdihs.phiA, iparams->pdihs.cpA);
1377 fprintf(fp, "kphi=%15.8e", iparams->cbtdihs.cbtcA[0]);
1378 for (i = 1; i < NR_CBTDIHS; i++)
1380 fprintf(fp, ", cbtcA[%d]=%15.8e", i-1, iparams->cbtdihs.cbtcA[i]);
1385 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1386 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1390 void pr_ilist(FILE *fp, int indent, const char *title,
1391 t_functype *functype, t_ilist *ilist, gmx_bool bShowNumbers)
1393 int i, j, k, type, ftype;
1396 if (available(fp, ilist, indent, title) && ilist->nr > 0)
1398 indent = pr_title(fp, indent, title);
1399 (void) pr_indent(fp, indent);
1400 fprintf(fp, "nr: %d\n", ilist->nr);
1403 (void) pr_indent(fp, indent);
1404 fprintf(fp, "iatoms:\n");
1405 iatoms = ilist->iatoms;
1406 for (i = j = 0; i < ilist->nr; )
1409 (void) pr_indent(fp, indent+INDENT);
1411 ftype = functype[type];
1412 (void) fprintf(fp, "%d type=%d (%s)",
1413 bShowNumbers ? j : -1, bShowNumbers ? type : -1,
1414 interaction_function[ftype].name);
1416 for (k = 0; k < interaction_function[ftype].nratoms; k++)
1418 (void) fprintf(fp, " %d", *(iatoms++));
1420 (void) fprintf(fp, "\n");
1421 i += 1+interaction_function[ftype].nratoms;
1423 fprintf(fp, "%5d%5d\n", i, iatoms[i]);
1431 static void pr_cmap(FILE *fp, int indent, const char *title,
1432 gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1437 dx = 360.0 / cmap_grid->grid_spacing;
1438 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1440 if (available(fp, cmap_grid, indent, title))
1442 fprintf(fp, "%s\n", title);
1444 for (i = 0; i < cmap_grid->ngrid; i++)
1447 fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
1449 fprintf(fp, "grid[%3d]={\n", bShowNumbers ? i : -1);
1451 for (j = 0; j < nelem; j++)
1453 if ( (j%cmap_grid->grid_spacing) == 0)
1455 fprintf(fp, "%8.1f\n", idx);
1459 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4]);
1460 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+1]);
1461 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+2]);
1462 fprintf(fp, "%8.3f\n", cmap_grid->cmapdata[i].cmap[j*4+3]);
1470 void pr_ffparams(FILE *fp, int indent, const char *title,
1471 gmx_ffparams_t *ffparams,
1472 gmx_bool bShowNumbers)
1476 indent = pr_title(fp, indent, title);
1477 (void) pr_indent(fp, indent);
1478 (void) fprintf(fp, "atnr=%d\n", ffparams->atnr);
1479 (void) pr_indent(fp, indent);
1480 (void) fprintf(fp, "ntypes=%d\n", ffparams->ntypes);
1481 for (i = 0; i < ffparams->ntypes; i++)
1483 (void) pr_indent(fp, indent+INDENT);
1484 (void) fprintf(fp, "functype[%d]=%s, ",
1485 bShowNumbers ? i : -1,
1486 interaction_function[ffparams->functype[i]].name);
1487 pr_iparams(fp, ffparams->functype[i], &ffparams->iparams[i]);
1489 (void) pr_double(fp, indent, "reppow", ffparams->reppow);
1490 (void) pr_real(fp, indent, "fudgeQQ", ffparams->fudgeQQ);
1491 pr_cmap(fp, indent, "cmap", &ffparams->cmap_grid, bShowNumbers);
1494 void pr_idef(FILE *fp, int indent, const char *title, t_idef *idef, gmx_bool bShowNumbers)
1498 if (available(fp, idef, indent, title))
1500 indent = pr_title(fp, indent, title);
1501 (void) pr_indent(fp, indent);
1502 (void) fprintf(fp, "atnr=%d\n", idef->atnr);
1503 (void) pr_indent(fp, indent);
1504 (void) fprintf(fp, "ntypes=%d\n", idef->ntypes);
1505 for (i = 0; i < idef->ntypes; i++)
1507 (void) pr_indent(fp, indent+INDENT);
1508 (void) fprintf(fp, "functype[%d]=%s, ",
1509 bShowNumbers ? i : -1,
1510 interaction_function[idef->functype[i]].name);
1511 pr_iparams(fp, idef->functype[i], &idef->iparams[i]);
1513 (void) pr_real(fp, indent, "fudgeQQ", idef->fudgeQQ);
1515 for (j = 0; (j < F_NRE); j++)
1517 pr_ilist(fp, indent, interaction_function[j].longname,
1518 idef->functype, &idef->il[j], bShowNumbers);
1523 static int pr_block_title(FILE *fp, int indent, const char *title, t_block *block)
1527 if (available(fp, block, indent, title))
1529 indent = pr_title(fp, indent, title);
1530 (void) pr_indent(fp, indent);
1531 (void) fprintf(fp, "nr=%d\n", block->nr);
1536 static int pr_blocka_title(FILE *fp, int indent, const char *title, t_blocka *block)
1540 if (available(fp, block, indent, title))
1542 indent = pr_title(fp, indent, title);
1543 (void) pr_indent(fp, indent);
1544 (void) fprintf(fp, "nr=%d\n", block->nr);
1545 (void) pr_indent(fp, indent);
1546 (void) fprintf(fp, "nra=%d\n", block->nra);
1551 static void low_pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1555 if (available(fp, block, indent, title))
1557 indent = pr_blocka_title(fp, indent, title, block);
1558 for (i = 0; i <= block->nr; i++)
1560 (void) pr_indent(fp, indent+INDENT);
1561 (void) fprintf(fp, "%s->index[%d]=%d\n",
1562 title, bShowNumbers ? i : -1, block->index[i]);
1564 for (i = 0; i < block->nra; i++)
1566 (void) pr_indent(fp, indent+INDENT);
1567 (void) fprintf(fp, "%s->a[%d]=%d\n",
1568 title, bShowNumbers ? i : -1, block->a[i]);
1573 void pr_block(FILE *fp, int indent, const char *title, t_block *block, gmx_bool bShowNumbers)
1575 int i, j, ok, size, start, end;
1577 if (available(fp, block, indent, title))
1579 indent = pr_block_title(fp, indent, title, block);
1582 if ((ok = (block->index[start] == 0)) == 0)
1584 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1588 for (i = 0; i < block->nr; i++)
1590 end = block->index[i+1];
1591 size = pr_indent(fp, indent);
1594 size += fprintf(fp, "%s[%d]={}\n", title, i);
1598 size += fprintf(fp, "%s[%d]={%d..%d}\n",
1599 title, bShowNumbers ? i : -1,
1600 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1608 void pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1610 int i, j, ok, size, start, end;
1612 if (available(fp, block, indent, title))
1614 indent = pr_blocka_title(fp, indent, title, block);
1617 if ((ok = (block->index[start] == 0)) == 0)
1619 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1623 for (i = 0; i < block->nr; i++)
1625 end = block->index[i+1];
1626 size = pr_indent(fp, indent);
1629 size += fprintf(fp, "%s[%d]={", title, i);
1633 size += fprintf(fp, "%s[%d][%d..%d]={",
1634 title, bShowNumbers ? i : -1,
1635 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1637 for (j = start; j < end; j++)
1641 size += fprintf(fp, ", ");
1643 if ((size) > (USE_WIDTH))
1645 (void) fprintf(fp, "\n");
1646 size = pr_indent(fp, indent+INDENT);
1648 size += fprintf(fp, "%d", block->a[j]);
1650 (void) fprintf(fp, "}\n");
1654 if ((end != block->nra) || (!ok))
1656 (void) pr_indent(fp, indent);
1657 (void) fprintf(fp, "tables inconsistent, dumping complete tables:\n");
1658 low_pr_blocka(fp, indent, title, block, bShowNumbers);
1663 static void pr_strings(FILE *fp, int indent, const char *title, char ***nm, int n, gmx_bool bShowNumbers)
1667 if (available(fp, nm, indent, title))
1669 indent = pr_title_n(fp, indent, title, n);
1670 for (i = 0; i < n; i++)
1672 (void) pr_indent(fp, indent);
1673 (void) fprintf(fp, "%s[%d]={name=\"%s\"}\n",
1674 title, bShowNumbers ? i : -1, *(nm[i]));
1679 static void pr_strings2(FILE *fp, int indent, const char *title,
1680 char ***nm, char ***nmB, int n, gmx_bool bShowNumbers)
1684 if (available(fp, nm, indent, title))
1686 indent = pr_title_n(fp, indent, title, n);
1687 for (i = 0; i < n; i++)
1689 (void) pr_indent(fp, indent);
1690 (void) fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1691 title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
1696 static void pr_resinfo(FILE *fp, int indent, const char *title, t_resinfo *resinfo, int n, gmx_bool bShowNumbers)
1700 if (available(fp, resinfo, indent, title))
1702 indent = pr_title_n(fp, indent, title, n);
1703 for (i = 0; i < n; i++)
1705 (void) pr_indent(fp, indent);
1706 (void) fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1707 title, bShowNumbers ? i : -1,
1708 *(resinfo[i].name), resinfo[i].nr,
1709 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1714 static void pr_atom(FILE *fp, int indent, const char *title, t_atom *atom, int n)
1718 if (available(fp, atom, indent, title))
1720 indent = pr_title_n(fp, indent, title, n);
1721 for (i = 0; i < n; i++)
1723 (void) pr_indent(fp, indent);
1724 fprintf(fp, "%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1725 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1726 title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype],
1727 atom[i].m, atom[i].q, atom[i].mB, atom[i].qB,
1728 atom[i].resind, atom[i].atomnumber);
1733 static void pr_grps(FILE *fp, const char *title, t_grps grps[], char **grpname[])
1737 for (i = 0; (i < egcNR); i++)
1739 fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
1740 for (j = 0; (j < grps[i].nr); j++)
1742 fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
1748 static void pr_groups(FILE *fp, int indent,
1749 gmx_groups_t *groups,
1750 gmx_bool bShowNumbers)
1755 pr_grps(fp, "grp", groups->grps, groups->grpname);
1756 pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
1758 (void) pr_indent(fp, indent);
1759 fprintf(fp, "groups ");
1760 for (g = 0; g < egcNR; g++)
1762 printf(" %5.5s", gtypes[g]);
1766 (void) pr_indent(fp, indent);
1767 fprintf(fp, "allocated ");
1769 for (g = 0; g < egcNR; g++)
1771 printf(" %5d", groups->ngrpnr[g]);
1772 nat_max = max(nat_max, groups->ngrpnr[g]);
1778 (void) pr_indent(fp, indent);
1779 fprintf(fp, "groupnr[%5s] =", "*");
1780 for (g = 0; g < egcNR; g++)
1782 fprintf(fp, " %3d ", 0);
1788 for (i = 0; i < nat_max; i++)
1790 (void) pr_indent(fp, indent);
1791 fprintf(fp, "groupnr[%5d] =", i);
1792 for (g = 0; g < egcNR; g++)
1794 fprintf(fp, " %3d ",
1795 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1802 void pr_atoms(FILE *fp, int indent, const char *title, t_atoms *atoms,
1803 gmx_bool bShownumbers)
1805 if (available(fp, atoms, indent, title))
1807 indent = pr_title(fp, indent, title);
1808 pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
1809 pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
1810 pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
1811 pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
1816 void pr_atomtypes(FILE *fp, int indent, const char *title, t_atomtypes *atomtypes,
1817 gmx_bool bShowNumbers)
1820 if (available(fp, atomtypes, indent, title))
1822 indent = pr_title(fp, indent, title);
1823 for (i = 0; i < atomtypes->nr; i++)
1825 pr_indent(fp, indent);
1827 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1828 bShowNumbers ? i : -1, atomtypes->radius[i], atomtypes->vol[i],
1829 atomtypes->gb_radius[i],
1830 atomtypes->surftens[i], atomtypes->atomnumber[i], atomtypes->S_hct[i]);
1835 static void pr_moltype(FILE *fp, int indent, const char *title,
1836 gmx_moltype_t *molt, int n,
1837 gmx_ffparams_t *ffparams,
1838 gmx_bool bShowNumbers)
1842 indent = pr_title_n(fp, indent, title, n);
1843 (void) pr_indent(fp, indent);
1844 (void) fprintf(fp, "name=\"%s\"\n", *(molt->name));
1845 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
1846 pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
1847 pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
1848 for (j = 0; (j < F_NRE); j++)
1850 pr_ilist(fp, indent, interaction_function[j].longname,
1851 ffparams->functype, &molt->ilist[j], bShowNumbers);
1855 static void pr_molblock(FILE *fp, int indent, const char *title,
1856 gmx_molblock_t *molb, int n,
1857 gmx_moltype_t *molt)
1859 indent = pr_title_n(fp, indent, title, n);
1860 (void) pr_indent(fp, indent);
1861 (void) fprintf(fp, "%-20s = %d \"%s\"\n",
1862 "moltype", molb->type, *(molt[molb->type].name));
1863 pr_int(fp, indent, "#molecules", molb->nmol);
1864 pr_int(fp, indent, "#atoms_mol", molb->natoms_mol);
1865 pr_int(fp, indent, "#posres_xA", molb->nposres_xA);
1866 if (molb->nposres_xA > 0)
1868 pr_rvecs(fp, indent, "posres_xA", molb->posres_xA, molb->nposres_xA);
1870 pr_int(fp, indent, "#posres_xB", molb->nposres_xB);
1871 if (molb->nposres_xB > 0)
1873 pr_rvecs(fp, indent, "posres_xB", molb->posres_xB, molb->nposres_xB);
1877 void pr_mtop(FILE *fp, int indent, const char *title, gmx_mtop_t *mtop,
1878 gmx_bool bShowNumbers)
1882 if (available(fp, mtop, indent, title))
1884 indent = pr_title(fp, indent, title);
1885 (void) pr_indent(fp, indent);
1886 (void) fprintf(fp, "name=\"%s\"\n", *(mtop->name));
1887 pr_int(fp, indent, "#atoms", mtop->natoms);
1888 pr_int(fp, indent, "#molblock", mtop->nmolblock);
1889 for (mb = 0; mb < mtop->nmolblock; mb++)
1891 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
1893 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
1894 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
1895 for (mt = 0; mt < mtop->nmoltype; mt++)
1897 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
1898 &mtop->ffparams, bShowNumbers);
1900 pr_groups(fp, indent, &mtop->groups, bShowNumbers);
1904 void pr_top(FILE *fp, int indent, const char *title, t_topology *top, gmx_bool bShowNumbers)
1906 if (available(fp, top, indent, title))
1908 indent = pr_title(fp, indent, title);
1909 (void) pr_indent(fp, indent);
1910 (void) fprintf(fp, "name=\"%s\"\n", *(top->name));
1911 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
1912 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
1913 pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
1914 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
1915 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
1916 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers);
1920 void pr_header(FILE *fp, int indent, const char *title, t_tpxheader *sh)
1924 if (available(fp, sh, indent, title))
1926 indent = pr_title(fp, indent, title);
1927 pr_indent(fp, indent);
1928 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
1929 pr_indent(fp, indent);
1930 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
1931 pr_indent(fp, indent);
1932 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
1933 pr_indent(fp, indent);
1934 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
1935 pr_indent(fp, indent);
1936 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
1937 pr_indent(fp, indent);
1938 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
1940 pr_indent(fp, indent);
1941 fprintf(fp, "natoms = %d\n", sh->natoms);
1942 pr_indent(fp, indent);
1943 fprintf(fp, "lambda = %e\n", sh->lambda);
1947 void pr_commrec(FILE *fp, int indent, t_commrec *cr)
1949 pr_indent(fp, indent);
1950 fprintf(fp, "commrec:\n");
1952 pr_indent(fp, indent);
1953 fprintf(fp, "rank = %d\n", cr->nodeid);
1954 pr_indent(fp, indent);
1955 fprintf(fp, "number of ranks = %d\n", cr->nnodes);
1956 pr_indent(fp, indent);
1957 fprintf(fp, "PME-only ranks = %d\n", cr->npmenodes);
1959 pr_indent(fp,indent);
1960 fprintf(fp,"threadid = %d\n",cr->threadid);
1961 pr_indent(fp,indent);
1962 fprintf(fp,"nthreads = %d\n",cr->nthreads);