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

This module contains procedures and generic interfaces for replacing patterns within arrays of various types.
More...

Data Types

interface  getReplaced
 Generate and return an arrayNew of the same type and kind as the input array, in which the requested instances of the input pattern have been replaced with the input replacement. More...
 
interface  setReplaced
 Replace the requested instances of the input pattern with the input replacement in the allocatable input/output array. More...
 

Functions/Subroutines

pure recursive character(:, SK) function, allocatable getReplaced_recurs_alloc (array, pattern, replacement)
 [LEGACY code]
Generate and return a copy of the input array where all instances of the input pattern in it are replaced with the input replacement.
More...
 

Variables

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

Detailed Description

This module contains procedures and generic interfaces for replacing patterns within arrays of various types.

Benchmarks:


Benchmark :: The runtime performance of getReplaced vs. setReplaced

1! Test the performance of `getReplaced()` vs. `setReplaced()`.
2program benchmark
3
4 use iso_fortran_env, only: error_unit
5 use pm_kind, only: IK, LK, RK, SK
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 = 12_IK
14 integer(IK) , parameter :: NBENCH = 3_IK
15 integer(IK) :: arraySize(NSIZE)
16 logical(LK) :: dummy = .true._LK
17 character(:, SK), allocatable :: array
18 character(*, SK), parameter :: pattern = "a"
19 character(*, SK), parameter :: Replacement = "**"
20 type(bench_type) :: bench(NBENCH)
21
22 bench(1) = bench_type(name = SK_"setReplaced", exec = setReplaced , overhead = setOverhead)
23 bench(2) = bench_type(name = SK_"getReplaced", exec = getReplaced , overhead = setOverhead)
24 bench(3) = bench_type(name = SK_"getReplaced_recurs_alloc", exec = getReplaced_recurs_alloc , overhead = setOverhead)
25
26 arraySize = [( 2_IK**isize, isize = 1_IK, NSIZE )]
27
28 write(*,"(*(g0,:,' '))")
29 write(*,"(*(g0,:,' '))") "setReplaced() vs. getReplaced() vs. getReplaced_recurs_alloc()"
30 write(*,"(*(g0,:,' '))")
31
32 open(newunit = fileUnit, file = "main.out", status = "replace")
33
34 write(fileUnit, "(*(g0,:,','))") "arraySize", (bench(i)%name, i = 1, NBENCH)
35
36 loopOverArraySize: do isize = 1, NSIZE
37
38 write(*,"(*(g0,:,' '))") "Benchmarking with size", arraySize(isize)
39
40 do i = 1, NBENCH
41 bench(i)%timing = bench(i)%getTiming(minsec = 0.05_RK)
42 end do
43
44 write(fileUnit,"(*(g0,:,','))") arraySize(isize), (bench(i)%timing%mean, i = 1, NBENCH)
45
46 end do loopOverArraySize
47 write(*,"(*(g0,:,' '))") dummy
48 write(*,"(*(g0,:,' '))")
49
50 close(fileUnit)
51
52contains
53
54 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
55 ! procedure wrappers.
56 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
57
58 subroutine setOverhead()
59 call initialize()
60 call finalize()
61 end subroutine
62
63 subroutine initialize()
64 allocate( character(arraySize(isize)) :: array )
65 !block
66 ! use pm_distUnif, only: setUnifRand
67 ! call setUnifRand(array)
68 !end block
69 block
70 array(:) = repeat(pattern, len(array))
71 end block
72 end subroutine
73
74 subroutine finalize()
75 dummy = dummy .and. array(1:1) == pattern
76 deallocate(array)
77 end subroutine
78
79 subroutine setReplaced()
80 block
82 call initialize()
83 call setReplaced(array, pattern, Replacement)
84 call finalize()
85 end block
86 end subroutine
87
88 subroutine getReplaced()
89 block
91 call initialize()
92 array = getReplaced(array, pattern, Replacement)
93 call finalize()
94 end block
95 end subroutine
96
97 subroutine getReplaced_recurs_alloc()
98 block
100 call initialize()
101 array = getReplaced_recurs_alloc(array, pattern, Replacement)
102 call finalize()
103 end block
104 end subroutine
105
106end program benchmark
Generate and return an arrayNew of the same type and kind as the input array, in which the requested ...
Replace the requested instances of the input pattern with the input replacement in the allocatable in...
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 replacing patterns within arrays of variou...
pure recursive character(:, SK) function, allocatable getReplaced_recurs_alloc(array, pattern, replacement)
[LEGACY code] Generate and return a copy of the input array where all instances of the input pattern...
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 = ["setReplaced", "getReplaced", "getReplaced_recurs_alloc"]
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("setReplaced() vs. getReplaced() vs. getReplaced_recurs_alloc()\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.getReplaced_vs_setReplaced.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 #, color = "black"
57 , linewidth = 2
58 )
59plt.plot( df["arraySize"].values
60 , df["getReplaced"].values / df["setReplaced"].values
61 , linewidth = 2
62 )
63plt.plot( df["arraySize"].values
64 , df["getReplaced_recurs_alloc"].values / df["setReplaced"].values
65 , linewidth = 2
66 )
67
68plt.xticks(fontsize = fontsize)
69plt.yticks(fontsize = fontsize)
70ax.set_xlabel("Array Size", fontsize = fontsize)
71ax.set_ylabel("Runtime compared to setReplaced()", fontsize = fontsize)
72ax.set_title("getReplaced / setReplaced vs. getReplaced_recurs_alloc / setReplaced\nLower means faster. Lower than 1 means faster than setReplaced().", fontsize = fontsize)
73ax.set_xscale("log")
74#ax.set_yscale("log")
75plt.minorticks_on()
76plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
77ax.tick_params(axis = "y", which = "minor")
78ax.tick_params(axis = "x", which = "minor")
79ax.legend ( ["setReplaced", "getReplaced", "getReplaced_recurs_alloc"]
80 #, bbox_to_anchor = (1, 0.5)
81 #, loc = "center left"
82 , fontsize = fontsize
83 )
84
85plt.tight_layout()
86plt.savefig("benchmark.getReplaced_vs_setReplaced.runtime.ratio.png")

