1 //c99 -o test_vswf_translations -ggdb -I .. test_vswf_translations.c -lqpms -lgsl -lm -lblas
2 #include <qpms/translations.h>
5 #include <gsl/gsl_rng.h>
6 #include <gsl/gsl_math.h>
7 #include <gsl/gsl_randist.h>
9 #include <qpms/vectors.h>
11 #include <qpms/indexing.h>
13 char *normstr(qpms_normalisation_t norm
) {
14 //int csphase = qpms_normalisation_t_csphase(norm);
15 norm
= norm
& QPMS_NORMALISATION_NORM_BITS
;
17 case QPMS_NORMALISATION_NORM_NONE
:
19 case QPMS_NORMALISATION_NORM_SPHARM
:
21 case QPMS_NORMALISATION_NORM_POWER
:
28 int test_sphwave_translation(const qpms_trans_calculator
*c
, qpms_bessel_t wavetype
,
29 cart3_t o2minuso1
, int npoints
, cart3_t
*o1points
);
30 //int test_planewave_decomposition(cart3_t k, ccart3_t E, qpms_l_t lMax, qpms_normalisation_t norm, int npoints, cart3_t *points);
31 //int test_planewave_decomposition_silent(cart3_t k, ccart3_t E, qpms_l_t lMax, qpms_normalisation_t norm, int npoints, cart3_t *points, double relerrthreshold, double *relerrs);
34 gsl_rng
*rng
= gsl_rng_alloc(gsl_rng_ranlxs0
);
35 gsl_rng_set(rng
, 666);
38 //qpms_l_t viewlMax = 2;
41 //double shiftsigma = 2.;
44 o2minuso1
.x
= 1; //gsl_ran_gaussian(rng, shiftsigma);
45 o2minuso1
.y
= 2; //gsl_ran_gaussian(rng, shiftsigma);
46 o2minuso1
.z
= 5; //gsl_ran_gaussian(rng, shiftsigma);
48 cart3_t points
[npoints
];
49 double relerrs
[npoints
];
50 memset(points
, 0, npoints
* sizeof(cart3_t
));
51 points
[0].x
= points
[1].y
= points
[2].z
= sigma
;
52 double relerrthreshold
= 1e-11;
53 for (unsigned i
= 3; i
< npoints
; ++i
) {
54 cart3_t
*w
= points
+i
;
55 w
->x
= gsl_ran_gaussian(rng
, sigma
);
56 w
->y
= gsl_ran_gaussian(rng
, sigma
);
57 w
->z
= gsl_ran_gaussian(rng
, sigma
);
60 for(int use_csbit
= 0; use_csbit
<= 1; ++use_csbit
) {
61 for(int i
= 0; i
< 3; ++i
){
62 qpms_normalisation_t norm
= ((qpms_normalisation_t
[])
63 { QPMS_NORMALISATION_NORM_SPHARM
,
64 QPMS_NORMALISATION_NORM_POWER
,
65 QPMS_NORMALISATION_NORM_NONE
67 | (use_csbit
? QPMS_NORMALISATION_CSPHASE
: 0);
68 qpms_trans_calculator
*c
= qpms_trans_calculator_init(lMax
, norm
);
69 for(int J
= 1; J
<= 4; ++J
)
70 test_sphwave_translation(c
, J
, o2minuso1
, npoints
, points
);
71 qpms_trans_calculator_free(c
);
78 int test_sphwave_translation(const qpms_trans_calculator
*c
, qpms_bessel_t wavetype
,
79 cart3_t sc
, int npoints
, cart3_t
*points
) {
80 puts("==============================================================");
81 printf("Test translation o2-o1 = %fx̂ + %fŷ + %fẑ", sc
.x
, sc
.y
, sc
.z
);
82 sph_t ss
= cart2sph(sc
);
83 printf(" = %fr̂ @ o1, θ = %f, φ = %f\n", ss
.r
, ss
.theta
, ss
.phi
);
84 printf("lMax = %d, norm: %s, csphase = %d\n",
85 (int)c
->lMax
, normstr(c
->normalisation
), qpms_normalisation_t_csphase(c
->normalisation
));
86 printf("wave type J = %d\n", wavetype
);
88 qpms_l_t lMax
= c
->lMax
;
89 qpms_y_t nelem
= c
->nelem
;
90 csphvec_t N1
[nelem
], /* N2[nelem], */ M1
[nelem
] /*, M2[nelem]*/;
92 for (int i
= 0; i
< npoints
; i
++) {
93 printf("-------- Point %d --------\n", i
);
94 cart3_t w1c
= points
[i
];
95 cart3_t w2c
= cart3_add(w1c
, cart3_scale(-1, sc
));
96 sph_t w1s
= cart2sph(w1c
);
97 sph_t w2s
= cart2sph(w2c
);
98 printf(" = %fx̂ + %fŷ + %fẑ @o1\n", w1c
.x
, w1c
.y
, w1c
.z
);
99 printf(" = %fx̂ + %fŷ + %fẑ @o2\n", w2c
.x
, w2c
.y
, w2c
.z
);
100 printf(" = %fr̂ @ o1, θ = %f, φ = %f\n", w1s
.r
, w1s
.theta
, w1s
.phi
);
101 printf(" = %fr̂ @ o2, θ = %f, φ = %f\n", w2s
.r
, w2s
.theta
, w2s
.phi
);
102 printf("Outside the sphere centered in o1 intersecting o2: %s; by %f\n", (w1s
.r
> ss
.r
) ? "true" : "false",
104 if(QPMS_SUCCESS
!= qpms_vswf_fill(NULL
, M1
, N1
, lMax
, w1s
, wavetype
, c
->normalisation
))
105 abort(); // original wave set
107 for(qpms_y_t y1
= 0; y1
< nelem
; ++y1
) { //index of the wave originating in o1 that will be reconstructed in o2
110 qpms_y2mn_p(y1
, &m1
, &l1
);
111 printf("*** wave l = %d, m = %d ***\n", l1
, m1
);
113 complex double A
[nelem
], B
[nelem
];
114 for(qpms_y_t y2
= 0; y2
< nelem
; ++y2
){
115 qpms_m_t m2
; qpms_l_t l2
;
116 qpms_y2mn_p(y2
, &m2
, &l2
);
117 if(qpms_trans_calculator_get_AB_p(c
, &(A
[y2
]), &(B
[y2
]), m2
, l2
, m1
, l1
, sph2csph(ss
), (w1s
.r
> ss
.r
) , wavetype
))
122 print_csphvec(M1
[y1
]);
123 printf(" @ o1\n = ");
124 ccart3_t M1c
= csphvec2ccart(M1
[y1
], w1s
);
127 csphvec_t M1s2
= ccart2csphvec(M1c
, w2s
);
130 csphvec_t M2s2
= qpms_eval_vswf(w2s
, NULL
, A
, B
, lMax
, wavetype
, c
->normalisation
);
136 print_csphvec(N1
[y1
]);
137 printf(" @ o1\n = ");
138 ccart3_t N1c
= csphvec2ccart(N1
[y1
], w1s
);
141 csphvec_t N1s2
= ccart2csphvec(N1c
, w2s
);
143 printf(" @o2\nNr= ");
144 csphvec_t N2s2
= qpms_eval_vswf(w2s
, NULL
, B
, A
, lMax
, wavetype
, c
->normalisation
);
150 return 0; // FIXME something more meaningful here...
159 int test_planewave_decomposition(cart3_t k
, ccart3_t E
, qpms_l_t lMax
, qpms_normalisation_t norm
, int npoints
, cart3_t
*points
){
160 qpms_y_t nelem
= qpms_lMax2nelem(lMax
);
161 complex double lc
[nelem
], mc
[nelem
], ec
[nelem
];
163 qpms_planewave2vswf_fill_cart(k
, E
, lc
, mc
, ec
, lMax
, norm
)) {
167 printf("==============================================================\n");
168 printf("Test wave k = %fx̂ + %fŷ + %fẑ", k
.x
, k
.y
, k
.z
);
169 sph_t k_sph
= cart2sph(k
);
170 printf(" = %fr̂ @ θ = %f, φ = %f\n", k_sph
.r
, k_sph
.theta
, k_sph
.phi
);
171 printf(" E_0 = (%f+%fj)x̂ + (%f+%fj)ŷ + (%f+%fj)ẑ",
172 creal(E
.x
),cimag(E
.x
),
173 creal(E
.y
),cimag(E
.y
),
174 creal(E
.z
),cimag(E
.z
));
175 csphvec_t E_s
= ccart2csphvec(E
, k_sph
);
176 printf(" = (%f+%fj)r̂ + (%f+%fj)θ̂ + (%f+%fj)φ̂ @ k\n",
177 creal(E_s
.rc
), cimag(E_s
.rc
), creal(E_s
.thetac
), cimag(E_s
.thetac
),
178 creal(E_s
.phic
), cimag(E_s
.phic
));
179 printf(" lMax = %d, norm: %s, csphase = %d\n",
180 (int)lMax
, normstr(norm
), qpms_normalisation_t_csphase(norm
));
182 for(qpms_y_t y
= 0; y
< nelem
; ++y
) printf("%g+%gj ", creal(lc
[y
]), cimag(lc
[y
]));
184 for(qpms_y_t y
= 0; y
< nelem
; ++y
) printf("%g+%gj ", creal(mc
[y
]), cimag(mc
[y
]));
186 for(qpms_y_t y
= 0; y
< nelem
; ++y
) printf("%g+%gj ", creal(ec
[y
]), cimag(ec
[y
]));
188 for (int i
= 0; i
< npoints
; i
++) {
189 cart3_t w
= points
[i
];
190 sph_t w_sph
= cart2sph(w
);
191 printf("Point %d: x = %f, y = %f, z = %f,\n", i
, w
.x
, w
.y
, w
.z
);
192 printf(" |r| = %f, θ = %f, φ = %f:\n", w_sph
.r
, w_sph
.theta
, w_sph
.phi
);
193 double phase
= cart3_dot(k
,w
);
194 printf(" k.r = %f\n", phase
);
195 complex double phfac
= cexp(phase
* I
);
196 ccart3_t Ew
= ccart3_scale(phfac
, E
);
197 printf(" pw E(r) = (%f+%fj)x̂ + (%f+%fj)ŷ + (%f+%fj)ẑ",
198 creal(Ew
.x
),cimag(Ew
.x
),
199 creal(Ew
.y
),cimag(Ew
.y
),
200 creal(Ew
.z
),cimag(Ew
.z
));
201 csphvec_t Ew_s
= ccart2csphvec(Ew
, w_sph
);
202 printf(" = (%f+%fj)r̂ + (%f+%fj)θ̂ + (%f+%fj)φ̂ @ r\n",
203 creal(Ew_s
.rc
), cimag(Ew_s
.rc
),
204 creal(Ew_s
.thetac
), cimag(Ew_s
.thetac
),
205 creal(Ew_s
.phic
), cimag(Ew_s
.phic
));
206 w_sph
.r
*= k_sph
.r
; /// NEVER FORGET THIS!!!
207 csphvec_t Ew_s_recomp
= qpms_eval_vswf(w_sph
,
208 lc
, mc
, ec
, lMax
, QPMS_BESSEL_REGULAR
, norm
);
209 ccart3_t Ew_recomp
= csphvec2ccart(Ew_s_recomp
, w_sph
);
210 printf(" rec E(r) = (%f+%fj)x̂ + (%f+%fj)ŷ + (%f+%fj)ẑ",
211 creal(Ew_recomp
.x
),cimag(Ew_recomp
.x
),
212 creal(Ew_recomp
.y
),cimag(Ew_recomp
.y
),
213 creal(Ew_recomp
.z
),cimag(Ew_recomp
.z
));
214 printf(" = (%f+%fj)r̂ + (%f+%fj)θ̂ + (%f+%fj)φ̂ @ r\n",
215 creal(Ew_s_recomp
.rc
), cimag(Ew_s_recomp
.rc
),
216 creal(Ew_s_recomp
.thetac
), cimag(Ew_s_recomp
.thetac
),
217 creal(Ew_s_recomp
.phic
), cimag(Ew_s_recomp
.phic
));
218 double relerrfac
= 2/(cabs(Ew_s_recomp
.rc
) + cabs(Ew_s
.rc
)
219 +cabs(Ew_s_recomp
.thetac
) + cabs(Ew_s
.thetac
)
220 +cabs(Ew_s_recomp
.phic
) + cabs(Ew_s
.phic
));
221 printf(" rel. err. magnitude: %g @ r̂, %g @ θ̂, %g @ φ̂\n",
222 cabs(Ew_s_recomp
.rc
- Ew_s
.rc
) * relerrfac
,
223 cabs(Ew_s_recomp
.thetac
- Ew_s
.thetac
) * relerrfac
,
224 cabs(Ew_s_recomp
.phic
- Ew_s
.phic
) * relerrfac
230 int test_planewave_decomposition_silent(cart3_t k
, ccart3_t E
, qpms_l_t lMax
, qpms_normalisation_t norm
, int npoints
, cart3_t
*points
, double relerrthreshold
, double *relerrs
) {
231 qpms_y_t nelem
= qpms_lMax2nelem(lMax
);
233 complex double lc
[nelem
], mc
[nelem
], ec
[nelem
];
235 qpms_planewave2vswf_fill_cart(k
, E
, lc
, mc
, ec
, lMax
, norm
)) {
239 sph_t k_sph
= cart2sph(k
);
240 csphvec_t E_s
= ccart2csphvec(E
, k_sph
);
241 for (int i
= 0; i
< npoints
; i
++) {
242 cart3_t w
= points
[i
];
243 sph_t w_sph
= cart2sph(w
);
245 double phase
= cart3_dot(k
,w
);
246 complex double phfac
= cexp(phase
* I
);
247 ccart3_t Ew
= ccart3_scale(phfac
, E
);
248 csphvec_t Ew_s
= ccart2csphvec(Ew
, w_sph
);
249 csphvec_t Ew_s_recomp
= qpms_eval_vswf(w_sph
,
250 lc
, mc
, ec
, lMax
, QPMS_BESSEL_REGULAR
, norm
);
251 ccart3_t Ew_recomp
= csphvec2ccart(Ew_s_recomp
, w_sph
);
252 double relerrfac
= 2/(cabs(Ew_s_recomp
.rc
) + cabs(Ew_s
.rc
)
253 +cabs(Ew_s_recomp
.thetac
) + cabs(Ew_s
.thetac
)
254 +cabs(Ew_s_recomp
.phic
) + cabs(Ew_s
.phic
));
256 double relerr
= (cabs(Ew_s_recomp
.rc
- Ew_s
.rc
)
257 + cabs(Ew_s_recomp
.thetac
- Ew_s
.thetac
)
258 + cabs(Ew_s_recomp
.phic
- Ew_s
.phic
)
260 if(relerrs
) relerrs
[i
] = relerr
;
261 if(relerr
> relerrthreshold
) ++failcount
;