Code beautification with uncrustify
[alexxy/gromacs.git] / src / gromacs / gmxlib / calch.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "macros.h"
40 #include "calch.h"
41 #include "maths.h"
42 #include "vec.h"
43 #include "physics.h"
44
45 #define xAI xa[0]
46 #define xAJ xa[1]
47 #define xAK xa[2]
48 #define xAL xa[3]
49 #define xH1 xh[0]
50 #define xH2 xh[1]
51 #define xH3 xh[2]
52 #define xH4 xh[3]
53
54 /* The source code in this file should be thread-safe.
55       Please keep it that way. */
56
57 static void gen_waterhydrogen(int nh, rvec xa[], rvec xh[], int *l)
58 {
59 #define AA 0.081649
60 #define BB 0.0
61 #define CC 0.0577350
62     const  rvec   matrix1[6] = {
63         { AA,     BB,     CC },
64         { AA,     BB,     CC },
65         { AA,     BB,     CC },
66         { -AA,    BB,     CC },
67         { -AA,    BB,     CC },
68         { BB,     AA,    -CC }
69     };
70     const  rvec   matrix2[6] = {
71         { -AA,   BB,   CC },
72         { BB,    AA,  -CC },
73         { BB,   -AA,  -CC },
74         { BB,    AA,  -CC },
75         { BB,   -AA,  -CC },
76         { BB,   -AA,  -CC }
77     };
78 #undef AA
79 #undef BB
80 #undef CC
81     int        m;
82     rvec       kkk;
83
84     /* This was copied from Gromos */
85     for (m = 0; (m < DIM); m++)
86     {
87         xH1[m] = xAI[m]+matrix1[*l][m];
88         xH2[m] = xAI[m]+matrix2[*l][m];
89     }
90     if (nh > 2)
91     {
92         copy_rvec(xAI, xH3);
93     }
94     if (nh > 3)
95     {
96         copy_rvec(xAI, xH4);
97     }
98
99     *l = (*l+1) % 6;
100 }
101
102 void calc_h_pos(int nht, rvec xa[], rvec xh[], int *l)
103 {
104 #define alfaH   (acos(-1/3.0)) /* 109.47 degrees */
105 #define alfaHpl (2*M_PI/3)     /* 120 degrees */
106 #define distH   0.1
107
108 #define alfaCOM (DEG2RAD*117)
109 #define alfaCO  (DEG2RAD*121)
110 #define alfaCOA (DEG2RAD*115)
111
112 #define distO   0.123
113 #define distOA  0.125
114 #define distOM  0.136
115
116     rvec sa, sb, sij;
117     real s6, rij, ra, rb, xd;
118     int  d;
119
120     s6 = 0.5*sqrt(3.e0);
121
122     /* common work for constructing one, two or three dihedral hydrogens */
123     switch (nht)
124     {
125         case 2:
126         case 3:
127         case 4:
128         case 8:
129         case 9:
130             rij = 0.e0;
131             for (d = 0; (d < DIM); d++)
132             {
133                 xd     = xAJ[d];
134                 sij[d] = xAI[d]-xd;
135                 sb[d]  = xd-xAK[d];
136                 rij   += sqr(sij[d]);
137             }
138             rij    = sqrt(rij);
139             sa[XX] = sij[YY]*sb[ZZ]-sij[ZZ]*sb[YY];
140             sa[YY] = sij[ZZ]*sb[XX]-sij[XX]*sb[ZZ];
141             sa[ZZ] = sij[XX]*sb[YY]-sij[YY]*sb[XX];
142             ra     = 0.e0;
143             for (d = 0; (d < DIM); d++)
144             {
145                 sij[d] = sij[d]/rij;
146                 ra    += sqr(sa[d]);
147             }
148             ra = sqrt(ra);
149             for (d = 0; (d < DIM); d++)
150             {
151                 sa[d] = sa[d]/ra;
152             }
153
154             sb[XX] = sa[YY]*sij[ZZ]-sa[ZZ]*sij[YY];
155             sb[YY] = sa[ZZ]*sij[XX]-sa[XX]*sij[ZZ];
156             sb[ZZ] = sa[XX]*sij[YY]-sa[YY]*sij[XX];
157             break;
158     } /* end switch */
159
160     switch (nht)
161     {
162         case 1: /* construct one planar hydrogen (peptide,rings) */
163             rij = 0.e0;
164             rb  = 0.e0;
165             for (d = 0; (d < DIM); d++)
166             {
167                 sij[d] = xAI[d]-xAJ[d];
168                 sb[d]  = xAI[d]-xAK[d];
169                 rij   += sqr(sij[d]);
170                 rb    += sqr(sb[d]);
171             }
172             rij = sqrt(rij);
173             rb  = sqrt(rb);
174             ra  = 0.e0;
175             for (d = 0; (d < DIM); d++)
176             {
177                 sa[d] = sij[d]/rij+sb[d]/rb;
178                 ra   += sqr(sa[d]);
179             }
180             ra = sqrt(ra);
181             for (d = 0; (d < DIM); d++)
182             {
183                 xH1[d] = xAI[d]+distH*sa[d]/ra;
184             }
185             break;
186         case 2: /* one single hydrogen, e.g. hydroxyl */
187             for (d = 0; (d < DIM); d++)
188             {
189                 xH1[d] = xAI[d]+distH*sin(alfaH)*sb[d]-distH*cos(alfaH)*sij[d];
190             }
191             break;
192         case 3: /* two planar hydrogens, e.g. -NH2 */
193             for (d = 0; (d < DIM); d++)
194             {
195                 xH1[d] = xAI[d]-distH*sin(alfaHpl)*sb[d]-distH*cos(alfaHpl)*sij[d];
196                 xH2[d] = xAI[d]+distH*sin(alfaHpl)*sb[d]-distH*cos(alfaHpl)*sij[d];
197             }
198             break;
199         case 4: /* two or three tetrahedral hydrogens, e.g. -CH3 */
200             for (d = 0; (d < DIM); d++)
201             {
202                 xH1[d] = xAI[d]+distH*sin(alfaH)*sb[d]-distH*cos(alfaH)*sij[d];
203                 xH2[d] = ( xAI[d]
204                            - distH*sin(alfaH)*0.5*sb[d]
205                            + distH*sin(alfaH)*s6*sa[d]
206                            - distH*cos(alfaH)*sij[d] );
207                 if (xH3[XX] != NOTSET && xH3[YY] != NOTSET && xH3[ZZ] != NOTSET)
208                 {
209                     xH3[d] = ( xAI[d]
210                                - distH*sin(alfaH)*0.5*sb[d]
211                                - distH*sin(alfaH)*s6*sa[d]
212                                - distH*cos(alfaH)*sij[d] );
213                 }
214             }
215             break;
216         case 5: /* one tetrahedral hydrogen, e.g. C3CH */
217         {
218             real center;
219             rvec dxc;
220
221             for (d = 0; (d < DIM); d++)
222             {
223                 center = (xAJ[d]+xAK[d]+xAL[d])/3.0;
224                 dxc[d] = xAI[d]-center;
225             }
226             center = norm(dxc);
227             for (d = 0; (d < DIM); d++)
228             {
229                 xH1[d] = xAI[d]+dxc[d]*distH/center;
230             }
231             break;
232         }
233         case 6: /* two tetrahedral hydrogens, e.g. C-CH2-C */
234         {
235             rvec rBB, rCC1, rCC2, rNN;
236             real bb, nn;
237
238             for (d = 0; (d < DIM); d++)
239             {
240                 rBB[d] = xAI[d]-0.5*(xAJ[d]+xAK[d]);
241             }
242             bb = norm(rBB);
243
244             rvec_sub(xAI, xAJ, rCC1);
245             rvec_sub(xAI, xAK, rCC2);
246             cprod(rCC1, rCC2, rNN);
247             nn = norm(rNN);
248
249             for (d = 0; (d < DIM); d++)
250             {
251                 xH1[d] = xAI[d]+distH*(cos(alfaH/2.0)*rBB[d]/bb+
252                                        sin(alfaH/2.0)*rNN[d]/nn);
253                 xH2[d] = xAI[d]+distH*(cos(alfaH/2.0)*rBB[d]/bb-
254                                        sin(alfaH/2.0)*rNN[d]/nn);
255             }
256             break;
257         }
258         case 7: /* two water hydrogens */
259             gen_waterhydrogen(2, xa, xh, l);
260             break;
261         case 10: /* three water hydrogens */
262             gen_waterhydrogen(3, xa, xh, l);
263             break;
264         case 11: /* four water hydrogens */
265             gen_waterhydrogen(4, xa, xh, l);
266             break;
267         case 8: /* two carboxyl oxygens, -COO- */
268             for (d = 0; (d < DIM); d++)
269             {
270                 xH1[d] = xAI[d]-distOM*sin(alfaCOM)*sb[d]-distOM*cos(alfaCOM)*sij[d];
271                 xH2[d] = xAI[d]+distOM*sin(alfaCOM)*sb[d]-distOM*cos(alfaCOM)*sij[d];
272             }
273             break;
274         case 9:          /* carboxyl oxygens and hydrogen, -COOH */
275         {
276             rvec xa2[4]; /* i,j,k,l   */
277
278             /* first add two oxygens */
279             for (d = 0; (d < DIM); d++)
280             {
281                 xH1[d] = xAI[d]-distO *sin(alfaCO )*sb[d]-distO *cos(alfaCO )*sij[d];
282                 xH2[d] = xAI[d]+distOA*sin(alfaCOA)*sb[d]-distOA*cos(alfaCOA)*sij[d];
283             }
284
285             /* now use rule 2 to add hydrogen to 2nd oxygen */
286             copy_rvec(xH2, xa2[0]); /* new i = n' */
287             copy_rvec(xAI, xa2[1]); /* new j = i  */
288             copy_rvec(xAJ, xa2[2]); /* new k = j  */
289             copy_rvec(xAK, xa2[3]); /* new l = k, not used */
290             calc_h_pos(2, xa2, (xh+2), l);
291
292             break;
293         }
294         default:
295             gmx_fatal(FARGS, "Invalid argument (%d) for nht in routine genh\n", nht);
296     } /* end switch */
297 }