Visualization of the benchmark output

Benchmark moral
  1. The procedures under the generic interface getReplaced are functions while the procedures under the generic interface setReplaced are subroutines.
    From the benchmark results, it appears that the functional interface performs slightly less efficiently than the subroutine interface.
    Note that this benchmark does not even include the cost of repeated reallcations, that is, the allocation of Replaced happen only once in all tests.
  2. Furthermore, the recursive getReplaced() implementation with recursive allocations appears to be 3-33 times slower than the subroutine implementation, depending on the size of the array within which the pattern is to be replaced.
  3. Note that this benchmark considers the worst-case scenario where all elements of the input array match the input pattern and must be therefore, replaced.
See also
pm_arrayInsert
pm_arrayRemove
Test:
test_pm_arrayReplace


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

Function/Subroutine Documentation

◆ getReplaced_recurs_alloc()

pure recursive character(:, SK) function, allocatable pm_arrayReplace::getReplaced_recurs_alloc ( character(*, SK), intent(in)  array,
character(*, SK), intent(in)  pattern,
character(*, SK), intent(in)  replacement 
)

[LEGACY code]
Generate and return a copy of the input array where all instances of the input pattern in it are replaced with the input replacement.

Parameters
[in]array: The input scalar character of default kind SK.
[in]pattern: The input scalar the same type and kind as the input array representing the pattern to be searched within array.
[in]replacement: The input scalar the same type and kind as the input array representing the replacement for instances of pattern in array.
Returns
arrayNew : The output copy of the input array where all instances of the input pattern in it are replaced with the input replacement.
Warning
This procedure is solely provided for benchmarking purposes to compare the performance difference between a naive recursive allocation method of replacement and the optimized method implemented under the generic interfaces in getReplaced and setReplaced.
Remarks
The procedures under discussion are pure.
The procedures under discussion are recursive.
See also
getReplaced
setReplaced

test_pm_arrayReplace


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, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin

Definition at line 14612 of file pm_arrayReplace.F90.

References getReplaced_recurs_alloc(), and pm_kind::IK.

Referenced by getReplaced_recurs_alloc(), test_pm_arrayReplace::test_getReplaced_recurs_alloc_1(), test_pm_arrayReplace::test_getReplaced_recurs_alloc_2(), and test_pm_arrayReplace::test_getReplaced_recurs_alloc_3().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ MODULE_NAME

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

Definition at line 60 of file pm_arrayReplace.F90.