diff --git a/Quda/.clang-format b/Quda/.clang-format new file mode 100644 index 0000000..9d54a25 --- /dev/null +++ b/Quda/.clang-format @@ -0,0 +1,14 @@ +{ + BasedOnStyle: LLVM, + UseTab: Never, + IndentWidth: 2, + TabWidth: 2, + BreakBeforeBraces: Allman, + AllowShortIfStatementsOnASingleLine: false, + IndentCaseLabels: false, + ColumnLimit: 90, + AccessModifierOffset: -4, + NamespaceIndentation: All, + FixNamespaceComments: false, + SortIncludes: true, +} diff --git a/Quda/Benchmark_Quda.cpp b/Quda/Benchmark_Quda.cpp new file mode 100644 index 0000000..4ce5e7b --- /dev/null +++ b/Quda/Benchmark_Quda.cpp @@ -0,0 +1,169 @@ +#include +#include +#include +#include +#include +// #include +#include +#include +#include + +#include +#include +#include + +using namespace quda; + +QudaPrecision smoother_halo_prec = QUDA_INVALID_PRECISION; + +std::array gridsize = {1, 1, 1, 4}; + +void initComms(int argc, char **argv, std::array const &commDims) +{ + // init MPI communication + MPI_Init(&argc, &argv); + + // this maps coordinates to rank number + auto lex_rank_from_coords = [](int const *coords, void *) + { + int rank = coords[0]; + for (int i = 1; i < 4; i++) + rank = gridsize[i] * rank + coords[i]; + return rank; + }; + initCommsGridQuda(4, commDims.data(), lex_rank_from_coords, nullptr); + + for (int d = 0; d < 4; d++) + if (gridsize[d] > 1) + commDimPartitionedSet(d); +} + +// creates a random gauge field +cudaGaugeField make_gauge_field(std::array const &geom) +{ + GaugeFieldParam param; + + // dimension and type of the lattice object + param.nDim = 4; + param.nColor = 3; + param.x[0] = geom[0]; + param.x[1] = geom[1]; + param.x[2] = geom[2]; + param.x[3] = geom[3]; + param.t_boundary = QUDA_PERIODIC_T; + param.siteSubset = QUDA_FULL_SITE_SUBSET; // no even/odd, just a full lattice + param.link_type = QUDA_SU3_LINKS; + param.setPrecision(QUDA_DOUBLE_PRECISION); + + param.create = QUDA_NULL_FIELD_CREATE; // do not (zero-) initilize the fields + param.location = QUDA_CUDA_FIELD_LOCATION; // field should live on the accelerator + + // turn off advanced features we dont care about for this benchmark + param.reconstruct = QUDA_RECONSTRUCT_NO; + param.ghostExchange = QUDA_GHOST_EXCHANGE_NO; + + // these control the physical data layout. Might be interesting to try out different + // settings + param.order = QUDA_FLOAT2_GAUGE_ORDER; + param.geometry = QUDA_SCALAR_GEOMETRY; + + // create the field and fill with random SU(3) matrices + auto U = cudaGaugeField(param); + quda::RNG rng(U, /*seed=*/1234); + gaugeGauss(U, rng, 1.0); + return U; +} + +// create a random source vector +ColorSpinorField make_source(std::array const &geom) +{ + ColorSpinorParam param; + param.nColor = 3; + param.nSpin = 4; + param.nVec = 1; // only a single vector + param.pad = 0; + param.siteSubset = QUDA_FULL_SITE_SUBSET; + param.nDim = 4; + param.x[0] = geom[0]; + param.x[1] = geom[1]; + param.x[2] = geom[2]; + param.x[3] = geom[3]; + param.x[4] = 1; // no fifth dimension + param.pc_type = QUDA_4D_PC; + param.siteOrder = QUDA_EVEN_ODD_SITE_ORDER; + param.gammaBasis = QUDA_DEGRAND_ROSSI_GAMMA_BASIS; + param.create = QUDA_NULL_FIELD_CREATE; // do not (zero-) initilize the field + param.setPrecision(QUDA_DOUBLE_PRECISION); + param.location = QUDA_CUDA_FIELD_LOCATION; + + // create the field and fill it with random values + auto src = ColorSpinorField(param); + quda::RNG rng(src, 1234); + spinorNoise(src, rng, QUDA_NOISE_GAUSS); + printfQuda( + "created src with norm = %f (sanity check: should be close to %f) and %f bytes\n", + blas::norm2(src), 2.0 * 12 * geom[0] * geom[1] * geom[2] * geom[3], + src.Bytes() * 1.0); + src.PrintDims(); + + return src; +} + +void benchmark(int L, int niter) +{ + std::array geom = {L, L, L, L}; + + auto U = make_gauge_field(geom); + auto src = make_source(geom); + + // create (Wilson) dirac operator + DiracParam param; + param.kappa = 0.10; + param.dagger = QUDA_DAG_NO; + param.matpcType = QUDA_MATPC_EVEN_EVEN; + auto dirac = DiracWilson(param); + + // insert gauge field into the dirac operator + // (the additional nullptr's are for smeared links and fancy preconditioners and such. + // Not used for simple Wilson fermions) + dirac.updateFields(&U, nullptr, nullptr, nullptr); + + auto tmp = ColorSpinorField(ColorSpinorParam(src)); + + printfQuda("benchmarking Dirac operator. geom=(%d,%d,%d,%d), niter=%d\n", geom[0], + geom[1], geom[2], geom[3], niter); + + // couple iterations without timing to warm up + for (int iter = 0; iter < 20; ++iter) + dirac.M(tmp, src); + + dirac.Flops(); // reset flops counter + device_timer_t device_timer; + device_timer.start(); + for (int iter = 0; iter < niter; ++iter) + dirac.M(tmp, src); + device_timer.stop(); + + double secs = device_timer.last(); + double gflops = (dirac.Flops() * 1e-9) / secs; + printfQuda("Gflops = %6.1f\n", gflops); +} + +int main(int argc, char **argv) +{ + initComms(argc, argv, gridsize); + + // -1 for multi-gpu. otherwise this selects the device to be used + initQuda(-1); + + // verbosity options are: + // SILENT, SUMMARIZE, VERBOSE, DEBUG_VERBOSE + setVerbosity(QUDA_SUMMARIZE); + + for (int L : {8, 16, 24, 32}) + benchmark(L, 1000); + + endQuda(); + quda::comm_finalize(); + MPI_Finalize(); +} diff --git a/Quda/build.sh b/Quda/build.sh new file mode 100755 index 0000000..ea3221f --- /dev/null +++ b/Quda/build.sh @@ -0,0 +1,11 @@ +#!/bin/bash + +set -e + + +COMPILE_FLAGS="-DMPI_COMMS -DMULTI_GPU -DQUDA_PRECISION=14 -DQUDA_RECONSTRUCT=7 -g -O3 -Wall -Wextra -pthread -std=c++17" +LINK_FLAGS="-g -O3 -Wl,-rpath -Wl,/mnt/lustre/tursafs1/home/dp207/dp207/shared/env/versions/220428/prefix/ompi_gpu/lib -Wl,--enable-new-dtags -L/mnt/lustre/tursafs1/home/dp207/dp207/shared/env/versions/220428/prefix/ompi_gpu/lib -pthread -Wl,-rpath,/home/dp207/dp207/dc-burg2/quda_build/tests:/home/dp207/dp207/dc-burg2/quda_build/lib:/home/dp207/dp207/shared/env/versions/220428/spack/opt/spack/linux-rhel8-zen2/gcc-9.4.0/cuda-11.4.0-etxow4jb23qdbs7j6txczy44cdatpj22/lib64/stubs:/home/dp207/dp207/shared/env/versions/220428/spack/opt/spack/linux-rhel8-zen2/gcc-9.4.0/cuda-11.4.0-etxow4jb23qdbs7j6txczy44cdatpj22/lib64: ../../quda_install/lib/libquda.so /home/dp207/dp207/shared/env/versions/220428/spack/opt/spack/linux-rhel8-zen2/gcc-9.4.0/cuda-11.4.0-etxow4jb23qdbs7j6txczy44cdatpj22/lib64/stubs/libcuda.so /home/dp207/dp207/shared/env/versions/220428/spack/opt/spack/linux-rhel8-zen2/gcc-9.4.0/cuda-11.4.0-etxow4jb23qdbs7j6txczy44cdatpj22/lib64/stubs/libnvidia-ml.so /home/dp207/dp207/shared/env/versions/220428/spack/opt/spack/linux-rhel8-zen2/gcc-9.4.0/cuda-11.4.0-etxow4jb23qdbs7j6txczy44cdatpj22/lib64/libcudart_static.a -ldl /usr/lib64/librt.so /home/dp207/dp207/shared/env/versions/220428/spack/opt/spack/linux-rhel8-zen2/gcc-9.4.0/cuda-11.4.0-etxow4jb23qdbs7j6txczy44cdatpj22/lib64/libcublas.so /home/dp207/dp207/shared/env/versions/220428/spack/opt/spack/linux-rhel8-zen2/gcc-9.4.0/cuda-11.4.0-etxow4jb23qdbs7j6txczy44cdatpj22/lib64/libcufft.so -lpthread /home/dp207/dp207/shared/env/versions/220428/prefix/ompi_gpu/lib/libmpi.so" + +$CXX $COMPILE_FLAGS -I/home/dp207/dp207/dc-burg2/quda_install/include/ -o Benchmark_Quda.o -c Benchmark_Quda.cpp +$CXX $LINK_FLAGS Benchmark_Quda.o -o Benchmark_Quda +