/** * compute the correction to the field and total energy */ FCSResult fcs_compute_dipole_correction(FCS handle, fcs_int local_particles, fcs_float* positions, fcs_float *charges, fcs_float epsilon, fcs_float *field_correction, fcs_float *energy_correction) { const char *fnc_name = "fcs_compute_dipole_correction"; CHECK_HANDLE_RETURN_RESULT(handle, fnc_name); /* Local dipole moment */ fcs_float local_dipole_moment[3] = {0.0, 0.0, 0.0}; /* Global dipole moment */ fcs_float dipole_moment[3]; fcs_int pid; fcs_int dim; if (fcs_float_is_zero(epsilon) || epsilon > 0.0) { /* Compute the global dipole moment */ for (pid = 0; pid < local_particles; pid++) for (dim = 0; dim < 3; dim++) local_dipole_moment[dim] += charges[pid]*positions[pid*3+dim]; MPI_Allreduce(local_dipole_moment, dipole_moment, 3, FCS_MPI_FLOAT, MPI_SUM, handle->communicator); const fcs_float *a = fcs_get_box_a(handle); const fcs_float *b = fcs_get_box_b(handle); const fcs_float *c = fcs_get_box_c(handle); /* Volume of the parallelepiped */ fcs_float volume = a[0] * (b[1]*c[2] - b[2]*c[1]) + a[1] * (b[2]*c[0] - b[0]*c[2]) + a[2] * (b[0]*c[1] - b[1]*c[0]); fcs_float pref = 4.0*3.14159265358979323846264338328 / (3.0*volume*(epsilon + 1.0)); if (energy_correction) *energy_correction = 0.5*pref*(dipole_moment[0]*dipole_moment[0] + dipole_moment[1]*dipole_moment[1] + dipole_moment[2]*dipole_moment[2]); if (field_correction) { field_correction[0] = -pref*dipole_moment[0]; field_correction[1] = -pref*dipole_moment[1]; field_correction[2] = -pref*dipole_moment[2]; } } else { /* metallic BC (epsilon=+infty) */ if (energy_correction) *energy_correction = 0.0; if (field_correction) { field_correction[0] = 0.0; field_correction[1] = 0.0; field_correction[2] = 0.0; } } return FCS_RESULT_SUCCESS; }
static void mmm_preparePolygammaOdd(fcs_int n, fcs_float binom, SizedList *series) { fcs_int order; fcs_float deriv; fcs_float maxx, x_order, coeff, pref; deriv = 2*n + 1; maxx = 0.5; // to get 1/(2n)! instead of 1/(2n+1)! pref = 2*deriv*(1 + deriv); mmm_init_doublelist(series); for (order = 0;; order++) { // only odd exponents of x x_order = 2*order + 1; coeff = pref*mmm_hzeta(1 + deriv + x_order, 2); if ((fcs_float_is_zero(fabs(maxx*coeff)*(4.0/3.0))) && (x_order > deriv)) break; mmm_realloc_doublelist(series, order + 1); series->e[order] = -binom*coeff; maxx *= 0.25; pref *= (1.0 + deriv/(x_order + 1)); pref *= (1.0 + deriv/(x_order + 2)); } series->n = order; }
/* compute the net charge of the system */ FCSResult mmm2d_check_system_charges(mmm2d_data_struct *d, fcs_int num_particles, fcs_float *charges) { /* fcs_int i; *totcharge=0.0; *anycharge=0; for (i = 0; i < num_particles; i++) { if (!fcs_float_is_zero(charges[i])) { *totcharge += charges[i]; *anycharge=1; } } return NULL; */ fcs_int i; fcs_float local_charge=0., total_charge=0.; for (i = 0; i < num_particles; i++) { if (!fcs_float_is_zero(charges[i])) { local_charge += charges[i]; } } fprintf(stderr,"check charges: rank %d\n",d->comm.rank); fprintf(stderr, "local_charges %d, net charge %e\n", num_particles, local_charge); MPI_Allreduce(&local_charge, &total_charge, 1, FCS_MPI_FLOAT, MPI_SUM, d->comm.mpicomm); d->total_charge=total_charge; fprintf(stderr, "total charge %e\n", d->total_charge); return NULL; }
static void mmm_preparePolygammaEven(fcs_int n, fcs_float binom, SizedList *series) { /* (-0.5 n) psi^2n/2n! (-0.5 n) and psi^(2n+1)/(2n)! series expansions note that BOTH carry 2n! */ fcs_int order; fcs_float deriv; fcs_float maxx, x_order, coeff, pref; deriv = 2*n; if (n == 0) { // psi^0 has a slightly different series expansion maxx = 0.25; mmm_alloc_doublelist(series, 1); series->e[0] = 2*(1 - MMM_COMMON_C_GAMMA); for (order = 1;; order += 1) { x_order = 2*order; coeff = -2*mmm_hzeta(x_order + 1, 2); if (fcs_float_is_zero(fabs(maxx*coeff)*(4.0/3.0))) break; mmm_realloc_doublelist(series, order + 1); series->e[order] = coeff; maxx *= 0.25; } series->n = order; } else { // even, n > 0 maxx = 1; pref = 2; mmm_init_doublelist(series); for (order = 0;; order++) { // only even exponents of x x_order = 2*order; coeff = pref*mmm_hzeta(1 + deriv + x_order, 2); if ((fcs_float_is_zero(fabs(maxx*coeff)*(4.0/3.0))) && (x_order > deriv)) break; mmm_realloc_doublelist(series, order + 1); series->e[order] = -binom*coeff; maxx *= 0.25; pref *= (1.0 + deriv/(x_order + 1)); pref *= (1.0 + deriv/(x_order + 2)); } series->n = order; } }
FCSResult mmm2d_tune(void* rd, fcs_int num_particles, fcs_float *positions, fcs_float *charges) { mmm2d_data_struct *d = (mmm2d_data_struct*)rd; const char* fnc_name = "mmm2d_tune"; FCSResult res; /* Check for charge existence and neutrality */ mmm2d_check_system_charges(d, num_particles, charges); if (!fcs_float_is_zero(d->total_charge)) return fcs_result_create(FCS_ERROR_LOGICAL_ERROR, fnc_name, "MMM2D requires a zero net charge."); d->my_bottom = d->comm.rank*d->box_l[2]/(fcs_float)(d->comm.size); /* Broadcast charges */ /* MPI_Bcast(&num_particles, 1, FCS_MPI_INT, 0, d->comm.mpicomm); MPI_Bcast(charges, num_particles, FCS_MPI_FLOAT, 0, d->comm.mpicomm); MPI_Bcast(positions, 3*num_particles, FCS_MPI_FLOAT, 0, d->comm.mpicomm); */ ///@TODO: skip next steps if d->needs_tuning==0. Check if this is flawless if(d->needs_tuning==0) return NULL; /* precalculate some constants */ mmm2d_setup_constants(d); /* tune near formula */ res=mmm2d_tune_near(d); if (res) return res; /* tune far formula */ if (d->comm.size*d->layers_per_node < 3) { d->far_cut = 0.0; if (d->dielectric_contrast_on) return fcs_result_create(FCS_ERROR_LOGICAL_ERROR, fnc_name, "Definition of dielectric contrasts requires more than 3 layers"); } else { res=mmm2d_tune_far(d); if (res) return res; } d->needs_tuning=0; /* printf("node %d, mmm2d_tune\n", d->comm.rank); printf("node %d: ux: %f, uy: %f\n", d->comm.rank, d->ux, d->uy); printf("node %d: min_far: %f, far_cut: %f\n", d->comm.rank, d->min_far, d->far_cut); printf("node %d: layer_h %f\n", d->comm.rank, d->layer_h); printf("node %d: n_layers %d\n", d->comm.rank, d->n_layers); */ return NULL; }
/** 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")); }
/*************************************************** **** K-SPACE CONTRIBUTION ***************************************************/ void ewald_compute_kspace(ewald_data_struct* d, fcs_int num_particles, fcs_float *positions, fcs_float *charges, fcs_float *fields, fcs_float *potentials) { FCS_INFO(fprintf(stderr, "ewald_compute_kspace started...\n")); /* DISTRIBUTE ALL PARTICLE DATA TO ALL NODES */ /* Gather all particle numbers */ int node_num_particles = num_particles; int node_particles[d->comm_size]; int node_particles3[d->comm_size]; int displs[d->comm_size]; int displs3[d->comm_size]; fcs_int total_particles; /* printf("%d: num_particles=%d\n", d->comm_rank, num_particles); */ MPI_Allgather(&node_num_particles, 1, MPI_INT, node_particles, 1, MPI_INT, d->comm); /* compute displacements for MPI_Gatherv */ total_particles = node_particles[0]; node_particles3[0] = node_particles[0]*3; displs[0] = 0; displs3[0] = 0; for (fcs_int i=1; i < d->comm_size; i++) { total_particles += node_particles[i]; node_particles3[i] = node_particles[i]*3; displs[i] = displs[i-1] + node_particles[i-1]; displs3[i] = displs3[i-1] + node_particles3[i-1]; } fcs_float *all_positions, *all_charges, *node_fields, *node_potentials; all_positions = malloc(sizeof(fcs_float) * 3 * total_particles); all_charges = malloc(sizeof(fcs_float) * total_particles); node_fields = malloc(sizeof(fcs_float) * 3 * total_particles); node_potentials = malloc(sizeof(fcs_float) * total_particles); /* gather all particle data at all nodes */ MPI_Allgatherv(positions, num_particles*3, FCS_MPI_FLOAT, all_positions, node_particles3, displs3, FCS_MPI_FLOAT, d->comm); MPI_Allgatherv(charges, num_particles, FCS_MPI_FLOAT, all_charges, node_particles, displs, FCS_MPI_FLOAT, d->comm); /* for (fcs_int i=0; i < total_particles; i++) { */ /* printf("%d: all_positions[%d]={%lf, %lf, %lf}\n", d->comm_rank, i, all_positions[i*3], all_positions[3*i+1], all_positions[3*i+2]); */ /* printf("%d: all_charges[%d]=%lf\n", d->comm_rank, i, all_charges[i]); */ /* } */ /* INIT ALGORITHM */ /* init fields and potentials */ if (fields != NULL) { for (fcs_int i=0; i < total_particles; i++) { node_fields[3*i] = 0.0; node_fields[3*i+1] = 0.0; node_fields[3*i+2] = 0.0; } } if (potentials != NULL) for (fcs_int i=0; i < total_particles; i++) node_potentials[i] = 0.0; /* COMPUTE FAR FIELDS */ /* evenly distribute the k-vectors onto all tasks */ fcs_int num_k_per_dir = 2*d->kmax+1; fcs_int num_k = num_k_per_dir * num_k_per_dir * num_k_per_dir; for (int k_ind = d->comm_rank; k_ind < num_k; k_ind += d->comm_size) { /* compute fields and potentials */ fcs_int nx = k_ind % num_k_per_dir - d->kmax; fcs_int ny = k_ind % (num_k_per_dir*num_k_per_dir) / num_k_per_dir - d->kmax; fcs_int nz = k_ind / (num_k_per_dir*num_k_per_dir) - d->kmax; if (nx*nx + ny*ny + nz*nz <= d->kmax*d->kmax) { /* system length vector L_vec */ const fcs_float Lx = d->box_l[0]; const fcs_float Ly = d->box_l[1]; const fcs_float Lz = d->box_l[2]; /* reciprocal vector k_vec */ const fcs_float kx = 2.0*M_PI*nx / Lx; const fcs_float ky = 2.0*M_PI*ny / Ly; const fcs_float kz = 2.0*M_PI*nz / Lz; /* reciprocal charge density rhohat */ fcs_float rhohat_re = 0.0; fcs_float rhohat_im = 0.0; /* compute Deserno, Holm (1998) eq. (8) */ for (fcs_int i=0; i < total_particles; i++) { /* charge q */ const fcs_float q = all_charges[i]; if (!fcs_float_is_zero(q)) { /* particle position r_vec */ const fcs_float rx = all_positions[3*i]; const fcs_float ry = all_positions[3*i+1]; const fcs_float rz = all_positions[3*i+2]; /* compute k_vec*r_vec */ fcs_float kr = kx*rx + ky*ry + kz*rz; /* rhohat = qi * exp(-i*k_vec*r_vec) */ rhohat_re += q * cos(kr); rhohat_im += q * -sin(kr); } } /* FCS_DEBUG(fprintf(stderr, " n_vec= (%d, %d, %d) rhohat_re=%e rhohat_im=%e\n", nx, ny, nz, rhohat_re, rhohat_im));*/ /* fetch influence function */ fcs_float g = d->G[linindex(abs(nx), abs(ny), abs(nz), d->kmax)]; for (fcs_int i=0; i < total_particles; i++) { /* particle position r_vec */ const fcs_float rx = all_positions[3*i]; const fcs_float ry = all_positions[3*i+1]; const fcs_float rz = all_positions[3*i+2]; /* compute k_vec*r_vec */ fcs_float kr = kx*rx + ky*ry + kz*rz; if (fields != NULL) { /* compute field at position of particle i compare to Deserno, Holm (1998) eq. (15) */ fcs_float fak1 = g * (rhohat_re*sin(kr) + rhohat_im*cos(kr)); node_fields[3*i] += kx * fak1; node_fields[3*i+1] += ky * fak1; node_fields[3*i+2] += kz * fak1; } if (potentials != NULL) { /* compute potential at position of particle i compare to Deserno, Holm (1998) eq. (9) */ node_potentials[i] += g * (rhohat_re*cos(kr) - rhohat_im*sin(kr)); } } } } /* printf("%d: node_fields[0]=%lf\n", d->comm_rank, node_fields[0]); */ /* REDISTRIBUTE COMPUTED FAR FIELDS AND POTENTIALS */ if (fields != NULL) { fcs_float all_fields[total_particles*3]; for (fcs_int pid=0; pid < total_particles; pid++) { all_fields[3*pid] = 0.0; all_fields[3*pid+1] = 0.0; all_fields[3*pid+2] = 0.0; } /* Combine all fields on master */ MPI_Reduce(node_fields, all_fields, total_particles*3, FCS_MPI_FLOAT, MPI_SUM, 0, d->comm); /* if (d->comm_rank == 0) */ /* printf("%d: all_fields[0]=%lf\n", d->comm_rank, all_fields[0]); */ /* Scatter the fields to the task that holds the particle */ MPI_Scatterv(all_fields, node_particles3, displs3, FCS_MPI_FLOAT, fields, num_particles*3, FCS_MPI_FLOAT, 0, d->comm); } if (potentials != NULL) { fcs_float all_potentials[total_particles]; for (fcs_int pid=0; pid < total_particles; pid++) all_potentials[pid] = 0.0; /* Combine all potentials on master */ MPI_Reduce(node_potentials, all_potentials, total_particles, FCS_MPI_FLOAT, MPI_SUM, 0, d->comm); /* Scatter the potentials to the task that holds the particle */ MPI_Scatterv(all_potentials, node_particles, displs, FCS_MPI_FLOAT, potentials, num_particles, FCS_MPI_FLOAT, 0, d->comm); /* subtract self potential */ FCS_INFO(fprintf(stderr, " subtracting self potential...\n")); for (fcs_int i=0; i < num_particles; i++) { potentials[i] -= charges[i] * M_2_SQRTPI * d->alpha; } } free(all_positions); free(all_charges); free(node_fields); free(node_potentials); /* now each task should have its far field components */ FCS_INFO(fprintf(stderr, "ewald_compute_kspace finished.\n")); }