void fcs_pair_int(real *pot, real *grad, real r2) { fcs_float r, fcs_pot, fcs_grad; FCSResult result; r = SQRT(r2); result = fcs_compute_near(handle, r, &fcs_pot, &fcs_grad); ASSERT_FCS(result); *pot = fcs_pot; /* return (1/r)*dV/dr as derivative */ *grad = fcs_grad / r; }
void calc_forces_fcs(int steps) { FCSResult result; fcs_float BoxX[3] = { box_x.x, box_x.y, box_x.z }; fcs_float BoxY[3] = { box_y.x, box_y.y, box_y.z }; fcs_float BoxZ[3] = { box_z.x, box_z.y, box_z.z }; #ifdef HOMDEF if ((lindef_int > 0) && (0 == steps % lindef_int)) { result = fcs_set_box_a(handle, BoxX); ASSERT_FCS(result); result = fcs_set_box_b(handle, BoxY); ASSERT_FCS(result); result = fcs_set_box_c(handle, BoxZ); ASSERT_FCS(result); } #endif pack_fcs(); /* dump_config_fcs( pos, chg, nloc, steps, myid ); */ result = fcs_run(handle, nloc, nloc_max, pos, chg, field, pot); ASSERT_FCS(result); unpack_fcs(); }
void fcs_cleanup(void) { FCSResult result; result = fcs_destroy(handle); ASSERT_FCS(result); }
void init_fcs(void) { FCSResult res; fcs_int srf = 1; char *method; fcs_int pbc [3] = { pbc_dirs.x, pbc_dirs.y, pbc_dirs.z }; fcs_float BoxX[3] = { box_x.x, box_x.y, box_x.z }; fcs_float BoxY[3] = { box_y.x, box_y.y, box_y.z }; fcs_float BoxZ[3] = { box_z.x, box_z.y, box_z.z }; fcs_float off [3] = { 0.0, 0.0, 0.0 }; /* subtract CM momentum */ if (0 == imdrestart) { int i, k; real ptot[4], ptot_2[4], px, py, pz; ptot[0] = 0.0; ptot[1] = 0.0; ptot[2] = 0.0, ptot[3] = 0.0; for (k=0; k<NCELLS; ++k) { /* loop over all cells */ cell *p = CELLPTR(k); for (i=0; i<p->n; i++) { ptot[0] += IMPULS(p,i,X); ptot[1] += IMPULS(p,i,Y); ptot[2] += IMPULS(p,i,Z); ptot[3] += MASSE(p,i); } } #ifdef MPI MPI_Allreduce( ptot, ptot_2, 4, REAL, MPI_SUM, cpugrid); ptot[0] = ptot_2[0]; ptot[1] = ptot_2[1]; ptot[2] = ptot_2[2]; ptot[3] = ptot_2[3]; #endif px = ptot[0]/ptot[3]; py = ptot[1]/ptot[3]; pz = ptot[2]/ptot[3]; for (k=0; k<NCELLS; ++k) { /* loop over all cells */ cell *p = CELLPTR(k); for (i=0; i<p->n; i++) { IMPULS(p,i,X) -= px * MASSE(p,i); IMPULS(p,i,Y) -= py * MASSE(p,i); IMPULS(p,i,Z) -= pz * MASSE(p,i); } } } switch (fcs_method) { case FCS_METH_DIRECT: method = "direct"; break; case FCS_METH_PEPC: method = "pepc"; break; case FCS_METH_FMM: method = "fmm"; break; case FCS_METH_P3M: method = "p3m"; srf = fcs_near_field_flag; break; case FCS_METH_P2NFFT: method = "p2nfft"; srf = fcs_near_field_flag; break; case FCS_METH_VMG: method = "vmg"; break; case FCS_METH_PP3MG: method = "pp3mg"; break; } /* initialize handle and set common parameters */ res = fcs_init(&handle, method, cpugrid); ASSERT_FCS(res); res = fcs_set_common(handle, srf, BoxX, BoxY, BoxZ, off, pbc, natoms); ASSERT_FCS(res); res = fcs_require_virial(handle, 1); ASSERT_FCS(res); /* set method specific parameters */ switch (fcs_method) { #ifdef FCS_ENABLE_DIRECT case FCS_METH_DIRECT: /* nothing to do */ break; #endif #ifdef FCS_ENABLE_PEPC case FCS_METH_PEPC: res = fcs_pepc_setup(handle, (fcs_float)fcs_pepc_eps, (fcs_float)fcs_pepc_theta ); ASSERT_FCS(res); res = fcs_pepc_set_num_walk_threads( handle, (fcs_int)fcs_pepc_nthreads ); ASSERT_FCS(res); break; #endif #ifdef FCS_ENABLE_FMM case FCS_METH_FMM: res = fcs_fmm_set_absrel(handle, (fcs_int)fcs_fmm_absrel); ASSERT_FCS(res); res = fcs_fmm_set_tolerance_energy(handle, (fcs_float)fcs_tolerance); ASSERT_FCS(res); break; #endif #ifdef FCS_ENABLE_P3M case FCS_METH_P3M: if (0==srf) { res = fcs_p3m_set_r_cut(handle, (fcs_float)fcs_rcut); ASSERT_FCS(res); } res = fcs_set_tolerance(handle, FCS_TOLERANCE_TYPE_FIELD, (fcs_float)fcs_tolerance); ASSERT_FCS(res); if (fcs_grid_dim.x) res = fcs_p3m_set_grid(handle, (fcs_int)fcs_grid_dim.x); ASSERT_FCS(res); break; #endif #ifdef FCS_ENABLE_P2NFFT case FCS_METH_P2NFFT: if (0==srf) { res = fcs_p2nfft_set_r_cut(handle, (fcs_float)fcs_rcut); ASSERT_FCS(res); } res = fcs_set_tolerance(handle, FCS_TOLERANCE_TYPE_FIELD, (fcs_float)fcs_tolerance); ASSERT_FCS(res); if (fcs_grid_dim.x) { res = fcs_p2nfft_set_grid(handle, (fcs_int)fcs_grid_dim.x, (fcs_int)fcs_grid_dim.y, (fcs_int)fcs_grid_dim.z); ASSERT_FCS(res); } if (fcs_p2nfft_intpol_order) { res = fcs_p2nfft_set_pnfft_interpolation_order(handle, (fcs_int)fcs_p2nfft_intpol_order); ASSERT_FCS(res); } if (fcs_p2nfft_epsI) { res = fcs_p2nfft_set_epsI(handle, (fcs_float)fcs_p2nfft_epsI); ASSERT_FCS(res); } //res = fcs_p2nfft_set_pnfft_window_by_name(handle, "bspline"); //ASSERT_FCS(res); break; #endif #ifdef FCS_ENABLE_VMG case FCS_METH_VMG: if (fcs_vmg_near_field_cells) { res = fcs_vmg_set_near_field_cells(handle, (fcs_int)fcs_vmg_near_field_cells); ASSERT_FCS(res); } if (fcs_vmg_interpol_order) { res = fcs_vmg_set_interpolation_order(handle, (fcs_int)fcs_vmg_interpol_order); ASSERT_FCS(res); } if (fcs_vmg_discr_order) { res = fcs_vmg_set_discretization_order(handle, (fcs_int)fcs_vmg_discr_order); ASSERT_FCS(res); } if (fcs_iter_tolerance > 0) { res = fcs_vmg_set_precision(handle, (fcs_float)fcs_iter_tolerance); ASSERT_FCS(res); } break; #endif #ifdef FCS_ENABLE_PP3MG case FCS_METH_PP3MG: /* use default values, if not specified otherwise */ if (fcs_grid_dim.x) { res = fcs_pp3mg_set_cells_x(handle, (fcs_int)fcs_grid_dim.x); ASSERT_FCS(res); res = fcs_pp3mg_set_cells_y(handle, (fcs_int)fcs_grid_dim.y); ASSERT_FCS(res); res = fcs_pp3mg_set_cells_z(handle, (fcs_int)fcs_grid_dim.z); ASSERT_FCS(res); } if (fcs_pp3mg_ghosts) { res = fcs_pp3mg_set_ghosts(handle, (fcs_int)fcs_pp3mg_ghosts); ASSERT_FCS(res); } if (fcs_pp3mg_degree) { res = fcs_pp3mg_set_degree(handle, (fcs_int)fcs_pp3mg_degree); ASSERT_FCS(res); } if (fcs_pp3mg_max_part) { res = fcs_pp3mg_set_max_particles(handle, (fcs_int)fcs_pp3mg_max_part); ASSERT_FCS(res); } if (fcs_max_iter) { res = fcs_pp3mg_set_max_iterations(handle,(fcs_int)fcs_max_iter); ASSERT_FCS(res); } if (fcs_iter_tolerance > 0) { res = fcs_pp3mg_set_tol(handle, (fcs_float)fcs_iter_tolerance); ASSERT_FCS(res); } break; #endif default: error("FCS method unknown or not implemented"); break; } pack_fcs(); res = fcs_tune(handle, nloc, nloc_max, pos, chg); ASSERT_FCS(res); /* inform about tuned parameters */ switch (fcs_method) { fcs_int grid_dim[3]; fcs_float r_cut; #ifdef FCS_ENABLE_P3M case FCS_METH_P3M: res = fcs_p3m_get_r_cut(handle, &r_cut); ASSERT_FCS(res); res = fcs_p3m_get_grid(handle, grid_dim); ASSERT_FCS(res); if (0==myid) printf("FCS: Tuned grid dimensions, cutoff: %d %d %d, %f\n", grid_dim[0], grid_dim[1], grid_dim[2], r_cut); break; #endif #ifdef FCS_ENABLE_P2NFFT case FCS_METH_P2NFFT: res = fcs_p2nfft_get_grid(handle, grid_dim, grid_dim+1, grid_dim+2); ASSERT_FCS(res); res = fcs_p2nfft_get_r_cut(handle, &r_cut); ASSERT_FCS(res); if (0==myid) printf("FCS: Tuned grid dimensions, cutoff: %d %d %d, %f\n", grid_dim[0], grid_dim[1], grid_dim[2], r_cut); break; #endif #ifdef FCS_ENABLE_PP3MG case FCS_METH_PP3MG: res = fcs_pp3mg_get_cells_x(handle, grid_dim ); ASSERT_FCS(res); res = fcs_pp3mg_get_cells_y(handle, grid_dim+1); ASSERT_FCS(res); res = fcs_pp3mg_get_cells_z(handle, grid_dim+2); if (0==myid) printf("FCS: Tuned grid dimensions: %d %d %d\n", grid_dim[0], grid_dim[1], grid_dim[2]); ASSERT_FCS(res); break; #endif default: break; } /* add near-field potential, after fcs_tune */ if (0==srf) fcs_update_pottab(); }
void unpack_fcs(void) { fcs_float vir[9] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 }; FCSResult result; real pot1, pot2, e, c, sum=0.0, fac=0.5; int n, m, k, i; /* extract output and distribute it to cell array */ n=0; m=0; pot1=0.0; for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { c = CHARGE(p,i) * coul_eng; KRAFT(p,i,X) += field[n++] * c; KRAFT(p,i,Y) += field[n++] * c; KRAFT(p,i,Z) += field[n++] * c; e = pot[m++] * c * fac; POTENG(p,i) += e; pot1 += e; } } /* unpack virial */ result = fcs_get_virial(handle, vir); ASSERT_FCS(result); #ifdef P_AXIAL vir_xx += vir[0]; vir_yy += vir[4]; vir_zz += vir[8]; #else virial += vir[0] + vir[4] + vir[8]; #endif #ifdef STRESS_TENS if (do_press_calc) { /* distribute virial tensor evenly on all atoms */ sym_tensor pp; pp.xx = vir[0] / natoms; pp.yy = vir[4] / natoms; pp.zz = vir[8] / natoms; pp.yz = (vir[5]+vir[7]) / (2*natoms); pp.zx = (vir[2]+vir[6]) / (2*natoms); pp.xy = (vir[1]+vir[3]) / (2*natoms); for (k=0; k<NCELLS; ++k) { cell *p = CELLPTR(k); for (i=0; i<p->n; ++i) { PRESSTENS(p,i,xx) += pp.xx; PRESSTENS(p,i,yy) += pp.yy; PRESSTENS(p,i,zz) += pp.zz; PRESSTENS(p,i,yz) += pp.yz; PRESSTENS(p,i,zx) += pp.zx; PRESSTENS(p,i,xy) += pp.xy; } } } #endif /* sum up potential energy */ #ifdef MPI MPI_Allreduce( &pot1, &pot2, 1, MPI_DOUBLE, MPI_SUM, cpugrid); tot_pot_energy += pot2; #else tot_pot_energy += pot1; #endif }
int main (int argc, char **argv) { char input_file_name[300], *conf_file_name; FILE *input_file, *conf_file; int i, j, help, p_local = 0; FCSResult fcs_result; FCS fcs_handle; char parameterstring[200]; int particles_i; fcs_int particles = 0; fcs_int local_particles = 0; fcs_float *pos = NULL; fcs_float *charges = NULL; fcs_float *field = NULL; fcs_float *potentials = NULL; int cells_x = 128, cells_y = 128, cells_z = 128; int periodic = 1; int ghosts = 4; int degree = 4; int max_particles_i; fcs_int max_particles = 2000000; int max_iterations; double err_bound = 1.e-3; int nu1, nu2; fcs_float box[DIM][DIM] = { {1.0, 0.0, 0.0}, {0.0, 1.0, 0.0}, {0.0, 0.0, 1.0}}; fcs_float offset[DIM] = {0.0, 0.0, 0.0}; fcs_int periodic_flags[DIM] = {1,1,1}; /* Variables for analyzing results */ int read_ret_value = 0; double f_sum_local[DIM]; double f_sum[DIM]; double f_sum_squared_local[DIM]; double f_sum_squared[DIM]; double f_max_abs_local[DIM]; double f_max_abs[DIM]; double e_sum_local = 0.0; double e_sum = 0.0; int min = 0; /* Variables for reading data */ double my_x, my_y, my_z, my_q; /* Size of local domain */ double x_start, y_start, z_start; double x_end, y_end, z_end; int my_rank; int mpi_size; MPI_Comm mpi_comm_cart; int mpi_dims_i[DIM]; fcs_int mpi_dims[DIM]; int mpi_periods[DIM]; int mpi_coords[DIM]; MPI_Status status; double starttime = 0.0, endtime = 0.0; MPI_Init(&argc, &argv); MPI_Comm_size(MPI_COMM_WORLD, &mpi_size); MPI_Comm_rank(MPI_COMM_WORLD, &my_rank); if (my_rank == 0) { fprintf(stderr, "----------------\n"); fprintf(stderr, "Running pp3mg test\n"); fprintf(stderr, "----------------\n"); fprintf(stderr, "Setting up MPI...\n"); fprintf(stderr, " Using %d tasks.\n", mpi_size); } if(argc == 1) { printf("No config file was specified!\n"); MPI_Finalize(); exit(-1); } /* read config and input files */ if(argc >= 2) { conf_file_name = argv[1]; read_ret_value = read_config_file(conf_file_name, input_file_name, &cells_x, &cells_y, &cells_z, &max_iterations, &max_particles_i, &nu1, &nu2, &ghosts, &err_bound); if(read_ret_value) { printf("Config file couldn't be read!\n"); MPI_Finalize(); exit(-1); } max_particles = max_particles_i; } if((input_file = fopen(input_file_name, "r")) == NULL) { printf("input file name: %s\n" , input_file_name); printf("No input file was found!\n"); MPI_Finalize(); exit(-1); } fscanf(input_file, " %d", &particles_i); particles = particles_i; /* create a cartesian communicator */ for(i = 0; i < DIM; ++i) { mpi_periods[i] = 1; mpi_dims_i[i] = 0; } assert(MPI_Dims_create(mpi_size, 3, mpi_dims_i) == MPI_SUCCESS); assert(MPI_Cart_create(MPI_COMM_WORLD, 3, mpi_dims_i, mpi_periods, 1, &mpi_comm_cart) == MPI_SUCCESS); assert(MPI_Comm_rank(mpi_comm_cart, &my_rank) == MPI_SUCCESS); assert(MPI_Cart_coords(mpi_comm_cart, my_rank, 3, mpi_coords) == MPI_SUCCESS); for (i = 0; i < DIM; ++i) mpi_dims[i] = mpi_dims_i[i]; if (mpi_comm_cart != MPI_COMM_NULL) { /* initialize parameters for particle grid */ x_start = ((double)1.)/mpi_dims[0] * (double)mpi_coords[0]; x_end = ((double)1.)/mpi_dims[0] * (double)(mpi_coords[0]+1.); y_start = ((double)1.)/mpi_dims[1] * (double)mpi_coords[1]; y_end = ((double)1.)/mpi_dims[1] * (double)(mpi_coords[1]+1.); z_start = ((double)1.)/mpi_dims[2] * (double)mpi_coords[2]; z_end = ((double)1.)/mpi_dims[2] * (double)(mpi_coords[2]+1.); local_particles = 0; /* count local partciles , domain decomposition */ for(i = 0; i < particles; ++i) { fscanf(input_file, " %lf", &my_x); fscanf(input_file, " %lf", &my_y); fscanf(input_file, " %lf", &my_z); fscanf(input_file, " %lf", &my_q); if(x_start <= my_x && my_x < x_end && y_start <= my_y && my_y < y_end && z_start <= my_z && my_z < z_end) local_particles ++; } } fclose(input_file); if(local_particles != 0) { pos = malloc(local_particles*DIM*sizeof(fcs_float)); assert(pos); charges = malloc(local_particles*sizeof(fcs_float)); assert(charges); field = malloc(local_particles*DIM*sizeof(fcs_float)); assert(field); potentials = malloc(local_particles*sizeof(fcs_float)); assert(potentials); /* open the input file again and read particle data */ if((input_file = fopen(input_file_name, "r")) == NULL) { printf("No input file found!\n"); MPI_Finalize(); exit(-1); } fscanf(input_file, " %d", &particles); p_local = 0; for(i = 0; i < particles; ++i ) { fscanf(input_file, " %lf", &my_x); fscanf(input_file, " %lf", &my_y); fscanf(input_file, " %lf", &my_z); fscanf(input_file, " %lf", &my_q); if(x_start <= my_x && my_x < x_end && y_start <= my_y && my_y < y_end && z_start <= my_z && my_z < z_end) { pos[p_local] = my_x; pos[p_local + 1] = my_y; pos[p_local + 2] = my_z; charges[p_local/3] = my_q; p_local += 3; } } } /* create FCSParameter object. IMPORTANT: pp3mg requires a cartesian MPI communicator which has to be created by the calling program */ fcs_result = fcs_init(&fcs_handle, "pp3mg", mpi_comm_cart); ASSERT_FCS(fcs_result); fcs_result = fcs_set_dimension (fcs_handle, DIM); ASSERT_FCS(fcs_result); fcs_result = fcs_set_common(fcs_handle, 1, box[0], box[1], box[2], offset, periodic_flags, particles); ASSERT_FCS(fcs_result); fcs_result = fcs_pp3mg_setup (fcs_handle, mpi_dims, cells_x, cells_y, cells_z, ghosts, max_iterations, max_particles, periodic_flags[0] /* FIXME */, degree, err_bound); ASSERT_FCS(fcs_result); fcs_result = fcs_tune(fcs_handle, local_particles, pos, charges); ASSERT_FCS(fcs_result); MPI_Barrier(mpi_comm_cart); if(my_rank == 0) { starttime = MPI_Wtime(); } /* 2. step: run pp3mg. IMPORTANT: domain decomposition, particles must be distributed */ /* according to the specified MPI communicator */ fcs_result = fcs_run(fcs_handle, local_particles, pos, charges, field, potentials); ASSERT_FCS(fcs_result); /* analyze results */ /* field = fcsOutput_getField(fcs_output); potentials = fcsOutput_getPotentials(fcs_output); */ MPI_Barrier(mpi_comm_cart); if(my_rank == 0) { endtime = MPI_Wtime(); printf("pp3mg runtime with %d processors: %f s\n", mpi_size, endtime - starttime); } e_sum_local = 0.0; for(i = 0; i < local_particles; ++i) e_sum_local += potentials[i]; e_sum = 0.0; MPI_Reduce(&e_sum_local, &e_sum, 1, MPI_DOUBLE, MPI_SUM, 0, mpi_comm_cart); if(my_rank == 0) { min = cells_x; if(cells_y < min || cells_z < min) min = (cells_y < cells_z) ? cells_y : cells_z; printf("Self energy: %e\n", 1.0/(4.0 * PI) * 14.0/(5.0*1./(2.*ghosts*min))); printf("Approx. Madelung's constant: %e\n", e_sum/particles *2.0 * PI); printf( "Total energy:: %e\n\n", e_sum); } for(i = 0; i < DIM; ++i) { f_sum_local[i] = 0.0; f_sum[i] = 0.0; f_sum_squared_local[i] = 0.0; f_sum_squared[i] = 0.0; f_max_abs_local[i] = 0.0; f_max_abs[i] = 0.0; } for(i = 0; i < local_particles; ++i) { for(j = 0; j < DIM; ++j) { f_sum_local[j] += charges[i] * field[i*DIM + j]; f_sum_squared_local[j] += charges[i] * charges[i] * field[i*DIM + j]*field[i*DIM + j]; if(fabs(field[i*DIM + j]) > f_max_abs_local[j]) { f_max_abs_local[j] = fabs(charges[i] * field[i*DIM + j]); }/* end if */ }/* end for-j */ }/* end for-i */ for(j = 0; j < DIM; ++j) { f_max_abs[j] = f_max_abs_local[j]; } MPI_Reduce(f_sum_local, f_sum, 3, MPI_DOUBLE, MPI_SUM, 0, mpi_comm_cart); MPI_Reduce(f_sum_squared_local, f_sum_squared, 3, MPI_DOUBLE, MPI_SUM, 0, mpi_comm_cart); MPI_Reduce(f_max_abs_local, f_max_abs, 3, MPI_DOUBLE, MPI_MAX, 0, mpi_comm_cart); if(my_rank == 0) { printf("Norm of sum of forces: %e\n", sqrt(f_sum[0]*f_sum[0] + f_sum[1]*f_sum[1] + f_sum[2]*f_sum[2]) ); printf("Sqrt. of sum of squares of all components of forces: %f\n", sqrt(f_sum_squared[0] + f_sum_squared[1] + f_sum_squared[2])); printf("Maximal absolute force components: \n"); for(i = 0; i < DIM; ++i) { printf("%f ", f_max_abs[i]); } printf("\n"); } if(local_particles != 0) fclose(input_file); /* 3. step: deallocate resources for pp3mg */ fcs_destroy(fcs_handle); free (pos); free (charges); free (field); free (potentials); MPI_Comm_free(&mpi_comm_cart); MPI_Finalize(); exit(0); }