/** Initializes the (inverse) grid constant \ref ifcs_p3m_struct::a (\ref ifcs_p3m_struct::ai) and the cutoff for charge assignment \ref ifcs_p3m_struct::cao_cut, which has to be done by \ref ifcs_p3m_init_charges once and by \ref ifcs_p3m_scaleby_box_l whenever the \ref box_l changed. */ static void ifcs_p3m_prepare_a_ai_cao_cut(ifcs_p3m_data_struct *d) { P3M_DEBUG(printf(" ifcs_p3m_prepare_a_ai_cao_cut() started... \n")); for (fcs_int i=0; i<3; i++) { d->ai[i] = (fcs_float)d->grid[i]/d->box_l[i]; d->a[i] = 1.0/d->ai[i]; d->cao_cut[i] = 0.5*d->a[i]*d->cao; } P3M_DEBUG(printf(" ifcs_p3m_prepare_a_ai_cao_cut() finished. \n")); }
/** Prepare the data structures and constants of the P3M algorithm. All parameters have to be set. */ void ifcs_p3m_prepare(ifcs_p3m_data_struct *d, fcs_int max_charges) { P3M_DEBUG(printf(" ifcs_p3m_prepare() started... \n")); P3M_DEBUG(printf(" system parameters: box_l=(%" FCS_LMOD_FLOAT "f, %" FCS_LMOD_FLOAT "f, %" FCS_LMOD_FLOAT "f)\n", \ d->box_l[0], d->box_l[1], d->box_l[2])); P3M_DEBUG(printf( \ " p3m params: r_cut=%" FCS_LMOD_FLOAT "f, grid=%d, cao=%d, alpha=%" FCS_LMOD_FLOAT "f" \ ", grid_off=(%" FCS_LMOD_FLOAT "f,%" FCS_LMOD_FLOAT "f,%" FCS_LMOD_FLOAT "f)\n", \ d->r_cut, d->grid[0], d->cao, d->alpha, \ d->grid_off[0], d->grid_off[1], d->grid_off[2])); /* initializes the (inverse) grid constant d->a (d->ai) and the cutoff for charge assignment d->cao_cut */ ifcs_p3m_prepare_a_ai_cao_cut(d); ifcs_p3m_calc_local_ca_grid(d); ifcs_p3m_calc_send_grid(d); P3M_DEBUG(ifcs_p3m_print_local_grid(d->local_grid)); P3M_DEBUG(ifcs_p3m_print_send_grid(d->sm)); d->send_grid = (fcs_float *) realloc(d->send_grid, sizeof(fcs_float)*d->sm.max); d->recv_grid = (fcs_float *) realloc(d->recv_grid, sizeof(fcs_float)*d->sm.max); P3M_INFO(printf(" Interpolating charge assignement function...\n")); ifcs_p3m_interpolate_charge_assignment_function(d); /* position offset for calc. of first gridpoint */ d->pos_shift = (fcs_float)((d->cao-1)/2) - (d->cao%2)/2.0; P3M_DEBUG(printf(" pos_shift = %" FCS_LMOD_FLOAT "f\n",d->pos_shift)); /* FFT */ P3M_INFO(printf(" Preparing FFTs...\n")); ifcs_fft_prepare(&d->fft, &d->comm, &d->rs_grid, &d->ks_grid, d->local_grid.dim,d->local_grid.margin, d->grid, d->grid_off, &d->ks_pnum); /* k-space part */ ifcs_p3m_calc_differential_operator(d); P3M_INFO(printf(" Calculating influence function...\n")); #if !defined(P3M_INTERLACE) && defined(P3M_IK) ifcs_p3m_calc_influence_function_ik(d); #elif defined(P3M_INTERLACE) && defined(P3M_IK) ifcs_p3m_calc_influence_function_iki(d); #else ifcs_p3m_calc_influence_function_adi(d); #endif P3M_DEBUG(printf(" ifcs_p3m_prepare() finished.\n")); }
/** Calculates the properties of the send/recv sub-grides of the local * FFT grid. In order to calculate the recv sub-grides there is a * communication of the margins between neighbouring nodes. */ static void ifcs_p3m_calc_send_grid(ifcs_p3m_data_struct *d) { fcs_int i,j, evenodd; fcs_int done[3]={0,0,0}; MPI_Status status; P3M_DEBUG(printf(" ifcs_p3m_calc_send_grid() started... \n")); /* send grids */ for (i=0; i<3; i++) { for (j=0; j<3; j++) { /* left */ d->sm.s_ld[i*2][j] = 0 + done[j]*d->local_grid.margin[j*2]; if (j==i) d->sm.s_ur[i*2][j] = d->local_grid.margin[j*2]; else d->sm.s_ur[i*2][j] = d->local_grid.dim[j]-done[j]*d->local_grid.margin[(j*2)+1]; /* right */ if (j==i) d->sm.s_ld[(i*2)+1][j] = d->local_grid.in_ur[j]; else d->sm.s_ld[(i*2)+1][j] = 0 + done[j]*d->local_grid.margin[j*2]; d->sm.s_ur[(i*2)+1][j] = d->local_grid.dim[j] - done[j]*d->local_grid.margin[(j*2)+1]; } done[i]=1; } d->sm.max=0; for (i=0; i<6; i++) { d->sm.s_size[i] = 1; for (j=0; j<3; j++) { d->sm.s_dim[i][j] = d->sm.s_ur[i][j]-d->sm.s_ld[i][j]; d->sm.s_size[i] *= d->sm.s_dim[i][j]; } if (d->sm.s_size[i]>d->sm.max) d->sm.max=d->sm.s_size[i]; } /* communication */ for (i=0; i<6; i++) { if (i%2 == 0) j = i+1; else j = i-1; if (d->comm.node_neighbors[i] != d->comm.rank) { /* two step communication: first all even positions than all odd */ for (evenodd=0; evenodd<2; evenodd++) { if ((d->comm.node_pos[i/2]+evenodd)%2 == 0) { P3M_DEBUG(printf(" %d: sending local_grid.margin to %d\n", \ d->comm.rank, d->comm.node_neighbors[i])); MPI_Send(&(d->local_grid.margin[i]), 1, FCS_MPI_INT, d->comm.node_neighbors[i], 0, d->comm.mpicomm); } else { P3M_DEBUG(printf(" %d: receiving local_grid.margin from %d\n", \ d->comm.rank, d->comm.node_neighbors[j])); MPI_Recv(&(d->local_grid.r_margin[j]), 1, FCS_MPI_INT, d->comm.node_neighbors[j], 0, d->comm.mpicomm, &status); } } } else { d->local_grid.r_margin[j] = d->local_grid.margin[i]; } } /* /\* communication *\/ */ /* for (i = 0; i < 3; i++) { */ /* /\* upshift *\/ */ /* MPI_Sendrecv(&(d->local_grid.margin[2*i]), 1, FCS_MPI_INT, */ /* d->comm.node_neighbors[2*i+1], 0, */ /* &(d->local_grid.r_margin[2*i]), 1, FCS_MPI_INT, */ /* d->comm.node_neighbors[2*i], 0, */ /* d->comm.mpicomm, &status); */ /* /\* downshift *\/ */ /* MPI_Sendrecv(&(d->local_grid.margin[2*i+1]), 1, FCS_MPI_INT, */ /* d->comm.node_neighbors[2*i], 0, */ /* &(d->local_grid.r_margin[2*i+1]), 1, FCS_MPI_INT, */ /* d->comm.node_neighbors[2*i+1], 0, */ /* d->comm.mpicomm, &status); */ /* } */ /* recv grids */ for (i=0; i<3; i++) for (j=0; j<3; j++) { if (j==i) { d->sm.r_ld[ i*2 ][j] = d->sm.s_ld[ i*2 ][j] + d->local_grid.margin[2*j]; d->sm.r_ur[ i*2 ][j] = d->sm.s_ur[ i*2 ][j] + d->local_grid.r_margin[2*j]; d->sm.r_ld[(i*2)+1][j] = d->sm.s_ld[(i*2)+1][j] - d->local_grid.r_margin[(2*j)+1]; d->sm.r_ur[(i*2)+1][j] = d->sm.s_ur[(i*2)+1][j] - d->local_grid.margin[(2*j)+1]; } else { d->sm.r_ld[ i*2 ][j] = d->sm.s_ld[ i*2 ][j]; d->sm.r_ur[ i*2 ][j] = d->sm.s_ur[ i*2 ][j]; d->sm.r_ld[(i*2)+1][j] = d->sm.s_ld[(i*2)+1][j]; d->sm.r_ur[(i*2)+1][j] = d->sm.s_ur[(i*2)+1][j]; } } for (i=0; i<6; i++) { d->sm.r_size[i] = 1; for (j=0;j<3;j++) { d->sm.r_dim[i][j] = d->sm.r_ur[i][j]-d->sm.r_ld[i][j]; d->sm.r_size[i] *= d->sm.r_dim[i][j]; } if (d->sm.r_size[i]>d->sm.max) d->sm.max=d->sm.r_size[i]; } P3M_DEBUG(printf(" ifcs_p3m_calc_send_grid() finished. \n")); }
/** Calculates properties of the local FFT grid for the charge assignment process. */ static void ifcs_p3m_calc_local_ca_grid(ifcs_p3m_data_struct *d) { fcs_int i; fcs_int ind[3]; /* total skin size */ fcs_float full_skin[3]; P3M_DEBUG(printf(" ifcs_p3m_calc_local_ca_grid() started... \n")); for(i=0;i<3;i++) full_skin[i]= d->cao_cut[i]+d->skin+d->additional_grid[i]; /* inner left down grid point (global index) */ for(i=0;i<3;i++) d->local_grid.in_ld[i] = (fcs_int)ceil(d->comm.my_left[i]*d->ai[i]-d->grid_off[i]); /* inner up right grid point (global index) */ for(i=0;i<3;i++) d->local_grid.in_ur[i] = (fcs_int)floor(d->comm.my_right[i]*d->ai[i]-d->grid_off[i]); /* correct roundof errors at boundary */ for(i=0;i<3;i++) { if (fcs_float_is_zero((d->comm.my_right[i] * d->ai[i] - d->grid_off[i]) - d->local_grid.in_ur[i])) d->local_grid.in_ur[i]--; if (fcs_float_is_zero(1.0+(d->comm.my_left[i] * d->ai[i] - d->grid_off[i]) - d->local_grid.in_ld[i])) d->local_grid.in_ld[i]--; } /* inner grid dimensions */ for(i=0; i<3; i++) d->local_grid.inner[i] = d->local_grid.in_ur[i] - d->local_grid.in_ld[i] + 1; /* index of left down grid point in global grid */ for(i=0; i<3; i++) d->local_grid.ld_ind[i] = (fcs_int)ceil((d->comm.my_left[i]-full_skin[i])*d->ai[i]-d->grid_off[i]); /* spatial position of left down grid point */ ifcs_p3m_calc_lm_ld_pos(d); /* left down margin */ for(i=0;i<3;i++) d->local_grid.margin[i*2] = d->local_grid.in_ld[i]-d->local_grid.ld_ind[i]; /* up right grid point */ for(i=0;i<3;i++) ind[i] = (fcs_int)floor((d->comm.my_right[i]+full_skin[i])*d->ai[i]-d->grid_off[i]); /* correct roundof errors at up right boundary */ for(i=0;i<3;i++) if (((d->comm.my_right[i]+full_skin[i])*d->ai[i]-d->grid_off[i])-ind[i]==0) ind[i]--; /* up right margin */ for(i=0;i<3;i++) d->local_grid.margin[(i*2)+1] = ind[i] - d->local_grid.in_ur[i]; /* grid dimension */ d->local_grid.size=1; for(i=0;i<3;i++) { d->local_grid.dim[i] = ind[i] - d->local_grid.ld_ind[i] + 1; d->local_grid.size *= d->local_grid.dim[i]; } /* reduce inner grid indices from global to local */ for(i=0;i<3;i++) d->local_grid.in_ld[i] = d->local_grid.margin[i*2]; for(i=0;i<3;i++) d->local_grid.in_ur[i] = d->local_grid.margin[i*2]+d->local_grid.inner[i]; d->local_grid.q_2_off = d->local_grid.dim[2] - d->cao; d->local_grid.q_21_off = d->local_grid.dim[2] * (d->local_grid.dim[1] - d->cao); P3M_DEBUG(printf(" ifcs_p3m_calc_local_ca_grid() finished. \n")); }
void Communication::prepare(p3m_float box_l[3]) { P3M_DEBUG(printf( " P3M::Communication::prepare() started...\n")); /* Test whether the communicator is cartesian and correct dimensionality */ bool comm_is_cart = false; int status; MPI_Topo_test(mpicomm_orig, &status); if (status == MPI_CART) { /* Communicator is cartesian, so test dimensionality */ int ndims; MPI_Cartdim_get(mpicomm_orig, &ndims); if (ndims == 3) { /* Correct dimensionality, so get grid and test periodicity */ int periodicity[3]; MPI_Cart_get(mpicomm_orig, 3, node_grid, periodicity, node_pos); if (periodicity[0] && periodicity[1] && periodicity[2]) { /* If periodicity is correct, we can just use this communicator */ mpicomm = mpicomm_orig; /* get the rank */ MPI_Comm_rank(mpicomm, &rank); comm_is_cart = true; } } } /* otherwise, we have to set up the cartesian communicator */ if (!comm_is_cart) { P3M_DEBUG(printf( " Setting up cartesian communicator...\n")); node_grid[0] = 0; node_grid[1] = 0; node_grid[2] = 0; /* compute node grid */ MPI_Dims_create(size, 3, node_grid); #ifdef P3M_ENABLE_INFO if (onMaster()) printf(" node_grid=%dx%dx%d\n", node_grid[0], node_grid[1], node_grid[2]); #endif /* create communicator */ int periodicity[3] = {1, 1, 1}; MPI_Cart_create(mpicomm_orig, 3, node_grid, periodicity, 1, &mpicomm); /* get the rank */ MPI_Comm_rank(mpicomm, &rank); /* get node pos */ MPI_Cart_coords(mpicomm, rank, 3, node_pos); } /* fetch neighborhood info */ for (int dir = 0; dir < 3; dir++) { MPI_Cart_shift(mpicomm, dir, 1, &node_neighbors[2*dir], &node_neighbors[2*dir+1]); P3M_DEBUG_LOCAL(printf( " %d: dir=%d: n1=%d n2=%d\n", rank, dir, \ node_neighbors[2*dir], \ node_neighbors[2*dir+1])); } /* init local points */ for (int i=0; i< 3; i++) { local_box_l[i] = 0.0; my_left[i] = 0.0; my_right[i] = 0.0; } /* compute box limits */ for(p3m_int i = 0; i < 3; i++) { local_box_l[i] = box_l[i]/(p3m_float)node_grid[i]; my_left[i] = node_pos[i] *local_box_l[i]; my_right[i] = (node_pos[i]+1)*local_box_l[i]; } P3M_DEBUG(printf(" local_box_l=" F3FLOAT "\n" \ " my_left=" F3FLOAT "\n" \ " my_right=" F3FLOAT "\n", \ local_box_l[0], \ local_box_l[1], \ local_box_l[2], \ my_left[0], my_left[1], my_left[2], \ my_right[0], my_right[1], my_right[2] \ )); P3M_DEBUG(printf(" P3M::Communication::prepare() finished.\n")); }