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(); }
int scafacos_p3m_tuning(){ int mesh[3] = {0, 0, 0}, tmp_mesh_points; int tmp_mesh[3]; double r_cut_iL_min, r_cut_iL_max, r_cut_iL = -1, tmp_r_cut_iL=0.0; int cao; double alpha_L = -1, tmp_alpha_L=0.0; double accuracy = -1, tmp_accuracy=0.0; double time_best=1e20, tmp_time; char b[3*ES_INTEGER_SPACE + 3*ES_DOUBLE_SPACE + 128]; int buffer_cao, buffer_grid; double buffer_alpha_L, buffer_r_cut_iL, buffer_tolerance_field; if (skin == -1) { // *log = strcat_alloc(*log, "p3m cannot be tuned, since the skin is not yet set"); return ES_ERROR; } /* preparation */ mpi_bcast_event(FCS_P3M_COUNT_CHARGES); /* Print Status */ sprintf(b, "P3M tune parameters: Accuracy goal = %.5e\n", scafacos_p3m.params.tolerance_field); // *log = strcat_alloc(*log, b); sprintf(b, "System: box_l = %.5e # charged part = %d Sum[q_i^2] = %.5e\n", box_l[0], scafacos_p3m.sum_qpart, scafacos_p3m.sum_q2); // *log = strcat_alloc(*log, b); if (scafacos_p3m.sum_qpart == 0) { // *log = strcat_alloc(*log, "no charged particles in the system, cannot tune P3M"); return ES_ERROR; } double density, r_cut; density = n_total_particles /(box_l[0]* box_l[1]* box_l[2]); int number = 100; r_cut= pow(number/ density/3.14 *3 /4 , 1/3.0); if(scafacos_p3m.params.r_cut_iL == 0.0) { // WARNING: r_cut_iL_min should not be small ! r_cut_iL_min = dmin(min_local_box_l, min_box_l/2); r_cut_iL_min *= 0.18; r_cut_iL_min = r_cut *0.8; r_cut_iL_max = dmin(min_local_box_l, min_box_l/2); r_cut_iL_max *= 0.28; r_cut_iL_max = r_cut *1.5; if(r_cut_iL_max >= 4.5) r_cut_iL_max = 4.5; // r_cut_iL_max = dmin(min_local_box_l, min_box_l/2);// - skin; // r_cut_iL_min *= box_l_i[0]; // r_cut_iL_max *= box_l_i[0]; fprintf(stderr, "r_cut_iL_max: %f r_cut_iL_min: %f \n", r_cut_iL_max, r_cut_iL_min); fprintf(stderr, "box_l is: %f \n", box_l[0]); } else{ fprintf(stderr, "Scafacos_p3m tuning varies the cutoff radius; you have to leave cutoff unspecified"); return ES_ERROR; } // *log = strcat_alloc(*log, "mesh cao r_cut_iL alpha_L err rs_err ks_err time [ms]\n"); for(tmp_r_cut_iL = r_cut_iL_min; tmp_r_cut_iL <= r_cut_iL_max; tmp_r_cut_iL += 0.1){ scafacos_p3m.params.cutoff = tmp_r_cut_iL; tmp_time = time_force_calc(1); // *log = strcat_alloc(*log, " %f ms "); fprintf(stderr, "time: %f \n", tmp_time); // *log = strcat_alloc(*log, b); if ((tmp_time < time_best) && (tmp_time > 0)) { time_best = tmp_time; fcs_p3m_get_grid(fcs_handle, &buffer_grid); fcs_p3m_get_cao(fcs_handle, &buffer_cao); fcs_p3m_get_r_cut(fcs_handle, &buffer_r_cut_iL); fcs_p3m_get_alpha(fcs_handle, &buffer_alpha_L); fcs_p3m_get_tolerance_field(fcs_handle, &buffer_tolerance_field); // target accuracy remains the same } /* no hope of further optimisation */ //else if (tmp_time > time_best + P3M_TIME_GRAN) { // break; //} } if(time_best == 1e20) { fprintf(stderr, "failed to tune P3M parameters to required accuracy\n"); return ES_ERROR; } /* prepare tuned p3m parameters for broadcast*/ scafacos_p3m.params.cao = buffer_cao; scafacos_p3m.params.alpha_L = buffer_alpha_L; scafacos_p3m.params.grid = buffer_grid; scafacos_p3m.params.cutoff = buffer_r_cut_iL; scafacos_p3m.params.alpha_L = buffer_alpha_L; /* Tell the user about the outcome */ fprintf(stderr, "\nresulting parameters:\n%-4d %-3d %.5e %.5e %.5e %-8d\n", buffer_grid, buffer_cao, buffer_r_cut_iL, buffer_alpha_L, buffer_tolerance_field, (int)time_best); return ES_OK; }