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

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: