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,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
47 #include "gromacs/listed_forces/bonded.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/topology/index.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/arraysize.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/smalloc.h"
57 real
ellipticity(int nres
, t_bb bb
[])
63 // Avoid warnings about narrowing conversions from double to real
65 # pragma warning(disable : 4838)
67 static const t_ppwstr ppw
[] = { { -67, -44, 0.31 }, { -66, -41, 0.31 }, { -59, -44, 0.44 },
68 { -57, -47, 0.56 }, { -53, -52, 0.78 }, { -48, -57, 1.00 },
69 { -70.5, -35.8, 0.15 }, { -57, -79, 0.23 }, { -38, -78, 1.20 },
70 { -60, -30, 0.24 }, { -54, -28, 0.46 }, { -44, -33, 0.68 } };
72 # pragma warning(default : 4838)
74 #define NPPW asize(ppw)
77 real ell
, pp2
, phi
, psi
;
80 for (i
= 0; (i
< nres
); i
++)
84 for (j
= 0; (j
< NPPW
); j
++)
86 pp2
= gmx::square(phi
- ppw
[j
].phi
) + gmx::square(psi
- ppw
[j
].psi
);
98 real
ahx_len(int gnx
, const int index
[], rvec x
[])
99 /* Assume we have a list of Calpha atoms only! */
103 rvec_sub(x
[index
[0]], x
[index
[gnx
- 1]], dx
);
108 real
radius(FILE* fp
, int nca
, const int ca_index
[], rvec x
[])
109 /* Assume we have all the backbone */
115 for (i
= 0; (i
< nca
); i
++)
118 dl2
= gmx::square(x
[ai
][XX
]) + gmx::square(x
[ai
][YY
]);
122 fprintf(fp
, " %10g", dl2
);
132 return std::sqrt(dlt
/ nca
);
135 static real
rot(rvec x1
, const rvec x2
)
137 real phi1
, dphi
, cp
, sp
;
140 phi1
= std::atan2(x1
[YY
], x1
[XX
]);
143 xx
= cp
* x2
[XX
] + sp
* x2
[YY
];
144 yy
= -sp
* x2
[XX
] + cp
* x2
[YY
];
146 dphi
= RAD2DEG
* std::atan2(yy
, xx
);
151 real
twist(int nca
, const int caindex
[], rvec x
[])
158 for (i
= 1; (i
< nca
); i
++)
162 dphi
= rot(x
[a0
], x
[a1
]);
171 return (pt
/ (nca
- 1));
174 real
ca_phi(int gnx
, const int index
[], rvec x
[])
175 /* Assume we have a list of Calpha atoms only! */
178 int i
, ai
, aj
, ak
, al
, t1
, t2
, t3
;
179 rvec r_ij
, r_kj
, r_kl
, m
, n
;
187 for (i
= 0; (i
< gnx
- 4); i
++)
194 * dih_angle(x
[ai
], x
[aj
], x
[ak
], x
[al
], nullptr, r_ij
, r_kj
, r_kl
, m
, n
, &t1
, &t2
, &t3
);
198 return (phitot
/ (gnx
- 4.0));
201 real
dip(int nbb
, int const bbind
[], const rvec x
[], const t_atom atom
[])
208 for (i
= 0; (i
< nbb
); i
++)
212 for (m
= 0; (m
< DIM
); m
++)
214 dipje
[m
] += x
[ai
][m
] * q
;
220 real
rise(int gnx
, const int index
[], rvec x
[])
221 /* Assume we have a list of Calpha atoms only! */
229 for (i
= 1; (i
< gnx
); i
++)
237 return (ztot
/ (gnx
- 1.0));
240 void av_hblen(FILE* fp3
, FILE* fp3a
, FILE* fp4
, FILE* fp4a
, FILE* fp5
, FILE* fp5a
, real t
, int nres
, t_bb bb
[])
242 int i
, n3
= 0, n4
= 0, n5
= 0;
243 real d3
= 0, d4
= 0, d5
= 0;
245 for (i
= 0; (i
< nres
- 3); i
++)
249 fprintf(fp3a
, "%10g", bb
[i
].d3
);
254 fprintf(fp4a
, "%10g", bb
[i
].d4
);
260 fprintf(fp5a
, "%10g", bb
[i
].d5
);
266 fprintf(fp3
, "%10g %10g\n", t
, d3
/ n3
);
267 fprintf(fp4
, "%10g %10g\n", t
, d4
/ n4
);
268 fprintf(fp5
, "%10g %10g\n", t
, d5
/ n5
);
274 void av_phipsi(FILE* fphi
, FILE* fpsi
, FILE* fphi2
, FILE* fpsi2
, real t
, int nres
, t_bb bb
[])
277 real phi
= 0, psi
= 0;
279 fprintf(fphi2
, "%10g", t
);
280 fprintf(fpsi2
, "%10g", t
);
281 for (i
= 0; (i
< nres
); i
++)
287 fprintf(fphi2
, " %10g", bb
[i
].phi
);
288 fprintf(fpsi2
, " %10g", bb
[i
].psi
);
292 fprintf(fphi
, "%10g %10g\n", t
, (phi
/ n
));
293 fprintf(fpsi
, "%10g %10g\n", t
, (psi
/ n
));
294 fprintf(fphi2
, "\n");
295 fprintf(fpsi2
, "\n");
298 static void set_ahcity(int nbb
, t_bb bb
[])
303 for (n
= 0; (n
< nbb
); n
++)
305 pp2
= gmx::square(bb
[n
].phi
- PHI_AHX
) + gmx::square(bb
[n
].psi
- PSI_AHX
);
307 bb
[n
].bHelix
= FALSE
;
310 if ((bb
[n
].d4
< 0.36) || ((n
> 0) && bb
[n
- 1].bHelix
))
318 t_bb
* mkbbind(const char* fn
,
328 static const char* bb_nm
[] = { "N", "H", "CA", "C", "O", "HN" };
329 #define NBB asize(bb_nm)
332 int ai
, i
, i0
, i1
, j
, k
, rnr
, gnx
, r0
, r1
;
334 fprintf(stderr
, "Please select a group containing the entire backbone\n");
335 rd_index(fn
, 1, &gnx
, index
, &grpname
);
337 fprintf(stderr
, "Checking group %s\n", grpname
);
338 r0
= r1
= atom
[(*index
)[0]].resind
;
339 for (i
= 1; (i
< gnx
); i
++)
341 r0
= std::min(r0
, atom
[(*index
)[i
]].resind
);
342 r1
= std::max(r1
, atom
[(*index
)[i
]].resind
);
345 fprintf(stderr
, "There are %d residues\n", rnr
);
347 for (i
= 0; (i
< rnr
); i
++)
349 bb
[i
].N
= bb
[i
].H
= bb
[i
].CA
= bb
[i
].C
= bb
[i
].O
= -1;
350 bb
[i
].resno
= res0
+ i
;
353 for (i
= j
= 0; (i
< gnx
); i
++)
356 // Create an index into the residue index for the topology.
357 int resindex
= atom
[ai
].resind
;
358 // Create an index into the residues present in the selected
360 int bbindex
= resindex
- r0
;
361 if (std::strcmp(*(resinfo
[resindex
].name
), "PRO") == 0)
363 // For PRO in a peptide, there is no H bound to backbone
364 // N, so use CD instead.
365 if (std::strcmp(*(atomname
[ai
]), "CD") == 0)
370 for (k
= 0; (k
< NBB
); k
++)
372 if (std::strcmp(bb_nm
[k
], *(atomname
[ai
])) == 0)
379 case 0: bb
[bbindex
].N
= ai
; break;
382 /* No attempt to address the case where some weird input has both H and HN atoms in the group */
385 case 2: bb
[bbindex
].CA
= ai
; break;
386 case 3: bb
[bbindex
].C
= ai
; break;
387 case 4: bb
[bbindex
].O
= ai
; break;
392 for (i0
= 0; (i0
< rnr
); i0
++)
394 if ((bb
[i0
].N
!= -1) && (bb
[i0
].H
!= -1) && (bb
[i0
].CA
!= -1) && (bb
[i0
].C
!= -1)
400 for (i1
= rnr
- 1; (i1
>= 0); i1
--)
402 if ((bb
[i1
].N
!= -1) && (bb
[i1
].H
!= -1) && (bb
[i1
].CA
!= -1) && (bb
[i1
].C
!= -1)
417 for (i
= i0
; (i
< i1
); i
++)
419 bb
[i
].Cprev
= bb
[i
- 1].C
;
420 bb
[i
].Nnext
= bb
[i
+ 1].N
;
422 rnr
= std::max(0, i1
- i0
+ 1);
423 fprintf(stderr
, "There are %d complete backbone residues (from %d to %d)\n", rnr
, bb
[i0
].resno
,
427 gmx_fatal(FARGS
, "Zero complete backbone residues were found, cannot proceed");
429 for (i
= 0; (i
< rnr
); i
++, i0
++)
435 for (i
= 0; (i
< rnr
); i
++)
437 int resindex
= atom
[bb
[i
].CA
].resind
;
438 sprintf(bb
[i
].label
, "%s%d", *(resinfo
[resindex
].name
), resinfo
[resindex
].nr
);
442 *nbb
= rnr
* asize(bb_nm
);
447 real
pprms(FILE* fp
, int nbb
, t_bb bb
[])
450 real rms
, rmst
, rms2
;
453 for (i
= n
= 0; (i
< nbb
); i
++)
457 rms
= std::sqrt(bb
[i
].pprms2
);
459 rms2
+= bb
[i
].pprms2
;
460 fprintf(fp
, "%10g ", rms
);
465 rms
= std::sqrt(rms2
/ n
- gmx::square(rmst
/ n
));
470 void calc_hxprops(int nres
, t_bb bb
[], const rvec x
[])
472 int i
, ao
, an
, t1
, t2
, t3
;
473 rvec dx
, r_ij
, r_kj
, r_kl
, m
, n
;
475 for (i
= 0; (i
< nres
); i
++)
478 bb
[i
].d4
= bb
[i
].d3
= bb
[i
].d5
= 0;
482 rvec_sub(x
[ao
], x
[an
], dx
);
488 rvec_sub(x
[ao
], x
[an
], dx
);
494 rvec_sub(x
[ao
], x
[an
], dx
);
499 * dih_angle(x
[bb
[i
].Cprev
], x
[bb
[i
].N
], x
[bb
[i
].CA
], x
[bb
[i
].C
], nullptr, r_ij
,
500 r_kj
, r_kl
, m
, n
, &t1
, &t2
, &t3
);
502 * dih_angle(x
[bb
[i
].N
], x
[bb
[i
].CA
], x
[bb
[i
].C
], x
[bb
[i
].Nnext
], nullptr, r_ij
,
503 r_kj
, r_kl
, m
, n
, &t1
, &t2
, &t3
);
504 bb
[i
].pprms2
= gmx::square(bb
[i
].phi
- PHI_AHX
) + gmx::square(bb
[i
].psi
- PSI_AHX
);
506 bb
[i
].jcaha
+= 1.4 * std::sin((bb
[i
].psi
+ 138.0) * DEG2RAD
)
507 - 4.1 * std::cos(2.0 * DEG2RAD
* (bb
[i
].psi
+ 138.0))
508 + 2.0 * std::cos(2.0 * DEG2RAD
* (bb
[i
].phi
+ 30.0));
512 static void check_ahx(int nres
, t_bb bb
[], int* hstart
, int* hend
)
514 int h0
, h1
, h0sav
, h1sav
;
516 set_ahcity(nres
, bb
);
517 h0
= h0sav
= h1sav
= 0;
520 for (; (!bb
[h0
].bHelix
) && (h0
< nres
- 4); h0
++) {}
521 for (h1
= h0
; bb
[h1
+ 1].bHelix
&& (h1
< nres
- 1); h1
++) {}
524 /*fprintf(stderr,"Helix from %d to %d\n",h0,h1);*/
525 if (h1
- h0
> h1sav
- h0sav
)
532 } while (h1
< nres
- 1);
537 void do_start_end(int nres
, t_bb bb
[], int* nbb
, int bbindex
[], int* nca
, int caindex
[], gmx_bool bRange
, int rStart
, int rEnd
)
539 int i
, j
, hstart
= 0, hend
= 0;
543 for (i
= 0; (i
< nres
); i
++)
545 if ((bb
[i
].resno
>= rStart
) && (bb
[i
].resno
<= rEnd
))
549 if (bb
[i
].resno
== rStart
)
553 if (bb
[i
].resno
== rEnd
)
561 /* Find start and end of longest helix fragment */
562 check_ahx(nres
, bb
, &hstart
, &hend
);
564 fprintf(stderr
, "helix from: %d through %d\n", bb
[hstart
].resno
, bb
[hend
].resno
);
566 for (j
= 0, i
= hstart
; (i
<= hend
); i
++)
568 bbindex
[j
++] = bb
[i
].N
;
569 bbindex
[j
++] = bb
[i
].H
;
570 bbindex
[j
++] = bb
[i
].CA
;
571 bbindex
[j
++] = bb
[i
].C
;
572 bbindex
[j
++] = bb
[i
].O
;
573 caindex
[i
- hstart
] = bb
[i
].CA
;
576 *nca
= (hend
- hstart
+ 1);
579 void pr_bb(FILE* fp
, int nres
, t_bb bb
[])
584 fprintf(fp
, "%3s %3s %3s %3s %3s %7s %7s %7s %7s %7s %3s\n", "AA", "N", "Ca", "C", "O", "Phi",
585 "Psi", "D3", "D4", "D5", "Hx?");
586 for (i
= 0; (i
< nres
); i
++)
588 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
,
589 bb
[i
].CA
, bb
[i
].C
, bb
[i
].O
, bb
[i
].phi
, bb
[i
].psi
, bb
[i
].d3
, bb
[i
].d4
, bb
[i
].d5
,
590 bb
[i
].bHelix
? "Yes" : "No");