Actually fix gmx helix segmentation fault
[gromacs.git] / src / gromacs / gmxana / hxprops.h
blob25f2e492691c847e69e17de42c117a09abdbdf26
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,2015,2016,2018, 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.
38 #ifndef _hxprops_h
39 #define _hxprops_h
41 #include <stdio.h>
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/topology/idef.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/real.h"
48 struct t_atom;
49 struct t_resinfo;
51 #ifdef __cplusplus
52 extern "C"
54 #endif
56 #define PHI_AHX (-55.0)
57 #define PSI_AHX (-45.0)
58 /* Canonical values of the helix phi/psi angles */
61 /*! \internal \brief Struct containing properties of a residue in a protein backbone. */
62 struct t_bb {
63 //! Protein backbone phi angle.
64 real phi;
65 //! Protein backbone psi angle.
66 real psi;
67 //! RMS distance of phi and psi angles from ideal helix
68 real pprms2;
69 //! Estimated J-coupling value
70 real jcaha;
71 //! Value of 3 turn helix?
72 real d3;
73 //! Value of 4 turn helix?
74 real d4;
75 //! Value of 5 turn?
76 real d5;
77 //! Average of RMS for analysis.
78 real rmsa;
79 //! If the structure is helical.
80 gmx_bool bHelix;
81 //! Number of elliptical elements
82 int nhx;
83 //! Average RMS Deviation when atoms of this residue are fitted to ideal helix
84 int nrms;
85 //! Residue index for output, relative to gmx_helix -r0 value
86 int resno;
87 //! Index for previous carbon.
88 int Cprev;
89 //! Index for backbone nitrogen.
90 int N;
91 //! Index for backbone NH hydrogen.
92 int H;
93 //! Index for alpha carbon.
94 int CA;
95 //! Index for carbonyl carbon.
96 int C;
97 //! Index for carbonyl oxygen.
98 int O;
99 //! Index for next backbone nitrogen.
100 int Nnext;
101 //! Name for this residue.
102 char label[32];
105 enum {
106 efhRAD, efhTWIST, efhRISE, efhLEN,
107 efhDIP, efhRMS, efhRMSA, efhCD222,
108 efhPPRMS, efhCPHI, efhPHI, efhPSI,
109 efhHB3, efhHB4, efhHB5, efhJCA,
110 efhAHX, efhNR
113 extern real ahx_len(int gnx, int index[], rvec x[]);
114 /* Assume we have a list of Calpha atoms only! */
116 extern real ellipticity(int nres, t_bb bb[]);
118 extern real radius(FILE *fp, int nca, int ca_index[], rvec x[]);
119 /* Assume we have calphas */
121 extern real twist(int nca, int caindex[], rvec x[]);
122 /* Calculate the twist of the helix */
124 extern real pprms(FILE *fp, int nbb, t_bb bb[]);
125 /* Calculate the average RMS from canonical phi/psi values
126 * and the distance per residue
129 extern real ca_phi(int gnx, int index[], rvec x[]);
130 /* Assume we have a list of Calpha atoms only! */
132 extern real dip(int nbb, const int bbind[], const rvec x[], const t_atom atom[]);
134 extern real rise(int gnx, int index[], rvec x[]);
135 /* Assume we have a list of Calpha atoms only! */
137 extern void av_hblen(FILE *fp3, FILE *fp3a,
138 FILE *fp4, FILE *fp4a,
139 FILE *fp5, FILE *fp5a,
140 real t, int nres, t_bb bb[]);
142 extern void av_phipsi(FILE *fphi, FILE *fpsi, FILE *fphi2, FILE *fpsi2,
143 real t, int nres, t_bb bb[]);
145 /*! \brief Allocate and fill an array of information about residues in a protein backbone.
147 * The user is propted for an index group of protein residues (little
148 * error checking occurs). For the number of residues found in the
149 * selected group, nbb entries are made in the returned array. Each
150 * entry contains the atom indices of the N, H, CA, C and O atoms (for
151 * PRO, H means CD), as well as the C of the previous residue and the
152 * N of the next (-1 if not found).
154 * In the output array, the first residue will be numbered starting
155 * from res0. */
156 extern t_bb *mkbbind(const char *fn, int *nres, int *nbb, int res0,
157 int *nall, int **index,
158 char ***atomname, t_atom atom[],
159 t_resinfo *resinfo);
161 extern void do_start_end(int nres, t_bb bb[], int *nbb,
162 int bbindex[], int *nca, int caindex[],
163 gmx_bool bRange, int rStart, int rEnd);
165 extern void calc_hxprops(int nres, t_bb bb[], const rvec x[]);
167 extern void pr_bb(FILE *fp, int nres, t_bb bb[]);
169 #ifdef __cplusplus
171 #endif
173 #endif