Ginkgo
Generated from pipelines/1330831941 branch based on master. Ginkgo version 1.8.0
A numerical linear algebra library targeting many-core architectures
|
The Ginkgo namespace. More...
Namespaces | |
accessor | |
The accessor namespace. | |
factorization | |
The Factorization namespace. | |
log | |
The logger namespace . | |
matrix | |
The matrix namespace. | |
multigrid | |
The multigrid components namespace. | |
name_demangling | |
The name demangling namespace. | |
preconditioner | |
The Preconditioner namespace. | |
reorder | |
The Reorder namespace. | |
solver | |
The ginkgo Solve namespace. | |
stop | |
The Stopping criterion namespace. | |
syn | |
The Synthesizer namespace. | |
xstd | |
The namespace for functionalities after C++14 standard. | |
Classes | |
class | AbsoluteComputable |
The AbsoluteComputable is an interface that allows to get the component wise absolute of a LinOp. More... | |
class | AbstractFactory |
The AbstractFactory is a generic interface template that enables easy implementation of the abstract factory design pattern. More... | |
class | AllocationError |
AllocationError is thrown if a memory allocation fails. More... | |
class | Allocator |
Provides generic allocation and deallocation functionality to be used by an Executor. More... | |
class | amd_device |
amd_device handles the number of executor on Amd devices and have the corresponding recursive_mutex. More... | |
struct | are_all_integral |
Evaluates if all template arguments Args fulfill std::is_integral. More... | |
struct | are_all_integral< First, Args... > |
class | array |
An array is a container which encapsulates fixed-sized arrays, stored on the Executor tied to the array. More... | |
class | BadDimension |
BadDimension is thrown if an operation is being applied to a LinOp with bad dimensions. More... | |
struct | batch_dim |
A type representing the dimensions of a multidimensional batch object. More... | |
class | BlockOperator |
A BlockOperator represents a linear operator that is partitioned into multiple blocks. More... | |
class | BlockSizeError |
Error that denotes issues between block sizes and matrix dimensions. More... | |
class | Combination |
The Combination class can be used to construct a linear combination of multiple linear operators c1 * op1 + c2 * op2 + ... More... | |
class | Composition |
The Composition class can be used to compose linear operators op1, op2, ..., opn and obtain the operator op1 * op2 * ... More... | |
class | ConvertibleTo |
ConvertibleTo interface is used to mark that the implementer can be converted to the object of ResultType. More... | |
class | CpuAllocator |
Allocator using new/delete. More... | |
class | CpuAllocatorBase |
Implement this interface to provide an allocator for OmpExecutor or ReferenceExecutor. More... | |
class | CpuTimer |
A timer using std::chrono::steady_clock for timing. More... | |
struct | cpx_real_type |
Access the underlying real type of a complex number. More... | |
struct | cpx_real_type< std::complex< T > > |
Specialization for complex types. | |
class | CublasError |
CublasError is thrown when a cuBLAS routine throws a non-zero error code. More... | |
class | cuda_stream |
An RAII wrapper for a custom CUDA stream. More... | |
class | CudaAllocator |
Allocator using cudaMalloc. More... | |
class | CudaAllocatorBase |
Implement this interface to provide an allocator for CudaExecutor. More... | |
class | CudaAsyncAllocator |
class | CudaError |
CudaError is thrown when a CUDA routine throws a non-zero error code. More... | |
class | CudaExecutor |
This is the Executor subclass which represents the CUDA device. More... | |
class | CudaHostAllocator |
class | CudaTimer |
A timer using events for timing on a CudaExecutor. More... | |
class | CudaUnifiedAllocator |
class | CufftError |
CufftError is thrown when a cuFFT routine throws a non-zero error code. More... | |
class | CurandError |
CurandError is thrown when a cuRAND routine throws a non-zero error code. More... | |
class | CusparseError |
CusparseError is thrown when a cuSPARSE routine throws a non-zero error code. More... | |
struct | default_converter |
Used to convert objects of type S to objects of type R using static_cast. More... | |
class | deferred_factory_parameter |
Represents a factory parameter of factory type that can either initialized by a pre-existing factory or by passing in a factory_parameters object whose .on(exec) will be called to instantiate a factory. More... | |
class | device_matrix_data |
This type is a device-side equivalent to matrix_data. More... | |
class | DiagonalExtractable |
The diagonal of a LinOp implementing this interface can be extracted. More... | |
class | DiagonalLinOpExtractable |
The diagonal of a LinOp can be extracted. More... | |
struct | dim |
A type representing the dimensions of a multidimensional object. More... | |
struct | dim< 1u, DimensionType > |
class | DimensionMismatch |
DimensionMismatch is thrown if an operation is being applied to LinOps of incompatible size. More... | |
class | DpcppExecutor |
This is the Executor subclass which represents a DPC++ enhanced device. More... | |
class | DpcppTimer |
A timer using kernels for timing on a DpcppExecutor in profiling mode. More... | |
class | enable_parameters_type |
The enable_parameters_type mixin is used to create a base implementation of the factory parameters structure. More... | |
class | EnableAbsoluteComputation |
The EnableAbsoluteComputation mixin provides the default implementations of compute_absolute_linop and the absolute interface. More... | |
class | EnableAbstractPolymorphicObject |
This mixin inherits from (a subclass of) PolymorphicObject and provides a base implementation of a new abstract object. More... | |
class | EnableCreateMethod |
This mixin implements a static create() method on ConcreteType that dynamically allocates the memory, uses the passed-in arguments to construct the object, and returns an std::unique_ptr to such an object. More... | |
class | EnableDefaultFactory |
This mixin provides a default implementation of a concrete factory. More... | |
class | EnableLinOp |
The EnableLinOp mixin can be used to provide sensible default implementations of the majority of the LinOp and PolymorphicObject interface. More... | |
class | EnablePolymorphicAssignment |
This mixin is used to enable a default PolymorphicObject::copy_from() implementation for objects that have implemented conversions between them. More... | |
class | EnablePolymorphicObject |
This mixin inherits from (a subclass of) PolymorphicObject and provides a base implementation of a new concrete polymorphic object. More... | |
struct | err |
class | Error |
The Error class is used to report exceptional behaviour in library functions. More... | |
class | Executor |
The first step in using the Ginkgo library consists of creating an executor. More... | |
class | executor_deleter |
This is a deleter that uses an executor's free method to deallocate the data. More... | |
class | executor_deleter< T[]> |
class | hip_stream |
An RAII wrapper for a custom HIP stream. More... | |
class | HipAllocator |
class | HipAllocatorBase |
Implement this interface to provide an allocator for HipExecutor. More... | |
class | HipAsyncAllocator |
class | HipblasError |
HipblasError is thrown when a hipBLAS routine throws a non-zero error code. More... | |
class | HipError |
HipError is thrown when a HIP routine throws a non-zero error code. More... | |
class | HipExecutor |
This is the Executor subclass which represents the HIP enhanced device. More... | |
class | HipfftError |
HipfftError is thrown when a hipFFT routine throws a non-zero error code. More... | |
class | HipHostAllocator |
class | HiprandError |
HiprandError is thrown when a hipRAND routine throws a non-zero error code. More... | |
class | HipsparseError |
HipsparseError is thrown when a hipSPARSE routine throws a non-zero error code. More... | |
class | HipTimer |
A timer using events for timing on a HipExecutor. More... | |
class | HipUnifiedAllocator |
class | index_set |
An index set class represents an ordered set of intervals. More... | |
class | InvalidStateError |
Exception thrown if an object is in an invalid state. More... | |
class | KernelNotFound |
KernelNotFound is thrown if Ginkgo cannot find a kernel which satisfies the criteria imposed by the input arguments. More... | |
class | LinOp |
class | LinOpFactory |
A LinOpFactory represents a higher order mapping which transforms one linear operator into another. More... | |
class | machine_topology |
The machine topology class represents the hierarchical topology of a machine, including NUMA nodes, cores and PCI Devices. More... | |
class | matrix_assembly_data |
This structure is used as an intermediate type to assemble a sparse matrix. More... | |
struct | matrix_data |
This structure is used as an intermediate data type to store a sparse matrix. More... | |
struct | matrix_data_entry |
Type used to store nonzeros. More... | |
class | MetisError |
MetisError is thrown when METIS routine throws an error code. More... | |
class | MpiError |
MpiError is thrown when a MPI routine throws a non-zero error code. More... | |
class | NotCompiled |
NotCompiled is thrown when attempting to call an operation which is a part of a module that was not compiled on the system. More... | |
class | NotImplemented |
NotImplemented is thrown in case an operation has not yet been implemented (but will be implemented in the future). More... | |
class | NotSupported |
NotSupported is thrown in case it is not possible to perform the requested operation on the given object type. More... | |
class | null_deleter |
This is a deleter that does not delete the object. More... | |
class | null_deleter< T[]> |
class | nvidia_device |
nvidia_device handles the number of executor on Nvidia devices and have the corresponding recursive_mutex. More... | |
class | OmpExecutor |
This is the Executor subclass which represents the OpenMP device (typically CPU). More... | |
class | Operation |
Operations can be used to define functionalities whose implementations differ among devices. More... | |
class | OutOfBoundsError |
OutOfBoundsError is thrown if a memory access is detected to be out-of-bounds. More... | |
class | OverflowError |
OverflowError is thrown when an index calculation for storage requirements overflows. More... | |
class | Permutable |
Linear operators which support permutation should implement the Permutable interface. More... | |
class | Perturbation |
The Perturbation class can be used to construct a LinOp to represent the operation (identity + scalar * basis * projector) . More... | |
class | PolymorphicObject |
A PolymorphicObject is the abstract base for all "heavy" objects in Ginkgo that behave polymorphically. More... | |
class | precision_reduction |
This class is used to encode storage precisions of low precision algorithms. More... | |
class | Preconditionable |
A LinOp implementing this interface can be preconditioned. More... | |
class | ptr_param |
This class is used for function parameters in the place of raw pointers. More... | |
class | range |
A range is a multidimensional view of the memory. More... | |
class | ReadableFromMatrixData |
A LinOp implementing this interface can read its data from a matrix_data structure. More... | |
class | ReferenceExecutor |
This is a specialization of the OmpExecutor, which runs the reference implementations of the kernels used for debugging purposes. More... | |
class | ScaledIdentityAddable |
Adds the operation M <- a I + b M for matrix M, identity operator I and scalars a and b, where M is the calling object. More... | |
class | scoped_device_id_guard |
This move-only class uses RAII to set the device id within a scoped block, if necessary. More... | |
struct | segmented_array |
A minimal interface for a segmented array. More... | |
struct | span |
A span is a lightweight structure used to create sub-ranges from other ranges. More... | |
class | stopping_status |
This class is used to keep track of the stopping status of one vector. More... | |
class | StreamError |
StreamError is thrown if accessing a stream failed. More... | |
class | time_point |
An opaque wrapper for a time point generated by a timer. More... | |
class | Timer |
Represents a generic timer that can be used to record time points and measure time differences on host or device streams. More... | |
class | Transposable |
Linear operators which support transposition should implement the Transposable interface. More... | |
class | truncated |
class | UnsupportedMatrixProperty |
Exception throws if a matrix does not have a property required by a numerical method. More... | |
class | UseComposition |
The UseComposition class can be used to store the composition information in LinOp. More... | |
class | ValueMismatch |
ValueMismatch is thrown if two values are not equal. More... | |
struct | version |
This structure is used to represent versions of various Ginkgo modules. More... | |
class | version_info |
Ginkgo uses version numbers to label new features and to communicate backward compatibility guarantees: More... | |
class | WritableToMatrixData |
A LinOp implementing this interface can write its data to a matrix_data structure. More... | |
Typedefs | |
template<typename ValueType > | |
using | Array = array< ValueType > |
template<typename ConcreteFactory , typename ConcreteLinOp , typename ParametersType , typename PolymorphicBase = LinOpFactory> | |
using | EnableDefaultLinOpFactory = EnableDefaultFactory< ConcreteFactory, ConcreteLinOp, ParametersType, PolymorphicBase > |
This is an alias for the EnableDefaultFactory mixin, which correctly sets the template parameters to enable a subclass of LinOpFactory. More... | |
using | MachineTopology = machine_topology |
template<typename T > | |
using | is_complex_s = detail::is_complex_impl< T > |
Allows to check if T is a complex value during compile time by accessing the value attribute of this struct. More... | |
template<typename T > | |
using | is_complex_or_scalar_s = detail::is_complex_or_scalar_impl< T > |
Allows to check if T is a complex or scalar value during compile time by accessing the value attribute of this struct. More... | |
template<typename T > | |
using | remove_complex = typename detail::remove_complex_s< T >::type |
Obtain the type which removed the complex of complex/scalar type or the template parameter of class by accessing the type attribute of this struct. More... | |
template<typename T > | |
using | to_complex = typename detail::to_complex_s< T >::type |
Obtain the type which adds the complex of complex/scalar type or the template parameter of class by accessing the type attribute of this struct. More... | |
template<typename T > | |
using | to_real = remove_complex< T > |
to_real is alias of remove_complex More... | |
template<typename T > | |
using | next_precision = typename detail::next_precision_impl< T >::type |
Obtains the next type in the singly-linked precision list. | |
template<typename T > | |
using | previous_precision = next_precision< T > |
Obtains the previous type in the singly-linked precision list. More... | |
template<typename T > | |
using | reduce_precision = typename detail::reduce_precision_impl< T >::type |
Obtains the next type in the hierarchy with lower precision than T. | |
template<typename T > | |
using | increase_precision = typename detail::increase_precision_impl< T >::type |
Obtains the next type in the hierarchy with higher precision than T. | |
template<typename... Ts> | |
using | highest_precision = typename detail::highest_precision_variadic< Ts... >::type |
Obtains the smallest arithmetic type that is able to store elements of all template parameter types exactly. More... | |
template<typename T , size_type Limit = sizeof(uint16) * byte_size> | |
using | truncate_type = std::conditional_t< detail::type_size_impl< T >::value >=2 *Limit, typename detail::truncate_type_impl< T >::type, T > |
Truncates the type by half (by dropping bits), but ensures that it is at least Limit bits wide. | |
using | size_type = std::size_t |
Integral type used for allocation quantities. | |
using | int8 = std::int8_t |
8-bit signed integral type. | |
using | int16 = std::int16_t |
16-bit signed integral type. | |
using | int32 = std::int32_t |
32-bit signed integral type. | |
using | int64 = std::int64_t |
64-bit signed integral type. | |
using | uint8 = std::uint8_t |
8-bit unsigned integral type. | |
using | uint16 = std::uint16_t |
16-bit unsigned integral type. | |
using | uint32 = std::uint32_t |
32-bit unsigned integral type. | |
using | uint64 = std::uint64_t |
64-bit unsigned integral type. | |
using | uintptr = std::uintptr_t |
Unsigned integer type capable of holding a pointer to void. | |
using | float16 = half |
Half precision floating point type. | |
using | float32 = float |
Single precision floating point type. | |
using | float64 = double |
Double precision floating point type. | |
using | full_precision = double |
The most precise floating-point type. | |
using | default_precision = double |
Precision used if no precision is explicitly specified. | |
Enumerations | |
enum | log_propagation_mode { log_propagation_mode::never, log_propagation_mode::automatic } |
How Logger events are propagated to their Executor. More... | |
enum | allocation_mode { device, unified_global, unified_host } |
Specify the mode of allocation for CUDA/HIP GPUs. More... | |
enum | layout_type { layout_type::array, layout_type::coordinate } |
Specifies the layout type when writing data in matrix market format. More... | |
Functions | |
template<typename ValueType > | |
ValueType | reduce_add (const array< ValueType > &input_arr, const ValueType init_val=0) |
Reduce (sum) the values in the array. More... | |
template<typename ValueType > | |
void | reduce_add (const array< ValueType > &input_arr, array< ValueType > &result) |
Reduce (sum) the values in the array. More... | |
template<typename ValueType > | |
array< ValueType > | make_array_view (std::shared_ptr< const Executor > exec, size_type size, ValueType *data) |
Helper function to create an array view deducing the value type. More... | |
template<typename ValueType > | |
detail::const_array_view< ValueType > | make_const_array_view (std::shared_ptr< const Executor > exec, size_type size, const ValueType *data) |
Helper function to create a const array view deducing the value type. More... | |
template<typename DimensionType > | |
batch_dim< 2, DimensionType > | transpose (const batch_dim< 2, DimensionType > &input) |
Returns a batch_dim object with its dimensions swapped for batched operators. More... | |
template<typename DimensionType > | |
constexpr dim< 2, DimensionType > | transpose (const dim< 2, DimensionType > &dimensions) noexcept |
Returns a dim<2> object with its dimensions swapped. More... | |
template<typename T > | |
constexpr bool | is_complex () |
Checks if T is a complex type. More... | |
template<typename T > | |
constexpr bool | is_complex_or_scalar () |
Checks if T is a complex/scalar type. More... | |
template<typename T > | |
constexpr reduce_precision< T > | round_down (T val) |
Reduces the precision of the input parameter. More... | |
template<typename T > | |
constexpr increase_precision< T > | round_up (T val) |
Increases the precision of the input parameter. More... | |
constexpr int64 | ceildiv (int64 num, int64 den) |
Performs integer division with rounding up. More... | |
template<typename T > | |
constexpr T | zero () |
Returns the additive identity for T. More... | |
template<typename T > | |
constexpr T | zero (const T &) |
Returns the additive identity for T. More... | |
template<typename T > | |
constexpr T | one () |
Returns the multiplicative identity for T. More... | |
template<typename T > | |
constexpr T | one (const T &) |
Returns the multiplicative identity for T. More... | |
template<typename T > | |
constexpr bool | is_zero (T value) |
Returns true if and only if the given value is zero. More... | |
template<typename T > | |
constexpr bool | is_nonzero (T value) |
Returns true if and only if the given value is not zero. More... | |
template<typename T > | |
constexpr T | max (const T &x, const T &y) |
Returns the larger of the arguments. More... | |
template<typename T > | |
constexpr T | min (const T &x, const T &y) |
Returns the smaller of the arguments. More... | |
template<typename T > | |
constexpr auto | real (const T &x) |
Returns the real part of the object. More... | |
template<typename T > | |
constexpr auto | imag (const T &x) |
Returns the imaginary part of the object. More... | |
template<typename T > | |
constexpr auto | conj (const T &x) |
Returns the conjugate of an object. More... | |
template<typename T > | |
constexpr auto | squared_norm (const T &x) -> decltype(real(conj(x) *x)) |
Returns the squared norm of the object. More... | |
template<typename T > | |
constexpr xstd::enable_if_t<!is_complex_s< T >::value, T > | abs (const T &x) |
Returns the absolute value of the object. More... | |
template<typename T > | |
constexpr xstd::enable_if_t< is_complex_s< T >::value, remove_complex< T > > | abs (const T &x) |
template<typename T > | |
constexpr T | pi () |
Returns the value of pi. More... | |
template<typename T > | |
constexpr std::complex< remove_complex< T > > | unit_root (int64 n, int64 k=1) |
Returns the value of exp(2 * pi * i * k / n), i.e. More... | |
template<typename T > | |
constexpr uint32 | get_significant_bit (const T &n, uint32 hint=0u) noexcept |
Returns the position of the most significant bit of the number. More... | |
template<typename T > | |
constexpr T | get_superior_power (const T &base, const T &limit, const T &hint=T{1}) noexcept |
Returns the smallest power of base not smaller than limit . More... | |
template<typename T > | |
std::enable_if_t<!is_complex_s< T >::value, bool > | is_finite (const T &value) |
Checks if a floating point number is finite, meaning it is neither +/- infinity nor NaN. More... | |
template<typename T > | |
std::enable_if_t< is_complex_s< T >::value, bool > | is_finite (const T &value) |
Checks if all components of a complex value are finite, meaning they are neither +/- infinity nor NaN. More... | |
template<typename T > | |
T | safe_divide (T a, T b) |
Computes the quotient of the given parameters, guarding against division by zero. More... | |
template<typename T > | |
std::enable_if_t<!is_complex_s< T >::value, bool > | is_nan (const T &value) |
Checks if a floating point number is NaN. More... | |
template<typename T > | |
std::enable_if_t< is_complex_s< T >::value, bool > | is_nan (const T &value) |
Checks if any component of a complex value is NaN. More... | |
template<typename T > | |
constexpr std::enable_if_t<!is_complex_s< T >::value, T > | nan () |
Returns a quiet NaN of the given type. More... | |
template<typename T > | |
constexpr std::enable_if_t< is_complex_s< T >::value, T > | nan () |
Returns a complex with both components quiet NaN. More... | |
template<typename ValueType = default_precision, typename IndexType = int32> | |
matrix_data< ValueType, IndexType > | read_raw (std::istream &is) |
Reads a matrix stored in matrix market format from an input stream. More... | |
template<typename ValueType = default_precision, typename IndexType = int32> | |
matrix_data< ValueType, IndexType > | read_binary_raw (std::istream &is) |
Reads a matrix stored in Ginkgo's binary matrix format from an input stream. More... | |
template<typename ValueType = default_precision, typename IndexType = int32> | |
matrix_data< ValueType, IndexType > | read_generic_raw (std::istream &is) |
Reads a matrix stored in either binary or matrix market format from an input stream. More... | |
template<typename ValueType , typename IndexType > | |
void | write_raw (std::ostream &os, const matrix_data< ValueType, IndexType > &data, layout_type layout=layout_type::coordinate) |
Writes a matrix_data structure to a stream in matrix market format. More... | |
template<typename ValueType , typename IndexType > | |
void | write_binary_raw (std::ostream &os, const matrix_data< ValueType, IndexType > &data) |
Writes a matrix_data structure to a stream in binary format. More... | |
template<typename MatrixType , typename StreamType , typename... MatrixArgs> | |
std::unique_ptr< MatrixType > | read (StreamType &&is, MatrixArgs &&... args) |
Reads a matrix stored in matrix market format from an input stream. More... | |
template<typename MatrixType , typename StreamType , typename... MatrixArgs> | |
std::unique_ptr< MatrixType > | read_binary (StreamType &&is, MatrixArgs &&... args) |
Reads a matrix stored in binary format from an input stream. More... | |
template<typename MatrixType , typename StreamType , typename... MatrixArgs> | |
std::unique_ptr< MatrixType > | read_generic (StreamType &&is, MatrixArgs &&... args) |
Reads a matrix stored either in binary or matrix market format from an input stream. More... | |
template<typename MatrixPtrType , typename StreamType > | |
void | write (StreamType &&os, MatrixPtrType &&matrix, layout_type layout=detail::mtx_io_traits< std::remove_cv_t< detail::pointee< MatrixPtrType >>>::default_layout) |
Writes a matrix into an output stream in matrix market format. More... | |
template<typename MatrixPtrType , typename StreamType > | |
void | write_binary (StreamType &&os, MatrixPtrType &&matrix) |
Writes a matrix into an output stream in binary format. More... | |
template<typename R , typename T > | |
std::unique_ptr< R, std::function< void(R *)> > | copy_and_convert_to (std::shared_ptr< const Executor > exec, T *obj) |
Converts the object to R and places it on Executor exec. More... | |
template<typename R , typename T > | |
std::unique_ptr< const R, std::function< void(const R *)> > | copy_and_convert_to (std::shared_ptr< const Executor > exec, const T *obj) |
Converts the object to R and places it on Executor exec. More... | |
template<typename R , typename T > | |
std::shared_ptr< R > | copy_and_convert_to (std::shared_ptr< const Executor > exec, std::shared_ptr< T > obj) |
Converts the object to R and places it on Executor exec. More... | |
template<typename R , typename T > | |
std::shared_ptr< const R > | copy_and_convert_to (std::shared_ptr< const Executor > exec, std::shared_ptr< const T > obj) |
template<typename ValueType , typename Ptr > | |
detail::temporary_conversion< std::conditional_t< std::is_const< detail::pointee< Ptr > >::value, const matrix::Dense< ValueType >, matrix::Dense< ValueType > > > | make_temporary_conversion (Ptr &&matrix) |
Convert the given LinOp from matrix::Dense<...> to matrix::Dense<ValueType>. More... | |
template<typename ValueType , typename Function , typename... Args> | |
void | precision_dispatch (Function fn, Args *... linops) |
Calls the given function with each given argument LinOp temporarily converted into matrix::Dense<ValueType> as parameters. More... | |
template<typename ValueType , typename Function > | |
void | precision_dispatch_real_complex (Function fn, const LinOp *in, LinOp *out) |
Calls the given function with the given LinOps temporarily converted to matrix::Dense<ValueType>* as parameters. More... | |
template<typename ValueType , typename Function > | |
void | precision_dispatch_real_complex (Function fn, const LinOp *alpha, const LinOp *in, LinOp *out) |
Calls the given function with the given LinOps temporarily converted to matrix::Dense<ValueType>* as parameters. More... | |
template<typename ValueType , typename Function > | |
void | precision_dispatch_real_complex (Function fn, const LinOp *alpha, const LinOp *in, const LinOp *beta, LinOp *out) |
Calls the given function with the given LinOps temporarily converted to matrix::Dense<ValueType>* as parameters. More... | |
template<typename ValueType , typename Function > | |
void | mixed_precision_dispatch (Function fn, const LinOp *in, LinOp *out) |
Calls the given function with each given argument LinOp converted into matrix::Dense<ValueType> as parameters. More... | |
template<typename ValueType , typename Function , std::enable_if_t< is_complex< ValueType >()> * = nullptr> | |
void | mixed_precision_dispatch_real_complex (Function fn, const LinOp *in, LinOp *out) |
Calls the given function with the given LinOps cast to their dynamic type matrix::Dense<ValueType>* as parameters. More... | |
constexpr bool | operator< (const span &first, const span &second) |
constexpr bool | operator<= (const span &first, const span &second) |
constexpr bool | operator> (const span &first, const span &second) |
constexpr bool | operator>= (const span &first, const span &second) |
constexpr bool | operator== (const span &first, const span &second) |
constexpr bool | operator!= (const span &first, const span &second) |
template<typename Accessor > | |
constexpr range< accessor::unary_plus< Accessor > > | operator+ (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::unary_minus< Accessor > > | operator- (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::logical_not< Accessor > > | operator! (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::bitwise_not< Accessor > > | operator~ (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::zero_operation< Accessor > > | zero (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::one_operation< Accessor > > | one (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::abs_operation< Accessor > > | abs (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::real_operation< Accessor > > | real (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::imag_operation< Accessor > > | imag (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::conj_operation< Accessor > > | conj (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::squared_norm_operation< Accessor > > | squared_norm (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::transpose_operation< Accessor > > | transpose (const range< Accessor > &operand) |
template<typename Accessor > | |
constexpr range< accessor::add< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator+ (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::add< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator+ (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::add< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator+ (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::add< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator+ (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::sub< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator- (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::sub< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator- (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::sub< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator- (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::sub< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator- (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::mul< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator* (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::mul< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator* (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::mul< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator* (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::mul< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator* (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::div< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator/ (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::div< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator/ (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::div< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator/ (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::div< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator/ (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::mod< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator% (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::mod< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator% (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::mod< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator% (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::mod< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator% (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::less< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator< (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::less< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator< (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::less< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator< (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::less< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator< (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::greater< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator> (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::greater< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator> (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::greater< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator> (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::greater< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator> (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::less_or_equal< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator<= (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::less_or_equal< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator<= (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::less_or_equal< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator<= (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::less_or_equal< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator<= (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::greater_or_equal< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator>= (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::greater_or_equal< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator>= (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::greater_or_equal< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator>= (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::greater_or_equal< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator>= (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::equal< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator== (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::equal< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator== (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::equal< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator== (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::equal< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator== (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::not_equal< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator!= (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::not_equal< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator!= (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::not_equal< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator!= (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::not_equal< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator!= (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::logical_or< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator|| (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::logical_or< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator|| (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::logical_or< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator|| (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::logical_or< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator|| (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::logical_and< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator&& (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::logical_and< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator&& (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::logical_and< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator&& (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::logical_and< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator&& (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::bitwise_or< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator| (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::bitwise_or< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator| (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::bitwise_or< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator| (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::bitwise_or< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator| (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::bitwise_and< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator& (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::bitwise_and< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator& (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::bitwise_and< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator& (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::bitwise_and< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator& (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::bitwise_xor< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator^ (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::bitwise_xor< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator^ (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::bitwise_xor< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator^ (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::bitwise_xor< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator^ (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::left_shift< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator<< (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::left_shift< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator<< (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::left_shift< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator<< (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::left_shift< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator<< (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::right_shift< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | operator>> (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::right_shift< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | operator>> (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::right_shift< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | operator>> (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::right_shift< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | operator>> (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::max_operation< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | max (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::max_operation< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | max (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::max_operation< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | max (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::max_operation< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | max (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::min_operation< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | min (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::min_operation< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | min (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::min_operation< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | min (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::min_operation< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | min (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Accessor > | |
constexpr range< accessor::mmul_operation< ::gko::detail::operation_kind::range_by_range, Accessor, Accessor > > | mmul (const range< Accessor > &first, const range< Accessor > &second) |
template<typename FirstAccessor , typename SecondAccessor > | |
constexpr range< accessor::mmul_operation< ::gko::detail::operation_kind::range_by_range, FirstAccessor, SecondAccessor > > | mmul (const range< FirstAccessor > &first, const range< SecondAccessor > &second) |
template<typename FirstAccessor , typename SecondOperand > | |
constexpr range< accessor::mmul_operation< ::gko::detail::operation_kind::range_by_scalar, FirstAccessor, SecondOperand > > | mmul (const range< FirstAccessor > &first, const SecondOperand &second) |
template<typename FirstOperand , typename SecondAccessor > | |
constexpr range< accessor::mmul_operation< ::gko::detail::operation_kind::scalar_by_range, FirstOperand, SecondAccessor > > | mmul (const FirstOperand &first, const range< SecondAccessor > &second) |
template<typename Ptr > | |
detail::temporary_clone< detail::pointee< Ptr > > | make_temporary_clone (std::shared_ptr< const Executor > exec, Ptr &&ptr) |
Creates a temporary_clone. More... | |
template<typename Ptr > | |
detail::temporary_clone< detail::pointee< Ptr > > | make_temporary_output_clone (std::shared_ptr< const Executor > exec, Ptr &&ptr) |
Creates a uninitialized temporary_clone that will be copied back to the input afterwards. More... | |
constexpr bool | operator== (precision_reduction x, precision_reduction y) noexcept |
Checks if two precision_reduction encodings are equal. More... | |
constexpr bool | operator!= (precision_reduction x, precision_reduction y) noexcept |
Checks if two precision_reduction encodings are different. More... | |
template<typename IndexType > | |
constexpr IndexType | invalid_index () |
Value for an invalid signed index type. | |
template<typename Pointer > | |
detail::cloned_type< Pointer > | clone (const Pointer &p) |
Creates a unique clone of the object pointed to by p . More... | |
template<typename Pointer > | |
detail::cloned_type< Pointer > | clone (std::shared_ptr< const Executor > exec, const Pointer &p) |
Creates a unique clone of the object pointed to by p on Executor exec . More... | |
template<typename OwningPointer > | |
detail::shared_type< OwningPointer > | share (OwningPointer &&p) |
Marks the object pointed to by p as shared. More... | |
template<typename OwningPointer > | |
std::remove_reference< OwningPointer >::type && | give (OwningPointer &&p) |
Marks that the object pointed to by p can be given to the callee. More... | |
template<typename Pointer > | |
std::enable_if< detail::have_ownership_s< Pointer >::value, detail::pointee< Pointer > * >::type | lend (const Pointer &p) |
Returns a non-owning (plain) pointer to the object pointed to by p . More... | |
template<typename Pointer > | |
std::enable_if<!detail::have_ownership_s< Pointer >::value, detail::pointee< Pointer > * >::type | lend (const Pointer &p) |
Returns a non-owning (plain) pointer to the object pointed to by p . More... | |
template<typename T , typename U > | |
std::decay_t< T > * | as (U *obj) |
Performs polymorphic type conversion. More... | |
template<typename T , typename U > | |
const std::decay_t< T > * | as (const U *obj) |
Performs polymorphic type conversion. More... | |
template<typename T , typename U > | |
std::decay_t< T > * | as (ptr_param< U > obj) |
Performs polymorphic type conversion on a ptr_param. More... | |
template<typename T , typename U > | |
const std::decay_t< T > * | as (ptr_param< const U > obj) |
Performs polymorphic type conversion. More... | |
template<typename T , typename U > | |
std::unique_ptr< std::decay_t< T > > | as (std::unique_ptr< U > &&obj) |
Performs polymorphic type conversion of a unique_ptr. More... | |
template<typename T , typename U > | |
std::shared_ptr< std::decay_t< T > > | as (std::shared_ptr< U > obj) |
Performs polymorphic type conversion of a shared_ptr. More... | |
template<typename T , typename U > | |
std::shared_ptr< const std::decay_t< T > > | as (std::shared_ptr< const U > obj) |
Performs polymorphic type conversion of a shared_ptr. More... | |
bool | operator== (const version &first, const version &second) |
bool | operator!= (const version &first, const version &second) |
bool | operator< (const version &first, const version &second) |
bool | operator<= (const version &first, const version &second) |
bool | operator> (const version &first, const version &second) |
bool | operator>= (const version &first, const version &second) |
std::ostream & | operator<< (std::ostream &os, const version &ver) |
Prints version information to a stream. More... | |
std::ostream & | operator<< (std::ostream &os, const version_info &ver_info) |
Prints library version information in human-readable format to a stream. More... | |
template<template< typename, typename > class MatrixType, typename... Args> | |
auto | with_matrix_type (Args &&... create_args) |
This function returns a type that delays a call to MatrixType::create. More... | |
template<typename VecPtr > | |
std::unique_ptr< matrix::Dense< typename detail::pointee< VecPtr >::value_type > > | make_dense_view (VecPtr &&vector) |
Creates a view of a given Dense vector. More... | |
template<typename VecPtr > | |
std::unique_ptr< const matrix::Dense< typename detail::pointee< VecPtr >::value_type > > | make_const_dense_view (VecPtr &&vector) |
Creates a view of a given Dense vector. More... | |
template<typename Matrix , typename... TArgs> | |
std::unique_ptr< Matrix > | initialize (size_type stride, std::initializer_list< typename Matrix::value_type > vals, std::shared_ptr< const Executor > exec, TArgs &&... create_args) |
Creates and initializes a column-vector. More... | |
template<typename Matrix , typename... TArgs> | |
std::unique_ptr< Matrix > | initialize (std::initializer_list< typename Matrix::value_type > vals, std::shared_ptr< const Executor > exec, TArgs &&... create_args) |
Creates and initializes a column-vector. More... | |
template<typename Matrix , typename... TArgs> | |
std::unique_ptr< Matrix > | initialize (size_type stride, std::initializer_list< std::initializer_list< typename Matrix::value_type >> vals, std::shared_ptr< const Executor > exec, TArgs &&... create_args) |
Creates and initializes a matrix. More... | |
template<typename Matrix , typename... TArgs> | |
std::unique_ptr< Matrix > | initialize (std::initializer_list< std::initializer_list< typename Matrix::value_type >> vals, std::shared_ptr< const Executor > exec, TArgs &&... create_args) |
Creates and initializes a matrix. More... | |
bool | operator== (const stopping_status &x, const stopping_status &y) noexcept |
Checks if two stopping statuses are equivalent. More... | |
bool | operator!= (const stopping_status &x, const stopping_status &y) noexcept |
Checks if two stopping statuses are different. More... | |
Variables | |
constexpr allocation_mode | default_cuda_alloc_mode |
constexpr allocation_mode | default_hip_alloc_mode |
constexpr size_type | byte_size = CHAR_BIT |
Number of bits in a byte. | |
The Ginkgo namespace.
using gko::highest_precision = typedef typename detail::highest_precision_variadic<Ts...>::type |
Obtains the smallest arithmetic type that is able to store elements of all template parameter types exactly.
All template type parameters need to be either real or complex types, mixing them is not possible.
Formally, it computes a right-fold over the type list, with the highest precision of a pair of real arithmetic types T1, T2 computed as decltype(T1{} + T2{})
, or std::complex<highest_precision<remove_complex<T1>, remove_complex<T2>>>
for complex types.
using gko::is_complex_or_scalar_s = typedef detail::is_complex_or_scalar_impl<T> |
Allows to check if T is a complex or scalar value during compile time by accessing the value
attribute of this struct.
If value
is true
, T is a complex/scalar type, if it is false
, T is not a complex/scalar type.
T | type to check |
using gko::is_complex_s = typedef detail::is_complex_impl<T> |
Allows to check if T is a complex value during compile time by accessing the value
attribute of this struct.
If value
is true
, T is a complex type, if it is false
, T is not a complex type.
T | type to check |
using gko::previous_precision = typedef next_precision<T> |
Obtains the previous type in the singly-linked precision list.
using gko::remove_complex = typedef typename detail::remove_complex_s<T>::type |
Obtain the type which removed the complex of complex/scalar type or the template parameter of class by accessing the type
attribute of this struct.
T | type to remove complex |
using gko::to_complex = typedef typename detail::to_complex_s<T>::type |
Obtain the type which adds the complex of complex/scalar type or the template parameter of class by accessing the type
attribute of this struct.
T | type to complex_type |
using gko::to_real = typedef remove_complex<T> |
to_real is alias of remove_complex
T | type to real |
|
strong |
Specify the mode of allocation for CUDA/HIP GPUs.
device
allocates memory on the device and Unified Memory model is not used.
unified_global
allocates memory on the device, but is accessible by the host through the Unified memory model.
unified_host
allocates memory on the host and it is not available on devices which do not have concurrent accesses switched on, but this access can be explicitly switched on, when necessary.
|
strong |
|
strong |
How Logger events are propagated to their Executor.
Enumerator | |
---|---|
never | Events only get reported at loggers attached to the triggering object. (Except for allocation/free, copy and Operations, since they happen at the Executor). |
automatic | Events get reported to loggers attached to the triggering object and propagating loggers (Logger::needs_propagation() return true) attached to its Executor. |
|
inlineconstexpr |
Returns the absolute value of the object.
T | the type of the object |
x | the object |
Referenced by is_finite().
|
inline |
Performs polymorphic type conversion.
This is the constant version of the function.
T | requested result type |
U | static type of the passed object |
obj | the object which should be converted |
|
inline |
Performs polymorphic type conversion.
This is the constant version of the function.
T | requested result type |
U | static type of the passed object |
obj | the object which should be converted |
References gko::ptr_param< T >::get().
|
inline |
Performs polymorphic type conversion on a ptr_param.
T | requested result type |
U | static type of the passed object |
obj | the object which should be converted |
References gko::ptr_param< T >::get().
|
inline |
Performs polymorphic type conversion of a shared_ptr.
This is the constant version of the function.
T | requested result type |
U | static type of the passed object |
obj | the shared_ptr to the object which should be converted. |
|
inline |
Performs polymorphic type conversion of a shared_ptr.
T | requested result type |
U | static type of the passed object |
obj | the shared_ptr to the object which should be converted. |
|
inline |
Performs polymorphic type conversion of a unique_ptr.
T | requested result type |
U | static type of the passed object |
obj | the unique_ptr to the object which should be converted. If successful, it will be reset to a nullptr. |
|
inline |
Performs polymorphic type conversion.
T | requested result type |
U | static type of the passed object |
obj | the object which should be converted |
Referenced by gko::preconditioner::Isai< IsaiType, ValueType, IndexType >::get_approximate_inverse().
Performs integer division with rounding up.
num | numerator |
den | denominator |
Referenced by gko::matrix::Csr< ValueType, IndexType >::load_balance::clac_size(), gko::preconditioner::block_interleaved_storage_scheme< index_type >::compute_storage_space(), and gko::matrix::Csr< ValueType, IndexType >::load_balance::process().
|
inline |
Creates a unique clone of the object pointed to by p
.
The pointee (i.e. *p
) needs to have a clone method that returns a std::unique_ptr in order for this method to work.
Pointer | type of pointer to the object (plain or smart pointer) |
p | a pointer to the object |
Referenced by gko::preconditioner::Ic< LSolverType, IndexType >::operator=(), gko::preconditioner::Ilu< LSolverType, USolverType, ReverseApply, IndexType >::operator=(), gko::solver::EnablePreconditionable< Bicg< ValueType > >::set_preconditioner(), and gko::solver::EnableIterativeBase< Bicg< ValueType > >::set_stop_criterion_factory().
|
inline |
Creates a unique clone of the object pointed to by p
on Executor exec
.
The pointee (i.e. *p
) needs to have a clone method that takes an executor and returns a std::unique_ptr in order for this method to work.
Pointer | type of pointer to the object (plain or smart pointer) |
exec | the executor where the cloned object should be stored |
p | a pointer to the object |
|
inlineconstexpr |
Returns the conjugate of an object.
x | the number to conjugate |
Referenced by squared_norm().
std::unique_ptr<const R, std::function<void(const R*)> > gko::copy_and_convert_to | ( | std::shared_ptr< const Executor > | exec, |
const T * | obj | ||
) |
Converts the object to R and places it on Executor exec.
If the object is already of the requested type and on the requested executor, the copy and conversion is avoided and a reference to the original object is returned instead.
R | the type to which the object should be converted |
T | the type of the input object |
exec | the executor where the result should be placed |
obj | the object that should be converted |
std::shared_ptr<const R> gko::copy_and_convert_to | ( | std::shared_ptr< const Executor > | exec, |
std::shared_ptr< const T > | obj | ||
) |
This is the version that takes in the std::shared_ptr and returns a std::shared_ptr
If the object is already of the requested type and on the requested executor, the copy and conversion is avoided and a reference to the original object is returned instead.
R | the type to which the object should be converted |
T | the type of the input object |
exec | the executor where the result should be placed |
obj | the object that should be converted |
std::shared_ptr<R> gko::copy_and_convert_to | ( | std::shared_ptr< const Executor > | exec, |
std::shared_ptr< T > | obj | ||
) |
Converts the object to R and places it on Executor exec.
This is the version that takes in the std::shared_ptr and returns a std::shared_ptr
If the object is already of the requested type and on the requested executor, the copy and conversion is avoided and a reference to the original object is returned instead.
R | the type to which the object should be converted |
T | the type of the input object |
exec | the executor where the result should be placed |
obj | the object that should be converted |
std::unique_ptr<R, std::function<void(R*)> > gko::copy_and_convert_to | ( | std::shared_ptr< const Executor > | exec, |
T * | obj | ||
) |
Converts the object to R and places it on Executor exec.
If the object is already of the requested type and on the requested executor, the copy and conversion is avoided and a reference to the original object is returned instead.
R | the type to which the object should be converted |
T | the type of the input object |
exec | the executor where the result should be placed |
obj | the object that should be converted |
|
constexprnoexcept |
Returns the position of the most significant bit of the number.
This is the same as the rounded down base-2 logarithm of the number.
T | a numeric type supporting bit shift and comparison |
n | a number |
hint | a lower bound for the position o the significant bit |
hint
and the significant bit position of n
|
constexprnoexcept |
Returns the smallest power of base
not smaller than limit
.
T | a numeric type supporting multiplication and comparison |
base | the base of the power to be returned |
limit | the lower limit on the size of the power returned |
hint | a lower bound on the result, has to be a power of base |
base
not smaller than limit
|
inline |
Marks that the object pointed to by p
can be given to the callee.
Effectively calls std::move(p)
.
OwningPointer | type of pointer with ownership to the object (has to be a smart pointer) |
p | a pointer to the object |
p
becomes invalid after this call.
|
inlineconstexpr |
Returns the imaginary part of the object.
T | type of the object |
x | the object |
|
inlineconstexpr |
Checks if T is a complex type.
T | type to check |
true
if T is a complex type, false
otherwise
|
inlineconstexpr |
Checks if T is a complex/scalar type.
T | type to check |
true
if T is a complex/scalar type, false
otherwise
|
inline |
Checks if a floating point number is finite, meaning it is neither +/- infinity nor NaN.
T | type of the value to check |
value | value to check |
true
if the value is finite, meaning it are neither +/- infinity nor NaN. References abs().
Referenced by is_finite().
|
inline |
Checks if all components of a complex value are finite, meaning they are neither +/- infinity nor NaN.
T | complex type of the value to check |
value | complex value to check |
true
if both components of the given value are finite, meaning they are neither +/- infinity nor NaN. References is_finite().
|
inline |
Checks if a floating point number is NaN.
T | type of the value to check |
value | value to check |
true
if the value is NaN.
|
inline |
Checks if any component of a complex value is NaN.
T | complex type of the value to check |
value | complex value to check |
true
if any component of the given value is NaN.
|
inlineconstexpr |
Returns true if and only if the given value is not zero.
T | the type of the value |
value | the given value |
value != zero<T>()
Referenced by gko::matrix_data< ValueType, IndexType >::diag().
|
inlineconstexpr |
Returns true if and only if the given value is zero.
T | the type of the value |
value | the given value |
value == zero<T>()
Referenced by gko::matrix_data< ValueType, IndexType >::remove_zeros().
|
inline |
Returns a non-owning (plain) pointer to the object pointed to by p
.
Pointer | type of pointer to the object (plain or smart pointer) |
p | a pointer to the object |
|
inline |
Returns a non-owning (plain) pointer to the object pointed to by p
.
Pointer | type of pointer to the object (plain or smart pointer) |
p | a pointer to the object |
p
. array<ValueType> gko::make_array_view | ( | std::shared_ptr< const Executor > | exec, |
size_type | size, | ||
ValueType * | data | ||
) |
Helper function to create an array view deducing the value type.
exec | the executor on which the array resides |
size | the number of elements for the array |
data | the pointer to the array we create a view on. |
ValueType | the type of the array elements |
array<ValueType>::view(exec, size, data)
References make_array_view().
Referenced by make_array_view().
detail::const_array_view<ValueType> gko::make_const_array_view | ( | std::shared_ptr< const Executor > | exec, |
size_type | size, | ||
const ValueType * | data | ||
) |
Helper function to create a const array view deducing the value type.
exec | the executor on which the array resides |
size | the number of elements for the array |
data | the pointer to the array we create a view on. |
ValueType | the type of the array elements |
array<ValueType>::const_view(exec, size, data)
References make_const_array_view().
Referenced by make_const_array_view().
std::unique_ptr< const matrix::Dense<typename detail::pointee<VecPtr>::value_type> > gko::make_const_dense_view | ( | VecPtr && | vector | ) |
Creates a view of a given Dense vector.
VecPtr | a (smart or raw) pointer to the vector. |
vector | the vector on which to create the view |
References gko::matrix::Dense< ValueType >::create_const_view_of().
std::unique_ptr<matrix::Dense<typename detail::pointee<VecPtr>::value_type> > gko::make_dense_view | ( | VecPtr && | vector | ) |
Creates a view of a given Dense vector.
VecPtr | a (smart or raw) pointer to the vector. |
vector | the vector on which to create the view |
References gko::matrix::Dense< ValueType >::create_view_of().
detail::temporary_clone<detail::pointee<Ptr> > gko::make_temporary_clone | ( | std::shared_ptr< const Executor > | exec, |
Ptr && | ptr | ||
) |
Creates a temporary_clone.
This is a helper function which avoids the need to explicitly specify the type of the object, as would be the case if using the constructor of temporary_clone.
exec | the executor where the clone will be created |
ptr | a pointer to the object of which the clone will be created |
Ptr | the (raw or smart) pointer type to be temporarily cloned |
Referenced by gko::ScaledIdentityAddable::add_scaled_identity(), gko::LinOp::apply(), gko::matrix::Csr< ValueType, IndexType >::inv_scale(), and gko::matrix::Csr< ValueType, IndexType >::scale().
detail::temporary_conversion<std::conditional_t< std::is_const<detail::pointee<Ptr> >::value, const matrix::Dense<ValueType>, matrix::Dense<ValueType> > > gko::make_temporary_conversion | ( | Ptr && | matrix | ) |
Convert the given LinOp from matrix::Dense<...> to matrix::Dense<ValueType>.
The conversion tries to convert the input LinOp to all Dense types with value type recursively reachable by next_precision<...> starting from the ValueType template parameter. This means that all real-to-real and complex-to-complex conversions for default precisions are being considered. If the input matrix is non-const, the contents of the modified converted object will be converted back to the input matrix when the returned object is destroyed. This may lead to a loss of precision!
matrix | the input matrix which is supposed to be converted. It is wrapped unchanged if it is already of type matrix::Dense<ValueType>, otherwise it will be converted to this type if possible. |
NotSupported | if the input matrix cannot be converted to matrix::Dense<ValueType> |
ValueType | the value type into whose associated matrix::Dense type to convert the input LinOp. |
detail::temporary_clone<detail::pointee<Ptr> > gko::make_temporary_output_clone | ( | std::shared_ptr< const Executor > | exec, |
Ptr && | ptr | ||
) |
Creates a uninitialized temporary_clone that will be copied back to the input afterwards.
It can be used for output parameters to avoid an unnecessary copy in make_temporary_clone.
This is a helper function which avoids the need to explicitly specify the type of the object, as would be the case if using the constructor of temporary_clone.
exec | the executor where the uninitialized clone will be created |
ptr | a pointer to the object of which the clone will be created |
Ptr | the (raw or smart) pointer type to be temporarily cloned |
|
inlineconstexpr |
Returns the larger of the arguments.
T | type of the arguments |
x | first argument |
y | second argument |
|
inlineconstexpr |
Returns the smaller of the arguments.
T | type of the arguments |
x | first argument |
y | second argument |
Referenced by gko::matrix::Csr< ValueType, IndexType >::load_balance::clac_size().
void gko::mixed_precision_dispatch | ( | Function | fn, |
const LinOp * | in, | ||
LinOp * | out | ||
) |
Calls the given function with each given argument LinOp converted into matrix::Dense<ValueType> as parameters.
If GINKGO_MIXED_PRECISION is defined, this means that the function will be called with its dynamic type as a static type, so the (templated/generic) function will be instantiated with all pairs of Dense<ValueType> and Dense<next_precision<ValueType>> parameter types, and the appropriate overload will be called based on the dynamic type of the parameter.
If GINKGO_MIXED_PRECISION is not defined, it will behave exactly like precision_dispatch.
fn | the given function. It will be called with one const and one non-const matrix::Dense<...> parameter based on the dynamic type of the inputs (GINKGO_MIXED_PRECISION) or of type matrix::Dense<ValueType> (no GINKGO_MIXED_PRECISION). |
in | The first parameter to be cast (GINKGO_MIXED_PRECISION) or converted (no GINKGO_MIXED_PRECISION) and used to call fn . |
out | The second parameter to be cast (GINKGO_MIXED_PRECISION) or converted (no GINKGO_MIXED_PRECISION) and used to call fn . |
ValueType | the value type to use for the parameters of fn (no GINKGO_MIXED_PRECISION). With GINKGO_MIXED_PRECISION enabled, it only matters whether this type is complex or real. |
Function | the function pointer, lambda or other functor type to call with the converted arguments. |
void gko::mixed_precision_dispatch_real_complex | ( | Function | fn, |
const LinOp * | in, | ||
LinOp * | out | ||
) |
Calls the given function with the given LinOps cast to their dynamic type matrix::Dense<ValueType>* as parameters.
If ValueType is real and both in
and out
are complex, uses matrix::Dense::get_real_view() to convert them into real matrices after precision conversion.
|
inlineconstexpr |
Returns a quiet NaN of the given type.
T | the type of the object |
Referenced by nan().
|
inlineconstexpr |
Returns a complex with both components quiet NaN.
T | the type of the object |
References nan().
|
inlineconstexpr |
Returns the multiplicative identity for T.
Referenced by unit_root().
|
inlineconstexpr |
Returns the multiplicative identity for T.
one<decltype(x)>()
. Instead, it allows one(x)
.
|
inlinenoexcept |
Checks if two stopping statuses are different.
x | a stopping status |
y | a stopping status |
!(x == y)
|
constexprnoexcept |
Checks if two precision_reduction encodings are different.
x | an encoding |
y | an encoding |
x
and y
are different encodings.
|
inline |
Prints version information to a stream.
os | output stream |
ver | version structure |
References gko::version::major, gko::version::minor, gko::version::patch, and gko::version::tag.
std::ostream& gko::operator<< | ( | std::ostream & | os, |
const version_info & | ver_info | ||
) |
Prints library version information in human-readable format to a stream.
os | output stream |
ver_info | version information |
|
inlinenoexcept |
Checks if two stopping statuses are equivalent.
x | a stopping status |
y | a stopping status |
x
and y
have the same mask and converged and finalized state
|
constexprnoexcept |
Checks if two precision_reduction encodings are equal.
x | an encoding |
y | an encoding |
x
and y
are the same encodings
|
inlineconstexpr |
Returns the value of pi.
T | the value type to return |
void gko::precision_dispatch | ( | Function | fn, |
Args *... | linops | ||
) |
Calls the given function with each given argument LinOp temporarily converted into matrix::Dense<ValueType> as parameters.
fn | the given function. It will be passed one (potentially const) matrix::Dense<ValueType>* parameter per parameter in the parameter pack linops . |
linops | the given arguments to be converted and passed on to fn. |
ValueType | the value type to use for the parameters of fn . |
Function | the function pointer, lambda or other functor type to call with the converted arguments. |
Args | the argument type list. |
void gko::precision_dispatch_real_complex | ( | Function | fn, |
const LinOp * | alpha, | ||
const LinOp * | in, | ||
const LinOp * | beta, | ||
LinOp * | out | ||
) |
Calls the given function with the given LinOps temporarily converted to matrix::Dense<ValueType>* as parameters.
If ValueType is real and both in
and out
are complex, uses matrix::Dense::get_real_view() to convert them into real matrices after precision conversion.
void gko::precision_dispatch_real_complex | ( | Function | fn, |
const LinOp * | alpha, | ||
const LinOp * | in, | ||
LinOp * | out | ||
) |
Calls the given function with the given LinOps temporarily converted to matrix::Dense<ValueType>* as parameters.
If ValueType is real and both in
and out
are complex, uses matrix::Dense::get_real_view() to convert them into real matrices after precision conversion.
void gko::precision_dispatch_real_complex | ( | Function | fn, |
const LinOp * | in, | ||
LinOp * | out | ||
) |
Calls the given function with the given LinOps temporarily converted to matrix::Dense<ValueType>* as parameters.
If ValueType is real and both input vectors are complex, uses matrix::Dense::get_real_view() to convert them into real matrices after precision conversion.
|
inline |
Reads a matrix stored in matrix market format from an input stream.
MatrixType | a ReadableFromMatrixData LinOp type used to store the matrix once it's been read from disk. |
StreamType | type of stream used to write the data to |
MatrixArgs | additional argument types passed to MatrixType constructor |
is | input stream from which to read the data |
args | additional arguments passed to MatrixType constructor |
References read_raw().
Referenced by gko::ReadableFromMatrixData< ValueType, int32 >::read().
|
inline |
Reads a matrix stored in binary format from an input stream.
MatrixType | a ReadableFromMatrixData LinOp type used to store the matrix once it's been read from disk. |
StreamType | type of stream used to write the data to |
MatrixArgs | additional argument types passed to MatrixType constructor |
is | input stream from which to read the data |
args | additional arguments passed to MatrixType constructor |
References read_binary_raw().
matrix_data<ValueType, IndexType> gko::read_binary_raw | ( | std::istream & | is | ) |
Reads a matrix stored in Ginkgo's binary matrix format from an input stream.
Note that this format depends on the processor's endianness, so files from a big endian processor can't be read from a little endian processor and vice-versa.
The binary format has the following structure (in system endianness):
ValueType | type of matrix values |
IndexType | type of matrix indexes |
is | input stream from which to read the data |
Referenced by read_binary().
|
inline |
Reads a matrix stored either in binary or matrix market format from an input stream.
MatrixType | a ReadableFromMatrixData LinOp type used to store the matrix once it's been read from disk. |
StreamType | type of stream used to write the data to |
MatrixArgs | additional argument types passed to MatrixType constructor |
is | input stream from which to read the data |
args | additional arguments passed to MatrixType constructor |
References read_generic_raw().
matrix_data<ValueType, IndexType> gko::read_generic_raw | ( | std::istream & | is | ) |
Reads a matrix stored in either binary or matrix market format from an input stream.
ValueType | type of matrix values |
IndexType | type of matrix indexes |
is | input stream from which to read the data |
Referenced by read_generic().
matrix_data<ValueType, IndexType> gko::read_raw | ( | std::istream & | is | ) |
Reads a matrix stored in matrix market format from an input stream.
ValueType | type of matrix values |
IndexType | type of matrix indexes |
is | input stream from which to read the data |
Referenced by read().
|
inlineconstexpr |
Returns the real part of the object.
T | type of the object |
x | the object |
Referenced by squared_norm().
void gko::reduce_add | ( | const array< ValueType > & | input_arr, |
array< ValueType > & | result | ||
) |
Reduce (sum) the values in the array.
The | type of the input data |
[in] | input_arr | the input array to be reduced |
[in,out] | result | the reduced value. The result is written into the first entry and the value in the first entry is used as the initial value for the reduce. |
ValueType gko::reduce_add | ( | const array< ValueType > & | input_arr, |
const ValueType | init_val = 0 |
||
) |
Reduce (sum) the values in the array.
The | type of the input data |
[in] | input_arr | the input array to be reduced |
[in] | init_val | the initial value |
|
inlineconstexpr |
Reduces the precision of the input parameter.
T | the original precision |
val | the value to round down |
|
inlineconstexpr |
Increases the precision of the input parameter.
T | the original precision |
val | the value to round up |
|
inline |
Computes the quotient of the given parameters, guarding against division by zero.
T | value type of the parameters |
a | the dividend |
b | the divisor |
a / b
if b is non-zero, zero otherwise.
|
inline |
Marks the object pointed to by p
as shared.
Effectively converts a pointer with ownership to std::shared_ptr.
OwningPointer | type of pointer with ownership to the object (has to be a smart pointer) |
p | a pointer to the object. It must be a temporary or explicitly marked movable (rvalue reference). |
p
becomes invalid after this call. Referenced by gko::preconditioner::Ic< LSolverType, IndexType >::conj_transpose(), gko::preconditioner::Ilu< LSolverType, USolverType, ReverseApply, IndexType >::conj_transpose(), gko::preconditioner::Ic< LSolverType, IndexType >::transpose(), and gko::preconditioner::Ilu< LSolverType, USolverType, ReverseApply, IndexType >::transpose().
|
inline |
Returns a batch_dim object with its dimensions swapped for batched operators.
DimensionType | datatype used to represent each dimension |
dimensions | original object |
References gko::batch_dim< Dimensionality, DimensionType >::get_common_size(), and gko::batch_dim< Dimensionality, DimensionType >::get_num_batch_items().
Referenced by gko::preconditioner::Ic< LSolverType, IndexType >::conj_transpose(), gko::preconditioner::Ilu< LSolverType, USolverType, ReverseApply, IndexType >::conj_transpose(), gko::preconditioner::Ic< LSolverType, IndexType >::transpose(), and gko::preconditioner::Ilu< LSolverType, USolverType, ReverseApply, IndexType >::transpose().
|
inlineconstexpr |
Returns the value of exp(2 * pi * i * k / n), i.e.
an nth root of unity.
n | the denominator of the argument |
k | the numerator of the argument. Defaults to 1. |
T | the corresponding real value type. |
References one().
auto gko::with_matrix_type | ( | Args &&... | create_args | ) |
This function returns a type that delays a call to MatrixType::create.
It can be used to set the used value and index type, as well as the executor at a later stage.
For example, the following code creates first a temporary object, which is then used later to construct an operator of the previously defined base type:
MatrixType | A template type that accepts two types, the first one will be set to the value type, the second one to the index type. |
Args | Types of the arguments passed to MatrixType::create. |
create_args | arguments that will be forwarded to MatrixType::create |
create<value_type, index_type>(executor)
.
|
inline |
Writes a matrix into an output stream in matrix market format.
MatrixPtrType | a (smart or raw) pointer to a WritableToMatrixData object providing data to be written. |
StreamType | type of stream used to write the data to |
os | output stream where the data is to be written |
matrix | the matrix to write |
layout | the layout used in the output |
References write_raw().
|
inline |
Writes a matrix into an output stream in binary format.
Note that this format depends on the processor's endianness, so files from a big endian processor can't be read from a little endian processor and vice-versa.
MatrixPtrType | a (smart or raw) pointer to a WritableToMatrixData object providing data to be written. |
StreamType | type of stream used to write the data to |
os | output stream where the data is to be written |
matrix | the matrix to write |
References write_binary_raw().
void gko::write_binary_raw | ( | std::ostream & | os, |
const matrix_data< ValueType, IndexType > & | data | ||
) |
Writes a matrix_data structure to a stream in binary format.
Note that this format depends on the processor's endianness, so files from a big endian processor can't be read from a little endian processor and vice-versa.
ValueType | type of matrix values |
IndexType | type of matrix indexes |
os | output stream where the data is to be written |
data | the matrix data to write |
Referenced by write_binary().
void gko::write_raw | ( | std::ostream & | os, |
const matrix_data< ValueType, IndexType > & | data, | ||
layout_type | layout = layout_type::coordinate |
||
) |
Writes a matrix_data structure to a stream in matrix market format.
ValueType | type of matrix values |
IndexType | type of matrix indexes |
os | output stream where the data is to be written |
data | the matrix data to write |
layout | the layout used in the output |
Referenced by write().
|
inlineconstexpr |
Returns the additive identity for T.
Referenced by gko::matrix::Hybrid< ValueType, IndexType >::strategy_type::strategy_type().
|
inlineconstexpr |
Returns the additive identity for T.
zero<decltype(x)>()
. Instead, it allows zero(x)
.
|
constexpr |
|
constexpr |