Line data Source code
1 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3 : !!!! !!!!
4 : !!!! ParaMonte: Parallel Monte Carlo and Machine Learning Library. !!!!
5 : !!!! !!!!
6 : !!!! Copyright (C) 2012-present, The Computational Data Science Lab !!!!
7 : !!!! !!!!
8 : !!!! This file is part of the ParaMonte library. !!!!
9 : !!!! !!!!
10 : !!!! LICENSE !!!!
11 : !!!! !!!!
12 : !!!! https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md !!!!
13 : !!!! !!!!
14 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16 :
17 : !> \brief
18 : !> This file contains the implementation details of the routines under the generic interfaces of [pm_sampleShift](@ref pm_sampleShift).
19 : !>
20 : !> \finmain
21 : !>
22 : !> \author
23 : !> \AmirShahmoradi, Saturday 2:33 AM, August 22, 2021, Dallas, TX
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : ! Define the runtime checks.
28 : #define CHECK_VAL_DIM \
29 : CHECK_ASSERTION(__LINE__, 1 <= dim .and. dim <= rank(sample), \
30 : SK_"@setShifted(): The condition `1 <= dim .and. dim <= rank(sample)` must hold. dim, rank(sample) = "//getStr([integer(IK) :: dim, rank(sample)]))
31 : #define CHECK_LEN_AMOUNT(DIM) \
32 : CHECK_ASSERTION(__LINE__, size(sample, DIM, IK) == size(amount, 1, IK), \
33 : SK_"@setShifted(): The condition `size(sample, dim, 1) == size(amount, 1)` must hold. dim, shape(sample), size(amount) = "\
34 : //getStr([DIM, shape(sample, IK), size(amount, 1, IK)]))
35 : ! Set the indexing rule.
36 : #if ONO_ENABLED
37 : #define GET_INDEX(I,J)I,J
38 : #elif OTH_ENABLED
39 : #define GET_INDEX(I,J)J,I
40 : #elif getShifted_ENABLED && D2_ENABLED
41 : #error "Unrecognized interface."
42 : #endif
43 : ! Set the conjugation rule.
44 : #if OTH_ENABLED && CK_ENABLED
45 : #define GET_CONJG(X)conjg(X)
46 : #elif ONO_ENABLED || (OTH_ENABLED && RK_ENABLED)
47 : #define GET_CONJG(X)X
48 : #elif getShifted_ENABLED && D2_ENABLED
49 : #error "Unrecognized interface."
50 : #endif
51 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
52 : #if getShifted_ENABLED && D1_ENABLED
53 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
54 :
55 46373 : sampleShifted = sample + amount
56 :
57 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
58 : #elif setShifted_ENABLED && (D1_ENABLED || (D2_ENABLED && ALL_ENABLED))
59 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
60 :
61 658 : sample = sample + amount
62 :
63 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
64 : #elif getShifted_ENABLED && D2_ENABLED && ALL_ENABLED
65 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
66 :
67 : integer(IK) :: idim, jdim, ndim, mdim
68 : ndim = size(sample, 1, IK)
69 : mdim = size(sample, 2, IK)
70 : do concurrent(idim = 1 : ndim, jdim = 1 : mdim)
71 : #if ONO_ENABLED
72 36132 : sampleShifted(idim, jdim) = GET_CONJG(sample(idim, jdim) + amount)
73 : #elif OTH_ENABLED
74 4128 : sampleShifted(jdim, idim) = GET_CONJG(sample(idim, jdim) + amount)
75 : #else
76 : #error "Unrecognized interface."
77 : #endif
78 : end do
79 :
80 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
81 : #elif getShifted_ENABLED && D2_ENABLED && DIM_ENABLED
82 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
83 :
84 : integer(IK) :: idim, isam, ndim, nsam
85 43728 : CHECK_VAL_DIM
86 102032 : CHECK_LEN_AMOUNT(3 - dim)
87 14576 : nsam = size(sample, dim, IK)
88 14576 : ndim = size(sample, 3 - dim, IK)
89 14576 : if (dim == 2_IK) then
90 : do concurrent(idim = 1 : ndim, isam = 1 : nsam)
91 135514 : sampleShifted(GET_INDEX(idim, isam)) = GET_CONJG(sample(idim, isam) + amount(idim))
92 : end do
93 : else
94 : do concurrent(idim = 1 : ndim, isam = 1 : nsam)
95 135541 : sampleShifted(GET_INDEX(isam, idim)) = GET_CONJG(sample(isam, idim) + amount(idim))
96 : end do
97 : end if
98 :
99 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
100 : #elif setShifted_ENABLED && D2_ENABLED
101 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
102 :
103 : integer(IK) :: ell, ndim, nsam
104 8028 : CHECK_VAL_DIM
105 18732 : CHECK_LEN_AMOUNT(3 - dim)
106 2676 : nsam = size(sample, dim, IK)
107 2676 : ndim = size(sample, 3 - dim, IK)
108 2676 : if (dim == 2_IK) then
109 : do concurrent(ell = 1 : nsam)
110 17590 : sample(1 : ndim, ell) = sample(1 : ndim, ell) + amount
111 : end do
112 : else
113 : do concurrent(ell = 1 : ndim)
114 18692 : sample(1 : nsam, ell) = sample(1 : nsam, ell) + amount(ell)
115 : end do
116 : end if
117 :
118 : #else
119 : !%%%%%%%%%%%%%%%%%%%%%%%%
120 : #error "Unrecognized interface."
121 : !%%%%%%%%%%%%%%%%%%%%%%%%
122 : #endif
123 : #undef CHECK_LEN_AMOUNT
124 : #undef CHECK_VAL_DIM
125 : #undef GET_CONJG
126 : #undef GET_INDEX
|