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

Generate and return an affine-transformation of the input sample of shape (1:ndim) or (1:ndim, 1:nsam) or (1:nsam, 1:ndim) based on the specified values for the translation tlate and tranformations tform along the specified axis dim.
More...

Detailed Description

Generate and return an affine-transformation of the input sample of shape (1:ndim) or (1:ndim, 1:nsam) or (1:nsam, 1:ndim) based on the specified values for the translation tlate and tranformations tform along the specified axis dim.

Here, ndim stands for the number of dimensions (data attributes) of the input sample and nsam represents the number of data points in the sample.
If the input tlate is the negative of the mean of the sample, then the returned sample will have a mean of zero.

Parameters
[in]sample: The input contiguous array of shape (nsam), (ndim, nsam), or (nsam, ndim) of,
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the sample to be affine-transformed.
[in]dim: The input scalar of type integer of default kind IK, whose value represents the dimension of the input sample containing different nsam observations:
  1. If dim = 1, the input sample is assumed to have the shape (1:nsam, 1:ndim).
  2. If dim = 2, the input sample is assumed to have the shape (1:ndim, 1:nsam).
(optional. It must be present if and only if the input arguments the condition rank(sample) > 1 holds.)
[in]tform: The input contiguous matrix of the same type and kind as sample of shape (1:ndim, 1:ndim), containing the transformation matrix of the affine-transformation.
[in]class: The input scalar constant that can be any of the following:
  1. The scalar constant genrecmat or equivalently an object of type genrecmat_type signifying the general rectangular matrix as the class of the input transformation matrix.
    All elements of the input matrix will be used in constructing the output sample.
  2. The scalar constant upperDiag or equivalently an object of type upperDiag_type signifying the upper-diagonal matrix as the class of the input transformation matrix.
    Only the upper-diagonal elements of tform will be used in constructing the output sample.
    The lower triangle elements of tform are assumed to be zero.
  3. The scalar constant lowerDiag or equivalently an object of type lowerDiag_type signifying the lower-diagonal matrix as the class of the input transformation matrix.
    Only the lower-diagonal elements of tform will be used in constructing the output sample.
    The upper triangle elements of tform are assumed to be zero.
  4. The scalar constant upperUnit or equivalently an object of type upperUnit_type signifying the unit-triangular matrix as the class of the input transformation matrix.
    Only the upper elements of tform will be used in constructing the output sample.
    The lower triangle elements of tform are assumed to be zero.
    The diagonal elements of tform are assumed to be one.
    The diagonal elements of tform are assumed to be zero.
  5. The scalar constant lowerUnit or equivalently an object of type lowerUnit_type signifying the unit-triangular matrix as the class of the input transformation matrix.
    Only the lower elements of tform will be used in constructing the output sample.
    The upper triangle elements of tform are assumed to be zero.
    The diagonal elements of tform are assumed to be zero.
(optional, default = genrecmat.)
[in]tlate: The input contiguous vector of the same type and kind as sample of shape (1:ndim), containing the amount by which the input sample must be translated (i.e., shifted).
  1. If the input argument dim = 1 then, size(tlate) == size(sample, 2) == ndim must hold.
  2. If the input argument dim = 2 then, size(tlate) == size(sample, 1) == ndim must hold.
(optional, default = [(0, i = 1, ndim)].)
Returns
affinity : The output object of the same type and kind and rank and shape as sample, containing the affine-transformed sample.


Possible calling interfaces

