Merge "Fix bug in selection subexpression handling."
[alexxy/gromacs.git] / src / gromacs / gmxana / pp2shift.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <stdlib.h>
40 #include <math.h>
41 #include "typedefs.h"
42 #include "futil.h"
43 #include "macros.h"
44 #include "physics.h"
45 #include "smalloc.h"
46 #include "statutil.h"
47 #include "gstat.h"
48 #include "matio.h"
49 #include "copyrite.h"
50 #include "gmx_fatal.h"
51
52 typedef struct {
53     int    nx, ny;
54     real   dx, dy;
55     real **data;
56 } t_shiftdata;
57
58 static real interpolate(real phi, real psi, t_shiftdata *sd)
59 {
60     int  iphi, ipsi, iphi1, ipsi1;
61     real fphi, fpsi, wx0, wx1, wy0, wy1;
62
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;
67      */
68     while (phi < 0)
69     {
70         phi += 2*M_PI;
71     }
72     while (psi < 0)
73     {
74         psi += 2*M_PI;
75     }
76     phi    = 2*M_PI-phi;
77
78     fphi  = phi*sd->dx;
79     fpsi  = psi*sd->dy;
80
81     iphi  = (int)fphi;
82     ipsi  = (int)fpsi;
83     fphi -= iphi; /* Fraction (offset from gridpoint) */
84     fpsi -= ipsi;
85
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;
94
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);
99 }
100
101 static void dump_sd(const char *fn, t_shiftdata *sd)
102 {
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 };
110
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++)
119     {
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++)
124         {
125             psi = j*2*M_PI/(nny-1);
126             if (i == 0)
127             {
128                 y_psi[j] = psi*RAD2DEG-180;
129             }
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            = min(lo, newdata[i][j]);
135             hi            = max(hi, newdata[i][j]);
136         }
137     }
138     sprintf(buf, "%s.xpm", fn);
139     fp = 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++)
143     {
144         sfree(newdata[i]);
145     }
146     sfree(newdata);
147     sfree(x_phi);
148     sfree(y_psi);
149 }
150
151 static t_shiftdata *read_shifts(const char *fn)
152 {
153     FILE        *fp;
154     double       xx;
155     int          i, j, nx, ny;
156     t_shiftdata *sd;
157
158     snew(sd, 1);
159     fp = libopen(fn);
160     if (2 != fscanf(fp, "%d%d", &nx, &ny))
161     {
162         gmx_fatal(FARGS, "Error reading from file %s", fn);
163     }
164
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++)
171     {
172         snew(sd->data[i], ny+1);
173         for (j = 0; (j < ny); j++)
174         {
175             if (i == nx)
176             {
177                 sd->data[i][j] = sd->data[0][j];
178             }
179             else
180             {
181                 if (1 != fscanf(fp, "%lf", &xx))
182                 {
183                     gmx_fatal(FARGS, "Error reading from file %s", fn);
184                 }
185                 sd->data[i][j] = xx;
186             }
187         }
188         sd->data[i][j] = sd->data[i][0];
189     }
190     ffclose(fp);
191
192     if (bDebugMode())
193     {
194         dump_sd(fn, sd);
195     }
196
197     return sd;
198 }
199
200
201 static void done_shifts(t_shiftdata *sd)
202 {
203     int i;
204
205     for (i = 0; (i <= sd->nx); i++)
206     {
207         sfree(sd->data[i]);
208     }
209     sfree(sd->data);
210     sfree(sd);
211 }
212
213 void do_pp2shifts(FILE *fp, int nf, int nlist, t_dlist dlist[], real **dih)
214 {
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;
219
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");
225
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++)
231     {
232         if ((has_dihedral(edPhi, &(dlist[i]))) &&
233             (has_dihedral(edPsi, &(dlist[i]))))
234         {
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++)
239             {
240                 phi = dih[Phi][j];
241                 psi = dih[Psi][j];
242
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);
247             }
248             fprintf(fp, "%12s  %10g  %10g  %10g  %10g\n",
249                     dlist[i].name, ca/nf, ha/nf, co/nf, cb/nf);
250         }
251     }
252     fprintf(fp, "\n");
253
254     /* Free memory */
255     done_shifts(ca_sd);
256     done_shifts(cb_sd);
257     done_shifts(co_sd);
258     done_shifts(ha_sd);
259 }