minor fixes in ditribution files
[gromacs/qmmm-gamess-us.git] / src / gmxlib / calch.c
blobac2348760bd6bfee4924b1744dbb7c060884f551
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GROningen Mixture of Alchemy and Childrens' Stories
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include "macros.h"
40 #include "calch.h"
41 #include "maths.h"
42 #include "vec.h"
43 #include "physics.h"
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]
54 /* The source code in this file should be thread-safe.
55 Please keep it that way. */
57 static void gen_waterhydrogen(int nh,rvec xa[], rvec xh[],int *l)
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 }
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 }
78 #undef AA
79 #undef BB
80 #undef CC
81 /*static int l=0; removed due to thread-safety issues*/
82 int m;
83 rvec kkk;
85 /* This was copied from Gromos */
86 for(m=0; (m<DIM); m++) {
87 xH1[m]=xAI[m]+matrix1[*l][m];
88 xH2[m]=xAI[m]+matrix2[*l][m];
90 if (nh > 2)
91 copy_rvec(xAI,xH3);
92 if (nh > 3)
93 copy_rvec(xAI,xH4);
95 *l=(*l+1) % 6;
98 void calc_h_pos(int nht, rvec xa[], rvec xh[], int *l)
100 #define alfaH (acos(-1/3.0)) /* 109.47 degrees */
101 #define alfaHpl (2*M_PI/3) /* 120 degrees */
102 #define distH 0.1
104 #define alfaCOM (DEG2RAD*117)
105 #define alfaCO (DEG2RAD*121)
106 #define alfaCOA (DEG2RAD*115)
108 #define distO 0.123
109 #define distOA 0.125
110 #define distOM 0.136
112 rvec sa,sb,sij;
113 real s6,rij,ra,rb,xd;
114 int d;
116 s6=0.5*sqrt(3.e0);
118 /* common work for constructing one, two or three dihedral hydrogens */
119 switch (nht) {
120 case 2:
121 case 3:
122 case 4:
123 case 8:
124 case 9:
125 rij = 0.e0;
126 for(d=0; (d<DIM); d++) {
127 xd = xAJ[d];
128 sij[d] = xAI[d]-xd;
129 sb[d] = xd-xAK[d];
130 rij += sqr(sij[d]);
132 rij = sqrt(rij);
133 sa[XX] = sij[YY]*sb[ZZ]-sij[ZZ]*sb[YY];
134 sa[YY] = sij[ZZ]*sb[XX]-sij[XX]*sb[ZZ];
135 sa[ZZ] = sij[XX]*sb[YY]-sij[YY]*sb[XX];
136 ra = 0.e0;
137 for(d=0; (d<DIM); d++) {
138 sij[d] = sij[d]/rij;
139 ra += sqr(sa[d]);
141 ra = sqrt(ra);
142 for(d=0; (d<DIM); d++)
143 sa[d] = sa[d]/ra;
145 sb[XX] = sa[YY]*sij[ZZ]-sa[ZZ]*sij[YY];
146 sb[YY] = sa[ZZ]*sij[XX]-sa[XX]*sij[ZZ];
147 sb[ZZ] = sa[XX]*sij[YY]-sa[YY]*sij[XX];
148 break;
149 }/* end switch */
151 switch (nht) {
152 case 1: /* construct one planar hydrogen (peptide,rings) */
153 rij = 0.e0;
154 rb = 0.e0;
155 for(d=0; (d<DIM); d++) {
156 sij[d] = xAI[d]-xAJ[d];
157 sb[d] = xAI[d]-xAK[d];
158 rij += sqr(sij[d]);
159 rb += sqr(sb[d]);
161 rij = sqrt(rij);
162 rb = sqrt(rb);
163 ra = 0.e0;
164 for(d=0; (d<DIM); d++) {
165 sa[d] = sij[d]/rij+sb[d]/rb;
166 ra += sqr(sa[d]);
168 ra = sqrt(ra);
169 for(d=0; (d<DIM); d++)
170 xH1[d] = xAI[d]+distH*sa[d]/ra;
171 break;
172 case 2: /* one single hydrogen, e.g. hydroxyl */
173 for(d=0; (d<DIM); d++) {
174 xH1[d] = xAI[d]+distH*sin(alfaH)*sb[d]-distH*cos(alfaH)*sij[d];
176 break;
177 case 3: /* two planar hydrogens, e.g. -NH2 */
178 for(d=0; (d<DIM); d++) {
179 xH1[d] = xAI[d]-distH*sin(alfaHpl)*sb[d]-distH*cos(alfaHpl)*sij[d];
180 xH2[d] = xAI[d]+distH*sin(alfaHpl)*sb[d]-distH*cos(alfaHpl)*sij[d];
182 break;
183 case 4: /* two or three tetrahedral hydrogens, e.g. -CH3 */
184 for(d=0; (d<DIM); d++) {
185 xH1[d] = xAI[d]+distH*sin(alfaH)*sb[d]-distH*cos(alfaH)*sij[d];
186 xH2[d] = ( xAI[d]
187 - distH*sin(alfaH)*0.5*sb[d]
188 + distH*sin(alfaH)*s6*sa[d]
189 - distH*cos(alfaH)*sij[d] );
190 if ( xH3[XX]!=NOTSET && xH3[YY]!=NOTSET && xH3[ZZ]!=NOTSET )
191 xH3[d] = ( xAI[d]
192 - distH*sin(alfaH)*0.5*sb[d]
193 - distH*sin(alfaH)*s6*sa[d]
194 - distH*cos(alfaH)*sij[d] );
196 break;
197 case 5: { /* one tetrahedral hydrogen, e.g. C3CH */
198 real center;
199 rvec dxc;
201 for(d=0; (d<DIM); d++) {
202 center=(xAJ[d]+xAK[d]+xAL[d])/3.0;
203 dxc[d]=xAI[d]-center;
205 center=norm(dxc);
206 for(d=0; (d<DIM); d++)
207 xH1[d]=xAI[d]+dxc[d]*distH/center;
208 break;
210 case 6: { /* two tetrahedral hydrogens, e.g. C-CH2-C */
211 rvec rBB,rCC1,rCC2,rNN;
212 real bb,nn;
214 for(d=0; (d<DIM); d++)
215 rBB[d]=xAI[d]-0.5*(xAJ[d]+xAK[d]);
216 bb=norm(rBB);
218 rvec_sub(xAI,xAJ,rCC1);
219 rvec_sub(xAI,xAK,rCC2);
220 cprod(rCC1,rCC2,rNN);
221 nn=norm(rNN);
223 for(d=0; (d<DIM); d++) {
224 xH1[d]=xAI[d]+distH*(cos(alfaH/2.0)*rBB[d]/bb+
225 sin(alfaH/2.0)*rNN[d]/nn);
226 xH2[d]=xAI[d]+distH*(cos(alfaH/2.0)*rBB[d]/bb-
227 sin(alfaH/2.0)*rNN[d]/nn);
229 break;
231 case 7: /* two water hydrogens */
232 gen_waterhydrogen(2, xa, xh, l);
233 break;
234 case 10: /* three water hydrogens */
235 gen_waterhydrogen(3, xa, xh, l);
236 break;
237 case 11: /* four water hydrogens */
238 gen_waterhydrogen(4, xa, xh, l);
239 break;
240 case 8: /* two carboxyl oxygens, -COO- */
241 for(d=0; (d<DIM); d++) {
242 xH1[d] = xAI[d]-distOM*sin(alfaCOM)*sb[d]-distOM*cos(alfaCOM)*sij[d];
243 xH2[d] = xAI[d]+distOM*sin(alfaCOM)*sb[d]-distOM*cos(alfaCOM)*sij[d];
245 break;
246 case 9: { /* carboxyl oxygens and hydrogen, -COOH */
247 rvec xa2[4]; /* i,j,k,l */
249 /* first add two oxygens */
250 for(d=0; (d<DIM); d++) {
251 xH1[d] = xAI[d]-distO *sin(alfaCO )*sb[d]-distO *cos(alfaCO )*sij[d];
252 xH2[d] = xAI[d]+distOA*sin(alfaCOA)*sb[d]-distOA*cos(alfaCOA)*sij[d];
255 /* now use rule 2 to add hydrogen to 2nd oxygen */
256 copy_rvec(xH2, xa2[0]); /* new i = n' */
257 copy_rvec(xAI, xa2[1]); /* new j = i */
258 copy_rvec(xAJ, xa2[2]); /* new k = j */
259 copy_rvec(xAK, xa2[3]); /* new l = k, not used */
260 calc_h_pos(2, xa2, (xh+2), l);
262 break;
264 default:
265 gmx_fatal(FARGS,"Invalid argument (%d) for nht in routine genh\n",nht);
266 } /* end switch */