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