Move symtab.* to topology/
[gromacs.git] / src / gromacs / gmxana / pp2shift.c
blobe95d80b2b7e24f86250c2b1200c39b69ca3fde73
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
41 #include <stdlib.h>
42 #include <math.h>
43 #include "typedefs.h"
44 #include "gromacs/utility/futil.h"
45 #include "macros.h"
46 #include "gromacs/math/units.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "gstat.h"
49 #include "gromacs/fileio/matio.h"
50 #include "copyrite.h"
51 #include "gromacs/utility/fatalerror.h"
53 typedef struct {
54 int nx, ny;
55 real dx, dy;
56 real **data;
57 } t_shiftdata;
59 static real interpolate(real phi, real psi, t_shiftdata *sd)
61 int iphi, ipsi, iphi1, ipsi1;
62 real fphi, fpsi, wx0, wx1, wy0, wy1;
64 /*phi += M_PI;
65 if (phi > 2*M_PI) phi -= 2*M_PI;
66 psi += M_PI;
67 if (psi > 2*M_PI) psi -= 2*M_PI;
69 while (phi < 0)
71 phi += 2*M_PI;
73 while (psi < 0)
75 psi += 2*M_PI;
77 phi = 2*M_PI-phi;
79 fphi = phi*sd->dx;
80 fpsi = psi*sd->dy;
82 iphi = (int)fphi;
83 ipsi = (int)fpsi;
84 fphi -= iphi; /* Fraction (offset from gridpoint) */
85 fpsi -= ipsi;
87 wx0 = 1.0-fphi;
88 wx1 = fphi;
89 wy0 = 1.0-fpsi;
90 wy1 = fpsi;
91 iphi = iphi % sd->nx;
92 ipsi = ipsi % sd->ny;
93 iphi1 = (iphi+1) % sd->nx;
94 ipsi1 = (ipsi+1) % sd->ny;
96 return (sd->data[iphi] [ipsi] * wx0*wy0 +
97 sd->data[iphi1] [ipsi] * wx1*wy0 +
98 sd->data[iphi] [ipsi1] * wx0*wy1 +
99 sd->data[iphi1] [ipsi1] * wx1*wy1);
102 static void dump_sd(const char *fn, t_shiftdata *sd)
104 FILE *fp;
105 int i, j;
106 char buf[256];
107 int nnx, nny, nfac = 4, nlevels = 20;
108 real phi, psi, *x_phi, *y_psi, **newdata;
109 real lo, hi;
110 t_rgb rlo = { 1, 0, 0 }, rhi = { 0, 0, 1 };
112 nnx = sd->nx*nfac+1;
113 nny = sd->ny*nfac+1;
114 snew(x_phi, nnx);
115 snew(y_psi, nny);
116 snew(newdata, nnx);
117 lo = 100000;
118 hi = -100000;
119 for (i = 0; (i < nnx); i++)
121 snew(newdata[i], nny);
122 phi = i*2*M_PI/(nnx-1);
123 x_phi[i] = phi*RAD2DEG-180;
124 for (j = 0; (j < nny); j++)
126 psi = j*2*M_PI/(nny-1);
127 if (i == 0)
129 y_psi[j] = psi*RAD2DEG-180;
131 /*if (((i % nfac) == 0) && ((j % nfac) == 0))
132 newdata[i][j] = sd->data[i/nfac][j/nfac];
133 else*/
134 newdata[i][j] = interpolate(phi, psi, sd);
135 lo = min(lo, newdata[i][j]);
136 hi = max(hi, newdata[i][j]);
139 sprintf(buf, "%s.xpm", fn);
140 fp = gmx_ffopen(buf, "w");
141 write_xpm(fp, 0, fn, fn, "Phi", "Psi", nnx, nny,
142 x_phi, y_psi, newdata, lo, hi, rlo, rhi, &nlevels);
143 for (i = 0; (i < nnx); i++)
145 sfree(newdata[i]);
147 sfree(newdata);
148 sfree(x_phi);
149 sfree(y_psi);
152 static t_shiftdata *read_shifts(const char *fn)
154 FILE *fp;
155 double xx;
156 int i, j, nx, ny;
157 t_shiftdata *sd;
159 snew(sd, 1);
160 fp = libopen(fn);
161 if (2 != fscanf(fp, "%d%d", &nx, &ny))
163 gmx_fatal(FARGS, "Error reading from file %s", fn);
166 sd->nx = nx;
167 sd->ny = ny;
168 sd->dx = nx/(2*M_PI);
169 sd->dy = ny/(2*M_PI);
170 snew(sd->data, nx+1);
171 for (i = 0; (i <= nx); i++)
173 snew(sd->data[i], ny+1);
174 for (j = 0; (j < ny); j++)
176 if (i == nx)
178 sd->data[i][j] = sd->data[0][j];
180 else
182 if (1 != fscanf(fp, "%lf", &xx))
184 gmx_fatal(FARGS, "Error reading from file %s", fn);
186 sd->data[i][j] = xx;
189 sd->data[i][j] = sd->data[i][0];
191 gmx_ffclose(fp);
193 if (bDebugMode())
195 dump_sd(fn, sd);
198 return sd;
202 static void done_shifts(t_shiftdata *sd)
204 int i;
206 for (i = 0; (i <= sd->nx); i++)
208 sfree(sd->data[i]);
210 sfree(sd->data);
211 sfree(sd);
214 void do_pp2shifts(FILE *fp, int nf, int nlist, t_dlist dlist[], real **dih)
216 t_shiftdata *ca_sd, *co_sd, *ha_sd, *cb_sd;
217 int i, j, Phi, Psi;
218 real phi, psi;
219 real ca, co, ha, cb;
221 /* Read the shift files */
222 ca_sd = read_shifts("ca-shift.dat");
223 cb_sd = read_shifts("cb-shift.dat");
224 ha_sd = read_shifts("ha-shift.dat");
225 co_sd = read_shifts("co-shift.dat");
227 fprintf(fp, "\n *** Chemical shifts from the chemical shift index ***\n");
228 please_cite(fp, "Wishart98a");
229 fprintf(fp, "%12s %10s %10s %10s %10s\n",
230 "Residue", "delta Ca", "delta Ha", "delta CO", "delta Cb");
231 for (i = 0; (i < nlist); i++)
233 if ((has_dihedral(edPhi, &(dlist[i]))) &&
234 (has_dihedral(edPsi, &(dlist[i]))))
236 Phi = dlist[i].j0[edPhi];
237 Psi = dlist[i].j0[edPsi];
238 ca = cb = co = ha = 0;
239 for (j = 0; (j < nf); j++)
241 phi = dih[Phi][j];
242 psi = dih[Psi][j];
244 ca += interpolate(phi, psi, ca_sd);
245 cb += interpolate(phi, psi, cb_sd);
246 co += interpolate(phi, psi, co_sd);
247 ha += interpolate(phi, psi, ha_sd);
249 fprintf(fp, "%12s %10g %10g %10g %10g\n",
250 dlist[i].name, ca/nf, ha/nf, co/nf, cb/nf);
253 fprintf(fp, "\n");
255 /* Free memory */
256 done_shifts(ca_sd);
257 done_shifts(cb_sd);
258 done_shifts(co_sd);
259 done_shifts(ha_sd);