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