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