#include #define STOKES extern "C" void InitDenColor(char *ID, double *Den, double *Phi, double das, double dbs, int Nx, int Ny, int Nz) { int i,j,k,n,N; N = Nx*Ny*Nz; for (n=0; n 0){ // Retrieve the color gradient nx = ColorGrad[n]; ny = ColorGrad[N+n]; nz = ColorGrad[2*N+n]; //...........Normalize the Color Gradient................................. C = sqrt(nx*nx+ny*ny+nz*nz); if (C==0.0) C=1.0; nx = nx/C; ny = ny/C; nz = nz/C; //......No color gradient at z-boundary if pressure BC are set............. // if (pBC && k==0) nx = ny = nz = 0.f; // if (pBC && k==Nz-1) nx = ny = nz = 0.f; //........................................................................ // READ THE DISTRIBUTIONS // (read from opposite array due to previous swap operation) //........................................................................ f2 = distodd[n]; f4 = distodd[N+n]; f6 = distodd[2*N+n]; f8 = distodd[3*N+n]; f10 = distodd[4*N+n]; f12 = distodd[5*N+n]; f14 = distodd[6*N+n]; f16 = distodd[7*N+n]; f18 = distodd[8*N+n]; //........................................................................ f0 = disteven[n]; f1 = disteven[N+n]; f3 = disteven[2*N+n]; f5 = disteven[3*N+n]; f7 = disteven[4*N+n]; f9 = disteven[5*N+n]; f11 = disteven[6*N+n]; f13 = disteven[7*N+n]; f15 = disteven[8*N+n]; f17 = disteven[9*N+n]; //........................................................................ // PERFORM RELAXATION PROCESS //........................................................................ //....................compute the moments............................................... rho = f0+f2+f1+f4+f3+f6+f5+f8+f7+f10+f9+f12+f11+f14+f13+f16+f15+f18+f17; m1 = -30*f0-11*(f2+f1+f4+f3+f6+f5)+8*(f8+f7+f10+f9+f12+f11+f14+f13+f16+f15+f18 +f17); m2 = 12*f0-4*(f2+f1 +f4+f3+f6 +f5)+f8+f7+f10+f9+f12+f11+f14+f13+f16+f15+f18+f17; jx = f1-f2+f7-f8+f9-f10+f11-f12+f13-f14; m4 = 4*(-f1+f2)+f7-f8+f9-f10+f11-f12+f13-f14; jy = f3-f4+f7-f8-f9+f10+f15-f16+f17-f18; m6 = -4*(f3-f4)+f7-f8-f9+f10+f15-f16+f17-f18; jz = f5-f6+f11-f12-f13+f14+f15-f16-f17+f18; m8 = -4*(f5-f6)+f11-f12-f13+f14+f15-f16-f17+f18; m9 = 2*(f1+f2)-f3-f4-f5-f6+f7+f8+f9+f10+f11+f12+f13+f14-2*(f15+f16+f17+f18); m10 = -4*(f1+f2)+2*(f4+f3+f6+f5)+f8+f7+f10+f9+f12+f11+f14+f13-2*(f16+f15+f18+f17); m11 = f4+f3-f6-f5+f8+f7+f10+f9-f12-f11-f14-f13; m12 = -2*(f4+f3-f6-f5)+f8+f7+f10+f9-f12-f11-f14-f13; m13 = f8+f7-f10-f9; m14 = f16+f15-f18-f17; m15 = f12+f11-f14-f13; m16 = f7-f8+f9-f10-f11+f12-f13+f14; m17 = -f7+f8+f9-f10+f15-f16+f17-f18; m18 = f11-f12-f13+f14-f15+f16+f17-f18; //..........Toelke, Fruediger et. al. 2006............... if (C == 0.0) nx = ny = nz = 1.0; #ifdef STOKES m1 = m1 + rlx_setA*(- 11*rho -alpha*C - m1); m2 = m2 + rlx_setA*(3*rho - m2); m4 = m4 + rlx_setB*((-0.6666666666666666*jx)- m4); m6 = m6 + rlx_setB*((-0.6666666666666666*jy)- m6); m8 = m8 + rlx_setB*((-0.6666666666666666*jz)- m8); m9 = m9 + rlx_setA*( 0.5*alpha*C*(2*nx*nx-ny*ny-nz*nz) - m9); m10 = m10 + rlx_setA*( - m10); m11 = m11 + rlx_setA*( 0.5*alpha*C*(ny*ny-nz*nz)- m11); m12 = m12 + rlx_setA*( - m12); m13 = m13 + rlx_setA*( 0.5*alpha*C*nx*ny - m13); m14 = m14 + rlx_setA*( 0.5*alpha*C*ny*nz - m14); m15 = m15 + rlx_setA*( 0.5*alpha*C*nx*nz - m15); m16 = m16 + rlx_setB*( - m16); m17 = m17 + rlx_setB*( - m17); m18 = m18 + rlx_setB*( - m18); #else m1 = m1 + rlx_setA*((19*(jx*jx+jy*jy+jz*jz)/rho - 11*rho) -alpha*C - m1); m2 = m2 + rlx_setA*((3*rho - 5.5*(jx*jx+jy*jy+jz*jz)/rho)- m2); m4 = m4 + rlx_setB*((-0.6666666666666666*jx)- m4); m6 = m6 + rlx_setB*((-0.6666666666666666*jy)- m6); m8 = m8 + rlx_setB*((-0.6666666666666666*jz)- m8); m9 = m9 + rlx_setA*(((2*jx*jx-jy*jy-jz*jz)/rho) + 0.5*alpha*C*(2*nx*nx-ny*ny-nz*nz) - m9); m10 = m10 + rlx_setA*( - m10); m11 = m11 + rlx_setA*(((jy*jy-jz*jz)/rho) + 0.5*alpha*C*(ny*ny-nz*nz)- m11); m12 = m12 + rlx_setA*( - m12); m13 = m13 + rlx_setA*( (jx*jy/rho) + 0.5*alpha*C*nx*ny - m13); m14 = m14 + rlx_setA*( (jy*jz/rho) + 0.5*alpha*C*ny*nz - m14); m15 = m15 + rlx_setA*( (jx*jz/rho) + 0.5*alpha*C*nx*nz - m15); m16 = m16 + rlx_setB*( - m16); m17 = m17 + rlx_setB*( - m17); m18 = m18 + rlx_setB*( - m18); #endif //.................inverse transformation...................................................... f0 = 0.05263157894736842*rho-0.012531328320802*m1+0.04761904761904762*m2; f1 = 0.05263157894736842*rho-0.004594820384294068*m1-0.01587301587301587*m2 +0.1*(jx-m4)+0.0555555555555555555555555*(m9-m10); f2 = 0.05263157894736842*rho-0.004594820384294068*m1-0.01587301587301587*m2 +0.1*(m4-jx)+0.0555555555555555555555555*(m9-m10); f3 = 0.05263157894736842*rho-0.004594820384294068*m1-0.01587301587301587*m2 +0.1*(jy-m6)+0.02777777777777778*(m10-m9)+0.08333333333333333*(m11-m12); f4 = 0.05263157894736842*rho-0.004594820384294068*m1-0.01587301587301587*m2 +0.1*(m6-jy)+0.02777777777777778*(m10-m9)+0.08333333333333333*(m11-m12); f5 = 0.05263157894736842*rho-0.004594820384294068*m1-0.01587301587301587*m2 +0.1*(jz-m8)+0.02777777777777778*(m10-m9)+0.08333333333333333*(m12-m11); f6 = 0.05263157894736842*rho-0.004594820384294068*m1-0.01587301587301587*m2 +0.1*(m8-jz)+0.02777777777777778*(m10-m9)+0.08333333333333333*(m12-m11); f7 = 0.05263157894736842*rho+0.003341687552213868*m1+0.003968253968253968*m2+0.1*(jx+jy)+0.025*(m4+m6) +0.02777777777777778*m9+0.01388888888888889*m10+0.08333333333333333*m11 +0.04166666666666666*m12+0.25*m13+0.125*(m16-m17); f8 = 0.05263157894736842*rho+0.003341687552213868*m1+0.003968253968253968*m2-0.1*(jx+jy)-0.025*(m4+m6) +0.02777777777777778*m9+0.01388888888888889*m10+0.08333333333333333*m11 +0.04166666666666666*m12+0.25*m13+0.125*(m17-m16); f9 = 0.05263157894736842*rho+0.003341687552213868*m1+0.003968253968253968*m2+0.1*(jx-jy)+0.025*(m4-m6) +0.02777777777777778*m9+0.01388888888888889*m10+0.08333333333333333*m11 +0.04166666666666666*m12-0.25*m13+0.125*(m16+m17); f10 = 0.05263157894736842*rho+0.003341687552213868*m1+0.003968253968253968*m2+0.1*(jy-jx)+0.025*(m6-m4) +0.02777777777777778*m9+0.01388888888888889*m10+0.08333333333333333*m11 +0.04166666666666666*m12-0.25*m13-0.125*(m16+m17); f11 = 0.05263157894736842*rho+0.003341687552213868*m1 +0.003968253968253968*m2+0.1*(jx+jz)+0.025*(m4+m8) +0.02777777777777778*m9+0.01388888888888889*m10-0.08333333333333333*m11 -0.04166666666666666*m12+0.25*m15+0.125*(m18-m16); f12 = 0.05263157894736842*rho+0.003341687552213868*m1 +0.003968253968253968*m2-0.1*(jx+jz)-0.025*(m4+m8) +0.02777777777777778*m9+0.01388888888888889*m10-0.08333333333333333*m11 -0.04166666666666666*m12+0.25*m15+0.125*(m16-m18); f13 = 0.05263157894736842*rho+0.003341687552213868*m1 +0.003968253968253968*m2+0.1*(jx-jz)+0.025*(m4-m8) +0.02777777777777778*m9+0.01388888888888889*m10-0.08333333333333333*m11 -0.04166666666666666*m12-0.25*m15-0.125*(m16+m18); f14 = 0.05263157894736842*rho+0.003341687552213868*m1 +0.003968253968253968*m2+0.1*(jz-jx)+0.025*(m8-m4) +0.02777777777777778*m9+0.01388888888888889*m10-0.08333333333333333*m11 -0.04166666666666666*m12-0.25*m15+0.125*(m16+m18); f15 = 0.05263157894736842*rho+0.003341687552213868*m1 +0.003968253968253968*m2+0.1*(jy+jz)+0.025*(m6+m8) -0.0555555555555555555555555*m9-0.02777777777777778*m10+0.25*m14+0.125*(m17-m18); f16 = 0.05263157894736842*rho+0.003341687552213868*m1 +0.003968253968253968*m2-0.1*(jy+jz)-0.025*(m6+m8) -0.0555555555555555555555555*m9-0.02777777777777778*m10+0.25*m14+0.125*(m18-m17); f17 = 0.05263157894736842*rho+0.003341687552213868*m1 +0.003968253968253968*m2+0.1*(jy-jz)+0.025*(m6-m8) -0.0555555555555555555555555*m9-0.02777777777777778*m10-0.25*m14+0.125*(m17+m18); f18 = 0.05263157894736842*rho+0.003341687552213868*m1 +0.003968253968253968*m2+0.1*(jz-jy)+0.025*(m8-m6) -0.0555555555555555555555555*m9-0.02777777777777778*m10-0.25*m14-0.125*(m17+m18); //....................................................................................................... // incorporate external force f1 += 0.16666666*Fx; f2 -= 0.16666666*Fx; f3 += 0.16666666*Fy; f4 -= 0.16666666*Fy; f5 += 0.16666666*Fz; f6 -= 0.16666666*Fz; f7 += 0.08333333333*(Fx+Fy); f8 -= 0.08333333333*(Fx+Fy); f9 += 0.08333333333*(Fx-Fy); f10 -= 0.08333333333*(Fx-Fy); f11 += 0.08333333333*(Fx+Fz); f12 -= 0.08333333333*(Fx+Fz); f13 += 0.08333333333*(Fx-Fz); f14 -= 0.08333333333*(Fx-Fz); f15 += 0.08333333333*(Fy+Fz); f16 -= 0.08333333333*(Fy+Fz); f17 += 0.08333333333*(Fy-Fz); f18 -= 0.08333333333*(Fy-Fz); //*********** WRITE UPDATED VALUES TO MEMORY ****************** // Write the updated distributions //....EVEN..................................... disteven[n] = f0; disteven[N+n] = f2; disteven[2*N+n] = f4; disteven[3*N+n] = f6; disteven[4*N+n] = f8; disteven[5*N+n] = f10; disteven[6*N+n] = f12; disteven[7*N+n] = f14; disteven[8*N+n] = f16; disteven[9*N+n] = f18; //....ODD...................................... distodd[n] = f1; distodd[N+n] = f3; distodd[2*N+n] = f5; distodd[3*N+n] = f7; distodd[4*N+n] = f9; distodd[5*N+n] = f11; distodd[6*N+n] = f13; distodd[7*N+n] = f15; distodd[8*N+n] = f17; //...Store the Velocity.......................... Velocity[n] = jx; Velocity[N+n] = jy; Velocity[2*N+n] = jz; /* Velocity[3*n] = jx; Velocity[3*n+1] = jy; Velocity[3*n+2] = jz; */ //...Store the Color Gradient.................... // ColorGrad[3*n] = nx*C; // ColorGrad[3*n+1] = ny*C; // ColorGrad[3*n+2] = nz*C; //............................................... //*************************************************************** } // check if n is in the solid } // loop over n } extern "C" void ColorCollideOpt( char *ID, double *disteven, double *distodd, double *phi, double *ColorGrad, double *Velocity, int Nx, int Ny, int Nz, double rlx_setA, double rlx_setB, double alpha, double beta, double Fx, double Fy, double Fz) { int i,j,k,n,nn,N; // distributions double f0,f1,f2,f3,f4,f5,f6,f7,f8,f9; double f10,f11,f12,f13,f14,f15,f16,f17,f18; // non-conserved moments double m1,m2,m4,m6,m8,m9,m10,m11,m12,m13,m14,m15,m16,m17,m18; // additional variables needed for computations double rho,jx,jy,jz,C,nx,ny,nz; N = Nx*Ny*Nz; char id; for (n=0; n 0){ //.......Back out the 3-D indices for node n.............. k = n/(Nx*Ny); j = (n-Nx*Ny*k)/Nx; i = n-Nx*Ny*k-Nx*j; //........................................................................ //........Get 1-D index for this thread.................... // n = S*blockIdx.x*blockDim.x + s*blockDim.x + threadIdx.x; //........................................................................ // COMPUTE THE COLOR GRADIENT //........................................................................ //.................Read Phase Indicator Values............................ //........................................................................ nn = n-1; // neighbor index (get convention) if (i-1<0) nn += Nx; // periodic BC along the x-boundary f1 = phi[nn]; // get neighbor for phi - 1 //........................................................................ nn = n+1; // neighbor index (get convention) if (!(i+10)) delta=0; a1 = na*(0.1111111111111111*(1+4.5*ux))+delta; b1 = nb*(0.1111111111111111*(1+4.5*ux))-delta; a2 = na*(0.1111111111111111*(1-4.5*ux))-delta; b2 = nb*(0.1111111111111111*(1-4.5*ux))+delta; A_odd[n] = a1; A_even[N+n] = a2; B_odd[n] = b1; B_even[N+n] = b2; //............................................... // q = 2 // Cq = {0,1,0} delta = beta*na*nb*nab*0.1111111111111111*ny; if (!(na*nb*nab>0)) delta=0; a1 = na*(0.1111111111111111*(1+4.5*uy))+delta; b1 = nb*(0.1111111111111111*(1+4.5*uy))-delta; a2 = na*(0.1111111111111111*(1-4.5*uy))-delta; b2 = nb*(0.1111111111111111*(1-4.5*uy))+delta; A_odd[N+n] = a1; A_even[2*N+n] = a2; B_odd[N+n] = b1; B_even[2*N+n] = b2; //............................................... // q = 4 // Cq = {0,0,1} delta = beta*na*nb*nab*0.1111111111111111*nz; if (!(na*nb*nab>0)) delta=0; a1 = na*(0.1111111111111111*(1+4.5*uz))+delta; b1 = nb*(0.1111111111111111*(1+4.5*uz))-delta; a2 = na*(0.1111111111111111*(1-4.5*uz))-delta; b2 = nb*(0.1111111111111111*(1-4.5*uz))+delta; A_odd[2*N+n] = a1; A_even[3*N+n] = a2; B_odd[2*N+n] = b1; B_even[3*N+n] = b2; //............................................... /* // Construction and streaming for the components for (idx=0; idx<3; idx++){ //............................................... // Distribution index q = 2*idx; // Associated discrete velocity Cqx = D3Q7[idx][0]; Cqy = D3Q7[idx][1]; Cqz = D3Q7[idx][2]; // Generate the Equilibrium Distribution a1 = na*feq[q]; b1 = nb*feq[q]; a2 = na*feq[q+1]; b2 = nb*feq[q+1]; // Recolor the distributions if (C > 0.0){ sp = nx*double(Cqx)+ny*double(Cqy)+nz*double(Cqz); //if (idx > 2) sp = 0.7071067811865475*sp; //delta = sp*min( min(a1,a2), min(b1,b2) ); delta = na*nb/(na+nb)*0.1111111111111111*sp; //if (a1>0 && b1>0){ a1 += beta*delta; a2 -= beta*delta; b1 -= beta*delta; b2 += beta*delta; } // Save the re-colored distributions A_odd[N*idx+n] = a1; A_even[N*(idx+1)+n] = a2; B_odd[N*idx+n] = b1; B_even[N*(idx+1)+n] = b2; //............................................... } */ } } } //************************************************************************* extern "C" void DensityStreamD3Q7(char *ID, double *Den, double *Copy, double *Phi, double *ColorGrad, double *Velocity, double beta, int Nx, int Ny, int Nz, bool pBC, int S) { char id; int idx; int in,jn,kn,n,nn,N; int q,Cqx,Cqy,Cqz; // int sendLoc; double na,nb; // density values double ux,uy,uz; // flow velocity double nx,ny,nz,C; // color gradient components double a1,a2,b1,b2; double sp,delta; double feq[6]; // equilibrium distributions // Set of Discrete velocities for the D3Q19 Model int D3Q7[3][3]={{1,0,0},{0,1,0},{0,0,1}}; N = Nx*Ny*Nz; for (n=0; n 0 && na+nb > 0.0){ //.......Back out the 3-D indices for node n.............. int k = n/(Nx*Ny); int j = (n-Nx*Ny*k)/Nx; int i = n-Nx*Ny*k-Nx*j; //.....Load the Color gradient......... nx = ColorGrad[n]; ny = ColorGrad[N+n]; nz = ColorGrad[2*N+n]; C = sqrt(nx*nx+ny*ny+nz*nz); nx = nx/C; ny = ny/C; nz = nz/C; //....Load the flow velocity........... ux = Velocity[n]; uy = Velocity[N+n]; uz = Velocity[2*N+n]; //....Instantiate the density distributions // Generate Equilibrium Distributions and stream // Stationary value - distribution 0 // Den[2*n] += 0.3333333333333333*na; // Den[2*n+1] += 0.3333333333333333*nb; Den[2*n] += 0.3333333333333333*na; Den[2*n+1] += 0.3333333333333333*nb; // Non-Stationary equilibrium distributions feq[0] = 0.1111111111111111*(1+3*ux); feq[1] = 0.1111111111111111*(1-3*ux); feq[2] = 0.1111111111111111*(1+3*uy); feq[3] = 0.1111111111111111*(1-3*uy); feq[4] = 0.1111111111111111*(1+3*uz); feq[5] = 0.1111111111111111*(1-3*uz); // Construction and streaming for the components for (idx=0; idx<3; idx++){ // Distribution index q = 2*idx; // Associated discrete velocity Cqx = D3Q7[idx][0]; Cqy = D3Q7[idx][1]; Cqz = D3Q7[idx][2]; // Generate the Equilibrium Distribution a1 = na*feq[q]; b1 = nb*feq[q]; a2 = na*feq[q+1]; b2 = nb*feq[q+1]; // Recolor the distributions if (C > 0.0){ sp = nx*double(Cqx)+ny*double(Cqy)+nz*double(Cqz); //if (idx > 2) sp = 0.7071067811865475*sp; //delta = sp*min( min(a1,a2), min(b1,b2) ); delta = na*nb/(na+nb)*0.1111111111111111*sp; //if (a1>0 && b1>0){ a1 += beta*delta; a2 -= beta*delta; b1 -= beta*delta; b2 += beta*delta; } // .......Get the neighbor node.............. //nn = n + Stride[idx]; in = i+Cqx; jn = j+Cqy; kn = k+Cqz; // Adjust for periodic BC, if necessary // if (in<0) in+= Nx; // if (jn<0) jn+= Ny; // if (kn<0) kn+= Nz; // if (!(in 0 ){ // Get the density value (Streaming already performed) Na = Den[n]; Nb = Den[N+n]; Phi[n] = (Na-Nb)/(Na+Nb); } } //................................................................... } extern "C" void SetPhiSlice_z(double *Phi, double value, int Nx, int Ny, int Nz, int Slice){ int n; for (n=Slice*Nx*Ny; n<(Slice+1)*Nx*Ny; n++){ Phi[n] = value; } } /* // ************************************************************************* extern "C" void InitDenColor( int nblocks, int nthreads, int S, char *ID, double *Den, double *Phi, double das, double dbs, int Nx, int Ny, int Nz) { InitDenColor <<>> (ID, Den, Phi, das, dbs, Nx, Ny, Nz, S); } // ************************************************************************* extern "C" void ComputeColorGradient(int nBlocks, int nthreads, int S, char *ID, double *Phi, double *ColorGrad, int Nx, int Ny, int Nz) { ComputeColorGradient<<>>(ID, Phi, ColorGrad, Nx, Ny, Nz, S); } // ************************************************************************* extern "C" void ColorCollide(int nBlocks, int nthreads, int S, char *ID, double *f_even, double *f_odd, double *ColorGrad, double *Velocity, double rlxA, double rlxB,double alpha, double beta, double Fx, double Fy, double Fz, int Nx, int Ny, int Nz, bool pBC) { ColorCollide<<>>(ID, f_even, f_odd, ColorGrad, Velocity, Nx, Ny, Nz, S, rlxA, rlxB, alpha, beta, Fx, Fy, Fz, pBC); } // ************************************************************************* extern "C" void ColorCollideOpt(int nBlocks, int nthreads, int S, char *ID, double *f_even, double *f_odd, double *Phi, double *ColorGrad, double *Velocity, int Nx, int Ny, int Nz,double rlxA, double rlxB, double alpha, double beta, double Fx, double Fy, double Fz) { ColorCollideOpt<<>>(ID, f_even, f_odd, Phi, ColorGrad, Velocity, Nx, Ny, Nz, S, rlxA, rlxB, alpha, beta, Fx, Fy, Fz); // bool pBC = false; // ColorCollide<<>>(ID, f_even, f_odd, ColorGrad, Velocity, Nx, Ny, Nz, S, // rlxA, rlxB, alpha, beta, Fx, Fy, Fz, pBC); } // ************************************************************************* extern "C" void DensityStreamD3Q7(int nBlocks, int nthreads, int S, char *ID, double *Den, double *Copy, double *Phi, double *ColorGrad, double *Velocity, double beta, int Nx, int Ny, int Nz, bool pBC) { DensityStreamD3Q7<<>>(ID,Den,Copy,Phi,ColorGrad,Velocity,beta,Nx,Ny,Nz,pBC,S); } // ************************************************************************* extern "C" void ComputePhi(int nBlocks, int nthreads, int S, char *ID, double *Phi, double *Copy, double *Den, int N) { ComputePhi<<>>(ID,Phi,Copy,Den,N,S); } // ************************************************************************* extern "C" void ComputePressure(int nBlocks, int nthreads, int S, char *ID, double *disteven, double *distodd, double *Pressure, int Nx, int Ny, int Nz) { ComputePressureD3Q19<<>>(ID,disteven,distodd,Pressure,Nx,Ny,Nz,S); } */