Redefine the default boolean type to gmx_bool.
[alexxy/gromacs.git] / src / mdlib / calcvir.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  * GROwing Monsters And Cloning Shrimps
34  */
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include "sysstuff.h"
41 #include "force.h"
42 #include "vec.h"
43 #include "mshift.h"
44 #include "macros.h"
45         
46 static void dprod1(rvec vir,real x,rvec f)
47 {
48   vir[XX]+=x*f[XX];
49   vir[YY]+=x*f[YY];
50   vir[ZZ]+=x*f[ZZ];
51 }
52
53 static void upd_vir(rvec vir,real dvx,real dvy,real dvz)
54 {
55   vir[XX]-=0.5*dvx;
56   vir[YY]-=0.5*dvy;
57   vir[ZZ]-=0.5*dvz;
58 }
59
60 void calc_vir(FILE *log,int nxf,rvec x[],rvec f[],tensor vir,
61               gmx_bool bScrewPBC,matrix box)
62 {
63   int      i,isx;
64   double   dvxx=0,dvxy=0,dvxz=0,dvyx=0,dvyy=0,dvyz=0,dvzx=0,dvzy=0,dvzz=0;
65     
66   for(i=0; (i<nxf); i++) {
67     dvxx+=x[i][XX]*f[i][XX];
68     dvxy+=x[i][XX]*f[i][YY];
69     dvxz+=x[i][XX]*f[i][ZZ];
70     
71     dvyx+=x[i][YY]*f[i][XX];
72     dvyy+=x[i][YY]*f[i][YY];
73     dvyz+=x[i][YY]*f[i][ZZ];
74     
75     dvzx+=x[i][ZZ]*f[i][XX];
76     dvzy+=x[i][ZZ]*f[i][YY];
77     dvzz+=x[i][ZZ]*f[i][ZZ];
78     
79     if (bScrewPBC) {
80       isx = IS2X(i);
81       /* We should correct all odd x-shifts, but the range of isx is -2 to 2 */
82       if (isx == 1 || isx == -1) {
83         dvyy += box[YY][YY]*f[i][YY];
84         dvyz += box[YY][YY]*f[i][ZZ];
85
86         dvzy += box[ZZ][ZZ]*f[i][YY];
87         dvzz += box[ZZ][ZZ]*f[i][ZZ];
88       }
89     }
90   }
91   
92   upd_vir(vir[XX],dvxx,dvxy,dvxz);
93   upd_vir(vir[YY],dvyx,dvyy,dvyz);
94   upd_vir(vir[ZZ],dvzx,dvzy,dvzz);
95 }
96
97
98 static void lo_fcv(int i0,int i1,int g0,
99                    real x[],real f[],tensor vir,
100                    int is[],real box[], gmx_bool bTriclinic)
101 {
102   int      i,i3,gg,g3,tx,ty,tz;
103   real     xx,yy,zz;
104   real     dvxx=0,dvxy=0,dvxz=0,dvyx=0,dvyy=0,dvyz=0,dvzx=0,dvzy=0,dvzz=0;
105
106   if(bTriclinic) {
107       for(i=i0,gg=g0; (i<i1); i++,gg++) {
108           i3=DIM*i;
109           g3=DIM*gg;
110           tx=is[g3+XX];
111           ty=is[g3+YY];
112           tz=is[g3+ZZ];
113           
114           xx=x[i3+XX]-tx*box[XXXX]-ty*box[YYXX]-tz*box[ZZXX];
115           dvxx+=xx*f[i3+XX];
116           dvxy+=xx*f[i3+YY];
117           dvxz+=xx*f[i3+ZZ];
118           
119           yy=x[i3+YY]-ty*box[YYYY]-tz*box[ZZYY];
120           dvyx+=yy*f[i3+XX];
121           dvyy+=yy*f[i3+YY];
122           dvyz+=yy*f[i3+ZZ];
123           
124           zz=x[i3+ZZ]-tz*box[ZZZZ]; 
125           dvzx+=zz*f[i3+XX];
126           dvzy+=zz*f[i3+YY];
127           dvzz+=zz*f[i3+ZZ];
128       }
129   } else {
130       for(i=i0,gg=g0; (i<i1); i++,gg++) {
131           i3=DIM*i;
132           g3=DIM*gg;
133           tx=is[g3+XX];
134           ty=is[g3+YY];
135           tz=is[g3+ZZ];
136           
137           xx=x[i3+XX]-tx*box[XXXX];
138           dvxx+=xx*f[i3+XX];
139           dvxy+=xx*f[i3+YY];
140           dvxz+=xx*f[i3+ZZ];
141           
142           yy=x[i3+YY]-ty*box[YYYY];
143           dvyx+=yy*f[i3+XX];
144           dvyy+=yy*f[i3+YY];
145           dvyz+=yy*f[i3+ZZ];
146           
147           zz=x[i3+ZZ]-tz*box[ZZZZ]; 
148           dvzx+=zz*f[i3+XX];
149           dvzy+=zz*f[i3+YY];
150           dvzz+=zz*f[i3+ZZ];
151       }
152   }
153   
154   upd_vir(vir[XX],dvxx,dvxy,dvxz);
155   upd_vir(vir[YY],dvyx,dvyy,dvyz);
156   upd_vir(vir[ZZ],dvzx,dvzy,dvzz);
157 }
158
159 static void lo_fcv2(int i0,int i1,
160                     rvec x[],rvec f[],tensor vir,
161                     ivec is[],matrix box, gmx_bool bTriclinic)
162 {
163   int      i,gg,tx,ty,tz;
164   real     xx,yy,zz;
165   real     dvxx=0,dvxy=0,dvxz=0,dvyx=0,dvyy=0,dvyz=0,dvzx=0,dvzy=0,dvzz=0;
166
167   if(bTriclinic) {
168       for(i=i0,gg=0; (i<i1); i++,gg++) {
169           tx=is[gg][XX];
170           ty=is[gg][YY];
171           tz=is[gg][ZZ];
172           
173           xx=x[i][XX]-tx*box[XX][XX]-ty*box[YY][XX]-tz*box[ZZ][XX];
174           dvxx+=xx*f[i][XX];
175           dvxy+=xx*f[i][YY];
176           dvxz+=xx*f[i][ZZ];
177           
178           yy=x[i][YY]-ty*box[YY][YY]-tz*box[ZZ][YY];
179           dvyx+=yy*f[i][XX];
180           dvyy+=yy*f[i][YY];
181           dvyz+=yy*f[i][ZZ];
182           
183           zz=x[i][ZZ]-tz*box[ZZ][ZZ];
184           dvzx+=zz*f[i][XX];
185           dvzy+=zz*f[i][YY];
186           dvzz+=zz*f[i][ZZ];
187       }
188   } else {
189       for(i=i0,gg=0; (i<i1); i++,gg++) {
190           tx=is[gg][XX];
191           ty=is[gg][YY];
192           tz=is[gg][ZZ];
193           
194           xx=x[i][XX]-tx*box[XX][XX];
195           dvxx+=xx*f[i][XX];
196           dvxy+=xx*f[i][YY];
197           dvxz+=xx*f[i][ZZ];
198           
199           yy=x[i][YY]-ty*box[YY][YY];
200           dvyx+=yy*f[i][XX];
201           dvyy+=yy*f[i][YY];
202           dvyz+=yy*f[i][ZZ];
203           
204           zz=x[i][ZZ]-tz*box[ZZ][ZZ];
205           dvzx+=zz*f[i][XX];
206           dvzy+=zz*f[i][YY];
207           dvzz+=zz*f[i][ZZ];
208       }
209   }
210   
211   upd_vir(vir[XX],dvxx,dvxy,dvxz);
212   upd_vir(vir[YY],dvyx,dvyy,dvyz);
213   upd_vir(vir[ZZ],dvzx,dvzy,dvzz);
214 }
215
216 void f_calc_vir(FILE *log,int i0,int i1,rvec x[],rvec f[],tensor vir,
217                 t_graph *g,matrix box)
218 {
219   int start,end;
220   
221   if (g && (g->nnodes > 0)) {
222     /* Calculate virial for bonded forces only when they belong to
223      * this node.
224      */
225     start = max(i0,g->start);
226     end   = min(i1,g->end+1);
227 #ifdef SAFE
228     lo_fcv2(start,end,x,f,vir,g->ishift,box,TRICLINIC(box));
229 #else
230     lo_fcv(start,end,0,x[0],f[0],vir,g->ishift[0],box[0],TRICLINIC(box));
231 #endif
232     
233     /* If not all atoms are bonded, calculate their virial contribution 
234      * anyway, without shifting back their coordinates.
235      * Note the nifty pointer arithmetic...
236      */
237     if (start > i0) 
238       calc_vir(log,start-i0,x + i0,f + i0,vir,FALSE,box);
239     if (end < i1)
240       calc_vir(log,i1-end,x + end,f + end,vir,FALSE,box);
241   }
242   else
243     calc_vir(log,i1-i0,x + i0,f + i0,vir,FALSE,box);
244 }