Voro++
v_compute.cc
Go to the documentation of this file.
00001 // Voro++, a 3D cell-based Voronoi library
00002 //
00003 // Author   : Chris H. Rycroft (LBL / UC Berkeley)
00004 // Email    : chr@alum.mit.edu
00005 // Date     : August 30th 2011
00006 
00007 /** \file v_compute.cc
00008  * \brief Function implementantions for the voro_compute template. */
00009 
00010 #include "worklist.hh"
00011 #include "v_compute.hh"
00012 #include "rad_option.hh"
00013 #include "container.hh"
00014 #include "container_prd.hh"
00015 
00016 namespace voro {
00017 
00018 /** The class constructor initializes constants from the container class, and
00019  * sets up the mask and queue used for Voronoi computations.
00020  * \param[in] con_ a reference to the container class to use.
00021  * \param[in] (hx_,hy_,hz_) the size of the mask to use. */
00022 template<class c_class>
00023 voro_compute<c_class>::voro_compute(c_class &con_,int hx_,int hy_,int hz_) :
00024         con(con_), boxx(con_.boxx), boxy(con_.boxy), boxz(con_.boxz),
00025         xsp(con_.xsp), ysp(con_.ysp), zsp(con_.zsp),
00026         hx(hx_), hy(hy_), hz(hz_), hxy(hx_*hy_), hxyz(hxy*hz_), ps(con_.ps),
00027         id(con_.id), p(con_.p), co(con_.co), bxsq(boxx*boxx+boxy*boxy+boxz*boxz),
00028         mv(0), qu_size(3*(3+hxy+hz*(hx+hy))), wl(con_.wl), mrad(con_.mrad),
00029         mask(new unsigned int[hxyz]), qu(new int[qu_size]), qu_l(qu+qu_size) {
00030         reset_mask();
00031 }
00032 
00033 /** Scans all of the particles within a block to see if any of them have a
00034  * smaller distance to the given test vector. If one is found, the routine
00035  * updates the minimum distance and store information about this particle.
00036  * \param[in] ijk the index of the block.
00037  * \param[in] (x,y,z) the test vector to consider (which may have already had a
00038  *                    periodic displacement applied to it).
00039  * \param[in] (di,dj,dk) the coordinates of the current block, to store if the
00040  *                       particle record is updated.
00041  * \param[in,out] w a reference to a particle record in which to store
00042  *                  information about the particle whose Voronoi cell the
00043  *                  vector is within.
00044  * \param[in,out] mrs the current minimum distance, that may be updated if a
00045  *                    closer particle is found. */
00046 template<class c_class>
00047 inline void voro_compute<c_class>::scan_all(int ijk,double x,double y,double z,int di,int dj,int dk,particle_record &w,double &mrs) {
00048         double x1,y1,z1,rs;bool in_block=false;
00049         for(int l=0;l<co[ijk];l++) {
00050                 x1=p[ijk][ps*l]-x;
00051                 y1=p[ijk][ps*l+1]-y;
00052                 z1=p[ijk][ps*l+2]-z;
00053                 rs=con.r_current_sub(x1*x1+y1*y1+z1*z1,ijk,l);
00054                 if(rs<mrs) {mrs=rs;w.l=l;in_block=true;}
00055         }
00056         if(in_block) {w.ijk=ijk;w.di=di;w.dj=dj,w.dk=dk;}
00057 }
00058 
00059 /** Finds the Voronoi cell that given vector is within. For containers that are
00060  * not radially dependent, this corresponds to findig the particle that is
00061  * closest to the vector; for the radical tessellation containers, this
00062  * corresponds to a finding the minimum weighted distance.
00063  * \param[in] (x,y,z) the vector to consider.
00064  * \param[in] (ci,cj,ck) the coordinates of the block that the test particle is
00065  *                       in relative to the container data structure.
00066  * \param[in] ijk the index of the block that the test particle is in.
00067  * \param[out] w a reference to a particle record in which to store information
00068  *               about the particle whose Voronoi cell the vector is within.
00069  * \param[out] mrs the minimum computed distance. */
00070 template<class c_class>
00071 void voro_compute<c_class>::find_voronoi_cell(double x,double y,double z,int ci,int cj,int ck,int ijk,particle_record &w,double &mrs) {
00072         double qx=0,qy=0,qz=0,rs;
00073         int i,j,k,di,dj,dk,ei,ej,ek,f,g,disp;
00074         double fx,fy,fz,mxs,mys,mzs,*radp;
00075         unsigned int q,*e,*mijk;
00076 
00077         // Init setup for parameters to return
00078         w.ijk=-1;mrs=large_number;
00079 
00080         con.initialize_search(ci,cj,ck,ijk,i,j,k,disp);
00081 
00082         // Test all particles in the particle's local region first
00083         scan_all(ijk,x,y,z,0,0,0,w,mrs);
00084 
00085         // Now compute the fractional position of the particle within its
00086         // region and store it in (fx,fy,fz). We use this to compute an index
00087         // (di,dj,dk) of which subregion the particle is within.
00088         unsigned int m1,m2;
00089         con.frac_pos(x,y,z,ci,cj,ck,fx,fy,fz);
00090         di=int(fx*xsp*wl_fgrid);dj=int(fy*ysp*wl_fgrid);dk=int(fz*zsp*wl_fgrid);
00091 
00092         // The indices (di,dj,dk) tell us which worklist to use, to test the
00093         // blocks in the optimal order. But we only store worklists for the
00094         // eighth of the region where di, dj, and dk are all less than half the
00095         // full grid. The rest of the cases are handled by symmetry. In this
00096         // section, we detect for these cases, by reflecting high values of di,
00097         // dj, and dk. For these cases, a mask is constructed in m1 and m2
00098         // which is used to flip the worklist information when it is loaded.
00099         if(di>=wl_hgrid) {
00100                 mxs=boxx-fx;
00101                 m1=127+(3<<21);m2=1+(1<<21);di=wl_fgrid-1-di;if(di<0) di=0;
00102         } else {m1=m2=0;mxs=fx;}
00103         if(dj>=wl_hgrid) {
00104                 mys=boxy-fy;
00105                 m1|=(127<<7)+(3<<24);m2|=(1<<7)+(1<<24);dj=wl_fgrid-1-dj;if(dj<0) dj=0;
00106         } else mys=fy;
00107         if(dk>=wl_hgrid) {
00108                 mzs=boxz-fz;
00109                 m1|=(127<<14)+(3<<27);m2|=(1<<14)+(1<<27);dk=wl_fgrid-1-dk;if(dk<0) dk=0;
00110         } else mzs=fz;
00111 
00112         // Do a quick test to account for the case when the minimum radius is
00113         // small enought that no other blocks need to be considered
00114         rs=con.r_max_add(mrs);
00115         if(mxs*mxs>rs&&mys*mys>rs&&mzs*mzs>rs) return;
00116 
00117         // Now compute which worklist we are going to use, and set radp and e to
00118         // point at the right offsets
00119         ijk=di+wl_hgrid*(dj+wl_hgrid*dk);
00120         radp=mrad+ijk*wl_seq_length;
00121         e=(const_cast<unsigned int*> (wl))+ijk*wl_seq_length;
00122 
00123         // Read in how many items in the worklist can be tested without having to
00124         // worry about writing to the mask
00125         f=e[0];g=0;
00126         do {
00127 
00128                 // If mrs is less than the minimum distance to any untested
00129                 // block, then we are done
00130                 if(con.r_max_add(mrs)<radp[g]) return;
00131                 g++;
00132 
00133                 // Load in a block off the worklist, permute it with the
00134                 // symmetry mask, and decode its position. These are all
00135                 // integer bit operations so they should run very fast.
00136                 q=e[g];q^=m1;q+=m2;
00137                 di=q&127;di-=64;
00138                 dj=(q>>7)&127;dj-=64;
00139                 dk=(q>>14)&127;dk-=64;
00140 
00141                 // Check that the worklist position is in range
00142                 ei=di+i;if(ei<0||ei>=hx) continue;
00143                 ej=dj+j;if(ej<0||ej>=hy) continue;
00144                 ek=dk+k;if(ek<0||ek>=hz) continue;
00145 
00146                 // Call the compute_min_max_radius() function. This returns
00147                 // true if the minimum distance to the block is bigger than the
00148                 // current mrs, in which case we skip this block and move on.
00149                 // Otherwise, it computes the maximum distance to the block and
00150                 // returns it in crs.
00151                 if(compute_min_radius(di,dj,dk,fx,fy,fz,mrs)) continue;
00152 
00153                 // Now compute which region we are going to loop over, adding a
00154                 // displacement for the periodic cases
00155                 ijk=con.region_index(ci,cj,ck,ei,ej,ek,qx,qy,qz,disp);
00156 
00157                 // If mrs is bigger than the maximum distance to the block,
00158                 // then we have to test all particles in the block for
00159                 // intersections. Otherwise, we do additional checks and skip
00160                 // those particles which can't possibly intersect the block.
00161                 scan_all(ijk,x-qx,y-qy,z-qz,di,dj,dk,w,mrs);
00162         } while(g<f);
00163 
00164         // Update mask value and initialize queue
00165         mv++;
00166         if(mv==0) {reset_mask();mv=1;}
00167         int *qu_s=qu,*qu_e=qu;
00168 
00169         while(g<wl_seq_length-1) {
00170 
00171                 // If mrs is less than the minimum distance to any untested
00172                 // block, then we are done
00173                 if(con.r_max_add(mrs)<radp[g]) return;
00174                 g++;
00175 
00176                 // Load in a block off the worklist, permute it with the
00177                 // symmetry mask, and decode its position. These are all
00178                 // integer bit operations so they should run very fast.
00179                 q=e[g];q^=m1;q+=m2;
00180                 di=q&127;di-=64;
00181                 dj=(q>>7)&127;dj-=64;
00182                 dk=(q>>14)&127;dk-=64;
00183 
00184                 // Compute the position in the mask of the current block. If
00185                 // this lies outside the mask, then skip it. Otherwise, mark
00186                 // it.
00187                 ei=di+i;if(ei<0||ei>=hx) continue;
00188                 ej=dj+j;if(ej<0||ej>=hy) continue;
00189                 ek=dk+k;if(ek<0||ek>=hz) continue;
00190                 mijk=mask+ei+hx*(ej+hy*ek);
00191                 *mijk=mv;
00192 
00193                 // Skip this block if it is further away than the current
00194                 // minimum radius
00195                 if(compute_min_radius(di,dj,dk,fx,fy,fz,mrs)) continue;
00196 
00197                 // Now compute which region we are going to loop over, adding a
00198                 // displacement for the periodic cases
00199                 ijk=con.region_index(ci,cj,ck,ei,ej,ek,qx,qy,qz,disp);
00200                 scan_all(ijk,x-qx,y-qy,z-qz,di,dj,dk,w,mrs);
00201 
00202                 if(qu_e>qu_l-18) add_list_memory(qu_s,qu_e);
00203                 scan_bits_mask_add(q,mijk,ei,ej,ek,qu_e);
00204         }
00205 
00206         // Do a check to see if we've reached the radius cutoff
00207         if(con.r_max_add(mrs)<radp[g]) return;
00208 
00209         // We were unable to completely compute the cell based on the blocks in
00210         // the worklist, so now we have to go block by block, reading in items
00211         // off the list
00212         while(qu_s!=qu_e) {
00213 
00214                 // Read the next entry of the queue
00215                 if(qu_s==qu_l) qu_s=qu;
00216                 ei=*(qu_s++);ej=*(qu_s++);ek=*(qu_s++);
00217                 di=ei-i;dj=ej-j;dk=ek-k;
00218                 if(compute_min_radius(di,dj,dk,fx,fy,fz,mrs)) continue;
00219 
00220                 ijk=con.region_index(ci,cj,ck,ei,ej,ek,qx,qy,qz,disp);
00221                 scan_all(ijk,x-qx,y-qy,z-qz,di,dj,dk,w,mrs);
00222 
00223                 // Test the neighbors of the current block, and add them to the
00224                 // block list if they haven't already been tested
00225                 if((qu_s<=qu_e?(qu_l-qu_e)+(qu_s-qu):qu_s-qu_e)<18) add_list_memory(qu_s,qu_e);
00226                 add_to_mask(ei,ej,ek,qu_e);
00227         }
00228 }
00229 
00230 /** Scans the six orthogonal neighbors of a given block and adds them to the
00231  * queue if they haven't been considered already. It assumes that the queue
00232  * will definitely have enough memory to add six entries at the end.
00233  * \param[in] (ei,ej,ek) the block to consider.
00234  * \param[in,out] qu_e a pointer to the end of the queue. */
00235 template<class c_class>
00236 inline void voro_compute<c_class>::add_to_mask(int ei,int ej,int ek,int *&qu_e) {
00237         unsigned int *mijk=mask+ei+hx*(ej+hy*ek);
00238         if(ek>0) if(*(mijk-hxy)!=mv) {if(qu_e==qu_l) qu_e=qu;*(mijk-hxy)=mv;*(qu_e++)=ei;*(qu_e++)=ej;*(qu_e++)=ek-1;}
00239         if(ej>0) if(*(mijk-hx)!=mv) {if(qu_e==qu_l) qu_e=qu;*(mijk-hx)=mv;*(qu_e++)=ei;*(qu_e++)=ej-1;*(qu_e++)=ek;}
00240         if(ei>0) if(*(mijk-1)!=mv) {if(qu_e==qu_l) qu_e=qu;*(mijk-1)=mv;*(qu_e++)=ei-1;*(qu_e++)=ej;*(qu_e++)=ek;}
00241         if(ei<hx-1) if(*(mijk+1)!=mv) {if(qu_e==qu_l) qu_e=qu;*(mijk+1)=mv;*(qu_e++)=ei+1;*(qu_e++)=ej;*(qu_e++)=ek;}
00242         if(ej<hy-1) if(*(mijk+hx)!=mv) {if(qu_e==qu_l) qu_e=qu;*(mijk+hx)=mv;*(qu_e++)=ei;*(qu_e++)=ej+1;*(qu_e++)=ek;}
00243         if(ek<hz-1) if(*(mijk+hxy)!=mv) {if(qu_e==qu_l) qu_e=qu;*(mijk+hxy)=mv;*(qu_e++)=ei;*(qu_e++)=ej;*(qu_e++)=ek+1;}
00244 }
00245 
00246 /** Scans a worklist entry and adds any blocks to the queue
00247  * \param[in] (ei,ej,ek) the block to consider.
00248  * \param[in,out] qu_e a pointer to the end of the queue. */
00249 template<class c_class>
00250 inline void voro_compute<c_class>::scan_bits_mask_add(unsigned int q,unsigned int *mijk,int ei,int ej,int ek,int *&qu_e) {
00251         const unsigned int b1=1<<21,b2=1<<22,b3=1<<24,b4=1<<25,b5=1<<27,b6=1<<28;
00252         if((q&b2)==b2) {
00253                 if(ei>0) {*(mijk-1)=mv;*(qu_e++)=ei-1;*(qu_e++)=ej;*(qu_e++)=ek;}
00254                 if((q&b1)==0&&ei<hx-1) {*(mijk+1)=mv;*(qu_e++)=ei+1;*(qu_e++)=ej;*(qu_e++)=ek;}
00255         } else if((q&b1)==b1&&ei<hx-1) {*(mijk+1)=mv;*(qu_e++)=ei+1;*(qu_e++)=ej;*(qu_e++)=ek;}
00256         if((q&b4)==b4) {
00257                 if(ej>0) {*(mijk-hx)=mv;*(qu_e++)=ei;*(qu_e++)=ej-1;*(qu_e++)=ek;}
00258                 if((q&b3)==0&&ej<hy-1) {*(mijk+hx)=mv;*(qu_e++)=ei;*(qu_e++)=ej+1;*(qu_e++)=ek;}
00259         } else if((q&b3)==b3&&ej<hy-1) {*(mijk+hx)=mv;*(qu_e++)=ei;*(qu_e++)=ej+1;*(qu_e++)=ek;}
00260         if((q&b6)==b6) {
00261                 if(ek>0) {*(mijk-hxy)=mv;*(qu_e++)=ei;*(qu_e++)=ej;*(qu_e++)=ek-1;}
00262                 if((q&b5)==0&&ek<hz-1) {*(mijk+hxy)=mv;*(qu_e++)=ei;*(qu_e++)=ej;*(qu_e++)=ek+1;}
00263         } else if((q&b5)==b5&&ek<hz-1) {*(mijk+hxy)=mv;*(qu_e++)=ei;*(qu_e++)=ej;*(qu_e++)=ek+1;}
00264 }
00265 
00266 /** This routine computes a Voronoi cell for a single particle in the
00267  * container. It can be called by the user, but is also forms the core part of
00268  * several of the main functions, such as store_cell_volumes(), print_all(),
00269  * and the drawing routines. The algorithm constructs the cell by testing over
00270  * the neighbors of the particle, working outwards until it reaches those
00271  * particles which could not possibly intersect the cell. For maximum
00272  * efficiency, this algorithm is divided into three parts. In the first
00273  * section, the algorithm tests over the blocks which are in the immediate
00274  * vicinity of the particle, by making use of one of the precomputed worklists.
00275  * The code then continues to test blocks on the worklist, but also begins to
00276  * construct a list of neighboring blocks outside the worklist which may need
00277  * to be test. In the third section, the routine starts testing these
00278  * neighboring blocks, evaluating whether or not a particle in them could
00279  * possibly intersect the cell. For blocks that intersect the cell, it tests
00280  * the particles in that block, and then adds the block neighbors to the list
00281  * of potential places to consider.
00282  * \param[in,out] c a reference to a voronoicell object.
00283  * \param[in] ijk the index of the block that the test particle is in.
00284  * \param[in] s the index of the particle within the test block.
00285  * \param[in] (ci,cj,ck) the coordinates of the block that the test particle is
00286  *                       in relative to the container data structure.
00287  * \return False if the Voronoi cell was completely removed during the
00288  *         computation and has zero volume, true otherwise. */
00289 template<class c_class>
00290 template<class v_cell>
00291 bool voro_compute<c_class>::compute_cell(v_cell &c,int ijk,int s,int ci,int cj,int ck) {
00292         static const int count_list[8]={7,11,15,19,26,35,45,59},*count_e=count_list+8;
00293         double x,y,z,x1,y1,z1,qx=0,qy=0,qz=0;
00294         double xlo,ylo,zlo,xhi,yhi,zhi,x2,y2,z2,rs;
00295         int i,j,k,di,dj,dk,ei,ej,ek,f,g,l,disp;
00296         double fx,fy,fz,gxs,gys,gzs,*radp;
00297         unsigned int q,*e,*mijk;
00298 
00299         if(!con.initialize_voronoicell(c,ijk,s,ci,cj,ck,i,j,k,x,y,z,disp)) return false;
00300         con.r_init(ijk,s);
00301 
00302         // Initialize the Voronoi cell to fill the entire container
00303         double crs,mrs;
00304 
00305         int next_count=3,*count_p=(const_cast<int*> (count_list));
00306 
00307         // Test all particles in the particle's local region first
00308         for(l=0;l<s;l++) {
00309                 x1=p[ijk][ps*l]-x;
00310                 y1=p[ijk][ps*l+1]-y;
00311                 z1=p[ijk][ps*l+2]-z;
00312                 rs=con.r_scale(x1*x1+y1*y1+z1*z1,ijk,l);
00313                 if(!c.nplane(x1,y1,z1,rs,id[ijk][l])) return false;
00314         }
00315         l++;
00316         while(l<co[ijk]) {
00317                 x1=p[ijk][ps*l]-x;
00318                 y1=p[ijk][ps*l+1]-y;
00319                 z1=p[ijk][ps*l+2]-z;
00320                 rs=con.r_scale(x1*x1+y1*y1+z1*z1,ijk,l);
00321                 if(!c.nplane(x1,y1,z1,rs,id[ijk][l])) return false;
00322                 l++;
00323         }
00324 
00325         // Now compute the maximum distance squared from the cell center to a
00326         // vertex. This is used to cut off the calculation since we only need
00327         // to test out to twice this range.
00328         mrs=c.max_radius_squared();
00329 
00330         // Now compute the fractional position of the particle within its
00331         // region and store it in (fx,fy,fz). We use this to compute an index
00332         // (di,dj,dk) of which subregion the particle is within.
00333         unsigned int m1,m2;
00334         con.frac_pos(x,y,z,ci,cj,ck,fx,fy,fz);
00335         di=int(fx*xsp*wl_fgrid);dj=int(fy*ysp*wl_fgrid);dk=int(fz*zsp*wl_fgrid);
00336 
00337         // The indices (di,dj,dk) tell us which worklist to use, to test the
00338         // blocks in the optimal order. But we only store worklists for the
00339         // eighth of the region where di, dj, and dk are all less than half the
00340         // full grid. The rest of the cases are handled by symmetry. In this
00341         // section, we detect for these cases, by reflecting high values of di,
00342         // dj, and dk. For these cases, a mask is constructed in m1 and m2
00343         // which is used to flip the worklist information when it is loaded.
00344         if(di>=wl_hgrid) {
00345                 gxs=fx;
00346                 m1=127+(3<<21);m2=1+(1<<21);di=wl_fgrid-1-di;if(di<0) di=0;
00347         } else {m1=m2=0;gxs=boxx-fx;}
00348         if(dj>=wl_hgrid) {
00349                 gys=fy;
00350                 m1|=(127<<7)+(3<<24);m2|=(1<<7)+(1<<24);dj=wl_fgrid-1-dj;if(dj<0) dj=0;
00351         } else gys=boxy-fy;
00352         if(dk>=wl_hgrid) {
00353                 gzs=fz;
00354                 m1|=(127<<14)+(3<<27);m2|=(1<<14)+(1<<27);dk=wl_fgrid-1-dk;if(dk<0) dk=0;
00355         } else gzs=boxz-fz;
00356         gxs*=gxs;gys*=gys;gzs*=gzs;
00357 
00358         // Now compute which worklist we are going to use, and set radp and e to
00359         // point at the right offsets
00360         ijk=di+wl_hgrid*(dj+wl_hgrid*dk);
00361         radp=mrad+ijk*wl_seq_length;
00362         e=(const_cast<unsigned int*> (wl))+ijk*wl_seq_length;
00363 
00364         // Read in how many items in the worklist can be tested without having to
00365         // worry about writing to the mask
00366         f=e[0];g=0;
00367         do {
00368 
00369                 // At the intervals specified by count_list, we recompute the
00370                 // maximum radius squared
00371                 if(g==next_count) {
00372                         mrs=c.max_radius_squared();
00373                         if(count_p!=count_e) next_count=*(count_p++);
00374                 }
00375 
00376                 // If mrs is less than the minimum distance to any untested
00377                 // block, then we are done
00378                 if(con.r_ctest(radp[g],mrs)) return true;
00379                 g++;
00380 
00381                 // Load in a block off the worklist, permute it with the
00382                 // symmetry mask, and decode its position. These are all
00383                 // integer bit operations so they should run very fast.
00384                 q=e[g];q^=m1;q+=m2;
00385                 di=q&127;di-=64;
00386                 dj=(q>>7)&127;dj-=64;
00387                 dk=(q>>14)&127;dk-=64;
00388 
00389                 // Check that the worklist position is in range
00390                 ei=di+i;if(ei<0||ei>=hx) continue;
00391                 ej=dj+j;if(ej<0||ej>=hy) continue;
00392                 ek=dk+k;if(ek<0||ek>=hz) continue;
00393 
00394                 // Call the compute_min_max_radius() function. This returns
00395                 // true if the minimum distance to the block is bigger than the
00396                 // current mrs, in which case we skip this block and move on.
00397                 // Otherwise, it computes the maximum distance to the block and
00398                 // returns it in crs.
00399                 if(compute_min_max_radius(di,dj,dk,fx,fy,fz,gxs,gys,gzs,crs,mrs)) continue;
00400 
00401                 // Now compute which region we are going to loop over, adding a
00402                 // displacement for the periodic cases
00403                 ijk=con.region_index(ci,cj,ck,ei,ej,ek,qx,qy,qz,disp);
00404 
00405                 // If mrs is bigger than the maximum distance to the block,
00406                 // then we have to test all particles in the block for
00407                 // intersections. Otherwise, we do additional checks and skip
00408                 // those particles which can't possibly intersect the block.
00409                 if(co[ijk]>0) {
00410                         l=0;x2=x-qx;y2=y-qy;z2=z-qz;
00411                         if(!con.r_ctest(crs,mrs)) {
00412                                 do {
00413                                         x1=p[ijk][ps*l]-x2;
00414                                         y1=p[ijk][ps*l+1]-y2;
00415                                         z1=p[ijk][ps*l+2]-z2;
00416                                         rs=con.r_scale(x1*x1+y1*y1+z1*z1,ijk,l);
00417                                         if(!c.nplane(x1,y1,z1,rs,id[ijk][l])) return false;
00418                                         l++;
00419                                 } while (l<co[ijk]);
00420                         } else {
00421                                 do {
00422                                         x1=p[ijk][ps*l]-x2;
00423                                         y1=p[ijk][ps*l+1]-y2;
00424                                         z1=p[ijk][ps*l+2]-z2;
00425                                         rs=x1*x1+y1*y1+z1*z1;
00426                                         if(con.r_scale_check(rs,mrs,ijk,l)&&!c.nplane(x1,y1,z1,rs,id[ijk][l])) return false;
00427                                         l++;
00428                                 } while (l<co[ijk]);
00429                         }
00430                 }
00431         } while(g<f);
00432 
00433         // If we reach here, we were unable to compute the entire cell using
00434         // the first part of the worklist. This section of the algorithm
00435         // continues the worklist, but it now starts preparing the mask that we
00436         // need if we end up going block by block. We do the same as before,
00437         // but we put a mark down on the mask for every block that's tested.
00438         // The worklist also contains information about which neighbors of each
00439         // block are not also on the worklist, and we start storing those
00440         // points in a list in case we have to go block by block. Update the
00441         // mask counter, and if it wraps around then reset the whole mask; that
00442         // will only happen once every 2^32 tries.
00443         mv++;
00444         if(mv==0) {reset_mask();mv=1;}
00445 
00446         // Set the queue pointers
00447         int *qu_s=qu,*qu_e=qu;
00448 
00449         while(g<wl_seq_length-1) {
00450 
00451                 // At the intervals specified by count_list, we recompute the
00452                 // maximum radius squared
00453                 if(g==next_count) {
00454                         mrs=c.max_radius_squared();
00455                         if(count_p!=count_e) next_count=*(count_p++);
00456                 }
00457 
00458                 // If mrs is less than the minimum distance to any untested
00459                 // block, then we are done
00460                 if(con.r_ctest(radp[g],mrs)) return true;
00461                 g++;
00462 
00463                 // Load in a block off the worklist, permute it with the
00464                 // symmetry mask, and decode its position. These are all
00465                 // integer bit operations so they should run very fast.
00466                 q=e[g];q^=m1;q+=m2;
00467                 di=q&127;di-=64;
00468                 dj=(q>>7)&127;dj-=64;
00469                 dk=(q>>14)&127;dk-=64;
00470 
00471                 // Compute the position in the mask of the current block. If
00472                 // this lies outside the mask, then skip it. Otherwise, mark
00473                 // it.
00474                 ei=di+i;if(ei<0||ei>=hx) continue;
00475                 ej=dj+j;if(ej<0||ej>=hy) continue;
00476                 ek=dk+k;if(ek<0||ek>=hz) continue;
00477                 mijk=mask+ei+hx*(ej+hy*ek);
00478                 *mijk=mv;
00479 
00480                 // Call the compute_min_max_radius() function. This returns
00481                 // true if the minimum distance to the block is bigger than the
00482                 // current mrs, in which case we skip this block and move on.
00483                 // Otherwise, it computes the maximum distance to the block and
00484                 // returns it in crs.
00485                 if(compute_min_max_radius(di,dj,dk,fx,fy,fz,gxs,gys,gzs,crs,mrs)) continue;
00486 
00487                 // Now compute which region we are going to loop over, adding a
00488                 // displacement for the periodic cases
00489                 ijk=con.region_index(ci,cj,ck,ei,ej,ek,qx,qy,qz,disp);
00490 
00491                 // If mrs is bigger than the maximum distance to the block,
00492                 // then we have to test all particles in the block for
00493                 // intersections. Otherwise, we do additional checks and skip
00494                 // those particles which can't possibly intersect the block.
00495                 if(co[ijk]>0) {
00496                         l=0;x2=x-qx;y2=y-qy;z2=z-qz;
00497                         if(!con.r_ctest(crs,mrs)) {
00498                                 do {
00499                                         x1=p[ijk][ps*l]-x2;
00500                                         y1=p[ijk][ps*l+1]-y2;
00501                                         z1=p[ijk][ps*l+2]-z2;
00502                                         rs=con.r_scale(x1*x1+y1*y1+z1*z1,ijk,l);
00503                                         if(!c.nplane(x1,y1,z1,rs,id[ijk][l])) return false;
00504                                         l++;
00505                                 } while (l<co[ijk]);
00506                         } else {
00507                                 do {
00508                                         x1=p[ijk][ps*l]-x2;
00509                                         y1=p[ijk][ps*l+1]-y2;
00510                                         z1=p[ijk][ps*l+2]-z2;
00511                                         rs=x1*x1+y1*y1+z1*z1;
00512                                         if(con.r_scale_check(rs,mrs,ijk,l)&&!c.nplane(x1,y1,z1,rs,id[ijk][l])) return false;
00513                                         l++;
00514                                 } while (l<co[ijk]);
00515                         }
00516                 }
00517 
00518                 // If there might not be enough memory on the list for these
00519                 // additions, then add more
00520                 if(qu_e>qu_l-18) add_list_memory(qu_s,qu_e);
00521 
00522                 // Test the parts of the worklist element which tell us what
00523                 // neighbors of this block are not on the worklist. Store them
00524                 // on the block list, and mark the mask.
00525                 scan_bits_mask_add(q,mijk,ei,ej,ek,qu_e);
00526         }
00527 
00528         // Do a check to see if we've reached the radius cutoff
00529         if(con.r_ctest(radp[g],mrs)) return true;
00530 
00531         // We were unable to completely compute the cell based on the blocks in
00532         // the worklist, so now we have to go block by block, reading in items
00533         // off the list
00534         while(qu_s!=qu_e) {
00535 
00536                 // If we reached the end of the list memory loop back to the
00537                 // start
00538                 if(qu_s==qu_l) qu_s=qu;
00539 
00540                 // Read in a block off the list, and compute the upper and lower
00541                 // coordinates in each of the three dimensions
00542                 ei=*(qu_s++);ej=*(qu_s++);ek=*(qu_s++);
00543                 xlo=(ei-i)*boxx-fx;xhi=xlo+boxx;
00544                 ylo=(ej-j)*boxy-fy;yhi=ylo+boxy;
00545                 zlo=(ek-k)*boxz-fz;zhi=zlo+boxz;
00546 
00547                 // Carry out plane tests to see if any particle in this block
00548                 // could possibly intersect the cell
00549                 if(ei>i) {
00550                         if(ej>j) {
00551                                 if(ek>k) {if(corner_test(c,xlo,ylo,zlo,xhi,yhi,zhi)) continue;}
00552                                 else if(ek<k) {if(corner_test(c,xlo,ylo,zhi,xhi,yhi,zlo)) continue;}
00553                                 else {if(edge_z_test(c,xlo,ylo,zlo,xhi,yhi,zhi)) continue;}
00554                         } else if(ej<j) {
00555                                 if(ek>k) {if(corner_test(c,xlo,yhi,zlo,xhi,ylo,zhi)) continue;}
00556                                 else if(ek<k) {if(corner_test(c,xlo,yhi,zhi,xhi,ylo,zlo)) continue;}
00557                                 else {if(edge_z_test(c,xlo,yhi,zlo,xhi,ylo,zhi)) continue;}
00558                         } else {
00559                                 if(ek>k) {if(edge_y_test(c,xlo,ylo,zlo,xhi,yhi,zhi)) continue;}
00560                                 else if(ek<k) {if(edge_y_test(c,xlo,ylo,zhi,xhi,yhi,zlo)) continue;}
00561                                 else {if(face_x_test(c,xlo,ylo,zlo,yhi,zhi)) continue;}
00562                         }
00563                 } else if(ei<i) {
00564                         if(ej>j) {
00565                                 if(ek>k) {if(corner_test(c,xhi,ylo,zlo,xlo,yhi,zhi)) continue;}
00566                                 else if(ek<k) {if(corner_test(c,xhi,ylo,zhi,xlo,yhi,zlo)) continue;}
00567                                 else {if(edge_z_test(c,xhi,ylo,zlo,xlo,yhi,zhi)) continue;}
00568                         } else if(ej<j) {
00569                                 if(ek>k) {if(corner_test(c,xhi,yhi,zlo,xlo,ylo,zhi)) continue;}
00570                                 else if(ek<k) {if(corner_test(c,xhi,yhi,zhi,xlo,ylo,zlo)) continue;}
00571                                 else {if(edge_z_test(c,xhi,yhi,zlo,xlo,ylo,zhi)) continue;}
00572                         } else {
00573                                 if(ek>k) {if(edge_y_test(c,xhi,ylo,zlo,xlo,yhi,zhi)) continue;}
00574                                 else if(ek<k) {if(edge_y_test(c,xhi,ylo,zhi,xlo,yhi,zlo)) continue;}
00575                                 else {if(face_x_test(c,xhi,ylo,zlo,yhi,zhi)) continue;}
00576                         }
00577                 } else {
00578                         if(ej>j) {
00579                                 if(ek>k) {if(edge_x_test(c,xlo,ylo,zlo,xhi,yhi,zhi)) continue;}
00580                                 else if(ek<k) {if(edge_x_test(c,xlo,ylo,zhi,xhi,yhi,zlo)) continue;}
00581                                 else {if(face_y_test(c,xlo,ylo,zlo,xhi,zhi)) continue;}
00582                         } else if(ej<j) {
00583                                 if(ek>k) {if(edge_x_test(c,xlo,yhi,zlo,xhi,ylo,zhi)) continue;}
00584                                 else if(ek<k) {if(edge_x_test(c,xlo,yhi,zhi,xhi,ylo,zlo)) continue;}
00585                                 else {if(face_y_test(c,xlo,yhi,zlo,xhi,zhi)) continue;}
00586                         } else {
00587                                 if(ek>k) {if(face_z_test(c,xlo,ylo,zlo,xhi,yhi)) continue;}
00588                                 else if(ek<k) {if(face_z_test(c,xlo,ylo,zhi,xhi,yhi)) continue;}
00589                                 else voro_fatal_error("Compute cell routine revisiting central block, which should never\nhappen.",VOROPP_INTERNAL_ERROR);
00590                         }
00591                 }
00592 
00593                 // Now compute the region that we are going to test over, and
00594                 // set a displacement vector for the periodic cases
00595                 ijk=con.region_index(ci,cj,ck,ei,ej,ek,qx,qy,qz,disp);
00596 
00597                 // Loop over all the elements in the block to test for cuts. It
00598                 // would be possible to exclude some of these cases by testing
00599                 // against mrs, but this will probably not save time.
00600                 if(co[ijk]>0) {
00601                         l=0;x2=x-qx;y2=y-qy;z2=z-qz;
00602                         do {
00603                                 x1=p[ijk][ps*l]-x2;
00604                                 y1=p[ijk][ps*l+1]-y2;
00605                                 z1=p[ijk][ps*l+2]-z2;
00606                                 rs=con.r_scale(x1*x1+y1*y1+z1*z1,ijk,l);
00607                                 if(!c.nplane(x1,y1,z1,rs,id[ijk][l])) return false;
00608                                 l++;
00609                         } while (l<co[ijk]);
00610                 }
00611 
00612                 // If there's not much memory on the block list then add more
00613                 if((qu_s<=qu_e?(qu_l-qu_e)+(qu_s-qu):qu_s-qu_e)<18) add_list_memory(qu_s,qu_e);
00614 
00615                 // Test the neighbors of the current block, and add them to the
00616                 // block list if they haven't already been tested
00617                 add_to_mask(ei,ej,ek,qu_e);
00618         }
00619 
00620         return true;
00621 }
00622 
00623 /** This function checks to see whether a particular block can possibly have
00624  * any intersection with a Voronoi cell, for the case when the closest point
00625  * from the cell center to the block is at a corner.
00626  * \param[in,out] c a reference to a Voronoi cell.
00627  * \param[in] (xl,yl,zl) the relative coordinates of the corner of the block
00628  *                       closest to the cell center.
00629  * \param[in] (xh,yh,zh) the relative coordinates of the corner of the block
00630  *                       furthest away from the cell center.
00631  * \return False if the block may intersect, true if does not. */
00632 template<class c_class>
00633 template<class v_cell>
00634 bool voro_compute<c_class>::corner_test(v_cell &c,double xl,double yl,double zl,double xh,double yh,double zh) {
00635         con.r_prime(xl*xl+yl*yl+zl*zl);
00636         if(c.plane_intersects_guess(xh,yl,zl,con.r_cutoff(xl*xh+yl*yl+zl*zl))) return false;
00637         if(c.plane_intersects(xh,yh,zl,con.r_cutoff(xl*xh+yl*yh+zl*zl))) return false;
00638         if(c.plane_intersects(xl,yh,zl,con.r_cutoff(xl*xl+yl*yh+zl*zl))) return false;
00639         if(c.plane_intersects(xl,yh,zh,con.r_cutoff(xl*xl+yl*yh+zl*zh))) return false;
00640         if(c.plane_intersects(xl,yl,zh,con.r_cutoff(xl*xl+yl*yl+zl*zh))) return false;
00641         if(c.plane_intersects(xh,yl,zh,con.r_cutoff(xl*xh+yl*yl+zl*zh))) return false;
00642         return true;
00643 }
00644 
00645 /** This function checks to see whether a particular block can possibly have
00646  * any intersection with a Voronoi cell, for the case when the closest point
00647  * from the cell center to the block is on an edge which points along the x
00648  * direction.
00649  * \param[in,out] c a reference to a Voronoi cell.
00650  * \param[in] (x0,x1) the minimum and maximum relative x coordinates of the
00651  *                    block.
00652  * \param[in] (yl,zl) the relative y and z coordinates of the corner of the
00653  *                    block closest to the cell center.
00654  * \param[in] (yh,zh) the relative y and z coordinates of the corner of the
00655  *                    block furthest away from the cell center.
00656  * \return False if the block may intersect, true if does not. */
00657 template<class c_class>
00658 template<class v_cell>
00659 inline bool voro_compute<c_class>::edge_x_test(v_cell &c,double x0,double yl,double zl,double x1,double yh,double zh) {
00660         con.r_prime(yl*yl+zl*zl);
00661         if(c.plane_intersects_guess(x0,yl,zh,con.r_cutoff(yl*yl+zl*zh))) return false;
00662         if(c.plane_intersects(x1,yl,zh,con.r_cutoff(yl*yl+zl*zh))) return false;
00663         if(c.plane_intersects(x1,yl,zl,con.r_cutoff(yl*yl+zl*zl))) return false;
00664         if(c.plane_intersects(x0,yl,zl,con.r_cutoff(yl*yl+zl*zl))) return false;
00665         if(c.plane_intersects(x0,yh,zl,con.r_cutoff(yl*yh+zl*zl))) return false;
00666         if(c.plane_intersects(x1,yh,zl,con.r_cutoff(yl*yh+zl*zl))) return false;
00667         return true;
00668 }
00669 
00670 /** This function checks to see whether a particular block can possibly have
00671  * any intersection with a Voronoi cell, for the case when the closest point
00672  * from the cell center to the block is on an edge which points along the y
00673  * direction.
00674  * \param[in,out] c a reference to a Voronoi cell.
00675  * \param[in] (y0,y1) the minimum and maximum relative y coordinates of the
00676  *                    block.
00677  * \param[in] (xl,zl) the relative x and z coordinates of the corner of the
00678  *                    block closest to the cell center.
00679  * \param[in] (xh,zh) the relative x and z coordinates of the corner of the
00680  *                    block furthest away from the cell center.
00681  * \return False if the block may intersect, true if does not. */
00682 template<class c_class>
00683 template<class v_cell>
00684 inline bool voro_compute<c_class>::edge_y_test(v_cell &c,double xl,double y0,double zl,double xh,double y1,double zh) {
00685         con.r_prime(xl*xl+zl*zl);
00686         if(c.plane_intersects_guess(xl,y0,zh,con.r_cutoff(xl*xl+zl*zh))) return false;
00687         if(c.plane_intersects(xl,y1,zh,con.r_cutoff(xl*xl+zl*zh))) return false;
00688         if(c.plane_intersects(xl,y1,zl,con.r_cutoff(xl*xl+zl*zl))) return false;
00689         if(c.plane_intersects(xl,y0,zl,con.r_cutoff(xl*xl+zl*zl))) return false;
00690         if(c.plane_intersects(xh,y0,zl,con.r_cutoff(xl*xh+zl*zl))) return false;
00691         if(c.plane_intersects(xh,y1,zl,con.r_cutoff(xl*xh+zl*zl))) return false;
00692         return true;
00693 }
00694 
00695 /** This function checks to see whether a particular block can possibly have
00696  * any intersection with a Voronoi cell, for the case when the closest point
00697  * from the cell center to the block is on an edge which points along the z
00698  * direction.
00699  * \param[in,out] c a reference to a Voronoi cell.
00700  * \param[in] (z0,z1) the minimum and maximum relative z coordinates of the block.
00701  * \param[in] (xl,yl) the relative x and y coordinates of the corner of the
00702  *                    block closest to the cell center.
00703  * \param[in] (xh,yh) the relative x and y coordinates of the corner of the
00704  *                    block furthest away from the cell center.
00705  * \return False if the block may intersect, true if does not. */
00706 template<class c_class>
00707 template<class v_cell>
00708 inline bool voro_compute<c_class>::edge_z_test(v_cell &c,double xl,double yl,double z0,double xh,double yh,double z1) {
00709         con.r_prime(xl*xl+yl*yl);
00710         if(c.plane_intersects_guess(xl,yh,z0,con.r_cutoff(xl*xl+yl*yh))) return false;
00711         if(c.plane_intersects(xl,yh,z1,con.r_cutoff(xl*xl+yl*yh))) return false;
00712         if(c.plane_intersects(xl,yl,z1,con.r_cutoff(xl*xl+yl*yl))) return false;
00713         if(c.plane_intersects(xl,yl,z0,con.r_cutoff(xl*xl+yl*yl))) return false;
00714         if(c.plane_intersects(xh,yl,z0,con.r_cutoff(xl*xh+yl*yl))) return false;
00715         if(c.plane_intersects(xh,yl,z1,con.r_cutoff(xl*xh+yl*yl))) return false;
00716         return true;
00717 }
00718 
00719 /** This function checks to see whether a particular block can possibly have
00720  * any intersection with a Voronoi cell, for the case when the closest point
00721  * from the cell center to the block is on a face aligned with the x direction.
00722  * \param[in,out] c a reference to a Voronoi cell.
00723  * \param[in] xl the minimum distance from the cell center to the face.
00724  * \param[in] (y0,y1) the minimum and maximum relative y coordinates of the
00725  *                    block.
00726  * \param[in] (z0,z1) the minimum and maximum relative z coordinates of the
00727  *                    block.
00728  * \return False if the block may intersect, true if does not. */
00729 template<class c_class>
00730 template<class v_cell>
00731 inline bool voro_compute<c_class>::face_x_test(v_cell &c,double xl,double y0,double z0,double y1,double z1) {
00732         con.r_prime(xl*xl);
00733         if(c.plane_intersects_guess(xl,y0,z0,con.r_cutoff(xl*xl))) return false;
00734         if(c.plane_intersects(xl,y0,z1,con.r_cutoff(xl*xl))) return false;
00735         if(c.plane_intersects(xl,y1,z1,con.r_cutoff(xl*xl))) return false;
00736         if(c.plane_intersects(xl,y1,z0,con.r_cutoff(xl*xl))) return false;
00737         return true;
00738 }
00739 
00740 /** This function checks to see whether a particular block can possibly have
00741  * any intersection with a Voronoi cell, for the case when the closest point
00742  * from the cell center to the block is on a face aligned with the y direction.
00743  * \param[in,out] c a reference to a Voronoi cell.
00744  * \param[in] yl the minimum distance from the cell center to the face.
00745  * \param[in] (x0,x1) the minimum and maximum relative x coordinates of the
00746  *                    block.
00747  * \param[in] (z0,z1) the minimum and maximum relative z coordinates of the
00748  *                    block.
00749  * \return False if the block may intersect, true if does not. */
00750 template<class c_class>
00751 template<class v_cell>
00752 inline bool voro_compute<c_class>::face_y_test(v_cell &c,double x0,double yl,double z0,double x1,double z1) {
00753         con.r_prime(yl*yl);
00754         if(c.plane_intersects_guess(x0,yl,z0,con.r_cutoff(yl*yl))) return false;
00755         if(c.plane_intersects(x0,yl,z1,con.r_cutoff(yl*yl))) return false;
00756         if(c.plane_intersects(x1,yl,z1,con.r_cutoff(yl*yl))) return false;
00757         if(c.plane_intersects(x1,yl,z0,con.r_cutoff(yl*yl))) return false;
00758         return true;
00759 }
00760 
00761 /** This function checks to see whether a particular block can possibly have
00762  * any intersection with a Voronoi cell, for the case when the closest point
00763  * from the cell center to the block is on a face aligned with the z direction.
00764  * \param[in,out] c a reference to a Voronoi cell.
00765  * \param[in] zl the minimum distance from the cell center to the face.
00766  * \param[in] (x0,x1) the minimum and maximum relative x coordinates of the
00767  *                    block.
00768  * \param[in] (y0,y1) the minimum and maximum relative y coordinates of the
00769  *                    block.
00770  * \return False if the block may intersect, true if does not. */
00771 template<class c_class>
00772 template<class v_cell>
00773 inline bool voro_compute<c_class>::face_z_test(v_cell &c,double x0,double y0,double zl,double x1,double y1) {
00774         con.r_prime(zl*zl);
00775         if(c.plane_intersects_guess(x0,y0,zl,con.r_cutoff(zl*zl))) return false;
00776         if(c.plane_intersects(x0,y1,zl,con.r_cutoff(zl*zl))) return false;
00777         if(c.plane_intersects(x1,y1,zl,con.r_cutoff(zl*zl))) return false;
00778         if(c.plane_intersects(x1,y0,zl,con.r_cutoff(zl*zl))) return false;
00779         return true;
00780 }
00781 
00782 
00783 /** This routine checks to see whether a point is within a particular distance
00784  * of a nearby region. If the point is within the distance of the region, then
00785  * the routine returns true, and computes the maximum distance from the point
00786  * to the region. Otherwise, the routine returns false.
00787  * \param[in] (di,dj,dk) the position of the nearby region to be tested,
00788  *                       relative to the region that the point is in.
00789  * \param[in] (fx,fy,fz) the displacement of the point within its region.
00790  * \param[in] (gxs,gys,gzs) the maximum squared distances from the point to the
00791  *                          sides of its region.
00792  * \param[out] crs a reference in which to return the maximum distance to the
00793  *                 region (only computed if the routine returns false).
00794  * \param[in] mrs the distance to be tested.
00795  * \return True if the region is further away than mrs, false if the region in
00796  *         within mrs. */
00797 template<class c_class>
00798 bool voro_compute<c_class>::compute_min_max_radius(int di,int dj,int dk,double fx,double fy,double fz,double gxs,double gys,double gzs,double &crs,double mrs) {
00799         double xlo,ylo,zlo;
00800         if(di>0) {
00801                 xlo=di*boxx-fx;
00802                 crs=xlo*xlo;
00803                 if(dj>0) {
00804                         ylo=dj*boxy-fy;
00805                         crs+=ylo*ylo;
00806                         if(dk>0) {
00807                                 zlo=dk*boxz-fz;
00808                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00809                                 crs+=bxsq+2*(boxx*xlo+boxy*ylo+boxz*zlo);
00810                         } else if(dk<0) {
00811                                 zlo=(dk+1)*boxz-fz;
00812                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00813                                 crs+=bxsq+2*(boxx*xlo+boxy*ylo-boxz*zlo);
00814                         } else {
00815                                 if(con.r_ctest(crs,mrs)) return true;
00816                                 crs+=boxx*(2*xlo+boxx)+boxy*(2*ylo+boxy)+gzs;
00817                         }
00818                 } else if(dj<0) {
00819                         ylo=(dj+1)*boxy-fy;
00820                         crs+=ylo*ylo;
00821                         if(dk>0) {
00822                                 zlo=dk*boxz-fz;
00823                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00824                                 crs+=bxsq+2*(boxx*xlo-boxy*ylo+boxz*zlo);
00825                         } else if(dk<0) {
00826                                 zlo=(dk+1)*boxz-fz;
00827                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00828                                 crs+=bxsq+2*(boxx*xlo-boxy*ylo-boxz*zlo);
00829                         } else {
00830                                 if(con.r_ctest(crs,mrs)) return true;
00831                                 crs+=boxx*(2*xlo+boxx)+boxy*(-2*ylo+boxy)+gzs;
00832                         }
00833                 } else {
00834                         if(dk>0) {
00835                                 zlo=dk*boxz-fz;
00836                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00837                                 crs+=boxz*(2*zlo+boxz);
00838                         } else if(dk<0) {
00839                                 zlo=(dk+1)*boxz-fz;
00840                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00841                                 crs+=boxz*(-2*zlo+boxz);
00842                         } else {
00843                                 if(con.r_ctest(crs,mrs)) return true;
00844                                 crs+=gzs;
00845                         }
00846                         crs+=gys+boxx*(2*xlo+boxx);
00847                 }
00848         } else if(di<0) {
00849                 xlo=(di+1)*boxx-fx;
00850                 crs=xlo*xlo;
00851                 if(dj>0) {
00852                         ylo=dj*boxy-fy;
00853                         crs+=ylo*ylo;
00854                         if(dk>0) {
00855                                 zlo=dk*boxz-fz;
00856                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00857                                 crs+=bxsq+2*(-boxx*xlo+boxy*ylo+boxz*zlo);
00858                         } else if(dk<0) {
00859                                 zlo=(dk+1)*boxz-fz;
00860                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00861                                 crs+=bxsq+2*(-boxx*xlo+boxy*ylo-boxz*zlo);
00862                         } else {
00863                                 if(con.r_ctest(crs,mrs)) return true;
00864                                 crs+=boxx*(-2*xlo+boxx)+boxy*(2*ylo+boxy)+gzs;
00865                         }
00866                 } else if(dj<0) {
00867                         ylo=(dj+1)*boxy-fy;
00868                         crs+=ylo*ylo;
00869                         if(dk>0) {
00870                                 zlo=dk*boxz-fz;
00871                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00872                                 crs+=bxsq+2*(-boxx*xlo-boxy*ylo+boxz*zlo);
00873                         } else if(dk<0) {
00874                                 zlo=(dk+1)*boxz-fz;
00875                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00876                                 crs+=bxsq+2*(-boxx*xlo-boxy*ylo-boxz*zlo);
00877                         } else {
00878                                 if(con.r_ctest(crs,mrs)) return true;
00879                                 crs+=boxx*(-2*xlo+boxx)+boxy*(-2*ylo+boxy)+gzs;
00880                         }
00881                 } else {
00882                         if(dk>0) {
00883                                 zlo=dk*boxz-fz;
00884                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00885                                 crs+=boxz*(2*zlo+boxz);
00886                         } else if(dk<0) {
00887                                 zlo=(dk+1)*boxz-fz;
00888                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00889                                 crs+=boxz*(-2*zlo+boxz);
00890                         } else {
00891                                 if(con.r_ctest(crs,mrs)) return true;
00892                                 crs+=gzs;
00893                         }
00894                         crs+=gys+boxx*(-2*xlo+boxx);
00895                 }
00896         } else {
00897                 if(dj>0) {
00898                         ylo=dj*boxy-fy;
00899                         crs=ylo*ylo;
00900                         if(dk>0) {
00901                                 zlo=dk*boxz-fz;
00902                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00903                                 crs+=boxz*(2*zlo+boxz);
00904                         } else if(dk<0) {
00905                                 zlo=(dk+1)*boxz-fz;
00906                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00907                                 crs+=boxz*(-2*zlo+boxz);
00908                         } else {
00909                                 if(con.r_ctest(crs,mrs)) return true;
00910                                 crs+=gzs;
00911                         }
00912                         crs+=boxy*(2*ylo+boxy);
00913                 } else if(dj<0) {
00914                         ylo=(dj+1)*boxy-fy;
00915                         crs=ylo*ylo;
00916                         if(dk>0) {
00917                                 zlo=dk*boxz-fz;
00918                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00919                                 crs+=boxz*(2*zlo+boxz);
00920                         } else if(dk<0) {
00921                                 zlo=(dk+1)*boxz-fz;
00922                                 crs+=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00923                                 crs+=boxz*(-2*zlo+boxz);
00924                         } else {
00925                                 if(con.r_ctest(crs,mrs)) return true;
00926                                 crs+=gzs;
00927                         }
00928                         crs+=boxy*(-2*ylo+boxy);
00929                 } else {
00930                         if(dk>0) {
00931                                 zlo=dk*boxz-fz;crs=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00932                                 crs+=boxz*(2*zlo+boxz);
00933                         } else if(dk<0) {
00934                                 zlo=(dk+1)*boxz-fz;crs=zlo*zlo;if(con.r_ctest(crs,mrs)) return true;
00935                                 crs+=boxz*(-2*zlo+boxz);
00936                         } else {
00937                                 crs=0;
00938                                 voro_fatal_error("Min/max radius function called for central block, which should never\nhappen.",VOROPP_INTERNAL_ERROR);
00939                         }
00940                         crs+=gys;
00941                 }
00942                 crs+=gxs;
00943         }
00944         return false;
00945 }
00946 
00947 template<class c_class>
00948 bool voro_compute<c_class>::compute_min_radius(int di,int dj,int dk,double fx,double fy,double fz,double mrs) {
00949         double t,crs;
00950 
00951         if(di>0) {t=di*boxx-fx;crs=t*t;}
00952         else if(di<0) {t=(di+1)*boxx-fx;crs=t*t;}
00953         else crs=0;
00954 
00955         if(dj>0) {t=dj*boxy-fy;crs+=t*t;}
00956         else if(dj<0) {t=(dj+1)*boxy-fy;crs+=t*t;}
00957 
00958         if(dk>0) {t=dk*boxz-fz;crs+=t*t;}
00959         else if(dk<0) {t=(dk+1)*boxz-fz;crs+=t*t;}
00960 
00961         return crs>con.r_max_add(mrs);
00962 }
00963 
00964 /** Adds memory to the queue.
00965  * \param[in,out] qu_s a reference to the queue start pointer.
00966  * \param[in,out] qu_e a reference to the queue end pointer. */
00967 template<class c_class>
00968 inline void voro_compute<c_class>::add_list_memory(int*& qu_s,int*& qu_e) {
00969         qu_size<<=1;
00970         int *qu_n=new int[qu_size],*qu_c=qu_n;
00971 #if VOROPP_VERBOSE >=2
00972         fprintf(stderr,"List memory scaled up to %d\n",qu_size);
00973 #endif
00974         if(qu_s<=qu_e) {
00975                 while(qu_s<qu_e) *(qu_c++)=*(qu_s++);
00976         } else {
00977                 while(qu_s<qu_l) *(qu_c++)=*(qu_s++);qu_s=qu;
00978                 while(qu_s<qu_e) *(qu_c++)=*(qu_s++);
00979         }
00980         delete [] qu;
00981         qu_s=qu=qu_n;
00982         qu_l=qu+qu_size;
00983         qu_e=qu_c;
00984 }
00985 
00986 // Explicit template instantiation
00987 template voro_compute<container>::voro_compute(container&,int,int,int);
00988 template voro_compute<container_poly>::voro_compute(container_poly&,int,int,int);
00989 template bool voro_compute<container>::compute_cell(voronoicell&,int,int,int,int,int);
00990 template bool voro_compute<container>::compute_cell(voronoicell_neighbor&,int,int,int,int,int);
00991 template void voro_compute<container>::find_voronoi_cell(double,double,double,int,int,int,int,particle_record&,double&);
00992 template bool voro_compute<container_poly>::compute_cell(voronoicell&,int,int,int,int,int);
00993 template bool voro_compute<container_poly>::compute_cell(voronoicell_neighbor&,int,int,int,int,int);
00994 template void voro_compute<container_poly>::find_voronoi_cell(double,double,double,int,int,int,int,particle_record&,double&);
00995 
00996 // Explicit template instantiation
00997 template voro_compute<container_periodic>::voro_compute(container_periodic&,int,int,int);
00998 template voro_compute<container_periodic_poly>::voro_compute(container_periodic_poly&,int,int,int);
00999 template bool voro_compute<container_periodic>::compute_cell(voronoicell&,int,int,int,int,int);
01000 template bool voro_compute<container_periodic>::compute_cell(voronoicell_neighbor&,int,int,int,int,int);
01001 template void voro_compute<container_periodic>::find_voronoi_cell(double,double,double,int,int,int,int,particle_record&,double&);
01002 template bool voro_compute<container_periodic_poly>::compute_cell(voronoicell&,int,int,int,int,int);
01003 template bool voro_compute<container_periodic_poly>::compute_cell(voronoicell_neighbor&,int,int,int,int,int);
01004 template void voro_compute<container_periodic_poly>::find_voronoi_cell(double,double,double,int,int,int,int,particle_record&,double&);
01005 
01006 }