This module contains procedures and generic interfaces relevant to generating and initializing matrices of arbitrary shapes (:, :)
.
More...
This module contains procedures and generic interfaces relevant to generating and initializing matrices of arbitrary shapes (:, :)
.
The following guidelines describe the various existing methods within the ParaMonte library for initializing matrices.
-
Use the interface getMatInit(..., subset = uppLowDia, ...) or setMatInit(..., subset = uppLowDia, ...) to initialize the upper, diagonal, and lower elements of a rectangular subset of a new or existing matrix.
-
Use the interface getMatInit(..., subset = uppLow, ...) or setMatInit(..., subset = uppLow, ...) to initialize the lower and upper elements of a rectangular subset of a new or existing matrix.
-
Use the interface getMatInit(..., subset = lowDia, ...) or setMatInit(..., subset = lowDia, ...) to initialize the diagonal and lower elements of a rectangular subset of a new or existing matrix.
-
Use the interface getMatInit(..., subset = uppDia, ...) or setMatInit(..., subset = uppDia, ...) to initialize the diagonal and upper elements of a rectangular subset of a new or existing matrix.
-
Use the interface getMatInit(..., subset = vdia, ...) or setMatInit(..., subset = vdia, ...) to initialize the diagonal elements of a rectangular subset of a new or existing matrix.
- BLAS/LAPACK equivalent:
- The procedures under discussion combine, modernize, and extend the interface and functionalities of Version 3.11 of BLAS/LAPACK routine(s):
SLASET
, DLASET
, CLASET
, ZLASET
.
- See also
- pm_arrayInit
pm_matrixCopy
pm_arrayCopy
pm_arrayCopy
- Benchmarks:
Benchmark :: The runtime performance of getMatInit vs. setMatInit ⛓
4 use iso_fortran_env,
only:
error_unit
11 integer(IK) :: fileUnit
12 integer(IK) :: rank, irank
13 integer(IK) ,
parameter :: NRANK
= 11_IK
14 integer(IK) ,
parameter :: NBENCH
= 3_IK
15 real(RKG) :: dummySum
= 0._RKG
16 real(RKG) ,
allocatable :: MatInit(:,:)
17 type(bench_type) :: bench(NBENCH)
19 bench(
1)
= bench_type(name
= SK_
"setMatInit", exec
= setMatInit , overhead
= setOverhead)
20 bench(
2)
= bench_type(name
= SK_
"getMatInit", exec
= getMatInit , overhead
= setOverhead)
21 bench(
3)
= bench_type(name
= SK_
"reshape_looping", exec
= getMatInit_D2_reshape , overhead
= setOverhead)
24 write(
*,
"(*(g0,:,' '))")
25 write(
*,
"(*(g0,:,' '))")
"setMatInit() vs. getMatInit() vs. getMatInit_D2_reshape()"
26 write(
*,
"(*(g0,:,' '))")
28 open(newunit
= fileUnit, file
= "main.out", status
= "replace")
30 write(fileUnit,
"(*(g0,:,', '))")
"MatrixRank", (bench(i)
%name, i
= 1, NBENCH)
32 loopOverMatrixRank:
do irank
= 1, NRANK
35 allocate(MatInit(rank, rank))
36 write(
*,
"(*(g0,:,' '))")
"Benchmarking with rank", rank
39 bench(i)
%timing
= bench(i)
%getTiming(minsec
= 0.07_RK)
42 write(fileUnit,
"(*(g0,:,', '))") rank, (bench(i)
%timing
%mean, i
= 1, NBENCH)
45 end do loopOverMatrixRank
46 write(
*,
"(*(g0,:,' '))") dummySum
47 write(
*,
"(*(g0,:,' '))")
57 subroutine setOverhead()
62 dummySum
= dummySum
+ MatInit(
1,
1)
65 subroutine setMatInit()
68 call setMatInit(MatInit, dia_type(),
1._RKG, rank,
0_IK,
0_IK)
73 subroutine getMatInit()
76 MatInit
= getMatInit([rank, rank], dia,
1._RKG)
81 subroutine getMatInit_D2_reshape()
82 MatInit
= reshape_looping(rank)
86 pure function reshape_looping(n)
result(MatInit)
87 integer(IK),
intent(in) :: n
88 real(RKG) :: MatInit(n,n)
90 MatInit
= reshape([
1._RKG, ([(
0._RKG, k
= 1, n)],
1._RKG, j
= 1, n
- 1)],
shape(MatInit))
Generate and return an object of type timing_type containing the benchmark timing information and sta...
Generate and return a matrix of shape (shape(1), shape(2)) with the upper/lower triangle and the diag...
Set the upper/lower triangle and the diagonal elements of the input matrix of arbitrary shape (:,...
This module contains abstract interfaces and types that facilitate benchmarking of different procedur...
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
This module contains procedures and generic interfaces relevant to generating and initializing matric...
This is the class for creating benchmark and performance-profiling objects.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Postprocessing of the benchmark output ⛓
3import matplotlib.pyplot
as plt
8dirname = os.path.basename(os.getcwd())
12df = pd.read_csv(
"main.out", delimiter =
", ")
13colnames = list(df.columns.values)
19ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
22for colname
in colnames[1:]:
23 plt.plot( df[colnames[0]].values
28plt.xticks(fontsize = fontsize)
29plt.yticks(fontsize = fontsize)
30ax.set_xlabel(colnames[0], fontsize = fontsize)
31ax.set_ylabel(
"Runtime [ seconds ]", fontsize = fontsize)
32ax.set_title(
" vs. ".join(colnames[1:])+
"\nLower is better.", fontsize = fontsize)
36plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
37ax.tick_params(axis =
"y", which =
"minor")
38ax.tick_params(axis =
"x", which =
"minor")
39ax.legend ( colnames[1:]
46plt.savefig(
"benchmark." + dirname +
".runtime.png")
52ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
55plt.plot( df[colnames[0]].values
56 , np.ones(len(df[colnames[0]].values))
61for colname
in colnames[2:]:
62 plt.plot( df[colnames[0]].values
63 , df[colname].values / df[colnames[1]].values
67plt.xticks(fontsize = fontsize)
68plt.yticks(fontsize = fontsize)
69ax.set_xlabel(colnames[0], fontsize = fontsize)
70ax.set_ylabel(
"Runtime compared to {}".format(colnames[1]), fontsize = fontsize)
71ax.set_title(
"Runtime Ratio Comparison. Lower means faster.\nLower than 1 means faster than {}().".format(colnames[1]), fontsize = fontsize)
75plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
76ax.tick_params(axis =
"y", which =
"minor")
77ax.tick_params(axis =
"x", which =
"minor")
78ax.legend ( colnames[1:]
85plt.savefig(
"benchmark." + dirname +
".runtime.ratio.png")
Visualization of the benchmark output ⛓
Benchmark moral ⛓
- The procedures under the generic interface getMatInit are functions while the procedures under the generic interface setMatInit are subroutines.
From the benchmark results, it appears that the functional interface performs less efficiently than the subroutine interface.
Note that this benchmark does not even include the cost of repeated reallcations, that is, the allocation of matrix happens only once in all tests prior to each benchmark.
- Furthermore, the
getMatInit()
implementation with Fortran-intrinsic reshape()
and implicit-loop appears to be significantly slower than both the subroutine and function implementations.
- Test:
- test_pm_matrixInit
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.
-
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.
-
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.
- Copyright
- Computational Data Science Lab
- Author:
- Amir Shahmoradi, Friday 1:54 AM, April 21, 2017, Institute for Computational Engineering and Sciences (ICES), The University of Texas, Austin, TX