float ngb_treefind(float xyz[3], int desngb, float hguess) { int numngb, iter; float h2max; if(hguess == 0) { hguess = Hmax; } iter = 0; do { iter++; /* printf("hguess= %g\n", hguess); */ numngb = ngb_treefind_variable(xyz, hguess); if(numngb < desngb) { hguess *= 1.26; continue; } if(numngb >= desngb) { qsort(R2list, numngb, sizeof(struct r2data), ngb_compare_key); h2max = R2list[desngb - 1].r2; break; } hguess *= 1.26; } while(1); return sqrt(h2max); }
void setup_nbr_sidm(void) { int i,j,k, numngb; float *r2list; int *ngblist; double tstart,tend; int ntot,ntotleft,npleft,nthis; int nstart,nbuffer,ncount,nchunk; int timelinecounter,startcounter; int *nrecv,*noffset; int level,sendTask,recvTask; int place, nexport; int num_collisionless; MPI_Status status; num_collisionless= NumForceUpdate-NumSphUpdate; MPI_Allreduce(&num_collisionless, &ntot, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD); /* compute maximum size of bunch that may come from this task */ if(num_collisionless==0) nthis=1; else nthis = ((double)(num_collisionless)*(All.BunchSizeSidm-NTask))/ntot + 1; if(num_collisionless>0) { nchunk= num_collisionless/nthis; /* number of bunches needed */ if((num_collisionless%nthis)>0) nchunk+=1; } else nchunk=0; nrecv= malloc(sizeof(int)*NTask); /* list of particle numbers that constituate current bunch */ noffset= malloc(sizeof(int)*NTask); /* offsets of bunches in common list */ ntotleft= ntot; /* particles left for all tasks together */ npleft= num_collisionless; /* particles left for this task */ nstart= IndFirstUpdate; /* first particle for this task */ startcounter=0; while(ntotleft > 0) { if(nthis > npleft) nthis=npleft; for(i=nstart, ncount=0, nexport=0, timelinecounter=startcounter; timelinecounter < NumForceUpdate && ncount<nthis; i=P[i].ForceFlag, timelinecounter++) { if(P[i].Type > 0) { ncount++; for(j=0; j<3; j++) { if(P[i].PosPred[j] < (DomainMin[P[i].Type][j] + P[i].HsmlVelDisp)) break; if(P[i].PosPred[j] > (DomainMax[P[i].Type][j] - P[i].HsmlVelDisp)) break; } if(j != 3) { /* particle lies NOT completely inside. needs to be sent to other processors */ nexport++; P[i].Type |= 8; } } } MPI_Allgather(&nexport, 1, MPI_INT, nrecv, 1, MPI_INT, MPI_COMM_WORLD); MPI_Allreduce(&ncount, &ntot, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD); for(i=nbuffer=0; i<NTask; i++) /* compute length of common list */ nbuffer += nrecv[i]; for(i=1, noffset[0]=0; i<NTask; i++) /* set-up offset-table */ noffset[i]= noffset[i-1] + nrecv[i-1]; /* fill in the own particles at the right place */ for(i=nstart, ncount=0, nexport=0, timelinecounter=startcounter; timelinecounter<NumForceUpdate && ncount<nthis; i=P[i].ForceFlag, timelinecounter++) { if(P[i].Type>0) { if(P[i].Type & 8) { place= noffset[ThisTask] + nexport; nexport++; } else place= nbuffer + (ncount-nexport); for(k=0; k<3; k++) { SidmDataIn[place].Pos[k]= P[i].PosPred[k]; } SidmDataIn[place].Hsml = P[i].HsmlVelDisp; SidmDataIn[place].Type = P[i].Type&7; ncount++; } } /* 4th big communication */ tstart= second(); for(level=1; level<NTask; level++) { sendTask = ThisTask; recvTask = ThisTask ^ level; MPI_Sendrecv(&SidmDataIn[noffset[sendTask]], nrecv[sendTask]*sizeof(struct sidmdata_in), MPI_BYTE, recvTask, TAG_ANY, &SidmDataIn[noffset[recvTask]], nrecv[recvTask]*sizeof(struct sidmdata_in), MPI_BYTE, recvTask, TAG_ANY, MPI_COMM_WORLD, &status); } tend=second(); All.CPU_CommSum+= timediff(tstart, tend); /* ok, all preparations are done. Now density evaluation and velocity stuff */ for(i=0; i<(nbuffer+ncount-nexport); i++) { numngb= ngb_treefind_variable(&SidmDataIn[i].Pos[0], SidmDataIn[i].Hsml, SidmDataIn[i].Type, &ngblist, &r2list); SidmDataResult[i].Ngb= numngb; } /* now communicate contributions, and sum up */ tstart=second(); for(level=1; level<NTask; level++) { sendTask = ThisTask; recvTask = ThisTask ^ level; MPI_Sendrecv(&SidmDataResult[noffset[recvTask]], nrecv[recvTask]*sizeof(struct sidmdata_out), MPI_BYTE, recvTask, TAG_ANY, &SidmDataPartialResult[0], nrecv[sendTask]*sizeof(struct sidmdata_out), MPI_BYTE, recvTask, TAG_ANY, MPI_COMM_WORLD, &status); for(i=0; i<nrecv[ThisTask]; i++) { SidmDataResult[noffset[ThisTask] + i].Ngb += SidmDataPartialResult[i].Ngb; } } tend=second(); All.CPU_CommSum+= timediff(tstart,tend); /* transfer the result to the particles */ for(i=nstart, ncount=0, nexport=0, timelinecounter=startcounter; timelinecounter<NumForceUpdate && ncount<nthis; i=P[i].ForceFlag, timelinecounter++) { if(P[i].Type>0) { if(P[i].Type & 8) { place= noffset[ThisTask] + nexport; nexport++; P[i].Type&= 7; } else place= nbuffer + (ncount-nexport); P[i].NgbVelDisp = SidmDataResult[place].Ngb; ncount++; } } nstart=i; /* continue at this point in the timeline in next iteration */ startcounter=timelinecounter; npleft-= ncount; ntotleft-= ntot; } free(nrecv); free(noffset); }
void sidm(void) { int i,ii,j,k,n, numngb; double h, hinv, hinv3, r, u, wk=0.0; float *r2list; int *ngblist; double tstart,tend; int ntot,ntotleft,npleft,nthis; int nstart,nbuffer,ncount,nchunk; int timelinecounter,startcounter; int *nrecv,*noffset; int level,sendTask,recvTask; int place, nexport; int num_collisionless; MPI_Status status; double dt_h0; double rand, C_Pmax, P_max, Prob; double rv, rvx, rvy, rvz; double rmass; double nx[3]; double dvx, dvy, dvz; double s_a, s_a_inverse; double CrossSectionCo; int nscat[3], nscatTot[3]; #if (CROSS_SECTION_TYPE == 2) double beta, v_dep; #elif (CROSS_SECTION_TYPE == 3) double n_cross_section; // sigma = sigma0*(dv/v_scale)^n double v_scale; #elif (CROSS_SECTION_TYPE == 4) double beta, cosO, sin22, sinO, denom; double rvv[3], nperp[3]; #endif /* file for scattering log */ #ifdef SCATTERLOG FILE *fp; char filename[128]; struct scatlog ScatLog; ScatLog.time= All.Time; sprintf(filename, "sct_%03d.%d", All.SnapshotFileCount, ThisTask); //sprintf(filename, "sct.%d", ThisTask); fp = fopen(filename, "ab"); #endif for(k=0; k<3; k++) nscat[k]= nscatTot[k]= 0; num_collisionless= NumForceUpdate-NumSphUpdate; MPI_Allreduce(&num_collisionless, &ntot, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD); /* compute maximum size of bunch that may come from this task */ if(num_collisionless==0) nthis=1; else nthis = ((double)(num_collisionless)*(All.BunchSizeSidm-NTask))/ntot + 1; if(num_collisionless>0) { nchunk= num_collisionless/nthis; /* number of bunches needed */ if((num_collisionless % nthis)>0) nchunk+=1; } else nchunk= 0; nrecv= malloc(sizeof(int)*NTask); /* list of particle numbers that constituate current bunch */ noffset= malloc(sizeof(int)*NTask); /* offsets of bunches in common list */ ntotleft= ntot; /* particles left for all tasks together */ npleft= num_collisionless; /* particles left for this task */ nstart= IndFirstUpdate; /* first particle for this task */ startcounter=0; while(ntotleft > 0){ if(nthis > npleft) nthis= npleft; for(i=nstart, ncount=0, nexport=0, timelinecounter=startcounter; timelinecounter<NumForceUpdate && ncount<nthis; i=P[i].ForceFlag, timelinecounter++) { if(P[i].Type>0) { ncount++; for(j=0; j<3; j++) { if(P[i].PosPred[j] < (DomainMin[P[i].Type][j] + P[i].HsmlVelDisp)) break; if(P[i].PosPred[j] > (DomainMax[P[i].Type][j] - P[i].HsmlVelDisp)) break; } if(j != 3) { /* particle lies NOT completely inside. needs to be sent to other processors */ nexport++; P[i].Type |= 8; } } } MPI_Allgather(&nexport, 1, MPI_INT, nrecv, 1, MPI_INT, MPI_COMM_WORLD); MPI_Allreduce(&ncount, &ntot, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD); for(i=nbuffer=0; i<NTask; i++) /* compute length of common list */ nbuffer += nrecv[i]; for(i=1, noffset[0]=0; i<NTask; i++) /* set-up offset-table */ noffset[i]= noffset[i-1] + nrecv[i-1]; /* fill in the own particles at the right place */ for(i=nstart, ncount=0, nexport=0, timelinecounter=startcounter; timelinecounter < NumForceUpdate && ncount < nthis; i=P[i].ForceFlag, timelinecounter++) { if(P[i].Type > 0) { if(P[i].Type & 8) { place= noffset[ThisTask] + nexport; nexport++; } else place= nbuffer + (ncount-nexport); for(k=0; k<3; k++) { SidmDataIn[place].Pos[k]= P[i].PosPred[k]; SidmDataIn[place].Vel[k]= P[i].Vel[k]; } SidmDataIn[place].Hsml= P[i].HsmlVelDisp; SidmDataIn[place].Type= P[i].Type & 7; SidmDataIn[place].Mass= P[i].Mass; if(P[i].dVel[0] != 0.0f) SidmDataIn[place].ID= 0; // Oct 9, 2005 else SidmDataIn[place].ID= P[i].ID; SidmDataIn[place].dt= 2*(All.Time - P[i].CurrentTime); SidmDataConfirm[place].Task= -1; ncount++; } } /* 1st big communication */ tstart= second(); for(level=1; level<NTask; level++) { sendTask = ThisTask; recvTask = ThisTask ^ level; MPI_Sendrecv(&SidmDataIn[noffset[sendTask]], nrecv[sendTask]*sizeof(struct sidmdata_in), MPI_BYTE, recvTask, TAG_ANY, &SidmDataIn[noffset[recvTask]], nrecv[recvTask]*sizeof(struct sidmdata_in), MPI_BYTE, recvTask, TAG_ANY, MPI_COMM_WORLD, &status); } tend=second(); All.CPU_CommSum += timediff(tstart, tend); /* ok, all preparations are done. */ /* SIDM: MAIN PART */ /* Common Factors */ /* s_a= a^{3/2} H ==> s_a_inverse da = a^{-1/2} dt ==> s_a_inverse da v_internal= a^{-1} v_phys dt = dx_comoving */ if(All.ComovingIntegrationOn) { s_a= All.Hubble*sqrt(All.Omega0 + All.Time* (1-All.Omega0-All.OmegaLambda)+pow(All.Time,3)*All.OmegaLambda); s_a_inverse= 1/s_a; #if (CROSS_SECTION_TYPE == 0) CrossSectionCo = All.CrossSectionInternal/pow(All.Time,2); C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)*2*vmax* CrossSectionCo; #elif (CROSS_SECTION_TYPE == 1) CrossSectionCo = All.CrossSectionInternal/pow(All.Time,2.5); // a^(1/2) factor considers that in s_a_inverse C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)* CrossSectionCo; #elif (CROSS_SECTION_TYPE == 2) CrossSectionCo = All.CrossSectionInternal/pow(All.Time,2); vc= All.YukawaVelocity/sqrt(All.Time); /* In internal unit */ if(2.0*vmax < vc/sqrt(3.0)) { beta= 2.0*vmax/vc; v_dep= 1.0/(1.0 + beta*beta); C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)* 2.0*vmax*v_dep*v_dep* CrossSectionCo; } else { C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)* (3.0*sqrt(3.0)/16.0)*vc* CrossSectionCo; } #elif (CROSS_SECTION_TYPE == 3) CrossSectionCo = All.CrossSectionInternal/pow(All.Time, 2); n_cross_section= All.CrossSectionPowLaw; v_scale= All.CrossSectionVelScale; C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)*2*v_scale* CrossSectionCo; #elif (CROSS_SECTION_TYPE == 4) vc= All.YukawaVelocity/sqrt(All.Time); CrossSectionCo = All.CrossSectionInternal/pow(All.Time,2); C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)*2*vmax* CrossSectionCo; #endif } else{ s_a_inverse= 1.0; #if (CROSS_SECTION_TYPE == 0) C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)*2*vmax* All.CrossSectionInternal; #elif (CROSS_SECTION_TYPE == 1) C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)* All.CrossSectionInternal; #elif (CROSS_SECTION_TYPE == 2) vc= All.YukawaVelocity; if(2.0*vmax < vc/sqrt(3.0)) { beta= 2.0*vmax/vc; v_dep= 1.0/(1.0 + beta*beta); C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)* 2.0*vmax*v_dep*v_dep* All.CrossSectionInternal; } else { C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)* (3.0*sqrt(3.0)/16.0)*vc* All.CrossSectionInternal; } #elif (CROSS_SECTION_TYPE == 3) n_cross_section= All.CrossSectionPowLaw; v_scale= All.CrossSectionVelScale; C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)*2*v_scale* All.CrossSectionInternal; #elif (CROSS_SECTION_TYPE == 4) vc= All.YukawaVelocity; C_Pmax= SAFEFACTOR * BALLINVERSE*(All.DesNumNgb+All.MaxNumNgbDeviation)*2*vmax* All.CrossSectionInternal; #endif CrossSectionCo = All.CrossSectionInternal; } /* For each active particles */ for(i=0; i<(nbuffer+ncount-nexport); i++) { numngb= ngb_treefind_variable(&SidmDataIn[i].Pos[0], SidmDataIn[i].Hsml, SidmDataIn[i].Type, &ngblist, &r2list); SidmDataResult[i].Ngb= numngb; dt_h0= SidmDataIn[i].dt * s_a_inverse; h = SCATKERNELFACTOR * SidmDataIn[i].Hsml; // Oct 31, 2005 hinv = 1.0/h; hinv3 = hinv*hinv*hinv; /* Initialize Variables */ SidmDataResult[i].dv[0] = 0.0; SidmDataResult[i].dv[1] = 0.0; SidmDataResult[i].dv[2] = 0.0; SidmTarget[i] = 0; /* for security */ P_max = C_Pmax*SidmDataIn[i].Mass*hinv3*dt_h0; // This P_max is not correct if Masses are different. rand = ran2(&iseed); if(P_max < rand) continue; else if(SidmDataIn[i].ID == 0) continue; // Already Scattered nscat[0]++; /* Passed first approximation */ Prob= 0.0; /* for each neighbors */ for(n=0; n<numngb; n++) { j = ngblist[n]+1; r = sqrt(r2list[n]); /* There was a BUG! */ if(P[j].dVel[0] != 0.0f) continue; // No double scattering. Oct 9, 2005 if(r < h) { u = r*hinv; ii = (int)(u*KERNEL_TABLE); wk =hinv3*( Kernel[ii] + (Kernel[ii+1]-Kernel[ii])*(u-KernelRad[ii])*KERNEL_TABLE); } /* Relative Velocities */ rvx = SidmDataIn[i].Vel[0] - P[j].Vel[0]; rvy = SidmDataIn[i].Vel[1] - P[j].Vel[1]; rvz = SidmDataIn[i].Vel[2] - P[j].Vel[2]; rv = sqrt(rvx*rvx + rvy*rvy + rvz*rvz); #if (CROSS_SECTION_TYPE == 0) Prob += 0.5*P[j].Mass*wk*rv*CrossSectionCo*dt_h0; #elif (CROSS_SECTION_TYPE == 1) Prob += 0.5*P[j].Mass*wk* CrossSectionCo*dt_h0; #elif (CROSS_SECTION_TYPE == 2) beta= rv/vc; v_dep= 1.0/(1.0+beta*beta); Prob += 0.5*P[j].Mass*wk*rv*v_dep*v_dep*CrossSectionCo*dt_h0; #elif (CROSS_SECTION_TYPE == 3) Prob += 0.5*P[j].Mass*wk*rv*pow(rv/v_scale, n_cross_section)*CrossSectionCo*dt_h0; #elif (CROSS_SECTION_TYPE == 4) Prob += 0.5*P[j].Mass*wk*rv*CrossSectionCo*dt_h0; #endif if(Prob < rand) continue; SidmTarget[i] = j; rmass = P[j].Mass / (SidmDataIn[i].Mass + P[j].Mass); #if (CROSS_SECTION_TYPE == 4) // Furture selection of angle beta= rv/vc; rand = ran2(&iseed); cosO = 2.0*ran2(&iseed) - 1.0; // cos(theta) of scatter direction sin22 = 0.5*(1.0 - cosO); // sin^2(theta/2) denom = 1.0 + beta*beta*sin22; if(rand >= 1/(denom*denom) || rv == 0.0) continue; if(rv == 0.0) { // debug //printf("Error: rv == 0. %e %e\n", rv, 0.5*P[j].Mass*wk*rv*CrossSectionCo*dt_h0); endrun(4006); } rvv[0]= rvx; rvv[1]= rvy; rvv[2]= rvz; random_direction(nx); perp(rvv, nx, nperp); //if(nperp[0] != nperp[0] || nperp[1] != nperp[1] || nperp[2] != nperp[2]) // endrun(4002); //assert(1.0 - cosO*cosO >= 0.0); sinO = 1.0 - cosO*cosO > 0.0 ? sqrt(1.0 - cosO*cosO) : 0.0; //if(sinO != sinO) endrun(4001); dvx= rmass*(-rvx + cosO*rvv[0] + sinO*rv*nperp[0]); dvy= rmass*(-rvy + cosO*rvv[1] + sinO*rv*nperp[1]); dvz= rmass*(-rvz + cosO*rvv[2] + sinO*rv*nperp[2]); if(dvx != dvx || dvy != dvy || dvz != dvz) { endrun(4000); } //debug /* rvv[0]= cosO*rvv[0] + sinO*rv*nperp[0]; rvv[1]= cosO*rvv[1] + sinO*rv*nperp[1]; rvv[2]= cosO*rvv[2] + sinO*rv*nperp[2]; rvv[0]= SidmDataIn[i].Vel[0] + dvx; rvv[1]= SidmDataIn[i].Vel[1] + dvy; rvv[2]= SidmDataIn[i].Vel[2] + dvz; fprintf(stderr, "%e %e\n", rv, norm(rvv)); */ #endif /* Scatter. Calc velocity change */ // Better way of producing spherical random numbers // Sep 19, 2005 #if (CROSS_SECTION_TYPE < 4) // isotropic scattering random_direction(nx); dvx= rmass*(-rvx+rv*nx[0]); dvy= rmass*(-rvy+rv*nx[1]); dvz= rmass*(-rvz+rv*nx[2]); #endif SidmDataResult[i].dv[0]= dvx; SidmDataResult[i].dv[1]= dvy; SidmDataResult[i].dv[2]= dvz; SidmDataConfirm[i].Task= ThisTask; break; } } /* 2nd communicate contributions, and sum up */ tstart=second(); for(level=1; level<NTask; level++) { sendTask= ThisTask; recvTask= ThisTask ^ level; MPI_Sendrecv(&SidmDataResult[noffset[recvTask]], nrecv[recvTask]*sizeof(struct sidmdata_out), MPI_BYTE, recvTask, TAG_ANY, &SidmDataPartialResult[0], nrecv[sendTask]*sizeof(struct sidmdata_out), MPI_BYTE, recvTask, TAG_ANY, MPI_COMM_WORLD, &status); for(i=0; i<nrecv[ThisTask]; i++) { SidmDataResult[noffset[ThisTask] + i].Ngb += SidmDataPartialResult[i].Ngb; // Keep First collision only if(SidmDataResult[noffset[ThisTask] + i].dv[0] == 0.0 && SidmDataPartialResult[i].dv[0] != 0.0 ) { for(k=0; k<3; k++) { SidmDataResult[noffset[ThisTask] + i].dv[k] = SidmDataPartialResult[i].dv[k]; } SidmDataConfirm[noffset[ThisTask] + i].Task = recvTask; } } } tend=second(); All.CPU_CommSum += timediff(tstart,tend); /* transfer the result to the particles */ for(i=nstart, ncount=0, nexport=0, timelinecounter=startcounter; timelinecounter<NumForceUpdate && ncount<nthis; i=P[i].ForceFlag, timelinecounter++) if(P[i].Type > 0) { if(P[i].Type & 8) { place= noffset[ThisTask] + nexport; nexport++; P[i].Type&= 7; } else place= nbuffer + (ncount-nexport); P[i].NgbVelDisp = SidmDataResult[place].Ngb; // Check # of neighbours if(P[i].NgbVelDisp < (All.DesNumNgb-All.MaxNumNgbDeviation) || (P[i].NgbVelDisp > (All.DesNumNgb+All.MaxNumNgbDeviation))) { // Scattering is invalid // Redo sidm in ensure_neighbours() function SidmDataConfirm[place].Task= -1; if(SidmDataResult[place].dv[0] != 0.0) nscat[2]++; // rejected; } else if(SidmDataResult[place].dv[0] != 0.0) { // OK. Update Velocity. /* Oct 09, 2005 P[i].Vel[0] += SidmDataResult[place].dv[0]; P[i].Vel[1] += SidmDataResult[place].dv[1]; P[i].Vel[2] += SidmDataResult[place].dv[2]; */ P[i].dVel[0] = SidmDataResult[place].dv[0]; P[i].dVel[1] = SidmDataResult[place].dv[1]; P[i].dVel[2] = SidmDataResult[place].dv[2]; nscat[1]++; // scattered; } else { // Feb 15,2006 Double scattering rejection (very rare) SidmDataConfirm[place].Task= -1; } ncount++; } /* Scattering Confirmation */ /* Communication; 3rd time */ tstart= second(); for(level=1; level<NTask; level++){ sendTask = ThisTask; recvTask = ThisTask ^ level; MPI_Sendrecv(&SidmDataConfirm[noffset[sendTask]], nrecv[sendTask]*sizeof(struct sidmdata_confirm), MPI_BYTE, recvTask, TAG_ANY, &SidmDataConfirm[noffset[recvTask]], nrecv[recvTask]*sizeof(struct sidmdata_confirm), MPI_BYTE, recvTask, TAG_ANY, MPI_COMM_WORLD, &status); } tend= second(); All.CPU_CommSum += timediff(tstart, tend); /* Update the targets of active particles if confirmed */ for(j=0; j<(nbuffer+ncount-nexport); j++) { if(SidmDataConfirm[j].Task == ThisTask) { if(SidmTarget[j] <= 0 || SidmTarget[j] > NumPart) { // for debug only printf("SIDM ERROR: SidmTarget is wrong\n"); endrun(0); } // Oct 9,2005 Vel[] => dVel[] -=-> =- P[SidmTarget[j]].dVel[0] = -SidmDataResult[j].dv[0]; P[SidmTarget[j]].dVel[1] = -SidmDataResult[j].dv[1]; P[SidmTarget[j]].dVel[2] = -SidmDataResult[j].dv[2]; #ifdef SCATTERLOG /* Scatter Log */ ScatLog.id1= SidmDataIn[j].ID; ScatLog.id2= P[SidmTarget[j]].ID; ScatLog.Hsml1= SidmDataIn[j].Hsml; ScatLog.Hsml2= P[SidmTarget[j]].HsmlVelDisp; ScatLog.x1[0]= SidmDataIn[j].Pos[0]; ScatLog.x1[1]= SidmDataIn[j].Pos[1]; ScatLog.x1[2]= SidmDataIn[j].Pos[2]; ScatLog.x2[0]= P[SidmTarget[j]].PosPred[0]; ScatLog.x2[1]= P[SidmTarget[j]].PosPred[1]; ScatLog.x2[2]= P[SidmTarget[j]].PosPred[2]; ScatLog.v1[0]= SidmDataIn[j].Vel[0]; ScatLog.v1[1]= SidmDataIn[j].Vel[1]; ScatLog.v1[2]= SidmDataIn[j].Vel[2]; ScatLog.v2[0]= P[SidmTarget[j]].Vel[0]; ScatLog.v2[1]= P[SidmTarget[j]].Vel[1]; ScatLog.v2[2]= P[SidmTarget[j]].Vel[2]; ScatLog.dv[0]= SidmDataResult[j].dv[0]; ScatLog.dv[1]= SidmDataResult[j].dv[1]; ScatLog.dv[2]= SidmDataResult[j].dv[2]; fwrite(&ScatLog, sizeof(ScatLog), 1, fp); #endif } } nstart= i; /* continue at this point in the timeline in next iteration */ startcounter=timelinecounter; npleft-= ncount; ntotleft-= ntot; } // Summerize the scattering statistics MPI_Reduce(nscat, nscatTot, 3, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD); #ifdef FINDNBRLOG if(ThisTask == 0){ //printf("SIDM SCT tot= %d 1st= %d scattered= %d rejected= %d\n", fprintf(stdout, "SCT %d %d %d %d\n", ntot, nscatTot[0], nscatTot[1], nscatTot[2]); } #endif #ifdef SCATTERLOG fclose(fp); #endif free(nrecv); free(noffset); }
/*! This function represents the core of the SPH density computation. The * target particle may either be local, or reside in the communication * buffer. */ void compute_smoothed_evaluate(int target, int mode) { int j, n; int startnode, numngb_inbox; double h, h2, hinv, hinv3, hinv4; double wk, dwk; double dx, dy, dz, r, r2, u, mass_j; FLOAT *pos; #ifdef CONDUCTION double smoothentr; #ifdef CONDUCTION_SATURATION double gradentr[3]; #endif #endif /* CONDUCTION */ #ifdef CR_DIFFUSION double rCR_SmoothE0 = 0.0; double rCR_Smoothn0 = 0.0; #endif /* CR_DIFFUSION */ #ifdef SMOOTH_PHI double SmoothPhi = 0.0; #endif /* SMOOTH_PHI */ #ifdef VOLUME_CORRECTION double graddensity[3],density,normdensity,dV; #endif /* VOLUME_CORRECTION */ #ifdef CONDUCTION smoothentr = 0; #ifdef CONDUCTION_SATURATION gradentr[0] = gradentr[1] = gradentr[2] = 0; #endif #endif /* CONDUCTION */ #ifdef VOLUME_CORRECTION graddensity[0] = graddensity[1] = graddensity[2] = normdensity = 0; #endif /* VOLUME_CORRECTION */ if(mode == 0) { pos = P[target].Pos; h = PPP[target].Hsml; #ifdef VOLUME_CORRECTION density = SphP[target].a2.Density; #endif } else { pos = DensDataGet[target].Pos; h = DensDataGet[target].Hsml; #ifdef VOLUME_CORRECTION density = DensDataGet[target].Density; #endif } #ifdef VOLUME_CORRECTION dV = P[target].Mass / density; #endif h2 = h * h; hinv = 1.0 / h; #ifndef TWODIMS hinv3 = hinv * hinv * hinv; #else hinv3 = hinv * hinv / boxSize_Z; #endif hinv4 = hinv3 * hinv; startnode = All.MaxPart; do { numngb_inbox = ngb_treefind_variable(&pos[0], h, &startnode); for(n = 0; n < numngb_inbox; n++) { j = Ngblist[n]; dx = pos[0] - P[j].Pos[0]; dy = pos[1] - P[j].Pos[1]; dz = pos[2] - P[j].Pos[2]; #ifdef PERIODIC /* now find the closest image in the given box size */ if(dx > boxHalf_X) dx -= boxSize_X; if(dx < -boxHalf_X) dx += boxSize_X; if(dy > boxHalf_Y) dy -= boxSize_Y; if(dy < -boxHalf_Y) dy += boxSize_Y; if(dz > boxHalf_Z) dz -= boxSize_Z; if(dz < -boxHalf_Z) dz += boxSize_Z; #endif r2 = dx * dx + dy * dy + dz * dz; if(r2 < h2) { r = sqrt(r2); u = r * hinv; if(u < 0.5) { wk = hinv3 * (KERNEL_COEFF_1 + KERNEL_COEFF_2 * (u - 1) * u * u); dwk = hinv4 * u * (KERNEL_COEFF_3 * u - KERNEL_COEFF_4); } else { wk = hinv3 * KERNEL_COEFF_5 * (1.0 - u) * (1.0 - u) * (1.0 - u); dwk = hinv4 * KERNEL_COEFF_6 * (1.0 - u) * (1.0 - u); } mass_j = P[j].Mass; #ifdef CONDUCTION smoothentr += mass_j * wk * pow(SphP[j].a2.Density, GAMMA_MINUS1) * SphP[j].Entropy; #ifdef CONDUCTION_SATURATION if(r > 0) { gradentr[0] += mass_j * dwk * dx / r * SphP[j].Entropy * pow(SphP[j].a2.Density, GAMMA_MINUS1); gradentr[1] += mass_j * dwk * dy / r * SphP[j].Entropy * pow(SphP[j].a2.Density, GAMMA_MINUS1); gradentr[2] += mass_j * dwk * dz / r * SphP[j].Entropy * pow(SphP[j].a2.Density, GAMMA_MINUS1); } #endif #endif /* CONDUCTION */ #ifdef CR_DIFFUSION rCR_SmoothE0 += mass_j * wk * SphP[j].CR_E0; rCR_Smoothn0 += mass_j * wk * SphP[j].CR_n0; #endif /* CR_DIFFUSION */ #ifdef SMOOTH_PHI SmoothPhi += mass_j * wk * SphP[j].PhiPred; #endif /* SMOOTH_PHI */ #ifdef VOLUME_CORRECTION if(r > 0) { graddensity[0] += mass_j * dwk * dx / r * (density - SphP[j].a2.Density); graddensity[1] += mass_j * dwk * dy / r * (density - SphP[j].a2.Density); graddensity[2] += mass_j * dwk * dz / r * (density - SphP[j].a2.Density); } #ifdef VOLUME_BIWEIGHT wk = NORM_COEFF_2 * (1 - u * u) * (1 - u * u); #endif #ifdef VOLUME_QUADRIC if(u < 0.2) wk = NORM_COEFF_4 * ( pow(5 * u + 5,4) - 5 * pow(5 * u + 3,4) + 10 * pow(5 * u + 1,4)); else if(u < 0.6) wk = NORM_COEFF_4 * ( pow(5 - 5 * u,4) - 5 * pow(3 - 5 * u,4)); else wk = NORM_COEFF_4 * ( pow(5 - 5 * u,4)); #endif #ifdef VOLUME_QUINTIC if(u < 1./3.) wk = NORM_COEFF_5 * ( pow(3 - 3 * u,5) - 6 * pow(2 - 3 * u,5) + 15 * pow(1 - 3 * u,5)); else if(u < 2./3) wk = NORM_COEFF_5 * ( pow(3 - 3 * u,5) - 6 * pow(2 - 3 * u,5)); else wk = NORM_COEFF_5 * ( pow(3 - 3 * u,5)); #endif #if defined(VOLUME_BIWEIGHT) || defined(VOLUME_QUADRIC) || defined(VOLUME_QUINTIC) normdensity += mass_j * hinv3 * pow(density/SphP[j].a2.Density,1.75) * wk; #else normdensity += mass_j * pow(density/SphP[j].a2.Density,1.75) * wk; #endif #endif /* VOLUME_CORRECTION */ } } } while(startnode >= 0); if(mode == 0) { #ifdef CONDUCTION SphP[target].SmoothedEntr = smoothentr; #ifdef CONDUCTION_SATURATION SphP[target].GradEntr[0] = gradentr[0]; SphP[target].GradEntr[1] = gradentr[1]; SphP[target].GradEntr[2] = gradentr[2]; #endif #endif /* CONDUCTION */ #ifdef CR_DIFFUSION SphP[target].CR_SmoothE0 = rCR_SmoothE0; SphP[target].CR_Smoothn0 = rCR_Smoothn0; #endif /* CR_DIFFUSION */ #ifdef SMOOTH_PHI SphP[target].SmoothPhi = SmoothPhi; #endif /* SMOOTH_PHI */ #ifdef VOLUME_CORRECTION SphP[target].GradDensity[0] = graddensity[0]; SphP[target].GradDensity[1] = graddensity[1]; SphP[target].GradDensity[2] = graddensity[2]; SphP[target].NormDensity = normdensity; #endif /* VOLUME_CORRECTION */ } else { #ifdef CONDUCTION DensDataResult[target].SmoothedEntr = smoothentr; #ifdef CONDUCTION_SATURATION DensDataResult[target].GradEntr[0] = gradentr[0]; DensDataResult[target].GradEntr[1] = gradentr[1]; DensDataResult[target].GradEntr[2] = gradentr[2]; #endif #endif /* CONDUCTION */ #ifdef CR_DIFFUSION DensDataResult[target].CR_SmoothE0 = rCR_SmoothE0; DensDataResult[target].CR_Smoothn0 = rCR_Smoothn0; #endif /* CR_DIFFUSION */ #ifdef SMOOTH_PHI DensDataResult[target].SmoothPhi = SmoothPhi; #endif /* SMOOTH_PHI */ #ifdef VOLUME_CORRECTION DensDataResult[target].GradDensity[0] = graddensity[0]; DensDataResult[target].GradDensity[1] = graddensity[1]; DensDataResult[target].GradDensity[2] = graddensity[2]; DensDataResult[target].NormDensity = normdensity; #endif /* VOLUME_CORRECTION */ } }
/*! This function represents the core of the SPH density computation. The * target particle may either be local, or reside in the communication * buffer. */ int star_density_evaluate(int target, int mode, int *nexport, int *nsend_local) { int j, n, numngb; int startnode, listindex = 0; double h, hinv, h2, mass_j, weight, hinv3; double wk, mass, density, lum; double dx, dy, dz, r, r2, u, a3inv; MyDouble *pos; #ifdef PERIODIC double boxsize, boxhalf; boxsize = All.BoxSize; boxhalf = 0.5 * All.BoxSize; #endif if(All.ComovingIntegrationOn) a3inv = 1.0 / (All.Time * All.Time * All.Time); else a3inv = 1.0; if(mode == 0) { pos = P[target].Pos; h = PPP[target].Hsml; density = P[target].DensAroundStar; mass = P[target].Mass; } else { pos = StarDataGet[target].Pos; h = StarDataGet[target].Hsml; density = StarDataGet[target].Density; mass = StarDataGet[target].Mass; } h2 = h * h; hinv = 1.0 / h; hinv3 = hinv * hinv * hinv; lum = mass * All.IonizingLumPerSolarMass * All.UnitTime_in_s / All.HubbleParam; if(mode == 0) { startnode = All.MaxPart; /* root node */ } else { startnode = StarDataGet[target].NodeList[0]; startnode = Nodes[startnode].u.d.nextnode; /* open it */ } while(startnode >= 0) { while(startnode >= 0) { numngb = ngb_treefind_variable(pos, h, target, &startnode, mode, nexport, nsend_local); if(numngb < 0) return -1; for(n = 0; n < numngb; n++) { j = Ngblist[n]; dx = pos[0] - P[j].Pos[0]; dy = pos[1] - P[j].Pos[1]; dz = pos[2] - P[j].Pos[2]; #ifdef PERIODIC if(dx > boxHalf_X) dx -= boxSize_X; if(dx < -boxHalf_X) dx += boxSize_X; if(dy > boxHalf_Y) dy -= boxSize_Y; if(dy < -boxHalf_Y) dy += boxSize_Y; if(dz > boxHalf_Z) dz -= boxSize_Z; if(dz < -boxHalf_Z) dz += boxSize_Z; #endif r2 = dx * dx + dy * dy + dz * dz; r = sqrt(r2); if(r2 < h2) { u = r * hinv; if(u < 0.5) wk = hinv3 * (KERNEL_COEFF_1 + KERNEL_COEFF_2 * (u - 1) * u * u); else wk = hinv3 * KERNEL_COEFF_5 * (1.0 - u) * (1.0 - u) * (1.0 - u); } else wk = 0; mass_j = P[j].Mass; weight = mass_j * wk / density; SphP[j].Je += lum * weight / mass_j * (PROTONMASS / SOLAR_MASS); } } if(mode == 1) { listindex++; if(listindex < NODELISTLENGTH) { startnode = StarDataGet[target].NodeList[listindex]; if(startnode >= 0) startnode = Nodes[startnode].u.d.nextnode; /* open it */ } } } return 0; }
void bsmooth_evaluate(int target, int mode) { int j, n; int startnode, numngb_inbox; double h, h2, hinv, hinv3; double wk; double dx, dy, dz, r, r2, u, mass_j; FLOAT pos[3]; FLOAT BSmooth[3], DensityNorm, mwk_d; BSmooth[0] = BSmooth[1] = BSmooth[2] = DensityNorm = 0; if(mode == 0) { pos[0] = P[target].Pos[0]; pos[1] = P[target].Pos[1]; pos[2] = P[target].Pos[2]; h = PPP[target].Hsml; } else { pos[0] = DensDataGet[target].Pos[0]; pos[1] = DensDataGet[target].Pos[1]; pos[2] = DensDataGet[target].Pos[2]; h = DensDataGet[target].Hsml; } h2 = h * h; hinv = 1.0 / h; #ifndef TWODIMS hinv3 = hinv * hinv * hinv; #else hinv3 = hinv * hinv / boxSize_Z; #endif startnode = All.MaxPart; do { numngb_inbox = ngb_treefind_variable(&pos[0], h, &startnode); for(n = 0; n < numngb_inbox; n++) { j = Ngblist[n]; /* MHD sould be solved as well for wind particles ! */ dx = pos[0] - P[j].Pos[0]; dy = pos[1] - P[j].Pos[1]; dz = pos[2] - P[j].Pos[2]; #ifdef PERIODIC /* now find the closest image in the given box size */ if(dx > boxHalf_X) dx -= boxSize_X; if(dx < -boxHalf_X) dx += boxSize_X; if(dy > boxHalf_Y) dy -= boxSize_Y; if(dy < -boxHalf_Y) dy += boxSize_Y; if(dz > boxHalf_Z) dz -= boxSize_Z; if(dz < -boxHalf_Z) dz += boxSize_Z; #endif r2 = dx * dx + dy * dy + dz * dz; if(r2 < h2) { r = sqrt(r2); u = r * hinv; if(u < 0.5) { wk = hinv3 * (KERNEL_COEFF_1 + KERNEL_COEFF_2 * (u - 1) * u * u); } else { wk = hinv3 * KERNEL_COEFF_5 * (1.0 - u) * (1.0 - u) * (1.0 - u); } mass_j = P[j].Mass; mwk_d = mass_j * wk / SphP[j].a2.Density; BSmooth[0] += SphP[j].BPred[0] * mwk_d; BSmooth[1] += SphP[j].BPred[1] * mwk_d; BSmooth[2] += SphP[j].BPred[2] * mwk_d; DensityNorm += mwk_d; /* mass_j * hinv3 * NORM_COEFF_2 * (1 - u * u) * (1 - u * u); */ } } } while(startnode >= 0); if(mode == 0) { SphP[target].BSmooth[0] = BSmooth[0]; SphP[target].BSmooth[1] = BSmooth[1]; SphP[target].BSmooth[2] = BSmooth[2]; SphP[target].DensityNorm = DensityNorm; } else { DensDataResult[target].BSmooth[0] = BSmooth[0]; DensDataResult[target].BSmooth[1] = BSmooth[1]; DensDataResult[target].BSmooth[2] = BSmooth[2]; DensDataResult[target].DensityNorm = DensityNorm; } }
void density_evaluate(int target, int mode) #endif { int ii, j, n, startnode, numngb, numngb_inbox; double h, h2, fac, hinv, hinv3, hinv4; double rho, divv, wk, dwk, divb; double dx, dy, dz, r, r2, u, mass_j; double dvx, dvy, dvz, rotv[3], gsci[3], vrel[3], sci, dBx, dBy, dBz; double weighted_numngb, dhsmlrho; FLOAT *pos, *vel, *bfd; #ifdef METALS_TG double a, hubble_param, dens, D_i, D_j, D_avg, dwk_r, K_ab; double sigma = 0.0; double const_A = 0.0; double const_B = 0.0; if(All.ComovingIntegrationOn) { a = All.Time; hubble_param = All.HubbleParam; } else a = hubble_param = 1.0; #endif if(mode == 0) { pos = P[target].Pos; vel = SphP[target].VelPred; bfd = SphP[target].bfield; sci = SphP[target].Sci; h = SphP[target].Hsml; #ifdef METALS_TG D_i = 2.0*h*SphP[target].Sigma; if(D_i < 1.0e-20) D_i = 1.0e-20; #endif } else { pos = DensDataGet[target].Pos; vel = DensDataGet[target].Vel; bfd = DensDataGet[target].bfield; sci = DensDataGet[target].Sci; h = DensDataGet[target].Hsml; #ifdef METALS_TG D_i = 2.0*h*DensDataGet[target].Sigma; if(D_i < 1.0e-20) D_i = 1.0e-20; #endif } h2 = h * h; hinv = 1.0 / h; #ifndef TWODIMS hinv3 = hinv * hinv * hinv; #else hinv3 = hinv * hinv / boxSize_Z; #endif hinv4 = hinv3 * hinv; rho = divv = rotv[0] = rotv[1] = rotv[2] = vrel[0] = vrel[1] = vrel[2] = 0; gsci[0] = gsci[1] = gsci[2] = sci = divb = 0; weighted_numngb = 0; dhsmlrho = 0; startnode = All.MaxPart; numngb = 0; do { numngb_inbox = ngb_treefind_variable(&pos[0], h, &startnode); for(n = 0; n < numngb_inbox; n++) { j = Ngblist[n]; dx = pos[0] - P[j].Pos[0]; dy = pos[1] - P[j].Pos[1]; dz = pos[2] - P[j].Pos[2]; #ifdef PERIODIC /* now find the closest image in the given box size */ if(dx > boxHalf_X) dx -= boxSize_X; if(dx < -boxHalf_X) dx += boxSize_X; if(dy > boxHalf_Y) dy -= boxSize_Y; if(dy < -boxHalf_Y) dy += boxSize_Y; if(dz > boxHalf_Z) dz -= boxSize_Z; if(dz < -boxHalf_Z) dz += boxSize_Z; #endif r2 = dx * dx + dy * dy + dz * dz; r = sqrt(r2); if(r2 < h2) { numngb++; u = r * hinv; if(u < 0.5) { wk = hinv3 * (KERNEL_COEFF_1 + KERNEL_COEFF_2 * (u - 1) * u * u); dwk = hinv4 * u * (KERNEL_COEFF_3 * u - KERNEL_COEFF_4); } else { wk = hinv3 * KERNEL_COEFF_5 * (1.0 - u) * (1.0 - u) * (1.0 - u); dwk = hinv4 * KERNEL_COEFF_6 * (1.0 - u) * (1.0 - u); } mass_j = P[j].Mass; if(P[j].ID < 0 || SphP[j].sink > 0.5) /*SINK*/ { printf("Skipping sink! sink = %g, ID = %d\n", SphP[j].sink, P[j].ID); mass_j = 0; } #ifdef METALS_TG if(metal_disperse == 1) { D_j = 2.0*SphP[j].Hsml*SphP[j].Sigma; if(D_j < 1.0e-20) D_j = 1.0e-20; D_avg = 4.0*D_i*D_j/(D_i+D_j); if(r == 0.0) { dwk_r = (double) fabs(hinv4/h*KERNEL_COEFF_4); } else dwk_r = (double) fabs(dwk/r); K_ab = mass_j/SphP[j].Density*D_avg*dwk_r*hubble_param/a; const_A += K_ab; const_B += K_ab*SphP[j].Metallicity; } else { sigma += (P[j].Vel[0]-vel[0])*(P[j].Vel[0]-vel[0])+(P[j].Vel[1]-vel[1])*(P[j].Vel[1]-vel[1])+(P[j].Vel[2]-vel[2])*(P[j].Vel[2]-vel[2]); #endif rho += mass_j * wk; weighted_numngb += NORM_COEFF * wk / hinv3; dhsmlrho += -mass_j * ((double)(NUMDIMS) * hinv * wk + u * dwk); if(r > 0) { fac = mass_j * dwk / r; dvx = vel[0] - SphP[j].VelPred[0]; dvy = vel[1] - SphP[j].VelPred[1]; dvz = vel[2] - SphP[j].VelPred[2]; dBx = bfd[0] - SphP[j].bfield[0]; dBy = bfd[1] - SphP[j].bfield[1]; dBz = bfd[2] - SphP[j].bfield[2]; divv -= fac * (dx * dvx + dy * dvy + dz * dvz); divb -= fac * (dx * dBx + dy * dBy + dz * dBz); gsci[0] -= fac * dx * (sci - SphP[j].Sci); gsci[1] -= fac * dy * (sci - SphP[j].Sci); gsci[2] -= fac * dz * (sci - SphP[j].Sci); rotv[0] += fac * (dz * dvy - dy * dvz); rotv[1] += fac * (dx * dvz - dz * dvx); rotv[2] += fac * (dy * dvx - dx * dvy); vrel[0] += mass_j * wk * dvx; vrel[1] += mass_j * wk * dvy; vrel[2] += mass_j * wk * dvz; } #ifdef METALS_TG } #endif } } } while(startnode >= 0); #ifdef METALS_TG if(metal_disperse == 0) { #endif if(mode == 0) { SphP[target].NumNgb = weighted_numngb; SphP[target].Density = rho; SphP[target].DivVel = divv; SphP[target].DivB = divb; SphP[target].DhsmlDensityFactor = dhsmlrho; SphP[target].Rot[0] = rotv[0]; SphP[target].Rot[1] = rotv[1]; SphP[target].Rot[2] = rotv[2]; SphP[target].GradSci[0] = gsci[0]; SphP[target].GradSci[1] = gsci[1]; SphP[target].GradSci[2] = gsci[2]; SphP[target].VelRel[0] = vrel[0]; SphP[target].VelRel[1] = vrel[1]; SphP[target].VelRel[2] = vrel[2]; #ifdef METALS_TG SphP[target].Sigma = sigma; #endif } else { DensDataResult[target].Rho = rho; DensDataResult[target].Div = divv; DensDataResult[target].Div_B = divb; DensDataResult[target].Ngb = weighted_numngb; DensDataResult[target].DhsmlDensity = dhsmlrho; DensDataResult[target].Rot[0] = rotv[0]; DensDataResult[target].Rot[1] = rotv[1]; DensDataResult[target].Rot[2] = rotv[2]; DensDataResult[target].GSci[0] = gsci[0]; DensDataResult[target].GSci[1] = gsci[1]; DensDataResult[target].GSci[2] = gsci[2]; DensDataResult[target].GSci[2] = gsci[2]; DensDataResult[target].VRel[0] = vrel[0]; DensDataResult[target].VRel[1] = vrel[1]; DensDataResult[target].VRel[2] = vrel[2]; #ifdef METALS_TG DensDataResult[target].Sigma = sigma; #endif } #ifdef METALS_TG } else { if(mode == 0) { SphP[target].const_A = const_A; SphP[target].const_B = const_B; } else { DensDataResult[target].const_A = const_A; DensDataResult[target].const_B = const_B; } } #endif }