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

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. More...

Data Types

interface  setRefilled
 Allocate or resize (shrink or expand) and refill an input allocatable scalar string or array of rank 1..3 to an arbitrary size while preserving the original contents or a subset of it.
More...
 

Variables

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

Detailed Description

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,

  1. 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.


  2. 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.


  3. 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.


  4. 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

1program benchmark
2
3 use iso_fortran_env, only: error_unit
4 use pm_bench, only: bench_type
5 use pm_kind, only: IK, RK, SK
6
7 implicit none
8
9 integer(IK) :: i
10 integer(IK) :: isize
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)
18
19 bench(1) = bench_type(name = SK_"setRefilled", exec = setRefilled , overhead = setOverhead)
20 bench(2) = bench_type(name = SK_"direct", exec = direct , overhead = setOverhead)
21
22 arraySize = [( 2_IK**isize, isize = 1_IK, NSIZE )]
23
24 write(*,"(*(g0,:,' '))")
25 write(*,"(*(g0,:,' '))") "setRefilled() vs. direct()"
26 write(*,"(*(g0,:,' '))")
27
28 open(newunit = fileUnit, file = "main.out", status = "replace")
29
30 write(fileUnit, "(*(g0,:,','))") "arraySize", (bench(i)%name, i = 1, NBENCH)
31
32 loopOverArraySize: do isize = 1, NSIZE - 1
33
34 write(*,"(*(g0,:,' '))") "Benchmarking with size", arraySize(isize)
35
36 allocate(array(arraySize(isize)))
37 call random_number(array)
38 do i = 1, NBENCH
39 bench(i)%timing = bench(i)%getTiming(minsec = 0.1_RK, miniter = 20_IK)
40 end do
41 deallocate(array)
42
43 write(fileUnit,"(*(g0,:,','))") arraySize(isize), (bench(i)%timing%mean, i = 1, NBENCH)
44
45 end do loopOverArraySize
46 write(*,"(*(g0,:,' '))") dummy
47 write(*,"(*(g0,:,' '))")
48
49 close(fileUnit)
50
51contains
52
53 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
54 ! procedure wrappers.
55 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
56
57 subroutine setOverhead()
58 call initialize()
59 call finalize()
60 end subroutine
61
62 subroutine initialize()
64 call setShuffled(array(1:arraySize(isize)))
65 end subroutine
66
67 subroutine finalize()
68 dummy = dummy + array(1)
69 end subroutine
70
71 subroutine setRefilled()
72 block
74 call initialize()
75 call setRefilled(array, 0._RK, arraySize(isize+1))
76 call finalize()
77 end block
78 end subroutine
79
80 subroutine direct()
81 real(RK), allocatable :: Temp(:)
82 call initialize()
83 Temp = array(1:arraySize(isize))
84 deallocate(array)
85 allocate(array(arraySize(isize+1)))
86 array(1:arraySize(isize)) = Temp
87 array(arraySize(isize)+1:arraySize(isize+1)) = 0._RK
88 call finalize()
89 deallocate(Temp)
90 end subroutine
91
92end program benchmark
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...
Definition: pm_bench.F90:574
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...
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 is the class for creating benchmark and performance-profiling objects.
Definition: pm_bench.F90:386

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 = ["setRefilled", "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("setRefilled() 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.setRefilled_vs_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["direct"].values / df["setRefilled"].values
60 , linewidth = 2
61 )
62
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)
68ax.set_xscale("log")
69#ax.set_yscale("log")
70plt.minorticks_on()
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"]
75 #, bbox_to_anchor = (1, 0.5)
76 #, loc = "center left"
77 , fontsize = fontsize
78 )
79
80plt.tight_layout()
81plt.savefig("benchmark.setRefilled_vs_direct.runtime.ratio.png")

Visualization of the benchmark output

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

  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_arrayRefill::MODULE_NAME = "@pm_arrayRefill"

Definition at line 85 of file pm_arrayRefill.F90.