Don't use POSIX fnmatch() for pattern matching.
[gromacs/qmmm-gamess-us.git] / src / gmxlib / txtdump.c
blob6a5b9b849848197abc3da7f00baa5b0a17e5b736
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
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
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 /* This file is completely threadsafe - please keep it that way! */
41 #include <thread_mpi.h>
44 #include <stdio.h>
45 #include "smalloc.h"
46 #include "typedefs.h"
47 #include "names.h"
48 #include "txtdump.h"
49 #include "string2.h"
50 #include "vec.h"
53 int pr_indent(FILE *fp,int n)
55 int i;
57 for (i=0; i<n; i++) (void) fprintf(fp," ");
58 return n;
61 int available(FILE *fp,void *p,int indent,const char *title)
63 if (!p) {
64 if (indent > 0)
65 pr_indent(fp,indent);
66 (void) fprintf(fp,"%s: not available\n",title);
68 return (p!=NULL);
71 int pr_title(FILE *fp,int indent,const char *title)
73 (void) pr_indent(fp,indent);
74 (void) fprintf(fp,"%s:\n",title);
75 return (indent+INDENT);
78 int pr_title_n(FILE *fp,int indent,const char *title,int n)
80 (void) pr_indent(fp,indent);
81 (void) fprintf(fp,"%s (%d):\n",title,n);
82 return (indent+INDENT);
85 int pr_title_nxn(FILE *fp,int indent,const char *title,int n1,int n2)
87 (void) pr_indent(fp,indent);
88 (void) fprintf(fp,"%s (%dx%d):\n",title,n1,n2);
89 return (indent+INDENT);
92 void pr_ivec(FILE *fp,int indent,const char *title,int vec[],int n, bool bShowNumbers)
94 int i;
96 if (available(fp,vec,indent,title))
98 indent=pr_title_n(fp,indent,title,n);
99 for (i=0; i<n; i++)
101 (void) pr_indent(fp,indent);
102 (void) fprintf(fp,"%s[%d]=%d\n",title,bShowNumbers?i:-1,vec[i]);
107 void pr_ivec_block(FILE *fp,int indent,const char *title,int vec[],int n, bool bShowNumbers)
109 int i,j;
111 if (available(fp,vec,indent,title))
113 indent=pr_title_n(fp,indent,title,n);
114 i = 0;
115 while (i < n)
117 j = i+1;
118 while (j < n && vec[j] == vec[j-1]+1)
120 j++;
122 /* Print consecutive groups of 3 or more as blocks */
123 if (j - i < 3)
125 while(i < j)
127 (void) pr_indent(fp,indent);
128 (void) fprintf(fp,"%s[%d]=%d\n",
129 title,bShowNumbers?i:-1,vec[i]);
130 i++;
133 else
135 (void) pr_indent(fp,indent);
136 (void) fprintf(fp,"%s[%d,...,%d] = {%d,...,%d}\n",
137 title,
138 bShowNumbers?i:-1,
139 bShowNumbers?j-1:-1,
140 vec[i],vec[j-1]);
141 i = j;
147 void pr_bvec(FILE *fp,int indent,const char *title,bool vec[],int n, bool bShowNumbers)
149 int i;
151 if (available(fp,vec,indent,title))
153 indent=pr_title_n(fp,indent,title,n);
154 for (i=0; i<n; i++)
156 (void) pr_indent(fp,indent);
157 (void) fprintf(fp,"%s[%d]=%s\n",title,bShowNumbers?i:-1,
158 BOOL(vec[i]));
163 void pr_ivecs(FILE *fp,int indent,const char *title,ivec vec[],int n, bool bShowNumbers)
165 int i,j;
167 if (available(fp,vec,indent,title))
169 indent=pr_title_nxn(fp,indent,title,n,DIM);
170 for (i=0; i<n; i++)
172 (void) pr_indent(fp,indent);
173 (void) fprintf(fp,"%s[%d]={",title,bShowNumbers?i:-1);
174 for (j=0; j<DIM; j++)
176 if (j!=0) (void) fprintf(fp,", ");
177 fprintf(fp,"%d",vec[i][j]);
179 (void) fprintf(fp,"}\n");
184 void pr_rvec(FILE *fp,int indent,const char *title,real vec[],int n, bool bShowNumbers)
186 int i;
188 if (available(fp,vec,indent,title))
190 indent=pr_title_n(fp,indent,title,n);
191 for (i=0; i<n; i++)
193 pr_indent(fp,indent);
194 fprintf(fp,"%s[%d]=%12.5e\n",title,bShowNumbers?i:-1,vec[i]);
199 void pr_dvec(FILE *fp,int indent,const char *title,double vec[],int n, bool bShowNumbers)
201 int i;
203 if (available(fp,vec,indent,title))
205 indent=pr_title_n(fp,indent,title,n);
206 for (i=0; i<n; i++)
208 pr_indent(fp,indent);
209 fprintf(fp,"%s[%d]=%12.5e\n",title,bShowNumbers?i:-1,vec[i]);
216 void pr_mat(FILE *fp,int indent,char *title,matrix m)
218 int i,j;
220 if (available(fp,m,indent,title)) {
221 indent=pr_title_n(fp,indent,title,n);
222 for(i=0; i<n; i++) {
223 pr_indent(fp,indent);
224 fprintf(fp,"%s[%d]=%12.5e %12.5e %12.5e\n",
225 title,bShowNumbers?i:-1,m[i][XX],m[i][YY],m[i][ZZ]);
231 void pr_rvecs_len(FILE *fp,int indent,const char *title,rvec vec[],int n)
233 int i,j;
235 if (available(fp,vec,indent,title)) {
236 indent=pr_title_nxn(fp,indent,title,n,DIM);
237 for (i=0; i<n; i++) {
238 (void) pr_indent(fp,indent);
239 (void) fprintf(fp,"%s[%5d]={",title,i);
240 for (j=0; j<DIM; j++) {
241 if (j != 0)
242 (void) fprintf(fp,", ");
243 (void) fprintf(fp,"%12.5e",vec[i][j]);
245 (void) fprintf(fp,"} len=%12.5e\n",norm(vec[i]));
250 void pr_rvecs(FILE *fp,int indent,const char *title,rvec vec[],int n)
252 const char *fshort = "%12.5e";
253 const char *flong = "%15.8e";
254 const char *format;
255 int i,j;
257 if (getenv("LONGFORMAT") != NULL)
258 format = flong;
259 else
260 format = fshort;
262 if (available(fp,vec,indent,title)) {
263 indent=pr_title_nxn(fp,indent,title,n,DIM);
264 for (i=0; i<n; i++) {
265 (void) pr_indent(fp,indent);
266 (void) fprintf(fp,"%s[%5d]={",title,i);
267 for (j=0; j<DIM; j++) {
268 if (j != 0)
269 (void) fprintf(fp,", ");
270 (void) fprintf(fp,format,vec[i][j]);
272 (void) fprintf(fp,"}\n");
278 void pr_reals(FILE *fp,int indent,const char *title,real *vec,int n)
280 int i;
282 if (available(fp,vec,indent,title)) {
283 (void) pr_indent(fp,indent);
284 (void) fprintf(fp,"%s:\t",title);
285 for(i=0; i<n; i++)
286 fprintf(fp," %10g",vec[i]);
287 (void) fprintf(fp,"\n");
291 void pr_doubles(FILE *fp,int indent,const char *title,double *vec,int n)
293 int i;
295 if (available(fp,vec,indent,title)) {
296 (void) pr_indent(fp,indent);
297 (void) fprintf(fp,"%s:\t",title);
298 for(i=0; i<n; i++)
299 fprintf(fp," %10g",vec[i]);
300 (void) fprintf(fp,"\n");
304 static void pr_int(FILE *fp,int indent,const char *title,int i)
306 pr_indent(fp,indent);
307 fprintf(fp,"%-20s = %d\n",title,i);
310 static void pr_gmx_large_int(FILE *fp,int indent,const char *title,gmx_large_int_t i)
312 char buf[22];
314 pr_indent(fp,indent);
315 fprintf(fp,"%-20s = %s\n",title,gmx_step_str(i,buf));
318 static void pr_real(FILE *fp,int indent,const char *title,real r)
320 pr_indent(fp,indent);
321 fprintf(fp,"%-20s = %g\n",title,r);
324 static void pr_double(FILE *fp,int indent,const char *title,double d)
326 pr_indent(fp,indent);
327 fprintf(fp,"%-20s = %g\n",title,d);
330 static void pr_str(FILE *fp,int indent,const char *title,const char *s)
332 pr_indent(fp,indent);
333 fprintf(fp,"%-20s = %s\n",title,s);
336 void pr_qm_opts(FILE *fp,int indent,const char *title,t_grpopts *opts)
338 int i,m,j;
340 fprintf(fp,"%s:\n",title);
342 pr_int(fp,indent,"ngQM",opts->ngQM);
343 if (opts->ngQM > 0) {
344 pr_ivec(fp,indent,"QMmethod",opts->QMmethod,opts->ngQM,FALSE);
345 pr_ivec(fp,indent,"QMbasis",opts->QMbasis,opts->ngQM,FALSE);
346 pr_ivec(fp,indent,"QMcharge",opts->QMcharge,opts->ngQM,FALSE);
347 pr_ivec(fp,indent,"QMmult",opts->QMmult,opts->ngQM,FALSE);
348 pr_bvec(fp,indent,"bSH",opts->bSH,opts->ngQM,FALSE);
349 pr_ivec(fp,indent,"CASorbitals",opts->CASorbitals,opts->ngQM,FALSE);
350 pr_ivec(fp,indent,"CASelectrons",opts->CASelectrons,opts->ngQM,FALSE);
351 pr_rvec(fp,indent,"SAon",opts->SAon,opts->ngQM,FALSE);
352 pr_rvec(fp,indent,"SAon",opts->SAon,opts->ngQM,FALSE);
353 pr_ivec(fp,indent,"SAsteps",opts->SAsteps,opts->ngQM,FALSE);
354 pr_bvec(fp,indent,"bOPT",opts->bOPT,opts->ngQM,FALSE);
355 pr_bvec(fp,indent,"bTS",opts->bTS,opts->ngQM,FALSE);
359 static void pr_grp_opts(FILE *out,int indent,const char *title,t_grpopts *opts,
360 bool bMDPformat)
362 int i,m,j;
364 if (!bMDPformat)
365 fprintf(out,"%s:\n",title);
367 pr_indent(out,indent);
368 fprintf(out,"nrdf%s",bMDPformat ? " = " : ":");
369 for(i=0; (i<opts->ngtc); i++)
370 fprintf(out," %10g",opts->nrdf[i]);
371 fprintf(out,"\n");
373 pr_indent(out,indent);
374 fprintf(out,"ref_t%s",bMDPformat ? " = " : ":");
375 for(i=0; (i<opts->ngtc); i++)
376 fprintf(out," %10g",opts->ref_t[i]);
377 fprintf(out,"\n");
379 pr_indent(out,indent);
380 fprintf(out,"tau_t%s",bMDPformat ? " = " : ":");
381 for(i=0; (i<opts->ngtc); i++)
382 fprintf(out," %10g",opts->tau_t[i]);
383 fprintf(out,"\n");
385 /* Pretty-print the simulated annealing info */
386 fprintf(out,"anneal%s",bMDPformat ? " = " : ":");
387 for(i=0; (i<opts->ngtc); i++)
388 fprintf(out," %10s",EANNEAL(opts->annealing[i]));
389 fprintf(out,"\n");
391 fprintf(out,"ann_npoints%s",bMDPformat ? " = " : ":");
392 for(i=0; (i<opts->ngtc); i++)
393 fprintf(out," %10d",opts->anneal_npoints[i]);
394 fprintf(out,"\n");
396 for(i=0; (i<opts->ngtc); i++) {
397 if(opts->anneal_npoints[i]>0) {
398 fprintf(out,"ann. times [%d]:\t",i);
399 for(j=0; (j<opts->anneal_npoints[i]); j++)
400 fprintf(out," %10.1f",opts->anneal_time[i][j]);
401 fprintf(out,"\n");
402 fprintf(out,"ann. temps [%d]:\t",i);
403 for(j=0; (j<opts->anneal_npoints[i]); j++)
404 fprintf(out," %10.1f",opts->anneal_temp[i][j]);
405 fprintf(out,"\n");
409 pr_indent(out,indent);
410 fprintf(out,"acc:\t");
411 for(i=0; (i<opts->ngacc); i++)
412 for(m=0; (m<DIM); m++)
413 fprintf(out," %10g",opts->acc[i][m]);
414 fprintf(out,"\n");
416 pr_indent(out,indent);
417 fprintf(out,"nfreeze:");
418 for(i=0; (i<opts->ngfrz); i++)
419 for(m=0; (m<DIM); m++)
420 fprintf(out," %10s",opts->nFreeze[i][m] ? "Y" : "N");
421 fprintf(out,"\n");
424 for(i=0; (i<opts->ngener); i++) {
425 pr_indent(out,indent);
426 fprintf(out,"energygrp_flags[%3d]:",i);
427 for(m=0; (m<opts->ngener); m++)
428 fprintf(out," %d",opts->egp_flags[opts->ngener*i+m]);
429 fprintf(out,"\n");
432 fflush(out);
435 static void pr_matrix(FILE *fp,int indent,const char *title,rvec *m,
436 bool bMDPformat)
438 if (bMDPformat)
439 fprintf(fp,"%-10s = %g %g %g %g %g %g\n",title,
440 m[XX][XX],m[YY][YY],m[ZZ][ZZ],m[XX][YY],m[XX][ZZ],m[YY][ZZ]);
441 else
442 pr_rvecs(fp,indent,title,m,DIM);
445 static void pr_cosine(FILE *fp,int indent,const char *title,t_cosines *cos,
446 bool bMDPformat)
448 int j;
450 if (bMDPformat) {
451 fprintf(fp,"%s = %d\n",title,cos->n);
453 else {
454 indent=pr_title(fp,indent,title);
455 (void) pr_indent(fp,indent);
456 fprintf(fp,"n = %d\n",cos->n);
457 if (cos->n > 0) {
458 (void) pr_indent(fp,indent+2);
459 fprintf(fp,"a =");
460 for(j=0; (j<cos->n); j++)
461 fprintf(fp," %e",cos->a[j]);
462 fprintf(fp,"\n");
463 (void) pr_indent(fp,indent+2);
464 fprintf(fp,"phi =");
465 for(j=0; (j<cos->n); j++)
466 fprintf(fp," %e",cos->phi[j]);
467 fprintf(fp,"\n");
472 #define PS(t,s) pr_str(fp,indent,t,s)
473 #define PI(t,s) pr_int(fp,indent,t,s)
474 #define PSTEP(t,s) pr_gmx_large_int(fp,indent,t,s)
475 #define PR(t,s) pr_real(fp,indent,t,s)
477 static void pr_pullgrp(FILE *fp,int indent,int g,t_pullgrp *pg)
479 pr_indent(fp,indent);
480 fprintf(fp,"pull_group %d:\n",g);
481 indent += 2;
482 pr_ivec_block(fp,indent,"atom",pg->ind,pg->nat,TRUE);
483 pr_rvec(fp,indent,"weight",pg->weight,pg->nweight,TRUE);
484 PI("pbcatom",pg->pbcatom);
485 pr_rvec(fp,indent,"vec",pg->vec,DIM,TRUE);
486 pr_rvec(fp,indent,"init",pg->init,DIM,TRUE);
487 PR("rate",pg->rate);
488 PR("k",pg->k);
489 PR("kB",pg->kB);
492 static void pr_pull(FILE *fp,int indent,t_pull *pull)
494 int g;
496 PS("pull_geometry",EPULLGEOM(pull->eGeom));
497 pr_ivec(fp,indent,"pull_dim",pull->dim,DIM,TRUE);
498 PR("pull_r1",pull->cyl_r1);
499 PR("pull_r0",pull->cyl_r0);
500 PR("pull_constr_tol",pull->constr_tol);
501 PI("pull_nstxout",pull->nstxout);
502 PI("pull_nstfout",pull->nstfout);
503 PI("pull_ngrp",pull->ngrp);
504 for(g=0; g<pull->ngrp+1; g++)
505 pr_pullgrp(fp,indent,g,&pull->grp[g]);
508 void pr_inputrec(FILE *fp,int indent,const char *title,t_inputrec *ir,
509 bool bMDPformat)
511 const char *infbuf="inf";
512 int i;
514 if (available(fp,ir,indent,title)) {
515 if (!bMDPformat)
516 indent=pr_title(fp,indent,title);
517 PS("integrator",EI(ir->eI));
518 PSTEP("nsteps",ir->nsteps);
519 PSTEP("init_step",ir->init_step);
520 PI("nstcalcenergy",ir->nstcalcenergy);
521 PS("ns_type",ENS(ir->ns_type));
522 PI("nstlist",ir->nstlist);
523 PI("ndelta",ir->ndelta);
524 PI("nstcomm",ir->nstcomm);
525 PS("comm_mode",ECOM(ir->comm_mode));
526 PI("nstlog",ir->nstlog);
527 PI("nstxout",ir->nstxout);
528 PI("nstvout",ir->nstvout);
529 PI("nstfout",ir->nstfout);
530 PI("nstenergy",ir->nstenergy);
531 PI("nstxtcout",ir->nstxtcout);
532 PR("init_t",ir->init_t);
533 PR("delta_t",ir->delta_t);
535 PR("xtcprec",ir->xtcprec);
536 PI("nkx",ir->nkx);
537 PI("nky",ir->nky);
538 PI("nkz",ir->nkz);
539 PI("pme_order",ir->pme_order);
540 PR("ewald_rtol",ir->ewald_rtol);
541 PR("ewald_geometry",ir->ewald_geometry);
542 PR("epsilon_surface",ir->epsilon_surface);
543 PS("optimize_fft",BOOL(ir->bOptFFT));
544 PS("ePBC",EPBC(ir->ePBC));
545 PS("bPeriodicMols",BOOL(ir->bPeriodicMols));
546 PS("bContinuation",BOOL(ir->bContinuation));
547 PS("bShakeSOR",BOOL(ir->bShakeSOR));
548 PS("etc",ETCOUPLTYPE(ir->etc));
549 PS("epc",EPCOUPLTYPE(ir->epc));
550 PS("epctype",EPCOUPLTYPETYPE(ir->epct));
551 PR("tau_p",ir->tau_p);
552 pr_matrix(fp,indent,"ref_p",ir->ref_p,bMDPformat);
553 pr_matrix(fp,indent,"compress",ir->compress,bMDPformat);
554 PS("refcoord_scaling",EREFSCALINGTYPE(ir->refcoord_scaling));
555 if (bMDPformat)
556 fprintf(fp,"posres_com = %g %g %g\n",ir->posres_com[XX],
557 ir->posres_com[YY],ir->posres_com[ZZ]);
558 else
559 pr_rvec(fp,indent,"posres_com",ir->posres_com,DIM,TRUE);
560 if (bMDPformat)
561 fprintf(fp,"posres_comB = %g %g %g\n",ir->posres_comB[XX],
562 ir->posres_comB[YY],ir->posres_comB[ZZ]);
563 else
564 pr_rvec(fp,indent,"posres_comB",ir->posres_comB,DIM,TRUE);
565 PI("andersen_seed",ir->andersen_seed);
566 PR("rlist",ir->rlist);
567 PR("rlistlong",ir->rlistlong);
568 PR("rtpi",ir->rtpi);
569 PS("coulombtype",EELTYPE(ir->coulombtype));
570 PR("rcoulomb_switch",ir->rcoulomb_switch);
571 PR("rcoulomb",ir->rcoulomb);
572 PS("vdwtype",EVDWTYPE(ir->vdwtype));
573 PR("rvdw_switch",ir->rvdw_switch);
574 PR("rvdw",ir->rvdw);
575 if (ir->epsilon_r != 0)
576 PR("epsilon_r",ir->epsilon_r);
577 else
578 PS("epsilon_r",infbuf);
579 if (ir->epsilon_rf != 0)
580 PR("epsilon_rf",ir->epsilon_rf);
581 else
582 PS("epsilon_rf",infbuf);
583 PR("tabext",ir->tabext);
584 PS("implicit_solvent",EIMPLICITSOL(ir->implicit_solvent));
585 PS("gb_algorithm",EGBALGORITHM(ir->gb_algorithm));
586 PR("gb_epsilon_solvent",ir->gb_epsilon_solvent);
587 PI("nstgbradii",ir->nstgbradii);
588 PR("rgbradii",ir->rgbradii);
589 PR("gb_saltconc",ir->gb_saltconc);
590 PR("gb_obc_alpha",ir->gb_obc_alpha);
591 PR("gb_obc_beta",ir->gb_obc_beta);
592 PR("gb_obc_gamma",ir->gb_obc_gamma);
593 PR("gb_dielectric_offset",ir->gb_dielectric_offset);
594 PS("sa_algorithm",ESAALGORITHM(ir->gb_algorithm));
595 PR("sa_surface_tension",ir->sa_surface_tension);
597 PS("DispCorr",EDISPCORR(ir->eDispCorr));
598 PS("free_energy",EFEPTYPE(ir->efep));
599 PR("init_lambda",ir->init_lambda);
600 PR("delta_lambda",ir->delta_lambda);
601 if (!bMDPformat)
603 PI("n_foreign_lambda",ir->n_flambda);
605 if (ir->n_flambda > 0)
607 pr_indent(fp,indent);
608 fprintf(fp,"foreign_lambda%s",bMDPformat ? " = " : ":");
609 for(i=0; i<ir->n_flambda; i++)
611 fprintf(fp," %10g",ir->flambda[i]);
613 fprintf(fp,"\n");
615 PR("sc_alpha",ir->sc_alpha);
616 PI("sc_power",ir->sc_power);
617 PR("sc_sigma",ir->sc_sigma);
618 PI("nstdhdl", ir->nstdhdl);
620 PI("nwall",ir->nwall);
621 PS("wall_type",EWALLTYPE(ir->wall_type));
622 PI("wall_atomtype[0]",ir->wall_atomtype[0]);
623 PI("wall_atomtype[1]",ir->wall_atomtype[1]);
624 PR("wall_density[0]",ir->wall_density[0]);
625 PR("wall_density[1]",ir->wall_density[1]);
626 PR("wall_ewald_zfac",ir->wall_ewald_zfac);
628 PS("pull",EPULLTYPE(ir->ePull));
629 if (ir->ePull != epullNO)
630 pr_pull(fp,indent,ir->pull);
632 PS("disre",EDISRETYPE(ir->eDisre));
633 PS("disre_weighting",EDISREWEIGHTING(ir->eDisreWeighting));
634 PS("disre_mixed",BOOL(ir->bDisreMixed));
635 PR("dr_fc",ir->dr_fc);
636 PR("dr_tau",ir->dr_tau);
637 PR("nstdisreout",ir->nstdisreout);
638 PR("orires_fc",ir->orires_fc);
639 PR("orires_tau",ir->orires_tau);
640 PR("nstorireout",ir->nstorireout);
642 PR("dihre-fc",ir->dihre_fc);
644 PR("em_stepsize",ir->em_stepsize);
645 PR("em_tol",ir->em_tol);
646 PI("niter",ir->niter);
647 PR("fc_stepsize",ir->fc_stepsize);
648 PI("nstcgsteep",ir->nstcgsteep);
649 PI("nbfgscorr",ir->nbfgscorr);
651 PS("ConstAlg",ECONSTRTYPE(ir->eConstrAlg));
652 PR("shake_tol",ir->shake_tol);
653 PI("lincs_order",ir->nProjOrder);
654 PR("lincs_warnangle",ir->LincsWarnAngle);
655 PI("lincs_iter",ir->nLincsIter);
656 PR("bd_fric",ir->bd_fric);
657 PI("ld_seed",ir->ld_seed);
658 PR("cos_accel",ir->cos_accel);
659 pr_matrix(fp,indent,"deform",ir->deform,bMDPformat);
660 PI("userint1",ir->userint1);
661 PI("userint2",ir->userint2);
662 PI("userint3",ir->userint3);
663 PI("userint4",ir->userint4);
664 PR("userreal1",ir->userreal1);
665 PR("userreal2",ir->userreal2);
666 PR("userreal3",ir->userreal3);
667 PR("userreal4",ir->userreal4);
668 pr_grp_opts(fp,indent,"grpopts",&(ir->opts),bMDPformat);
669 pr_cosine(fp,indent,"efield-x",&(ir->ex[XX]),bMDPformat);
670 pr_cosine(fp,indent,"efield-xt",&(ir->et[XX]),bMDPformat);
671 pr_cosine(fp,indent,"efield-y",&(ir->ex[YY]),bMDPformat);
672 pr_cosine(fp,indent,"efield-yt",&(ir->et[YY]),bMDPformat);
673 pr_cosine(fp,indent,"efield-z",&(ir->ex[ZZ]),bMDPformat);
674 pr_cosine(fp,indent,"efield-zt",&(ir->et[ZZ]),bMDPformat);
675 PS("bQMMM",BOOL(ir->bQMMM));
676 PI("QMconstraints",ir->QMconstraints);
677 PI("QMMMscheme",ir->QMMMscheme);
678 PR("scalefactor",ir->scalefactor);
679 pr_qm_opts(fp,indent,"qm_opts",&(ir->opts));
682 #undef PS
683 #undef PR
684 #undef PI
686 static void pr_harm(FILE *fp,t_iparams *iparams,const char *r,const char *kr)
688 fprintf(fp,"%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e\n",
689 r,iparams->harmonic.rA,kr,iparams->harmonic.krA,
690 r,iparams->harmonic.rB,kr,iparams->harmonic.krB);
693 void pr_iparams(FILE *fp,t_functype ftype,t_iparams *iparams)
695 int i;
696 real VA[4],VB[4],*rbcA,*rbcB;
698 switch (ftype) {
699 case F_ANGLES:
700 case F_G96ANGLES:
701 pr_harm(fp,iparams,"th","ct");
702 break;
703 case F_CROSS_BOND_BONDS:
704 fprintf(fp,"r1e=%15.8e, r2e=%15.8e, krr=%15.8e\n",
705 iparams->cross_bb.r1e,iparams->cross_bb.r2e,
706 iparams->cross_bb.krr);
707 break;
708 case F_CROSS_BOND_ANGLES:
709 fprintf(fp,"r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e\n",
710 iparams->cross_ba.r1e,iparams->cross_ba.r2e,
711 iparams->cross_ba.r3e,iparams->cross_ba.krt);
712 break;
713 case F_UREY_BRADLEY:
714 fprintf(fp,"theta=%15.8e, ktheta=%15.8e, r13=%15.8e, kUB=%15.8e\n",
715 iparams->u_b.theta,iparams->u_b.ktheta,iparams->u_b.r13,iparams->u_b.kUB);
716 break;
717 case F_QUARTIC_ANGLES:
718 fprintf(fp,"theta=%15.8e",iparams->qangle.theta);
719 for(i=0; i<5; i++)
720 fprintf(fp,", c%c=%15.8e",'0'+i,iparams->qangle.c[i]);
721 fprintf(fp,"\n");
722 break;
723 case F_BHAM:
724 fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
725 iparams->bham.a,iparams->bham.b,iparams->bham.c);
726 break;
727 case F_BONDS:
728 case F_G96BONDS:
729 case F_HARMONIC:
730 pr_harm(fp,iparams,"b0","cb");
731 break;
732 case F_IDIHS:
733 pr_harm(fp,iparams,"xi","cx");
734 break;
735 case F_MORSE:
736 fprintf(fp,"b0=%15.8e, cb=%15.8e, beta=%15.8e\n",
737 iparams->morse.b0,iparams->morse.cb,iparams->morse.beta);
738 break;
739 case F_CUBICBONDS:
740 fprintf(fp,"b0=%15.8e, kb=%15.8e, kcub=%15.8e\n",
741 iparams->cubic.b0,iparams->cubic.kb,iparams->cubic.kcub);
742 break;
743 case F_CONNBONDS:
744 fprintf(fp,"\n");
745 break;
746 case F_FENEBONDS:
747 fprintf(fp,"bm=%15.8e, kb=%15.8e\n",iparams->fene.bm,iparams->fene.kb);
748 break;
749 case F_TABBONDS:
750 case F_TABBONDSNC:
751 case F_TABANGLES:
752 case F_TABDIHS:
753 fprintf(fp,"tab=%d, kA=%15.8e, kB=%15.8e\n",
754 iparams->tab.table,iparams->tab.kA,iparams->tab.kB);
755 break;
756 case F_POLARIZATION:
757 fprintf(fp,"alpha=%15.8e\n",iparams->polarize.alpha);
758 break;
759 case F_THOLE_POL:
760 fprintf(fp,"a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e\n",
761 iparams->thole.a,iparams->thole.alpha1,iparams->thole.alpha2,
762 iparams->thole.rfac);
763 break;
764 case F_WATER_POL:
765 fprintf(fp,"al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f\n",
766 iparams->wpol.al_x,iparams->wpol.al_y,iparams->wpol.al_z,
767 iparams->wpol.rOH,iparams->wpol.rHH,iparams->wpol.rOD);
768 break;
769 case F_LJ:
770 fprintf(fp,"c6=%15.8e, c12=%15.8e\n",iparams->lj.c6,iparams->lj.c12);
771 break;
772 case F_LJ14:
773 fprintf(fp,"c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e\n",
774 iparams->lj14.c6A,iparams->lj14.c12A,
775 iparams->lj14.c6B,iparams->lj14.c12B);
776 break;
777 case F_LJC14_Q:
778 fprintf(fp,"fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
779 iparams->ljc14.fqq,
780 iparams->ljc14.qi,iparams->ljc14.qj,
781 iparams->ljc14.c6,iparams->ljc14.c12);
782 break;
783 case F_LJC_PAIRS_NB:
784 fprintf(fp,"qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e\n",
785 iparams->ljcnb.qi,iparams->ljcnb.qj,
786 iparams->ljcnb.c6,iparams->ljcnb.c12);
787 break;
788 case F_PDIHS:
789 case F_ANGRES:
790 case F_ANGRESZ:
791 fprintf(fp,"phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d\n",
792 iparams->pdihs.phiA,iparams->pdihs.cpA,
793 iparams->pdihs.phiB,iparams->pdihs.cpB,
794 iparams->pdihs.mult);
795 break;
796 case F_DISRES:
797 fprintf(fp,"label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)\n",
798 iparams->disres.label,iparams->disres.type,
799 iparams->disres.low,iparams->disres.up1,
800 iparams->disres.up2,iparams->disres.kfac);
801 break;
802 case F_ORIRES:
803 fprintf(fp,"ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)\n",
804 iparams->orires.ex,iparams->orires.label,iparams->orires.power,
805 iparams->orires.c,iparams->orires.obs,iparams->orires.kfac);
806 break;
807 case F_DIHRES:
808 fprintf(fp,"label=%d, power=%4d phi=%15.8e, dphi=%15.8e, kfac=%15.8e)\n",
809 iparams->dihres.label,iparams->dihres.power,
810 iparams->dihres.phi,iparams->dihres.dphi,iparams->dihres.kfac);
811 break;
812 case F_POSRES:
813 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",
814 iparams->posres.pos0A[XX],iparams->posres.pos0A[YY],
815 iparams->posres.pos0A[ZZ],iparams->posres.fcA[XX],
816 iparams->posres.fcA[YY],iparams->posres.fcA[ZZ],
817 iparams->posres.pos0B[XX],iparams->posres.pos0B[YY],
818 iparams->posres.pos0B[ZZ],iparams->posres.fcB[XX],
819 iparams->posres.fcB[YY],iparams->posres.fcB[ZZ]);
820 break;
821 case F_RBDIHS:
822 for (i=0; i<NR_RBDIHS; i++)
823 fprintf(fp,"%srbcA[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcA[i]);
824 fprintf(fp,"\n");
825 for (i=0; i<NR_RBDIHS; i++)
826 fprintf(fp,"%srbcB[%d]=%15.8e",i==0?"":", ",i,iparams->rbdihs.rbcB[i]);
827 fprintf(fp,"\n");
828 break;
829 case F_FOURDIHS:
830 /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get the
831 * OPLS potential constants back.
833 rbcA = iparams->rbdihs.rbcA;
834 rbcB = iparams->rbdihs.rbcB;
836 VA[3] = -0.25*rbcA[4];
837 VA[2] = -0.5*rbcA[3];
838 VA[1] = 4.0*VA[3]-rbcA[2];
839 VA[0] = 3.0*VA[2]-2.0*rbcA[1];
841 VB[3] = -0.25*rbcB[4];
842 VB[2] = -0.5*rbcB[3];
843 VB[1] = 4.0*VB[3]-rbcB[2];
844 VB[0] = 3.0*VB[2]-2.0*rbcB[1];
846 for (i=0; i<NR_FOURDIHS; i++)
847 fprintf(fp,"%sFourA[%d]=%15.8e",i==0?"":", ",i,VA[i]);
848 fprintf(fp,"\n");
849 for (i=0; i<NR_FOURDIHS; i++)
850 fprintf(fp,"%sFourB[%d]=%15.8e",i==0?"":", ",i,VB[i]);
851 fprintf(fp,"\n");
852 break;
854 case F_CONSTR:
855 case F_CONSTRNC:
856 fprintf(fp,"dA=%15.8e, dB=%15.8e\n",iparams->constr.dA,iparams->constr.dB);
857 break;
858 case F_SETTLE:
859 fprintf(fp,"doh=%15.8e, dhh=%15.8e\n",iparams->settle.doh,
860 iparams->settle.dhh);
861 break;
862 case F_VSITE2:
863 fprintf(fp,"a=%15.8e\n",iparams->vsite.a);
864 break;
865 case F_VSITE3:
866 case F_VSITE3FD:
867 case F_VSITE3FAD:
868 fprintf(fp,"a=%15.8e, b=%15.8e\n",iparams->vsite.a,iparams->vsite.b);
869 break;
870 case F_VSITE3OUT:
871 case F_VSITE4FD:
872 case F_VSITE4FDN:
873 fprintf(fp,"a=%15.8e, b=%15.8e, c=%15.8e\n",
874 iparams->vsite.a,iparams->vsite.b,iparams->vsite.c);
875 break;
876 case F_VSITEN:
877 fprintf(fp,"n=%2d, a=%15.8e\n",iparams->vsiten.n,iparams->vsiten.a);
878 break;
879 case F_GB12:
880 case F_GB13:
881 case F_GB14:
882 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);
883 break;
884 case F_CMAP:
885 fprintf(fp, "cmapA=%1d, cmapB=%1d\n",iparams->cmap.cmapA, iparams->cmap.cmapB);
886 break;
887 default:
888 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
889 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
893 void pr_ilist(FILE *fp,int indent,const char *title,
894 t_functype *functype,t_ilist *ilist, bool bShowNumbers)
896 int i,j,k,type,ftype;
897 t_iatom *iatoms;
899 if (available(fp,ilist,indent,title) && ilist->nr > 0)
901 indent=pr_title(fp,indent,title);
902 (void) pr_indent(fp,indent);
903 fprintf(fp,"nr: %d\n",ilist->nr);
904 if (ilist->nr > 0) {
905 (void) pr_indent(fp,indent);
906 fprintf(fp,"iatoms:\n");
907 iatoms=ilist->iatoms;
908 for (i=j=0; i<ilist->nr;) {
909 #ifndef DEBUG
910 (void) pr_indent(fp,indent+INDENT);
911 type=*(iatoms++);
912 ftype=functype[type];
913 (void) fprintf(fp,"%d type=%d (%s)",
914 bShowNumbers?j:-1,bShowNumbers?type:-1,
915 interaction_function[ftype].name);
916 j++;
917 for (k=0; k<interaction_function[ftype].nratoms; k++)
918 (void) fprintf(fp," %u",*(iatoms++));
919 (void) fprintf(fp,"\n");
920 i+=1+interaction_function[ftype].nratoms;
921 #else
922 fprintf(fp,"%5d%5d\n",i,iatoms[i]);
923 i++;
924 #endif
930 void pr_ffparams(FILE *fp,int indent,const char *title,
931 gmx_ffparams_t *ffparams,
932 bool bShowNumbers)
934 int i,j;
936 indent=pr_title(fp,indent,title);
937 (void) pr_indent(fp,indent);
938 (void) fprintf(fp,"atnr=%d\n",ffparams->atnr);
939 (void) pr_indent(fp,indent);
940 (void) fprintf(fp,"ntypes=%d\n",ffparams->ntypes);
941 for (i=0; i<ffparams->ntypes; i++) {
942 (void) pr_indent(fp,indent+INDENT);
943 (void) fprintf(fp,"functype[%d]=%s, ",
944 bShowNumbers?i:-1,
945 interaction_function[ffparams->functype[i]].name);
946 pr_iparams(fp,ffparams->functype[i],&ffparams->iparams[i]);
948 (void) pr_double(fp,indent,"reppow",ffparams->reppow);
949 (void) pr_real(fp,indent,"fudgeQQ",ffparams->fudgeQQ);
952 void pr_idef(FILE *fp,int indent,const char *title,t_idef *idef, bool bShowNumbers)
954 int i,j;
956 if (available(fp,idef,indent,title)) {
957 indent=pr_title(fp,indent,title);
958 (void) pr_indent(fp,indent);
959 (void) fprintf(fp,"atnr=%d\n",idef->atnr);
960 (void) pr_indent(fp,indent);
961 (void) fprintf(fp,"ntypes=%d\n",idef->ntypes);
962 for (i=0; i<idef->ntypes; i++) {
963 (void) pr_indent(fp,indent+INDENT);
964 (void) fprintf(fp,"functype[%d]=%s, ",
965 bShowNumbers?i:-1,
966 interaction_function[idef->functype[i]].name);
967 pr_iparams(fp,idef->functype[i],&idef->iparams[i]);
969 (void) pr_real(fp,indent,"fudgeQQ",idef->fudgeQQ);
971 for(j=0; (j<F_NRE); j++)
972 pr_ilist(fp,indent,interaction_function[j].longname,
973 idef->functype,&idef->il[j],bShowNumbers);
977 static int pr_block_title(FILE *fp,int indent,const char *title,t_block *block)
979 int i;
981 if (available(fp,block,indent,title))
983 indent=pr_title(fp,indent,title);
984 (void) pr_indent(fp,indent);
985 (void) fprintf(fp,"nr=%d\n",block->nr);
987 return indent;
990 static int pr_blocka_title(FILE *fp,int indent,const char *title,t_blocka *block)
992 int i;
994 if (available(fp,block,indent,title))
996 indent=pr_title(fp,indent,title);
997 (void) pr_indent(fp,indent);
998 (void) fprintf(fp,"nr=%d\n",block->nr);
999 (void) pr_indent(fp,indent);
1000 (void) fprintf(fp,"nra=%d\n",block->nra);
1002 return indent;
1005 static void low_pr_block(FILE *fp,int indent,const char *title,t_block *block, bool bShowNumbers)
1007 int i;
1009 if (available(fp,block,indent,title))
1011 indent=pr_block_title(fp,indent,title,block);
1012 for (i=0; i<=block->nr; i++)
1014 (void) pr_indent(fp,indent+INDENT);
1015 (void) fprintf(fp,"%s->index[%d]=%u\n",
1016 title,bShowNumbers?i:-1,block->index[i]);
1021 static void low_pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block, bool bShowNumbers)
1023 int i;
1025 if (available(fp,block,indent,title))
1027 indent=pr_blocka_title(fp,indent,title,block);
1028 for (i=0; i<=block->nr; i++)
1030 (void) pr_indent(fp,indent+INDENT);
1031 (void) fprintf(fp,"%s->index[%d]=%u\n",
1032 title,bShowNumbers?i:-1,block->index[i]);
1034 for (i=0; i<block->nra; i++)
1036 (void) pr_indent(fp,indent+INDENT);
1037 (void) fprintf(fp,"%s->a[%d]=%u\n",
1038 title,bShowNumbers?i:-1,block->a[i]);
1043 void pr_block(FILE *fp,int indent,const char *title,t_block *block,bool bShowNumbers)
1045 int i,j,ok,size,start,end;
1047 if (available(fp,block,indent,title))
1049 indent=pr_block_title(fp,indent,title,block);
1050 start=0;
1051 end=start;
1052 if ((ok=(block->index[start]==0))==0)
1053 (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1054 else
1055 for (i=0; i<block->nr; i++)
1057 end=block->index[i+1];
1058 size=pr_indent(fp,indent);
1059 if (end<=start)
1060 size+=fprintf(fp,"%s[%d]={}\n",title,i);
1061 else
1062 size+=fprintf(fp,"%s[%d]={%d..%d}\n",
1063 title,bShowNumbers?i:-1,
1064 bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1065 start=end;
1070 void pr_blocka(FILE *fp,int indent,const char *title,t_blocka *block,bool bShowNumbers)
1072 int i,j,ok,size,start,end;
1074 if (available(fp,block,indent,title))
1076 indent=pr_blocka_title(fp,indent,title,block);
1077 start=0;
1078 end=start;
1079 if ((ok=(block->index[start]==0))==0)
1080 (void) fprintf(fp,"block->index[%d] should be 0\n",start);
1081 else
1082 for (i=0; i<block->nr; i++)
1084 end=block->index[i+1];
1085 size=pr_indent(fp,indent);
1086 if (end<=start)
1087 size+=fprintf(fp,"%s[%d]={",title,i);
1088 else
1089 size+=fprintf(fp,"%s[%d][%d..%d]={",
1090 title,bShowNumbers?i:-1,
1091 bShowNumbers?start:-1,bShowNumbers?end-1:-1);
1092 for (j=start; j<end; j++)
1094 if (j>start) size+=fprintf(fp,", ");
1095 if ((size)>(USE_WIDTH))
1097 (void) fprintf(fp,"\n");
1098 size=pr_indent(fp,indent+INDENT);
1100 size+=fprintf(fp,"%u",block->a[j]);
1102 (void) fprintf(fp,"}\n");
1103 start=end;
1105 if ((end!=block->nra)||(!ok))
1107 (void) pr_indent(fp,indent);
1108 (void) fprintf(fp,"tables inconsistent, dumping complete tables:\n");
1109 low_pr_blocka(fp,indent,title,block,bShowNumbers);
1114 static void pr_strings(FILE *fp,int indent,const char *title,char ***nm,int n, bool bShowNumbers)
1116 int i;
1118 if (available(fp,nm,indent,title))
1120 indent=pr_title_n(fp,indent,title,n);
1121 for (i=0; i<n; i++)
1123 (void) pr_indent(fp,indent);
1124 (void) fprintf(fp,"%s[%d]={name=\"%s\"}\n",
1125 title,bShowNumbers?i:-1,*(nm[i]));
1130 static void pr_strings2(FILE *fp,int indent,const char *title,
1131 char ***nm,char ***nmB,int n, bool bShowNumbers)
1133 int i;
1135 if (available(fp,nm,indent,title))
1137 indent=pr_title_n(fp,indent,title,n);
1138 for (i=0; i<n; i++)
1140 (void) pr_indent(fp,indent);
1141 (void) fprintf(fp,"%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
1142 title,bShowNumbers?i:-1,*(nm[i]),*(nmB[i]));
1147 static void pr_resinfo(FILE *fp,int indent,const char *title,t_resinfo *resinfo,int n, bool bShowNumbers)
1149 int i;
1151 if (available(fp,resinfo,indent,title))
1153 indent=pr_title_n(fp,indent,title,n);
1154 for (i=0; i<n; i++)
1156 (void) pr_indent(fp,indent);
1157 (void) fprintf(fp,"%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
1158 title,bShowNumbers?i:-1,
1159 *(resinfo[i].name),resinfo[i].nr,
1160 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
1165 static void pr_atom(FILE *fp,int indent,const char *title,t_atom *atom,int n)
1167 int i,j;
1169 if (available(fp,atom,indent,title)) {
1170 indent=pr_title_n(fp,indent,title,n);
1171 for (i=0; i<n; i++) {
1172 (void) pr_indent(fp,indent);
1173 fprintf(fp,"%s[%6d]={type=%3d, typeB=%3d, ptype=%8s, m=%12.5e, "
1174 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
1175 title,i,atom[i].type,atom[i].typeB,ptype_str[atom[i].ptype],
1176 atom[i].m,atom[i].q,atom[i].mB,atom[i].qB,
1177 atom[i].resind,atom[i].atomnumber);
1182 static void pr_grps(FILE *fp,int indent,const char *title,t_grps grps[],int ngrp,
1183 char **grpname[], bool bShowNumbers)
1185 int i,j;
1187 for(i=0; (i<ngrp); i++) {
1188 fprintf(fp,"%s[%d] nr=%d, name=[",title,bShowNumbers?i:-1,grps[i].nr);
1189 for(j=0; (j<grps[i].nr); j++)
1190 fprintf(fp," %s",*(grpname[grps[i].nm_ind[j]]));
1191 fprintf(fp,"]\n");
1195 static void pr_groups(FILE *fp,int indent,const char *title,
1196 gmx_groups_t *groups,
1197 bool bShowNumbers)
1199 int grpnr[egcNR];
1200 int nat_max,i,g;
1202 pr_grps(fp,indent,"grp",groups->grps,egcNR,groups->grpname,bShowNumbers);
1203 pr_strings(fp,indent,"grpname",groups->grpname,groups->ngrpname,bShowNumbers);
1205 (void) pr_indent(fp,indent);
1206 fprintf(fp,"groups ");
1207 for(g=0; g<egcNR; g++)
1209 printf(" %5.5s",gtypes[g]);
1211 printf("\n");
1213 (void) pr_indent(fp,indent);
1214 fprintf(fp,"allocated ");
1215 nat_max = 0;
1216 for(g=0; g<egcNR; g++)
1218 printf(" %5d",groups->ngrpnr[g]);
1219 nat_max = max(nat_max,groups->ngrpnr[g]);
1221 printf("\n");
1223 if (nat_max == 0)
1225 (void) pr_indent(fp,indent);
1226 fprintf(fp,"groupnr[%5s] =","*");
1227 for(g=0; g<egcNR; g++)
1229 fprintf(fp," %3d ",0);
1231 fprintf(fp,"\n");
1233 else
1235 for(i=0; i<nat_max; i++)
1237 (void) pr_indent(fp,indent);
1238 fprintf(fp,"groupnr[%5d] =",i);
1239 for(g=0; g<egcNR; g++)
1241 fprintf(fp," %3d ",
1242 groups->grpnr[g] ? groups->grpnr[g][i] : 0);
1244 fprintf(fp,"\n");
1249 void pr_atoms(FILE *fp,int indent,const char *title,t_atoms *atoms,
1250 bool bShownumbers)
1252 if (available(fp,atoms,indent,title))
1254 indent=pr_title(fp,indent,title);
1255 pr_atom(fp,indent,"atom",atoms->atom,atoms->nr);
1256 pr_strings(fp,indent,"atom",atoms->atomname,atoms->nr,bShownumbers);
1257 pr_strings2(fp,indent,"type",atoms->atomtype,atoms->atomtypeB,atoms->nr,bShownumbers);
1258 pr_resinfo(fp,indent,"residue",atoms->resinfo,atoms->nres,bShownumbers);
1262 void pr_cmap(FILE *fp, int indent, const char *title, gmx_cmap_t *cmap_grid,
1263 bool bShowNumbers)
1265 int i,j,nelem;
1266 real dx,idx;
1268 dx = 360.0 / cmap_grid->grid_spacing;
1269 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1271 if(available(fp,cmap_grid,indent,title))
1273 fprintf(fp,"%s\n",title);
1275 for(i=0;i<cmap_grid->ngrid;i++)
1277 idx = -180.0;
1278 fprintf(fp,"%8s %8s %8s %8s\n","V","dVdx","dVdy","d2dV");
1280 fprintf(fp,"grid[%3d]={\n",bShowNumbers?i:-1);
1282 for(j=0;j<nelem;j++)
1284 if( (j%cmap_grid->grid_spacing)==0)
1286 fprintf(fp,"%8.1f\n",idx);
1287 idx+=dx;
1290 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4]);
1291 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+1]);
1292 fprintf(fp,"%8.3f ",cmap_grid->cmapdata[i].cmap[j*4+2]);
1293 fprintf(fp,"%8.3f\n",cmap_grid->cmapdata[i].cmap[j*4+3]);
1295 fprintf(fp,"\n");
1302 void pr_atomtypes(FILE *fp,int indent,const char *title,t_atomtypes *atomtypes,
1303 bool bShowNumbers)
1305 int i;
1306 if (available(fp,atomtypes,indent,title))
1308 indent=pr_title(fp,indent,title);
1309 for(i=0;i<atomtypes->nr;i++) {
1310 pr_indent(fp,indent);
1311 fprintf(fp,
1312 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
1313 bShowNumbers?i:-1,atomtypes->radius[i],atomtypes->vol[i],
1314 atomtypes->gb_radius[i],
1315 atomtypes->surftens[i],atomtypes->atomnumber[i],atomtypes->S_hct[i]);
1320 static void pr_moltype(FILE *fp,int indent,const char *title,
1321 gmx_moltype_t *molt,int n,
1322 gmx_ffparams_t *ffparams,
1323 bool bShowNumbers)
1325 int j;
1327 indent = pr_title_n(fp,indent,title,n);
1328 (void) pr_indent(fp,indent);
1329 (void) fprintf(fp,"name=\"%s\"\n",*(molt->name));
1330 pr_atoms(fp,indent,"atoms",&(molt->atoms),bShowNumbers);
1331 pr_block(fp,indent,"cgs",&molt->cgs, bShowNumbers);
1332 pr_blocka(fp,indent,"excls",&molt->excls, bShowNumbers);
1333 for(j=0; (j<F_NRE); j++) {
1334 pr_ilist(fp,indent,interaction_function[j].longname,
1335 ffparams->functype,&molt->ilist[j],bShowNumbers);
1339 static void pr_molblock(FILE *fp,int indent,const char *title,
1340 gmx_molblock_t *molb,int n,
1341 gmx_moltype_t *molt,
1342 bool bShowNumbers)
1344 indent = pr_title_n(fp,indent,title,n);
1345 (void) pr_indent(fp,indent);
1346 (void) fprintf(fp,"%-20s = %d \"%s\"\n",
1347 "moltype",molb->type,*(molt[molb->type].name));
1348 pr_int(fp,indent,"#molecules",molb->nmol);
1349 pr_int(fp,indent,"#atoms_mol",molb->natoms_mol);
1350 pr_int(fp,indent,"#posres_xA",molb->nposres_xA);
1351 if (molb->nposres_xA > 0) {
1352 pr_rvecs(fp,indent,"posres_xA",molb->posres_xA,molb->nposres_xA);
1354 pr_int(fp,indent,"#posres_xB",molb->nposres_xB);
1355 if (molb->nposres_xB > 0) {
1356 pr_rvecs(fp,indent,"posres_xB",molb->posres_xB,molb->nposres_xB);
1360 void pr_mtop(FILE *fp,int indent,const char *title,gmx_mtop_t *mtop,
1361 bool bShowNumbers)
1363 int mt,mb;
1365 if (available(fp,mtop,indent,title)) {
1366 indent=pr_title(fp,indent,title);
1367 (void) pr_indent(fp,indent);
1368 (void) fprintf(fp,"name=\"%s\"\n",*(mtop->name));
1369 pr_int(fp,indent,"#atoms",mtop->natoms);
1370 for(mb=0; mb<mtop->nmolblock; mb++) {
1371 pr_molblock(fp,indent,"molblock",&mtop->molblock[mb],mb,
1372 mtop->moltype,bShowNumbers);
1374 pr_ffparams(fp,indent,"ffparams",&(mtop->ffparams),bShowNumbers);
1375 pr_atomtypes(fp,indent,"atomtypes",&(mtop->atomtypes),bShowNumbers);
1376 for(mt=0; mt<mtop->nmoltype; mt++) {
1377 pr_moltype(fp,indent,"moltype",&mtop->moltype[mt],mt,
1378 &mtop->ffparams,bShowNumbers);
1380 pr_groups(fp,indent,"groups",&mtop->groups,bShowNumbers);
1381 pr_cmap(fp,indent,"cmap",&mtop->cmap_grid,bShowNumbers);
1385 void pr_top(FILE *fp,int indent,const char *title,t_topology *top, bool bShowNumbers)
1387 if (available(fp,top,indent,title)) {
1388 indent=pr_title(fp,indent,title);
1389 (void) pr_indent(fp,indent);
1390 (void) fprintf(fp,"name=\"%s\"\n",*(top->name));
1391 pr_atoms(fp,indent,"atoms",&(top->atoms),bShowNumbers);
1392 pr_atomtypes(fp,indent,"atomtypes",&(top->atomtypes),bShowNumbers);
1393 pr_block(fp,indent,"cgs",&top->cgs, bShowNumbers);
1394 pr_block(fp,indent,"mols",&top->mols, bShowNumbers);
1395 pr_blocka(fp,indent,"excls",&top->excls, bShowNumbers);
1396 pr_idef(fp,indent,"idef",&top->idef,bShowNumbers);
1400 void pr_header(FILE *fp,int indent,const char *title,t_tpxheader *sh)
1402 char buf[22];
1404 if (available(fp,sh,indent,title))
1406 indent=pr_title(fp,indent,title);
1407 pr_indent(fp,indent);
1408 fprintf(fp,"bIr = %spresent\n",sh->bIr?"":"not ");
1409 pr_indent(fp,indent);
1410 fprintf(fp,"bBox = %spresent\n",sh->bBox?"":"not ");
1411 pr_indent(fp,indent);
1412 fprintf(fp,"bTop = %spresent\n",sh->bTop?"":"not ");
1413 pr_indent(fp,indent);
1414 fprintf(fp,"bX = %spresent\n",sh->bX?"":"not ");
1415 pr_indent(fp,indent);
1416 fprintf(fp,"bV = %spresent\n",sh->bV?"":"not ");
1417 pr_indent(fp,indent);
1418 fprintf(fp,"bF = %spresent\n",sh->bF?"":"not ");
1420 pr_indent(fp,indent);
1421 fprintf(fp,"natoms = %d\n",sh->natoms);
1422 pr_indent(fp,indent);
1423 fprintf(fp,"lambda = %e\n",sh->lambda);
1427 void pr_commrec(FILE *fp,int indent,t_commrec *cr)
1429 pr_indent(fp,indent);
1430 fprintf(fp,"commrec:\n");
1431 indent+=2;
1432 pr_indent(fp,indent);
1433 fprintf(fp,"nodeid = %d\n",cr->nodeid);
1434 pr_indent(fp,indent);
1435 fprintf(fp,"nnodes = %d\n",cr->nnodes);
1436 pr_indent(fp,indent);
1437 fprintf(fp,"npmenodes = %d\n",cr->npmenodes);
1438 pr_indent(fp,indent);
1439 fprintf(fp,"threadid = %d\n",cr->threadid);
1440 pr_indent(fp,indent);
1441 fprintf(fp,"nthreads = %d\n",cr->nthreads);