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_sampleNorm](@ref pm_sampleNorm).
19 : !>
20 : !> \finmain
21 : !>
22 : !> \author
23 : !> \AmirShahmoradi, Saturday 2:33 AM, August 22, 2021, Dallas, TX
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : ! Define the runtime checks.
28 : #define CHECK_VAL_DIM \
29 : CHECK_ASSERTION(__LINE__, 1 <= dim .and. dim <= rank(sample), \
30 : SK_"@setNormed(): The condition `1 <= dim .and. dim <= rank(sample)` must hold. dim, rank(sample) = "//getStr([integer(IK) :: dim, rank(sample)]))
31 : #define CHECK_LEN_SCALE(DIM) \
32 : CHECK_ASSERTION(__LINE__, size(sample, DIM, IK) == size(scale, 1, IK), \
33 : SK_"@setNormed(): The condition `size(sample, dim, 1) == size(scale, 1)` must hold. dim, shape(sample), size(scale) = "\
34 : //getStr([DIM, shape(sample, IK), size(scale, 1, IK)]))
35 : #define CHECK_LEN_SHIFT(DIM) \
36 : CHECK_ASSERTION(__LINE__, size(sample, DIM, IK) == size(shift, 1, IK), \
37 : SK_"@setNormed(): The condition `size(sample, dim, 1) == size(shift, 1)` must hold. dim, shape(sample), size(shift) = "\
38 : //getStr([DIM, shape(sample, IK), size(shift, 1, IK)]))
39 :
40 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
41 : #if getNormed_ENABLED && D1_ENABLED
42 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
43 :
44 108 : sampleNormed = (sample + shift) * scale
45 :
46 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
47 : #elif setNormed_ENABLED && D1_ENABLED
48 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
49 :
50 24 : sample = (sample + shift) * scale
51 :
52 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
53 : #elif getNormed_ENABLED && D2_ENABLED
54 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
55 :
56 :
57 : ! Set the indexing rule.
58 : #if ONO_ENABLED
59 : #define GET_INDEX(I,J)I,J
60 : #elif OTH_ENABLED
61 : #define GET_INDEX(I,J)J,I
62 : #else
63 : #error "Unrecognized interface."
64 : #endif
65 : ! Set the conjugation rule.
66 : #if OTH_ENABLED && CK_ENABLED
67 : #define GET_CONJG(X)conjg(X)
68 : #elif ONO_ENABLED || (OTH_ENABLED && RK_ENABLED)
69 : #define GET_CONJG(X)X
70 : #else
71 : #error "Unrecognized interface."
72 : #endif
73 : integer(IK) :: idim, isam, ndim, nsam
74 12 : CHECK_VAL_DIM
75 28 : CHECK_LEN_SCALE(3 - dim)
76 28 : CHECK_LEN_SHIFT(3 - dim)
77 4 : nsam = size(sample, dim, IK)
78 4 : ndim = size(sample, 3 - dim, IK)
79 4 : if (dim == 2_IK) then
80 : do concurrent(idim = 1 : ndim, isam = 1 : nsam)
81 58 : sampleNormed(GET_INDEX(idim, isam)) = GET_CONJG((sample(idim, isam) + shift(idim)) * scale(idim))
82 : end do
83 : else
84 : do concurrent(idim = 1 : ndim, isam = 1 : nsam)
85 0 : sampleNormed(GET_INDEX(isam, idim)) = GET_CONJG((sample(isam, idim) + shift(idim)) * scale(idim))
86 : end do
87 : end if
88 : #undef GET_CONJG
89 : #undef GET_INDEX
90 :
91 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
92 : #elif setNormed_ENABLED && D2_ENABLED
93 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
94 :
95 : integer(IK) :: idim, isam, ndim, nsam
96 6 : CHECK_VAL_DIM
97 14 : CHECK_LEN_SCALE(3 - dim)
98 14 : CHECK_LEN_SHIFT(3 - dim)
99 2 : nsam = size(sample, dim, IK)
100 2 : ndim = size(sample, 3 - dim, IK)
101 2 : if (dim == 2_IK) then
102 : do concurrent(isam = 1 : nsam, idim = 1 : ndim)
103 32 : sample(idim, isam) = (sample(idim, isam) + shift(idim)) * scale(idim)
104 : end do
105 : else
106 : do concurrent(idim = 1 : ndim, isam = 1 : nsam)
107 0 : sample(isam, idim) = (sample(isam, idim) + shift(idim)) * scale(idim)
108 : end do
109 : end if
110 :
111 : #else
112 : !%%%%%%%%%%%%%%%%%%%%%%%%
113 : #error "Unrecognized interface."
114 : !%%%%%%%%%%%%%%%%%%%%%%%%
115 : #endif
116 : #undef CHECK_LEN_SCALE
117 : #undef CHECK_LEN_SHIFT
118 : #undef CHECK_VAL_DIM
|