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, +7
14ub
15+6, +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.437240243, -0.509490132, +0.282073259, +0.287244081, +0.876468897, +0.808184624, -0.120214701, +0.498818636
23-0.461711287, +0.660466790, -0.166719437, +0.446940184, -0.578876972, +0.581213236, +0.968951106, -0.628086090
24cov(:,:,0) = getCov(sample, dim)
25cov(:,:,0) ! reference
26+0.184757173, -0.148447484
27-0.148447484, +0.349839181
28do isam = 1, nsam
29 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
30 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
31end do
32call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
33covMerged
34+0.184757188, -0.148447484
35+0.307599026E-40, +0.349839121
36call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
37covMerged
38+0.184757188, -0.148447484
39-0.148447484, +0.349839121
40call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
41cov(:,:,2)
42+0.184757188, -0.148447484
43-0.247154817, +0.349839121
44cov(:,:,0) ! reference
45+0.184757173, -0.148447484
46-0.148447484, +0.349839181
47
48
49dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
50do isam = 2, nsam
51 lb(isam) = ub(isam - 1) + 1
52 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
53end do
54lb
55+1, +8
56ub
57+7, +11
58ndim = getUnifRand(1, minval(ub - lb + 1, 1))
59call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
60call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
61call setResized(covMerged, [ndim, ndim])
62sample = getUnifRand(-1., +1., ndim, ub(nsam))
63sample
64-0.802463293E-1, +0.328312397, -0.492934108, +0.639440894, +0.472739220, +0.960628271, -0.764507771, -0.748263717, +0.865222931, +0.976664305, +0.926007509
65cov(:,:,0) = getCov(sample, dim)
66cov(:,:,0) ! reference
67+0.432055891
68do isam = 1, nsam
69 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
70 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
71end do
72call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
73covMerged
74+0.432055861
75call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
76covMerged
77+0.432055861
78call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
79cov(:,:,2)
80+0.432055861
81cov(:,:,0) ! reference
82+0.432055891
83
84
85dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
86do isam = 2, nsam
87 lb(isam) = ub(isam - 1) + 1
88 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
89end do
90lb
91+1, +6
92ub
93+5, +11
94ndim = getUnifRand(1, minval(ub - lb + 1, 1))
95call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
96call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
97call setResized(covMerged, [ndim, ndim])
98sample = getUnifRand(-1., +1., ndim, ub(nsam))
99sample
100-0.786222816, +0.398960471, -0.736645579, -0.108120918, +0.429905772, +0.605038047, -0.997594476, +0.561249256E-2, -0.821478009, -0.784788966, +0.800313354
101cov(:,:,0) = getCov(sample, dim)
102cov(:,:,0) ! reference
103+0.404287428
104do isam = 1, nsam
105 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
106 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
107end do
108call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
109covMerged
110+0.404287457
111call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
112covMerged
113+0.404287457
114call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
115cov(:,:,2)
116+0.404287457
117cov(:,:,0) ! reference
118+0.404287428
119
120
121dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
122do isam = 2, nsam
123 lb(isam) = ub(isam - 1) + 1
124 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
125end do
126lb
127+1, +7
128ub
129+6, +12
130ndim = getUnifRand(1, minval(ub - lb + 1, 1))
131call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
132call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
133call setResized(covMerged, [ndim, ndim])
134sample = getUnifRand(-1., +1., ndim, ub(nsam))
135sample
136+0.856189251, -0.841922164, -0.602826476, +0.140571594E-2, -0.660218477, -0.180970430, +0.657186985, +0.422808647, -0.627171636, +0.997480154E-1, -0.671745300, +0.604941130
137-0.257167578, -0.403482080, +0.588209867, -0.273119211, -0.470341086, -0.605908990, +0.726447582, +0.508324385, +0.312461972, +0.853343844, +0.980333686, -0.209288597
138-0.992402196, -0.260337234, +0.833070993, +0.593615770, -0.399094224, +0.394656539, -0.451255083, +0.289817572, +0.766197801, +0.301555872, +0.276217937, -0.868330956
139-0.505458951, -0.670212388, -0.757353663, +0.949979186, -0.744159222E-1, +0.915979505, -0.393143177, +0.975389481E-1, +0.299537539, -0.148490667E-1, +0.194208622, -0.990006328
140cov(:,:,0) = getCov(sample, dim)
141cov(:,:,0) ! reference
142+0.335920572, +0.987617206E-2, -0.186426491, -0.614627227E-1
143+0.987617206E-2, +0.299959540, +0.111246869, -0.369814336E-1
144-0.186426491, +0.111246869, +0.349206060, +0.198542640
145-0.614627227E-1, -0.369814336E-1, +0.198542640, +0.351849973
146do isam = 1, nsam
147 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
148 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
149end do
150call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
151covMerged
152+0.335920483, +0.987614691E-2, -0.186426342, -0.614626482E-1
153+0.00000000, +0.299959481, +0.111246891, -0.369814187E-1
154+0.00000000, +0.00000000, +0.349205971, +0.198542565
155+0.00000000, +0.00000000, +0.00000000, +0.351849973
156call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
157covMerged
158+0.335920483, +0.987614691E-2, -0.186426342, -0.614626482E-1
159+0.987614691E-2, +0.299959481, +0.111246891, -0.369814187E-1
160-0.186426342, +0.111246891, +0.349205971, +0.198542565
161-0.614626482E-1, -0.369814187E-1, +0.198542565, +0.351849973
162call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
163cov(:,:,2)
164+0.335920483, +0.987614691E-2, -0.186426342, -0.614626482E-1
165-0.767387748E-1, +0.299959481, +0.111246891, -0.369814187E-1
166-0.226792455, +0.100445889, +0.349205971, +0.198542565
167-0.178418741, +0.111487411, +0.231040940, +0.351849973
168cov(:,:,0) ! reference
169+0.335920572, +0.987617206E-2, -0.186426491, -0.614627227E-1
170+0.987617206E-2, +0.299959540, +0.111246869, -0.369814336E-1
171-0.186426491, +0.111246869, +0.349206060, +0.198542640
172-0.614627227E-1, -0.369814336E-1, +0.198542640, +0.351849973
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, +8
182ub
183+7, +10
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.565130711E-1, +0.604866743E-1, +0.764050245, -0.767296433, -0.171962023, -0.552497029, -0.634489059E-1, -0.398125887, -0.557760119, -0.535832763
191+0.421899557, +0.646954656, -0.734021664, +0.600680470, -0.283736110, -0.728355527, -0.996326685, +0.172049880, +0.669707537, -0.255620480E-1
192cov(:,:,0) = getCov(sample, dim)
193cov(:,:,0) ! reference
194+0.180583522, -0.928269327E-1
195-0.928269327E-1, +0.357200235
196do isam = 1, nsam
197 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
198 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
199end do
200call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
201covMerged
202+0.180583537, -0.928269327E-1
203+0.987614691E-2, +0.357200205
204call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
205covMerged
206+0.180583537, -0.928269327E-1
207-0.928269327E-1, +0.357200205
208call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
209cov(:,:,2)
210+0.180583537, -0.928269327E-1
211-0.749734556E-2, +0.357200205
212cov(:,:,0) ! reference
213+0.180583522, -0.928269327E-1
214-0.928269327E-1, +0.357200235
215
216
217dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
218do isam = 2, nsam
219 lb(isam) = ub(isam - 1) + 1
220 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
221end do
222lb
223+1, +8
224ub
225+7, +9
226ndim = getUnifRand(1, minval(ub - lb + 1, 1))
227call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
228call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
229call setResized(covMerged, [ndim, ndim])
230sample = getUnifRand(-1., +1., ndim, ub(nsam))
231sample
232-0.336373925, -0.731604457, +0.813004136, -0.851869583E-1, -0.438123941E-1, -0.288273335, +0.697150588, -0.770933867, -0.663641095
233cov(:,:,0) = getCov(sample, dim)
234cov(:,:,0) ! reference
235+0.300180703
236do isam = 1, nsam
237 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
238 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
239end do
240call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
241covMerged
242+0.300180733
243call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
244covMerged
245+0.300180733
246call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
247cov(:,:,2)
248+0.300180733
249cov(:,:,0) ! reference
250+0.300180703
251
252
253dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
254do isam = 2, nsam
255 lb(isam) = ub(isam - 1) + 1
256 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
257end do
258lb
259+1, +7
260ub
261+6, +9
262ndim = getUnifRand(1, minval(ub - lb + 1, 1))
263call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
264call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
265call setResized(covMerged, [ndim, ndim])
266sample = getUnifRand(-1., +1., ndim, ub(nsam))
267sample
268-0.487814307, +0.335594654, +0.834491849, +0.350853205E-1, -0.405078053, +0.653033733, -0.649943352, -0.417102695, +0.738690853
269-0.232359409, -0.494017243, +0.143267870, +0.554116249, +0.960776806, -0.754070282, -0.216568351, -0.888762832, +0.314500332E-1
270cov(:,:,0) = getCov(sample, dim)
271cov(:,:,0) ! reference
272+0.303968996, -0.218831245E-1
273-0.218831245E-1, +0.318431526
274do isam = 1, nsam
275 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
276 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
277end do
278call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
279covMerged
280+0.303968996, -0.218831189E-1
281+0.307599026E-40, +0.318431526
282call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
283covMerged
284+0.303968996, -0.218831189E-1
285-0.218831189E-1, +0.318431526
286call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
287cov(:,:,2)
288+0.303968996, -0.218831189E-1
289+0.139052093, +0.318431526
290cov(:,:,0) ! reference
291+0.303968996, -0.218831245E-1
292-0.218831245E-1, +0.318431526
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, +5
302ub
303+4, +9
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.229614854, -0.957527518, -0.910255432, +0.563752651E-2, +0.842261910, +0.400862336, +0.592784524, +0.164602995, +0.342426538
311cov(:,:,0) = getCov(sample, dim)
312cov(:,:,0) ! reference
313+0.350778669
314do isam = 1, nsam
315 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
316 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
317end do
318call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
319covMerged
320+0.350778639
321call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
322covMerged
323+0.350778639
324call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
325cov(:,:,2)
326+0.350778639
327cov(:,:,0) ! reference
328+0.350778669
329
330
331dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
332do isam = 2, nsam
333 lb(isam) = ub(isam - 1) + 1
334 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
335end do
336lb
337+1, +7
338ub
339+6, +11
340ndim = getUnifRand(1, minval(ub - lb + 1, 1))
341call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
342call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
343call setResized(covMerged, [ndim, ndim])
344sample = getUnifRand(-1., +1., ndim, ub(nsam))
345sample
346-0.446139693, +0.376419902, +0.832920074E-1, -0.425267816, +0.877103925, -0.871086717, +0.119835019, -0.987196565, -0.783539057, -0.875635743, +0.865844965
347-0.208117366, +0.686280966, +0.408572674, -0.336873889, -0.712618589, +0.626461506, +0.531866789, -0.828013539, +0.856381536, -0.386170149, -0.860902429
348+0.183855057, -0.133424163, +0.732053280, -0.559751630, -0.524893522, +0.788945675, -0.643589973, +0.368193388E-1, +0.220641494, +0.555992723, -0.431228876
349cov(:,:,0) = getCov(sample, dim)
350cov(:,:,0) ! reference
351+0.435248315, -0.801191479E-1, -0.189553574
352-0.801191479E-1, +0.389329225, +0.109399885
353-0.189553574, +0.109399885, +0.250317305
354do isam = 1, nsam
355 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
356 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
357end do
358call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
359covMerged
360+0.435248315, -0.801191553E-1, -0.189553559
361+0.307599026E-40, +0.389329195, +0.109399930
362+0.00000000, -0.288273335, +0.250317305
363call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
364covMerged
365+0.435248315, -0.801191553E-1, -0.189553559
366-0.801191553E-1, +0.389329195, +0.109399930
367-0.189553559, +0.109399930, +0.250317305
368call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
369cov(:,:,2)
370+0.435248315, -0.801191553E-1, -0.189553559
371-0.850497708E-1, +0.389329195, +0.109399930
372-0.246677548, -0.126404762E-1, +0.250317305
373cov(:,:,0) ! reference
374+0.435248315, -0.801191479E-1, -0.189553574
375-0.801191479E-1, +0.389329225, +0.109399885
376-0.189553574, +0.109399885, +0.250317305
377
378
379dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
380do isam = 2, nsam
381 lb(isam) = ub(isam - 1) + 1
382 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
383end do
384lb
385+1, +8
386ub
387+7, +12
388ndim = getUnifRand(1, minval(ub - lb + 1, 1))
389call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
390call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
391call setResized(covMerged, [ndim, ndim])
392sample = getUnifRand(-1., +1., ndim, ub(nsam))
393sample
394-0.177401304, -0.170461059, +0.584437490, -0.170793414, +0.596334219, -0.847934127, +0.935771823, +0.772306442, +0.565444231E-1, -0.447962999, +0.312679410, -0.882995129E-2
395+0.817470074, -0.257990956, +0.795777440, +0.741047859E-1, -0.100444317, +0.946513772, +0.974837542, -0.337901115, -0.805449009, -0.870478392, +0.421383977, +0.848143101E-1
396-0.362956166, +0.127947330E-1, +0.999591231, -0.102045536E-1, +0.814240456, -0.927634478, -0.738519788, -0.198680997, +0.330262542, +0.416713595, +0.683423519, -0.421222091
397-0.356921077, +0.460367441, +0.637933493, -0.863015056, +0.359622598, +0.225723147, -0.813788176E-1, +0.137403965, +0.113751769, -0.439707041, +0.830465674, +0.365170717
398cov(:,:,0) = getCov(sample, dim)
399cov(:,:,0) ! reference
400+0.259015411, +0.339356363E-1, +0.878078789E-1, +0.709858090E-1
401+0.339356363E-1, +0.390183687, -0.131768554, +0.462670922E-1
402+0.878078789E-1, -0.131768554, +0.344750345, +0.914173350E-1
403+0.709858090E-1, +0.462670922E-1, +0.914173350E-1, +0.213774800
404do isam = 1, nsam
405 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
406 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
407end do
408call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
409covMerged
410+0.259015411, +0.339356288E-1, +0.878078863E-1, +0.709858313E-1
411-0.801191553E-1, +0.390183806, -0.131768554, +0.462670811E-1
412-0.189553559, +0.109399930, +0.344750375, +0.914173350E-1
413+0.00000000, -0.369814187E-1, +0.198542565, +0.213774815
414call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
415covMerged
416+0.259015411, +0.339356288E-1, +0.878078863E-1, +0.709858313E-1
417+0.339356288E-1, +0.390183806, -0.131768554, +0.462670811E-1
418+0.878078863E-1, -0.131768554, +0.344750375, +0.914173350E-1
419+0.709858313E-1, +0.462670811E-1, +0.914173350E-1, +0.213774815
420call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
421cov(:,:,2)
422+0.259015411, +0.339356288E-1, +0.878078863E-1, +0.709858313E-1
423+0.841821283E-1, +0.390183806, -0.131768554, +0.462670811E-1
424-0.430048220E-1, -0.129942792E-1, +0.344750375, +0.914173350E-1
425+0.856099650E-1, +0.185856864, +0.155060915E-1, +0.213774815
426cov(:,:,0) ! reference
427+0.259015411, +0.339356363E-1, +0.878078789E-1, +0.709858090E-1
428+0.339356363E-1, +0.390183687, -0.131768554, +0.462670922E-1
429+0.878078789E-1, -0.131768554, +0.344750345, +0.914173350E-1
430+0.709858090E-1, +0.462670922E-1, +0.914173350E-1, +0.213774800
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, +5
445ub
446+4, +7
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.969636440E-2, +0.305743933, -0.974210262, +0.312156320, +0.533624530, +0.416775107, -0.825224757
454-0.487743378, +0.640795588, -0.165283918, +0.638889194, -0.344024777, +0.954592109, -0.519934893E-1
455-0.403031707, -0.293452263, -0.655549407, +0.498218775, +0.715052962, -0.845716357, +0.642480016
456iweight = getUnifRand(1, 10, size(sample, dim, IK))
457iweight
458+3, +10, +8, +4, +5, +9, +2
459cov(:,:,0) = getCov(sample, dim, iweight)
460cov(:,:,0) ! reference
461+0.322554708, +0.155263871, +0.688224584E-1
462+0.155263871, +0.277615517, -0.109992489
463+0.688224584E-1, -0.109992489, +0.319208324
464do isam = 1, nsam
465 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
466 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
467end do
468call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
469covMerged
470+0.322554708, +0.155263871, +0.688224584E-1
471+0.307599026E-40, +0.277615428, -0.109992467
472+0.00000000, -0.301526040, +0.319208294
473call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
474covMerged
475+0.322554708, +0.155263871, +0.688224584E-1
476+0.155263871, +0.277615428, -0.109992467
477+0.688224584E-1, -0.109992467, +0.319208294
478call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
479cov(:,:,2)
480+0.322554708, +0.155263871, +0.688224584E-1
481+0.457286835E-1, +0.277615428, -0.109992467
482-0.940511078E-1, -0.462436438, +0.319208294
483cov(:,:,0) ! reference
484+0.322554708, +0.155263871, +0.688224584E-1
485+0.155263871, +0.277615517, -0.109992489
486+0.688224584E-1, -0.109992489, +0.319208324
487
488
489dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
490do isam = 2, nsam
491 lb(isam) = ub(isam - 1) + 1
492 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
493end do
494lb
495+1, +8
496ub
497+7, +11
498ndim = getUnifRand(1, minval(ub - lb + 1, 1))
499call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
500call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
501call setResized(covMerged, [ndim, ndim])
502sample = getUnifRand(-1., +1., ndim, ub(nsam))
503sample
504-0.646307588, +0.698085308, -0.109918118, +0.894554853, -0.115206599, +0.521912694, +0.516319036, +0.661313534E-2, +0.136246085, +0.310346484, +0.895836711
505iweight = getUnifRand(1, 10, size(sample, dim, IK))
506iweight
507+3, +2, +4, +1, +7, +4, +6, +5, +8, +8, +9
508cov(:,:,0) = getCov(sample, dim, iweight)
509cov(:,:,0) ! reference
510+0.167223081
511do isam = 1, nsam
512 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
513 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
514end do
515call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
516covMerged
517+0.167223111
518call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
519covMerged
520+0.167223111
521call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
522cov(:,:,2)
523+0.167223111
524cov(:,:,0) ! reference
525+0.167223081
526
527
528dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
529do isam = 2, nsam
530 lb(isam) = ub(isam - 1) + 1
531 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
532end do
533lb
534+1, +3
535ub
536+2, +5
537ndim = getUnifRand(1, minval(ub - lb + 1, 1))
538call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
539call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
540call setResized(covMerged, [ndim, ndim])
541sample = getUnifRand(-1., +1., ndim, ub(nsam))
542sample
543+0.873802304, +0.986678004, +0.835249424, -0.175638080, -0.922498345
544+0.361487627, +0.491217613, +0.301380634, +0.241978765, +0.506595731
545iweight = getUnifRand(1, 10, size(sample, dim, IK))
546iweight
547+2, +2, +10, +4, +8
548cov(:,:,0) = getCov(sample, dim, iweight)
549cov(:,:,0) ! reference
550+0.644975901, -0.494365655E-1
551-0.494365655E-1, +0.111880973E-1
552do isam = 1, nsam
553 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
554 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
555end do
556call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
557covMerged
558+0.644975841, -0.494365580E-1
559+0.307599026E-40, +0.111880973E-1
560call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
561covMerged
562+0.644975841, -0.494365580E-1
563-0.494365580E-1, +0.111880973E-1
564call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
565cov(:,:,2)
566+0.644975841, -0.494365580E-1
567-0.677264184E-1, +0.111880973E-1
568cov(:,:,0) ! reference
569+0.644975901, -0.494365655E-1
570-0.494365655E-1, +0.111880973E-1
571
572
573dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
574do isam = 2, nsam
575 lb(isam) = ub(isam - 1) + 1
576 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
577end do
578lb
579+1, +6
580ub
581+5, +12
582ndim = getUnifRand(1, minval(ub - lb + 1, 1))
583call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
584call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
585call setResized(covMerged, [ndim, ndim])
586sample = getUnifRand(-1., +1., ndim, ub(nsam))
587sample
588+0.397652388E-1, -0.259557009, -0.538977504, +0.320966244E-1, +0.780536413, +0.433959723, +0.307669759, -0.715632915, +0.822037697, +0.969168305, +0.443890929, -0.176961303
589iweight = getUnifRand(1, 10, size(sample, dim, IK))
590iweight
591+3, +7, +7, +8, +1, +8, +3, +10, +6, +3, +4, +1
592cov(:,:,0) = getCov(sample, dim, iweight)
593cov(:,:,0) ! reference
594+0.289318621
595do isam = 1, nsam
596 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
597 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
598end do
599call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
600covMerged
601+0.289318651
602call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
603covMerged
604+0.289318651
605call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
606cov(:,:,2)
607+0.289318651
608cov(:,:,0) ! reference
609+0.289318621
610
611
612dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
613do isam = 2, nsam
614 lb(isam) = ub(isam - 1) + 1
615 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
616end do
617lb
618+1, +3
619ub
620+2, +7
621ndim = getUnifRand(1, minval(ub - lb + 1, 1))
622call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
623call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
624call setResized(covMerged, [ndim, ndim])
625sample = getUnifRand(-1., +1., ndim, ub(nsam))
626sample
627+0.471642971, -0.430481195, -0.156780124, -0.595144033, -0.633973241, -0.350999475, +0.400300384
628-0.967359185, -0.380933285E-1, -0.540009975, +0.150626898, +0.860162973E-1, +0.920194387E-1, -0.232621789
629iweight = getUnifRand(1, 10, size(sample, dim, IK))
630iweight
631+7, +9, +5, +3, +8, +9, +7
632cov(:,:,0) = getCov(sample, dim, iweight)
633cov(:,:,0) ! reference
634+0.173456475, -0.124214932
635-0.124214932, +0.140285343
636do isam = 1, nsam
637 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
638 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
639end do
640call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
641covMerged
642+0.173456430, -0.124214955
643+0.307599026E-40, +0.140285373
644call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
645covMerged
646+0.173456430, -0.124214955
647-0.124214955, +0.140285373
648call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
649cov(:,:,2)
650+0.173456430, -0.124214955
651-0.566969700E-1, +0.140285373
652cov(:,:,0) ! reference
653+0.173456475, -0.124214932
654-0.124214932, +0.140285343
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, +7
664ub
665+6, +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.453070283, +0.220144272, -0.394856095, -0.229476452, -0.156873703, +0.703125119, -0.996817589, -0.388846636, +0.624248624, +0.717509151, -0.415325999
673-0.296703219, -0.643605947, -0.183574677, +0.541148305, -0.802192450, +0.700549603, -0.670517325, +0.420852542, +0.674514771E-1, -0.244485974, +0.663476706
674iweight = getUnifRand(1, 10, size(sample, dim, IK))
675iweight
676+8, +9, +1, +6, +6, +2, +3, +3, +8, +1, +5
677cov(:,:,0) = getCov(sample, dim, iweight)
678cov(:,:,0) ! reference
679+0.218925819, +0.151605839E-1
680+0.151605839E-1, +0.281323403
681do isam = 1, nsam
682 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
683 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
684end do
685call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
686covMerged
687+0.218925774, +0.151605960E-1
688-0.124214955, +0.281323373
689call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
690covMerged
691+0.218925774, +0.151605960E-1
692+0.151605960E-1, +0.281323373
693call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
694cov(:,:,2)
695+0.218925774, +0.151605960E-1
696+0.186300278E-1, +0.281323373
697cov(:,:,0) ! reference
698+0.218925819, +0.151605839E-1
699+0.151605839E-1, +0.281323403
700
701
702dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
703do isam = 2, nsam
704 lb(isam) = ub(isam - 1) + 1
705 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
706end do
707lb
708+1, +6
709ub
710+5, +9
711ndim = getUnifRand(1, minval(ub - lb + 1, 1))
712call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
713call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
714call setResized(covMerged, [ndim, ndim])
715sample = getUnifRand(-1., +1., ndim, ub(nsam))
716sample
717+0.813121796, -0.871680737, -0.403814435, +0.733884811, -0.839575648, +0.579057455, +0.495270014, +0.595306754, -0.375130177E-1
718+0.891590714, -0.694035053, -0.260246515, -0.332813263E-1, -0.380794048, +0.445414901, +0.301734328, +0.470428467, -0.493112326
719-0.684803247, +0.108013272, -0.448408961, +0.942111254, +0.321065426, +0.743355513, +0.965480208, +0.108347058, +0.844238758
720+0.311780453, +0.358275175, +0.182612538, +0.582923055, -0.117170811E-2, -0.508075953, -0.279668689, -0.155853152, +0.285991669
721iweight = getUnifRand(1, 10, size(sample, dim, IK))
722iweight
723+6, +2, +1, +4, +6, +2, +9, +6, +7
724cov(:,:,0) = getCov(sample, dim, iweight)
725cov(:,:,0) ! reference
726+0.347090065, +0.237425908, -0.206882115E-1, -0.118817948E-1
727+0.237425908, +0.244165689, -0.131911606, -0.430569761E-1
728-0.206882115E-1, -0.131911606, +0.329096228, -0.417070612E-1
729-0.118817948E-1, -0.430569761E-1, -0.417070612E-1, +0.932646021E-1
730do isam = 1, nsam
731 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
732 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
733end do
734call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
735covMerged
736+0.347090095, +0.237425894, -0.206881836E-1, -0.118818153E-1
737+0.151605960E-1, +0.244165689, -0.131911442, -0.430569984E-1
738+0.00000000, -0.131768554, +0.329096049, -0.417070016E-1
739+0.00000000, +0.462670811E-1, +0.914173350E-1, +0.932645798E-1
740call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
741covMerged
742+0.347090095, +0.237425894, -0.206881836E-1, -0.118818153E-1
743+0.237425894, +0.244165689, -0.131911442, -0.430569984E-1
744-0.206881836E-1, -0.131911442, +0.329096049, -0.417070016E-1
745-0.118818153E-1, -0.430569984E-1, -0.417070016E-1, +0.932645798E-1
746call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
747cov(:,:,2)
748+0.347090095, +0.237425894, -0.206881836E-1, -0.118818153E-1
749+0.106810600, +0.244165689, -0.131911442, -0.430569984E-1
750-0.372382030E-1, -0.583472960E-1, +0.329096049, -0.417070016E-1
751-0.645684823E-1, -0.972173214E-1, +0.513126468E-2, +0.932645798E-1
752cov(:,:,0) ! reference
753+0.347090065, +0.237425908, -0.206882115E-1, -0.118817948E-1
754+0.237425908, +0.244165689, -0.131911606, -0.430569761E-1
755-0.206882115E-1, -0.131911606, +0.329096228, -0.417070612E-1
756-0.118817948E-1, -0.430569761E-1, -0.417070612E-1, +0.932646021E-1
757
758
759dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
760do isam = 2, nsam
761 lb(isam) = ub(isam - 1) + 1
762 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
763end do
764lb
765+1, +8
766ub
767+7, +9
768ndim = getUnifRand(1, minval(ub - lb + 1, 1))
769call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
770call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
771call setResized(covMerged, [ndim, ndim])
772sample = getUnifRand(-1., +1., ndim, ub(nsam))
773sample
774-0.805975676, -0.981281757, -0.581297874, -0.953873515, -0.902015209, +0.805794239, +0.663238406, +0.580360770, -0.507354140
775+0.886942506, +0.310795426, -0.892228484, -0.785974741, +0.615517974, +0.161004782, +0.538974643, +0.581898332, +0.507910609
776iweight = getUnifRand(1, 10, size(sample, dim, IK))
777iweight
778+1, +2, +8, +7, +2, +7, +2, +4, +4
779cov(:,:,0) = getCov(sample, dim, iweight)
780cov(:,:,0) ! reference
781+0.515104175, +0.219625473
782+0.219625473, +0.412931085
783do isam = 1, nsam
784 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
785 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
786end do
787call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
788covMerged
789+0.515104115, +0.219625413
790+0.237425894, +0.412931263
791call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
792covMerged
793+0.515104115, +0.219625413
794+0.219625413, +0.412931263
795call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
796cov(:,:,2)
797+0.515104115, +0.219625413
798+0.201193877E-1, +0.412931263
799cov(:,:,0) ! reference
800+0.515104175, +0.219625473
801+0.219625473, +0.412931085
802
803
804dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
805do isam = 2, nsam
806 lb(isam) = ub(isam - 1) + 1
807 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
808end do
809lb
810+1, +6
811ub
812+5, +10
813ndim = getUnifRand(1, minval(ub - lb + 1, 1))
814call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
815call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
816call setResized(covMerged, [ndim, ndim])
817sample = getUnifRand(-1., +1., ndim, ub(nsam))
818sample
819-0.561125278, -0.284564376, +0.125218987, -0.168395400, +0.877732277, -0.291383505, +0.238629818, +0.834685206, +0.994071007, +0.623575330
820+0.389709592, -0.178529859, -0.740751147, -0.123763084E-2, -0.568816662, -0.390716553, +0.512215972, -0.558581471, +0.576901078, -0.770619035
821+0.912296176, +0.962348700, -0.239195824, -0.153777003, +0.506079197, +0.838832855E-1, -0.987587094, +0.803172708, +0.613986969, +0.260746717
822-0.928229332, -0.141882420, -0.370680332, -0.627827644, +0.218896389, -0.537421703, -0.103929043E-1, +0.823149323, -0.109847426, -0.598375916
823iweight = getUnifRand(1, 10, size(sample, dim, IK))
824iweight
825+7, +7, +3, +4, +2, +1, +2, +9, +1, +6
826cov(:,:,0) = getCov(sample, dim, iweight)
827cov(:,:,0) ! reference
828+0.303562313, -0.158064619, -0.316820145E-1, +0.238386020
829-0.158064619, +0.200517938, +0.725471415E-2, -0.115787551
830-0.316820145E-1, +0.725471415E-2, +0.281115592, +0.607889742E-1
831+0.238386020, -0.115787551, +0.607889742E-1, +0.367995232
832do isam = 1, nsam
833 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
834 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
835end do
836call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
837covMerged
838+0.303562373, -0.158064634, -0.316820256E-1, +0.238386154
839+0.219625413, +0.200517878, +0.725473464E-2, -0.115787603
840+0.00000000, +0.965480208, +0.281115592, +0.607889481E-1
841+0.00000000, -0.279668689, -0.155853152, +0.367995232
842call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
843covMerged
844+0.303562373, -0.158064634, -0.316820256E-1, +0.238386154
845-0.158064634, +0.200517878, +0.725473464E-2, -0.115787603
846-0.316820256E-1, +0.725473464E-2, +0.281115592, +0.607889481E-1
847+0.238386154, -0.115787603, +0.607889481E-1, +0.367995232
848call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
849cov(:,:,2)
850+0.303562373, -0.158064634, -0.316820256E-1, +0.238386154
851-0.328560434E-1, +0.200517878, +0.725473464E-2, -0.115787603
852+0.116075113, -0.136575684, +0.281115592, +0.607889481E-1
853+0.101363383, +0.862368476E-2, +0.193473294, +0.367995232
854cov(:,:,0) ! reference
855+0.303562313, -0.158064619, -0.316820145E-1, +0.238386020
856-0.158064619, +0.200517938, +0.725471415E-2, -0.115787551
857-0.316820145E-1, +0.725471415E-2, +0.281115592, +0.607889742E-1
858+0.238386020, -0.115787551, +0.607889742E-1, +0.367995232
859
860
861dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
862do isam = 2, nsam
863 lb(isam) = ub(isam - 1) + 1
864 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
865end do
866lb
867+1, +8
868ub
869+7, +11
870ndim = getUnifRand(1, minval(ub - lb + 1, 1))
871call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
872call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
873call setResized(covMerged, [ndim, ndim])
874sample = getUnifRand(-1., +1., ndim, ub(nsam))
875sample
876+0.951781154, +0.972033024, +0.783034205, -0.831316352, -0.881828308, -0.589258313, -0.606879115, -0.801753283, +0.247774124, +0.429144502, -0.692235827
877iweight = getUnifRand(1, 10, size(sample, dim, IK))
878iweight
879+6, +5, +2, +3, +10, +7, +9, +7, +10, +8, +1
880cov(:,:,0) = getCov(sample, dim, iweight)
881cov(:,:,0) ! reference
882+0.482558310
883do isam = 1, nsam
884 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
885 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
886end do
887call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
888covMerged
889+0.482558399
890call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
891covMerged
892+0.482558399
893call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
894cov(:,:,2)
895+0.482558399
896cov(:,:,0) ! reference
897+0.482558310
898
899
900!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
901!Compute the biased merged covariance of a reliability weighted multivariate sample.
902!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
903
904
905dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
906do isam = 2, nsam
907 lb(isam) = ub(isam - 1) + 1
908 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
909end do
910lb
911+1, +3
912ub
913+2, +8
914ndim = getUnifRand(1, minval(ub - lb + 1, 1))
915call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
916call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
917call setResized(covMerged, [ndim, ndim])
918sample = getUnifRand(-1., +1., ndim, ub(nsam))
919sample
920-0.208812952, +0.172437429E-1, -0.536489725, -0.195361137, +0.886425614, -0.363358021, +0.619984269, +0.300908208
921-0.966960430, +0.943702579, -0.338876247, -0.536715031, +0.778634548E-1, +0.361042261, -0.916269183, -0.496102929
922rweight = getUnifRand(1., 2., size(sample, dim, IK))
923rweight
924+1.02934909, +1.17991400, +1.48440814, +1.33803344, +1.04135633, +1.32981753, +1.53234267, +1.07241571
925cov(:,:,0) = getCov(sample, dim, rweight)
926cov(:,:,0) ! reference
927+0.218063891, -0.381696895E-1
928-0.381696895E-1, +0.368289858
929do isam = 1, nsam
930 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
931 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
932end do
933call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
934covMerged
935+0.218063816, -0.381697342E-1
936+0.307599026E-40, +0.368289679
937call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
938covMerged
939+0.218063816, -0.381697342E-1
940-0.381697342E-1, +0.368289679
941call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
942cov(:,:,2)
943+0.218063816, -0.647561923E-1
944-0.381697342E-1, +0.368289679
945cov(:,:,0) ! reference
946+0.218063891, -0.381696895E-1
947-0.381696895E-1, +0.368289858
948
949
950dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
951do isam = 2, nsam
952 lb(isam) = ub(isam - 1) + 1
953 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
954end do
955lb
956+1, +6
957ub
958+5, +10
959ndim = getUnifRand(1, minval(ub - lb + 1, 1))
960call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
961call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
962call setResized(covMerged, [ndim, ndim])
963sample = getUnifRand(-1., +1., ndim, ub(nsam))
964sample
965+0.915432692, -0.653598905, -0.838334560E-1, -0.743202448, +0.295577884, -0.487685323, -0.239836812, -0.215451598, +0.351986527, -0.319292068
966rweight = getUnifRand(1., 2., size(sample, dim, IK))
967rweight
968+1.72926068, +1.93589008, +1.66682434, +1.99141932, +1.31718850, +1.56028974, +1.65035868, +1.44783390, +1.06330442, +1.88710332
969cov(:,:,0) = getCov(sample, dim, rweight)
970cov(:,:,0) ! reference
971+0.242539853
972do isam = 1, nsam
973 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
974 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
975end do
976call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
977covMerged
978+0.242539778
979call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
980covMerged
981+0.242539778
982call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
983cov(:,:,2)
984+0.242539778
985cov(:,:,0) ! reference
986+0.242539853
987
988
989dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
990do isam = 2, nsam
991 lb(isam) = ub(isam - 1) + 1
992 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
993end do
994lb
995+1, +6
996ub
997+5, +11
998ndim = getUnifRand(1, minval(ub - lb + 1, 1))
999call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1000call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1001call setResized(covMerged, [ndim, ndim])
1002sample = getUnifRand(-1., +1., ndim, ub(nsam))
1003sample
1004+0.817428350, +0.223605990, -0.603304029, -0.190736890, -0.217534661, -0.581344962, +0.888407230E-1, +0.106315613, -0.874715805, +0.599988341, -0.968065143
1005rweight = getUnifRand(1., 2., size(sample, dim, IK))
1006rweight
1007+1.92081749, +1.37140512, +1.20661354, +1.55333233, +1.64777613, +1.64123273, +1.91226554, +1.91070175, +1.30830121, +1.18622780, +1.07210815
1008cov(:,:,0) = getCov(sample, dim, rweight)
1009cov(:,:,0) ! reference
1010+0.287185639
1011do isam = 1, nsam
1012 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1013 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1014end do
1015call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1016covMerged
1017+0.287185580
1018call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1019covMerged
1020+0.287185580
1021call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1022cov(:,:,2)
1023+0.287185580
1024cov(:,:,0) ! reference
1025+0.287185639
1026
1027
1028dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1029do isam = 2, nsam
1030 lb(isam) = ub(isam - 1) + 1
1031 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1032end do
1033lb
1034+1, +5
1035ub
1036+4, +9
1037ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1038call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1039call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1040call setResized(covMerged, [ndim, ndim])
1041sample = getUnifRand(-1., +1., ndim, ub(nsam))
1042sample
1043+0.729969978, +0.325574040, +0.633107424, +0.875135303, +0.231574774E-1, +0.510561466E-2, +0.188273549, -0.825449944, +0.657613039
1044rweight = getUnifRand(1., 2., size(sample, dim, IK))
1045rweight
1046+1.38936448, +1.35191381, +1.46252060, +1.14231896, +1.25327516, +1.87033701, +1.43848753, +1.69921613, +1.81293392
1047cov(:,:,0) = getCov(sample, dim, rweight)
1048cov(:,:,0) ! reference
1049+0.256056577
1050do isam = 1, nsam
1051 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1052 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1053end do
1054call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1055covMerged
1056+0.256056607
1057call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1058covMerged
1059+0.256056607
1060call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1061cov(:,:,2)
1062+0.256056607
1063cov(:,:,0) ! reference
1064+0.256056577
1065
1066
1067dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1068do isam = 2, nsam
1069 lb(isam) = ub(isam - 1) + 1
1070 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1071end do
1072lb
1073+1, +3
1074ub
1075+2, +9
1076ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1077call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1078call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1079call setResized(covMerged, [ndim, ndim])
1080sample = getUnifRand(-1., +1., ndim, ub(nsam))
1081sample
1082+0.310665727, -0.784798861E-1, -0.408750772E-1, +0.236700773E-1, -0.589460135, +0.605108380, -0.586218119, -0.736391544, +0.923646688E-1
1083rweight = getUnifRand(1., 2., size(sample, dim, IK))
1084rweight
1085+1.51354599, +1.73482478, +1.80212164, +1.13209867, +1.94563532, +1.34622145, +1.78995442, +1.61077118, +1.24866509
1086cov(:,:,0) = getCov(sample, dim, rweight)
1087cov(:,:,0) ! reference
1088+0.176771611
1089do isam = 1, nsam
1090 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1091 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1092end do
1093call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1094covMerged
1095+0.176771656
1096call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1097covMerged
1098+0.176771656
1099call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1100cov(:,:,2)
1101+0.176771656
1102cov(:,:,0) ! reference
1103+0.176771611
1104
1105
1106dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1107do isam = 2, nsam
1108 lb(isam) = ub(isam - 1) + 1
1109 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1110end do
1111lb
1112+1, +7
1113ub
1114+6, +10
1115ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1116call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1117call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1118call setResized(covMerged, [ndim, ndim])
1119sample = getUnifRand(-1., +1., ndim, ub(nsam))
1120sample
1121-0.408083797, -0.320973396, +0.925886035, -0.809636116E-1, +0.596170425, -0.375377297, +0.506871700, -0.556127310, -0.181728721, -0.940762162
1122+0.178442121, -0.687025309, -0.238163710, -0.775015831, -0.543785214, -0.252047300, -0.203243732, -0.916440487E-1, -0.658980966, -0.124249458E-1
1123+0.566037536, -0.245252132, +0.358351231, +0.737571239, +0.125814557, +0.725347996E-1, -0.446675777, +0.603518128, -0.241560578, -0.698436499
1124rweight = getUnifRand(1., 2., size(sample, dim, IK))
1125rweight
1126+1.18939352, +1.53044820, +1.71274900, +1.76704192, +1.05411363, +1.15521646, +1.17673028, +1.05407143, +1.36943102, +1.37501192
1127cov(:,:,0) = getCov(sample, dim, rweight)
1128cov(:,:,0) ! reference
1129+0.314311415, -0.331808701E-1, +0.583453812E-1
1130-0.331808701E-1, +0.956496894E-1, -0.112399152E-1
1131+0.583453812E-1, -0.112399152E-1, +0.218885809
1132do isam = 1, nsam
1133 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1134 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1135end do
1136call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1137covMerged
1138+0.314311445, -0.331808589E-1, +0.583453923E-1
1139+0.307599026E-40, +0.956496224E-1, -0.112399347E-1
1140+0.00000000, +1.56028974, +0.218885839
1141call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1142covMerged
1143+0.314311445, -0.331808589E-1, +0.583453923E-1
1144-0.331808589E-1, +0.956496224E-1, -0.112399347E-1
1145+0.583453923E-1, -0.112399347E-1, +0.218885839
1146call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1147cov(:,:,2)
1148+0.314311445, -0.550740920E-1, -0.406017667E-2
1149-0.331808589E-1, +0.956496224E-1, -0.390417594E-2
1150+0.583453923E-1, -0.112399347E-1, +0.218885839
1151cov(:,:,0) ! reference
1152+0.314311415, -0.331808701E-1, +0.583453812E-1
1153-0.331808701E-1, +0.956496894E-1, -0.112399152E-1
1154+0.583453812E-1, -0.112399152E-1, +0.218885809
1155
1156
1157dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1158do isam = 2, nsam
1159 lb(isam) = ub(isam - 1) + 1
1160 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1161end do
1162lb
1163+1, +5
1164ub
1165+4, +11
1166ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1167call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1168call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1169call setResized(covMerged, [ndim, ndim])
1170sample = getUnifRand(-1., +1., ndim, ub(nsam))
1171sample
1172+0.341432214, +0.881614089, -0.399983644, +0.121986508, -0.629821777, -0.989700913, -0.313008666, +0.485669732, -0.426841259, -0.730681419E-1, -0.632340789
1173+0.270296335E-1, +0.960321307, +0.123806000, -0.832738996, +0.585676551, -0.332803130, +0.632914662, -0.709569335, -0.768916965, +0.896968365, +0.700377822
1174rweight = getUnifRand(1., 2., size(sample, dim, IK))
1175rweight
1176+1.36839724, +1.22569716, +1.53228807, +1.43287218, +1.19800758, +1.80557775, +1.30822086, +1.96927047, +1.65572071, +1.94618654, +1.35516381
1177cov(:,:,0) = getCov(sample, dim, rweight)
1178cov(:,:,0) ! reference
1179+0.277987570, +0.232558488E-2
1180+0.232558488E-2, +0.440771341
1181do isam = 1, nsam
1182 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1183 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1184end do
1185call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1186covMerged
1187+0.277987510, +0.232557347E-2
1188-0.331808589E-1, +0.440771312
1189call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1190covMerged
1191+0.277987510, +0.232557347E-2
1192+0.232557347E-2, +0.440771312
1193call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1194cov(:,:,2)
1195+0.277987510, -0.583272949E-1
1196+0.232557347E-2, +0.440771312
1197cov(:,:,0) ! reference
1198+0.277987570, +0.232558488E-2
1199+0.232558488E-2, +0.440771341
1200
1201
1202dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1203do isam = 2, nsam
1204 lb(isam) = ub(isam - 1) + 1
1205 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1206end do
1207lb
1208+1, +4
1209ub
1210+3, +9
1211ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1212call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1213call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1214call setResized(covMerged, [ndim, ndim])
1215sample = getUnifRand(-1., +1., ndim, ub(nsam))
1216sample
1217-0.703363776, +0.365755200, +0.485166788, +0.877479315, -0.228544474E-1, +0.341466427, +0.713805676, -0.425514817, -0.471305132
1218+0.161477685, -0.383538008E-1, -0.361167669, -0.686302781, +0.998305082, -0.851873279, -0.922629476, +0.188710928, +0.393612862
1219rweight = getUnifRand(1., 2., size(sample, dim, IK))
1220rweight
1221+1.54010212, +1.89975560, +1.66566694, +1.38662136, +1.46789765, +1.68640089, +1.33119249, +1.50009692, +1.93372250
1222cov(:,:,0) = getCov(sample, dim, rweight)
1223cov(:,:,0) ! reference
1224+0.269198686, -0.210975319
1225-0.210975319, +0.341622353
1226do isam = 1, nsam
1227 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1228 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1229end do
1230call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1231covMerged
1232+0.269198686, -0.210975319
1233+0.232557347E-2, +0.341622353
1234call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1235covMerged
1236+0.269198686, -0.210975319
1237-0.210975319, +0.341622353
1238call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1239cov(:,:,2)
1240+0.269198686, -0.276925981
1241-0.210975319, +0.341622353
1242cov(:,:,0) ! reference
1243+0.269198686, -0.210975319
1244-0.210975319, +0.341622353
1245
1246
1247dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1248do isam = 2, nsam
1249 lb(isam) = ub(isam - 1) + 1
1250 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1251end do
1252lb
1253+1, +6
1254ub
1255+5, +11
1256ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1257call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1258call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1259call setResized(covMerged, [ndim, ndim])
1260sample = getUnifRand(-1., +1., ndim, ub(nsam))
1261sample
1262+0.231803775, -0.405511141, +0.100138426, +0.110578656, -0.815388680, +0.529161692E-1, -0.383244634, +0.933436871, +0.768813372, +0.413139462, +0.917377234
1263-0.494561315, +0.677074194E-1, +0.495207191, -0.518746138, -0.309982181, +0.412847281, +0.829150081, +0.160709620E-1, +0.425965190, -0.868157625, -0.439377069
1264-0.226817369, +0.459925175, -0.234521627, +0.150696158, -0.294090867, +0.708761692, +0.540903807E-1, +0.898070693, -0.831350088, -0.343743086, -0.402327776E-1
1265+0.698737860, +0.154460907, -0.666509628, -0.666931868E-1, -0.802560449, +0.647636533, +0.121201873, -0.332906127, +0.910137892, +0.242470622, +0.903353572
1266rweight = getUnifRand(1., 2., size(sample, dim, IK))
1267rweight
1268+1.43745375, +1.66492128, +1.20702338, +1.92775202, +1.99444473, +1.73931384, +1.15162110, +1.35957360, +1.12056088, +1.65248513, +1.72152197
1269cov(:,:,0) = getCov(sample, dim, rweight)
1270cov(:,:,0) ! reference
1271+0.301746637, -0.442910679E-1, +0.871652924E-2, +0.173276588
1272-0.442910679E-1, +0.240171909, +0.523611270E-1, -0.561171537E-2
1273+0.871652924E-2, +0.523611270E-1, +0.213028520, -0.140174637E-1
1274+0.173276588, -0.561171537E-2, -0.140174637E-1, +0.325453520
1275do isam = 1, nsam
1276 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1277 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1278end do
1279call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1280covMerged
1281+0.301746666, -0.442910492E-1, +0.871652924E-2, +0.173276573
1282-0.210975319, +0.240171894, +0.523611009E-1, -0.561170280E-2
1283+0.00000000, -0.953873515, +0.213028461, -0.140174627E-1
1284+0.00000000, -0.785974741, +0.161004782, +0.325453460
1285call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1286covMerged
1287+0.301746666, -0.442910492E-1, +0.871652924E-2, +0.173276573
1288-0.442910492E-1, +0.240171894, +0.523611009E-1, -0.561170280E-2
1289+0.871652924E-2, +0.523611009E-1, +0.213028461, -0.140174627E-1
1290+0.173276573, -0.561170280E-2, -0.140174627E-1, +0.325453460
1291call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1292cov(:,:,2)
1293+0.301746666, -0.139328852, -0.310193896E-1, +0.237882007E-1
1294-0.442910492E-1, +0.240171894, +0.810756236E-1, -0.171762111E-2
1295+0.871652924E-2, +0.523611009E-1, +0.213028461, -0.121108644
1296+0.173276573, -0.561170280E-2, -0.140174627E-1, +0.325453460
1297cov(:,:,0) ! reference
1298+0.301746637, -0.442910679E-1, +0.871652924E-2, +0.173276588
1299-0.442910679E-1, +0.240171909, +0.523611270E-1, -0.561171537E-2
1300+0.871652924E-2, +0.523611270E-1, +0.213028520, -0.140174637E-1
1301+0.173276588, -0.561171537E-2, -0.140174637E-1, +0.325453520
1302
1303
1304dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1305do isam = 2, nsam
1306 lb(isam) = ub(isam - 1) + 1
1307 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1308end do
1309lb
1310+1, +3
1311ub
1312+2, +8
1313ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1314call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1315call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1316call setResized(covMerged, [ndim, ndim])
1317sample = getUnifRand(-1., +1., ndim, ub(nsam))
1318sample
1319-0.299010873, +0.176291347, -0.858146787, -0.699359655, -0.857713223E-1, -0.647481441, -0.197677374, +0.367065191
1320rweight = getUnifRand(1., 2., size(sample, dim, IK))
1321rweight
1322+1.16289353, +1.44778609, +1.33882809, +1.10557842, +1.73289824, +1.95168519, +1.53807294, +1.21041763
1323cov(:,:,0) = getCov(sample, dim, rweight)
1324cov(:,:,0) ! reference
1325+0.155770376
1326do isam = 1, nsam
1327 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1328 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1329end do
1330call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1331covMerged
1332+0.155770332
1333call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1334covMerged
1335+0.155770332
1336call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1337cov(:,:,2)
1338+0.155770332
1339cov(:,:,0) ! reference
1340+0.155770376
1341
1342
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 Austin

Definition at line 7197 of file pm_sampleCov.F90.


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