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,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.
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.
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"
56 real
ellipticity(int nres
, t_bb bb
[])
61 // Avoid warnings about narrowing conversions from double to real
63 #pragma warning(disable: 4838)
65 static const t_ppwstr ppw
[] = {
72 { -70.5, -35.8, 0.15 },
80 #pragma warning(default: 4838)
82 #define NPPW asize(ppw)
85 real ell
, pp2
, phi
, psi
;
88 for (i
= 0; (i
< nres
); i
++)
92 for (j
= 0; (j
< NPPW
); j
++)
94 pp2
= gmx::square(phi
-ppw
[j
].phi
)+gmx::square(psi
-ppw
[j
].psi
);
106 real
ahx_len(int gnx
, const int index
[], rvec x
[])
107 /* Assume we have a list of Calpha atoms only! */
111 rvec_sub(x
[index
[0]], x
[index
[gnx
-1]], dx
);
116 real
radius(FILE *fp
, int nca
, const int ca_index
[], rvec x
[])
117 /* Assume we have all the backbone */
123 for (i
= 0; (i
< nca
); i
++)
126 dl2
= gmx::square(x
[ai
][XX
])+gmx::square(x
[ai
][YY
]);
130 fprintf(fp
, " %10g", dl2
);
140 return std::sqrt(dlt
/nca
);
143 static real
rot(rvec x1
, const rvec x2
)
145 real phi1
, dphi
, cp
, sp
;
148 phi1
= std::atan2(x1
[YY
], x1
[XX
]);
151 xx
= cp
*x2
[XX
]+sp
*x2
[YY
];
152 yy
= -sp
*x2
[XX
]+cp
*x2
[YY
];
154 dphi
= RAD2DEG
*std::atan2(yy
, xx
);
159 real
twist(int nca
, const int caindex
[], rvec x
[])
166 for (i
= 1; (i
< nca
); i
++)
170 dphi
= rot(x
[a0
], x
[a1
]);
182 real
ca_phi(int gnx
, const int index
[], rvec x
[])
183 /* Assume we have a list of Calpha atoms only! */
186 int i
, ai
, aj
, ak
, al
, t1
, t2
, t3
;
187 rvec r_ij
, r_kj
, r_kl
, m
, n
;
195 for (i
= 0; (i
< gnx
-4); i
++)
202 dih_angle(x
[ai
], x
[aj
], x
[ak
], x
[al
], nullptr,
203 r_ij
, r_kj
, r_kl
, m
, n
,
208 return (phitot
/(gnx
-4.0));
211 real
dip(int nbb
, int const bbind
[], const rvec x
[], const t_atom atom
[])
218 for (i
= 0; (i
< nbb
); i
++)
222 for (m
= 0; (m
< DIM
); m
++)
224 dipje
[m
] += x
[ai
][m
]*q
;
230 real
rise(int gnx
, const int index
[], rvec x
[])
231 /* Assume we have a list of Calpha atoms only! */
239 for (i
= 1; (i
< gnx
); i
++)
247 return (ztot
/(gnx
-1.0));
250 void av_hblen(FILE *fp3
, FILE *fp3a
,
251 FILE *fp4
, FILE *fp4a
,
252 FILE *fp5
, FILE *fp5a
,
253 real t
, int nres
, t_bb bb
[])
255 int i
, n3
= 0, n4
= 0, n5
= 0;
256 real d3
= 0, d4
= 0, d5
= 0;
258 for (i
= 0; (i
< nres
-3); i
++)
262 fprintf(fp3a
, "%10g", bb
[i
].d3
);
267 fprintf(fp4a
, "%10g", bb
[i
].d4
);
273 fprintf(fp5a
, "%10g", bb
[i
].d5
);
279 fprintf(fp3
, "%10g %10g\n", t
, d3
/n3
);
280 fprintf(fp4
, "%10g %10g\n", t
, d4
/n4
);
281 fprintf(fp5
, "%10g %10g\n", t
, d5
/n5
);
288 void av_phipsi(FILE *fphi
, FILE *fpsi
, FILE *fphi2
, FILE *fpsi2
,
289 real t
, int nres
, t_bb bb
[])
292 real phi
= 0, psi
= 0;
294 fprintf(fphi2
, "%10g", t
);
295 fprintf(fpsi2
, "%10g", t
);
296 for (i
= 0; (i
< nres
); i
++)
302 fprintf(fphi2
, " %10g", bb
[i
].phi
);
303 fprintf(fpsi2
, " %10g", bb
[i
].psi
);
307 fprintf(fphi
, "%10g %10g\n", t
, (phi
/n
));
308 fprintf(fpsi
, "%10g %10g\n", t
, (psi
/n
));
309 fprintf(fphi2
, "\n");
310 fprintf(fpsi2
, "\n");
313 static void set_ahcity(int nbb
, t_bb bb
[])
318 for (n
= 0; (n
< nbb
); n
++)
320 pp2
= gmx::square(bb
[n
].phi
-PHI_AHX
)+gmx::square(bb
[n
].psi
-PSI_AHX
);
322 bb
[n
].bHelix
= FALSE
;
325 if ((bb
[n
].d4
< 0.36) || ((n
> 0) && bb
[n
-1].bHelix
))
333 t_bb
*mkbbind(const char *fn
, int *nres
, int *nbb
, int res0
,
334 int *nall
, int **index
,
335 char ***atomname
, t_atom atom
[],
338 static const char * bb_nm
[] = { "N", "H", "CA", "C", "O", "HN" };
339 #define NBB asize(bb_nm)
342 int ai
, i
, i0
, i1
, j
, k
, rnr
, gnx
, r0
, r1
;
344 fprintf(stderr
, "Please select a group containing the entire backbone\n");
345 rd_index(fn
, 1, &gnx
, index
, &grpname
);
347 fprintf(stderr
, "Checking group %s\n", grpname
);
348 r0
= r1
= atom
[(*index
)[0]].resind
;
349 for (i
= 1; (i
< gnx
); i
++)
351 r0
= std::min(r0
, atom
[(*index
)[i
]].resind
);
352 r1
= std::max(r1
, atom
[(*index
)[i
]].resind
);
355 fprintf(stderr
, "There are %d residues\n", rnr
);
357 for (i
= 0; (i
< rnr
); i
++)
359 bb
[i
].N
= bb
[i
].H
= bb
[i
].CA
= bb
[i
].C
= bb
[i
].O
= -1; bb
[i
].resno
= res0
+i
;
362 for (i
= j
= 0; (i
< gnx
); i
++)
365 // Create an index into the residue index for the topology.
366 int resindex
= atom
[ai
].resind
;
367 // Create an index into the residues present in the selected
369 int bbindex
= resindex
-r0
;
370 if (std::strcmp(*(resinfo
[resindex
].name
), "PRO") == 0)
372 // For PRO in a peptide, there is no H bound to backbone
373 // N, so use CD instead.
374 if (std::strcmp(*(atomname
[ai
]), "CD") == 0)
379 for (k
= 0; (k
< NBB
); k
++)
381 if (std::strcmp(bb_nm
[k
], *(atomname
[ai
])) == 0)
393 /* No attempt to address the case where some weird input has both H and HN atoms in the group */
410 for (i0
= 0; (i0
< rnr
); i0
++)
412 if ((bb
[i0
].N
!= -1) && (bb
[i0
].H
!= -1) &&
414 (bb
[i0
].C
!= -1) && (bb
[i0
].O
!= -1))
419 for (i1
= rnr
-1; (i1
>= 0); i1
--)
421 if ((bb
[i1
].N
!= -1) && (bb
[i1
].H
!= -1) &&
423 (bb
[i1
].C
!= -1) && (bb
[i1
].O
!= -1))
437 for (i
= i0
; (i
< i1
); i
++)
439 bb
[i
].Cprev
= bb
[i
-1].C
;
440 bb
[i
].Nnext
= bb
[i
+1].N
;
442 rnr
= std::max(0, i1
-i0
+1);
443 fprintf(stderr
, "There are %d complete backbone residues (from %d to %d)\n",
444 rnr
, bb
[i0
].resno
, bb
[i1
].resno
);
447 gmx_fatal(FARGS
, "Zero complete backbone residues were found, cannot proceed");
449 for (i
= 0; (i
< rnr
); i
++, i0
++)
455 for (i
= 0; (i
< rnr
); i
++)
457 int resindex
= atom
[bb
[i
].CA
].resind
;
458 sprintf(bb
[i
].label
, "%s%d", *(resinfo
[resindex
].name
), resinfo
[resindex
].nr
);
462 *nbb
= rnr
*asize(bb_nm
);
467 real
pprms(FILE *fp
, int nbb
, t_bb bb
[])
470 real rms
, rmst
, rms2
;
473 for (i
= n
= 0; (i
< nbb
); i
++)
477 rms
= std::sqrt(bb
[i
].pprms2
);
479 rms2
+= bb
[i
].pprms2
;
480 fprintf(fp
, "%10g ", rms
);
485 rms
= std::sqrt(rms2
/n
-gmx::square(rmst
/n
));
490 void calc_hxprops(int nres
, t_bb bb
[], const rvec x
[])
492 int i
, ao
, an
, t1
, t2
, t3
;
493 rvec dx
, r_ij
, r_kj
, r_kl
, m
, n
;
495 for (i
= 0; (i
< nres
); i
++)
498 bb
[i
].d4
= bb
[i
].d3
= bb
[i
].d5
= 0;
502 rvec_sub(x
[ao
], x
[an
], dx
);
508 rvec_sub(x
[ao
], x
[an
], dx
);
514 rvec_sub(x
[ao
], x
[an
], dx
);
519 dih_angle(x
[bb
[i
].Cprev
], x
[bb
[i
].N
], x
[bb
[i
].CA
], x
[bb
[i
].C
], nullptr,
520 r_ij
, r_kj
, r_kl
, m
, n
,
523 dih_angle(x
[bb
[i
].N
], x
[bb
[i
].CA
], x
[bb
[i
].C
], x
[bb
[i
].Nnext
], nullptr,
524 r_ij
, r_kj
, r_kl
, m
, n
,
526 bb
[i
].pprms2
= gmx::square(bb
[i
].phi
-PHI_AHX
)+gmx::square(bb
[i
].psi
-PSI_AHX
);
529 1.4*std::sin((bb
[i
].psi
+138.0)*DEG2RAD
) -
530 4.1*std::cos(2.0*DEG2RAD
*(bb
[i
].psi
+138.0)) +
531 2.0*std::cos(2.0*DEG2RAD
*(bb
[i
].phi
+30.0));
535 static void check_ahx(int nres
, t_bb bb
[],
536 int *hstart
, int *hend
)
538 int h0
, h1
, h0sav
, h1sav
;
540 set_ahcity(nres
, bb
);
541 h0
= h0sav
= h1sav
= 0;
544 for (; (!bb
[h0
].bHelix
) && (h0
< nres
-4); h0
++)
548 for (h1
= h0
; bb
[h1
+1].bHelix
&& (h1
< nres
-1); h1
++)
554 /*fprintf(stderr,"Helix from %d to %d\n",h0,h1);*/
555 if (h1
-h0
> h1sav
-h0sav
)
568 void do_start_end(int nres
, t_bb bb
[], int *nbb
, int bbindex
[],
569 int *nca
, int caindex
[],
570 gmx_bool bRange
, int rStart
, int rEnd
)
572 int i
, j
, hstart
= 0, hend
= 0;
576 for (i
= 0; (i
< nres
); i
++)
578 if ((bb
[i
].resno
>= rStart
) && (bb
[i
].resno
<= rEnd
))
582 if (bb
[i
].resno
== rStart
)
586 if (bb
[i
].resno
== rEnd
)
594 /* Find start and end of longest helix fragment */
595 check_ahx(nres
, bb
, &hstart
, &hend
);
597 fprintf(stderr
, "helix from: %d through %d\n",
598 bb
[hstart
].resno
, bb
[hend
].resno
);
600 for (j
= 0, i
= hstart
; (i
<= hend
); i
++)
602 bbindex
[j
++] = bb
[i
].N
;
603 bbindex
[j
++] = bb
[i
].H
;
604 bbindex
[j
++] = bb
[i
].CA
;
605 bbindex
[j
++] = bb
[i
].C
;
606 bbindex
[j
++] = bb
[i
].O
;
607 caindex
[i
-hstart
] = bb
[i
].CA
;
610 *nca
= (hend
-hstart
+1);
613 void pr_bb(FILE *fp
, int nres
, t_bb bb
[])
618 fprintf(fp
, "%3s %3s %3s %3s %3s %7s %7s %7s %7s %7s %3s\n",
619 "AA", "N", "Ca", "C", "O", "Phi", "Psi", "D3", "D4", "D5", "Hx?");
620 for (i
= 0; (i
< nres
); i
++)
622 fprintf(fp
, "%3d %3d %3d %3d %3d %7.2f %7.2f %7.3f %7.3f %7.3f %3s\n",
623 bb
[i
].resno
, bb
[i
].N
, bb
[i
].CA
, bb
[i
].C
, bb
[i
].O
,
624 bb
[i
].phi
, bb
[i
].psi
, bb
[i
].d3
, bb
[i
].d4
, bb
[i
].d5
,
625 bb
[i
].bHelix
? "Yes" : "No");