RunSparseMatrixBenchmark {RHPCBenchmark}R Documentation

Runs all of the sparse matrix microbenchmarks

Description

RunSparseMatrixBenchmark runs all of the microbenchmarks for performance testing the sparse matrix linear algebra kernels

Usage

RunSparseMatrixBenchmark(runIdentifier, resultsDirectory,
  matrixVectorMicrobenchmarks = GetSparseMatrixVectorDefaultMicrobenchmarks(),
  choleskyMicrobenchmarks = GetSparseCholeskyDefaultMicrobenchmarks(),
  luMicrobenchmarks = GetSparseLuDefaultMicrobenchmarks(),
  qrMicrobenchmarks = GetSparseQrDefaultMicrobenchmarks())

Arguments

runIdentifier

a character string specifying the suffix to be appended to the base of the file name of the output CSV format files

resultsDirectory

a character string specifying the directory where all of the CSV performance results files will be saved

matrixVectorMicrobenchmarks

a list of SparseMatrixMicrobenchmark objects defining the matrix-vector multiplication microbenchmarks to execute as part of the sparse matrix benchmark. Default values are provided by the function GetSparseMatrixVectorDefaultMicrobenchmarks. If the value is NULL, then all of the matrix-vector multiplication microbenchmarks will be skipped.

choleskyMicrobenchmarks

a list of SparseMatrixMicrobenchmark objects defining the Cholesky factorization microbenchmarks to execute as part of the sparse matrix benchmark. Default values are provided by the function GetSparseCholeskyDefaultMicrobenchmarks. If the value is NULL, then all of the Cholesky factorization microbenchmarks will be skipped.

luMicrobenchmarks

a list of SparseMatrixMicrobenchmark objects defining the LU factorization microbenchmarks to execute as part of the sparse matrix benchmark. Default values are provided by the function GetSparseLuDefaultMicrobenchmarks. If the value is NULL, then all of the LU factorization microbenchmarks will be skipped.

qrMicrobenchmarks

a list of SparseMatrixMicrobenchmark objects defining the QR factorization microbenchmarks to execute as part of the sparse matrix benchmark. Default values are provided by the function GetSparseQrDefaultMicrobenchmarks. If the value is NULL, then all of the QR factorization microbenchmarks will be skipped.

Details

This function runs the sparse matrix microbenchmarks, which are divided into four categories supported by this benchmark, defined in the matrixVectorMicrobenchmarks, choleskyMicrobenchmarks, luMicrobenchmarks, and qrMicrobenchmarks input lists For each microbenchmark, it attempts to create a separate output file in CSV format containing the performance results for each matrix tested by the microbenchmark. The names of the output files follow the format benchmarkName_runIdentifier.csv, where benchmarkName is specified in the SparseMatrixMicrobenchmark object of each microbenchmark and runIdentifier is an input parameter to this function. If the file, already exists, the results will be appended to the existing file. Each input lists contains instances of the SparseMatrixMicrobenchmark class defining each microbenchmark. Each microbenchmark object with the active field set to TRUE will be executed. The lists of default microbenchmarks are generated by the functions GetSparseMatrixVectorDefaultMicrobenchmarks, GetSparseCholeskyDefaultMicrobenchmarks, GetSparseLuDefaultMicrobenchmarks, and GetSparseQrDefaultMicrobenchmarks. Each SparseMatrixMicrobenchmark specifies an R data file which contains the sparse matrix object needed by the microbenchmark. The needed R data files should either be given in an attached R package or given in the data subdirectory of the current working directory, and they should have the extension .RData. If the linear algebra kernels are multithreaded, by linking to multithreaded BLAS or LAPACK libraries for example, then the number of threads must be retrievable from an environment variable which is set before execution of the R programming environment. The name of the environment variable specifying the number of threads must be provided in the R HPC benchmark environment variable R_BENCH_NUM_THREADS_VARIABLE. This function will retrieve the number of threads through R_BENCH_NUM_THREADS_VARIABLE so that the number of threads can be printed to the results files and recorded in data frames for reporting purposes. This function utilizes the number of threads only for reporting purposes and is not used by the benchmark to effect the actual number of threads utilized by the kernels, as that is assumed to be controlled by the numerical library. An error exception will be thrown if the environment variable R_BENCH_NUM_THREADS_VARIABLE and the variable it is set to are not both set.

Value

a data frame containing the benchmark name, user, system, and elapsed (wall clock) times of each performance trial for each microbenchmark

See Also

GetSparseMatrixVectorDefaultMicrobenchmarks GetSparseCholeskyDefaultMicrobenchmarks GetSparseLuDefaultMicrobenchmarks GetSparseQrDefaultMicrobenchmarks GetSparseMatrixVectorExampleMicrobenchmarks GetSparseCholeskyExampleMicrobenchmarks

Examples

## Not run: 
# Set needed environment variables for multithreading.  Only single threading
# is used in the example.
#
# Note: The environment variables are usually set by the user before starting
#       the R programming environment; they are set here only to facilitate
#       a working example.  See the section on multithreading in the vignette
#       for further details.
Sys.setenv(R_BENCH_NUM_THREADS_VARIABLE="MKL_NUM_THREADS")
Sys.setenv(MKL_NUM_THREADS="1")
#
# Generate example microbenchmarks that can be run in a few minutes; see
# the vignette for more involved examples.  The matvec_laplacian7pt_100
# and cholesky_ct20stif microbenchmarks are defined in the examples.
#
# Note: The example microbenchmarks are different than the microbenchmarks
#       generated by
#       \code{\link{GetSparseMatrixVectorDefaultMicrobenchmarks}},
#       \code{\link{GetSparseCholeskyDefaultMicrobenchmarks}},
#       \code{\link{GetSparseLuDefaultMicrobenchmarks}}, and
#       \code{\link{GetSparseQrDefaultMicrobenchmarks}};
#       they were chosen for their short run times and suitability for
#       example code. 
exampleMatrixVectorMicrobenchmarks <- GetSparseMatrixVectorExampleMicrobenchmarks()
exampleCholeskyMicrobenchmarks <- GetSparseCholeskyExampleMicrobenchmarks()
# Set the output directory of the CSV summary results files
resultsDirectory <- "./SparseMatrixExampleOutput"
# Create the output directory
dir.create(resultsDirectory)
# Set an appropriate run identifier
runIdentifier <- "example"
# Run only the matrix-vector and Cholesky factorization microbenchmarks, as
# the others take a long time
resultsFrame <- RunSparseMatrixBenchmark(runIdentifier, resultsDirectory,
   matrixVectorMicrobenchmarks=exampleMatrixVectorMicrobenchmarks,
   choleskyMicrobenchmarks=exampleCholeskyMicrobenchmarks,
   luMicrobenchmarks=NULL,
   qrMicrobenchmarks=NULL)

# This example runs only the Cholesky factorization microbenchmarks.
runIdentifier <- "choleksy_only"
# Run only the sparse Choleksy factorization microbenchmarks
choleskyResults <- RunSparseMatrixBenchmark(runIdentifier, resultsDirectory,
   matrixVectorMicrobenchmarks=NULL, luMicrobenchmarks=NULL,
   qrMicrobenchmarks=NULL)


## End(Not run)


[Package RHPCBenchmark version 0.1.0 Index]