ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_matrixMulAdd::setMatMulAdd Interface Reference

Return the result of the multiplication of the input matrices/vector matA and matB in the user-specified form.
More...

Detailed Description

Return the result of the multiplication of the input matrices/vector matA and matB in the user-specified form.

Return the result of the multiplication of the input matrices matA and matB in one of the following forms,

\begin{align*} & \ms{matC} \leftarrow \alpha \ms{matA}~~ \ms{matB} + \beta \ms{matC} && \ms{matC} \leftarrow \alpha \ms{matA}~~ \ms{matB}^T + \beta \ms{matC} && \ms{matC} \leftarrow \alpha \ms{matA}~~ \ms{matB}^H + \beta \ms{matC} \\ & \ms{matC} \leftarrow \alpha \ms{matA}^T \ms{matB} + \beta \ms{matC} && \ms{matC} \leftarrow \alpha \ms{matA}^T \ms{matB}^T + \beta \ms{matC} && \ms{matC} \leftarrow \alpha \ms{matA}^T \ms{matB}^H + \beta \ms{matC} \\ & \ms{matC} \leftarrow \alpha \ms{matA}^H \ms{matB} + \beta \ms{matC} && \ms{matC} \leftarrow \alpha \ms{matA}^H \ms{matB}^T + \beta \ms{matC} && \ms{matC} \leftarrow \alpha \ms{matA}^H \ms{matB}^H + \beta \ms{matC} \end{align*}

where \(\cdot^T\) represents a Symmetric transpose, and \(\cdot^H\) represents a Hermitian transpose.

The following figure illustrates the form of the general matrix-matrix-matrix or matrix-vector multiplication depending on the input values.

  • Default Multiplication (no transposition involved).

  • Same as the default but with operationA set to transSymm or transHerm.

  • Same as the default but with operationB set to transSymm or transHerm.
