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,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.
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.
46 #include "gromacs/listed_forces/bonded.h"
47 #include "gromacs/math/functions.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/topology/index.h"
51 #include "gromacs/topology/topology.h"
52 #include "gromacs/utility/arraysize.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/smalloc.h"
56 real ellipticity(int nres, t_bb bb[])
62 // Avoid warnings about narrowing conversions from double to real
64 # pragma warning(disable : 4838)
66 static const t_ppwstr ppw[] = { { -67, -44, 0.31 }, { -66, -41, 0.31 }, { -59, -44, 0.44 },
67 { -57, -47, 0.56 }, { -53, -52, 0.78 }, { -48, -57, 1.00 },
68 { -70.5, -35.8, 0.15 }, { -57, -79, 0.23 }, { -38, -78, 1.20 },
69 { -60, -30, 0.24 }, { -54, -28, 0.46 }, { -44, -33, 0.68 } };
71 # pragma warning(default : 4838)
73 #define NPPW asize(ppw)
76 real ell, pp2, phi, psi;
79 for (i = 0; (i < nres); i++)
83 for (j = 0; (j < NPPW); j++)
85 pp2 = gmx::square(phi - ppw[j].phi) + gmx::square(psi - ppw[j].psi);
97 real ahx_len(int gnx, const int index[], rvec x[])
98 /* Assume we have a list of Calpha atoms only! */
102 rvec_sub(x[index[0]], x[index[gnx - 1]], dx);
107 real radius(FILE* fp, int nca, const int ca_index[], rvec x[])
108 /* Assume we have all the backbone */
114 for (i = 0; (i < nca); i++)
117 dl2 = gmx::square(x[ai][XX]) + gmx::square(x[ai][YY]);
121 fprintf(fp, " %10g", dl2);
131 return std::sqrt(dlt / nca);
134 static real rot(rvec x1, const rvec x2)
136 real phi1, dphi, cp, sp;
139 phi1 = std::atan2(x1[YY], x1[XX]);
142 xx = cp * x2[XX] + sp * x2[YY];
143 yy = -sp * x2[XX] + cp * x2[YY];
145 dphi = RAD2DEG * std::atan2(yy, xx);
150 real twist(int nca, const int caindex[], rvec x[])
157 for (i = 1; (i < nca); i++)
161 dphi = rot(x[a0], x[a1]);
170 return (pt / (nca - 1));
173 real ca_phi(int gnx, const int index[], rvec x[])
174 /* Assume we have a list of Calpha atoms only! */
177 int i, ai, aj, ak, al, t1, t2, t3;
178 rvec r_ij, r_kj, r_kl, m, n;
186 for (i = 0; (i < gnx - 4); i++)
193 * dih_angle(x[ai], x[aj], x[ak], x[al], nullptr, r_ij, r_kj, r_kl, m, n, &t1, &t2, &t3);
197 return (phitot / (gnx - 4.0));
200 real dip(int nbb, int const bbind[], const rvec x[], const t_atom atom[])
207 for (i = 0; (i < nbb); i++)
211 for (m = 0; (m < DIM); m++)
213 dipje[m] += x[ai][m] * q;
219 real rise(int gnx, const int index[], rvec x[])
220 /* Assume we have a list of Calpha atoms only! */
228 for (i = 1; (i < gnx); i++)
236 return (ztot / (gnx - 1.0));
239 void av_hblen(FILE* fp3, FILE* fp3a, FILE* fp4, FILE* fp4a, FILE* fp5, FILE* fp5a, real t, int nres, t_bb bb[])
241 int i, n3 = 0, n4 = 0, n5 = 0;
242 real d3 = 0, d4 = 0, d5 = 0;
244 for (i = 0; (i < nres - 3); i++)
248 fprintf(fp3a, "%10g", bb[i].d3);
253 fprintf(fp4a, "%10g", bb[i].d4);
259 fprintf(fp5a, "%10g", bb[i].d5);
265 fprintf(fp3, "%10g %10g\n", t, d3 / n3);
266 fprintf(fp4, "%10g %10g\n", t, d4 / n4);
267 fprintf(fp5, "%10g %10g\n", t, d5 / n5);
273 void av_phipsi(FILE* fphi, FILE* fpsi, FILE* fphi2, FILE* fpsi2, real t, int nres, t_bb bb[])
276 real phi = 0, psi = 0;
278 fprintf(fphi2, "%10g", t);
279 fprintf(fpsi2, "%10g", t);
280 for (i = 0; (i < nres); i++)
286 fprintf(fphi2, " %10g", bb[i].phi);
287 fprintf(fpsi2, " %10g", bb[i].psi);
291 fprintf(fphi, "%10g %10g\n", t, (phi / n));
292 fprintf(fpsi, "%10g %10g\n", t, (psi / n));
293 fprintf(fphi2, "\n");
294 fprintf(fpsi2, "\n");
297 static void set_ahcity(int nbb, t_bb bb[])
302 for (n = 0; (n < nbb); n++)
304 pp2 = gmx::square(bb[n].phi - PHI_AHX) + gmx::square(bb[n].psi - PSI_AHX);
306 bb[n].bHelix = FALSE;
309 if ((bb[n].d4 < 0.36) || ((n > 0) && bb[n - 1].bHelix))
317 t_bb* mkbbind(const char* fn,
327 static const char* bb_nm[] = { "N", "H", "CA", "C", "O", "HN" };
328 #define NBB asize(bb_nm)
331 int ai, i, i0, i1, j, k, rnr, gnx, r0, r1;
333 fprintf(stderr, "Please select a group containing the entire backbone\n");
334 rd_index(fn, 1, &gnx, index, &grpname);
336 fprintf(stderr, "Checking group %s\n", grpname);
337 r0 = r1 = atom[(*index)[0]].resind;
338 for (i = 1; (i < gnx); i++)
340 r0 = std::min(r0, atom[(*index)[i]].resind);
341 r1 = std::max(r1, atom[(*index)[i]].resind);
344 fprintf(stderr, "There are %d residues\n", rnr);
346 for (i = 0; (i < rnr); i++)
348 bb[i].N = bb[i].H = bb[i].CA = bb[i].C = bb[i].O = -1;
349 bb[i].resno = res0 + i;
352 for (i = j = 0; (i < gnx); i++)
355 // Create an index into the residue index for the topology.
356 int resindex = atom[ai].resind;
357 // Create an index into the residues present in the selected
359 int bbindex = resindex - r0;
360 if (std::strcmp(*(resinfo[resindex].name), "PRO") == 0)
362 // For PRO in a peptide, there is no H bound to backbone
363 // N, so use CD instead.
364 if (std::strcmp(*(atomname[ai]), "CD") == 0)
369 for (k = 0; (k < NBB); k++)
371 if (std::strcmp(bb_nm[k], *(atomname[ai])) == 0)
378 case 0: bb[bbindex].N = ai; break;
381 /* No attempt to address the case where some weird input has both H and HN atoms in the group */
384 case 2: bb[bbindex].CA = ai; break;
385 case 3: bb[bbindex].C = ai; break;
386 case 4: bb[bbindex].O = ai; break;
391 for (i0 = 0; (i0 < rnr); i0++)
393 if ((bb[i0].N != -1) && (bb[i0].H != -1) && (bb[i0].CA != -1) && (bb[i0].C != -1)
399 for (i1 = rnr - 1; (i1 >= 0); i1--)
401 if ((bb[i1].N != -1) && (bb[i1].H != -1) && (bb[i1].CA != -1) && (bb[i1].C != -1)
416 for (i = i0; (i < i1); i++)
418 bb[i].Cprev = bb[i - 1].C;
419 bb[i].Nnext = bb[i + 1].N;
421 rnr = std::max(0, i1 - i0 + 1);
422 fprintf(stderr, "There are %d complete backbone residues (from %d to %d)\n", rnr, bb[i0].resno,
426 gmx_fatal(FARGS, "Zero complete backbone residues were found, cannot proceed");
428 for (i = 0; (i < rnr); i++, i0++)
434 for (i = 0; (i < rnr); i++)
436 int resindex = atom[bb[i].CA].resind;
437 sprintf(bb[i].label, "%s%d", *(resinfo[resindex].name), resinfo[resindex].nr);
441 *nbb = rnr * asize(bb_nm);
446 real pprms(FILE* fp, int nbb, t_bb bb[])
449 real rms, rmst, rms2;
452 for (i = n = 0; (i < nbb); i++)
456 rms = std::sqrt(bb[i].pprms2);
458 rms2 += bb[i].pprms2;
459 fprintf(fp, "%10g ", rms);
464 rms = std::sqrt(rms2 / n - gmx::square(rmst / n));
469 void calc_hxprops(int nres, t_bb bb[], const rvec x[])
471 int i, ao, an, t1, t2, t3;
472 rvec dx, r_ij, r_kj, r_kl, m, n;
474 for (i = 0; (i < nres); i++)
477 bb[i].d4 = bb[i].d3 = bb[i].d5 = 0;
481 rvec_sub(x[ao], x[an], dx);
487 rvec_sub(x[ao], x[an], dx);
493 rvec_sub(x[ao], x[an], dx);
498 * dih_angle(x[bb[i].Cprev], x[bb[i].N], x[bb[i].CA], x[bb[i].C], nullptr, r_ij,
499 r_kj, r_kl, m, n, &t1, &t2, &t3);
501 * dih_angle(x[bb[i].N], x[bb[i].CA], x[bb[i].C], x[bb[i].Nnext], nullptr, r_ij,
502 r_kj, r_kl, m, n, &t1, &t2, &t3);
503 bb[i].pprms2 = gmx::square(bb[i].phi - PHI_AHX) + gmx::square(bb[i].psi - PSI_AHX);
505 bb[i].jcaha += 1.4 * std::sin((bb[i].psi + 138.0) * DEG2RAD)
506 - 4.1 * std::cos(2.0 * DEG2RAD * (bb[i].psi + 138.0))
507 + 2.0 * std::cos(2.0 * DEG2RAD * (bb[i].phi + 30.0));
511 static void check_ahx(int nres, t_bb bb[], int* hstart, int* hend)
513 int h0, h1, h0sav, h1sav;
515 set_ahcity(nres, bb);
516 h0 = h0sav = h1sav = 0;
519 for (; (!bb[h0].bHelix) && (h0 < nres - 4); h0++) {}
520 for (h1 = h0; bb[h1 + 1].bHelix && (h1 < nres - 1); h1++) {}
523 /*fprintf(stderr,"Helix from %d to %d\n",h0,h1);*/
524 if (h1 - h0 > h1sav - h0sav)
531 } while (h1 < nres - 1);
536 void do_start_end(int nres, t_bb bb[], int* nbb, int bbindex[], int* nca, int caindex[], gmx_bool bRange, int rStart, int rEnd)
538 int i, j, hstart = 0, hend = 0;
542 for (i = 0; (i < nres); i++)
544 if ((bb[i].resno >= rStart) && (bb[i].resno <= rEnd))
548 if (bb[i].resno == rStart)
552 if (bb[i].resno == rEnd)
560 /* Find start and end of longest helix fragment */
561 check_ahx(nres, bb, &hstart, &hend);
563 fprintf(stderr, "helix from: %d through %d\n", bb[hstart].resno, bb[hend].resno);
565 for (j = 0, i = hstart; (i <= hend); i++)
567 bbindex[j++] = bb[i].N;
568 bbindex[j++] = bb[i].H;
569 bbindex[j++] = bb[i].CA;
570 bbindex[j++] = bb[i].C;
571 bbindex[j++] = bb[i].O;
572 caindex[i - hstart] = bb[i].CA;
575 *nca = (hend - hstart + 1);
578 void pr_bb(FILE* fp, int nres, t_bb bb[])
583 fprintf(fp, "%3s %3s %3s %3s %3s %7s %7s %7s %7s %7s %3s\n", "AA", "N", "Ca", "C", "O", "Phi",
584 "Psi", "D3", "D4", "D5", "Hx?");
585 for (i = 0; (i < nres); i++)
587 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,
588 bb[i].CA, bb[i].C, bb[i].O, bb[i].phi, bb[i].psi, bb[i].d3, bb[i].d4, bb[i].d5,
589 bb[i].bHelix ? "Yes" : "No");