Moved first batch of analysis tool source to C++
[alexxy/gromacs.git] / src / gromacs / gmxana / pp2shift.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include <cmath>
40 #include <cstdlib>
41
42 #include <algorithm>
43
44 #include "gromacs/fileio/matio.h"
45 #include "gromacs/gmxana/gstat.h"
46 #include "gromacs/legacyheaders/copyrite.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/legacyheaders/typedefs.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/futil.h"
52 #include "gromacs/utility/smalloc.h"
53
54 typedef struct {
55     int    nx, ny;
56     real   dx, dy;
57     real **data;
58 } t_shiftdata;
59
60 static real interpolate(real phi, real psi, t_shiftdata *sd)
61 {
62     int  iphi, ipsi, iphi1, ipsi1;
63     real fphi, fpsi, wx0, wx1, wy0, wy1;
64
65     /*phi  += M_PI;
66        if (phi > 2*M_PI) phi -= 2*M_PI;
67        psi  += M_PI;
68        if (psi > 2*M_PI) psi -= 2*M_PI;
69      */
70     while (phi < 0)
71     {
72         phi += 2*M_PI;
73     }
74     while (psi < 0)
75     {
76         psi += 2*M_PI;
77     }
78     phi    = 2*M_PI-phi;
79
80     fphi  = phi*sd->dx;
81     fpsi  = psi*sd->dy;
82
83     iphi  = static_cast<int>(fphi);
84     ipsi  = static_cast<int>(fpsi);
85     fphi -= iphi; /* Fraction (offset from gridpoint) */
86     fpsi -= ipsi;
87
88     wx0   = 1.0-fphi;
89     wx1   = fphi;
90     wy0   = 1.0-fpsi;
91     wy1   = fpsi;
92     iphi  = iphi % sd->nx;
93     ipsi  = ipsi % sd->ny;
94     iphi1 = (iphi+1) % sd->nx;
95     ipsi1 = (ipsi+1) % sd->ny;
96
97     return (sd->data[iphi]  [ipsi]  * wx0*wy0 +
98             sd->data[iphi1] [ipsi]  * wx1*wy0 +
99             sd->data[iphi]  [ipsi1] * wx0*wy1 +
100             sd->data[iphi1] [ipsi1] * wx1*wy1);
101 }
102
103 static void dump_sd(const char *fn, t_shiftdata *sd)
104 {
105     FILE  *fp;
106     int    i, j;
107     char   buf[256];
108     int    nnx, nny, nfac = 4, nlevels = 20;
109     real   phi, psi, *x_phi, *y_psi, **newdata;
110     real   lo, hi;
111     t_rgb  rlo = { 1, 0, 0 }, rhi = { 0, 0, 1 };
112
113     nnx = sd->nx*nfac+1;
114     nny = sd->ny*nfac+1;
115     snew(x_phi, nnx);
116     snew(y_psi, nny);
117     snew(newdata, nnx);
118     lo =  100000;
119     hi = -100000;
120     for (i = 0; (i < nnx); i++)
121     {
122         snew(newdata[i], nny);
123         phi      = i*2*M_PI/(nnx-1);
124         x_phi[i] = phi*RAD2DEG-180;
125         for (j = 0; (j < nny); j++)
126         {
127             psi = j*2*M_PI/(nny-1);
128             if (i == 0)
129             {
130                 y_psi[j] = psi*RAD2DEG-180;
131             }
132             /*if (((i % nfac) == 0) && ((j % nfac) == 0))
133                newdata[i][j] = sd->data[i/nfac][j/nfac];
134                else*/
135             newdata[i][j] = interpolate(phi, psi, sd);
136             lo            = std::min(lo, newdata[i][j]);
137             hi            = std::max(hi, newdata[i][j]);
138         }
139     }
140     sprintf(buf, "%s.xpm", fn);
141     fp = gmx_ffopen(buf, "w");
142     write_xpm(fp, 0, fn, fn, "Phi", "Psi", nnx, nny,
143               x_phi, y_psi, newdata, lo, hi, rlo, rhi, &nlevels);
144     for (i = 0; (i < nnx); i++)
145     {
146         sfree(newdata[i]);
147     }
148     sfree(newdata);
149     sfree(x_phi);
150     sfree(y_psi);
151 }
152
153 static t_shiftdata *read_shifts(const char *fn)
154 {
155     FILE        *fp;
156     double       xx;
157     int          i, j, nx, ny;
158     t_shiftdata *sd;
159
160     snew(sd, 1);
161     fp = libopen(fn);
162     if (2 != fscanf(fp, "%d%d", &nx, &ny))
163     {
164         gmx_fatal(FARGS, "Error reading from file %s", fn);
165     }
166
167     sd->nx = nx;
168     sd->ny = ny;
169     sd->dx = nx/(2*M_PI);
170     sd->dy = ny/(2*M_PI);
171     snew(sd->data, nx+1);
172     for (i = 0; (i <= nx); i++)
173     {
174         snew(sd->data[i], ny+1);
175         for (j = 0; (j < ny); j++)
176         {
177             if (i == nx)
178             {
179                 sd->data[i][j] = sd->data[0][j];
180             }
181             else
182             {
183                 if (1 != fscanf(fp, "%lf", &xx))
184                 {
185                     gmx_fatal(FARGS, "Error reading from file %s", fn);
186                 }
187                 sd->data[i][j] = xx;
188             }
189         }
190         sd->data[i][j] = sd->data[i][0];
191     }
192     gmx_ffclose(fp);
193
194     if (bDebugMode())
195     {
196         dump_sd(fn, sd);
197     }
198
199     return sd;
200 }
201
202
203 static void done_shifts(t_shiftdata *sd)
204 {
205     int i;
206
207     for (i = 0; (i <= sd->nx); i++)
208     {
209         sfree(sd->data[i]);
210     }
211     sfree(sd->data);
212     sfree(sd);
213 }
214
215 void do_pp2shifts(FILE *fp, int nf, int nlist, t_dlist dlist[], real **dih)
216 {
217     t_shiftdata *ca_sd, *co_sd, *ha_sd, *cb_sd;
218     int          i, j, Phi, Psi;
219     real         phi, psi;
220     real         ca, co, ha, cb;
221
222     /* Read the shift files */
223     ca_sd = read_shifts("ca-shift.dat");
224     cb_sd = read_shifts("cb-shift.dat");
225     ha_sd = read_shifts("ha-shift.dat");
226     co_sd = read_shifts("co-shift.dat");
227
228     fprintf(fp, "\n *** Chemical shifts from the chemical shift index ***\n");
229     please_cite(fp, "Wishart98a");
230     fprintf(fp, "%12s  %10s  %10s  %10s  %10s\n",
231             "Residue", "delta Ca", "delta Ha", "delta CO", "delta Cb");
232     for (i = 0; (i < nlist); i++)
233     {
234         if ((has_dihedral(edPhi, &(dlist[i]))) &&
235             (has_dihedral(edPsi, &(dlist[i]))))
236         {
237             Phi = dlist[i].j0[edPhi];
238             Psi = dlist[i].j0[edPsi];
239             ca  = cb = co = ha = 0;
240             for (j = 0; (j < nf); j++)
241             {
242                 phi = dih[Phi][j];
243                 psi = dih[Psi][j];
244
245                 ca += interpolate(phi, psi, ca_sd);
246                 cb += interpolate(phi, psi, cb_sd);
247                 co += interpolate(phi, psi, co_sd);
248                 ha += interpolate(phi, psi, ha_sd);
249             }
250             fprintf(fp, "%12s  %10g  %10g  %10g  %10g\n",
251                     dlist[i].name, ca/nf, ha/nf, co/nf, cb/nf);
252         }
253     }
254     fprintf(fp, "\n");
255
256     /* Free memory */
257     done_shifts(ca_sd);
258     done_shifts(cb_sd);
259     done_shifts(co_sd);
260     done_shifts(ha_sd);
261 }