cbb4a0d58a2ba9dcb260bfa992a8432afda264bc
[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,2018,2019,2020, 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/math/units.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/futil.h"
49 #include "gromacs/utility/gmxassert.h"
50 #include "gromacs/utility/pleasecite.h"
51 #include "gromacs/utility/smalloc.h"
52
53 typedef struct
54 {
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 + sd->data[iphi1][ipsi] * wx1 * wy0
98             + sd->data[iphi][ipsi1] * wx0 * wy1 + 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            = std::min(lo, newdata[i][j]);
135             hi            = std::max(hi, newdata[i][j]);
136         }
137     }
138     sprintf(buf, "%s.xpm", fn);
139     fp = gmx_ffopen(buf, "w");
140     write_xpm(fp, 0, fn, fn, "Phi", "Psi", nnx, nny, x_phi, y_psi, newdata, lo, hi, rlo, rhi, &nlevels);
141     for (i = 0; (i < nnx); i++)
142     {
143         sfree(newdata[i]);
144     }
145     sfree(newdata);
146     sfree(x_phi);
147     sfree(y_psi);
148 }
149
150 static t_shiftdata* read_shifts(const char* fn)
151 {
152     double       xx;
153     int          i, j, nx, ny;
154     t_shiftdata* sd;
155
156     snew(sd, 1);
157     gmx::FilePtr fp = gmx::openLibraryFile(fn);
158     if (2 != fscanf(fp.get(), "%d%d", &nx, &ny))
159     {
160         gmx_fatal(FARGS, "Error reading from file %s", fn);
161     }
162     GMX_ASSERT(nx > 0, "");
163     sd->nx = nx;
164     sd->ny = ny;
165     sd->dx = nx / (2 * M_PI);
166     sd->dy = ny / (2 * M_PI);
167     snew(sd->data, nx + 1);
168     for (i = 0; (i <= nx); i++)
169     {
170         snew(sd->data[i], ny + 1);
171         for (j = 0; (j < ny); j++)
172         {
173             if (i == nx)
174             {
175                 sd->data[i][j] = sd->data[0][j];
176             }
177             else
178             {
179                 if (1 != fscanf(fp.get(), "%lf", &xx))
180                 {
181                     gmx_fatal(FARGS, "Error reading from file %s", fn);
182                 }
183                 sd->data[i][j] = xx;
184             }
185         }
186         sd->data[i][j] = sd->data[i][0];
187     }
188
189     if (bDebugMode())
190     {
191         dump_sd(fn, sd);
192     }
193
194     return sd;
195 }
196
197
198 static void done_shifts(t_shiftdata* sd)
199 {
200     int i;
201
202     for (i = 0; (i <= sd->nx); i++)
203     {
204         sfree(sd->data[i]);
205     }
206     sfree(sd->data);
207     sfree(sd);
208 }
209
210 void do_pp2shifts(FILE* fp, int nf, int nlist, t_dlist dlist[], real** dih)
211 {
212     t_shiftdata *ca_sd, *co_sd, *ha_sd, *cb_sd;
213     int          i, j, Phi, Psi;
214     real         phi, psi;
215     real         ca, co, ha, cb;
216
217     /* Read the shift files */
218     ca_sd = read_shifts("ca-shift.dat");
219     cb_sd = read_shifts("cb-shift.dat");
220     ha_sd = read_shifts("ha-shift.dat");
221     co_sd = read_shifts("co-shift.dat");
222
223     fprintf(fp, "\n *** Chemical shifts from the chemical shift index ***\n");
224     please_cite(fp, "Wishart98a");
225     fprintf(fp, "%12s  %10s  %10s  %10s  %10s\n", "Residue", "delta Ca", "delta Ha", "delta CO",
226             "delta Cb");
227     for (i = 0; (i < nlist); i++)
228     {
229         if ((has_dihedral(edPhi, &(dlist[i]))) && (has_dihedral(edPsi, &(dlist[i]))))
230         {
231             Phi = dlist[i].j0[edPhi];
232             Psi = dlist[i].j0[edPsi];
233             ca = cb = co = ha = 0;
234             for (j = 0; (j < nf); j++)
235             {
236                 phi = dih[Phi][j];
237                 psi = dih[Psi][j];
238
239                 ca += interpolate(phi, psi, ca_sd);
240                 cb += interpolate(phi, psi, cb_sd);
241                 co += interpolate(phi, psi, co_sd);
242                 ha += interpolate(phi, psi, ha_sd);
243             }
244             fprintf(fp, "%12s  %10g  %10g  %10g  %10g\n", dlist[i].name, ca / nf, ha / nf, co / nf,
245                     cb / nf);
246         }
247     }
248     fprintf(fp, "\n");
249
250     /* Free memory */
251     done_shifts(ca_sd);
252     done_shifts(cb_sd);
253     done_shifts(co_sd);
254     done_shifts(ha_sd);
255 }