Avoid writing xvgr formatting with -xvg none
authorErik Lindahl <erik@kth.se>
Fri, 6 Jun 2014 06:57:43 +0000 (08:57 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Wed, 11 Jun 2014 18:43:49 +0000 (20:43 +0200)
Several tools were writing xvgr formatting code directly
to output files, even when the users selects -xvg none as
a command line option.

Fixes #1407, #1479.

Change-Id: I1db10c2ad7455332e2a415a922885cddaec5efd1

18 files changed:
src/tools/anadih.c
src/tools/gmx_anaeig.c
src/tools/gmx_analyze.c
src/tools/gmx_angle.c
src/tools/gmx_chi.c
src/tools/gmx_cluster.c
src/tools/gmx_dos.c
src/tools/gmx_enemat.c
src/tools/gmx_energy.c
src/tools/gmx_kinetics.c
src/tools/gmx_mdmat.c
src/tools/gmx_mindist.c
src/tools/gmx_polystat.c
src/tools/gmx_rama.c
src/tools/gmx_rms.c
src/tools/gmx_tcaf.c
src/tools/gmx_vanhove.c
src/tools/gmx_wham.c

index af9e6eb0be33703a1584951ed3b5dcef011610e5..5f5da0eefba94486a9cf09ebfc4f653ca0e5e3bc 100644 (file)
@@ -568,9 +568,12 @@ void get_chi_product_traj (real **dih, int nframes, int nangles, int nlist,
                 sprintf(histitle, "cumulative rotamer distribution for %s", dlist[i].name);
                 fprintf(stderr, "  and %s  ", hisfile);
                 fp = xvgropen(hisfile, histitle, "number", "", oenv);
-                fprintf(fp, "@ xaxis tick on\n");
-                fprintf(fp, "@ xaxis tick major 1\n");
-                fprintf(fp, "@ type xy\n");
+                if(output_env_get_print_xvgr_codes(oenv))
+                {
+                    fprintf(fp, "@ xaxis tick on\n");
+                    fprintf(fp, "@ xaxis tick major 1\n");
+                    fprintf(fp, "@ type xy\n");
+                }
                 for (k = 0; (k < nbin); k++)
                 {
                     if (bNormalize)
@@ -582,7 +585,7 @@ void get_chi_product_traj (real **dih, int nframes, int nangles, int nlist,
                         fprintf(fp, "%5d  %10d\n", k, chi_prhist[k]);
                     }
                 }
-                fprintf(fp, "&\n");
+                fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
                 ffclose(fp);
             }
 
index b01c5754baad4f2cde943cbaa65f874eba034257..54f2dae555016870a05dca26f98bb74985c3fe30 100644 (file)
@@ -260,18 +260,12 @@ static void write_xvgr_graphs(const char *file, int ngraphs, int nsetspergraph,
             {
                 if (bSplit && i > 0 && fabs(x[i]) < 1e-5)
                 {
-                    if (output_env_get_print_xvgr_codes(oenv))
-                    {
-                        fprintf(out, "&\n");
-                    }
+                    fprintf(out, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
                 }
                 fprintf(out, "%10.4f %10.5f\n",
                         x[i]*scale_x, y ? y[g][i] : sy[g][s][i]);
             }
-            if (output_env_get_print_xvgr_codes(oenv))
-            {
-                fprintf(out, "&\n");
-            }
+            fprintf(out, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         }
     }
     ffclose(out);
@@ -473,8 +467,10 @@ static void overlap(const char *outfile, int natoms,
 
     out = xvgropen(outfile, "Subspace overlap",
                    "Eigenvectors of trajectory 2", "Overlap", oenv);
-    fprintf(out, "@ subtitle \"using %d eigenvectors of trajectory 1\"\n",
-            noutvec);
+    if(output_env_get_print_xvgr_codes(oenv))
+    {
+        fprintf(out, "@ subtitle \"using %d eigenvectors of trajectory 1\"\n",noutvec);
+    }
     overlap = 0;
     for (x = 0; x < nvec2; x++)
     {
@@ -673,7 +669,7 @@ static void project(const char *trajfile, t_topology *top, int ePBC, matrix topb
         {
             if (bSplit && i > 0 && fabs(inprod[noutvec][i]) < 1e-5)
             {
-                fprintf(xvgrout, "&\n");
+                fprintf(xvgrout, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             }
             fprintf(xvgrout, "%10.5f %10.5f\n", inprod[0][i], inprod[noutvec-1][i]);
         }
index aef05b964b6e928d337ab7e2e63d0c8c4274df2b..2b8f2827b3bf3c95c2e14fefe43874bdc2319a7f 100644 (file)
@@ -295,7 +295,7 @@ void histogram(const char *distfile, real binwidth, int n, int nset, real **val,
         }
         if (s < nset-1)
         {
-            fprintf(fp, "&\n");
+            fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         }
     }
     ffclose(fp);
@@ -615,9 +615,11 @@ static void estimate_error(const char *eefile, int nb_min, int resol, int n,
         }
         fprintf(stdout, "Set %3d:  err.est. %g  a %g  tau1 %g  tau2 %g\n",
                 s+1, ee, a, tau1, tau2);
-        fprintf(fp, "@ legend string %d \"av %f\"\n", 2*s, av[s]);
-        fprintf(fp, "@ legend string %d \"ee %6g\"\n",
-                2*s+1, sig[s]*anal_ee_inf(fitparm, n*dt));
+        if(output_env_get_print_xvgr_codes(oenv))
+        {
+            fprintf(fp, "@ legend string %d \"av %f\"\n", 2*s, av[s]);
+            fprintf(fp, "@ legend string %d \"ee %6g\"\n",2*s+1, sig[s]*anal_ee_inf(fitparm, n*dt));
+        }
         for (i = 0; i < nbs; i++)
         {
             fprintf(fp, "%g %g %g\n", tbs[i], sig[s]*sqrt(ybs[i]/(n*dt)),
@@ -673,7 +675,7 @@ static void estimate_error(const char *eefile, int nb_min, int resol, int n,
                     s+1, sig[s]*anal_ee_inf(ac_fit, n*dt),
                     ac_fit[1], ac_fit[0], ac_fit[2]);
 
-            fprintf(fp, "&\n");
+            fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             for (i = 0; i < nbs; i++)
             {
                 fprintf(fp, "%g %g\n", tbs[i],
@@ -684,7 +686,7 @@ static void estimate_error(const char *eefile, int nb_min, int resol, int n,
         }
         if (s < nset-1)
         {
-            fprintf(fp, "&\n");
+            fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         }
     }
     sfree(fitsig);
@@ -1354,7 +1356,7 @@ int gmx_analyze(int argc, char *argv[])
             }
             if (s < nset-1)
             {
-                fprintf(out, "&\n");
+                fprintf(out, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             }
         }
         ffclose(out);
index f3b6321987b055eafb186b16ddecf6677e1b3b56..da7f4e5865419655fd405a77e46f0f096d58154b 100644 (file)
@@ -442,15 +442,18 @@ int gmx_g_angle(int argc, char *argv[])
         {
             maxstat = max(maxstat, angstat[i]*norm_fac);
         }
-        fprintf(out, "@with g0\n");
-        fprintf(out, "@    world xmin -180\n");
-        fprintf(out, "@    world xmax  180\n");
-        fprintf(out, "@    world ymin 0\n");
-        fprintf(out, "@    world ymax %g\n", maxstat*1.05);
-        fprintf(out, "@    xaxis  tick major 60\n");
-        fprintf(out, "@    xaxis  tick minor 30\n");
-        fprintf(out, "@    yaxis  tick major 0.005\n");
-        fprintf(out, "@    yaxis  tick minor 0.0025\n");
+        if (output_env_get_print_xvgr_codes(oenv))
+        {
+            fprintf(out, "@with g0\n");
+            fprintf(out, "@    world xmin -180\n");
+            fprintf(out, "@    world xmax  180\n");
+            fprintf(out, "@    world ymin 0\n");
+            fprintf(out, "@    world ymax %g\n", maxstat*1.05);
+            fprintf(out, "@    xaxis  tick major 60\n");
+            fprintf(out, "@    xaxis  tick minor 30\n");
+            fprintf(out, "@    yaxis  tick major 0.005\n");
+            fprintf(out, "@    yaxis  tick minor 0.0025\n");
+        }
     }
     for (i = first; (i <= last); i++)
     {
index 40fc9dc255600b018b55a6a524feb41c68815dce..ffc3d95ce6438d376edc3f0b03ea8ddf141bf754 100644 (file)
@@ -766,16 +766,22 @@ static void histogramming(FILE *log, int nbin, gmx_residuetype_t rt,
                 strcpy(hhisfile, hisfile);
                 strcat(hhisfile, ".xvg");
                 fp = xvgropen(hhisfile, title, "Degrees", "", oenv);
-                fprintf(fp, "@ with g0\n");
+                if (output_env_get_print_xvgr_codes(oenv))
+                {
+                    fprintf(fp, "@ with g0\n");
+                }
                 xvgr_world(fp, -180, 0, 180, 0.1, oenv);
-                fprintf(fp, "# this effort to set graph size fails unless you run with -autoscale none or -autoscale y flags\n");
-                fprintf(fp, "@ xaxis tick on\n");
-                fprintf(fp, "@ xaxis tick major 90\n");
-                fprintf(fp, "@ xaxis tick minor 30\n");
-                fprintf(fp, "@ xaxis ticklabel prec 0\n");
-                fprintf(fp, "@ yaxis tick off\n");
-                fprintf(fp, "@ yaxis ticklabel off\n");
-                fprintf(fp, "@ type xy\n");
+                if (output_env_get_print_xvgr_codes(oenv))
+                {
+                    fprintf(fp, "# this effort to set graph size fails unless you run with -autoscale none or -autoscale y flags\n");
+                    fprintf(fp, "@ xaxis tick on\n");
+                    fprintf(fp, "@ xaxis tick major 90\n");
+                    fprintf(fp, "@ xaxis tick minor 30\n");
+                    fprintf(fp, "@ xaxis ticklabel prec 0\n");
+                    fprintf(fp, "@ yaxis tick off\n");
+                    fprintf(fp, "@ yaxis ticklabel off\n");
+                    fprintf(fp, "@ type xy\n");
+                }
                 if (bSSHisto)
                 {
                     for (k = 0; (k < 3); k++)
@@ -804,13 +810,13 @@ static void histogramming(FILE *log, int nbin, gmx_residuetype_t rt,
                         }
                     }
                 }
-                fprintf(fp, "&\n");
+                fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
                 ffclose(fp);
                 if (bSSHisto)
                 {
                     for (k = 0; (k < 3); k++)
                     {
-                        fprintf(ssfp[k], "&\n");
+                        fprintf(ssfp[k], "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
                         ffclose(ssfp[k]);
                     }
                 }
@@ -845,30 +851,35 @@ static FILE *rama_file(const char *fn, const char *title, const char *xaxis,
     FILE *fp;
 
     fp = xvgropen(fn, title, xaxis, yaxis, oenv);
-    fprintf(fp, "@ with g0\n");
+    if (output_env_get_print_xvgr_codes(oenv))
+    {
+        fprintf(fp, "@ with g0\n");
+    }
     xvgr_world(fp, -180, -180, 180, 180, oenv);
-    fprintf(fp, "@ xaxis tick on\n");
-    fprintf(fp, "@ xaxis tick major 90\n");
-    fprintf(fp, "@ xaxis tick minor 30\n");
-    fprintf(fp, "@ xaxis ticklabel prec 0\n");
-    fprintf(fp, "@ yaxis tick on\n");
-    fprintf(fp, "@ yaxis tick major 90\n");
-    fprintf(fp, "@ yaxis tick minor 30\n");
-    fprintf(fp, "@ yaxis ticklabel prec 0\n");
-    fprintf(fp, "@    s0 type xy\n");
-    fprintf(fp, "@    s0 symbol 2\n");
-    fprintf(fp, "@    s0 symbol size 0.410000\n");
-    fprintf(fp, "@    s0 symbol fill 1\n");
-    fprintf(fp, "@    s0 symbol color 1\n");
-    fprintf(fp, "@    s0 symbol linewidth 1\n");
-    fprintf(fp, "@    s0 symbol linestyle 1\n");
-    fprintf(fp, "@    s0 symbol center false\n");
-    fprintf(fp, "@    s0 symbol char 0\n");
-    fprintf(fp, "@    s0 skip 0\n");
-    fprintf(fp, "@    s0 linestyle 0\n");
-    fprintf(fp, "@    s0 linewidth 1\n");
-    fprintf(fp, "@ type xy\n");
-
+    if (output_env_get_print_xvgr_codes(oenv))
+    {
+        fprintf(fp, "@ xaxis tick on\n");
+        fprintf(fp, "@ xaxis tick major 90\n");
+        fprintf(fp, "@ xaxis tick minor 30\n");
+        fprintf(fp, "@ xaxis ticklabel prec 0\n");
+        fprintf(fp, "@ yaxis tick on\n");
+        fprintf(fp, "@ yaxis tick major 90\n");
+        fprintf(fp, "@ yaxis tick minor 30\n");
+        fprintf(fp, "@ yaxis ticklabel prec 0\n");
+        fprintf(fp, "@    s0 type xy\n");
+        fprintf(fp, "@    s0 symbol 2\n");
+        fprintf(fp, "@    s0 symbol size 0.410000\n");
+        fprintf(fp, "@    s0 symbol fill 1\n");
+        fprintf(fp, "@    s0 symbol color 1\n");
+        fprintf(fp, "@    s0 symbol linewidth 1\n");
+        fprintf(fp, "@    s0 symbol linestyle 1\n");
+        fprintf(fp, "@    s0 symbol center false\n");
+        fprintf(fp, "@    s0 symbol char 0\n");
+        fprintf(fp, "@    s0 skip 0\n");
+        fprintf(fp, "@    s0 linestyle 0\n");
+        fprintf(fp, "@    s0 linewidth 1\n");
+        fprintf(fp, "@ type xy\n");
+    }
     return fp;
 }
 
index 7658137242a1b8657c6308c14d3d1a7c5dec3420..c106ea68817862aebe764c791155a4095183a8bb 100644 (file)
@@ -1039,9 +1039,12 @@ static void analyze_clusters(int nf, t_clusters *clust, real **rmsd,
     if (clustidfn)
     {
         fp = xvgropen(clustidfn, "Clusters", output_env_get_xvgr_tlabel(oenv), "Cluster #", oenv);
-        fprintf(fp, "@    s0 symbol 2\n");
-        fprintf(fp, "@    s0 symbol size 0.2\n");
-        fprintf(fp, "@    s0 linestyle 0\n");
+        if(output_env_get_print_xvgr_codes(oenv))
+        {
+            fprintf(fp, "@    s0 symbol 2\n");
+            fprintf(fp, "@    s0 symbol size 0.2\n");
+            fprintf(fp, "@    s0 linestyle 0\n");
+        }
         for (i = 0; i < nf; i++)
         {
             fprintf(fp, "%8g %8d\n", time[i], clust->cl[i]);
@@ -1051,7 +1054,10 @@ static void analyze_clusters(int nf, t_clusters *clust, real **rmsd,
     if (sizefn)
     {
         fp = xvgropen(sizefn, "Cluster Sizes", "Cluster #", "# Structures", oenv);
-        fprintf(fp, "@g%d type %s\n", 0, "bar");
+        if(output_env_get_print_xvgr_codes(oenv))
+        {
+            fprintf(fp, "@g%d type %s\n", 0, "bar");
+        }
     }
     snew(structure, nf);
     fprintf(log, "\n%3s | %3s  %4s | %6s %4s | cluster members\n",
index f5432dc7b1a949b5cb65b1dc42f740eabea3ea51..9c9fcf284c6429ca860b2d95b19ac92e380e61a9 100644 (file)
@@ -225,8 +225,11 @@ static void dump_fy(output_env_t oenv, real toler)
     DD = pow(10.0, 0.125);
     fp = xvgropen("fy.xvg", "Fig. 2, Lin2003a", "Delta", "y or fy", oenv);
     xvgr_legend(fp, asize(leg), leg, oenv);
-    fprintf(fp, "@    world 1e-05, 0, 1000, 1\n");
-    fprintf(fp, "@    xaxes scale Logarithmic\n");
+    if (output_env_get_print_xvgr_codes(oenv))
+    {
+        fprintf(fp, "@    world 1e-05, 0, 1000, 1\n");
+        fprintf(fp, "@    xaxes scale Logarithmic\n");
+    }
     for (Delta = 1e-5; (Delta <= 1000); Delta *= DD)
     {
         f = calc_fluidicity(Delta, toler);
index f9e76c79ea208b64a97ebe4737741d9ec4e4d547..f60972502b42b6edbd0c4e14d27e3b0a48b70da9 100644 (file)
@@ -512,39 +512,43 @@ int gmx_enemat(int argc, char *argv[])
                        oenv);
         xvgr_legend(out, 0, NULL, oenv);
         j = 0;
-        for (m = 0; (m < egNR+egSP); m++)
+        if (output_env_get_print_xvgr_codes(oenv))
         {
-            if (egrp_use[m])
+            for (m = 0; (m < egNR+egSP); m++)
             {
-                fprintf(out, "@ legend string %d \"%s\"\n", j++, egrp_nm[m]);
+                if (egrp_use[m])
+                {
+                    fprintf(out, "@ legend string %d \"%s\"\n", j++, egrp_nm[m]);
+                }
             }
-        }
-        if (bFree)
-        {
-            fprintf(out, "@ legend string %d \"%s\"\n", j++, "Free");
-        }
-        if (bFree)
-        {
-            fprintf(out, "@ legend string %d \"%s\"\n", j++, "Diff");
-        }
-        fprintf(out, "@TYPE xy\n");
-        fprintf(out, "#%3s", "grp");
-        for (m = 0; (m < egNR+egSP); m++)
-        {
-            if (egrp_use[m])
+            if (bFree)
             {
-                fprintf(out, " %9s", egrp_nm[m]);
+                fprintf(out, "@ legend string %d \"%s\"\n", j++, "Free");
             }
+            if (bFree)
+            {
+                fprintf(out, "@ legend string %d \"%s\"\n", j++, "Diff");
+            }
+            fprintf(out, "@TYPE xy\n");
+            fprintf(out, "#%3s", "grp");
+
+            for (m = 0; (m < egNR+egSP); m++)
+            {
+                if (egrp_use[m])
+                {
+                    fprintf(out, " %9s", egrp_nm[m]);
+                }
+            }
+            if (bFree)
+            {
+                fprintf(out, " %9s", "Free");
+            }
+            if (bFree)
+            {
+                fprintf(out, " %9s", "Diff");
+            }
+            fprintf(out, "\n");
         }
-        if (bFree)
-        {
-            fprintf(out, " %9s", "Free");
-        }
-        if (bFree)
-        {
-            fprintf(out, " %9s", "Diff");
-        }
-        fprintf(out, "\n");
         for (i = 0; (i < ngroups); i++)
         {
             fprintf(out, "%3.0f", groupnr[i]);
index 3323f35786b0e755d16dfabbfe5c436e1f278ddd..7341ddbd34376e77bfc9aebe97bfb422d37e5d72 100644 (file)
@@ -2253,7 +2253,7 @@ int gmx_energy(int argc, char *argv[])
                 {
                     fort = xvgropen(opt2fn("-ort", NFILE, fnm), "Calculated orientations",
                                     "Time (ps)", "", oenv);
-                    if (bOrinst)
+                    if (bOrinst && output_env_get_print_xvgr_codes(oenv))
                     {
                         fprintf(fort, "%s", orinst_sub);
                     }
@@ -2264,7 +2264,7 @@ int gmx_energy(int argc, char *argv[])
                     fodt = xvgropen(opt2fn("-odt", NFILE, fnm),
                                     "Orientation restraint deviation",
                                     "Time (ps)", "", oenv);
-                    if (bOrinst)
+                    if (bOrinst && output_env_get_print_xvgr_codes(oenv))
                     {
                         fprintf(fodt, "%s", orinst_sub);
                     }
@@ -2730,7 +2730,7 @@ int gmx_energy(int argc, char *argv[])
         out = xvgropen(opt2fn("-ora", NFILE, fnm),
                        "Average calculated orientations",
                        "Restraint label", "", oenv);
-        if (bOrinst)
+        if (bOrinst && output_env_get_print_xvgr_codes(oenv))
         {
             fprintf(out, "%s", orinst_sub);
         }
@@ -2745,7 +2745,7 @@ int gmx_energy(int argc, char *argv[])
         out = xvgropen(opt2fn("-oda", NFILE, fnm),
                        "Average restraint deviation",
                        "Restraint label", "", oenv);
-        if (bOrinst)
+        if (bOrinst && output_env_get_print_xvgr_codes(oenv))
         {
             fprintf(out, "%s", orinst_sub);
         }
@@ -2760,7 +2760,7 @@ int gmx_energy(int argc, char *argv[])
         out = xvgropen(opt2fn("-odr", NFILE, fnm),
                        "RMS orientation restraint deviations",
                        "Restraint label", "", oenv);
-        if (bOrinst)
+        if (bOrinst && output_env_get_print_xvgr_codes(oenv))
         {
             fprintf(out, "%s", orinst_sub);
         }
index dce1cb7bec0068f6672b7089803cda297135c232..06fa11329d7a926c57108d6a0268aa1c24a70f50 100644 (file)
@@ -803,7 +803,10 @@ static void dump_remd_parameters(FILE *gp, t_remd_data *d, const char *fn,
         for (j = 0; (j < d->nparams); j++)
         {
             /* Reset all parameters to optimized values */
-            fprintf(hp, "@type xy\n");
+            if(output_env_get_print_xvgr_codes(oenv))
+            {
+                fprintf(hp, "@type xy\n");
+            }
             for (i = 0; (i < d->nparams); i++)
             {
                 d->params[i] = params[i];
@@ -816,7 +819,7 @@ static void dump_remd_parameters(FILE *gp, t_remd_data *d, const char *fn,
                 fprintf(gp, "%s = %12g  d2 = %12g\n", epnm(np, j), d->params[j], d2[i]);
                 fprintf(hp, "%12g  %12g\n", fac[i], d2[i]);
             }
-            fprintf(hp, "&\n");
+            fprintf(hp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         }
         ffclose(hp);
         for (i = 0; (i < d->nparams); i++)
index 5e73faca2ca187a457ed74b40e8959322414330f..a5df7ef22de39ab761c40b195c4ed449f1e2aa89 100644 (file)
@@ -371,17 +371,20 @@ int gmx_mdmat(int argc, char *argv[])
         tot_nmat(nres, natoms, nframes, totnmat, tot_n, mean_n);
         fp = xvgropen(ftp2fn(efXVG, NFILE, fnm),
                       "Increase in number of contacts", "Residue", "Ratio", oenv);
-        fprintf(fp, "@ legend on\n");
-        fprintf(fp, "@ legend box on\n");
-        fprintf(fp, "@ legend loctype view\n");
-        fprintf(fp, "@ legend 0.75, 0.8\n");
-        fprintf(fp, "@ legend string 0 \"Total/mean\"\n");
-        fprintf(fp, "@ legend string 1 \"Total\"\n");
-        fprintf(fp, "@ legend string 2 \"Mean\"\n");
-        fprintf(fp, "@ legend string 3 \"# atoms\"\n");
-        fprintf(fp, "@ legend string 4 \"Mean/# atoms\"\n");
-        fprintf(fp, "#%3s %8s  %3s  %8s  %3s  %8s\n",
-                "res", "ratio", "tot", "mean", "natm", "mean/atm");
+        if(output_env_get_print_xvgr_codes(oenv))
+        {
+            fprintf(fp, "@ legend on\n");
+            fprintf(fp, "@ legend box on\n");
+            fprintf(fp, "@ legend loctype view\n");
+            fprintf(fp, "@ legend 0.75, 0.8\n");
+            fprintf(fp, "@ legend string 0 \"Total/mean\"\n");
+            fprintf(fp, "@ legend string 1 \"Total\"\n");
+            fprintf(fp, "@ legend string 2 \"Mean\"\n");
+            fprintf(fp, "@ legend string 3 \"# atoms\"\n");
+            fprintf(fp, "@ legend string 4 \"Mean/# atoms\"\n");
+            fprintf(fp, "#%3s %8s  %3s  %8s  %3s  %8s\n",
+                    "res", "ratio", "tot", "mean", "natm", "mean/atm");
+        }
         for (i = 0; (i < nres); i++)
         {
             if (mean_n[i] == 0)
index 30d5abbff38811151395b16a2313c6d92cb5b248..94718f20108ef3def64c60649419e61c2dbc43d8 100644 (file)
@@ -174,7 +174,7 @@ static void periodic_mindist_plot(const char *trxfn, const char *outfn,
         }
         if (bSplit && !bFirst && fabs(t/output_env_get_time_factor(oenv)) < 1e-5)
         {
-            fprintf(out, "&\n");
+            fprintf(out, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         }
         fprintf(out, "\t%g\t%6.3f %6.3f %6.3f %6.3f %6.3f\n",
                 output_env_conv_time(oenv, t), rmin, rmax, norm(box[0]), norm(box[1]), norm(box[2]));
@@ -425,14 +425,14 @@ void dist_plot(const char *fn, const char *afile, const char *dfile,
     {
         if (bSplit && !bFirst && fabs(t/output_env_get_time_factor(oenv)) < 1e-5)
         {
-            fprintf(dist, "&\n");
+            fprintf(dist, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             if (num)
             {
-                fprintf(num, "&\n");
+                fprintf(num, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             }
             if (atm)
             {
-                fprintf(atm, "&\n");
+                fprintf(atm, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             }
         }
         fprintf(dist, "%12e", output_env_conv_time(oenv, t));
index 064bbc1bad9f476cd6af6b6752e9e19176a87a95..1b3e5f17cef443f91062768a9a270dd87ac23192 100644 (file)
@@ -493,7 +493,10 @@ int gmx_polystat(int argc, char *argv[])
     /* Handle printing of internal distances. */
     if (outi)
     {
-        fprintf(outi, "@    xaxes scale Logarithmic\n");
+        if(output_env_get_print_xvgr_codes(oenv))
+        {
+            fprintf(outi, "@    xaxes scale Logarithmic\n");
+        }
         ymax = -1;
         ymin = 1e300;
         j    = index[molind[1]-1] - index[molind[0]]; /* Polymer length -1. */
index 044bfd439328d7a151d1b829583167c7d22ade56..89c268c0bd6f6e9f626e8665d93ff80fafd13363 100644 (file)
@@ -102,10 +102,12 @@ int gmx_rama(int argc, char *argv[])
     xvgr_line_props(out, 0, elNone, ecFrank, oenv);
     xvgr_view(out, 0.2, 0.2, 0.8, 0.8, oenv);
     xvgr_world(out, -180, -180, 180, 180, oenv);
-    fprintf(out, "@    xaxis  tick on\n@    xaxis  tick major 60\n@    xaxis  tick minor 30\n");
-    fprintf(out, "@    yaxis  tick on\n@    yaxis  tick major 60\n@    yaxis  tick minor 30\n");
-    fprintf(out, "@ s0 symbol 2\n@ s0 symbol size 0.4\n@ s0 symbol fill 1\n");
-
+    if(output_env_get_print_xvgr_codes(oenv))
+    {
+        fprintf(out, "@    xaxis  tick on\n@    xaxis  tick major 60\n@    xaxis  tick minor 30\n");
+        fprintf(out, "@    yaxis  tick on\n@    yaxis  tick major 60\n@    yaxis  tick minor 30\n");
+        fprintf(out, "@ s0 symbol 2\n@ s0 symbol size 0.4\n@ s0 symbol fill 1\n");
+    }
     j = 0;
     do
     {
index 895bdb836e7719e1b0fde0a2418f1fb85ac74097..488f2a36d8663a830c02d34d765442461c71abf5 100644 (file)
@@ -1139,7 +1139,7 @@ int gmx_rms(int argc, char *argv[])
         if (bSplit && i > 0 &&
             fabs(time[bPrev ? freq*i : i]/output_env_get_time_factor(oenv)) < 1e-5)
         {
-            fprintf(fp, "&\n");
+            fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         }
         fprintf(fp, "%12.7f", time[bPrev ? freq*i : i]);
         for (j = 0; (j < nrms); j++)
@@ -1181,7 +1181,7 @@ int gmx_rms(int argc, char *argv[])
         {
             if (bSplit && i > 0 && fabs(time[i]) < 1e-5)
             {
-                fprintf(fp, "&\n");
+                fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             }
             fprintf(fp, "%12.7f", time[i]);
             for (j = 0; (j < nrms); j++)
index 87d303d9b738560f8a8a08f32280e5c813208ef1..074e53c12407824919de7e2189eacd1ec3c8d2db 100644 (file)
@@ -185,20 +185,23 @@ static void process_tcaf(int nframes, real dt, int nkc, real **tc, rvec *kfac,
                 tcafc[kc][i] /= tcafc[kc][0];
                 fprintf(fp_cub, "%g %g\n", i*dt, tcafc[kc][i]);
             }
-            fprintf(fp_cub, "&\n");
+            fprintf(fp_cub, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             tcafc[kc][0] = 1.0;
         }
     }
 
     fp_vk = xvgropen(fn_vk, "Fits", "k (nm\\S-1\\N)",
                      "\\8h\\4 (10\\S-3\\N kg m\\S-1\\N s\\S-1\\N)", oenv);
-    fprintf(fp_vk, "@    s0 symbol 2\n");
-    fprintf(fp_vk, "@    s0 symbol color 1\n");
-    fprintf(fp_vk, "@    s0 linestyle 0\n");
-    if (fn_cub)
+    if(output_env_get_print_xvgr_codes(oenv))
     {
-        fprintf(fp_vk, "@    s1 symbol 3\n");
-        fprintf(fp_vk, "@    s1 symbol color 2\n");
+        fprintf(fp_vk, "@    s0 symbol 2\n");
+        fprintf(fp_vk, "@    s0 symbol color 1\n");
+        fprintf(fp_vk, "@    s0 linestyle 0\n");
+        if (fn_cub)
+        {
+            fprintf(fp_vk, "@    s1 symbol 3\n");
+            fprintf(fp_vk, "@    s1 symbol color 2\n");
+        }
     }
     fp = xvgropen(fn_tcf, "TCAF Fits", "Time (ps)", "", oenv);
     for (k = 0; k < nk; k++)
@@ -217,7 +220,7 @@ static void process_tcaf(int nframes, real dt, int nkc, real **tc, rvec *kfac,
         {
             fprintf(fp, "%g %g\n", i*dt, fit_function(effnVAC, fitparms, i*dt));
         }
-        fprintf(fp, "&\n");
+        fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
     }
     ffclose(fp);
     do_view(oenv, fn_tcf, "-nxy");
@@ -225,7 +228,7 @@ static void process_tcaf(int nframes, real dt, int nkc, real **tc, rvec *kfac,
     if (fn_cub)
     {
         fprintf(stdout, "Averaged over k-vectors:\n");
-        fprintf(fp_vk, "&\n");
+        fprintf(fp_vk, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         for (k = 0; k < nkc; k++)
         {
             tcafc[k][0]  = 1.0;
@@ -243,9 +246,9 @@ static void process_tcaf(int nframes, real dt, int nkc, real **tc, rvec *kfac,
             {
                 fprintf(fp_cub, "%g %g\n", i*dt, fit_function(effnVAC, fitparms, i*dt));
             }
-            fprintf(fp_cub, "&\n");
+            fprintf(fp_cub, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         }
-        fprintf(fp_vk, "&\n");
+        fprintf(fp_vk, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         ffclose(fp_cub);
         do_view(oenv, fn_cub, "-nxy");
     }
index aea369f3f26a82672656aa39daa939b6af810ab1..93ec6b85c65e0331d1bc0e420595d1172bc04b31 100644 (file)
@@ -435,7 +435,10 @@ int gmx_vanhove(int argc, char *argv[])
     if (orfile)
     {
         fp = xvgropen(orfile, "Van Hove function", "r (nm)", "G (nm\\S-1\\N)", oenv);
-        fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
+        if(output_env_get_print_xvgr_codes(oenv))
+        {
+            fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
+        }
         snew(legend, nr);
         for (fbin = 0; fbin < nr; fbin++)
         {
@@ -460,7 +463,10 @@ int gmx_vanhove(int argc, char *argv[])
     {
         sprintf(buf, "Probability of moving less than %g nm", rint);
         fp = xvgropen(otfile, buf, "t (ps)", "", oenv);
-        fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
+        if(output_env_get_print_xvgr_codes(oenv))
+        {
+            fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
+        }
         for (f = 0; f <= ftmax; f++)
         {
             fprintf(fp, "%g %g\n", f*dt, (real)pt[f]/(tcount[f]*isize));
index 22257ead8599450caffb771f6326af1d1e4ffb80..38fd567f6ee37c9274fe6fe0afd4458ca8fbcea8 100644 (file)
@@ -1596,13 +1596,16 @@ void do_bootstrapping(const char *fnres, const char* fnprof, const char *fnhist,
             bsProfiles_av2[i] += tmp*tmp;
             fprintf(fp, "%e\t%e\n", (i+0.5)*opt->dz+opt->min, tmp);
         }
-        fprintf(fp, "&\n");
+        fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(opt->oenv) ? "&" : "");
     }
     ffclose(fp);
 
     /* write average and stddev */
     fp = xvgropen(fnres, "Average and stddev from bootstrapping", "z", ylabel, opt->oenv);
-    fprintf(fp, "@TYPE xydy\n");
+    if(output_env_get_print_xvgr_codes(opt->oenv))
+    {
+        fprintf(fp, "@TYPE xydy\n");
+    }
     for (i = 0; i < opt->bins; i++)
     {
         bsProfiles_av [i] /= opt->nBootStrap;
@@ -2523,7 +2526,7 @@ void calcIntegratedAutocorrelationTimes(t_UmbrellaWindow *window, int nwins,
                 {
                     fprintf(fpcorr, "%g  %g\n", k*dt, corr[k]);
                 }
-                fprintf(fpcorr, "&\n");
+                fprintf(fpcorr, "%s\n", output_env_get_print_xvgr_codes(opt->oenv) ? "&" : "");
             }
 
             /* esimate integrated correlation time, fitting is too unstable */
@@ -2551,16 +2554,19 @@ void calcIntegratedAutocorrelationTimes(t_UmbrellaWindow *window, int nwins,
 
     /* plot IACT along reaction coordinate */
     fp = xvgropen(fn, "Integrated autocorrelation times", "z", "IACT [ps]", opt->oenv);
-    fprintf(fp, "@    s0 symbol 1\n@    s0 symbol size 0.5\n@    s0 line linestyle 0\n");
-    fprintf(fp, "#  WIN   tau(gr1)  tau(gr2) ...\n");
-    for (i = 0; i < nwins; i++)
+    if(output_env_get_print_xvgr_codes(opt->oenv))
     {
-        fprintf(fp, "# %3d   ", i);
-        for (ig = 0; ig < window[i].nPull; ig++)
+        fprintf(fp, "@    s0 symbol 1\n@    s0 symbol size 0.5\n@    s0 line linestyle 0\n");
+        fprintf(fp, "#  WIN   tau(gr1)  tau(gr2) ...\n");
+        for (i = 0; i < nwins; i++)
         {
-            fprintf(fp, " %11g", window[i].tau[ig]);
+            fprintf(fp, "# %3d   ", i);
+            for (ig = 0; ig < window[i].nPull; ig++)
+            {
+                fprintf(fp, " %11g", window[i].tau[ig]);
+            }
+            fprintf(fp, "\n");
         }
-        fprintf(fp, "\n");
     }
     for (i = 0; i < nwins; i++)
     {
@@ -2575,8 +2581,12 @@ void calcIntegratedAutocorrelationTimes(t_UmbrellaWindow *window, int nwins,
                opt->sigSmoothIact);
         /* smooth IACT along reaction coordinate and overwrite g=1+2tau */
         smoothIact(window, nwins, opt);
-        fprintf(fp, "&\n@    s1 symbol 1\n@    s1 symbol size 0.5\n@    s1 line linestyle 0\n");
-        fprintf(fp, "@    s1 symbol color 2\n");
+        fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(opt->oenv) ? "&" : "");
+        if(output_env_get_print_xvgr_codes(opt->oenv))
+        {
+            fprintf(fp, "@    s1 symbol 1\n@    s1 symbol size 0.5\n@    s1 line linestyle 0\n");
+            fprintf(fp, "@    s1 symbol color 2\n");
+        }
         for (i = 0; i < nwins; i++)
         {
             for (ig = 0; ig < window[i].nPull; ig++)