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

Return the merged covariance of a sample resulting from the merger of two separate (potentially weighted) samples \(A\) and \(B\).
More...

Detailed Description

Return the merged covariance of a sample resulting from the merger of two separate (potentially weighted) samples \(A\) and \(B\).

See the documentation of pm_sampleCov for more information and definition online updating of sample covariance.

Note
The input and output variances of this generic interface are all biased variances.
A biased covariance can be readily unbiased by multiplying it with the appropriate bias-correction factors detailed in the documentation of pm_sampleCov.
Parameters
[out]cov: The output object of the same type and kind and rank and shape as covA, containing the biased covariance of the sample resulting form the merger of the two samples.
(optional. If missing, the resulting merged covariance will be written to the argument covB.)
[in,out]covB: The input or input/output object of the same type and kind and rank and shape as the input argument covA, containing the biased covariance of the second sample that must be merged with the first sample.
If the input argument cov is missing, then covB contains the updated biased covariance of the merged sample on return.
Otherwise, the contents of covB remain intact upon return.
[in]covA: The input contiguous matrix of shape (1 : ndim, 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 biased covariance of the first sample that must be merged with the second sample.
[in]meanDiff: The input object of the same type and kind as the input argument covA of size size(covA, 1), containing the difference of mean of the two samples meanDiff = meanA - meanB.
The subtraction order is immaterial.
[in]fracA: The input scalar of type real of the same kind as covA, 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)).
[in]subset: The input scalar constant that can be:
  1. The constant lowDia, implying that only the lower-diagonal subset of the input cov, covB, and covA matrices must be accessed and/or manipulated.
  2. The constant uppDia, implying that only the upper-diagonal subset of the input cov, covB, and covA matrices must be accessed and/or manipulated.
This input argument is merely serves to resolve the different procedures of this generic interface from each other at compile-time.


Possible calling interfaces

