Disable "hard" examples in CI
[qpms.git] / qpms / tmatrices.c
blob5ae6beb8da95d32be81ca6035f901577ddecd05b
1 #define _POSIX_C_SOURCE 200809L // for getline()
2 #define lapack_int int
3 #define lapack_complex_double complex double
4 #define lapack_complex_double_real(z) (creal(z))
5 #define lapack_complex_double_imag(z) (cimag(z))
6 #include <lapacke.h>
7 #include <stdio.h>
8 #include <stdlib.h>
9 #include <cblas.h>
10 #include <unistd.h>
11 #include "scatsystem.h"
12 #include "indexing.h"
13 #include "vswf.h"
14 #include "groups.h"
15 #include "symmetries.h"
16 #include <assert.h>
17 #include "vectors.h"
18 #include "quaternions.h"
19 #include <string.h>
20 #include "qpms_error.h"
21 #include "tmatrices.h"
22 #include "qpms_specfunc.h"
23 #include "normalisation.h"
24 #include <errno.h>
25 #include <pthread.h>
26 #include <gsl/gsl_integration.h>
27 #include <gsl/gsl_errno.h>
28 #include "optim.h"
29 #include "oshacks.h"
31 // These are quite arbitrarily chosen constants for the quadrature in qpms_tmatrix_axialsym_fill()
32 #define TMATRIX_AXIALSYM_INTEGRAL_EPSREL (1e-5)
33 #define TMATRIX_AXIALSYM_INTEGRAL_EPSABS (1e-8)
35 #define HBAR (1.05457162825e-34)
36 #define ELECTRONVOLT (1.602176487e-19)
37 #define SCUFF_OMEGAUNIT (3e14)
39 qpms_tmatrix_t *qpms_tmatrix_init(const qpms_vswf_set_spec_t *bspec) {
40 qpms_tmatrix_t *t;
41 QPMS_CRASHING_MALLOC(t, sizeof(qpms_tmatrix_t));
42 t->spec = bspec;
43 size_t n = bspec->n;
44 QPMS_CRASHING_CALLOC(t->m, n*n, sizeof(complex double));
45 t->owns_m = true;
46 return t;
49 qpms_tmatrix_t *qpms_tmatrix_init_from_generator(
50 const qpms_vswf_set_spec_t *bspec,
51 qpms_tmatrix_generator_t gen,
52 complex double omega) {
53 qpms_tmatrix_t *t = qpms_tmatrix_init(bspec);
54 QPMS_ENSURE_SUCCESS(gen.function(t, omega, gen.params));
55 return t;
58 qpms_tmatrix_t *qpms_tmatrix_copy(const qpms_tmatrix_t *T) {
59 qpms_tmatrix_t *t = qpms_tmatrix_init(T->spec);
60 size_t n = T->spec->n;
61 for(size_t i = 0; i < n*n; ++i)
62 t->m[i] = T->m[i];
63 return t;
66 void qpms_tmatrix_free(qpms_tmatrix_t *t){
67 if(t && t->owns_m) free(t->m);
68 free(t);
71 qpms_tmatrix_t *qpms_tmatrix_apply_symop_inplace(
72 qpms_tmatrix_t *T,
73 const complex double *M
76 //qpms_tmatrix_t *t = qpms_tmatrix_init(T->spec);
77 const size_t n = T->spec->n;
78 complex double tmp[n][n];
79 // tmp = M T
80 const complex double one = 1, zero = 0;
81 cblas_zgemm(CblasRowMajor,
82 CblasNoTrans,
83 CblasNoTrans,
84 n, n, n, &one, M, n, T->m, n, &zero, tmp, n);
85 // t->m = tmp M* = M T M*
86 cblas_zgemm(CblasRowMajor,
87 CblasNoTrans,
88 CblasConjTrans,
89 n, n, n, &one, tmp, n, M, n, &zero, T->m, n);
90 return T;
93 qpms_tmatrix_t *qpms_tmatrix_apply_symop(
94 const qpms_tmatrix_t *T,
95 const complex double *M
98 qpms_tmatrix_t *t = qpms_tmatrix_init(T->spec);
99 const size_t n = T->spec->n;
100 complex double tmp[n][n];
101 // tmp = M T
102 const complex double one = 1, zero = 0;
103 cblas_zgemm(CblasRowMajor,
104 CblasNoTrans,
105 CblasNoTrans,
106 n, n, n, &one, M, n, T->m, n, &zero, tmp, n);
107 // t->m = tmp M* = M T M*
108 cblas_zgemm(CblasRowMajor,
109 CblasNoTrans,
110 CblasConjTrans,
111 n, n, n, &one, tmp, n, M, n, &zero, t->m, n);
112 return t;
115 qpms_errno_t qpms_symmetrise_tmdata_irot3arr(
116 complex double *tmdata, const size_t tmcount,
117 const qpms_vswf_set_spec_t *bspec,
118 const size_t n_symops, const qpms_irot3_t *symops) {
119 const size_t n = bspec->n;
120 qpms_tmatrix_t *tmcopy = qpms_tmatrix_init(bspec);
121 complex double *symop_matrices = malloc(n*n*sizeof(complex double) * n_symops);
122 if(!symop_matrices) qpms_pr_error_at_flf(__FILE__, __LINE__, __func__,
123 "malloc() failed.");
124 for (size_t i = 0; i < n_symops; ++i)
125 qpms_irot3_uvswfi_dense(symop_matrices + i*n*n, bspec, symops[i]);
126 complex double tmp[n][n];
127 const complex double one = 1, zero = 0;
128 for (size_t tmi = 0; tmi < tmcount; ++tmi) {
129 // Move the data in tmcopy; we will then write the sum directly into tmdata.
130 memcpy(tmcopy->m, tmdata+n*n*tmi, n*n*sizeof(complex double));
131 memset(tmdata+n*n*tmi, 0, n*n*sizeof(complex double));
132 for (size_t i = 0; i < n_symops; ++i) {
133 const complex double *const M = symop_matrices + i*n*n;
134 // tmp = M T
135 cblas_zgemm(CblasRowMajor, CblasNoTrans, CblasNoTrans,
136 n, n, n, &one, M, n, tmcopy->m, n, &zero, tmp, n);
137 // tmdata[...] += tmp M* = M T M*
138 cblas_zgemm(CblasRowMajor, CblasNoTrans, CblasConjTrans,
139 n, n, n, &one, tmp, n, M, n, &one, tmdata + tmi*n*n, n);
141 for (size_t ii = 0; ii < n*n; ++ii)
142 tmdata[n*n*tmi+ii] /= n_symops;
144 free(symop_matrices);
145 qpms_tmatrix_free(tmcopy);
146 return QPMS_SUCCESS;
149 qpms_errno_t qpms_symmetrise_tmdata_finite_group(
150 complex double *tmdata, const size_t tmcount,
151 const qpms_vswf_set_spec_t *bspec,
152 const qpms_finite_group_t *pointgroup) {
153 if (!(pointgroup->rep3d)) qpms_pr_error_at_flf(__FILE__, __LINE__, __func__,
154 "This function requires pointgroup->rep3d to be set correctly!");
155 return qpms_symmetrise_tmdata_irot3arr(tmdata, tmcount, bspec,
156 pointgroup->order, pointgroup->rep3d);
159 qpms_tmatrix_t *qpms_tmatrix_symmetrise_irot3arr_inplace(
160 qpms_tmatrix_t *T,
161 size_t n_symops,
162 const qpms_irot3_t *symops
164 if(qpms_symmetrise_tmdata_irot3arr(T->m, 1,
165 T->spec, n_symops, symops) != QPMS_SUCCESS)
166 return NULL;
167 else return T;
170 qpms_tmatrix_t *qpms_tmatrix_symmetrise_finite_group_inplace(
171 qpms_tmatrix_t *T,
172 const qpms_finite_group_t *pointgroup
174 if(qpms_symmetrise_tmdata_finite_group(T->m, 1,
175 T->spec, pointgroup) != QPMS_SUCCESS)
176 return NULL;
177 else return T;
180 qpms_tmatrix_t *qpms_tmatrix_symmetrise_involution_inplace(
181 qpms_tmatrix_t *T,
182 const complex double *M
185 qpms_tmatrix_t *t = qpms_tmatrix_apply_symop(T, M);
186 const size_t n = T->spec->n;
187 for(size_t i = 0; i < n*n; ++i)
188 T->m[i] = 0.5 * (t->m[i] + T->m[i]);
189 qpms_tmatrix_free(t);
190 return T;
193 qpms_tmatrix_t *qpms_tmatrix_symmetrise_involution(
194 const qpms_tmatrix_t *T,
195 const complex double *M
198 qpms_tmatrix_t *t = qpms_tmatrix_init(T->spec);
199 const size_t n = T->spec->n;
200 complex double tmp[n][n];
201 // tmp = M T
202 const complex double one = 1, zero = 0;
203 cblas_zgemm(CblasRowMajor,
204 CblasNoTrans,
205 CblasNoTrans,
206 n, n, n, &one, M, n, T->m, n, &zero, tmp, n);
207 // t->m = tmp M* = M T M*
208 cblas_zgemm(CblasRowMajor,
209 CblasNoTrans,
210 CblasConjTrans,
211 n, n, n, &one, tmp, n, M, n, &zero, t->m, n);
212 for(size_t i = 0; i < n*n; ++i)
213 t->m[i] = 0.5 * (t->m[i] + T->m[i]);
214 return t;
217 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_inf(const qpms_tmatrix_t *T) {
218 qpms_tmatrix_t *t = qpms_tmatrix_copy(T);
219 return qpms_tmatrix_symmetrise_C_inf_inplace(t);
222 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_inf_inplace(qpms_tmatrix_t *T) {
223 const size_t n = T->spec->n;
224 for (size_t row = 0; row < n; row++) {
225 qpms_m_t rm = qpms_uvswfi2m(T->spec->ilist[row]);
226 for (size_t col = 0; col < n; col++) {
227 qpms_m_t cm = qpms_uvswfi2m(T->spec->ilist[col]);
228 if (rm == cm)
229 ;// No-op // t->m[n*row + col] = T->m[n*row + col];
230 else
231 T->m[n*row + col] = 0;
234 return T;
237 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_N(const qpms_tmatrix_t *T, int N) {
238 qpms_tmatrix_t *t = qpms_tmatrix_copy(T);
239 return qpms_tmatrix_symmetrise_C_N_inplace(t, N);
242 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_N_inplace(qpms_tmatrix_t *T, int N) {
243 const size_t n = T->spec->n;
244 for (size_t row = 0; row < n; row++) {
245 qpms_m_t rm = qpms_uvswfi2m(T->spec->ilist[row]);
246 for (size_t col = 0; col < n; col++) {
247 qpms_m_t cm = qpms_uvswfi2m(T->spec->ilist[col]);
248 if (((rm - cm) % N) == 0)
249 ; // T->m[n*row + col] = T->m[n*row + col];
250 else
251 T->m[n*row + col] = 0;
254 return T;
257 bool qpms_tmatrix_isclose(const qpms_tmatrix_t *A, const qpms_tmatrix_t *B,
258 const double rtol, const double atol)
260 if (!qpms_vswf_set_spec_isidentical(A->spec, B->spec))
261 return false;
262 if (A->m == B->m)
263 return true;
264 const size_t n = A->spec->n;
265 for (size_t i = 0; i < n*n; ++i) {
266 const double tol = atol + rtol * (cabs(B->m[i]));
267 if ( cabs(B->m[i] - A->m[i]) > tol )
268 return false;
270 return true;
273 qpms_tmatrix_interpolator_t *qpms_tmatrix_interpolator_create(const size_t incount,
274 const double *freqs, const qpms_tmatrix_t *ta, const gsl_interp_type *iptype//, const bool copy_bspec
276 if (incount <= 0) return NULL;
277 qpms_tmatrix_interpolator_t *ip = malloc(sizeof(qpms_tmatrix_interpolator_t));
279 if (copy_bspec) {
280 ip->bspec = qpms_vswf_set_spec_copy(ta[0].spec);
281 ip->owns_bspec = true;
283 else {
285 ip->bspec = ta[0].spec;
286 // ip->owns_bspec = false;
288 const size_t n = ip->bspec->n;
290 // check if all matrices have the same bspec
291 for (size_t i = 0; i < incount; ++i)
292 QPMS_ENSURE(qpms_vswf_set_spec_isidentical(ip->bspec, ta[i].spec),
293 "All T-matrices for interpolation must have the same basis!");
295 QPMS_CRASHING_CALLOC(ip->splines_real, n*n, sizeof(gsl_spline *));
296 QPMS_CRASHING_CALLOC(ip->splines_imag, n*n,sizeof(gsl_spline *));
297 for (size_t row = 0; row < n; ++row)
298 for (size_t col = 0; col < n; ++col) {
299 double y_real[incount], y_imag[incount];
300 bool n0_real = false, n0_imag = false;
301 for (size_t i = 0; i < incount; ++i) {
302 complex double telem = ta[i].m[n * row + col];
303 if ((y_real[i] = creal(telem))) n0_real = true;
304 if ((y_imag[i] = cimag(telem))) n0_imag = true;
306 if (n0_real) {
307 gsl_spline *s =
308 ip->splines_real[n * row + col] = gsl_spline_alloc(iptype, incount);
309 QPMS_ENSURE_SUCCESS(gsl_spline_init(s, freqs, y_real, incount));
311 else ip->splines_real[n * row + col] = NULL;
312 if (n0_imag) {
313 gsl_spline *s =
314 ip->splines_imag[n * row + col] = gsl_spline_alloc(iptype, incount);
315 QPMS_ENSURE_SUCCESS(gsl_spline_init(s, freqs, y_imag, incount));
317 else ip->splines_imag[n * row + col] = NULL;
319 return ip;
322 void qpms_tmatrix_interpolator_free(qpms_tmatrix_interpolator_t *ip) {
323 if (ip) {
324 const size_t n = ip->bspec->n;
325 for (size_t i = 0; i < n*n; ++i) {
326 if (ip->splines_real[i]) gsl_spline_free(ip->splines_real[i]);
327 if (ip->splines_imag[i]) gsl_spline_free(ip->splines_imag[i]);
329 //if (ip->owns_bspec)
330 // qpms_vswf_set_spec_free(ip->bspec);
331 free(ip);
335 qpms_tmatrix_t *qpms_tmatrix_interpolator_eval(const qpms_tmatrix_interpolator_t *ip, double freq) {
336 qpms_tmatrix_t *t = qpms_tmatrix_init(ip->bspec);
337 QPMS_ENSURE_SUCCESS(qpms_tmatrix_interpolator_eval_fill(t, ip, freq));
338 return t;
341 qpms_errno_t qpms_tmatrix_interpolator_eval_fill(qpms_tmatrix_t *t,
342 const qpms_tmatrix_interpolator_t *ip, double freq) {
343 QPMS_ENSURE(qpms_vswf_set_spec_isidentical(t->spec, ip->bspec), "Tried to fill a T-matrix with an incompatible interpolator!");
344 const size_t n = ip->bspec->n;
345 for (size_t i = 0; i < n*n; ++i){
346 if (ip->splines_real[i]) t->m[i] = gsl_spline_eval(ip->splines_real[i], freq, NULL /*does this work?*/);
347 if (ip->splines_imag[i]) t->m[i] += I* gsl_spline_eval(ip->splines_imag[i], freq, NULL /*does this work?*/);
349 return QPMS_SUCCESS;
352 qpms_errno_t qpms_tmatrix_generator_interpolator(qpms_tmatrix_t *t,
353 complex double omega, const void *interp)
355 return qpms_tmatrix_interpolator_eval_fill(t,
356 (const qpms_tmatrix_interpolator_t *) interp, omega);
359 double qpms_SU2eV(double e_SU) {
360 return e_SU * SCUFF_OMEGAUNIT / (ELECTRONVOLT / HBAR);
363 double qpms_SU2SI(double e_SU) {
364 return e_SU * SCUFF_OMEGAUNIT;
367 /// TODO doc and more checks
368 qpms_errno_t qpms_read_scuff_tmatrix(
369 FILE *f, ///< file handle
370 const qpms_vswf_set_spec_t * bs, ///< VSWF set spec
371 size_t *const n, ///< Number of successfully loaded t-matrices
372 double* *const freqs, ///< Frequencies in SI units
373 double* *const freqs_su, ///< Frequencies in SCUFF units (optional)
374 qpms_tmatrix_t* *const tmatrices_array, ///< The resulting T-matrices (optional).
375 complex double* *const tmdata
377 if (!(freqs && n && tmdata))
378 qpms_pr_error_at_flf(__FILE__, __LINE__, __func__,
379 "freqs, n, and tmdata are mandatory arguments and must not be NULL.");
380 if(bs->norm & (QPMS_NORMALISATION_REVERSE_AZIMUTHAL_PHASE
381 | QPMS_NORMALISATION_SPHARM_REAL))
382 QPMS_NOT_IMPLEMENTED("Sorry, only standard complex-spherical harmonic based waves are supported right now");
383 int n_alloc = 128; // First chunk to allocate
384 *n = 0;
385 *freqs = malloc(n_alloc * sizeof(double));
386 if (freqs_su) *freqs_su = malloc(n_alloc * sizeof(double));
387 *tmdata = malloc(sizeof(complex double) * bs->n * bs->n * n_alloc);
388 if (!*freqs || (!freqs_su != !*freqs_su) || !*tmdata)
389 qpms_pr_error_at_flf(__FILE__, __LINE__, __func__,
390 "malloc() failed.");
391 size_t linebufsz = 256;
392 char *linebuf = malloc(linebufsz);
393 ssize_t readchars;
394 double lastfreq_su = NAN;
395 while((readchars = getline(&linebuf, &linebufsz, f)) != -1) {
396 if (linebuf[0] == '#') continue;
397 int Alpha, LAlpha, MAlpha, PAlpha, Beta, LBeta, MBeta, PBeta;
398 double currentfreq_su, tr, ti;
399 QPMS_ENSURE(11 == sscanf(linebuf, "%lf %d %d %d %d %d %d %d %d %lf %lf",
400 &currentfreq_su, &Alpha, &LAlpha, &MAlpha, &PAlpha,
401 &Beta, &LBeta, &MBeta, &PBeta, &tr, &ti),
402 "Malformed T-matrix file");
403 if (currentfreq_su != lastfreq_su) { // New frequency -> new T-matrix
404 ++*n;
405 lastfreq_su = currentfreq_su;
406 if(*n > n_alloc) {
407 n_alloc *= 2;
408 QPMS_CRASHING_REALLOC(*freqs, n_alloc * sizeof(double));
409 if (freqs_su) {QPMS_CRASHING_REALLOC(*freqs_su, n_alloc * sizeof(double));}
410 QPMS_CRASHING_REALLOC(*tmdata, sizeof(complex double) * bs->n * bs->n * n_alloc);
412 if (freqs_su) (*freqs_su)[*n-1] = currentfreq_su;
413 (*freqs)[*n-1] = qpms_SU2SI(currentfreq_su);
415 for(size_t i = 0; i < bs->n * bs->n; ++i)
416 (*tmdata)[(*n-1)*bs->n*bs->n + i] = NAN + I*NAN;
418 qpms_vswf_type_t TAlpha, TBeta;
419 switch(PAlpha) {
420 case 0: TAlpha = QPMS_VSWF_MAGNETIC; break;
421 case 1: TAlpha = QPMS_VSWF_ELECTRIC; break;
422 default: QPMS_WTF;
424 switch(PBeta) {
425 case 0: TBeta = QPMS_VSWF_MAGNETIC; break;
426 case 1: TBeta = QPMS_VSWF_ELECTRIC; break;
427 default: QPMS_WTF;
429 qpms_uvswfi_t srcui = qpms_tmn2uvswfi(TAlpha, MAlpha, LAlpha),
430 destui = qpms_tmn2uvswfi(TBeta, MBeta, LBeta);
431 ssize_t srci = qpms_vswf_set_spec_find_uvswfi(bs, srcui),
432 desti = qpms_vswf_set_spec_find_uvswfi(bs, destui);
433 if (srci == -1 || desti == -1)
434 /* This element has not been requested in bs->ilist. */
435 continue;
436 else
437 (*tmdata)[(*n-1)*bs->n*bs->n + desti*bs->n + srci] = (tr + I*ti)
438 * qpms_normalisation_factor_uvswfi(bs->norm, srcui)
439 / qpms_normalisation_factor_uvswfi(bs->norm, destui)
440 * qpms_normalisation_factor_uvswfi(QPMS_NORMALISATION_CONVENTION_SCUFF, destui)
441 / qpms_normalisation_factor_uvswfi(QPMS_NORMALISATION_CONVENTION_SCUFF, srcui);
443 free(linebuf);
444 // free some more memory
445 n_alloc = *n;
446 *freqs = realloc(*freqs, n_alloc * sizeof(double));
447 if (freqs_su) *freqs_su = realloc(*freqs_su, n_alloc * sizeof(double));
448 if (tmatrices_array) *tmatrices_array = realloc(*tmatrices_array, n_alloc * sizeof(qpms_tmatrix_t));
449 *tmdata = realloc(*tmdata, sizeof(complex double) * bs->n * bs->n * n_alloc);
450 if (!*freqs || (!freqs_su != !*freqs_su) || !*tmdata)
451 qpms_pr_error_at_flf(__FILE__, __LINE__, __func__,
452 "realloc() failed.");
453 if (tmatrices_array) {
454 *tmatrices_array = malloc(n_alloc * sizeof(qpms_tmatrix_t));
455 if (!*tmatrices_array)
456 qpms_pr_error_at_flf(__FILE__, __LINE__, __func__,
457 "malloc() failed.");
458 for (size_t i = 0; i < *n; ++i) {
459 (*tmatrices_array)[i].spec = bs;
460 (*tmatrices_array)[i].m = (*tmdata) + i * bs->n * bs->n;
461 (*tmatrices_array)[i].owns_m = false;
464 return QPMS_SUCCESS;
467 bool qpms_load_scuff_tmatrix_crash_on_failure = true;
469 qpms_errno_t qpms_load_scuff_tmatrix(
470 const char *path, ///< file path
471 const qpms_vswf_set_spec_t * bs, ///< VSWF set spec
472 size_t *const n, ///< Number of successfully loaded t-matrices
473 double **const freqs, ///< Frequencies in SI units
474 double ** const freqs_su, ///< Frequencies in SCUFF units (optional)
475 qpms_tmatrix_t ** const tmatrices_array, ///< The resulting T-matrices (optional).
476 complex double ** const tmdata
478 FILE *f = fopen(path, "r");
479 if (!f) {
480 if (qpms_load_scuff_tmatrix_crash_on_failure) {
481 QPMS_PR_ERROR("Could not open T-matrix file %s", path);
482 } else return errno;
484 qpms_errno_t retval =
485 qpms_read_scuff_tmatrix(f, bs, n, freqs, freqs_su, tmatrices_array, tmdata);
487 for (size_t i = 0; i < *n * bs->n * bs->n; ++i)
488 if(isnan(creal((*tmdata)[i])) || isnan(cimag((*tmdata)[i]))) {
489 QPMS_WARN("Encountered NAN in a loaded T-matrix");
490 retval |= QPMS_NAN_ENCOUNTERED;
491 break;
494 if(fclose(f)) { QPMS_PR_ERROR("Could not close the T-matrix file %s "
495 "(well, that's weird, since it's read only).", path); }
497 return retval;
500 complex double *qpms_mie_coefficients_reflection(
501 complex double *target, ///< Target array of length bspec->n. If NULL, a new one will be allocated.
502 const qpms_vswf_set_spec_t *bspec, ///< Defines which of the coefficients are calculated.
503 double a, ///< Radius of the sphere.
504 complex double k_i, ///< Wave number of the internal material of the sphere.
505 complex double k_e, ///< Wave number of the surrounding medium.
506 complex double mu_i, ///< Relative permeability of the sphere material.
507 complex double mu_e, ///< Relative permeability of the surrounding medium.
508 qpms_bessel_t J_ext,
509 qpms_bessel_t J_scat
510 // TODO J_ext, J_scat?
513 * This implementation pretty much copies mie_coefficients() from qpms_p.py, so
514 * any bugs there should affect this function as well and perhaps vice versa.
516 * Compared to mie_coefficients() from qpms_p.py, this has opposite sign.
518 QPMS_ENSURE(J_ext != J_scat, "J_ext and J_scat must not be equal. Perhaps you want J_ext = 1 and J_scat = 3 anyways.");
519 if (!target) QPMS_CRASHING_MALLOC(target, bspec->n * sizeof(complex double));
520 qpms_l_t lMax = bspec->lMax;
521 complex double x_i = k_i * a;
522 complex double x_e = k_e * a;
523 // complex double m = k_i / k_e;
524 complex double eta_inv_i = k_i / mu_i;
525 complex double eta_inv_e = k_e / mu_e;
527 complex double zi[lMax + 2];
528 complex double ze[lMax + 2];
529 complex double zs[lMax + 2];
530 complex double RH[lMax + 1] /* MAGNETIC */, RV[lMax+1] /* ELECTRIC */;
531 QPMS_ENSURE_SUCCESS(qpms_sph_bessel_fill(QPMS_BESSEL_REGULAR, lMax+1, x_i, zi));
532 QPMS_ENSURE_SUCCESS(qpms_sph_bessel_fill(J_ext, lMax+1, x_e, ze));
533 QPMS_ENSURE_SUCCESS(qpms_sph_bessel_fill(J_scat, lMax+1, x_e, zs));
534 for (qpms_l_t l = 0; l <= lMax; ++l) {
535 // Bessel function derivatives as in DLMF 10.51.2
536 complex double dzi = -zi[l+1] + l/x_i*zi[l];
537 complex double dze = -ze[l+1] + l/x_e*ze[l];
538 complex double dzs = -zs[l+1] + l/x_e*zs[l];
539 complex double fi = zi[l]/x_i+dzi;
540 complex double fs = zs[l]/x_e+dzs;
541 complex double fe = ze[l]/x_e+dze;
542 RV[l] = (-eta_inv_i * fe * zi[l] + eta_inv_e * ze[l] * fi)/(-eta_inv_e * fi * zs[l] + eta_inv_i * zi[l] * fs);
543 RH[l] = (-eta_inv_e * fe * zi[l] + eta_inv_i * ze[l] * fi)/(-eta_inv_i * fi * zs[l] + eta_inv_e * zi[l] * fs);
546 for (size_t i = 0; i < bspec->n; ++i) {
547 qpms_l_t l; qpms_m_t m; qpms_vswf_type_t t;
548 QPMS_ENSURE_SUCCESS(qpms_uvswfi2tmn(bspec->ilist[i], &t, &m, &l));
549 assert(l <= lMax);
550 switch(t) {
551 case QPMS_VSWF_ELECTRIC:
552 target[i] = RV[l];
553 break;
554 case QPMS_VSWF_MAGNETIC:
555 target[i] = RH[l];
556 break;
557 default: QPMS_WTF;
560 return target;
563 /// Replaces the contents of an existing T-matrix with that of a spherical nanoparticle calculated using the Lorentz-mie theory.
564 qpms_errno_t qpms_tmatrix_spherical_fill(qpms_tmatrix_t *t, ///< T-matrix whose contents are to be replaced. Not NULL.
565 double a, ///< Radius of the sphere.
566 complex double k_i, ///< Wave number of the internal material of the sphere.
567 complex double k_e, ///< Wave number of the surrounding medium.
568 complex double mu_i, ///< Relative permeability of the sphere material.
569 complex double mu_e ///< Relative permeability of the surrounding medium.
571 complex double *miecoeffs = qpms_mie_coefficients_reflection(NULL, t->spec, a, k_i, k_e, mu_i, mu_e,
572 QPMS_BESSEL_REGULAR, QPMS_HANKEL_PLUS);
573 memset(t->m, 0, SQ(t->spec->n)*sizeof(complex double));
574 for(size_t i = 0; i < t->spec->n; ++i)
575 t->m[t->spec->n*i + i] = miecoeffs[i];
576 free(miecoeffs);
577 return QPMS_SUCCESS;
580 qpms_errno_t qpms_tmatrix_generator_sphere(qpms_tmatrix_t *t,
581 complex double omega, const void *params)
583 const qpms_tmatrix_generator_sphere_param_t *p = params;
584 qpms_epsmu_t out = p->outside.function(omega, p->outside.params);
585 qpms_epsmu_t in = p->inside.function(omega, p->inside.params);
586 return qpms_tmatrix_spherical_fill(t,
587 p->radius,
588 qpms_wavenumber(omega, in),
589 qpms_wavenumber(omega, out),
590 in.mu, out.mu);
593 /// Convenience function to calculate T-matrix of a non-magnetic spherical \
594 particle using the permittivity values, replacing existing T-matrix data.
595 qpms_errno_t qpms_tmatrix_spherical_mu0_fill(
596 qpms_tmatrix_t *t, ///< T-matrix whose contents are to be replaced. Not NULL.
597 double a, ///< Radius of the sphere.
598 double omega, ///< Angular frequency.
599 complex double epsilon_fg, ///< Permittivity of the sphere.
600 complex double epsilon_bg ///< Permittivity of the background medium.
603 complex double k_i = csqrt(epsilon_fg) * omega / SPEED_OF_LIGHT;
604 complex double k_e = csqrt(epsilon_bg) * omega / SPEED_OF_LIGHT;
605 return qpms_tmatrix_spherical_fill(t, a, k_i, k_e, 1, 1);
608 complex double *qpms_apply_tmatrix(
609 complex double *f, const complex double *a, const qpms_tmatrix_t *T)
611 const size_t n = T->spec->n;
612 if(!f)
613 QPMS_CRASHING_CALLOC(f, n, sizeof(complex double));
614 const complex double one = 1;
615 const complex double zero = 0;
616 cblas_zgemv(CblasRowMajor, CblasNoTrans, n, n, &one, T->m, n, a, 1, &zero, f, 1);
617 return f;
620 qpms_arc_function_retval_t qpms_arc_sphere(double theta, const void *R) {
621 qpms_arc_function_retval_t retval = {*(const double*)R, 0};
622 return retval;
625 qpms_arc_function_retval_t qpms_arc_cylinder(double theta, const void *param) {
626 const qpms_arc_cylinder_params_t *p = param;
627 double thresh = atan(2 * p->R / p->h);
628 QPMS_ENSURE(theta >= 0 && theta <= M_PI,
629 "theta = %g, but it must lie in interval [0, M_PI]", theta);
630 qpms_arc_function_retval_t res;
631 if (theta < thresh) {
632 // Upper base
633 res.r = 0.5 * p->h / cos(theta);
634 res.beta = -theta;
635 } else if (theta <= M_PI - thresh) {
636 // Side
637 res.r = p->R / cos(theta - M_PI_2);
638 res.beta = -theta + M_PI_2;
639 } else {
640 // Lower base
641 res.r = 0.5 * p->h / cos(theta - M_PI);
642 res.beta = -theta + M_PI;
644 return res;
648 struct tmatrix_axialsym_integral_param_t {
649 const qpms_vswf_set_spec_t *bspec;
650 qpms_l_t l, l_in;
651 qpms_m_t m; // m_in = -m
652 qpms_vswf_type_t t, t_in;
653 qpms_arc_function_t f;
654 complex double k_in, k, z_in, z;
655 bool realpart; // Otherwise imaginary part
656 qpms_bessel_t btype; // For Q QPMS_HANKEL_PLUS, for R QPMS_BESSEL_REGULAR
659 static double tmatrix_axialsym_integrand(double theta, void *param) {
660 // Pretty inefficient; either real or imaginary part is thrown away etc.
661 struct tmatrix_axialsym_integral_param_t *p = param;
662 qpms_arc_function_retval_t rb = p->f.function(theta, p->f.params);
663 csph_t kr = {rb.r * p->k, theta, 0},
664 kr_in = {rb.r * p->k_in, theta, 0};
666 csphvec_t y_el = qpms_vswf_single_el_csph(p->m, p->l, kr, p->btype, p->bspec->norm);
667 csphvec_t y_mg = qpms_vswf_single_mg_csph(p->m, p->l, kr, p->btype, p->bspec->norm);
668 csphvec_t v_in_el = qpms_vswf_single_el_csph(-p->m, p->l_in, kr_in, QPMS_BESSEL_REGULAR, p->bspec->norm);
669 csphvec_t v_in_mg = qpms_vswf_single_mg_csph(-p->m, p->l_in, kr_in, QPMS_BESSEL_REGULAR, p->bspec->norm);
670 csphvec_t y1, y2, v_in1, v_in2;
671 switch(p->t) {
672 case QPMS_VSWF_ELECTRIC: y1 = y_el; y2 = y_mg; break;
673 case QPMS_VSWF_MAGNETIC: y1 = y_mg; y2 = y_el; break;
674 default: QPMS_WTF;
676 switch(p->t_in) {
677 case QPMS_VSWF_ELECTRIC: v_in1 = v_in_mg; v_in2 = v_in_el; break;
678 case QPMS_VSWF_MAGNETIC: v_in1 = v_in_el; v_in2 = v_in_mg; break;
679 default: QPMS_WTF;
681 // Normal vector components
682 double nrc = cos(rb.beta), nthetac = sin(rb.beta);
683 // First triple product
684 complex double tp1 = nrc * (y1.thetac * v_in1.phic - y1.phic * v_in1.thetac)
685 + nthetac * (y1.phic * v_in1.rc - y1.rc * v_in1.phic);
686 // Second thiple product
687 complex double tp2 = nrc * (y2.thetac * v_in2.phic - y2.phic * v_in2.thetac)
688 + nthetac * (y2.phic * v_in2.rc - y2.rc * v_in2.phic);
689 double jac = SQ(rb.r) * sin(theta) / nrc; // Jacobian
690 complex double res = jac * (p->z/p->z_in * tp1 + tp2);
691 return p->realpart ? creal(res) : cimag(res);
696 long qpms_axsym_tmatrix_integration_nthreads_default = 4;
697 long qpms_axsym_tmatrix_integration_nthreads_override = 0;
699 void qpms_axsym_tmatrix_integration_set_nthreads(long n) {
700 qpms_axsym_tmatrix_integration_nthreads_override = n;
703 struct qpms_tmatrix_axialsym_fill_integration_thread_arg {
704 complex double *Q, *R;
705 const qpms_vswf_set_spec_t *bspecQR;
706 double epsrel, epsabs;
707 size_t intlimit;
708 qpms_arc_function_t f;
709 const struct tmatrix_axialsym_integral_param_t *p_template; // Thread must copy this and set its own realpart and btype fields.
710 size_t *i2start_ptr;
711 pthread_mutex_t *i2start_mutex;
714 // This wrapper retries to calculate the T-matrix integral and if it fails, reduces torelance.
715 static inline int axint_wrapper(const gsl_function *f, double epsabs, double epsrel,
716 size_t intlimit, gsl_integration_workspace *w, double *result, double *abserr) {
717 static int tol_exceeded = 0;
718 double errfac;
719 int retval;
720 for (errfac = 1; epsabs*errfac < 0.01 && epsrel*errfac < 0.01; errfac *= 8) {
721 retval = gsl_integration_qag(f, 0, M_PI, epsabs*errfac, epsrel*errfac,
722 intlimit, 2, w, result, abserr);
723 if(QPMS_LIKELY(!retval)) {
724 if(QPMS_LIKELY(errfac == 1))
725 return 0;
726 else {
727 QPMS_DEBUG(QPMS_DBGMSG_INTEGRATION,
728 "T-matrix quadrature succeeded only after %g-fold"
729 "tolerance reduction."); // TODO more details here?
730 return GSL_EROUND;
733 else if(retval == GSL_EROUND) {
734 if (!tol_exceeded) {
735 QPMS_WARN("T-matrix quadrature failed; retrying with worse tolerances."
736 " (this warning will be shown only once).");
738 tol_exceeded += 1;
740 else QPMS_ENSURE_SUCCESS(retval);
742 QPMS_PR_ERROR("T-matrix quadrature failed even after harsh tolerance reduction.");
746 static void * qpms_tmatrix_axialsym_fill_integration_thread(void *arg) {
747 const struct qpms_tmatrix_axialsym_fill_integration_thread_arg *a = arg;
748 struct tmatrix_axialsym_integral_param_t p = *a->p_template;
749 const gsl_function f = {tmatrix_axialsym_integrand, (void *) &p};
750 gsl_integration_workspace *w = gsl_integration_workspace_alloc(a->intlimit);
752 while(1) {
753 QPMS_ENSURE_SUCCESS(pthread_mutex_lock(a->i2start_mutex));
754 if(*(a->i2start_ptr) >= a->bspecQR->n) {// Everything already done, leave thread
755 QPMS_ENSURE_SUCCESS(pthread_mutex_unlock(a->i2start_mutex));
756 break;
758 const size_t i2 = *(a->i2start_ptr);
759 ++*(a->i2start_ptr);
760 QPMS_ENSURE_SUCCESS(pthread_mutex_unlock(a->i2start_mutex));
761 for(size_t i1 = 0; i1 < a->bspecQR->n; ++i1) {
762 // NOTE that the Q', R' matrices are !TRANSPOSED! here because of zgetrs
763 const size_t iQR = i1 + i2 * a->bspecQR->n;
764 qpms_m_t m1, m2;
765 qpms_l_t l1, l2;
766 qpms_vswf_type_t t1, t2;
767 const qpms_uvswfi_t u1 = a->bspecQR->ilist[i1], u2 = a->bspecQR->ilist[i2];
768 QPMS_ENSURE_SUCCESS(qpms_uvswfi2tmn(u1, &t1, &m1, &l1));
769 QPMS_ENSURE_SUCCESS(qpms_uvswfi2tmn(u2, &t2, &m2, &l2));
770 if (m1 + m2) {
771 a->Q[iQR] = 0;
772 a->R[iQR] = 0;
773 } else {
774 p.m = m1;
775 p.l = l1; p.l_in = l2;
776 p.t = t1; p.t_in = t2;
778 double result;
779 double abserr; // gsl_integration_qag() needs this; thrown away for now
780 // Re(R')
781 p.btype = QPMS_BESSEL_REGULAR;
782 p.realpart = true;
783 axint_wrapper(&f, a->epsabs, a->epsrel, a->intlimit, w, &result, &abserr);
784 a->R[iQR] = result;
786 // Im(R')
787 p.realpart = false;
788 axint_wrapper(&f, a->epsabs, a->epsrel, a->intlimit, w, &result, &abserr);
789 a->R[iQR] += I*result;
791 // Re(Q')
792 p.btype = QPMS_HANKEL_PLUS;
793 p.realpart = true;
794 axint_wrapper(&f, a->epsabs, a->epsrel, a->intlimit, w, &result, &abserr);
795 a->Q[iQR] = result;
797 // Im(Q')
798 p.realpart = false;
799 axint_wrapper(&f, a->epsabs, a->epsrel, a->intlimit, w, &result, &abserr);
800 a->Q[iQR] += I*result;
804 gsl_integration_workspace_free(w);
807 qpms_errno_t qpms_tmatrix_axialsym_fill(
808 qpms_tmatrix_t *t, complex double omega, qpms_epsmu_t outside,
809 qpms_epsmu_t inside,qpms_arc_function_t shape, qpms_l_t lMaxQR)
811 QPMS_UNTESTED; // TODO also check whether this is valid for all norms.
812 const qpms_vswf_set_spec_t *bspec = t->spec;
813 struct tmatrix_axialsym_integral_param_t p;
814 p.k = qpms_wavenumber(omega, outside);
815 p.k_in = qpms_wavenumber(omega, inside);
816 p.z = qpms_waveimpedance(outside);
817 p.z_in = qpms_waveimpedance(inside);
818 p.f = shape;
819 p.bspec = bspec;
821 if (lMaxQR < bspec->lMax) lMaxQR = bspec->lMax;
822 qpms_vswf_set_spec_t *bspecQR = qpms_vswf_set_spec_from_lMax(lMaxQR, bspec->norm);
823 size_t *reindex = qpms_vswf_set_reindex(bspec, bspecQR);
825 // Q' and R' matrices
826 complex double *Q, *R;
827 QPMS_CRASHING_MALLOC(Q, sizeof(complex double) * SQ(bspecQR->n));
828 QPMS_CRASHING_MALLOC(R, sizeof(complex double) * SQ(bspecQR->n));
830 // Not sure what the size should be, but this should be more than enough.
831 const size_t intlimit = 1024;
832 const double epsrel = TMATRIX_AXIALSYM_INTEGRAL_EPSREL;
833 const double epsabs = TMATRIX_AXIALSYM_INTEGRAL_EPSABS;
835 //==== multi-threaded integration begin
836 size_t i2start = 0;
837 pthread_mutex_t i2start_mutex;
838 QPMS_ENSURE_SUCCESS(pthread_mutex_init(&i2start_mutex, NULL));
839 const struct qpms_tmatrix_axialsym_fill_integration_thread_arg
840 arg = {Q, R, bspecQR, epsrel, epsabs, intlimit, shape, &p, &i2start, &i2start_mutex};
842 // FIXME THIS IS NOT PORTABLE (_SC_NPROCESSORS_ONLN is not a standard POSIX thing)
843 long nthreads;
844 if (qpms_axsym_tmatrix_integration_nthreads_override > 0) {
845 nthreads = qpms_axsym_tmatrix_integration_nthreads_override;
846 QPMS_DEBUG(QPMS_DBGMSG_THREADS, "Using overriding value of %ld thread(s).",
847 nthreads);
848 } else {
849 nthreads = get_ncpus();
850 if (nthreads < 1) {
851 QPMS_DEBUG(QPMS_DBGMSG_THREADS, "_SC_NPROCESSORS_ONLN returned %ld, using %ld thread(s) instead.",
852 nthreads, qpms_axsym_tmatrix_integration_nthreads_default);
853 nthreads = qpms_axsym_tmatrix_integration_nthreads_default;
854 } else {
855 QPMS_DEBUG(QPMS_DBGMSG_THREADS, "_SC_NRPOCESSORS_ONLN returned %ld.", nthreads);
858 nthreads = MIN(nthreads, bspec->n); // don't make redundant threads
859 pthread_t thread_ids[nthreads];
860 for(long thi = 0; thi < nthreads; ++thi)
861 QPMS_ENSURE_SUCCESS(pthread_create(thread_ids + thi, NULL,
862 qpms_tmatrix_axialsym_fill_integration_thread,
863 (void *) &arg));
864 for(long thi = 0; thi < nthreads; ++thi) {
865 void *retval;
866 QPMS_ENSURE_SUCCESS(pthread_join(thread_ids[thi], &retval));
868 QPMS_ENSURE_SUCCESS(pthread_mutex_destroy(&i2start_mutex));
869 //===== multi-threaded integration end
872 // Compute T-matrix; maybe it would be better solved with some sparse matrix mechanism,
873 // but fukkit.
874 const size_t n = bspecQR->n;
875 lapack_int *pivot;
876 QPMS_CRASHING_MALLOC(pivot, sizeof(lapack_int) * n);
877 QPMS_ENSURE_SUCCESS(LAPACKE_zgetrf(LAPACK_ROW_MAJOR, n, n, Q, n, pivot));
878 // Solve Q'^T X = R'^T where X will be -T^T
879 // Note that Q'^T, R'^T are already (transposed) in Q, R.
880 QPMS_ENSURE_SUCCESS(LAPACKE_zgetrs(LAPACK_ROW_MAJOR, 'N', n, n /*nrhs*/,
881 Q, n, pivot, R, n));
882 // R now contains -T^T.
883 for(size_t i1 = 0; i1 < bspec->n; ++i1)
884 for(size_t i2 = 0; i2 < bspec->n; ++i2) {
885 if (reindex[i1] == ~(size_t) 0 && reindex[i2] == ~(size_t) 0) QPMS_WTF;
886 const size_t it = i1 * bspec->n + i2;
887 const size_t iQR = reindex[i1] + reindex[i2] * bspecQR->n;
888 t->m[it] = -R[iQR];
891 free(pivot); free(R); free(Q); free(reindex);
892 qpms_vswf_set_spec_free(bspecQR);
893 return QPMS_SUCCESS;
896 qpms_errno_t qpms_tmatrix_axialsym_RQ_transposed_fill(complex double *target,
897 complex double omega, qpms_epsmu_t outside, qpms_epsmu_t inside,
898 qpms_arc_function_t shape, qpms_l_t lMaxQR, qpms_normalisation_t norm,qpms_bessel_t J)
900 QPMS_UNTESTED; // TODO also check whether this is valid for all norms.
901 struct tmatrix_axialsym_integral_param_t p;
903 qpms_vswf_set_spec_t *bspecQR = qpms_vswf_set_spec_from_lMax(lMaxQR, norm);
905 p.k = qpms_wavenumber(omega, outside);
906 p.k_in = qpms_wavenumber(omega, inside);
907 p.z = qpms_waveimpedance(outside);
908 p.z_in = qpms_waveimpedance(inside);
909 p.f = shape;
910 p.bspec = bspecQR;
911 p.btype = J;
912 const gsl_function f = {tmatrix_axialsym_integrand, (void *) &p};
915 // Not sure what the size should be, but this should be more than enough.
916 const size_t intlimit = 1024;
917 const double epsrel = TMATRIX_AXIALSYM_INTEGRAL_EPSREL;
918 const double epsabs = TMATRIX_AXIALSYM_INTEGRAL_EPSABS;
919 gsl_integration_workspace *w = gsl_integration_workspace_alloc(intlimit);
920 for(size_t i1 = 0; i1 < bspecQR->n; ++i1)
921 for(size_t i2 = 0; i2 < bspecQR->n; ++i2) {
922 // NOTE that the Q', R' matrices are !TRANSPOSED! here because of zgetrs
923 const size_t iQR = i1 + i2 * bspecQR->n;
924 qpms_m_t m1, m2;
925 qpms_l_t l1, l2;
926 qpms_vswf_type_t t1, t2;
927 const qpms_uvswfi_t u1 = bspecQR->ilist[i1], u2 = bspecQR->ilist[i2];
928 QPMS_ENSURE_SUCCESS(qpms_uvswfi2tmn(u1, &t1, &m1, &l1));
929 QPMS_ENSURE_SUCCESS(qpms_uvswfi2tmn(u2, &t2, &m2, &l2));
930 if (m1 + m2) {
931 target[iQR] = 0;
932 } else {
933 p.m = m1;
934 p.l = l1; p.l_in = l2;
935 p.t = t1; p.t_in = t2;
937 double result;
938 double abserr; // gsl_integration_qag() needs this; thrown away for now
939 // Re(R' or Q')
940 p.realpart = true;
941 QPMS_ENSURE_SUCCESS(gsl_integration_qag(&f, 0, M_PI, epsabs, epsrel,
942 intlimit, 2, w, &result, &abserr));
943 target[iQR] = result;
945 // Im(R' or Q')
946 p.realpart = false;
947 QPMS_ENSURE_SUCCESS(gsl_integration_qag(&f, 0, M_PI, epsabs, epsrel,
948 intlimit, 2, w, &result, &abserr));
949 target[iQR] += I*result;
952 gsl_integration_workspace_free(w);
954 qpms_vswf_set_spec_free(bspecQR);
955 return QPMS_SUCCESS;
958 qpms_errno_t qpms_tmatrix_generator_axialsym_RQ_transposed_fill(complex double *target,
959 complex double omega, const qpms_tmatrix_generator_axialsym_param_t *param,
960 qpms_normalisation_t norm, qpms_bessel_t J)
962 const qpms_tmatrix_generator_axialsym_param_t *p = param;
963 return qpms_tmatrix_axialsym_RQ_transposed_fill(target, omega,
964 p->outside.function(omega, p->outside.params),
965 p->inside.function(omega, p->inside.params),
966 p->shape,
967 p->lMax_extend, norm, J);
970 qpms_errno_t qpms_tmatrix_generator_axialsym(qpms_tmatrix_t *t, complex double omega, const void *param)
972 const qpms_tmatrix_generator_axialsym_param_t *p = param;
973 return qpms_tmatrix_axialsym_fill(t, omega,
974 p->outside.function(omega, p->outside.params),
975 p->inside.function(omega, p->inside.params),
976 p->shape,
977 p->lMax_extend);
980 qpms_errno_t qpms_tmatrix_generator_constant(qpms_tmatrix_t *t,
981 complex double omega, const void *tmorig) {
982 QPMS_UNTESTED;
983 // TODO more thorough check on bspec!
984 const qpms_tmatrix_t *orig = tmorig;
985 const size_t tocopy = SQ(MIN(t->spec->n, orig->spec->n));
986 memcpy(t->m, orig->m, tocopy*sizeof(complex double));
987 memset(t->m+tocopy, 0, (SQ(t->spec->n)-tocopy)*sizeof(complex double));
988 return QPMS_SUCCESS;
991 void qpms_tmatrix_operation_clear(qpms_tmatrix_operation_t *f) {
992 switch(f->typ) {
993 case QPMS_TMATRIX_OPERATION_NOOP:
994 break;
995 case QPMS_TMATRIX_OPERATION_LRMATRIX:
996 if(f->op.lrmatrix.owns_m)
997 free(f->op.lrmatrix.m);
998 break;
999 case QPMS_TMATRIX_OPERATION_SCMULZ:
1000 if(f->op.scmulz.owns_m)
1001 free(f->op.scmulz.m);
1002 break;
1003 case QPMS_TMATRIX_OPERATION_IROT3:
1004 break;
1005 case QPMS_TMATRIX_OPERATION_IROT3ARR:
1006 if(f->op.irot3arr.owns_ops)
1007 free(f->op.irot3arr.ops);
1008 break;
1009 case QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE: // Group never owned
1010 break;
1011 case QPMS_TMATRIX_OPERATION_COMPOSE_SUM:
1013 struct qpms_tmatrix_operation_compose_sum * const o =
1014 &(f->op.compose_sum);
1015 if(o->opmem) {
1016 for(size_t i = 0; i < o->n; ++i)
1017 qpms_tmatrix_operation_clear(&(o->opmem[i]));
1018 free(o->opmem);
1020 free(o->ops);
1022 break;
1023 case QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN:
1025 struct qpms_tmatrix_operation_compose_chain * const o =
1026 &(f->op.compose_chain);
1027 if(o->opmem) {
1028 for(size_t i = 0; i < o->n; ++i)
1029 if(o->ops_owned[i]) {
1030 // A bit complicated yet imperfect sanity/bound check,
1031 // but this way we at least get rid of the const discard warning.
1032 ptrdiff_t d = o->ops[i] - o->opmem;
1033 QPMS_ENSURE(d >= 0 && d < o->opmem_size,
1034 "Bound check failed; is there a bug in the"
1035 " construction of compose_chain operation?\n"
1036 "opmem_size = %zu, o->ops[%zu] - o->opmem = %td.",
1037 o->opmem_size, i, d);
1038 qpms_tmatrix_operation_clear(o->opmem + d);
1040 free(o->opmem);
1041 free(o->ops_owned);
1043 free(o->ops);
1045 break;
1046 default:
1047 QPMS_WTF;
1051 void qpms_tmatrix_operation_copy(qpms_tmatrix_operation_t *dest, const qpms_tmatrix_operation_t *src) {
1052 memcpy(dest, src, sizeof(qpms_tmatrix_operation_t));
1053 switch(dest->typ) {
1054 case QPMS_TMATRIX_OPERATION_NOOP:
1055 break;
1056 case QPMS_TMATRIX_OPERATION_LRMATRIX:
1057 QPMS_CRASHING_MALLOCPY(dest->op.lrmatrix.m, src->op.lrmatrix.m,
1058 sizeof(complex double) * src->op.lrmatrix.m_size);
1059 dest->op.lrmatrix.owns_m = true;
1060 break;
1061 case QPMS_TMATRIX_OPERATION_SCMULZ:
1062 QPMS_CRASHING_MALLOCPY(dest->op.scmulz.m, src->op.scmulz.m,
1063 sizeof(complex double) * src->op.scmulz.m_size);
1064 dest->op.scmulz.owns_m = true;
1065 break;
1066 case QPMS_TMATRIX_OPERATION_IROT3:
1067 break;
1068 case QPMS_TMATRIX_OPERATION_IROT3ARR:
1069 QPMS_CRASHING_MALLOCPY(dest->op.irot3arr.ops, src->op.irot3arr.ops,
1070 src->op.irot3arr.n * sizeof(qpms_irot3_t));
1071 dest->op.irot3arr.owns_ops = true;
1072 break;
1073 case QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE: // Group never owned
1074 break;
1075 case QPMS_TMATRIX_OPERATION_COMPOSE_SUM:
1077 struct qpms_tmatrix_operation_compose_sum * const o =
1078 &(dest->op.compose_sum);
1079 QPMS_CRASHING_MALLOC(o->ops, o->n * sizeof(*(o->ops)));
1080 QPMS_CRASHING_MALLOC(o->opmem, o->n * sizeof(*(o->opmem)));
1081 for(size_t i = 0; i < o->n; ++i) {
1082 qpms_tmatrix_operation_copy(o->opmem + i, src->op.compose_sum.ops[i]);
1083 o->ops[i] = o->opmem + i;
1086 break;
1087 case QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN:
1089 struct qpms_tmatrix_operation_compose_chain * const o =
1090 &(dest->op.compose_chain);
1091 QPMS_CRASHING_MALLOC(o->ops, o->n * sizeof(*(o->ops)));
1092 QPMS_CRASHING_MALLOC(o->ops_owned, o->n * sizeof(*(o->ops_owned)));
1093 QPMS_CRASHING_MALLOC(o->opmem, o->n * sizeof(*(o->opmem)));
1094 for(size_t i = 0; i < o->n; ++i) {
1095 qpms_tmatrix_operation_copy(o->opmem + i, src->op.compose_chain.ops[i]);
1096 o->ops[i] = o->opmem + i;
1097 o->ops_owned[i] = true;
1100 break;
1101 default:
1102 QPMS_WTF;
1104 dest->typ = src->typ;
1107 void qpms_tmatrix_operation_compose_chain_init(qpms_tmatrix_operation_t *dest,
1108 size_t nops, size_t opmem_size) {
1109 if (nops == 0) {
1110 QPMS_WARN("Tried to create a composed (chain) operation of zero size;"
1111 "setting to no-op instead.");
1112 *dest = qpms_tmatrix_operation_noop;
1114 if (nops < opmem_size)
1115 QPMS_WARN("Allocating buffer for %zu operations, in a chained operation of"
1116 " only %zu elemens, that does not seem to make sense.", opmem_size, nops);
1117 dest->typ = QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN;
1118 struct qpms_tmatrix_operation_compose_chain * const o =
1119 &(dest->op.compose_chain);
1120 o->n = nops;
1121 QPMS_CRASHING_MALLOC(o->ops, nops * sizeof(*(o->ops)));
1122 if (opmem_size != 0) {
1123 QPMS_CRASHING_CALLOC(o->ops_owned, o->n, sizeof(_Bool));
1124 QPMS_CRASHING_MALLOC(o->opmem, opmem_size * sizeof(*(o->opmem)));
1125 } else {
1126 o->ops_owned = NULL;
1127 o->opmem = NULL;
1129 o->opmem_size = opmem_size;
1132 qpms_tmatrix_t *qpms_tmatrix_mv(qpms_tmatrix_t *dest,
1133 const qpms_tmatrix_t *orig) {
1134 QPMS_ENSURE(qpms_vswf_set_spec_isidentical(dest->spec, orig->spec),
1135 "Basis specifications must be identical!");
1136 memcpy(dest->m, orig->m, SQ(orig->spec->n)*sizeof(complex double));
1137 return dest;
1140 qpms_tmatrix_t *qpms_tmatrix_apply_operation_replace(qpms_tmatrix_t *dest,
1141 const qpms_tmatrix_operation_t *op, const qpms_tmatrix_t *orig) {
1142 //QPMS_ENSURE(qpms_vswf_set_spec_isidentical(dest->spec, orig->spec),
1143 // "Basis specifications must be identical!");
1144 // TODO should I check also for dest->owns_m?
1145 // FIXME this is highly inoptimal; the hierarchy should be such
1146 // that _operation() and operation_inplace() call this, and not the
1147 // other way around
1148 qpms_tmatrix_mv(dest, orig);
1149 return qpms_tmatrix_apply_operation_inplace(op, dest);
1152 qpms_tmatrix_t *qpms_tmatrix_apply_operation(
1153 const qpms_tmatrix_operation_t *f, const qpms_tmatrix_t *orig) {
1154 // Certain operations could be optimized, but the effect would be marginal.
1155 qpms_tmatrix_t *res = qpms_tmatrix_copy(orig);
1156 return qpms_tmatrix_apply_operation_inplace(f, res);
1159 static qpms_tmatrix_t *qtao_compose_sum_inplace(qpms_tmatrix_t *T,
1160 const struct qpms_tmatrix_operation_compose_sum *cs) {
1161 qpms_tmatrix_t *tmp_target = qpms_tmatrix_init(T->spec);
1162 qpms_tmatrix_t *sum = qpms_tmatrix_init(T->spec);
1163 for (size_t i = 0; i < cs->n; ++i) {
1164 memcpy(tmp_target->m, T->m, SQ(T->spec->n) * sizeof(complex double));
1165 QPMS_ENSURE(qpms_tmatrix_apply_operation_inplace(cs->ops[i] , tmp_target),
1166 "Got NULL pointer from qpms_tmatrix_apply_operation_inplace, hupsis!");
1167 for (size_t j = 0; j < SQ(T->spec->n); ++j)
1168 sum->m[j] += tmp_target->m[j];
1170 for(size_t j = 0; j < SQ(T->spec->n); ++j)
1171 T->m[j] = sum->m[j] * cs->factor;
1172 qpms_tmatrix_free(sum);
1173 qpms_tmatrix_free(tmp_target);
1174 return T;
1177 static qpms_tmatrix_t *qtao_compose_chain_inplace(qpms_tmatrix_t *T,
1178 const struct qpms_tmatrix_operation_compose_chain *cc) {
1179 for(size_t i = 0; i < cc->n; ++i)
1180 qpms_tmatrix_apply_operation_inplace(cc->ops[i], T);
1181 return T;
1184 static qpms_tmatrix_t *qtao_scmulz_inplace(qpms_tmatrix_t *T,
1185 const struct qpms_tmatrix_operation_scmulz *s) {
1186 for(size_t i = 0; i < SQ(T->spec->n); ++i)
1187 T->m[i] *= s->m[i];
1188 return T;
1191 qpms_tmatrix_t *qpms_tmatrix_apply_operation_inplace(
1192 const qpms_tmatrix_operation_t *f, qpms_tmatrix_t *T) {
1193 switch(f->typ) {
1194 case QPMS_TMATRIX_OPERATION_NOOP:
1195 return T;
1196 case QPMS_TMATRIX_OPERATION_LRMATRIX:
1197 return qpms_tmatrix_apply_symop_inplace(T, f->op.lrmatrix.m);
1198 case QPMS_TMATRIX_OPERATION_IROT3:
1199 return qpms_tmatrix_symmetrise_irot3arr_inplace(T, 1, &(f->op.irot3));
1200 case QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE:
1201 return qpms_tmatrix_symmetrise_finite_group_inplace(T, f->op.finite_group);
1202 case QPMS_TMATRIX_OPERATION_COMPOSE_SUM:
1203 return qtao_compose_sum_inplace(T, &(f->op.compose_sum));
1204 case QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN:
1205 return qtao_compose_chain_inplace(T, &(f->op.compose_chain));
1206 case QPMS_TMATRIX_OPERATION_SCMULZ:
1207 return qtao_scmulz_inplace(T, &(f->op.scmulz));
1208 default:
1209 QPMS_WTF;