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

This module contains procedures and generic interfaces for remapping arrays of various types. More...

Data Types

interface  getRemapped
 Generate a copy of the input array whose elements are reordered according to the input index array such that,
  arrayNew = array(index) or,
  arrayNew = array(index(size(index):1:-1)) if action = reverse or,
holds for the output. More...
 
interface  setRemapped
 Reorder the elements of the input array according to the input index array, such that
  array = array(index) or,
  array = array(index(size(index):1:-1)) if action = reverse or,
  arrayNew = array(index) if arrayNew is specified as input argument or,
  arrayNew = array(index(size(index):1:-1)) if arrayNew and action = reverse are specified as input arguments,
More...
 

Variables

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

Detailed Description

This module contains procedures and generic interfaces for remapping arrays of various types.

Benchmarks:


Benchmark :: The runtime performance of setRemapped vs. direct remapping

1program benchmark
2
3 use pm_kind, only: IK, LK, RK, SK
4 use iso_fortran_env, only: error_unit
6 use pm_bench, only: bench_type
7
8 implicit none
9
10 integer(IK) :: i
11 integer(IK) :: isize
12 integer(IK) :: fileUnit
13 integer(IK) , parameter :: NSIZE = 15_IK
14 integer(IK) , parameter :: NBENCH = 3_IK
15 integer(IK) :: arraySize(NSIZE)
16 real(RK) :: dummy = 0._RK
17 real(RK) , allocatable :: array(:)
18 integer(IK) , allocatable :: index(:)
19 type(bench_type) :: bench(NBENCH)
20
21 bench(1) = bench_type(name = SK_"setRemapped", exec = setRemapped, overhead = setOverhead)
22 bench(2) = bench_type(name = SK_"getRemapped", exec = getRemapped, overhead = setOverhead)
23 bench(3) = bench_type(name = SK_"direct", exec = direct , overhead = setOverhead)
24
25 arraySize = [( 2_IK**isize, isize = 1_IK, NSIZE )]
26
27 write(*,"(*(g0,:,' '))")
28 write(*,"(*(g0,:,' '))") "setRemapped() vs. getRemapped() vs. direct()"
29 write(*,"(*(g0,:,' '))")
30
31 open(newunit = fileUnit, file = "main.out", status = "replace")
32
33 write(fileUnit, "(*(g0,:,','))") "arraySize", (bench(i)%name, i = 1, NBENCH)
34
35 loopOverArraySize: do isize = 1, NSIZE
36
37 write(*,"(*(g0,:,' '))") "Benchmarking with size", arraySize(isize)
38
39 index = [( i, i = 1, arraySize(isize) )]
40 allocate(array(arraySize(isize)))
41 call random_number(array)
42 call setShuffled(index)
43 do i = 1, NBENCH
44 bench(i)%timing = bench(i)%getTiming(minsec = 0.07_RK)
45 end do
46 deallocate(array)
47
48 write(fileUnit,"(*(g0,:,','))") arraySize(isize), (bench(i)%timing%mean, i = 1, NBENCH)
49
50 end do loopOverArraySize
51 write(*,"(*(g0,:,' '))") dummy
52 write(*,"(*(g0,:,' '))")
53
54 close(fileUnit)
55
56contains
57
58 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
59 ! procedure wrappers.
60 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
61
62 subroutine setOverhead()
63 call finalize()
64 end subroutine
65
66 subroutine finalize()
67 dummy = dummy + array(1)
68 end subroutine
69
70 subroutine setRemapped()
71 block
72 use pm_arrayRemap, only: setRemapped, reverse
73 call setRemapped(array, index, reverse)
74 call finalize()
75 end block
76 end subroutine
77
78 subroutine getRemapped()
79 block
80 use pm_arrayRemap, only: getRemapped, reverse
81 array = getRemapped(array, index, reverse)
82 call finalize()
83 end block
84 end subroutine
85
86 subroutine direct()
87 array = array(index(arraySize(isize):1:-1))
88 call finalize()
89 end subroutine
90
91end program benchmark
Generate a copy of the input array whose elements are reordered according to the input index array su...
Reorder the elements of the input array according to the input index array, such that   array = ar...
Perform an unbiased random shuffling of the input array, known as the Knuth or Fisher-Yates shuffle.
Generate and return an object of type timing_type containing the benchmark timing information and sta...
Definition: pm_bench.F90:574
This module contains procedures and generic interfaces for remapping arrays of various types.
This module contains procedures and generic interfaces for shuffling arrays of various types.
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 LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
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 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
7fontsize = 14
8
9methods = ["setRemapped", "getRemapped", "direct"]
10
11df = pd.read_csv("main.out")
12
13
16
17ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
18ax = plt.subplot()
19
20for method in methods:
21 plt.plot( df["arraySize"].values
22 , df[method].values
23 , linewidth = 2
24 )
25
26plt.xticks(fontsize = fontsize)
27plt.yticks(fontsize = fontsize)
28ax.set_xlabel("Array Size", fontsize = fontsize)
29ax.set_ylabel("Runtime [ seconds ]", fontsize = fontsize)
30ax.set_title("setRemapped() vs. getRemapped() vs. direct method.\nLower is better.", fontsize = fontsize)
31ax.set_xscale("log")
32ax.set_yscale("log")
33plt.minorticks_on()
34plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
35ax.tick_params(axis = "y", which = "minor")
36ax.tick_params(axis = "x", which = "minor")
37ax.legend ( methods
38 #, loc='center left'
39 #, bbox_to_anchor=(1, 0.5)
40 , fontsize = fontsize
41 )
42
43plt.tight_layout()
44plt.savefig("benchmark.setRemapped_getRemapped_direct.runtime.png")
45
46
49
50ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
51ax = plt.subplot()
52
53plt.plot( df["arraySize"].values
54 , np.ones(len(df["arraySize"].values))
55 #, linestyle = "--"
56 , linewidth = 2
57 )
58plt.plot( df["arraySize"].values
59 , df["getRemapped"].values / df["setRemapped"].values
60 , linewidth = 2
61 )
62plt.plot( df["arraySize"].values
63 , df["direct"].values / df["setRemapped"].values
64 , linewidth = 2
65 )
66
67plt.xticks(fontsize = fontsize)
68plt.yticks(fontsize = fontsize)
69ax.set_xlabel("Array Size", fontsize = fontsize)
70ax.set_ylabel("Runtime compared to setRemapped()", fontsize = fontsize)
71ax.set_title("direct vs. getRemapped() to setRemapped() Runtime Ratio.\nLower means faster. Lower than 1 means faster than setRemapped().", 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 ( ["setRemapped()", "getRemapped()", "Direct Method"]
79 #, bbox_to_anchor = (1, 0.5)
80 #, loc = "center left"
81 , fontsize = fontsize
82 )
83
84plt.tight_layout()
85plt.savefig("benchmark.setRemapped_getRemapped_direct.runtime.ratio.png")

Visualization of the benchmark output

Benchmark moral
  1. The procedures under the generic interface setRemapped tend to be significantly faster than directly remapping arrays. This likely only true for remapping of allocatable arrays.
    The primary reason for the better performance of setRemapped is that setRemapped avoids a final data copy from a dummy array to the original array by copying the allocation descriptor instead of the whole remapped array to the original array.
  2. Note that the observed performance benefit slightly diminishes if the remapping is not in action = reverse mode.
  3. The other benefit of setRemapped is that it provides a unified seamless generic interface for reversing all intrinsic types and kinds of arrays as well as assumed-length and assumed-shape characters.
Test:
test_pm_arrayRemap


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:
Fatemeh Bagheri, Wednesday 12:20 AM, October 13, 2021, Dallas, TX

Variable Documentation

◆ MODULE_NAME

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

Definition at line 54 of file pm_arrayRemap.F90.