This module contains procedures and generic interfaces for resizing allocatable
arrays of various types, relocating their contents and filling the newly added elements with specific values.
In brief,
-
Use setResized to resize strings/arrays without changing the lower bound.
The following figure illustrates example resizing of a 1D array and transferal of its contents.
-
Use setRefilled to resize strings/arrays without changing the lower bound and to initialize the new elements.
The following figure illustrates example refilling of a 1D array and transferal of its contents.
-
Use setRebound to resize arrays by changing their lower and/or upper bounds.
The following figure illustrates example rebinding of a 1D array and transferal of its contents.
-
Use setRebilled to resize arrays by changing their lower and/or upper bounds and to initialize the new elements.
The following figure illustrates example rebinding and refilling of a 1D array and transferal of its contents.
- Benchmarks:
Benchmark :: The runtime performance of setRefilled vs. direct reversing ⛓
3 use iso_fortran_env,
only:
error_unit
11 integer(IK) :: fileUnit
12 integer(IK) ,
parameter :: NSIZE
= 21_IK
13 integer(IK) ,
parameter :: NBENCH
= 2_IK
14 integer(IK) :: arraySize(NSIZE)
15 real(RK) :: dummy
= 0._RK
16 real(RK) ,
allocatable :: array(:)
17 type(bench_type) :: bench(NBENCH)
19 bench(
1)
= bench_type(name
= SK_
"setRefilled", exec
= setRefilled , overhead
= setOverhead)
20 bench(
2)
= bench_type(name
= SK_
"direct", exec
= direct , overhead
= setOverhead)
22 arraySize
= [(
2_IK**isize, isize
= 1_IK, NSIZE )]
24 write(
*,
"(*(g0,:,' '))")
25 write(
*,
"(*(g0,:,' '))")
"setRefilled() vs. direct()"
26 write(
*,
"(*(g0,:,' '))")
28 open(newunit
= fileUnit, file
= "main.out", status
= "replace")
30 write(fileUnit,
"(*(g0,:,','))")
"arraySize", (bench(i)
%name, i
= 1, NBENCH)
32 loopOverArraySize:
do isize
= 1, NSIZE
- 1
34 write(
*,
"(*(g0,:,' '))")
"Benchmarking with size", arraySize(isize)
36 allocate(array(arraySize(isize)))
37 call random_number(array)
39 bench(i)
%timing
= bench(i)
%getTiming(minsec
= 0.1_RK, miniter
= 20_IK)
43 write(fileUnit,
"(*(g0,:,','))") arraySize(isize), (bench(i)
%timing
%mean, i
= 1, NBENCH)
45 end do loopOverArraySize
46 write(
*,
"(*(g0,:,' '))") dummy
47 write(
*,
"(*(g0,:,' '))")
57 subroutine setOverhead()
62 subroutine initialize()
68 dummy
= dummy
+ array(
1)
71 subroutine setRefilled()
81 real(RK),
allocatable :: Temp(:)
83 Temp
= array(
1:arraySize(isize))
85 allocate(array(arraySize(isize
+1)))
86 array(
1:arraySize(isize))
= Temp
87 array(arraySize(isize)
+1:arraySize(isize
+1))
= 0._RK
Allocate or resize (shrink or expand) and refill an input allocatable scalar string or array of rank ...
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...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
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...
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 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
9methods = [
"setRefilled",
"direct"]
11df = pd.read_csv(
"main.out")
17ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
21 plt.plot( df[
"arraySize"].values
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(
"setRefilled() vs. direct method.\nLower is better.", fontsize = fontsize)
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")
44plt.savefig(
"benchmark.setRefilled_vs_direct.runtime.png")
50ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
53plt.plot( df[
"arraySize"].values
54 , np.ones(len(df[
"arraySize"].values))
58plt.plot( df[
"arraySize"].values
59 , df[
"direct"].values / df[
"setRefilled"].values
63plt.xticks(fontsize = fontsize)
64plt.yticks(fontsize = fontsize)
65ax.set_xlabel(
"Array Size", fontsize = fontsize)
66ax.set_ylabel(
"Runtime compared to setRefilled()", fontsize = fontsize)
67ax.set_title(
"direct method to setRefilled() Runtime Ratio.\nLower means faster. Lower than 1 means faster than setRefilled().", fontsize = fontsize)
71plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
72ax.tick_params(axis =
"y", which =
"minor")
73ax.tick_params(axis =
"x", which =
"minor")
74ax.legend ( [
"setRefilled()",
"Direct Method"]
81plt.savefig(
"benchmark.setRefilled_vs_direct.runtime.ratio.png")
Visualization of the benchmark output ⛓
Benchmark moral ⛓
- The procedures under the generic interface setRefilled are about \(\ms{10-20%}\) faster than naively copying and resizing via temporary arrays and initializing the values, simply because setRefilled avoids an extra unnecessary reallocation and data copy.
- Test:
- test_pm_arrayRefill
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:
- Fatemeh Bagheri, Wednesday 12:20 AM, October 13, 2021, Dallas, TX