Note
For triangular matrix-matrix or matrix-vector multiplications, see pm_matrixMulTri.
Parameters
[in]matA: The input contiguous matrix of either shape (:,:) or (:) (in Linear Full Packing (LFP) known as the LAPACK packed storage format) of either,
  1. type integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64),
  2. type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128),
  3. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the subset of matrix \(\ms{matA}\) to be used in the matrix-matrix or matrix-vector multiplication.
  1. If the input argument packA is present, then matA must be a vector of shape (1 : ndim * (ndim + 1) / 2) containing the upper/lower triangle of a Symmetric/Hermitian matrix in Linear Full Packing (LFP) format (also known as the LAPACK packed storage format) to be used in the multiplication.
  2. If matA is a matrix (of shape (:,:)), then matA is assumed to contain a general full matrix or a subset of a Symmetric/Hermitian matrix in regular (Rectangular Default Packing (RDP) format with the optionally-specified subsetA format.
[in]classA: The input scalar that can be either,
  1. the constant symmetric for matA of any type of any kind, implying that matA is a Symmetric matrix with the specified subsetA.
  2. the constant hermitian for matA of any type of any kind, implying that matA is a Hermitian (conjugate) matrix with the specified subsetA.
This argument is merely a convenience to resolve the different procedure functionalities within this generic interface.
(optional, default = general matrix. It must be present if and only if the input argument subsetA is also present.)
[in]subsetA: The input scalar that can be either,
  1. the constant uppDia implying that the upper-diagonal triangular subset of the input (Symmetric/Hermitian) matA should be used in the multiplication without referencing the lower triangle.
  2. the constant lowDia implying that the lower-diagonal triangular subset of the input (Symmetric/Hermitian) matA should be used in the multiplication without referencing the upper triangle.
This optional argument is merely a convenience to resolve the different procedure functionalities within this generic interface.
If missing, the input matA is assumed to a generic matrix whose full specified subset will be used in multiplication.
(optional. It can must be present if and only if the input argument classA is present.)
[in]packA: The input scalar that can be:
  1. The constant lfpack signifying the Linear Full Packing (LFP) format (also known as the LAPACK packed storage format of the input subset of the matrix represented by matA.
(optional, default = rdpack (Rectangular Default Packing). It must be present if and only if the input arguments classA and subsetA are both present.)
[in]operationA: The input scalar that can be either,
  1. the constant transSymm for matA of any type of any kind, implying that a Symmetric transpose of the specified subset of matA is to be used in the matrix-matrix or matrix-vector multiplication.
  2. the constant transHerm for matA of any type of any kind, implying that a Hermitian transpose of the specified subset of matA is to be used in the matrix-matrix or matrix-vector multiplication.
Specifying this argument changes the shape of the subset of matA used in the matrix-matrix or matrix-vector multiplication.
See the description of the input argument ndum.
This argument is merely a convenience to resolve the different procedure functionalities within this generic interface.
(optional. If missing, the specified subset of matA will be used as is, without any operations performed on it.)
[in]matB: The input contiguous vector of shape (:) or matrix of shape (:,:) of the same type and kind as matA, containing the subset of vector/matrix \(\ms{matB}\) in the matrix-matrix or matrix-vector multiplication.
[in]classB: The input scalar that can be either,
  1. the constant symmetric for matB of any type of any kind, implying that matB is a Symmetric matrix with the specified subsetB.
  2. the constant hermitian for matB of any type of any kind, implying that matB is a Hermitian (conjugate) matrix with the specified subsetB.
This argument is merely a convenience to resolve the different procedure functionalities within this generic interface.
(optional. It must be present if and only if the input argument subsetB is also present.)
[in]subsetB: The input scalar that can be either,
  1. the constant uppDia implying that the upper-diagonal triangular block of the input Symmetric/Hermitian matB should be used in the multiplication without referencing the lower triangle.
  2. the constant lowDia implying that the lower-diagonal triangular block of the input Symmetric/Hermitian matB should be used in the multiplication without referencing the upper triangle.
This optional argument is merely a convenience to resolve the different procedure functionalities within this generic interface.
If missing, the input matB is assumed to a generic matrix whose full specified subset will be used in multiplication.
(optional. It can be present only if either classB is present.)
[in]operationB: The input scalar that can be either,
  1. the constant transSymm for matB of any type of any kind, implying that a Symmetric transpose of the specified subset of matB is to be used in the matrix-matrix or matrix-vector multiplication.
  2. the constant transHerm for matB of any type of any kind, implying that a Hermitian transpose of the specified subset of matB is to be used in the matrix-matrix or matrix-vector multiplication.
Specifying this argument changes the shape of the subset of matB used in the matrix-matrix or matrix-vector multiplication.
See the description of the input argument ndum.
This argument is merely a convenience to resolve the different procedure functionalities within this generic interface.
(optional. If missing, the specified subset of matB will be used as is, without any operations performed on it.)
[in,out]matC: The input/output contiguous vector of shape (:) or matrix of shape (:,:) of the same type and kind as matA, containing the subset of vector/matrix \(\ms{matC}\) in the matrix-matrix or matrix-vector multiplication.
[in]alpha: The input scalar of the same type and kind as matA, containing the coefficient \(\alpha\) in the matrix-matrix or matrix-vector multiplication.
(optional, default = 1. It must be present if any of the input arguments nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC, incB, incC are also present.)
[in]beta: The input scalar of the same type and kind as matA, containing the coefficient \(\beta\) in the matrix-matrix or matrix-vector multiplication.
(optional, default = 1. It must be present if any of the input arguments nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC, incB, incC are also present.)
[in]nrow: The input non-negative scalar of type integer of default kind IK, containing the number of rows of matC used in the matrix-matrix or matrix-vector multiplication matC(roffC + 1 : roffC + nrow, coffC + 1 : coffC + ncol).
(optional, default = size(matC, 1). It must be present if any of the input arguments nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC, incB, incC are also present.)
[in]ncol: The input non-negative scalar of type integer of default kind IK, containing the number of columns of matC used in the matrix-matrix or matrix-vector multiplication matC(roffC + 1 : roffC + nrow, coffC + 1 : coffC + ncol).
(optional, default = size(matC, 2). It must be present if any of the input arguments nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC, incB, incC are also present.)
[in]ndum: The input non-negative scalar of type integer of default kind IK, containing the number of dummy rows or columns of matA and matB used in the matrix-matrix or matrix-vector multiplication.
  1. For matrix matA,
    1. If the input argument operationA is missing, then the subset matA(roffA + 1 : roffA + nrow, coffA + 1 : coffA + ndum) is used in the matrix-matrix or matrix-vector multiplication.
    2. If the input argument operationA is transSymm, then the subset transpose(matA(roffA + 1 : roffA + ndum, coffA + 1 : coffA + nrow) is used in the matrix-matrix or matrix-vector multiplication.
    3. If the input argument operationA is transHerm, then the subset conjg(transpose(matA(roffA + 1 : roffA + ndum, coffA + 1 : coffA + nrow)) is used in the matrix-matrix or matrix-vector multiplication.
  2. For matrix matB,
    1. If the input argument operationB is missing, then the subset matA(roffB + 1 : roffB + nrow, coffB + 1 : coffB + ndum) is used in the matrix-matrix or matrix-vector multiplication.
    2. If the input argument operationB is transSymm, then the subset transpose(matA(roffB + 1 : roffB + ndum, coffB + 1 : coffB + ncol) is used in the matrix-matrix or matrix-vector multiplication.
    3. If the input argument operationB is transHerm, then the subset conjg(transpose(matA(roffB + 1 : roffB + ncol, coffB + 1 : coffB + ndum)) is used in the matrix-matrix or matrix-vector multiplication.
(optional. It must be missing if the input arguments subsetA or subsetB are present, that is, when either matA or matB is a square upper/lower-triangular Symmetric matrix.)
[in]ndim: The input non-negative scalar of type integer of default kind IK, containing the number of rows and columns of the square Symmetric/Hermitian subset of matA used in the matrix-matrix or matrix-vector multiplication matA(roffA + 1 : roffA + ndim, coffA + 1 : coffA + ndim).
(optional, default = size(matA, 1). It can be present only if the arguments matB and matC are of rank 1 and the input arguments incB and incC are present.)
[in]roffA: The input non-negative scalar of type integer of default kind IK containing the offset from the first row of the input matrix matA, such that matA(1 + roffA, 1 + coffA) marks the top-left corner of the block of matA used in the matrix-matrix or matrix-vector multiplication.
(optional, default = 0)
[in]coffA: The input non-negative scalar of type integer of default kind IK containing the offset from the first column of the input matrix matA, such that matA(1 + roffA, 1 + coffA) marks the top-left corner of the block of matA used in the matrix-matrix or matrix-vector multiplication.
(optional, default = 0)
[in]roffB: The input non-negative scalar of type integer of default kind IK containing the offset from the first row of the input matrix matB, such that matB(1 + roffB, 1 + coffB) marks the top-left corner of the block of matB used in the matrix-matrix or matrix-vector multiplication.
(optional, default = 0)
[in]coffB: The input non-negative scalar of type integer of default kind IK containing the offset from the first column of the input matrix matB, such that matB(1 + roffB, 1 + coffB) marks the top-left corner of the block of matB used in the matrix-matrix or matrix-vector multiplication.
(optional, default = 0)
[in]roffC: The input non-negative scalar of type integer of default kind IK containing the offset from the first row of the input matrix matC, such that matC(1 + roffC, 1 + coffC) marks the top-left corner of the block of matC used in the matrix-matrix or matrix-vector multiplication.
(optional, default = 0)
[in]coffC: The input non-negative scalar of type integer of default kind IK containing the offset from the first column of the input matrix matC, such that matC(1 + roffC, 1 + coffC) marks the top-left corner of the block of matC used in the matrix-matrix or matrix-vector multiplication.
(optional, default = 0)
[in]incB: The input non-zero scalar integer of default kind IK, containing the stride of the input argument matB(:).
  1. A positive value implies the multiplication to be performed on the subset matB(1 : 1 + (ndim - 1) * incB : incB).
  2. A negative value implies the multiplication to be performed on the subset matB(1 + (1 - ndim) * incB : 1 : incB).
(optional, default = 1. It can be present only if matB and matC are of rank 1.)
[in]incC: The input non-zero scalar integer of default kind IK, containing the stride of the input/output argument matC(:).
  1. A positive value implies the multiplication to be performed on the subset matC(1 : 1 + (ndim - 1) * incC : incC).
  2. A negative value implies the multiplication to be performed on the subset matC(1 + (1 - ndim) * incC : 1 : incC).
(optional, default = 1. It can be present only if matB and matC are of rank 1.)


Possible calling interfaces

use pm_matrixMulAdd, only: uppDiaA, lowDiaA
use pm_matrixMulAdd, only: uppDiaB, lowDiaB
use pm_matrixMulAdd, only: transSymm, transHerm
use pm_matrixMulAdd, only: symmetric, hermitian
! BLAS - LEVEL 2: ?GEMV - SIMPLIFIED INTERFACE.
call setMatMulAdd(matA(:,:), matB(:), matC(:), alpha = alpha, beta = beta)
call setMatMulAdd(matA(:,:), operationA, matB(:), matC(:), alpha = alpha, beta = beta) ! operationA = transSymm/transHerm
! BLAS - LEVEL 2: ?GEMV - contiguous interface.
call setMatMulAdd(matA(:,:), matB(:), matC(:), alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
call setMatMulAdd(matA(:,:), operationA, matB(:), matC(:), alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! operationA = transSymm/transHerm
! BLAS - LEVEL 3: ?GEMM - SIMPLIFIED INTERFACE.
call setMatMulAdd(matA(:,:), matB(:,:), matC(:,:), alpha = alpha, beta = beta)
call setMatMulAdd(matA(:,:), operationA, matB(:,:), matC(:,:), alpha = alpha, beta = beta) ! operationA = transSymm/transHerm
call setMatMulAdd(matA(:,:), matB(:,:), operationB, matC(:,:), alpha = alpha, beta = beta) ! operationB = transSymm/transHerm
call setMatMulAdd(matA(:,:), operationA, matB(:,:), operationB, matC(:,:), alpha = alpha, beta = beta) ! operationA/operationB = transSymm/transHerm
! BLAS - LEVEL 3: ?GEMM - contiguous interface.
call setMatMulAdd(matA(:,:), matB(:,:), matC(:,:), alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
call setMatMulAdd(matA(:,:), operationA, matB(:,:), matC(:,:), alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC) ! operationA = transSymm/transHerm
call setMatMulAdd(matA(:,:), matB(:,:), operationB, matC(:,:), alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC) ! operationB = transSymm/transHerm
call setMatMulAdd(matA(:,:), operationA, matB(:,:), operationB, matC(:,:), alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC) ! operationA/operationB = transSymm/transHerm
! BLAS - LEVEL 2: ?SYMV / ?HEMV - SIMPLIFIED INTERFACE.
call setMatMulAdd(matA(:,:), classA, subsetA, matB(:), matC(:), alpha = alpha, beta = beta) ! classA = symmetric/hermitian, subsetA = uppDia/lowDia
! BLAS - LEVEL 2: ?SYMV / ?HEMV - contiguous interface.
call setMatMulAdd(matA(:,:), classA, subsetA, matB(:), matC(:), alpha, beta, ndim, roffA, coffA, incB, incC) ! classA = symmetric/hermitian, subsetA = uppDia/lowDia
! BLAS - LEVEL 3: ?SYMM / ?HEMM - SIMPLIFIED INTERFACE.
call setMatMulAdd(matA(:,:), classA, subsetA, matB(:,:), matC(:,:), alpha = alpha, beta = beta) ! classA = symmetric/hermitian, subsetA = uppDia/lowDia
call setMatMulAdd(matA(:,:), matB(:,:), classB, subsetB, matC(:,:), alpha = alpha, beta = beta) ! classB = symmetric/hermitian, subsetB = uppDia/lowDia
! BLAS - LEVEL 3: ?SYMM / ?HEMM - contiguous interface.
call setMatMulAdd(matA(:,:), classA, subsetA, matB(:,:), matC(:,:), alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC) ! classA = symmetric/hermitian, subsetA = uppDia/lowDia
call setMatMulAdd(matA(:,:), matB(:,:), classB, subsetB, matC(:,:), alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC) ! classB = symmetric/hermitian, subsetB = uppDia/lowDia
! BLAS - LEVEL 2: ?SPMV / ?HPMV - SIMPLIFIED INTERFACE.
call setMatMulAdd(matA(:), classA, subsetA, packA, matB(:), matC(:), alpha = alpha, beta = beta) ! classA = symmetric/hermitian, subsetA = uppDia/lowDia, packA = lfpack
! BLAS - LEVEL 2: ?SPMV / ?HPMV - contiguous interface.
call setMatMulAdd(matA(:), classA, subsetA, packA, matB(:), matC(:), alpha, beta, ndim, incB, incC) ! classA = symmetric/hermitian, subsetA = uppDia/lowDia, packA = lfpack
!
Return the result of the multiplication of the input matrices/vector matA and matB in the user-specif...
This module contains procedures and generic interfaces relevant to combined matrix-matrix or matrix-v...
Warning
The condition 0 <= ndum must hold for the corresponding input arguments.
The condition all(0 /= [incB, incC]) must hold for the corresponding input arguments.
The condition all(0 <= [nrow, ncol]) must hold for the corresponding input arguments.
The condition all(0 <= [roffA, coffA]) must hold for the corresponding input arguments.
The condition all(0 <= [roffB, coffB]) must hold for the corresponding input arguments.
The condition all(0 <= [roffC, coffC]) must hold for the corresponding input arguments.
The condition all(0 <= [ndim, roffA, coffA]) must hold for the corresponding input arguments.
The condition abs(incB) * (ndim - 1) < size(matB(:)) must hold for the corresponding input arguments.
The condition abs(incC) * (ndim - 1) < size(matC(:)) must hold for the corresponding input arguments.
The condition all(ndim + [roffA, coffA] <= shape(matA)) must hold for the corresponding input arguments.
The condition all(ndim == [shape(matA(:,:)), shape(matB(:)), shape(matC(:))]) must hold for the corresponding input arguments.
The condition all([roffA + nrow, coffA + ndum] <= shape(matA)) must hold for the corresponding input arguments when the input argument operationA is missing.
The condition all([roffA + ndum, coffA + nrow] <= shape(matA)) must hold for the corresponding input arguments when the input argument operationA is present.
The condition all([roffB + ndum, coffB + ncol] <= shape(matB)) must hold for the corresponding input arguments when the input argument operationB is missing.
The condition all([roffB + ncol, coffB + ndum] <= shape(matB)) must hold for the corresponding input arguments when the input argument operationB is present.
All inequalities in the above conditions become equalities when the optional shape arguments ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC are missing.
There are many more runtime bound and sanity checks that are performed within the procedures, but not listed above.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
BLAS/LAPACK equivalent:
The procedures under discussion combine, modernize, and extend the interface and functionalities of Version 3.11 of BLAS/LAPACK routine(s): SGEMV, DGEMV, CGEMV, ZGEMV, SSPMV, DSPMV, CSPMV, ZSPMV, SHPMV, DHPMV, CHPMV, ZHPMV, SSYMV, DSYMV, CSYMV, ZSYMV, SHEMV, DHEMV, CHEMV, ZHEMV, SSYMM, DSYMM, CSYMM, ZSYMM, SHEMM, DHEMM, CHEMM, ZHEMM, SGEMM, DGEMM, CGEMM, ZGEMM.
In particular multiplications of matrices of type integer of arbitrary kinds are also possible.
BLAS/LAPACK interface arguments setMatMulAdd interface arguments
side = "L" classA argument is present and classB argument is missing.
side = "B" classA argument is missing and classB argument is present.
uplo = "U" subsetA = uppDiaA or subsetB = uppDiaB.
uplo = "L" subsetA = lowDiaA or subsetB = lowDiaB.
transa = "N" operationA argument is missing.
transa = "T" operationA = transSymm
transb = "C" operationA = transHerm
transb = "N" operationB argument is missing.
transb = "T" operationB = transSymm
transb = "C" operationB = transHerm
m nrow
n ncol
k ndum
alpha alpha
a = A(i,j) matA = A, roffA = i - 1, coffA = j - 1
lda NONE (passed implicitly).
b = B(i,j) matB = B, roffB = i - 1, coffB = j - 1
ldb NONE (passed implicitly).
beta beta
c = C(i,j) matC = C, roffC = i - 1, coffC = j - 1
ldc NONE (passed implicitly).
See also
lowDia
uppDia
symmetric
hermitian
transSymm
transHerm
setMatMulTri


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, TKG => RKS
4 use pm_matrixMulAdd, only: symmetric, hermitian
5 use pm_matrixMulAdd, only: transSymm, transHerm
6 use pm_matrixMulAdd, only: uppDia, lowDia
9 use pm_io, only: display_type
10 use pm_io, only: getFormat
11
12 implicit none
13
14 type(display_type) :: disp
15 character(:, SK), allocatable :: cform, rform, iform
16 integer(IK) :: nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC, incB, incC
17 cform = getFormat([cmplx(0., 0., TKG)], ed = SK_'f', signed = .true.)
18 rform = getFormat([real(0., TKG)], ed = SK_'f', signed = .true.)
19 iform = getFormat([0_IK], ed = SK_'i', signed = .true.)
20
21 disp = display_type(file = "main.out.F90")
22
23 call disp%skip()
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%show("! BLAS 2 - GEMV: General matrix-vector multiplications - integer")
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%skip()
28
29 block
30
31 use pm_kind, only: TKG => IKS
32 integer(TKG) :: alpha, beta
33 integer(TKG), parameter :: DUM = huge(DUM)
34 integer(TKG), allocatable :: matA(:,:), matB(:), matC(:), refC(:)
35
36 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
37
38 matA = reshape([ integer(TKG) :: DUM, DUM, DUM, DUM, DUM &
39 , DUM, DUM, 1.0, 2.0, 3.0 &
40 , DUM, DUM, 2.0, 2.0, 4.0 &
41 , DUM, DUM, 3.0, 2.0, 2.0 &
42 , DUM, DUM, 4.0, 2.0, 1.0 &
43 , DUM, DUM, DUM, DUM, DUM &
44 , DUM, DUM, DUM, DUM, DUM &
45 , DUM, DUM, DUM, DUM, DUM &
46 , DUM, DUM, DUM, DUM, DUM &
47 , DUM, DUM, DUM, DUM, DUM &
48 , DUM, DUM, DUM, DUM, DUM ], shape = [11, 5], order = [2, 1])
49 matB = [ integer(TKG) :: DUM, DUM, 3.0, 2.0, 1.0 ]
50 matC = [ integer(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
51 refC = [ integer(TKG) :: 14.0, DUM, 19.0, DUM, 17.0, DUM, 20.0 ]
52 call disp%skip()
53 call disp%show("matA")
54 call disp%show( matA , format = iform )
55 call disp%show("matB")
56 call disp%show( matB , format = iform )
57 call disp%show("matC")
58 call disp%show( matC , format = iform )
59 call disp%show("alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;")
60 alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
61 call disp%skip()
62 call disp%show("call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
63 call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
64 call disp%show("matC")
65 call disp%show( matC , format = iform )
66 call disp%show("matC - refC")
67 call disp%show( matC - refC , format = iform )
68 call disp%skip()
69 call disp%show("matC = [ integer(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]")
70 matC = [ integer(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
71 call disp%show("call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.")
72 call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC))
73 call disp%show("matC")
74 call disp%show( matC , format = iform )
75 call disp%show("matC - refC")
76 call disp%show( matC - refC , format = iform )
77 call disp%skip()
78
79 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
80
81 matA = reshape([ integer(TKG) :: DUM, DUM, DUM, DUM, DUM &
82 , DUM, DUM, 1.0, 2.0, 3.0 &
83 , DUM, DUM, 2.0, 2.0, 4.0 &
84 , DUM, DUM, 3.0, 2.0, 2.0 &
85 , DUM, DUM, 4.0, 2.0, 1.0 &
86 , DUM, DUM, DUM, DUM, DUM &
87 , DUM, DUM, DUM, DUM, DUM &
88 , DUM, DUM, DUM, DUM, DUM &
89 , DUM, DUM, DUM, DUM, DUM &
90 , DUM, DUM, DUM, DUM, DUM &
91 , DUM, DUM, DUM, DUM, DUM ], shape = [11, 5], order = [2, 1])
92 matB = [ integer(TKG) :: DUM, DUM, 3.0, 2.0, 1.0, 4.0 ]
93 matC = [ integer(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]
94 refC = [ integer(TKG) :: 28.0, DUM, 24.0, DUM, 29.0 ]
95 call disp%skip()
96 call disp%show("matA")
97 call disp%show( matA , format = iform )
98 call disp%show("matB")
99 call disp%show( matB , format = iform )
100 call disp%show("matC")
101 call disp%show( matC , format = iform )
102 call disp%show("alpha = 1._TKG; beta = 2._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;")
103 alpha = 1._TKG; beta = 2._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
104 call disp%skip()
105 call disp%show("call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
106 call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
107 call disp%show("matC")
108 call disp%show( matC , format = iform )
109 call disp%show("matC - refC")
110 call disp%show( matC - refC , format = iform )
111 call disp%skip()
112 call disp%show("matC = [ integer(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]")
113 matC = [ integer(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]
114 call disp%show("call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta) ! simplified interface.")
115 call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta)
116 call disp%show("matC")
117 call disp%show( matC , format = iform )
118 call disp%show("matC - refC")
119 call disp%show( matC - refC , format = iform )
120 call disp%skip()
121
122 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
123
124 matA = reshape([ integer(TKG) :: DUM, DUM, DUM, DUM, DUM &
125 , DUM, DUM, 1.0, 2.0, 3.0 &
126 , DUM, DUM, 2.0, 2.0, 4.0 &
127 , DUM, DUM, 3.0, 2.0, 2.0 &
128 , DUM, DUM, 4.0, 2.0, 1.0 &
129 , DUM, DUM, DUM, DUM, DUM &
130 , DUM, DUM, DUM, DUM, DUM &
131 , DUM, DUM, DUM, DUM, DUM &
132 , DUM, DUM, DUM, DUM, DUM &
133 , DUM, DUM, DUM, DUM, DUM &
134 , DUM, DUM, DUM, DUM, DUM ], shape = [11, 5], order = [2, 1])
135 matB = [ integer(TKG) :: DUM, DUM, 3.0, 2.0, 1.0 ]
136 matC = [ integer(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
137 refC = [ integer(TKG) :: 14.0, DUM, 19.0, DUM, 17.0, DUM, 20.0 ]
138 call disp%skip()
139 call disp%show("matA")
140 call disp%show( matA , format = iform )
141 call disp%show("matB")
142 call disp%show( matB , format = iform )
143 call disp%show("matC")
144 call disp%show( matC , format = iform )
145 call disp%show("alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;")
146 alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
147 call disp%skip()
148 call disp%show("call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
149 call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
150 call disp%show("matC")
151 call disp%show( matC , format = iform )
152 call disp%show("matC - refC")
153 call disp%show( matC - refC , format = iform )
154 call disp%skip()
155 call disp%show("matC = [ integer(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]")
156 matC = [ integer(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
157 call disp%show("call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.")
158 call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC))
159 call disp%show("matC")
160 call disp%show( matC , format = iform )
161 call disp%show("matC - refC")
162 call disp%show( matC - refC , format = iform )
163 call disp%skip()
164
165 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
166
167 end block
168
169 call disp%skip()
170 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
171 call disp%show("! BLAS 2 - GEMV: General matrix-vector multiplications - complex")
172 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
173 call disp%skip()
174
175 block
176
177 use pm_kind, only: TKG => RKS
178 complex(TKG) :: alpha, beta
179 complex(TKG), parameter :: COMPLEXDUM = cmplx(huge(0._TKG), huge(0._TKG), TKG)
180 complex(TKG), allocatable :: matA(:,:), matB(:), matC(:), refC(:)
181
182 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
183
184 matA = reshape([ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
185 , COMPLEXDUM, COMPLEXDUM, (1.0, 2.0), (3.0, 5.0), (2.0, 0.0) &
186 , COMPLEXDUM, COMPLEXDUM, (2.0, 3.0), (7.0, 9.0), (4.0, 8.0) &
187 , COMPLEXDUM, COMPLEXDUM, (7.0, 4.0), (1.0, 4.0), (6.0, 0.0) &
188 , COMPLEXDUM, COMPLEXDUM, (8.0, 2.0), (2.0, 5.0), (8.0, 0.0) &
189 , COMPLEXDUM, COMPLEXDUM, (9.0, 1.0), (3.0, 6.0), (1.0, 0.0) &
190 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
191 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
192 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
193 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
194 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ], shape = [11, 5], order = [2, 1])
195 matB = [ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, (1.0, 2.0), (4.0, 0.0), (1.0, 1.0) ]
196 matC = [ complex(TKG) :: (1.0, 2.0), (4.0, 0.0), (1.0, -1.0), (3.0, 4.0), (2.0, 0.0) ]
197 refC = [ complex(TKG) :: (12.0, 28.0), (24.0, 55.0), (10.0, 39.0), (23.0, 50.0), (22.0, 44.0) ]
198 call disp%skip()
199 call disp%show("matA")
200 call disp%show( matA , format = cform )
201 call disp%show("matB")
202 call disp%show( matB , format = cform )
203 call disp%show("matC")
204 call disp%show( matC , format = cform )
205 call disp%show("alpha = 1._TKG; beta = 1._TKG; nrow = 5; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 1;")
206 alpha = 1._TKG; beta = 1._TKG; nrow = 5; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 1;
207 call disp%skip()
208 call disp%show("call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
209 call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
210 call disp%show("matC")
211 call disp%show( matC , format = cform )
212 call disp%show("matC - refC")
213 call disp%show( matC - refC , format = cform )
214 call disp%skip()
215 call disp%show("matC = [ complex(TKG) :: (1.0, 2.0), (4.0, 0.0), (1.0, -1.0), (3.0, 4.0), (2.0, 0.0) ]")
216 matC = [ complex(TKG) :: (1.0, 2.0), (4.0, 0.0), (1.0, -1.0), (3.0, 4.0), (2.0, 0.0) ]
217 call disp%show("call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.")
218 call setMatMulAdd(matA(2:6, 3:5), matB(3:), matC(1::incC))
219 call disp%show("matC")
220 call disp%show( matC , format = cform )
221 call disp%show("matC - refC")
222 call disp%show( matC - refC , format = cform )
223 call disp%skip()
224
225 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
226
227 matA = reshape([ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
228 , COMPLEXDUM, COMPLEXDUM, (1.0, 2.0), (3.0, 5.0), (2.0, 0.0) &
229 , COMPLEXDUM, COMPLEXDUM, (2.0, 3.0), (7.0, 9.0), (4.0, 8.0) &
230 , COMPLEXDUM, COMPLEXDUM, (7.0, 4.0), (1.0, 4.0), (6.0, 0.0) &
231 , COMPLEXDUM, COMPLEXDUM, (8.0, 2.0), (2.0, 5.0), (8.0, 0.0) &
232 , COMPLEXDUM, COMPLEXDUM, (9.0, 1.0), (3.0, 6.0), (1.0, 0.0) &
233 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
234 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
235 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
236 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
237 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ], shape = [11, 5], order = [2, 1])
238 matB = [ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, (1.0, 2.0), (4.0, 0.0), (1.0, 1.0), (3.0, 4.0), (2.0, 0.0) ]
239 matC = [ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ]
240 refC = [ complex(TKG) :: (42.0, 67.0), (10.0, 87.0), (50.0, 74.0) ]
241 call disp%skip()
242 call disp%show("matA")
243 call disp%show( matA , format = cform )
244 call disp%show("matB")
245 call disp%show( matB , format = cform )
246 call disp%show("matC")
247 call disp%show( matC , format = cform )
248 call disp%show("alpha = 1._TKG; beta = 0._TKG; nrow = 5; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 1;")
249 alpha = 1._TKG; beta = 0._TKG; nrow = 5; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 1;
250 call disp%skip()
251 call disp%show("call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
252 call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
253 call disp%show("matC")
254 call disp%show( matC , format = cform )
255 call disp%show("matC - refC")
256 call disp%show( matC - refC , format = cform )
257 call disp%skip()
258 call disp%show("matC = [ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ]")
259 matC = [ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ]
260 call disp%show("call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta) ! simplified interface.")
261 call setMatMulAdd(matA(2:6, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta)
262 call disp%show("matC")
263 call disp%show( matC , format = cform )
264 call disp%show("matC - refC")
265 call disp%show( matC - refC , format = cform )
266 call disp%skip()
267
268 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
269
270 matA = reshape([ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
271 , COMPLEXDUM, COMPLEXDUM, (1.0, 2.0), (3.0, 5.0), (2.0, 0.0) &
272 , COMPLEXDUM, COMPLEXDUM, (2.0, 3.0), (7.0, 9.0), (4.0, 8.0) &
273 , COMPLEXDUM, COMPLEXDUM, (7.0, 4.0), (1.0, 4.0), (6.0, 0.0) &
274 , COMPLEXDUM, COMPLEXDUM, (8.0, 2.0), (2.0, 5.0), (8.0, 0.0) &
275 , COMPLEXDUM, COMPLEXDUM, (9.0, 1.0), (3.0, 6.0), (1.0, 0.0) &
276 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
277 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
278 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
279 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
280 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ], shape = [11, 5], order = [2, 1])
281 matB = [ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, (1.0, 2.0), (4.0, 0.0), (1.0, 1.0), (3.0, 4.0), (2.0, 0.0) ]
282 matC = [ complex(TKG) :: (1.0, 2.0), (4.0, 0.0), (1.0, -1.0) ]
283 refC = [ complex(TKG) :: (-73.0, -13.0), (-74.0, 57.0), (-49.0, -11.0) ]
284 call disp%skip()
285 call disp%show("matA")
286 call disp%show( matA , format = cform )
287 call disp%show("matB")
288 call disp%show( matB , format = cform )
289 call disp%show("matC")
290 call disp%show( matC , format = cform )
291 call disp%show("alpha = -1._TKG; beta = 1._TKG; nrow = 5; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 1;")
292 alpha = -1._TKG; beta = 1._TKG; nrow = 5; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 1;
293 call disp%skip()
294 call disp%show("call setMatMulAdd(matA, transHerm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
295 call setMatMulAdd(matA, transHerm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
296 call disp%show("matC")
297 call disp%show( matC , format = cform )
298 call disp%show("matC - refC")
299 call disp%show( matC - refC , format = cform )
300 call disp%skip()
301 call disp%show("matC = [ complex(TKG) :: (1.0, 2.0), (4.0, 0.0), (1.0, -1.0) ]")
302 matC = [ complex(TKG) :: (1.0, 2.0), (4.0, 0.0), (1.0, -1.0) ]
303 call disp%show("call setMatMulAdd(matA(2:5, 3:5), transHerm, matB(3:), matC(1::incC), alpha) ! simplified interface.")
304 call setMatMulAdd(matA(2:6, 3:5), transHerm, matB(3:), matC(1::incC), alpha)
305 call disp%show("matC")
306 call disp%show( matC , format = cform )
307 call disp%show("matC - refC")
308 call disp%show( matC - refC , format = cform )
309 call disp%skip()
310
311 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
312
313 end block
314
315 call disp%skip()
316 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
317 call disp%show("! BLAS 2 - GEMV: General matrix-vector multiplications - real")
318 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
319 call disp%skip()
320
321 block
322
323 use pm_kind, only: TKG => RKS
324 real(TKG) :: alpha, beta
325 real(TKG), parameter :: DUM = huge(DUM)
326 real(TKG), allocatable :: matA(:,:), matB(:), matC(:), refC(:)
327
328 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
329
330 matA = reshape([ real(TKG) :: DUM, DUM, DUM, DUM, DUM &
331 , DUM, DUM, 1.0, 2.0, 3.0 &
332 , DUM, DUM, 2.0, 2.0, 4.0 &
333 , DUM, DUM, 3.0, 2.0, 2.0 &
334 , DUM, DUM, 4.0, 2.0, 1.0 &
335 , DUM, DUM, DUM, DUM, DUM &
336 , DUM, DUM, DUM, DUM, DUM &
337 , DUM, DUM, DUM, DUM, DUM &
338 , DUM, DUM, DUM, DUM, DUM &
339 , DUM, DUM, DUM, DUM, DUM &
340 , DUM, DUM, DUM, DUM, DUM ], shape = [11, 5], order = [2, 1])
341 matB = [ real(TKG) :: DUM, DUM, 3.0, 2.0, 1.0 ]
342 matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
343 refC = [ real(TKG) :: 14.0, DUM, 19.0, DUM, 17.0, DUM, 20.0 ]
344 call disp%skip()
345 call disp%show("matA")
346 call disp%show( matA , format = rform )
347 call disp%show("matB")
348 call disp%show( matB , format = rform )
349 call disp%show("matC")
350 call disp%show( matC , format = rform )
351 call disp%show("alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;")
352 alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
353 call disp%skip()
354 call disp%show("call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
355 call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
356 call disp%show("matC")
357 call disp%show( matC , format = rform )
358 call disp%show("matC - refC")
359 call disp%show( matC - refC , format = rform )
360 call disp%skip()
361 call disp%show("matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]")
362 matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
363 call disp%show("call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.")
364 call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC))
365 call disp%show("matC")
366 call disp%show( matC , format = rform )
367 call disp%show("matC - refC")
368 call disp%show( matC - refC , format = rform )
369 call disp%skip()
370
371 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
372
373 matA = reshape([ real(TKG) :: DUM, DUM, DUM, DUM, DUM &
374 , DUM, DUM, 1.0, 2.0, 3.0 &
375 , DUM, DUM, 2.0, 2.0, 4.0 &
376 , DUM, DUM, 3.0, 2.0, 2.0 &
377 , DUM, DUM, 4.0, 2.0, 1.0 &
378 , DUM, DUM, DUM, DUM, DUM &
379 , DUM, DUM, DUM, DUM, DUM &
380 , DUM, DUM, DUM, DUM, DUM &
381 , DUM, DUM, DUM, DUM, DUM &
382 , DUM, DUM, DUM, DUM, DUM &
383 , DUM, DUM, DUM, DUM, DUM ], shape = [11, 5], order = [2, 1])
384 matB = [ real(TKG) :: DUM, DUM, 3.0, 2.0, 1.0, 4.0 ]
385 matC = [ real(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]
386 refC = [ real(TKG) :: 28.0, DUM, 24.0, DUM, 29.0 ]
387 call disp%skip()
388 call disp%show("matA")
389 call disp%show( matA , format = rform )
390 call disp%show("matB")
391 call disp%show( matB , format = rform )
392 call disp%show("matC")
393 call disp%show( matC , format = rform )
394 call disp%show("alpha = 1._TKG; beta = 2._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;")
395 alpha = 1._TKG; beta = 2._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
396 call disp%skip()
397 call disp%show("call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
398 call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
399 call disp%show("matC")
400 call disp%show( matC , format = rform )
401 call disp%show("matC - refC")
402 call disp%show( matC - refC , format = rform )
403 call disp%skip()
404 call disp%show("matC = [ real(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]")
405 matC = [ real(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]
406 call disp%show("call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta) ! simplified interface.")
407 call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta)
408 call disp%show("matC")
409 call disp%show( matC , format = rform )
410 call disp%show("matC - refC")
411 call disp%show( matC - refC , format = rform )
412 call disp%skip()
413
414 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
415
416 matA = reshape([ real(TKG) :: DUM, DUM, DUM, DUM, DUM &
417 , DUM, DUM, 1.0, 2.0, 3.0 &
418 , DUM, DUM, 2.0, 2.0, 4.0 &
419 , DUM, DUM, 3.0, 2.0, 2.0 &
420 , DUM, DUM, 4.0, 2.0, 1.0 &
421 , DUM, DUM, DUM, DUM, DUM &
422 , DUM, DUM, DUM, DUM, DUM &
423 , DUM, DUM, DUM, DUM, DUM &
424 , DUM, DUM, DUM, DUM, DUM &
425 , DUM, DUM, DUM, DUM, DUM &
426 , DUM, DUM, DUM, DUM, DUM ], shape = [11, 5], order = [2, 1])
427 matB = [ real(TKG) :: DUM, DUM, 3.0, 2.0, 1.0 ]
428 matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
429 refC = [ real(TKG) :: 14.0, DUM, 19.0, DUM, 17.0, DUM, 20.0 ]
430 call disp%skip()
431 call disp%show("matA")
432 call disp%show( matA , format = rform )
433 call disp%show("matB")
434 call disp%show( matB , format = rform )
435 call disp%show("matC")
436 call disp%show( matC , format = rform )
437 call disp%show("alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;")
438 alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
439 call disp%skip()
440 call disp%show("call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
441 call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
442 call disp%show("matC")
443 call disp%show( matC , format = rform )
444 call disp%show("matC - refC")
445 call disp%show( matC - refC , format = rform )
446 call disp%skip()
447 call disp%show("matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]")
448 matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
449 call disp%show("call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.")
450 call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC))
451 call disp%show("matC")
452 call disp%show( matC , format = rform )
453 call disp%show("matC - refC")
454 call disp%show( matC - refC , format = rform )
455 call disp%skip()
456
457 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
458
459 end block
460
461 block
462
463 use pm_kind, only: TKG => RKH
464 real(TKG) :: alpha, beta
465 real(TKG), parameter :: DUM = huge(DUM)
466 real(TKG), allocatable :: matA(:,:), matB(:), matC(:), refC(:)
467
468 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
469
470 matA = reshape([ real(TKG) :: DUM, DUM, DUM, DUM, DUM &
471 , DUM, DUM, 1.0, 2.0, 3.0 &
472 , DUM, DUM, 2.0, 2.0, 4.0 &
473 , DUM, DUM, 3.0, 2.0, 2.0 &
474 , DUM, DUM, 4.0, 2.0, 1.0 &
475 , DUM, DUM, DUM, DUM, DUM &
476 , DUM, DUM, DUM, DUM, DUM &
477 , DUM, DUM, DUM, DUM, DUM &
478 , DUM, DUM, DUM, DUM, DUM &
479 , DUM, DUM, DUM, DUM, DUM &
480 , DUM, DUM, DUM, DUM, DUM ], shape = [11, 5], order = [2, 1])
481 matB = [ real(TKG) :: DUM, DUM, 3.0, 2.0, 1.0 ]
482 matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
483 refC = [ real(TKG) :: 14.0, DUM, 19.0, DUM, 17.0, DUM, 20.0 ]
484 call disp%skip()
485 call disp%show("matA")
486 call disp%show( matA , format = rform )
487 call disp%show("matB")
488 call disp%show( matB , format = rform )
489 call disp%show("matC")
490 call disp%show( matC , format = rform )
491 call disp%show("alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;")
492 alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
493 call disp%skip()
494 call disp%show("call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
495 call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
496 call disp%show("matC")
497 call disp%show( matC , format = rform )
498 call disp%show("matC - refC")
499 call disp%show( matC - refC , format = rform )
500 call disp%skip()
501 call disp%show("matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]")
502 matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
503 call disp%show("call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.")
504 call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC))
505 call disp%show("matC")
506 call disp%show( matC , format = rform )
507 call disp%show("matC - refC")
508 call disp%show( matC - refC , format = rform )
509 call disp%skip()
510
511 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
512
513 matA = reshape([ real(TKG) :: DUM, DUM, DUM, DUM, DUM &
514 , DUM, DUM, 1.0, 2.0, 3.0 &
515 , DUM, DUM, 2.0, 2.0, 4.0 &
516 , DUM, DUM, 3.0, 2.0, 2.0 &
517 , DUM, DUM, 4.0, 2.0, 1.0 &
518 , DUM, DUM, DUM, DUM, DUM &
519 , DUM, DUM, DUM, DUM, DUM &
520 , DUM, DUM, DUM, DUM, DUM &
521 , DUM, DUM, DUM, DUM, DUM &
522 , DUM, DUM, DUM, DUM, DUM &
523 , DUM, DUM, DUM, DUM, DUM ], shape = [11, 5], order = [2, 1])
524 matB = [ real(TKG) :: DUM, DUM, 3.0, 2.0, 1.0, 4.0 ]
525 matC = [ real(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]
526 refC = [ real(TKG) :: 28.0, DUM, 24.0, DUM, 29.0 ]
527 call disp%skip()
528 call disp%show("matA")
529 call disp%show( matA , format = rform )
530 call disp%show("matB")
531 call disp%show( matB , format = rform )
532 call disp%show("matC")
533 call disp%show( matC , format = rform )
534 call disp%show("alpha = 1._TKG; beta = 2._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;")
535 alpha = 1._TKG; beta = 2._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
536 call disp%skip()
537 call disp%show("call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
538 call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
539 call disp%show("matC")
540 call disp%show( matC , format = rform )
541 call disp%show("matC - refC")
542 call disp%show( matC - refC , format = rform )
543 call disp%skip()
544 call disp%show("matC = [ real(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]")
545 matC = [ real(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]
546 call disp%show("call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta) ! simplified interface.")
547 call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta)
548 call disp%show("matC")
549 call disp%show( matC , format = rform )
550 call disp%show("matC - refC")
551 call disp%show( matC - refC , format = rform )
552 call disp%skip()
553
554 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
555
556 matA = reshape([ real(TKG) :: DUM, DUM, DUM, DUM, DUM &
557 , DUM, DUM, 1.0, 2.0, 3.0 &
558 , DUM, DUM, 2.0, 2.0, 4.0 &
559 , DUM, DUM, 3.0, 2.0, 2.0 &
560 , DUM, DUM, 4.0, 2.0, 1.0 &
561 , DUM, DUM, DUM, DUM, DUM &
562 , DUM, DUM, DUM, DUM, DUM &
563 , DUM, DUM, DUM, DUM, DUM &
564 , DUM, DUM, DUM, DUM, DUM &
565 , DUM, DUM, DUM, DUM, DUM &
566 , DUM, DUM, DUM, DUM, DUM ], shape = [11, 5], order = [2, 1])
567 matB = [ real(TKG) :: DUM, DUM, 3.0, 2.0, 1.0 ]
568 matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
569 refC = [ real(TKG) :: 14.0, DUM, 19.0, DUM, 17.0, DUM, 20.0 ]
570 call disp%skip()
571 call disp%show("matA")
572 call disp%show( matA , format = rform )
573 call disp%show("matB")
574 call disp%show( matB , format = rform )
575 call disp%show("matC")
576 call disp%show( matC , format = rform )
577 call disp%show("alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;")
578 alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
579 call disp%skip()
580 call disp%show("call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.")
581 call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC)
582 call disp%show("matC")
583 call disp%show( matC , format = rform )
584 call disp%show("matC - refC")
585 call disp%show( matC - refC , format = rform )
586 call disp%skip()
587 call disp%show("matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]")
588 matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
589 call disp%show("call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.")
590 call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC))
591 call disp%show("matC")
592 call disp%show( matC , format = rform )
593 call disp%show("matC - refC")
594 call disp%show( matC - refC , format = rform )
595 call disp%skip()
596
597 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
598
599 end block
600
601 call disp%skip()
602 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
603 call disp%show("! Complete general integer matrix-matrix multiplications.")
604 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
605 call disp%skip()
606
607 block
608
609 use pm_kind, only: TKG => IK
610 use pm_distUnif, only: setUnifRand
611 integer(TKG), allocatable, dimension(:,:) :: matA, matB, matC, refC
612 integer(TKG) :: alpha, beta
613 integer(IK) :: nrow, ncol
614
615 alpha = 2_TKG; beta = 3_TKG
616 nrow = 2; ncol = 2; ndum = 3
617 allocate(matA(nrow, ndum), matB(ndum, ncol), matC(nrow, ncol))
618
619 call disp%skip()
620 call disp%show("call setUnifRand(matA, lb = -10_TKG, ub = +10_TKG)")
621 call setUnifRand(matA, lb = -10_TKG, ub = +10_TKG)
622 call disp%show("call setUnifRand(matB, lb = -10_TKG, ub = +10_TKG)")
623 call setUnifRand(matB, lb = -10_TKG, ub = +10_TKG)
624 call disp%show("call setUnifRand(matC, lb = -10_TKG, ub = +10_TKG)")
625 call setUnifRand(matC, lb = -10_TKG, ub = +10_TKG)
626 call disp%show("matA")
627 call disp%show( matA , format = iform )
628 call disp%show("matB")
629 call disp%show( matB , format = iform )
630 call disp%show("matC")
631 call disp%show( matC , format = iform )
632 call disp%show("[alpha, beta]")
633 call disp%show( [alpha, beta] )
634 call disp%show("[nrow, ncol, ndum]")
635 call disp%show( [nrow, ncol, ndum] )
636 call disp%show("refC = matmul(alpha * matA, matB) + beta * matC ! reference value.")
637 refC = matmul(alpha * matA, matB) + beta * matC
638 call disp%show("call setMatMulAdd(matA, matB, matC, alpha, beta)")
639 call setMatMulAdd(matA, matB, matC, alpha, beta)
640 call disp%show("matC - refC")
641 call disp%show( matC - refC , format = iform )
642 call disp%skip()
643 call disp%show("call setUnifRand(matC, lb = -10_TKG, ub = +10_TKG) ! reset for new multiplication.")
644 call setUnifRand(matC, lb = -10_TKG, ub = +10_TKG)
645 call disp%show("matA = transpose(matA)")
646 matA = transpose(matA)
647 call disp%show("matA")
648 call disp%show( matA , format = iform )
649 call disp%show("matC")
650 call disp%show( matC , format = iform )
651 call disp%show("refC = matmul(alpha * transpose(matA), matB) + beta * matC ! reference value.")
652 refC = matmul(alpha * transpose(matA), matB) + beta * matC
653 call disp%show("call setMatMulAdd(matA, transSymm, matB, matC, alpha, beta)")
654 call setMatMulAdd(matA, transSymm, matB, matC, alpha, beta)
655 call disp%show("matC - refC")
656 call disp%show( matC - refC , format = iform )
657 call disp%skip()
658
659 end block
660
661 call disp%skip()
662 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
663 call disp%show("! Complete general complex matrix-matrix multiplications.")
664 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
665 call disp%skip()
666
667 block
668
669 use pm_kind, only: TKG => CKS
670 use pm_distUnif, only: setUnifRand
671 complex(TKG), allocatable, dimension(:,:) :: matA, matB, matC, refC
672 complex(TKG) :: alpha, beta
673 integer(IK) :: nrow, ncol
674
675 nrow = 2; ncol = 2; ndum = 3
676 alpha = (1._TKG, 0._TKG); beta = (0._TKG, 0._TKG)
677 allocate(matA(nrow, ndum), matB(ndum, ncol), matC(nrow, ncol))
678
679 call disp%skip()
680 call disp%show("call setUnifRand(matA)")
681 call setUnifRand(matA)
682 call disp%show("call setUnifRand(matB)")
683 call setUnifRand(matB)
684 call disp%show("call setUnifRand(matC)")
685 call setUnifRand(matC)
686 call disp%show("matA")
687 call disp%show( matA , format = cform )
688 call disp%show("matB")
689 call disp%show( matB , format = cform )
690 call disp%show("matC")
691 call disp%show( matC , format = cform )
692 call disp%show("[alpha, beta]")
693 call disp%show( [alpha, beta] )
694 call disp%show("[nrow, ncol, ndum]")
695 call disp%show( [nrow, ncol, ndum] )
696 call disp%show("refC = matmul(alpha * matA, matB) + beta * matC ! reference value.")
697 refC = matmul(alpha * matA, matB) + beta * matC
698 call disp%show("call setMatMulAdd(matA, matB, matC, alpha, beta)")
699 call setMatMulAdd(matA, matB, matC, alpha, beta)
700 call disp%show("matC - refC")
701 call disp%show( matC - refC , format = cform )
702 call disp%skip()
703 call disp%show("call setUnifRand(matC) ! reset for new multiplication.")
704 call setUnifRand(matC)
705 call disp%show("matA = conjg(transpose(matA))")
706 matA = conjg(transpose(matA))
707 call disp%show("matA")
708 call disp%show( matA , format = cform )
709 call disp%show("matC")
710 call disp%show( matC , format = cform )
711 call disp%show("refC = matmul(alpha * conjg(transpose(matA)), matB) + beta * matC ! reference value.")
712 refC = matmul(alpha * conjg(transpose(matA)), matB) + beta * matC
713 call disp%show("call setMatMulAdd(matA, transHerm, matB, matC, alpha, beta)")
714 call setMatMulAdd(matA, transHerm, matB, matC, alpha, beta)
715 call disp%show("matC - refC")
716 call disp%show( matC - refC , format = cform )
717 call disp%skip()
718
719 end block
720
721 call disp%skip()
722 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
723 call disp%show("! Subset general integer matrix-matrix multiplication.")
724 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
725 call disp%skip()
726
727 block
728
729 use pm_kind, only: IKG => IK
730 integer(TKG) :: alpha, beta
731 integer(TKG), parameter :: DUM = huge(DUM)
732 integer(TKG), allocatable, dimension(:,:) :: matA, matB, matC, refC
733
734 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
735
736 matA = reshape( [ integer(TKG) :: 1.0, 2.0, -1.0, -1.0, 4.0 &
737 , 2.0, 0.0, 1.0, 1.0, -1.0 &
738 , 1.0, -1.0, -1.0, 1.0, 2.0 &
739 , -3.0, 2.0, 2.0, 2.0, 0.0 &
740 , 4.0, 0.0, -2.0, 1.0, -1.0 &
741 , -1.0, -1.0, 1.0, -3.0, 2.0 &
742 , DUM, DUM, DUM, DUM, DUM &
743 , DUM, DUM, DUM, DUM, DUM ], shape = [8, 5], order = [2, 1])
744 matB = reshape( [ integer(TKG) :: 1.0, -1.0, 0.0, 2.0 &
745 , 2.0, 2.0, -1.0, -2.0 &
746 , 1.0, 0.0, -1.0, 1.0 &
747 , -3.0, -1.0, 1.0, -1.0 &
748 , 4.0, 2.0, -1.0, 1.0 &
749 , DUM, DUM, DUM, DUM ], shape = [6, 4], order = [2, 1])
750 matC = reshape( [ integer(TKG) :: 1.0, 1.0, 1.0, 1.0 &
751 , 1.0, 1.0, 1.0, 1.0 &
752 , 1.0, 1.0, 1.0, 1.0 &
753 , 1.0, 1.0, 1.0, 1.0 &
754 , 1.0, 1.0, 1.0, 1.0 &
755 , 1.0, 1.0, 1.0, 1.0 &
756 , DUM, DUM, DUM, DUM ], shape = [6, 4], order = [2, 1])
757 refC = reshape( [ integer(TKG) :: 24.0, 13.0, -5.0, 3.0 &
758 , -3.0, -4.0, 2.0, 4.0 &
759 , 4.0, 1.0, 2.0, 5.0 &
760 , -2.0, 6.0, -1.0, -9.0 &
761 , -4.0, -6.0, 5.0, 5.0 &
762 , 16.0, 7.0, -4.0, 7.0 &
763 , DUM, DUM, DUM, DUM ], shape = [6, 4], order = [2, 1])
764 call disp%skip()
765 call disp%show("matA")
766 call disp%show( matA , format = iform )
767 call disp%show("matB")
768 call disp%show( matB , format = iform )
769 call disp%show("matC")
770 call disp%show( matC , format = iform )
771 call disp%show("alpha = 1._TKG; beta = 1._TKG; nrow = 6; ncol = 4; ndum = 5; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
772 alpha = 1._TKG; beta = 1._TKG; nrow = 6; ncol = 4; ndum = 5; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
773 call disp%show("call setMatMulAdd(matA, matB, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)")
774 call setMatMulAdd(matA, matB, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
775 call disp%show("matC")
776 call disp%show( matC , format = iform )
777 call disp%show("matC - refC")
778 call disp%show( matC - refC , format = iform )
779 call disp%skip()
780
781 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
782
783 matA = reshape( [ integer(TKG) :: 1.0, -3.0 &
784 , 2.0, 4.0 &
785 , 1.0, -1.0 &
786 , DUM, DUM ], shape = [4, 2], order = [2, 1])
787 matB = reshape( [ integer(TKG) :: 1.0, -3.0 &
788 , 2.0, 4.0 &
789 , 1.0, -1.0 ], shape = [3, 2], order = [2, 1])
790 matC = reshape( [ integer(TKG) :: 1.0, 1.0, 1.0 &
791 , 1.0, 1.0, 1.0 &
792 , 1.0, 1.0, 1.0 &
793 , DUM, DUM, DUM &
794 , DUM, DUM, DUM ], shape = [5, 3], order = [2, 1])
795 refC = reshape( [ integer(TKG) :: 11.0, -9., 5.0 &
796 , -9.0, 21., -1.0 &
797 , 5.0, -1., 3.0 &
798 , DUM, DUM, DUM &
799 , DUM, DUM, DUM ], shape = [5, 3], order = [2, 1])
800 call disp%skip()
801 call disp%show("matA")
802 call disp%show( matA , format = iform )
803 call disp%show("matB")
804 call disp%show( matB , format = iform )
805 call disp%show("matC")
806 call disp%show( matC , format = iform )
807 call disp%show("alpha = 1._TKG; beta = 1._TKG; nrow = 3; ncol = 3; ndum = 2; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
808 alpha = 1._TKG; beta = 1._TKG; nrow = 3; ncol = 3; ndum = 2; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
809 call disp%show("call setMatMulAdd(matA, matB, transSymm, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)")
810 call setMatMulAdd(matA, matB, transSymm, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
811 call disp%show("matC")
812 call disp%show( matC , format = iform )
813 call disp%show("matC - refC")
814 call disp%show( matC - refC , format = iform )
815 call disp%skip()
816
817 end block
818
819 call disp%skip()
820 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
821 call disp%show("! Subset general complex matrix-matrix multiplication.")
822 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
823 call disp%skip()
824
825 block
826
827 use pm_kind, only: TKG => CKS
828 complex(TKG) :: alpha, beta
829 complex(TKG), parameter :: COMPLEXDUM = cmplx(huge(0._TKG), huge(0._TKG), TKG)
830 complex(TKG), allocatable, dimension(:,:) :: matA, matB, matC, refC
831
832 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
833
834 matA = reshape( [ complex(TKG) :: (1.0, 5.0), (9.0, 2.0), (1.0, 9.0) &
835 , (2.0, 4.0), (8.0, 3.0), (1.0, 8.0) &
836 , (3.0, 3.0), (7.0, 5.0), (1.0, 7.0) &
837 , (4.0, 2.0), (4.0, 7.0), (1.0, 5.0) &
838 , (5.0, 1.0), (5.0, 1.0), (1.0, 6.0) &
839 , (6.0, 6.0), (3.0, 6.0), (1.0, 4.0) &
840 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
841 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ], shape = [8, 3], order = [2, 1])
842 matB = reshape( [ complex(TKG) :: (1.0, 8.0), (2.0, 7.0) &
843 , (4.0, 4.0), (6.0, 8.0) &
844 , (6.0, 2.0), (4.0, 5.0) &
845 , COMPLEXDUM, COMPLEXDUM ], shape = [4, 2], order = [2, 1])
846 matC = reshape( [ complex(TKG) :: (0.5, 0.0), (0.5, 0.0) &
847 , (0.5, 0.0), (0.5, 0.0) &
848 , (0.5, 0.0), (0.5, 0.0) &
849 , (0.5, 0.0), (0.5, 0.0) &
850 , (0.5, 0.0), (0.5, 0.0) &
851 , (0.5, 0.0), (0.5, 0.0) &
852 , COMPLEXDUM, COMPLEXDUM &
853 , COMPLEXDUM, COMPLEXDUM ], shape = [8, 2], order = [2, 1])
854 refC = reshape( [ complex(TKG) :: (-22.0, 113.0), (-35.0, 142.0) &
855 , (-19.0, 114.0), (-35.0, 141.0) &
856 , (-20.0, 119.0), (-43.0, 146.0) &
857 , (-27.0, 110.0), (-58.0, 131.0) &
858 , (8.0, 103.0), (0.0, 112.0) &
859 , (-55.0, 116.0), (-75.0, 135.0) &
860 , COMPLEXDUM, COMPLEXDUM &
861 , COMPLEXDUM, COMPLEXDUM ], shape = [8, 2], order = [2, 1])
862 call disp%skip()
863 call disp%show("matA")
864 call disp%show( matA , format = cform )
865 call disp%show("matB")
866 call disp%show( matB , format = cform )
867 call disp%show("matC")
868 call disp%show( matC , format = cform )
869 call disp%show("alpha = (1._TKG, 0._TKG); beta = (2._TKG, 0._TKG); nrow = 6; ncol = 2; ndum = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
870 alpha = (1._TKG, 0._TKG); beta = (2._TKG, 0._TKG); nrow = 6; ncol = 2; ndum = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
871 call disp%show("call setMatMulAdd(matA, matB, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)")
872 call setMatMulAdd(matA, matB, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
873 call disp%show("matC")
874 call disp%show( matC , format = cform )
875 call disp%show("matC - refC")
876 call disp%show( matC - refC , format = cform )
877 call disp%skip()
878
879 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
880
881 matA = reshape( [ complex(TKG) :: (1.0, 3.0), (-3.0, 2.0) &
882 , (2.0, 5.0), (4.0, 6.0) &
883 , (1.0, 1.0), (-1.0, 9.0) ], shape = [3, 2], order = [2, 1])
884 matB = reshape( [ complex(TKG) :: (1.0, 2.0), (-3.0, 2.0) &
885 , (2.0, 6.0), (4.0, 5.0) &
886 , (1.0, 2.0), (-1.0, 8.0) &
887 , COMPLEXDUM, COMPLEXDUM ], shape = [4, 2], order = [2, 1])
888 matC = reshape( [ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
889 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
890 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
891 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ], shape = [4, 3], order = [2, 1])
892 refC = reshape( [ complex(TKG) :: (20.0, 1.0), (18.0, 23.0), (26.0, 23.0) &
893 , (12.0, -25.0), (80.0, 2.0), (56.0, -37.0) &
894 , (24.0, -26.0), (49.0, 37.0), (76.0, -2.0) &
895 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ], shape = [4, 3], order = [2, 1])
896 call disp%skip()
897 call disp%show("matA")
898 call disp%show( matA , format = cform )
899 call disp%show("matB")
900 call disp%show( matB , format = cform )
901 call disp%show("matC ! Note that the initialization is irrelevant because `beta = (0., 0.)`.")
902 call disp%show( matC , format = cform )
903 call disp%show("alpha = (1._TKG, 0._TKG); beta = (0._TKG, 0._TKG); nrow = 3; ncol = 3; ndum = 2; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
904 alpha = (1._TKG, 0._TKG); beta = (0._TKG, 0._TKG); nrow = 3; ncol = 3; ndum = 2; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
905 call disp%show("call setMatMulAdd(matA, matB, transHerm, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)")
906 call setMatMulAdd(matA, matB, transHerm, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
907 call disp%show("matC")
908 call disp%show( matC , format = cform )
909 call disp%show("matC - refC")
910 call disp%show( matC - refC , format = cform )
911 call disp%skip()
912
913 end block
914
915 call disp%skip()
916 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
917 call disp%show("! Subset general real matrix-matrix multiplication.")
918 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
919 call disp%skip()
920
921 block
922
923 real(TKG) :: alpha, beta
924 real(TKG), parameter :: DUM = huge(DUM)
925 real(TKG), allocatable, dimension(:,:) :: matA, matB, matC, refC
926
927 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
928
929 matA = reshape([ real(TKG) :: 1.0, 2.0, -1.0, -1.0, 4.0 &
930 , 2.0, 0.0, 1.0, 1.0, -1.0 &
931 , 1.0, -1.0, -1.0, 1.0, 2.0 &
932 , -3.0, 2.0, 2.0, 2.0, 0.0 &
933 , 4.0, 0.0, -2.0, 1.0, -1.0 &
934 , -1.0, -1.0, 1.0, -3.0, 2.0 &
935 , DUM, DUM, DUM, DUM, DUM &
936 , DUM, DUM, DUM, DUM, DUM ], shape = [8, 5], order = [2, 1])
937 matB = reshape([ real(TKG) :: 1.0, -1.0, 0.0, 2.0 &
938 , 2.0, 2.0, -1.0, -2.0 &
939 , 1.0, 0.0, -1.0, 1.0 &
940 , -3.0, -1.0, 1.0, -1.0 &
941 , 4.0, 2.0, -1.0, 1.0 &
942 , DUM, DUM, DUM, DUM ], shape = [6, 4], order = [2, 1])
943 matC = reshape([ real(TKG) :: 0.5, 0.5, 0.5, 0.5 &
944 , 0.5, 0.5, 0.5, 0.5 &
945 , 0.5, 0.5, 0.5, 0.5 &
946 , 0.5, 0.5, 0.5, 0.5 &
947 , 0.5, 0.5, 0.5, 0.5 &
948 , 0.5, 0.5, 0.5, 0.5 &
949 , DUM, DUM, DUM, DUM ], shape = [6, 4], order = [2, 1])
950 refC = reshape([ real(TKG) :: 24.0, 13.0, -5.0, 3.0 &
951 , -3.0, -4.0, 2.0, 4.0 &
952 , 4.0, 1.0, 2.0, 5.0 &
953 , -2.0, 6.0, -1.0, -9.0 &
954 , -4.0, -6.0, 5.0, 5.0 &
955 , 16.0, 7.0, -4.0, 7.0 &
956 , DUM, DUM, DUM, DUM ], shape = [6, 4], order = [2, 1])
957 call disp%skip()
958 call disp%show("matA")
959 call disp%show( matA , format = rform )
960 call disp%show("matB")
961 call disp%show( matB , format = rform )
962 call disp%show("matC")
963 call disp%show( matC , format = rform )
964 call disp%show("alpha = 1._TKG; beta = 2._TKG; nrow = 6; ncol = 4; ndum = 5; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
965 alpha = 1._TKG; beta = 2._TKG; nrow = 6; ncol = 4; ndum = 5; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
966 call disp%show("call setMatMulAdd(matA, matB, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)")
967 call setMatMulAdd(matA, matB, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
968 call disp%show("matC")
969 call disp%show( matC , format = rform )
970 call disp%show("matC - refC")
971 call disp%show( matC - refC , format = rform )
972 call disp%skip()
973
974 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
975
976 matA = reshape([ real(TKG) :: 1.0, -3.0 &
977 , 2.0, 4.0 &
978 , 1.0, -1.0 &
979 , DUM, DUM ], shape = [4, 2], order = [2, 1])
980 matB = reshape([ real(TKG) :: 1.0, -3.0 &
981 , 2.0, 4.0 &
982 , 1.0, -1.0 ], shape = [3, 2], order = [2, 1])
983 matC = reshape([ real(TKG) :: 0.5, 0.5, 0.5 &
984 , 0.5, 0.5, 0.5 &
985 , 0.5, 0.5, 0.5 &
986 , DUM, DUM, DUM &
987 , DUM, DUM, DUM ], shape = [5, 3], order = [2, 1])
988 refC = reshape([ real(TKG) :: 11.0, -9., 5.0 &
989 , -9.0, 21., -1.0 &
990 , 5.0, -1., 3.0 &
991 , DUM, DUM, DUM &
992 , DUM, DUM, DUM ], shape = [5, 3], order = [2, 1])
993 call disp%skip()
994 call disp%show("matA")
995 call disp%show( matA , format = rform )
996 call disp%show("matB")
997 call disp%show( matB , format = rform )
998 call disp%show("matC")
999 call disp%show( matC , format = rform )
1000 call disp%show("alpha = 1._TKG; beta = 2._TKG; nrow = 3; ncol = 3; ndum = 2; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
1001 alpha = 1._TKG; beta = 2._TKG; nrow = 3; ncol = 3; ndum = 2; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
1002 call disp%show("call setMatMulAdd(matA, matB, transSymm, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)")
1003 call setMatMulAdd(matA, matB, transSymm, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
1004 call disp%show("matC")
1005 call disp%show( matC , format = rform )
1006 call disp%show("matC - refC")
1007 call disp%show( matC - refC , format = rform )
1008 call disp%skip()
1009
1010 end block
1011
1012 call disp%skip()
1013 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1014 call disp%show("! Subset symmetric complex matrix-matrix multiplication.")
1015 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1016 call disp%skip()
1017
1018 block
1019
1020 use pm_kind, only: TKG => CKS
1021 complex(TKG) :: alpha, beta
1022 real(TKG), parameter :: DUM = huge(DUM)
1023 complex(TKG), parameter :: COMPLEXDUM = cmplx(huge(0._TKG), huge(0._TKG), TKG)
1024 complex(TKG), allocatable, dimension(:,:) :: matA, matB, matC, refC
1025
1026 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1027
1028 matB = reshape( [ complex(TKG) :: (1.0, 5.0), (-3.0, 2.0), (1.0, 6.0) &
1029 , COMPLEXDUM, (4.0, 5.0), (-1.0, 4.0) &
1030 , COMPLEXDUM, COMPLEXDUM, (2.0, 5.0) &
1031 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1032 ], shape = [4, 3], order = [2, 1])
1033 matA = reshape( [ complex(TKG) :: (1.0, 1.0), (-3.0, 2.0), (3.0, 3.0) &
1034 , (2.0, 6.0), (4.0, 5.0), (-1.0, 4.0) &
1035 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1036 ], shape = [3, 3], order = [2, 1])
1037 matC = reshape( [ complex(TKG) :: (13.0, 6.0), (-18.0, 6.0), (10.0, 7.0) &
1038 , (-11.0, 8.0), (11.0, 1.0), (-4.0, 2.0) &
1039 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1040 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1041 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1042 ], shape = [5, 3], order = [2, 1])
1043 refC = reshape( [ complex(TKG) :: (-96.0, 72.0), (-141.0, -226.0), (-112.0, 38.0) &
1044 , (-230.0, -269.0), (-133.0, -23.0), (-272.0, -198.0) &
1045 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1046 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1047 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1048 ], shape = [5, 3], order = [2, 1])
1049 call disp%skip()
1050 call disp%show("matA")
1051 call disp%show( matA , format = cform )
1052 call disp%show("matB")
1053 call disp%show( matB , format = cform )
1054 call disp%show("matC")
1055 call disp%show( matC , format = cform )
1056 call disp%show("alpha = (2._TKG, 3._TKG); beta = (1._TKG, 6._TKG); nrow = 2; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
1057 alpha = (2._TKG, 3._TKG); beta = (1._TKG, 6._TKG); nrow = 2; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
1058 call disp%show("call setMatMulAdd(matA, matB, symmetric, uppDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)")
1059 call setMatMulAdd(matA, matB, symmetric, uppDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
1060 call disp%show("matC")
1061 call disp%show( matC , format = cform )
1062 call disp%show("matC - refC")
1063 call disp%show( matC - refC , format = cform )
1064 call disp%skip()
1065
1066 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1067
1068 matB = reshape( [ complex(TKG) :: (1.0, DUM), COMPLEXDUM, COMPLEXDUM &
1069 , (3.0, 2.0), (4.0, DUM), COMPLEXDUM &
1070 , (-1.0, 6.0), (1.0, 4.0), (2.0, DUM) &
1071 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1072 ], shape = [4, 3], order = [2, 1])
1073 matA = reshape( [ complex(TKG) :: (1.0, 1.0), (-3.0, 2.0), (3.0, 3.0) &
1074 , (2.0, 6.0), (4.0, 5.0), (-1.0, 4.0) &
1075 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1076 ], shape = [3, 3], order = [2, 1])
1077 matC = reshape( [ complex(TKG) :: (13.0, 6.0), (-18.0, 6.0), (10.0, 7.0) &
1078 , (-11.0, 8.0), (11.0, 1.0), (-4.0, 2.0) &
1079 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1080 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1081 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1082 ], shape = [5, 3], order = [2, 1])
1083 refC = reshape( [ complex(TKG) :: (-137.0, 17.0), (-158.0, -102.0), (-39.0, 141.0) &
1084 , (-154.0, -77.0), (-63.0, 186.0), (159.0, 104.0) &
1085 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1086 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1087 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1088 ], shape = [5, 3], order = [2, 1])
1089 call disp%skip()
1090 call disp%show("matA")
1091 call disp%show( matA , format = cform )
1092 call disp%show("matB")
1093 call disp%show( matB , format = cform )
1094 call disp%show("matC")
1095 call disp%show( matC , format = cform )
1096 call disp%show("alpha = (2._TKG, 3._TKG); beta = (1._TKG, 6._TKG); nrow = 2; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
1097 alpha = (2._TKG, 3._TKG); beta = (1._TKG, 6._TKG); nrow = 2; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
1098 call disp%show("call setMatMulAdd(matA, matB, hermitian, lowDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)")
1099 call setMatMulAdd(matA, matB, hermitian, lowDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
1100 call disp%show("matC")
1101 call disp%show( matC , format = cform )
1102 call disp%show("matC - refC")
1103 call disp%show( matC - refC , format = cform )
1104 call disp%skip()
1105
1106 end block
1107
1108 call disp%skip()
1109 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1110 call disp%show("! Subset symmetric real matrix-matrix multiplication.")
1111 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1112 call disp%skip()
1113
1114 block
1115
1116 real(TKG) :: alpha, beta
1117 real(TKG), parameter :: DUM = huge(DUM)
1118 real(TKG), allocatable, dimension(:,:) :: matA, matB, matC, refC
1119
1120 matA = reshape([ real(TKG) :: 1.0, 2.0, -1.0, -1.0, 4.0 &
1121 , DUM, 0.0, 1.0, 1.0, -1.0 &
1122 , DUM, DUM, -1.0, 1.0, 2.0 &
1123 , DUM, DUM, DUM, 2.0, 0.0 &
1124 , DUM, DUM, DUM, DUM, -1.0 &
1125 , DUM, DUM, DUM, DUM, DUM &
1126 , DUM, DUM, DUM, DUM, DUM &
1127 , DUM, DUM, DUM, DUM, DUM &
1128 ], shape = [8, 5], order = [2, 1])
1129 matB = reshape([ real(TKG) :: 1.0, -1.0, 0.0, 2.0 &
1130 , 2.0, 2.0, -1.0, -2.0 &
1131 , 1.0, 0.0, -1.0, 1.0 &
1132 , -3.0, -1.0, 1.0, -1.0 &
1133 , 4.0, 2.0, -1.0, 1.0 &
1134 , DUM, DUM, DUM, DUM &
1135 ], shape = [6, 4], order = [2, 1])
1136 matC = reshape([ real(TKG) :: 23.0, 12.0, -6.0, 2.0 &
1137 , -4.0, -5.0, 1.0, 3.0 &
1138 , 5.0, 6.0, -1.0, -4.0 &
1139 , -4.0, 1.0, 0.0, -5.0 &
1140 , 8.0, -4.0, -2.0, 13.0 &
1141 ], shape = [5, 4], order = [2, 1])
1142 refC = reshape([ real(TKG) :: 69.0, 36.0, -18.0, 6.0 &
1143 , -12.0, -15.0, 3.0, 9.0 &
1144 , 15.0, 18.0, -3.0, -12.0 &
1145 , -12.0, 3.0, 0.0, -15.0 &
1146 , 8.0, -20.0, -2.0, 35.0 &
1147 ], shape = [5, 4], order = [2, 1])
1148 call disp%skip()
1149 call disp%show("matA")
1150 call disp%show( matA , format = rform )
1151 call disp%show("matB")
1152 call disp%show( matB , format = rform )
1153 call disp%show("matC")
1154 call disp%show( matC , format = rform )
1155 call disp%show("alpha = 2._TKG; beta = 1._TKG; nrow = 5; ncol = 4; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
1156 alpha = 2._TKG; beta = 1._TKG; nrow = 5; ncol = 4; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
1157 call disp%show("call setMatMulAdd(matA, symmetric, uppDia, matB, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)")
1158 call setMatMulAdd(matA, symmetric, uppDia, matB, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
1159 call disp%show("matC")
1160 call disp%show( matC , format = rform )
1161 call disp%show("matC - refC")
1162 call disp%show( matC - refC , format = rform )
1163 call disp%skip()
1164
1165 matC = reshape([ real(TKG) :: 23.0, 12.0, -6.0, 2.0 &
1166 , -4.0, -5.0, 1.0, 3.0 &
1167 , 5.0, 6.0, -1.0, -4.0 &
1168 , -4.0, 1.0, 0.0, -5.0 &
1169 , 8.0, -4.0, -2.0, 13.0 &
1170 ], shape = [5, 4], order = [2, 1])
1171 call disp%show("matC")
1172 call disp%show( matC , format = rform )
1173 call disp%show("call setMatMulAdd(matA(1:5, 1:5), symmetric, uppDia, matB(1:5, 1:4), matC, alpha, beta)")
1174 call setMatMulAdd(matA(1:5, 1:5), symmetric, uppDia, matB(1:5, 1:4), matC, alpha, beta)
1175 call disp%show("matC")
1176 call disp%show( matC , format = rform )
1177 call disp%show("matC - refC")
1178 call disp%show( matC - refC , format = rform )
1179 call disp%skip()
1180
1181 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1182
1183 matA = reshape([ real(TKG) :: 1.0, DUM, DUM &
1184 , 2.0, 4.0, DUM &
1185 , 1.0, -1.0, -1.0 &
1186 , DUM, DUM, DUM &
1187 ], shape = [4, 3], order = [2, 1])
1188 matB = reshape([ real(TKG) :: 1.0, -3.0, 2.0, 2.0, -1.0, 2.0 &
1189 , 2.0, 4.0, 0.0, 0.0, 1.0, -2.0 &
1190 , 1.0, -1.0, -1.0, -1.0, -1.0, 1.0 &
1191 ], shape = [3, 6], order = [2, 1])
1192 matC = reshape([ real(TKG) :: 6.0, 4.0, 1.0, 1.0, 0.0, -1.0 &
1193 , 9.0, 11.0, 5.0, 5.0, 3.0, -5.0 &
1194 , -2.0, -6.0, 3.0, 3.0, -1.0, 32.0 &
1195 , DUM, DUM, DUM, DUM, DUM, DUM &
1196 , DUM, DUM, DUM, DUM, DUM, DUM &
1197 ], shape = [5, 6], order = [2, 1])
1198 refC = reshape([ real(TKG) :: 24.0, 16.0, 4.0, 4.0, 0.0, -4.0 &
1199 , 36.0, 44.0, 20.0, 20.0, 12.0, -20.0 &
1200 , -8.0, -24.0, 12.0, 12.0, -4.0, 70.0 &
1201 , DUM, DUM, DUM, DUM, DUM, DUM &
1202 , DUM, DUM, DUM, DUM, DUM, DUM &
1203 ], shape = [5, 6], order = [2, 1])
1204 call disp%skip()
1205 call disp%show("matA")
1206 call disp%show( matA , format = rform )
1207 call disp%show("matB")
1208 call disp%show( matB , format = rform )
1209 call disp%show("matC")
1210 call disp%show( matC , format = rform )
1211 call disp%show("alpha = 2._TKG; beta = 2._TKG; nrow = 3; ncol = 6; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
1212 alpha = 2._TKG; beta = 2._TKG; nrow = 3; ncol = 6; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
1213 call disp%show("call setMatMulAdd(matA, symmetric, lowDia, matB, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)")
1214 call setMatMulAdd(matA, symmetric, lowDia, matB, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
1215 call disp%show("matC")
1216 call disp%show( matC , format = rform )
1217 call disp%show("matC - refC")
1218 call disp%show( matC - refC , format = rform )
1219 call disp%skip()
1220
1221 matC = reshape([ real(TKG) :: 6.0, 4.0, 1.0, 1.0, 0.0, -1.0 &
1222 , 9.0, 11.0, 5.0, 5.0, 3.0, -5.0 &
1223 , -2.0, -6.0, 3.0, 3.0, -1.0, 32.0 &
1224 , DUM, DUM, DUM, DUM, DUM, DUM &
1225 , DUM, DUM, DUM, DUM, DUM, DUM &
1226 ], shape = [5, 6], order = [2, 1])
1227 call disp%show("matC")
1228 call disp%show( matC , format = rform )
1229 call disp%show("call setMatMulAdd(matA(1:3, 1:3), symmetric, lowDia, matB(1:3, 1:6), matC(1:3, 1:6), alpha, beta)")
1230 call setMatMulAdd(matA(1:3, 1:3), symmetric, lowDia, matB(1:3, 1:6), matC(1:3, 1:6), alpha, beta)
1231 call disp%show("matC")
1232 call disp%show( matC , format = rform )
1233 call disp%show("matC - refC")
1234 call disp%show( matC - refC , format = rform )
1235 call disp%skip()
1236
1237 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1238
1239 matB = reshape([ real(TKG) :: 1.0, -3.0, 1.0 &
1240 , DUM, 4.0, -1.0 &
1241 , DUM, DUM, 2.0 &
1242 , DUM, DUM, DUM &
1243 ], shape = [4, 3], order = [2, 1])
1244 matA = reshape([ real(TKG) :: 1.0, -3.0, 3.0 &
1245 , 2.0, 4.0, -1.0 &
1246 , DUM, DUM, DUM &
1247 ], shape = [3, 3], order = [2, 1])
1248 matC = reshape([ real(TKG) :: 13.0, -18.0, 10.0 &
1249 , -11.0, 11.0, -4.0 &
1250 , DUM, DUM, DUM &
1251 , DUM, DUM, DUM &
1252 , DUM, DUM, DUM &
1253 ], shape = [5, 3], order = [2, 1])
1254 refC = reshape([ real(TKG) :: 39.0, -54.0, 30.0 &
1255 , -33.0, 33.0, -12.0 &
1256 , DUM, DUM, DUM &
1257 , DUM, DUM, DUM &
1258 , DUM, DUM, DUM &
1259 ], shape = [5, 3], order = [2, 1])
1260 call disp%skip()
1261 call disp%show("matA")
1262 call disp%show( matA , format = rform )
1263 call disp%show("matB")
1264 call disp%show( matB , format = rform )
1265 call disp%show("matC")
1266 call disp%show( matC , format = rform )
1267 call disp%show("alpha = 2._TKG; beta = 1._TKG; nrow = 2; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
1268 alpha = 2._TKG; beta = 1._TKG; nrow = 2; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
1269 call disp%show("call setMatMulAdd(matA, matB, symmetric, uppDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)")
1270 call setMatMulAdd(matA, matB, symmetric, uppDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
1271 call disp%show("matC")
1272 call disp%show( matC , format = rform )
1273 call disp%show("matC - refC")
1274 call disp%show( matC - refC , format = rform )
1275 call disp%skip()
1276
1277 matC = reshape([ real(TKG) :: 13.0, -18.0, 10.0 &
1278 , -11.0, 11.0, -4.0 &
1279 , DUM, DUM, DUM &
1280 , DUM, DUM, DUM &
1281 , DUM, DUM, DUM &
1282 ], shape = [5, 3], order = [2, 1])
1283 call disp%show("matC")
1284 call disp%show( matC , format = rform )
1285 call disp%show("call setMatMulAdd(matA(1:2, 1:3), matB(1:3, 1:3), symmetric, uppDia, matC(1:2, 1:3), alpha, beta)")
1286 call setMatMulAdd(matA(1:2, 1:3), matB(1:3, 1:3), symmetric, uppDia, matC(1:2, 1:3), alpha, beta)
1287 call disp%show("matC")
1288 call disp%show( matC , format = rform )
1289 call disp%show("matC - refC")
1290 call disp%show( matC - refC , format = rform )
1291 call disp%skip()
1292
1293 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1294
1295 matB = reshape([ real(TKG) :: 1.0, DUM, DUM &
1296 , 2.0, 10.0, DUM &
1297 , 1.0, 11.0, 4.0 &
1298 ], shape = [3, 3], order = [2, 1])
1299 matA = reshape([ real(TKG) :: 1.0, -3.0, 2.0 &
1300 , 2.0, 4.0, 0.0 &
1301 , 1.0, -1.0, -1.0 &
1302 ], shape = [3, 3], order = [2, 1])
1303 matC = reshape([ real(TKG) :: 1.0, 5.0, -9.0 &
1304 , -3.0, 10.0, -2.0 &
1305 , -2.0, 8.0, 0.0 &
1306 ], shape = [3, 3], order = [2, 1])
1307 refC = reshape([ real(TKG) :: 4.0, 11.0, 15.0 &
1308 , -13.0, -34.0, -48.0 &
1309 , 0.0, 27.0, 14.0 &
1310 ], shape = [3, 3], order = [2, 1])
1311 call disp%skip()
1312 call disp%show("matA")
1313 call disp%show( matA , format = rform )
1314 call disp%show("matB")
1315 call disp%show( matB , format = rform )
1316 call disp%show("matC")
1317 call disp%show( matC , format = rform )
1318 call disp%show("alpha = -1._TKG; beta = 1._TKG; nrow = 3; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;")
1319 alpha = -1._TKG; beta = 1._TKG; nrow = 3; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
1320 call disp%show("call setMatMulAdd(matA, matB, symmetric, lowDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)")
1321 call setMatMulAdd(matA, matB, symmetric, lowDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
1322 call disp%show("matC")
1323 call disp%show( matC , format = rform )
1324 call disp%show("matC - refC")
1325 call disp%show( matC - refC , format = rform )
1326 call disp%skip()
1327
1328 matC = reshape([ real(TKG) :: 1.0, 5.0, -9.0 &
1329 , -3.0, 10.0, -2.0 &
1330 , -2.0, 8.0, 0.0 &
1331 ], shape = [3, 3], order = [2, 1])
1332 call disp%show("matC")
1333 call disp%show( matC , format = rform )
1334 call disp%show("call setMatMulAdd(matA(1:3, 1:3), matB(1:3, 1:3), symmetric, lowDia, matC(1:3, 1:3), alpha, beta)")
1335 call setMatMulAdd(matA(1:3, 1:3), matB(1:3, 1:3), symmetric, lowDia, matC(1:3, 1:3), alpha, beta)
1336 call disp%show("matC")
1337 call disp%show( matC , format = rform )
1338 call disp%show("matC - refC")
1339 call disp%show( matC - refC , format = rform )
1340 call disp%skip()
1341
1342 end block
1343
1344 call disp%skip()
1345 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1346 call disp%show("! Subset symmetric real matrix-vector multiplication.")
1347 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1348 call disp%skip()
1349
1350 block
1351
1352 real(TKG) :: alpha, beta
1353 real(TKG) , parameter :: DUM = huge(DUM)
1354 real(TKG) , allocatable :: matA(:,:), matB(:), matC(:), refC(:)
1355 integer(IK) :: ndim, incB, incC
1356
1357 matA = reshape([ real(TKG) :: 8.0, DUM, DUM &
1358 , 4.0, 6.0, DUM &
1359 , 2.0, 7.0, 3.0 &
1360 ], shape = [3, 3], order = [2, 1])
1361 matB = [3.00, 2.0, 1.00]
1362 matC = [5.00, DUM, 3.00, DUM, 2.00]
1363 refC = [39.0, DUM, 34.0, DUM, 25.0]
1364 call disp%skip()
1365 call disp%show("matA")
1366 call disp%show( matA , format = rform )
1367 call disp%show("matB")
1368 call disp%show( matB , format = rform )
1369 call disp%show("matC")
1370 call disp%show( matC , format = rform )
1371 call disp%show("alpha = 1._TKG; beta = 1._TKG; ndim = 3; roffA = 0; coffA = 0; incB = 1; incC = 2;")
1372 alpha = 1._TKG; beta = 1._TKG; ndim = 3; roffA = 0; coffA = 0; incB = 1; incC = 2;
1373 call disp%show("call setMatMulAdd(matA, symmetric, lowDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)")
1374 call setMatMulAdd(matA, symmetric, lowDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)
1375 call disp%show("matC")
1376 call disp%show( matC , format = rform )
1377 call disp%show("matC - refC")
1378 call disp%show( matC - refC , format = rform )
1379 call disp%skip()
1380
1381 matB = [3.00, 2.0, 1.00]
1382 matC = [5.00, 3.00, 2.00]
1383 refC = [39.0, 34.0, 25.0]
1384 call disp%show("matB")
1385 call disp%show( matB , format = rform )
1386 call disp%show("matC")
1387 call disp%show( matC , format = rform )
1388 call disp%show("alpha = 1._TKG; beta = 1._TKG;")
1389 alpha = 1._TKG; beta = 1._TKG;
1390 call disp%show("call setMatMulAdd(matA, symmetric, lowDia, matB, matC, alpha, beta)")
1391 call setMatMulAdd(matA, symmetric, lowDia, matB, matC, alpha, beta)
1392 call disp%show("matC")
1393 call disp%show( matC , format = rform )
1394 call disp%show("matC - refC")
1395 call disp%show( matC - refC , format = rform )
1396 call disp%skip()
1397
1398 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1399
1400 matA = reshape([ real(TKG) :: 8.0, 4.0, 2.0 &
1401 , DUM, 6.0, 7.0 &
1402 , DUM, DUM, 3.0 &
1403 , DUM, DUM, DUM &
1404 ], shape = [3, 3], order = [2, 1])
1405 matB = [4.0, DUM, 2.0, DUM, 1.0]
1406 matC = [6.0, 5.0, 4.0]
1407 refC = [36.0, 54.0, 36.0]
1408 call disp%skip()
1409 call disp%show("matA")
1410 call disp%show( matA , format = rform )
1411 call disp%show("matB")
1412 call disp%show( matB , format = rform )
1413 call disp%show("matC")
1414 call disp%show( matC , format = rform )
1415 call disp%show("alpha = 1._TKG; beta = 2._TKG; ndim = 3; roffA = 0; coffA = 0; incB = -2; incC = 1;")
1416 alpha = 1._TKG; beta = 2._TKG; ndim = 3; roffA = 0; coffA = 0; incB = -2; incC = 1;
1417 call disp%show("call setMatMulAdd(matA, symmetric, uppDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)")
1418 call setMatMulAdd(matA, symmetric, uppDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)
1419 call disp%show("matC")
1420 call disp%show( matC , format = rform )
1421 call disp%show("matC - refC")
1422 call disp%show( matC - refC , format = rform )
1423 call disp%skip()
1424
1425 matB = matB(size(matB):1:incB)
1426 matC = [6.0, 5.0, 4.0]
1427 call disp%show("matB")
1428 call disp%show( matB , format = rform )
1429 call disp%show("matC")
1430 call disp%show( matC , format = rform )
1431 call disp%show("alpha = 1._TKG; beta = 2._TKG;")
1432 alpha = 1._TKG; beta = 2._TKG;
1433 call disp%show("call setMatMulAdd(matA, symmetric, uppDia, matB, matC, alpha, beta)")
1434 call setMatMulAdd(matA, symmetric, uppDia, matB, matC, alpha, beta)
1435 call disp%show("matC")
1436 call disp%show( matC , format = rform )
1437 call disp%show("matC - refC")
1438 call disp%show( matC - refC , format = rform )
1439 call disp%skip()
1440
1441 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1442
1443 end block
1444
1445 call disp%skip()
1446 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1447 call disp%show("! Subset Hermitian complex matrix-vector multiplication.")
1448 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1449 call disp%skip()
1450
1451 block
1452
1453 use pm_kind, only: TKG => CKS
1454 complex(TKG) :: alpha, beta
1455 real(TKG) , parameter :: DUM = huge(DUM)
1456 complex(TKG) , parameter :: COMPLEXDUM = cmplx(huge(0._TKG), huge(0._TKG), TKG)
1457 complex(TKG) , allocatable :: matA(:,:), matB(:), matC(:), refC(:)
1458 integer(IK) :: ndim, incB, incC
1459
1460 matA = reshape([ real(TKG) :: (1.0, DUM), COMPLEXDUM, COMPLEXDUM &
1461 , (3.0, -5.0), (7.0, DUM), COMPLEXDUM &
1462 , (2.0, 3.0), (4.0, 8.0), (6.0, DUM) &
1463 ], shape = [3, 3], order = [2, 1])
1464 matB = [(1.0, 2.0), (4.0, 0.0), (3.0, 4.0)]
1465 matC = [(1.0, 0.0), COMPLEXDUM, (2.0, -1.0), COMPLEXDUM, (2.0, 1.0)]
1466 refC = [(20., 10.), COMPLEXDUM, (45., 21.0), COMPLEXDUM, (38., 29.)]
1467 call disp%skip()
1468 call disp%show("matA")
1469 call disp%show( matA , format = cform )
1470 call disp%show("matB")
1471 call disp%show( matB , format = cform )
1472 call disp%show("matC")
1473 call disp%show( matC , format = cform )
1474 call disp%show("alpha = (1._TKG, 0._TKG); beta = (1._TKG, 0._TKG); ndim = 3; roffA = 0; coffA = 0; incB = 1; incC = 2;")
1475 alpha = (1._TKG, 0._TKG); beta = (1._TKG, 0._TKG); ndim = 3; roffA = 0; coffA = 0; incB = 1; incC = 2;
1476 call disp%show("call setMatMulAdd(matA, hermitian, lowDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)")
1477 call setMatMulAdd(matA, hermitian, lowDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)
1478 call disp%show("matC")
1479 call disp%show( matC , format = cform )
1480 call disp%show("matC - refC")
1481 call disp%show( matC - refC , format = cform )
1482 call disp%skip()
1483
1484 matC = [(1.0, 0.0), (2.0, -1.0), (2.0, 1.0)]
1485 refC = refC(1::incC)
1486 call disp%show("matB")
1487 call disp%show( matB , format = cform )
1488 call disp%show("matC")
1489 call disp%show( matC , format = cform )
1490 call disp%show("call setMatMulAdd(matA, hermitian, lowDia, matB, matC, alpha, beta)")
1491 call setMatMulAdd(matA, hermitian, lowDia, matB, matC, alpha, beta)
1492 call disp%show("matC")
1493 call disp%show( matC , format = cform )
1494 call disp%show("matC - refC")
1495 call disp%show( matC - refC , format = cform )
1496 call disp%skip()
1497
1498 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1499
1500 matA = reshape([ real(TKG) :: COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM &
1501 , COMPLEXDUM, COMPLEXDUM, (1.0, DUM), (3.0, 5.0), (2.0, -3.0) &
1502 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, (7.0, DUM), (4.0, -8.0) &
1503 , COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, (6.0, DUM) &
1504 ], shape = [4, 5], order = [2, 1])
1505 matB = [(3.0, 4.0), COMPLEXDUM, (4.0, 0.0), COMPLEXDUM, (1.0, 2.0)]
1506 matC = [COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM, COMPLEXDUM]
1507 refC = [(19., 10.), COMPLEXDUM, (43., 22.), COMPLEXDUM, (36., 28.)]
1508 call disp%skip()
1509 call disp%show("matA")
1510 call disp%show( matA , format = cform )
1511 call disp%show("matB")
1512 call disp%show( matB , format = cform )
1513 call disp%show("matC")
1514 call disp%show( matC , format = cform )
1515 call disp%show("alpha = (1._TKG, 0._TKG); beta = (0._TKG, 0._TKG); ndim = 3; roffA = 1; coffA = 2; incB = -2; incC = 2;")
1516 alpha = (1._TKG, 0._TKG); beta = (0._TKG, 0._TKG); ndim = 3; roffA = 1; coffA = 2; incB = -2; incC = 2;
1517 call disp%show("call setMatMulAdd(matA, hermitian, uppDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)")
1518 call setMatMulAdd(matA, hermitian, uppDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)
1519 call disp%show("matC")
1520 call disp%show( matC , format = cform )
1521 call disp%show("matC - refC")
1522 call disp%show( matC - refC , format = cform )
1523 call disp%skip()
1524
1525 matA = matA(2:, 3:)
1526 matB = matB(size(matB):1:incB)
1527 matC = [COMPLEXDUM, COMPLEXDUM, COMPLEXDUM]
1528 refC = refC(1::-incB)
1529 call disp%show("matA")
1530 call disp%show( matA , format = cform )
1531 call disp%show("matB")
1532 call disp%show( matB , format = cform )
1533 call disp%show("matC")
1534 call disp%show( matC , format = cform )
1535 call disp%show("call setMatMulAdd(matA, hermitian, uppDia, matB, matC, alpha, beta)")
1536 call setMatMulAdd(matA, hermitian, uppDia, matB, matC, alpha, beta)
1537 call disp%show("matC")
1538 call disp%show( matC , format = cform )
1539 call disp%show("matC - refC")
1540 call disp%show( matC - refC , format = cform )
1541 call disp%skip()
1542
1543 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1544
1545 end block
1546
1547end program example
Generate and return an output array whose elements are the reversed-order elements of the input array...
Return a uniform random scalar or contiguous array of arbitrary rank of randomly uniformly distribute...
Generate and return a generic or type/kind-specific IO format with the requested specifications that ...
Definition: pm_io.F90:18485
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
This module contains procedures and generic interfaces for reversing the order of elements in arrays ...
This module contains classes and procedures for computing various statistical quantities related to t...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter CKS
The single-precision complex kind in Fortran mode. On most platforms, this is a 32-bit real kind.
Definition: pm_kind.F90:570
integer, parameter IKS
The single-precision integer kind in Fortran mode. On most platforms, this is a 32-bit integer kind.
Definition: pm_kind.F90:563
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3! BLAS 2 - GEMV: General matrix-vector multiplications - integer
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7matA
8*, *, *, *, *
9*, *, *, *, *
10*, *, *, *, *
11*, *, *, *, *
12*, *, *, *, *
13*, *, *, *, *
14*, *, *, *, *
15*, *, *, *, *
16*, *, *, *, *
17*, *, *, *, *
18*, *, *, *, *
19matB
20*, *, *, *, *
21matC
22*, *, *, *, *, *, *
23alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
24
25call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
26matC
27*, *, *, *, *, *, *
28matC - refC
29*, *, *, *, *, *, *
30
31matC = [ integer(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
32call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.
33matC
34*, *, *, *, *, *, *
35matC - refC
36*, *, *, *, *, *, *
37
38
39matA
40*, *, *, *, *
41*, *, *, *, *
42*, *, *, *, *
43*, *, *, *, *
44*, *, *, *, *
45*, *, *, *, *
46*, *, *, *, *
47*, *, *, *, *
48*, *, *, *, *
49*, *, *, *, *
50*, *, *, *, *
51matB
52*, *, *, *, *, *
53matC
54*, *, *, *, *
55alpha = 1._TKG; beta = 2._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
56
57call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
58matC
59*, *, *, *, *
60matC - refC
61*, *, *, *, *
62
63matC = [ integer(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]
64call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta) ! simplified interface.
65matC
66*, *, *, *, *
67matC - refC
68*, *, *, *, *
69
70
71matA
72*, *, *, *, *
73*, *, *, *, *
74*, *, *, *, *
75*, *, *, *, *
76*, *, *, *, *
77*, *, *, *, *
78*, *, *, *, *
79*, *, *, *, *
80*, *, *, *, *
81*, *, *, *, *
82*, *, *, *, *
83matB
84*, *, *, *, *
85matC
86*, *, *, *, *, *, *
87alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
88
89call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
90matC
91*, *, *, *, *, *, *
92matC - refC
93*, *, *, *, *, *, *
94
95matC = [ integer(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
96call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.
97matC
98*, *, *, *, *, *, *
99matC - refC
100*, *, *, *, *, *, *
101
102
103!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
104! BLAS 2 - GEMV: General matrix-vector multiplications - complex
105!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
106
107
108matA
109(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
110(***************, ***************), (***************, ***************), ( +1.000000, +2.000000), ( +3.000000, +5.000000), ( +2.000000, +0.000000)
111(***************, ***************), (***************, ***************), ( +2.000000, +3.000000), ( +7.000000, +9.000000), ( +4.000000, +8.000000)
112(***************, ***************), (***************, ***************), ( +7.000000, +4.000000), ( +1.000000, +4.000000), ( +6.000000, +0.000000)
113(***************, ***************), (***************, ***************), ( +8.000000, +2.000000), ( +2.000000, +5.000000), ( +8.000000, +0.000000)
114(***************, ***************), (***************, ***************), ( +9.000000, +1.000000), ( +3.000000, +6.000000), ( +1.000000, +0.000000)
115(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
116(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
117(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
118(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
119(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
120matB
121(***************, ***************), (***************, ***************), ( +1.000000, +2.000000), ( +4.000000, +0.000000), ( +1.000000, +1.000000)
122matC
123( +1.000000, +2.000000), ( +4.000000, +0.000000), ( +1.000000, -1.000000), ( +3.000000, +4.000000), ( +2.000000, +0.000000)
124alpha = 1._TKG; beta = 1._TKG; nrow = 5; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 1;
125
126call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
127matC
128( +12.000000, +28.000000), ( +24.000000, +55.000000), ( +10.000000, +39.000000), ( +23.000000, +50.000000), ( +22.000000, +44.000000)
129matC - refC
130( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
131
132matC = [ complex(TKG) :: (1.0, 2.0), (4.0, 0.0), (1.0, -1.0), (3.0, 4.0), (2.0, 0.0) ]
133call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.
134matC
135( +12.000000, +28.000000), ( +24.000000, +55.000000), ( +10.000000, +39.000000), ( +23.000000, +50.000000), ( +22.000000, +44.000000)
136matC - refC
137( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
138
139
140matA
141(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
142(***************, ***************), (***************, ***************), ( +1.000000, +2.000000), ( +3.000000, +5.000000), ( +2.000000, +0.000000)
143(***************, ***************), (***************, ***************), ( +2.000000, +3.000000), ( +7.000000, +9.000000), ( +4.000000, +8.000000)
144(***************, ***************), (***************, ***************), ( +7.000000, +4.000000), ( +1.000000, +4.000000), ( +6.000000, +0.000000)
145(***************, ***************), (***************, ***************), ( +8.000000, +2.000000), ( +2.000000, +5.000000), ( +8.000000, +0.000000)
146(***************, ***************), (***************, ***************), ( +9.000000, +1.000000), ( +3.000000, +6.000000), ( +1.000000, +0.000000)
147(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
148(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
149(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
150(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
151(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
152matB
153(***************, ***************), (***************, ***************), ( +1.000000, +2.000000), ( +4.000000, +0.000000), ( +1.000000, +1.000000), ( +3.000000, +4.000000), ( +2.000000, +0.000000)
154matC
155(***************, ***************), (***************, ***************), (***************, ***************)
156alpha = 1._TKG; beta = 0._TKG; nrow = 5; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 1;
157
158call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
159matC
160( +42.000000, +67.000000), ( +10.000000, +87.000000), ( +50.000000, +74.000000)
161matC - refC
162( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
163
164matC = [ complex(TKG) :: COMPLEXDUM, COMPLEXDUM, COMPLEXDUM ]
165call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta) ! simplified interface.
166matC
167( +42.000000, +67.000000), ( +10.000000, +87.000000), ( +50.000000, +74.000000)
168matC - refC
169( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
170
171
172matA
173(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
174(***************, ***************), (***************, ***************), ( +1.000000, +2.000000), ( +3.000000, +5.000000), ( +2.000000, +0.000000)
175(***************, ***************), (***************, ***************), ( +2.000000, +3.000000), ( +7.000000, +9.000000), ( +4.000000, +8.000000)
176(***************, ***************), (***************, ***************), ( +7.000000, +4.000000), ( +1.000000, +4.000000), ( +6.000000, +0.000000)
177(***************, ***************), (***************, ***************), ( +8.000000, +2.000000), ( +2.000000, +5.000000), ( +8.000000, +0.000000)
178(***************, ***************), (***************, ***************), ( +9.000000, +1.000000), ( +3.000000, +6.000000), ( +1.000000, +0.000000)
179(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
180(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
181(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
182(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
183(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
184matB
185(***************, ***************), (***************, ***************), ( +1.000000, +2.000000), ( +4.000000, +0.000000), ( +1.000000, +1.000000), ( +3.000000, +4.000000), ( +2.000000, +0.000000)
186matC
187( +1.000000, +2.000000), ( +4.000000, +0.000000), ( +1.000000, -1.000000)
188alpha = -1._TKG; beta = 1._TKG; nrow = 5; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 1;
189
190call setMatMulAdd(matA, transHerm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
191matC
192( -73.000000, -13.000000), ( -74.000000, +57.000000), ( -49.000000, -11.000000)
193matC - refC
194( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
195
196matC = [ complex(TKG) :: (1.0, 2.0), (4.0, 0.0), (1.0, -1.0) ]
197call setMatMulAdd(matA(2:5, 3:5), transHerm, matB(3:), matC(1::incC), alpha) ! simplified interface.
198matC
199( -73.000000, -13.000000), ( -74.000000, +57.000000), ( -49.000000, -11.000000)
200matC - refC
201( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
202
203
204!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
205! BLAS 2 - GEMV: General matrix-vector multiplications - real
206!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
207
208
209matA
210***************, ***************, ***************, ***************, ***************
211***************, ***************, +1.000000, +2.000000, +3.000000
212***************, ***************, +2.000000, +2.000000, +4.000000
213***************, ***************, +3.000000, +2.000000, +2.000000
214***************, ***************, +4.000000, +2.000000, +1.000000
215***************, ***************, ***************, ***************, ***************
216***************, ***************, ***************, ***************, ***************
217***************, ***************, ***************, ***************, ***************
218***************, ***************, ***************, ***************, ***************
219***************, ***************, ***************, ***************, ***************
220***************, ***************, ***************, ***************, ***************
221matB
222***************, ***************, +3.000000, +2.000000, +1.000000
223matC
224 +4.000000, ***************, +5.000000, ***************, +2.000000, ***************, +3.000000
225alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
226
227call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
228matC
229 +14.000000, ***************, +19.000000, ***************, +17.000000, ***************, +20.000000
230matC - refC
231 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
232
233matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
234call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.
235matC
236 +14.000000, ***************, +19.000000, ***************, +17.000000, ***************, +20.000000
237matC - refC
238 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
239
240
241matA
242***************, ***************, ***************, ***************, ***************
243***************, ***************, +1.000000, +2.000000, +3.000000
244***************, ***************, +2.000000, +2.000000, +4.000000
245***************, ***************, +3.000000, +2.000000, +2.000000
246***************, ***************, +4.000000, +2.000000, +1.000000
247***************, ***************, ***************, ***************, ***************
248***************, ***************, ***************, ***************, ***************
249***************, ***************, ***************, ***************, ***************
250***************, ***************, ***************, ***************, ***************
251***************, ***************, ***************, ***************, ***************
252***************, ***************, ***************, ***************, ***************
253matB
254***************, ***************, +3.000000, +2.000000, +1.000000, +4.000000
255matC
256 +1.000000, ***************, +2.000000, ***************, +3.000000
257alpha = 1._TKG; beta = 2._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
258
259call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
260matC
261 +28.000000, ***************, +24.000000, ***************, +29.000000
262matC - refC
263 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
264
265matC = [ real(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]
266call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta) ! simplified interface.
267matC
268 +28.000000, ***************, +24.000000, ***************, +29.000000
269matC - refC
270 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
271
272
273matA
274***************, ***************, ***************, ***************, ***************
275***************, ***************, +1.000000, +2.000000, +3.000000
276***************, ***************, +2.000000, +2.000000, +4.000000
277***************, ***************, +3.000000, +2.000000, +2.000000
278***************, ***************, +4.000000, +2.000000, +1.000000
279***************, ***************, ***************, ***************, ***************
280***************, ***************, ***************, ***************, ***************
281***************, ***************, ***************, ***************, ***************
282***************, ***************, ***************, ***************, ***************
283***************, ***************, ***************, ***************, ***************
284***************, ***************, ***************, ***************, ***************
285matB
286***************, ***************, +3.000000, +2.000000, +1.000000
287matC
288 +4.000000, ***************, +5.000000, ***************, +2.000000, ***************, +3.000000
289alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
290
291call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
292matC
293 +14.000000, ***************, +19.000000, ***************, +17.000000, ***************, +20.000000
294matC - refC
295 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
296
297matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
298call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.
299matC
300 +14.000000, ***************, +19.000000, ***************, +17.000000, ***************, +20.000000
301matC - refC
302 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
303
304
305matA
306***************, ***************, ***************, ***************, ***************
307***************, ***************, +1.000000, +2.000000, +3.000000
308***************, ***************, +2.000000, +2.000000, +4.000000
309***************, ***************, +3.000000, +2.000000, +2.000000
310***************, ***************, +4.000000, +2.000000, +1.000000
311***************, ***************, ***************, ***************, ***************
312***************, ***************, ***************, ***************, ***************
313***************, ***************, ***************, ***************, ***************
314***************, ***************, ***************, ***************, ***************
315***************, ***************, ***************, ***************, ***************
316***************, ***************, ***************, ***************, ***************
317matB
318***************, ***************, +3.000000, +2.000000, +1.000000
319matC
320 +4.000000, ***************, +5.000000, ***************, +2.000000, ***************, +3.000000
321alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
322
323call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
324matC
325 +14.000000, ***************, +19.000000, ***************, +17.000000, ***************, +20.000000
326matC - refC
327 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
328
329matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
330call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.
331matC
332 +14.000000, ***************, +19.000000, ***************, +17.000000, ***************, +20.000000
333matC - refC
334 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
335
336
337matA
338***************, ***************, ***************, ***************, ***************
339***************, ***************, +1.000000, +2.000000, +3.000000
340***************, ***************, +2.000000, +2.000000, +4.000000
341***************, ***************, +3.000000, +2.000000, +2.000000
342***************, ***************, +4.000000, +2.000000, +1.000000
343***************, ***************, ***************, ***************, ***************
344***************, ***************, ***************, ***************, ***************
345***************, ***************, ***************, ***************, ***************
346***************, ***************, ***************, ***************, ***************
347***************, ***************, ***************, ***************, ***************
348***************, ***************, ***************, ***************, ***************
349matB
350***************, ***************, +3.000000, +2.000000, +1.000000, +4.000000
351matC
352 +1.000000, ***************, +2.000000, ***************, +3.000000
353alpha = 1._TKG; beta = 2._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
354
355call setMatMulAdd(matA, transSymm, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
356matC
357 +28.000000, ***************, +24.000000, ***************, +29.000000
358matC - refC
359 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
360
361matC = [ real(TKG) :: 1.0, DUM, 2.0, DUM, 3.0 ]
362call setMatMulAdd(matA(2:5, 3:5), transSymm, matB(3:), matC(1::incC), beta = beta) ! simplified interface.
363matC
364 +28.000000, ***************, +24.000000, ***************, +29.000000
365matC - refC
366 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
367
368
369matA
370***************, ***************, ***************, ***************, ***************
371***************, ***************, +1.000000, +2.000000, +3.000000
372***************, ***************, +2.000000, +2.000000, +4.000000
373***************, ***************, +3.000000, +2.000000, +2.000000
374***************, ***************, +4.000000, +2.000000, +1.000000
375***************, ***************, ***************, ***************, ***************
376***************, ***************, ***************, ***************, ***************
377***************, ***************, ***************, ***************, ***************
378***************, ***************, ***************, ***************, ***************
379***************, ***************, ***************, ***************, ***************
380***************, ***************, ***************, ***************, ***************
381matB
382***************, ***************, +3.000000, +2.000000, +1.000000
383matC
384 +4.000000, ***************, +5.000000, ***************, +2.000000, ***************, +3.000000
385alpha = 1._TKG; beta = 1._TKG; nrow = 4; ncol = 3; roffA = 1; coffA = 2; incB = 1; incC = 2;
386
387call setMatMulAdd(matA, matB(3:), matC, alpha, beta, nrow, ncol, roffA, coffA, incB, incC) ! full contiguous interface.
388matC
389 +14.000000, ***************, +19.000000, ***************, +17.000000, ***************, +20.000000
390matC - refC
391 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
392
393matC = [ real(TKG) :: 4.0, DUM, 5.0, DUM, 2.0, DUM, 3.0 ]
394call setMatMulAdd(matA(2:5, 3:5), matB(3:), matC(1::incC)) ! simplified interface.
395matC
396 +14.000000, ***************, +19.000000, ***************, +17.000000, ***************, +20.000000
397matC - refC
398 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
399
400
401!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
402! Complete general integer matrix-matrix multiplications.
403!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
404
405
406call setUnifRand(matA, lb = -10_TKG, ub = +10_TKG)
407call setUnifRand(matB, lb = -10_TKG, ub = +10_TKG)
408call setUnifRand(matC, lb = -10_TKG, ub = +10_TKG)
409matA
410*, *, *
411*, *, *
412matB
413*, *
414*, *
415*, *
416matC
417*, *
418*, *
419[alpha, beta]
420+2, +3
421[nrow, ncol, ndum]
422+2, +2, +3
423refC = matmul(alpha * matA, matB) + beta * matC ! reference value.
424call setMatMulAdd(matA, matB, matC, alpha, beta)
425matC - refC
426*, *
427*, *
428
429call setUnifRand(matC, lb = -10_TKG, ub = +10_TKG) ! reset for new multiplication.
430matA = transpose(matA)
431matA
432*, *
433*, *
434*, *
435matC
436*, *
437*, *
438refC = matmul(alpha * transpose(matA), matB) + beta * matC ! reference value.
439call setMatMulAdd(matA, transSymm, matB, matC, alpha, beta)
440matC - refC
441*, *
442*, *
443
444
445!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
446! Complete general complex matrix-matrix multiplications.
447!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
448
449
450call setUnifRand(matA)
451call setUnifRand(matB)
452call setUnifRand(matC)
453matA
454( +0.877938, +0.316958), ( +0.118660, +0.973250), ( +0.991527, +0.255652)
455( +0.728736, +0.980935), ( +0.489920, +0.331403), ( +0.032310, +0.608260)
456matB
457( +0.542904, +0.140401), ( +0.289301, +0.383171)
458( +0.379832, +0.537845), ( +0.231047, +0.546758)
459( +0.447297, +0.754162), ( +0.255238, +0.112539)
460matC
461( +0.236097, +0.905751), ( +0.859790, +0.629443)
462( +0.123540, +0.043919), ( +0.410959, +0.742840)
463[alpha, beta]
464(+1.00000000, +0.00000000), (+0.00000000, +0.00000000)
465[nrow, ncol, ndum]
466+2, +2, +3
467refC = matmul(alpha * matA, matB) + beta * matC ! reference value.
468call setMatMulAdd(matA, matB, matC, alpha, beta)
469matC - refC
470( +0.000000, +0.000000), ( +0.000000, +0.000000)
471( +0.000000, +0.000000), ( +0.000000, +0.000000)
472
473call setUnifRand(matC) ! reset for new multiplication.
474matA = conjg(transpose(matA))
475matA
476( +0.877938, -0.316958), ( +0.728736, -0.980935)
477( +0.118660, -0.973250), ( +0.489920, -0.331403)
478( +0.991527, -0.255652), ( +0.032310, -0.608260)
479matC
480( +0.994759, +0.675531), ( +0.123368, +0.815181)
481( +0.632780, +0.332137), ( +0.198606, +0.457706)
482refC = matmul(alpha * conjg(transpose(matA)), matB) + beta * matC ! reference value.
483call setMatMulAdd(matA, transHerm, matB, matC, alpha, beta)
484matC - refC
485( +0.000000, +0.000000), ( +0.000000, +0.000000)
486( +0.000000, +0.000000), ( +0.000000, +0.000000)
487
488
489!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
490! Subset general integer matrix-matrix multiplication.
491!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
492
493
494matA
495*, *, *, *, *
496*, *, *, *, *
497*, *, *, *, *
498*, *, *, *, *
499*, *, *, *, *
500*, *, *, *, *
501*, *, *, *, *
502*, *, *, *, *
503matB
504*, *, *, *
505*, *, *, *
506*, *, *, *
507*, *, *, *
508*, *, *, *
509*, *, *, *
510matC
511*, *, *, *
512*, *, *, *
513*, *, *, *
514*, *, *, *
515*, *, *, *
516*, *, *, *
517alpha = 1._TKG; beta = 1._TKG; nrow = 6; ncol = 4; ndum = 5; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
518call setMatMulAdd(matA, matB, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
519matC
520*, *, *, *
521*, *, *, *
522*, *, *, *
523*, *, *, *
524*, *, *, *
525*, *, *, *
526matC - refC
527*, *, *, *
528*, *, *, *
529*, *, *, *
530*, *, *, *
531*, *, *, *
532*, *, *, *
533
534
535matA
536*, *
537*, *
538*, *
539*, *
540matB
541*, *
542*, *
543*, *
544matC
545*, *, *
546*, *, *
547*, *, *
548*, *, *
549*, *, *
550alpha = 1._TKG; beta = 1._TKG; nrow = 3; ncol = 3; ndum = 2; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
551call setMatMulAdd(matA, matB, transSymm, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
552matC
553*, *, *
554*, *, *
555*, *, *
556*, *, *
557*, *, *
558matC - refC
559*, *, *
560*, *, *
561*, *, *
562*, *, *
563*, *, *
564
565
566!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
567! Subset general complex matrix-matrix multiplication.
568!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
569
570
571matA
572( +1.000000, +5.000000), ( +9.000000, +2.000000), ( +1.000000, +9.000000)
573( +2.000000, +4.000000), ( +8.000000, +3.000000), ( +1.000000, +8.000000)
574( +3.000000, +3.000000), ( +7.000000, +5.000000), ( +1.000000, +7.000000)
575( +4.000000, +2.000000), ( +4.000000, +7.000000), ( +1.000000, +5.000000)
576( +5.000000, +1.000000), ( +5.000000, +1.000000), ( +1.000000, +6.000000)
577( +6.000000, +6.000000), ( +3.000000, +6.000000), ( +1.000000, +4.000000)
578(***************, ***************), (***************, ***************), (***************, ***************)
579(***************, ***************), (***************, ***************), (***************, ***************)
580matB
581( +1.000000, +8.000000), ( +2.000000, +7.000000)
582( +4.000000, +4.000000), ( +6.000000, +8.000000)
583( +6.000000, +2.000000), ( +4.000000, +5.000000)
584(***************, ***************), (***************, ***************)
585matC
586( +0.500000, +0.000000), ( +0.500000, +0.000000)
587( +0.500000, +0.000000), ( +0.500000, +0.000000)
588( +0.500000, +0.000000), ( +0.500000, +0.000000)
589( +0.500000, +0.000000), ( +0.500000, +0.000000)
590( +0.500000, +0.000000), ( +0.500000, +0.000000)
591( +0.500000, +0.000000), ( +0.500000, +0.000000)
592(***************, ***************), (***************, ***************)
593(***************, ***************), (***************, ***************)
594alpha = (1._TKG, 0._TKG); beta = (2._TKG, 0._TKG); nrow = 6; ncol = 2; ndum = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
595call setMatMulAdd(matA, matB, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
596matC
597( -22.000000, +113.000000), ( -35.000000, +142.000000)
598( -19.000000, +114.000000), ( -35.000000, +141.000000)
599( -20.000000, +119.000000), ( -43.000000, +146.000000)
600( -27.000000, +110.000000), ( -58.000000, +131.000000)
601( +8.000000, +103.000000), ( +0.000000, +112.000000)
602( -55.000000, +116.000000), ( -75.000000, +135.000000)
603(***************, ***************), (***************, ***************)
604(***************, ***************), (***************, ***************)
605matC - refC
606( +0.000000, +0.000000), ( +0.000000, +0.000000)
607( +0.000000, +0.000000), ( +0.000000, +0.000000)
608( +0.000000, +0.000000), ( +0.000000, +0.000000)
609( +0.000000, +0.000000), ( +0.000000, +0.000000)
610( +0.000000, +0.000000), ( +0.000000, +0.000000)
611( +0.000000, +0.000000), ( +0.000000, +0.000000)
612( +0.000000, +0.000000), ( +0.000000, +0.000000)
613( +0.000000, +0.000000), ( +0.000000, +0.000000)
614
615
616matA
617( +1.000000, +3.000000), ( -3.000000, +2.000000)
618( +2.000000, +5.000000), ( +4.000000, +6.000000)
619( +1.000000, +1.000000), ( -1.000000, +9.000000)
620matB
621( +1.000000, +2.000000), ( -3.000000, +2.000000)
622( +2.000000, +6.000000), ( +4.000000, +5.000000)
623( +1.000000, +2.000000), ( -1.000000, +8.000000)
624(***************, ***************), (***************, ***************)
625matC ! Note that the initialization is irrelevant because `beta = (0., 0.)`.
626(***************, ***************), (***************, ***************), (***************, ***************)
627(***************, ***************), (***************, ***************), (***************, ***************)
628(***************, ***************), (***************, ***************), (***************, ***************)
629(***************, ***************), (***************, ***************), (***************, ***************)
630alpha = (1._TKG, 0._TKG); beta = (0._TKG, 0._TKG); nrow = 3; ncol = 3; ndum = 2; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
631call setMatMulAdd(matA, matB, transHerm, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
632matC
633( +20.000000, +1.000000), ( +18.000000, +23.000000), ( +26.000000, +23.000000)
634( +12.000000, -25.000000), ( +80.000000, +2.000000), ( +56.000000, -37.000000)
635( +24.000000, -26.000000), ( +49.000000, +37.000000), ( +76.000000, -2.000000)
636(***************, ***************), (***************, ***************), (***************, ***************)
637matC - refC
638( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
639( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
640( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
641( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
642
643
644!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
645! Subset general real matrix-matrix multiplication.
646!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
647
648
649matA
650 +1.000000, +2.000000, -1.000000, -1.000000, +4.000000
651 +2.000000, +0.000000, +1.000000, +1.000000, -1.000000
652 +1.000000, -1.000000, -1.000000, +1.000000, +2.000000
653 -3.000000, +2.000000, +2.000000, +2.000000, +0.000000
654 +4.000000, +0.000000, -2.000000, +1.000000, -1.000000
655 -1.000000, -1.000000, +1.000000, -3.000000, +2.000000
656***************, ***************, ***************, ***************, ***************
657***************, ***************, ***************, ***************, ***************
658matB
659 +1.000000, -1.000000, +0.000000, +2.000000
660 +2.000000, +2.000000, -1.000000, -2.000000
661 +1.000000, +0.000000, -1.000000, +1.000000
662 -3.000000, -1.000000, +1.000000, -1.000000
663 +4.000000, +2.000000, -1.000000, +1.000000
664***************, ***************, ***************, ***************
665matC
666 +0.500000, +0.500000, +0.500000, +0.500000
667 +0.500000, +0.500000, +0.500000, +0.500000
668 +0.500000, +0.500000, +0.500000, +0.500000
669 +0.500000, +0.500000, +0.500000, +0.500000
670 +0.500000, +0.500000, +0.500000, +0.500000
671 +0.500000, +0.500000, +0.500000, +0.500000
672alpha = 1._TKG; beta = 2._TKG; nrow = 6; ncol = 4; ndum = 5; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
673call setMatMulAdd(matA, matB, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
674matC
675 +24.000000, +13.000000, -5.000000, +3.000000
676 -3.000000, -4.000000, +2.000000, +4.000000
677 +4.000000, +1.000000, +2.000000, +5.000000
678 -2.000000, +6.000000, -1.000000, -9.000000
679 -4.000000, -6.000000, +5.000000, +5.000000
680 +16.000000, +7.000000, -4.000000, +7.000000
681matC - refC
682 +0.000000, +0.000000, +0.000000, +0.000000
683 +0.000000, +0.000000, +0.000000, +0.000000
684 +0.000000, +0.000000, +0.000000, +0.000000
685 +0.000000, +0.000000, +0.000000, +0.000000
686 +0.000000, +0.000000, +0.000000, +0.000000
687 +0.000000, +0.000000, +0.000000, +0.000000
688
689
690matA
691 +1.000000, -3.000000
692 +2.000000, +4.000000
693 +1.000000, -1.000000
694***************, ***************
695matB
696 +1.000000, -3.000000
697 +2.000000, +4.000000
698 +1.000000, -1.000000
699matC
700 +0.500000, +0.500000, +0.500000
701 +0.500000, +0.500000, +0.500000
702 +0.500000, +0.500000, +0.500000
703***************, ***************, ***************
704***************, ***************, ***************
705alpha = 1._TKG; beta = 2._TKG; nrow = 3; ncol = 3; ndum = 2; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
706call setMatMulAdd(matA, matB, transSymm, matC, alpha, beta, nrow, ncol, ndum, roffA, coffA, roffB, coffB, roffC, coffC)
707matC
708 +11.000000, -9.000000, +5.000000
709 -9.000000, +21.000000, -1.000000
710 +5.000000, -1.000000, +3.000000
711***************, ***************, ***************
712***************, ***************, ***************
713matC - refC
714 +0.000000, +0.000000, +0.000000
715 +0.000000, +0.000000, +0.000000
716 +0.000000, +0.000000, +0.000000
717 +0.000000, +0.000000, +0.000000
718 +0.000000, +0.000000, +0.000000
719
720
721!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
722! Subset symmetric complex matrix-matrix multiplication.
723!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
724
725
726matA
727( +1.000000, +1.000000), ( -3.000000, +2.000000), ( +3.000000, +3.000000)
728( +2.000000, +6.000000), ( +4.000000, +5.000000), ( -1.000000, +4.000000)
729(***************, ***************), (***************, ***************), (***************, ***************)
730matB
731( +1.000000, +5.000000), ( -3.000000, +2.000000), ( +1.000000, +6.000000)
732(***************, ***************), ( +4.000000, +5.000000), ( -1.000000, +4.000000)
733(***************, ***************), (***************, ***************), ( +2.000000, +5.000000)
734(***************, ***************), (***************, ***************), (***************, ***************)
735matC
736( +13.000000, +6.000000), ( -18.000000, +6.000000), ( +10.000000, +7.000000)
737( -11.000000, +8.000000), ( +11.000000, +1.000000), ( -4.000000, +2.000000)
738(***************, ***************), (***************, ***************), (***************, ***************)
739(***************, ***************), (***************, ***************), (***************, ***************)
740(***************, ***************), (***************, ***************), (***************, ***************)
741alpha = (2._TKG, 3._TKG); beta = (1._TKG, 6._TKG); nrow = 2; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
742call setMatMulAdd(matA, matB, symmetric, uppDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
743matC
744( -96.000000, +72.000000), ( -141.000000, -226.000000), ( -112.000000, +38.000000)
745( -230.000000, -269.000000), ( -133.000000, -23.000000), ( -272.000000, -198.000000)
746(***************, ***************), (***************, ***************), (***************, ***************)
747(***************, ***************), (***************, ***************), (***************, ***************)
748(***************, ***************), (***************, ***************), (***************, ***************)
749matC - refC
750( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
751( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
752( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
753( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
754( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
755
756
757matA
758( +1.000000, +1.000000), ( -3.000000, +2.000000), ( +3.000000, +3.000000)
759( +2.000000, +6.000000), ( +4.000000, +5.000000), ( -1.000000, +4.000000)
760(***************, ***************), (***************, ***************), (***************, ***************)
761matB
762( +1.000000, ***************), (***************, ***************), (***************, ***************)
763( +3.000000, +2.000000), ( +4.000000, ***************), (***************, ***************)
764( -1.000000, +6.000000), ( +1.000000, +4.000000), ( +2.000000, ***************)
765(***************, ***************), (***************, ***************), (***************, ***************)
766matC
767( +13.000000, +6.000000), ( -18.000000, +6.000000), ( +10.000000, +7.000000)
768( -11.000000, +8.000000), ( +11.000000, +1.000000), ( -4.000000, +2.000000)
769(***************, ***************), (***************, ***************), (***************, ***************)
770(***************, ***************), (***************, ***************), (***************, ***************)
771(***************, ***************), (***************, ***************), (***************, ***************)
772alpha = (2._TKG, 3._TKG); beta = (1._TKG, 6._TKG); nrow = 2; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
773call setMatMulAdd(matA, matB, hermitian, lowDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
774matC
775( -137.000000, +17.000000), ( -158.000000, -102.000000), ( -39.000000, +141.000000)
776( -154.000000, -77.000000), ( -63.000000, +186.000000), ( +159.000000, +104.000000)
777(***************, ***************), (***************, ***************), (***************, ***************)
778(***************, ***************), (***************, ***************), (***************, ***************)
779(***************, ***************), (***************, ***************), (***************, ***************)
780matC - refC
781( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
782( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
783( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
784( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
785( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
786
787
788!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
789! Subset symmetric real matrix-matrix multiplication.
790!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
791
792
793matA
794 +1.000000, +2.000000, -1.000000, -1.000000, +4.000000
795***************, +0.000000, +1.000000, +1.000000, -1.000000
796***************, ***************, -1.000000, +1.000000, +2.000000
797***************, ***************, ***************, +2.000000, +0.000000
798***************, ***************, ***************, ***************, -1.000000
799***************, ***************, ***************, ***************, ***************
800***************, ***************, ***************, ***************, ***************
801***************, ***************, ***************, ***************, ***************
802matB
803 +1.000000, -1.000000, +0.000000, +2.000000
804 +2.000000, +2.000000, -1.000000, -2.000000
805 +1.000000, +0.000000, -1.000000, +1.000000
806 -3.000000, -1.000000, +1.000000, -1.000000
807 +4.000000, +2.000000, -1.000000, +1.000000
808***************, ***************, ***************, ***************
809matC
810 +23.000000, +12.000000, -6.000000, +2.000000
811 -4.000000, -5.000000, +1.000000, +3.000000
812 +5.000000, +6.000000, -1.000000, -4.000000
813 -4.000000, +1.000000, +0.000000, -5.000000
814 +8.000000, -4.000000, -2.000000, +13.000000
815alpha = 2._TKG; beta = 1._TKG; nrow = 5; ncol = 4; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
816call setMatMulAdd(matA, symmetric, uppDia, matB, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
817matC
818 +69.000000, +36.000000, -18.000000, +6.000000
819 -12.000000, -15.000000, +3.000000, +9.000000
820 +15.000000, +18.000000, -3.000000, -12.000000
821 -12.000000, +3.000000, +0.000000, -15.000000
822 +8.000000, -20.000000, -2.000000, +35.000000
823matC - refC
824 +0.000000, +0.000000, +0.000000, +0.000000
825 +0.000000, +0.000000, +0.000000, +0.000000
826 +0.000000, +0.000000, +0.000000, +0.000000
827 +0.000000, +0.000000, +0.000000, +0.000000
828 +0.000000, +0.000000, +0.000000, +0.000000
829
830matC
831 +23.000000, +12.000000, -6.000000, +2.000000
832 -4.000000, -5.000000, +1.000000, +3.000000
833 +5.000000, +6.000000, -1.000000, -4.000000
834 -4.000000, +1.000000, +0.000000, -5.000000
835 +8.000000, -4.000000, -2.000000, +13.000000
836call setMatMulAdd(matA(1:5, 1:5), symmetric, uppDia, matB(1:5, 1:4), matC, alpha, beta)
837matC
838 +69.000000, +36.000000, -18.000000, +6.000000
839 -12.000000, -15.000000, +3.000000, +9.000000
840 +15.000000, +18.000000, -3.000000, -12.000000
841 -12.000000, +3.000000, +0.000000, -15.000000
842 +8.000000, -20.000000, -2.000000, +35.000000
843matC - refC
844 +0.000000, +0.000000, +0.000000, +0.000000
845 +0.000000, +0.000000, +0.000000, +0.000000
846 +0.000000, +0.000000, +0.000000, +0.000000
847 +0.000000, +0.000000, +0.000000, +0.000000
848 +0.000000, +0.000000, +0.000000, +0.000000
849
850
851matA
852 +1.000000, ***************, ***************
853 +2.000000, +4.000000, ***************
854 +1.000000, -1.000000, -1.000000
855***************, ***************, ***************
856matB
857 +1.000000, -3.000000, +2.000000, +2.000000, -1.000000, +2.000000
858 +2.000000, +4.000000, +0.000000, +0.000000, +1.000000, -2.000000
859 +1.000000, -1.000000, -1.000000, -1.000000, -1.000000, +1.000000
860matC
861 +6.000000, +4.000000, +1.000000, +1.000000, +0.000000, -1.000000
862 +9.000000, +11.000000, +5.000000, +5.000000, +3.000000, -5.000000
863 -2.000000, -6.000000, +3.000000, +3.000000, -1.000000, +32.000000
864***************, ***************, ***************, ***************, ***************, ***************
865***************, ***************, ***************, ***************, ***************, ***************
866alpha = 2._TKG; beta = 2._TKG; nrow = 3; ncol = 6; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
867call setMatMulAdd(matA, symmetric, lowDia, matB, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
868matC
869 +24.000000, +16.000000, +4.000000, +4.000000, +0.000000, -4.000000
870 +36.000000, +44.000000, +20.000000, +20.000000, +12.000000, -20.000000
871 -8.000000, -24.000000, +12.000000, +12.000000, -4.000000, +70.000000
872***************, ***************, ***************, ***************, ***************, ***************
873***************, ***************, ***************, ***************, ***************, ***************
874matC - refC
875 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
876 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
877 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
878 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
879 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
880
881matC
882 +6.000000, +4.000000, +1.000000, +1.000000, +0.000000, -1.000000
883 +9.000000, +11.000000, +5.000000, +5.000000, +3.000000, -5.000000
884 -2.000000, -6.000000, +3.000000, +3.000000, -1.000000, +32.000000
885***************, ***************, ***************, ***************, ***************, ***************
886***************, ***************, ***************, ***************, ***************, ***************
887call setMatMulAdd(matA(1:3, 1:3), symmetric, lowDia, matB(1:3, 1:6), matC(1:3, 1:6), alpha, beta)
888matC
889 +24.000000, +16.000000, +4.000000, +4.000000, +0.000000, -4.000000
890 +36.000000, +44.000000, +20.000000, +20.000000, +12.000000, -20.000000
891 -8.000000, -24.000000, +12.000000, +12.000000, -4.000000, +70.000000
892***************, ***************, ***************, ***************, ***************, ***************
893***************, ***************, ***************, ***************, ***************, ***************
894matC - refC
895 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
896 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
897 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
898 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
899 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
900
901
902matA
903 +1.000000, -3.000000, +3.000000
904 +2.000000, +4.000000, -1.000000
905***************, ***************, ***************
906matB
907 +1.000000, -3.000000, +1.000000
908***************, +4.000000, -1.000000
909***************, ***************, +2.000000
910***************, ***************, ***************
911matC
912 +13.000000, -18.000000, +10.000000
913 -11.000000, +11.000000, -4.000000
914***************, ***************, ***************
915***************, ***************, ***************
916***************, ***************, ***************
917alpha = 2._TKG; beta = 1._TKG; nrow = 2; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
918call setMatMulAdd(matA, matB, symmetric, uppDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
919matC
920 +39.000000, -54.000000, +30.000000
921 -33.000000, +33.000000, -12.000000
922***************, ***************, ***************
923***************, ***************, ***************
924***************, ***************, ***************
925matC - refC
926 +0.000000, +0.000000, +0.000000
927 +0.000000, +0.000000, +0.000000
928 +0.000000, +0.000000, +0.000000
929 +0.000000, +0.000000, +0.000000
930 +0.000000, +0.000000, +0.000000
931
932matC
933 +13.000000, -18.000000, +10.000000
934 -11.000000, +11.000000, -4.000000
935***************, ***************, ***************
936***************, ***************, ***************
937***************, ***************, ***************
938call setMatMulAdd(matA(1:2, 1:3), matB(1:3, 1:3), symmetric, uppDia, matC(1:2, 1:3), alpha, beta)
939matC
940 +39.000000, -54.000000, +30.000000
941 -33.000000, +33.000000, -12.000000
942***************, ***************, ***************
943***************, ***************, ***************
944***************, ***************, ***************
945matC - refC
946 +0.000000, +0.000000, +0.000000
947 +0.000000, +0.000000, +0.000000
948 +0.000000, +0.000000, +0.000000
949 +0.000000, +0.000000, +0.000000
950 +0.000000, +0.000000, +0.000000
951
952
953matA
954 +1.000000, -3.000000, +2.000000
955 +2.000000, +4.000000, +0.000000
956 +1.000000, -1.000000, -1.000000
957matB
958 +1.000000, ***************, ***************
959 +2.000000, +10.000000, ***************
960 +1.000000, +11.000000, +4.000000
961matC
962 +1.000000, +5.000000, -9.000000
963 -3.000000, +10.000000, -2.000000
964 -2.000000, +8.000000, +0.000000
965alpha = -1._TKG; beta = 1._TKG; nrow = 3; ncol = 3; roffA = 0; coffA = 0; roffB = 0; coffB = 0; roffC = 0; coffC = 0;
966call setMatMulAdd(matA, matB, symmetric, lowDia, matC, alpha, beta, nrow, ncol, roffA, coffA, roffB, coffB, roffC, coffC)
967matC
968 +4.000000, +11.000000, +15.000000
969 -13.000000, -34.000000, -48.000000
970 +0.000000, +27.000000, +14.000000
971matC - refC
972 +0.000000, +0.000000, +0.000000
973 +0.000000, +0.000000, +0.000000
974 +0.000000, +0.000000, +0.000000
975
976matC
977 +1.000000, +5.000000, -9.000000
978 -3.000000, +10.000000, -2.000000
979 -2.000000, +8.000000, +0.000000
980call setMatMulAdd(matA(1:3, 1:3), matB(1:3, 1:3), symmetric, lowDia, matC(1:3, 1:3), alpha, beta)
981matC
982 +4.000000, +11.000000, +15.000000
983 -13.000000, -34.000000, -48.000000
984 +0.000000, +27.000000, +14.000000
985matC - refC
986 +0.000000, +0.000000, +0.000000
987 +0.000000, +0.000000, +0.000000
988 +0.000000, +0.000000, +0.000000
989
990
991!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
992! Subset symmetric real matrix-vector multiplication.
993!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
994
995
996matA
997 +8.000000, ***************, ***************
998 +4.000000, +6.000000, ***************
999 +2.000000, +7.000000, +3.000000
1000matB
1001 +3.000000, +2.000000, +1.000000
1002matC
1003 +5.000000, ***************, +3.000000, ***************, +2.000000
1004alpha = 1._TKG; beta = 1._TKG; ndim = 3; roffA = 0; coffA = 0; incB = 1; incC = 2;
1005call setMatMulAdd(matA, symmetric, lowDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)
1006matC
1007 +39.000000, ***************, +34.000000, ***************, +25.000000
1008matC - refC
1009 +0.000000, +0.000000, +0.000000, +0.000000, +0.000000
1010
1011matB
1012 +3.000000, +2.000000, +1.000000
1013matC
1014 +5.000000, +3.000000, +2.000000
1015alpha = 1._TKG; beta = 1._TKG;
1016call setMatMulAdd(matA, symmetric, lowDia, matB, matC, alpha, beta)
1017matC
1018 +39.000000, +34.000000, +25.000000
1019matC - refC
1020 +0.000000, +0.000000, +0.000000
1021
1022
1023matA
1024 +8.000000, +4.000000, +2.000000
1025***************, +6.000000, +7.000000
1026***************, ***************, +3.000000
1027matB
1028 +4.000000, ***************, +2.000000, ***************, +1.000000
1029matC
1030 +6.000000, +5.000000, +4.000000
1031alpha = 1._TKG; beta = 2._TKG; ndim = 3; roffA = 0; coffA = 0; incB = -2; incC = 1;
1032call setMatMulAdd(matA, symmetric, uppDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)
1033matC
1034 +36.000000, +54.000000, +36.000000
1035matC - refC
1036 +0.000000, +0.000000, +0.000000
1037
1038matB
1039 +1.000000, +2.000000, +4.000000
1040matC
1041 +6.000000, +5.000000, +4.000000
1042alpha = 1._TKG; beta = 2._TKG;
1043call setMatMulAdd(matA, symmetric, uppDia, matB, matC, alpha, beta)
1044matC
1045 +36.000000, +54.000000, +36.000000
1046matC - refC
1047 +0.000000, +0.000000, +0.000000
1048
1049
1050!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1051! Subset Hermitian complex matrix-vector multiplication.
1052!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1053
1054
1055matA
1056( +1.000000, +0.000000), (***************, +0.000000), (***************, +0.000000)
1057( +3.000000, +0.000000), ( +7.000000, +0.000000), (***************, +0.000000)
1058( +2.000000, +0.000000), ( +4.000000, +0.000000), ( +6.000000, +0.000000)
1059matB
1060( +1.000000, +2.000000), ( +4.000000, +0.000000), ( +3.000000, +4.000000)
1061matC
1062( +1.000000, +0.000000), (***************, ***************), ( +2.000000, -1.000000), (***************, ***************), ( +2.000000, +1.000000)
1063alpha = (1._TKG, 0._TKG); beta = (1._TKG, 0._TKG); ndim = 3; roffA = 0; coffA = 0; incB = 1; incC = 2;
1064call setMatMulAdd(matA, hermitian, lowDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)
1065matC
1066( +20.000000, +10.000000), (***************, ***************), ( +45.000000, +21.000000), (***************, ***************), ( +38.000000, +29.000000)
1067matC - refC
1068( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
1069
1070matB
1071( +1.000000, +2.000000), ( +4.000000, +0.000000), ( +3.000000, +4.000000)
1072matC
1073( +1.000000, +0.000000), ( +2.000000, -1.000000), ( +2.000000, +1.000000)
1074call setMatMulAdd(matA, hermitian, lowDia, matB, matC, alpha, beta)
1075matC
1076( +20.000000, +10.000000), ( +45.000000, +21.000000), ( +38.000000, +29.000000)
1077matC - refC
1078( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
1079
1080
1081matA
1082(***************, +0.000000), (***************, +0.000000), (***************, +0.000000), (***************, +0.000000), (***************, +0.000000)
1083(***************, +0.000000), (***************, +0.000000), ( +1.000000, +0.000000), ( +3.000000, +0.000000), ( +2.000000, +0.000000)
1084(***************, +0.000000), (***************, +0.000000), (***************, +0.000000), ( +7.000000, +0.000000), ( +4.000000, +0.000000)
1085(***************, +0.000000), (***************, +0.000000), (***************, +0.000000), (***************, +0.000000), ( +6.000000, +0.000000)
1086matB
1087( +3.000000, +4.000000), (***************, ***************), ( +4.000000, +0.000000), (***************, ***************), ( +1.000000, +2.000000)
1088matC
1089(***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************), (***************, ***************)
1090alpha = (1._TKG, 0._TKG); beta = (0._TKG, 0._TKG); ndim = 3; roffA = 1; coffA = 2; incB = -2; incC = 2;
1091call setMatMulAdd(matA, hermitian, uppDia, matB, matC, alpha, beta, ndim, roffA, coffA, incB, incC)
1092matC
1093( +19.000000, +10.000000), (***************, ***************), ( +43.000000, +22.000000), (***************, ***************), ( +36.000000, +28.000000)
1094matC - refC
1095( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
1096
1097matA
1098( +1.000000, +0.000000), ( +3.000000, +0.000000), ( +2.000000, +0.000000)
1099(***************, +0.000000), ( +7.000000, +0.000000), ( +4.000000, +0.000000)
1100(***************, +0.000000), (***************, +0.000000), ( +6.000000, +0.000000)
1101matB
1102( +1.000000, +2.000000), ( +4.000000, +0.000000), ( +3.000000, +4.000000)
1103matC
1104(***************, ***************), (***************, ***************), (***************, ***************)
1105call setMatMulAdd(matA, hermitian, uppDia, matB, matC, alpha, beta)
1106matC
1107( +19.000000, +10.000000), ( +43.000000, +22.000000), ( +36.000000, +28.000000)
1108matC - refC
1109( +0.000000, +0.000000), ( +0.000000, +0.000000), ( +0.000000, +0.000000)
1110
1111
Test:
test_pm_matrixMulAdd
Internal naming convention:
The following illustrates the internal naming convention used for the procedures within this generic interface.
XXXX_ASS_CNA_CNB_SFA_SFB_TNA_TNB_IK5()
|||| ||| ||| ||| ||| ||| ||| ||| |||
|||| ||| ||| ||| ||| ||| ||| ||| The type and kind of the input matrix arguments.
|||| ||| ||| ||| ||| ||| ||| The transposition of `matB`: TNB/TSB/THB => Transposition None / Transposition Symmetric / Transposition Hermitian
|||| ||| ||| ||| ||| ||| The transposition of `matA`: TNA/TSA/THA => Transposition None / Transposition Symmetric / Transposition Hermitian
|||| ||| ||| ||| ||| The subset of `matB` used: SFB/SUB/SLB => Subset Full / Subset Upper / Subset Lower
|||| ||| ||| ||| The subset of `matA` used: SFA/SUA/SLA => Subset Full / Subset Upper / Subset Lower
|||| ||| ||| The class of `matB`: CNB/CSB/CHB => General (None) / Symmetric / Hermitian
|||| ||| The class of `matA`: CNA/CSA/CHA => General (None) / Symmetric / Hermitian
|||| The explicitness of the array bounds: ASS/EXP => ASSumed bounds passed / EXPlicit bounds passed
The BLAS routine implemented: gemv, symv, hemv, gemm, symm, hemm, spmv, hpmv


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, Friday 1:54 AM, April 21, 2017, Institute for Computational Engineering and Sciences (ICES), The University of Texas, Austin, TX

Definition at line 485 of file pm_matrixMulAdd.F90.


The documentation for this interface was generated from the following file: