Remove all unnecessary HAVE_CONFIG_H
[gromacs.git] / src / gromacs / gmxana / pp2shift.c
blob53011c92dd495f635f81b4a84ca282fd3dc901c4
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 #include "config.h"
39 #include <stdlib.h>
40 #include <math.h>
41 #include "typedefs.h"
42 #include "gromacs/utility/futil.h"
43 #include "macros.h"
44 #include "gromacs/math/units.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "gstat.h"
47 #include "gromacs/fileio/matio.h"
48 #include "copyrite.h"
49 #include "gromacs/utility/fatalerror.h"
51 typedef struct {
52 int nx, ny;
53 real dx, dy;
54 real **data;
55 } t_shiftdata;
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;
62 /*phi += M_PI;
63 if (phi > 2*M_PI) phi -= 2*M_PI;
64 psi += M_PI;
65 if (psi > 2*M_PI) psi -= 2*M_PI;
67 while (phi < 0)
69 phi += 2*M_PI;
71 while (psi < 0)
73 psi += 2*M_PI;
75 phi = 2*M_PI-phi;
77 fphi = phi*sd->dx;
78 fpsi = psi*sd->dy;
80 iphi = (int)fphi;
81 ipsi = (int)fpsi;
82 fphi -= iphi; /* Fraction (offset from gridpoint) */
83 fpsi -= ipsi;
85 wx0 = 1.0-fphi;
86 wx1 = fphi;
87 wy0 = 1.0-fpsi;
88 wy1 = fpsi;
89 iphi = iphi % sd->nx;
90 ipsi = ipsi % sd->ny;
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)
102 FILE *fp;
103 int i, j;
104 char buf[256];
105 int nnx, nny, nfac = 4, nlevels = 20;
106 real phi, psi, *x_phi, *y_psi, **newdata;
107 real lo, hi;
108 t_rgb rlo = { 1, 0, 0 }, rhi = { 0, 0, 1 };
110 nnx = sd->nx*nfac+1;
111 nny = sd->ny*nfac+1;
112 snew(x_phi, nnx);
113 snew(y_psi, nny);
114 snew(newdata, nnx);
115 lo = 100000;
116 hi = -100000;
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);
125 if (i == 0)
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];
131 else*/
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 = gmx_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++)
143 sfree(newdata[i]);
145 sfree(newdata);
146 sfree(x_phi);
147 sfree(y_psi);
150 static t_shiftdata *read_shifts(const char *fn)
152 FILE *fp;
153 double xx;
154 int i, j, nx, ny;
155 t_shiftdata *sd;
157 snew(sd, 1);
158 fp = libopen(fn);
159 if (2 != fscanf(fp, "%d%d", &nx, &ny))
161 gmx_fatal(FARGS, "Error reading from file %s", fn);
164 sd->nx = nx;
165 sd->ny = ny;
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++)
174 if (i == nx)
176 sd->data[i][j] = sd->data[0][j];
178 else
180 if (1 != fscanf(fp, "%lf", &xx))
182 gmx_fatal(FARGS, "Error reading from file %s", fn);
184 sd->data[i][j] = xx;
187 sd->data[i][j] = sd->data[i][0];
189 gmx_ffclose(fp);
191 if (bDebugMode())
193 dump_sd(fn, sd);
196 return sd;
200 static void done_shifts(t_shiftdata *sd)
202 int i;
204 for (i = 0; (i <= sd->nx); i++)
206 sfree(sd->data[i]);
208 sfree(sd->data);
209 sfree(sd);
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;
215 int i, j, Phi, Psi;
216 real phi, psi;
217 real ca, co, ha, cb;
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++)
239 phi = dih[Phi][j];
240 psi = dih[Psi][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);
251 fprintf(fp, "\n");
253 /* Free memory */
254 done_shifts(ca_sd);
255 done_shifts(cb_sd);
256 done_shifts(co_sd);
257 done_shifts(ha_sd);