Simplify gmx chi
[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,2021, 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/math/utilities.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/futil.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/pleasecite.h"
52 #include "gromacs/utility/smalloc.h"
53
54 typedef struct
55 {
56     int    nx, ny;
57     real   dx, dy;
58     real** data;
59 } t_shiftdata;
60
61 static real interpolate(real phi, real psi, t_shiftdata* sd)
62 {
63     int  iphi, ipsi, iphi1, ipsi1;
64     real fphi, fpsi, wx0, wx1, wy0, wy1;
65
66     /*phi  += M_PI;
67        if (phi > 2*M_PI) phi -= 2*M_PI;
68        psi  += M_PI;
69        if (psi > 2*M_PI) psi -= 2*M_PI;
70      */
71     while (phi < 0)
72     {
73         phi += 2 * M_PI;
74     }
75     while (psi < 0)
76     {
77         psi += 2 * M_PI;
78     }
79     phi = 2 * M_PI - phi;
80
81     fphi = phi * sd->dx;
82     fpsi = psi * sd->dy;
83
84     iphi = static_cast<int>(fphi);
85     ipsi = static_cast<int>(fpsi);
86     fphi -= iphi; /* Fraction (offset from gridpoint) */
87     fpsi -= ipsi;
88
89     wx0   = 1.0 - fphi;
90     wx1   = fphi;
91     wy0   = 1.0 - fpsi;
92     wy1   = fpsi;
93     iphi  = iphi % sd->nx;
94     ipsi  = ipsi % sd->ny;
95     iphi1 = (iphi + 1) % sd->nx;
96     ipsi1 = (ipsi + 1) % sd->ny;
97
98     return (sd->data[iphi][ipsi] * wx0 * wy0 + sd->data[iphi1][ipsi] * wx1 * wy0
99             + sd->data[iphi][ipsi1] * wx0 * wy1 + sd->data[iphi1][ipsi1] * wx1 * wy1);
100 }
101
102 static void dump_sd(const char* fn, t_shiftdata* sd)
103 {
104     FILE* fp;
105     int   i, j;
106     char  buf[256];
107     int   nnx, nny, nfac = 4, nlevels = 20;
108     real  phi, psi, *x_phi, *y_psi, **newdata;
109     real  lo, hi;
110     t_rgb rlo = { 1, 0, 0 }, rhi = { 0, 0, 1 };
111
112     nnx = sd->nx * nfac + 1;
113     nny = sd->ny * nfac + 1;
114     snew(x_phi, nnx);
115     snew(y_psi, nny);
116     snew(newdata, nnx);
117     lo = 100000;
118     hi = -100000;
119     for (i = 0; (i < nnx); i++)
120     {
121         snew(newdata[i], nny);
122         phi      = i * 2 * M_PI / (nnx - 1);
123         x_phi[i] = phi * gmx::c_rad2Deg - 180;
124         for (j = 0; (j < nny); j++)
125         {
126             psi = j * 2 * M_PI / (nny - 1);
127             if (i == 0)
128             {
129                 y_psi[j] = psi * gmx::c_rad2Deg - 180;
130             }
131             /*if (((i % nfac) == 0) && ((j % nfac) == 0))
132                newdata[i][j] = sd->data[i/nfac][j/nfac];
133                else*/
134             newdata[i][j] = interpolate(phi, psi, sd);
135             lo            = std::min(lo, newdata[i][j]);
136             hi            = std::max(hi, newdata[i][j]);
137         }
138     }
139     sprintf(buf, "%s.xpm", fn);
140     fp = gmx_ffopen(buf, "w");
141     write_xpm(fp, 0, fn, fn, "Phi", "Psi", nnx, nny, 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     double       xx;
154     int          i, j, nx, ny;
155     t_shiftdata* sd;
156
157     snew(sd, 1);
158     gmx::FilePtr fp = gmx::openLibraryFile(fn);
159     if (2 != fscanf(fp.get(), "%d%d", &nx, &ny))
160     {
161         gmx_fatal(FARGS, "Error reading from file %s", fn);
162     }
163     GMX_ASSERT(nx > 0, "");
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++)
170     {
171         snew(sd->data[i], ny + 1);
172         for (j = 0; (j < ny); j++)
173         {
174             if (i == nx)
175             {
176                 sd->data[i][j] = sd->data[0][j];
177             }
178             else
179             {
180                 if (1 != fscanf(fp.get(), "%lf", &xx))
181                 {
182                     gmx_fatal(FARGS, "Error reading from file %s", fn);
183                 }
184                 sd->data[i][j] = xx;
185             }
186         }
187         sd->data[i][j] = sd->data[i][0];
188     }
189
190     if (bDebugMode())
191     {
192         dump_sd(fn, sd);
193     }
194
195     return sd;
196 }
197
198
199 static void done_shifts(t_shiftdata* sd)
200 {
201     int i;
202
203     for (i = 0; (i <= sd->nx); i++)
204     {
205         sfree(sd->data[i]);
206     }
207     sfree(sd->data);
208     sfree(sd);
209 }
210
211 void do_pp2shifts(FILE* fp, int nf, gmx::ArrayRef<const t_dlist> dlist, real** dih)
212 {
213     /* Read the shift files */
214     t_shiftdata* ca_sd = read_shifts("ca-shift.dat");
215     t_shiftdata* cb_sd = read_shifts("cb-shift.dat");
216     t_shiftdata* ha_sd = read_shifts("ha-shift.dat");
217     t_shiftdata* co_sd = read_shifts("co-shift.dat");
218
219     fprintf(fp, "\n *** Chemical shifts from the chemical shift index ***\n");
220     please_cite(fp, "Wishart98a");
221     fprintf(fp,
222             "%12s  %10s  %10s  %10s  %10s\n",
223             "Residue",
224             "delta Ca",
225             "delta Ha",
226             "delta CO",
227             "delta Cb");
228     for (const auto& dihedral : dlist)
229     {
230         if ((has_dihedral(edPhi, dihedral)) && (has_dihedral(edPsi, dihedral)))
231         {
232             int  Phi = dihedral.j0[edPhi];
233             int  Psi = dihedral.j0[edPsi];
234             real ca = 0, cb = 0, co = 0, ha = 0;
235             for (int j = 0; (j < nf); j++)
236             {
237                 real phi = dih[Phi][j];
238                 real psi = dih[Psi][j];
239
240                 ca += interpolate(phi, psi, ca_sd);
241                 cb += interpolate(phi, psi, cb_sd);
242                 co += interpolate(phi, psi, co_sd);
243                 ha += interpolate(phi, psi, ha_sd);
244             }
245             fprintf(fp, "%12s  %10g  %10g  %10g  %10g\n", dihedral.name, ca / nf, ha / nf, co / nf, 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 }