changed reading hint
[gromacs/adressmacs.git] / src / tools / hxprops.h
blob56fe22495d66ebc6a251bd2984cad89261b21f99
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 2.0
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
16 * Please refer to:
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
23 * or e-mail to:
24 * gromacs@chem.rug.nl
26 * And Hey:
27 * Great Red Oystrich Makes All Chemists Sane
30 #ifndef _hxprops_h
31 #define _hxprops_h
33 static char *SRCID_hxprops_h = "$Id$";
35 #include <stdio.h>
36 #include "typedefs.h"
38 #define PHI_AHX (-55.0)
39 #define PSI_AHX (-45.0)
40 /* Canonical values of the helix phi/psi angles */
43 typedef struct {
44 real phi,psi,pprms2;
45 real jcaha;
46 real d3,d4,d5,rmsa;
47 bool bHelix;
48 int nhx;
49 int nrms,resno;
50 int Cprev,N,H,CA,C,O,Nnext;
51 char label[32];
52 } t_bb;
54 enum {
55 efhRAD, efhTWIST, efhRISE, efhLEN,
56 efhDIP, efhRMS, efhRMSA, efhCD222,
57 efhPPRMS,efhCPHI, efhPHI, efhPSI,
58 efhHB3, efhHB4, efhHB5, efhJCA,
59 efhAHX, efhNR
62 extern real ahx_len(int gnx,atom_id index[],rvec x[],matrix box);
63 /* Assume we have a list of Calpha atoms only! */
65 extern real ellipticity(int nres,t_bb bb[]);
67 extern real radius(FILE *fp,int nca,atom_id ca_index[],rvec x[]);
68 /* Assume we have calphas */
70 extern real twist(FILE *fp,int nca,atom_id caindex[],rvec x[]);
71 /* Calculate the twist of the helix */
73 extern real pprms(FILE *fp,int nbb,t_bb bb[]);
74 /* Calculate the average RMS from canonical phi/psi values
75 * and the distance per residue
78 extern real ca_phi(int gnx,atom_id index[],rvec x[],matrix box);
79 /* Assume we have a list of Calpha atoms only! */
81 extern real dip(int nbb,atom_id bbind[],rvec x[],t_atom atom[]);
83 extern real rise(int gnx,atom_id index[],rvec x[]);
84 /* Assume we have a list of Calpha atoms only! */
86 extern void av_hblen(FILE *fp3,FILE *fp3a,
87 FILE *fp4,FILE *fp4a,
88 FILE *fp5,FILE *fp5a,
89 real t,int nres,t_bb bb[]);
91 extern void av_phipsi(FILE *fphi,FILE *fpsi,FILE *fphi2,FILE *fpsi2,
92 real t,int nres,t_bb bb[]);
94 extern t_bb *mkbbind(char *fn,int *nres,int *nbb,int res0,
95 int *nall,atom_id **index,
96 char ***atomname,t_atom atom[],
97 char ***resname);
99 extern void do_start_end(int nres,t_bb bb[],rvec x[],int *nbb,
100 atom_id bbindex[],int *nca,atom_id caindex[],
101 bool bRange,int rStart,int rEnd);
103 extern void calc_hxprops(int nres,t_bb bb[],rvec x[],matrix box);
105 extern void pr_bb(FILE *fp,int nres,t_bb bb[]);
107 #endif