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

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

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
[out]affinity: The output object of the same type and kind and rank and shape as sample, containing the affine-transformed sample.
[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.
[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)].)


Possible calling interfaces

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

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
Todo:
High Priority: The performance of the algorithm for the case dim = 1 can be further improved by looping over the sample in the innermost layer.


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, Thursday 2:45 AM, August 19, 2021, Dallas, TX
Fatemeh Bagheri, Wednesday 00:01 AM, August 25, 2021, Dallas, TX

Definition at line 502 of file pm_sampleAffinity.F90.


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