Move symtab.* to topology/
[gromacs.git] / src / gromacs / gmxana / hxprops.c
blob5a999cec6655af0c7cac25bfd6a7590f65be86dc
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
41 #include <math.h>
42 #include <string.h>
44 #include "typedefs.h"
45 #include "macros.h"
46 #include "gromacs/math/units.h"
47 #include "gromacs/math/vec.h"
48 #include "index.h"
49 #include "hxprops.h"
50 #include "gromacs/utility/smalloc.h"
51 #include "bondf.h"
53 #include "gromacs/utility/fatalerror.h"
55 real ellipticity(int nres, t_bb bb[])
57 typedef struct {
58 real phi, psi, w;
59 } t_ppwstr;
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 }
75 #define NPPW asize(ppw)
77 int i, j;
78 const real Xi = 1.0;
79 real ell, pp2, phi, psi;
81 ell = 0;
82 for (i = 0; (i < nres); i++)
84 phi = bb[i].phi;
85 psi = bb[i].psi;
86 for (j = 0; (j < NPPW); j++)
88 pp2 = sqr(phi-ppw[j].phi)+sqr(psi-ppw[j].psi);
89 if (pp2 < 64)
91 bb[i].nhx++;
92 ell += ppw[j].w;
93 break;
97 return ell;
100 real ahx_len(int gnx, atom_id index[], rvec x[])
101 /* Assume we have a list of Calpha atoms only! */
103 rvec dx;
105 rvec_sub(x[index[0]], x[index[gnx-1]], dx);
107 return norm(dx);
110 real radius(FILE *fp, int nca, atom_id ca_index[], rvec x[])
111 /* Assume we have all the backbone */
113 real dl2, dlt;
114 int i, ai;
116 dlt = 0;
117 for (i = 0; (i < nca); i++)
119 ai = ca_index[i];
120 dl2 = sqr(x[ai][XX])+sqr(x[ai][YY]);
122 if (fp)
124 fprintf(fp, " %10g", dl2);
127 dlt += dl2;
129 if (fp)
131 fprintf(fp, "\n");
134 return sqrt(dlt/nca);
137 static real rot(rvec x1, rvec x2)
139 real phi1, dphi, cp, sp;
140 real xx, yy;
142 phi1 = atan2(x1[YY], x1[XX]);
143 cp = cos(phi1);
144 sp = sin(phi1);
145 xx = cp*x2[XX]+sp*x2[YY];
146 yy = -sp*x2[XX]+cp*x2[YY];
148 dphi = RAD2DEG*atan2(yy, xx);
150 return dphi;
153 real twist(int nca, atom_id caindex[], rvec x[])
155 real pt, dphi;
156 int i, a0, a1;
158 pt = 0;
159 a0 = caindex[0];
160 for (i = 1; (i < nca); i++)
162 a1 = caindex[i];
164 dphi = rot(x[a0], x[a1]);
165 if (dphi < -90)
167 dphi += 360;
169 pt += dphi;
170 a0 = a1;
173 return (pt/(nca-1));
176 real ca_phi(int gnx, atom_id index[], rvec x[])
177 /* Assume we have a list of Calpha atoms only! */
179 real phi, phitot;
180 int i, ai, aj, ak, al, t1, t2, t3;
181 rvec r_ij, r_kj, r_kl, m, n;
182 real sign;
184 if (gnx <= 4)
186 return 0;
189 phitot = 0;
190 for (i = 0; (i < gnx-4); i++)
192 ai = index[i+0];
193 aj = index[i+1];
194 ak = index[i+2];
195 al = index[i+3];
196 phi = RAD2DEG*
197 dih_angle(x[ai], x[aj], x[ak], x[al], NULL,
198 r_ij, r_kj, r_kl, m, n,
199 &sign, &t1, &t2, &t3);
200 phitot += phi;
203 return (phitot/(gnx-4.0));
206 real dip(int nbb, atom_id bbind[], rvec x[], t_atom atom[])
208 int i, m, ai;
209 rvec dipje;
210 real q;
212 clear_rvec(dipje);
213 for (i = 0; (i < nbb); i++)
215 ai = bbind[i];
216 q = atom[ai].q;
217 for (m = 0; (m < DIM); m++)
219 dipje[m] += x[ai][m]*q;
222 return norm(dipje);
225 real rise(int gnx, atom_id index[], rvec x[])
226 /* Assume we have a list of Calpha atoms only! */
228 real z, z0, ztot;
229 int i, ai;
231 ai = index[0];
232 z0 = x[ai][ZZ];
233 ztot = 0;
234 for (i = 1; (i < gnx); i++)
236 ai = index[i];
237 z = x[ai][ZZ];
238 ztot += (z-z0);
239 z0 = z;
242 return (ztot/(gnx-1.0));
245 void av_hblen(FILE *fp3, FILE *fp3a,
246 FILE *fp4, FILE *fp4a,
247 FILE *fp5, FILE *fp5a,
248 real t, int nres, t_bb bb[])
250 int i, n3 = 0, n4 = 0, n5 = 0;
251 real d3 = 0, d4 = 0, d5 = 0;
253 for (i = 0; (i < nres-3); i++)
255 if (bb[i].bHelix)
257 fprintf(fp3a, "%10g", bb[i].d3);
258 n3++;
259 d3 += bb[i].d3;
260 if (i < nres-4)
262 fprintf(fp4a, "%10g", bb[i].d4);
263 n4++;
264 d4 += bb[i].d4;
266 if (i < nres-5)
268 fprintf(fp5a, "%10g", bb[i].d5);
269 n5++;
270 d5 += bb[i].d5;
274 fprintf(fp3, "%10g %10g\n", t, d3/n3);
275 fprintf(fp4, "%10g %10g\n", t, d4/n4);
276 fprintf(fp5, "%10g %10g\n", t, d5/n5);
277 fprintf(fp3a, "\n");
278 fprintf(fp4a, "\n");
279 fprintf(fp5a, "\n");
283 void av_phipsi(FILE *fphi, FILE *fpsi, FILE *fphi2, FILE *fpsi2,
284 real t, int nres, t_bb bb[])
286 int i, n = 0;
287 real phi = 0, psi = 0;
289 fprintf(fphi2, "%10g", t);
290 fprintf(fpsi2, "%10g", t);
291 for (i = 0; (i < nres); i++)
293 if (bb[i].bHelix)
295 phi += bb[i].phi;
296 psi += bb[i].psi;
297 fprintf(fphi2, " %10g", bb[i].phi);
298 fprintf(fpsi2, " %10g", bb[i].psi);
299 n++;
302 fprintf(fphi, "%10g %10g\n", t, (phi/n));
303 fprintf(fpsi, "%10g %10g\n", t, (psi/n));
304 fprintf(fphi2, "\n");
305 fprintf(fpsi2, "\n");
308 static void set_ahcity(int nbb, t_bb bb[])
310 real pp2;
311 int n;
313 for (n = 0; (n < nbb); n++)
315 pp2 = sqr(bb[n].phi-PHI_AHX)+sqr(bb[n].psi-PSI_AHX);
317 bb[n].bHelix = FALSE;
318 if (pp2 < 2500)
320 if ((bb[n].d4 < 0.36) || ((n > 0) && bb[n-1].bHelix))
322 bb[n].bHelix = TRUE;
328 t_bb *mkbbind(const char *fn, int *nres, int *nbb, int res0,
329 int *nall, atom_id **index,
330 char ***atomname, t_atom atom[],
331 t_resinfo *resinfo)
333 static const char * bb_nm[] = { "N", "H", "CA", "C", "O", "HN" };
334 #define NBB asize(bb_nm)
335 t_bb *bb;
336 char *grpname;
337 int ai, i, i0, i1, j, k, ri, rnr, gnx, r0, r1;
339 fprintf(stderr, "Please select a group containing the entire backbone\n");
340 rd_index(fn, 1, &gnx, index, &grpname);
341 *nall = gnx;
342 fprintf(stderr, "Checking group %s\n", grpname);
343 r0 = r1 = atom[(*index)[0]].resind;
344 for (i = 1; (i < gnx); i++)
346 r0 = min(r0, atom[(*index)[i]].resind);
347 r1 = max(r1, atom[(*index)[i]].resind);
349 rnr = r1-r0+1;
350 fprintf(stderr, "There are %d residues\n", rnr);
351 snew(bb, rnr);
352 for (i = 0; (i < rnr); i++)
354 bb[i].N = bb[i].H = bb[i].CA = bb[i].C = bb[i].O = -1, bb[i].resno = res0+i;
357 for (i = j = 0; (i < gnx); i++)
359 ai = (*index)[i];
360 ri = atom[ai].resind-r0;
361 if (strcmp(*(resinfo[ri].name), "PRO") == 0)
363 if (strcmp(*(atomname[ai]), "CD") == 0)
365 bb[ri].H = ai;
368 for (k = 0; (k < NBB); k++)
370 if (strcmp(bb_nm[k], *(atomname[ai])) == 0)
372 j++;
373 break;
376 switch (k)
378 case 0:
379 bb[ri].N = ai;
380 break;
381 case 1:
382 case 5:
383 /* No attempt to address the case where some weird input has both H and HN atoms in the group */
384 bb[ri].H = ai;
385 break;
386 case 2:
387 bb[ri].CA = ai;
388 break;
389 case 3:
390 bb[ri].C = ai;
391 break;
392 case 4:
393 bb[ri].O = ai;
394 break;
395 default:
396 break;
400 for (i0 = 0; (i0 < rnr); i0++)
402 if ((bb[i0].N != -1) && (bb[i0].H != -1) &&
403 (bb[i0].CA != -1) &&
404 (bb[i0].C != -1) && (bb[i0].O != -1))
406 break;
409 for (i1 = rnr-1; (i1 >= 0); i1--)
411 if ((bb[i1].N != -1) && (bb[i1].H != -1) &&
412 (bb[i1].CA != -1) &&
413 (bb[i1].C != -1) && (bb[i1].O != -1))
415 break;
418 if (i0 == 0)
420 i0++;
422 if (i1 == rnr-1)
424 i1--;
427 for (i = i0; (i < i1); i++)
429 bb[i].Cprev = bb[i-1].C;
430 bb[i].Nnext = bb[i+1].N;
432 rnr = max(0, i1-i0+1);
433 fprintf(stderr, "There are %d complete backbone residues (from %d to %d)\n",
434 rnr, bb[i0].resno, bb[i1].resno);
435 if (rnr == 0)
437 gmx_fatal(FARGS, "Zero complete backbone residues were found, cannot proceed");
439 for (i = 0; (i < rnr); i++, i0++)
441 bb[i] = bb[i0];
444 /* Set the labels */
445 for (i = 0; (i < rnr); i++)
447 ri = atom[bb[i].CA].resind;
448 sprintf(bb[i].label, "%s%d", *(resinfo[ri].name), ri+res0);
451 *nres = rnr;
452 *nbb = rnr*asize(bb_nm);
454 return bb;
457 real pprms(FILE *fp, int nbb, t_bb bb[])
459 int i, n;
460 real rms, rmst, rms2;
462 rmst = rms2 = 0;
463 for (i = n = 0; (i < nbb); i++)
465 if (bb[i].bHelix)
467 rms = sqrt(bb[i].pprms2);
468 rmst += rms;
469 rms2 += bb[i].pprms2;
470 fprintf(fp, "%10g ", rms);
471 n++;
474 fprintf(fp, "\n");
475 rms = sqrt(rms2/n-sqr(rmst/n));
477 return rms;
480 void calc_hxprops(int nres, t_bb bb[], rvec x[])
482 int i, ao, an, t1, t2, t3;
483 rvec dx, r_ij, r_kj, r_kl, m, n;
484 real sign;
486 for (i = 0; (i < nres); i++)
488 ao = bb[i].O;
489 bb[i].d4 = bb[i].d3 = bb[i].d5 = 0;
490 if (i < nres-3)
492 an = bb[i+3].N;
493 rvec_sub(x[ao], x[an], dx);
494 bb[i].d3 = norm(dx);
496 if (i < nres-4)
498 an = bb[i+4].N;
499 rvec_sub(x[ao], x[an], dx);
500 bb[i].d4 = norm(dx);
502 if (i < nres-5)
504 an = bb[i+5].N;
505 rvec_sub(x[ao], x[an], dx);
506 bb[i].d5 = norm(dx);
509 bb[i].phi = RAD2DEG*
510 dih_angle(x[bb[i].Cprev], x[bb[i].N], x[bb[i].CA], x[bb[i].C], NULL,
511 r_ij, r_kj, r_kl, m, n,
512 &sign, &t1, &t2, &t3);
513 bb[i].psi = RAD2DEG*
514 dih_angle(x[bb[i].N], x[bb[i].CA], x[bb[i].C], x[bb[i].Nnext], NULL,
515 r_ij, r_kj, r_kl, m, n,
516 &sign, &t1, &t2, &t3);
517 bb[i].pprms2 = sqr(bb[i].phi-PHI_AHX)+sqr(bb[i].psi-PSI_AHX);
519 bb[i].jcaha +=
520 1.4*sin((bb[i].psi+138.0)*DEG2RAD) -
521 4.1*cos(2.0*DEG2RAD*(bb[i].psi+138.0)) +
522 2.0*cos(2.0*DEG2RAD*(bb[i].phi+30.0));
526 static void check_ahx(int nres, t_bb bb[],
527 int *hstart, int *hend)
529 int h0, h1, h0sav, h1sav;
531 set_ahcity(nres, bb);
532 h0 = h0sav = h1sav = 0;
535 for (; (!bb[h0].bHelix) && (h0 < nres-4); h0++)
539 for (h1 = h0; bb[h1+1].bHelix && (h1 < nres-1); h1++)
543 if (h1 > h0)
545 /*fprintf(stderr,"Helix from %d to %d\n",h0,h1);*/
546 if (h1-h0 > h1sav-h0sav)
548 h0sav = h0;
549 h1sav = h1;
552 h0 = h1+1;
554 while (h1 < nres-1);
555 *hstart = h0sav;
556 *hend = h1sav;
559 void do_start_end(int nres, t_bb bb[], int *nbb, atom_id bbindex[],
560 int *nca, atom_id caindex[],
561 gmx_bool bRange, int rStart, int rEnd)
563 int i, j, hstart = 0, hend = 0;
565 if (bRange)
567 for (i = 0; (i < nres); i++)
569 if ((bb[i].resno >= rStart) && (bb[i].resno <= rEnd))
571 bb[i].bHelix = TRUE;
573 if (bb[i].resno == rStart)
575 hstart = i;
577 if (bb[i].resno == rEnd)
579 hend = i;
583 else
585 /* Find start and end of longest helix fragment */
586 check_ahx(nres, bb, &hstart, &hend);
588 fprintf(stderr, "helix from: %d through %d\n",
589 bb[hstart].resno, bb[hend].resno);
591 for (j = 0, i = hstart; (i <= hend); i++)
593 bbindex[j++] = bb[i].N;
594 bbindex[j++] = bb[i].H;
595 bbindex[j++] = bb[i].CA;
596 bbindex[j++] = bb[i].C;
597 bbindex[j++] = bb[i].O;
598 caindex[i-hstart] = bb[i].CA;
600 *nbb = j;
601 *nca = (hend-hstart+1);
604 void pr_bb(FILE *fp, int nres, t_bb bb[])
606 int i;
608 fprintf(fp, "\n");
609 fprintf(fp, "%3s %3s %3s %3s %3s %7s %7s %7s %7s %7s %3s\n",
610 "AA", "N", "Ca", "C", "O", "Phi", "Psi", "D3", "D4", "D5", "Hx?");
611 for (i = 0; (i < nres); i++)
613 fprintf(fp, "%3d %3d %3d %3d %3d %7.2f %7.2f %7.3f %7.3f %7.3f %3s\n",
614 bb[i].resno, bb[i].N, bb[i].CA, bb[i].C, bb[i].O,
615 bb[i].phi, bb[i].psi, bb[i].d3, bb[i].d4, bb[i].d5,
616 bb[i].bHelix ? "Yes" : "No");
618 fprintf(fp, "\n");