Moved copyrite.* to fileio from gmxlib and legacyheaders.
[gromacs.git] / src / gromacs / gmxana / pp2shift.cpp
blobdc6cbdcbbbc879562c501d63693425e8376d3304
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,2015, 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 "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
42 #include <algorithm>
44 #include "gromacs/fileio/copyrite.h"
45 #include "gromacs/fileio/matio.h"
46 #include "gromacs/gmxana/gstat.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/futil.h"
50 #include "gromacs/utility/smalloc.h"
52 typedef struct {
53 int nx, ny;
54 real dx, dy;
55 real **data;
56 } t_shiftdata;
58 static real interpolate(real phi, real psi, t_shiftdata *sd)
60 int iphi, ipsi, iphi1, ipsi1;
61 real fphi, fpsi, wx0, wx1, wy0, wy1;
63 /*phi += M_PI;
64 if (phi > 2*M_PI) phi -= 2*M_PI;
65 psi += M_PI;
66 if (psi > 2*M_PI) psi -= 2*M_PI;
68 while (phi < 0)
70 phi += 2*M_PI;
72 while (psi < 0)
74 psi += 2*M_PI;
76 phi = 2*M_PI-phi;
78 fphi = phi*sd->dx;
79 fpsi = psi*sd->dy;
81 iphi = static_cast<int>(fphi);
82 ipsi = static_cast<int>(fpsi);
83 fphi -= iphi; /* Fraction (offset from gridpoint) */
84 fpsi -= ipsi;
86 wx0 = 1.0-fphi;
87 wx1 = fphi;
88 wy0 = 1.0-fpsi;
89 wy1 = fpsi;
90 iphi = iphi % sd->nx;
91 ipsi = ipsi % sd->ny;
92 iphi1 = (iphi+1) % sd->nx;
93 ipsi1 = (ipsi+1) % sd->ny;
95 return (sd->data[iphi] [ipsi] * wx0*wy0 +
96 sd->data[iphi1] [ipsi] * wx1*wy0 +
97 sd->data[iphi] [ipsi1] * wx0*wy1 +
98 sd->data[iphi1] [ipsi1] * wx1*wy1);
101 static void dump_sd(const char *fn, t_shiftdata *sd)
103 FILE *fp;
104 int i, j;
105 char buf[256];
106 int nnx, nny, nfac = 4, nlevels = 20;
107 real phi, psi, *x_phi, *y_psi, **newdata;
108 real lo, hi;
109 t_rgb rlo = { 1, 0, 0 }, rhi = { 0, 0, 1 };
111 nnx = sd->nx*nfac+1;
112 nny = sd->ny*nfac+1;
113 snew(x_phi, nnx);
114 snew(y_psi, nny);
115 snew(newdata, nnx);
116 lo = 100000;
117 hi = -100000;
118 for (i = 0; (i < nnx); i++)
120 snew(newdata[i], nny);
121 phi = i*2*M_PI/(nnx-1);
122 x_phi[i] = phi*RAD2DEG-180;
123 for (j = 0; (j < nny); j++)
125 psi = j*2*M_PI/(nny-1);
126 if (i == 0)
128 y_psi[j] = psi*RAD2DEG-180;
130 /*if (((i % nfac) == 0) && ((j % nfac) == 0))
131 newdata[i][j] = sd->data[i/nfac][j/nfac];
132 else*/
133 newdata[i][j] = interpolate(phi, psi, sd);
134 lo = std::min(lo, newdata[i][j]);
135 hi = std::max(hi, newdata[i][j]);
138 sprintf(buf, "%s.xpm", fn);
139 fp = gmx_ffopen(buf, "w");
140 write_xpm(fp, 0, fn, fn, "Phi", "Psi", nnx, nny,
141 x_phi, y_psi, newdata, lo, hi, rlo, rhi, &nlevels);
142 for (i = 0; (i < nnx); i++)
144 sfree(newdata[i]);
146 sfree(newdata);
147 sfree(x_phi);
148 sfree(y_psi);
151 static t_shiftdata *read_shifts(const char *fn)
153 FILE *fp;
154 double xx;
155 int i, j, nx, ny;
156 t_shiftdata *sd;
158 snew(sd, 1);
159 fp = libopen(fn);
160 if (2 != fscanf(fp, "%d%d", &nx, &ny))
162 gmx_fatal(FARGS, "Error reading from file %s", fn);
165 sd->nx = nx;
166 sd->ny = ny;
167 sd->dx = nx/(2*M_PI);
168 sd->dy = ny/(2*M_PI);
169 snew(sd->data, nx+1);
170 for (i = 0; (i <= nx); i++)
172 snew(sd->data[i], ny+1);
173 for (j = 0; (j < ny); j++)
175 if (i == nx)
177 sd->data[i][j] = sd->data[0][j];
179 else
181 if (1 != fscanf(fp, "%lf", &xx))
183 gmx_fatal(FARGS, "Error reading from file %s", fn);
185 sd->data[i][j] = xx;
188 sd->data[i][j] = sd->data[i][0];
190 gmx_ffclose(fp);
192 if (bDebugMode())
194 dump_sd(fn, sd);
197 return sd;
201 static void done_shifts(t_shiftdata *sd)
203 int i;
205 for (i = 0; (i <= sd->nx); i++)
207 sfree(sd->data[i]);
209 sfree(sd->data);
210 sfree(sd);
213 void do_pp2shifts(FILE *fp, int nf, int nlist, t_dlist dlist[], real **dih)
215 t_shiftdata *ca_sd, *co_sd, *ha_sd, *cb_sd;
216 int i, j, Phi, Psi;
217 real phi, psi;
218 real ca, co, ha, cb;
220 /* Read the shift files */
221 ca_sd = read_shifts("ca-shift.dat");
222 cb_sd = read_shifts("cb-shift.dat");
223 ha_sd = read_shifts("ha-shift.dat");
224 co_sd = read_shifts("co-shift.dat");
226 fprintf(fp, "\n *** Chemical shifts from the chemical shift index ***\n");
227 please_cite(fp, "Wishart98a");
228 fprintf(fp, "%12s %10s %10s %10s %10s\n",
229 "Residue", "delta Ca", "delta Ha", "delta CO", "delta Cb");
230 for (i = 0; (i < nlist); i++)
232 if ((has_dihedral(edPhi, &(dlist[i]))) &&
233 (has_dihedral(edPsi, &(dlist[i]))))
235 Phi = dlist[i].j0[edPhi];
236 Psi = dlist[i].j0[edPsi];
237 ca = cb = co = ha = 0;
238 for (j = 0; (j < nf); j++)
240 phi = dih[Phi][j];
241 psi = dih[Psi][j];
243 ca += interpolate(phi, psi, ca_sd);
244 cb += interpolate(phi, psi, cb_sd);
245 co += interpolate(phi, psi, co_sd);
246 ha += interpolate(phi, psi, ha_sd);
248 fprintf(fp, "%12s %10g %10g %10g %10g\n",
249 dlist[i].name, ca/nf, ha/nf, co/nf, cb/nf);
252 fprintf(fp, "\n");
254 /* Free memory */
255 done_shifts(ca_sd);
256 done_shifts(cb_sd);
257 done_shifts(co_sd);
258 done_shifts(ha_sd);