697777768e6b2e6ca3464c8dc7e08c0f48bf8423
[alexxy/gromacs.git] / src / gromacs / gmxana / hxprops.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,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "hxprops.h"
41
42 #include <cmath>
43 #include <cstring>
44
45 #include <algorithm>
46
47 #include "gromacs/listed_forces/bonded.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/math/utilities.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"
56
57 real ellipticity(int nres, t_bb bb[])
58 {
59     typedef struct
60     {
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)
75
76     int  i, j;
77     real ell, pp2, phi, psi;
78
79     ell = 0;
80     for (i = 0; (i < nres); i++)
81     {
82         phi = bb[i].phi;
83         psi = bb[i].psi;
84         for (j = 0; (j < NPPW); j++)
85         {
86             pp2 = gmx::square(phi - ppw[j].phi) + gmx::square(psi - ppw[j].psi);
87             if (pp2 < 64)
88             {
89                 bb[i].nhx++;
90                 ell += ppw[j].w;
91                 break;
92             }
93         }
94     }
95     return ell;
96 }
97
98 real ahx_len(int gnx, const int index[], rvec x[])
99 /* Assume we have a list of Calpha atoms only! */
100 {
101     rvec dx;
102
103     rvec_sub(x[index[0]], x[index[gnx - 1]], dx);
104
105     return norm(dx);
106 }
107
108 real radius(FILE* fp, int nca, const int ca_index[], rvec x[])
109 /* Assume we have all the backbone */
110 {
111     real dl2, dlt;
112     int  i, ai;
113
114     dlt = 0;
115     for (i = 0; (i < nca); i++)
116     {
117         ai  = ca_index[i];
118         dl2 = gmx::square(x[ai][XX]) + gmx::square(x[ai][YY]);
119
120         if (fp)
121         {
122             fprintf(fp, "  %10g", dl2);
123         }
124
125         dlt += dl2;
126     }
127     if (fp)
128     {
129         fprintf(fp, "\n");
130     }
131
132     return std::sqrt(dlt / nca);
133 }
134
135 static real rot(rvec x1, const rvec x2)
136 {
137     real phi1, dphi, cp, sp;
138     real xx, yy;
139
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];
145
146     dphi = gmx::c_rad2Deg * std::atan2(yy, xx);
147
148     return dphi;
149 }
150
151 real twist(int nca, const int caindex[], rvec x[])
152 {
153     real pt, dphi;
154     int  i, a0, a1;
155
156     pt = 0;
157     a0 = caindex[0];
158     for (i = 1; (i < nca); i++)
159     {
160         a1 = caindex[i];
161
162         dphi = rot(x[a0], x[a1]);
163         if (dphi < -90)
164         {
165             dphi += 360;
166         }
167         pt += dphi;
168         a0 = a1;
169     }
170
171     return (pt / (nca - 1));
172 }
173
174 real ca_phi(int gnx, const int index[], rvec x[])
175 /* Assume we have a list of Calpha atoms only! */
176 {
177     real phi, phitot;
178     int  i, ai, aj, ak, al, t1, t2, t3;
179     rvec r_ij, r_kj, r_kl, m, n;
180
181     if (gnx <= 4)
182     {
183         return 0;
184     }
185
186     phitot = 0;
187     for (i = 0; (i < gnx - 4); i++)
188     {
189         ai  = index[i + 0];
190         aj  = index[i + 1];
191         ak  = index[i + 2];
192         al  = index[i + 3];
193         phi = gmx::c_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;
196     }
197
198     return (phitot / (gnx - 4.0));
199 }
200
201 real dip(int nbb, int const bbind[], const rvec x[], const t_atom atom[])
202 {
203     int  i, m, ai;
204     rvec dipje;
205     real q;
206
207     clear_rvec(dipje);
208     for (i = 0; (i < nbb); i++)
209     {
210         ai = bbind[i];
211         q  = atom[ai].q;
212         for (m = 0; (m < DIM); m++)
213         {
214             dipje[m] += x[ai][m] * q;
215         }
216     }
217     return norm(dipje);
218 }
219
220 real rise(int gnx, const int index[], rvec x[])
221 /* Assume we have a list of Calpha atoms only! */
222 {
223     real z, z0, ztot;
224     int  i, ai;
225
226     ai   = index[0];
227     z0   = x[ai][ZZ];
228     ztot = 0;
229     for (i = 1; (i < gnx); i++)
230     {
231         ai = index[i];
232         z  = x[ai][ZZ];
233         ztot += (z - z0);
234         z0 = z;
235     }
236
237     return (ztot / (gnx - 1.0));
238 }
239
240 void av_hblen(FILE* fp3, FILE* fp3a, FILE* fp4, FILE* fp4a, FILE* fp5, FILE* fp5a, real t, int nres, t_bb bb[])
241 {
242     int  i, n3 = 0, n4 = 0, n5 = 0;
243     real d3 = 0, d4 = 0, d5 = 0;
244
245     for (i = 0; (i < nres - 3); i++)
246     {
247         if (bb[i].bHelix)
248         {
249             fprintf(fp3a, "%10g", bb[i].d3);
250             n3++;
251             d3 += bb[i].d3;
252             if (i < nres - 4)
253             {
254                 fprintf(fp4a, "%10g", bb[i].d4);
255                 n4++;
256                 d4 += bb[i].d4;
257             }
258             if (i < nres - 5)
259             {
260                 fprintf(fp5a, "%10g", bb[i].d5);
261                 n5++;
262                 d5 += bb[i].d5;
263             }
264         }
265     }
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");
272 }
273
274 void av_phipsi(FILE* fphi, FILE* fpsi, FILE* fphi2, FILE* fpsi2, real t, int nres, t_bb bb[])
275 {
276     int  i, n = 0;
277     real phi = 0, psi = 0;
278
279     fprintf(fphi2, "%10g", t);
280     fprintf(fpsi2, "%10g", t);
281     for (i = 0; (i < nres); i++)
282     {
283         if (bb[i].bHelix)
284         {
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++;
290         }
291     }
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");
296 }
297
298 static void set_ahcity(int nbb, t_bb bb[])
299 {
300     real pp2;
301     int  n;
302
303     for (n = 0; (n < nbb); n++)
304     {
305         pp2 = gmx::square(bb[n].phi - PHI_AHX) + gmx::square(bb[n].psi - PSI_AHX);
306
307         bb[n].bHelix = FALSE;
308         if (pp2 < 2500)
309         {
310             if ((bb[n].d4 < 0.36) || ((n > 0) && bb[n - 1].bHelix))
311             {
312                 bb[n].bHelix = TRUE;
313             }
314         }
315     }
316 }
317
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)
327 {
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   gnx, r0, r1;
333
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 (int i = 1; (i < gnx); i++)
340     {
341         r0 = std::min(r0, atom[(*index)[i]].resind);
342         r1 = std::max(r1, atom[(*index)[i]].resind);
343     }
344     int rnr = r1 - r0 + 1;
345     fprintf(stderr, "There are %d residues\n", rnr);
346     snew(bb, rnr);
347     for (int i = 0; (i < rnr); i++)
348     {
349         bb[i].N = bb[i].H = bb[i].CA = bb[i].C = bb[i].O = -1;
350         bb[i].resno                                      = res0 + i;
351     }
352
353     for (int i = 0; (i < gnx); i++)
354     {
355         int 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)
362         {
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)
366             {
367                 bb[bbindex].H = ai;
368             }
369         }
370         int k = 0;
371         for (; (k < NBB); k++)
372         {
373             if (std::strcmp(bb_nm[k], *(atomname[ai])) == 0)
374             {
375                 break;
376             }
377         }
378         switch (k)
379         {
380             case 0: bb[bbindex].N = ai; break;
381             case 1:
382             case 5:
383                 /* No attempt to address the case where some weird input has both H and HN atoms in the group */
384                 bb[bbindex].H = ai;
385                 break;
386             case 2: bb[bbindex].CA = ai; break;
387             case 3: bb[bbindex].C = ai; break;
388             case 4: bb[bbindex].O = ai; break;
389             default: break;
390         }
391     }
392
393     int i0 = 0;
394     for (; (i0 < rnr); i0++)
395     {
396         if ((bb[i0].N != -1) && (bb[i0].H != -1) && (bb[i0].CA != -1) && (bb[i0].C != -1)
397             && (bb[i0].O != -1))
398         {
399             break;
400         }
401     }
402     int i1 = rnr - 1;
403     for (; (i1 >= 0); i1--)
404     {
405         if ((bb[i1].N != -1) && (bb[i1].H != -1) && (bb[i1].CA != -1) && (bb[i1].C != -1)
406             && (bb[i1].O != -1))
407         {
408             break;
409         }
410     }
411     if (i0 == 0)
412     {
413         i0++;
414     }
415     if (i1 == rnr - 1)
416     {
417         i1--;
418     }
419
420     for (int i = i0; (i < i1); i++)
421     {
422         bb[i].Cprev = bb[i - 1].C;
423         bb[i].Nnext = bb[i + 1].N;
424     }
425     rnr = std::max(0, i1 - i0 + 1);
426     fprintf(stderr, "There are %d complete backbone residues (from %d to %d)\n", rnr, bb[i0].resno, bb[i1].resno);
427     if (rnr == 0)
428     {
429         gmx_fatal(FARGS, "Zero complete backbone residues were found, cannot proceed");
430     }
431     for (int i = 0; (i < rnr); i++, i0++)
432     {
433         bb[i] = bb[i0];
434     }
435
436     /* Set the labels */
437     for (int i = 0; (i < rnr); i++)
438     {
439         int resindex = atom[bb[i].CA].resind;
440         sprintf(bb[i].label, "%s%d", *(resinfo[resindex].name), resinfo[resindex].nr);
441     }
442
443     *nres = rnr;
444     *nbb  = rnr * asize(bb_nm);
445
446     return bb;
447 }
448
449 real pprms(FILE* fp, int nbb, t_bb bb[])
450 {
451     int  i, n;
452     real rms, rmst, rms2;
453
454     rmst = rms2 = 0;
455     for (i = n = 0; (i < nbb); i++)
456     {
457         if (bb[i].bHelix)
458         {
459             rms = std::sqrt(bb[i].pprms2);
460             rmst += rms;
461             rms2 += bb[i].pprms2;
462             fprintf(fp, "%10g  ", rms);
463             n++;
464         }
465     }
466     fprintf(fp, "\n");
467     rms = std::sqrt(rms2 / n - gmx::square(rmst / n));
468
469     return rms;
470 }
471
472 void calc_hxprops(int nres, t_bb bb[], const rvec x[])
473 {
474     int  i, ao, an, t1, t2, t3;
475     rvec dx, r_ij, r_kj, r_kl, m, n;
476
477     for (i = 0; (i < nres); i++)
478     {
479         ao       = bb[i].O;
480         bb[i].d4 = bb[i].d3 = bb[i].d5 = 0;
481         if (i < nres - 3)
482         {
483             an = bb[i + 3].N;
484             rvec_sub(x[ao], x[an], dx);
485             bb[i].d3 = norm(dx);
486         }
487         if (i < nres - 4)
488         {
489             an = bb[i + 4].N;
490             rvec_sub(x[ao], x[an], dx);
491             bb[i].d4 = norm(dx);
492         }
493         if (i < nres - 5)
494         {
495             an = bb[i + 5].N;
496             rvec_sub(x[ao], x[an], dx);
497             bb[i].d5 = norm(dx);
498         }
499
500         bb[i].phi = gmx::c_rad2Deg
501                     * dih_angle(x[bb[i].Cprev],
502                                 x[bb[i].N],
503                                 x[bb[i].CA],
504                                 x[bb[i].C],
505                                 nullptr,
506                                 r_ij,
507                                 r_kj,
508                                 r_kl,
509                                 m,
510                                 n,
511                                 &t1,
512                                 &t2,
513                                 &t3);
514         bb[i].psi = gmx::c_rad2Deg
515                     * dih_angle(x[bb[i].N],
516                                 x[bb[i].CA],
517                                 x[bb[i].C],
518                                 x[bb[i].Nnext],
519                                 nullptr,
520                                 r_ij,
521                                 r_kj,
522                                 r_kl,
523                                 m,
524                                 n,
525                                 &t1,
526                                 &t2,
527                                 &t3);
528         bb[i].pprms2 = gmx::square(bb[i].phi - PHI_AHX) + gmx::square(bb[i].psi - PSI_AHX);
529
530         bb[i].jcaha += 1.4 * std::sin((bb[i].psi + 138.0) * gmx::c_deg2Rad)
531                        - 4.1 * std::cos(2.0 * gmx::c_deg2Rad * (bb[i].psi + 138.0))
532                        + 2.0 * std::cos(2.0 * gmx::c_deg2Rad * (bb[i].phi + 30.0));
533     }
534 }
535
536 static void check_ahx(int nres, t_bb bb[], int* hstart, int* hend)
537 {
538     int h0, h1, h0sav, h1sav;
539
540     set_ahcity(nres, bb);
541     h0 = h0sav = h1sav = 0;
542     do
543     {
544         for (; (!bb[h0].bHelix) && (h0 < nres - 4); h0++) {}
545         for (h1 = h0; bb[h1 + 1].bHelix && (h1 < nres - 1); h1++) {}
546         if (h1 > h0)
547         {
548             /*fprintf(stderr,"Helix from %d to %d\n",h0,h1);*/
549             if (h1 - h0 > h1sav - h0sav)
550             {
551                 h0sav = h0;
552                 h1sav = h1;
553             }
554         }
555         h0 = h1 + 1;
556     } while (h1 < nres - 1);
557     *hstart = h0sav;
558     *hend   = h1sav;
559 }
560
561 void do_start_end(int nres, t_bb bb[], int* nbb, int bbindex[], int* nca, int caindex[], gmx_bool bRange, int rStart, int rEnd)
562 {
563     int i, j, hstart = 0, hend = 0;
564
565     if (bRange)
566     {
567         for (i = 0; (i < nres); i++)
568         {
569             if ((bb[i].resno >= rStart) && (bb[i].resno <= rEnd))
570             {
571                 bb[i].bHelix = TRUE;
572             }
573             if (bb[i].resno == rStart)
574             {
575                 hstart = i;
576             }
577             if (bb[i].resno == rEnd)
578             {
579                 hend = i;
580             }
581         }
582     }
583     else
584     {
585         /* Find start and end of longest helix fragment */
586         check_ahx(nres, bb, &hstart, &hend);
587     }
588     fprintf(stderr, "helix from: %d through %d\n", bb[hstart].resno, bb[hend].resno);
589
590     for (j = 0, i = hstart; (i <= hend); i++)
591     {
592         bbindex[j++]        = bb[i].N;
593         bbindex[j++]        = bb[i].H;
594         bbindex[j++]        = bb[i].CA;
595         bbindex[j++]        = bb[i].C;
596         bbindex[j++]        = bb[i].O;
597         caindex[i - hstart] = bb[i].CA;
598     }
599     *nbb = j;
600     *nca = (hend - hstart + 1);
601 }
602
603 void pr_bb(FILE* fp, int nres, t_bb bb[])
604 {
605     int i;
606
607     fprintf(fp, "\n");
608     fprintf(fp,
609             "%3s %3s %3s %3s %3s %7s %7s %7s %7s %7s %3s\n",
610             "AA",
611             "N",
612             "Ca",
613             "C",
614             "O",
615             "Phi",
616             "Psi",
617             "D3",
618             "D4",
619             "D5",
620             "Hx?");
621     for (i = 0; (i < nres); i++)
622     {
623         fprintf(fp,
624                 "%3d %3d %3d %3d %3d %7.2f %7.2f %7.3f %7.3f %7.3f %3s\n",
625                 bb[i].resno,
626                 bb[i].N,
627                 bb[i].CA,
628                 bb[i].C,
629                 bb[i].O,
630                 bb[i].phi,
631                 bb[i].psi,
632                 bb[i].d3,
633                 bb[i].d4,
634                 bb[i].d5,
635                 bb[i].bHelix ? "Yes" : "No");
636     }
637     fprintf(fp, "\n");
638 }