Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / pp2shift.cpp
blobcbb4a0d58a2ba9dcb260bfa992a8432afda264bc
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
42 #include <algorithm>
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"
53 typedef struct
55 int nx, ny;
56 real dx, dy;
57 real** data;
58 } t_shiftdata;
60 static real interpolate(real phi, real psi, t_shiftdata* sd)
62 int iphi, ipsi, iphi1, ipsi1;
63 real fphi, fpsi, wx0, wx1, wy0, wy1;
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;
70 while (phi < 0)
72 phi += 2 * M_PI;
74 while (psi < 0)
76 psi += 2 * M_PI;
78 phi = 2 * M_PI - phi;
80 fphi = phi * sd->dx;
81 fpsi = psi * sd->dy;
83 iphi = static_cast<int>(fphi);
84 ipsi = static_cast<int>(fpsi);
85 fphi -= iphi; /* Fraction (offset from gridpoint) */
86 fpsi -= ipsi;
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;
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);
101 static void dump_sd(const char* fn, t_shiftdata* sd)
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 };
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++)
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++)
125 psi = j * 2 * M_PI / (nny - 1);
126 if (i == 0)
128 y_psi[j] = psi * RAD2DEG - 180;
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]);
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++)
143 sfree(newdata[i]);
145 sfree(newdata);
146 sfree(x_phi);
147 sfree(y_psi);
150 static t_shiftdata* read_shifts(const char* fn)
152 double xx;
153 int i, j, nx, ny;
154 t_shiftdata* sd;
156 snew(sd, 1);
157 gmx::FilePtr fp = gmx::openLibraryFile(fn);
158 if (2 != fscanf(fp.get(), "%d%d", &nx, &ny))
160 gmx_fatal(FARGS, "Error reading from file %s", fn);
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++)
170 snew(sd->data[i], ny + 1);
171 for (j = 0; (j < ny); j++)
173 if (i == nx)
175 sd->data[i][j] = sd->data[0][j];
177 else
179 if (1 != fscanf(fp.get(), "%lf", &xx))
181 gmx_fatal(FARGS, "Error reading from file %s", fn);
183 sd->data[i][j] = xx;
186 sd->data[i][j] = sd->data[i][0];
189 if (bDebugMode())
191 dump_sd(fn, sd);
194 return sd;
198 static void done_shifts(t_shiftdata* sd)
200 int i;
202 for (i = 0; (i <= sd->nx); i++)
204 sfree(sd->data[i]);
206 sfree(sd->data);
207 sfree(sd);
210 void do_pp2shifts(FILE* fp, int nf, int nlist, t_dlist dlist[], real** dih)
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;
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");
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++)
229 if ((has_dihedral(edPhi, &(dlist[i]))) && (has_dihedral(edPsi, &(dlist[i]))))
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++)
236 phi = dih[Phi][j];
237 psi = dih[Psi][j];
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);
244 fprintf(fp, "%12s %10g %10g %10g %10g\n", dlist[i].name, ca / nf, ha / nf, co / nf,
245 cb / nf);
248 fprintf(fp, "\n");
250 /* Free memory */
251 done_shifts(ca_sd);
252 done_shifts(cb_sd);
253 done_shifts(co_sd);
254 done_shifts(ha_sd);