double reb_communication_distance2_of_proc_to_node(struct reb_simulation* const r, int proc_id, struct reb_treecell* node){ int nghostxcol = (r->nghostx>0?1:0); int nghostycol = (r->nghosty>0?1:0); int nghostzcol = (r->nghostz>0?1:0); double distance2 = r->root_size*(double)r->root_n; // A conservative estimate for the minimum distance. distance2 *= distance2; for (int gbx=-nghostxcol; gbx<=nghostxcol; gbx++){ for (int gby=-nghostycol; gby<=nghostycol; gby++){ for (int gbz=-nghostzcol; gbz<=nghostzcol; gbz++){ struct reb_ghostbox gb = reb_boundary_get_ghostbox(r, gbx,gby,gbz); struct reb_aabb boundingbox = reb_communication_boundingbox_for_proc(r, proc_id); boundingbox.xmin+=gb.shiftx; boundingbox.xmax+=gb.shiftx; boundingbox.ymin+=gb.shifty; boundingbox.ymax+=gb.shifty; boundingbox.zmin+=gb.shiftz; boundingbox.zmax+=gb.shiftz; // calculate distance double distance2new = reb_communication_distance2_of_aabb_to_cell(boundingbox,node); if (distance2 > distance2new) distance2 = distance2new; } } } return distance2; }
void reb_display(void){ if (reb_dc.pause){ return; } sem_wait(reb_dc.mutex); const struct reb_particle* particles = reb_dc.r->particles; if (reb_dc.clear){ glClear(GL_DEPTH_BUFFER_BIT | GL_COLOR_BUFFER_BIT | GL_STENCIL_BUFFER_BIT); } glEnable(GL_POINT_SMOOTH); glVertexPointer(3, GL_DOUBLE, sizeof(struct reb_particle), particles); //int _N_active = ((N_active==-1)?N:N_active); if (reb_dc.reference>=0){ glTranslatef(-particles[reb_dc.reference].x,-particles[reb_dc.reference].y,-particles[reb_dc.reference].z); } for (int i=-reb_dc.ghostboxes*reb_dc.r->nghostx;i<=reb_dc.ghostboxes*reb_dc.r->nghostx;i++){ for (int j=-reb_dc.ghostboxes*reb_dc.r->nghosty;j<=reb_dc.ghostboxes*reb_dc.r->nghosty;j++){ for (int k=-reb_dc.ghostboxes*reb_dc.r->nghostz;k<=reb_dc.ghostboxes*reb_dc.r->nghostz;k++){ struct reb_ghostbox gb = reb_boundary_get_ghostbox(reb_dc.r, i,j,k); glTranslatef(gb.shiftx,gb.shifty,gb.shiftz); if (!(!reb_dc.clear&&reb_dc.wire)){ if (reb_dc.spheres==0 || reb_dc.spheres==2){ // Drawing Points glEnableClientState(GL_VERTEX_ARRAY); glPointSize(3.); glColor4f(1.0,1.0,1.0,0.5); //glDrawArrays(GL_POINTS, _N_active, N-_N_active); glColor4f(1.0,1.0,0.0,0.9); glPointSize(5.); glDrawArrays(GL_POINTS, 0, reb_dc.r->N-reb_dc.r->N_var); glDisableClientState(GL_VERTEX_ARRAY); } if (reb_dc.spheres){ glDisable(GL_BLEND); glEnable(GL_DEPTH_TEST); glEnable(GL_LIGHTING); glEnable(GL_LIGHT0); GLfloat lightpos[] = {0, reb_dc.r->boxsize_max, reb_dc.r->boxsize_max, 0.f}; glLightfv(GL_LIGHT0, GL_POSITION, lightpos); // Drawing Spheres glColor4f(1.0,1.0,1.0,1.0); for (int i=0;i<reb_dc.r->N-reb_dc.r->N_var;i++){ struct reb_particle p = particles[i]; if (p.r>0){ glTranslatef(p.x,p.y,p.z); glScalef(p.r,p.r,p.r); #ifdef _APPLE glCallList(reb_dc.dlist_sphere); #else //_APPLE glutSolidSphere(1,40,10); #endif //_APPLE glScalef(1./p.r,1./p.r,1./p.r); glTranslatef(-p.x,-p.y,-p.z); } } glEnable(GL_BLEND); glDisable(GL_DEPTH_TEST); glDisable(GL_LIGHTING); glDisable(GL_LIGHT0); } } // Drawing wires if (reb_dc.wire){ if(reb_dc.r->integrator!=REB_INTEGRATOR_SEI){ double radius = 0; struct reb_particle com = particles[0]; for (int i=1;i<reb_dc.r->N-reb_dc.r->N_var;i++){ struct reb_particle p = particles[i]; if (reb_dc.r->N_active>0){ // Different colors for active/test particles if (i>=reb_dc.r->N_active){ glColor4f(0.9,1.0,0.9,0.9); }else{ glColor4f(1.0,0.9,0.0,0.9); } }else{ // Alternating colors if (i%2 == 1){ glColor4f(0.0,1.0,0.0,0.9); }else{ glColor4f(0.0,0.0,1.0,0.9); } } //if (reb_dc.r->integrator==REB_INTEGRATOR_WHFAST && reb_dc.r->ri_whfast.is_synchronized==0){ // double m = p.m; // p = reb_dc.r->ri_whfast.p_j[i]; // p.m = m; //} struct reb_orbit o = reb_tools_particle_to_orbit(reb_dc.r->G, p,com); glPushMatrix(); glTranslatef(com.x,com.y,com.z); glRotatef(o.Omega/DEG2RAD,0,0,1); glRotatef(o.inc/DEG2RAD,1,0,0); glRotatef(o.omega/DEG2RAD,0,0,1); glBegin(GL_LINE_LOOP); for (double trueAnom=0; trueAnom < 2.*M_PI; trueAnom+=M_PI/100.) { //convert degrees into radians radius = o.a * (1. - o.e*o.e) / (1. + o.e*cos(trueAnom)); glVertex3f(radius*cos(trueAnom),radius*sin(trueAnom),0); } glEnd(); glPopMatrix(); com = reb_get_com_of_pair(p,com); } }else{ for (int i=1;i<reb_dc.r->N;i++){ struct reb_particle p = particles[i]; glBegin(GL_LINE_LOOP); for (double _t=-100.*reb_dc.r->dt;_t<=100.*reb_dc.r->dt;_t+=20.*reb_dc.r->dt){ double frac = 1.-fabs(_t/(120.*reb_dc.r->dt)); glColor4f(1.0,(_t+100.*reb_dc.r->dt)/(200.*reb_dc.r->dt),0.0,frac); glVertex3f(p.x+p.vx*_t, p.y+p.vy*_t, p.z+p.vz*_t); } glEnd(); } } } glTranslatef(-gb.shiftx,-gb.shifty,-gb.shiftz); } } } glColor4f(1.0,0.0,0.0,0.4); glScalef(reb_dc.r->boxsize.x,reb_dc.r->boxsize.y,reb_dc.r->boxsize.z); if (reb_dc.r->boundary == REB_BOUNDARY_NONE){ glBegin(GL_LINES); glVertex3f(0,0,0.04); glVertex3f(0,0,-0.04); glVertex3f(0,0.04,0); glVertex3f(0,-0.04,0); glVertex3f(0.04,0,0); glVertex3f(-0.04,0,0); glEnd(); }else{ glutWireCube(1); } glScalef(1./reb_dc.r->boxsize.x,1./reb_dc.r->boxsize.y,1./reb_dc.r->boxsize.z); if (reb_dc.reference>=0){ glTranslatef(particles[reb_dc.reference].x,particles[reb_dc.reference].y,particles[reb_dc.reference].z); } glFlush(); sem_post(reb_dc.mutex); }
/** * Main Gravity Routine */ void reb_calculate_acceleration(struct reb_simulation* r){ struct reb_particle* const particles = r->particles; const int N = r->N; const int N_active = r->N_active; const double G = r->G; const double softening2 = r->softening*r->softening; const unsigned int _gravity_ignore_10 = r->gravity_ignore_10; const int _N_start = (r->integrator==REB_INTEGRATOR_WH?1:0); const int _N_active = ((N_active==-1)?N:N_active) - r->N_var; const int _N_real = N - r->N_var; switch (r->gravity){ case REB_GRAVITY_NONE: // Do nothing. break; case REB_GRAVITY_BASIC: { const int nghostx = r->nghostx; const int nghosty = r->nghosty; const int nghostz = r->nghostz; #pragma omp parallel for schedule(guided) for (int i=0; i<N; i++){ particles[i].ax = 0; particles[i].ay = 0; particles[i].az = 0; } // Summing over all Ghost Boxes for (int gbx=-nghostx; gbx<=nghostx; gbx++){ for (int gby=-nghosty; gby<=nghosty; gby++){ for (int gbz=-nghostz; gbz<=nghostz; gbz++){ struct reb_ghostbox gb = reb_boundary_get_ghostbox(r, gbx,gby,gbz); // Summing over all particle pairs #pragma omp parallel for schedule(guided) for (int i=_N_start; i<_N_real; i++){ for (int j=_N_start; j<_N_active; j++){ if (_gravity_ignore_10 && j==1 && i==0 ) continue; if (i==j) continue; const double dx = (gb.shiftx+particles[i].x) - particles[j].x; const double dy = (gb.shifty+particles[i].y) - particles[j].y; const double dz = (gb.shiftz+particles[i].z) - particles[j].z; const double _r = sqrt(dx*dx + dy*dy + dz*dz + softening2); const double prefact = -G/(_r*_r*_r)*particles[j].m; particles[i].ax += prefact*dx; particles[i].ay += prefact*dy; particles[i].az += prefact*dz; } } } } } } break; case REB_GRAVITY_COMPENSATED: { if (r->gravity_cs_allocatedN<_N_real){ r->gravity_cs = realloc(r->gravity_cs,_N_real*sizeof(struct reb_vec3d)); r->gravity_cs_allocatedN = _N_real; } struct reb_vec3d* restrict const cs = r->gravity_cs; #pragma omp parallel for schedule(guided) for (int i=0; i<_N_real; i++){ particles[i].ax = 0.; particles[i].ay = 0.; particles[i].az = 0.; cs[i].x = 0.; cs[i].y = 0.; cs[i].z = 0.; } // Summing over all massive particle pairs #pragma omp parallel for schedule(guided) for (int i=_N_start; i<_N_active; i++){ for (int j=i+1; j<_N_active; j++){ if (_gravity_ignore_10 && j==1 && i==0 ) continue; const double dx = particles[i].x - particles[j].x; const double dy = particles[i].y - particles[j].y; const double dz = particles[i].z - particles[j].z; const double r2 = dx*dx + dy*dy + dz*dz + softening2; const double r = sqrt(r2); const double prefact = -G/(r2*r); const double prefacti = prefact*particles[i].m; const double prefactj = prefact*particles[j].m; { double ax = particles[i].ax; cs[i].x += prefactj*dx; particles[i].ax = ax + cs[i].x; cs[i].x += ax - particles[i].ax; double ay = particles[i].ay; cs[i].y += prefactj*dy; particles[i].ay = ay + cs[i].y; cs[i].y += ay - particles[i].ay; double az = particles[i].az; cs[i].z += prefactj*dz; particles[i].az = az + cs[i].z; cs[i].z += az - particles[i].az; } { double ax = particles[j].ax; cs[j].x -= prefacti*dx; particles[j].ax = ax + cs[j].x; cs[j].x += ax - particles[j].ax; double ay = particles[j].ay; cs[j].y -= prefacti*dy; particles[j].ay = ay + cs[j].y; cs[j].y += ay - particles[j].ay; double az = particles[j].az; cs[j].z -= prefacti*dz; particles[j].az = az + cs[j].z; cs[j].z += az - particles[j].az; } } } // Testparticles #pragma omp parallel for schedule(guided) for (int i=_N_active; i<_N_real; i++){ for (int j=_N_start; j<_N_active; j++){ if (_gravity_ignore_10 && i==1 && j==0 ) continue; const double dx = particles[i].x - particles[j].x; const double dy = particles[i].y - particles[j].y; const double dz = particles[i].z - particles[j].z; const double r2 = dx*dx + dy*dy + dz*dz + softening2; const double r = sqrt(r2); const double prefact = -G/(r2*r)*particles[j].m; double ax = particles[i].ax; cs[i].x += prefact*dx; particles[i].ax = ax + cs[i].x; cs[i].x += ax - particles[i].ax; double ay = particles[i].ay; cs[i].y += prefact*dy; particles[i].ay = ay + cs[i].y; cs[i].y += ay - particles[i].ay; double az = particles[i].az; cs[i].z += prefact*dz; particles[i].az = az + cs[i].z; cs[i].z += az - particles[i].az; } } } break; case REB_GRAVITY_TREE: { #pragma omp parallel for schedule(guided) for (int i=0; i<N; i++){ particles[i].ax = 0; particles[i].ay = 0; particles[i].az = 0; } // Summing over all Ghost Boxes for (int gbx=-r->nghostx; gbx<=r->nghostx; gbx++){ for (int gby=-r->nghosty; gby<=r->nghosty; gby++){ for (int gbz=-r->nghostz; gbz<=r->nghostz; gbz++){ // Summing over all particle pairs #pragma omp parallel for schedule(guided) for (int i=0; i<N; i++){ struct reb_ghostbox gb = reb_boundary_get_ghostbox(r, gbx,gby,gbz); // Precalculated shifted position gb.shiftx += particles[i].x; gb.shifty += particles[i].y; gb.shiftz += particles[i].z; reb_calculate_acceleration_for_particle(r, i, gb); } } } } } break; default: fprintf(stderr,"\n\033[1mError!\033[0m Gravity calculation not yet implemented.\n"); exit(EXIT_FAILURE); break; break; } #ifdef MPI // Distribute active particles from root to all other nodes. // This assures that round-off errors do not accumulate and // the copies of active particles do not diverge. MPI_Bcast(particles, N_active, mpi_particle, 0, MPI_COMM_WORLD); #endif }