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

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: