Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / hxprops.cpp
blob59206dd2ea61e1f5169738804c72cbcd72d4a1a1
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,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "gmxpre.h"
40 #include "hxprops.h"
42 #include <cmath>
43 #include <cstring>
45 #include <algorithm>
47 #include "gromacs/listed_forces/bonded.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/topology/index.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/arraysize.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/smalloc.h"
57 real ellipticity(int nres, t_bb bb[])
59 typedef struct
61 real phi, psi, w;
62 } t_ppwstr;
63 // Avoid warnings about narrowing conversions from double to real
64 #ifdef _MSC_VER
65 # pragma warning(disable : 4838)
66 #endif
67 static const t_ppwstr ppw[] = { { -67, -44, 0.31 }, { -66, -41, 0.31 }, { -59, -44, 0.44 },
68 { -57, -47, 0.56 }, { -53, -52, 0.78 }, { -48, -57, 1.00 },
69 { -70.5, -35.8, 0.15 }, { -57, -79, 0.23 }, { -38, -78, 1.20 },
70 { -60, -30, 0.24 }, { -54, -28, 0.46 }, { -44, -33, 0.68 } };
71 #ifdef _MSC_VER
72 # pragma warning(default : 4838)
73 #endif
74 #define NPPW asize(ppw)
76 int i, j;
77 real ell, pp2, phi, psi;
79 ell = 0;
80 for (i = 0; (i < nres); i++)
82 phi = bb[i].phi;
83 psi = bb[i].psi;
84 for (j = 0; (j < NPPW); j++)
86 pp2 = gmx::square(phi - ppw[j].phi) + gmx::square(psi - ppw[j].psi);
87 if (pp2 < 64)
89 bb[i].nhx++;
90 ell += ppw[j].w;
91 break;
95 return ell;
98 real ahx_len(int gnx, const int index[], rvec x[])
99 /* Assume we have a list of Calpha atoms only! */
101 rvec dx;
103 rvec_sub(x[index[0]], x[index[gnx - 1]], dx);
105 return norm(dx);
108 real radius(FILE* fp, int nca, const int ca_index[], rvec x[])
109 /* Assume we have all the backbone */
111 real dl2, dlt;
112 int i, ai;
114 dlt = 0;
115 for (i = 0; (i < nca); i++)
117 ai = ca_index[i];
118 dl2 = gmx::square(x[ai][XX]) + gmx::square(x[ai][YY]);
120 if (fp)
122 fprintf(fp, " %10g", dl2);
125 dlt += dl2;
127 if (fp)
129 fprintf(fp, "\n");
132 return std::sqrt(dlt / nca);
135 static real rot(rvec x1, const rvec x2)
137 real phi1, dphi, cp, sp;
138 real xx, yy;
140 phi1 = std::atan2(x1[YY], x1[XX]);
141 cp = std::cos(phi1);
142 sp = std::sin(phi1);
143 xx = cp * x2[XX] + sp * x2[YY];
144 yy = -sp * x2[XX] + cp * x2[YY];
146 dphi = RAD2DEG * std::atan2(yy, xx);
148 return dphi;
151 real twist(int nca, const int caindex[], rvec x[])
153 real pt, dphi;
154 int i, a0, a1;
156 pt = 0;
157 a0 = caindex[0];
158 for (i = 1; (i < nca); i++)
160 a1 = caindex[i];
162 dphi = rot(x[a0], x[a1]);
163 if (dphi < -90)
165 dphi += 360;
167 pt += dphi;
168 a0 = a1;
171 return (pt / (nca - 1));
174 real ca_phi(int gnx, const int index[], rvec x[])
175 /* Assume we have a list of Calpha atoms only! */
177 real phi, phitot;
178 int i, ai, aj, ak, al, t1, t2, t3;
179 rvec r_ij, r_kj, r_kl, m, n;
181 if (gnx <= 4)
183 return 0;
186 phitot = 0;
187 for (i = 0; (i < gnx - 4); i++)
189 ai = index[i + 0];
190 aj = index[i + 1];
191 ak = index[i + 2];
192 al = index[i + 3];
193 phi = RAD2DEG
194 * dih_angle(x[ai], x[aj], x[ak], x[al], nullptr, r_ij, r_kj, r_kl, m, n, &t1, &t2, &t3);
195 phitot += phi;
198 return (phitot / (gnx - 4.0));
201 real dip(int nbb, int const bbind[], const rvec x[], const t_atom atom[])
203 int i, m, ai;
204 rvec dipje;
205 real q;
207 clear_rvec(dipje);
208 for (i = 0; (i < nbb); i++)
210 ai = bbind[i];
211 q = atom[ai].q;
212 for (m = 0; (m < DIM); m++)
214 dipje[m] += x[ai][m] * q;
217 return norm(dipje);
220 real rise(int gnx, const int index[], rvec x[])
221 /* Assume we have a list of Calpha atoms only! */
223 real z, z0, ztot;
224 int i, ai;
226 ai = index[0];
227 z0 = x[ai][ZZ];
228 ztot = 0;
229 for (i = 1; (i < gnx); i++)
231 ai = index[i];
232 z = x[ai][ZZ];
233 ztot += (z - z0);
234 z0 = z;
237 return (ztot / (gnx - 1.0));
240 void av_hblen(FILE* fp3, FILE* fp3a, FILE* fp4, FILE* fp4a, FILE* fp5, FILE* fp5a, real t, int nres, t_bb bb[])
242 int i, n3 = 0, n4 = 0, n5 = 0;
243 real d3 = 0, d4 = 0, d5 = 0;
245 for (i = 0; (i < nres - 3); i++)
247 if (bb[i].bHelix)
249 fprintf(fp3a, "%10g", bb[i].d3);
250 n3++;
251 d3 += bb[i].d3;
252 if (i < nres - 4)
254 fprintf(fp4a, "%10g", bb[i].d4);
255 n4++;
256 d4 += bb[i].d4;
258 if (i < nres - 5)
260 fprintf(fp5a, "%10g", bb[i].d5);
261 n5++;
262 d5 += bb[i].d5;
266 fprintf(fp3, "%10g %10g\n", t, d3 / n3);
267 fprintf(fp4, "%10g %10g\n", t, d4 / n4);
268 fprintf(fp5, "%10g %10g\n", t, d5 / n5);
269 fprintf(fp3a, "\n");
270 fprintf(fp4a, "\n");
271 fprintf(fp5a, "\n");
274 void av_phipsi(FILE* fphi, FILE* fpsi, FILE* fphi2, FILE* fpsi2, real t, int nres, t_bb bb[])
276 int i, n = 0;
277 real phi = 0, psi = 0;
279 fprintf(fphi2, "%10g", t);
280 fprintf(fpsi2, "%10g", t);
281 for (i = 0; (i < nres); i++)
283 if (bb[i].bHelix)
285 phi += bb[i].phi;
286 psi += bb[i].psi;
287 fprintf(fphi2, " %10g", bb[i].phi);
288 fprintf(fpsi2, " %10g", bb[i].psi);
289 n++;
292 fprintf(fphi, "%10g %10g\n", t, (phi / n));
293 fprintf(fpsi, "%10g %10g\n", t, (psi / n));
294 fprintf(fphi2, "\n");
295 fprintf(fpsi2, "\n");
298 static void set_ahcity(int nbb, t_bb bb[])
300 real pp2;
301 int n;
303 for (n = 0; (n < nbb); n++)
305 pp2 = gmx::square(bb[n].phi - PHI_AHX) + gmx::square(bb[n].psi - PSI_AHX);
307 bb[n].bHelix = FALSE;
308 if (pp2 < 2500)
310 if ((bb[n].d4 < 0.36) || ((n > 0) && bb[n - 1].bHelix))
312 bb[n].bHelix = TRUE;
318 t_bb* mkbbind(const char* fn,
319 int* nres,
320 int* nbb,
321 int res0,
322 int* nall,
323 int** index,
324 char*** atomname,
325 t_atom atom[],
326 t_resinfo* resinfo)
328 static const char* bb_nm[] = { "N", "H", "CA", "C", "O", "HN" };
329 #define NBB asize(bb_nm)
330 t_bb* bb;
331 char* grpname;
332 int ai, i, i0, i1, j, k, rnr, gnx, r0, r1;
334 fprintf(stderr, "Please select a group containing the entire backbone\n");
335 rd_index(fn, 1, &gnx, index, &grpname);
336 *nall = gnx;
337 fprintf(stderr, "Checking group %s\n", grpname);
338 r0 = r1 = atom[(*index)[0]].resind;
339 for (i = 1; (i < gnx); i++)
341 r0 = std::min(r0, atom[(*index)[i]].resind);
342 r1 = std::max(r1, atom[(*index)[i]].resind);
344 rnr = r1 - r0 + 1;
345 fprintf(stderr, "There are %d residues\n", rnr);
346 snew(bb, rnr);
347 for (i = 0; (i < rnr); i++)
349 bb[i].N = bb[i].H = bb[i].CA = bb[i].C = bb[i].O = -1;
350 bb[i].resno = res0 + i;
353 for (i = j = 0; (i < gnx); i++)
355 ai = (*index)[i];
356 // Create an index into the residue index for the topology.
357 int resindex = atom[ai].resind;
358 // Create an index into the residues present in the selected
359 // index group.
360 int bbindex = resindex - r0;
361 if (std::strcmp(*(resinfo[resindex].name), "PRO") == 0)
363 // For PRO in a peptide, there is no H bound to backbone
364 // N, so use CD instead.
365 if (std::strcmp(*(atomname[ai]), "CD") == 0)
367 bb[bbindex].H = ai;
370 for (k = 0; (k < NBB); k++)
372 if (std::strcmp(bb_nm[k], *(atomname[ai])) == 0)
374 break;
377 switch (k)
379 case 0: bb[bbindex].N = ai; break;
380 case 1:
381 case 5:
382 /* No attempt to address the case where some weird input has both H and HN atoms in the group */
383 bb[bbindex].H = ai;
384 break;
385 case 2: bb[bbindex].CA = ai; break;
386 case 3: bb[bbindex].C = ai; break;
387 case 4: bb[bbindex].O = ai; break;
388 default: break;
392 for (i0 = 0; (i0 < rnr); i0++)
394 if ((bb[i0].N != -1) && (bb[i0].H != -1) && (bb[i0].CA != -1) && (bb[i0].C != -1)
395 && (bb[i0].O != -1))
397 break;
400 for (i1 = rnr - 1; (i1 >= 0); i1--)
402 if ((bb[i1].N != -1) && (bb[i1].H != -1) && (bb[i1].CA != -1) && (bb[i1].C != -1)
403 && (bb[i1].O != -1))
405 break;
408 if (i0 == 0)
410 i0++;
412 if (i1 == rnr - 1)
414 i1--;
417 for (i = i0; (i < i1); i++)
419 bb[i].Cprev = bb[i - 1].C;
420 bb[i].Nnext = bb[i + 1].N;
422 rnr = std::max(0, i1 - i0 + 1);
423 fprintf(stderr, "There are %d complete backbone residues (from %d to %d)\n", rnr, bb[i0].resno,
424 bb[i1].resno);
425 if (rnr == 0)
427 gmx_fatal(FARGS, "Zero complete backbone residues were found, cannot proceed");
429 for (i = 0; (i < rnr); i++, i0++)
431 bb[i] = bb[i0];
434 /* Set the labels */
435 for (i = 0; (i < rnr); i++)
437 int resindex = atom[bb[i].CA].resind;
438 sprintf(bb[i].label, "%s%d", *(resinfo[resindex].name), resinfo[resindex].nr);
441 *nres = rnr;
442 *nbb = rnr * asize(bb_nm);
444 return bb;
447 real pprms(FILE* fp, int nbb, t_bb bb[])
449 int i, n;
450 real rms, rmst, rms2;
452 rmst = rms2 = 0;
453 for (i = n = 0; (i < nbb); i++)
455 if (bb[i].bHelix)
457 rms = std::sqrt(bb[i].pprms2);
458 rmst += rms;
459 rms2 += bb[i].pprms2;
460 fprintf(fp, "%10g ", rms);
461 n++;
464 fprintf(fp, "\n");
465 rms = std::sqrt(rms2 / n - gmx::square(rmst / n));
467 return rms;
470 void calc_hxprops(int nres, t_bb bb[], const rvec x[])
472 int i, ao, an, t1, t2, t3;
473 rvec dx, r_ij, r_kj, r_kl, m, n;
475 for (i = 0; (i < nres); i++)
477 ao = bb[i].O;
478 bb[i].d4 = bb[i].d3 = bb[i].d5 = 0;
479 if (i < nres - 3)
481 an = bb[i + 3].N;
482 rvec_sub(x[ao], x[an], dx);
483 bb[i].d3 = norm(dx);
485 if (i < nres - 4)
487 an = bb[i + 4].N;
488 rvec_sub(x[ao], x[an], dx);
489 bb[i].d4 = norm(dx);
491 if (i < nres - 5)
493 an = bb[i + 5].N;
494 rvec_sub(x[ao], x[an], dx);
495 bb[i].d5 = norm(dx);
498 bb[i].phi = RAD2DEG
499 * dih_angle(x[bb[i].Cprev], x[bb[i].N], x[bb[i].CA], x[bb[i].C], nullptr, r_ij,
500 r_kj, r_kl, m, n, &t1, &t2, &t3);
501 bb[i].psi = RAD2DEG
502 * dih_angle(x[bb[i].N], x[bb[i].CA], x[bb[i].C], x[bb[i].Nnext], nullptr, r_ij,
503 r_kj, r_kl, m, n, &t1, &t2, &t3);
504 bb[i].pprms2 = gmx::square(bb[i].phi - PHI_AHX) + gmx::square(bb[i].psi - PSI_AHX);
506 bb[i].jcaha += 1.4 * std::sin((bb[i].psi + 138.0) * DEG2RAD)
507 - 4.1 * std::cos(2.0 * DEG2RAD * (bb[i].psi + 138.0))
508 + 2.0 * std::cos(2.0 * DEG2RAD * (bb[i].phi + 30.0));
512 static void check_ahx(int nres, t_bb bb[], int* hstart, int* hend)
514 int h0, h1, h0sav, h1sav;
516 set_ahcity(nres, bb);
517 h0 = h0sav = h1sav = 0;
520 for (; (!bb[h0].bHelix) && (h0 < nres - 4); h0++) {}
521 for (h1 = h0; bb[h1 + 1].bHelix && (h1 < nres - 1); h1++) {}
522 if (h1 > h0)
524 /*fprintf(stderr,"Helix from %d to %d\n",h0,h1);*/
525 if (h1 - h0 > h1sav - h0sav)
527 h0sav = h0;
528 h1sav = h1;
531 h0 = h1 + 1;
532 } while (h1 < nres - 1);
533 *hstart = h0sav;
534 *hend = h1sav;
537 void do_start_end(int nres, t_bb bb[], int* nbb, int bbindex[], int* nca, int caindex[], gmx_bool bRange, int rStart, int rEnd)
539 int i, j, hstart = 0, hend = 0;
541 if (bRange)
543 for (i = 0; (i < nres); i++)
545 if ((bb[i].resno >= rStart) && (bb[i].resno <= rEnd))
547 bb[i].bHelix = TRUE;
549 if (bb[i].resno == rStart)
551 hstart = i;
553 if (bb[i].resno == rEnd)
555 hend = i;
559 else
561 /* Find start and end of longest helix fragment */
562 check_ahx(nres, bb, &hstart, &hend);
564 fprintf(stderr, "helix from: %d through %d\n", bb[hstart].resno, bb[hend].resno);
566 for (j = 0, i = hstart; (i <= hend); i++)
568 bbindex[j++] = bb[i].N;
569 bbindex[j++] = bb[i].H;
570 bbindex[j++] = bb[i].CA;
571 bbindex[j++] = bb[i].C;
572 bbindex[j++] = bb[i].O;
573 caindex[i - hstart] = bb[i].CA;
575 *nbb = j;
576 *nca = (hend - hstart + 1);
579 void pr_bb(FILE* fp, int nres, t_bb bb[])
581 int i;
583 fprintf(fp, "\n");
584 fprintf(fp, "%3s %3s %3s %3s %3s %7s %7s %7s %7s %7s %3s\n", "AA", "N", "Ca", "C", "O", "Phi",
585 "Psi", "D3", "D4", "D5", "Hx?");
586 for (i = 0; (i < nres); i++)
588 fprintf(fp, "%3d %3d %3d %3d %3d %7.2f %7.2f %7.3f %7.3f %7.3f %3s\n", bb[i].resno, bb[i].N,
589 bb[i].CA, bb[i].C, bb[i].O, bb[i].phi, bb[i].psi, bb[i].d3, bb[i].d4, bb[i].d5,
590 bb[i].bHelix ? "Yes" : "No");
592 fprintf(fp, "\n");