affinity(1:ndim) = getAffinity(sample(1:ndim), tform(1:ndim, 1:ndim), class = class, tlate = tlate(1:ndim))
affinity(@shape(sample)) = getAffinity(sample(:,:), dim, tform(1:size(sample, 3 - dim), 1:size(sample, 3 - dim)), class = class, tlate = tlate(1:size(sample, 3 - dim)))
Generate and return an affine-transformation of the input sample of shape (1:ndim) or (1:ndim,...
This module contains classes and procedures for affine transformation of multivariate samples.
Warning
The condition 1 <= dim .and. dim <= rank(sample) must hold for the corresponding input arguments.
The condition size(tlate) == size(sample, 3 - dim) must hold for the corresponding input arguments.
The condition all(shape(tform) == size(sample, 3 - dim)) 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
getShifted
setShifted
setAffinity
setAffinity


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_matrixChol, only: getMatChol, uppDia, lowDia
5 use pm_sampleAffinity, only: upperDiag, upperUnit
6 use pm_sampleAffinity, only: lowerDiag, lowerUnit
7 use pm_matrixInit, only: setMatInit, lowDia
9 use pm_distCov, only: getCovRand
10 use pm_arrayResize, only: setResized
11 use pm_distUnif, only: getUnifRand
12 use pm_arrayFill, only: getFilled
13 use pm_matrixInv, only: getMatInv
14 use pm_io, only: display_type
15
16 implicit none
17
18 type(display_type) :: disp
19 integer(IK) :: itry, ntry = 5
20 integer(IK) :: dim, isam, ndim, nsam
21 disp = display_type(file = "main.out.F90")
22
23 call disp%skip()
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%show("!Transform 2D sample along the second dimension.")
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%skip()
28
29 block
30 use pm_kind, only: RKG => RKS
31 real(RKG), allocatable :: tlate(:), tform(:,:), sample(:,:), affinInv(:,:), affinity(:,:)
32 do itry = 1, ntry
33 call disp%skip()
34 call disp%show("dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)")
35 dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
36 call disp%show("[dim, ndim, nsam]")
37 call disp%show( [dim, ndim, nsam] )
38 call disp%show("tlate = getUnifRand(-5, +5, ndim)")
39 tlate = getUnifRand(-5, +5, ndim)
40 call disp%show("tlate")
41 call disp%show( tlate )
42 call disp%show("sample = getUnifRand(-9, +9, ndim, nsam)")
43 sample = getUnifRand(-9, +9, ndim, nsam)
44 call disp%show("sample")
45 call disp%show( sample )
46 call disp%show("call setResized(affinity, shape(sample, IK))")
47 call setResized(affinity, shape(sample, IK))
48 call disp%show("call setResized(affinInv, shape(sample, IK))")
49 call setResized(affinInv, shape(sample, IK))
50 call disp%show("tform = getCovRand(1., ndim)")
51 tform = getCovRand(1., ndim)
52 call disp%show("tform")
53 call disp%show( tform )
54 call disp%show("affinity = getAffinity(sample, dim, tform) ! whole sample transformation.")
55 affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
56 call disp%show("transpose(affinity)")
57 call disp%show( transpose(affinity) )
58 call disp%show("do isam = 1, nsam")
59 call disp%show(" affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation")
60 call disp%show("end do")
61 do isam = 1, nsam
62 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
63 end do
64 call disp%show("transpose(affinity)")
65 call disp%show( transpose(affinity) )
66 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.")
67 affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
68 call disp%show("transpose(affinInv)")
69 call disp%show( transpose(affinInv) )
70 call disp%show("transpose(sample) ! for comparison with affinInv.")
71 call disp%show( transpose(sample) )
72 call disp%skip()
73 call disp%show("call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.")
74 call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
75 call disp%show("tform")
76 call disp%show( tform )
77 call disp%show("affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.")
78 affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
79 call disp%show("transpose(affinity)")
80 call disp%show( transpose(affinity) )
81 call disp%show("do isam = 1, nsam")
82 call disp%show(" affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation")
83 call disp%show("end do")
84 do isam = 1, nsam
85 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
86 end do
87 call disp%show("transpose(affinity)")
88 call disp%show( transpose(affinity) )
89 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.")
90 affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
91 call disp%show("transpose(affinInv)")
92 call disp%show( transpose(affinInv) )
93 call disp%show("transpose(sample) ! for comparison with affinInv.")
94 call disp%show( transpose(sample) )
95 call disp%skip()
96 call disp%show("affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.")
97 affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
98 call disp%show("transpose(affinity)")
99 call disp%show( transpose(affinity) )
100 call disp%show("do isam = 1, nsam")
101 call disp%show(" affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation")
102 call disp%show("end do")
103 do isam = 1, nsam
104 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
105 end do
106 call disp%show("transpose(affinity)")
107 call disp%show( transpose(affinity) )
108 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.")
109 affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
110 call disp%show("transpose(affinInv)")
111 call disp%show( transpose(affinInv) )
112 call disp%show("transpose(sample) ! for comparison with affinInv.")
113 call disp%show( transpose(sample) )
114 call disp%skip()
115 call disp%show("affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.")
116 affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
117 call disp%show("transpose(affinity)")
118 call disp%show( transpose(affinity) )
119 call disp%show("do isam = 1, nsam")
120 call disp%show(" affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation")
121 call disp%show("end do")
122 do isam = 1, nsam
123 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
124 end do
125 call disp%show("transpose(affinity)")
126 call disp%show( transpose(affinity) )
127 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.")
128 affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
129 call disp%show("transpose(affinInv)")
130 call disp%show( transpose(affinInv) )
131 call disp%show("transpose(sample) ! for comparison with affinInv.")
132 call disp%show( transpose(sample) )
133 call disp%skip()
134 call disp%show("affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.")
135 affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
136 call disp%show("transpose(affinity)")
137 call disp%show( transpose(affinity) )
138 call disp%show("do isam = 1, nsam")
139 call disp%show(" affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation")
140 call disp%show("end do")
141 do isam = 1, nsam
142 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
143 end do
144 call disp%show("transpose(affinity)")
145 call disp%show( transpose(affinity) )
146 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.")
147 affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
148 call disp%show("transpose(affinInv)")
149 call disp%show( transpose(affinInv) )
150 call disp%show("transpose(sample) ! for comparison with affinInv.")
151 call disp%show( transpose(sample) )
152 end do
153 end block
154
155 call disp%skip()
156 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
157 call disp%show("!Transform 2D sample along the first dimension.")
158 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
159 call disp%skip()
160
161 block
162 use pm_kind, only: RKG => RKS
163 real(RKG), allocatable :: tlate(:), tform(:,:), sample(:,:), affinInv(:,:), affinity(:,:)
164 do itry = 1, ntry
165 call disp%skip()
166 call disp%show("dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)")
167 dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
168 call disp%show("[dim, ndim, nsam]")
169 call disp%show( [dim, ndim, nsam] )
170 call disp%show("tlate = getUnifRand(-5, +5, ndim)")
171 tlate = getUnifRand(-5, +5, ndim)
172 call disp%show("tlate")
173 call disp%show( tlate )
174 call disp%show("sample = getUnifRand(-9, +9, nsam, ndim)")
175 sample = getUnifRand(-9, +9, nsam, ndim)
176 call disp%show("sample")
177 call disp%show( sample )
178 call disp%show("call setResized(affinity, shape(sample, IK))")
179 call setResized(affinity, shape(sample, IK))
180 call disp%show("call setResized(affinInv, shape(sample, IK))")
181 call setResized(affinInv, shape(sample, IK))
182 call disp%show("tform = getCovRand(1., ndim)")
183 tform = getCovRand(1., ndim)
184 call disp%show("tform")
185 call disp%show( tform )
186 call disp%show("affinity = getAffinity(sample, dim, tform) ! whole sample transformation.")
187 affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
188 call disp%show("affinity")
189 call disp%show( affinity )
190 call disp%show("do isam = 1, nsam")
191 call disp%show(" affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation")
192 call disp%show("end do")
193 do isam = 1, nsam
194 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
195 end do
196 call disp%show("affinity")
197 call disp%show( affinity )
198 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.")
199 affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
200 call disp%show("affinInv")
201 call disp%show( affinInv )
202 call disp%show("sample ! for comparison with affinInv.")
203 call disp%show( sample )
204 call disp%skip()
205 call disp%show("call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.")
206 call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
207 call disp%show("tform")
208 call disp%show( tform )
209 call disp%show("affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.")
210 affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
211 call disp%show("affinity")
212 call disp%show( affinity )
213 call disp%show("do isam = 1, nsam")
214 call disp%show(" affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation")
215 call disp%show("end do")
216 do isam = 1, nsam
217 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
218 end do
219 call disp%show("affinity")
220 call disp%show( affinity )
221 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.")
222 affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
223 call disp%show("affinInv")
224 call disp%show( affinInv )
225 call disp%show("sample ! for comparison with affinInv.")
226 call disp%show( sample )
227 call disp%skip()
228 call disp%show("affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.")
229 affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
230 call disp%show("affinity")
231 call disp%show( affinity )
232 call disp%show("do isam = 1, nsam")
233 call disp%show(" affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation")
234 call disp%show("end do")
235 do isam = 1, nsam
236 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
237 end do
238 call disp%show("affinity")
239 call disp%show( affinity )
240 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.")
241 affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
242 call disp%show("affinInv")
243 call disp%show( affinInv )
244 call disp%show("sample ! for comparison with affinInv.")
245 call disp%show( sample )
246 call disp%skip()
247 call disp%show("affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.")
248 affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
249 call disp%show("affinity")
250 call disp%show( affinity )
251 call disp%show("do isam = 1, nsam")
252 call disp%show(" affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation")
253 call disp%show("end do")
254 do isam = 1, nsam
255 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
256 end do
257 call disp%show("affinity")
258 call disp%show( affinity )
259 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.")
260 affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
261 call disp%show("affinInv")
262 call disp%show( affinInv )
263 call disp%show("sample ! for comparison with affinInv.")
264 call disp%show( sample )
265 call disp%skip()
266 call disp%show("affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.")
267 affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
268 call disp%show("affinity")
269 call disp%show( affinity )
270 call disp%show("do isam = 1, nsam")
271 call disp%show(" affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation")
272 call disp%show("end do")
273 do isam = 1, nsam
274 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
275 end do
276 call disp%show("affinity")
277 call disp%show( affinity )
278 call disp%show("affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.")
279 affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
280 call disp%show("affinInv")
281 call disp%show( affinInv )
282 call disp%show("sample ! for comparison with affinInv.")
283 call disp%show( sample )
284 end do
285 end block
286
287 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
288 ! Output an example rand array for visualization.
289 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
290
291 block
292 use pm_mathConst, only: PI
293 use pm_container, only: css_type
294 use pm_matrixInit, only: uppLowDia
295 use pm_matrixInit, only: getMatInit
296 use pm_arraySpace, only: getLinSpace
297 use pm_sampleAffinity, only: genrecmat
298 use pm_io, only: getErrTableWrite, trans
299 use pm_distUnifEll, only: getUnifEllRand, uppDia
300 use pm_sampleShift, only: getShifted
301 use pm_matrixInv, only: getMatInv
302 use pm_kind, only: RKG => RKS
303 type(css_type), allocatable :: shapes(:), tforms(:)
304 integer(IK) :: ishape, itform
305 integer(IK), parameter :: nsam = 1000, ndim = 2, dim = 2
306 real(RKG), allocatable :: angle, tlate(:), cov(:,:), tform(:,:), sample(:,:), affinity(:,:), affinInv(:,:)
307 tforms = [css_type(SK_"warp"), css_type(SK_"rotation")]
308 shapes = [css_type(SK_"circle"), css_type(SK_"square")]
309 tlate = getUnifRand(0, 5, ndim)
310 do ishape = 1, size(shapes)
311 if (shapes(ishape)%val == "circle") then
312 sample = getUnifEllRand(getMatInit([ndim, ndim], uppLowDia, 0., 0., getLinSpace(1., real(ndim), ndim)), uppDia, nsam)
313 elseif (shapes(ishape)%val == "square") then
314 sample = getUnifRand(0., 1., ndim, nsam)
315 else
316 error stop "Unrecognized shape."
317 end if
318 do itform = 1, size(tforms)
319 if (tforms(itform)%val == "rotation") then
320 angle = PI / 4
321 tform = reshape([cos(angle), -sin(angle), sin(angle), cos(angle)], [ndim, ndim])
322 affinity = getAffinity(sample, dim, tform, genrecmat, tlate)
323 elseif (tforms(itform)%val == "warp") then
324 cov = reshape([1., .5, .5, 1.], shape = [ndim, ndim]) ! getCovRand(1., ndim)
325 tform = getMatChol(cov, lowDia)
326 affinity = getAffinity(sample, dim, tform, lowerDiag, tlate)
327 else
328 error stop "Unrecognized tform."
329 end if
330 affinInv = getAffinity(getShifted(affinity, dim, -tlate), dim, getMatInv(tform))
331 if (0 /= getErrTableWrite(SK_"getAffinity."//tforms(itform)%val//SK_"."//shapes(ishape)%val//SK_".sample.txt", sample, trans)) error stop "Failed table-write."
332 if (0 /= getErrTableWrite(SK_"getAffinity."//tforms(itform)%val//SK_"."//shapes(ishape)%val//SK_".affinity.txt", affinity, trans)) error stop "Failed table-write."
333 if (0 /= getErrTableWrite(SK_"getAffinity."//tforms(itform)%val//SK_"."//shapes(ishape)%val//SK_".affinInv.txt", affinInv, trans)) error stop "Failed table-write."
334 end do
335 end do
336 end block
337
338end program example
Generate and return an array of the specified rank and shape of arbitrary intrinsic type and kind wit...
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Generate and return a random positive-definite (correlation or covariance) matrix using the Gram meth...
Definition: pm_distCov.F90:394
Generate and return a (collection) of random vector(s) of size ndim from the ndim-dimensional MultiVa...
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
Generate and return the iostat code resulting from writing the input table of rank 1 or 2 to the spec...
Definition: pm_io.F90:5940
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 upper or the lower Cholesky factorization of the input symmetric positive-def...
Generate and return a matrix of shape (shape(1), shape(2)) with the upper/lower triangle and the diag...
Set the upper/lower triangle and the diagonal elements of the input matrix of arbitrary shape (:,...
Generate and return the full inverse of an input matrix of general or triangular form directly or thr...
Generate a sample of shape (nsam), or (ndim, nsam) or (nsam, ndim) that is shifted by the specified i...
This module contains procedures and generic interfaces for convenient allocation and filling of array...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains procedures and generic interfaces for generating arrays with linear or logarithm...
This module contains the derived types for generating allocatable containers of scalar,...
This module contains classes and procedures for generating random matrices distributed on the space o...
Definition: pm_distCov.F90:72
This module contains classes and procedures for computing various statistical quantities related to t...
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 LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
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 relevant mathematical constants.
real(RKB), parameter PI
The scalar real constant of kind with highest available precision RKB representing the irrational num...
This module contains procedures and generic interfaces for computing the Cholesky factorization of po...
This module contains procedures and generic interfaces relevant to generating and initializing matric...
This module contains abstract and concrete derived types and procedures related to the inversion of s...
This module contains classes and procedures for shifting univariate or multivariate samples by arbitr...
This is the css_type type for generating instances of container of scalar of string objects.
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!Transform 2D sample along the second dimension.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
8[dim, ndim, nsam]
9+2, +2, +4
10tlate = getUnifRand(-5, +5, ndim)
11tlate
12-4.00000000, +2.00000000
13sample = getUnifRand(-9, +9, ndim, nsam)
14sample
15+9.00000000, +8.00000000, +1.00000000, +6.00000000
16+5.00000000, +1.00000000, -6.00000000, +3.00000000
17call setResized(affinity, shape(sample, IK))
18call setResized(affinInv, shape(sample, IK))
19tform = getCovRand(1., ndim)
20tform
21+1.00000000, +0.895324886
22+0.895324886, +1.00000000
23affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
24transpose(affinity)
25+13.4766245, +13.0579243
26+8.89532471, +8.16259956
27-4.37194920, -5.10467529
28+8.68597412, +8.37194920
29do isam = 1, nsam
30 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
31end do
32transpose(affinity)
33+13.4766245, +13.0579243
34+8.89532471, +8.16259956
35-4.37194920, -5.10467529
36+8.68597412, +8.37194920
37affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
38transpose(affinInv)
39+8.99999619, +5.00000381
40+8.00000000, +1.00000381
41+1.00000000, -6.00000191
42+5.99999619, +3.00000000
43transpose(sample) ! for comparison with affinInv.
44+9.00000000, +5.00000000
45+8.00000000, +1.00000000
46+1.00000000, -6.00000000
47+6.00000000, +3.00000000
48
49call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
50tform
51+1.00000000, +0.895324886
52+0.00000000, +1.00000000
53affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
54transpose(affinity)
55+13.4766245, +5.00000000
56+8.89532471, +1.00000000
57-4.37194920, -6.00000000
58+8.68597412, +3.00000000
59do isam = 1, nsam
60 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
61end do
62transpose(affinity)
63+13.4766245, +5.00000000
64+8.89532471, +1.00000000
65-4.37194920, -6.00000000
66+8.68597412, +3.00000000
67affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
68transpose(affinInv)
69+9.00000000, +5.00000000
70+8.00000000, +1.00000000
71+1.00000000, -6.00000000
72+5.99999952, +3.00000000
73transpose(sample) ! for comparison with affinInv.
74+9.00000000, +5.00000000
75+8.00000000, +1.00000000
76+1.00000000, -6.00000000
77+6.00000000, +3.00000000
78
79affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
80transpose(affinity)
81+13.4766245, +5.00000000
82+8.89532471, +1.00000000
83-4.37194920, -6.00000000
84+8.68597412, +3.00000000
85do isam = 1, nsam
86 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
87end do
88transpose(affinity)
89+13.4766245, +5.00000000
90+8.89532471, +1.00000000
91-4.37194920, -6.00000000
92+8.68597412, +3.00000000
93affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
94transpose(affinInv)
95+9.00000000, +5.00000000
96+8.00000000, +1.00000000
97+1.00000000, -6.00000000
98+5.99999952, +3.00000000
99transpose(sample) ! for comparison with affinInv.
100+9.00000000, +5.00000000
101+8.00000000, +1.00000000
102+1.00000000, -6.00000000
103+6.00000000, +3.00000000
104
105affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
106transpose(affinity)
107+9.00000000, +5.00000000
108+8.00000000, +1.00000000
109+1.00000000, -6.00000000
110+6.00000000, +3.00000000
111do isam = 1, nsam
112 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
113end do
114transpose(affinity)
115+9.00000000, +5.00000000
116+8.00000000, +1.00000000
117+1.00000000, -6.00000000
118+6.00000000, +3.00000000
119affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
120transpose(affinInv)
121+9.00000000, +5.00000000
122+8.00000000, +1.00000000
123+1.00000000, -6.00000000
124+6.00000000, +3.00000000
125transpose(sample) ! for comparison with affinInv.
126+9.00000000, +5.00000000
127+8.00000000, +1.00000000
128+1.00000000, -6.00000000
129+6.00000000, +3.00000000
130
131affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
132transpose(affinity)
133+9.00000000, +5.00000000
134+8.00000000, +1.00000000
135+1.00000000, -6.00000000
136+6.00000000, +3.00000000
137do isam = 1, nsam
138 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
139end do
140transpose(affinity)
141+9.00000000, +5.00000000
142+8.00000000, +1.00000000
143+1.00000000, -6.00000000
144+6.00000000, +3.00000000
145affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
146transpose(affinInv)
147+9.00000000, +5.00000000
148+8.00000000, +1.00000000
149+1.00000000, -6.00000000
150+6.00000000, +3.00000000
151transpose(sample) ! for comparison with affinInv.
152+9.00000000, +5.00000000
153+8.00000000, +1.00000000
154+1.00000000, -6.00000000
155+6.00000000, +3.00000000
156
157dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
158[dim, ndim, nsam]
159+2, +1, +4
160tlate = getUnifRand(-5, +5, ndim)
161tlate
162+5.00000000
163sample = getUnifRand(-9, +9, ndim, nsam)
164sample
165+4.00000000, +2.00000000, +9.00000000, +1.00000000
166call setResized(affinity, shape(sample, IK))
167call setResized(affinInv, shape(sample, IK))
168tform = getCovRand(1., ndim)
169tform
170+1.00000000
171affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
172transpose(affinity)
173+4.00000000
174+2.00000000
175+9.00000000
176+1.00000000
177do isam = 1, nsam
178 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
179end do
180transpose(affinity)
181+4.00000000
182+2.00000000
183+9.00000000
184+1.00000000
185affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
186transpose(affinInv)
187+4.00000000
188+2.00000000
189+9.00000000
190+1.00000000
191transpose(sample) ! for comparison with affinInv.
192+4.00000000
193+2.00000000
194+9.00000000
195+1.00000000
196
197call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
198tform
199+1.00000000
200affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
201transpose(affinity)
202+4.00000000
203+2.00000000
204+9.00000000
205+1.00000000
206do isam = 1, nsam
207 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
208end do
209transpose(affinity)
210+4.00000000
211+2.00000000
212+9.00000000
213+1.00000000
214affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
215transpose(affinInv)
216+4.00000000
217+2.00000000
218+9.00000000
219+1.00000000
220transpose(sample) ! for comparison with affinInv.
221+4.00000000
222+2.00000000
223+9.00000000
224+1.00000000
225
226affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
227transpose(affinity)
228+4.00000000
229+2.00000000
230+9.00000000
231+1.00000000
232do isam = 1, nsam
233 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
234end do
235transpose(affinity)
236+4.00000000
237+2.00000000
238+9.00000000
239+1.00000000
240affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
241transpose(affinInv)
242+4.00000000
243+2.00000000
244+9.00000000
245+1.00000000
246transpose(sample) ! for comparison with affinInv.
247+4.00000000
248+2.00000000
249+9.00000000
250+1.00000000
251
252affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
253transpose(affinity)
254+4.00000000
255+2.00000000
256+9.00000000
257+1.00000000
258do isam = 1, nsam
259 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
260end do
261transpose(affinity)
262+4.00000000
263+2.00000000
264+9.00000000
265+1.00000000
266affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
267transpose(affinInv)
268+4.00000000
269+2.00000000
270+9.00000000
271+1.00000000
272transpose(sample) ! for comparison with affinInv.
273+4.00000000
274+2.00000000
275+9.00000000
276+1.00000000
277
278affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
279transpose(affinity)
280+4.00000000
281+2.00000000
282+9.00000000
283+1.00000000
284do isam = 1, nsam
285 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
286end do
287transpose(affinity)
288+4.00000000
289+2.00000000
290+9.00000000
291+1.00000000
292affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
293transpose(affinInv)
294+4.00000000
295+2.00000000
296+9.00000000
297+1.00000000
298transpose(sample) ! for comparison with affinInv.
299+4.00000000
300+2.00000000
301+9.00000000
302+1.00000000
303
304dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
305[dim, ndim, nsam]
306+2, +1, +4
307tlate = getUnifRand(-5, +5, ndim)
308tlate
309-3.00000000
310sample = getUnifRand(-9, +9, ndim, nsam)
311sample
312+5.00000000, +5.00000000, -8.00000000, -7.00000000
313call setResized(affinity, shape(sample, IK))
314call setResized(affinInv, shape(sample, IK))
315tform = getCovRand(1., ndim)
316tform
317+1.00000000
318affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
319transpose(affinity)
320+5.00000000
321+5.00000000
322-8.00000000
323-7.00000000
324do isam = 1, nsam
325 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
326end do
327transpose(affinity)
328+5.00000000
329+5.00000000
330-8.00000000
331-7.00000000
332affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
333transpose(affinInv)
334+5.00000000
335+5.00000000
336-8.00000000
337-7.00000000
338transpose(sample) ! for comparison with affinInv.
339+5.00000000
340+5.00000000
341-8.00000000
342-7.00000000
343
344call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
345tform
346+1.00000000
347affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
348transpose(affinity)
349+5.00000000
350+5.00000000
351-8.00000000
352-7.00000000
353do isam = 1, nsam
354 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
355end do
356transpose(affinity)
357+5.00000000
358+5.00000000
359-8.00000000
360-7.00000000
361affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
362transpose(affinInv)
363+5.00000000
364+5.00000000
365-8.00000000
366-7.00000000
367transpose(sample) ! for comparison with affinInv.
368+5.00000000
369+5.00000000
370-8.00000000
371-7.00000000
372
373affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
374transpose(affinity)
375+5.00000000
376+5.00000000
377-8.00000000
378-7.00000000
379do isam = 1, nsam
380 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
381end do
382transpose(affinity)
383+5.00000000
384+5.00000000
385-8.00000000
386-7.00000000
387affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
388transpose(affinInv)
389+5.00000000
390+5.00000000
391-8.00000000
392-7.00000000
393transpose(sample) ! for comparison with affinInv.
394+5.00000000
395+5.00000000
396-8.00000000
397-7.00000000
398
399affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
400transpose(affinity)
401+5.00000000
402+5.00000000
403-8.00000000
404-7.00000000
405do isam = 1, nsam
406 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
407end do
408transpose(affinity)
409+5.00000000
410+5.00000000
411-8.00000000
412-7.00000000
413affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
414transpose(affinInv)
415+5.00000000
416+5.00000000
417-8.00000000
418-7.00000000
419transpose(sample) ! for comparison with affinInv.
420+5.00000000
421+5.00000000
422-8.00000000
423-7.00000000
424
425affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
426transpose(affinity)
427+5.00000000
428+5.00000000
429-8.00000000
430-7.00000000
431do isam = 1, nsam
432 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
433end do
434transpose(affinity)
435+5.00000000
436+5.00000000
437-8.00000000
438-7.00000000
439affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
440transpose(affinInv)
441+5.00000000
442+5.00000000
443-8.00000000
444-7.00000000
445transpose(sample) ! for comparison with affinInv.
446+5.00000000
447+5.00000000
448-8.00000000
449-7.00000000
450
451dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
452[dim, ndim, nsam]
453+2, +3, +5
454tlate = getUnifRand(-5, +5, ndim)
455tlate
456-4.00000000, +5.00000000, +5.00000000
457sample = getUnifRand(-9, +9, ndim, nsam)
458sample
459-9.00000000, -4.00000000, +0.00000000, +3.00000000, +4.00000000
460+3.00000000, +7.00000000, +6.00000000, -3.00000000, +3.00000000
461+6.00000000, +4.00000000, -3.00000000, +8.00000000, -5.00000000
462call setResized(affinity, shape(sample, IK))
463call setResized(affinInv, shape(sample, IK))
464tform = getCovRand(1., ndim)
465tform
466+1.00000000, +0.992524028, -0.746748388
467+0.992524028, +1.00000000, -0.820488811
468-0.746748388, -0.820488811, +1.00000000
469affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
470transpose(affinity)
471-10.5029182, -10.8556490, +10.2592697
472-0.393254757E-1, -0.252051353, +1.24357200
473+8.19538879, +8.46146584, -7.92293262
474-5.95155907, -6.58633852, +8.22122097
475+10.7113142, +11.0725403, -10.4484596
476do isam = 1, nsam
477 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
478end do
479transpose(affinity)
480-10.5029182, -10.8556490, +10.2592697
481-0.393254757E-1, -0.252051353, +1.24357200
482+8.19538879, +8.46146584, -7.92293262
483-5.95155907, -6.58633852, +8.22122097
484+10.7113142, +11.0725403, -10.4484596
485affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
486transpose(affinInv)
487-9.00146484, +3.00122070, +6.00006104
488-4.00045776, +7.00057983, +4.00009918
489+0.610351562E-3, +5.99877930, -2.99990845
490+2.99926758, -2.99975586, +7.99990845
491+4.00097656, +2.99877930, -5.00000000
492transpose(sample) ! for comparison with affinInv.
493-9.00000000, +3.00000000, +6.00000000
494-4.00000000, +7.00000000, +4.00000000
495+0.00000000, +6.00000000, -3.00000000
496+3.00000000, -3.00000000, +8.00000000
497+4.00000000, +3.00000000, -5.00000000
498
499call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
500tform
501+1.00000000, +0.992524028, -0.746748388
502+0.00000000, +1.00000000, -0.820488811
503+0.00000000, +0.00000000, +1.00000000
504affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
505transpose(affinity)
506-10.5029182, -1.92293262, +6.00000000
507-0.393254757E-1, +3.71804476, +4.00000000
508+8.19538879, +8.46146584, -3.00000000
509-5.95155907, -9.56391048, +8.00000000
510+10.7113142, +7.10244417, -5.00000000
511do isam = 1, nsam
512 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
513end do
514transpose(affinity)
515-10.5029182, -1.92293262, +6.00000000
516-0.393254757E-1, +3.71804476, +4.00000000
517+8.19538879, +8.46146584, -3.00000000
518-5.95155907, -9.56391048, +8.00000000
519+10.7113142, +7.10244417, -5.00000000
520affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
521transpose(affinInv)
522-9.00000000, +3.00000000, +6.00000000
523-4.00000000, +7.00000000, +4.00000000
524+0.00000000, +5.99999952, -3.00000000
525+3.00000000, -3.00000000, +8.00000000
526+4.00000000, +3.00000000, -5.00000000
527transpose(sample) ! for comparison with affinInv.
528-9.00000000, +3.00000000, +6.00000000
529-4.00000000, +7.00000000, +4.00000000
530+0.00000000, +6.00000000, -3.00000000
531+3.00000000, -3.00000000, +8.00000000
532+4.00000000, +3.00000000, -5.00000000
533
534affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
535transpose(affinity)
536-10.5029182, -1.92293262, +6.00000000
537-0.393254757E-1, +3.71804476, +4.00000000
538+8.19538879, +8.46146584, -3.00000000
539-5.95155907, -9.56391048, +8.00000000
540+10.7113142, +7.10244417, -5.00000000
541do isam = 1, nsam
542 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
543end do
544transpose(affinity)
545-10.5029182, -1.92293262, +6.00000000
546-0.393254757E-1, +3.71804476, +4.00000000
547+8.19538879, +8.46146584, -3.00000000
548-5.95155907, -9.56391048, +8.00000000
549+10.7113142, +7.10244417, -5.00000000
550affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
551transpose(affinInv)
552-9.00000000, +3.00000000, +6.00000000
553-4.00000000, +7.00000000, +4.00000000
554+0.00000000, +5.99999952, -3.00000000
555+3.00000000, -3.00000000, +8.00000000
556+4.00000000, +3.00000000, -5.00000000
557transpose(sample) ! for comparison with affinInv.
558-9.00000000, +3.00000000, +6.00000000
559-4.00000000, +7.00000000, +4.00000000
560+0.00000000, +6.00000000, -3.00000000
561+3.00000000, -3.00000000, +8.00000000
562+4.00000000, +3.00000000, -5.00000000
563
564affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
565transpose(affinity)
566-9.00000000, +3.00000000, +6.00000000
567-4.00000000, +7.00000000, +4.00000000
568+0.00000000, +6.00000000, -3.00000000
569+3.00000000, -3.00000000, +8.00000000
570+4.00000000, +3.00000000, -5.00000000
571do isam = 1, nsam
572 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
573end do
574transpose(affinity)
575-9.00000000, +3.00000000, +6.00000000
576-4.00000000, +7.00000000, +4.00000000
577+0.00000000, +6.00000000, -3.00000000
578+3.00000000, -3.00000000, +8.00000000
579+4.00000000, +3.00000000, -5.00000000
580affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
581transpose(affinInv)
582-9.00000000, +3.00000000, +6.00000000
583-4.00000000, +7.00000000, +4.00000000
584+0.00000000, +6.00000000, -3.00000000
585+3.00000000, -3.00000000, +8.00000000
586+4.00000000, +3.00000000, -5.00000000
587transpose(sample) ! for comparison with affinInv.
588-9.00000000, +3.00000000, +6.00000000
589-4.00000000, +7.00000000, +4.00000000
590+0.00000000, +6.00000000, -3.00000000
591+3.00000000, -3.00000000, +8.00000000
592+4.00000000, +3.00000000, -5.00000000
593
594affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
595transpose(affinity)
596-9.00000000, +3.00000000, +6.00000000
597-4.00000000, +7.00000000, +4.00000000
598+0.00000000, +6.00000000, -3.00000000
599+3.00000000, -3.00000000, +8.00000000
600+4.00000000, +3.00000000, -5.00000000
601do isam = 1, nsam
602 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
603end do
604transpose(affinity)
605-9.00000000, +3.00000000, +6.00000000
606-4.00000000, +7.00000000, +4.00000000
607+0.00000000, +6.00000000, -3.00000000
608+3.00000000, -3.00000000, +8.00000000
609+4.00000000, +3.00000000, -5.00000000
610affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
611transpose(affinInv)
612-9.00000000, +3.00000000, +6.00000000
613-4.00000000, +7.00000000, +4.00000000
614+0.00000000, +6.00000000, -3.00000000
615+3.00000000, -3.00000000, +8.00000000
616+4.00000000, +3.00000000, -5.00000000
617transpose(sample) ! for comparison with affinInv.
618-9.00000000, +3.00000000, +6.00000000
619-4.00000000, +7.00000000, +4.00000000
620+0.00000000, +6.00000000, -3.00000000
621+3.00000000, -3.00000000, +8.00000000
622+4.00000000, +3.00000000, -5.00000000
623
624dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
625[dim, ndim, nsam]
626+2, +2, +5
627tlate = getUnifRand(-5, +5, ndim)
628tlate
629+2.00000000, +1.00000000
630sample = getUnifRand(-9, +9, ndim, nsam)
631sample
632+4.00000000, +0.00000000, +4.00000000, -5.00000000, +5.00000000
633+2.00000000, +5.00000000, -9.00000000, +1.00000000, +9.00000000
634call setResized(affinity, shape(sample, IK))
635call setResized(affinInv, shape(sample, IK))
636tform = getCovRand(1., ndim)
637tform
638+1.00000000, -0.592477202
639-0.592477202, +1.00000000
640affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
641transpose(affinity)
642+2.81504560, -0.369908810
643-2.96238613, +5.00000000
644+9.33229446, -11.3699093
645-5.59247732, +3.96238613
646-0.332294941, +6.03761387
647do isam = 1, nsam
648 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
649end do
650transpose(affinity)
651+2.81504560, -0.369908810
652-2.96238613, +5.00000000
653+9.33229446, -11.3699093
654-5.59247732, +3.96238613
655-0.332294941, +6.03761387
656affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
657transpose(affinInv)
658+4.00000000, +2.00000000
659+0.00000000, +5.00000000
660+3.99999809, -9.00000191
661-5.00000000, +1.00000048
662+4.99999952, +9.00000000
663transpose(sample) ! for comparison with affinInv.
664+4.00000000, +2.00000000
665+0.00000000, +5.00000000
666+4.00000000, -9.00000000
667-5.00000000, +1.00000000
668+5.00000000, +9.00000000
669
670call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
671tform
672+1.00000000, -0.592477202
673+0.00000000, +1.00000000
674affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
675transpose(affinity)
676+2.81504560, +2.00000000
677-2.96238613, +5.00000000
678+9.33229446, -9.00000000
679-5.59247732, +1.00000000
680-0.332294941, +9.00000000
681do isam = 1, nsam
682 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
683end do
684transpose(affinity)
685+2.81504560, +2.00000000
686-2.96238613, +5.00000000
687+9.33229446, -9.00000000
688-5.59247732, +1.00000000
689-0.332294941, +9.00000000
690affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
691transpose(affinInv)
692+4.00000000, +2.00000000
693+0.00000000, +5.00000000
694+3.99999952, -9.00000000
695-5.00000000, +1.00000000
696+5.00000000, +9.00000000
697transpose(sample) ! for comparison with affinInv.
698+4.00000000, +2.00000000
699+0.00000000, +5.00000000
700+4.00000000, -9.00000000
701-5.00000000, +1.00000000
702+5.00000000, +9.00000000
703
704affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
705transpose(affinity)
706+2.81504560, +2.00000000
707-2.96238613, +5.00000000
708+9.33229446, -9.00000000
709-5.59247732, +1.00000000
710-0.332294941, +9.00000000
711do isam = 1, nsam
712 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
713end do
714transpose(affinity)
715+2.81504560, +2.00000000
716-2.96238613, +5.00000000
717+9.33229446, -9.00000000
718-5.59247732, +1.00000000
719-0.332294941, +9.00000000
720affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
721transpose(affinInv)
722+4.00000000, +2.00000000
723+0.00000000, +5.00000000
724+3.99999952, -9.00000000
725-5.00000000, +1.00000000
726+5.00000000, +9.00000000
727transpose(sample) ! for comparison with affinInv.
728+4.00000000, +2.00000000
729+0.00000000, +5.00000000
730+4.00000000, -9.00000000
731-5.00000000, +1.00000000
732+5.00000000, +9.00000000
733
734affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
735transpose(affinity)
736+4.00000000, +2.00000000
737+0.00000000, +5.00000000
738+4.00000000, -9.00000000
739-5.00000000, +1.00000000
740+5.00000000, +9.00000000
741do isam = 1, nsam
742 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
743end do
744transpose(affinity)
745+4.00000000, +2.00000000
746+0.00000000, +5.00000000
747+4.00000000, -9.00000000
748-5.00000000, +1.00000000
749+5.00000000, +9.00000000
750affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
751transpose(affinInv)
752+4.00000000, +2.00000000
753+0.00000000, +5.00000000
754+4.00000000, -9.00000000
755-5.00000000, +1.00000000
756+5.00000000, +9.00000000
757transpose(sample) ! for comparison with affinInv.
758+4.00000000, +2.00000000
759+0.00000000, +5.00000000
760+4.00000000, -9.00000000
761-5.00000000, +1.00000000
762+5.00000000, +9.00000000
763
764affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
765transpose(affinity)
766+4.00000000, +2.00000000
767+0.00000000, +5.00000000
768+4.00000000, -9.00000000
769-5.00000000, +1.00000000
770+5.00000000, +9.00000000
771do isam = 1, nsam
772 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
773end do
774transpose(affinity)
775+4.00000000, +2.00000000
776+0.00000000, +5.00000000
777+4.00000000, -9.00000000
778-5.00000000, +1.00000000
779+5.00000000, +9.00000000
780affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
781transpose(affinInv)
782+4.00000000, +2.00000000
783+0.00000000, +5.00000000
784+4.00000000, -9.00000000
785-5.00000000, +1.00000000
786+5.00000000, +9.00000000
787transpose(sample) ! for comparison with affinInv.
788+4.00000000, +2.00000000
789+0.00000000, +5.00000000
790+4.00000000, -9.00000000
791-5.00000000, +1.00000000
792+5.00000000, +9.00000000
793
794!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
795!Transform 2D sample along the first dimension.
796!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
797
798
799dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
800[dim, ndim, nsam]
801+1, +3, +4
802tlate = getUnifRand(-5, +5, ndim)
803tlate
804-1.00000000, +3.00000000, +4.00000000
805sample = getUnifRand(-9, +9, nsam, ndim)
806sample
807+9.00000000, -6.00000000, +6.00000000
808-9.00000000, -3.00000000, +6.00000000
809+6.00000000, +4.00000000, +7.00000000
810+0.00000000, +4.00000000, -1.00000000
811call setResized(affinity, shape(sample, IK))
812call setResized(affinInv, shape(sample, IK))
813tform = getCovRand(1., ndim)
814tform
815+1.00000000, +0.840149701, +0.516399145
816+0.840149701, +1.00000000, -0.167174637E-1
817+0.516399145, -0.167174637E-1, +1.00000000
818affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
819affinity
820+7.05749655, +1.46104264, +10.7478971
821-8.42205429, -10.6616526, +1.40255976
822+12.9753923, +8.92387581, +10.0315247
823+2.84419966, +4.01671743, -1.06686985
824do isam = 1, nsam
825 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
826end do
827affinity
828+7.05749655, +1.46104264, +10.7478971
829-8.42205429, -10.6616526, +1.40255976
830+12.9753923, +8.92387581, +10.0315247
831+2.84419966, +4.01671743, -1.06686985
832affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
833affinInv
834+8.99993896, -5.99993896, +6.00004578
835-8.99985886, -3.00009155, +5.99995041
836+5.99987793, +4.00009155, +7.00004578
837-0.152587891E-4, +4.00001526, -0.999992371
838sample ! for comparison with affinInv.
839+9.00000000, -6.00000000, +6.00000000
840-9.00000000, -3.00000000, +6.00000000
841+6.00000000, +4.00000000, +7.00000000
842+0.00000000, +4.00000000, -1.00000000
843
844call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
845tform
846+1.00000000, +0.840149701, +0.516399145
847+0.00000000, +1.00000000, -0.167174637E-1
848+0.00000000, +0.00000000, +1.00000000
849affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
850affinity
851+7.05749655, -6.10030460, +6.00000000
852-8.42205429, -3.10030484, +6.00000000
853+12.9753933, +3.88297772, +7.00000000
854+2.84419966, +4.01671743, -1.00000000
855do isam = 1, nsam
856 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
857end do
858affinity
859+7.05749655, -6.10030460, +6.00000000
860-8.42205429, -3.10030484, +6.00000000
861+12.9753933, +3.88297772, +7.00000000
862+2.84419966, +4.01671743, -1.00000000
863affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
864affinInv
865+9.00000000, -6.00000000, +6.00000000
866-9.00000000, -3.00000000, +6.00000000
867+6.00000048, +4.00000000, +7.00000000
868+0.00000000, +4.00000000, -1.00000000
869sample ! for comparison with affinInv.
870+9.00000000, -6.00000000, +6.00000000
871-9.00000000, -3.00000000, +6.00000000
872+6.00000000, +4.00000000, +7.00000000
873+0.00000000, +4.00000000, -1.00000000
874
875affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
876affinity
877+7.05749655, -6.10030460, +6.00000000
878-8.42205429, -3.10030484, +6.00000000
879+12.9753933, +3.88297772, +7.00000000
880+2.84419966, +4.01671743, -1.00000000
881do isam = 1, nsam
882 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
883end do
884affinity
885+7.05749655, -6.10030460, +6.00000000
886-8.42205429, -3.10030484, +6.00000000
887+12.9753933, +3.88297772, +7.00000000
888+2.84419966, +4.01671743, -1.00000000
889affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
890affinInv
891+9.00000000, -6.00000000, +6.00000000
892-9.00000000, -3.00000000, +6.00000000
893+6.00000048, +4.00000000, +7.00000000
894+0.00000000, +4.00000000, -1.00000000
895sample ! for comparison with affinInv.
896+9.00000000, -6.00000000, +6.00000000
897-9.00000000, -3.00000000, +6.00000000
898+6.00000000, +4.00000000, +7.00000000
899+0.00000000, +4.00000000, -1.00000000
900
901affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
902affinity
903+9.00000000, -6.00000000, +6.00000000
904-9.00000000, -3.00000000, +6.00000000
905+6.00000000, +4.00000000, +7.00000000
906+0.00000000, +4.00000000, -1.00000000
907do isam = 1, nsam
908 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
909end do
910affinity
911+9.00000000, -6.00000000, +6.00000000
912-9.00000000, -3.00000000, +6.00000000
913+6.00000000, +4.00000000, +7.00000000
914+0.00000000, +4.00000000, -1.00000000
915affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
916affinInv
917+9.00000000, -6.00000000, +6.00000000
918-9.00000000, -3.00000000, +6.00000000
919+6.00000000, +4.00000000, +7.00000000
920+0.00000000, +4.00000000, -1.00000000
921sample ! for comparison with affinInv.
922+9.00000000, -6.00000000, +6.00000000
923-9.00000000, -3.00000000, +6.00000000
924+6.00000000, +4.00000000, +7.00000000
925+0.00000000, +4.00000000, -1.00000000
926
927affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
928affinity
929+9.00000000, -6.00000000, +6.00000000
930-9.00000000, -3.00000000, +6.00000000
931+6.00000000, +4.00000000, +7.00000000
932+0.00000000, +4.00000000, -1.00000000
933do isam = 1, nsam
934 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
935end do
936affinity
937+9.00000000, -6.00000000, +6.00000000
938-9.00000000, -3.00000000, +6.00000000
939+6.00000000, +4.00000000, +7.00000000
940+0.00000000, +4.00000000, -1.00000000
941affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
942affinInv
943+9.00000000, -6.00000000, +6.00000000
944-9.00000000, -3.00000000, +6.00000000
945+6.00000000, +4.00000000, +7.00000000
946+0.00000000, +4.00000000, -1.00000000
947sample ! for comparison with affinInv.
948+9.00000000, -6.00000000, +6.00000000
949-9.00000000, -3.00000000, +6.00000000
950+6.00000000, +4.00000000, +7.00000000
951+0.00000000, +4.00000000, -1.00000000
952
953dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
954[dim, ndim, nsam]
955+1, +3, +5
956tlate = getUnifRand(-5, +5, ndim)
957tlate
958+3.00000000, -2.00000000, -2.00000000
959sample = getUnifRand(-9, +9, nsam, ndim)
960sample
961-4.00000000, -5.00000000, -2.00000000
962-7.00000000, +1.00000000, -6.00000000
963+5.00000000, +5.00000000, +5.00000000
964+7.00000000, -9.00000000, -3.00000000
965-1.00000000, -8.00000000, +9.00000000
966call setResized(affinity, shape(sample, IK))
967call setResized(affinInv, shape(sample, IK))
968tform = getCovRand(1., ndim)
969tform
970+1.00000000, -0.961070478, +0.609573305
971-0.961070478, +1.00000000, -0.767626822
972+0.609573305, -0.767626822, +1.00000000
973affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
974affinity
975-0.413794398, +0.379535556, -0.600159168
976-11.6185102, +12.3332539, -11.0346394
977+3.24251437, -3.64348626, +4.20973253
978+13.8209143, -13.4246130, +8.17565441
979+12.1747236, -13.9475708, +14.5314407
980do isam = 1, nsam
981 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
982end do
983affinity
984-0.413794398, +0.379535556, -0.600159168
985-11.6185102, +12.3332539, -11.0346394
986+3.24251437, -3.64348626, +4.20973253
987+13.8209143, -13.4246130, +8.17565441
988+12.1747236, -13.9475708, +14.5314407
989affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
990affinInv
991-4.00000477, -5.00000668, -2.00000238
992-6.99996948, +1.00003052, -6.00000763
993+5.00001144, +5.00000763, +5.00000763
994+6.99996185, -9.00006866, -2.99999619
995-1.00003052, -8.00006104, +9.00000000
996sample ! for comparison with affinInv.
997-4.00000000, -5.00000000, -2.00000000
998-7.00000000, +1.00000000, -6.00000000
999+5.00000000, +5.00000000, +5.00000000
1000+7.00000000, -9.00000000, -3.00000000
1001-1.00000000, -8.00000000, +9.00000000
1002
1003call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
1004tform
1005+1.00000000, -0.961070478, +0.609573305
1006+0.00000000, +1.00000000, -0.767626822
1007+0.00000000, +0.00000000, +1.00000000
1008affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
1009affinity
1010-0.413794518, -3.46474648, -2.00000000
1011-11.6185102, +5.60576105, -6.00000000
1012+3.24251437, +1.16186595, +5.00000000
1013+13.8209143, -6.69711971, -3.00000000
1014+12.1747236, -14.9086418, +9.00000000
1015do isam = 1, nsam
1016 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
1017end do
1018affinity
1019-0.413794518, -3.46474648, -2.00000000
1020-11.6185102, +5.60576105, -6.00000000
1021+3.24251437, +1.16186595, +5.00000000
1022+13.8209143, -6.69711971, -3.00000000
1023+12.1747236, -14.9086418, +9.00000000
1024affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
1025affinInv
1026-4.00000048, -5.00000000, -2.00000000
1027-7.00000000, +1.00000000, -6.00000000
1028+5.00000048, +5.00000000, +5.00000000
1029+6.99999952, -9.00000000, -3.00000000
1030-1.00000000, -8.00000000, +9.00000000
1031sample ! for comparison with affinInv.
1032-4.00000000, -5.00000000, -2.00000000
1033-7.00000000, +1.00000000, -6.00000000
1034+5.00000000, +5.00000000, +5.00000000
1035+7.00000000, -9.00000000, -3.00000000
1036-1.00000000, -8.00000000, +9.00000000
1037
1038affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
1039affinity
1040-0.413794518, -3.46474648, -2.00000000
1041-11.6185102, +5.60576105, -6.00000000
1042+3.24251437, +1.16186595, +5.00000000
1043+13.8209143, -6.69711971, -3.00000000
1044+12.1747236, -14.9086418, +9.00000000
1045do isam = 1, nsam
1046 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
1047end do
1048affinity
1049-0.413794518, -3.46474648, -2.00000000
1050-11.6185102, +5.60576105, -6.00000000
1051+3.24251437, +1.16186595, +5.00000000
1052+13.8209143, -6.69711971, -3.00000000
1053+12.1747236, -14.9086418, +9.00000000
1054affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
1055affinInv
1056-4.00000048, -5.00000000, -2.00000000
1057-7.00000000, +1.00000000, -6.00000000
1058+5.00000048, +5.00000000, +5.00000000
1059+6.99999952, -9.00000000, -3.00000000
1060-1.00000000, -8.00000000, +9.00000000
1061sample ! for comparison with affinInv.
1062-4.00000000, -5.00000000, -2.00000000
1063-7.00000000, +1.00000000, -6.00000000
1064+5.00000000, +5.00000000, +5.00000000
1065+7.00000000, -9.00000000, -3.00000000
1066-1.00000000, -8.00000000, +9.00000000
1067
1068affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
1069affinity
1070-4.00000000, -5.00000000, -2.00000000
1071-7.00000000, +1.00000000, -6.00000000
1072+5.00000000, +5.00000000, +5.00000000
1073+7.00000000, -9.00000000, -3.00000000
1074-1.00000000, -8.00000000, +9.00000000
1075do isam = 1, nsam
1076 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
1077end do
1078affinity
1079-4.00000000, -5.00000000, -2.00000000
1080-7.00000000, +1.00000000, -6.00000000
1081+5.00000000, +5.00000000, +5.00000000
1082+7.00000000, -9.00000000, -3.00000000
1083-1.00000000, -8.00000000, +9.00000000
1084affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
1085affinInv
1086-4.00000000, -5.00000000, -2.00000000
1087-7.00000000, +1.00000000, -6.00000000
1088+5.00000000, +5.00000000, +5.00000000
1089+7.00000000, -9.00000000, -3.00000000
1090-1.00000000, -8.00000000, +9.00000000
1091sample ! for comparison with affinInv.
1092-4.00000000, -5.00000000, -2.00000000
1093-7.00000000, +1.00000000, -6.00000000
1094+5.00000000, +5.00000000, +5.00000000
1095+7.00000000, -9.00000000, -3.00000000
1096-1.00000000, -8.00000000, +9.00000000
1097
1098affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
1099affinity
1100-4.00000000, -5.00000000, -2.00000000
1101-7.00000000, +1.00000000, -6.00000000
1102+5.00000000, +5.00000000, +5.00000000
1103+7.00000000, -9.00000000, -3.00000000
1104-1.00000000, -8.00000000, +9.00000000
1105do isam = 1, nsam
1106 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
1107end do
1108affinity
1109-4.00000000, -5.00000000, -2.00000000
1110-7.00000000, +1.00000000, -6.00000000
1111+5.00000000, +5.00000000, +5.00000000
1112+7.00000000, -9.00000000, -3.00000000
1113-1.00000000, -8.00000000, +9.00000000
1114affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
1115affinInv
1116-4.00000000, -5.00000000, -2.00000000
1117-7.00000000, +1.00000000, -6.00000000
1118+5.00000000, +5.00000000, +5.00000000
1119+7.00000000, -9.00000000, -3.00000000
1120-1.00000000, -8.00000000, +9.00000000
1121sample ! for comparison with affinInv.
1122-4.00000000, -5.00000000, -2.00000000
1123-7.00000000, +1.00000000, -6.00000000
1124+5.00000000, +5.00000000, +5.00000000
1125+7.00000000, -9.00000000, -3.00000000
1126-1.00000000, -8.00000000, +9.00000000
1127
1128dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
1129[dim, ndim, nsam]
1130+1, +3, +5
1131tlate = getUnifRand(-5, +5, ndim)
1132tlate
1133-4.00000000, +5.00000000, -5.00000000
1134sample = getUnifRand(-9, +9, nsam, ndim)
1135sample
1136-7.00000000, -4.00000000, +7.00000000
1137+4.00000000, -9.00000000, -7.00000000
1138-6.00000000, -1.00000000, +2.00000000
1139+0.00000000, +5.00000000, -8.00000000
1140-6.00000000, -2.00000000, -7.00000000
1141call setResized(affinity, shape(sample, IK))
1142call setResized(affinInv, shape(sample, IK))
1143tform = getCovRand(1., ndim)
1144tform
1145+1.00000000, +0.983967423, +0.568975270
1146+0.983967423, +1.00000000, +0.631563604
1147+0.568975270, +0.631563604, +1.00000000
1148affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
1149affinity
1150-6.95304203, -6.46682644, +0.490918636
1151-8.83853436, -9.48507500, -10.4081717
1152-5.84601688, -5.64067745, -2.04541540
1153+0.368034840, -0.525088310E-1, -4.84218216
1154-11.9507618, -12.3247499, -11.6769791
1155do isam = 1, nsam
1156 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
1157end do
1158affinity
1159-6.95304203, -6.46682644, +0.490918636
1160-8.83853436, -9.48507500, -10.4081717
1161-5.84601688, -5.64067745, -2.04541540
1162+0.368034840, -0.525088310E-1, -4.84218216
1163-11.9507618, -12.3247499, -11.6769791
1164affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
1165affinInv
1166-6.99997234, -4.00003242, +6.99999905
1167+3.99996567, -8.99995804, -7.00000572
1168-5.99998522, -1.00000668, +1.99999881
1169+0.476837158E-5, +4.99999237, -7.99999952
1170-5.99997330, -2.00002289, -7.00000191
1171sample ! for comparison with affinInv.
1172-7.00000000, -4.00000000, +7.00000000
1173+4.00000000, -9.00000000, -7.00000000
1174-6.00000000, -1.00000000, +2.00000000
1175+0.00000000, +5.00000000, -8.00000000
1176-6.00000000, -2.00000000, -7.00000000
1177
1178call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
1179tform
1180+1.00000000, +0.983967423, +0.568975270
1181+0.00000000, +1.00000000, +0.631563604
1182+0.00000000, +0.00000000, +1.00000000
1183affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
1184affinity
1185-6.95304298, +0.420945168, +7.00000000
1186-8.83853436, -13.4209452, -7.00000000
1187-5.84601688, +0.263127208, +2.00000000
1188+0.368034840, -0.525088310E-1, -8.00000000
1189-11.9507618, -6.42094517, -7.00000000
1190do isam = 1, nsam
1191 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
1192end do
1193affinity
1194-6.95304298, +0.420945168, +7.00000000
1195-8.83853436, -13.4209452, -7.00000000
1196-5.84601688, +0.263127208, +2.00000000
1197+0.368034840, -0.525088310E-1, -8.00000000
1198-11.9507618, -6.42094517, -7.00000000
1199affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
1200affinInv
1201-7.00000000, -4.00000000, +7.00000000
1202+3.99999905, -9.00000000, -7.00000000
1203-6.00000000, -1.00000000, +2.00000000
1204-0.238418579E-6, +5.00000000, -8.00000000
1205-6.00000048, -2.00000000, -7.00000000
1206sample ! for comparison with affinInv.
1207-7.00000000, -4.00000000, +7.00000000
1208+4.00000000, -9.00000000, -7.00000000
1209-6.00000000, -1.00000000, +2.00000000
1210+0.00000000, +5.00000000, -8.00000000
1211-6.00000000, -2.00000000, -7.00000000
1212
1213affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
1214affinity
1215-6.95304298, +0.420945168, +7.00000000
1216-8.83853436, -13.4209452, -7.00000000
1217-5.84601688, +0.263127208, +2.00000000
1218+0.368034840, -0.525088310E-1, -8.00000000
1219-11.9507618, -6.42094517, -7.00000000
1220do isam = 1, nsam
1221 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
1222end do
1223affinity
1224-6.95304298, +0.420945168, +7.00000000
1225-8.83853436, -13.4209452, -7.00000000
1226-5.84601688, +0.263127208, +2.00000000
1227+0.368034840, -0.525088310E-1, -8.00000000
1228-11.9507618, -6.42094517, -7.00000000
1229affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
1230affinInv
1231-7.00000000, -4.00000000, +7.00000000
1232+3.99999905, -9.00000000, -7.00000000
1233-6.00000000, -1.00000000, +2.00000000
1234-0.238418579E-6, +5.00000000, -8.00000000
1235-6.00000048, -2.00000000, -7.00000000
1236sample ! for comparison with affinInv.
1237-7.00000000, -4.00000000, +7.00000000
1238+4.00000000, -9.00000000, -7.00000000
1239-6.00000000, -1.00000000, +2.00000000
1240+0.00000000, +5.00000000, -8.00000000
1241-6.00000000, -2.00000000, -7.00000000
1242
1243affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
1244affinity
1245-7.00000000, -4.00000000, +7.00000000
1246+4.00000000, -9.00000000, -7.00000000
1247-6.00000000, -1.00000000, +2.00000000
1248+0.00000000, +5.00000000, -8.00000000
1249-6.00000000, -2.00000000, -7.00000000
1250do isam = 1, nsam
1251 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
1252end do
1253affinity
1254-7.00000000, -4.00000000, +7.00000000
1255+4.00000000, -9.00000000, -7.00000000
1256-6.00000000, -1.00000000, +2.00000000
1257+0.00000000, +5.00000000, -8.00000000
1258-6.00000000, -2.00000000, -7.00000000
1259affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
1260affinInv
1261-7.00000000, -4.00000000, +7.00000000
1262+4.00000000, -9.00000000, -7.00000000
1263-6.00000000, -1.00000000, +2.00000000
1264+0.00000000, +5.00000000, -8.00000000
1265-6.00000000, -2.00000000, -7.00000000
1266sample ! for comparison with affinInv.
1267-7.00000000, -4.00000000, +7.00000000
1268+4.00000000, -9.00000000, -7.00000000
1269-6.00000000, -1.00000000, +2.00000000
1270+0.00000000, +5.00000000, -8.00000000
1271-6.00000000, -2.00000000, -7.00000000
1272
1273affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
1274affinity
1275-7.00000000, -4.00000000, +7.00000000
1276+4.00000000, -9.00000000, -7.00000000
1277-6.00000000, -1.00000000, +2.00000000
1278+0.00000000, +5.00000000, -8.00000000
1279-6.00000000, -2.00000000, -7.00000000
1280do isam = 1, nsam
1281 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
1282end do
1283affinity
1284-7.00000000, -4.00000000, +7.00000000
1285+4.00000000, -9.00000000, -7.00000000
1286-6.00000000, -1.00000000, +2.00000000
1287+0.00000000, +5.00000000, -8.00000000
1288-6.00000000, -2.00000000, -7.00000000
1289affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
1290affinInv
1291-7.00000000, -4.00000000, +7.00000000
1292+4.00000000, -9.00000000, -7.00000000
1293-6.00000000, -1.00000000, +2.00000000
1294+0.00000000, +5.00000000, -8.00000000
1295-6.00000000, -2.00000000, -7.00000000
1296sample ! for comparison with affinInv.
1297-7.00000000, -4.00000000, +7.00000000
1298+4.00000000, -9.00000000, -7.00000000
1299-6.00000000, -1.00000000, +2.00000000
1300+0.00000000, +5.00000000, -8.00000000
1301-6.00000000, -2.00000000, -7.00000000
1302
1303dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
1304[dim, ndim, nsam]
1305+1, +1, +4
1306tlate = getUnifRand(-5, +5, ndim)
1307tlate
1308-1.00000000
1309sample = getUnifRand(-9, +9, nsam, ndim)
1310sample
1311-3.00000000
1312-5.00000000
1313+9.00000000
1314+4.00000000
1315call setResized(affinity, shape(sample, IK))
1316call setResized(affinInv, shape(sample, IK))
1317tform = getCovRand(1., ndim)
1318tform
1319+1.00000000
1320affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
1321affinity
1322-3.00000000
1323-5.00000000
1324+9.00000000
1325+4.00000000
1326do isam = 1, nsam
1327 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
1328end do
1329affinity
1330-3.00000000
1331-5.00000000
1332+9.00000000
1333+4.00000000
1334affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
1335affinInv
1336-3.00000000
1337-5.00000000
1338+9.00000000
1339+4.00000000
1340sample ! for comparison with affinInv.
1341-3.00000000
1342-5.00000000
1343+9.00000000
1344+4.00000000
1345
1346call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
1347tform
1348+1.00000000
1349affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
1350affinity
1351-3.00000000
1352-5.00000000
1353+9.00000000
1354+4.00000000
1355do isam = 1, nsam
1356 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
1357end do
1358affinity
1359-3.00000000
1360-5.00000000
1361+9.00000000
1362+4.00000000
1363affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
1364affinInv
1365-3.00000000
1366-5.00000000
1367+9.00000000
1368+4.00000000
1369sample ! for comparison with affinInv.
1370-3.00000000
1371-5.00000000
1372+9.00000000
1373+4.00000000
1374
1375affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
1376affinity
1377-3.00000000
1378-5.00000000
1379+9.00000000
1380+4.00000000
1381do isam = 1, nsam
1382 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
1383end do
1384affinity
1385-3.00000000
1386-5.00000000
1387+9.00000000
1388+4.00000000
1389affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
1390affinInv
1391-3.00000000
1392-5.00000000
1393+9.00000000
1394+4.00000000
1395sample ! for comparison with affinInv.
1396-3.00000000
1397-5.00000000
1398+9.00000000
1399+4.00000000
1400
1401affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
1402affinity
1403-3.00000000
1404-5.00000000
1405+9.00000000
1406+4.00000000
1407do isam = 1, nsam
1408 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
1409end do
1410affinity
1411-3.00000000
1412-5.00000000
1413+9.00000000
1414+4.00000000
1415affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
1416affinInv
1417-3.00000000
1418-5.00000000
1419+9.00000000
1420+4.00000000
1421sample ! for comparison with affinInv.
1422-3.00000000
1423-5.00000000
1424+9.00000000
1425+4.00000000
1426
1427affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
1428affinity
1429-3.00000000
1430-5.00000000
1431+9.00000000
1432+4.00000000
1433do isam = 1, nsam
1434 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
1435end do
1436affinity
1437-3.00000000
1438-5.00000000
1439+9.00000000
1440+4.00000000
1441affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
1442affinInv
1443-3.00000000
1444-5.00000000
1445+9.00000000
1446+4.00000000
1447sample ! for comparison with affinInv.
1448-3.00000000
1449-5.00000000
1450+9.00000000
1451+4.00000000
1452
1453dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
1454[dim, ndim, nsam]
1455+1, +2, +4
1456tlate = getUnifRand(-5, +5, ndim)
1457tlate
1458-1.00000000, -4.00000000
1459sample = getUnifRand(-9, +9, nsam, ndim)
1460sample
1461+1.00000000, -1.00000000
1462+5.00000000, +4.00000000
1463+7.00000000, -2.00000000
1464+3.00000000, +9.00000000
1465call setResized(affinity, shape(sample, IK))
1466call setResized(affinInv, shape(sample, IK))
1467tform = getCovRand(1., ndim)
1468tform
1469+1.00000000, +0.244880095E-1
1470+0.244880095E-1, +1.00000000
1471affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
1472affinity
1473+0.975511968, -0.975511968
1474+5.09795189, +4.12243986
1475+6.95102406, -1.82858396
1476+3.22039199, +9.07346439
1477do isam = 1, nsam
1478 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
1479end do
1480affinity
1481+0.975511968, -0.975511968
1482+5.09795189, +4.12243986
1483+6.95102406, -1.82858396
1484+3.22039199, +9.07346439
1485affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
1486affinInv
1487+0.999999940, -1.00000000
1488+4.99999952, +4.00000000
1489+6.99999952, -2.00000024
1490+2.99999976, +9.00000095
1491sample ! for comparison with affinInv.
1492+1.00000000, -1.00000000
1493+5.00000000, +4.00000000
1494+7.00000000, -2.00000000
1495+3.00000000, +9.00000000
1496
1497call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
1498tform
1499+1.00000000, +0.244880095E-1
1500+0.00000000, +1.00000000
1501affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
1502affinity
1503+0.975511968, -1.00000000
1504+5.09795189, +4.00000000
1505+6.95102406, -2.00000000
1506+3.22039199, +9.00000000
1507do isam = 1, nsam
1508 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
1509end do
1510affinity
1511+0.975511968, -1.00000000
1512+5.09795189, +4.00000000
1513+6.95102406, -2.00000000
1514+3.22039199, +9.00000000
1515affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
1516affinInv
1517+1.00000000, -1.00000000
1518+5.00000000, +4.00000000
1519+7.00000000, -2.00000000
1520+3.00000000, +9.00000000
1521sample ! for comparison with affinInv.
1522+1.00000000, -1.00000000
1523+5.00000000, +4.00000000
1524+7.00000000, -2.00000000
1525+3.00000000, +9.00000000
1526
1527affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
1528affinity
1529+0.975511968, -1.00000000
1530+5.09795189, +4.00000000
1531+6.95102406, -2.00000000
1532+3.22039199, +9.00000000
1533do isam = 1, nsam
1534 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
1535end do
1536affinity
1537+0.975511968, -1.00000000
1538+5.09795189, +4.00000000
1539+6.95102406, -2.00000000
1540+3.22039199, +9.00000000
1541affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
1542affinInv
1543+1.00000000, -1.00000000
1544+5.00000000, +4.00000000
1545+7.00000000, -2.00000000
1546+3.00000000, +9.00000000
1547sample ! for comparison with affinInv.
1548+1.00000000, -1.00000000
1549+5.00000000, +4.00000000
1550+7.00000000, -2.00000000
1551+3.00000000, +9.00000000
1552
1553affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
1554affinity
1555+1.00000000, -1.00000000
1556+5.00000000, +4.00000000
1557+7.00000000, -2.00000000
1558+3.00000000, +9.00000000
1559do isam = 1, nsam
1560 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
1561end do
1562affinity
1563+1.00000000, -1.00000000
1564+5.00000000, +4.00000000
1565+7.00000000, -2.00000000
1566+3.00000000, +9.00000000
1567affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
1568affinInv
1569+1.00000000, -1.00000000
1570+5.00000000, +4.00000000
1571+7.00000000, -2.00000000
1572+3.00000000, +9.00000000
1573sample ! for comparison with affinInv.
1574+1.00000000, -1.00000000
1575+5.00000000, +4.00000000
1576+7.00000000, -2.00000000
1577+3.00000000, +9.00000000
1578
1579affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
1580affinity
1581+1.00000000, -1.00000000
1582+5.00000000, +4.00000000
1583+7.00000000, -2.00000000
1584+3.00000000, +9.00000000
1585do isam = 1, nsam
1586 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
1587end do
1588affinity
1589+1.00000000, -1.00000000
1590+5.00000000, +4.00000000
1591+7.00000000, -2.00000000
1592+3.00000000, +9.00000000
1593affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
1594affinInv
1595+1.00000000, -1.00000000
1596+5.00000000, +4.00000000
1597+7.00000000, -2.00000000
1598+3.00000000, +9.00000000
1599sample ! for comparison with affinInv.
1600+1.00000000, -1.00000000
1601+5.00000000, +4.00000000
1602+7.00000000, -2.00000000
1603+3.00000000, +9.00000000
1604

Postprocessing of the example output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6import glob
7import sys
8
9linewidth = 2
10fontsize = 17
11
12for kind in ["warp", "rotation"]:
13
14 pattern = "*."+kind+"*.txt"
15 fileList = glob.glob(pattern)
16
17 for file in fileList:
18
19 df = pd.read_csv(file, delimiter = ",", header = None)
20
21 # definitions for the axes
22 left, width = 0.1, 0.65
23 bottom, height = 0.1, 0.65
24 spacing = 0.015
25
26 # start with a square Figure
27 fig = plt.figure(figsize = (8, 8))
28
29 plt.rcParams.update({'font.size': fontsize - 2})
30 ax = fig.add_axes([left, bottom, width, height]) # scatter plot
31 ax_histx = fig.add_axes([left, bottom + height + spacing, width, 0.2], sharex = ax) # histx
32 ax_histy = fig.add_axes([left + width + spacing, bottom, 0.2, height], sharey = ax) # histy
33
34 for axes in [ax, ax_histx, ax_histy]:
35 axes.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
36 axes.tick_params(axis = "y", which = "minor")
37 axes.tick_params(axis = "x", which = "minor")
38
39 # no labels
40 ax_histy.tick_params(axis = "y", labelleft = False)
41 ax_histx.tick_params(axis = "x", labelbottom = False)
42
43 # the scatter plot:
44 ax.scatter ( df.values[:, 0]
45 , df.values[:, 1]
46 , s = 8
47 , zorder = 1000
48 )
49
50 ax_histx.hist(df.values[:, 0], bins = 50, zorder = 1000)
51 ax_histy.hist(df.values[:, 1], bins = 50, orientation = "horizontal", zorder = 1000)
52
53 ax.set_xlabel("X", fontsize = 17)
54 ax.set_ylabel("Y", fontsize = 17)
55 ax.legend([file.split(".")[-2]], fontsize = fontsize)
56
57 #plt.tight_layout()
58 plt.savefig(file.replace(".txt",".png"))

Visualization of the example output
Test:
test_pm_sampleAffinity


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, Saturday 2:48 AM, August 22, 2021, Dallas, TX

Definition at line 225 of file pm_sampleAffinity.F90.


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