1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
40 /* This file is completely threadsafe - please keep it that way! */
52 int pr_indent(FILE *fp, int n)
56 for (i = 0; i < n; i++)
58 (void) fprintf(fp, " ");
63 int available(FILE *fp, void *p, int indent, const char *title)
69 pr_indent(fp, indent);
71 (void) fprintf(fp, "%s: not available\n", title);
76 int pr_title(FILE *fp, int indent, const char *title)
78 (void) pr_indent(fp, indent);
79 (void) fprintf(fp, "%s:\n", title);
80 return (indent+INDENT);
83 int pr_title_n(FILE *fp, int indent, const char *title, int n)
85 (void) pr_indent(fp, indent);
86 (void) fprintf(fp, "%s (%d):\n", title, n);
87 return (indent+INDENT);
90 int pr_title_nxn(FILE *fp, int indent, const char *title, int n1, int n2)
92 (void) pr_indent(fp, indent);
93 (void) fprintf(fp, "%s (%dx%d):\n", title, n1, n2);
94 return (indent+INDENT);
97 void pr_ivec(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
101 if (available(fp, vec, indent, title))
103 indent = pr_title_n(fp, indent, title, n);
104 for (i = 0; i < n; i++)
106 (void) pr_indent(fp, indent);
107 (void) fprintf(fp, "%s[%d]=%d\n", title, bShowNumbers ? i : -1, vec[i]);
112 void pr_ivec_block(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
116 if (available(fp, vec, indent, title))
118 indent = pr_title_n(fp, indent, title, n);
123 while (j < n && vec[j] == vec[j-1]+1)
127 /* Print consecutive groups of 3 or more as blocks */
132 (void) pr_indent(fp, indent);
133 (void) fprintf(fp, "%s[%d]=%d\n",
134 title, bShowNumbers ? i : -1, vec[i]);
140 (void) pr_indent(fp, indent);
141 (void) fprintf(fp, "%s[%d,...,%d] = {%d,...,%d}\n",
143 bShowNumbers ? i : -1,
144 bShowNumbers ? j-1 : -1,
152 void pr_bvec(FILE *fp, int indent, const char *title, gmx_bool vec[], int n, gmx_bool bShowNumbers)
156 if (available(fp, vec, indent, title))
158 indent = pr_title_n(fp, indent, title, n);
159 for (i = 0; i < n; i++)
161 (void) pr_indent(fp, indent);
162 (void) fprintf(fp, "%s[%d]=%s\n", title, bShowNumbers ? i : -1,
168 void pr_ivecs(FILE *fp, int indent, const char *title, ivec vec[], int n, gmx_bool bShowNumbers)
172 if (available(fp, vec, indent, title))
174 indent = pr_title_nxn(fp, indent, title, n, DIM);
175 for (i = 0; i < n; i++)
177 (void) pr_indent(fp, indent);
178 (void) fprintf(fp, "%s[%d]={", title, bShowNumbers ? i : -1);
179 for (j = 0; j < DIM; j++)
183 (void) fprintf(fp, ", ");
185 fprintf(fp, "%d", vec[i][j]);
187 (void) fprintf(fp, "}\n");
192 void pr_rvec(FILE *fp, int indent, const char *title, real vec[], int n, gmx_bool bShowNumbers)
196 if (available(fp, vec, indent, title))
198 indent = pr_title_n(fp, indent, title, n);
199 for (i = 0; i < n; i++)
201 pr_indent(fp, indent);
202 fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
207 void pr_dvec(FILE *fp, int indent, const char *title, double vec[], int n, gmx_bool bShowNumbers)
211 if (available(fp, vec, indent, title))
213 indent = pr_title_n(fp, indent, title, n);
214 for (i = 0; i < n; i++)
216 pr_indent(fp, indent);
217 fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
224 void pr_mat(FILE *fp,int indent,char *title,matrix m)
228 if (available(fp,m,indent,title)) {
229 indent=pr_title_n(fp,indent,title,n);
231 pr_indent(fp,indent);
232 fprintf(fp,"%s[%d]=%12.5e %12.5e %12.5e\n",
233 title,bShowNumbers?i:-1,m[i][XX],m[i][YY],m[i][ZZ]);
239 void pr_rvecs_len(FILE *fp, int indent, const char *title, rvec vec[], int n)
243 if (available(fp, vec, indent, title))
245 indent = pr_title_nxn(fp, indent, title, n, DIM);
246 for (i = 0; i < n; i++)
248 (void) pr_indent(fp, indent);
249 (void) fprintf(fp, "%s[%5d]={", title, i);
250 for (j = 0; j < DIM; j++)
254 (void) fprintf(fp, ", ");
256 (void) fprintf(fp, "%12.5e", vec[i][j]);
258 (void) fprintf(fp, "} len=%12.5e\n", norm(vec[i]));
263 void pr_rvecs(FILE *fp, int indent, const char *title, rvec vec[], int n)
265 const char *fshort = "%12.5e";
266 const char *flong = "%15.8e";
270 if (getenv("LONGFORMAT") != NULL)
279 if (available(fp, vec, indent, title))
281 indent = pr_title_nxn(fp, indent, title, n, DIM);
282 for (i = 0; i < n; i++)
284 (void) pr_indent(fp, indent);
285 (void) fprintf(fp, "%s[%5d]={", title, i);
286 for (j = 0; j < DIM; j++)
290 (void) fprintf(fp, ", ");
292 (void) fprintf(fp, format, vec[i][j]);
294 (void) fprintf(fp, "}\n");
300 void pr_reals(FILE *fp, int indent, const char *title, real *vec, int n)
304 if (available(fp, vec, indent, title))
306 (void) pr_indent(fp, indent);
307 (void) fprintf(fp, "%s:\t", title);
308 for (i = 0; i < n; i++)
310 fprintf(fp, " %10g", vec[i]);
312 (void) fprintf(fp, "\n");
316 void pr_doubles(FILE *fp, int indent, const char *title, double *vec, int n)
320 if (available(fp, vec, indent, title))
322 (void) pr_indent(fp, indent);
323 (void) fprintf(fp, "%s:\t", title);
324 for (i = 0; i < n; i++)
326 fprintf(fp, " %10g", vec[i]);
328 (void) fprintf(fp, "\n");
332 static void pr_int(FILE *fp, int indent, const char *title, int i)
334 pr_indent(fp, indent);
335 fprintf(fp, "%-20s = %d\n", title, i);
338 static void pr_gmx_large_int(FILE *fp, int indent, const char *title, gmx_large_int_t i)
340 char buf[STEPSTRSIZE];
342 pr_indent(fp, indent);
343 fprintf(fp, "%-20s = %s\n", title, gmx_step_str(i, buf));
346 static void pr_real(FILE *fp, int indent, const char *title, real r)
348 pr_indent(fp, indent);
349 fprintf(fp, "%-20s = %g\n", title, r);
352 static void pr_double(FILE *fp, int indent, const char *title, double d)
354 pr_indent(fp, indent);
355 fprintf(fp, "%-20s = %g\n", title, d);
358 static void pr_str(FILE *fp, int indent, const char *title, const char *s)
360 pr_indent(fp, indent);
361 fprintf(fp, "%-20s = %s\n", title, s);
364 void pr_qm_opts(FILE *fp, int indent, const char *title, t_grpopts *opts)
368 fprintf(fp, "%s:\n", title);
370 pr_int(fp, indent, "ngQM", opts->ngQM);
373 pr_ivec(fp, indent, "QMmethod", opts->QMmethod, opts->ngQM, FALSE);
374 pr_ivec(fp, indent, "QMbasis", opts->QMbasis, opts->ngQM, FALSE);
375 pr_ivec(fp, indent, "QMcharge", opts->QMcharge, opts->ngQM, FALSE);
376 pr_ivec(fp, indent, "QMmult", opts->QMmult, opts->ngQM, FALSE);
377 pr_bvec(fp, indent, "bSH", opts->bSH, opts->ngQM, FALSE);
378 pr_ivec(fp, indent, "CASorbitals", opts->CASorbitals, opts->ngQM, FALSE);
379 pr_ivec(fp, indent, "CASelectrons", opts->CASelectrons, opts->ngQM, FALSE);
380 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
381 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
382 pr_ivec(fp, indent, "SAsteps", opts->SAsteps, opts->ngQM, FALSE);
383 pr_bvec(fp, indent, "bOPT", opts->bOPT, opts->ngQM, FALSE);
384 pr_bvec(fp, indent, "bTS", opts->bTS, opts->ngQM, FALSE);
388 static void pr_grp_opts(FILE *out, int indent, const char *title, t_grpopts *opts,
395 fprintf(out, "%s:\n", title);
398 pr_indent(out, indent);
399 fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
400 for (i = 0; (i < opts->ngtc); i++)
402 fprintf(out, " %10g", opts->nrdf[i]);
406 pr_indent(out, indent);
407 fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
408 for (i = 0; (i < opts->ngtc); i++)
410 fprintf(out, " %10g", opts->ref_t[i]);
414 pr_indent(out, indent);
415 fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
416 for (i = 0; (i < opts->ngtc); i++)
418 fprintf(out, " %10g", opts->tau_t[i]);
422 /* Pretty-print the simulated annealing info */
423 fprintf(out, "anneal%s", bMDPformat ? " = " : ":");
424 for (i = 0; (i < opts->ngtc); i++)
426 fprintf(out, " %10s", EANNEAL(opts->annealing[i]));
430 fprintf(out, "ann-npoints%s", bMDPformat ? " = " : ":");
431 for (i = 0; (i < opts->ngtc); i++)
433 fprintf(out, " %10d", opts->anneal_npoints[i]);
437 for (i = 0; (i < opts->ngtc); i++)
439 if (opts->anneal_npoints[i] > 0)
441 fprintf(out, "ann. times [%d]:\t", i);
442 for (j = 0; (j < opts->anneal_npoints[i]); j++)
444 fprintf(out, " %10.1f", opts->anneal_time[i][j]);
447 fprintf(out, "ann. temps [%d]:\t", i);
448 for (j = 0; (j < opts->anneal_npoints[i]); j++)
450 fprintf(out, " %10.1f", opts->anneal_temp[i][j]);
456 pr_indent(out, indent);
457 fprintf(out, "acc:\t");
458 for (i = 0; (i < opts->ngacc); i++)
460 for (m = 0; (m < DIM); m++)
462 fprintf(out, " %10g", opts->acc[i][m]);
467 pr_indent(out, indent);
468 fprintf(out, "nfreeze:");
469 for (i = 0; (i < opts->ngfrz); i++)
471 for (m = 0; (m < DIM); m++)
473 fprintf(out, " %10s", opts->nFreeze[i][m] ? "Y" : "N");
479 for (i = 0; (i < opts->ngener); i++)
481 pr_indent(out, indent);
482 fprintf(out, "energygrp-flags[%3d]:", i);
483 for (m = 0; (m < opts->ngener); m++)
485 fprintf(out, " %d", opts->egp_flags[opts->ngener*i+m]);
493 static void pr_matrix(FILE *fp, int indent, const char *title, rvec *m,
498 fprintf(fp, "%-10s = %g %g %g %g %g %g\n", title,
499 m[XX][XX], m[YY][YY], m[ZZ][ZZ], m[XX][YY], m[XX][ZZ], m[YY][ZZ]);
503 pr_rvecs(fp, indent, title, m, DIM);
507 static void pr_cosine(FILE *fp, int indent, const char *title, t_cosines *cos,
514 fprintf(fp, "%s = %d\n", title, cos->n);
518 indent = pr_title(fp, indent, title);
519 (void) pr_indent(fp, indent);
520 fprintf(fp, "n = %d\n", cos->n);
523 (void) pr_indent(fp, indent+2);
525 for (j = 0; (j < cos->n); j++)
527 fprintf(fp, " %e", cos->a[j]);
530 (void) pr_indent(fp, indent+2);
531 fprintf(fp, "phi =");
532 for (j = 0; (j < cos->n); j++)
534 fprintf(fp, " %e", cos->phi[j]);
541 #define PS(t, s) pr_str(fp, indent, t, s)
542 #define PI(t, s) pr_int(fp, indent, t, s)
543 #define PSTEP(t, s) pr_gmx_large_int(fp, indent, t, s)
544 #define PR(t, s) pr_real(fp, indent, t, s)
545 #define PD(t, s) pr_double(fp, indent, t, s)
547 static void pr_pullgrp(FILE *fp, int indent, int g, t_pullgrp *pg)
549 pr_indent(fp, indent);
550 fprintf(fp, "pull-group %d:\n", g);
552 pr_ivec_block(fp, indent, "atom", pg->ind, pg->nat, TRUE);
553 pr_rvec(fp, indent, "weight", pg->weight, pg->nweight, TRUE);
554 PI("pbcatom", pg->pbcatom);
555 pr_rvec(fp, indent, "vec", pg->vec, DIM, TRUE);
556 pr_rvec(fp, indent, "init", pg->init, DIM, TRUE);
557 PR("rate", pg->rate);
562 static void pr_simtempvals(FILE *fp, int indent, t_simtemp *simtemp, int n_lambda)
564 PR("simtemp_low", simtemp->simtemp_low);
565 PR("simtemp_high", simtemp->simtemp_high);
566 PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
567 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
570 static void pr_expandedvals(FILE *fp, int indent, t_expanded *expand, int n_lambda)
573 PI("nstexpanded", expand->nstexpanded);
574 PS("lambda-stats", elamstats_names[expand->elamstats]);
575 PS("lambda-mc-move", elmcmove_names[expand->elmcmove]);
576 PI("lmc-repeats", expand->lmc_repeats);
577 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
578 PI("lmc-nstart", expand->lmc_forced_nstart);
579 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
580 PI("nst-transition-matrix", expand->nstTij);
581 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
582 PI("weight-c-range", expand->c_range); /* default is just C=0 */
583 PR("wl-scale", expand->wl_scale);
584 PR("init-wl-delta", expand->init_wl_delta);
585 PR("wl-ratio", expand->wl_ratio);
586 PS("bWLoneovert", EBOOL(expand->bWLoneovert));
587 PI("lmc-seed", expand->lmc_seed);
588 PR("mc-temperature", expand->mc_temp);
589 PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
590 if (expand->elmceq == elmceqNUMATLAM)
592 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
594 if (expand->elmceq == elmceqSAMPLES)
596 PI("weight-equil-number-samples", expand->equil_samples);
598 if (expand->elmceq == elmceqSTEPS)
600 PI("weight-equil-number-steps", expand->equil_steps);
602 if (expand->elmceq == elmceqWLDELTA)
604 PR("weight-equil-wl-delta", expand->equil_wl_delta);
606 if (expand->elmceq == elmceqRATIO)
608 PR("weight-equil-count-ratio", expand->equil_ratio);
611 pr_indent(fp, indent);
612 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
613 PS("init-weights", EBOOL(expand->bInit_weights));
616 static void pr_fepvals(FILE *fp, int indent, t_lambda *fep, gmx_bool bMDPformat)
620 PI("nstdhdl", fep->nstdhdl);
621 PI("init-lambda-state", fep->init_fep_state);
622 PR("init-lambda", fep->init_lambda);
623 PR("delta-lambda", fep->delta_lambda);
626 PI("n-lambdas", fep->n_lambda);
628 if (fep->n_lambda > 0)
630 pr_indent(fp, indent);
631 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
632 for (i = 0; i < efptNR; i++)
634 fprintf(fp, "%18s = ", efpt_names[i]);
635 if (fep->separate_dvdl[i])
637 fprintf(fp, " TRUE");
641 fprintf(fp, " FALSE");
645 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
646 for (i = 0; i < efptNR; i++)
648 fprintf(fp, "%18s = ", efpt_names[i]);
649 for (j = 0; j < fep->n_lambda; j++)
651 fprintf(fp, " %10g", fep->all_lambda[i][j]);
656 PI("calc-lambda-neighbors", fep->lambda_neighbors);
658 PR("sc-alpha", fep->sc_alpha);
659 PS("bScCoul", EBOOL(fep->bScCoul));
660 PS("bScPrintEnergy", EBOOL(fep->bPrintEnergy));
661 PI("sc-power", fep->sc_power);
662 PR("sc-r-power", fep->sc_r_power);
663 PR("sc-sigma", fep->sc_sigma);
664 PR("sc-sigma-min", fep->sc_sigma_min);
665 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
666 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
667 PI("dh-hist-size", fep->dh_hist_size);
668 PD("dh-hist-spacing", fep->dh_hist_spacing);
671 static void pr_pull(FILE *fp, int indent, t_pull *pull)
675 PS("pull-geometry", EPULLGEOM(pull->eGeom));
676 pr_ivec(fp, indent, "pull-dim", pull->dim, DIM, TRUE);
677 PR("pull-r1", pull->cyl_r1);
678 PR("pull-r0", pull->cyl_r0);
679 PR("pull-constr-tol", pull->constr_tol);
680 PI("pull-nstxout", pull->nstxout);
681 PI("pull-nstfout", pull->nstfout);
682 PI("pull-ngrp", pull->ngrp);
683 for (g = 0; g < pull->ngrp+1; g++)
685 pr_pullgrp(fp, indent, g, &pull->grp[g]);
689 static void pr_rotgrp(FILE *fp, int indent, int g, t_rotgrp *rotg)
691 pr_indent(fp, indent);
692 fprintf(fp, "rotation_group %d:\n", g);
694 PS("type", EROTGEOM(rotg->eType));
695 PS("massw", EBOOL(rotg->bMassW));
696 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
697 pr_rvecs(fp, indent, "x_ref", rotg->x_ref, rotg->nat);
698 pr_rvec(fp, indent, "vec", rotg->vec, DIM, TRUE);
699 pr_rvec(fp, indent, "pivot", rotg->pivot, DIM, TRUE);
700 PR("rate", rotg->rate);
702 PR("slab_dist", rotg->slab_dist);
703 PR("min_gaussian", rotg->min_gaussian);
704 PR("epsilon", rotg->eps);
705 PS("fit_method", EROTFIT(rotg->eFittype));
706 PI("potfitangle_nstep", rotg->PotAngle_nstep);
707 PR("potfitangle_step", rotg->PotAngle_step);
710 static void pr_rot(FILE *fp, int indent, t_rot *rot)
714 PI("rot_nstrout", rot->nstrout);
715 PI("rot_nstsout", rot->nstsout);
716 PI("rot_ngrp", rot->ngrp);
717 for (g = 0; g < rot->ngrp; g++)
719 pr_rotgrp(fp, indent, g, &rot->grp[g]);
723 void pr_inputrec(FILE *fp, int indent, const char *title, t_inputrec *ir,
726 const char *infbuf = "inf";
729 if (available(fp, ir, indent, title))
733 indent = pr_title(fp, indent, title);
735 PS("integrator", EI(ir->eI));
736 PSTEP("nsteps", ir->nsteps);
737 PSTEP("init-step", ir->init_step);
738 PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
739 PS("ns_type", ENS(ir->ns_type));
740 PI("nstlist", ir->nstlist);
741 PI("ndelta", ir->ndelta);
742 PI("nstcomm", ir->nstcomm);
743 PS("comm-mode", ECOM(ir->comm_mode));
744 PI("nstlog", ir->nstlog);
745 PI("nstxout", ir->nstxout);
746 PI("nstvout", ir->nstvout);
747 PI("nstfout", ir->nstfout);
748 PI("nstcalcenergy", ir->nstcalcenergy);
749 PI("nstenergy", ir->nstenergy);
750 PI("nstxtcout", ir->nstxtcout);
751 PR("init-t", ir->init_t);
752 PR("delta-t", ir->delta_t);
754 PR("xtcprec", ir->xtcprec);
755 PR("fourierspacing", ir->fourier_spacing);
759 PI("pme-order", ir->pme_order);
760 PR("ewald-rtol", ir->ewald_rtol);
761 PR("ewald-geometry", ir->ewald_geometry);
762 PR("epsilon-surface", ir->epsilon_surface);
763 PS("optimize-fft", EBOOL(ir->bOptFFT));
764 PS("ePBC", EPBC(ir->ePBC));
765 PS("bPeriodicMols", EBOOL(ir->bPeriodicMols));
766 PS("bContinuation", EBOOL(ir->bContinuation));
767 PS("bShakeSOR", EBOOL(ir->bShakeSOR));
768 PS("etc", ETCOUPLTYPE(ir->etc));
769 PS("bPrintNHChains", EBOOL(ir->bPrintNHChains));
770 PI("nsttcouple", ir->nsttcouple);
771 PS("epc", EPCOUPLTYPE(ir->epc));
772 PS("epctype", EPCOUPLTYPETYPE(ir->epct));
773 PI("nstpcouple", ir->nstpcouple);
774 PR("tau-p", ir->tau_p);
775 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
776 pr_matrix(fp, indent, "compress", ir->compress, bMDPformat);
777 PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
780 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX],
781 ir->posres_com[YY], ir->posres_com[ZZ]);
785 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
789 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX],
790 ir->posres_comB[YY], ir->posres_comB[ZZ]);
794 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
796 PR("verlet-buffer-drift", ir->verletbuf_drift);
797 PR("rlist", ir->rlist);
798 PR("rlistlong", ir->rlistlong);
799 PR("nstcalclr", ir->nstcalclr);
800 PR("rtpi", ir->rtpi);
801 PS("coulombtype", EELTYPE(ir->coulombtype));
802 PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
803 PR("rcoulomb-switch", ir->rcoulomb_switch);
804 PR("rcoulomb", ir->rcoulomb);
805 PS("vdwtype", EVDWTYPE(ir->vdwtype));
806 PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
807 PR("rvdw-switch", ir->rvdw_switch);
808 PR("rvdw", ir->rvdw);
809 if (ir->epsilon_r != 0)
811 PR("epsilon-r", ir->epsilon_r);
815 PS("epsilon-r", infbuf);
817 if (ir->epsilon_rf != 0)
819 PR("epsilon-rf", ir->epsilon_rf);
823 PS("epsilon-rf", infbuf);
825 PR("tabext", ir->tabext);
826 PS("implicit-solvent", EIMPLICITSOL(ir->implicit_solvent));
827 PS("gb-algorithm", EGBALGORITHM(ir->gb_algorithm));
828 PR("gb-epsilon-solvent", ir->gb_epsilon_solvent);
829 PI("nstgbradii", ir->nstgbradii);
830 PR("rgbradii", ir->rgbradii);
831 PR("gb-saltconc", ir->gb_saltconc);
832 PR("gb-obc-alpha", ir->gb_obc_alpha);
833 PR("gb-obc-beta", ir->gb_obc_beta);
834 PR("gb-obc-gamma", ir->gb_obc_gamma);
835 PR("gb-dielectric-offset", ir->gb_dielectric_offset);
836 PS("sa-algorithm", ESAALGORITHM(ir->gb_algorithm));
837 PR("sa-surface-tension", ir->sa_surface_tension);
838 PS("DispCorr", EDISPCORR(ir->eDispCorr));
839 PS("bSimTemp", EBOOL(ir->bSimTemp));
842 pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
844 PS("free-energy", EFEPTYPE(ir->efep));
845 if (ir->efep != efepNO || ir->bSimTemp)
847 pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
851 pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda);
854 PI("nwall", ir->nwall);
855 PS("wall-type", EWALLTYPE(ir->wall_type));
856 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
857 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
858 PR("wall-density[0]", ir->wall_density[0]);
859 PR("wall-density[1]", ir->wall_density[1]);
860 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
862 PS("pull", EPULLTYPE(ir->ePull));
863 if (ir->ePull != epullNO)
865 pr_pull(fp, indent, ir->pull);
868 PS("rotation", EBOOL(ir->bRot));
871 pr_rot(fp, indent, ir->rot);
874 PS("disre", EDISRETYPE(ir->eDisre));
875 PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
876 PS("disre-mixed", EBOOL(ir->bDisreMixed));
877 PR("dr-fc", ir->dr_fc);
878 PR("dr-tau", ir->dr_tau);
879 PR("nstdisreout", ir->nstdisreout);
880 PR("orires-fc", ir->orires_fc);
881 PR("orires-tau", ir->orires_tau);
882 PR("nstorireout", ir->nstorireout);
884 PR("dihre-fc", ir->dihre_fc);
886 PR("em-stepsize", ir->em_stepsize);
887 PR("em-tol", ir->em_tol);
888 PI("niter", ir->niter);
889 PR("fc-stepsize", ir->fc_stepsize);
890 PI("nstcgsteep", ir->nstcgsteep);
891 PI("nbfgscorr", ir->nbfgscorr);
893 PS("ConstAlg", ECONSTRTYPE(ir->eConstrAlg));
894 PR("shake-tol", ir->shake_tol);
895 PI("lincs-order", ir->nProjOrder);
896 PR("lincs-warnangle", ir->LincsWarnAngle);
897 PI("lincs-iter", ir->nLincsIter);
898 PR("bd-fric", ir->bd_fric);
899 PI("ld-seed", ir->ld_seed);
900 PR("cos-accel", ir->cos_accel);
901 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
903 PS("adress", EBOOL(ir->bAdress));
906 PS("adress_type", EADRESSTYPE(ir->adress->type));
907 PR("adress_const_wf", ir->adress->const_wf);
908 PR("adress_ex_width", ir->adress->ex_width);
909 PR("adress_hy_width", ir->adress->hy_width);
910 PS("adress_interface_correction", EADRESSICTYPE(ir->adress->icor));
911 PS("adress_site", EADRESSSITETYPE(ir->adress->site));
912 PR("adress_ex_force_cap", ir->adress->ex_forcecap);
913 PS("adress_do_hybridpairs", EBOOL(ir->adress->do_hybridpairs));
915 pr_rvec(fp, indent, "adress_reference_coords", ir->adress->refs, DIM, TRUE);
917 PI("userint1", ir->userint1);
918 PI("userint2", ir->userint2);
919 PI("userint3", ir->userint3);
920 PI("userint4", ir->userint4);
921 PR("userreal1", ir->userreal1);
922 PR("userreal2", ir->userreal2);
923 PR("userreal3", ir->userreal3);
924 PR("userreal4", ir->userreal4);
925 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
926 pr_cosine(fp, indent, "efield-x", &(ir->ex[XX]), bMDPformat);
927 pr_cosine(fp, indent, "efield-xt", &(ir->et[XX]), bMDPformat);
928 pr_cosine(fp, indent, "efield-y", &(ir->ex[YY]), bMDPformat);
929 pr_cosine(fp, indent, "efield-yt", &(ir->et[YY]), bMDPformat);
930 pr_cosine(fp, indent, "efield-z", &(ir->ex[ZZ]), bMDPformat);
931 pr_cosine(fp, indent, "efield-zt", &(ir->et[ZZ]), bMDPformat);
932 PS("bQMMM", EBOOL(ir->bQMMM));
933 PI("QMconstraints", ir->QMconstraints);
934 PI("QMMMscheme", ir->QMMMscheme);
935 PR("scalefactor", ir->scalefactor);
936 pr_qm_opts(fp, indent, "qm-opts", &(ir->opts));
943 static void pr_harm(FILE *fp, t_iparams *iparams, const char *r, const char *kr)
945 fprintf(fp, "%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
946 r, iparams->harmonic.rA, kr, iparams->harmonic.krA,
947 r, iparams->harmonic.rB, kr, iparams->harmonic.krB);
950 void pr_iparams(FILE *fp, t_functype ftype, t_iparams *iparams)
953 real VA[4], VB[4], *rbcA, *rbcB;
959 pr_harm(fp, iparams, "th", "ct");
961 case F_CROSS_BOND_BONDS:
962 fprintf(fp, "r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
963 iparams->cross_bb.r1e, iparams->cross_bb.r2e,
964 iparams->cross_bb.krr);
966 case F_CROSS_BOND_ANGLES:
967 fprintf(fp, "r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
968 iparams->cross_ba.r1e, iparams->cross_ba.r2e,
969 iparams->cross_ba.r3e, iparams->cross_ba.krt);
971 case F_LINEAR_ANGLES:
972 fprintf(fp, "klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
973 iparams->linangle.klinA, iparams->linangle.aA,
974 iparams->linangle.klinB, iparams->linangle.aB);
977 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);
979 case F_QUARTIC_ANGLES:
980 fprintf(fp, "theta=%15.8e", iparams->qangle.theta);
981 for (i = 0; i < 5; i++)
983 fprintf(fp, ", c%c=%15.8e", '0'+i, iparams->qangle.c[i]);
988 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
989 iparams->bham.a, iparams->bham.b, iparams->bham.c);
994 pr_harm(fp, iparams, "b0", "cb");
997 pr_harm(fp, iparams, "xi", "cx");
1000 fprintf(fp, "b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
1001 iparams->morse.b0A, iparams->morse.cbA, iparams->morse.betaA,
1002 iparams->morse.b0B, iparams->morse.cbB, iparams->morse.betaB);
1005 fprintf(fp, "b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
1006 iparams->cubic.b0, iparams->cubic.kb, iparams->cubic.kcub);
1012 fprintf(fp, "bm=%15.8e, kb=%15.8e\n", iparams->fene.bm, iparams->fene.kb);
1015 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",
1016 iparams->restraint.lowA, iparams->restraint.up1A,
1017 iparams->restraint.up2A, iparams->restraint.kA,
1018 iparams->restraint.lowB, iparams->restraint.up1B,
1019 iparams->restraint.up2B, iparams->restraint.kB);
1025 fprintf(fp, "tab=%d, kA=%15.8e, kB=%15.8e\n",
1026 iparams->tab.table, iparams->tab.kA, iparams->tab.kB);
1028 case F_POLARIZATION:
1029 fprintf(fp, "alpha=%15.8e\n", iparams->polarize.alpha);
1032 fprintf(fp, "alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
1033 iparams->anharm_polarize.alpha,
1034 iparams->anharm_polarize.drcut,
1035 iparams->anharm_polarize.khyp);
1038 fprintf(fp, "a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
1039 iparams->thole.a, iparams->thole.alpha1, iparams->thole.alpha2,
1040 iparams->thole.rfac);
1043 fprintf(fp, "al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
1044 iparams->wpol.al_x, iparams->wpol.al_y, iparams->wpol.al_z,
1045 iparams->wpol.rOH, iparams->wpol.rHH, iparams->wpol.rOD);
1048 fprintf(fp, "c6=%15.8e, c12=%15.8e\n", iparams->lj.c6, iparams->lj.c12);
1051 fprintf(fp, "c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
1052 iparams->lj14.c6A, iparams->lj14.c12A,
1053 iparams->lj14.c6B, iparams->lj14.c12B);
1056 fprintf(fp, "fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1058 iparams->ljc14.qi, iparams->ljc14.qj,
1059 iparams->ljc14.c6, iparams->ljc14.c12);
1061 case F_LJC_PAIRS_NB:
1062 fprintf(fp, "qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1063 iparams->ljcnb.qi, iparams->ljcnb.qj,
1064 iparams->ljcnb.c6, iparams->ljcnb.c12);
1070 fprintf(fp, "phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
1071 iparams->pdihs.phiA, iparams->pdihs.cpA,
1072 iparams->pdihs.phiB, iparams->pdihs.cpB,
1073 iparams->pdihs.mult);
1076 fprintf(fp, "label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
1077 iparams->disres.label, iparams->disres.type,
1078 iparams->disres.low, iparams->disres.up1,
1079 iparams->disres.up2, iparams->disres.kfac);
1082 fprintf(fp, "ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
1083 iparams->orires.ex, iparams->orires.label, iparams->orires.power,
1084 iparams->orires.c, iparams->orires.obs, iparams->orires.kfac);
1087 fprintf(fp, "phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
1088 iparams->dihres.phiA, iparams->dihres.dphiA, iparams->dihres.kfacA,
1089 iparams->dihres.phiB, iparams->dihres.dphiB, iparams->dihres.kfacB);
1092 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",
1093 iparams->posres.pos0A[XX], iparams->posres.pos0A[YY],
1094 iparams->posres.pos0A[ZZ], iparams->posres.fcA[XX],
1095 iparams->posres.fcA[YY], iparams->posres.fcA[ZZ],
1096 iparams->posres.pos0B[XX], iparams->posres.pos0B[YY],
1097 iparams->posres.pos0B[ZZ], iparams->posres.fcB[XX],
1098 iparams->posres.fcB[YY], iparams->posres.fcB[ZZ]);
1101 fprintf(fp, "pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e\n",
1102 iparams->fbposres.pos0[XX], iparams->fbposres.pos0[YY],
1103 iparams->fbposres.pos0[ZZ], iparams->fbposres.geom,
1104 iparams->fbposres.r, iparams->fbposres.k);
1107 for (i = 0; i < NR_RBDIHS; i++)
1109 fprintf(fp, "%srbcA[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcA[i]);
1112 for (i = 0; i < NR_RBDIHS; i++)
1114 fprintf(fp, "%srbcB[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcB[i]);
1119 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
1120 * OPLS potential constants back.
1122 rbcA = iparams->rbdihs.rbcA;
1123 rbcB = iparams->rbdihs.rbcB;
1125 VA[3] = -0.25*rbcA[4];
1126 VA[2] = -0.5*rbcA[3];
1127 VA[1] = 4.0*VA[3]-rbcA[2];
1128 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1130 VB[3] = -0.25*rbcB[4];
1131 VB[2] = -0.5*rbcB[3];
1132 VB[1] = 4.0*VB[3]-rbcB[2];
1133 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1135 for (i = 0; i < NR_FOURDIHS; i++)
1137 fprintf(fp, "%sFourA[%d]=%15.8e", i == 0 ? "" : ", ", i, VA[i]);
1140 for (i = 0; i < NR_FOURDIHS; i++)
1142 fprintf(fp, "%sFourB[%d]=%15.8e", i == 0 ? "" : ", ", i, VB[i]);
1149 fprintf(fp, "dA=%15.8e, dB=%15.8e\n", iparams->constr.dA, iparams->constr.dB);
1152 fprintf(fp, "doh=%15.8e, dhh=%15.8e\n", iparams->settle.doh,
1153 iparams->settle.dhh);
1156 fprintf(fp, "a=%15.8e\n", iparams->vsite.a);
1161 fprintf(fp, "a=%15.8e, b=%15.8e\n", iparams->vsite.a, iparams->vsite.b);
1166 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1167 iparams->vsite.a, iparams->vsite.b, iparams->vsite.c);
1170 fprintf(fp, "n=%2d, a=%15.8e\n", iparams->vsiten.n, iparams->vsiten.a);
1175 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);
1178 fprintf(fp, "cmapA=%1d, cmapB=%1d\n", iparams->cmap.cmapA, iparams->cmap.cmapB);
1181 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1182 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1186 void pr_ilist(FILE *fp, int indent, const char *title,
1187 t_functype *functype, t_ilist *ilist, gmx_bool bShowNumbers)
1189 int i, j, k, type, ftype;
1192 if (available(fp, ilist, indent, title) && ilist->nr > 0)
1194 indent = pr_title(fp, indent, title);
1195 (void) pr_indent(fp, indent);
1196 fprintf(fp, "nr: %d\n", ilist->nr);
1199 (void) pr_indent(fp, indent);
1200 fprintf(fp, "iatoms:\n");
1201 iatoms = ilist->iatoms;
1202 for (i = j = 0; i < ilist->nr; )
1205 (void) pr_indent(fp, indent+INDENT);
1207 ftype = functype[type];
1208 (void) fprintf(fp, "%d type=%d (%s)",
1209 bShowNumbers ? j : -1, bShowNumbers ? type : -1,
1210 interaction_function[ftype].name);
1212 for (k = 0; k < interaction_function[ftype].nratoms; k++)
1214 (void) fprintf(fp, " %u", *(iatoms++));
1216 (void) fprintf(fp, "\n");
1217 i += 1+interaction_function[ftype].nratoms;
1219 fprintf(fp, "%5d%5d\n", i, iatoms[i]);
1227 static void pr_cmap(FILE *fp, int indent, const char *title,
1228 gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1233 dx = 360.0 / cmap_grid->grid_spacing;
1234 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1236 if (available(fp, cmap_grid, indent, title))
1238 fprintf(fp, "%s\n", title);
1240 for (i = 0; i < cmap_grid->ngrid; i++)
1243 fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
1245 fprintf(fp, "grid[%3d]={\n", bShowNumbers ? i : -1);
1247 for (j = 0; j < nelem; j++)
1249 if ( (j%cmap_grid->grid_spacing) == 0)
1251 fprintf(fp, "%8.1f\n", idx);
1255 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4]);
1256 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+1]);
1257 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+2]);
1258 fprintf(fp, "%8.3f\n", cmap_grid->cmapdata[i].cmap[j*4+3]);
1266 void pr_ffparams(FILE *fp, int indent, const char *title,
1267 gmx_ffparams_t *ffparams,
1268 gmx_bool bShowNumbers)
1272 indent = pr_title(fp, indent, title);
1273 (void) pr_indent(fp, indent);
1274 (void) fprintf(fp, "atnr=%d\n", ffparams->atnr);
1275 (void) pr_indent(fp, indent);
1276 (void) fprintf(fp, "ntypes=%d\n", ffparams->ntypes);
1277 for (i = 0; i < ffparams->ntypes; i++)
1279 (void) pr_indent(fp, indent+INDENT);
1280 (void) fprintf(fp, "functype[%d]=%s, ",
1281 bShowNumbers ? i : -1,
1282 interaction_function[ffparams->functype[i]].name);
1283 pr_iparams(fp, ffparams->functype[i], &ffparams->iparams[i]);
1285 (void) pr_double(fp, indent, "reppow", ffparams->reppow);
1286 (void) pr_real(fp, indent, "fudgeQQ", ffparams->fudgeQQ);
1287 pr_cmap(fp, indent, "cmap", &ffparams->cmap_grid, bShowNumbers);
1290 void pr_idef(FILE *fp, int indent, const char *title, t_idef *idef, gmx_bool bShowNumbers)
1294 if (available(fp, idef, indent, title))
1296 indent = pr_title(fp, indent, title);
1297 (void) pr_indent(fp, indent);
1298 (void) fprintf(fp, "atnr=%d\n", idef->atnr);
1299 (void) pr_indent(fp, indent);
1300 (void) fprintf(fp, "ntypes=%d\n", idef->ntypes);
1301 for (i = 0; i < idef->ntypes; i++)
1303 (void) pr_indent(fp, indent+INDENT);
1304 (void) fprintf(fp, "functype[%d]=%s, ",
1305 bShowNumbers ? i : -1,
1306 interaction_function[idef->functype[i]].name);
1307 pr_iparams(fp, idef->functype[i], &idef->iparams[i]);
1309 (void) pr_real(fp, indent, "fudgeQQ", idef->fudgeQQ);
1311 for (j = 0; (j < F_NRE); j++)
1313 pr_ilist(fp, indent, interaction_function[j].longname,
1314 idef->functype, &idef->il[j], bShowNumbers);
1319 static int pr_block_title(FILE *fp, int indent, const char *title, t_block *block)
1323 if (available(fp, block, indent, title))
1325 indent = pr_title(fp, indent, title);
1326 (void) pr_indent(fp, indent);
1327 (void) fprintf(fp, "nr=%d\n", block->nr);
1332 static int pr_blocka_title(FILE *fp, int indent, const char *title, t_blocka *block)
1336 if (available(fp, block, indent, title))
1338 indent = pr_title(fp, indent, title);
1339 (void) pr_indent(fp, indent);
1340 (void) fprintf(fp, "nr=%d\n", block->nr);
1341 (void) pr_indent(fp, indent);
1342 (void) fprintf(fp, "nra=%d\n", block->nra);
1347 static void low_pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1351 if (available(fp, block, indent, title))
1353 indent = pr_blocka_title(fp, indent, title, block);
1354 for (i = 0; i <= block->nr; i++)
1356 (void) pr_indent(fp, indent+INDENT);
1357 (void) fprintf(fp, "%s->index[%d]=%u\n",
1358 title, bShowNumbers ? i : -1, block->index[i]);
1360 for (i = 0; i < block->nra; i++)
1362 (void) pr_indent(fp, indent+INDENT);
1363 (void) fprintf(fp, "%s->a[%d]=%u\n",
1364 title, bShowNumbers ? i : -1, block->a[i]);
1369 void pr_block(FILE *fp, int indent, const char *title, t_block *block, gmx_bool bShowNumbers)
1371 int i, j, ok, size, start, end;
1373 if (available(fp, block, indent, title))
1375 indent = pr_block_title(fp, indent, title, block);
1378 if ((ok = (block->index[start] == 0)) == 0)
1380 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1384 for (i = 0; i < block->nr; i++)
1386 end = block->index[i+1];
1387 size = pr_indent(fp, indent);
1390 size += fprintf(fp, "%s[%d]={}\n", title, i);
1394 size += fprintf(fp, "%s[%d]={%d..%d}\n",
1395 title, bShowNumbers ? i : -1,
1396 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1404 void pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1406 int i, j, ok, size, start, end;
1408 if (available(fp, block, indent, title))
1410 indent = pr_blocka_title(fp, indent, title, block);
1413 if ((ok = (block->index[start] == 0)) == 0)
1415 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1419 for (i = 0; i < block->nr; i++)
1421 end = block->index[i+1];
1422 size = pr_indent(fp, indent);
1425 size += fprintf(fp, "%s[%d]={", title, i);
1429 size += fprintf(fp, "%s[%d][%d..%d]={",
1430 title, bShowNumbers ? i : -1,
1431 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1433 for (j = start; j < end; j++)
1437 size += fprintf(fp, ", ");
1439 if ((size) > (USE_WIDTH))
1441 (void) fprintf(fp, "\n");
1442 size = pr_indent(fp, indent+INDENT);
1444 size += fprintf(fp, "%u", block->a[j]);
1446 (void) fprintf(fp, "}\n");
1450 if ((end != block->nra) || (!ok))
1452 (void) pr_indent(fp, indent);
1453 (void) fprintf(fp, "tables inconsistent, dumping complete tables:\n");
1454 low_pr_blocka(fp, indent, title, block, bShowNumbers);
1459 static void pr_strings(FILE *fp, int indent, const char *title, char ***nm, int n, gmx_bool bShowNumbers)
1463 if (available(fp, nm, indent, title))
1465 indent = pr_title_n(fp, indent, title, n);
1466 for (i = 0; i < n; i++)
1468 (void) pr_indent(fp, indent);
1469 (void) fprintf(fp, "%s[%d]={name=\"%s\"}\n",
1470 title, bShowNumbers ? i : -1, *(nm[i]));
1475 static void pr_strings2(FILE *fp, int indent, const char *title,
1476 char ***nm, char ***nmB, int n, gmx_bool bShowNumbers)
1480 if (available(fp, nm, indent, title))
1482 indent = pr_title_n(fp, indent, title, n);
1483 for (i = 0; i < n; i++)
1485 (void) pr_indent(fp, indent);
1486 (void) fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1487 title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
1492 static void pr_resinfo(FILE *fp, int indent, const char *title, t_resinfo *resinfo, int n, gmx_bool bShowNumbers)
1496 if (available(fp, resinfo, indent, title))
1498 indent = pr_title_n(fp, indent, title, n);
1499 for (i = 0; i < n; i++)
1501 (void) pr_indent(fp, indent);
1502 (void) fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1503 title, bShowNumbers ? i : -1,
1504 *(resinfo[i].name), resinfo[i].nr,
1505 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1510 static void pr_atom(FILE *fp, int indent, const char *title, t_atom *atom, int n)
1514 if (available(fp, atom, indent, title))
1516 indent = pr_title_n(fp, indent, title, n);
1517 for (i = 0; i < n; i++)
1519 (void) pr_indent(fp, indent);
1520 fprintf(fp, "%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1521 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1522 title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype],
1523 atom[i].m, atom[i].q, atom[i].mB, atom[i].qB,
1524 atom[i].resind, atom[i].atomnumber);
1529 static void pr_grps(FILE *fp, const char *title, t_grps grps[], char **grpname[])
1533 for (i = 0; (i < egcNR); i++)
1535 fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
1536 for (j = 0; (j < grps[i].nr); j++)
1538 fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
1544 static void pr_groups(FILE *fp, int indent,
1545 gmx_groups_t *groups,
1546 gmx_bool bShowNumbers)
1551 pr_grps(fp, "grp", groups->grps, groups->grpname);
1552 pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
1554 (void) pr_indent(fp, indent);
1555 fprintf(fp, "groups ");
1556 for (g = 0; g < egcNR; g++)
1558 printf(" %5.5s", gtypes[g]);
1562 (void) pr_indent(fp, indent);
1563 fprintf(fp, "allocated ");
1565 for (g = 0; g < egcNR; g++)
1567 printf(" %5d", groups->ngrpnr[g]);
1568 nat_max = max(nat_max, groups->ngrpnr[g]);
1574 (void) pr_indent(fp, indent);
1575 fprintf(fp, "groupnr[%5s] =", "*");
1576 for (g = 0; g < egcNR; g++)
1578 fprintf(fp, " %3d ", 0);
1584 for (i = 0; i < nat_max; i++)
1586 (void) pr_indent(fp, indent);
1587 fprintf(fp, "groupnr[%5d] =", i);
1588 for (g = 0; g < egcNR; g++)
1590 fprintf(fp, " %3d ",
1591 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1598 void pr_atoms(FILE *fp, int indent, const char *title, t_atoms *atoms,
1599 gmx_bool bShownumbers)
1601 if (available(fp, atoms, indent, title))
1603 indent = pr_title(fp, indent, title);
1604 pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
1605 pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
1606 pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
1607 pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
1612 void pr_atomtypes(FILE *fp, int indent, const char *title, t_atomtypes *atomtypes,
1613 gmx_bool bShowNumbers)
1616 if (available(fp, atomtypes, indent, title))
1618 indent = pr_title(fp, indent, title);
1619 for (i = 0; i < atomtypes->nr; i++)
1621 pr_indent(fp, indent);
1623 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1624 bShowNumbers ? i : -1, atomtypes->radius[i], atomtypes->vol[i],
1625 atomtypes->gb_radius[i],
1626 atomtypes->surftens[i], atomtypes->atomnumber[i], atomtypes->S_hct[i]);
1631 static void pr_moltype(FILE *fp, int indent, const char *title,
1632 gmx_moltype_t *molt, int n,
1633 gmx_ffparams_t *ffparams,
1634 gmx_bool bShowNumbers)
1638 indent = pr_title_n(fp, indent, title, n);
1639 (void) pr_indent(fp, indent);
1640 (void) fprintf(fp, "name=\"%s\"\n", *(molt->name));
1641 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
1642 pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
1643 pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
1644 for (j = 0; (j < F_NRE); j++)
1646 pr_ilist(fp, indent, interaction_function[j].longname,
1647 ffparams->functype, &molt->ilist[j], bShowNumbers);
1651 static void pr_molblock(FILE *fp, int indent, const char *title,
1652 gmx_molblock_t *molb, int n,
1653 gmx_moltype_t *molt)
1655 indent = pr_title_n(fp, indent, title, n);
1656 (void) pr_indent(fp, indent);
1657 (void) fprintf(fp, "%-20s = %d \"%s\"\n",
1658 "moltype", molb->type, *(molt[molb->type].name));
1659 pr_int(fp, indent, "#molecules", molb->nmol);
1660 pr_int(fp, indent, "#atoms_mol", molb->natoms_mol);
1661 pr_int(fp, indent, "#posres_xA", molb->nposres_xA);
1662 if (molb->nposres_xA > 0)
1664 pr_rvecs(fp, indent, "posres_xA", molb->posres_xA, molb->nposres_xA);
1666 pr_int(fp, indent, "#posres_xB", molb->nposres_xB);
1667 if (molb->nposres_xB > 0)
1669 pr_rvecs(fp, indent, "posres_xB", molb->posres_xB, molb->nposres_xB);
1673 void pr_mtop(FILE *fp, int indent, const char *title, gmx_mtop_t *mtop,
1674 gmx_bool bShowNumbers)
1678 if (available(fp, mtop, indent, title))
1680 indent = pr_title(fp, indent, title);
1681 (void) pr_indent(fp, indent);
1682 (void) fprintf(fp, "name=\"%s\"\n", *(mtop->name));
1683 pr_int(fp, indent, "#atoms", mtop->natoms);
1684 pr_int(fp, indent, "#molblock", mtop->nmolblock);
1685 for (mb = 0; mb < mtop->nmolblock; mb++)
1687 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb, mtop->moltype);
1689 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
1690 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
1691 for (mt = 0; mt < mtop->nmoltype; mt++)
1693 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
1694 &mtop->ffparams, bShowNumbers);
1696 pr_groups(fp, indent, &mtop->groups, bShowNumbers);
1700 void pr_top(FILE *fp, int indent, const char *title, t_topology *top, gmx_bool bShowNumbers)
1702 if (available(fp, top, indent, title))
1704 indent = pr_title(fp, indent, title);
1705 (void) pr_indent(fp, indent);
1706 (void) fprintf(fp, "name=\"%s\"\n", *(top->name));
1707 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
1708 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
1709 pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
1710 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
1711 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
1712 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers);
1716 void pr_header(FILE *fp, int indent, const char *title, t_tpxheader *sh)
1720 if (available(fp, sh, indent, title))
1722 indent = pr_title(fp, indent, title);
1723 pr_indent(fp, indent);
1724 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
1725 pr_indent(fp, indent);
1726 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
1727 pr_indent(fp, indent);
1728 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
1729 pr_indent(fp, indent);
1730 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
1731 pr_indent(fp, indent);
1732 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
1733 pr_indent(fp, indent);
1734 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
1736 pr_indent(fp, indent);
1737 fprintf(fp, "natoms = %d\n", sh->natoms);
1738 pr_indent(fp, indent);
1739 fprintf(fp, "lambda = %e\n", sh->lambda);
1743 void pr_commrec(FILE *fp, int indent, t_commrec *cr)
1745 pr_indent(fp, indent);
1746 fprintf(fp, "commrec:\n");
1748 pr_indent(fp, indent);
1749 fprintf(fp, "nodeid = %d\n", cr->nodeid);
1750 pr_indent(fp, indent);
1751 fprintf(fp, "nnodes = %d\n", cr->nnodes);
1752 pr_indent(fp, indent);
1753 fprintf(fp, "npmenodes = %d\n", cr->npmenodes);
1755 pr_indent(fp,indent);
1756 fprintf(fp,"threadid = %d\n",cr->threadid);
1757 pr_indent(fp,indent);
1758 fprintf(fp,"nthreads = %d\n",cr->nthreads);