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