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>
55 int pr_indent(FILE *fp,int n)
59 for (i=0; i<n; i++) (void) fprintf(fp," ");
63 int available(FILE *fp,void *p,int indent,const char *title)
68 (void) fprintf(fp,"%s: not available\n",title);
73 int pr_title(FILE *fp,int indent,const char *title)
75 (void) pr_indent(fp,indent);
76 (void) fprintf(fp,"%s:\n",title);
77 return (indent+INDENT);
80 int pr_title_n(FILE *fp,int indent,const char *title,int n)
82 (void) pr_indent(fp,indent);
83 (void) fprintf(fp,"%s (%d):\n",title,n);
84 return (indent+INDENT);
87 int pr_title_nxn(FILE *fp,int indent,const char *title,int n1,int n2)
89 (void) pr_indent(fp,indent);
90 (void) fprintf(fp,"%s (%dx%d):\n",title,n1,n2);
91 return (indent+INDENT);
94 void pr_ivec(FILE *fp,int indent,const char *title,int vec[],int n, gmx_bool bShowNumbers)
98 if (available(fp,vec,indent,title))
100 indent=pr_title_n(fp,indent,title,n);
103 (void) pr_indent(fp,indent);
104 (void) fprintf(fp,"%s[%d]=%d\n",title,bShowNumbers?i:-1,vec[i]);
109 void pr_ivec_block(FILE *fp,int indent,const char *title,int vec[],int n, gmx_bool bShowNumbers)
113 if (available(fp,vec,indent,title))
115 indent=pr_title_n(fp,indent,title,n);
120 while (j < n && vec[j] == vec[j-1]+1)
124 /* Print consecutive groups of 3 or more as blocks */
129 (void) pr_indent(fp,indent);
130 (void) fprintf(fp,"%s[%d]=%d\n",
131 title,bShowNumbers?i:-1,vec[i]);
137 (void) pr_indent(fp,indent);
138 (void) fprintf(fp,"%s[%d,...,%d] = {%d,...,%d}\n",
149 void pr_bvec(FILE *fp,int indent,const char *title,gmx_bool vec[],int n, gmx_bool bShowNumbers)
153 if (available(fp,vec,indent,title))
155 indent=pr_title_n(fp,indent,title,n);
158 (void) pr_indent(fp,indent);
159 (void) fprintf(fp,"%s[%d]=%s\n",title,bShowNumbers?i:-1,
165 void pr_ivecs(FILE *fp,int indent,const char *title,ivec vec[],int n, gmx_bool bShowNumbers)
169 if (available(fp,vec,indent,title))
171 indent=pr_title_nxn(fp,indent,title,n,DIM);
174 (void) pr_indent(fp,indent);
175 (void) fprintf(fp,"%s[%d]={",title,bShowNumbers?i:-1);
176 for (j=0; j<DIM; j++)
178 if (j!=0) (void) fprintf(fp,", ");
179 fprintf(fp,"%d",vec[i][j]);
181 (void) fprintf(fp,"}\n");
186 void pr_rvec(FILE *fp,int indent,const char *title,real vec[],int n, gmx_bool bShowNumbers)
190 if (available(fp,vec,indent,title))
192 indent=pr_title_n(fp,indent,title,n);
195 pr_indent(fp,indent);
196 fprintf(fp,"%s[%d]=%12.5e\n",title,bShowNumbers?i:-1,vec[i]);
201 void pr_dvec(FILE *fp,int indent,const char *title,double vec[],int n, gmx_bool bShowNumbers)
205 if (available(fp,vec,indent,title))
207 indent=pr_title_n(fp,indent,title,n);
210 pr_indent(fp,indent);
211 fprintf(fp,"%s[%d]=%12.5e\n",title,bShowNumbers?i:-1,vec[i]);
218 void pr_mat(FILE *fp,int indent,char *title,matrix m)
222 if (available(fp,m,indent,title)) {
223 indent=pr_title_n(fp,indent,title,n);
225 pr_indent(fp,indent);
226 fprintf(fp,"%s[%d]=%12.5e %12.5e %12.5e\n",
227 title,bShowNumbers?i:-1,m[i][XX],m[i][YY],m[i][ZZ]);
233 void pr_rvecs_len(FILE *fp,int indent,const char *title,rvec vec[],int n)
237 if (available(fp,vec,indent,title)) {
238 indent=pr_title_nxn(fp,indent,title,n,DIM);
239 for (i=0; i<n; i++) {
240 (void) pr_indent(fp,indent);
241 (void) fprintf(fp,"%s[%5d]={",title,i);
242 for (j=0; j<DIM; j++) {
244 (void) fprintf(fp,", ");
245 (void) fprintf(fp,"%12.5e",vec[i][j]);
247 (void) fprintf(fp,"} len=%12.5e\n",norm(vec[i]));
252 void pr_rvecs(FILE *fp,int indent,const char *title,rvec vec[],int n)
254 const char *fshort = "%12.5e";
255 const char *flong = "%15.8e";
259 if (getenv("LONGFORMAT") != NULL)
264 if (available(fp,vec,indent,title)) {
265 indent=pr_title_nxn(fp,indent,title,n,DIM);
266 for (i=0; i<n; i++) {
267 (void) pr_indent(fp,indent);
268 (void) fprintf(fp,"%s[%5d]={",title,i);
269 for (j=0; j<DIM; j++) {
271 (void) fprintf(fp,", ");
272 (void) fprintf(fp,format,vec[i][j]);
274 (void) fprintf(fp,"}\n");
280 void pr_reals(FILE *fp,int indent,const char *title,real *vec,int n)
284 if (available(fp,vec,indent,title)) {
285 (void) pr_indent(fp,indent);
286 (void) fprintf(fp,"%s:\t",title);
288 fprintf(fp," %10g",vec[i]);
289 (void) fprintf(fp,"\n");
293 void pr_doubles(FILE *fp,int indent,const char *title,double *vec,int n)
297 if (available(fp,vec,indent,title)) {
298 (void) pr_indent(fp,indent);
299 (void) fprintf(fp,"%s:\t",title);
301 fprintf(fp," %10g",vec[i]);
302 (void) fprintf(fp,"\n");
306 static void pr_int(FILE *fp,int indent,const char *title,int i)
308 pr_indent(fp,indent);
309 fprintf(fp,"%-20s = %d\n",title,i);
312 static void pr_gmx_large_int(FILE *fp,int indent,const char *title,gmx_large_int_t i)
314 char buf[STEPSTRSIZE];
316 pr_indent(fp,indent);
317 fprintf(fp,"%-20s = %s\n",title,gmx_step_str(i,buf));
320 static void pr_real(FILE *fp,int indent,const char *title,real r)
322 pr_indent(fp,indent);
323 fprintf(fp,"%-20s = %g\n",title,r);
326 static void pr_double(FILE *fp,int indent,const char *title,double d)
328 pr_indent(fp,indent);
329 fprintf(fp,"%-20s = %g\n",title,d);
332 static void pr_str(FILE *fp,int indent,const char *title,const char *s)
334 pr_indent(fp,indent);
335 fprintf(fp,"%-20s = %s\n",title,s);
338 void pr_qm_opts(FILE *fp,int indent,const char *title,t_grpopts *opts)
342 fprintf(fp,"%s:\n",title);
344 pr_int(fp,indent,"ngQM",opts->ngQM);
345 if (opts->ngQM > 0) {
346 pr_ivec(fp,indent,"QMmethod",opts->QMmethod,opts->ngQM,FALSE);
347 pr_ivec(fp,indent,"QMbasis",opts->QMbasis,opts->ngQM,FALSE);
348 pr_ivec(fp,indent,"QMcharge",opts->QMcharge,opts->ngQM,FALSE);
349 pr_ivec(fp,indent,"QMmult",opts->QMmult,opts->ngQM,FALSE);
350 pr_bvec(fp,indent,"bSH",opts->bSH,opts->ngQM,FALSE);
351 pr_ivec(fp,indent,"CASorbitals",opts->CASorbitals,opts->ngQM,FALSE);
352 pr_ivec(fp,indent,"CASelectrons",opts->CASelectrons,opts->ngQM,FALSE);
353 pr_rvec(fp,indent,"SAon",opts->SAon,opts->ngQM,FALSE);
354 pr_rvec(fp,indent,"SAon",opts->SAon,opts->ngQM,FALSE);
355 pr_ivec(fp,indent,"SAsteps",opts->SAsteps,opts->ngQM,FALSE);
356 pr_bvec(fp,indent,"bOPT",opts->bOPT,opts->ngQM,FALSE);
357 pr_bvec(fp,indent,"bTS",opts->bTS,opts->ngQM,FALSE);
361 static void pr_grp_opts(FILE *out,int indent,const char *title,t_grpopts *opts,
367 fprintf(out,"%s:\n",title);
369 pr_indent(out,indent);
370 fprintf(out,"nrdf%s",bMDPformat ? " = " : ":");
371 for(i=0; (i<opts->ngtc); i++)
372 fprintf(out," %10g",opts->nrdf[i]);
375 pr_indent(out,indent);
376 fprintf(out,"ref_t%s",bMDPformat ? " = " : ":");
377 for(i=0; (i<opts->ngtc); i++)
378 fprintf(out," %10g",opts->ref_t[i]);
381 pr_indent(out,indent);
382 fprintf(out,"tau_t%s",bMDPformat ? " = " : ":");
383 for(i=0; (i<opts->ngtc); i++)
384 fprintf(out," %10g",opts->tau_t[i]);
387 /* Pretty-print the simulated annealing info */
388 fprintf(out,"anneal%s",bMDPformat ? " = " : ":");
389 for(i=0; (i<opts->ngtc); i++)
390 fprintf(out," %10s",EANNEAL(opts->annealing[i]));
393 fprintf(out,"ann_npoints%s",bMDPformat ? " = " : ":");
394 for(i=0; (i<opts->ngtc); i++)
395 fprintf(out," %10d",opts->anneal_npoints[i]);
398 for(i=0; (i<opts->ngtc); i++) {
399 if(opts->anneal_npoints[i]>0) {
400 fprintf(out,"ann. times [%d]:\t",i);
401 for(j=0; (j<opts->anneal_npoints[i]); j++)
402 fprintf(out," %10.1f",opts->anneal_time[i][j]);
404 fprintf(out,"ann. temps [%d]:\t",i);
405 for(j=0; (j<opts->anneal_npoints[i]); j++)
406 fprintf(out," %10.1f",opts->anneal_temp[i][j]);
411 pr_indent(out,indent);
412 fprintf(out,"acc:\t");
413 for(i=0; (i<opts->ngacc); i++)
414 for(m=0; (m<DIM); m++)
415 fprintf(out," %10g",opts->acc[i][m]);
418 pr_indent(out,indent);
419 fprintf(out,"nfreeze:");
420 for(i=0; (i<opts->ngfrz); i++)
421 for(m=0; (m<DIM); m++)
422 fprintf(out," %10s",opts->nFreeze[i][m] ? "Y" : "N");
426 for(i=0; (i<opts->ngener); i++) {
427 pr_indent(out,indent);
428 fprintf(out,"energygrp_flags[%3d]:",i);
429 for(m=0; (m<opts->ngener); m++)
430 fprintf(out," %d",opts->egp_flags[opts->ngener*i+m]);
437 static void pr_matrix(FILE *fp,int indent,const char *title,rvec *m,
441 fprintf(fp,"%-10s = %g %g %g %g %g %g\n",title,
442 m[XX][XX],m[YY][YY],m[ZZ][ZZ],m[XX][YY],m[XX][ZZ],m[YY][ZZ]);
444 pr_rvecs(fp,indent,title,m,DIM);
447 static void pr_cosine(FILE *fp,int indent,const char *title,t_cosines *cos,
453 fprintf(fp,"%s = %d\n",title,cos->n);
456 indent=pr_title(fp,indent,title);
457 (void) pr_indent(fp,indent);
458 fprintf(fp,"n = %d\n",cos->n);
460 (void) pr_indent(fp,indent+2);
462 for(j=0; (j<cos->n); j++)
463 fprintf(fp," %e",cos->a[j]);
465 (void) pr_indent(fp,indent+2);
467 for(j=0; (j<cos->n); j++)
468 fprintf(fp," %e",cos->phi[j]);
474 #define PS(t,s) pr_str(fp,indent,t,s)
475 #define PI(t,s) pr_int(fp,indent,t,s)
476 #define PSTEP(t,s) pr_gmx_large_int(fp,indent,t,s)
477 #define PR(t,s) pr_real(fp,indent,t,s)
478 #define PD(t,s) pr_double(fp,indent,t,s)
480 static void pr_pullgrp(FILE *fp,int indent,int g,t_pullgrp *pg)
482 pr_indent(fp,indent);
483 fprintf(fp,"pull_group %d:\n",g);
485 pr_ivec_block(fp,indent,"atom",pg->ind,pg->nat,TRUE);
486 pr_rvec(fp,indent,"weight",pg->weight,pg->nweight,TRUE);
487 PI("pbcatom",pg->pbcatom);
488 pr_rvec(fp,indent,"vec",pg->vec,DIM,TRUE);
489 pr_rvec(fp,indent,"init",pg->init,DIM,TRUE);
495 static void pr_pull(FILE *fp,int indent,t_pull *pull)
499 PS("pull_geometry",EPULLGEOM(pull->eGeom));
500 pr_ivec(fp,indent,"pull_dim",pull->dim,DIM,TRUE);
501 PR("pull_r1",pull->cyl_r1);
502 PR("pull_r0",pull->cyl_r0);
503 PR("pull_constr_tol",pull->constr_tol);
504 PI("pull_nstxout",pull->nstxout);
505 PI("pull_nstfout",pull->nstfout);
506 PI("pull_ngrp",pull->ngrp);
507 for(g=0; g<pull->ngrp+1; g++)
508 pr_pullgrp(fp,indent,g,&pull->grp[g]);
511 static void pr_rotgrp(FILE *fp,int indent,int g,t_rotgrp *rotg)
513 pr_indent(fp,indent);
514 fprintf(fp,"rotation_group %d:\n",g);
516 PS("type",EROTGEOM(rotg->eType));
517 PS("massw",BOOL(rotg->bMassW));
518 pr_ivec_block(fp,indent,"atom",rotg->ind,rotg->nat,TRUE);
519 pr_rvecs(fp,indent,"x_ref",rotg->x_ref,rotg->nat);
520 pr_rvec(fp,indent,"vec",rotg->vec,DIM,TRUE);
521 pr_rvec(fp,indent,"pivot",rotg->pivot,DIM,TRUE);
522 PR("rate",rotg->rate);
524 PR("slab_dist",rotg->slab_dist);
525 PR("min_gaussian",rotg->min_gaussian);
526 PR("epsilon",rotg->eps);
527 PS("fit_method",EROTFIT(rotg->eFittype));
528 PI("potfitangle_nstep",rotg->PotAngle_nstep);
529 PR("potfitangle_step",rotg->PotAngle_step);
532 static void pr_rot(FILE *fp,int indent,t_rot *rot)
536 PI("rot_nstrout",rot->nstrout);
537 PI("rot_nstsout",rot->nstsout);
538 PI("rot_ngrp",rot->ngrp);
539 for(g=0; g<rot->ngrp; g++)
540 pr_rotgrp(fp,indent,g,&rot->grp[g]);
543 void pr_inputrec(FILE *fp,int indent,const char *title,t_inputrec *ir,
546 const char *infbuf="inf";
549 if (available(fp,ir,indent,title)) {
551 indent=pr_title(fp,indent,title);
552 PS("integrator",EI(ir->eI));
553 PSTEP("nsteps",ir->nsteps);
554 PSTEP("init_step",ir->init_step);
555 PS("ns_type",ENS(ir->ns_type));
556 PI("nstlist",ir->nstlist);
557 PI("ndelta",ir->ndelta);
558 PI("nstcomm",ir->nstcomm);
559 PS("comm_mode",ECOM(ir->comm_mode));
560 PI("nstlog",ir->nstlog);
561 PI("nstxout",ir->nstxout);
562 PI("nstvout",ir->nstvout);
563 PI("nstfout",ir->nstfout);
564 PI("nstcalcenergy",ir->nstcalcenergy);
565 PI("nstenergy",ir->nstenergy);
566 PI("nstxtcout",ir->nstxtcout);
567 PR("init_t",ir->init_t);
568 PR("delta_t",ir->delta_t);
570 PR("xtcprec",ir->xtcprec);
574 PI("pme_order",ir->pme_order);
575 PR("ewald_rtol",ir->ewald_rtol);
576 PR("ewald_geometry",ir->ewald_geometry);
577 PR("epsilon_surface",ir->epsilon_surface);
578 PS("optimize_fft",BOOL(ir->bOptFFT));
579 PS("ePBC",EPBC(ir->ePBC));
580 PS("bPeriodicMols",BOOL(ir->bPeriodicMols));
581 PS("bContinuation",BOOL(ir->bContinuation));
582 PS("bShakeSOR",BOOL(ir->bShakeSOR));
583 PS("etc",ETCOUPLTYPE(ir->etc));
584 PI("nsttcouple",ir->nsttcouple);
585 PS("epc",EPCOUPLTYPE(ir->epc));
586 PS("epctype",EPCOUPLTYPETYPE(ir->epct));
587 PI("nstpcouple",ir->nstpcouple);
588 PR("tau_p",ir->tau_p);
589 pr_matrix(fp,indent,"ref_p",ir->ref_p,bMDPformat);
590 pr_matrix(fp,indent,"compress",ir->compress,bMDPformat);
591 PS("refcoord_scaling",EREFSCALINGTYPE(ir->refcoord_scaling));
593 fprintf(fp,"posres_com = %g %g %g\n",ir->posres_com[XX],
594 ir->posres_com[YY],ir->posres_com[ZZ]);
596 pr_rvec(fp,indent,"posres_com",ir->posres_com,DIM,TRUE);
598 fprintf(fp,"posres_comB = %g %g %g\n",ir->posres_comB[XX],
599 ir->posres_comB[YY],ir->posres_comB[ZZ]);
601 pr_rvec(fp,indent,"posres_comB",ir->posres_comB,DIM,TRUE);
602 PI("andersen_seed",ir->andersen_seed);
603 PR("rlist",ir->rlist);
604 PR("rlistlong",ir->rlistlong);
606 PS("coulombtype",EELTYPE(ir->coulombtype));
607 PR("rcoulomb_switch",ir->rcoulomb_switch);
608 PR("rcoulomb",ir->rcoulomb);
609 PS("vdwtype",EVDWTYPE(ir->vdwtype));
610 PR("rvdw_switch",ir->rvdw_switch);
612 if (ir->epsilon_r != 0)
613 PR("epsilon_r",ir->epsilon_r);
615 PS("epsilon_r",infbuf);
616 if (ir->epsilon_rf != 0)
617 PR("epsilon_rf",ir->epsilon_rf);
619 PS("epsilon_rf",infbuf);
620 PR("tabext",ir->tabext);
621 PS("implicit_solvent",EIMPLICITSOL(ir->implicit_solvent));
622 PS("gb_algorithm",EGBALGORITHM(ir->gb_algorithm));
623 PR("gb_epsilon_solvent",ir->gb_epsilon_solvent);
624 PI("nstgbradii",ir->nstgbradii);
625 PR("rgbradii",ir->rgbradii);
626 PR("gb_saltconc",ir->gb_saltconc);
627 PR("gb_obc_alpha",ir->gb_obc_alpha);
628 PR("gb_obc_beta",ir->gb_obc_beta);
629 PR("gb_obc_gamma",ir->gb_obc_gamma);
630 PR("gb_dielectric_offset",ir->gb_dielectric_offset);
631 PS("sa_algorithm",ESAALGORITHM(ir->gb_algorithm));
632 PR("sa_surface_tension",ir->sa_surface_tension);
634 PS("DispCorr",EDISPCORR(ir->eDispCorr));
635 PS("free_energy",EFEPTYPE(ir->efep));
636 PR("init_lambda",ir->init_lambda);
637 PR("delta_lambda",ir->delta_lambda);
640 PI("n_foreign_lambda",ir->n_flambda);
642 if (ir->n_flambda > 0)
644 pr_indent(fp,indent);
645 fprintf(fp,"foreign_lambda%s",bMDPformat ? " = " : ":");
646 for(i=0; i<ir->n_flambda; i++)
648 fprintf(fp," %10g",ir->flambda[i]);
652 PR("sc_alpha",ir->sc_alpha);
653 PI("sc_power",ir->sc_power);
654 PR("sc_sigma",ir->sc_sigma);
655 PR("sc_sigma_min",ir->sc_sigma_min);
656 PI("nstdhdl", ir->nstdhdl);
657 PS("separate_dhdl_file", SEPDHDLFILETYPE(ir->separate_dhdl_file));
658 PS("dhdl_derivatives", DHDLDERIVATIVESTYPE(ir->dhdl_derivatives));
659 PI("dh_hist_size", ir->dh_hist_size);
660 PD("dh_hist_spacing", ir->dh_hist_spacing);
662 PI("nwall",ir->nwall);
663 PS("wall_type",EWALLTYPE(ir->wall_type));
664 PI("wall_atomtype[0]",ir->wall_atomtype[0]);
665 PI("wall_atomtype[1]",ir->wall_atomtype[1]);
666 PR("wall_density[0]",ir->wall_density[0]);
667 PR("wall_density[1]",ir->wall_density[1]);
668 PR("wall_ewald_zfac",ir->wall_ewald_zfac);
670 PS("pull",EPULLTYPE(ir->ePull));
671 if (ir->ePull != epullNO)
672 pr_pull(fp,indent,ir->pull);
674 PS("rotation",BOOL(ir->bRot));
676 pr_rot(fp,indent,ir->rot);
678 PS("disre",EDISRETYPE(ir->eDisre));
679 PS("disre_weighting",EDISREWEIGHTING(ir->eDisreWeighting));
680 PS("disre_mixed",BOOL(ir->bDisreMixed));
681 PR("dr_fc",ir->dr_fc);
682 PR("dr_tau",ir->dr_tau);
683 PR("nstdisreout",ir->nstdisreout);
684 PR("orires_fc",ir->orires_fc);
685 PR("orires_tau",ir->orires_tau);
686 PR("nstorireout",ir->nstorireout);
688 PR("dihre-fc",ir->dihre_fc);
690 PR("em_stepsize",ir->em_stepsize);
691 PR("em_tol",ir->em_tol);
692 PI("niter",ir->niter);
693 PR("fc_stepsize",ir->fc_stepsize);
694 PI("nstcgsteep",ir->nstcgsteep);
695 PI("nbfgscorr",ir->nbfgscorr);
697 PS("ConstAlg",ECONSTRTYPE(ir->eConstrAlg));
698 PR("shake_tol",ir->shake_tol);
699 PI("lincs_order",ir->nProjOrder);
700 PR("lincs_warnangle",ir->LincsWarnAngle);
701 PI("lincs_iter",ir->nLincsIter);
702 PR("bd_fric",ir->bd_fric);
703 PI("ld_seed",ir->ld_seed);
704 PR("cos_accel",ir->cos_accel);
705 pr_matrix(fp,indent,"deform",ir->deform,bMDPformat);
706 PI("userint1",ir->userint1);
707 PI("userint2",ir->userint2);
708 PI("userint3",ir->userint3);
709 PI("userint4",ir->userint4);
710 PR("userreal1",ir->userreal1);
711 PR("userreal2",ir->userreal2);
712 PR("userreal3",ir->userreal3);
713 PR("userreal4",ir->userreal4);
714 pr_grp_opts(fp,indent,"grpopts",&(ir->opts),bMDPformat);
715 pr_cosine(fp,indent,"efield-x",&(ir->ex[XX]),bMDPformat);
716 pr_cosine(fp,indent,"efield-xt",&(ir->et[XX]),bMDPformat);
717 pr_cosine(fp,indent,"efield-y",&(ir->ex[YY]),bMDPformat);
718 pr_cosine(fp,indent,"efield-yt",&(ir->et[YY]),bMDPformat);
719 pr_cosine(fp,indent,"efield-z",&(ir->ex[ZZ]),bMDPformat);
720 pr_cosine(fp,indent,"efield-zt",&(ir->et[ZZ]),bMDPformat);
721 PS("bQMMM",BOOL(ir->bQMMM));
722 PI("QMconstraints",ir->QMconstraints);
723 PI("QMMMscheme",ir->QMMMscheme);
724 PR("scalefactor",ir->scalefactor);
725 pr_qm_opts(fp,indent,"qm_opts",&(ir->opts));
732 static void pr_harm(FILE *fp,t_iparams *iparams,const char *r,const char *kr)
734 fprintf(fp,"%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
735 r,iparams->harmonic.rA,kr,iparams->harmonic.krA,
736 r,iparams->harmonic.rB,kr,iparams->harmonic.krB);
739 void pr_iparams(FILE *fp,t_functype ftype,t_iparams *iparams)
742 real VA[4],VB[4],*rbcA,*rbcB;
747 pr_harm(fp,iparams,"th","ct");
749 case F_CROSS_BOND_BONDS:
750 fprintf(fp,"r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
751 iparams->cross_bb.r1e,iparams->cross_bb.r2e,
752 iparams->cross_bb.krr);
754 case F_CROSS_BOND_ANGLES:
755 fprintf(fp,"r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
756 iparams->cross_ba.r1e,iparams->cross_ba.r2e,
757 iparams->cross_ba.r3e,iparams->cross_ba.krt);
760 fprintf(fp,"theta=%15.8e, ktheta=%15.8e, r13=%15.8e, kUB=%15.8e\n",
761 iparams->u_b.theta,iparams->u_b.ktheta,iparams->u_b.r13,iparams->u_b.kUB);
763 case F_QUARTIC_ANGLES:
764 fprintf(fp,"theta=%15.8e",iparams->qangle.theta);
766 fprintf(fp,", c%c=%15.8e",'0'+i,iparams->qangle.c[i]);
770 fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
771 iparams->bham.a,iparams->bham.b,iparams->bham.c);
776 pr_harm(fp,iparams,"b0","cb");
779 pr_harm(fp,iparams,"xi","cx");
782 fprintf(fp,"b0=%15.8e, cb=%15.8e, beta=%15.8e\n",
783 iparams->morse.b0,iparams->morse.cb,iparams->morse.beta);
786 fprintf(fp,"b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
787 iparams->cubic.b0,iparams->cubic.kb,iparams->cubic.kcub);
793 fprintf(fp,"bm=%15.8e, kb=%15.8e\n",iparams->fene.bm,iparams->fene.kb);
796 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",
797 iparams->restraint.lowA,iparams->restraint.up1A,
798 iparams->restraint.up2A,iparams->restraint.kA,
799 iparams->restraint.lowB,iparams->restraint.up1B,
800 iparams->restraint.up2B,iparams->restraint.kB);
806 fprintf(fp,"tab=%d, kA=%15.8e, kB=%15.8e\n",
807 iparams->tab.table,iparams->tab.kA,iparams->tab.kB);
810 fprintf(fp,"alpha=%15.8e\n",iparams->polarize.alpha);
813 fprintf(fp,"a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
814 iparams->thole.a,iparams->thole.alpha1,iparams->thole.alpha2,
815 iparams->thole.rfac);
818 fprintf(fp,"al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
819 iparams->wpol.al_x,iparams->wpol.al_y,iparams->wpol.al_z,
820 iparams->wpol.rOH,iparams->wpol.rHH,iparams->wpol.rOD);
823 fprintf(fp,"c6=%15.8e, c12=%15.8e\n",iparams->lj.c6,iparams->lj.c12);
826 fprintf(fp,"c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
827 iparams->lj14.c6A,iparams->lj14.c12A,
828 iparams->lj14.c6B,iparams->lj14.c12B);
831 fprintf(fp,"fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
833 iparams->ljc14.qi,iparams->ljc14.qj,
834 iparams->ljc14.c6,iparams->ljc14.c12);
837 fprintf(fp,"qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
838 iparams->ljcnb.qi,iparams->ljcnb.qj,
839 iparams->ljcnb.c6,iparams->ljcnb.c12);
845 fprintf(fp,"phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
846 iparams->pdihs.phiA,iparams->pdihs.cpA,
847 iparams->pdihs.phiB,iparams->pdihs.cpB,
848 iparams->pdihs.mult);
851 fprintf(fp,"label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
852 iparams->disres.label,iparams->disres.type,
853 iparams->disres.low,iparams->disres.up1,
854 iparams->disres.up2,iparams->disres.kfac);
857 fprintf(fp,"ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
858 iparams->orires.ex,iparams->orires.label,iparams->orires.power,
859 iparams->orires.c,iparams->orires.obs,iparams->orires.kfac);
862 fprintf(fp,"label=%d, power=%4d phi=%15.8e, dphi=%15.8e, kfac=%15.8e)\n",
863 iparams->dihres.label,iparams->dihres.power,
864 iparams->dihres.phi,iparams->dihres.dphi,iparams->dihres.kfac);
867 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",
868 iparams->posres.pos0A[XX],iparams->posres.pos0A[YY],
869 iparams->posres.pos0A[ZZ],iparams->posres.fcA[XX],
870 iparams->posres.fcA[YY],iparams->posres.fcA[ZZ],
871 iparams->posres.pos0B[XX],iparams->posres.pos0B[YY],
872 iparams->posres.pos0B[ZZ],iparams->posres.fcB[XX],
873 iparams->posres.fcB[YY],iparams->posres.fcB[ZZ]);
876 for (i=0; i<NR_RBDIHS; i++)
877 fprintf(fp,"%srbcA[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcA[i]);
879 for (i=0; i<NR_RBDIHS; i++)
880 fprintf(fp,"%srbcB[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcB[i]);
884 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
885 * OPLS potential constants back.
887 rbcA = iparams->rbdihs.rbcA;
888 rbcB = iparams->rbdihs.rbcB;
890 VA[3] = -0.25*rbcA[4];
891 VA[2] = -0.5*rbcA[3];
892 VA[1] = 4.0*VA[3]-rbcA[2];
893 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
895 VB[3] = -0.25*rbcB[4];
896 VB[2] = -0.5*rbcB[3];
897 VB[1] = 4.0*VB[3]-rbcB[2];
898 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
900 for (i=0; i<NR_FOURDIHS; i++)
901 fprintf(fp,"%sFourA[%d]=%15.8e",i==0?"":", ",i,VA[i]);
903 for (i=0; i<NR_FOURDIHS; i++)
904 fprintf(fp,"%sFourB[%d]=%15.8e",i==0?"":", ",i,VB[i]);
910 fprintf(fp,"dA=%15.8e, dB=%15.8e\n",iparams->constr.dA,iparams->constr.dB);
913 fprintf(fp,"doh=%15.8e, dhh=%15.8e\n",iparams->settle.doh,
914 iparams->settle.dhh);
917 fprintf(fp,"a=%15.8e\n",iparams->vsite.a);
922 fprintf(fp,"a=%15.8e, b=%15.8e\n",iparams->vsite.a,iparams->vsite.b);
927 fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
928 iparams->vsite.a,iparams->vsite.b,iparams->vsite.c);
931 fprintf(fp,"n=%2d, a=%15.8e\n",iparams->vsiten.n,iparams->vsiten.a);
936 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);
939 fprintf(fp, "cmapA=%1d, cmapB=%1d\n",iparams->cmap.cmapA, iparams->cmap.cmapB);
942 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
943 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
947 void pr_ilist(FILE *fp,int indent,const char *title,
948 t_functype *functype,t_ilist *ilist, gmx_bool bShowNumbers)
950 int i,j,k,type,ftype;
953 if (available(fp,ilist,indent,title) && ilist->nr > 0)
955 indent=pr_title(fp,indent,title);
956 (void) pr_indent(fp,indent);
957 fprintf(fp,"nr: %d\n",ilist->nr);
959 (void) pr_indent(fp,indent);
960 fprintf(fp,"iatoms:\n");
961 iatoms=ilist->iatoms;
962 for (i=j=0; i<ilist->nr;) {
964 (void) pr_indent(fp,indent+INDENT);
966 ftype=functype[type];
967 (void) fprintf(fp,"%d type=%d (%s)",
968 bShowNumbers?j:-1,bShowNumbers?type:-1,
969 interaction_function[ftype].name);
971 for (k=0; k<interaction_function[ftype].nratoms; k++)
972 (void) fprintf(fp," %u",*(iatoms++));
973 (void) fprintf(fp,"\n");
974 i+=1+interaction_function[ftype].nratoms;
976 fprintf(fp,"%5d%5d\n",i,iatoms[i]);
984 static void pr_cmap(FILE *fp, int indent, const char *title,
985 gmx_cmap_t *cmap_grid, gmx_bool bShowNumbers)
990 dx = 360.0 / cmap_grid->grid_spacing;
991 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
993 if(available(fp,cmap_grid,indent,title))
995 fprintf(fp,"%s\n",title);
997 for(i=0;i<cmap_grid->ngrid;i++)
1000 fprintf(fp,"%8s %8s %8s %8s\n","V","dVdx","dVdy","d2dV");
1002 fprintf(fp,"grid[%3d]={\n",bShowNumbers?i:-1);
1004 for(j=0;j<nelem;j++)
1006 if( (j%cmap_grid->grid_spacing)==0)
1008 fprintf(fp,"%8.1f\n",idx);
1012 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4]);
1013 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+1]);
1014 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+2]);
1015 fprintf(fp,"%8.3f\n",cmap_grid->cmapdata[i].cmap[j*4+3]);
1023 void pr_ffparams(FILE *fp,int indent,const char *title,
1024 gmx_ffparams_t *ffparams,
1025 gmx_bool bShowNumbers)
1029 indent=pr_title(fp,indent,title);
1030 (void) pr_indent(fp,indent);
1031 (void) fprintf(fp,"atnr=%d\n",ffparams->atnr);
1032 (void) pr_indent(fp,indent);
1033 (void) fprintf(fp,"ntypes=%d\n",ffparams->ntypes);
1034 for (i=0; i<ffparams->ntypes; i++) {
1035 (void) pr_indent(fp,indent+INDENT);
1036 (void) fprintf(fp,"functype[%d]=%s, ",
1038 interaction_function[ffparams->functype[i]].name);
1039 pr_iparams(fp,ffparams->functype[i],&ffparams->iparams[i]);
1041 (void) pr_double(fp,indent,"reppow",ffparams->reppow);
1042 (void) pr_real(fp,indent,"fudgeQQ",ffparams->fudgeQQ);
1043 pr_cmap(fp,indent,"cmap",&ffparams->cmap_grid,bShowNumbers);
1046 void pr_idef(FILE *fp,int indent,const char *title,t_idef *idef, gmx_bool bShowNumbers)
1050 if (available(fp,idef,indent,title)) {
1051 indent=pr_title(fp,indent,title);
1052 (void) pr_indent(fp,indent);
1053 (void) fprintf(fp,"atnr=%d\n",idef->atnr);
1054 (void) pr_indent(fp,indent);
1055 (void) fprintf(fp,"ntypes=%d\n",idef->ntypes);
1056 for (i=0; i<idef->ntypes; i++) {
1057 (void) pr_indent(fp,indent+INDENT);
1058 (void) fprintf(fp,"functype[%d]=%s, ",
1060 interaction_function[idef->functype[i]].name);
1061 pr_iparams(fp,idef->functype[i],&idef->iparams[i]);
1063 (void) pr_real(fp,indent,"fudgeQQ",idef->fudgeQQ);
1065 for(j=0; (j<F_NRE); j++)
1066 pr_ilist(fp,indent,interaction_function[j].longname,
1067 idef->functype,&idef->il[j],bShowNumbers);
1071 static int pr_block_title(FILE *fp,int indent,const char *title,t_block *block)
1075 if (available(fp,block,indent,title))
1077 indent=pr_title(fp,indent,title);
1078 (void) pr_indent(fp,indent);
1079 (void) fprintf(fp,"nr=%d\n",block->nr);
1084 static int pr_blocka_title(FILE *fp,int indent,const char *title,t_blocka *block)
1088 if (available(fp,block,indent,title))
1090 indent=pr_title(fp,indent,title);
1091 (void) pr_indent(fp,indent);
1092 (void) fprintf(fp,"nr=%d\n",block->nr);
1093 (void) pr_indent(fp,indent);
1094 (void) fprintf(fp,"nra=%d\n",block->nra);
1099 static void low_pr_block(FILE *fp,int indent,const char *title,t_block *block, gmx_bool bShowNumbers)
1103 if (available(fp,block,indent,title))
1105 indent=pr_block_title(fp,indent,title,block);
1106 for (i=0; i<=block->nr; i++)
1108 (void) pr_indent(fp,indent+INDENT);
1109 (void) fprintf(fp,"%s->index[%d]=%u\n",
1110 title,bShowNumbers?i:-1,block->index[i]);
1115 static void low_pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block, gmx_bool bShowNumbers)
1119 if (available(fp,block,indent,title))
1121 indent=pr_blocka_title(fp,indent,title,block);
1122 for (i=0; i<=block->nr; i++)
1124 (void) pr_indent(fp,indent+INDENT);
1125 (void) fprintf(fp,"%s->index[%d]=%u\n",
1126 title,bShowNumbers?i:-1,block->index[i]);
1128 for (i=0; i<block->nra; i++)
1130 (void) pr_indent(fp,indent+INDENT);
1131 (void) fprintf(fp,"%s->a[%d]=%u\n",
1132 title,bShowNumbers?i:-1,block->a[i]);
1137 void pr_block(FILE *fp,int indent,const char *title,t_block *block,gmx_bool bShowNumbers)
1139 int i,j,ok,size,start,end;
1141 if (available(fp,block,indent,title))
1143 indent=pr_block_title(fp,indent,title,block);
1146 if ((ok=(block->index[start]==0))==0)
1147 (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1149 for (i=0; i<block->nr; i++)
1151 end=block->index[i+1];
1152 size=pr_indent(fp,indent);
1154 size+=fprintf(fp,"%s[%d]={}\n",title,i);
1156 size+=fprintf(fp,"%s[%d]={%d..%d}\n",
1157 title,bShowNumbers?i:-1,
1158 bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1164 void pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block,gmx_bool bShowNumbers)
1166 int i,j,ok,size,start,end;
1168 if (available(fp,block,indent,title))
1170 indent=pr_blocka_title(fp,indent,title,block);
1173 if ((ok=(block->index[start]==0))==0)
1174 (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1176 for (i=0; i<block->nr; i++)
1178 end=block->index[i+1];
1179 size=pr_indent(fp,indent);
1181 size+=fprintf(fp,"%s[%d]={",title,i);
1183 size+=fprintf(fp,"%s[%d][%d..%d]={",
1184 title,bShowNumbers?i:-1,
1185 bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1186 for (j=start; j<end; j++)
1188 if (j>start) size+=fprintf(fp,", ");
1189 if ((size)>(USE_WIDTH))
1191 (void) fprintf(fp,"\n");
1192 size=pr_indent(fp,indent+INDENT);
1194 size+=fprintf(fp,"%u",block->a[j]);
1196 (void) fprintf(fp,"}\n");
1199 if ((end!=block->nra)||(!ok))
1201 (void) pr_indent(fp,indent);
1202 (void) fprintf(fp,"tables inconsistent, dumping complete tables:\n");
1203 low_pr_blocka(fp,indent,title,block,bShowNumbers);
1208 static void pr_strings(FILE *fp,int indent,const char *title,char ***nm,int n, gmx_bool bShowNumbers)
1212 if (available(fp,nm,indent,title))
1214 indent=pr_title_n(fp,indent,title,n);
1217 (void) pr_indent(fp,indent);
1218 (void) fprintf(fp,"%s[%d]={name=\"%s\"}\n",
1219 title,bShowNumbers?i:-1,*(nm[i]));
1224 static void pr_strings2(FILE *fp,int indent,const char *title,
1225 char ***nm,char ***nmB,int n, gmx_bool bShowNumbers)
1229 if (available(fp,nm,indent,title))
1231 indent=pr_title_n(fp,indent,title,n);
1234 (void) pr_indent(fp,indent);
1235 (void) fprintf(fp,"%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1236 title,bShowNumbers?i:-1,*(nm[i]),*(nmB[i]));
1241 static void pr_resinfo(FILE *fp,int indent,const char *title,t_resinfo *resinfo,int n, gmx_bool bShowNumbers)
1245 if (available(fp,resinfo,indent,title))
1247 indent=pr_title_n(fp,indent,title,n);
1250 (void) pr_indent(fp,indent);
1251 (void) fprintf(fp,"%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1252 title,bShowNumbers?i:-1,
1253 *(resinfo[i].name),resinfo[i].nr,
1254 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1259 static void pr_atom(FILE *fp,int indent,const char *title,t_atom *atom,int n)
1263 if (available(fp,atom,indent,title)) {
1264 indent=pr_title_n(fp,indent,title,n);
1265 for (i=0; i<n; i++) {
1266 (void) pr_indent(fp,indent);
1267 fprintf(fp,"%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1268 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1269 title,i,atom[i].type,atom[i].typeB,ptype_str[atom[i].ptype],
1270 atom[i].m,atom[i].q,atom[i].mB,atom[i].qB,
1271 atom[i].resind,atom[i].atomnumber);
1276 static void pr_grps(FILE *fp,int indent,const char *title,t_grps grps[],
1277 char **grpname[], gmx_bool bShowNumbers)
1281 for(i=0; (i<egcNR); i++)
1283 fprintf(fp,"%s[%-12s] nr=%d, name=[",title,gtypes[i],grps[i].nr);
1284 for(j=0; (j<grps[i].nr); j++)
1286 fprintf(fp," %s",*(grpname[grps[i].nm_ind[j]]));
1292 static void pr_groups(FILE *fp,int indent,const char *title,
1293 gmx_groups_t *groups,
1294 gmx_bool bShowNumbers)
1299 pr_grps(fp,indent,"grp",groups->grps,groups->grpname,bShowNumbers);
1300 pr_strings(fp,indent,"grpname",groups->grpname,groups->ngrpname,bShowNumbers);
1302 (void) pr_indent(fp,indent);
1303 fprintf(fp,"groups ");
1304 for(g=0; g<egcNR; g++)
1306 printf(" %5.5s",gtypes[g]);
1310 (void) pr_indent(fp,indent);
1311 fprintf(fp,"allocated ");
1313 for(g=0; g<egcNR; g++)
1315 printf(" %5d",groups->ngrpnr[g]);
1316 nat_max = max(nat_max,groups->ngrpnr[g]);
1322 (void) pr_indent(fp,indent);
1323 fprintf(fp,"groupnr[%5s] =","*");
1324 for(g=0; g<egcNR; g++)
1326 fprintf(fp," %3d ",0);
1332 for(i=0; i<nat_max; i++)
1334 (void) pr_indent(fp,indent);
1335 fprintf(fp,"groupnr[%5d] =",i);
1336 for(g=0; g<egcNR; g++)
1339 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1346 void pr_atoms(FILE *fp,int indent,const char *title,t_atoms *atoms,
1347 gmx_bool bShownumbers)
1349 if (available(fp,atoms,indent,title))
1351 indent=pr_title(fp,indent,title);
1352 pr_atom(fp,indent,"atom",atoms->atom,atoms->nr);
1353 pr_strings(fp,indent,"atom",atoms->atomname,atoms->nr,bShownumbers);
1354 pr_strings2(fp,indent,"type",atoms->atomtype,atoms->atomtypeB,atoms->nr,bShownumbers);
1355 pr_resinfo(fp,indent,"residue",atoms->resinfo,atoms->nres,bShownumbers);
1360 void pr_atomtypes(FILE *fp,int indent,const char *title,t_atomtypes *atomtypes,
1361 gmx_bool bShowNumbers)
1364 if (available(fp,atomtypes,indent,title))
1366 indent=pr_title(fp,indent,title);
1367 for(i=0;i<atomtypes->nr;i++) {
1368 pr_indent(fp,indent);
1370 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1371 bShowNumbers?i:-1,atomtypes->radius[i],atomtypes->vol[i],
1372 atomtypes->gb_radius[i],
1373 atomtypes->surftens[i],atomtypes->atomnumber[i],atomtypes->S_hct[i]);
1378 static void pr_moltype(FILE *fp,int indent,const char *title,
1379 gmx_moltype_t *molt,int n,
1380 gmx_ffparams_t *ffparams,
1381 gmx_bool bShowNumbers)
1385 indent = pr_title_n(fp,indent,title,n);
1386 (void) pr_indent(fp,indent);
1387 (void) fprintf(fp,"name=\"%s\"\n",*(molt->name));
1388 pr_atoms(fp,indent,"atoms",&(molt->atoms),bShowNumbers);
1389 pr_block(fp,indent,"cgs",&molt->cgs, bShowNumbers);
1390 pr_blocka(fp,indent,"excls",&molt->excls, bShowNumbers);
1391 for(j=0; (j<F_NRE); j++) {
1392 pr_ilist(fp,indent,interaction_function[j].longname,
1393 ffparams->functype,&molt->ilist[j],bShowNumbers);
1397 static void pr_molblock(FILE *fp,int indent,const char *title,
1398 gmx_molblock_t *molb,int n,
1399 gmx_moltype_t *molt,
1400 gmx_bool bShowNumbers)
1402 indent = pr_title_n(fp,indent,title,n);
1403 (void) pr_indent(fp,indent);
1404 (void) fprintf(fp,"%-20s = %d \"%s\"\n",
1405 "moltype",molb->type,*(molt[molb->type].name));
1406 pr_int(fp,indent,"#molecules",molb->nmol);
1407 pr_int(fp,indent,"#atoms_mol",molb->natoms_mol);
1408 pr_int(fp,indent,"#posres_xA",molb->nposres_xA);
1409 if (molb->nposres_xA > 0) {
1410 pr_rvecs(fp,indent,"posres_xA",molb->posres_xA,molb->nposres_xA);
1412 pr_int(fp,indent,"#posres_xB",molb->nposres_xB);
1413 if (molb->nposres_xB > 0) {
1414 pr_rvecs(fp,indent,"posres_xB",molb->posres_xB,molb->nposres_xB);
1418 void pr_mtop(FILE *fp,int indent,const char *title,gmx_mtop_t *mtop,
1419 gmx_bool bShowNumbers)
1423 if (available(fp,mtop,indent,title)) {
1424 indent=pr_title(fp,indent,title);
1425 (void) pr_indent(fp,indent);
1426 (void) fprintf(fp,"name=\"%s\"\n",*(mtop->name));
1427 pr_int(fp,indent,"#atoms",mtop->natoms);
1428 for(mb=0; mb<mtop->nmolblock; mb++) {
1429 pr_molblock(fp,indent,"molblock",&mtop->molblock[mb],mb,
1430 mtop->moltype,bShowNumbers);
1432 pr_ffparams(fp,indent,"ffparams",&(mtop->ffparams),bShowNumbers);
1433 pr_atomtypes(fp,indent,"atomtypes",&(mtop->atomtypes),bShowNumbers);
1434 for(mt=0; mt<mtop->nmoltype; mt++) {
1435 pr_moltype(fp,indent,"moltype",&mtop->moltype[mt],mt,
1436 &mtop->ffparams,bShowNumbers);
1438 pr_groups(fp,indent,"groups",&mtop->groups,bShowNumbers);
1442 void pr_top(FILE *fp,int indent,const char *title,t_topology *top, gmx_bool bShowNumbers)
1444 if (available(fp,top,indent,title)) {
1445 indent=pr_title(fp,indent,title);
1446 (void) pr_indent(fp,indent);
1447 (void) fprintf(fp,"name=\"%s\"\n",*(top->name));
1448 pr_atoms(fp,indent,"atoms",&(top->atoms),bShowNumbers);
1449 pr_atomtypes(fp,indent,"atomtypes",&(top->atomtypes),bShowNumbers);
1450 pr_block(fp,indent,"cgs",&top->cgs, bShowNumbers);
1451 pr_block(fp,indent,"mols",&top->mols, bShowNumbers);
1452 pr_blocka(fp,indent,"excls",&top->excls, bShowNumbers);
1453 pr_idef(fp,indent,"idef",&top->idef,bShowNumbers);
1457 void pr_header(FILE *fp,int indent,const char *title,t_tpxheader *sh)
1461 if (available(fp,sh,indent,title))
1463 indent=pr_title(fp,indent,title);
1464 pr_indent(fp,indent);
1465 fprintf(fp,"bIr = %spresent\n",sh->bIr?"":"not ");
1466 pr_indent(fp,indent);
1467 fprintf(fp,"bBox = %spresent\n",sh->bBox?"":"not ");
1468 pr_indent(fp,indent);
1469 fprintf(fp,"bTop = %spresent\n",sh->bTop?"":"not ");
1470 pr_indent(fp,indent);
1471 fprintf(fp,"bX = %spresent\n",sh->bX?"":"not ");
1472 pr_indent(fp,indent);
1473 fprintf(fp,"bV = %spresent\n",sh->bV?"":"not ");
1474 pr_indent(fp,indent);
1475 fprintf(fp,"bF = %spresent\n",sh->bF?"":"not ");
1477 pr_indent(fp,indent);
1478 fprintf(fp,"natoms = %d\n",sh->natoms);
1479 pr_indent(fp,indent);
1480 fprintf(fp,"lambda = %e\n",sh->lambda);
1484 void pr_commrec(FILE *fp,int indent,t_commrec *cr)
1486 pr_indent(fp,indent);
1487 fprintf(fp,"commrec:\n");
1489 pr_indent(fp,indent);
1490 fprintf(fp,"nodeid = %d\n",cr->nodeid);
1491 pr_indent(fp,indent);
1492 fprintf(fp,"nnodes = %d\n",cr->nnodes);
1493 pr_indent(fp,indent);
1494 fprintf(fp,"npmenodes = %d\n",cr->npmenodes);
1496 pr_indent(fp,indent);
1497 fprintf(fp,"threadid = %d\n",cr->threadid);
1498 pr_indent(fp,indent);
1499 fprintf(fp,"nthreads = %d\n",cr->nthreads);