Line data Source code
1 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3 : !!!! !!!!
4 : !!!! ParaMonte: Parallel Monte Carlo and Machine Learning Library. !!!!
5 : !!!! !!!!
6 : !!!! Copyright (C) 2012-present, The Computational Data Science Lab !!!!
7 : !!!! !!!!
8 : !!!! This file is part of the ParaMonte library. !!!!
9 : !!!! !!!!
10 : !!!! LICENSE !!!!
11 : !!!! !!!!
12 : !!!! https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md !!!!
13 : !!!! !!!!
14 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16 :
17 : !> \brief
18 : !> This file contains the implementation details of the routines under the generic interfaces of [pm_arrayCenter](@ref pm_arrayCenter).
19 : !>
20 : !> \finmain
21 : !>
22 : !> \author
23 : !> \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : !> \bug
28 : !> gfortran as of version 10.3 cannot handle regular allocation for assumed-length allocatable character types and returns the following error:
29 : !> Fortran runtime error: Integer overflow when calculating the amount of memory to allocate
30 : !> The following preprocessor condition bypasses gfortran's bug.
31 : !#if getCenteredAsis_D0_SK_ENABLED || getCenteredMarg_D0_SK_ENABLED || setCenteredAsis_D0_SK_ENABLED || setCenteredMarg_D0_SK_ENABLED
32 : !#define ALLOCATE_NEW_WITH_STAT allocate(character(size + lmsize + rmsize, SKC) :: arrayCentered, stat = stat)
33 : !#define ALLOCATE_NEW allocate(character(size + lmsize + rmsize, SKC) :: arrayCentered)
34 : !#elif getCenteredAsis_D1_SK_ENABLED || getCenteredMarg_D1_SK_ENABLED || setCenteredAsis_D1_SK_ENABLED || setCenteredMarg_D1_SK_ENABLED
35 : !#define ALLOCATE_NEW_WITH_STAT allocate(character(len(array),SKC) :: arrayCentered(lb : lb + size + lmsize + rmsize - 1_IK), stat = stat)
36 : !#define ALLOCATE_NEW allocate(character(len(array),SKC) :: arrayCentered(lb : lb + size + lmsize + rmsize - 1_IK))
37 : !#else
38 : !#define ALLOCATE_NEW_WITH_STAT allocate(arrayCentered(lb : lb + size + lmsize + rmsize - 1_IK), stat = stat)
39 : !#define ALLOCATE_NEW allocate(arrayCentered(lb : lb + size + lmsize + rmsize - 1_IK))
40 : !#endif
41 : integer(IK) :: i, lenArray, cbeg, cend ! content beginning location in the new array.
42 : ! Set the default margins.
43 : #if Asis_ENABLED
44 : integer(IK), parameter :: lmsize = 0_IK, rmsize = 0_IK
45 : #elif !Marg_ENABLED
46 : #error "Unrecognized interface."
47 : #endif
48 : ! Set the array bounds.
49 : #if D0_ENABLED && SK_ENABLED
50 : #define GET_SIZE(ARR) len(ARR, kind = IK)
51 : #define GET_INDEX(i) i:i
52 : #elif D1_ENABLED
53 : #define GET_SIZE(ARR) ubound(ARR, 1, IK)
54 : #define GET_INDEX(i) i
55 : #else
56 : #error "Unrecognized interface."
57 : #endif
58 : integer(IK) , parameter :: lb = 1_IK
59 : #if setCentered_ENABLED
60 : integer(IK) :: size ! size of arrayCentered
61 782 : size = GET_SIZE(arrayCentered) - lmsize - rmsize ! fpp
62 : #endif
63 1562 : lenArray = GET_SIZE(array) ! fpp
64 : ! Verify the validity of the input.
65 3250 : CHECK_ASSERTION(__LINE__, size >= 0_IK, SK_"@getCentered()/setCentered(): The condition `size >= 0_IK` must hold. size = "//getStr(size)) ! fpp
66 : #if Marg_ENABLED
67 2846 : CHECK_ASSERTION(__LINE__, lmsize >= 0_IK, SK_"@getCentered()/setCentered(): The condition `lmsize >= 0_IK` must hold. lmsize = "//getStr(lmsize)) ! fpp
68 2846 : CHECK_ASSERTION(__LINE__, rmsize >= 0_IK, SK_"@getCentered()/setCentered(): The condition `rmsize >= 0_IK` must hold. rmsize = "//getStr(rmsize)) ! fpp
69 : #endif
70 : #if setCentered_ENABLED
71 9592 : CHECK_ASSERTION(__LINE__, size >= 0_IK, \
72 : SK_"@setCentered(): `GET_SIZE(arrayCentered) - lmsize - rmsize >= 0` must hold. GET_SIZE(arrayCentered), lmsize, rmsize = "\
73 : //getStr([GET_SIZE(arrayCentered), lmsize, rmsize])) ! fpp
74 : #endif
75 : #if setCentered_ENABLED && D1_ENABLED && SK_ENABLED
76 132 : CHECK_ASSERTION(__LINE__, len(array,IK) == len(arrayCentered,IK), \
77 : SK_"@setCentered(): `len(array) == len(arrayCentered)` must hold. len(array), len(arrayCentered) = "//\
78 : getStr([len(array,IK), len(arrayCentered,IK)])) ! fpp
79 : #endif
80 : ! Fill the left margin, if any.
81 : #if Marg_ENABLED
82 2846 : if (present(lmfill)) then
83 : do concurrent(i = 1 : lb + lmsize - 1_IK)
84 8615 : arrayCentered(GET_INDEX(i)) = lmfill
85 : end do
86 : end if
87 : #endif
88 : ! Center the array contents in the new array.
89 3250 : if (size > lenArray) then
90 2406 : cbeg = lb + lmsize + (size - lenArray) / 2_IK
91 2406 : cend = cbeg + lenArray - 1_IK
92 2406 : if (present(fill)) then
93 1878 : do concurrent(i = lb + lmsize : cbeg - 1_IK)
94 81791 : arrayCentered(GET_INDEX(i)) = fill
95 : end do
96 3278 : arrayCentered(cbeg:cend) = array
97 1894 : do concurrent(i = cend + 1_IK : lb + lmsize + size - 1_IK)
98 82523 : arrayCentered(GET_INDEX(i)) = fill
99 : end do
100 : else
101 1540 : arrayCentered(cbeg:cend) = array
102 : end if
103 : else
104 844 : cbeg = lb + (lenArray - size) / 2_IK
105 672 : cend = cbeg + size - 1_IK
106 2804 : arrayCentered(lb + lmsize : lb + lmsize + size - 1_IK) = array(cbeg:cend)
107 : end if
108 : ! Fill the right margin, if any.
109 : #if Marg_ENABLED
110 2846 : if (present(rmfill)) then
111 2206 : do concurrent(i = lb + lmsize + size : lb + lmsize + size + rmsize - 1_IK)
112 9310 : arrayCentered(GET_INDEX(i)) = rmfill
113 : end do
114 : end if
115 : #endif
116 :
117 : #undef GET_INDEX
118 : #undef GET_SIZE
|