3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
42 #include "gromacs/fileio/futil.h"
47 #include "gromacs/fileio/matio.h"
49 #include "gmx_fatal.h"
57 static real interpolate(real phi, real psi, t_shiftdata *sd)
59 int iphi, ipsi, iphi1, ipsi1;
60 real fphi, fpsi, wx0, wx1, wy0, wy1;
63 if (phi > 2*M_PI) phi -= 2*M_PI;
65 if (psi > 2*M_PI) psi -= 2*M_PI;
82 fphi -= iphi; /* Fraction (offset from gridpoint) */
91 iphi1 = (iphi+1) % sd->nx;
92 ipsi1 = (ipsi+1) % sd->ny;
94 return (sd->data[iphi] [ipsi] * wx0*wy0 +
95 sd->data[iphi1] [ipsi] * wx1*wy0 +
96 sd->data[iphi] [ipsi1] * wx0*wy1 +
97 sd->data[iphi1] [ipsi1] * wx1*wy1);
100 static void dump_sd(const char *fn, t_shiftdata *sd)
105 int nnx, nny, nfac = 4, nlevels = 20;
106 real phi, psi, *x_phi, *y_psi, **newdata;
108 t_rgb rlo = { 1, 0, 0 }, rhi = { 0, 0, 1 };
117 for (i = 0; (i < nnx); i++)
119 snew(newdata[i], nny);
120 phi = i*2*M_PI/(nnx-1);
121 x_phi[i] = phi*RAD2DEG-180;
122 for (j = 0; (j < nny); j++)
124 psi = j*2*M_PI/(nny-1);
127 y_psi[j] = psi*RAD2DEG-180;
129 /*if (((i % nfac) == 0) && ((j % nfac) == 0))
130 newdata[i][j] = sd->data[i/nfac][j/nfac];
132 newdata[i][j] = interpolate(phi, psi, sd);
133 lo = min(lo, newdata[i][j]);
134 hi = max(hi, newdata[i][j]);
137 sprintf(buf, "%s.xpm", fn);
138 fp = ffopen(buf, "w");
139 write_xpm(fp, 0, fn, fn, "Phi", "Psi", nnx, nny,
140 x_phi, y_psi, newdata, lo, hi, rlo, rhi, &nlevels);
141 for (i = 0; (i < nnx); i++)
150 static t_shiftdata *read_shifts(const char *fn)
159 if (2 != fscanf(fp, "%d%d", &nx, &ny))
161 gmx_fatal(FARGS, "Error reading from file %s", fn);
166 sd->dx = nx/(2*M_PI);
167 sd->dy = ny/(2*M_PI);
168 snew(sd->data, nx+1);
169 for (i = 0; (i <= nx); i++)
171 snew(sd->data[i], ny+1);
172 for (j = 0; (j < ny); j++)
176 sd->data[i][j] = sd->data[0][j];
180 if (1 != fscanf(fp, "%lf", &xx))
182 gmx_fatal(FARGS, "Error reading from file %s", fn);
187 sd->data[i][j] = sd->data[i][0];
200 static void done_shifts(t_shiftdata *sd)
204 for (i = 0; (i <= sd->nx); i++)
212 void do_pp2shifts(FILE *fp, int nf, int nlist, t_dlist dlist[], real **dih)
214 t_shiftdata *ca_sd, *co_sd, *ha_sd, *cb_sd;
219 /* Read the shift files */
220 ca_sd = read_shifts("ca-shift.dat");
221 cb_sd = read_shifts("cb-shift.dat");
222 ha_sd = read_shifts("ha-shift.dat");
223 co_sd = read_shifts("co-shift.dat");
225 fprintf(fp, "\n *** Chemical shifts from the chemical shift index ***\n");
226 please_cite(fp, "Wishart98a");
227 fprintf(fp, "%12s %10s %10s %10s %10s\n",
228 "Residue", "delta Ca", "delta Ha", "delta CO", "delta Cb");
229 for (i = 0; (i < nlist); i++)
231 if ((has_dihedral(edPhi, &(dlist[i]))) &&
232 (has_dihedral(edPsi, &(dlist[i]))))
234 Phi = dlist[i].j0[edPhi];
235 Psi = dlist[i].j0[edPsi];
236 ca = cb = co = ha = 0;
237 for (j = 0; (j < nf); j++)
242 ca += interpolate(phi, psi, ca_sd);
243 cb += interpolate(phi, psi, cb_sd);
244 co += interpolate(phi, psi, co_sd);
245 ha += interpolate(phi, psi, ha_sd);
247 fprintf(fp, "%12s %10g %10g %10g %10g\n",
248 dlist[i].name, ca/nf, ha/nf, co/nf, cb/nf);