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! */
42 #include <thread_mpi.h>
56 int pr_indent(FILE *fp, int n)
60 for (i = 0; i < n; i++)
62 (void) fprintf(fp, " ");
67 int available(FILE *fp, void *p, int indent, const char *title)
73 pr_indent(fp, indent);
75 (void) fprintf(fp, "%s: not available\n", title);
80 int pr_title(FILE *fp, int indent, const char *title)
82 (void) pr_indent(fp, indent);
83 (void) fprintf(fp, "%s:\n", title);
84 return (indent+INDENT);
87 int pr_title_n(FILE *fp, int indent, const char *title, int n)
89 (void) pr_indent(fp, indent);
90 (void) fprintf(fp, "%s (%d):\n", title, n);
91 return (indent+INDENT);
94 int pr_title_nxn(FILE *fp, int indent, const char *title, int n1, int n2)
96 (void) pr_indent(fp, indent);
97 (void) fprintf(fp, "%s (%dx%d):\n", title, n1, n2);
98 return (indent+INDENT);
101 void pr_ivec(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
105 if (available(fp, vec, indent, title))
107 indent = pr_title_n(fp, indent, title, n);
108 for (i = 0; i < n; i++)
110 (void) pr_indent(fp, indent);
111 (void) fprintf(fp, "%s[%d]=%d\n", title, bShowNumbers ? i : -1, vec[i]);
116 void pr_ivec_block(FILE *fp, int indent, const char *title, int vec[], int n, gmx_bool bShowNumbers)
120 if (available(fp, vec, indent, title))
122 indent = pr_title_n(fp, indent, title, n);
127 while (j < n && vec[j] == vec[j-1]+1)
131 /* Print consecutive groups of 3 or more as blocks */
136 (void) pr_indent(fp, indent);
137 (void) fprintf(fp, "%s[%d]=%d\n",
138 title, bShowNumbers ? i : -1, vec[i]);
144 (void) pr_indent(fp, indent);
145 (void) fprintf(fp, "%s[%d,...,%d] = {%d,...,%d}\n",
147 bShowNumbers ? i : -1,
148 bShowNumbers ? j-1 : -1,
156 void pr_bvec(FILE *fp, int indent, const char *title, gmx_bool vec[], int n, gmx_bool bShowNumbers)
160 if (available(fp, vec, indent, title))
162 indent = pr_title_n(fp, indent, title, n);
163 for (i = 0; i < n; i++)
165 (void) pr_indent(fp, indent);
166 (void) fprintf(fp, "%s[%d]=%s\n", title, bShowNumbers ? i : -1,
172 void pr_ivecs(FILE *fp, int indent, const char *title, ivec vec[], int n, gmx_bool bShowNumbers)
176 if (available(fp, vec, indent, title))
178 indent = pr_title_nxn(fp, indent, title, n, DIM);
179 for (i = 0; i < n; i++)
181 (void) pr_indent(fp, indent);
182 (void) fprintf(fp, "%s[%d]={", title, bShowNumbers ? i : -1);
183 for (j = 0; j < DIM; j++)
187 (void) fprintf(fp, ", ");
189 fprintf(fp, "%d", vec[i][j]);
191 (void) fprintf(fp, "}\n");
196 void pr_rvec(FILE *fp, int indent, const char *title, real vec[], int n, gmx_bool bShowNumbers)
200 if (available(fp, vec, indent, title))
202 indent = pr_title_n(fp, indent, title, n);
203 for (i = 0; i < n; i++)
205 pr_indent(fp, indent);
206 fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
211 void pr_dvec(FILE *fp, int indent, const char *title, double vec[], int n, gmx_bool bShowNumbers)
215 if (available(fp, vec, indent, title))
217 indent = pr_title_n(fp, indent, title, n);
218 for (i = 0; i < n; i++)
220 pr_indent(fp, indent);
221 fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
228 void pr_mat(FILE *fp,int indent,char *title,matrix m)
232 if (available(fp,m,indent,title)) {
233 indent=pr_title_n(fp,indent,title,n);
235 pr_indent(fp,indent);
236 fprintf(fp,"%s[%d]=%12.5e %12.5e %12.5e\n",
237 title,bShowNumbers?i:-1,m[i][XX],m[i][YY],m[i][ZZ]);
243 void pr_rvecs_len(FILE *fp, int indent, const char *title, rvec vec[], int n)
247 if (available(fp, vec, indent, title))
249 indent = pr_title_nxn(fp, indent, title, n, DIM);
250 for (i = 0; i < n; i++)
252 (void) pr_indent(fp, indent);
253 (void) fprintf(fp, "%s[%5d]={", title, i);
254 for (j = 0; j < DIM; j++)
258 (void) fprintf(fp, ", ");
260 (void) fprintf(fp, "%12.5e", vec[i][j]);
262 (void) fprintf(fp, "} len=%12.5e\n", norm(vec[i]));
267 void pr_rvecs(FILE *fp, int indent, const char *title, rvec vec[], int n)
269 const char *fshort = "%12.5e";
270 const char *flong = "%15.8e";
274 if (getenv("LONGFORMAT") != NULL)
283 if (available(fp, vec, indent, title))
285 indent = pr_title_nxn(fp, indent, title, n, DIM);
286 for (i = 0; i < n; i++)
288 (void) pr_indent(fp, indent);
289 (void) fprintf(fp, "%s[%5d]={", title, i);
290 for (j = 0; j < DIM; j++)
294 (void) fprintf(fp, ", ");
296 (void) fprintf(fp, format, vec[i][j]);
298 (void) fprintf(fp, "}\n");
304 void pr_reals(FILE *fp, int indent, const char *title, real *vec, int n)
308 if (available(fp, vec, indent, title))
310 (void) pr_indent(fp, indent);
311 (void) fprintf(fp, "%s:\t", title);
312 for (i = 0; i < n; i++)
314 fprintf(fp, " %10g", vec[i]);
316 (void) fprintf(fp, "\n");
320 void pr_doubles(FILE *fp, int indent, const char *title, double *vec, int n)
324 if (available(fp, vec, indent, title))
326 (void) pr_indent(fp, indent);
327 (void) fprintf(fp, "%s:\t", title);
328 for (i = 0; i < n; i++)
330 fprintf(fp, " %10g", vec[i]);
332 (void) fprintf(fp, "\n");
336 static void pr_int(FILE *fp, int indent, const char *title, int i)
338 pr_indent(fp, indent);
339 fprintf(fp, "%-20s = %d\n", title, i);
342 static void pr_gmx_large_int(FILE *fp, int indent, const char *title, gmx_large_int_t i)
344 char buf[STEPSTRSIZE];
346 pr_indent(fp, indent);
347 fprintf(fp, "%-20s = %s\n", title, gmx_step_str(i, buf));
350 static void pr_real(FILE *fp, int indent, const char *title, real r)
352 pr_indent(fp, indent);
353 fprintf(fp, "%-20s = %g\n", title, r);
356 static void pr_double(FILE *fp, int indent, const char *title, double d)
358 pr_indent(fp, indent);
359 fprintf(fp, "%-20s = %g\n", title, d);
362 static void pr_str(FILE *fp, int indent, const char *title, const char *s)
364 pr_indent(fp, indent);
365 fprintf(fp, "%-20s = %s\n", title, s);
368 void pr_qm_opts(FILE *fp, int indent, const char *title, t_grpopts *opts)
372 fprintf(fp, "%s:\n", title);
374 pr_int(fp, indent, "ngQM", opts->ngQM);
377 pr_ivec(fp, indent, "QMmethod", opts->QMmethod, opts->ngQM, FALSE);
378 pr_ivec(fp, indent, "QMbasis", opts->QMbasis, opts->ngQM, FALSE);
379 pr_ivec(fp, indent, "QMcharge", opts->QMcharge, opts->ngQM, FALSE);
380 pr_ivec(fp, indent, "QMmult", opts->QMmult, opts->ngQM, FALSE);
381 pr_bvec(fp, indent, "bSH", opts->bSH, opts->ngQM, FALSE);
382 pr_ivec(fp, indent, "CASorbitals", opts->CASorbitals, opts->ngQM, FALSE);
383 pr_ivec(fp, indent, "CASelectrons", opts->CASelectrons, opts->ngQM, FALSE);
384 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
385 pr_rvec(fp, indent, "SAon", opts->SAon, opts->ngQM, FALSE);
386 pr_ivec(fp, indent, "SAsteps", opts->SAsteps, opts->ngQM, FALSE);
387 pr_bvec(fp, indent, "bOPT", opts->bOPT, opts->ngQM, FALSE);
388 pr_bvec(fp, indent, "bTS", opts->bTS, opts->ngQM, FALSE);
392 static void pr_grp_opts(FILE *out, int indent, const char *title, t_grpopts *opts,
399 fprintf(out, "%s:\n", title);
402 pr_indent(out, indent);
403 fprintf(out, "nrdf%s", bMDPformat ? " = " : ":");
404 for (i = 0; (i < opts->ngtc); i++)
406 fprintf(out, " %10g", opts->nrdf[i]);
410 pr_indent(out, indent);
411 fprintf(out, "ref-t%s", bMDPformat ? " = " : ":");
412 for (i = 0; (i < opts->ngtc); i++)
414 fprintf(out, " %10g", opts->ref_t[i]);
418 pr_indent(out, indent);
419 fprintf(out, "tau-t%s", bMDPformat ? " = " : ":");
420 for (i = 0; (i < opts->ngtc); i++)
422 fprintf(out, " %10g", opts->tau_t[i]);
426 /* Pretty-print the simulated annealing info */
427 fprintf(out, "anneal%s", bMDPformat ? " = " : ":");
428 for (i = 0; (i < opts->ngtc); i++)
430 fprintf(out, " %10s", EANNEAL(opts->annealing[i]));
434 fprintf(out, "ann-npoints%s", bMDPformat ? " = " : ":");
435 for (i = 0; (i < opts->ngtc); i++)
437 fprintf(out, " %10d", opts->anneal_npoints[i]);
441 for (i = 0; (i < opts->ngtc); i++)
443 if (opts->anneal_npoints[i] > 0)
445 fprintf(out, "ann. times [%d]:\t", i);
446 for (j = 0; (j < opts->anneal_npoints[i]); j++)
448 fprintf(out, " %10.1f", opts->anneal_time[i][j]);
451 fprintf(out, "ann. temps [%d]:\t", i);
452 for (j = 0; (j < opts->anneal_npoints[i]); j++)
454 fprintf(out, " %10.1f", opts->anneal_temp[i][j]);
460 pr_indent(out, indent);
461 fprintf(out, "acc:\t");
462 for (i = 0; (i < opts->ngacc); i++)
464 for (m = 0; (m < DIM); m++)
466 fprintf(out, " %10g", opts->acc[i][m]);
471 pr_indent(out, indent);
472 fprintf(out, "nfreeze:");
473 for (i = 0; (i < opts->ngfrz); i++)
475 for (m = 0; (m < DIM); m++)
477 fprintf(out, " %10s", opts->nFreeze[i][m] ? "Y" : "N");
483 for (i = 0; (i < opts->ngener); i++)
485 pr_indent(out, indent);
486 fprintf(out, "energygrp-flags[%3d]:", i);
487 for (m = 0; (m < opts->ngener); m++)
489 fprintf(out, " %d", opts->egp_flags[opts->ngener*i+m]);
497 static void pr_matrix(FILE *fp, int indent, const char *title, rvec *m,
502 fprintf(fp, "%-10s = %g %g %g %g %g %g\n", title,
503 m[XX][XX], m[YY][YY], m[ZZ][ZZ], m[XX][YY], m[XX][ZZ], m[YY][ZZ]);
507 pr_rvecs(fp, indent, title, m, DIM);
511 static void pr_cosine(FILE *fp, int indent, const char *title, t_cosines *cos,
518 fprintf(fp, "%s = %d\n", title, cos->n);
522 indent = pr_title(fp, indent, title);
523 (void) pr_indent(fp, indent);
524 fprintf(fp, "n = %d\n", cos->n);
527 (void) pr_indent(fp, indent+2);
529 for (j = 0; (j < cos->n); j++)
531 fprintf(fp, " %e", cos->a[j]);
534 (void) pr_indent(fp, indent+2);
535 fprintf(fp, "phi =");
536 for (j = 0; (j < cos->n); j++)
538 fprintf(fp, " %e", cos->phi[j]);
545 #define PS(t, s) pr_str(fp, indent, t, s)
546 #define PI(t, s) pr_int(fp, indent, t, s)
547 #define PSTEP(t, s) pr_gmx_large_int(fp, indent, t, s)
548 #define PR(t, s) pr_real(fp, indent, t, s)
549 #define PD(t, s) pr_double(fp, indent, t, s)
551 static void pr_pullgrp(FILE *fp, int indent, int g, t_pullgrp *pg)
553 pr_indent(fp, indent);
554 fprintf(fp, "pull-group %d:\n", g);
556 pr_ivec_block(fp, indent, "atom", pg->ind, pg->nat, TRUE);
557 pr_rvec(fp, indent, "weight", pg->weight, pg->nweight, TRUE);
558 PI("pbcatom", pg->pbcatom);
559 pr_rvec(fp, indent, "vec", pg->vec, DIM, TRUE);
560 pr_rvec(fp, indent, "init", pg->init, DIM, TRUE);
561 PR("rate", pg->rate);
566 static void pr_simtempvals(FILE *fp, int indent, t_simtemp *simtemp, int n_lambda, gmx_bool bMDPformat)
568 PR("simtemp_low", simtemp->simtemp_low);
569 PR("simtemp_high", simtemp->simtemp_high);
570 PS("simulated-tempering-scaling", ESIMTEMP(simtemp->eSimTempScale));
571 pr_rvec(fp, indent, "simulated tempering temperatures", simtemp->temperatures, n_lambda, TRUE);
574 static void pr_expandedvals(FILE *fp, int indent, t_expanded *expand, int n_lambda, gmx_bool bMDPformat)
577 PI("nstexpanded", expand->nstexpanded);
578 PS("lambda-stats", elamstats_names[expand->elamstats]);
579 PS("lambda-mc-move", elmcmove_names[expand->elmcmove]);
580 PI("lmc-repeats", expand->lmc_repeats);
581 PI("lmc-gibbsdelta", expand->gibbsdeltalam);
582 PI("lmc-nstart", expand->lmc_forced_nstart);
583 PS("symmetrized-transition-matrix", EBOOL(expand->bSymmetrizedTMatrix));
584 PI("nst-transition-matrix", expand->nstTij);
585 PI("mininum-var-min", expand->minvarmin); /*default is reasonable */
586 PI("weight-c-range", expand->c_range); /* default is just C=0 */
587 PR("wl-scale", expand->wl_scale);
588 PR("init-wl-delta", expand->init_wl_delta);
589 PR("wl-ratio", expand->wl_ratio);
590 PS("bWLoneovert", EBOOL(expand->bWLoneovert));
591 PI("lmc-seed", expand->lmc_seed);
592 PR("mc-temperature", expand->mc_temp);
593 PS("lmc-weights-equil", elmceq_names[expand->elmceq]);
594 if (expand->elmceq == elmceqNUMATLAM)
596 PI("weight-equil-number-all-lambda", expand->equil_n_at_lam);
598 if (expand->elmceq == elmceqSAMPLES)
600 PI("weight-equil-number-samples", expand->equil_samples);
602 if (expand->elmceq == elmceqSTEPS)
604 PI("weight-equil-number-steps", expand->equil_steps);
606 if (expand->elmceq == elmceqWLDELTA)
608 PR("weight-equil-wl-delta", expand->equil_wl_delta);
610 if (expand->elmceq == elmceqRATIO)
612 PR("weight-equil-count-ratio", expand->equil_ratio);
615 pr_indent(fp, indent);
616 pr_rvec(fp, indent, "init-lambda-weights", expand->init_lambda_weights, n_lambda, TRUE);
617 PS("init-weights", EBOOL(expand->bInit_weights));
620 static void pr_fepvals(FILE *fp, int indent, t_lambda *fep, gmx_bool bMDPformat)
624 PI("nstdhdl", fep->nstdhdl);
625 PI("init-lambda-state", fep->init_fep_state);
626 PR("init-lambda", fep->init_lambda);
627 PR("delta-lambda", fep->delta_lambda);
630 PI("n-lambdas", fep->n_lambda);
632 if (fep->n_lambda > 0)
634 pr_indent(fp, indent);
635 fprintf(fp, "separate-dvdl%s\n", bMDPformat ? " = " : ":");
636 for (i = 0; i < efptNR; i++)
638 fprintf(fp, "%18s = ", efpt_names[i]);
639 if (fep->separate_dvdl[i])
641 fprintf(fp, " TRUE");
645 fprintf(fp, " FALSE");
649 fprintf(fp, "all-lambdas%s\n", bMDPformat ? " = " : ":");
650 for (i = 0; i < efptNR; i++)
652 fprintf(fp, "%18s = ", efpt_names[i]);
653 for (j = 0; j < fep->n_lambda; j++)
655 fprintf(fp, " %10g", fep->all_lambda[i][j]);
660 PI("calc-lambda-neighbors", fep->lambda_neighbors);
662 PR("sc-alpha", fep->sc_alpha);
663 PS("bScCoul", EBOOL(fep->bScCoul));
664 PS("bScPrintEnergy", EBOOL(fep->bPrintEnergy));
665 PI("sc-power", fep->sc_power);
666 PR("sc-r-power", fep->sc_r_power);
667 PR("sc-sigma", fep->sc_sigma);
668 PR("sc-sigma-min", fep->sc_sigma_min);
669 PS("separate-dhdl-file", SEPDHDLFILETYPE(fep->separate_dhdl_file));
670 PS("dhdl-derivatives", DHDLDERIVATIVESTYPE(fep->dhdl_derivatives));
671 PI("dh-hist-size", fep->dh_hist_size);
672 PD("dh-hist-spacing", fep->dh_hist_spacing);
675 static void pr_pull(FILE *fp, int indent, t_pull *pull)
679 PS("pull-geometry", EPULLGEOM(pull->eGeom));
680 pr_ivec(fp, indent, "pull-dim", pull->dim, DIM, TRUE);
681 PR("pull-r1", pull->cyl_r1);
682 PR("pull-r0", pull->cyl_r0);
683 PR("pull-constr-tol", pull->constr_tol);
684 PI("pull-nstxout", pull->nstxout);
685 PI("pull-nstfout", pull->nstfout);
686 PI("pull-ngrp", pull->ngrp);
687 for (g = 0; g < pull->ngrp+1; g++)
689 pr_pullgrp(fp, indent, g, &pull->grp[g]);
693 static void pr_rotgrp(FILE *fp, int indent, int g, t_rotgrp *rotg)
695 pr_indent(fp, indent);
696 fprintf(fp, "rotation_group %d:\n", g);
698 PS("type", EROTGEOM(rotg->eType));
699 PS("massw", EBOOL(rotg->bMassW));
700 pr_ivec_block(fp, indent, "atom", rotg->ind, rotg->nat, TRUE);
701 pr_rvecs(fp, indent, "x_ref", rotg->x_ref, rotg->nat);
702 pr_rvec(fp, indent, "vec", rotg->vec, DIM, TRUE);
703 pr_rvec(fp, indent, "pivot", rotg->pivot, DIM, TRUE);
704 PR("rate", rotg->rate);
706 PR("slab_dist", rotg->slab_dist);
707 PR("min_gaussian", rotg->min_gaussian);
708 PR("epsilon", rotg->eps);
709 PS("fit_method", EROTFIT(rotg->eFittype));
710 PI("potfitangle_nstep", rotg->PotAngle_nstep);
711 PR("potfitangle_step", rotg->PotAngle_step);
714 static void pr_rot(FILE *fp, int indent, t_rot *rot)
718 PI("rot_nstrout", rot->nstrout);
719 PI("rot_nstsout", rot->nstsout);
720 PI("rot_ngrp", rot->ngrp);
721 for (g = 0; g < rot->ngrp; g++)
723 pr_rotgrp(fp, indent, g, &rot->grp[g]);
727 void pr_inputrec(FILE *fp, int indent, const char *title, t_inputrec *ir,
730 const char *infbuf = "inf";
733 if (available(fp, ir, indent, title))
737 indent = pr_title(fp, indent, title);
739 PS("integrator", EI(ir->eI));
740 PSTEP("nsteps", ir->nsteps);
741 PSTEP("init-step", ir->init_step);
742 PS("cutoff-scheme", ECUTSCHEME(ir->cutoff_scheme));
743 PS("ns_type", ENS(ir->ns_type));
744 PI("nstlist", ir->nstlist);
745 PI("ndelta", ir->ndelta);
746 PI("nstcomm", ir->nstcomm);
747 PS("comm-mode", ECOM(ir->comm_mode));
748 PI("nstlog", ir->nstlog);
749 PI("nstxout", ir->nstxout);
750 PI("nstvout", ir->nstvout);
751 PI("nstfout", ir->nstfout);
752 PI("nstcalcenergy", ir->nstcalcenergy);
753 PI("nstenergy", ir->nstenergy);
754 PI("nstxtcout", ir->nstxtcout);
755 PR("init-t", ir->init_t);
756 PR("delta-t", ir->delta_t);
758 PR("xtcprec", ir->xtcprec);
759 PR("fourierspacing", ir->fourier_spacing);
763 PI("pme-order", ir->pme_order);
764 PR("ewald-rtol", ir->ewald_rtol);
765 PR("ewald-geometry", ir->ewald_geometry);
766 PR("epsilon-surface", ir->epsilon_surface);
767 PS("optimize-fft", EBOOL(ir->bOptFFT));
768 PS("ePBC", EPBC(ir->ePBC));
769 PS("bPeriodicMols", EBOOL(ir->bPeriodicMols));
770 PS("bContinuation", EBOOL(ir->bContinuation));
771 PS("bShakeSOR", EBOOL(ir->bShakeSOR));
772 PS("etc", ETCOUPLTYPE(ir->etc));
773 PS("bPrintNHChains", EBOOL(ir->bPrintNHChains));
774 PI("nsttcouple", ir->nsttcouple);
775 PS("epc", EPCOUPLTYPE(ir->epc));
776 PS("epctype", EPCOUPLTYPETYPE(ir->epct));
777 PI("nstpcouple", ir->nstpcouple);
778 PR("tau-p", ir->tau_p);
779 pr_matrix(fp, indent, "ref-p", ir->ref_p, bMDPformat);
780 pr_matrix(fp, indent, "compress", ir->compress, bMDPformat);
781 PS("refcoord-scaling", EREFSCALINGTYPE(ir->refcoord_scaling));
784 fprintf(fp, "posres-com = %g %g %g\n", ir->posres_com[XX],
785 ir->posres_com[YY], ir->posres_com[ZZ]);
789 pr_rvec(fp, indent, "posres-com", ir->posres_com, DIM, TRUE);
793 fprintf(fp, "posres-comB = %g %g %g\n", ir->posres_comB[XX],
794 ir->posres_comB[YY], ir->posres_comB[ZZ]);
798 pr_rvec(fp, indent, "posres-comB", ir->posres_comB, DIM, TRUE);
800 PR("verlet-buffer-drift", ir->verletbuf_drift);
801 PR("rlist", ir->rlist);
802 PR("rlistlong", ir->rlistlong);
803 PR("nstcalclr", ir->nstcalclr);
804 PR("rtpi", ir->rtpi);
805 PS("coulombtype", EELTYPE(ir->coulombtype));
806 PS("coulomb-modifier", INTMODIFIER(ir->coulomb_modifier));
807 PR("rcoulomb-switch", ir->rcoulomb_switch);
808 PR("rcoulomb", ir->rcoulomb);
809 PS("vdwtype", EVDWTYPE(ir->vdwtype));
810 PS("vdw-modifier", INTMODIFIER(ir->vdw_modifier));
811 PR("rvdw-switch", ir->rvdw_switch);
812 PR("rvdw", ir->rvdw);
813 if (ir->epsilon_r != 0)
815 PR("epsilon-r", ir->epsilon_r);
819 PS("epsilon-r", infbuf);
821 if (ir->epsilon_rf != 0)
823 PR("epsilon-rf", ir->epsilon_rf);
827 PS("epsilon-rf", infbuf);
829 PR("tabext", ir->tabext);
830 PS("implicit-solvent", EIMPLICITSOL(ir->implicit_solvent));
831 PS("gb-algorithm", EGBALGORITHM(ir->gb_algorithm));
832 PR("gb-epsilon-solvent", ir->gb_epsilon_solvent);
833 PI("nstgbradii", ir->nstgbradii);
834 PR("rgbradii", ir->rgbradii);
835 PR("gb-saltconc", ir->gb_saltconc);
836 PR("gb-obc-alpha", ir->gb_obc_alpha);
837 PR("gb-obc-beta", ir->gb_obc_beta);
838 PR("gb-obc-gamma", ir->gb_obc_gamma);
839 PR("gb-dielectric-offset", ir->gb_dielectric_offset);
840 PS("sa-algorithm", ESAALGORITHM(ir->gb_algorithm));
841 PR("sa-surface-tension", ir->sa_surface_tension);
842 PS("DispCorr", EDISPCORR(ir->eDispCorr));
843 PS("bSimTemp", EBOOL(ir->bSimTemp));
846 pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda, bMDPformat);
848 PS("free-energy", EFEPTYPE(ir->efep));
849 if (ir->efep != efepNO || ir->bSimTemp)
851 pr_fepvals(fp, indent, ir->fepvals, bMDPformat);
855 pr_expandedvals(fp, indent, ir->expandedvals, ir->fepvals->n_lambda, bMDPformat);
858 PI("nwall", ir->nwall);
859 PS("wall-type", EWALLTYPE(ir->wall_type));
860 PI("wall-atomtype[0]", ir->wall_atomtype[0]);
861 PI("wall-atomtype[1]", ir->wall_atomtype[1]);
862 PR("wall-density[0]", ir->wall_density[0]);
863 PR("wall-density[1]", ir->wall_density[1]);
864 PR("wall-ewald-zfac", ir->wall_ewald_zfac);
866 PS("pull", EPULLTYPE(ir->ePull));
867 if (ir->ePull != epullNO)
869 pr_pull(fp, indent, ir->pull);
872 PS("rotation", EBOOL(ir->bRot));
875 pr_rot(fp, indent, ir->rot);
878 PS("disre", EDISRETYPE(ir->eDisre));
879 PS("disre-weighting", EDISREWEIGHTING(ir->eDisreWeighting));
880 PS("disre-mixed", EBOOL(ir->bDisreMixed));
881 PR("dr-fc", ir->dr_fc);
882 PR("dr-tau", ir->dr_tau);
883 PR("nstdisreout", ir->nstdisreout);
884 PR("orires-fc", ir->orires_fc);
885 PR("orires-tau", ir->orires_tau);
886 PR("nstorireout", ir->nstorireout);
888 PR("dihre-fc", ir->dihre_fc);
890 PR("em-stepsize", ir->em_stepsize);
891 PR("em-tol", ir->em_tol);
892 PI("niter", ir->niter);
893 PR("fc-stepsize", ir->fc_stepsize);
894 PI("nstcgsteep", ir->nstcgsteep);
895 PI("nbfgscorr", ir->nbfgscorr);
897 PS("ConstAlg", ECONSTRTYPE(ir->eConstrAlg));
898 PR("shake-tol", ir->shake_tol);
899 PI("lincs-order", ir->nProjOrder);
900 PR("lincs-warnangle", ir->LincsWarnAngle);
901 PI("lincs-iter", ir->nLincsIter);
902 PR("bd-fric", ir->bd_fric);
903 PI("ld-seed", ir->ld_seed);
904 PR("cos-accel", ir->cos_accel);
905 pr_matrix(fp, indent, "deform", ir->deform, bMDPformat);
907 PS("adress", EBOOL(ir->bAdress));
910 PS("adress_type", EADRESSTYPE(ir->adress->type));
911 PR("adress_const_wf", ir->adress->const_wf);
912 PR("adress_ex_width", ir->adress->ex_width);
913 PR("adress_hy_width", ir->adress->hy_width);
914 PS("adress_interface_correction", EADRESSICTYPE(ir->adress->icor));
915 PS("adress_site", EADRESSSITETYPE(ir->adress->site));
916 PR("adress_ex_force_cap", ir->adress->ex_forcecap);
917 PS("adress_do_hybridpairs", EBOOL(ir->adress->do_hybridpairs));
919 pr_rvec(fp, indent, "adress_reference_coords", ir->adress->refs, DIM, TRUE);
921 PI("userint1", ir->userint1);
922 PI("userint2", ir->userint2);
923 PI("userint3", ir->userint3);
924 PI("userint4", ir->userint4);
925 PR("userreal1", ir->userreal1);
926 PR("userreal2", ir->userreal2);
927 PR("userreal3", ir->userreal3);
928 PR("userreal4", ir->userreal4);
929 pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
930 pr_cosine(fp, indent, "efield-x", &(ir->ex[XX]), bMDPformat);
931 pr_cosine(fp, indent, "efield-xt", &(ir->et[XX]), bMDPformat);
932 pr_cosine(fp, indent, "efield-y", &(ir->ex[YY]), bMDPformat);
933 pr_cosine(fp, indent, "efield-yt", &(ir->et[YY]), bMDPformat);
934 pr_cosine(fp, indent, "efield-z", &(ir->ex[ZZ]), bMDPformat);
935 pr_cosine(fp, indent, "efield-zt", &(ir->et[ZZ]), bMDPformat);
936 PS("bQMMM", EBOOL(ir->bQMMM));
937 PI("QMconstraints", ir->QMconstraints);
938 PI("QMMMscheme", ir->QMMMscheme);
939 PR("scalefactor", ir->scalefactor);
940 pr_qm_opts(fp, indent, "qm-opts", &(ir->opts));
947 static void pr_harm(FILE *fp, t_iparams *iparams, const char *r, const char *kr)
949 fprintf(fp, "%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
950 r, iparams->harmonic.rA, kr, iparams->harmonic.krA,
951 r, iparams->harmonic.rB, kr, iparams->harmonic.krB);
954 void pr_iparams(FILE *fp, t_functype ftype, t_iparams *iparams)
957 real VA[4], VB[4], *rbcA, *rbcB;
963 pr_harm(fp, iparams, "th", "ct");
965 case F_CROSS_BOND_BONDS:
966 fprintf(fp, "r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
967 iparams->cross_bb.r1e, iparams->cross_bb.r2e,
968 iparams->cross_bb.krr);
970 case F_CROSS_BOND_ANGLES:
971 fprintf(fp, "r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
972 iparams->cross_ba.r1e, iparams->cross_ba.r2e,
973 iparams->cross_ba.r3e, iparams->cross_ba.krt);
975 case F_LINEAR_ANGLES:
976 fprintf(fp, "klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e\n",
977 iparams->linangle.klinA, iparams->linangle.aA,
978 iparams->linangle.klinB, iparams->linangle.aB);
981 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);
983 case F_QUARTIC_ANGLES:
984 fprintf(fp, "theta=%15.8e", iparams->qangle.theta);
985 for (i = 0; i < 5; i++)
987 fprintf(fp, ", c%c=%15.8e", '0'+i, iparams->qangle.c[i]);
992 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
993 iparams->bham.a, iparams->bham.b, iparams->bham.c);
998 pr_harm(fp, iparams, "b0", "cb");
1001 pr_harm(fp, iparams, "xi", "cx");
1004 fprintf(fp, "b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e\n",
1005 iparams->morse.b0A, iparams->morse.cbA, iparams->morse.betaA,
1006 iparams->morse.b0B, iparams->morse.cbB, iparams->morse.betaB);
1009 fprintf(fp, "b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
1010 iparams->cubic.b0, iparams->cubic.kb, iparams->cubic.kcub);
1016 fprintf(fp, "bm=%15.8e, kb=%15.8e\n", iparams->fene.bm, iparams->fene.kb);
1019 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",
1020 iparams->restraint.lowA, iparams->restraint.up1A,
1021 iparams->restraint.up2A, iparams->restraint.kA,
1022 iparams->restraint.lowB, iparams->restraint.up1B,
1023 iparams->restraint.up2B, iparams->restraint.kB);
1029 fprintf(fp, "tab=%d, kA=%15.8e, kB=%15.8e\n",
1030 iparams->tab.table, iparams->tab.kA, iparams->tab.kB);
1032 case F_POLARIZATION:
1033 fprintf(fp, "alpha=%15.8e\n", iparams->polarize.alpha);
1036 fprintf(fp, "alpha=%15.8e drcut=%15.8e khyp=%15.8e\n",
1037 iparams->anharm_polarize.alpha,
1038 iparams->anharm_polarize.drcut,
1039 iparams->anharm_polarize.khyp);
1042 fprintf(fp, "a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
1043 iparams->thole.a, iparams->thole.alpha1, iparams->thole.alpha2,
1044 iparams->thole.rfac);
1047 fprintf(fp, "al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
1048 iparams->wpol.al_x, iparams->wpol.al_y, iparams->wpol.al_z,
1049 iparams->wpol.rOH, iparams->wpol.rHH, iparams->wpol.rOD);
1052 fprintf(fp, "c6=%15.8e, c12=%15.8e\n", iparams->lj.c6, iparams->lj.c12);
1055 fprintf(fp, "c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
1056 iparams->lj14.c6A, iparams->lj14.c12A,
1057 iparams->lj14.c6B, iparams->lj14.c12B);
1060 fprintf(fp, "fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1062 iparams->ljc14.qi, iparams->ljc14.qj,
1063 iparams->ljc14.c6, iparams->ljc14.c12);
1065 case F_LJC_PAIRS_NB:
1066 fprintf(fp, "qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
1067 iparams->ljcnb.qi, iparams->ljcnb.qj,
1068 iparams->ljcnb.c6, iparams->ljcnb.c12);
1074 fprintf(fp, "phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
1075 iparams->pdihs.phiA, iparams->pdihs.cpA,
1076 iparams->pdihs.phiB, iparams->pdihs.cpB,
1077 iparams->pdihs.mult);
1080 fprintf(fp, "label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
1081 iparams->disres.label, iparams->disres.type,
1082 iparams->disres.low, iparams->disres.up1,
1083 iparams->disres.up2, iparams->disres.kfac);
1086 fprintf(fp, "ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
1087 iparams->orires.ex, iparams->orires.label, iparams->orires.power,
1088 iparams->orires.c, iparams->orires.obs, iparams->orires.kfac);
1091 fprintf(fp, "phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, kfacB=%15.8e\n",
1092 iparams->dihres.phiA, iparams->dihres.dphiA, iparams->dihres.kfacA,
1093 iparams->dihres.phiB, iparams->dihres.dphiB, iparams->dihres.kfacB);
1096 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",
1097 iparams->posres.pos0A[XX], iparams->posres.pos0A[YY],
1098 iparams->posres.pos0A[ZZ], iparams->posres.fcA[XX],
1099 iparams->posres.fcA[YY], iparams->posres.fcA[ZZ],
1100 iparams->posres.pos0B[XX], iparams->posres.pos0B[YY],
1101 iparams->posres.pos0B[ZZ], iparams->posres.fcB[XX],
1102 iparams->posres.fcB[YY], iparams->posres.fcB[ZZ]);
1105 fprintf(fp, "pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e\n",
1106 iparams->fbposres.pos0[XX], iparams->fbposres.pos0[YY],
1107 iparams->fbposres.pos0[ZZ], iparams->fbposres.geom,
1108 iparams->fbposres.r, iparams->fbposres.k);
1111 for (i = 0; i < NR_RBDIHS; i++)
1113 fprintf(fp, "%srbcA[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcA[i]);
1116 for (i = 0; i < NR_RBDIHS; i++)
1118 fprintf(fp, "%srbcB[%d]=%15.8e", i == 0 ? "" : ", ", i, iparams->rbdihs.rbcB[i]);
1123 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
1124 * OPLS potential constants back.
1126 rbcA = iparams->rbdihs.rbcA;
1127 rbcB = iparams->rbdihs.rbcB;
1129 VA[3] = -0.25*rbcA[4];
1130 VA[2] = -0.5*rbcA[3];
1131 VA[1] = 4.0*VA[3]-rbcA[2];
1132 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
1134 VB[3] = -0.25*rbcB[4];
1135 VB[2] = -0.5*rbcB[3];
1136 VB[1] = 4.0*VB[3]-rbcB[2];
1137 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
1139 for (i = 0; i < NR_FOURDIHS; i++)
1141 fprintf(fp, "%sFourA[%d]=%15.8e", i == 0 ? "" : ", ", i, VA[i]);
1144 for (i = 0; i < NR_FOURDIHS; i++)
1146 fprintf(fp, "%sFourB[%d]=%15.8e", i == 0 ? "" : ", ", i, VB[i]);
1153 fprintf(fp, "dA=%15.8e, dB=%15.8e\n", iparams->constr.dA, iparams->constr.dB);
1156 fprintf(fp, "doh=%15.8e, dhh=%15.8e\n", iparams->settle.doh,
1157 iparams->settle.dhh);
1160 fprintf(fp, "a=%15.8e\n", iparams->vsite.a);
1165 fprintf(fp, "a=%15.8e, b=%15.8e\n", iparams->vsite.a, iparams->vsite.b);
1170 fprintf(fp, "a=%15.8e, b=%15.8e, c=%15.8e\n",
1171 iparams->vsite.a, iparams->vsite.b, iparams->vsite.c);
1174 fprintf(fp, "n=%2d, a=%15.8e\n", iparams->vsiten.n, iparams->vsiten.a);
1179 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);
1182 fprintf(fp, "cmapA=%1d, cmapB=%1d\n", iparams->cmap.cmapA, iparams->cmap.cmapB);
1185 gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d",
1186 ftype, interaction_function[ftype].name, __FILE__, __LINE__);
1190 void pr_ilist(FILE *fp, int indent, const char *title,
1191 t_functype *functype, t_ilist *ilist, gmx_bool bShowNumbers)
1193 int i, j, k, type, ftype;
1196 if (available(fp, ilist, indent, title) && ilist->nr > 0)
1198 indent = pr_title(fp, indent, title);
1199 (void) pr_indent(fp, indent);
1200 fprintf(fp, "nr: %d\n", ilist->nr);
1203 (void) pr_indent(fp, indent);
1204 fprintf(fp, "iatoms:\n");
1205 iatoms = ilist->iatoms;
1206 for (i = j = 0; i < ilist->nr; )
1209 (void) pr_indent(fp, indent+INDENT);
1211 ftype = functype[type];
1212 (void) fprintf(fp, "%d type=%d (%s)",
1213 bShowNumbers ? j : -1, bShowNumbers ? type : -1,
1214 interaction_function[ftype].name);
1216 for (k = 0; k < interaction_function[ftype].nratoms; k++)
1218 (void) fprintf(fp, " %u", *(iatoms++));
1220 (void) fprintf(fp, "\n");
1221 i += 1+interaction_function[ftype].nratoms;
1223 fprintf(fp, "%5d%5d\n", i, iatoms[i]);
1231 static void pr_cmap(FILE *fp, int indent, const char *title,
1232 gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
1237 dx = 360.0 / cmap_grid->grid_spacing;
1238 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1240 if (available(fp, cmap_grid, indent, title))
1242 fprintf(fp, "%s\n", title);
1244 for (i = 0; i < cmap_grid->ngrid; i++)
1247 fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
1249 fprintf(fp, "grid[%3d]={\n", bShowNumbers ? i : -1);
1251 for (j = 0; j < nelem; j++)
1253 if ( (j%cmap_grid->grid_spacing) == 0)
1255 fprintf(fp, "%8.1f\n", idx);
1259 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4]);
1260 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+1]);
1261 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j*4+2]);
1262 fprintf(fp, "%8.3f\n", cmap_grid->cmapdata[i].cmap[j*4+3]);
1270 void pr_ffparams(FILE *fp, int indent, const char *title,
1271 gmx_ffparams_t *ffparams,
1272 gmx_bool bShowNumbers)
1276 indent = pr_title(fp, indent, title);
1277 (void) pr_indent(fp, indent);
1278 (void) fprintf(fp, "atnr=%d\n", ffparams->atnr);
1279 (void) pr_indent(fp, indent);
1280 (void) fprintf(fp, "ntypes=%d\n", ffparams->ntypes);
1281 for (i = 0; i < ffparams->ntypes; i++)
1283 (void) pr_indent(fp, indent+INDENT);
1284 (void) fprintf(fp, "functype[%d]=%s, ",
1285 bShowNumbers ? i : -1,
1286 interaction_function[ffparams->functype[i]].name);
1287 pr_iparams(fp, ffparams->functype[i], &ffparams->iparams[i]);
1289 (void) pr_double(fp, indent, "reppow", ffparams->reppow);
1290 (void) pr_real(fp, indent, "fudgeQQ", ffparams->fudgeQQ);
1291 pr_cmap(fp, indent, "cmap", &ffparams->cmap_grid, bShowNumbers);
1294 void pr_idef(FILE *fp, int indent, const char *title, t_idef *idef, gmx_bool bShowNumbers)
1298 if (available(fp, idef, indent, title))
1300 indent = pr_title(fp, indent, title);
1301 (void) pr_indent(fp, indent);
1302 (void) fprintf(fp, "atnr=%d\n", idef->atnr);
1303 (void) pr_indent(fp, indent);
1304 (void) fprintf(fp, "ntypes=%d\n", idef->ntypes);
1305 for (i = 0; i < idef->ntypes; i++)
1307 (void) pr_indent(fp, indent+INDENT);
1308 (void) fprintf(fp, "functype[%d]=%s, ",
1309 bShowNumbers ? i : -1,
1310 interaction_function[idef->functype[i]].name);
1311 pr_iparams(fp, idef->functype[i], &idef->iparams[i]);
1313 (void) pr_real(fp, indent, "fudgeQQ", idef->fudgeQQ);
1315 for (j = 0; (j < F_NRE); j++)
1317 pr_ilist(fp, indent, interaction_function[j].longname,
1318 idef->functype, &idef->il[j], bShowNumbers);
1323 static int pr_block_title(FILE *fp, int indent, const char *title, t_block *block)
1327 if (available(fp, block, indent, title))
1329 indent = pr_title(fp, indent, title);
1330 (void) pr_indent(fp, indent);
1331 (void) fprintf(fp, "nr=%d\n", block->nr);
1336 static int pr_blocka_title(FILE *fp, int indent, const char *title, t_blocka *block)
1340 if (available(fp, block, indent, title))
1342 indent = pr_title(fp, indent, title);
1343 (void) pr_indent(fp, indent);
1344 (void) fprintf(fp, "nr=%d\n", block->nr);
1345 (void) pr_indent(fp, indent);
1346 (void) fprintf(fp, "nra=%d\n", block->nra);
1351 static void low_pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1355 if (available(fp, block, indent, title))
1357 indent = pr_blocka_title(fp, indent, title, block);
1358 for (i = 0; i <= block->nr; i++)
1360 (void) pr_indent(fp, indent+INDENT);
1361 (void) fprintf(fp, "%s->index[%d]=%u\n",
1362 title, bShowNumbers ? i : -1, block->index[i]);
1364 for (i = 0; i < block->nra; i++)
1366 (void) pr_indent(fp, indent+INDENT);
1367 (void) fprintf(fp, "%s->a[%d]=%u\n",
1368 title, bShowNumbers ? i : -1, block->a[i]);
1373 void pr_block(FILE *fp, int indent, const char *title, t_block *block, gmx_bool bShowNumbers)
1375 int i, j, ok, size, start, end;
1377 if (available(fp, block, indent, title))
1379 indent = pr_block_title(fp, indent, title, block);
1382 if ((ok = (block->index[start] == 0)) == 0)
1384 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1388 for (i = 0; i < block->nr; i++)
1390 end = block->index[i+1];
1391 size = pr_indent(fp, indent);
1394 size += fprintf(fp, "%s[%d]={}\n", title, i);
1398 size += fprintf(fp, "%s[%d]={%d..%d}\n",
1399 title, bShowNumbers ? i : -1,
1400 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1408 void pr_blocka(FILE *fp, int indent, const char *title, t_blocka *block, gmx_bool bShowNumbers)
1410 int i, j, ok, size, start, end;
1412 if (available(fp, block, indent, title))
1414 indent = pr_blocka_title(fp, indent, title, block);
1417 if ((ok = (block->index[start] == 0)) == 0)
1419 (void) fprintf(fp, "block->index[%d] should be 0\n", start);
1423 for (i = 0; i < block->nr; i++)
1425 end = block->index[i+1];
1426 size = pr_indent(fp, indent);
1429 size += fprintf(fp, "%s[%d]={", title, i);
1433 size += fprintf(fp, "%s[%d][%d..%d]={",
1434 title, bShowNumbers ? i : -1,
1435 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
1437 for (j = start; j < end; j++)
1441 size += fprintf(fp, ", ");
1443 if ((size) > (USE_WIDTH))
1445 (void) fprintf(fp, "\n");
1446 size = pr_indent(fp, indent+INDENT);
1448 size += fprintf(fp, "%u", block->a[j]);
1450 (void) fprintf(fp, "}\n");
1454 if ((end != block->nra) || (!ok))
1456 (void) pr_indent(fp, indent);
1457 (void) fprintf(fp, "tables inconsistent, dumping complete tables:\n");
1458 low_pr_blocka(fp, indent, title, block, bShowNumbers);
1463 static void pr_strings(FILE *fp, int indent, const char *title, char ***nm, int n, gmx_bool bShowNumbers)
1467 if (available(fp, nm, indent, title))
1469 indent = pr_title_n(fp, indent, title, n);
1470 for (i = 0; i < n; i++)
1472 (void) pr_indent(fp, indent);
1473 (void) fprintf(fp, "%s[%d]={name=\"%s\"}\n",
1474 title, bShowNumbers ? i : -1, *(nm[i]));
1479 static void pr_strings2(FILE *fp, int indent, const char *title,
1480 char ***nm, char ***nmB, int n, gmx_bool bShowNumbers)
1484 if (available(fp, nm, indent, title))
1486 indent = pr_title_n(fp, indent, title, n);
1487 for (i = 0; i < n; i++)
1489 (void) pr_indent(fp, indent);
1490 (void) fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1491 title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
1496 static void pr_resinfo(FILE *fp, int indent, const char *title, t_resinfo *resinfo, int n, gmx_bool bShowNumbers)
1500 if (available(fp, resinfo, indent, title))
1502 indent = pr_title_n(fp, indent, title, n);
1503 for (i = 0; i < n; i++)
1505 (void) pr_indent(fp, indent);
1506 (void) fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1507 title, bShowNumbers ? i : -1,
1508 *(resinfo[i].name), resinfo[i].nr,
1509 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1514 static void pr_atom(FILE *fp, int indent, const char *title, t_atom *atom, int n)
1518 if (available(fp, atom, indent, title))
1520 indent = pr_title_n(fp, indent, title, n);
1521 for (i = 0; i < n; i++)
1523 (void) pr_indent(fp, indent);
1524 fprintf(fp, "%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1525 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1526 title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype],
1527 atom[i].m, atom[i].q, atom[i].mB, atom[i].qB,
1528 atom[i].resind, atom[i].atomnumber);
1533 static void pr_grps(FILE *fp, int indent, const char *title, t_grps grps[],
1534 char **grpname[], gmx_bool bShowNumbers)
1538 for (i = 0; (i < egcNR); i++)
1540 fprintf(fp, "%s[%-12s] nr=%d, name=[", title, gtypes[i], grps[i].nr);
1541 for (j = 0; (j < grps[i].nr); j++)
1543 fprintf(fp, " %s", *(grpname[grps[i].nm_ind[j]]));
1549 static void pr_groups(FILE *fp, int indent, const char *title,
1550 gmx_groups_t *groups,
1551 gmx_bool bShowNumbers)
1556 pr_grps(fp, indent, "grp", groups->grps, groups->grpname, bShowNumbers);
1557 pr_strings(fp, indent, "grpname", groups->grpname, groups->ngrpname, bShowNumbers);
1559 (void) pr_indent(fp, indent);
1560 fprintf(fp, "groups ");
1561 for (g = 0; g < egcNR; g++)
1563 printf(" %5.5s", gtypes[g]);
1567 (void) pr_indent(fp, indent);
1568 fprintf(fp, "allocated ");
1570 for (g = 0; g < egcNR; g++)
1572 printf(" %5d", groups->ngrpnr[g]);
1573 nat_max = max(nat_max, groups->ngrpnr[g]);
1579 (void) pr_indent(fp, indent);
1580 fprintf(fp, "groupnr[%5s] =", "*");
1581 for (g = 0; g < egcNR; g++)
1583 fprintf(fp, " %3d ", 0);
1589 for (i = 0; i < nat_max; i++)
1591 (void) pr_indent(fp, indent);
1592 fprintf(fp, "groupnr[%5d] =", i);
1593 for (g = 0; g < egcNR; g++)
1595 fprintf(fp, " %3d ",
1596 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1603 void pr_atoms(FILE *fp, int indent, const char *title, t_atoms *atoms,
1604 gmx_bool bShownumbers)
1606 if (available(fp, atoms, indent, title))
1608 indent = pr_title(fp, indent, title);
1609 pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
1610 pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
1611 pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
1612 pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
1617 void pr_atomtypes(FILE *fp, int indent, const char *title, t_atomtypes *atomtypes,
1618 gmx_bool bShowNumbers)
1621 if (available(fp, atomtypes, indent, title))
1623 indent = pr_title(fp, indent, title);
1624 for (i = 0; i < atomtypes->nr; i++)
1626 pr_indent(fp, indent);
1628 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1629 bShowNumbers ? i : -1, atomtypes->radius[i], atomtypes->vol[i],
1630 atomtypes->gb_radius[i],
1631 atomtypes->surftens[i], atomtypes->atomnumber[i], atomtypes->S_hct[i]);
1636 static void pr_moltype(FILE *fp, int indent, const char *title,
1637 gmx_moltype_t *molt, int n,
1638 gmx_ffparams_t *ffparams,
1639 gmx_bool bShowNumbers)
1643 indent = pr_title_n(fp, indent, title, n);
1644 (void) pr_indent(fp, indent);
1645 (void) fprintf(fp, "name=\"%s\"\n", *(molt->name));
1646 pr_atoms(fp, indent, "atoms", &(molt->atoms), bShowNumbers);
1647 pr_block(fp, indent, "cgs", &molt->cgs, bShowNumbers);
1648 pr_blocka(fp, indent, "excls", &molt->excls, bShowNumbers);
1649 for (j = 0; (j < F_NRE); j++)
1651 pr_ilist(fp, indent, interaction_function[j].longname,
1652 ffparams->functype, &molt->ilist[j], bShowNumbers);
1656 static void pr_molblock(FILE *fp, int indent, const char *title,
1657 gmx_molblock_t *molb, int n,
1658 gmx_moltype_t *molt,
1659 gmx_bool bShowNumbers)
1661 indent = pr_title_n(fp, indent, title, n);
1662 (void) pr_indent(fp, indent);
1663 (void) fprintf(fp, "%-20s = %d \"%s\"\n",
1664 "moltype", molb->type, *(molt[molb->type].name));
1665 pr_int(fp, indent, "#molecules", molb->nmol);
1666 pr_int(fp, indent, "#atoms_mol", molb->natoms_mol);
1667 pr_int(fp, indent, "#posres_xA", molb->nposres_xA);
1668 if (molb->nposres_xA > 0)
1670 pr_rvecs(fp, indent, "posres_xA", molb->posres_xA, molb->nposres_xA);
1672 pr_int(fp, indent, "#posres_xB", molb->nposres_xB);
1673 if (molb->nposres_xB > 0)
1675 pr_rvecs(fp, indent, "posres_xB", molb->posres_xB, molb->nposres_xB);
1679 void pr_mtop(FILE *fp, int indent, const char *title, gmx_mtop_t *mtop,
1680 gmx_bool bShowNumbers)
1684 if (available(fp, mtop, indent, title))
1686 indent = pr_title(fp, indent, title);
1687 (void) pr_indent(fp, indent);
1688 (void) fprintf(fp, "name=\"%s\"\n", *(mtop->name));
1689 pr_int(fp, indent, "#atoms", mtop->natoms);
1690 pr_int(fp, indent, "#molblock", mtop->nmolblock);
1691 for (mb = 0; mb < mtop->nmolblock; mb++)
1693 pr_molblock(fp, indent, "molblock", &mtop->molblock[mb], mb,
1694 mtop->moltype, bShowNumbers);
1696 pr_ffparams(fp, indent, "ffparams", &(mtop->ffparams), bShowNumbers);
1697 pr_atomtypes(fp, indent, "atomtypes", &(mtop->atomtypes), bShowNumbers);
1698 for (mt = 0; mt < mtop->nmoltype; mt++)
1700 pr_moltype(fp, indent, "moltype", &mtop->moltype[mt], mt,
1701 &mtop->ffparams, bShowNumbers);
1703 pr_groups(fp, indent, "groups", &mtop->groups, bShowNumbers);
1707 void pr_top(FILE *fp, int indent, const char *title, t_topology *top, gmx_bool bShowNumbers)
1709 if (available(fp, top, indent, title))
1711 indent = pr_title(fp, indent, title);
1712 (void) pr_indent(fp, indent);
1713 (void) fprintf(fp, "name=\"%s\"\n", *(top->name));
1714 pr_atoms(fp, indent, "atoms", &(top->atoms), bShowNumbers);
1715 pr_atomtypes(fp, indent, "atomtypes", &(top->atomtypes), bShowNumbers);
1716 pr_block(fp, indent, "cgs", &top->cgs, bShowNumbers);
1717 pr_block(fp, indent, "mols", &top->mols, bShowNumbers);
1718 pr_blocka(fp, indent, "excls", &top->excls, bShowNumbers);
1719 pr_idef(fp, indent, "idef", &top->idef, bShowNumbers);
1723 void pr_header(FILE *fp, int indent, const char *title, t_tpxheader *sh)
1727 if (available(fp, sh, indent, title))
1729 indent = pr_title(fp, indent, title);
1730 pr_indent(fp, indent);
1731 fprintf(fp, "bIr = %spresent\n", sh->bIr ? "" : "not ");
1732 pr_indent(fp, indent);
1733 fprintf(fp, "bBox = %spresent\n", sh->bBox ? "" : "not ");
1734 pr_indent(fp, indent);
1735 fprintf(fp, "bTop = %spresent\n", sh->bTop ? "" : "not ");
1736 pr_indent(fp, indent);
1737 fprintf(fp, "bX = %spresent\n", sh->bX ? "" : "not ");
1738 pr_indent(fp, indent);
1739 fprintf(fp, "bV = %spresent\n", sh->bV ? "" : "not ");
1740 pr_indent(fp, indent);
1741 fprintf(fp, "bF = %spresent\n", sh->bF ? "" : "not ");
1743 pr_indent(fp, indent);
1744 fprintf(fp, "natoms = %d\n", sh->natoms);
1745 pr_indent(fp, indent);
1746 fprintf(fp, "lambda = %e\n", sh->lambda);
1750 void pr_commrec(FILE *fp, int indent, t_commrec *cr)
1752 pr_indent(fp, indent);
1753 fprintf(fp, "commrec:\n");
1755 pr_indent(fp, indent);
1756 fprintf(fp, "nodeid = %d\n", cr->nodeid);
1757 pr_indent(fp, indent);
1758 fprintf(fp, "nnodes = %d\n", cr->nnodes);
1759 pr_indent(fp, indent);
1760 fprintf(fp, "npmenodes = %d\n", cr->npmenodes);
1762 pr_indent(fp,indent);
1763 fprintf(fp,"threadid = %d\n",cr->threadid);
1764 pr_indent(fp,indent);
1765 fprintf(fp,"nthreads = %d\n",cr->nthreads);