/* * ======================================================================================= * * Author: Jan Eitzinger (je), jan.eitzinger@fau.de * Copyright (c) 2020 RRZE, University Erlangen-Nuremberg * * This file is part of MD-Bench. * * MD-Bench is free software: you can redistribute it and/or modify it * under the terms of the GNU Lesser General Public License as published * by the Free Software Foundation, either version 3 of the License, or * (at your option) any later version. * * MD-Bench is distributed in the hope that it will be useful, but WITHOUT ANY * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A * PARTICULAR PURPOSE. See the GNU Lesser General Public License for more * details. * * You should have received a copy of the GNU Lesser General Public License along * with MD-Bench. If not, see . * ======================================================================================= */ #include #include //-- #include //-- #include #include #include #include #include #include #include #include #include #include #include #include #include #define HLINE "----------------------------------------------------------------------------\n" extern double computeForceLJ_ref(Parameter*, Atom*, Neighbor*, Stats*); extern double computeForceLJ_4xn(Parameter*, Atom*, Neighbor*, Stats*); extern double computeForceLJ_2xnn(Parameter*, Atom*, Neighbor*, Stats*); extern double computeForceEam(Eam*, Parameter*, Atom*, Neighbor*, Stats*); double setup(Parameter *param, Eam *eam, Atom *atom, Neighbor *neighbor, Stats *stats) { if(param->force_field == FF_EAM) { initEam(eam, param); } double S, E; param->lattice = pow((4.0 / param->rho), (1.0 / 3.0)); param->xprd = param->nx * param->lattice; param->yprd = param->ny * param->lattice; param->zprd = param->nz * param->lattice; S = getTimeStamp(); initAtom(atom); initPbc(atom); initStats(stats); initNeighbor(neighbor, param); if(param->input_file == NULL) { createAtom(atom, param); } else { readAtom(atom, param); } setupNeighbor(param, atom); setupThermo(param, atom->Natoms); if(param->input_file == NULL) { adjustThermo(param, atom); } buildClusters(atom); defineJClusters(atom); setupPbc(atom, param); binClusters(atom); buildNeighbor(atom, neighbor); E = getTimeStamp(); return E-S; } double reneighbour(Parameter *param, Atom *atom, Neighbor *neighbor) { double S, E; S = getTimeStamp(); LIKWID_MARKER_START("reneighbour"); updateSingleAtoms(atom); updateAtomsPbc(atom, param); buildClusters(atom); defineJClusters(atom); setupPbc(atom, param); binClusters(atom); buildNeighbor(atom, neighbor); LIKWID_MARKER_STOP("reneighbour"); E = getTimeStamp(); return E-S; } void initialIntegrate(Parameter *param, Atom *atom) { DEBUG_MESSAGE("initialIntegrate start\n"); for(int ci = 0; ci < atom->Nclusters_local; ci++) { int ci_vec_base = CI_VECTOR_BASE_INDEX(ci); MD_FLOAT *ci_x = &atom->cl_x[ci_vec_base]; MD_FLOAT *ci_v = &atom->cl_v[ci_vec_base]; MD_FLOAT *ci_f = &atom->cl_f[ci_vec_base]; for(int cii = 0; cii < atom->iclusters[ci].natoms; cii++) { ci_v[CL_X_OFFSET + cii] += param->dtforce * ci_f[CL_X_OFFSET + cii]; ci_v[CL_Y_OFFSET + cii] += param->dtforce * ci_f[CL_Y_OFFSET + cii]; ci_v[CL_Z_OFFSET + cii] += param->dtforce * ci_f[CL_Z_OFFSET + cii]; ci_x[CL_X_OFFSET + cii] += param->dt * ci_v[CL_X_OFFSET + cii]; ci_x[CL_Y_OFFSET + cii] += param->dt * ci_v[CL_Y_OFFSET + cii]; ci_x[CL_Z_OFFSET + cii] += param->dt * ci_v[CL_Z_OFFSET + cii]; } } DEBUG_MESSAGE("initialIntegrate end\n"); } void finalIntegrate(Parameter *param, Atom *atom) { DEBUG_MESSAGE("finalIntegrate start\n"); for(int ci = 0; ci < atom->Nclusters_local; ci++) { int ci_vec_base = CI_VECTOR_BASE_INDEX(ci); MD_FLOAT *ci_v = &atom->cl_v[ci_vec_base]; MD_FLOAT *ci_f = &atom->cl_f[ci_vec_base]; for(int cii = 0; cii < atom->iclusters[ci].natoms; cii++) { ci_v[CL_X_OFFSET + cii] += param->dtforce * ci_f[CL_X_OFFSET + cii]; ci_v[CL_Y_OFFSET + cii] += param->dtforce * ci_f[CL_Y_OFFSET + cii]; ci_v[CL_Z_OFFSET + cii] += param->dtforce * ci_f[CL_Z_OFFSET + cii]; } } DEBUG_MESSAGE("finalIntegrate end\n"); } void printAtomState(Atom *atom) { printf("Atom counts: Natoms=%d Nlocal=%d Nghost=%d Nmax=%d\n", atom->Natoms, atom->Nlocal, atom->Nghost, atom->Nmax); /* int nall = atom->Nlocal + atom->Nghost; */ /* for (int i=0; ix[i], atom->y[i], atom->z[i]); */ /* } */ } int main(int argc, char** argv) { double timer[NUMTIMER]; Eam eam; Atom atom; Neighbor neighbor; Stats stats; Parameter param; LIKWID_MARKER_INIT; #pragma omp parallel { LIKWID_MARKER_REGISTER("force"); //LIKWID_MARKER_REGISTER("reneighbour"); //LIKWID_MARKER_REGISTER("pbc"); } initParameter(¶m); for(int i = 0; i < argc; i++) { if((strcmp(argv[i], "-p") == 0)) { readParameter(¶m, argv[++i]); continue; } if((strcmp(argv[i], "-f") == 0)) { if((param.force_field = str2ff(argv[++i])) < 0) { fprintf(stderr, "Invalid force field!\n"); exit(-1); } continue; } if((strcmp(argv[i], "-i") == 0)) { param.input_file = strdup(argv[++i]); continue; } if((strcmp(argv[i], "-e") == 0)) { param.eam_file = strdup(argv[++i]); continue; } if((strcmp(argv[i], "-n") == 0) || (strcmp(argv[i], "--nsteps") == 0)) { param.ntimes = atoi(argv[++i]); continue; } if((strcmp(argv[i], "-nx") == 0)) { param.nx = atoi(argv[++i]); continue; } if((strcmp(argv[i], "-ny") == 0)) { param.ny = atoi(argv[++i]); continue; } if((strcmp(argv[i], "-nz") == 0)) { param.nz = atoi(argv[++i]); continue; } if((strcmp(argv[i], "-half") == 0)) { param.half_neigh = atoi(argv[++i]); continue; } if((strcmp(argv[i], "-m") == 0) || (strcmp(argv[i], "--mass") == 0)) { param.mass = atof(argv[++i]); continue; } if((strcmp(argv[i], "-r") == 0) || (strcmp(argv[i], "--radius") == 0)) { param.cutforce = atof(argv[++i]); continue; } if((strcmp(argv[i], "-s") == 0) || (strcmp(argv[i], "--skin") == 0)) { param.skin = atof(argv[++i]); continue; } if((strcmp(argv[i], "--freq") == 0)) { param.proc_freq = atof(argv[++i]); continue; } if((strcmp(argv[i], "--vtk") == 0)) { param.vtk_file = strdup(argv[++i]); continue; } if((strcmp(argv[i], "--xtc") == 0)) { #ifndef XTC_OUTPUT fprintf(stderr, "XTC not available, set XTC_OUTPUT option in config.mk file and recompile MD-Bench!"); exit(-1); #else param.xtc_file = strdup(argv[++i]); #endif continue; } if((strcmp(argv[i], "-h") == 0) || (strcmp(argv[i], "--help") == 0)) { printf("MD Bench: A minimalistic re-implementation of miniMD\n"); printf(HLINE); printf("-p : file to read parameters from (can be specified more than once)\n"); printf("-f : force field (lj or eam), default lj\n"); printf("-i : input file with atom positions (dump)\n"); printf("-e : input file for EAM\n"); printf("-n / --nsteps : set number of timesteps for simulation\n"); printf("-nx/-ny/-nz : set linear dimension of systembox in x/y/z direction\n"); printf("-r / --radius : set cutoff radius\n"); printf("-s / --skin : set skin (verlet buffer)\n"); printf("--freq : processor frequency (GHz)\n"); printf("--vtk : VTK file for visualization\n"); printf("--xtc : XTC file for visualization\n"); printf(HLINE); exit(EXIT_SUCCESS); } } param.cutneigh = param.cutforce + param.skin; setup(¶m, &eam, &atom, &neighbor, &stats); printParameter(¶m); printf("step\ttemp\t\tpressure\n"); computeThermo(0, ¶m, &atom); #if defined(MEM_TRACER) || defined(INDEX_TRACER) traceAddresses(¶m, &atom, &neighbor, n + 1); #endif if(param.force_field == FF_EAM) { timer[FORCE] = computeForceEam(&eam, ¶m, &atom, &neighbor, &stats); } else { timer[FORCE] = computeForceLJ(¶m, &atom, &neighbor, &stats); } timer[NEIGH] = 0.0; timer[TOTAL] = getTimeStamp(); if(param.vtk_file != NULL) { write_data_to_vtk_file(param.vtk_file, &atom, 0); } if(param.xtc_file != NULL) { xtc_init(param.xtc_file, &atom, 0); } for(int n = 0; n < param.ntimes; n++) { initialIntegrate(¶m, &atom); if((n + 1) % param.reneigh_every) { if(!((n + 1) % param.prune_every)) { pruneNeighbor(¶m, &atom, &neighbor); } updatePbc(&atom, ¶m, 0); } else { timer[NEIGH] += reneighbour(¶m, &atom, &neighbor); } #if defined(MEM_TRACER) || defined(INDEX_TRACER) traceAddresses(¶m, &atom, &neighbor, n + 1); #endif if(param.force_field == FF_EAM) { timer[FORCE] += computeForceEam(&eam, ¶m, &atom, &neighbor, &stats); } else { timer[FORCE] += computeForceLJ(¶m, &atom, &neighbor, &stats); } finalIntegrate(¶m, &atom); if(!((n + 1) % param.nstat) && (n+1) < param.ntimes) { computeThermo(n + 1, ¶m, &atom); } int write_pos = !((n + 1) % param.x_out_every); int write_vel = !((n + 1) % param.v_out_every); if(write_pos || write_vel) { if(param.vtk_file != NULL) { write_data_to_vtk_file(param.vtk_file, &atom, n + 1); } if(param.xtc_file != NULL) { xtc_write(&atom, n + 1, write_pos, write_vel); } } } timer[TOTAL] = getTimeStamp() - timer[TOTAL]; updateSingleAtoms(&atom); computeThermo(-1, ¶m, &atom); if(param.xtc_file != NULL) { xtc_end(); } printf(HLINE); printf("Kernel: %s, MxN: %dx%d, Vector width: %d\n", KERNEL_NAME, CLUSTER_M, CLUSTER_N, VECTOR_WIDTH); printf("Data layout for positions: %s\n", POS_DATA_LAYOUT); #if PRECISION == 1 printf("Using single precision floating point.\n"); #else printf("Using double precision floating point.\n"); #endif printf(HLINE); printf("System: %d atoms %d ghost atoms, Steps: %d\n", atom.Natoms, atom.Nghost, param.ntimes); printf("TOTAL %.2fs FORCE %.2fs NEIGH %.2fs REST %.2fs\n", timer[TOTAL], timer[FORCE], timer[NEIGH], timer[TOTAL]-timer[FORCE]-timer[NEIGH]); printf(HLINE); printf("Performance: %.2f million atom updates per second\n", 1e-6 * (double) atom.Natoms * param.ntimes / timer[TOTAL]); #ifdef COMPUTE_STATS displayStatistics(&atom, ¶m, &stats, timer); #endif LIKWID_MARKER_CLOSE; return EXIT_SUCCESS; }