lattice-benchmarks/Quda/Benchmark_Quda.cpp

448 lines
14 KiB
C++
Raw Normal View History

2023-03-31 18:03:39 +01:00
#include <algorithm>
#include <array>
#include <blas_quda.h>
2023-04-24 11:35:47 +01:00
#include <cassert>
2023-06-19 18:22:24 +01:00
#include <chrono>
2023-03-31 18:03:39 +01:00
#include <color_spinor_field.h>
2023-04-24 11:35:47 +01:00
#include <dirac_quda.h>
2023-06-19 18:08:50 +01:00
#include <fstream>
2023-04-24 11:35:47 +01:00
#include <gauge_tools.h>
2023-03-31 18:03:39 +01:00
#include <memory>
2023-04-24 11:35:47 +01:00
#include <mpi.h>
2023-03-31 18:03:39 +01:00
#include <stdio.h>
#include <stdlib.h>
2023-06-05 17:07:07 +01:00
// remove to use QUDA's own flop counting instead of Grid's convention
#define FLOP_COUNTING_GRID
2023-06-19 18:08:50 +01:00
#include "json.hpp"
using nlohmann::json;
json json_results;
using namespace quda;
2023-06-19 18:22:24 +01:00
// timestamp = seconds since program start.
// these are written to the json output with the goal of later matching them against
// power-measurments to determine energy efficiency.
using Clock = std::chrono::steady_clock;
Clock::time_point program_start_time = Clock::now();
double get_timestamp()
{
auto dur = Clock::now() - program_start_time;
return std::chrono::duration_cast<std::chrono::microseconds>(dur).count() * 1.0e-6;
}
2023-04-24 11:35:47 +01:00
// This is the MPI grid, i.e. the layout of ranks
int nranks = -1;
std::array<int, 4> mpi_grid = {1, 1, 1, 1};
2023-03-31 18:03:39 +01:00
2023-04-24 11:35:47 +01:00
void initComms(int argc, char **argv)
2023-03-31 18:03:39 +01:00
{
// init MPI communication
MPI_Init(&argc, &argv);
2023-04-24 11:35:47 +01:00
MPI_Comm_size(MPI_COMM_WORLD, &nranks);
assert(1 <= nranks && nranks <= 100000);
mpi_grid[3] = nranks;
2023-03-31 18:03:39 +01:00
// 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++)
2023-04-24 11:35:47 +01:00
rank = mpi_grid[i] * rank + coords[i];
2023-03-31 18:03:39 +01:00
return rank;
};
2023-04-21 10:38:28 +01:00
2023-04-24 11:35:47 +01:00
initCommsGridQuda(4, mpi_grid.data(), lex_rank_from_coords, nullptr);
2023-03-31 18:03:39 +01:00
for (int d = 0; d < 4; d++)
2023-04-24 11:35:47 +01:00
if (mpi_grid[d] > 1)
2023-03-31 18:03:39 +01:00
commDimPartitionedSet(d);
2023-06-19 18:08:50 +01:00
json_results["geometry"]["ranks"] = nranks;
json_results["geometry"]["mpi"] = mpi_grid;
2023-03-31 18:03:39 +01:00
}
2023-04-24 11:35:47 +01:00
// creates a random gauge field. L = local(!) size
cudaGaugeField make_gauge_field(int L)
2023-03-31 18:03:39 +01:00
{
GaugeFieldParam param;
// dimension and type of the lattice object
param.nDim = 4;
2023-04-24 11:35:47 +01:00
param.x[0] = L;
param.x[1] = L;
param.x[2] = L;
param.x[3] = L;
2023-04-21 10:38:28 +01:00
// number of colors. potentially confusingly, QUDA sometimes uses the word "color" to
// things unrelated with physical color. things like "nColor=32" do pop up in deflation
// solvers where it (to my understanding) refers to the number of (parallely processed)
// deflation vectors.
param.nColor = 3;
// boundary conditions (dont really care for benchmark)
2023-03-31 18:03:39 +01:00
param.t_boundary = QUDA_PERIODIC_T;
2023-04-21 10:38:28 +01:00
// for this benchmark we only need "SINGLE" and/or "DOUBLE" precision. But smaller
// precisions are available in QUDA too
2023-06-05 17:07:07 +01:00
param.setPrecision(QUDA_SINGLE_PRECISION);
2023-03-31 18:03:39 +01:00
2023-04-21 10:38:28 +01:00
// no even/odd subset, we want a full lattice
param.siteSubset = QUDA_FULL_SITE_SUBSET;
// what kind of 3x3 matrices the field contains. A proper gauge field has SU(3)
// matrices, but (for example) smeared/thick links could have non-unitary links.
param.link_type = QUDA_SU3_LINKS;
// "NULL" does not initialize the field upon creation, "ZERO" would set everything to 0
param.create = QUDA_NULL_FIELD_CREATE;
2023-03-31 18:03:39 +01:00
2023-04-21 10:38:28 +01:00
// field should be allocated directly on the accelerator/GPU
param.location = QUDA_CUDA_FIELD_LOCATION;
// "reconstruct" here means reconstructing a SU(3) matrix from fewer than 18 real
// numbers (=3x3 complex numbers). Great feature in production (saving
// memory/cache/network bandwidth), not used for this benchmark.
2023-03-31 18:03:39 +01:00
param.reconstruct = QUDA_RECONSTRUCT_NO;
2023-04-21 10:38:28 +01:00
// "ghostExchange" would often be called "halo exchange" outside of Quda. This has
// nothing to do with ghost fields from continuum/perturbative qcd.
2023-03-31 18:03:39 +01:00
param.ghostExchange = QUDA_GHOST_EXCHANGE_NO;
2023-04-21 10:38:28 +01:00
// This controls the physical order of elements. "float2" is the the default
2023-03-31 18:03:39 +01:00
param.order = QUDA_FLOAT2_GAUGE_ORDER;
2023-04-21 10:38:28 +01:00
// this means the field is a LORENTZ vector (which a gauge field must be). Has nothing
// to do with spin.
param.geometry = QUDA_VECTOR_GEOMETRY;
2023-03-31 18:03:39 +01:00
// create the field and fill with random SU(3) matrices
2023-04-21 10:38:28 +01:00
// std::cout << param << std::endl; // double-check parameters
2023-03-31 18:03:39 +01:00
auto U = cudaGaugeField(param);
2023-04-21 10:38:28 +01:00
gaugeGauss(U, /*seed=*/1234, 1.0);
2023-03-31 18:03:39 +01:00
return U;
}
2023-04-24 11:35:47 +01:00
// create a random source vector (L = local size)
2023-06-05 17:07:07 +01:00
ColorSpinorField make_source(int L, int Ls = 1)
2023-03-31 18:03:39 +01:00
{
// NOTE: `param.x` directly determines the size of the (local, per rank) memory
// allocation. Thus for checkerboarding, we have to specifly x=(L/2,L,L,L) to get a
// physical local volume of L^4, thus implicity choosing a dimension for the
// checkerboarding (shouldnt really matter of course which one).
2023-03-31 18:03:39 +01:00
ColorSpinorParam param;
param.nColor = 3;
param.nSpin = 4;
param.nVec = 1; // only a single vector
param.pad = 0;
param.siteSubset = QUDA_PARITY_SITE_SUBSET;
2023-06-05 17:07:07 +01:00
param.nDim = Ls == 1 ? 4 : 5;
param.x[0] = L / 2;
2023-04-24 11:35:47 +01:00
param.x[1] = L;
param.x[2] = L;
param.x[3] = L;
2023-06-05 17:07:07 +01:00
param.x[4] = Ls;
2023-03-31 18:03:39 +01:00
param.pc_type = QUDA_4D_PC;
param.siteOrder = QUDA_EVEN_ODD_SITE_ORDER;
// somewhat surprisingly, the DiracWilson::Dslash(...) function only works with the
// UKQCD_GAMMA_BASIS
param.gammaBasis = QUDA_UKQCD_GAMMA_BASIS;
2023-03-31 18:03:39 +01:00
param.create = QUDA_NULL_FIELD_CREATE; // do not (zero-) initilize the field
2023-06-05 17:07:07 +01:00
param.setPrecision(QUDA_SINGLE_PRECISION);
2023-03-31 18:03:39 +01:00
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);
2023-04-21 10:38:28 +01:00
/*printfQuda(
2023-03-31 18:03:39 +01:00
"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],
2023-04-21 10:38:28 +01:00
src.Bytes() * 1.0);*/
// src.PrintDims();
2023-03-31 18:03:39 +01:00
return src;
}
2023-06-19 18:08:50 +01:00
void benchmark_wilson(std::vector<int> const &L_list, int niter)
2023-03-31 18:03:39 +01:00
{
2023-04-24 11:35:47 +01:00
int niter_warmup = 10;
printfQuda("==================== wilson dirac operator ====================\n");
2023-06-05 17:07:07 +01:00
#ifdef FLOP_COUNTING_GRID
printfQuda("IMPORTANT: flop counting as in Benchmark_Grid\n");
#else
printfQuda("IMPORTANT: flop counting by QUDA's own convention (different from "
"Benchmark_Grid)\n");
#endif
2023-04-24 11:35:47 +01:00
printfQuda("%5s %15s %15s\n", "L", "time (usec)", "Gflop/s/rank");
2023-06-19 18:08:50 +01:00
for (int L : L_list)
2023-04-24 11:35:47 +01:00
{
auto U = make_gauge_field(L);
auto src = make_source(L);
// 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);
2023-06-19 18:08:50 +01:00
auto res = ColorSpinorField(ColorSpinorParam(src));
2023-04-24 11:35:47 +01:00
// couple iterations without timing to warm up
for (int iter = 0; iter < niter_warmup; ++iter)
2023-06-19 18:08:50 +01:00
dirac.Dslash(res, src, QUDA_EVEN_PARITY);
2023-04-24 11:35:47 +01:00
// actual benchmark with timings
dirac.Flops(); // reset flops counter
device_timer_t device_timer;
device_timer.start();
2023-06-19 18:22:24 +01:00
double start_time = get_timestamp();
2023-04-24 11:35:47 +01:00
for (int iter = 0; iter < niter; ++iter)
2023-06-19 18:08:50 +01:00
dirac.Dslash(res, src, QUDA_EVEN_PARITY);
2023-06-19 18:22:24 +01:00
double end_time = get_timestamp();
2023-04-24 11:35:47 +01:00
device_timer.stop();
double secs = device_timer.last() / niter;
2023-06-05 17:07:07 +01:00
#ifdef FLOP_COUNTING_GRID
// this is the flop counting from Benchmark_Grid
double Nc = 3;
double Nd = 4;
double Ns = 4;
double flops =
(Nc * (6 + (Nc - 1) * 8) * Ns * Nd + 2 * Nd * Nc * Ns + 2 * Nd * Nc * Ns * 2);
flops *= L * L * L * L / 2.0;
#else
double flops = 1.0 * dirac.Flops() / niter;
#endif
printfQuda("%5d %15.2f %15.2f\n", L, secs * 1e6, flops / secs * 1e-9);
2023-06-19 18:08:50 +01:00
json tmp;
tmp["L"] = L;
tmp["Gflops_wilson"] = flops / secs * 1e-9;
2023-06-19 18:22:24 +01:00
tmp["start_time"] = start_time;
tmp["end_time"] = end_time;
2023-06-19 18:08:50 +01:00
json_results["flops"]["results"].push_back(tmp);
2023-06-05 17:07:07 +01:00
}
}
2023-06-19 18:08:50 +01:00
void benchmark_dwf(std::vector<int> const &L_list, int niter)
2023-06-05 17:07:07 +01:00
{
int niter_warmup = 10;
printfQuda("==================== domain wall dirac operator ====================\n");
#ifdef FLOP_COUNTING_GRID
printfQuda("IMPORTANT: flop counting as in Benchmark_Grid\n");
#else
printfQuda("IMPORTANT: flop counting by QUDA's own convention (different from "
"Benchmark_Grid)\n");
#endif
printfQuda("%5s %15s %15s\n", "L", "time (usec)", "Gflop/s/rank");
int Ls = 12;
2023-06-19 18:08:50 +01:00
for (int L : L_list)
2023-06-05 17:07:07 +01:00
{
auto U = make_gauge_field(L);
auto src = make_source(L, Ls);
// create dirac operator
DiracParam param;
param.kappa = 0.10;
param.Ls = Ls;
param.m5 = 0.1;
param.dagger = QUDA_DAG_NO;
param.matpcType = QUDA_MATPC_EVEN_EVEN;
auto dirac = DiracDomainWall(param);
// insert gauge field into the dirac operator
// (the additional nullptr's are for smeared links and fancy preconditioners and such)
dirac.updateFields(&U, nullptr, nullptr, nullptr);
2023-06-19 18:08:50 +01:00
auto res = ColorSpinorField(ColorSpinorParam(src));
2023-06-05 17:07:07 +01:00
// couple iterations without timing to warm up
for (int iter = 0; iter < niter_warmup; ++iter)
2023-06-19 18:08:50 +01:00
dirac.Dslash(res, src, QUDA_EVEN_PARITY);
2023-06-05 17:07:07 +01:00
// actual benchmark with timings
dirac.Flops(); // reset flops counter
device_timer_t device_timer;
device_timer.start();
2023-06-19 18:22:24 +01:00
double start_time = get_timestamp();
2023-06-05 17:07:07 +01:00
for (int iter = 0; iter < niter; ++iter)
2023-06-19 18:08:50 +01:00
dirac.Dslash(res, src, QUDA_EVEN_PARITY);
2023-06-19 18:22:24 +01:00
double end_time = get_timestamp();
2023-06-05 17:07:07 +01:00
device_timer.stop();
double secs = device_timer.last() / niter;
#ifdef FLOP_COUNTING_GRID
// this is the flop counting from Benchmark_Grid
double Nc = 3;
double Nd = 4;
double Ns = 4;
double flops =
(Nc * (6 + (Nc - 1) * 8) * Ns * Nd + 2 * Nd * Nc * Ns + 2 * Nd * Nc * Ns * 2);
flops *= L * L * L * L * Ls / 2.0;
#else
2023-04-24 11:35:47 +01:00
double flops = 1.0 * dirac.Flops() / niter;
2023-06-05 17:07:07 +01:00
#endif
2023-04-24 11:35:47 +01:00
printfQuda("%5d %15.2f %15.2f\n", L, secs * 1e6, flops / secs * 1e-9);
2023-06-19 18:08:50 +01:00
json tmp;
tmp["L"] = L;
tmp["Gflops_dwf4"] = flops / secs * 1e-9;
2023-06-19 18:22:24 +01:00
tmp["start_time"] = start_time;
tmp["end_time"] = end_time;
2023-06-19 18:08:50 +01:00
json_results["flops"]["results"].push_back(tmp);
2023-04-24 11:35:47 +01:00
}
2023-04-21 10:38:28 +01:00
}
2023-06-19 18:08:50 +01:00
void benchmark_axpy(std::vector<int> const &L_list, int niter)
2023-04-21 10:38:28 +01:00
{
2023-04-24 11:35:47 +01:00
// number of iterations for warmup / measurement
// (feel free to change for noise/time tradeoff)
constexpr int niter_warmup = 10;
printfQuda("==================== axpy / memory ====================\n");
2023-04-21 10:38:28 +01:00
ColorSpinorParam param;
2023-04-24 11:35:47 +01:00
param.nDim = 4; // 4-dimensional lattice
param.x[4] = 1; // no fifth dimension
param.nColor = 3; // supported values for nSpin/nColor are configured when compiling
// QUDA. "3*4" will probably always be enabled, so we stick with this
2023-04-21 10:38:28 +01:00
param.nSpin = 4;
2023-04-24 11:35:47 +01:00
param.nVec = 1; // just a single vector
param.siteSubset = QUDA_FULL_SITE_SUBSET; // full lattice = no odd/even
param.pad = 0; // no padding
param.create = QUDA_NULL_FIELD_CREATE; // do not (zero-) initilize the field
param.location = QUDA_CUDA_FIELD_LOCATION; // field should reside on GPU
2023-06-05 17:07:07 +01:00
param.setPrecision(QUDA_SINGLE_PRECISION);
2023-04-24 11:35:47 +01:00
// the following dont matter for an axpy benchmark, but need to choose something
2023-04-21 10:38:28 +01:00
param.pc_type = QUDA_4D_PC;
param.siteOrder = QUDA_EVEN_ODD_SITE_ORDER;
param.gammaBasis = QUDA_DEGRAND_ROSSI_GAMMA_BASIS;
2023-04-24 11:35:47 +01:00
printfQuda("%5s %15s %15s %15s %15s\n", "L", "size (MiB/rank)", "time (usec)",
"GiB/s/rank", "Gflop/s/rank");
for (int L : L_list)
{
// IMPORTANT: all of `param.x`, `field_elements`, `field.Bytes()`
// are LOCAL, i.e. per rank / per GPU
param.x[0] = L;
param.x[1] = L;
param.x[2] = L;
param.x[3] = L;
// number of (real) elements in one (local) field
size_t field_elements = 2 * param.x[0] * param.x[1] * param.x[2] * param.x[3] *
param.nColor * param.nSpin;
// create the field(s)
auto fieldA = ColorSpinorField(param);
auto fieldB = ColorSpinorField(param);
2023-06-05 17:07:07 +01:00
assert(fieldA.Bytes() == sizeof(float) * field_elements); // sanity check
assert(fieldB.Bytes() == sizeof(float) * field_elements); // sanity check
2023-04-24 11:35:47 +01:00
// fill fields with random values
quda::RNG rng(fieldA, 1234);
spinorNoise(fieldA, rng, QUDA_NOISE_GAUSS);
spinorNoise(fieldB, rng, QUDA_NOISE_GAUSS);
// number of operations / bytes per iteration
// axpy is one addition, one multiplication, two read, one write
double flops = 2 * field_elements;
2023-06-05 17:07:07 +01:00
double memory = 3 * sizeof(float) * field_elements;
2023-04-24 11:35:47 +01:00
// do some iterations to to let QUDA do its internal tuning and also stabilize cache
// behaviour and such
for (int iter = 0; iter < niter_warmup; ++iter)
blas::axpy(1.234, fieldA, fieldB);
// running the actual benchmark
device_timer_t device_timer;
device_timer.start();
2023-06-19 18:22:24 +01:00
double start_time = get_timestamp();
2023-04-24 11:35:47 +01:00
for (int iter = 0; iter < niter; ++iter)
blas::axpy(1.234, fieldA, fieldB);
2023-06-19 18:22:24 +01:00
double end_time = get_timestamp();
2023-04-24 11:35:47 +01:00
device_timer.stop();
double secs = device_timer.last() / niter; // seconds per iteration
2023-06-19 18:08:50 +01:00
double mem_MiB = memory / 1024. / 1024.;
double GBps = mem_MiB / 1024 / secs;
printfQuda("%5d %15.2f %15.2f %15.2f %15.2f\n", L, mem_MiB, secs * 1e6, GBps,
flops / secs * 1e-9);
json tmp;
tmp["L"] = L;
tmp["size_MB"] = mem_MiB;
tmp["GBps"] = GBps;
tmp["GFlops"] = flops / secs * 1e-9;
2023-06-19 18:22:24 +01:00
tmp["start_time"] = start_time;
tmp["end_time"] = end_time;
2023-06-19 18:08:50 +01:00
json_results["axpy"].push_back(tmp);
2023-04-24 11:35:47 +01:00
}
2023-03-31 18:03:39 +01:00
}
int main(int argc, char **argv)
{
2023-06-19 18:08:50 +01:00
std::string json_filename = ""; // empty indicates no json output
for (int i = 0; i < argc; i++)
{
if (std::string(argv[i]) == "--json-out")
json_filename = argv[i + 1];
}
2023-04-24 11:35:47 +01:00
initComms(argc, argv);
2023-03-31 18:03:39 +01:00
2023-04-21 10:38:28 +01:00
initQuda(-1); // -1 for multi-gpu. otherwise this selects the device to be used
2023-03-31 18:03:39 +01:00
2023-04-21 10:38:28 +01:00
// verbosity options are:
// SILENT, SUMMARIZE, VERBOSE, DEBUG_VERBOSE
2023-04-24 11:35:47 +01:00
setVerbosity(QUDA_SUMMARIZE);
2023-03-31 18:03:39 +01:00
2023-04-24 11:35:47 +01:00
printfQuda("MPI layout = %d %d %d %d\n", mpi_grid[0], mpi_grid[1], mpi_grid[2],
mpi_grid[3]);
2023-06-19 18:08:50 +01:00
benchmark_axpy({8, 12, 16, 24, 32, 48}, 20);
2023-04-24 11:35:47 +01:00
setVerbosity(QUDA_SILENT);
2023-06-19 18:08:50 +01:00
benchmark_wilson({8, 12, 16, 24, 32, 48}, 20);
benchmark_dwf({8, 12, 16, 24, 32}, 20);
2023-04-24 11:35:47 +01:00
setVerbosity(QUDA_SUMMARIZE);
2023-03-31 18:03:39 +01:00
2023-04-24 11:35:47 +01:00
printfQuda("==================== done with all benchmarks ====================\n");
2023-06-19 18:08:50 +01:00
if (!json_filename.empty())
{
printfQuda("writing benchmark results to %s\n", json_filename.c_str());
int me = 0;
MPI_Comm_rank(MPI_COMM_WORLD, &me);
if (me == 0)
{
std::ofstream json_file(json_filename);
json_file << std::setw(2) << json_results;
}
}
2023-03-31 18:03:39 +01:00
endQuda();
quda::comm_finalize();
MPI_Finalize();
}