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, +3, +5
10tlate = getUnifRand(-5, +5, ndim)
11tlate
12-2.00000000, +1.00000000, +3.00000000
13sample = getUnifRand(-9, +9, ndim, nsam)
14sample
15-3.00000000, -2.00000000, +0.00000000, +0.00000000, -7.00000000
16-9.00000000, -8.00000000, -4.00000000, +8.00000000, +8.00000000
17+2.00000000, +1.00000000, +8.00000000, +4.00000000, -4.00000000
18call setResized(affinity, shape(sample, IK))
19call setResized(affinInv, shape(sample, IK))
20tform = getCovRand(1., ndim)
21tform
22+1.00000000, +0.225927606E-1, +0.975619927E-1
23+0.225927606E-1, +1.00000000, +0.936427712
24+0.975619927E-1, +0.936427712, +1.00000000
25affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
26transpose(affinity)
27-3.00821090, -7.19492340, -6.72053528
28-2.08317995, -7.10875750, -6.68654585
29+0.690124869, +3.49142170, +4.25428915
30+0.570990086, +11.7457104, +11.4914217
31-7.20950556, +4.09613991, +2.80848789
32do isam = 1, nsam
33 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
34end do
35transpose(affinity)
36-3.00821090, -7.19492340, -6.72053528
37-2.08317995, -7.10875750, -6.68654585
38+0.690124869, +3.49142170, +4.25428915
39+0.570990086, +11.7457104, +11.4914217
40-7.20950556, +4.09613991, +2.80848789
41affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
42transpose(affinInv)
43-3.00000048, -9.00000000, +2.00000381
44-1.99999952, -7.99999619, +0.999996185
45+0.00000000, -4.00000381, +8.00000000
46+0.00000000, +7.99999237, +4.00000763
47-7.00000000, +7.99999809, -4.00000191
48transpose(sample) ! for comparison with affinInv.
49-3.00000000, -9.00000000, +2.00000000
50-2.00000000, -8.00000000, +1.00000000
51+0.00000000, -4.00000000, +8.00000000
52+0.00000000, +8.00000000, +4.00000000
53-7.00000000, +8.00000000, -4.00000000
54
55call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
56tform
57+1.00000000, +0.225927606E-1, +0.975619927E-1
58+0.00000000, +1.00000000, +0.936427712
59+0.00000000, +0.00000000, +1.00000000
60affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
61transpose(affinity)
62-3.00821090, -7.12714481, +2.00000000
63-2.08318019, -7.06357241, +1.00000000
64+0.690124869, +3.49142170, +8.00000000
65+0.570990086, +11.7457104, +4.00000000
66-7.20950603, +4.25428915, -4.00000000
67do isam = 1, nsam
68 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
69end do
70transpose(affinity)
71-3.00821090, -7.12714481, +2.00000000
72-2.08318019, -7.06357241, +1.00000000
73+0.690124869, +3.49142170, +8.00000000
74+0.570990086, +11.7457104, +4.00000000
75-7.20950603, +4.25428915, -4.00000000
76affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
77transpose(affinInv)
78-3.00000000, -9.00000000, +2.00000000
79-2.00000000, -8.00000000, +1.00000000
80+0.00000000, -4.00000000, +8.00000000
81+0.596046448E-7, +7.99999952, +4.00000000
82-7.00000000, +8.00000000, -4.00000000
83transpose(sample) ! for comparison with affinInv.
84-3.00000000, -9.00000000, +2.00000000
85-2.00000000, -8.00000000, +1.00000000
86+0.00000000, -4.00000000, +8.00000000
87+0.00000000, +8.00000000, +4.00000000
88-7.00000000, +8.00000000, -4.00000000
89
90affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
91transpose(affinity)
92-3.00821090, -7.12714481, +2.00000000
93-2.08318019, -7.06357241, +1.00000000
94+0.690124869, +3.49142170, +8.00000000
95+0.570990086, +11.7457104, +4.00000000
96-7.20950603, +4.25428915, -4.00000000
97do isam = 1, nsam
98 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
99end do
100transpose(affinity)
101-3.00821090, -7.12714481, +2.00000000
102-2.08318019, -7.06357241, +1.00000000
103+0.690124869, +3.49142170, +8.00000000
104+0.570990086, +11.7457104, +4.00000000
105-7.20950603, +4.25428915, -4.00000000
106affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
107transpose(affinInv)
108-3.00000000, -9.00000000, +2.00000000
109-2.00000000, -8.00000000, +1.00000000
110+0.00000000, -4.00000000, +8.00000000
111+0.596046448E-7, +7.99999952, +4.00000000
112-7.00000000, +8.00000000, -4.00000000
113transpose(sample) ! for comparison with affinInv.
114-3.00000000, -9.00000000, +2.00000000
115-2.00000000, -8.00000000, +1.00000000
116+0.00000000, -4.00000000, +8.00000000
117+0.00000000, +8.00000000, +4.00000000
118-7.00000000, +8.00000000, -4.00000000
119
120affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
121transpose(affinity)
122-3.00000000, -9.00000000, +2.00000000
123-2.00000000, -8.00000000, +1.00000000
124+0.00000000, -4.00000000, +8.00000000
125+0.00000000, +8.00000000, +4.00000000
126-7.00000000, +8.00000000, -4.00000000
127do isam = 1, nsam
128 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
129end do
130transpose(affinity)
131-3.00000000, -9.00000000, +2.00000000
132-2.00000000, -8.00000000, +1.00000000
133+0.00000000, -4.00000000, +8.00000000
134+0.00000000, +8.00000000, +4.00000000
135-7.00000000, +8.00000000, -4.00000000
136affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
137transpose(affinInv)
138-3.00000000, -9.00000000, +2.00000000
139-2.00000000, -8.00000000, +1.00000000
140+0.00000000, -4.00000000, +8.00000000
141+0.00000000, +8.00000000, +4.00000000
142-7.00000000, +8.00000000, -4.00000000
143transpose(sample) ! for comparison with affinInv.
144-3.00000000, -9.00000000, +2.00000000
145-2.00000000, -8.00000000, +1.00000000
146+0.00000000, -4.00000000, +8.00000000
147+0.00000000, +8.00000000, +4.00000000
148-7.00000000, +8.00000000, -4.00000000
149
150affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
151transpose(affinity)
152-3.00000000, -9.00000000, +2.00000000
153-2.00000000, -8.00000000, +1.00000000
154+0.00000000, -4.00000000, +8.00000000
155+0.00000000, +8.00000000, +4.00000000
156-7.00000000, +8.00000000, -4.00000000
157do isam = 1, nsam
158 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
159end do
160transpose(affinity)
161-3.00000000, -9.00000000, +2.00000000
162-2.00000000, -8.00000000, +1.00000000
163+0.00000000, -4.00000000, +8.00000000
164+0.00000000, +8.00000000, +4.00000000
165-7.00000000, +8.00000000, -4.00000000
166affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
167transpose(affinInv)
168-3.00000000, -9.00000000, +2.00000000
169-2.00000000, -8.00000000, +1.00000000
170+0.00000000, -4.00000000, +8.00000000
171+0.00000000, +8.00000000, +4.00000000
172-7.00000000, +8.00000000, -4.00000000
173transpose(sample) ! for comparison with affinInv.
174-3.00000000, -9.00000000, +2.00000000
175-2.00000000, -8.00000000, +1.00000000
176+0.00000000, -4.00000000, +8.00000000
177+0.00000000, +8.00000000, +4.00000000
178-7.00000000, +8.00000000, -4.00000000
179
180dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
181[dim, ndim, nsam]
182+2, +1, +4
183tlate = getUnifRand(-5, +5, ndim)
184tlate
185+1.00000000
186sample = getUnifRand(-9, +9, ndim, nsam)
187sample
188+4.00000000, -3.00000000, -1.00000000, -9.00000000
189call setResized(affinity, shape(sample, IK))
190call setResized(affinInv, shape(sample, IK))
191tform = getCovRand(1., ndim)
192tform
193+1.00000000
194affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
195transpose(affinity)
196+4.00000000
197-3.00000000
198-1.00000000
199-9.00000000
200do isam = 1, nsam
201 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
202end do
203transpose(affinity)
204+4.00000000
205-3.00000000
206-1.00000000
207-9.00000000
208affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
209transpose(affinInv)
210+4.00000000
211-3.00000000
212-1.00000000
213-9.00000000
214transpose(sample) ! for comparison with affinInv.
215+4.00000000
216-3.00000000
217-1.00000000
218-9.00000000
219
220call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
221tform
222+1.00000000
223affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
224transpose(affinity)
225+4.00000000
226-3.00000000
227-1.00000000
228-9.00000000
229do isam = 1, nsam
230 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
231end do
232transpose(affinity)
233+4.00000000
234-3.00000000
235-1.00000000
236-9.00000000
237affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
238transpose(affinInv)
239+4.00000000
240-3.00000000
241-1.00000000
242-9.00000000
243transpose(sample) ! for comparison with affinInv.
244+4.00000000
245-3.00000000
246-1.00000000
247-9.00000000
248
249affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
250transpose(affinity)
251+4.00000000
252-3.00000000
253-1.00000000
254-9.00000000
255do isam = 1, nsam
256 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
257end do
258transpose(affinity)
259+4.00000000
260-3.00000000
261-1.00000000
262-9.00000000
263affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
264transpose(affinInv)
265+4.00000000
266-3.00000000
267-1.00000000
268-9.00000000
269transpose(sample) ! for comparison with affinInv.
270+4.00000000
271-3.00000000
272-1.00000000
273-9.00000000
274
275affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
276transpose(affinity)
277+4.00000000
278-3.00000000
279-1.00000000
280-9.00000000
281do isam = 1, nsam
282 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
283end do
284transpose(affinity)
285+4.00000000
286-3.00000000
287-1.00000000
288-9.00000000
289affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
290transpose(affinInv)
291+4.00000000
292-3.00000000
293-1.00000000
294-9.00000000
295transpose(sample) ! for comparison with affinInv.
296+4.00000000
297-3.00000000
298-1.00000000
299-9.00000000
300
301affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
302transpose(affinity)
303+4.00000000
304-3.00000000
305-1.00000000
306-9.00000000
307do isam = 1, nsam
308 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
309end do
310transpose(affinity)
311+4.00000000
312-3.00000000
313-1.00000000
314-9.00000000
315affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
316transpose(affinInv)
317+4.00000000
318-3.00000000
319-1.00000000
320-9.00000000
321transpose(sample) ! for comparison with affinInv.
322+4.00000000
323-3.00000000
324-1.00000000
325-9.00000000
326
327dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
328[dim, ndim, nsam]
329+2, +3, +5
330tlate = getUnifRand(-5, +5, ndim)
331tlate
332+3.00000000, +5.00000000, +0.00000000
333sample = getUnifRand(-9, +9, ndim, nsam)
334sample
335+2.00000000, +7.00000000, +7.00000000, +6.00000000, -1.00000000
336-6.00000000, -3.00000000, -4.00000000, +8.00000000, -5.00000000
337+5.00000000, -8.00000000, +7.00000000, -9.00000000, -4.00000000
338call setResized(affinity, shape(sample, IK))
339call setResized(affinInv, shape(sample, IK))
340tform = getCovRand(1., ndim)
341tform
342+1.00000000, +0.845762074, -0.792924222E-2
343+0.845762074, +1.00000000, -0.386259466
344-0.792924222E-2, -0.386259466, +1.00000000
345affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
346transpose(affinity)
347-3.11421871, -6.23977327, +7.30169868
348+4.52614784, +6.01041031, -6.89672613
349+3.56144691, -0.783481836, +8.48953342
350+12.8374596, +16.5509071, -12.1376514
351-5.19709349, -4.30072451, -2.06077337
352do isam = 1, nsam
353 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
354end do
355transpose(affinity)
356-3.11421871, -6.23977327, +7.30169868
357+4.52614784, +6.01041031, -6.89672613
358+3.56144691, -0.783481836, +8.48953342
359+12.8374596, +16.5509071, -12.1376514
360-5.19709349, -4.30072451, -2.06077337
361affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
362transpose(affinInv)
363+1.99999809, -5.99999619, +5.00000095
364+7.00000286, -3.00000191, -8.00000000
365+7.00000000, -4.00000000, +7.00000000
366+6.00000572, +7.99999237, -9.00000381
367-1.00000000, -5.00000095, -4.00000048
368transpose(sample) ! for comparison with affinInv.
369+2.00000000, -6.00000000, +5.00000000
370+7.00000000, -3.00000000, -8.00000000
371+7.00000000, -4.00000000, +7.00000000
372+6.00000000, +8.00000000, -9.00000000
373-1.00000000, -5.00000000, -4.00000000
374
375call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
376tform
377+1.00000000, +0.845762074, -0.792924222E-2
378+0.00000000, +1.00000000, -0.386259466
379+0.00000000, +0.00000000, +1.00000000
380affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
381transpose(affinity)
382-3.11421871, -7.93129730, +5.00000000
383+4.52614784, +0.900757313E-1, -8.00000000
384+3.56144691, -6.70381641, +7.00000000
385+12.8374596, +11.4763355, -9.00000000
386-5.19709349, -3.45496225, -4.00000000
387do isam = 1, nsam
388 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
389end do
390transpose(affinity)
391-3.11421871, -7.93129730, +5.00000000
392+4.52614784, +0.900757313E-1, -8.00000000
393+3.56144691, -6.70381641, +7.00000000
394+12.8374596, +11.4763355, -9.00000000
395-5.19709349, -3.45496225, -4.00000000
396affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
397transpose(affinInv)
398+2.00000000, -6.00000000, +5.00000000
399+7.00000000, -3.00000000, -8.00000000
400+7.00000000, -4.00000000, +7.00000000
401+6.00000000, +8.00000000, -9.00000000
402-1.00000000, -5.00000000, -4.00000000
403transpose(sample) ! for comparison with affinInv.
404+2.00000000, -6.00000000, +5.00000000
405+7.00000000, -3.00000000, -8.00000000
406+7.00000000, -4.00000000, +7.00000000
407+6.00000000, +8.00000000, -9.00000000
408-1.00000000, -5.00000000, -4.00000000
409
410affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
411transpose(affinity)
412-3.11421871, -7.93129730, +5.00000000
413+4.52614784, +0.900757313E-1, -8.00000000
414+3.56144691, -6.70381641, +7.00000000
415+12.8374596, +11.4763355, -9.00000000
416-5.19709349, -3.45496225, -4.00000000
417do isam = 1, nsam
418 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
419end do
420transpose(affinity)
421-3.11421871, -7.93129730, +5.00000000
422+4.52614784, +0.900757313E-1, -8.00000000
423+3.56144691, -6.70381641, +7.00000000
424+12.8374596, +11.4763355, -9.00000000
425-5.19709349, -3.45496225, -4.00000000
426affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
427transpose(affinInv)
428+2.00000000, -6.00000000, +5.00000000
429+7.00000000, -3.00000000, -8.00000000
430+7.00000000, -4.00000000, +7.00000000
431+6.00000000, +8.00000000, -9.00000000
432-1.00000000, -5.00000000, -4.00000000
433transpose(sample) ! for comparison with affinInv.
434+2.00000000, -6.00000000, +5.00000000
435+7.00000000, -3.00000000, -8.00000000
436+7.00000000, -4.00000000, +7.00000000
437+6.00000000, +8.00000000, -9.00000000
438-1.00000000, -5.00000000, -4.00000000
439
440affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
441transpose(affinity)
442+2.00000000, -6.00000000, +5.00000000
443+7.00000000, -3.00000000, -8.00000000
444+7.00000000, -4.00000000, +7.00000000
445+6.00000000, +8.00000000, -9.00000000
446-1.00000000, -5.00000000, -4.00000000
447do isam = 1, nsam
448 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
449end do
450transpose(affinity)
451+2.00000000, -6.00000000, +5.00000000
452+7.00000000, -3.00000000, -8.00000000
453+7.00000000, -4.00000000, +7.00000000
454+6.00000000, +8.00000000, -9.00000000
455-1.00000000, -5.00000000, -4.00000000
456affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
457transpose(affinInv)
458+2.00000000, -6.00000000, +5.00000000
459+7.00000000, -3.00000000, -8.00000000
460+7.00000000, -4.00000000, +7.00000000
461+6.00000000, +8.00000000, -9.00000000
462-1.00000000, -5.00000000, -4.00000000
463transpose(sample) ! for comparison with affinInv.
464+2.00000000, -6.00000000, +5.00000000
465+7.00000000, -3.00000000, -8.00000000
466+7.00000000, -4.00000000, +7.00000000
467+6.00000000, +8.00000000, -9.00000000
468-1.00000000, -5.00000000, -4.00000000
469
470affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
471transpose(affinity)
472+2.00000000, -6.00000000, +5.00000000
473+7.00000000, -3.00000000, -8.00000000
474+7.00000000, -4.00000000, +7.00000000
475+6.00000000, +8.00000000, -9.00000000
476-1.00000000, -5.00000000, -4.00000000
477do isam = 1, nsam
478 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
479end do
480transpose(affinity)
481+2.00000000, -6.00000000, +5.00000000
482+7.00000000, -3.00000000, -8.00000000
483+7.00000000, -4.00000000, +7.00000000
484+6.00000000, +8.00000000, -9.00000000
485-1.00000000, -5.00000000, -4.00000000
486affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
487transpose(affinInv)
488+2.00000000, -6.00000000, +5.00000000
489+7.00000000, -3.00000000, -8.00000000
490+7.00000000, -4.00000000, +7.00000000
491+6.00000000, +8.00000000, -9.00000000
492-1.00000000, -5.00000000, -4.00000000
493transpose(sample) ! for comparison with affinInv.
494+2.00000000, -6.00000000, +5.00000000
495+7.00000000, -3.00000000, -8.00000000
496+7.00000000, -4.00000000, +7.00000000
497+6.00000000, +8.00000000, -9.00000000
498-1.00000000, -5.00000000, -4.00000000
499
500dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
501[dim, ndim, nsam]
502+2, +3, +4
503tlate = getUnifRand(-5, +5, ndim)
504tlate
505+4.00000000, +2.00000000, +0.00000000
506sample = getUnifRand(-9, +9, ndim, nsam)
507sample
508+0.00000000, +4.00000000, +8.00000000, +5.00000000
509+0.00000000, -1.00000000, -2.00000000, -5.00000000
510-9.00000000, -6.00000000, +2.00000000, +7.00000000
511call setResized(affinity, shape(sample, IK))
512call setResized(affinInv, shape(sample, IK))
513tform = getCovRand(1., ndim)
514tform
515+1.00000000, +0.966987073, +0.293903768
516+0.966987073, +1.00000000, +0.118419081
517+0.293903768, +0.118419081, +1.00000000
518affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
519transpose(affinity)
520-2.64513397, -1.06577170, -9.00000000
521+1.26959026, +2.15743375, -4.94280386
522+6.65383339, +5.97273493, +4.11439228
523+2.22239113, +0.663868785, +7.87742329
524do isam = 1, nsam
525 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
526end do
527transpose(affinity)
528-2.64513397, -1.06577170, -9.00000000
529+1.26959026, +2.15743375, -4.94280386
530+6.65383339, +5.97273493, +4.11439228
531+2.22239113, +0.663868785, +7.87742329
532affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
533transpose(affinInv)
534+0.00000000, -0.762939453E-5, -9.00000000
535+3.99999237, -0.999992371, -5.99999952
536+7.99995613, -1.99997330, +2.00000381
537+5.00000381, -5.00000381, +6.99999905
538transpose(sample) ! for comparison with affinInv.
539+0.00000000, +0.00000000, -9.00000000
540+4.00000000, -1.00000000, -6.00000000
541+8.00000000, -2.00000000, +2.00000000
542+5.00000000, -5.00000000, +7.00000000
543
544call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
545tform
546+1.00000000, +0.966987073, +0.293903768
547+0.00000000, +1.00000000, +0.118419081
548+0.00000000, +0.00000000, +1.00000000
549affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
550transpose(affinity)
551-2.64513397, -1.06577170, -9.00000000
552+1.26959038, -1.71051455, -6.00000000
553+6.65383339, -1.76316190, +2.00000000
554+2.22239113, -4.17106628, +7.00000000
555do isam = 1, nsam
556 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
557end do
558transpose(affinity)
559-2.64513397, -1.06577170, -9.00000000
560+1.26959038, -1.71051455, -6.00000000
561+6.65383339, -1.76316190, +2.00000000
562+2.22239113, -4.17106628, +7.00000000
563affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
564transpose(affinInv)
565-0.238418579E-6, +0.00000000, -9.00000000
566+4.00000000, -1.00000000, -6.00000000
567+8.00000000, -2.00000000, +2.00000000
568+5.00000000, -5.00000000, +7.00000000
569transpose(sample) ! for comparison with affinInv.
570+0.00000000, +0.00000000, -9.00000000
571+4.00000000, -1.00000000, -6.00000000
572+8.00000000, -2.00000000, +2.00000000
573+5.00000000, -5.00000000, +7.00000000
574
575affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
576transpose(affinity)
577-2.64513397, -1.06577170, -9.00000000
578+1.26959038, -1.71051455, -6.00000000
579+6.65383339, -1.76316190, +2.00000000
580+2.22239113, -4.17106628, +7.00000000
581do isam = 1, nsam
582 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
583end do
584transpose(affinity)
585-2.64513397, -1.06577170, -9.00000000
586+1.26959038, -1.71051455, -6.00000000
587+6.65383339, -1.76316190, +2.00000000
588+2.22239113, -4.17106628, +7.00000000
589affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
590transpose(affinInv)
591-0.238418579E-6, +0.00000000, -9.00000000
592+4.00000000, -1.00000000, -6.00000000
593+8.00000000, -2.00000000, +2.00000000
594+5.00000000, -5.00000000, +7.00000000
595transpose(sample) ! for comparison with affinInv.
596+0.00000000, +0.00000000, -9.00000000
597+4.00000000, -1.00000000, -6.00000000
598+8.00000000, -2.00000000, +2.00000000
599+5.00000000, -5.00000000, +7.00000000
600
601affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
602transpose(affinity)
603+0.00000000, +0.00000000, -9.00000000
604+4.00000000, -1.00000000, -6.00000000
605+8.00000000, -2.00000000, +2.00000000
606+5.00000000, -5.00000000, +7.00000000
607do isam = 1, nsam
608 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
609end do
610transpose(affinity)
611+0.00000000, +0.00000000, -9.00000000
612+4.00000000, -1.00000000, -6.00000000
613+8.00000000, -2.00000000, +2.00000000
614+5.00000000, -5.00000000, +7.00000000
615affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
616transpose(affinInv)
617+0.00000000, +0.00000000, -9.00000000
618+4.00000000, -1.00000000, -6.00000000
619+8.00000000, -2.00000000, +2.00000000
620+5.00000000, -5.00000000, +7.00000000
621transpose(sample) ! for comparison with affinInv.
622+0.00000000, +0.00000000, -9.00000000
623+4.00000000, -1.00000000, -6.00000000
624+8.00000000, -2.00000000, +2.00000000
625+5.00000000, -5.00000000, +7.00000000
626
627affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
628transpose(affinity)
629+0.00000000, +0.00000000, -9.00000000
630+4.00000000, -1.00000000, -6.00000000
631+8.00000000, -2.00000000, +2.00000000
632+5.00000000, -5.00000000, +7.00000000
633do isam = 1, nsam
634 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
635end do
636transpose(affinity)
637+0.00000000, +0.00000000, -9.00000000
638+4.00000000, -1.00000000, -6.00000000
639+8.00000000, -2.00000000, +2.00000000
640+5.00000000, -5.00000000, +7.00000000
641affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
642transpose(affinInv)
643+0.00000000, +0.00000000, -9.00000000
644+4.00000000, -1.00000000, -6.00000000
645+8.00000000, -2.00000000, +2.00000000
646+5.00000000, -5.00000000, +7.00000000
647transpose(sample) ! for comparison with affinInv.
648+0.00000000, +0.00000000, -9.00000000
649+4.00000000, -1.00000000, -6.00000000
650+8.00000000, -2.00000000, +2.00000000
651+5.00000000, -5.00000000, +7.00000000
652
653dim = 2; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
654[dim, ndim, nsam]
655+2, +1, +5
656tlate = getUnifRand(-5, +5, ndim)
657tlate
658+1.00000000
659sample = getUnifRand(-9, +9, ndim, nsam)
660sample
661+7.00000000, +2.00000000, +8.00000000, -2.00000000, +6.00000000
662call setResized(affinity, shape(sample, IK))
663call setResized(affinInv, shape(sample, IK))
664tform = getCovRand(1., ndim)
665tform
666+1.00000000
667affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
668transpose(affinity)
669+7.00000000
670+2.00000000
671+8.00000000
672-2.00000000
673+6.00000000
674do isam = 1, nsam
675 affinity(:,isam) = getAffinity(sample(:,isam), tform) ! point-wise transformation
676end do
677transpose(affinity)
678+7.00000000
679+2.00000000
680+8.00000000
681-2.00000000
682+6.00000000
683affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
684transpose(affinInv)
685+7.00000000
686+2.00000000
687+8.00000000
688-2.00000000
689+6.00000000
690transpose(sample) ! for comparison with affinInv.
691+7.00000000
692+2.00000000
693+8.00000000
694-2.00000000
695+6.00000000
696
697call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
698tform
699+1.00000000
700affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
701transpose(affinity)
702+7.00000000
703+2.00000000
704+8.00000000
705-2.00000000
706+6.00000000
707do isam = 1, nsam
708 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperDiag) ! point-wise transformation
709end do
710transpose(affinity)
711+7.00000000
712+2.00000000
713+8.00000000
714-2.00000000
715+6.00000000
716affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
717transpose(affinInv)
718+7.00000000
719+2.00000000
720+8.00000000
721-2.00000000
722+6.00000000
723transpose(sample) ! for comparison with affinInv.
724+7.00000000
725+2.00000000
726+8.00000000
727-2.00000000
728+6.00000000
729
730affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
731transpose(affinity)
732+7.00000000
733+2.00000000
734+8.00000000
735-2.00000000
736+6.00000000
737do isam = 1, nsam
738 affinity(:,isam) = getAffinity(sample(:,isam), tform, upperUnit) ! point-wise transformation
739end do
740transpose(affinity)
741+7.00000000
742+2.00000000
743+8.00000000
744-2.00000000
745+6.00000000
746affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
747transpose(affinInv)
748+7.00000000
749+2.00000000
750+8.00000000
751-2.00000000
752+6.00000000
753transpose(sample) ! for comparison with affinInv.
754+7.00000000
755+2.00000000
756+8.00000000
757-2.00000000
758+6.00000000
759
760affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
761transpose(affinity)
762+7.00000000
763+2.00000000
764+8.00000000
765-2.00000000
766+6.00000000
767do isam = 1, nsam
768 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerDiag) ! point-wise transformation
769end do
770transpose(affinity)
771+7.00000000
772+2.00000000
773+8.00000000
774-2.00000000
775+6.00000000
776affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
777transpose(affinInv)
778+7.00000000
779+2.00000000
780+8.00000000
781-2.00000000
782+6.00000000
783transpose(sample) ! for comparison with affinInv.
784+7.00000000
785+2.00000000
786+8.00000000
787-2.00000000
788+6.00000000
789
790affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
791transpose(affinity)
792+7.00000000
793+2.00000000
794+8.00000000
795-2.00000000
796+6.00000000
797do isam = 1, nsam
798 affinity(:,isam) = getAffinity(sample(:,isam), tform, lowerUnit) ! point-wise transformation
799end do
800transpose(affinity)
801+7.00000000
802+2.00000000
803+8.00000000
804-2.00000000
805+6.00000000
806affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
807transpose(affinInv)
808+7.00000000
809+2.00000000
810+8.00000000
811-2.00000000
812+6.00000000
813transpose(sample) ! for comparison with affinInv.
814+7.00000000
815+2.00000000
816+8.00000000
817-2.00000000
818+6.00000000
819
820!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
821!Transform 2D sample along the first dimension.
822!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
823
824
825dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
826[dim, ndim, nsam]
827+1, +1, +5
828tlate = getUnifRand(-5, +5, ndim)
829tlate
830-3.00000000
831sample = getUnifRand(-9, +9, nsam, ndim)
832sample
833+0.00000000
834+7.00000000
835+3.00000000
836-6.00000000
837-1.00000000
838call setResized(affinity, shape(sample, IK))
839call setResized(affinInv, shape(sample, IK))
840tform = getCovRand(1., ndim)
841tform
842+1.00000000
843affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
844affinity
845+0.00000000
846+7.00000000
847+3.00000000
848-6.00000000
849-1.00000000
850do isam = 1, nsam
851 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
852end do
853affinity
854+0.00000000
855+7.00000000
856+3.00000000
857-6.00000000
858-1.00000000
859affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
860affinInv
861+0.00000000
862+7.00000000
863+3.00000000
864-6.00000000
865-1.00000000
866sample ! for comparison with affinInv.
867+0.00000000
868+7.00000000
869+3.00000000
870-6.00000000
871-1.00000000
872
873call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
874tform
875+1.00000000
876affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
877affinity
878+0.00000000
879+7.00000000
880+3.00000000
881-6.00000000
882-1.00000000
883do isam = 1, nsam
884 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
885end do
886affinity
887+0.00000000
888+7.00000000
889+3.00000000
890-6.00000000
891-1.00000000
892affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
893affinInv
894+0.00000000
895+7.00000000
896+3.00000000
897-6.00000000
898-1.00000000
899sample ! for comparison with affinInv.
900+0.00000000
901+7.00000000
902+3.00000000
903-6.00000000
904-1.00000000
905
906affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
907affinity
908+0.00000000
909+7.00000000
910+3.00000000
911-6.00000000
912-1.00000000
913do isam = 1, nsam
914 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
915end do
916affinity
917+0.00000000
918+7.00000000
919+3.00000000
920-6.00000000
921-1.00000000
922affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
923affinInv
924+0.00000000
925+7.00000000
926+3.00000000
927-6.00000000
928-1.00000000
929sample ! for comparison with affinInv.
930+0.00000000
931+7.00000000
932+3.00000000
933-6.00000000
934-1.00000000
935
936affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
937affinity
938+0.00000000
939+7.00000000
940+3.00000000
941-6.00000000
942-1.00000000
943do isam = 1, nsam
944 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
945end do
946affinity
947+0.00000000
948+7.00000000
949+3.00000000
950-6.00000000
951-1.00000000
952affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
953affinInv
954+0.00000000
955+7.00000000
956+3.00000000
957-6.00000000
958-1.00000000
959sample ! for comparison with affinInv.
960+0.00000000
961+7.00000000
962+3.00000000
963-6.00000000
964-1.00000000
965
966affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
967affinity
968+0.00000000
969+7.00000000
970+3.00000000
971-6.00000000
972-1.00000000
973do isam = 1, nsam
974 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
975end do
976affinity
977+0.00000000
978+7.00000000
979+3.00000000
980-6.00000000
981-1.00000000
982affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
983affinInv
984+0.00000000
985+7.00000000
986+3.00000000
987-6.00000000
988-1.00000000
989sample ! for comparison with affinInv.
990+0.00000000
991+7.00000000
992+3.00000000
993-6.00000000
994-1.00000000
995
996dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
997[dim, ndim, nsam]
998+1, +1, +5
999tlate = getUnifRand(-5, +5, ndim)
1000tlate
1001+2.00000000
1002sample = getUnifRand(-9, +9, nsam, ndim)
1003sample
1004+4.00000000
1005+9.00000000
1006+3.00000000
1007+2.00000000
1008+6.00000000
1009call setResized(affinity, shape(sample, IK))
1010call setResized(affinInv, shape(sample, IK))
1011tform = getCovRand(1., ndim)
1012tform
1013+1.00000000
1014affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
1015affinity
1016+4.00000000
1017+9.00000000
1018+3.00000000
1019+2.00000000
1020+6.00000000
1021do isam = 1, nsam
1022 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
1023end do
1024affinity
1025+4.00000000
1026+9.00000000
1027+3.00000000
1028+2.00000000
1029+6.00000000
1030affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
1031affinInv
1032+4.00000000
1033+9.00000000
1034+3.00000000
1035+2.00000000
1036+6.00000000
1037sample ! for comparison with affinInv.
1038+4.00000000
1039+9.00000000
1040+3.00000000
1041+2.00000000
1042+6.00000000
1043
1044call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
1045tform
1046+1.00000000
1047affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
1048affinity
1049+4.00000000
1050+9.00000000
1051+3.00000000
1052+2.00000000
1053+6.00000000
1054do isam = 1, nsam
1055 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
1056end do
1057affinity
1058+4.00000000
1059+9.00000000
1060+3.00000000
1061+2.00000000
1062+6.00000000
1063affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
1064affinInv
1065+4.00000000
1066+9.00000000
1067+3.00000000
1068+2.00000000
1069+6.00000000
1070sample ! for comparison with affinInv.
1071+4.00000000
1072+9.00000000
1073+3.00000000
1074+2.00000000
1075+6.00000000
1076
1077affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
1078affinity
1079+4.00000000
1080+9.00000000
1081+3.00000000
1082+2.00000000
1083+6.00000000
1084do isam = 1, nsam
1085 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
1086end do
1087affinity
1088+4.00000000
1089+9.00000000
1090+3.00000000
1091+2.00000000
1092+6.00000000
1093affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
1094affinInv
1095+4.00000000
1096+9.00000000
1097+3.00000000
1098+2.00000000
1099+6.00000000
1100sample ! for comparison with affinInv.
1101+4.00000000
1102+9.00000000
1103+3.00000000
1104+2.00000000
1105+6.00000000
1106
1107affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
1108affinity
1109+4.00000000
1110+9.00000000
1111+3.00000000
1112+2.00000000
1113+6.00000000
1114do isam = 1, nsam
1115 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
1116end do
1117affinity
1118+4.00000000
1119+9.00000000
1120+3.00000000
1121+2.00000000
1122+6.00000000
1123affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
1124affinInv
1125+4.00000000
1126+9.00000000
1127+3.00000000
1128+2.00000000
1129+6.00000000
1130sample ! for comparison with affinInv.
1131+4.00000000
1132+9.00000000
1133+3.00000000
1134+2.00000000
1135+6.00000000
1136
1137affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
1138affinity
1139+4.00000000
1140+9.00000000
1141+3.00000000
1142+2.00000000
1143+6.00000000
1144do isam = 1, nsam
1145 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
1146end do
1147affinity
1148+4.00000000
1149+9.00000000
1150+3.00000000
1151+2.00000000
1152+6.00000000
1153affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
1154affinInv
1155+4.00000000
1156+9.00000000
1157+3.00000000
1158+2.00000000
1159+6.00000000
1160sample ! for comparison with affinInv.
1161+4.00000000
1162+9.00000000
1163+3.00000000
1164+2.00000000
1165+6.00000000
1166
1167dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
1168[dim, ndim, nsam]
1169+1, +3, +4
1170tlate = getUnifRand(-5, +5, ndim)
1171tlate
1172+2.00000000, +3.00000000, +0.00000000
1173sample = getUnifRand(-9, +9, nsam, ndim)
1174sample
1175-5.00000000, +2.00000000, -3.00000000
1176+1.00000000, +8.00000000, -7.00000000
1177+9.00000000, -8.00000000, +9.00000000
1178+2.00000000, +7.00000000, -7.00000000
1179call setResized(affinity, shape(sample, IK))
1180call setResized(affinInv, shape(sample, IK))
1181tform = getCovRand(1., ndim)
1182tform
1183+1.00000000, -0.815945625, -0.683251798
1184-0.815945625, +1.00000000, +0.962275386
1185-0.683251798, +0.962275386, +1.00000000
1186affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
1187affinity
1188-4.58213615, +3.19290209, +2.34080982
1189-0.744802475, +0.448126793, +0.149512291E-1
1190+9.37829876, -6.68303204, -4.84746933
1191+1.07114315, -1.36781883, -1.63057613
1192do isam = 1, nsam
1193 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
1194end do
1195affinity
1196-4.58213615, +3.19290209, +2.34080982
1197-0.744802475, +0.448126793, +0.149512291E-1
1198+9.37829876, -6.68303204, -4.84746933
1199+1.07114315, -1.36781883, -1.63057613
1200affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
1201affinInv
1202-5.00000763, +1.99998474, -3.00000000
1203+0.999997795, +7.99999571, -6.99999666
1204+9.00001907, -7.99996948, +8.99999237
1205+2.00000286, +7.00000000, -7.00000191
1206sample ! for comparison with affinInv.
1207-5.00000000, +2.00000000, -3.00000000
1208+1.00000000, +8.00000000, -7.00000000
1209+9.00000000, -8.00000000, +9.00000000
1210+2.00000000, +7.00000000, -7.00000000
1211
1212call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
1213tform
1214+1.00000000, -0.815945625, -0.683251798
1215+0.00000000, +1.00000000, +0.962275386
1216+0.00000000, +0.00000000, +1.00000000
1217affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
1218affinity
1219-4.58213615, -0.886826038, -3.00000000
1220-0.744802475, +1.26407242, -7.00000000
1221+9.37829876, +0.660478592, +9.00000000
1222+1.07114315, +0.264072418, -7.00000000
1223do isam = 1, nsam
1224 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
1225end do
1226affinity
1227-4.58213615, -0.886826038, -3.00000000
1228-0.744802475, +1.26407242, -7.00000000
1229+9.37829876, +0.660478592, +9.00000000
1230+1.07114315, +0.264072418, -7.00000000
1231affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
1232affinInv
1233-5.00000000, +2.00000000, -3.00000000
1234+1.00000024, +8.00000000, -7.00000000
1235+9.00000000, -8.00000000, +9.00000000
1236+2.00000024, +7.00000000, -7.00000000
1237sample ! for comparison with affinInv.
1238-5.00000000, +2.00000000, -3.00000000
1239+1.00000000, +8.00000000, -7.00000000
1240+9.00000000, -8.00000000, +9.00000000
1241+2.00000000, +7.00000000, -7.00000000
1242
1243affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
1244affinity
1245-4.58213615, -0.886826038, -3.00000000
1246-0.744802475, +1.26407242, -7.00000000
1247+9.37829876, +0.660478592, +9.00000000
1248+1.07114315, +0.264072418, -7.00000000
1249do isam = 1, nsam
1250 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
1251end do
1252affinity
1253-4.58213615, -0.886826038, -3.00000000
1254-0.744802475, +1.26407242, -7.00000000
1255+9.37829876, +0.660478592, +9.00000000
1256+1.07114315, +0.264072418, -7.00000000
1257affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
1258affinInv
1259-5.00000000, +2.00000000, -3.00000000
1260+1.00000024, +8.00000000, -7.00000000
1261+9.00000000, -8.00000000, +9.00000000
1262+2.00000024, +7.00000000, -7.00000000
1263sample ! for comparison with affinInv.
1264-5.00000000, +2.00000000, -3.00000000
1265+1.00000000, +8.00000000, -7.00000000
1266+9.00000000, -8.00000000, +9.00000000
1267+2.00000000, +7.00000000, -7.00000000
1268
1269affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
1270affinity
1271-5.00000000, +2.00000000, -3.00000000
1272+1.00000000, +8.00000000, -7.00000000
1273+9.00000000, -8.00000000, +9.00000000
1274+2.00000000, +7.00000000, -7.00000000
1275do isam = 1, nsam
1276 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
1277end do
1278affinity
1279-5.00000000, +2.00000000, -3.00000000
1280+1.00000000, +8.00000000, -7.00000000
1281+9.00000000, -8.00000000, +9.00000000
1282+2.00000000, +7.00000000, -7.00000000
1283affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
1284affinInv
1285-5.00000000, +2.00000000, -3.00000000
1286+1.00000000, +8.00000000, -7.00000000
1287+9.00000000, -8.00000000, +9.00000000
1288+2.00000000, +7.00000000, -7.00000000
1289sample ! for comparison with affinInv.
1290-5.00000000, +2.00000000, -3.00000000
1291+1.00000000, +8.00000000, -7.00000000
1292+9.00000000, -8.00000000, +9.00000000
1293+2.00000000, +7.00000000, -7.00000000
1294
1295affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
1296affinity
1297-5.00000000, +2.00000000, -3.00000000
1298+1.00000000, +8.00000000, -7.00000000
1299+9.00000000, -8.00000000, +9.00000000
1300+2.00000000, +7.00000000, -7.00000000
1301do isam = 1, nsam
1302 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
1303end do
1304affinity
1305-5.00000000, +2.00000000, -3.00000000
1306+1.00000000, +8.00000000, -7.00000000
1307+9.00000000, -8.00000000, +9.00000000
1308+2.00000000, +7.00000000, -7.00000000
1309affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
1310affinInv
1311-5.00000000, +2.00000000, -3.00000000
1312+1.00000000, +8.00000000, -7.00000000
1313+9.00000000, -8.00000000, +9.00000000
1314+2.00000000, +7.00000000, -7.00000000
1315sample ! for comparison with affinInv.
1316-5.00000000, +2.00000000, -3.00000000
1317+1.00000000, +8.00000000, -7.00000000
1318+9.00000000, -8.00000000, +9.00000000
1319+2.00000000, +7.00000000, -7.00000000
1320
1321dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
1322[dim, ndim, nsam]
1323+1, +2, +5
1324tlate = getUnifRand(-5, +5, ndim)
1325tlate
1326+2.00000000, +2.00000000
1327sample = getUnifRand(-9, +9, nsam, ndim)
1328sample
1329+5.00000000, -9.00000000
1330+0.00000000, +9.00000000
1331+6.00000000, -7.00000000
1332+1.00000000, -2.00000000
1333+8.00000000, +1.00000000
1334call setResized(affinity, shape(sample, IK))
1335call setResized(affinInv, shape(sample, IK))
1336tform = getCovRand(1., ndim)
1337tform
1338+1.00000000, +0.999979317
1339+0.999979317, +1.00000000
1340affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
1341affinity
1342-3.99981403, -4.00010347
1343+8.99981403, +9.00000000
1344-0.999855042, -1.00012398
1345-0.999958634, -1.00002074
1346+8.99997902, +8.99983406
1347do isam = 1, nsam
1348 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
1349end do
1350affinity
1351-3.99981403, -4.00010347
1352+8.99981403, +9.00000000
1353-0.999855042, -1.00012398
1354-0.999958634, -1.00002074
1355+8.99997902, +8.99983406
1356affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
1357affinInv
1358+4.99218750, -8.99218750
1359+0.156250000E-1, +8.98437500
1360+6.00000000, -7.00195312
1361+1.00195312, -2.00195312
1362+8.00000000, +1.00000000
1363sample ! for comparison with affinInv.
1364+5.00000000, -9.00000000
1365+0.00000000, +9.00000000
1366+6.00000000, -7.00000000
1367+1.00000000, -2.00000000
1368+8.00000000, +1.00000000
1369
1370call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
1371tform
1372+1.00000000, +0.999979317
1373+0.00000000, +1.00000000
1374affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
1375affinity
1376-3.99981403, -9.00000000
1377+8.99981403, +9.00000000
1378-0.999855042, -7.00000000
1379-0.999958634, -2.00000000
1380+8.99997902, +1.00000000
1381do isam = 1, nsam
1382 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
1383end do
1384affinity
1385-3.99981403, -9.00000000
1386+8.99981403, +9.00000000
1387-0.999855042, -7.00000000
1388-0.999958634, -2.00000000
1389+8.99997902, +1.00000000
1390affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
1391affinInv
1392+5.00000000, -9.00000000
1393+0.00000000, +9.00000000
1394+6.00000000, -7.00000000
1395+1.00000000, -2.00000000
1396+7.99999952, +1.00000000
1397sample ! for comparison with affinInv.
1398+5.00000000, -9.00000000
1399+0.00000000, +9.00000000
1400+6.00000000, -7.00000000
1401+1.00000000, -2.00000000
1402+8.00000000, +1.00000000
1403
1404affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
1405affinity
1406-3.99981403, -9.00000000
1407+8.99981403, +9.00000000
1408-0.999855042, -7.00000000
1409-0.999958634, -2.00000000
1410+8.99997902, +1.00000000
1411do isam = 1, nsam
1412 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
1413end do
1414affinity
1415-3.99981403, -9.00000000
1416+8.99981403, +9.00000000
1417-0.999855042, -7.00000000
1418-0.999958634, -2.00000000
1419+8.99997902, +1.00000000
1420affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
1421affinInv
1422+5.00000000, -9.00000000
1423+0.00000000, +9.00000000
1424+6.00000000, -7.00000000
1425+1.00000000, -2.00000000
1426+7.99999952, +1.00000000
1427sample ! for comparison with affinInv.
1428+5.00000000, -9.00000000
1429+0.00000000, +9.00000000
1430+6.00000000, -7.00000000
1431+1.00000000, -2.00000000
1432+8.00000000, +1.00000000
1433
1434affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
1435affinity
1436+5.00000000, -9.00000000
1437+0.00000000, +9.00000000
1438+6.00000000, -7.00000000
1439+1.00000000, -2.00000000
1440+8.00000000, +1.00000000
1441do isam = 1, nsam
1442 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
1443end do
1444affinity
1445+5.00000000, -9.00000000
1446+0.00000000, +9.00000000
1447+6.00000000, -7.00000000
1448+1.00000000, -2.00000000
1449+8.00000000, +1.00000000
1450affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
1451affinInv
1452+5.00000000, -9.00000000
1453+0.00000000, +9.00000000
1454+6.00000000, -7.00000000
1455+1.00000000, -2.00000000
1456+8.00000000, +1.00000000
1457sample ! for comparison with affinInv.
1458+5.00000000, -9.00000000
1459+0.00000000, +9.00000000
1460+6.00000000, -7.00000000
1461+1.00000000, -2.00000000
1462+8.00000000, +1.00000000
1463
1464affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
1465affinity
1466+5.00000000, -9.00000000
1467+0.00000000, +9.00000000
1468+6.00000000, -7.00000000
1469+1.00000000, -2.00000000
1470+8.00000000, +1.00000000
1471do isam = 1, nsam
1472 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
1473end do
1474affinity
1475+5.00000000, -9.00000000
1476+0.00000000, +9.00000000
1477+6.00000000, -7.00000000
1478+1.00000000, -2.00000000
1479+8.00000000, +1.00000000
1480affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
1481affinInv
1482+5.00000000, -9.00000000
1483+0.00000000, +9.00000000
1484+6.00000000, -7.00000000
1485+1.00000000, -2.00000000
1486+8.00000000, +1.00000000
1487sample ! for comparison with affinInv.
1488+5.00000000, -9.00000000
1489+0.00000000, +9.00000000
1490+6.00000000, -7.00000000
1491+1.00000000, -2.00000000
1492+8.00000000, +1.00000000
1493
1494dim = 1; ndim = getUnifRand(1, 3); nsam = getUnifRand(4, 5)
1495[dim, ndim, nsam]
1496+1, +3, +5
1497tlate = getUnifRand(-5, +5, ndim)
1498tlate
1499+5.00000000, -2.00000000, -2.00000000
1500sample = getUnifRand(-9, +9, nsam, ndim)
1501sample
1502-5.00000000, +8.00000000, +5.00000000
1503-4.00000000, -7.00000000, -5.00000000
1504+6.00000000, +1.00000000, +7.00000000
1505-6.00000000, +2.00000000, -2.00000000
1506+6.00000000, -2.00000000, +1.00000000
1507call setResized(affinity, shape(sample, IK))
1508call setResized(affinInv, shape(sample, IK))
1509tform = getCovRand(1., ndim)
1510tform
1511+1.00000000, +0.992879331, +0.809138238
1512+0.992879331, +1.00000000, +0.850124955
1513+0.809138238, +0.850124955, +1.00000000
1514affinity = getAffinity(sample, dim, tform) ! whole sample transformation.
1515affinity
1516+6.98872566, +7.28622818, +7.75530863
1517-14.9958458, -15.2221422, -14.1874275
1518+12.6568470, +12.9081507, +12.7049541
1519-5.63251781, -5.65752602, -5.15457916
1520+4.82337952, +4.80740070, +4.15457916
1521do isam = 1, nsam
1522 affinity(isam,:) = getAffinity(sample(isam,:), tform) ! point-wise transformation
1523end do
1524affinity
1525+6.98872566, +7.28622818, +7.75530863
1526-14.9958458, -15.2221422, -14.1874275
1527+12.6568470, +12.9081507, +12.7049541
1528-5.63251781, -5.65752602, -5.15457916
1529+4.82337952, +4.80740070, +4.15457916
1530affinInv = getAffinity(affinity, dim, getMatInv(tform)) ! inverse transformation.
1531affinInv
1532-5.00003815, +8.00009155, +4.99999237
1533-3.99984741, -7.00028992, -4.99996185
1534+5.99987793, +1.00019836, +6.99998474
1535-5.99993134, +1.99996948, -1.99999619
1536+6.00001907, -1.99997711, +0.999994278
1537sample ! for comparison with affinInv.
1538-5.00000000, +8.00000000, +5.00000000
1539-4.00000000, -7.00000000, -5.00000000
1540+6.00000000, +1.00000000, +7.00000000
1541-6.00000000, +2.00000000, -2.00000000
1542+6.00000000, -2.00000000, +1.00000000
1543
1544call setMatInit(tform(2:,1:ndim-1), lowDia, vlow = 0._RKG, vdia = 0._RKG) ! make tform an upper-diagonal matrix.
1545tform
1546+1.00000000, +0.992879331, +0.809138238
1547+0.00000000, +1.00000000, +0.850124955
1548+0.00000000, +0.00000000, +1.00000000
1549affinity = getAffinity(sample, dim, tform, class = upperDiag) ! whole sample transformation.
1550affinity
1551+6.98872566, +12.2506247, +5.00000000
1552-14.9958458, -11.2506247, -5.00000000
1553+12.6568470, +6.95087481, +7.00000000
1554-5.63251781, +0.299750090, -2.00000000
1555+4.82337952, -1.14987504, +1.00000000
1556do isam = 1, nsam
1557 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperDiag) ! point-wise transformation
1558end do
1559affinity
1560+6.98872566, +12.2506247, +5.00000000
1561-14.9958458, -11.2506247, -5.00000000
1562+12.6568470, +6.95087481, +7.00000000
1563-5.63251781, +0.299750090, -2.00000000
1564+4.82337952, -1.14987504, +1.00000000
1565affinInv = getAffinity(affinity, dim, getMatInv(tform, upperDiag), class = upperDiag) ! inverse transformation.
1566affinInv
1567-5.00000000, +8.00000000, +5.00000000
1568-3.99999905, -7.00000000, -5.00000000
1569+6.00000000, +1.00000000, +7.00000000
1570-6.00000000, +2.00000000, -2.00000000
1571+6.00000000, -2.00000000, +1.00000000
1572sample ! for comparison with affinInv.
1573-5.00000000, +8.00000000, +5.00000000
1574-4.00000000, -7.00000000, -5.00000000
1575+6.00000000, +1.00000000, +7.00000000
1576-6.00000000, +2.00000000, -2.00000000
1577+6.00000000, -2.00000000, +1.00000000
1578
1579affinity = getAffinity(sample, dim, tform, class = upperUnit) ! whole sample transformation.
1580affinity
1581+6.98872566, +12.2506247, +5.00000000
1582-14.9958458, -11.2506247, -5.00000000
1583+12.6568470, +6.95087481, +7.00000000
1584-5.63251781, +0.299750090, -2.00000000
1585+4.82337952, -1.14987504, +1.00000000
1586do isam = 1, nsam
1587 affinity(isam,:) = getAffinity(sample(isam,:), tform, upperUnit) ! point-wise transformation
1588end do
1589affinity
1590+6.98872566, +12.2506247, +5.00000000
1591-14.9958458, -11.2506247, -5.00000000
1592+12.6568470, +6.95087481, +7.00000000
1593-5.63251781, +0.299750090, -2.00000000
1594+4.82337952, -1.14987504, +1.00000000
1595affinInv = getAffinity(affinity, dim, getMatInv(tform), class = upperUnit) ! inverse transformation.
1596affinInv
1597-5.00000000, +8.00000000, +5.00000000
1598-3.99999905, -7.00000000, -5.00000000
1599+6.00000000, +1.00000000, +7.00000000
1600-6.00000000, +2.00000000, -2.00000000
1601+6.00000000, -2.00000000, +1.00000000
1602sample ! for comparison with affinInv.
1603-5.00000000, +8.00000000, +5.00000000
1604-4.00000000, -7.00000000, -5.00000000
1605+6.00000000, +1.00000000, +7.00000000
1606-6.00000000, +2.00000000, -2.00000000
1607+6.00000000, -2.00000000, +1.00000000
1608
1609affinity = getAffinity(sample, dim, tform, class = lowerDiag) ! whole sample transformation.
1610affinity
1611-5.00000000, +8.00000000, +5.00000000
1612-4.00000000, -7.00000000, -5.00000000
1613+6.00000000, +1.00000000, +7.00000000
1614-6.00000000, +2.00000000, -2.00000000
1615+6.00000000, -2.00000000, +1.00000000
1616do isam = 1, nsam
1617 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerDiag) ! point-wise transformation
1618end do
1619affinity
1620-5.00000000, +8.00000000, +5.00000000
1621-4.00000000, -7.00000000, -5.00000000
1622+6.00000000, +1.00000000, +7.00000000
1623-6.00000000, +2.00000000, -2.00000000
1624+6.00000000, -2.00000000, +1.00000000
1625affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerDiag) ! inverse transformation.
1626affinInv
1627-5.00000000, +8.00000000, +5.00000000
1628-4.00000000, -7.00000000, -5.00000000
1629+6.00000000, +1.00000000, +7.00000000
1630-6.00000000, +2.00000000, -2.00000000
1631+6.00000000, -2.00000000, +1.00000000
1632sample ! for comparison with affinInv.
1633-5.00000000, +8.00000000, +5.00000000
1634-4.00000000, -7.00000000, -5.00000000
1635+6.00000000, +1.00000000, +7.00000000
1636-6.00000000, +2.00000000, -2.00000000
1637+6.00000000, -2.00000000, +1.00000000
1638
1639affinity = getAffinity(sample, dim, tform, class = lowerUnit) ! whole sample transformation.
1640affinity
1641-5.00000000, +8.00000000, +5.00000000
1642-4.00000000, -7.00000000, -5.00000000
1643+6.00000000, +1.00000000, +7.00000000
1644-6.00000000, +2.00000000, -2.00000000
1645+6.00000000, -2.00000000, +1.00000000
1646do isam = 1, nsam
1647 affinity(isam,:) = getAffinity(sample(isam,:), tform, lowerUnit) ! point-wise transformation
1648end do
1649affinity
1650-5.00000000, +8.00000000, +5.00000000
1651-4.00000000, -7.00000000, -5.00000000
1652+6.00000000, +1.00000000, +7.00000000
1653-6.00000000, +2.00000000, -2.00000000
1654+6.00000000, -2.00000000, +1.00000000
1655affinInv = getAffinity(affinity, dim, getMatInv(tform), class = lowerUnit) ! inverse transformation.
1656affinInv
1657-5.00000000, +8.00000000, +5.00000000
1658-4.00000000, -7.00000000, -5.00000000
1659+6.00000000, +1.00000000, +7.00000000
1660-6.00000000, +2.00000000, -2.00000000
1661+6.00000000, -2.00000000, +1.00000000
1662sample ! for comparison with affinInv.
1663-5.00000000, +8.00000000, +5.00000000
1664-4.00000000, -7.00000000, -5.00000000
1665+6.00000000, +1.00000000, +7.00000000
1666-6.00000000, +2.00000000, -2.00000000
1667+6.00000000, -2.00000000, +1.00000000
1668

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: