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, +6
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.887790442, +0.194713950, -0.644258261, +0.611379266, +0.394568563, -0.408532143
23+0.457165003, -0.532902479, -0.113757491, +0.269741297, -0.705339909E-1, -0.124546051
24cov(:,:,0) = getCov(sample, dim)
25cov(:,:,0) ! reference
26+0.293125987, +0.971965045E-1
27+0.971965045E-1, +0.994958282E-1
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.293125957, +0.971964225E-1
35+0.307122584E-40, +0.994958282E-1
36call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
37covMerged
38+0.293125957, +0.971964225E-1
39+0.971964225E-1, +0.994958282E-1
40call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
41cov(:,:,2)
42+0.293125957, +0.971964225E-1
43+0.108442809E-1, +0.994958282E-1
44cov(:,:,0) ! reference
45+0.293125987, +0.971965045E-1
46+0.971965045E-1, +0.994958282E-1
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, +4
56ub
57+3, +5
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.770795345E-2, +0.325453281E-2, +0.461736202, +0.205103278, -0.273051500
65cov(:,:,0) = getCov(sample, dim)
66cov(:,:,0) ! reference
67+0.599157028E-1
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.599157065E-1
75call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
76covMerged
77+0.599157065E-1
78call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
79cov(:,:,2)
80+0.599157065E-1
81cov(:,:,0) ! reference
82+0.599157028E-1
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, +8
92ub
93+7, +10
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.938814998, +0.831867456, +0.454279900, +0.947416425, +0.297788739, -0.765782595E-1, +0.273837209, +0.123464346, +0.500893593E-1, -0.198791265
101cov(:,:,0) = getCov(sample, dim)
102cov(:,:,0) ! reference
103+0.259277642
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.259277701
111call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
112covMerged
113+0.259277701
114call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
115cov(:,:,2)
116+0.259277701
117cov(:,:,0) ! reference
118+0.259277642
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, +4
128ub
129+3, +8
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.892275929, +0.385918498, -0.911497712, +0.930299640, -0.469763637, -0.489342690, -0.774952888, +0.938227415
137cov(:,:,0) = getCov(sample, dim)
138cov(:,:,0) ! reference
139+0.547055900
140do isam = 1, nsam
141 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
142 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
143end do
144call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
145covMerged
146+0.547055781
147call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
148covMerged
149+0.547055781
150call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
151cov(:,:,2)
152+0.547055781
153cov(:,:,0) ! reference
154+0.547055900
155
156
157dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
158do isam = 2, nsam
159 lb(isam) = ub(isam - 1) + 1
160 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
161end do
162lb
163+1, +7
164ub
165+6, +8
166ndim = getUnifRand(1, minval(ub - lb + 1, 1))
167call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
168call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
169call setResized(covMerged, [ndim, ndim])
170sample = getUnifRand(-1., +1., ndim, ub(nsam))
171sample
172-0.976064444, +0.517642498E-1, -0.891297221, -0.750822067, -0.120664597, +0.366348028, -0.726908088, +0.354174733
173+0.214452744E-1, -0.728033781E-1, -0.223361135, +0.867613792, +0.751181364, +0.429640412, -0.823048830, -0.716594815
174cov(:,:,0) = getCov(sample, dim)
175cov(:,:,0) ! reference
176+0.276160538, +0.162542704E-2
177+0.162542704E-2, +0.342667162
178do isam = 1, nsam
179 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
180 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
181end do
182call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
183covMerged
184+0.276160568, +0.162541121E-2
185+0.307122584E-40, +0.342667162
186call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
187covMerged
188+0.276160568, +0.162541121E-2
189+0.162541121E-2, +0.342667162
190call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
191cov(:,:,2)
192+0.276160568, +0.162541121E-2
193+0.287714023E-1, +0.342667162
194cov(:,:,0) ! reference
195+0.276160538, +0.162542704E-2
196+0.162542704E-2, +0.342667162
197
198
199dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
200do isam = 2, nsam
201 lb(isam) = ub(isam - 1) + 1
202 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
203end do
204lb
205+1, +8
206ub
207+7, +10
208ndim = getUnifRand(1, minval(ub - lb + 1, 1))
209call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
210call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
211call setResized(covMerged, [ndim, ndim])
212sample = getUnifRand(-1., +1., ndim, ub(nsam))
213sample
214-0.625729561E-2, -0.194524646, +0.496241212, -0.253749251, -0.371255517, -0.158991814E-1, -0.619093418, -0.135685563, +0.299021125, -0.513912082
215cov(:,:,0) = getCov(sample, dim)
216cov(:,:,0) ! reference
217+0.106886022
218do isam = 1, nsam
219 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
220 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
221end do
222call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
223covMerged
224+0.106886029
225call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
226covMerged
227+0.106886029
228call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
229cov(:,:,2)
230+0.106886029
231cov(:,:,0) ! reference
232+0.106886022
233
234
235dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
236do isam = 2, nsam
237 lb(isam) = ub(isam - 1) + 1
238 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
239end do
240lb
241+1, +5
242ub
243+4, +6
244ndim = getUnifRand(1, minval(ub - lb + 1, 1))
245call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
246call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
247call setResized(covMerged, [ndim, ndim])
248sample = getUnifRand(-1., +1., ndim, ub(nsam))
249sample
250-0.572700143, -0.344327211, +0.314110160, -0.148440957, -0.705322266, -0.657382011
251-0.343778849, +0.439948678, -0.135738850, -0.480763316, +0.951527119, +0.380480289E-2
252cov(:,:,0) = getCov(sample, dim)
253cov(:,:,0) ! reference
254+0.125333428, -0.743735284E-1
255-0.743735284E-1, +0.239196345
256do isam = 1, nsam
257 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
258 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
259end do
260call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
261covMerged
262+0.125333428, -0.743735284E-1
263+0.307122584E-40, +0.239196345
264call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
265covMerged
266+0.125333428, -0.743735284E-1
267-0.743735284E-1, +0.239196345
268call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
269cov(:,:,2)
270+0.125333428, -0.743735284E-1
271-0.113585126E-1, +0.239196345
272cov(:,:,0) ! reference
273+0.125333428, -0.743735284E-1
274-0.743735284E-1, +0.239196345
275
276
277dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
278do isam = 2, nsam
279 lb(isam) = ub(isam - 1) + 1
280 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
281end do
282lb
283+1, +6
284ub
285+5, +10
286ndim = getUnifRand(1, minval(ub - lb + 1, 1))
287call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
288call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
289call setResized(covMerged, [ndim, ndim])
290sample = getUnifRand(-1., +1., ndim, ub(nsam))
291sample
292-0.893119574, -0.196503043, -0.172145009, -0.490762353, -0.700281620, +0.604099751, -0.737175226, +0.232497334, -0.267077446, +0.507488251E-1
293-0.999871016, -0.641601324, +0.495658755, +0.360296845, +0.632520676, -0.282853007, +0.850572348, -0.753382087, -0.820460796, +0.730622649
294+0.742694616, -0.654029846E-1, +0.542109847, -0.466190577E-1, -0.997512817, +0.887774229, -0.436052084, -0.305809617, -0.918288589, -0.404472470
295cov(:,:,0) = getCov(sample, dim)
296cov(:,:,0) ! reference
297+0.197313070, -0.512966625E-1, +0.731636062E-1
298-0.512966625E-1, +0.474663645, -0.105630495
299+0.731636062E-1, -0.105630495, +0.382530123
300do isam = 1, nsam
301 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
302 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
303end do
304call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
305covMerged
306+0.197313115, -0.512965582E-1, +0.731634647E-1
307-0.743735284E-1, +0.474663615, -0.105630465
308+0.00000000, -0.158991814E-1, +0.382529974
309call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
310covMerged
311+0.197313115, -0.512965582E-1, +0.731634647E-1
312-0.512965582E-1, +0.474663615, -0.105630465
313+0.731634647E-1, -0.105630465, +0.382529974
314call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
315cov(:,:,2)
316+0.197313115, -0.512965582E-1, +0.731634647E-1
317-0.144657686, +0.474663615, -0.105630465
318+0.196772382, +0.289344782E-3, +0.382529974
319cov(:,:,0) ! reference
320+0.197313070, -0.512966625E-1, +0.731636062E-1
321-0.512966625E-1, +0.474663645, -0.105630495
322+0.731636062E-1, -0.105630495, +0.382530123
323
324
325dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
326do isam = 2, nsam
327 lb(isam) = ub(isam - 1) + 1
328 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
329end do
330lb
331+1, +6
332ub
333+5, +11
334ndim = getUnifRand(1, minval(ub - lb + 1, 1))
335call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
336call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
337call setResized(covMerged, [ndim, ndim])
338sample = getUnifRand(-1., +1., ndim, ub(nsam))
339sample
340+0.811735034, +0.638141155, -0.112208962, +0.953350425, +0.801234245, +0.169725299, +0.711966515, +0.350052714, +0.654049397, -0.618703723, -0.861874104
341cov(:,:,0) = getCov(sample, dim)
342cov(:,:,0) ! reference
343+0.339018017
344do isam = 1, nsam
345 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
346 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
347end do
348call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
349covMerged
350+0.339017987
351call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
352covMerged
353+0.339017987
354call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
355cov(:,:,2)
356+0.339017987
357cov(:,:,0) ! reference
358+0.339018017
359
360
361dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
362do isam = 2, nsam
363 lb(isam) = ub(isam - 1) + 1
364 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
365end do
366lb
367+1, +5
368ub
369+4, +8
370ndim = getUnifRand(1, minval(ub - lb + 1, 1))
371call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
372call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
373call setResized(covMerged, [ndim, ndim])
374sample = getUnifRand(-1., +1., ndim, ub(nsam))
375sample
376+0.451006889E-1, +0.669973969, +0.322528005, -0.842590928, -0.508468390, +0.557825089, +0.427691817, -0.475879073
377+0.561920881, +0.805896997, -0.514683008, -0.515792370E-1, +0.963925838, -0.593674660, +0.269054413, +0.195440173
378cov(:,:,0) = getCov(sample, dim)
379cov(:,:,0) ! reference
380+0.279895276, -0.495773517E-1
381-0.495773517E-1, +0.286286086
382do isam = 1, nsam
383 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim)
384 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim)
385end do
386call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
387covMerged
388+0.279895246, -0.495773293E-1
389+0.307122584E-40, +0.286286086
390call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), lowDia)
391covMerged
392+0.279895246, -0.495773293E-1
393-0.495773293E-1, +0.286286086
394call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(ub(1), TKG) / real(ub(2), TKG), uppDia)
395cov(:,:,2)
396+0.279895246, -0.495773293E-1
397-0.199867457, +0.286286086
398cov(:,:,0) ! reference
399+0.279895276, -0.495773517E-1
400-0.495773517E-1, +0.286286086
401
402
403!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
404!Compute the biased merged covariance of a frequency weighted multivariate sample.
405!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
406
407
408dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
409do isam = 2, nsam
410 lb(isam) = ub(isam - 1) + 1
411 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
412end do
413lb
414+1, +6
415ub
416+5, +9
417ndim = getUnifRand(1, minval(ub - lb + 1, 1))
418call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
419call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
420call setResized(covMerged, [ndim, ndim])
421sample = getUnifRand(-1., +1., ndim, ub(nsam))
422sample
423-0.485566378, -0.353631377, +0.514807105, -0.625597596, +0.457923651, +0.109889507, -0.981246233E-1, -0.294900775, +0.325146079
424iweight = getUnifRand(1, 10, size(sample, dim, IK))
425iweight
426+5, +5, +1, +10, +8, +9, +10, +10, +2
427cov(:,:,0) = getCov(sample, dim, iweight)
428cov(:,:,0) ! reference
429+0.128732443
430do isam = 1, nsam
431 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
432 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
433end do
434call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
435covMerged
436+0.128732443
437call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
438covMerged
439+0.128732443
440call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
441cov(:,:,2)
442+0.128732443
443cov(:,:,0) ! reference
444+0.128732443
445
446
447dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
448do isam = 2, nsam
449 lb(isam) = ub(isam - 1) + 1
450 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
451end do
452lb
453+1, +3
454ub
455+2, +5
456ndim = getUnifRand(1, minval(ub - lb + 1, 1))
457call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
458call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
459call setResized(covMerged, [ndim, ndim])
460sample = getUnifRand(-1., +1., ndim, ub(nsam))
461sample
462-0.671662211, +0.370629191, +0.730410576, +0.309888601, -0.650067806
463-0.261559367, +0.306724072, -0.566862702, -0.531929970, +0.925978422
464iweight = getUnifRand(1, 10, size(sample, dim, IK))
465iweight
466+1, +7, +5, +2, +9
467cov(:,:,0) = getCov(sample, dim, iweight)
468cov(:,:,0) ! reference
469+0.336275518, -0.289024085
470-0.289024085, +0.372983068
471do isam = 1, nsam
472 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
473 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
474end do
475call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
476covMerged
477+0.336275667, -0.289024115
478+0.307122584E-40, +0.372983098
479call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
480covMerged
481+0.336275667, -0.289024115
482-0.289024115, +0.372983098
483call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
484cov(:,:,2)
485+0.336275667, -0.289024115
486-0.461234450, +0.372983098
487cov(:,:,0) ! reference
488+0.336275518, -0.289024085
489-0.289024085, +0.372983068
490
491
492dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
493do isam = 2, nsam
494 lb(isam) = ub(isam - 1) + 1
495 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
496end do
497lb
498+1, +4
499ub
500+3, +8
501ndim = getUnifRand(1, minval(ub - lb + 1, 1))
502call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
503call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
504call setResized(covMerged, [ndim, ndim])
505sample = getUnifRand(-1., +1., ndim, ub(nsam))
506sample
507+0.114958525, +0.681636333E-1, -0.666265607, -0.406709075, +0.473294973, -0.220571995, +0.385011554, +0.571779609
508+0.580688119, +0.251309872, -0.475486636, +0.138509274, -0.739216089, +0.805498362E-1, +0.599527478, -0.505210280
509iweight = getUnifRand(1, 10, size(sample, dim, IK))
510iweight
511+1, +3, +7, +9, +1, +5, +2, +7
512cov(:,:,0) = getCov(sample, dim, iweight)
513cov(:,:,0) ! reference
514+0.207125053, -0.156532694E-1
515-0.156532694E-1, +0.143770218
516do isam = 1, nsam
517 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
518 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
519end do
520call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
521covMerged
522+0.207125038, -0.156532433E-1
523-0.289024115, +0.143770248
524call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
525covMerged
526+0.207125038, -0.156532433E-1
527-0.156532433E-1, +0.143770248
528call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
529cov(:,:,2)
530+0.207125038, -0.156532433E-1
531-0.103226453, +0.143770248
532cov(:,:,0) ! reference
533+0.207125053, -0.156532694E-1
534-0.156532694E-1, +0.143770218
535
536
537dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
538do isam = 2, nsam
539 lb(isam) = ub(isam - 1) + 1
540 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
541end do
542lb
543+1, +8
544ub
545+7, +10
546ndim = getUnifRand(1, minval(ub - lb + 1, 1))
547call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
548call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
549call setResized(covMerged, [ndim, ndim])
550sample = getUnifRand(-1., +1., ndim, ub(nsam))
551sample
552+0.241144776, +0.597155333, +0.955941796, -0.879346609, +0.387665033E-1, +0.991892576, +0.845819831, +0.961139202E-1, -0.233176112, -0.929979324
553+0.615392923, +0.667261481, +0.412907243, +0.176240206E-1, +0.630695462, -0.544753075E-1, +0.247068524, -0.679912090, +0.764704466, -0.613764048
554iweight = getUnifRand(1, 10, size(sample, dim, IK))
555iweight
556+7, +8, +7, +3, +5, +7, +4, +3, +4, +9
557cov(:,:,0) = getCov(sample, dim, iweight)
558cov(:,:,0) ! reference
559+0.477623731, +0.162313640
560+0.162313640, +0.251158118
561do isam = 1, nsam
562 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
563 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
564end do
565call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
566covMerged
567+0.477623641, +0.162313625
568-0.156532433E-1, +0.251158059
569call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
570covMerged
571+0.477623641, +0.162313625
572+0.162313625, +0.251158059
573call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
574cov(:,:,2)
575+0.477623641, +0.162313625
576+0.105616361, +0.251158059
577cov(:,:,0) ! reference
578+0.477623731, +0.162313640
579+0.162313640, +0.251158118
580
581
582dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
583do isam = 2, nsam
584 lb(isam) = ub(isam - 1) + 1
585 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
586end do
587lb
588+1, +5
589ub
590+4, +7
591ndim = getUnifRand(1, minval(ub - lb + 1, 1))
592call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
593call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
594call setResized(covMerged, [ndim, ndim])
595sample = getUnifRand(-1., +1., ndim, ub(nsam))
596sample
597-0.813448429, +0.475425363, +0.377520204, -0.527169704E-1, +0.649902225, +0.681423664, -0.584281683
598iweight = getUnifRand(1, 10, size(sample, dim, IK))
599iweight
600+6, +6, +5, +7, +10, +2, +6
601cov(:,:,0) = getCov(sample, dim, iweight)
602cov(:,:,0) ! reference
603+0.307305664
604do isam = 1, nsam
605 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
606 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
607end do
608call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
609covMerged
610+0.307305753
611call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
612covMerged
613+0.307305753
614call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
615cov(:,:,2)
616+0.307305753
617cov(:,:,0) ! reference
618+0.307305664
619
620
621dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
622do isam = 2, nsam
623 lb(isam) = ub(isam - 1) + 1
624 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
625end do
626lb
627+1, +7
628ub
629+6, +8
630ndim = getUnifRand(1, minval(ub - lb + 1, 1))
631call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
632call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
633call setResized(covMerged, [ndim, ndim])
634sample = getUnifRand(-1., +1., ndim, ub(nsam))
635sample
636+0.556881428E-1, -0.863156438, -0.562242985, +0.397259951, +0.313982964E-1, -0.797294378E-1, +0.191023707, +0.799912930
637iweight = getUnifRand(1, 10, size(sample, dim, IK))
638iweight
639+9, +10, +9, +3, +2, +5, +9, +2
640cov(:,:,0) = getCov(sample, dim, iweight)
641cov(:,:,0) ! reference
642+0.219997913
643do isam = 1, nsam
644 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
645 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
646end do
647call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
648covMerged
649+0.219997883
650call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
651covMerged
652+0.219997883
653call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
654cov(:,:,2)
655+0.219997883
656cov(:,:,0) ! reference
657+0.219997913
658
659
660dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
661do isam = 2, nsam
662 lb(isam) = ub(isam - 1) + 1
663 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
664end do
665lb
666+1, +4
667ub
668+3, +7
669ndim = getUnifRand(1, minval(ub - lb + 1, 1))
670call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
671call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
672call setResized(covMerged, [ndim, ndim])
673sample = getUnifRand(-1., +1., ndim, ub(nsam))
674sample
675+0.109876513, -0.432837605, -0.398776531, +0.501464605E-1, +0.385011435, -0.849950075, -0.776914477
676iweight = getUnifRand(1, 10, size(sample, dim, IK))
677iweight
678+3, +2, +2, +5, +8, +8, +4
679cov(:,:,0) = getCov(sample, dim, iweight)
680cov(:,:,0) ! reference
681+0.255182505
682do isam = 1, nsam
683 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
684 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
685end do
686call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
687covMerged
688+0.255182505
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.255182505
692call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
693cov(:,:,2)
694+0.255182505
695cov(:,:,0) ! reference
696+0.255182505
697
698
699dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
700do isam = 2, nsam
701 lb(isam) = ub(isam - 1) + 1
702 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
703end do
704lb
705+1, +7
706ub
707+6, +13
708ndim = getUnifRand(1, minval(ub - lb + 1, 1))
709call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
710call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
711call setResized(covMerged, [ndim, ndim])
712sample = getUnifRand(-1., +1., ndim, ub(nsam))
713sample
714-0.128016353, -0.324343443E-1, -0.941679239, -0.464350820, +0.354151130, +0.266821980, -0.477802038, +0.552601814E-1, +0.413681269, +0.235629559, +0.227665663, -0.145036817, +0.203377366
715-0.182116270, -0.352576494, +0.929882526, -0.115843177, -0.538827300, +0.937554479, +0.515086651, -0.853407383, +0.722092271, +0.295398116, -0.471053123, -0.247512460, -0.999714375
716iweight = getUnifRand(1, 10, size(sample, dim, IK))
717iweight
718+10, +7, +2, +9, +2, +6, +4, +10, +4, +1, +5, +10, +6
719cov(:,:,0) = getCov(sample, dim, iweight)
720cov(:,:,0) ! reference
721+0.891153738E-1, -0.273560192E-1
722-0.273560192E-1, +0.327686310
723do isam = 1, nsam
724 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
725 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
726end do
727call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
728covMerged
729+0.891153738E-1, -0.273560360E-1
730+0.307122584E-40, +0.327686280
731call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
732covMerged
733+0.891153738E-1, -0.273560360E-1
734-0.273560360E-1, +0.327686280
735call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
736cov(:,:,2)
737+0.891153738E-1, -0.273560360E-1
738-0.270112045E-1, +0.327686280
739cov(:,:,0) ! reference
740+0.891153738E-1, -0.273560192E-1
741-0.273560192E-1, +0.327686310
742
743
744dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
745do isam = 2, nsam
746 lb(isam) = ub(isam - 1) + 1
747 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
748end do
749lb
750+1, +3
751ub
752+2, +7
753ndim = getUnifRand(1, minval(ub - lb + 1, 1))
754call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
755call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
756call setResized(covMerged, [ndim, ndim])
757sample = getUnifRand(-1., +1., ndim, ub(nsam))
758sample
759-0.314625382, +0.445162773, +0.146602392E-1, +0.776642680, +0.803261995, +0.951674938, +0.242754102
760-0.603415370, +0.112859488, -0.342069149, -0.687069654, -0.140752554, +0.699705720, -0.246074080
761iweight = getUnifRand(1, 10, size(sample, dim, IK))
762iweight
763+10, +4, +10, +9, +3, +4, +5
764cov(:,:,0) = getCov(sample, dim, iweight)
765cov(:,:,0) ! reference
766+0.204305604, +0.753671974E-1
767+0.753671974E-1, +0.156673118
768do isam = 1, nsam
769 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
770 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
771end do
772call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
773covMerged
774+0.204305619, +0.753672123E-1
775-0.273560360E-1, +0.156673118
776call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
777covMerged
778+0.204305619, +0.753672123E-1
779+0.753672123E-1, +0.156673118
780call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
781cov(:,:,2)
782+0.204305619, +0.753672123E-1
783+0.370159410E-1, +0.156673118
784cov(:,:,0) ! reference
785+0.204305604, +0.753671974E-1
786+0.753671974E-1, +0.156673118
787
788
789dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
790do isam = 2, nsam
791 lb(isam) = ub(isam - 1) + 1
792 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
793end do
794lb
795+1, +7
796ub
797+6, +8
798ndim = getUnifRand(1, minval(ub - lb + 1, 1))
799call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
800call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
801call setResized(covMerged, [ndim, ndim])
802sample = getUnifRand(-1., +1., ndim, ub(nsam))
803sample
804-0.129903913, +0.713821888, +0.700427771, +0.229957104E-1, +0.803664088, +0.953495502E-2, -0.549821496, -0.874510884
805iweight = getUnifRand(1, 10, size(sample, dim, IK))
806iweight
807+9, +6, +9, +10, +4, +2, +2, +3
808cov(:,:,0) = getCov(sample, dim, iweight)
809cov(:,:,0) ! reference
810+0.249968603
811do isam = 1, nsam
812 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
813 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, iweight(lb(isam):ub(isam)))
814end do
815call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
816covMerged
817+0.249968603
818call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), lowDia)
819covMerged
820+0.249968603
821call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(iweight(:ub(1))), TKG) / real(sum(iweight), TKG), uppDia)
822cov(:,:,2)
823+0.249968603
824cov(:,:,0) ! reference
825+0.249968603
826
827
828!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
829!Compute the biased merged covariance of a reliability weighted multivariate sample.
830!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
831
832
833dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
834do isam = 2, nsam
835 lb(isam) = ub(isam - 1) + 1
836 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
837end do
838lb
839+1, +4
840ub
841+3, +9
842ndim = getUnifRand(1, minval(ub - lb + 1, 1))
843call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
844call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
845call setResized(covMerged, [ndim, ndim])
846sample = getUnifRand(-1., +1., ndim, ub(nsam))
847sample
848+0.336129785, +0.399799824, -0.655887246, -0.276351333, +0.302923799, +0.299154282, +0.593482137, -0.267784119, -0.287802339
849-0.982821107, +0.674128175, +0.112987995, +0.414405942, -0.783872247, +0.393260717, +0.413521290, -0.542125583, -0.192581534
850+0.809373140, -0.470817685, +0.182070017, -0.686879992, -0.659692526, +0.284913898, +0.322795749, -0.741516232, +0.567436218E-4
851rweight = getUnifRand(1., 2., size(sample, dim, IK))
852rweight
853+1.73429620, +1.44691420, +1.69528043, +1.98914814, +1.05251694, +1.83924091, +1.90870798, +1.17813289, +1.75417948
854cov(:,:,0) = getCov(sample, dim, rweight)
855cov(:,:,0) ! reference
856+0.167132929, +0.114290174E-1, +0.655004904E-1
857+0.114290174E-1, +0.298997104, -0.610715635E-1
858+0.655004904E-1, -0.610715635E-1, +0.264802158
859do isam = 1, nsam
860 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
861 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
862end do
863call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
864covMerged
865+0.167132914, +0.114290435E-1, +0.655005127E-1
866+0.307122584E-40, +0.298997045, -0.610715747E-1
867+0.00000000, +0.280259693E-44, +0.264802158
868call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
869covMerged
870+0.167132914, +0.114290435E-1, +0.655005127E-1
871+0.114290435E-1, +0.298997045, -0.610715747E-1
872+0.655005127E-1, -0.610715747E-1, +0.264802158
873call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
874cov(:,:,2)
875+0.167132914, +0.453154147E-1, +0.106991664
876+0.114290435E-1, +0.298997045, +0.104022406
877+0.655005127E-1, -0.610715747E-1, +0.264802158
878cov(:,:,0) ! reference
879+0.167132929, +0.114290174E-1, +0.655004904E-1
880+0.114290174E-1, +0.298997104, -0.610715635E-1
881+0.655004904E-1, -0.610715635E-1, +0.264802158
882
883
884dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
885do isam = 2, nsam
886 lb(isam) = ub(isam - 1) + 1
887 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
888end do
889lb
890+1, +8
891ub
892+7, +11
893ndim = getUnifRand(1, minval(ub - lb + 1, 1))
894call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
895call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
896call setResized(covMerged, [ndim, ndim])
897sample = getUnifRand(-1., +1., ndim, ub(nsam))
898sample
899-0.246406794, +0.280203819E-1, -0.689879060, +0.791810870, +0.545506358, -0.484165907, +0.362807035, +0.538294315, -0.527532578, -0.997618198, +0.349747896
900-0.582575798E-2, -0.572107911, +0.747431397, -0.835165620, -0.663823485, +0.989633679, -0.974120378, +0.165681124, -0.800786614, +0.405744195, -0.267635345
901-0.207938790, -0.571386576, +0.696033478, -0.542981744, -0.764543653, +0.568241596, -0.758426785, -0.361368299, -0.770800710, +0.444702625, +0.403569937
902rweight = getUnifRand(1., 2., size(sample, dim, IK))
903rweight
904+1.21100819, +1.30155742, +1.70866656, +1.27854586, +1.09147906, +1.24509311, +1.30848384, +1.09948349, +1.51119423, +1.47240877, +1.32550883
905cov(:,:,0) = getCov(sample, dim, rweight)
906cov(:,:,0) ! reference
907+0.326425940, -0.230338082, -0.186394781
908-0.230338082, +0.429419011, +0.331138372
909-0.186394781, +0.331138372, +0.326418370
910do isam = 1, nsam
911 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
912 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
913end do
914call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
915covMerged
916+0.326425910, -0.230338082, -0.186394796
917+0.114290435E-1, +0.429418981, +0.331138372
918+0.655005127E-1, -0.610715747E-1, +0.326418400
919call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
920covMerged
921+0.326425910, -0.230338082, -0.186394796
922-0.230338082, +0.429418981, +0.331138372
923-0.186394796, +0.331138372, +0.326418400
924call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
925cov(:,:,2)
926+0.326425910, -0.294733923E-1, -0.275275931E-1
927-0.230338082, +0.429418981, +0.172943875
928-0.186394796, +0.331138372, +0.326418400
929cov(:,:,0) ! reference
930+0.326425940, -0.230338082, -0.186394781
931-0.230338082, +0.429419011, +0.331138372
932-0.186394781, +0.331138372, +0.326418370
933
934
935dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
936do isam = 2, nsam
937 lb(isam) = ub(isam - 1) + 1
938 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
939end do
940lb
941+1, +8
942ub
943+7, +9
944ndim = getUnifRand(1, minval(ub - lb + 1, 1))
945call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
946call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
947call setResized(covMerged, [ndim, ndim])
948sample = getUnifRand(-1., +1., ndim, ub(nsam))
949sample
950+0.642199993, +0.106531382, -0.528338194, -0.496698022, -0.176626563, +0.261250615, -0.830453634, -0.751735210, -0.710292697
951rweight = getUnifRand(1., 2., size(sample, dim, IK))
952rweight
953+1.23513579, +1.96217811, +1.29759443, +1.17209673, +1.68872261, +1.55529404, +1.72909141, +1.64058971, +1.79977822
954cov(:,:,0) = getCov(sample, dim, rweight)
955cov(:,:,0) ! reference
956+0.225324735
957do isam = 1, nsam
958 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
959 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
960end do
961call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
962covMerged
963+0.225324705
964call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
965covMerged
966+0.225324705
967call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
968cov(:,:,2)
969+0.225324705
970cov(:,:,0) ! reference
971+0.225324735
972
973
974dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
975do isam = 2, nsam
976 lb(isam) = ub(isam - 1) + 1
977 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
978end do
979lb
980+1, +6
981ub
982+5, +11
983ndim = getUnifRand(1, minval(ub - lb + 1, 1))
984call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
985call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
986call setResized(covMerged, [ndim, ndim])
987sample = getUnifRand(-1., +1., ndim, ub(nsam))
988sample
989-0.843564272E-1, +0.897690654, -0.388406396, +0.262057185, +0.799131393E-2, +0.859550357, +0.541457772, +0.322756529, -0.100277185, +0.448757887, -0.269633532E-1
990+0.621384382E-1, +0.910715818, +0.747956038E-1, -0.741100073, -0.886915445, +0.791664720, -0.952938795E-1, -0.995840073, +0.618917942E-1, +0.355942965, +0.193558931E-1
991rweight = getUnifRand(1., 2., size(sample, dim, IK))
992rweight
993+1.50359178, +1.21140695, +1.89611542, +1.84549332, +1.19793797, +1.88025129, +1.64584541, +1.33082783, +1.54656208, +1.29058266, +1.94891465
994cov(:,:,0) = getCov(sample, dim, rweight)
995cov(:,:,0) ! reference
996+0.158984810, +0.946634933E-1
997+0.946634933E-1, +0.326425970
998do isam = 1, nsam
999 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1000 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1001end do
1002call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1003covMerged
1004+0.158984780, +0.946634635E-1
1005+0.307122584E-40, +0.326426029
1006call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1007covMerged
1008+0.158984780, +0.946634635E-1
1009+0.946634635E-1, +0.326426029
1010call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1011cov(:,:,2)
1012+0.158984780, +0.784763321E-1
1013+0.946634635E-1, +0.326426029
1014cov(:,:,0) ! reference
1015+0.158984810, +0.946634933E-1
1016+0.946634933E-1, +0.326425970
1017
1018
1019dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1020do isam = 2, nsam
1021 lb(isam) = ub(isam - 1) + 1
1022 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1023end do
1024lb
1025+1, +5
1026ub
1027+4, +7
1028ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1029call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1030call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1031call setResized(covMerged, [ndim, ndim])
1032sample = getUnifRand(-1., +1., ndim, ub(nsam))
1033sample
1034-0.970712543, -0.779816985, -0.814021826E-1, +0.633818150, +0.635938406, -0.616514444, +0.904412866
1035+0.248281360, +0.656293035, -0.695680261, -0.273201942, +0.750689983, -0.472333431, -0.898378968
1036-0.412211895, -0.325180769, +0.821596265, +0.831295609, +0.252224565, +0.696952820, -0.506856203
1037rweight = getUnifRand(1., 2., size(sample, dim, IK))
1038rweight
1039+1.43145823, +1.52432418, +1.98789155, +1.28443980, +1.05062294, +1.57613015, +1.03030109
1040cov(:,:,0) = getCov(sample, dim, rweight)
1041cov(:,:,0) ! reference
1042+0.455360591, -0.129172876, +0.854417384E-1
1043-0.129172876, +0.345723301, -0.145424947
1044+0.854417384E-1, -0.145424947, +0.315681368
1045do isam = 1, nsam
1046 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1047 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1048end do
1049call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1050covMerged
1051+0.455360651, -0.129172847, +0.854417682E-1
1052+0.946634635E-1, +0.345723242, -0.145424977
1053+0.00000000, +0.261250615, +0.315681398
1054call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1055covMerged
1056+0.455360651, -0.129172847, +0.854417682E-1
1057-0.129172847, +0.345723242, -0.145424977
1058+0.854417682E-1, -0.145424977, +0.315681398
1059call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1060cov(:,:,2)
1061+0.455360651, +0.751461089E-1, -0.307770282
1062-0.129172847, +0.345723242, +0.962442905E-1
1063+0.854417682E-1, -0.145424977, +0.315681398
1064cov(:,:,0) ! reference
1065+0.455360591, -0.129172876, +0.854417384E-1
1066-0.129172876, +0.345723301, -0.145424947
1067+0.854417384E-1, -0.145424947, +0.315681368
1068
1069
1070dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1071do isam = 2, nsam
1072 lb(isam) = ub(isam - 1) + 1
1073 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1074end do
1075lb
1076+1, +7
1077ub
1078+6, +13
1079ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1080call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1081call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1082call setResized(covMerged, [ndim, ndim])
1083sample = getUnifRand(-1., +1., ndim, ub(nsam))
1084sample
1085+0.287596345, +0.609653115, +0.701835155E-1, -0.371857166, -0.107286572, -0.796721578, -0.667690516, +0.959987044, +0.710416555, -0.663437843E-1, +0.794315577, +0.904552698, -0.558568716
1086-0.532123685, +0.220215201, +0.521948576, -0.712247491, +0.296623945, -0.522906542, +0.473205566, +0.910298586, +0.557248473, -0.341901779E-1, +0.374021411, -0.978827119, -0.857307076
1087+0.893969059, +0.618004084, -0.998070955, -0.156439543E-1, +0.805911303, +0.570536613, -0.393446326, +0.548493028, +0.750734687, -0.743613362, +0.466376781, -0.271179199, -0.352984548
1088+0.123402238, +0.764948130, -0.232205153, -0.446975231E-2, +0.957462549, -0.688403010, +0.484785318, -0.360543370, +0.765505314, +0.701770067, -0.694312096, +0.514655232, +0.398856401E-1
1089-0.809070110, +0.100823522, -0.598138213, -0.366575718, +0.419979572, -0.235639691, +0.961288691, -0.813308120, +0.502279401, -0.950531125, +0.222866297, +0.666248441, -0.831619859
1090+0.407885432, -0.583597898, -0.709315062, -0.914039850, -0.640670061, -0.136670709, -0.617813826, +0.784289956, -0.370497823, +0.428459525, +0.288957834, +0.792332768, -0.716278672
1091rweight = getUnifRand(1., 2., size(sample, dim, IK))
1092rweight
1093+1.40707326, +1.90542054, +1.16595840, +1.78542459, +1.44870067, +1.56057048, +1.37476301, +1.77938473, +1.31904650, +1.94967294, +1.81528306, +1.46870184, +1.37535346
1094cov(:,:,0) = getCov(sample, dim, rweight)
1095cov(:,:,0) ! reference
1096+0.358092040, +0.135784149, +0.115340963, +0.138025507E-1, +0.427269153E-1, +0.207665399
1097+0.135784149, +0.350155115, +0.638591871E-1, +0.600517448E-2, +0.481187701E-1, +0.215528999E-1
1098+0.115340963, +0.638591871E-1, +0.353281200, -0.211410671E-1, +0.724910200E-1, +0.294857826E-1
1099+0.138025507E-1, +0.600517448E-2, -0.211410671E-1, +0.303482622, +0.103677869, -0.676821247E-1
1100+0.427269153E-1, +0.481187701E-1, +0.724910200E-1, +0.103677869, +0.377885282, -0.734018087E-1
1101+0.207665399, +0.215528999E-1, +0.294857826E-1, -0.676821247E-1, -0.734018087E-1, +0.353962362
1102do isam = 1, nsam
1103 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1104 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1105end do
1106call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1107covMerged
1108+0.358092099, +0.135784075, +0.115340956, +0.138025153E-1, +0.427267551E-1, +0.207665354
1109-0.129172847, +0.350155115, +0.638591498E-1, +0.600520056E-2, +0.481188744E-1, +0.215528421E-1
1110+0.854417682E-1, -0.145424977, +0.353281111, -0.211410373E-1, +0.724909082E-1, +0.294857956E-1
1111+0.00000000, +0.00000000, +0.00000000, +0.303482562, +0.103677869, -0.676821470E-1
1112+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.377885222, -0.734018385E-1
1113+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.353962332
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.358092099, +0.135784075, +0.115340956, +0.138025153E-1, +0.427267551E-1, +0.207665354
1117+0.135784075, +0.350155115, +0.638591498E-1, +0.600520056E-2, +0.481188744E-1, +0.215528421E-1
1118+0.115340956, +0.638591498E-1, +0.353281111, -0.211410373E-1, +0.724909082E-1, +0.294857956E-1
1119+0.138025153E-1, +0.600520056E-2, -0.211410373E-1, +0.303482562, +0.103677869, -0.676821470E-1
1120+0.427267551E-1, +0.481188744E-1, +0.724909082E-1, +0.103677869, +0.377885222, -0.734018385E-1
1121+0.207665354, +0.215528421E-1, +0.294857956E-1, -0.676821470E-1, -0.734018385E-1, +0.353962332
1122call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1123cov(:,:,2)
1124+0.358092099, +0.110046245, +0.243737906, -0.127852440, +0.466685891E-1, +0.271181405
1125+0.135784075, +0.350155115, +0.213146374, -0.107149333, +0.107361667E-2, +0.293361265E-1
1126+0.115340956, +0.638591498E-1, +0.353281111, -0.151122853, +0.796599686E-1, +0.358438082E-1
1127+0.138025153E-1, +0.600520056E-2, -0.211410373E-1, +0.303482562, +0.664112493E-1, -0.711643174E-1
1128+0.427267551E-1, +0.481188744E-1, +0.724909082E-1, +0.103677869, +0.377885222, -0.101460040
1129+0.207665354, +0.215528421E-1, +0.294857956E-1, -0.676821470E-1, -0.734018385E-1, +0.353962332
1130cov(:,:,0) ! reference
1131+0.358092040, +0.135784149, +0.115340963, +0.138025507E-1, +0.427269153E-1, +0.207665399
1132+0.135784149, +0.350155115, +0.638591871E-1, +0.600517448E-2, +0.481187701E-1, +0.215528999E-1
1133+0.115340963, +0.638591871E-1, +0.353281200, -0.211410671E-1, +0.724910200E-1, +0.294857826E-1
1134+0.138025507E-1, +0.600517448E-2, -0.211410671E-1, +0.303482622, +0.103677869, -0.676821247E-1
1135+0.427269153E-1, +0.481187701E-1, +0.724910200E-1, +0.103677869, +0.377885282, -0.734018087E-1
1136+0.207665399, +0.215528999E-1, +0.294857826E-1, -0.676821247E-1, -0.734018087E-1, +0.353962362
1137
1138
1139dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1140do isam = 2, nsam
1141 lb(isam) = ub(isam - 1) + 1
1142 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1143end do
1144lb
1145+1, +3
1146ub
1147+2, +5
1148ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1149call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1150call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1151call setResized(covMerged, [ndim, ndim])
1152sample = getUnifRand(-1., +1., ndim, ub(nsam))
1153sample
1154+0.393426657, +0.871836424, +0.168076754E-1, +0.878648400, -0.516014934
1155+0.248538852, -0.486588717, -0.339158535, +0.417438984, +0.158906937
1156rweight = getUnifRand(1., 2., size(sample, dim, IK))
1157rweight
1158+1.03938174, +1.47683001, +1.78322101, +1.49937773, +1.23638535
1159cov(:,:,0) = getCov(sample, dim, rweight)
1160cov(:,:,0) ! reference
1161+0.276883870, -0.494360633E-3
1162-0.494360633E-3, +0.128370881
1163do isam = 1, nsam
1164 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1165 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1166end do
1167call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1168covMerged
1169+0.276883841, -0.494353473E-3
1170+0.135784075, +0.128370881
1171call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1172covMerged
1173+0.276883841, -0.494353473E-3
1174-0.494353473E-3, +0.128370881
1175call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1176cov(:,:,2)
1177+0.276883841, +0.894543231E-1
1178-0.494353473E-3, +0.128370881
1179cov(:,:,0) ! reference
1180+0.276883870, -0.494360633E-3
1181-0.494360633E-3, +0.128370881
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, +6
1191ub
1192+5, +11
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.424062014E-1, +0.391265750, -0.286736250, +0.133393049, +0.605205536, +0.753888011, +0.387332678, -0.505940676, +0.582209587, -0.846098661E-1, -0.301306248
1200-0.324899912, +0.393074512, -0.626505375, +0.762512803, -0.765664458, -0.953210950, +0.712095737, +0.473800540, +0.740413666E-1, +0.361572266, +0.912321210
1201rweight = getUnifRand(1., 2., size(sample, dim, IK))
1202rweight
1203+1.65867400, +1.02932429, +1.27839708, +1.60448480, +1.00452054, +1.07561648, +1.78529668, +1.62345636, +1.59856486, +1.91399944, +1.84600949
1204cov(:,:,0) = getCov(sample, dim, rweight)
1205cov(:,:,0) ! reference
1206+0.150918290, -0.877128020E-1
1207-0.877128020E-1, +0.355788380
1208do isam = 1, nsam
1209 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1210 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1211end do
1212call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1213covMerged
1214+0.150918335, -0.877128243E-1
1215-0.494353473E-3, +0.355788350
1216call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1217covMerged
1218+0.150918335, -0.877128243E-1
1219-0.877128243E-1, +0.355788350
1220call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1221cov(:,:,2)
1222+0.150918335, -0.150912106
1223-0.877128243E-1, +0.355788350
1224cov(:,:,0) ! reference
1225+0.150918290, -0.877128020E-1
1226-0.877128020E-1, +0.355788380
1227
1228
1229dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1230do isam = 2, nsam
1231 lb(isam) = ub(isam - 1) + 1
1232 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1233end do
1234lb
1235+1, +4
1236ub
1237+3, +5
1238ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1239call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1240call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1241call setResized(covMerged, [ndim, ndim])
1242sample = getUnifRand(-1., +1., ndim, ub(nsam))
1243sample
1244-0.625773668, -0.648405671, -0.711457729, +0.795756817, +0.507063031
1245-0.330471039, -0.336011291, -0.214211464, -0.713448167, +0.304769158
1246rweight = getUnifRand(1., 2., size(sample, dim, IK))
1247rweight
1248+1.13904023, +1.30430484, +1.27402830, +1.18407238, +1.05133414
1249cov(:,:,0) = getCov(sample, dim, rweight)
1250cov(:,:,0) ! reference
1251+0.419111252, -0.107254609E-1
1252-0.107254609E-1, +0.997806042E-1
1253do isam = 1, nsam
1254 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1255 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1256end do
1257call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1258covMerged
1259+0.419111252, -0.107254628E-1
1260-0.877128243E-1, +0.997805968E-1
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.419111252, -0.107254628E-1
1264-0.107254628E-1, +0.997805968E-1
1265call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1266cov(:,:,2)
1267+0.419111252, -0.732291341E-1
1268-0.107254628E-1, +0.997805968E-1
1269cov(:,:,0) ! reference
1270+0.419111252, -0.107254609E-1
1271-0.107254609E-1, +0.997806042E-1
1272
1273
1274dim = 2; lb(1) = 1; ub(1) = getUnifRand(2, 7)
1275do isam = 2, nsam
1276 lb(isam) = ub(isam - 1) + 1
1277 ub(isam) = ub(isam - 1) + getUnifRand(2, 7)
1278end do
1279lb
1280+1, +7
1281ub
1282+6, +12
1283ndim = getUnifRand(1, minval(ub - lb + 1, 1))
1284call setRebound(cov, [1_IK, 1_IK, 0_IK], [ndim, ndim, nsam])
1285call setRebound(mean, [1_IK, 1_IK], [ndim, nsam])
1286call setResized(covMerged, [ndim, ndim])
1287sample = getUnifRand(-1., +1., ndim, ub(nsam))
1288sample
1289+0.881685853, -0.114970922, +0.905977368, -0.865088582, +0.705984712, +0.140058994E-2, -0.228908420, -0.851893663, -0.607634664, +0.988011360E-1, +0.543577433, +0.792512059
1290-0.197809339, +0.469179630, -0.515428782E-1, +0.588172317, +0.150875449, -0.849586487, +0.369192362E-1, -0.248943806, -0.572045803, -0.945229530, -0.528795123, +0.846990824
1291-0.416366220, +0.428619385E-1, +0.267181993, -0.380258322, -0.583978891, -0.520030022, -0.378844142, +0.721287489, +0.871469498, -0.480708718, -0.542719603, -0.486914635
1292+0.859933376, +0.893869400, +0.981066704, +0.941281319E-1, -0.642131567E-1, -0.985565186, -0.391182423, -0.627261281, -0.631874084, +0.165917635, +0.882730365, -0.625619292
1293+0.855932117, -0.866293907, -0.684256673, +0.727418184, -0.686436057, -0.469346404, +0.223051310, +0.369931221, +0.895434737, -0.262888551, -0.773062825, -0.906424880
1294rweight = getUnifRand(1., 2., size(sample, dim, IK))
1295rweight
1296+1.40264380, +1.72922564, +1.92329550, +1.21741176, +1.12335014, +1.28854942, +1.24030256, +1.91249955, +1.11259890, +1.40713239, +1.54322839, +1.20640409
1297cov(:,:,0) = getCov(sample, dim, rweight)
1298cov(:,:,0) ! reference
1299+0.412137538, +0.176595803E-1, -0.148879126, +0.231813774, -0.221702054
1300+0.176595803E-1, +0.270661592, -0.139986966E-1, +0.491984598E-1, -0.558295771E-1
1301-0.148879126, -0.139986966E-1, +0.241976827, -0.486526527E-1, +0.106253289
1302+0.231813774, +0.491984598E-1, -0.486526527E-1, +0.493014127, -0.127302885
1303-0.221702054, -0.558295771E-1, +0.106253289, -0.127302885, +0.435951948
1304do isam = 1, nsam
1305 cov(:,:,isam) = getCov(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1306 mean(:,isam) = getMean(sample(:,lb(isam):ub(isam)), dim, rweight(lb(isam):ub(isam)))
1307end do
1308call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), uppDia)
1309covMerged
1310+0.412137419, +0.176596157E-1, -0.148879141, +0.231813729, -0.221702233
1311-0.107254628E-1, +0.270661622, -0.139986994E-1, +0.491984524E-1, -0.558295585E-1
1312+0.00000000, +0.00000000, +0.241976887, -0.486527160E-1, +0.106253326
1313+0.00000000, +0.00000000, +0.00000000, +0.493014157, -0.127302989
1314+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.435951889
1315call setCovMerged(covMerged, cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1316covMerged
1317+0.412137419, +0.176596157E-1, -0.148879141, +0.231813729, -0.221702233
1318+0.176596157E-1, +0.270661622, -0.139986994E-1, +0.491984524E-1, -0.558295585E-1
1319-0.148879141, -0.139986994E-1, +0.241976887, -0.486527160E-1, +0.106253326
1320+0.231813729, +0.491984524E-1, -0.486527160E-1, +0.493014157, -0.127302989
1321-0.221702233, -0.558295585E-1, +0.106253326, -0.127302989, +0.435951889
1322call setCovMerged(cov(:,:,2), cov(:,:,1), mean(:,1) - mean(:,2), real(sum(rweight(:ub(1))), TKG) / real(sum(rweight), TKG), lowDia)
1323cov(:,:,2)
1324+0.412137419, +0.100470573, -0.316032201, +0.190600798, -0.341439635
1325+0.176596157E-1, +0.270661622, -0.461273007E-1, -0.154335245, -0.101799726
1326-0.148879141, -0.139986994E-1, +0.241976887, -0.215977818, +0.310325861
1327+0.231813729, +0.491984524E-1, -0.486527160E-1, +0.493014157, -0.206127733
1328-0.221702233, -0.558295585E-1, +0.106253326, -0.127302989, +0.435951889
1329cov(:,:,0) ! reference
1330+0.412137538, +0.176595803E-1, -0.148879126, +0.231813774, -0.221702054
1331+0.176595803E-1, +0.270661592, -0.139986966E-1, +0.491984598E-1, -0.558295771E-1
1332-0.148879126, -0.139986966E-1, +0.241976827, -0.486526527E-1, +0.106253289
1333+0.231813774, +0.491984598E-1, -0.486526527E-1, +0.493014127, -0.127302885
1334-0.221702054, -0.558295771E-1, +0.106253289, -0.127302885, +0.435951948
1335
1336
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: