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