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