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