1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This file is part of Gromacs Copyright (c) 1991-2008
5 * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
7 * This program is free software; you can redistribute it and/or
8 * modify it under the terms of the GNU General Public License
9 * as published by the Free Software Foundation; either version 2
10 * of the License, or (at your option) any later version.
12 * To help us fund GROMACS development, we humbly ask that you cite
13 * the research papers on the package. Check out http://www.gromacs.org
16 * Gnomes, ROck Monsters And Chili Sauce
27 #include "domdec_network.h"
31 static void calc_cgcm_av_stddev(t_block *cgs,int n,rvec *x,rvec av,rvec stddev,
37 int cg,d,k0,k1,k,nrcg;
52 copy_rvec(x[k0],cg_cm);
59 for(k=k0; (k<k1); k++)
63 for(d=0; (d<DIM); d++)
71 s2[d] += cg_cm[d]*cg_cm[d];
83 gmx_sumd(7,buf,cr_sum);
89 n = (int)(buf[6] + 0.5);
98 stddev[d] = sqrt(s2[d] - s1[d]*s1[d]);
102 static void set_tric_dir(ivec *dd_nc,gmx_ddbox_t *ddbox,matrix box)
106 real dep,inv_skew_fac2;
108 npbcdim = ddbox->npbcdim;
109 normal = ddbox->normal;
112 ddbox->tric_dir[d] = 0;
113 for(j=d+1; j<npbcdim; j++)
117 ddbox->tric_dir[d] = 1;
118 if (dd_nc != NULL && (*dd_nc)[j] > 1 && (*dd_nc)[d] == 1)
120 gmx_fatal(FARGS,"Domain decomposition has not been implemented for box vectors that have non-zero components in directions that do not use domain decomposition: ncells = %d %d %d, box vector[%d] = %f %f %f",
121 dd_nc[XX],dd_nc[YY],dd_nc[ZZ],
122 j+1,box[j][XX],box[j][YY],box[j][ZZ]);
127 /* Convert box vectors to orthogonal vectors for this dimension,
128 * for use in distance calculations.
129 * Set the trilinic skewing factor that translates
130 * the thickness of a slab perpendicular to this dimension
131 * into the real thickness of the slab.
133 if (ddbox->tric_dir[d])
137 if (d == XX || d == YY)
139 /* Normalize such that the "diagonal" is 1 */
140 svmul(1/box[d+1][d+1],box[d+1],v[d+1]);
145 inv_skew_fac2 += sqr(v[d+1][d]);
148 /* Normalize such that the "diagonal" is 1 */
149 svmul(1/box[d+2][d+2],box[d+2],v[d+2]);
154 /* Make vector [d+2] perpendicular to vector [d+1],
155 * this does not affect the normalization.
157 dep = iprod(v[d+1],v[d+2])/norm2(v[d+1]);
160 v[d+2][i] -= dep*v[d+1][i];
162 inv_skew_fac2 += sqr(v[d+2][d]);
164 cprod(v[d+1],v[d+2],normal[d]);
168 /* cross product with (1,0,0) */
170 normal[d][YY] = v[d+1][ZZ];
171 normal[d][ZZ] = -v[d+1][YY];
175 fprintf(debug,"box[%d] %.3f %.3f %.3f\n",
176 d,box[d][XX],box[d][YY],box[d][ZZ]);
177 for(i=d+1; i<DIM; i++)
179 fprintf(debug," v[%d] %.3f %.3f %.3f\n",
180 i,v[i][XX],v[i][YY],v[i][ZZ]);
184 ddbox->skew_fac[d] = 1.0/sqrt(inv_skew_fac2);
185 /* Set the normal vector length to skew_fac */
186 dep = ddbox->skew_fac[d]/norm(normal[d]);
187 svmul(dep,normal[d],normal[d]);
191 fprintf(debug,"skew_fac[%d] = %f\n",d,ddbox->skew_fac[d]);
192 fprintf(debug,"normal[%d] %.3f %.3f %.3f\n",
193 d,normal[d][XX],normal[d][YY],normal[d][ZZ]);
198 ddbox->skew_fac[d] = 1;
202 clear_rvec(ddbox->v[d][i]);
203 ddbox->v[d][i][i] = 1;
205 clear_rvec(normal[d]);
211 static void low_set_ddbox(t_inputrec *ir,ivec *dd_nc,matrix box,
212 gmx_bool bCalcUnboundedSize,int ncg,t_block *cgs,rvec *x,
220 ddbox->npbcdim = ePBC2npbcdim(ir->ePBC);
221 ddbox->nboundeddim = inputrec2nboundeddim(ir);
223 for(d=0; d<ddbox->nboundeddim; d++)
226 ddbox->box_size[d] = box[d][d];
229 if (ddbox->nboundeddim < DIM && bCalcUnboundedSize)
231 calc_cgcm_av_stddev(cgs,ncg,x,av,stddev,cr_sum);
233 /* GRID_STDDEV_FAC * stddev
234 * gives a uniform load for a rectangular block of cg's.
235 * For a sphere it is not a bad approximation for 4x1x1 up to 4x2x2.
237 for(d=ddbox->nboundeddim; d<DIM; d++)
239 b0 = av[d] - GRID_STDDEV_FAC*stddev[d];
240 b1 = av[d] + GRID_STDDEV_FAC*stddev[d];
243 fprintf(debug,"Setting global DD grid boundaries to %f - %f\n",
247 ddbox->box_size[d] = b1 - b0;
251 set_tric_dir(dd_nc,ddbox,box);
254 void set_ddbox(gmx_domdec_t *dd,gmx_bool bMasterState,t_commrec *cr_sum,
255 t_inputrec *ir,matrix box,
256 gmx_bool bCalcUnboundedSize,t_block *cgs,rvec *x,
259 if (!bMasterState || DDMASTER(dd))
261 low_set_ddbox(ir,&dd->nc,box,bCalcUnboundedSize,
262 bMasterState ? cgs->nr : dd->ncg_home,cgs,x,
263 bMasterState ? NULL : cr_sum,
269 dd_bcast(dd,sizeof(gmx_ddbox_t),ddbox);
273 void set_ddbox_cr(t_commrec *cr,ivec *dd_nc,
274 t_inputrec *ir,matrix box,t_block *cgs,rvec *x,
279 low_set_ddbox(ir,dd_nc,box,TRUE,cgs->nr,cgs,x,NULL,ddbox);
282 gmx_bcast(sizeof(gmx_ddbox_t),ddbox,cr);