ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_sampleMean::getMeanMerged Interface Reference

Generate and return the (weighted) merged mean of a sample resulting from the merger of two separate (weighted) samples \(A\) and \(B\).
More...

Detailed Description

Generate and return the (weighted) merged mean of a sample resulting from the merger of two separate (weighted) samples \(A\) and \(B\).

See the documentation of pm_sampelMean for more information and definition online updating of sample mean.

Parameters
[in]meanB: The input object of the same type and kind and rank and shape as the input argument meanA, containing the mean of the second sample that must be merged with the first sample.
[in]meanA: The input scalar or contiguous vector of shape (1:ndim) of,
  1. type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128),
  2. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the mean of the first sample.
[in]fracA: The input scalar of type real of the same kind as kind of meanA, containing the ratio of the sum of the weights of all points in sample \(A\) to sum of weights of all points in the merged sample.
If the sample is unweighted, then fracA is simply size(sampleA) / (size(sampleA) + size(sampleB)).
Returns
meanMerged : The output object of the same type and kind and rank and shape as meanA, containing the mean of the sample resulting form the merger of the two samples.


Possible calling interfaces

! univariate sample.
meanMerged = getMeanMerged(meanB, meanA, fracA)
! ndim-dimensional sample.
meanMerged(1:ndim) = getMeanMerged(meanB(1:ndim), meanA(1:ndim), fracA)
Generate and return the (weighted) merged mean of a sample resulting from the merger of two separate ...
This module contains classes and procedures for computing the first moment (i.e., the statistical mea...
Warning
The condition 0 < fracA .and. fracA < 1 must hold for the corresponding input arguments.
The condition size(meanB) == size(meanA) must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
See also
getCor
setCor
getCov
setCov
getVar
setVar
getMean
setMean
getCovMerged
setCovMerged
getVarMerged
setVarMerged
getMeanMerged
setMeanMerged


Example usage

1program example
2
3 use pm_kind, only: SK, IK
4 use pm_kind, only: TKG => RKS ! All other real types are also supported.
5 use pm_sampleMean, only: getMean
9 use pm_distUnif, only: getUnifRand
10 use pm_arrayRange, only: getRange
11 use pm_io, only: display_type
12
13 implicit none
14
15 integer(IK) , parameter :: nsam = 2
16 integer(IK) , allocatable :: iweight(:)
17 integer(IK) :: isam, ndim, lb(nsam), ub(nsam)
18 integer(IK) :: dim, itry, ntry = 10
19 type(display_type) :: disp
20 disp = display_type(file = "main.out.F90")
21
22 call disp%skip()
23 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
24 call disp%show("!Compute the merged mean of a univariate sample.")
25 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%skip()
27
28 block
29 real(TKG) :: mean(0:nsam), meanMerged
30 real(TKG), allocatable :: sample(:)
31 do itry = 1, ntry
32 call disp%skip()
33 call disp%show("lb(1) = 1; ub(1) = getUnifRand(1, 7)")
34 lb(1) = 1; ub(1) = getUnifRand(1, 7)
35 call disp%show("do isam = 2, nsam")
36 call disp%show(" lb(isam) = ub(isam - 1) + 1")
37 call disp%show(" ub(isam) = ub(isam - 1) + getUnifRand(1, 7)")
38 call disp%show("end do")
39 do isam = 2, nsam
40 lb(isam) = ub(isam - 1) + 1
41 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
42 end do
43 call disp%show("lb")
44 call disp%show( lb )
45 call disp%show("ub")
46 call disp%show( ub )
47 call disp%show("sample = getUnifRand(0., 1., ub(nsam))")
48 sample = getUnifRand(0., 1., ub(nsam))
49 call disp%show("sample")
50 call disp%show( sample )
51 call disp%show("mean(0) = getMean(sample)")
52 mean(0) = getMean(sample)
53 call disp%show("mean(0) ! reference")
54 call disp%show( mean(0) )
55 call disp%show("do isam = 1, nsam")
56 call disp%show(" mean(isam) = getMean(sample(lb(isam):ub(isam)))")
57 call disp%show("end do")
58 do isam = 1, nsam
59 mean(isam) = getMean(sample(lb(isam):ub(isam)))
60 end do
61 call disp%show("call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))")
62 call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
63 call disp%show("meanMerged")
64 call disp%show( meanMerged )
65 call disp%show("call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))")
66 call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
67 call disp%show("mean(2)")
68 call disp%show( mean(2) )
69 call disp%show("mean(0) ! reference")
70 call disp%show( mean(0) )
71 call disp%skip()
72 end do
73 end block
74
75 call disp%skip()
76 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
77 call disp%show("!Compute the merged mean of a frequency weighted univariate sample.")
78 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
79 call disp%skip()
80
81 block
82 real(TKG) :: mean(0:nsam), meanMerged
83 real(TKG), allocatable :: sample(:)
84 do itry = 1, ntry
85 call disp%skip()
86 call disp%show("lb(1) = 1; ub(1) = getUnifRand(1, 7)")
87 lb(1) = 1; ub(1) = getUnifRand(1, 7)
88 call disp%show("do isam = 2, nsam")
89 call disp%show(" lb(isam) = ub(isam - 1) + 1")
90 call disp%show(" ub(isam) = ub(isam - 1) + getUnifRand(1, 7)")
91 call disp%show("end do")
92 do isam = 2, nsam
93 lb(isam) = ub(isam - 1) + 1
94 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
95 end do
96 call disp%show("lb")
97 call disp%show( lb )
98 call disp%show("ub")
99 call disp%show( ub )
100 call disp%show("sample = getUnifRand(0., 1., ub(nsam))")
101 sample = getUnifRand(0., 1., ub(nsam))
102 call disp%show("sample")
103 call disp%show( sample )
104 call disp%show("iweight = getUnifRand(1, 10, size(sample, 1, IK))")
105 iweight = getUnifRand(1, 10, size(sample, 1, IK))
106 call disp%show("iweight")
107 call disp%show( iweight )
108 call disp%show("mean(0) = getMean(sample, iweight)")
109 mean(0) = getMean(sample, iweight)
110 call disp%show("mean(0) ! reference")
111 call disp%show( mean(0) )
112 call disp%show("do isam = 1, nsam")
113 call disp%show(" mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))")
114 call disp%show("end do")
115 do isam = 1, nsam
116 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
117 end do
118 call disp%show("call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))")
119 call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
120 call disp%show("meanMerged")
121 call disp%show( meanMerged )
122 call disp%show("call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))")
123 call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
124 call disp%show("mean(2)")
125 call disp%show( mean(2) )
126 call disp%show("mean(0) ! reference")
127 call disp%show( mean(0) )
128 call disp%skip()
129 end do
130 end block
131
132 call disp%skip()
133 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
134 call disp%show("!Compute the merged mean of a reliability weighted univariate sample.")
135 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
136 call disp%skip()
137
138 block
139 real(TKG) :: mean(0:nsam), meanMerged
140 real(TKG), allocatable :: sample(:)
141 real(TKG), allocatable :: rweight(:)
142 do itry = 1, ntry
143 call disp%skip()
144 call disp%show("lb(1) = 1; ub(1) = getUnifRand(1, 7)")
145 lb(1) = 1; ub(1) = getUnifRand(1, 7)
146 call disp%show("do isam = 2, nsam")
147 call disp%show(" lb(isam) = ub(isam - 1) + 1")
148 call disp%show(" ub(isam) = ub(isam - 1) + getUnifRand(1, 7)")
149 call disp%show("end do")
150 do isam = 2, nsam
151 lb(isam) = ub(isam - 1) + 1
152 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
153 end do
154 call disp%show("lb")
155 call disp%show( lb )
156 call disp%show("ub")
157 call disp%show( ub )
158 call disp%show("sample = getUnifRand(0., 1., ub(nsam))")
159 sample = getUnifRand(0., 1., ub(nsam))
160 call disp%show("sample")
161 call disp%show( sample )
162 call disp%show("rweight = getUnifRand(1., 2., size(sample, 1, IK))")
163 rweight = getUnifRand(1., 2., size(sample, 1, IK))
164 call disp%show("rweight")
165 call disp%show( rweight )
166 call disp%show("mean(0) = getMean(sample, rweight)")
167 mean(0) = getMean(sample, rweight)
168 call disp%show("mean(0) ! reference")
169 call disp%show( mean(0) )
170 call disp%show("do isam = 1, nsam")
171 call disp%show(" mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))")
172 call disp%show("end do")
173 do isam = 1, nsam
174 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
175 end do
176 call disp%show("call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))")
177 call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
178 call disp%show("meanMerged")
179 call disp%show( meanMerged )
180 call disp%show("call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))")
181 call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
182 call disp%show("mean(2)")
183 call disp%show( mean(2) )
184 call disp%show("mean(0) ! reference")
185 call disp%show( mean(0) )
186 call disp%skip()
187 end do
188 end block
189
190 call disp%skip()
191 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
192 call disp%show("!Compute the merged mean of a multivariate sample.")
193 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
194 call disp%skip()
195
196 block
197 real(TKG), allocatable :: mean(:,:), meanMerged(:)
198 real(TKG), allocatable :: sample(:,:)
199 do itry = 1, ntry
200 call disp%skip()
201 call disp%show("dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)")
202 dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
203 call disp%show("do isam = 2, nsam")
204 call disp%show(" lb(isam) = ub(isam - 1) + 1")
205 call disp%show(" ub(isam) = ub(isam - 1) + getUnifRand(1, 7)")
206 call disp%show("end do")
207 do isam = 2, nsam
208 lb(isam) = ub(isam - 1) + 1
209 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
210 end do
211 call disp%show("lb")
212 call disp%show( lb )
213 call disp%show("ub")
214 call disp%show( ub )
215 call disp%show("ndim = getUnifRand(1, minval(ub - lb + 1, 1))")
216 ndim = getUnifRand(1, minval(ub - lb + 1, 1))
217 call disp%show("call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])")
218 call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
219 call disp%show("call setResized(meanMerged, ndim)")
220 call setResized(meanMerged, ndim)
221 call disp%show("sample = getUnifRand(-1., +1., ndim, ub(nsam))")
222 sample = getUnifRand(-1., +1., ndim, ub(nsam))
223 call disp%show("sample")
224 call disp%show( sample )
225 call disp%show("mean(:,0) = getMean(sample, dim)")
226 mean(:,0) = getMean(sample, dim)
227 call disp%show("mean(:,0) ! reference")
228 call disp%show( mean(:,0) )
229 call disp%show("do isam = 1, nsam")
230 call disp%show(" mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)")
231 call disp%show("end do")
232 do isam = 1, nsam
233 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
234 end do
235 call disp%show("call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))")
236 call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
237 call disp%show("meanMerged")
238 call disp%show( meanMerged )
239 call disp%show("call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))")
240 call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
241 call disp%show("mean(:,2)")
242 call disp%show( mean(:,2) )
243 call disp%show("mean(:,0) ! reference")
244 call disp%show( mean(:,0) )
245 call disp%skip()
246 end do
247 end block
248
249 call disp%skip()
250 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
251 call disp%show("!Compute the merged mean of a frequency weighted multivariate sample.")
252 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
253 call disp%skip()
254
255 block
256 real(TKG), allocatable :: mean(:,:), meanMerged(:)
257 real(TKG), allocatable :: sample(:,:)
258 do itry = 1, ntry
259 call disp%skip()
260 call disp%show("dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)")
261 dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
262 call disp%show("do isam = 2, nsam")
263 call disp%show(" lb(isam) = ub(isam - 1) + 1")
264 call disp%show(" ub(isam) = ub(isam - 1) + getUnifRand(1, 7)")
265 call disp%show("end do")
266 do isam = 2, nsam
267 lb(isam) = ub(isam - 1) + 1
268 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
269 end do
270 call disp%show("lb")
271 call disp%show( lb )
272 call disp%show("ub")
273 call disp%show( ub )
274 call disp%show("ndim = getUnifRand(1, minval(ub - lb + 1, 1))")
275 ndim = getUnifRand(1, minval(ub - lb + 1, 1))
276 call disp%show("call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])")
277 call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
278 call disp%show("call setResized(meanMerged, ndim)")
279 call setResized(meanMerged, ndim)
280 call disp%show("sample = getUnifRand(-1., +1., ndim, ub(nsam))")
281 sample = getUnifRand(-1., +1., ndim, ub(nsam))
282 call disp%show("sample")
283 call disp%show( sample )
284 call disp%show("iweight = getUnifRand(1, 10, size(sample, dim, IK))")
285 iweight = getUnifRand(1, 10, size(sample, dim, IK))
286 call disp%show("iweight")
287 call disp%show( iweight )
288 call disp%show("mean(:,0) = getMean(sample, 2_IK, iweight)")
289 mean(:,0) = getMean(sample, 2_IK, iweight)
290 call disp%show("mean(:,0) ! reference")
291 call disp%show( mean(:,0) )
292 call disp%show("do isam = 1, nsam")
293 call disp%show(" mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))")
294 call disp%show("end do")
295 do isam = 1, nsam
296 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
297 end do
298 call disp%show("call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))")
299 call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
300 call disp%show("meanMerged")
301 call disp%show( meanMerged )
302 call disp%show("call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))")
303 call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
304 call disp%show("mean(:,2)")
305 call disp%show( mean(:,2) )
306 call disp%show("mean(:,0) ! reference")
307 call disp%show( mean(:,0) )
308 call disp%skip()
309 end do
310 end block
311
312 call disp%skip()
313 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
314 call disp%show("!Compute the merged mean of a reliability weighted multivariate sample.")
315 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
316 call disp%skip()
317
318 block
319 real(TKG), allocatable :: mean(:,:), meanMerged(:)
320 real(TKG), allocatable :: sample(:,:)
321 real(TKG), allocatable :: rweight(:)
322 do itry = 1, ntry
323 call disp%skip()
324 call disp%show("lb(1) = 1; ub(1) = getUnifRand(1, 7)")
325 lb(1) = 1; ub(1) = getUnifRand(1, 7)
326 call disp%show("do isam = 2, nsam")
327 call disp%show(" lb(isam) = ub(isam - 1) + 1")
328 call disp%show(" ub(isam) = ub(isam - 1) + getUnifRand(1, 7)")
329 call disp%show("end do")
330 do isam = 2, nsam
331 lb(isam) = ub(isam - 1) + 1
332 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
333 end do
334 call disp%show("lb")
335 call disp%show( lb )
336 call disp%show("ub")
337 call disp%show( ub )
338 call disp%show("ndim = getUnifRand(1, minval(ub - lb + 1, 1))")
339 ndim = getUnifRand(1, minval(ub - lb + 1, 1))
340 call disp%show("call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])")
341 call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
342 call disp%show("call setResized(meanMerged, ndim)")
343 call setResized(meanMerged, ndim)
344 call disp%show("sample = getUnifRand(-1., +1., ndim, ub(nsam))")
345 sample = getUnifRand(-1., +1., ndim, ub(nsam))
346 call disp%show("sample")
347 call disp%show( sample )
348 call disp%show("rweight = getUnifRand(1, 10, size(sample, dim, IK))")
349 rweight = getUnifRand(1, 10, size(sample, dim, IK))
350 call disp%show("rweight")
351 call disp%show( rweight )
352 call disp%show("mean(:,0) = getMean(sample, dim, rweight)")
353 mean(:,0) = getMean(sample, dim, rweight)
354 call disp%show("mean(:,0) ! reference")
355 call disp%show( mean(:,0) )
356 call disp%show("do isam = 1, nsam")
357 call disp%show(" mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))")
358 call disp%show("end do")
359 do isam = 1, nsam
360 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
361 end do
362 call disp%show("call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))")
363 call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
364 call disp%show("meanMerged")
365 call disp%show( meanMerged )
366 call disp%show("call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))")
367 call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
368 call disp%show("mean(:,2)")
369 call disp%show( mean(:,2) )
370 call disp%show("mean(:,0) ! reference")
371 call disp%show( mean(:,0) )
372 call disp%skip()
373 end do
374 end block
375
376end program example
Generate minimally-spaced character, integer, real sequences or sequences at fixed intervals of size ...
Resize (shrink or expand) an input allocatable array of rank 1..3 to arbitrary new lower and upper bo...
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
Generate and return the (weighted) mean of an input sample of nsam observations with ndim = 1 or 2 at...
Return the (weighted) merged mean of a sample resulting from the merger of two separate (weighted) sa...
This module contains procedures and generic interfaces for generating ranges of discrete character,...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains classes and procedures for computing various statistical quantities related to t...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
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
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

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

Example output
1
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3!Compute the merged mean of a univariate sample.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7lb(1) = 1; ub(1) = getUnifRand(1, 7)
8do isam = 2, nsam
9 lb(isam) = ub(isam - 1) + 1
10 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
11end do
12lb
13+1, +7
14ub
15+6, +13
16sample = getUnifRand(0., 1., ub(nsam))
17sample
18+0.450917125, +0.647481143, +0.579088211, +0.919775665, +0.116213977, +0.626923084, +0.252830148, +0.812254667, +0.764650822, +0.402873576, +0.611293614, +0.274069428, +0.698515773
19mean(0) = getMean(sample)
20mean(0) ! reference
21+0.550529838
22do isam = 1, nsam
23 mean(isam) = getMean(sample(lb(isam):ub(isam)))
24end do
25call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
26meanMerged
27+0.550529838
28call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
29mean(2)
30+0.550529838
31mean(0) ! reference
32+0.550529838
33
34
35lb(1) = 1; ub(1) = getUnifRand(1, 7)
36do isam = 2, nsam
37 lb(isam) = ub(isam - 1) + 1
38 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
39end do
40lb
41+1, +6
42ub
43+5, +7
44sample = getUnifRand(0., 1., ub(nsam))
45sample
46+0.597962201, +0.783627629, +0.778012991, +0.839513838, +0.665288568E-1, +0.574622750E-1, +0.287638009
47mean(0) = getMean(sample)
48mean(0) ! reference
49+0.487249374
50do isam = 1, nsam
51 mean(isam) = getMean(sample(lb(isam):ub(isam)))
52end do
53call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
54meanMerged
55+0.487249404
56call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
57mean(2)
58+0.487249404
59mean(0) ! reference
60+0.487249374
61
62
63lb(1) = 1; ub(1) = getUnifRand(1, 7)
64do isam = 2, nsam
65 lb(isam) = ub(isam - 1) + 1
66 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
67end do
68lb
69+1, +3
70ub
71+2, +8
72sample = getUnifRand(0., 1., ub(nsam))
73sample
74+0.522541583, +0.871412575, +0.480284035, +0.971569538, +0.653817952, +0.366228819E-1, +0.211793125, +0.485796511
75mean(0) = getMean(sample)
76mean(0) ! reference
77+0.529229820
78do isam = 1, nsam
79 mean(isam) = getMean(sample(lb(isam):ub(isam)))
80end do
81call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
82meanMerged
83+0.529229760
84call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
85mean(2)
86+0.529229760
87mean(0) ! reference
88+0.529229820
89
90
91lb(1) = 1; ub(1) = getUnifRand(1, 7)
92do isam = 2, nsam
93 lb(isam) = ub(isam - 1) + 1
94 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
95end do
96lb
97+1, +7
98ub
99+6, +10
100sample = getUnifRand(0., 1., ub(nsam))
101sample
102+0.701328456, +0.717610896, +0.981206298, +0.235992670, +0.513632417, +0.684459209, +0.128934205, +0.514072537, +0.752650261, +0.330836773
103mean(0) = getMean(sample)
104mean(0) ! reference
105+0.556072354
106do isam = 1, nsam
107 mean(isam) = getMean(sample(lb(isam):ub(isam)))
108end do
109call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
110meanMerged
111+0.556072354
112call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
113mean(2)
114+0.556072354
115mean(0) ! reference
116+0.556072354
117
118
119lb(1) = 1; ub(1) = getUnifRand(1, 7)
120do isam = 2, nsam
121 lb(isam) = ub(isam - 1) + 1
122 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
123end do
124lb
125+1, +4
126ub
127+3, +8
128sample = getUnifRand(0., 1., ub(nsam))
129sample
130+0.641783416, +0.943830848, +0.969381273, +0.338236094E-1, +0.402671814, +0.632590532, +0.337394714, +0.189478397E-1
131mean(0) = getMean(sample)
132mean(0) ! reference
133+0.497552991
134do isam = 1, nsam
135 mean(isam) = getMean(sample(lb(isam):ub(isam)))
136end do
137call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
138meanMerged
139+0.497552991
140call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
141mean(2)
142+0.497552991
143mean(0) ! reference
144+0.497552991
145
146
147lb(1) = 1; ub(1) = getUnifRand(1, 7)
148do isam = 2, nsam
149 lb(isam) = ub(isam - 1) + 1
150 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
151end do
152lb
153+1, +7
154ub
155+6, +9
156sample = getUnifRand(0., 1., ub(nsam))
157sample
158+0.926578522, +0.492532432, +0.584347427, +0.594964683, +0.367180228, +0.699692965E-2, +0.717345178, +0.390479803, +0.529461861
159mean(0) = getMean(sample)
160mean(0) ! reference
161+0.512209654
162do isam = 1, nsam
163 mean(isam) = getMean(sample(lb(isam):ub(isam)))
164end do
165call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
166meanMerged
167+0.512209654
168call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
169mean(2)
170+0.512209654
171mean(0) ! reference
172+0.512209654
173
174
175lb(1) = 1; ub(1) = getUnifRand(1, 7)
176do isam = 2, nsam
177 lb(isam) = ub(isam - 1) + 1
178 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
179end do
180lb
181+1, +3
182ub
183+2, +3
184sample = getUnifRand(0., 1., ub(nsam))
185sample
186+0.338360667, +0.140468717, +0.724098504
187mean(0) = getMean(sample)
188mean(0) ! reference
189+0.400975943
190do isam = 1, nsam
191 mean(isam) = getMean(sample(lb(isam):ub(isam)))
192end do
193call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
194meanMerged
195+0.400975943
196call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
197mean(2)
198+0.400975943
199mean(0) ! reference
200+0.400975943
201
202
203lb(1) = 1; ub(1) = getUnifRand(1, 7)
204do isam = 2, nsam
205 lb(isam) = ub(isam - 1) + 1
206 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
207end do
208lb
209+1, +5
210ub
211+4, +11
212sample = getUnifRand(0., 1., ub(nsam))
213sample
214+0.310075819, +0.423070192E-1, +0.579325616, +0.569845974, +0.194519162, +0.624657273, +0.625143647E-1, +0.435734630, +0.460889161, +0.947877705, +0.317531228E-1
215mean(0) = getMean(sample)
216mean(0) ! reference
217+0.387227237
218do isam = 1, nsam
219 mean(isam) = getMean(sample(lb(isam):ub(isam)))
220end do
221call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
222meanMerged
223+0.387227237
224call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
225mean(2)
226+0.387227237
227mean(0) ! reference
228+0.387227237
229
230
231lb(1) = 1; ub(1) = getUnifRand(1, 7)
232do isam = 2, nsam
233 lb(isam) = ub(isam - 1) + 1
234 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
235end do
236lb
237+1, +8
238ub
239+7, +9
240sample = getUnifRand(0., 1., ub(nsam))
241sample
242+0.701725006, +0.641202450, +0.524549246, +0.442632616, +0.168878436, +0.173875213, +0.876675248, +0.943157434, +0.165064931E-1
243mean(0) = getMean(sample)
244mean(0) ! reference
245+0.498800218
246do isam = 1, nsam
247 mean(isam) = getMean(sample(lb(isam):ub(isam)))
248end do
249call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
250meanMerged
251+0.498800218
252call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
253mean(2)
254+0.498800218
255mean(0) ! reference
256+0.498800218
257
258
259lb(1) = 1; ub(1) = getUnifRand(1, 7)
260do isam = 2, nsam
261 lb(isam) = ub(isam - 1) + 1
262 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
263end do
264lb
265+1, +5
266ub
267+4, +5
268sample = getUnifRand(0., 1., ub(nsam))
269sample
270+0.526680708, +0.296144009, +0.329744816E-2, +0.148903370, +0.488615394
271mean(0) = getMean(sample)
272mean(0) ! reference
273+0.292728186
274do isam = 1, nsam
275 mean(isam) = getMean(sample(lb(isam):ub(isam)))
276end do
277call setMeanMerged(meanMerged, mean(2), mean(1), ub(1) / real(ub(2), TKG))
278meanMerged
279+0.292728186
280call setMeanMerged(mean(2), mean(1), ub(1) / real(ub(2), TKG))
281mean(2)
282+0.292728186
283mean(0) ! reference
284+0.292728186
285
286
287!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
288!Compute the merged mean of a frequency weighted univariate sample.
289!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
290
291
292lb(1) = 1; ub(1) = getUnifRand(1, 7)
293do isam = 2, nsam
294 lb(isam) = ub(isam - 1) + 1
295 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
296end do
297lb
298+1, +8
299ub
300+7, +9
301sample = getUnifRand(0., 1., ub(nsam))
302sample
303+0.762919188E-1, +0.999307156, +0.328115582, +0.229218841, +0.201227903, +0.557934225, +0.675776422, +0.897380888, +0.985947430
304iweight = getUnifRand(1, 10, size(sample, 1, IK))
305iweight
306+8, +1, +5, +3, +9, +3, +3, +8, +9
307mean(0) = getMean(sample, iweight)
308mean(0) ! reference
309+0.520461857
310do isam = 1, nsam
311 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
312end do
313call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
314meanMerged
315+0.520461917
316call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
317mean(2)
318+0.520461917
319mean(0) ! reference
320+0.520461857
321
322
323lb(1) = 1; ub(1) = getUnifRand(1, 7)
324do isam = 2, nsam
325 lb(isam) = ub(isam - 1) + 1
326 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
327end do
328lb
329+1, +7
330ub
331+6, +8
332sample = getUnifRand(0., 1., ub(nsam))
333sample
334+0.222939193, +0.285459101, +0.305185556, +0.301230967, +0.495533168, +0.634538293, +0.813965201, +0.473072708
335iweight = getUnifRand(1, 10, size(sample, 1, IK))
336iweight
337+1, +7, +2, +2, +10, +3, +8, +6
338mean(0) = getMean(sample, iweight)
339mean(0) ! reference
340+0.503669024
341do isam = 1, nsam
342 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
343end do
344call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
345meanMerged
346+0.503669024
347call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
348mean(2)
349+0.503669024
350mean(0) ! reference
351+0.503669024
352
353
354lb(1) = 1; ub(1) = getUnifRand(1, 7)
355do isam = 2, nsam
356 lb(isam) = ub(isam - 1) + 1
357 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
358end do
359lb
360+1, +6
361ub
362+5, +7
363sample = getUnifRand(0., 1., ub(nsam))
364sample
365+0.148988485, +0.875638485, +0.153747797E-1, +0.803034306, +0.470391035, +0.685037196, +0.842763603
366iweight = getUnifRand(1, 10, size(sample, 1, IK))
367iweight
368+9, +1, +2, +2, +6, +1, +9
369mean(0) = getMean(sample, iweight)
370mean(0) ! reference
371+0.498186946
372do isam = 1, nsam
373 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
374end do
375call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
376meanMerged
377+0.498186946
378call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
379mean(2)
380+0.498186946
381mean(0) ! reference
382+0.498186946
383
384
385lb(1) = 1; ub(1) = getUnifRand(1, 7)
386do isam = 2, nsam
387 lb(isam) = ub(isam - 1) + 1
388 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
389end do
390lb
391+1, +3
392ub
393+2, +7
394sample = getUnifRand(0., 1., ub(nsam))
395sample
396+0.826362610, +0.792909563, +0.241511822, +0.954819739, +0.388998926, +0.910269499, +0.960105479
397iweight = getUnifRand(1, 10, size(sample, 1, IK))
398iweight
399+5, +1, +1, +6, +2, +3, +7
400mean(0) = getMean(sample, iweight)
401mean(0) ! reference
402+0.844987929
403do isam = 1, nsam
404 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
405end do
406call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
407meanMerged
408+0.844987869
409call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
410mean(2)
411+0.844987869
412mean(0) ! reference
413+0.844987929
414
415
416lb(1) = 1; ub(1) = getUnifRand(1, 7)
417do isam = 2, nsam
418 lb(isam) = ub(isam - 1) + 1
419 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
420end do
421lb
422+1, +3
423ub
424+2, +3
425sample = getUnifRand(0., 1., ub(nsam))
426sample
427+0.280467212, +0.837615788, +0.817623377
428iweight = getUnifRand(1, 10, size(sample, 1, IK))
429iweight
430+10, +9, +1
431mean(0) = getMean(sample, iweight)
432mean(0) ! reference
433+0.558041871
434do isam = 1, nsam
435 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
436end do
437call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
438meanMerged
439+0.558041871
440call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
441mean(2)
442+0.558041871
443mean(0) ! reference
444+0.558041871
445
446
447lb(1) = 1; ub(1) = getUnifRand(1, 7)
448do isam = 2, nsam
449 lb(isam) = ub(isam - 1) + 1
450 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
451end do
452lb
453+1, +3
454ub
455+2, +7
456sample = getUnifRand(0., 1., ub(nsam))
457sample
458+0.234112203, +0.263211846, +0.483875990, +0.937600076, +0.293082535, +0.730139375, +0.170546353
459iweight = getUnifRand(1, 10, size(sample, 1, IK))
460iweight
461+2, +3, +4, +9, +10, +6, +7
462mean(0) = getMean(sample, iweight)
463mean(0) ! reference
464+0.491152465
465do isam = 1, nsam
466 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
467end do
468call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
469meanMerged
470+0.491152465
471call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
472mean(2)
473+0.491152465
474mean(0) ! reference
475+0.491152465
476
477
478lb(1) = 1; ub(1) = getUnifRand(1, 7)
479do isam = 2, nsam
480 lb(isam) = ub(isam - 1) + 1
481 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
482end do
483lb
484+1, +8
485ub
486+7, +14
487sample = getUnifRand(0., 1., ub(nsam))
488sample
489+0.372631550, +0.318311393, +0.158778369, +0.263324797, +0.575749218, +0.199925721, +0.569555163, +0.341664970, +0.218133628, +0.784770250E-1, +0.595189154, +0.561443925, +0.283682048, +0.743308663E-1
490iweight = getUnifRand(1, 10, size(sample, 1, IK))
491iweight
492+7, +9, +10, +7, +10, +3, +1, +6, +8, +5, +1, +5, +10, +6
493mean(0) = getMean(sample, iweight)
494mean(0) ! reference
495+0.303451896
496do isam = 1, nsam
497 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
498end do
499call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
500meanMerged
501+0.303451866
502call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
503mean(2)
504+0.303451866
505mean(0) ! reference
506+0.303451896
507
508
509lb(1) = 1; ub(1) = getUnifRand(1, 7)
510do isam = 2, nsam
511 lb(isam) = ub(isam - 1) + 1
512 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
513end do
514lb
515+1, +2
516ub
517+1, +8
518sample = getUnifRand(0., 1., ub(nsam))
519sample
520+0.964238048, +0.451485991, +0.515642166, +0.842446625, +0.720165610, +0.527779520, +0.763204515, +0.412703097
521iweight = getUnifRand(1, 10, size(sample, 1, IK))
522iweight
523+5, +1, +3, +1, +1, +10, +8, +3
524mean(0) = getMean(sample, iweight)
525mean(0) ! reference
526+0.656367421
527do isam = 1, nsam
528 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
529end do
530call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
531meanMerged
532+0.656367362
533call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
534mean(2)
535+0.656367362
536mean(0) ! reference
537+0.656367421
538
539
540lb(1) = 1; ub(1) = getUnifRand(1, 7)
541do isam = 2, nsam
542 lb(isam) = ub(isam - 1) + 1
543 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
544end do
545lb
546+1, +7
547ub
548+6, +10
549sample = getUnifRand(0., 1., ub(nsam))
550sample
551+0.818193436, +0.719758928, +0.379511297, +0.931423962, +0.775915027, +0.996842086, +0.627090812, +0.884779632, +0.324178994, +0.762720704
552iweight = getUnifRand(1, 10, size(sample, 1, IK))
553iweight
554+8, +10, +10, +1, +5, +8, +6, +10, +10, +4
555mean(0) = getMean(sample, iweight)
556mean(0) ! reference
557+0.683708429
558do isam = 1, nsam
559 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
560end do
561call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
562meanMerged
563+0.683708310
564call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
565mean(2)
566+0.683708310
567mean(0) ! reference
568+0.683708429
569
570
571lb(1) = 1; ub(1) = getUnifRand(1, 7)
572do isam = 2, nsam
573 lb(isam) = ub(isam - 1) + 1
574 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
575end do
576lb
577+1, +5
578ub
579+4, +11
580sample = getUnifRand(0., 1., ub(nsam))
581sample
582+0.152686119, +0.208874285, +0.506744385E-1, +0.617766142, +0.680997014, +0.365606904, +0.453252792E-1, +0.304934323, +0.695051908, +0.777131855, +0.759560645
583iweight = getUnifRand(1, 10, size(sample, 1, IK))
584iweight
585+3, +1, +2, +10, +8, +8, +9, +2, +10, +1, +4
586mean(0) = getMean(sample, iweight)
587mean(0) ! reference
588+0.467283875
589do isam = 1, nsam
590 mean(isam) = getMean(sample(lb(isam):ub(isam)), iweight(lb(isam):ub(isam)))
591end do
592call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
593meanMerged
594+0.467283845
595call setMeanMerged(mean(2), mean(1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
596mean(2)
597+0.467283845
598mean(0) ! reference
599+0.467283875
600
601
602!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
603!Compute the merged mean of a reliability weighted univariate sample.
604!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
605
606
607lb(1) = 1; ub(1) = getUnifRand(1, 7)
608do isam = 2, nsam
609 lb(isam) = ub(isam - 1) + 1
610 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
611end do
612lb
613+1, +2
614ub
615+1, +8
616sample = getUnifRand(0., 1., ub(nsam))
617sample
618+0.219523549, +0.530630350E-1, +0.591719151, +0.554355145, +0.855650306E-1, +0.390973747, +0.686549723, +0.876149595
619rweight = getUnifRand(1., 2., size(sample, 1, IK))
620rweight
621+1.98177838, +1.93096590, +1.66918647, +1.06585813, +1.89013302, +1.56982970, +1.27188182, +1.20056570
622mean(0) = getMean(sample, rweight)
623mean(0) ! reference
624+0.382874727
625do isam = 1, nsam
626 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
627end do
628call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
629meanMerged
630+0.382874727
631call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
632mean(2)
633+0.382874727
634mean(0) ! reference
635+0.382874727
636
637
638lb(1) = 1; ub(1) = getUnifRand(1, 7)
639do isam = 2, nsam
640 lb(isam) = ub(isam - 1) + 1
641 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
642end do
643lb
644+1, +6
645ub
646+5, +9
647sample = getUnifRand(0., 1., ub(nsam))
648sample
649+0.198383570, +0.804668903, +0.708835721E-1, +0.592108309, +0.514955163, +0.497140467, +0.950913489, +0.890209377, +0.194909513
650rweight = getUnifRand(1., 2., size(sample, 1, IK))
651rweight
652+1.54535902, +1.55188823, +1.93018126, +1.36929965, +1.75334907, +1.84446430, +1.91420722, +1.23851562, +1.71659327
653mean(0) = getMean(sample, rweight)
654mean(0) ! reference
655+0.509971201
656do isam = 1, nsam
657 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
658end do
659call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
660meanMerged
661+0.509971261
662call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
663mean(2)
664+0.509971261
665mean(0) ! reference
666+0.509971201
667
668
669lb(1) = 1; ub(1) = getUnifRand(1, 7)
670do isam = 2, nsam
671 lb(isam) = ub(isam - 1) + 1
672 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
673end do
674lb
675+1, +6
676ub
677+5, +6
678sample = getUnifRand(0., 1., ub(nsam))
679sample
680+0.336120367, +0.197943091, +0.902010739, +0.243265152, +0.763902426, +0.316251516
681rweight = getUnifRand(1., 2., size(sample, 1, IK))
682rweight
683+1.50913501, +1.32110560, +1.80416369, +1.94302225, +1.03593516, +1.43990326
684mean(0) = getMean(sample, rweight)
685mean(0) ! reference
686+0.454590082
687do isam = 1, nsam
688 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
689end do
690call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
691meanMerged
692+0.454590082
693call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
694mean(2)
695+0.454590082
696mean(0) ! reference
697+0.454590082
698
699
700lb(1) = 1; ub(1) = getUnifRand(1, 7)
701do isam = 2, nsam
702 lb(isam) = ub(isam - 1) + 1
703 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
704end do
705lb
706+1, +3
707ub
708+2, +7
709sample = getUnifRand(0., 1., ub(nsam))
710sample
711+0.202634156, +0.247774780, +0.121375084, +0.754495919, +0.339238048E-1, +0.470085084, +0.981331944
712rweight = getUnifRand(1., 2., size(sample, 1, IK))
713rweight
714+1.59052038, +1.59817326, +1.90501332, +1.15679741, +1.14449167, +1.57645488, +1.53720415
715mean(0) = getMean(sample, rweight)
716mean(0) ! reference
717+0.391172975
718do isam = 1, nsam
719 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
720end do
721call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
722meanMerged
723+0.391172945
724call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
725mean(2)
726+0.391172945
727mean(0) ! reference
728+0.391172975
729
730
731lb(1) = 1; ub(1) = getUnifRand(1, 7)
732do isam = 2, nsam
733 lb(isam) = ub(isam - 1) + 1
734 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
735end do
736lb
737+1, +4
738ub
739+3, +6
740sample = getUnifRand(0., 1., ub(nsam))
741sample
742+0.790228426, +0.782313526, +0.158574581E-1, +0.461433649, +0.472879410E-2, +0.663071275
743rweight = getUnifRand(1., 2., size(sample, 1, IK))
744rweight
745+1.52090764, +1.49448633, +1.40140331, +1.67496121, +1.09953547, +1.62354589
746mean(0) = getMean(sample, rweight)
747mean(0) ! reference
748+0.481897980
749do isam = 1, nsam
750 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
751end do
752call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
753meanMerged
754+0.481898010
755call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
756mean(2)
757+0.481898010
758mean(0) ! reference
759+0.481897980
760
761
762lb(1) = 1; ub(1) = getUnifRand(1, 7)
763do isam = 2, nsam
764 lb(isam) = ub(isam - 1) + 1
765 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
766end do
767lb
768+1, +5
769ub
770+4, +7
771sample = getUnifRand(0., 1., ub(nsam))
772sample
773+0.134169936, +0.822881997, +0.938389301, +0.721104383, +0.373824298, +0.668917716, +0.153163910
774rweight = getUnifRand(1., 2., size(sample, 1, IK))
775rweight
776+1.52584398, +1.71873546, +1.60745549, +1.77679586, +1.63114715, +1.35978794, +1.74717808
777mean(0) = getMean(sample, rweight)
778mean(0) ! reference
779+0.545060039
780do isam = 1, nsam
781 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
782end do
783call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
784meanMerged
785+0.545060039
786call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
787mean(2)
788+0.545060039
789mean(0) ! reference
790+0.545060039
791
792
793lb(1) = 1; ub(1) = getUnifRand(1, 7)
794do isam = 2, nsam
795 lb(isam) = ub(isam - 1) + 1
796 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
797end do
798lb
799+1, +5
800ub
801+4, +8
802sample = getUnifRand(0., 1., ub(nsam))
803sample
804+0.212155700, +0.303383946, +0.101725161, +0.863918900, +0.623905778, +0.472624063, +0.570023000, +0.147542298
805rweight = getUnifRand(1., 2., size(sample, 1, IK))
806rweight
807+1.75116062, +1.94414186, +1.92987323, +1.17084336, +1.24150038, +1.29195631, +1.89036167, +1.71409440
808mean(0) = getMean(sample, rweight)
809mean(0) ! reference
810+0.377674073
811do isam = 1, nsam
812 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
813end do
814call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
815meanMerged
816+0.377674103
817call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
818mean(2)
819+0.377674103
820mean(0) ! reference
821+0.377674073
822
823
824lb(1) = 1; ub(1) = getUnifRand(1, 7)
825do isam = 2, nsam
826 lb(isam) = ub(isam - 1) + 1
827 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
828end do
829lb
830+1, +4
831ub
832+3, +10
833sample = getUnifRand(0., 1., ub(nsam))
834sample
835+0.925814509, +0.818524837, +0.663796067, +0.374792159, +0.661942065, +0.123193860E-1, +0.794063985, +0.938345492, +0.216213584, +0.858668685
836rweight = getUnifRand(1., 2., size(sample, 1, IK))
837rweight
838+1.06552732, +1.54301429, +1.46615005, +1.15717053, +1.21852922, +1.94094324, +1.47867775, +1.58476901, +1.02172172, +1.46851611
839mean(0) = getMean(sample, rweight)
840mean(0) ! reference
841+0.618859708
842do isam = 1, nsam
843 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
844end do
845call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
846meanMerged
847+0.618859708
848call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
849mean(2)
850+0.618859708
851mean(0) ! reference
852+0.618859708
853
854
855lb(1) = 1; ub(1) = getUnifRand(1, 7)
856do isam = 2, nsam
857 lb(isam) = ub(isam - 1) + 1
858 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
859end do
860lb
861+1, +7
862ub
863+6, +8
864sample = getUnifRand(0., 1., ub(nsam))
865sample
866+0.302460909, +0.918179691, +0.773340881, +0.581627488E-1, +0.236778378, +0.314706564, +0.416233838, +0.892307341
867rweight = getUnifRand(1., 2., size(sample, 1, IK))
868rweight
869+1.12129414, +1.42488420, +1.04970360, +1.67097211, +1.34633660, +1.71282721, +1.61148119, +1.92837989
870mean(0) = getMean(sample, rweight)
871mean(0) ! reference
872+0.489276558
873do isam = 1, nsam
874 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
875end do
876call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
877meanMerged
878+0.489276558
879call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
880mean(2)
881+0.489276558
882mean(0) ! reference
883+0.489276558
884
885
886lb(1) = 1; ub(1) = getUnifRand(1, 7)
887do isam = 2, nsam
888 lb(isam) = ub(isam - 1) + 1
889 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
890end do
891lb
892+1, +5
893ub
894+4, +9
895sample = getUnifRand(0., 1., ub(nsam))
896sample
897+0.282088637, +0.648042798, +0.295889616, +0.151668191, +0.713053942, +0.620828569, +0.523952007, +0.509549379, +0.329731762
898rweight = getUnifRand(1., 2., size(sample, 1, IK))
899rweight
900+1.43753850, +1.78597784, +1.91361976, +1.66635203, +1.79687297, +1.18335116, +1.43050218, +1.29144895, +1.57350898
901mean(0) = getMean(sample, rweight)
902mean(0) ! reference
903+0.449187577
904do isam = 1, nsam
905 mean(isam) = getMean(sample(lb(isam):ub(isam)), rweight(lb(isam):ub(isam)))
906end do
907call setMeanMerged(meanMerged, mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
908meanMerged
909+0.449187547
910call setMeanMerged(mean(2), mean(1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
911mean(2)
912+0.449187547
913mean(0) ! reference
914+0.449187577
915
916
917!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
918!Compute the merged mean of a multivariate sample.
919!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
920
921
922dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
923do isam = 2, nsam
924 lb(isam) = ub(isam - 1) + 1
925 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
926end do
927lb
928+1, +2
929ub
930+1, +7
931ndim = getUnifRand(1, minval(ub - lb + 1, 1))
932call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
933call setResized(meanMerged, ndim)
934sample = getUnifRand(-1., +1., ndim, ub(nsam))
935sample
936+0.792492628E-1, -0.814815283, -0.527313471, -0.769526958, -0.952939153, +0.106447458, -0.987016201
937mean(:,0) = getMean(sample, dim)
938mean(:,0) ! reference
939-0.552273452
940do isam = 1, nsam
941 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
942end do
943call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
944meanMerged
945-0.552273512
946call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
947mean(:,2)
948-0.552273512
949mean(:,0) ! reference
950-0.552273452
951
952
953dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
954do isam = 2, nsam
955 lb(isam) = ub(isam - 1) + 1
956 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
957end do
958lb
959+1, +3
960ub
961+2, +3
962ndim = getUnifRand(1, minval(ub - lb + 1, 1))
963call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
964call setResized(meanMerged, ndim)
965sample = getUnifRand(-1., +1., ndim, ub(nsam))
966sample
967-0.981962085, -0.651926517, -0.735267639
968mean(:,0) = getMean(sample, dim)
969mean(:,0) ! reference
970-0.789718807
971do isam = 1, nsam
972 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
973end do
974call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
975meanMerged
976-0.789718747
977call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
978mean(:,2)
979-0.789718747
980mean(:,0) ! reference
981-0.789718807
982
983
984dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
985do isam = 2, nsam
986 lb(isam) = ub(isam - 1) + 1
987 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
988end do
989lb
990+1, +5
991ub
992+4, +6
993ndim = getUnifRand(1, minval(ub - lb + 1, 1))
994call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
995call setResized(meanMerged, ndim)
996sample = getUnifRand(-1., +1., ndim, ub(nsam))
997sample
998+0.863795877, -0.360506058, +0.924755096, +0.939900517, +0.948748112, +0.819427133
999-0.364199758, +0.722237945, -0.718837142, +0.304505348, -0.886546612, +0.672012687
1000mean(:,0) = getMean(sample, dim)
1001mean(:,0) ! reference
1002+0.689353466, -0.451379232E-1
1003do isam = 1, nsam
1004 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
1005end do
1006call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1007meanMerged
1008+0.689353466, -0.451379195E-1
1009call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1010mean(:,2)
1011+0.689353466, -0.451379195E-1
1012mean(:,0) ! reference
1013+0.689353466, -0.451379232E-1
1014
1015
1016dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1017do isam = 2, nsam
1018 lb(isam) = ub(isam - 1) + 1
1019 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1020end do
1021lb
1022+1, +2
1023ub
1024+1, +7
1025ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1026call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1027call setResized(meanMerged, ndim)
1028sample = getUnifRand(-1., +1., ndim, ub(nsam))
1029sample
1030-0.830671430, -0.171727538, -0.417463899, -0.731832266, +0.671187639, -0.894712806, +0.508095384
1031mean(:,0) = getMean(sample, dim)
1032mean(:,0) ! reference
1033-0.266732156
1034do isam = 1, nsam
1035 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
1036end do
1037call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1038meanMerged
1039-0.266732156
1040call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1041mean(:,2)
1042-0.266732156
1043mean(:,0) ! reference
1044-0.266732156
1045
1046
1047dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1048do isam = 2, nsam
1049 lb(isam) = ub(isam - 1) + 1
1050 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1051end do
1052lb
1053+1, +3
1054ub
1055+2, +6
1056ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1057call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1058call setResized(meanMerged, ndim)
1059sample = getUnifRand(-1., +1., ndim, ub(nsam))
1060sample
1061-0.650840878, -0.889423847, +0.224046469, -0.217802882, -0.441380978, -0.530807972
1062mean(:,0) = getMean(sample, dim)
1063mean(:,0) ! reference
1064-0.417701691
1065do isam = 1, nsam
1066 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
1067end do
1068call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1069meanMerged
1070-0.417701662
1071call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1072mean(:,2)
1073-0.417701662
1074mean(:,0) ! reference
1075-0.417701691
1076
1077
1078dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1079do isam = 2, nsam
1080 lb(isam) = ub(isam - 1) + 1
1081 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1082end do
1083lb
1084+1, +3
1085ub
1086+2, +5
1087ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1088call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1089call setResized(meanMerged, ndim)
1090sample = getUnifRand(-1., +1., ndim, ub(nsam))
1091sample
1092-0.348052502, +0.466094136, -0.484915733, -0.759090066, -0.997612596
1093+0.772763848, -0.393468142E-1, +0.243729353, -0.548236370E-1, -0.484283328
1094mean(:,0) = getMean(sample, dim)
1095mean(:,0) ! reference
1096-0.424715340, +0.876078829E-1
1097do isam = 1, nsam
1098 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
1099end do
1100call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1101meanMerged
1102-0.424715370, +0.876078829E-1
1103call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1104mean(:,2)
1105-0.424715370, +0.876078829E-1
1106mean(:,0) ! reference
1107-0.424715340, +0.876078829E-1
1108
1109
1110dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1111do isam = 2, nsam
1112 lb(isam) = ub(isam - 1) + 1
1113 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1114end do
1115lb
1116+1, +6
1117ub
1118+5, +10
1119ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1120call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1121call setResized(meanMerged, ndim)
1122sample = getUnifRand(-1., +1., ndim, ub(nsam))
1123sample
1124+0.396198750, -0.762101054, +0.969689369, +0.204674721, +0.795172453, -0.854734540, +0.329230547, -0.234186053, -0.215780854, +0.314898849
1125-0.185340166, -0.274718761, -0.306277514, +0.301264405, -0.267217159E-1, +0.811058521, +0.658574343, -0.489765525, -0.956857443, -0.833419800
1126-0.399742723, +0.608554840, -0.359128714E-1, +0.669092059, -0.329301000, +0.851661325, -0.249085426, +0.185821414, +0.263315082, +0.261289597
1127-0.332344055, -0.396055222, +0.368460655, -0.634359360, -0.381998777, +0.901958227, +0.560599923, -0.559601068, -0.351397157, +0.666177034
1128-0.135547638, +0.815477252, +0.891591311E-1, +0.271032333, -0.780617595, +0.978300214, -0.767898560E-2, -0.629989862, +0.253776312, -0.507682443
1129mean(:,0) = getMean(sample, dim)
1130mean(:,0) ! reference
1131+0.943062231E-1, -0.130220369, +0.182569236, -0.158559810E-1, +0.346228741E-1
1132do isam = 1, nsam
1133 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
1134end do
1135call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1136meanMerged
1137+0.943062156E-1, -0.130220369, +0.182569236, -0.158559754E-1, +0.346228704E-1
1138call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1139mean(:,2)
1140+0.943062156E-1, -0.130220369, +0.182569236, -0.158559754E-1, +0.346228704E-1
1141mean(:,0) ! reference
1142+0.943062231E-1, -0.130220369, +0.182569236, -0.158559810E-1, +0.346228741E-1
1143
1144
1145dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1146do isam = 2, nsam
1147 lb(isam) = ub(isam - 1) + 1
1148 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1149end do
1150lb
1151+1, +4
1152ub
1153+3, +8
1154ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1155call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1156call setResized(meanMerged, ndim)
1157sample = getUnifRand(-1., +1., ndim, ub(nsam))
1158sample
1159-0.438054562, +0.977368116, +0.262950063, -0.437987804, +0.531090498E-1, +0.335477114, -0.878040314, -0.496594548
1160-0.887449384, -0.255854487, +0.749666452, -0.704798341, +0.632519722E-1, +0.500726104, +0.627201676, -0.743068576
1161mean(:,0) = getMean(sample, dim)
1162mean(:,0) ! reference
1163-0.777216107E-1, -0.812905729E-1
1164do isam = 1, nsam
1165 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
1166end do
1167call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1168meanMerged
1169-0.777216107E-1, -0.812905729E-1
1170call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1171mean(:,2)
1172-0.777216107E-1, -0.812905729E-1
1173mean(:,0) ! reference
1174-0.777216107E-1, -0.812905729E-1
1175
1176
1177dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1178do isam = 2, nsam
1179 lb(isam) = ub(isam - 1) + 1
1180 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1181end do
1182lb
1183+1, +6
1184ub
1185+5, +11
1186ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1187call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1188call setResized(meanMerged, ndim)
1189sample = getUnifRand(-1., +1., ndim, ub(nsam))
1190sample
1191+0.279461145, +0.356027722, +0.428345919, -0.444567561, -0.917892456, -0.291271448, -0.481053829, -0.439249873, -0.753668666, +0.574824929, +0.199058056
1192-0.777057290, -0.511084795, +0.236571193, -0.738955975, +0.980914474, +0.650825620, -0.168111682, +0.216277838, -0.635682106, +0.419185996, +0.214123487
1193+0.944226980E-1, -0.260097980, -0.148499846, -0.759849072, -0.103393078, -0.822131991, +0.996027231, +0.794895291, +0.807016015, -0.417787433, -0.657907724E-1
1194-0.917671561, -0.996422768E-2, -0.553924561, +0.124514103E-2, +0.434498429, +0.597463131, +0.996397018, +0.294252753, -0.905705690, +0.755450368, +0.667023301
1195mean(:,0) = getMean(sample, dim)
1196mean(:,0) ! reference
1197-0.135453284, -0.102721127E-1, +0.104373693E-1, +0.123551287
1198do isam = 1, nsam
1199 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
1200end do
1201call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1202meanMerged
1203-0.135453284, -0.102721229E-1, +0.104373619E-1, +0.123551294
1204call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1205mean(:,2)
1206-0.135453284, -0.102721229E-1, +0.104373619E-1, +0.123551294
1207mean(:,0) ! reference
1208-0.135453284, -0.102721127E-1, +0.104373693E-1, +0.123551287
1209
1210
1211dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1212do isam = 2, nsam
1213 lb(isam) = ub(isam - 1) + 1
1214 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1215end do
1216lb
1217+1, +7
1218ub
1219+6, +7
1220ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1221call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1222call setResized(meanMerged, ndim)
1223sample = getUnifRand(-1., +1., ndim, ub(nsam))
1224sample
1225-0.684290051, +0.782229781, +0.281648278, -0.645682335, -0.726508856, +0.632212996, +0.282883763
1226mean(:,0) = getMean(sample, dim)
1227mean(:,0) ! reference
1228-0.110723469E-1
1229do isam = 1, nsam
1230 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
1231end do
1232call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1233meanMerged
1234-0.110723488E-1
1235call setMeanMerged(mean(:,2), mean(:,1), real(ub(1), TKG) / real(ub(2), TKG))
1236mean(:,2)
1237-0.110723488E-1
1238mean(:,0) ! reference
1239-0.110723469E-1
1240
1241
1242!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1243!Compute the merged mean of a frequency weighted multivariate sample.
1244!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1245
1246
1247dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1248do isam = 2, nsam
1249 lb(isam) = ub(isam - 1) + 1
1250 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1251end do
1252lb
1253+1, +2
1254ub
1255+1, +2
1256ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1257call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1258call setResized(meanMerged, ndim)
1259sample = getUnifRand(-1., +1., ndim, ub(nsam))
1260sample
1261+0.482810855, +0.937024593
1262iweight = getUnifRand(1, 10, size(sample, dim, IK))
1263iweight
1264+5, +4
1265mean(:,0) = getMean(sample, 2_IK, iweight)
1266mean(:,0) ! reference
1267+0.684683621
1268do isam = 1, nsam
1269 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1270end do
1271call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1272meanMerged
1273+0.684683681
1274call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1275mean(:,2)
1276+0.684683681
1277mean(:,0) ! reference
1278+0.684683621
1279
1280
1281dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1282do isam = 2, nsam
1283 lb(isam) = ub(isam - 1) + 1
1284 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1285end do
1286lb
1287+1, +6
1288ub
1289+5, +7
1290ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1291call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1292call setResized(meanMerged, ndim)
1293sample = getUnifRand(-1., +1., ndim, ub(nsam))
1294sample
1295-0.576741457, -0.139735341, -0.573289156, +0.326125622E-1, -0.897829294, +0.283338547, -0.511916757
1296iweight = getUnifRand(1, 10, size(sample, dim, IK))
1297iweight
1298+4, +6, +1, +2, +6, +7, +10
1299mean(:,0) = getMean(sample, 2_IK, iweight)
1300mean(:,0) ! reference
1301-0.338228196
1302do isam = 1, nsam
1303 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1304end do
1305call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1306meanMerged
1307-0.338228196
1308call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1309mean(:,2)
1310-0.338228196
1311mean(:,0) ! reference
1312-0.338228196
1313
1314
1315dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1316do isam = 2, nsam
1317 lb(isam) = ub(isam - 1) + 1
1318 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1319end do
1320lb
1321+1, +3
1322ub
1323+2, +4
1324ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1325call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1326call setResized(meanMerged, ndim)
1327sample = getUnifRand(-1., +1., ndim, ub(nsam))
1328sample
1329-0.841204166, +0.985905647, -0.919012427, +0.287443757
1330iweight = getUnifRand(1, 10, size(sample, dim, IK))
1331iweight
1332+10, +9, +7, +4
1333mean(:,0) = getMean(sample, 2_IK, iweight)
1334mean(:,0) ! reference
1335-0.160740092
1336do isam = 1, nsam
1337 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1338end do
1339call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1340meanMerged
1341-0.160740092
1342call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1343mean(:,2)
1344-0.160740092
1345mean(:,0) ! reference
1346-0.160740092
1347
1348
1349dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1350do isam = 2, nsam
1351 lb(isam) = ub(isam - 1) + 1
1352 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1353end do
1354lb
1355+1, +4
1356ub
1357+3, +6
1358ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1359call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1360call setResized(meanMerged, ndim)
1361sample = getUnifRand(-1., +1., ndim, ub(nsam))
1362sample
1363-0.644922614, +0.565256953, -0.828196526, +0.732427835, +0.527273297, -0.257063031
1364iweight = getUnifRand(1, 10, size(sample, dim, IK))
1365iweight
1366+2, +10, +8, +10, +3, +1
1367mean(:,0) = getMean(sample, 2_IK, iweight)
1368mean(:,0) ! reference
1369+0.187829047
1370do isam = 1, nsam
1371 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1372end do
1373call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1374meanMerged
1375+0.187829033
1376call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1377mean(:,2)
1378+0.187829033
1379mean(:,0) ! reference
1380+0.187829047
1381
1382
1383dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1384do isam = 2, nsam
1385 lb(isam) = ub(isam - 1) + 1
1386 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1387end do
1388lb
1389+1, +8
1390ub
1391+7, +11
1392ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1393call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1394call setResized(meanMerged, ndim)
1395sample = getUnifRand(-1., +1., ndim, ub(nsam))
1396sample
1397+0.484885097, +0.110787153E-1, -0.757265687, +0.447370410, -0.790321827, -0.423355579, +0.737119436, -0.442779064E-2, -0.373531342, +0.897284985, +0.957667589
1398iweight = getUnifRand(1, 10, size(sample, dim, IK))
1399iweight
1400+6, +10, +6, +1, +2, +10, +10, +9, +8, +1, +7
1401mean(:,0) = getMean(sample, 2_IK, iweight)
1402mean(:,0) ! reference
1403+0.721960887E-1
1404do isam = 1, nsam
1405 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1406end do
1407call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1408meanMerged
1409+0.721960887E-1
1410call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1411mean(:,2)
1412+0.721960887E-1
1413mean(:,0) ! reference
1414+0.721960887E-1
1415
1416
1417dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1418do isam = 2, nsam
1419 lb(isam) = ub(isam - 1) + 1
1420 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1421end do
1422lb
1423+1, +2
1424ub
1425+1, +6
1426ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1427call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1428call setResized(meanMerged, ndim)
1429sample = getUnifRand(-1., +1., ndim, ub(nsam))
1430sample
1431-0.331791997, +0.163888931E-2, +0.600777745, -0.913177252, +0.160655975E-1, -0.749979854
1432iweight = getUnifRand(1, 10, size(sample, dim, IK))
1433iweight
1434+2, +1, +8, +7, +3, +1
1435mean(:,0) = getMean(sample, 2_IK, iweight)
1436mean(:,0) ! reference
1437-0.134079412
1438do isam = 1, nsam
1439 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1440end do
1441call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1442meanMerged
1443-0.134079397
1444call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1445mean(:,2)
1446-0.134079397
1447mean(:,0) ! reference
1448-0.134079412
1449
1450
1451dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1452do isam = 2, nsam
1453 lb(isam) = ub(isam - 1) + 1
1454 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1455end do
1456lb
1457+1, +5
1458ub
1459+4, +11
1460ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1461call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1462call setResized(meanMerged, ndim)
1463sample = getUnifRand(-1., +1., ndim, ub(nsam))
1464sample
1465-0.487235069, +0.339401484, -0.905389786, -0.504239202, +0.351963758, +0.431632638, +0.205532193, -0.898764491, -0.145130992, +0.679128885, +0.553917170
1466+0.788933158, -0.761665106E-1, +0.873067021, +0.236579061, -0.158794761, -0.958598852E-1, +0.153317928, +0.974743128, -0.469558358, +0.768145561, -0.303305864
1467-0.758252621, -0.166301608, +0.885466933, -0.781515956, +0.444861531, -0.961319327, +0.109414339, -0.751347542E-1, +0.428465486, -0.640737653, -0.348250866E-1
1468+0.340287924, +0.560660362, +0.126956105, +0.890405178, -0.931557894, +0.238565683, +0.883358121, -0.776199460, -0.485056162, -0.979454756, -0.906894207
1469iweight = getUnifRand(1, 10, size(sample, dim, IK))
1470iweight
1471+7, +7, +1, +5, +10, +3, +6, +5, +6, +1, +3
1472mean(:,0) = getMean(sample, 2_IK, iweight)
1473mean(:,0) ! reference
1474-0.266188998E-1, +0.148227811, -0.107832253, -0.538071059E-1
1475do isam = 1, nsam
1476 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1477end do
1478call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1479meanMerged
1480-0.266189016E-1, +0.148227811, -0.107832253, -0.538071245E-1
1481call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1482mean(:,2)
1483-0.266189016E-1, +0.148227811, -0.107832253, -0.538071245E-1
1484mean(:,0) ! reference
1485-0.266188998E-1, +0.148227811, -0.107832253, -0.538071059E-1
1486
1487
1488dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1489do isam = 2, nsam
1490 lb(isam) = ub(isam - 1) + 1
1491 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1492end do
1493lb
1494+1, +8
1495ub
1496+7, +14
1497ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1498call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1499call setResized(meanMerged, ndim)
1500sample = getUnifRand(-1., +1., ndim, ub(nsam))
1501sample
1502-0.557179093, -0.528724432, -0.443182349, -0.707323909, +0.916378140, +0.123061895, -0.873217344, -0.624743938, +0.630951285, +0.595556021, +0.743969679E-1, -0.967899799, +0.823269606, +0.250709176
1503+0.674735785, +0.910176396, +0.586703658, -0.718786836, +0.308602452, +0.693881512, -0.855715632, +0.352848887, -0.156742334, +0.991984129, +0.867736578, -0.141610742, +0.650489926, +0.775660276
1504+0.300323248, +0.238235950, +0.433225513, -0.355881214, +0.819830060, +0.811681986, +0.479061723, -0.901137233, +0.111945748, +0.397492528, -0.279687285, +0.882044435, -0.918330908, +0.407699585
1505-0.393614769, +0.168538332, -0.803650618E-1, +0.872630239, -0.610213757, +0.281566262, +0.812224746, +0.147985220, -0.317014456, -0.440087318, +0.190901160, -0.300594568E-1, -0.749153733, +0.874656796
1506iweight = getUnifRand(1, 10, size(sample, dim, IK))
1507iweight
1508+2, +2, +1, +2, +8, +7, +10, +9, +3, +1, +2, +6, +4, +1
1509mean(:,0) = getMean(sample, 2_IK, iweight)
1510mean(:,0) ! reference
1511-0.169283718, +0.156010613, +0.205516815, +0.766632408E-1
1512do isam = 1, nsam
1513 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1514end do
1515call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1516meanMerged
1517-0.169283733, +0.156010598, +0.205516815, +0.766632408E-1
1518call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1519mean(:,2)
1520-0.169283733, +0.156010598, +0.205516815, +0.766632408E-1
1521mean(:,0) ! reference
1522-0.169283718, +0.156010613, +0.205516815, +0.766632408E-1
1523
1524
1525dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1526do isam = 2, nsam
1527 lb(isam) = ub(isam - 1) + 1
1528 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1529end do
1530lb
1531+1, +3
1532ub
1533+2, +6
1534ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1535call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1536call setResized(meanMerged, ndim)
1537sample = getUnifRand(-1., +1., ndim, ub(nsam))
1538sample
1539+0.876397729, +0.836336255, +0.976139665, -0.693070889E-2, -0.476712465, -0.982724071
1540-0.257339001, +0.708578825E-1, +0.947160125, -0.254223228, -0.348530769, -0.551466703
1541iweight = getUnifRand(1, 10, size(sample, dim, IK))
1542iweight
1543+7, +6, +4, +8, +5, +1
1544mean(:,0) = getMean(sample, 2_IK, iweight)
1545mean(:,0) ! reference
1546+0.375342816, -0.617900491E-1
1547do isam = 1, nsam
1548 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1549end do
1550call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1551meanMerged
1552+0.375342846, -0.617900528E-1
1553call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1554mean(:,2)
1555+0.375342846, -0.617900528E-1
1556mean(:,0) ! reference
1557+0.375342816, -0.617900491E-1
1558
1559
1560dim = 2; lb(1) = 1; ub(1) = getUnifRand(1, 7)
1561do isam = 2, nsam
1562 lb(isam) = ub(isam - 1) + 1
1563 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1564end do
1565lb
1566+1, +3
1567ub
1568+2, +5
1569ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1570call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1571call setResized(meanMerged, ndim)
1572sample = getUnifRand(-1., +1., ndim, ub(nsam))
1573sample
1574-0.646990180, -0.668978691E-1, +0.583444357, +0.650681853, -0.343840957
1575iweight = getUnifRand(1, 10, size(sample, dim, IK))
1576iweight
1577+4, +2, +3, +10, +7
1578mean(:,0) = getMean(sample, 2_IK, iweight)
1579mean(:,0) ! reference
1580+0.120327279
1581do isam = 1, nsam
1582 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), 2_IK, iweight(lb(isam):ub(isam)))
1583end do
1584call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1585meanMerged
1586+0.120327257
1587call setMeanMerged(mean(:,2), mean(:,1), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG))
1588mean(:,2)
1589+0.120327257
1590mean(:,0) ! reference
1591+0.120327279
1592
1593
1594!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1595!Compute the merged mean of a reliability weighted multivariate sample.
1596!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1597
1598
1599lb(1) = 1; ub(1) = getUnifRand(1, 7)
1600do isam = 2, nsam
1601 lb(isam) = ub(isam - 1) + 1
1602 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1603end do
1604lb
1605+1, +7
1606ub
1607+6, +11
1608ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1609call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1610call setResized(meanMerged, ndim)
1611sample = getUnifRand(-1., +1., ndim, ub(nsam))
1612sample
1613+0.490081310E-2, -0.618352771, -0.959584713E-1, +0.265206099E-1, +0.385996222, -0.763190150, +0.227920651, +0.899359584, +0.703910232, -0.857494235, -0.138871908
1614-0.589534402, +0.805036664, +0.458170056, -0.920660496E-1, -0.412250161, +0.467877626, -0.306022406, +0.731227398E-1, +0.243213415, -0.406769276, +0.761890411
1615-0.390493512, -0.837315321, +0.380594492, -0.259415269, +0.214852452, -0.805881739, +0.429843307, -0.336398005, +0.998563766E-1, +0.351866484E-1, +0.867916346
1616-0.613202572, +0.735453367E-1, +0.917518139E-1, -0.957813740, -0.390845656, +0.644869804, +0.616778135, -0.392547965, +0.766606331, +0.914812326, +0.414704919
1617+0.212546587E-1, -0.720114708E-1, +0.494683504, +0.113169074, +0.904285908E-1, -0.741189837, -0.825785637, +0.310242414, -0.897150755, +0.789320827, +0.988754988
1618rweight = getUnifRand(1, 10, size(sample, dim, IK))
1619rweight
1620+3.00000000, +5.00000000, +4.00000000, +10.0000000, +1.00000000, +10.0000000, +2.00000000, +7.00000000, +10.0000000, +10.0000000, +3.00000000
1621mean(:,0) = getMean(sample, dim, rweight)
1622mean(:,0) ! reference
1623-0.868106559E-1, +0.122846968, -0.181764871, +0.183367133, -0.322961174E-1
1624do isam = 1, nsam
1625 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1626end do
1627call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1628meanMerged
1629-0.868106633E-1, +0.122846983, -0.181764901, +0.183367103, -0.322961286E-1
1630call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1631mean(:,2)
1632-0.868106633E-1, +0.122846983, -0.181764901, +0.183367103, -0.322961286E-1
1633mean(:,0) ! reference
1634-0.868106559E-1, +0.122846968, -0.181764871, +0.183367133, -0.322961174E-1
1635
1636
1637lb(1) = 1; ub(1) = getUnifRand(1, 7)
1638do isam = 2, nsam
1639 lb(isam) = ub(isam - 1) + 1
1640 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1641end do
1642lb
1643+1, +6
1644ub
1645+5, +11
1646ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1647call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1648call setResized(meanMerged, ndim)
1649sample = getUnifRand(-1., +1., ndim, ub(nsam))
1650sample
1651-0.918369651, +0.993519902, +0.234226704, +0.768038392, +0.150652885, -0.555616617E-1, -0.226068377, -0.328015089E-1, -0.721464515, +0.480678082E-1, -0.709615827
1652+0.482911587, -0.518410206, -0.843890429, +0.131964087, +0.945477605, +0.479867339, -0.501396298, -0.435482144, +0.971004367, -0.346634507, +0.269029379
1653+0.436396599E-1, +0.813761592, +0.332093716, +0.910766244, +0.840726852, +0.358553767, +0.573832393, +0.269068956, -0.693935513, -0.681281567, -0.384454727E-1
1654-0.370262742, -0.192218900, -0.286040187, +0.647199154E-1, +0.864447355E-1, +0.885421872, +0.265064597, -0.850031972, -0.779427052, -0.957997680, -0.629420280
1655-0.195411921, -0.383368015, -0.438025475, +0.485347867, -0.124623775E-1, +0.823004246E-1, +0.737084270, -0.126240253E-1, -0.840724230, -0.191241026, +0.241984487
1656rweight = getUnifRand(1, 10, size(sample, dim, IK))
1657rweight
1658+10.0000000, +7.00000000, +5.00000000, +1.00000000, +3.00000000, +6.00000000, +8.00000000, +8.00000000, +2.00000000, +9.00000000, +7.00000000
1659mean(:,0) = getMean(sample, dim, rweight)
1660mean(:,0) ! reference
1661-0.124543712, -0.600171201E-1, +0.186850801, -0.304671049, -0.272591468E-1
1662do isam = 1, nsam
1663 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1664end do
1665call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1666meanMerged
1667-0.124543726, -0.600171313E-1, +0.186850816, -0.304671049, -0.272591487E-1
1668call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1669mean(:,2)
1670-0.124543726, -0.600171313E-1, +0.186850816, -0.304671049, -0.272591487E-1
1671mean(:,0) ! reference
1672-0.124543712, -0.600171201E-1, +0.186850801, -0.304671049, -0.272591468E-1
1673
1674
1675lb(1) = 1; ub(1) = getUnifRand(1, 7)
1676do isam = 2, nsam
1677 lb(isam) = ub(isam - 1) + 1
1678 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1679end do
1680lb
1681+1, +6
1682ub
1683+5, +11
1684ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1685call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1686call setResized(meanMerged, ndim)
1687sample = getUnifRand(-1., +1., ndim, ub(nsam))
1688sample
1689+0.693571806, +0.548965931, +0.298013926, -0.640588284, -0.775781751, +0.853347182, +0.903656006, +0.770961881, -0.834916234, -0.784127593, +0.929877520
1690+0.444394469, -0.820051312, +0.254682660, +0.775581360, +0.545632005, +0.397066236, +0.730875611, +0.822177291, -0.130782485, -0.693093181, -0.394552827
1691-0.213444471, -0.549933910E-1, -0.169571161, -0.758736014, +0.875785708, -0.469015718, +0.931631804, +0.374847651, +0.321420312, -0.103204966, +0.725391030
1692rweight = getUnifRand(1, 10, size(sample, dim, IK))
1693rweight
1694+2.00000000, +7.00000000, +1.00000000, +9.00000000, +1.00000000, +6.00000000, +5.00000000, +2.00000000, +7.00000000, +7.00000000, +4.00000000
1695mean(:,0) = getMean(sample, dim, rweight)
1696mean(:,0) ! reference
1697+0.500653498E-1, +0.640150085E-1, +0.173686072E-2
1698do isam = 1, nsam
1699 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1700end do
1701call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1702meanMerged
1703+0.500653386E-1, +0.640150309E-1, +0.173684955E-2
1704call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1705mean(:,2)
1706+0.500653386E-1, +0.640150309E-1, +0.173684955E-2
1707mean(:,0) ! reference
1708+0.500653498E-1, +0.640150085E-1, +0.173686072E-2
1709
1710
1711lb(1) = 1; ub(1) = getUnifRand(1, 7)
1712do isam = 2, nsam
1713 lb(isam) = ub(isam - 1) + 1
1714 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1715end do
1716lb
1717+1, +6
1718ub
1719+5, +12
1720ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1721call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1722call setResized(meanMerged, ndim)
1723sample = getUnifRand(-1., +1., ndim, ub(nsam))
1724sample
1725-0.729637146, +0.176713705, -0.665539503, -0.943921566, +0.374110699, -0.523375511, +0.643920183, +0.239128709, +0.351717472, -0.780999422, -0.730462551, +0.254852891
1726+0.355911136, +0.771650076E-1, -0.928746462, +0.724209189, +0.811487675, +0.618055701, +0.679419041E-1, +0.495584726, +0.204825401, +0.276764631, -0.900089622, -0.473653674
1727-0.236621380, -0.528571606E-1, +0.187513947, -0.233240604, -0.515419960, -0.522814870, +0.264549375, -0.564676166, +0.343389630, +0.785636306, +0.825565696, +0.772434235
1728+0.611832142E-1, +0.535980105, -0.179605961, +0.979906440, +0.502887845, +0.772213697, +0.240775347, -0.134008050, +0.996972084, -0.982038498, -0.192457914, -0.917604566
1729-0.168083191, -0.218716383, -0.166893363, +0.283057809, -0.477798223, -0.436976910, -0.819481611, +0.772327900, -0.885001421E-1, -0.920504808, +0.243986726, +0.228114843
1730rweight = getUnifRand(1, 10, size(sample, dim, IK))
1731rweight
1732+8.00000000, +7.00000000, +2.00000000, +9.00000000, +2.00000000, +6.00000000, +7.00000000, +4.00000000, +7.00000000, +9.00000000, +5.00000000, +4.00000000
1733mean(:,0) = getMean(sample, dim, rweight)
1734mean(:,0) ! reference
1735-0.265063196, +0.190952465, +0.116147481, +0.185678214, -0.195115939
1736do isam = 1, nsam
1737 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1738end do
1739call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1740meanMerged
1741-0.265063226, +0.190952480, +0.116147473, +0.185678229, -0.195115939
1742call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1743mean(:,2)
1744-0.265063226, +0.190952480, +0.116147473, +0.185678229, -0.195115939
1745mean(:,0) ! reference
1746-0.265063196, +0.190952465, +0.116147481, +0.185678214, -0.195115939
1747
1748
1749lb(1) = 1; ub(1) = getUnifRand(1, 7)
1750do isam = 2, nsam
1751 lb(isam) = ub(isam - 1) + 1
1752 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1753end do
1754lb
1755+1, +5
1756ub
1757+4, +7
1758ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1759call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1760call setResized(meanMerged, ndim)
1761sample = getUnifRand(-1., +1., ndim, ub(nsam))
1762sample
1763-0.914240122, -0.505501628, -0.715151191, +0.546790838, +0.200616121, +0.304987907, -0.711024404
1764-0.852367759, +0.144938827, -0.208925128, -0.993572235, -0.845995426, -0.106132984, +0.403549671E-1
1765rweight = getUnifRand(1, 10, size(sample, dim, IK))
1766rweight
1767+8.00000000, +8.00000000, +2.00000000, +4.00000000, +9.00000000, +10.0000000, +1.00000000
1768mean(:,0) = getMean(sample, dim, rweight)
1769mean(:,0) ! reference
1770-0.153730333, -0.444916785
1771do isam = 1, nsam
1772 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1773end do
1774call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1775meanMerged
1776-0.153730363, -0.444916785
1777call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1778mean(:,2)
1779-0.153730363, -0.444916785
1780mean(:,0) ! reference
1781-0.153730333, -0.444916785
1782
1783
1784lb(1) = 1; ub(1) = getUnifRand(1, 7)
1785do isam = 2, nsam
1786 lb(isam) = ub(isam - 1) + 1
1787 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1788end do
1789lb
1790+1, +4
1791ub
1792+3, +10
1793ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1794call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1795call setResized(meanMerged, ndim)
1796sample = getUnifRand(-1., +1., ndim, ub(nsam))
1797sample
1798+0.166830301, +0.615802288, +0.511212111, +0.590934753E-1, +0.824588180, -0.316116333, +0.905315399, +0.723817110, +0.189445853, -0.876321673
1799-0.617155790, +0.163237691, -0.169392943, +0.494131565, +0.346638441, -0.311914444, -0.737410188, +0.248980045, -0.206802964, +0.824367285
1800-0.977749705, +0.717214823, -0.839931726, -0.608612180, +0.586789489, -0.769068360, +0.738113284, +0.424382091, +0.961120963, +0.898550510
1801rweight = getUnifRand(1, 10, size(sample, dim, IK))
1802rweight
1803+5.00000000, +6.00000000, +6.00000000, +3.00000000, +3.00000000, +9.00000000, +7.00000000, +7.00000000, +3.00000000, +5.00000000
1804mean(:,0) = getMean(sample, dim, rweight)
1805mean(:,0) ! reference
1806+0.277646184, -0.615780354E-1, +0.537305512E-1
1807do isam = 1, nsam
1808 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1809end do
1810call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1811meanMerged
1812+0.277646124, -0.615780354E-1, +0.537305400E-1
1813call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1814mean(:,2)
1815+0.277646124, -0.615780354E-1, +0.537305400E-1
1816mean(:,0) ! reference
1817+0.277646184, -0.615780354E-1, +0.537305512E-1
1818
1819
1820lb(1) = 1; ub(1) = getUnifRand(1, 7)
1821do isam = 2, nsam
1822 lb(isam) = ub(isam - 1) + 1
1823 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1824end do
1825lb
1826+1, +8
1827ub
1828+7, +12
1829ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1830call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1831call setResized(meanMerged, ndim)
1832sample = getUnifRand(-1., +1., ndim, ub(nsam))
1833sample
1834+0.353289008, -0.817571044, +0.190581203, +0.488667130, +0.619613647, +0.248415589, -0.932450056, +0.510061979E-1, -0.357687473, +0.239421487, -0.563956857, -0.292816162E-1
1835+0.569867373, -0.604049683, -0.900990963, -0.488433838E-1, +0.518765092, -0.883710980, -0.196419954, +0.637945056, -0.988142252, -0.407198548, -0.556697845E-1, +0.585686922
1836-0.543030977, +0.920944810, -0.563188910, -0.890938163, +0.976711392, -0.826065660, -0.225739002, +0.859814405, -0.247413158, -0.636211634, +0.180789709, -0.167532086
1837rweight = getUnifRand(1, 10, size(sample, dim, IK))
1838rweight
1839+7.00000000, +10.0000000, +7.00000000, +2.00000000, +4.00000000, +7.00000000, +4.00000000, +8.00000000, +9.00000000, +1.00000000, +10.0000000, +2.00000000
1840mean(:,0) = getMean(sample, dim, rweight)
1841mean(:,0) ! reference
1842-0.157374218, -0.238519534, +0.337166227E-1
1843do isam = 1, nsam
1844 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1845end do
1846call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1847meanMerged
1848-0.157374233, -0.238519520, +0.337166265E-1
1849call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1850mean(:,2)
1851-0.157374233, -0.238519520, +0.337166265E-1
1852mean(:,0) ! reference
1853-0.157374218, -0.238519534, +0.337166227E-1
1854
1855
1856lb(1) = 1; ub(1) = getUnifRand(1, 7)
1857do isam = 2, nsam
1858 lb(isam) = ub(isam - 1) + 1
1859 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1860end do
1861lb
1862+1, +2
1863ub
1864+1, +4
1865ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1866call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1867call setResized(meanMerged, ndim)
1868sample = getUnifRand(-1., +1., ndim, ub(nsam))
1869sample
1870-0.376463890, +0.658104539, -0.252119660, +0.448137283
1871rweight = getUnifRand(1, 10, size(sample, dim, IK))
1872rweight
1873+7.00000000, +6.00000000, +3.00000000, +4.00000000
1874mean(:,0) = getMean(sample, dim, rweight)
1875mean(:,0) ! reference
1876+0.117478512
1877do isam = 1, nsam
1878 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1879end do
1880call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1881meanMerged
1882+0.117478520
1883call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1884mean(:,2)
1885+0.117478520
1886mean(:,0) ! reference
1887+0.117478512
1888
1889
1890lb(1) = 1; ub(1) = getUnifRand(1, 7)
1891do isam = 2, nsam
1892 lb(isam) = ub(isam - 1) + 1
1893 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1894end do
1895lb
1896+1, +7
1897ub
1898+6, +11
1899ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1900call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1901call setResized(meanMerged, ndim)
1902sample = getUnifRand(-1., +1., ndim, ub(nsam))
1903sample
1904+0.125908852E-1, +0.772338748, -0.360266924, -0.378568053, +0.121783853, -0.377302051, -0.209380031, -0.317876220, -0.827109694, -0.731416821, +0.821354270
1905+0.551213026E-1, +0.304964662, -0.839131594, +0.226134181, -0.183019638E-1, +0.570032716, -0.306092978, +0.602171898, -0.242266417, -0.673467159, -0.452617407
1906-0.284770489, +0.566995502, +0.837630868, -0.285176039, -0.772216916, -0.207923055, -0.351421952, -0.948767543, -0.510708451, -0.507664800, -0.107165933
1907+0.560374737, -0.782676935E-1, +0.229095101, +0.491057873, +0.456467986, -0.725031018, +0.308146477E-1, +0.487562537, +0.599642754, -0.892860889, +0.724265218
1908-0.173839331E-1, +0.546600103, -0.461957335, -0.693558574, +0.982072353E-1, +0.544248343, -0.483799458, -0.772729278, +0.534604311, +0.716113567, +0.668600321
1909rweight = getUnifRand(1, 10, size(sample, dim, IK))
1910rweight
1911+8.00000000, +8.00000000, +7.00000000, +7.00000000, +6.00000000, +9.00000000, +5.00000000, +8.00000000, +8.00000000, +5.00000000, +6.00000000
1912mean(:,0) = getMean(sample, dim, rweight)
1913mean(:,0) ! reference
1914-0.136274934, -0.145983510E-1, -0.220699504, +0.179798156, +0.636463612E-1
1915do isam = 1, nsam
1916 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1917end do
1918call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1919meanMerged
1920-0.136274919, -0.145983659E-1, -0.220699519, +0.179798156, +0.636463538E-1
1921call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1922mean(:,2)
1923-0.136274919, -0.145983659E-1, -0.220699519, +0.179798156, +0.636463538E-1
1924mean(:,0) ! reference
1925-0.136274934, -0.145983510E-1, -0.220699504, +0.179798156, +0.636463612E-1
1926
1927
1928lb(1) = 1; ub(1) = getUnifRand(1, 7)
1929do isam = 2, nsam
1930 lb(isam) = ub(isam - 1) + 1
1931 ub(isam) = ub(isam - 1) + getUnifRand(1, 7)
1932end do
1933lb
1934+1, +8
1935ub
1936+7, +13
1937ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1938call setRebound(mean, [1_IK, 0_IK], [ndim, nsam])
1939call setResized(meanMerged, ndim)
1940sample = getUnifRand(-1., +1., ndim, ub(nsam))
1941sample
1942-0.395796537, -0.259393215, +0.869132876, +0.925002456, +0.695981741, -0.735191226, +0.750091314, +0.119754672, +0.130768061, -0.449383140, -0.798324466, +0.110868692, -0.595369935
1943-0.552621722, -0.912419081, -0.169597030, -0.452745795, +0.252257705, +0.743482471, +0.844059587, +0.516170859, +0.379201770, +0.983427048, +0.234758735, +0.520820022, -0.182270646
1944+0.340086102, +0.165648460, -0.461425543, -0.568191528, +0.576981306E-1, +0.314194560, +0.547724128, -0.896787643E-1, -0.390197039, -0.308937430, -0.910579562, -0.952359676, -0.898418546
1945+0.769450665E-1, -0.592040062, +0.292707086, -0.850330591E-1, -0.104155540, -0.867553711, +0.931960821, +0.904703259, +0.324153781, +0.275048494, +0.264284372, +0.414369464, +0.121279716
1946rweight = getUnifRand(1, 10, size(sample, dim, IK))
1947rweight
1948+6.00000000, +6.00000000, +4.00000000, +9.00000000, +9.00000000, +3.00000000, +9.00000000, +6.00000000, +10.0000000, +5.00000000, +9.00000000, +9.00000000, +7.00000000
1949mean(:,0) = getMean(sample, dim, rweight)
1950mean(:,0) ! reference
1951+0.880880356E-1, +0.172657058, -0.288845122, +0.208308592
1952do isam = 1, nsam
1953 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1954end do
1955call setMeanMerged(meanMerged, mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1956meanMerged
1957+0.880880430E-1, +0.172657058, -0.288845122, +0.208308592
1958call setMeanMerged(mean(:,2), mean(:,1), real(sum(rweight(:ub(1))) / sum(rweight), TKG))
1959mean(:,2)
1960+0.880880430E-1, +0.172657058, -0.288845122, +0.208308592
1961mean(:,0) ! reference
1962+0.880880356E-1, +0.172657058, -0.288845122, +0.208308592
1963
1964
Test:
test_pm_sampleMean


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, April 21, 2017, 1:54 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin

Definition at line 4194 of file pm_sampleMean.F90.


The documentation for this interface was generated from the following file: