ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_matrixInit Module Reference

This module contains procedures and generic interfaces relevant to generating and initializing matrices of arbitrary shapes (:, :).
More...

Data Types

interface  getMatInit
 Generate and return a matrix of shape (shape(1), shape(2)) with the upper/lower triangle and the diagonal elements of the matrix set to the corresponding requested input values.
More...
 
interface  setMatInit
 Set the upper/lower triangle and the diagonal elements of the input matrix of arbitrary shape (:,:) to the requested input values.
More...
 

Variables

character(*, SK), parameter MODULE_NAME = "@pm_matrixInit"
 

Detailed Description

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.

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

1! Test the performance of `getMatInit()` vs. `setMatInit()`.
2program benchmark
3
4 use iso_fortran_env, only: error_unit
5 use pm_bench, only: bench_type
6 use pm_kind, only: IK, RKG => RK, RK, SK
7
8 implicit none
9
10 integer(IK) :: i
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)
18
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)
22
23
24 write(*,"(*(g0,:,' '))")
25 write(*,"(*(g0,:,' '))") "setMatInit() vs. getMatInit() vs. getMatInit_D2_reshape()"
26 write(*,"(*(g0,:,' '))")
27
28 open(newunit = fileUnit, file = "main.out", status = "replace")
29
30 write(fileUnit, "(*(g0,:,', '))") "MatrixRank", (bench(i)%name, i = 1, NBENCH)
31
32 loopOverMatrixRank: do irank = 1, NRANK
33
34 rank = 2_IK**irank
35 allocate(MatInit(rank, rank))
36 write(*,"(*(g0,:,' '))") "Benchmarking with rank", rank
37
38 do i = 1, NBENCH
39 bench(i)%timing = bench(i)%getTiming(minsec = 0.07_RK)
40 end do
41
42 write(fileUnit,"(*(g0,:,', '))") rank, (bench(i)%timing%mean, i = 1, NBENCH)
43 deallocate(MatInit)
44
45 end do loopOverMatrixRank
46 write(*,"(*(g0,:,' '))") dummySum
47 write(*,"(*(g0,:,' '))")
48
49 close(fileUnit)
50
51contains
52
53 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
54 ! procedure wrappers.
55 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
56
57 subroutine setOverhead()
58 call getDummy()
59 end subroutine
60
61 subroutine getDummy()
62 dummySum = dummySum + MatInit(1,1)
63 end subroutine
64
65 subroutine setMatInit()
66 block
67 use pm_matrixInit, only: setMatInit, dia_type
68 call setMatInit(MatInit, dia_type(), 1._RKG, rank, 0_IK, 0_IK)
69 call getDummy()
70 end block
71 end subroutine
72
73 subroutine getMatInit()
74 block
75 use pm_matrixInit, only: getMatInit, dia
76 MatInit = getMatInit([rank, rank], dia, 1._RKG)
77 call getDummy()
78 end block
79 end subroutine
80
81 subroutine getMatInit_D2_reshape()
82 MatInit = reshape_looping(rank)
83 call getDummy()
84 end subroutine
85
86 pure function reshape_looping(n) result(MatInit)
87 integer(IK), intent(in) :: n
88 real(RKG) :: MatInit(n,n)
89 integer(IK) :: k, j
90 MatInit = reshape([1._RKG, ([(0._RKG, k = 1, n)], 1._RKG, j = 1, n - 1)], shape(MatInit))
91 end function
92
93end program benchmark
Generate and return an object of type timing_type containing the benchmark timing information and sta...
Definition: pm_bench.F90:574
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...
Definition: pm_bench.F90:41
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
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
This module contains procedures and generic interfaces relevant to generating and initializing matric...
This is the class for creating benchmark and performance-profiling objects.
Definition: pm_bench.F90:386
subroutine bench(sort, arraySize)

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

Postprocessing of the benchmark output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6
7import os
8dirname = os.path.basename(os.getcwd())
9
10fontsize = 14
11
12df = pd.read_csv("main.out", delimiter = ", ")
13colnames = list(df.columns.values)
14
15
18
19ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
20ax = plt.subplot()
21
22for colname in colnames[1:]:
23 plt.plot( df[colnames[0]].values
24 , df[colname].values
25 , linewidth = 2
26 )
27
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)
33ax.set_xscale("log")
34ax.set_yscale("log")
35plt.minorticks_on()
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:]
40 #, loc='center left'
41 #, bbox_to_anchor=(1, 0.5)
42 , fontsize = fontsize
43 )
44
45plt.tight_layout()
46plt.savefig("benchmark." + dirname + ".runtime.png")
47
48
51
52ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
53ax = plt.subplot()
54
55plt.plot( df[colnames[0]].values
56 , np.ones(len(df[colnames[0]].values))
57 , linestyle = "--"
58 #, color = "black"
59 , linewidth = 2
60 )
61for colname in colnames[2:]:
62 plt.plot( df[colnames[0]].values
63 , df[colname].values / df[colnames[1]].values
64 , linewidth = 2
65 )
66
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)
72ax.set_xscale("log")
73#ax.set_yscale("log")
74plt.minorticks_on()
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:]
79 #, bbox_to_anchor = (1, 0.5)
80 #, loc = "center left"
81 , fontsize = fontsize
82 )
83
84plt.tight_layout()
85plt.savefig("benchmark." + dirname + ".runtime.ratio.png")

Visualization of the benchmark output

Benchmark moral
  1. 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.
  2. 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.

  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

Variable Documentation

◆ MODULE_NAME

character(*, SK), parameter pm_matrixInit::MODULE_NAME = "@pm_matrixInit"

Definition at line 129 of file pm_matrixInit.F90.