call setCovMerged(covB(1:ndim, 1:ndim), covA(1:ndim, 1:ndim), meanDiff(1:ndim), fracA, subset)
call setCovMerged(cov(1:ndim, 1:ndim), covB(1:ndim, 1:ndim), covA(1:ndim, 1:ndim), meanDiff(1:ndim), fracA, subset)
Return the merged covariance of a sample resulting from the merger of two separate (potentially weigh...
This module contains classes and procedures for computing the properties related to the covariance ma...
Warning
The condition 0 < fracA .and. fracA < 1 must hold for the corresponding input arguments.
The condition all(shape(covB) == shape(covA)) must hold for the corresponding input arguments.
The condition all(size(meanDiff) == shape(covA)) must hold for the corresponding input arguments.
The condition all(shape(cov) == shape(covA)) 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
setCovMeanMerged
getVarCorrection
getMeanMerged
setMeanMerged
getVarMerged
setVarMerged


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_sampleCov, only: getCov
6 use pm_sampleMean, only: getMean
7 use pm_sampleCov, only: uppDia, lowDia
10 use pm_arrayResize, only: setResized
11 use pm_distUnif, only: getUnifRand
12 use pm_arrayRange, only: getRange
13 use pm_io, only: display_type
14
15 implicit none
16
17 integer(IK) , parameter :: nsam = 2
18 integer(IK) , allocatable :: iweight(:)
19 integer(IK) :: isam, ndim, lb(nsam), ub(nsam)
20 integer(IK) :: dim, itry, ntry = 10
21 type(display_type) :: disp
22 disp = display_type(file = "main.out.F90")
23
24 call disp%skip()
25 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%show("!Compute the biased merged covariance of a multivariate sample.")
27 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%skip()
29
30 block
31 real(TKG), allocatable :: mean(:,:), cov(:,:,:), covMerged(:,:)
32 real(TKG), allocatable :: sample(:,:)
33 do itry = 1, ntry
34 call disp%skip()
35 call disp%show("dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)")
36 dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
37 call disp%show("do isam = 2, nsam")
38 call disp%show(" lb(isam) = ub(isam - 1) + 1")
39 call disp%show(" ub(isam) = ub(isam - 1) + getUnifRand(2, 7)")
40 call disp%show("end do")
41 do isam = 2, nsam
42 lb(isam) = ub(isam - 1) + 1
43 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
44 end do
45 call disp%show("lb")
46 call disp%show( lb )
47 call disp%show("ub")
48 call disp%show( ub )
49 call disp%show("ndim = getUnifRand(1, minval(ub - lb + 1, 1))")
50 ndim = getUnifRand(1, minval(ub - lb + 1, 1))
51 call disp%show("call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])")
52 call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
53 call disp%show("call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])")
54 call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
55 call disp%show("call setResized(covMerged, [ndim, ndim])")
56 call setResized(covMerged, [ndim, ndim])
57 call disp%show("sample = getUnifRand(-1., +1., ndim, ub(nsam))")
58 sample = getUnifRand(-1., +1., ndim, ub(nsam))
59 call disp%show("sample")
60 call disp%show( sample )
61 call disp%show("cov(:,:,0) = getCov(sample, dim)")
62 cov(:,:,0) = getCov(sample, dim)
63 call disp%show("cov(:,:,0) ! reference")
64 call disp%show( cov(:,:,0) )
65 call disp%show("do isam = 1, nsam")
66 call disp%show(" cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)")
67 call disp%show(" mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)")
68 call disp%show("end do")
69 do isam = 1, nsam
70 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
71 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
72 end do
73 call disp%show("call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)")
74 call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
75 call disp%show("covMerged")
76 call disp%show( covMerged )
77 call disp%show("call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)")
78 call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
79 call disp%show("covMerged")
80 call disp%show( covMerged )
81 call disp%show("call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)")
82 call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
83 call disp%show("cov(:,:,2)")
84 call disp%show( cov(:,:,2) )
85 call disp%show("cov(:,:,0) ! reference")
86 call disp%show( cov(:,:,0) )
87 call disp%skip()
88 end do
89 end block
90
91 call disp%skip()
92 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
93 call disp%show("!Compute the biased merged covariance of a frequency weighted multivariate sample.")
94 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
95 call disp%skip()
96
97 block
98 real(TKG), allocatable :: mean(:,:), cov(:,:,:), covMerged(:,:)
99 real(TKG), allocatable :: sample(:,:)
100 do itry = 1, ntry
101 call disp%skip()
102 call disp%show("dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)")
103 dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
104 call disp%show("do isam = 2, nsam")
105 call disp%show(" lb(isam) = ub(isam - 1) + 1")
106 call disp%show(" ub(isam) = ub(isam - 1) + getUnifRand(2, 7)")
107 call disp%show("end do")
108 do isam = 2, nsam
109 lb(isam) = ub(isam - 1) + 1
110 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
111 end do
112 call disp%show("lb")
113 call disp%show( lb )
114 call disp%show("ub")
115 call disp%show( ub )
116 call disp%show("ndim = getUnifRand(1, minval(ub - lb + 1, 1))")
117 ndim = getUnifRand(1, minval(ub - lb + 1, 1))
118 call disp%show("call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])")
119 call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
120 call disp%show("call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])")
121 call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
122 call disp%show("call setResized(covMerged, [ndim, ndim])")
123 call setResized(covMerged, [ndim, ndim])
124 call disp%show("sample = getUnifRand(-1., +1., ndim, ub(nsam))")
125 sample = getUnifRand(-1., +1., ndim, ub(nsam))
126 call disp%show("sample")
127 call disp%show( sample )
128 call disp%show("iweight = getUnifRand(1, 10, size(sample, dim, IK))")
129 iweight = getUnifRand(1, 10, size(sample, dim, IK))
130 call disp%show("iweight")
131 call disp%show( iweight )
132 call disp%show("cov(:,:,0) = getCov(sample, dim, iweight)")
133 cov(:,:,0) = getCov(sample, dim, iweight)
134 call disp%show("cov(:,:,0) ! reference")
135 call disp%show( cov(:,:,0) )
136 call disp%show("do isam = 1, nsam")
137 call disp%show(" cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))")
138 call disp%show(" mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))")
139 call disp%show("end do")
140 do isam = 1, nsam
141 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
142 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
143 end do
144 call disp%show("call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)")
145 call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
146 call disp%show("covMerged")
147 call disp%show( covMerged )
148 call disp%show("call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)")
149 call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
150 call disp%show("covMerged")
151 call disp%show( covMerged )
152 call disp%show("call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)")
153 call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
154 call disp%show("cov(:,:,2)")
155 call disp%show( cov(:,:,2) )
156 call disp%show("cov(:,:,0) ! reference")
157 call disp%show( cov(:,:,0) )
158 call disp%skip()
159 end do
160 end block
161
162 call disp%skip()
163 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
164 call disp%show("!Compute the biased merged covariance of a reliability weighted multivariate sample.")
165 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
166 call disp%skip()
167
168 block
169 real(TKG), allocatable :: mean(:,:), cov(:,:,:), covMerged(:,:)
170 real(TKG), allocatable :: sample(:,:)
171 real(TKG), allocatable :: rweight(:)
172 do itry = 1, ntry
173 call disp%skip()
174 call disp%show("dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)")
175 dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
176 call disp%show("do isam = 2, nsam")
177 call disp%show(" lb(isam) = ub(isam - 1) + 1")
178 call disp%show(" ub(isam) = ub(isam - 1) + getUnifRand(2, 7)")
179 call disp%show("end do")
180 do isam = 2, nsam
181 lb(isam) = ub(isam - 1) + 1
182 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
183 end do
184 call disp%show("lb")
185 call disp%show( lb )
186 call disp%show("ub")
187 call disp%show( ub )
188 call disp%show("ndim = getUnifRand(1, minval(ub - lb + 1, 1))")
189 ndim = getUnifRand(1, minval(ub - lb + 1, 1))
190 call disp%show("call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])")
191 call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
192 call disp%show("call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])")
193 call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
194 call disp%show("call setResized(covMerged, [ndim, ndim])")
195 call setResized(covMerged, [ndim, ndim])
196 call disp%show("sample = getUnifRand(-1., +1., ndim, ub(nsam))")
197 sample = getUnifRand(-1., +1., ndim, ub(nsam))
198 call disp%show("sample")
199 call disp%show( sample )
200 call disp%show("rweight = getUnifRand(1., 2., size(sample, dim, IK))")
201 rweight = getUnifRand(1., 2., size(sample, dim, IK))
202 call disp%show("rweight")
203 call disp%show( rweight )
204 call disp%show("cov(:,:,0) = getCov(sample, dim, rweight)")
205 cov(:,:,0) = getCov(sample, dim, rweight)
206 call disp%show("cov(:,:,0) ! reference")
207 call disp%show( cov(:,:,0) )
208 call disp%show("do isam = 1, nsam")
209 call disp%show(" cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))")
210 call disp%show(" mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))")
211 call disp%show("end do")
212 do isam = 1, nsam
213 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
214 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
215 end do
216 call disp%show("call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)")
217 call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
218 call disp%show("covMerged")
219 call disp%show( covMerged )
220 call disp%show("call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)")
221 call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
222 call disp%show("covMerged")
223 call disp%show( covMerged )
224 call disp%show("call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)")
225 call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
226 call disp%show("cov(:,:,2)")
227 call disp%show( cov(:,:,2) )
228 call disp%show("cov(:,:,0) ! reference")
229 call disp%show( cov(:,:,0) )
230 call disp%skip()
231 end do
232 end block
233
234end 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 (optionally unbiased) covariance matrix of a pair of (potentially weighted) t...
Generate and return the (weighted) mean of an input sample of nsam observations with ndim = 1 or 2 at...
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
This module contains classes and procedures for computing the first moment (i.e., the statistical mea...
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 biased merged covariance of a multivariate sample.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
8do isam = 2, nsam
9 lb(isam) = ub(isam - 1) + 1
10 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
11end do
12lb
13+1, +5
14ub
15+4, +8
16ndim = getUnifRand(1, minval(ub - lb + 1, 1))
17call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
18call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
19call setResized(covMerged, [ndim, ndim])
20sample = getUnifRand(-1., +1., ndim, ub(nsam))
21sample
22-0.730689526, +0.174534321E-1, +0.183695674, +0.683780670, -0.559766293, +0.301911831, +0.747577667, +0.578095913E-1
23+0.477809548, -0.452203751, +0.482423067, +0.439996362, -0.726205826, +0.746573091, -0.845496297, -0.195052624E-1
24-0.131714344E-1, -0.830839396, -0.650135636, +0.553827286, -0.140310526E-1, +0.615285873, +0.131123543, -0.559346437
25cov(:,:,0) = getCov(sample, dim)
26cov(:,:,0) ! reference
27+0.242581904, +0.276127458E-2, +0.726257265E-1
28+0.276127458E-2, +0.332221717, +0.848705173E-1
29+0.726257265E-1, +0.848705173E-1, +0.256889105
30do isam = 1, nsam
31 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
32 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
33end do
34call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
35covMerged
36+0.242581859, +0.276129507E-2, +0.726257265E-1
37+0.700649232E-44, +0.332221687, +0.848705173E-1
38+0.00000000, +0.229218613E-17, +0.256889105
39call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
40covMerged
41+0.242581859, +0.276129507E-2, +0.726257265E-1
42+0.276129507E-2, +0.332221687, +0.848705173E-1
43+0.726257265E-1, +0.848705173E-1, +0.256889105
44call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
45cov(:,:,2)
46+0.242581859, +0.276129507E-2, +0.726257265E-1
47+0.285798013E-1, +0.332221687, +0.848705173E-1
48+0.589051247E-1, +0.101532042, +0.256889105
49cov(:,:,0) ! reference
50+0.242581904, +0.276127458E-2, +0.726257265E-1
51+0.276127458E-2, +0.332221717, +0.848705173E-1
52+0.726257265E-1, +0.848705173E-1, +0.256889105
53
54
55dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
56do isam = 2, nsam
57 lb(isam) = ub(isam - 1) + 1
58 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
59end do
60lb
61+1, +3
62ub
63+2, +4
64ndim = getUnifRand(1, minval(ub - lb + 1, 1))
65call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
66call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
67call setResized(covMerged, [ndim, ndim])
68sample = getUnifRand(-1., +1., ndim, ub(nsam))
69sample
70-0.946599126, -0.990906954, +0.977254272, -0.387032866
71+0.384290814, -0.765278220, -0.679953933, +0.679383636
72cov(:,:,0) = getCov(sample, dim)
73cov(:,:,0) ! reference
74+0.632243276, -0.165349573
75-0.165349573, +0.405208290
76do isam = 1, nsam
77 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
78 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
79end do
80call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
81covMerged
82+0.632243216, -0.165349558
83+0.276129507E-2, +0.405208290
84call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
85covMerged
86+0.632243216, -0.165349558
87-0.165349558, +0.405208290
88call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
89cov(:,:,2)
90+0.632243216, -0.165349558
91-0.463631690, +0.405208290
92cov(:,:,0) ! reference
93+0.632243276, -0.165349573
94-0.165349573, +0.405208290
95
96
97dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
98do isam = 2, nsam
99 lb(isam) = ub(isam - 1) + 1
100 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
101end do
102lb
103+1, +6
104ub
105+5, +10
106ndim = getUnifRand(1, minval(ub - lb + 1, 1))
107call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
108call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
109call setResized(covMerged, [ndim, ndim])
110sample = getUnifRand(-1., +1., ndim, ub(nsam))
111sample
112-0.462055445, -0.787469745, -0.821069121, +0.628360391, -0.953035951, -0.536895037, +0.833436251E-1, +0.980564356E-1, -0.853996992, -0.755120277
113cov(:,:,0) = getCov(sample, dim)
114cov(:,:,0) ! reference
115+0.251435012
116do isam = 1, nsam
117 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
118 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
119end do
120call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
121covMerged
122+0.251435041
123call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
124covMerged
125+0.251435041
126call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
127cov(:,:,2)
128+0.251435041
129cov(:,:,0) ! reference
130+0.251435012
131
132
133dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
134do isam = 2, nsam
135 lb(isam) = ub(isam - 1) + 1
136 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
137end do
138lb
139+1, +4
140ub
141+3, +9
142ndim = getUnifRand(1, minval(ub - lb + 1, 1))
143call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
144call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
145call setResized(covMerged, [ndim, ndim])
146sample = getUnifRand(-1., +1., ndim, ub(nsam))
147sample
148-0.379431725, -0.836435914, -0.912897825, +0.529910326, +0.443402410, -0.593599200, -0.311179161, +0.589435935, +0.177296877
149-0.901885629, +0.936819196, -0.620117545, -0.370423198, +0.496287465, +0.361477137E-1, +0.895470738, +0.405763268, +0.410896540
150cov(:,:,0) = getCov(sample, dim)
151cov(:,:,0) ! reference
152+0.310727090, +0.384052582E-1
153+0.384052582E-1, +0.379015505
154do isam = 1, nsam
155 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
156 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
157end do
158call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
159covMerged
160+0.310727060, +0.384052694E-1
161+0.306926403E-40, +0.379015714
162call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
163covMerged
164+0.310727060, +0.384052694E-1
165+0.384052694E-1, +0.379015714
166call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
167cov(:,:,2)
168+0.310727060, +0.384052694E-1
169-0.375374369E-1, +0.379015714
170cov(:,:,0) ! reference
171+0.310727090, +0.384052582E-1
172+0.384052582E-1, +0.379015505
173
174
175dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
176do isam = 2, nsam
177 lb(isam) = ub(isam - 1) + 1
178 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
179end do
180lb
181+1, +4
182ub
183+3, +8
184ndim = getUnifRand(1, minval(ub - lb + 1, 1))
185call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
186call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
187call setResized(covMerged, [ndim, ndim])
188sample = getUnifRand(-1., +1., ndim, ub(nsam))
189sample
190+0.374810696E-1, +0.411958218, -0.492504120, -0.254053473, +0.359842300, -0.803278089, -0.631706715, +0.978754759E-1
191+0.911535025, -0.708054304, -0.819551587, -0.214745879, +0.444084287, -0.847451687, +0.599342346, +0.779874086
192+0.434346080, +0.815168619, -0.727621913, -0.450674295E-1, +0.637009263, -0.814595938, +0.837425828, -0.515295625
193cov(:,:,0) = getCov(sample, dim)
194cov(:,:,0) ! reference
195+0.182323232, +0.952538848E-1, +0.140626520
196+0.952538848E-1, +0.491273522, +0.185841322
197+0.140626520, +0.185841322, +0.421564698
198do isam = 1, nsam
199 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
200 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
201end do
202call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
203covMerged
204+0.182323217, +0.952538773E-1, +0.140626505
205+0.384052694E-1, +0.491273493, +0.185841322
206+0.00000000, -0.536895037, +0.421564728
207call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
208covMerged
209+0.182323217, +0.952538773E-1, +0.140626505
210+0.952538773E-1, +0.491273493, +0.185841322
211+0.140626505, +0.185841322, +0.421564728
212call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
213cov(:,:,2)
214+0.182323217, +0.952538773E-1, +0.140626505
215+0.156050190, +0.491273493, +0.185841322
216+0.680148974E-1, +0.213558272, +0.421564728
217cov(:,:,0) ! reference
218+0.182323232, +0.952538848E-1, +0.140626520
219+0.952538848E-1, +0.491273522, +0.185841322
220+0.140626520, +0.185841322, +0.421564698
221
222
223dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
224do isam = 2, nsam
225 lb(isam) = ub(isam - 1) + 1
226 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
227end do
228lb
229+1, +6
230ub
231+5, +7
232ndim = getUnifRand(1, minval(ub - lb + 1, 1))
233call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
234call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
235call setResized(covMerged, [ndim, ndim])
236sample = getUnifRand(-1., +1., ndim, ub(nsam))
237sample
238+0.364428043, -0.515398860, +0.652534008, -0.669225335, -0.665750146, -0.984044671, -0.577415824
239cov(:,:,0) = getCov(sample, dim)
240cov(:,:,0) ! reference
241+0.313962549
242do isam = 1, nsam
243 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
244 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
245end do
246call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
247covMerged
248+0.313962609
249call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
250covMerged
251+0.313962609
252call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
253cov(:,:,2)
254+0.313962609
255cov(:,:,0) ! reference
256+0.313962549
257
258
259dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
260do isam = 2, nsam
261 lb(isam) = ub(isam - 1) + 1
262 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
263end do
264lb
265+1, +4
266ub
267+3, +5
268ndim = getUnifRand(1, minval(ub - lb + 1, 1))
269call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
270call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
271call setResized(covMerged, [ndim, ndim])
272sample = getUnifRand(-1., +1., ndim, ub(nsam))
273sample
274-0.346350789, +0.476465344, -0.276171088, +0.197638154, +0.397109389
275cov(:,:,0) = getCov(sample, dim)
276cov(:,:,0) ! reference
277+0.115948103
278do isam = 1, nsam
279 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
280 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
281end do
282call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
283covMerged
284+0.115948111
285call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
286covMerged
287+0.115948111
288call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
289cov(:,:,2)
290+0.115948111
291cov(:,:,0) ! reference
292+0.115948103
293
294
295dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
296do isam = 2, nsam
297 lb(isam) = ub(isam - 1) + 1
298 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
299end do
300lb
301+1, +8
302ub
303+7, +14
304ndim = getUnifRand(1, minval(ub - lb + 1, 1))
305call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
306call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
307call setResized(covMerged, [ndim, ndim])
308sample = getUnifRand(-1., +1., ndim, ub(nsam))
309sample
310-0.859471560E-1, -0.539392233E-1, +0.373067975, -0.930675507, +0.681902170, +0.910912275, +0.373710155, -0.789715886, +0.948826909, +0.799502254, +0.781024218, +0.891802073, +0.666453838E-1, -0.902711034
311-0.828424335, +0.524303436, -0.131525040, -0.627688050, -0.260316730, +0.430219531, -0.911553264, -0.646403670, -0.803606153, -0.187167764, -0.556549788, -0.831809998, +0.767171383, +0.171775818E-1
312+0.175639510, +0.520481825, -0.715081453, -0.355767012, +0.539922595, -0.254513741, +0.178725600, -0.904102325E-1, +0.352688670, +0.256346583, +0.454981327, +0.583103538, -0.531791806, -0.446020246
313+0.115087390, -0.164973497, +0.244309187, -0.192922950, +0.712866783E-1, -0.646193981, -0.452042103, -0.457638502E-1, -0.818601370, -0.877954125, -0.377118587, -0.947598934, +0.339508891, -0.437502861
314cov(:,:,0) = getCov(sample, dim)
315cov(:,:,0) ! reference
316+0.440503895, -0.146509930E-1, +0.137967303, -0.120467961
317-0.146509930E-1, +0.284001112, -0.858512968E-1, +0.653083324E-1
318+0.137967303, -0.858512968E-1, +0.180216253, -0.763353407E-1
319-0.120467961, +0.653083324E-1, -0.763353407E-1, +0.164279029
320do isam = 1, nsam
321 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
322 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
323end do
324call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
325covMerged
326+0.440503955, -0.146510061E-1, +0.137967288, -0.120467931
327+0.00000000, +0.284001112, -0.858513117E-1, +0.653083324E-1
328+0.00000000, +0.00000000, +0.180216268, -0.763353258E-1
329+0.00000000, +0.00000000, +0.00000000, +0.164279044
330call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
331covMerged
332+0.440503955, -0.146510061E-1, +0.137967288, -0.120467931
333-0.146510061E-1, +0.284001112, -0.858513117E-1, +0.653083324E-1
334+0.137967288, -0.858513117E-1, +0.180216268, -0.763353258E-1
335-0.120467931, +0.653083324E-1, -0.763353258E-1, +0.164279044
336call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
337cov(:,:,2)
338+0.440503955, -0.146510061E-1, +0.137967288, -0.120467931
339-0.138241649, +0.284001112, -0.858513117E-1, +0.653083324E-1
340+0.243580267, -0.177357689, +0.180216268, -0.763353258E-1
341-0.193252400, +0.155635208, -0.136775777, +0.164279044
342cov(:,:,0) ! reference
343+0.440503895, -0.146509930E-1, +0.137967303, -0.120467961
344-0.146509930E-1, +0.284001112, -0.858512968E-1, +0.653083324E-1
345+0.137967303, -0.858512968E-1, +0.180216253, -0.763353407E-1
346-0.120467961, +0.653083324E-1, -0.763353407E-1, +0.164279029
347
348
349dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
350do isam = 2, nsam
351 lb(isam) = ub(isam - 1) + 1
352 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
353end do
354lb
355+1, +6
356ub
357+5, +7
358ndim = getUnifRand(1, minval(ub - lb + 1, 1))
359call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
360call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
361call setResized(covMerged, [ndim, ndim])
362sample = getUnifRand(-1., +1., ndim, ub(nsam))
363sample
364-0.658514977, +0.568775535, +0.121245623, +0.843169928, -0.488507986, +0.115345001, -0.210928202
365+0.851544261, +0.754020572, +0.840786695E-1, +0.820228815, -0.905382276, +0.559113622, -0.138116002
366cov(:,:,0) = getCov(sample, dim)
367cov(:,:,0) ! reference
368+0.252450913, +0.145961121
369+0.145961121, +0.362690657
370do isam = 1, nsam
371 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
372 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
373end do
374call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
375covMerged
376+0.252450913, +0.145961136
377-0.146510061E-1, +0.362690657
378call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
379covMerged
380+0.252450913, +0.145961136
381+0.145961136, +0.362690657
382call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
383cov(:,:,2)
384+0.252450913, +0.145961136
385+0.568718351E-1, +0.362690657
386cov(:,:,0) ! reference
387+0.252450913, +0.145961121
388+0.145961121, +0.362690657
389
390
391dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
392do isam = 2, nsam
393 lb(isam) = ub(isam - 1) + 1
394 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
395end do
396lb
397+1, +3
398ub
399+2, +7
400ndim = getUnifRand(1, minval(ub - lb + 1, 1))
401call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
402call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
403call setResized(covMerged, [ndim, ndim])
404sample = getUnifRand(-1., +1., ndim, ub(nsam))
405sample
406+0.233716965E-1, -0.427911639, -0.718568563E-1, +0.105924726, +0.804329872, -0.567386031, +0.342302322
407+0.142705083, +0.978607059, +0.813134789, -0.265102863, -0.338086009, +0.600153804, -0.785182476
408cov(:,:,0) = getCov(sample, dim)
409cov(:,:,0) ! reference
410+0.182836547, -0.202477112
411-0.202477112, +0.373258680
412do isam = 1, nsam
413 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
414 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
415end do
416call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
417covMerged
418+0.182836577, -0.202477127
419+0.145961136, +0.373258680
420call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
421covMerged
422+0.182836577, -0.202477127
423-0.202477127, +0.373258680
424call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
425cov(:,:,2)
426+0.182836577, -0.202477127
427-0.194157600, +0.373258680
428cov(:,:,0) ! reference
429+0.182836547, -0.202477112
430-0.202477112, +0.373258680
431
432
433!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
434!Compute the biased merged covariance of a frequency weighted multivariate sample.
435!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
436
437
438dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
439do isam = 2, nsam
440 lb(isam) = ub(isam - 1) + 1
441 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
442end do
443lb
444+1, +6
445ub
446+5, +12
447ndim = getUnifRand(1, minval(ub - lb + 1, 1))
448call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
449call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
450call setResized(covMerged, [ndim, ndim])
451sample = getUnifRand(-1., +1., ndim, ub(nsam))
452sample
453+0.663856387, +0.370827556, -0.847152233, -0.812152386, -0.829524755, +0.984963179, +0.806908131, -0.675281048, -0.852650881, -0.291886568, -0.352904320, +0.327542543
454+0.928904533, -0.202757239, -0.771390915, +0.441682816, -0.559954405, -0.913024783, +0.421713829, -0.983847380E-1, +0.350212097, -0.725432396, -0.642382145, -0.732728243E-1
455iweight = getUnifRand(1, 10, size(sample, dim, IK))
456iweight
457+7, +6, +7, +9, +2, +1, +7, +10, +8, +4, +7, +2
458cov(:,:,0) = getCov(sample, dim, iweight)
459cov(:,:,0) ! reference
460+0.413921237, +0.114881843
461+0.114881843, +0.299755603
462do isam = 1, nsam
463 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
464 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
465end do
466call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
467covMerged
468+0.413921297, +0.114881814
469+0.306926403E-40, +0.299755484
470call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
471covMerged
472+0.413921297, +0.114881814
473+0.114881814, +0.299755484
474call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
475cov(:,:,2)
476+0.413921297, +0.114881814
477+0.309136230E-1, +0.299755484
478cov(:,:,0) ! reference
479+0.413921237, +0.114881843
480+0.114881843, +0.299755603
481
482
483dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
484do isam = 2, nsam
485 lb(isam) = ub(isam - 1) + 1
486 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
487end do
488lb
489+1, +5
490ub
491+4, +9
492ndim = getUnifRand(1, minval(ub - lb + 1, 1))
493call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
494call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
495call setResized(covMerged, [ndim, ndim])
496sample = getUnifRand(-1., +1., ndim, ub(nsam))
497sample
498-0.407470942, -0.713191986, +0.470277667, +0.375116706, +0.545417070E-1, +0.243642926, +0.625978470, -0.471876025, -0.357685685
499iweight = getUnifRand(1, 10, size(sample, dim, IK))
500iweight
501+5, +10, +6, +2, +10, +4, +1, +10, +2
502cov(:,:,0) = getCov(sample, dim, iweight)
503cov(:,:,0) ! reference
504+0.181750387
505do isam = 1, nsam
506 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
507 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
508end do
509call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
510covMerged
511+0.181750387
512call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
513covMerged
514+0.181750387
515call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
516cov(:,:,2)
517+0.181750387
518cov(:,:,0) ! reference
519+0.181750387
520
521
522dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
523do isam = 2, nsam
524 lb(isam) = ub(isam - 1) + 1
525 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
526end do
527lb
528+1, +4
529ub
530+3, +8
531ndim = getUnifRand(1, minval(ub - lb + 1, 1))
532call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
533call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
534call setResized(covMerged, [ndim, ndim])
535sample = getUnifRand(-1., +1., ndim, ub(nsam))
536sample
537-0.700062752, +0.183719397, -0.225917935, -0.130042434, +0.173640966, -0.456657410, -0.664977551, +0.388418317
538-0.832275152, -0.783651352, +0.533238053, +0.136787891E-1, +0.270397305, +0.264082909, -0.162866950, +0.855092525
539iweight = getUnifRand(1, 10, size(sample, dim, IK))
540iweight
541+4, +1, +8, +2, +2, +10, +8, +4
542cov(:,:,0) = getCov(sample, dim, iweight)
543cov(:,:,0) ! reference
544+0.118133374, +0.104674026
545+0.104674026, +0.227452010
546do isam = 1, nsam
547 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
548 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
549end do
550call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
551covMerged
552+0.118133366, +0.104673967
553+0.306926403E-40, +0.227452025
554call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
555covMerged
556+0.118133366, +0.104673967
557+0.104673967, +0.227452025
558call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
559cov(:,:,2)
560+0.118133366, +0.104673967
561+0.106482521, +0.227452025
562cov(:,:,0) ! reference
563+0.118133374, +0.104674026
564+0.104674026, +0.227452010
565
566
567dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
568do isam = 2, nsam
569 lb(isam) = ub(isam - 1) + 1
570 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
571end do
572lb
573+1, +5
574ub
575+4, +8
576ndim = getUnifRand(1, minval(ub - lb + 1, 1))
577call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
578call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
579call setResized(covMerged, [ndim, ndim])
580sample = getUnifRand(-1., +1., ndim, ub(nsam))
581sample
582-0.504294038, +0.469903946E-1, +0.690548897, +0.807892323, -0.184420228, +0.149602890E-1, -0.276815057, +0.196006060
583iweight = getUnifRand(1, 10, size(sample, dim, IK))
584iweight
585+1, +4, +7, +6, +6, +6, +6, +10
586cov(:,:,0) = getCov(sample, dim, iweight)
587cov(:,:,0) ! reference
588+0.150893211
589do isam = 1, nsam
590 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
591 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
592end do
593call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
594covMerged
595+0.150893107
596call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
597covMerged
598+0.150893107
599call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
600cov(:,:,2)
601+0.150893107
602cov(:,:,0) ! reference
603+0.150893211
604
605
606dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
607do isam = 2, nsam
608 lb(isam) = ub(isam - 1) + 1
609 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
610end do
611lb
612+1, +5
613ub
614+4, +7
615ndim = getUnifRand(1, minval(ub - lb + 1, 1))
616call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
617call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
618call setResized(covMerged, [ndim, ndim])
619sample = getUnifRand(-1., +1., ndim, ub(nsam))
620sample
621-0.182912111, +0.121653080E-2, -0.447187424, -0.544398546, +0.342447400, +0.698184967E-2, +0.182927251
622+0.770518064, -0.689107895, -0.174644470, +0.484755874, +0.726112604, -0.731326580, +0.174746037
623-0.981148481, -0.731261611, -0.839406490, -0.645702958, -0.536799312, +0.631785750, -0.418595314
624iweight = getUnifRand(1, 10, size(sample, dim, IK))
625iweight
626+4, +3, +9, +5, +7, +1, +2
627cov(:,:,0) = getCov(sample, dim, iweight)
628cov(:,:,0) ! reference
629+0.115784690, +0.520413034E-1, +0.402940586E-1
630+0.520413034E-1, +0.262709767, -0.217508189E-1
631+0.402940586E-1, -0.217508189E-1, +0.840230733E-1
632do isam = 1, nsam
633 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
634 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
635end do
636call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
637covMerged
638+0.115784682, +0.520412959E-1, +0.402940437E-1
639+0.306926403E-40, +0.262709796, -0.217507929E-1
640+0.00000000, +0.840779079E-44, +0.840230435E-1
641call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
642covMerged
643+0.115784682, +0.520412959E-1, +0.402940437E-1
644+0.520412959E-1, +0.262709796, -0.217507929E-1
645+0.402940437E-1, -0.217507929E-1, +0.840230435E-1
646call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
647cov(:,:,2)
648+0.115784682, +0.520412959E-1, +0.402940437E-1
649+0.497264005E-1, +0.262709796, -0.217507929E-1
650-0.337774232E-1, -0.147378668, +0.840230435E-1
651cov(:,:,0) ! reference
652+0.115784690, +0.520413034E-1, +0.402940586E-1
653+0.520413034E-1, +0.262709767, -0.217508189E-1
654+0.402940586E-1, -0.217508189E-1, +0.840230733E-1
655
656
657dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
658do isam = 2, nsam
659 lb(isam) = ub(isam - 1) + 1
660 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
661end do
662lb
663+1, +5
664ub
665+4, +11
666ndim = getUnifRand(1, minval(ub - lb + 1, 1))
667call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
668call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
669call setResized(covMerged, [ndim, ndim])
670sample = getUnifRand(-1., +1., ndim, ub(nsam))
671sample
672-0.802071691, +0.828577161, -0.238057494, -0.982773066, -0.387693882, +0.437667966, +0.105431795, +0.973488092, +0.410717010, +0.643182755, -0.126897931
673iweight = getUnifRand(1, 10, size(sample, dim, IK))
674iweight
675+7, +1, +4, +1, +8, +5, +2, +10, +7, +9, +6
676cov(:,:,0) = getCov(sample, dim, iweight)
677cov(:,:,0) ! reference
678+0.355057031
679do isam = 1, nsam
680 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
681 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
682end do
683call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
684covMerged
685+0.355056882
686call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
687covMerged
688+0.355056882
689call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
690cov(:,:,2)
691+0.355056882
692cov(:,:,0) ! reference
693+0.355057031
694
695
696dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
697do isam = 2, nsam
698 lb(isam) = ub(isam - 1) + 1
699 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
700end do
701lb
702+1, +6
703ub
704+5, +7
705ndim = getUnifRand(1, minval(ub - lb + 1, 1))
706call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
707call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
708call setResized(covMerged, [ndim, ndim])
709sample = getUnifRand(-1., +1., ndim, ub(nsam))
710sample
711+0.869203568, -0.175490379, +0.161351442, +0.192652464, +0.454697609, -0.140027046, +0.888376713
712iweight = getUnifRand(1, 10, size(sample, dim, IK))
713iweight
714+3, +3, +3, +7, +7, +3, +4
715cov(:,:,0) = getCov(sample, dim, iweight)
716cov(:,:,0) ! reference
717+0.129042953
718do isam = 1, nsam
719 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
720 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
721end do
722call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
723covMerged
724+0.129043013
725call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
726covMerged
727+0.129043013
728call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
729cov(:,:,2)
730+0.129043013
731cov(:,:,0) ! reference
732+0.129042953
733
734
735dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
736do isam = 2, nsam
737 lb(isam) = ub(isam - 1) + 1
738 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
739end do
740lb
741+1, +7
742ub
743+6, +10
744ndim = getUnifRand(1, minval(ub - lb + 1, 1))
745call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
746call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
747call setResized(covMerged, [ndim, ndim])
748sample = getUnifRand(-1., +1., ndim, ub(nsam))
749sample
750-0.312891006E-1, -0.643205047, +0.785211802, +0.927758455, -0.985159755, +0.881144166, +0.162091374, +0.515345335E-1, +0.675864816, +0.680277228
751-0.465568423, -0.680073023, +0.890941620, +0.264859319, +0.829611778, +0.123730898, +0.706754684, -0.701587677, +0.966189981, -0.619274139
752iweight = getUnifRand(1, 10, size(sample, dim, IK))
753iweight
754+7, +9, +7, +10, +9, +4, +5, +5, +9, +3
755cov(:,:,0) = getCov(sample, dim, iweight)
756cov(:,:,0) ! reference
757+0.467610121, +0.100294508
758+0.100294508, +0.436529011
759do isam = 1, nsam
760 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
761 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
762end do
763call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
764covMerged
765+0.467610151, +0.100294486
766+0.306926403E-40, +0.436528951
767call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
768covMerged
769+0.467610151, +0.100294486
770+0.100294486, +0.436528951
771call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
772cov(:,:,2)
773+0.467610151, +0.100294486
774+0.971634537E-1, +0.436528951
775cov(:,:,0) ! reference
776+0.467610121, +0.100294508
777+0.100294508, +0.436529011
778
779
780dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
781do isam = 2, nsam
782 lb(isam) = ub(isam - 1) + 1
783 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
784end do
785lb
786+1, +3
787ub
788+2, +8
789ndim = getUnifRand(1, minval(ub - lb + 1, 1))
790call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
791call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
792call setResized(covMerged, [ndim, ndim])
793sample = getUnifRand(-1., +1., ndim, ub(nsam))
794sample
795+0.668229699, -0.351782441, -0.265810490E-1, -0.441505671, -0.972831249, +0.875248909E-1, +0.550184011, +0.733505964
796+0.907635450, -0.304288387, -0.650396466, -0.265776157, -0.367276669E-1, -0.729820967, -0.492727995, +0.612545013E-1
797iweight = getUnifRand(1, 10, size(sample, dim, IK))
798iweight
799+6, +6, +1, +7, +8, +4, +1, +7
800cov(:,:,0) = getCov(sample, dim, iweight)
801cov(:,:,0) ! reference
802+0.405942261, +0.125824735
803+0.125824735, +0.217718601
804do isam = 1, nsam
805 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
806 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
807end do
808call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
809covMerged
810+0.405942172, +0.125824720
811+0.100294486, +0.217718601
812call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
813covMerged
814+0.405942172, +0.125824720
815+0.125824720, +0.217718601
816call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
817cov(:,:,2)
818+0.405942172, +0.125824720
819-0.332863000E-2, +0.217718601
820cov(:,:,0) ! reference
821+0.405942261, +0.125824735
822+0.125824735, +0.217718601
823
824
825dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
826do isam = 2, nsam
827 lb(isam) = ub(isam - 1) + 1
828 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
829end do
830lb
831+1, +4
832ub
833+3, +5
834ndim = getUnifRand(1, minval(ub - lb + 1, 1))
835call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
836call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
837call setResized(covMerged, [ndim, ndim])
838sample = getUnifRand(-1., +1., ndim, ub(nsam))
839sample
840+0.388753414E-1, -0.322981477, -0.633426547, -0.999809146, -0.717864394
841-0.377255678, -0.483701944, +0.146544456, -0.766735077E-1, -0.565060496
842iweight = getUnifRand(1, 10, size(sample, dim, IK))
843iweight
844+8, +8, +8, +9, +10
845cov(:,:,0) = getCov(sample, dim, iweight)
846cov(:,:,0) ! reference
847+0.124290712, -0.338921919E-1
848-0.338921919E-1, +0.708793178E-1
849do isam = 1, nsam
850 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
851 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
852end do
853call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
854covMerged
855+0.124290705, -0.338921882E-1
856+0.125824720, +0.708793178E-1
857call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
858covMerged
859+0.124290705, -0.338921882E-1
860-0.338921882E-1, +0.708793178E-1
861call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
862cov(:,:,2)
863+0.124290705, -0.338921882E-1
864-0.343291759E-1, +0.708793178E-1
865cov(:,:,0) ! reference
866+0.124290712, -0.338921919E-1
867-0.338921919E-1, +0.708793178E-1
868
869
870!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
871!Compute the biased merged covariance of a reliability weighted multivariate sample.
872!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
873
874
875dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
876do isam = 2, nsam
877 lb(isam) = ub(isam - 1) + 1
878 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
879end do
880lb
881+1, +5
882ub
883+4, +6
884ndim = getUnifRand(1, minval(ub - lb + 1, 1))
885call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
886call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
887call setResized(covMerged, [ndim, ndim])
888sample = getUnifRand(-1., +1., ndim, ub(nsam))
889sample
890+0.870274186, -0.991506100, +0.368181825, -0.547364950, +0.999474883, -0.154445529
891+0.511259794, +0.408414721, +0.820556760, +0.940234423, +0.733433366, -0.939947009
892rweight = getUnifRand(1., 2., size(sample, dim, IK))
893rweight
894+1.03623295, +1.26095629, +1.63207519, +1.67217135, +1.34095991, +1.18148041
895cov(:,:,0) = getCov(sample, dim, rweight)
896cov(:,:,0) ! reference
897+0.502747834, +0.620557740E-1
898+0.620557740E-1, +0.371575356
899do isam = 1, nsam
900 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
901 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
902end do
903call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
904covMerged
905+0.502747834, +0.620558113E-1
906+0.306926403E-40, +0.371575415
907call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
908covMerged
909+0.502747834, +0.620558113E-1
910+0.620558113E-1, +0.371575415
911call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
912cov(:,:,2)
913+0.502747834, +0.480807334
914+0.620558113E-1, +0.371575415
915cov(:,:,0) ! reference
916+0.502747834, +0.620557740E-1
917+0.620557740E-1, +0.371575356
918
919
920dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
921do isam = 2, nsam
922 lb(isam) = ub(isam - 1) + 1
923 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
924end do
925lb
926+1, +3
927ub
928+2, +6
929ndim = getUnifRand(1, minval(ub - lb + 1, 1))
930call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
931call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
932call setResized(covMerged, [ndim, ndim])
933sample = getUnifRand(-1., +1., ndim, ub(nsam))
934sample
935+0.691048265, +0.740772843, -0.840807080, -0.682481170, -0.177993774, +0.467064857
936rweight = getUnifRand(1., 2., size(sample, dim, IK))
937rweight
938+1.15253329, +1.55860722, +1.55412745, +1.76245379, +1.66519260, +1.57567096
939cov(:,:,0) = getCov(sample, dim, rweight)
940cov(:,:,0) ! reference
941+0.401382148
942do isam = 1, nsam
943 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
944 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
945end do
946call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
947covMerged
948+0.401382148
949call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
950covMerged
951+0.401382148
952call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
953cov(:,:,2)
954+0.401382148
955cov(:,:,0) ! reference
956+0.401382148
957
958
959dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
960do isam = 2, nsam
961 lb(isam) = ub(isam - 1) + 1
962 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
963end do
964lb
965+1, +6
966ub
967+5, +7
968ndim = getUnifRand(1, minval(ub - lb + 1, 1))
969call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
970call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
971call setResized(covMerged, [ndim, ndim])
972sample = getUnifRand(-1., +1., ndim, ub(nsam))
973sample
974+0.398148179, -0.527587295, -0.583539009, +0.976370931, +0.653176665, -0.795626044, +0.250953317
975-0.837774038, -0.343427896, +0.202022076, -0.647808313, +0.799802780, -0.298324823, +0.417882442
976rweight = getUnifRand(1., 2., size(sample, dim, IK))
977rweight
978+1.27997518, +1.20901752, +1.17209601, +1.52410865, +1.29827309, +1.04914808, +1.79405916
979cov(:,:,0) = getCov(sample, dim, rweight)
980cov(:,:,0) ! reference
981+0.384199560, -0.108779976E-1
982-0.108779976E-1, +0.311271399
983do isam = 1, nsam
984 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
985 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
986end do
987call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
988covMerged
989+0.384199560, -0.108780172E-1
990+0.306926403E-40, +0.311271489
991call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
992covMerged
993+0.384199560, -0.108780172E-1
994-0.108780172E-1, +0.311271489
995call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
996cov(:,:,2)
997+0.384199560, +0.174528942
998-0.108780172E-1, +0.311271489
999cov(:,:,0) ! reference
1000+0.384199560, -0.108779976E-1
1001-0.108779976E-1, +0.311271399
1002
1003
1004dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1005do isam = 2, nsam
1006 lb(isam) = ub(isam - 1) + 1
1007 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1008end do
1009lb
1010+1, +3
1011ub
1012+2, +8
1013ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1014call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1015call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1016call setResized(covMerged, [ndim, ndim])
1017sample = getUnifRand(-1., +1., ndim, ub(nsam))
1018sample
1019-0.960645795, -0.627111554, +0.796313286, +0.771232009, -0.603791952, +0.174230933, +0.487685919, -0.849461555E-2
1020rweight = getUnifRand(1., 2., size(sample, dim, IK))
1021rweight
1022+1.80196762, +1.04828525, +1.67848074, +1.49191523, +1.41728687, +1.75760233, +1.41154253, +1.15368915
1023cov(:,:,0) = getCov(sample, dim, rweight)
1024cov(:,:,0) ! reference
1025+0.419050932
1026do isam = 1, nsam
1027 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1028 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1029end do
1030call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1031covMerged
1032+0.419050813
1033call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1034covMerged
1035+0.419050813
1036call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1037cov(:,:,2)
1038+0.419050813
1039cov(:,:,0) ! reference
1040+0.419050932
1041
1042
1043dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1044do isam = 2, nsam
1045 lb(isam) = ub(isam - 1) + 1
1046 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1047end do
1048lb
1049+1, +6
1050ub
1051+5, +11
1052ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1053call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1054call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1055call setResized(covMerged, [ndim, ndim])
1056sample = getUnifRand(-1., +1., ndim, ub(nsam))
1057sample
1058+0.974749804, -0.615403652, -0.426107287, -0.308581352, -0.417941689, +0.298369527, -0.159184217, +0.648955822, +0.950057149, +0.478833795, +0.916819692
1059rweight = getUnifRand(1., 2., size(sample, dim, IK))
1060rweight
1061+1.86715579, +1.94725442, +1.97527015, +1.89408898, +1.03421712, +1.84572613, +1.73083353, +1.41008782, +1.79386687, +1.02597880, +1.48380542
1062cov(:,:,0) = getCov(sample, dim, rweight)
1063cov(:,:,0) ! reference
1064+0.358218163
1065do isam = 1, nsam
1066 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1067 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1068end do
1069call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1070covMerged
1071+0.358218253
1072call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1073covMerged
1074+0.358218253
1075call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1076cov(:,:,2)
1077+0.358218253
1078cov(:,:,0) ! reference
1079+0.358218163
1080
1081
1082dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1083do isam = 2, nsam
1084 lb(isam) = ub(isam - 1) + 1
1085 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1086end do
1087lb
1088+1, +6
1089ub
1090+5, +9
1091ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1092call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1093call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1094call setResized(covMerged, [ndim, ndim])
1095sample = getUnifRand(-1., +1., ndim, ub(nsam))
1096sample
1097+0.418940544, -0.238521576, +0.577281237, -0.907191515, +0.229027987, -0.419914961, +0.664691925E-1, +0.370019555, +0.940028429
1098+0.336863637, -0.952550173E-1, +0.309643149, -0.930635571, +0.242483974, +0.281979561, +0.719059706E-1, +0.711679459, -0.552076578
1099rweight = getUnifRand(1., 2., size(sample, dim, IK))
1100rweight
1101+1.27995753, +1.51105499, +1.35377562, +1.80416691, +1.46107578, +1.49078465, +1.66627026, +1.00531673, +1.08269715
1102cov(:,:,0) = getCov(sample, dim, rweight)
1103cov(:,:,0) ! reference
1104+0.289658338, +0.125951514
1105+0.125951514, +0.229427144
1106do isam = 1, nsam
1107 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1108 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1109end do
1110call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1111covMerged
1112+0.289658368, +0.125951514
1113+0.306926403E-40, +0.229427144
1114call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1115covMerged
1116+0.289658368, +0.125951514
1117+0.125951514, +0.229427144
1118call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1119cov(:,:,2)
1120+0.289658368, -0.109702945
1121+0.125951514, +0.229427144
1122cov(:,:,0) ! reference
1123+0.289658338, +0.125951514
1124+0.125951514, +0.229427144
1125
1126
1127dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1128do isam = 2, nsam
1129 lb(isam) = ub(isam - 1) + 1
1130 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1131end do
1132lb
1133+1, +7
1134ub
1135+6, +13
1136ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1137call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1138call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1139call setResized(covMerged, [ndim, ndim])
1140sample = getUnifRand(-1., +1., ndim, ub(nsam))
1141sample
1142-0.116342306E-1, +0.122706294, -0.108284235, -0.845807314, +0.435151815, +0.247197628, +0.976198673, +0.488863587, +0.270626068, -0.777511120, -0.244551063, +0.132064700, +0.472542763
1143-0.246724486, -0.301316977, +0.627304196, -0.510698676, -0.692843914, -0.265240669E-1, +0.625134349, +0.279241681, -0.310252547, +0.229469538, +0.637716413, -0.969700933, -0.936755300
1144+0.944617867, +0.142844796, +0.239712238, +0.274327993E-1, -0.611779690, -0.647225857, -0.930213690, +0.360064149, -0.620508194, -0.534427166E-2, -0.294285178, +0.717393756, -0.443050265
1145+0.592299819, +0.553686380, -0.252180696, +0.519427180, +0.879992962, -0.513559222, -0.963591814, -0.139496446, -0.697885871, +0.527045250, +0.909265280E-1, +0.525159001, +0.547252655
1146rweight = getUnifRand(1., 2., size(sample, dim, IK))
1147rweight
1148+1.56536078, +1.55579054, +1.60155630, +1.32733440, +1.66204166, +1.08505321, +1.27863312, +1.58737767, +1.92929327, +1.82076478, +1.32020164, +1.03257108, +1.37034059
1149cov(:,:,0) = getCov(sample, dim, rweight)
1150cov(:,:,0) ! reference
1151+0.235387415, -0.194309521E-1, -0.962856039E-1, -0.113151498
1152-0.194309521E-1, +0.282119334, -0.147216897E-1, -0.165110275
1153-0.962856039E-1, -0.147216897E-1, +0.283279181, +0.132138297
1154-0.113151498, -0.165110275, +0.132138297, +0.318178713
1155do isam = 1, nsam
1156 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1157 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1158end do
1159call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1160covMerged
1161+0.235387519, -0.194309540E-1, -0.962856039E-1, -0.113151543
1162+0.125951514, +0.282119364, -0.147216981E-1, -0.165110230
1163+0.00000000, -0.858513117E-1, +0.283279151, +0.132138208
1164+0.00000000, +0.653083324E-1, -0.763353258E-1, +0.318178564
1165call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1166covMerged
1167+0.235387519, -0.194309540E-1, -0.962856039E-1, -0.113151543
1168-0.194309540E-1, +0.282119364, -0.147216981E-1, -0.165110230
1169-0.962856039E-1, -0.147216981E-1, +0.283279151, +0.132138208
1170-0.113151543, -0.165110230, +0.132138208, +0.318178564
1171call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1172cov(:,:,2)
1173+0.235387519, -0.363548845E-1, -0.940104872E-1, -0.189470693
1174-0.194309540E-1, +0.282119364, -0.666179433E-1, -0.131700113
1175-0.962856039E-1, -0.147216981E-1, +0.283279151, +0.177193880
1176-0.113151543, -0.165110230, +0.132138208, +0.318178564
1177cov(:,:,0) ! reference
1178+0.235387415, -0.194309521E-1, -0.962856039E-1, -0.113151498
1179-0.194309521E-1, +0.282119334, -0.147216897E-1, -0.165110275
1180-0.962856039E-1, -0.147216897E-1, +0.283279181, +0.132138297
1181-0.113151498, -0.165110275, +0.132138297, +0.318178713
1182
1183
1184dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1185do isam = 2, nsam
1186 lb(isam) = ub(isam - 1) + 1
1187 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1188end do
1189lb
1190+1, +3
1191ub
1192+2, +7
1193ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1194call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1195call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1196call setResized(covMerged, [ndim, ndim])
1197sample = getUnifRand(-1., +1., ndim, ub(nsam))
1198sample
1199-0.695487022, -0.232330918, -0.518701434, -0.469264507, -0.903761029, -0.573341250, -0.690165997
1200rweight = getUnifRand(1., 2., size(sample, dim, IK))
1201rweight
1202+1.69031453, +1.28638148, +1.11450469, +1.47378707, +1.89527369, +1.28675365, +1.86891818
1203cov(:,:,0) = getCov(sample, dim, rweight)
1204cov(:,:,0) ! reference
1205+0.387720130E-1
1206do isam = 1, nsam
1207 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1208 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1209end do
1210call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1211covMerged
1212+0.387720168E-1
1213call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1214covMerged
1215+0.387720168E-1
1216call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1217cov(:,:,2)
1218+0.387720168E-1
1219cov(:,:,0) ! reference
1220+0.387720130E-1
1221
1222
1223dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1224do isam = 2, nsam
1225 lb(isam) = ub(isam - 1) + 1
1226 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1227end do
1228lb
1229+1, +5
1230ub
1231+4, +11
1232ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1233call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1234call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1235call setResized(covMerged, [ndim, ndim])
1236sample = getUnifRand(-1., +1., ndim, ub(nsam))
1237sample
1238+0.813830495, -0.697535753, +0.531675339, -0.663321972, -0.224798441, -0.639760375, -0.740989566, -0.628053784, +0.370065093, -0.474365711, +0.643649459
1239+0.641063213, +0.108627081E-1, +0.620541096, +0.226497650E-1, +0.727605820E-2, +0.940457702, +0.477403045, -0.963360071, +0.332710505, -0.612947583, +0.322023869
1240+0.463461876, -0.269112587E-1, +0.536152124, +0.935295582, -0.501708150, +0.165737152, -0.618223429, -0.639381528, -0.693918228, +0.472055316, +0.761637092
1241+0.976543784, -0.595544338, -0.214109421E-1, -0.295390129, +0.558245182E-2, -0.423438191, -0.586519599, -0.966666698, +0.111133933, -0.193339229, -0.206089258
1242rweight = getUnifRand(1., 2., size(sample, dim, IK))
1243rweight
1244+1.33977342, +1.56714368, +1.69980264, +1.48273349, +1.15075803, +1.30249929, +1.80330729, +1.73553944, +1.04150343, +1.63476944, +1.40844989
1245cov(:,:,0) = getCov(sample, dim, rweight)
1246cov(:,:,0) ! reference
1247+0.343572348, +0.133891404, +0.120037317, +0.212991893
1248+0.133891404, +0.306821913, +0.724486113E-1, +0.126839072
1249+0.120037317, +0.724486113E-1, +0.331439912, +0.106509455
1250+0.212991893, +0.126839072, +0.106509455, +0.225562170
1251do isam = 1, nsam
1252 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1253 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1254end do
1255call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1256covMerged
1257+0.343572110, +0.133891329, +0.120037332, +0.212991804
1258+0.700649232E-44, +0.306821913, +0.724486262E-1, +0.126839086
1259+0.00000000, -0.147216981E-1, +0.331439883, +0.106509343
1260+0.00000000, -0.165110230, +0.132138208, +0.225562274
1261call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1262covMerged
1263+0.343572110, +0.133891329, +0.120037332, +0.212991804
1264+0.133891329, +0.306821913, +0.724486262E-1, +0.126839086
1265+0.120037332, +0.724486262E-1, +0.331439883, +0.106509343
1266+0.212991804, +0.126839086, +0.106509343, +0.225562274
1267call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1268cov(:,:,2)
1269+0.343572110, +0.583240241E-1, +0.112679943, +0.106323555
1270+0.133891329, +0.306821913, +0.396544114E-1, +0.807652250E-1
1271+0.120037332, +0.724486262E-1, +0.331439883, +0.628527701E-1
1272+0.212991804, +0.126839086, +0.106509343, +0.225562274
1273cov(:,:,0) ! reference
1274+0.343572348, +0.133891404, +0.120037317, +0.212991893
1275+0.133891404, +0.306821913, +0.724486113E-1, +0.126839072
1276+0.120037317, +0.724486113E-1, +0.331439912, +0.106509455
1277+0.212991893, +0.126839072, +0.106509455, +0.225562170
1278
1279
1280dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1281do isam = 2, nsam
1282 lb(isam) = ub(isam - 1) + 1
1283 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1284end do
1285lb
1286+1, +5
1287ub
1288+4, +7
1289ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1290call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1291call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1292call setResized(covMerged, [ndim, ndim])
1293sample = getUnifRand(-1., +1., ndim, ub(nsam))
1294sample
1295+0.582108974, +0.301891685, -0.126112819, +0.113933206, +0.951265812, -0.782470226, +0.366373420
1296rweight = getUnifRand(1., 2., size(sample, dim, IK))
1297rweight
1298+1.24149013, +1.79682541, +1.76411724, +1.70669603, +1.27832305, +1.02916932, +1.75521207
1299cov(:,:,0) = getCov(sample, dim, rweight)
1300cov(:,:,0) ! reference
1301+0.204382941
1302do isam = 1, nsam
1303 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1304 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1305end do
1306call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1307covMerged
1308+0.204382926
1309call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1310covMerged
1311+0.204382926
1312call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1313cov(:,:,2)
1314+0.204382926
1315cov(:,:,0) ! reference
1316+0.204382941
1317
1318
Test:
test_pm_sampleCov


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 at Austin

Definition at line 7197 of file pm_sampleCov.F90.


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