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_sampleAffinity](@ref pm_sampleAffinity).
19 : !>
20 : !> \finmain
21 : !>
22 : !> \author
23 : !> \AmirShahmoradi, Saturday 2:33 AM, August 22, 2021, Dallas, TX
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : ! Define the default mean.
28 : #if ATL_ENABLED
29 : #define CHECK_LEN_TLATE(DIM) \
30 : CHECK_ASSERTION(__LINE__, size(sample, DIM, IK) == size(tlate, 1, IK), \
31 : SK_"@setAffinity(): The condition `size(sample, dim, 1) == size(tlate, 1)` must hold. dim, shape(sample), size(tlate) = "\
32 : //getStr([DIM, shape(sample, IK), size(tlate, 1, IK)]))
33 : #define PLUS(X)+ X
34 : #elif DTL_ENABLED
35 : #define CHECK_LEN_TLATE(DIM)
36 : #define PLUS(X)
37 : #elif setAffinity_ENABLED
38 : #error "Unrecognized interface."
39 : #endif
40 : ! Define the runtime checks.
41 : #define CHECK_VAL_DIM \
42 : CHECK_ASSERTION(__LINE__, 1 <= dim .and. dim <= rank(sample), \
43 : SK_"@setAffinity(): The condition `1 <= dim .and. dim <= rank(sample)` must hold. dim, rank(sample) = "//getStr([integer(IK) :: dim, rank(sample)]))
44 : #define CHECK_SHAPE_SAMAFF \
45 : CHECK_ASSERTION(__LINE__, all(shape(affinity, IK) == shape(sample, IK)), \
46 : SK_"@setAffinity(): The condition `all(shape(affinity) == shape(sample))` must hold. shape(affinity), shape(sample) = "\
47 : //getStr([shape(affinity, IK), shape(sample, IK)]))
48 :
49 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
50 : #if getAffinity_ENABLED && D1_ENABLED
51 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
52 :
53 : #define CALL_SET_AFFINE(CLASS) \
54 : if (present(tlate)) then; \
55 : call setAffinity(affinity, sample, tform, CLASS, tlate); \
56 : else; \
57 : call setAffinity(affinity, sample, tform, CLASS); \
58 : end if;
59 230 : if (present(class)) then
60 : select type (class)
61 : type is (genrecmat_type)
62 0 : CALL_SET_AFFINE(class)
63 : type is (upperDiag_type)
64 46 : CALL_SET_AFFINE(class)
65 : type is (lowerDiag_type)
66 46 : CALL_SET_AFFINE(class)
67 : type is (upperUnit_type)
68 46 : CALL_SET_AFFINE(class)
69 : type is (lowerUnit_type)
70 46 : CALL_SET_AFFINE(class)
71 : class default
72 0 : error stop "Unrecognized `class`."
73 : end select
74 : return
75 : end if
76 46 : CALL_SET_AFFINE(genrecmat)
77 : #undef CALL_SET_AFFINE
78 :
79 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
80 : #elif getAffinity_ENABLED && D2_ENABLED
81 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
82 :
83 : #define CALL_SET_AFFINE(CLASS) \
84 : if (present(tlate)) then; \
85 : call setAffinity(affinity, sample, dim, tform, CLASS, tlate); \
86 : else; \
87 : call setAffinity(affinity, sample, dim, tform, CLASS); \
88 : end if;
89 108 : if (present(class)) then
90 : select type (class)
91 : type is (genrecmat_type)
92 2 : CALL_SET_AFFINE(class)
93 : type is (upperDiag_type)
94 20 : CALL_SET_AFFINE(class)
95 : type is (lowerDiag_type)
96 22 : CALL_SET_AFFINE(class)
97 : type is (upperUnit_type)
98 20 : CALL_SET_AFFINE(class)
99 : type is (lowerUnit_type)
100 20 : CALL_SET_AFFINE(class)
101 : class default
102 0 : error stop "Unrecognized `class`."
103 : end select
104 : return
105 : end if
106 24 : CALL_SET_AFFINE(genrecmat)
107 : #undef CALL_SET_AFFINE
108 :
109 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
110 : #elif setAffinity_ENABLED && D1_ENABLED
111 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
112 :
113 : #if !CGR_ENABLED
114 : integer(IK) :: idim
115 : #endif
116 : integer(IK) :: ndim
117 460 : ndim = size(sample, 1, IK)
118 0 : CHECK_LEN_TLATE(1)
119 2760 : CHECK_SHAPE_SAMAFF
120 460 : if (ndim < 1_IK) return
121 : #if CGR_ENABLED
122 870 : affinity = matmul(tform, sample) PLUS(tlate)
123 : #elif CUD_ENABLED
124 : !do idim = 1, ndim
125 : ! affinity(idim) = dot_product(sample(idim : ndim), tform(idim, idim : ndim)) PLUS(tlate(idim))
126 : !end do
127 274 : affinity(1 : ndim) = sample(ndim) * tform(1 : ndim, ndim) PLUS(tlate)
128 182 : do idim = ndim - 1, 1, -1
129 298 : affinity(1 : idim) = affinity(1 : idim) + sample(idim) * tform(1 : idim, idim)
130 : end do
131 : #elif CUU_ENABLED
132 : !do idim = 1, ndim
133 : ! affinity(idim) = dot_product(sample(idim : ndim), tform(idim, idim : ndim)) PLUS(tlate(idim))
134 : !end do
135 182 : affinity(1 : ndim - 1) = sample(ndim) * tform(1 : ndim - 1, ndim) PLUS(tlate(1 : ndim - 1))
136 92 : affinity(ndim) = sample(ndim) PLUS(tlate(ndim))
137 182 : do idim = ndim - 1, 1, -1
138 116 : affinity(1 : idim - 1) = affinity(1 : idim - 1) + sample(idim) * tform(1 : idim - 1, idim)
139 182 : affinity(idim) = affinity(idim) + sample(idim) * tform(idim, idim)
140 : end do
141 : #elif CLD_ENABLED
142 274 : affinity(1 : ndim) = sample(1) * tform(1 : ndim, 1) PLUS(tlate)
143 182 : do idim = 2, ndim
144 298 : affinity(idim : ndim) = affinity(idim : ndim) + sample(idim) * tform(idim : ndim, idim)
145 : end do
146 : #elif CLU_ENABLED
147 92 : affinity(1) = sample(1) PLUS(tlate(1))
148 182 : affinity(2 : ndim) = sample(1) * tform(2 : ndim, 1) PLUS(tlate(2 : ndim))
149 182 : do idim = 2, ndim
150 90 : affinity(idim) = affinity(idim) + sample(idim)
151 208 : affinity(idim + 1 : ndim) = affinity(idim + 1 : ndim) + sample(idim) * tform(idim + 1 : ndim, idim)
152 : end do
153 : #else
154 : #error "Unrecognized interface."
155 : #endif
156 :
157 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
158 : #elif setAffinity_ENABLED && D2_ENABLED
159 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
160 :
161 : #if !CGR_ENABLED
162 : integer(IK) :: idim
163 : #endif
164 : integer(IK) :: ndim, isam, nsam
165 648 : CHECK_VAL_DIM
166 2376 : CHECK_SHAPE_SAMAFF
167 56 : CHECK_LEN_TLATE(3 - dim)
168 216 : nsam = size(sample, dim, IK)
169 216 : ndim = size(sample, 3 - dim, IK)
170 216 : if (ndim < 1_IK) return
171 216 : if (dim == 2_IK) then
172 : do concurrent(isam = 1 : nsam)
173 : #if CGR_ENABLED
174 116840 : affinity(1 : ndim, isam) = matmul(tform, sample(1 : ndim, isam)) PLUS(tlate)
175 : #elif CUD_ENABLED
176 : !do idim = 1, ndim
177 : ! affinity(idim, isam) = dot_product(sample(idim : ndim, isam), tform(idim, idim : ndim)) PLUS(tlate(idim))
178 : !end do
179 268 : affinity(1 : ndim, isam) = sample(ndim, isam) * tform(1 : ndim, ndim) PLUS(tlate)
180 288 : do idim = ndim - 1, 1, -1
181 270 : affinity(1 : idim, isam) = affinity(1 : idim, isam) + sample(idim, isam) * tform(1 : idim, idim)
182 : end do
183 : #elif CUU_ENABLED
184 : !do idim = 1, ndim
185 : ! affinity(idim, isam) = dot_product(sample(idim : ndim, isam), tform(idim, idim : ndim)) PLUS(tlate(idim))
186 : !end do
187 174 : affinity(1 : ndim - 1, isam) = sample(ndim, isam) * tform(1 : ndim - 1, ndim) PLUS(tlate(1 : ndim - 1))
188 94 : affinity(ndim, isam) = sample(ndim, isam) PLUS(tlate(ndim))
189 288 : do idim = ndim - 1, 1, -1
190 96 : affinity(1 : idim - 1, isam) = affinity(1 : idim - 1, isam) + sample(idim, isam) * tform(1 : idim - 1, idim)
191 174 : affinity(idim, isam) = affinity(idim, isam) + sample(idim, isam) * tform(idim, idim)
192 : end do
193 : #elif CLD_ENABLED
194 12268 : affinity(1 : ndim, isam) = sample(1, isam) * tform(1 : ndim, 1) PLUS(tlate)
195 12292 : do idim = 2, ndim
196 12270 : affinity(idim : ndim, isam) = affinity(idim : ndim, isam) + sample(idim, isam) * tform(idim : ndim, idim)
197 : end do
198 : #elif CLU_ENABLED
199 94 : affinity(1, isam) = sample(1, isam) PLUS(tlate(1))
200 174 : affinity(2 : ndim, isam) = sample(1, isam) * tform(2 : ndim, 1) PLUS(tlate(2 : ndim))
201 288 : do idim = 2, ndim
202 80 : affinity(idim, isam) = affinity(idim, isam) + sample(idim, isam)
203 190 : affinity(idim + 1 : ndim, isam) = affinity(idim + 1 : ndim, isam) + sample(idim, isam) * tform(idim + 1 : ndim, idim)
204 : end do
205 : #else
206 : #error "Unrecognized interface."
207 : #endif
208 : end do
209 : else
210 : ! \todo:
211 : ! The performance of this block can be improved by looping over `nsam` in the innermost layer.
212 : do concurrent(isam = 1 : nsam)
213 : #if CGR_ENABLED
214 952 : affinity(isam, 1 : ndim) = matmul(tform, sample(isam, 1 : ndim)) PLUS(tlate)
215 : #elif CUD_ENABLED
216 280 : affinity(isam, 1 : ndim) = sample(isam, ndim) * tform(1 : ndim, ndim) PLUS(tlate)
217 300 : do idim = ndim - 1, 1, -1
218 326 : affinity(isam, 1 : idim) = affinity(isam, 1 : idim) + sample(isam, idim) * tform(1 : idim, idim)
219 : end do
220 : #elif CUU_ENABLED
221 190 : affinity(isam, 1 : ndim - 1) = sample(isam, ndim) * tform(1 : ndim - 1, ndim) PLUS(tlate(1 : ndim - 1))
222 90 : affinity(isam, ndim) = sample(isam, ndim) PLUS(tlate(ndim))
223 300 : do idim = ndim - 1, 1, -1
224 136 : affinity(isam, 1 : idim - 1) = affinity(isam, 1 : idim - 1) + sample(isam, idim) * tform(1 : idim - 1, idim)
225 190 : affinity(isam, idim) = affinity(isam, idim) + sample(isam, idim) * tform(idim, idim)
226 : end do
227 : #elif CLD_ENABLED
228 280 : affinity(isam, 1 : ndim) = sample(isam, 1) * tform(1 : ndim, 1) PLUS(tlate)
229 300 : do idim = 2, ndim
230 326 : affinity(isam, idim : ndim) = affinity(isam, idim : ndim) + sample(isam, idim) * tform(idim : ndim, idim)
231 : end do
232 : #elif CLU_ENABLED
233 90 : affinity(isam, 1) = sample(isam, 1) PLUS(tlate(1))
234 190 : affinity(isam, 2 : ndim) = sample(isam, 1) * tform(2 : ndim, 1) PLUS(tlate(2 : ndim))
235 300 : do idim = 2, ndim
236 100 : affinity(isam, idim) = affinity(isam, idim) + sample(isam, idim)
237 226 : affinity(isam, idim + 1 : ndim) = affinity(isam, idim + 1 : ndim) + sample(isam, idim) * tform(idim + 1 : ndim, idim)
238 : end do
239 : #else
240 : #error "Unrecognized interface."
241 : #endif
242 : end do
243 : end if
244 :
245 : #else
246 : !%%%%%%%%%%%%%%%%%%%%%%%%
247 : #error "Unrecognized interface."
248 : !%%%%%%%%%%%%%%%%%%%%%%%%
249 : #endif
250 : #undef CHECK_LEN_TLATE
251 : #undef CHECK_VAL_DIM
252 : #undef PLUS
|