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 interface [pm_sampleMean](@ref pm_sampleMean).
19 : !>
20 : !> \author
21 : !> \FatemehBagheri, Saturday 4:40 PM, August 21, 2021, Dallas, TX
22 :
23 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
24 :
25 : ! Define weight rules.
26 : #if WNO_ENABLED
27 : #define INCREMENT(X,Y)
28 : #define GET_WEIGHTED(X,W)X
29 : #elif WTI_ENABLED || WTR_ENABLED
30 : #define GET_WEIGHTED(X,W)X * W
31 : #define INCREMENT(X,Y)X = X + Y
32 : #elif !(getMeanMerged_ENABLED || setMeanMerged_ENABLED)
33 : #error "Unrecognized interface."
34 : #endif
35 : ! Define weight type and kind.
36 : #if WTI_ENABLED
37 : #define TYPE_OF_WEIGHT integer(IK)
38 : #elif WTR_ENABLED
39 : #define TYPE_OF_WEIGHT real(TKC)
40 : #elif !(WNO_ENABLED || getMeanMerged_ENABLED || setMeanMerged_ENABLED)
41 : #error "Unrecognized interface."
42 : #endif
43 : ! Define the runtime checks.
44 : #define CHECK_LEN_MEAN(DIM) \
45 : CHECK_ASSERTION(__LINE__, size(sample, DIM, IK) == size(mean, 1, IK), \
46 : SK_"@setMean(): The condition `size(sample, dim, 1) == size(mean, 1)` must hold. dim, shape(sample), size(mean) = "\
47 : //getStr([DIM, shape(sample, IK), size(mean, 1, IK)]))
48 : #define CHECK_LEN_WEI(DIM) \
49 : CHECK_ASSERTION(__LINE__, size(sample, DIM, IK) == size(weight, 1, IK), \
50 : SK_"@setMean(): The condition `size(sample, dim, 1) == size(weight, 1)` must hold. dim, shape(sample), size(weight) = "\
51 : //getStr([DIM, shape(sample, IK), size(weight, 1, IK)]))
52 : #define CHECK_VAL_WEI \
53 : CHECK_ASSERTION(__LINE__, all(0 <= weight), \
54 : SK_"@setMean(): The condition `all(0. <= weight)` must hold. weight = "//getStr(weight))
55 : #define CHECK_SUM_WEI \
56 : CHECK_ASSERTION(__LINE__, 0 < sum(weight), \
57 : SK_"@setMean(): The condition `0 < sum(weight)` must hold. weight = "//getStr(weight))
58 : #define CHECK_SHAPE_SAMPLE \
59 : CHECK_ASSERTION(__LINE__, all(0_IK < shape(sample, IK)), \
60 : SK_"@setMean(): The condition `all(0 < shape(sample))` must hold. shape(sample) = "//getStr(shape(sample, IK)))
61 : #if setMean_ENABLED && DIM_ENABLED
62 : #define CHECK_VAL_DIM \
63 : CHECK_ASSERTION(__LINE__, 1 <= dim .and. dim <= rank(sample), \
64 : SK_"@setMean(): The condition `1 <= dim .and. dim <= rank(sample)` must hold. dim, rank(sample) = "//getStr([integer(IK) :: dim, rank(sample)]))
65 : #elif setMean_ENABLED && ALL_ENABLED
66 : #define CHECK_VAL_DIM
67 : #endif
68 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
69 : #if getMeanMerged_ENABLED && New_ENABLED
70 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
71 :
72 860 : call setMeanMerged(meanMerged, meanB, meanA, fracA)
73 :
74 : !%%%%%%%%%%%%%%%%%%%%
75 : #elif setMeanMerged_ENABLED
76 : !%%%%%%%%%%%%%%%%%%%%
77 :
78 : ! Define the output value for setMeanMerged.
79 : #if Old_ENABLED
80 : #define TARGET_AVG meanB
81 : #define CHECK_LEN_TARGET_MEAN
82 : #elif New_ENABLED
83 : #define CHECK_LEN_TARGET_MEAN \
84 : CHECK_ASSERTION(__LINE__, size(meanMerged, 1, IK) == size(meanB, 1, IK), SK_"@setMeanMerged(): The condition `size(meanMerged) == size(meanB)` must hold. size(meanMerged), size(meanB) = "//getStr([size(meanMerged, 1, IK), size(meanB, 1, IK)]))
85 : #define TARGET_AVG meanMerged
86 : #else
87 : #error "Unrecognized interface."
88 : #endif
89 : #if D1_ENABLED
90 : real(TKC) :: fracB
91 : integer(IK) :: idim
92 : #endif
93 2580 : CHECK_ASSERTION(__LINE__, 0._TKC < fracA .and. fracA < 1._TKC, SK_"@setMeanMerged(): The condition `0 < fracA .and. fracA < 1` must hold. fracA = "//getStr(fracA))
94 : #if D0_ENABLED
95 1290 : TARGET_AVG = fracA * meanA + (1._TKC - fracA) * meanB
96 : #elif D1_ENABLED
97 2580 : CHECK_LEN_TARGET_MEAN
98 3870 : CHECK_ASSERTION(__LINE__, size(meanA, 1, IK) == size(meanB, 1, IK), SK_"@setMeanMerged(): The condition `size(meanA) == size(meanB)` must hold. size(meanA), size(meanB) = "//getStr([size(meanA, 1, IK), size(meanB, 1, IK)]))
99 1290 : fracB = 1._TKC - fracA
100 : do concurrent(idim = 1 : size(meanB, 1, IK))
101 5118 : TARGET_AVG(idim) = fracA * meanA(idim) + fracB * meanB(idim)
102 : end do
103 : #else
104 : #error "Unrecognized interface."
105 : #endif
106 : #undef TARGET_AVG
107 :
108 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%
109 : #elif getMean_ENABLED && XY_ENABLED
110 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%
111 :
112 : #if ALL_ENABLED && WNO_ENABLED
113 : call setMean(mean, x, y)
114 : #elif ALL_ENABLED && (WTI_ENABLED || WTR_ENABLED)
115 : TYPE_OF_WEIGHT :: weisum
116 : call setMean(mean, x, y, weight, weisum)
117 : #else
118 : #error "Unrecognized interface."
119 : #endif
120 :
121 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
122 : #elif getMean_ENABLED && (D1_ENABLED || D2_ENABLED)
123 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
124 :
125 : #if ALL_ENABLED && WNO_ENABLED
126 40116 : call setMean(mean, sample)
127 : #elif DIM_ENABLED && WNO_ENABLED
128 30412 : call setMean(mean, sample, dim)
129 : #elif ALL_ENABLED && (WTI_ENABLED || WTR_ENABLED)
130 : TYPE_OF_WEIGHT :: weisum
131 12057 : call setMean(mean, sample, weight, weisum)
132 : #elif DIM_ENABLED && (WTI_ENABLED || WTR_ENABLED)
133 : TYPE_OF_WEIGHT :: weisum
134 20201 : call setMean(mean, sample, dim, weight, weisum)
135 : #else
136 : #error "Unrecognized interface."
137 : #endif
138 :
139 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%
140 : #elif setMean_ENABLED && XY_ENABLED
141 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%
142 :
143 : integer(IK) :: isam
144 : #if WNO_ENABLED
145 : real(TKC) :: weisum
146 3207 : weisum = size(x, 1, IK)
147 : #elif WTI_ENABLED || WTR_ENABLED
148 19224 : CHECK_ASSERTION(__LINE__, size(x, 1, IK) == size(weight, 1, IK), SK_"@setMean(): The condition `size(x) == size(weight)` must hold. size(x), size(weight) = "//getStr([size(x, 1, IK), size(weight, 1, IK)]))
149 132362 : CHECK_SUM_WEI
150 132362 : CHECK_VAL_WEI
151 6408 : weisum = 0
152 : #else
153 : #error "Unrecognized interface."
154 : #endif
155 9615 : mean = 0._TKC
156 28845 : CHECK_ASSERTION(__LINE__, size(x, 1, IK) == size(y, 1, IK), SK_"@setMean(): The condition `size(x) == size(y)` must hold. size(x), size(y) = "//getStr([size(x, 1, IK), size(y, 1, IK)]))
157 198576 : do isam = 1_IK, size(x, 1, IK)
158 125954 : INCREMENT(weisum,weight(isam))
159 188961 : mean(1) = mean(1) + GET_WEIGHTED(x(isam),weight(isam))
160 198576 : mean(2) = mean(2) + GET_WEIGHTED(y(isam),weight(isam))
161 : end do
162 28845 : mean = mean / weisum
163 :
164 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%
165 : #elif setMean_ENABLED && D1_ENABLED
166 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%
167 :
168 : integer(IK) :: isam
169 : #if WNO_ENABLED
170 : real(TKC) :: weisum
171 41320 : weisum = size(sample, 1, IK)
172 : #elif WTI_ENABLED || WTR_ENABLED
173 74310 : CHECK_LEN_WEI(1_IK)
174 74975 : CHECK_SUM_WEI
175 74975 : CHECK_VAL_WEI
176 14862 : weisum = 0
177 : #else
178 : #error "Unrecognized interface."
179 : #endif
180 56182 : mean = 0._TKC
181 30069 : CHECK_VAL_DIM
182 112364 : CHECK_SHAPE_SAMPLE
183 3939828 : do isam = 1_IK, size(sample)
184 60113 : INCREMENT(weisum,weight(isam))
185 3939828 : mean = mean + GET_WEIGHTED(sample(isam),weight(isam))
186 : end do
187 56182 : mean = mean / weisum
188 :
189 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
190 : #elif setMean_ENABLED && D2_ENABLED && ALL_ENABLED
191 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
192 :
193 : integer(IK) :: idim, jdim
194 : integer(IK) :: mdim, ndim
195 : #if WNO_ENABLED
196 : real(TKC) :: weisum
197 3614 : weisum = size(sample, kind = IK)
198 : #elif WTI_ENABLED || WTR_ENABLED
199 : integer(IK) :: iwei
200 50463 : CHECK_ASSERTION(__LINE__, size(sample, kind = IK) == size(weight, 1, IK), SK_"@setMean(): The condition `size(sample) == size(weight)` must hold. shape(sample), size(weight) = "//getStr([shape(sample, IK), shape(weight, IK)]))
201 93208 : CHECK_SUM_WEI
202 93208 : CHECK_VAL_WEI
203 7209 : weisum = 0
204 : iwei = 0
205 : #else
206 : #error "Unrecognized interface."
207 : #endif
208 10823 : mdim = size(sample, 1, IK)
209 10823 : ndim = size(sample, 2, IK)
210 10823 : mean = 0._TKC
211 48773 : do jdim = 1, ndim
212 177940 : do idim = 1, mdim
213 85999 : INCREMENT(iwei,1_IK)
214 85999 : INCREMENT(weisum,weight(iwei))
215 167117 : mean = mean + GET_WEIGHTED(sample(idim, jdim),weight(iwei))
216 : end do
217 : end do
218 10823 : mean = mean / weisum
219 :
220 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
221 : #elif setMean_ENABLED && D2_ENABLED && WNO_ENABLED && DIM_ENABLED
222 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
223 :
224 : real(TKC) :: normfac
225 : integer(IK) :: ndim, nsam, idim, isam
226 97551 : CHECK_VAL_DIM
227 97551 : CHECK_SHAPE_SAMPLE
228 227619 : CHECK_LEN_MEAN(3 - dim)
229 32517 : nsam = size(sample, dim, IK)
230 32517 : ndim = size(sample, 3 - dim, IK)
231 32517 : normfac = 1._TKC / nsam
232 32517 : if (dim == 2_IK) then
233 : ! attributes are along the columns.
234 86296 : mean = 0._TKC
235 170468 : do isam = 1, nsam
236 530996 : mean = mean + sample(1 : ndim, isam)
237 : end do
238 86296 : mean = mean * normfac
239 : else
240 : ! attributes are along the rows.
241 : do concurrent(idim = 1 : ndim)
242 208997 : mean(idim) = sum(sample(1 : nsam, idim)) * normfac
243 : end do
244 : end if
245 :
246 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
247 : #elif setMean_ENABLED && D2_ENABLED && DIM_ENABLED
248 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
249 :
250 : real(TKC) :: normfac
251 : integer(IK) :: ndim, nsam, idim, isam
252 74873 : ndim = size(sample, 3 - dim, IK)
253 74873 : nsam = size(sample, dim, IK)
254 524111 : CHECK_LEN_MEAN(3 - dim)
255 224619 : CHECK_SHAPE_SAMPLE
256 224619 : CHECK_VAL_DIM
257 : #if WNO_ENABLED
258 : normfac = 1._TKC / nsam
259 : #elif WTI_ENABLED || WTR_ENABLED
260 524111 : CHECK_LEN_WEI(dim)
261 1020890 : CHECK_SUM_WEI
262 1020890 : CHECK_VAL_WEI
263 74873 : weisum = 0
264 : #else
265 : #error "Unrecognized interface."
266 : #endif
267 74873 : if (dim == 2_IK) then ! attributes are along the columns.
268 254048 : mean = 0._TKC
269 881931 : do isam = 1, nsam
270 820299 : INCREMENT(weisum,weight(isam))
271 3780261 : mean = mean + GET_WEIGHTED(sample(1 : ndim, isam),weight(isam))
272 : end do
273 : #if WTI_ENABLED || WTR_ENABLED
274 58106 : normfac = 1._TKC / weisum
275 : #endif
276 254048 : mean = mean * normfac
277 : else ! attributes are along the rows.
278 : #if WNO_ENABLED
279 : do concurrent(idim = 1 : ndim)
280 : mean(idim) = sum(sample(1 : nsam, idim)) * normfac
281 : end do
282 : #elif WTI_ENABLED || WTR_ENABLED
283 13241 : mean(1) = 0._TKC
284 138959 : do isam = 1, nsam
285 125718 : INCREMENT(weisum,weight(isam))
286 138959 : mean(1) = mean(1) + GET_WEIGHTED(sample(isam, 1),weight(isam))
287 : end do
288 13241 : normfac = 1._TKC / weisum
289 13241 : mean(1) = mean(1) * normfac
290 : do concurrent(idim = 2 : ndim)
291 : ! warning: dot_product() complex-conjugates its first argument.
292 275289 : mean(idim) = dot_product(weight, sample(1 : nsam, idim)) * normfac
293 : end do
294 : end if
295 : #endif
296 :
297 : #else
298 : !%%%%%%%%%%%%%%%%%%%%%%%%
299 : #error "Unrecognized interface."
300 : !%%%%%%%%%%%%%%%%%%%%%%%%
301 : #endif
302 : #undef CHECK_LEN_TARGET_MEAN
303 : #undef CHECK_SHAPE_SAMPLE
304 : #undef CHECK_LEN_MEAN
305 : #undef TYPE_OF_WEIGHT
306 : #undef CHECK_LEN_WEI
307 : #undef CHECK_SUM_WEI
308 : #undef CHECK_VAL_WEI
309 : #undef CHECK_VAL_DIM
310 : #undef GET_WEIGHTED
311 : #undef INCREMENT
|