SLATE 2024.05.31
Software for Linear Algebra Targeting Exascale
Loading...
Searching...
No Matches

Functions

template<typename scalar_t >
void slate::scale (blas::real_type< scalar_t > numer, blas::real_type< scalar_t > denom, Matrix< scalar_t > &A, Options const &opts)
 Scale matrix entries by the real scalar numer/denom.
 
template<typename scalar_t >
void slate::scale (blas::real_type< scalar_t > numer, blas::real_type< scalar_t > denom, BaseTrapezoidMatrix< scalar_t > &A, Options const &opts)
 Scale matrix entries by the real scalar numer/denom.
 
template<typename scalar_t >
void slate::set (scalar_t offdiag_value, scalar_t diag_value, Matrix< scalar_t > &A, Options const &opts)
 Set matrix entries.
 
template<typename scalar_t >
void slate::set (scalar_t offdiag_value, scalar_t diag_value, BaseTrapezoidMatrix< scalar_t > &A, Options const &opts)
 Set matrix entries.
 
template<typename scalar_t >
void slate::set (std::function< scalar_t(int64_t i, int64_t j) > const &value, Matrix< scalar_t > &A, Options const &opts)
 Set matrix entries.
 
template<typename scalar_t >
void slate::set (std::function< scalar_t(int64_t i, int64_t j) > const &value, BaseTrapezoidMatrix< scalar_t > &A, Options const &opts)
 Set matrix entries.
 

Detailed Description

Function Documentation

◆ scale() [1/2]

template<typename scalar_t >
void slate::scale ( blas::real_type< scalar_t >  numer,
blas::real_type< scalar_t >  denom,
BaseTrapezoidMatrix< scalar_t > &  A,
Options const &  opts 
)

Scale matrix entries by the real scalar numer/denom.

Transposition is currently ignored. TODO: Inspect transposition?

Template Parameters
scalar_tOne of float, double, std::complex<float>, std::complex<double>.
Parameters
[in]AThe m-by-n matrix A.
[in]optsAdditional options, as map of name = value pairs. Possible options:
  • Option::Target: Implementation to target. Possible values:
    • HostTask: OpenMP tasks on CPU host [default].
    • HostNest: nested OpenMP parallel for loop on CPU host.
    • HostBatch: batched BLAS on CPU host.
    • Devices: batched BLAS on GPU device.

◆ scale() [2/2]

template<typename scalar_t >
void slate::scale ( blas::real_type< scalar_t >  numer,
blas::real_type< scalar_t >  denom,
Matrix< scalar_t > &  A,
Options const &  opts 
)

Scale matrix entries by the real scalar numer/denom.

Transposition is currently ignored. TODO: Inspect transposition?

Template Parameters
scalar_tOne of float, double, std::complex<float>, std::complex<double>.
Parameters
[in]AThe m-by-n matrix A.
[in]optsAdditional options, as map of name = value pairs. Possible options:
  • Option::Target: Implementation to target. Possible values:
    • HostTask: OpenMP tasks on CPU host [default].
    • HostNest: nested OpenMP parallel for loop on CPU host.
    • HostBatch: batched BLAS on CPU host.
    • Devices: batched BLAS on GPU device.

◆ set() [1/4]

template<typename scalar_t >
void slate::set ( scalar_t  offdiag_value,
scalar_t  diag_value,
BaseTrapezoidMatrix< scalar_t > &  A,
Options const &  opts 
)

Set matrix entries.

Transposition is currently ignored. TODO: Inspect transposition?

Template Parameters
scalar_tOne of float, double, std::complex<float>, std::complex<double>.
Parameters
[in]offdiag_valueValue to set off-diagonal elements to.
[in]diag_valueValue to set diagonal elements to.
[in,out]AThe m-by-n matrix A.
[in]optsAdditional options, as map of name = value pairs. Possible options:
  • Option::Target: Implementation to target. Possible values:
    • HostTask: OpenMP tasks on CPU host [default].
    • HostNest: nested OpenMP parallel for loop on CPU host.
    • HostBatch: batched BLAS on CPU host.
    • Devices: batched BLAS on GPU device.

◆ set() [2/4]

template<typename scalar_t >
void slate::set ( scalar_t  offdiag_value,
scalar_t  diag_value,
Matrix< scalar_t > &  A,
Options const &  opts 
)

Set matrix entries.

Transposition is currently ignored. TODO: Inspect transposition?

Template Parameters
scalar_tOne of float, double, std::complex<float>, std::complex<double>.
Parameters
[in]offdiag_valueValue to set off-diagonal elements to.
[in]diag_valueValue to set diagonal elements to.
[in,out]AThe m-by-n matrix A.
[in]optsAdditional options, as map of name = value pairs. Possible options:
  • Option::Target: Implementation to target. Possible values:
    • HostTask: OpenMP tasks on CPU host [default].
    • HostNest: nested OpenMP parallel for loop on CPU host.
    • HostBatch: batched BLAS on CPU host.
    • Devices: batched BLAS on GPU device.

◆ set() [3/4]

template<typename scalar_t >
void slate::set ( std::function< scalar_t(int64_t i, int64_t j) > const &  value,
BaseTrapezoidMatrix< scalar_t > &  A,
Options const &  opts 
)

Set matrix entries.

Transposition is automatically handled.

Template Parameters
scalar_tOne of float, double, std::complex<float>, std::complex<double>.
Parameters
[in]valueA function that takes global row and column indices i and j, and returns the scalar value for entry Aij.
[in,out]AThe m-by-n matrix A.
[in]optsAdditional options, as map of name = value pairs. Currently no options. It always uses Target = HostTask, since lambda is a CPU function.

◆ set() [4/4]

template<typename scalar_t >
void slate::set ( std::function< scalar_t(int64_t i, int64_t j) > const &  value,
Matrix< scalar_t > &  A,
Options const &  opts 
)

Set matrix entries.

Transposition is automatically handled.

Template Parameters
scalar_tOne of float, double, std::complex<float>, std::complex<double>.
Parameters
[in]valueA function that takes global row and column indices i and j, and returns the scalar value for entry Aij.
[in,out]AThe m-by-n matrix A.
[in]optsAdditional options, as map of name = value pairs. Currently no options. It always uses Target = HostTask, since lambda is a CPU function.