2d09699aee67f60879ef9617ac4d4c7b6b86bcf6
[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, 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, 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, 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, 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, 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, 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     real sign;
189
190     if (gnx <= 4)
191     {
192         return 0;
193     }
194
195     phitot = 0;
196     for (i = 0; (i < gnx-4); i++)
197     {
198         ai  = index[i+0];
199         aj  = index[i+1];
200         ak  = index[i+2];
201         al  = index[i+3];
202         phi = RAD2DEG*
203             dih_angle(x[ai], x[aj], x[ak], x[al], NULL,
204                       r_ij, r_kj, r_kl, m, n,
205                       &sign, &t1, &t2, &t3);
206         phitot += phi;
207     }
208
209     return (phitot/(gnx-4.0));
210 }
211
212 real dip(int nbb, int const bbind[], const rvec x[], const t_atom atom[])
213 {
214     int  i, m, ai;
215     rvec dipje;
216     real q;
217
218     clear_rvec(dipje);
219     for (i = 0; (i < nbb); i++)
220     {
221         ai = bbind[i];
222         q  = atom[ai].q;
223         for (m = 0; (m < DIM); m++)
224         {
225             dipje[m] += x[ai][m]*q;
226         }
227     }
228     return norm(dipje);
229 }
230
231 real rise(int gnx, int index[], rvec x[])
232 /* Assume we have a list of Calpha atoms only! */
233 {
234     real z, z0, ztot;
235     int  i, ai;
236
237     ai   = index[0];
238     z0   = x[ai][ZZ];
239     ztot = 0;
240     for (i = 1; (i < gnx); i++)
241     {
242         ai    = index[i];
243         z     = x[ai][ZZ];
244         ztot += (z-z0);
245         z0    = z;
246     }
247
248     return (ztot/(gnx-1.0));
249 }
250
251 void av_hblen(FILE *fp3, FILE *fp3a,
252               FILE *fp4, FILE *fp4a,
253               FILE *fp5, FILE *fp5a,
254               real t, int nres, t_bb bb[])
255 {
256     int  i, n3 = 0, n4 = 0, n5 = 0;
257     real d3 = 0, d4 = 0, d5 = 0;
258
259     for (i = 0; (i < nres-3); i++)
260     {
261         if (bb[i].bHelix)
262         {
263             fprintf(fp3a,  "%10g", bb[i].d3);
264             n3++;
265             d3 += bb[i].d3;
266             if (i < nres-4)
267             {
268                 fprintf(fp4a, "%10g", bb[i].d4);
269                 n4++;
270                 d4 += bb[i].d4;
271             }
272             if (i < nres-5)
273             {
274                 fprintf(fp5a, "%10g", bb[i].d5);
275                 n5++;
276                 d5 += bb[i].d5;
277             }
278         }
279     }
280     fprintf(fp3, "%10g  %10g\n", t, d3/n3);
281     fprintf(fp4, "%10g  %10g\n", t, d4/n4);
282     fprintf(fp5, "%10g  %10g\n", t, d5/n5);
283     fprintf(fp3a, "\n");
284     fprintf(fp4a, "\n");
285     fprintf(fp5a, "\n");
286
287 }
288
289 void av_phipsi(FILE *fphi, FILE *fpsi, FILE *fphi2, FILE *fpsi2,
290                real t, int nres, t_bb bb[])
291 {
292     int  i, n = 0;
293     real phi = 0, psi = 0;
294
295     fprintf(fphi2, "%10g", t);
296     fprintf(fpsi2, "%10g", t);
297     for (i = 0; (i < nres); i++)
298     {
299         if (bb[i].bHelix)
300         {
301             phi += bb[i].phi;
302             psi += bb[i].psi;
303             fprintf(fphi2, "  %10g", bb[i].phi);
304             fprintf(fpsi2, "  %10g", bb[i].psi);
305             n++;
306         }
307     }
308     fprintf(fphi, "%10g  %10g\n", t, (phi/n));
309     fprintf(fpsi, "%10g  %10g\n", t, (psi/n));
310     fprintf(fphi2, "\n");
311     fprintf(fpsi2, "\n");
312 }
313
314 static void set_ahcity(int nbb, t_bb bb[])
315 {
316     real pp2;
317     int  n;
318
319     for (n = 0; (n < nbb); n++)
320     {
321         pp2 = gmx::square(bb[n].phi-PHI_AHX)+gmx::square(bb[n].psi-PSI_AHX);
322
323         bb[n].bHelix = FALSE;
324         if (pp2 < 2500)
325         {
326             if ((bb[n].d4 < 0.36) || ((n > 0) && bb[n-1].bHelix))
327             {
328                 bb[n].bHelix = TRUE;
329             }
330         }
331     }
332 }
333
334 t_bb *mkbbind(const char *fn, int *nres, int *nbb, int res0,
335               int *nall, int **index,
336               char ***atomname, t_atom atom[],
337               t_resinfo *resinfo)
338 {
339     static const char * bb_nm[] = { "N", "H", "CA", "C", "O", "HN" };
340 #define NBB asize(bb_nm)
341     t_bb               *bb;
342     char               *grpname;
343     int                 ai, i, i0, i1, j, k, ri, rnr, gnx, r0, r1;
344
345     fprintf(stderr, "Please select a group containing the entire backbone\n");
346     rd_index(fn, 1, &gnx, index, &grpname);
347     *nall = gnx;
348     fprintf(stderr, "Checking group %s\n", grpname);
349     r0 = r1 = atom[(*index)[0]].resind;
350     for (i = 1; (i < gnx); i++)
351     {
352         r0 = std::min(r0, atom[(*index)[i]].resind);
353         r1 = std::max(r1, atom[(*index)[i]].resind);
354     }
355     rnr = r1-r0+1;
356     fprintf(stderr, "There are %d residues\n", rnr);
357     snew(bb, rnr);
358     for (i = 0; (i < rnr); i++)
359     {
360         bb[i].N = bb[i].H = bb[i].CA = bb[i].C = bb[i].O = -1, bb[i].resno = res0+i;
361     }
362
363     for (i = j = 0; (i < gnx); i++)
364     {
365         ai = (*index)[i];
366         ri = atom[ai].resind-r0;
367         if (std::strcmp(*(resinfo[ri].name), "PRO") == 0)
368         {
369             if (std::strcmp(*(atomname[ai]), "CD") == 0)
370             {
371                 bb[ri].H = ai;
372             }
373         }
374         for (k = 0; (k < NBB); k++)
375         {
376             if (std::strcmp(bb_nm[k], *(atomname[ai])) == 0)
377             {
378                 j++;
379                 break;
380             }
381         }
382         switch (k)
383         {
384             case 0:
385                 bb[ri].N = ai;
386                 break;
387             case 1:
388             case 5:
389                 /* No attempt to address the case where some weird input has both H and HN atoms in the group */
390                 bb[ri].H = ai;
391                 break;
392             case 2:
393                 bb[ri].CA = ai;
394                 break;
395             case 3:
396                 bb[ri].C = ai;
397                 break;
398             case 4:
399                 bb[ri].O = ai;
400                 break;
401             default:
402                 break;
403         }
404     }
405
406     for (i0 = 0; (i0 < rnr); i0++)
407     {
408         if ((bb[i0].N != -1) && (bb[i0].H != -1) &&
409             (bb[i0].CA != -1) &&
410             (bb[i0].C != -1) && (bb[i0].O != -1))
411         {
412             break;
413         }
414     }
415     for (i1 = rnr-1; (i1 >= 0); i1--)
416     {
417         if ((bb[i1].N != -1) && (bb[i1].H != -1) &&
418             (bb[i1].CA != -1) &&
419             (bb[i1].C != -1) && (bb[i1].O != -1))
420         {
421             break;
422         }
423     }
424     if (i0 == 0)
425     {
426         i0++;
427     }
428     if (i1 == rnr-1)
429     {
430         i1--;
431     }
432
433     for (i = i0; (i < i1); i++)
434     {
435         bb[i].Cprev = bb[i-1].C;
436         bb[i].Nnext = bb[i+1].N;
437     }
438     rnr = std::max(0, i1-i0+1);
439     fprintf(stderr, "There are %d complete backbone residues (from %d to %d)\n",
440             rnr, bb[i0].resno, bb[i1].resno);
441     if (rnr == 0)
442     {
443         gmx_fatal(FARGS, "Zero complete backbone residues were found, cannot proceed");
444     }
445     for (i = 0; (i < rnr); i++, i0++)
446     {
447         bb[i] = bb[i0];
448     }
449
450     /* Set the labels */
451     for (i = 0; (i < rnr); i++)
452     {
453         ri = atom[bb[i].CA].resind;
454         sprintf(bb[i].label, "%s%d", *(resinfo[ri].name), ri+res0);
455     }
456
457     *nres = rnr;
458     *nbb  = rnr*asize(bb_nm);
459
460     return bb;
461 }
462
463 real pprms(FILE *fp, int nbb, t_bb bb[])
464 {
465     int  i, n;
466     real rms, rmst, rms2;
467
468     rmst = rms2 = 0;
469     for (i = n = 0; (i < nbb); i++)
470     {
471         if (bb[i].bHelix)
472         {
473             rms   = std::sqrt(bb[i].pprms2);
474             rmst += rms;
475             rms2 += bb[i].pprms2;
476             fprintf(fp, "%10g  ", rms);
477             n++;
478         }
479     }
480     fprintf(fp, "\n");
481     rms = std::sqrt(rms2/n-gmx::square(rmst/n));
482
483     return rms;
484 }
485
486 void calc_hxprops(int nres, t_bb bb[], const rvec x[])
487 {
488     int  i, ao, an, t1, t2, t3;
489     rvec dx, r_ij, r_kj, r_kl, m, n;
490     real sign;
491
492     for (i = 0; (i < nres); i++)
493     {
494         ao       = bb[i].O;
495         bb[i].d4 = bb[i].d3 = bb[i].d5 = 0;
496         if (i < nres-3)
497         {
498             an = bb[i+3].N;
499             rvec_sub(x[ao], x[an], dx);
500             bb[i].d3 = norm(dx);
501         }
502         if (i < nres-4)
503         {
504             an = bb[i+4].N;
505             rvec_sub(x[ao], x[an], dx);
506             bb[i].d4 = norm(dx);
507         }
508         if (i < nres-5)
509         {
510             an = bb[i+5].N;
511             rvec_sub(x[ao], x[an], dx);
512             bb[i].d5 = norm(dx);
513         }
514
515         bb[i].phi = RAD2DEG*
516             dih_angle(x[bb[i].Cprev], x[bb[i].N], x[bb[i].CA], x[bb[i].C], NULL,
517                       r_ij, r_kj, r_kl, m, n,
518                       &sign, &t1, &t2, &t3);
519         bb[i].psi = RAD2DEG*
520             dih_angle(x[bb[i].N], x[bb[i].CA], x[bb[i].C], x[bb[i].Nnext], NULL,
521                       r_ij, r_kj, r_kl, m, n,
522                       &sign, &t1, &t2, &t3);
523         bb[i].pprms2 = gmx::square(bb[i].phi-PHI_AHX)+gmx::square(bb[i].psi-PSI_AHX);
524
525         bb[i].jcaha +=
526             1.4*std::sin((bb[i].psi+138.0)*DEG2RAD) -
527             4.1*std::cos(2.0*DEG2RAD*(bb[i].psi+138.0)) +
528             2.0*std::cos(2.0*DEG2RAD*(bb[i].phi+30.0));
529     }
530 }
531
532 static void check_ahx(int nres, t_bb bb[],
533                       int *hstart, int *hend)
534 {
535     int  h0, h1, h0sav, h1sav;
536
537     set_ahcity(nres, bb);
538     h0 = h0sav = h1sav = 0;
539     do
540     {
541         for (; (!bb[h0].bHelix) && (h0 < nres-4); h0++)
542         {
543             ;
544         }
545         for (h1 = h0; bb[h1+1].bHelix && (h1 < nres-1); h1++)
546         {
547             ;
548         }
549         if (h1 > h0)
550         {
551             /*fprintf(stderr,"Helix from %d to %d\n",h0,h1);*/
552             if (h1-h0 > h1sav-h0sav)
553             {
554                 h0sav = h0;
555                 h1sav = h1;
556             }
557         }
558         h0 = h1+1;
559     }
560     while (h1 < nres-1);
561     *hstart = h0sav;
562     *hend   = h1sav;
563 }
564
565 void do_start_end(int nres, t_bb bb[], int *nbb, int bbindex[],
566                   int *nca, int caindex[],
567                   gmx_bool bRange, int rStart, int rEnd)
568 {
569     int    i, j, hstart = 0, hend = 0;
570
571     if (bRange)
572     {
573         for (i = 0; (i < nres); i++)
574         {
575             if ((bb[i].resno >= rStart) && (bb[i].resno <= rEnd))
576             {
577                 bb[i].bHelix = TRUE;
578             }
579             if (bb[i].resno == rStart)
580             {
581                 hstart = i;
582             }
583             if (bb[i].resno == rEnd)
584             {
585                 hend = i;
586             }
587         }
588     }
589     else
590     {
591         /* Find start and end of longest helix fragment */
592         check_ahx(nres, bb, &hstart, &hend);
593     }
594     fprintf(stderr, "helix from: %d through %d\n",
595             bb[hstart].resno, bb[hend].resno);
596
597     for (j = 0, i = hstart; (i <= hend); i++)
598     {
599         bbindex[j++]      = bb[i].N;
600         bbindex[j++]      = bb[i].H;
601         bbindex[j++]      = bb[i].CA;
602         bbindex[j++]      = bb[i].C;
603         bbindex[j++]      = bb[i].O;
604         caindex[i-hstart] = bb[i].CA;
605     }
606     *nbb = j;
607     *nca = (hend-hstart+1);
608 }
609
610 void pr_bb(FILE *fp, int nres, t_bb bb[])
611 {
612     int i;
613
614     fprintf(fp, "\n");
615     fprintf(fp, "%3s %3s %3s %3s %3s %7s %7s %7s %7s %7s %3s\n",
616             "AA", "N", "Ca", "C", "O", "Phi", "Psi", "D3", "D4", "D5", "Hx?");
617     for (i = 0; (i < nres); i++)
618     {
619         fprintf(fp, "%3d %3d %3d %3d %3d %7.2f %7.2f %7.3f %7.3f %7.3f %3s\n",
620                 bb[i].resno, bb[i].N, bb[i].CA, bb[i].C, bb[i].O,
621                 bb[i].phi, bb[i].psi, bb[i].d3, bb[i].d4, bb[i].d5,
622                 bb[i].bHelix ? "Yes" : "No");
623     }
624     fprintf(fp, "\n");
625 }