Fix saving lists of arrays with recent versions of numpy
[qpms.git] / qpms / tmatrices.h
blobf9e3fd8c6f17dc701ff35d9e0b599f03a6c14bfc
1 /*! \file tmatrices.h
2 * \brief T-matrices for scattering systems.
3 */
4 #ifndef TMATRICES_H
5 #define TMATRICES_H
6 // #include "qpms_types.h" // included via materials.h
7 // #include <gsl/gsl_spline.h> // included via materials.h
8 #include "materials.h"
9 #include <stdio.h>
13 struct qpms_finite_group_t;
14 typedef struct qpms_finite_group_t qpms_finite_group_t;
16 /// Returns a pointer to the beginning of the T-matrix row \a rowno.
17 static inline complex double *qpms_tmatrix_row(qpms_tmatrix_t *t, size_t rowno){
18 return t->m + (t->spec->n * rowno);
21 /// Initialises a zero T-matrix.
22 /** \sa qpms_tmatrix_init_from_generator()
23 * \sa qpms_tmatrix_init_from_function() */
24 qpms_tmatrix_t *qpms_tmatrix_init(const qpms_vswf_set_spec_t *bspec);
27 /// Copies a T-matrix, allocating new array for the T-matrix data.
28 qpms_tmatrix_t *qpms_tmatrix_copy(const qpms_tmatrix_t *T);
30 /// Copies a T-matrix contents between two pre-allocated T-matrices.
31 /** orig->spec and dest->spec must be identical.
33 * \returns \a dest
34 */
35 qpms_tmatrix_t *qpms_tmatrix_mv(qpms_tmatrix_t *dest, const qpms_tmatrix_t *orig);
37 /// Destroys a T-matrix.
38 void qpms_tmatrix_free(qpms_tmatrix_t *t);
40 /// Check T-matrix equality/similarity.
41 /**
42 * This function actually checks for identical vswf specs.
43 * TODO define constants with "default" atol, rtol for this function.
45 bool qpms_tmatrix_isclose(const qpms_tmatrix_t *T1, const qpms_tmatrix_t *T2,
46 const double rtol, const double atol);
48 /// Creates a T-matrix from another matrix and a symmetry operation.
49 /** The symmetry operation is expected to be a unitary (square)
50 * matrix \a M with the same
51 * VSWF basis spec as the T-matrix (i.e. \a t->spec). The new T-matrix will then
52 * correspond to CHECKME \f[ T' = MTM^\dagger \f]
54 qpms_tmatrix_t *qpms_tmatrix_apply_symop(
55 const qpms_tmatrix_t *T, ///< the original T-matrix
56 const complex double *M ///< the symmetry op matrix in row-major format
59 /// Applies a symmetry operation onto a T-matrix, rewriting the original T-matrix data.
60 /** The symmetry operation is expected to be a unitary (square)
61 * matrix \a M with the same
62 * VSWF basis spec as the T-matrix (i.e. \a t->spec). The new T-matrix will then
63 * correspond to CHECKME \f[ T' = MTM^\dagger \f]
65 qpms_tmatrix_t *qpms_tmatrix_apply_symop_inplace(
66 qpms_tmatrix_t *T, ///< the original T-matrix
67 const complex double *M ///< the symmetry op matrix in row-major format
70 /// Symmetrizes a T-matrix with an involution symmetry operation.
71 /** The symmetry operation is expected to be a unitary (square)
72 * matrix \a M with the same
73 * VSWF basis spec as the T-matrix (i.e. \a t->spec). The new T-matrix will then
74 * correspond to CHECKME \f[ T' = \frac{T + MTM^\dagger}{2} \f]
76 qpms_tmatrix_t *qpms_tmatrix_symmetrise_involution(
77 const qpms_tmatrix_t *T, ///< the original T-matrix
78 const complex double *M ///< the symmetry op matrix in row-major format
81 /// Creates a \f$ C_\infty \f$ -symmetrized version of a T-matrix.
82 /**
83 * \f[ {T'}_{tlm}^{\tau\lambda\mu} = T_{tlm}^{\tau\lambda\mu} \delta_{m\mu} \f]
85 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_inf(
86 const qpms_tmatrix_t *T ///< the original T-matrix
88 /// Creates a \f$ C_N \f$ -symmetrized version of a T-matrix.
89 /**
90 * \f[ {T'}_{tlm}^{\tau\lambda\mu} = \begin{cases}
91 * T{}_{lm}^{\lambda\mu} & (m-\mu)\mod N=0\\
92 * 0 & (m-\mu)\mod N\ne0
93 * \end{cases} . \f]
95 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_N(
96 const qpms_tmatrix_t *T, ///< the original T-matrix
97 int N ///< number of z-axis rotations in the group
100 /// Symmetrizes a T-matrix with an involution symmetry operation, rewriting the original one.
101 /** The symmetry operation is expected to be a unitary (square)
102 * matrix \a M with the same
103 * VSWF basis spec as the T-matrix (i.e. \a t->spec). The new T-matrix will then
104 * correspond to CHECKME \f[ T' = \frac{T + MTM^\dagger}{2} \f]
106 qpms_tmatrix_t *qpms_tmatrix_symmetrise_involution_inplace(
107 qpms_tmatrix_t *T, ///< the original T-matrix
108 const complex double *M ///< the symmetry op matrix in row-major format
111 /// Creates a \f$ C_\infty \f$ -symmetrized version of a T-matrix, rewriting the original one.
113 * \f[ {T'}_{tlm}^{\tau\lambda\mu} = T_{tlm}^{\tau\lambda\mu} \delta_{m\mu} \f]
115 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_inf_inplace(
116 qpms_tmatrix_t *T ///< the original T-matrix
118 /// Creates a \f$ C_N \f$ -symmetrized version of a T-matrix, rewriting the original one.
120 * \f[ {T'}_{tlm}^{\tau\lambda\mu} = \begin{cases}
121 * T{}_{lm}^{\lambda\mu} & (m-\mu)\mod N=0\\
122 * 0 & (m-\mu)\mod N\ne0
123 * \end{cases} . \f]
125 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_N_inplace(
126 qpms_tmatrix_t *T, ///< the original T-matrix
127 int N ///< number of z-axis rotations in the group
130 /// Reads an open scuff-tmatrix generated file.
131 /**
132 * \a *freqs, \a *freqs_su, \a *tmatrices_array and \a *tmdata
133 * arrays are allocated by this function
134 * and have to be freed by the caller after use.
135 * \a freqs_su and \a tmatrices_array can be NULL, in that case
136 * the respective arrays are not filled nor allocated.
138 * The contents of tmatrices_array is NOT
139 * supposed to be freed element per element.
141 * TODO more checks and options regarding NANs etc.
144 qpms_errno_t qpms_load_scuff_tmatrix(
145 const char *path, ///< Path to the TMatrix file
146 const qpms_vswf_set_spec_t *bspec, ///< VSWF set spec
147 size_t *n, ///< Number of successfully loaded t-matrices
148 double **freqs, ///< Frequencies in SI units..
149 double **freqs_su, ///< Frequencies in SCUFF units (optional).
150 /// The resulting T-matrices (optional).
151 qpms_tmatrix_t **tmatrices_array,
152 complex double **tmdata ///< The T-matrices raw contents
155 /// Tells whether qpms_load_scuff_tmatrix should crash if fopen() fails.
156 /** If true (default), the function causes the program
157 * die e.g. when the tmatrix file
158 * does not exists.
160 * If false, it does nothing and returns an appropriate error value instead.
161 * This is desirable e.g. when used in Python (so that proper exception can
162 * be thrown).
164 extern bool qpms_load_scuff_tmatrix_crash_on_failure;
166 /// Loads a scuff-tmatrix generated file.
167 /** A simple wrapper over qpms_read_scuff_tmatrix() that needs a
168 * path instead of open FILE.
170 * The T-matrix is transformed from the VSWF basis defined by
171 * QPMS_NORMALISATION_CONVENTION_SCUFF into the basis defined
172 * by convention bspec->norm.
174 * Right now, bspec->norm with real or "reversed complex" spherical
175 * harmonics are not supported.
177 qpms_errno_t qpms_read_scuff_tmatrix(
178 FILE *f, ///< An open stream with the T-matrix data.
179 const qpms_vswf_set_spec_t *bspec, ///< VSWF set spec
180 size_t *n, ///< Number of successfully loaded t-matrices
181 double **freqs, ///< Frequencies in SI units.
182 double **freqs_su, ///< Frequencies in SCUFF units (optional).
183 /// The resulting T-matrices (optional).
184 qpms_tmatrix_t **tmatrices_array,
185 /// The T-matrices raw contents.
186 /** The coefficient of outgoing wave defined by
187 * \a bspec->ilist[desti] as a result of incoming wave
188 * \a bspec->ilist[srci] at frequency \a (*freqs)[fi]
189 * is accessed via
190 * (*tmdata)[bspec->n*bspec->n*fi + desti*bspec->n + srci].
192 complex double ** tmdata
195 /// In-place application of point group elements on raw T-matrix data.
196 /** \a tmdata can be e.g. obtained by qpms_load_scuff_tmatrix().
197 * The \a symops array should always contain all elements of a finite
198 * point (sub)group, including the identity operation.
200 * TODO more doc.
202 qpms_errno_t qpms_symmetrise_tmdata_irot3arr(
203 complex double *tmdata, const size_t tmcount,
204 const qpms_vswf_set_spec_t *bspec,
205 size_t n_symops,
206 const qpms_irot3_t *symops
209 /// In-place application of a point group on raw T-matrix data.
210 /** This does the same as qpms_symmetrise_tmdata_irot3arr(),
211 * but takes a valid finite point group as an argument.
213 * TODO more doc.
215 qpms_errno_t qpms_symmetrise_tmdata_finite_group(
216 complex double *tmdata, const size_t tmcount,
217 const qpms_vswf_set_spec_t *bspec,
218 const qpms_finite_group_t *pointgroup
221 /// In-place application of point group elements on a T-matrix.
222 /** The \a symops array should always contain all elements of a finite
223 * point (sub)group, including the identity operation.
225 * TODO more doc.
227 qpms_tmatrix_t *qpms_tmatrix_symmetrise_irot3arr_inplace(
228 qpms_tmatrix_t *T,
229 size_t n_symops,
230 const qpms_irot3_t *symops
233 /// In-place application of point group elements on a T-matrix.
234 /** This does the same as qpms_tmatrix_symmetrise_irot3arr(),
235 * but takes a valid finite point group as an argument.
237 * TODO more doc.
239 qpms_tmatrix_t *qpms_tmatrix_symmetrise_finite_group_inplace(
240 qpms_tmatrix_t *T,
241 const qpms_finite_group_t *pointgroup
244 /// Application of T-matrix on a vector of incident field coefficients, \f$ f = Ta \f$.
245 complex double *qpms_apply_tmatrix(
246 complex double *f_target, ///< Scattered field coefficient array of size T->spec->n; if NULL, a new one is allocated.
247 const complex double *a, ///< Incident field coefficient array of size T->spec->n.
248 const qpms_tmatrix_t *T ///< T-matrix \a T to apply.
251 /// Generic T-matrix generator function that fills a pre-initialised qpms_tmatrix_t given a frequency.
252 /** Implemented by:
253 * qpms_tmatrix_generator_axialsym()
254 * qpms_tmatrix_generator_interpolator()
255 * qpms_tmatrix_generator_sphere()
256 * qpms_tmatrix_generator_constant()
258 typedef struct qpms_tmatrix_generator_t {
259 qpms_errno_t (*function) (qpms_tmatrix_t *t, ///< T-matrix to fill.
260 complex double omega, ///< Angular frequency.
261 const void *params ///< Implementation dependent parameters.
263 const void *params; ///< Parameter pointer passed to the function.
264 } qpms_tmatrix_generator_t;
266 /// Initialises and evaluates a new T-matrix using a generator.
267 qpms_tmatrix_t *qpms_tmatrix_init_from_generator(
268 const qpms_vswf_set_spec_t *bspec,
269 qpms_tmatrix_generator_t gen,
270 complex double omega);
273 /// Implementation of qpms_matrix_generator_t that just copies a constant matrix.
274 /** N.B. this does almost no checks at all, so use it only for t-matrices with
275 * the same base spec.
277 qpms_errno_t qpms_tmatrix_generator_constant(qpms_tmatrix_t *t,
278 complex double omega,
279 /// Source T-matrix, real type is (const qpms_tmatrix_t*).
280 const void *tmatrix_orig
283 /* Fuck this, include the whole <gsl/gsl_spline.h>
284 typedef struct gsl_spline gsl_spline; // Forward declaration for the interpolator struct
285 typedef struct gsl_interp_type gsl_interp_type;
286 extern const gsl_interp_type * gsl_interp_linear;
287 extern const gsl_interp_type * gsl_interp_polynomial;
288 extern const gsl_interp_type * gsl_interp_cspline;
289 extern const gsl_interp_type * gsl_interp_cspline_periodic;
290 extern const gsl_interp_type * gsl_interp_akima;
291 extern const gsl_interp_type * gsl_interp_akima_periodic;
292 extern const gsl_interp_type * gsl_interp_steffen;
295 // struct gsl_interp_accel; // use if lookup proves to be too slow
296 typedef struct qpms_tmatrix_interpolator_t {
297 const qpms_vswf_set_spec_t *bspec;
298 //bool owns_bspec;
299 gsl_spline **splines_real; ///< There will be a spline object for each nonzero element
300 gsl_spline **splines_imag; ///< There will be a spline object for each nonzero element
301 // gsl_interp_accel **accel_real;
302 // gsl_interp_accel **accel_imag;
303 } qpms_tmatrix_interpolator_t;
305 /// Free a T-matrix interpolator.
306 void qpms_tmatrix_interpolator_free(qpms_tmatrix_interpolator_t *interp);
308 /// Fills an existing T-matrix with new interpolated values.
309 qpms_errno_t qpms_tmatrix_interpolator_eval_fill(qpms_tmatrix_t *target, ///< T-matrix to be updated, not NULL.
310 const qpms_tmatrix_interpolator_t *interp, double freq);
312 /// Evaluate a T-matrix interpolated value.
313 /** The result is to be freed using qpms_tmatrix_free().*/
314 qpms_tmatrix_t *qpms_tmatrix_interpolator_eval(const qpms_tmatrix_interpolator_t *interp, double freq);
316 /// Create a T-matrix interpolator from frequency and T-matrix arrays.
317 qpms_tmatrix_interpolator_t *qpms_tmatrix_interpolator_create(size_t n, ///< Number of freqs and T-matrices provided.
318 const double *freqs, const qpms_tmatrix_t *tmatrices_array, ///< N.B. array of qpms_tmatrix_t, not pointers!
319 const gsl_interp_type *iptype
320 //, bool copy_bspec ///< if true, copies its own copy of basis spec from the first T-matrix.
321 /*, ...? */);
324 /// qpms_tmatrix_interpolator for qpms_tmatrix_generator_t.
325 /** As in qpms_tmatrix_interpolator_eval(), the imaginary part of frequency is discarded!
327 qpms_errno_t qpms_tmatrix_generator_interpolator(qpms_tmatrix_t *t, ///< T-matrix to fill.
328 complex double omega, ///< Angular frequency.
329 const void *interpolator ///< Parameter of type qpms_tmatrix_interpolator_t *.
332 /// Calculates the reflection Mie-Lorentz coefficients for a spherical particle.
334 * This function is based on the previous python implementation mie_coefficients() from qpms_p.py,
335 * so any bugs therein should affect this function as well and perhaps vice versa.
337 * Most importantly, the case of magnetic material, \a mu_i != 0 or \a mu_e != 0 has never been tested
338 * and might give wrong results.
340 * \return Array with the Mie-Lorentz reflection coefficients in the order determined by bspec.
341 * If \a target was not NULL, this is target, otherwise a newly allocated array.
343 * TODO better doc.
345 complex double *qpms_mie_coefficients_reflection(
346 complex double *target, ///< Target array of length bspec->n. If NULL, a new one will be allocated.
347 const qpms_vswf_set_spec_t *bspec, ///< Defines which of the coefficients are calculated.
348 double a, ///< Radius of the sphere.
349 complex double k_i, ///< Wave number of the internal material of the sphere.
350 complex double k_e, ///< Wave number of the surrounding medium.
351 complex double mu_i, ///< Relative permeability of the sphere material.
352 complex double mu_e, ///< Relative permeability of the surrounding medium.
353 qpms_bessel_t J_ext, ///< Kind of the "incoming" waves. Most likely QPMS_BESSEL_REGULAR.
354 qpms_bessel_t J_scat ///< Kind of the "scattered" waves. Most likely QPMS_HANKEL_PLUS.
357 /// Replaces the contents of an existing T-matrix with that of a spherical nanoparticle calculated using the Lorentz-mie theory.
358 qpms_errno_t qpms_tmatrix_spherical_fill(qpms_tmatrix_t *t, ///< T-matrix whose contents are to be replaced. Not NULL.
359 double a, ///< Radius of the sphere.
360 complex double k_i, ///< Wave number of the internal material of the sphere.
361 complex double k_e, ///< Wave number of the surrounding medium.
362 complex double mu_i, ///< Relative permeability of the sphere material.
363 complex double mu_e ///< Relative permeability of the surrounding medium.
366 /// Parameter structure for qpms_tmatrix_generator_sphere().
367 typedef struct qpms_tmatrix_generator_sphere_param_t {
368 qpms_epsmu_generator_t outside;
369 qpms_epsmu_generator_t inside;
370 double radius;
371 } qpms_tmatrix_generator_sphere_param_t;
373 /// T-matrix generator for spherical particles using Lorentz-Mie solution.
374 qpms_errno_t qpms_tmatrix_generator_sphere(qpms_tmatrix_t *t,
375 complex double omega,
376 const void *params ///< Of type qpms_tmatrix_generator_sphere_param_t.
379 /// Creates a new T-matrix of a spherical particle using the Lorentz-Mie theory.
380 static inline qpms_tmatrix_t *qpms_tmatrix_spherical(
381 const qpms_vswf_set_spec_t *bspec,
382 double a, ///< Radius of the sphere.
383 complex double k_i, ///< Wave number of the internal material of the sphere.
384 complex double k_e, ///< Wave number of the surrounding medium.
385 complex double mu_i, ///< Relative permeability of the sphere material.
386 complex double mu_e ///< Relative permeability of the surrounding medium.
388 qpms_tmatrix_t *t = qpms_tmatrix_init(bspec);
389 qpms_tmatrix_spherical_fill(t, a, k_i, k_e, mu_i, mu_e);
390 return t;
393 /// Convenience function to calculate T-matrix of a non-magnetic spherical particle using the permittivity values, replacing existing T-matrix data.
394 qpms_errno_t qpms_tmatrix_spherical_mu0_fill(
395 qpms_tmatrix_t *t, ///< T-matrix whose contents are to be replaced. Not NULL.
396 double a, ///< Radius of the sphere.
397 double omega, ///< Angular frequency.
398 complex double epsilon_fg, ///< Relative permittivity of the sphere.
399 complex double epsilon_bg ///< Relative permittivity of the background medium.
402 /// Convenience function to calculate T-matrix of a non-magnetic spherical particle using the permittivity values.
403 static inline qpms_tmatrix_t *qpms_tmatrix_spherical_mu0(
404 const qpms_vswf_set_spec_t *bspec,
405 double a, ///< Radius of the sphere.
406 double omega, ///< Angular frequency.
407 complex double epsilon_fg, ///< Relative permittivity of the sphere.
408 complex double epsilon_bg ///< Relative permittivity of the background medium.
410 qpms_tmatrix_t *t = qpms_tmatrix_init(bspec);
411 qpms_tmatrix_spherical_mu0_fill(t, a, omega, epsilon_fg, epsilon_bg);
412 return t;
415 /// Return value type for qpms_arc_function_t.
416 typedef struct qpms_arc_function_retval_t {
417 double r; ///< Distance from the origin.
418 double beta; ///< Angle between surface normal and radial direction.
419 } qpms_arc_function_retval_t;
421 /// Prototype for general parametrisation of \f$ C_\infty \f$-symmetric particle's surface.
422 typedef struct qpms_arc_function_t {
423 /// Arc parametrisation function.
424 /** TODO link to notes.
426 * Implemented by:
427 * qpms_arc_cylinder(),
428 * qpms_arc_sphere().
430 qpms_arc_function_retval_t (*function) (
431 double theta, ///< Polar angle from interval \f$ [0, \pi] \f$.
432 const void *params ///< Pointer to implementation specific parameters.
434 const void *params;
435 } qpms_arc_function_t;
437 /// Parameter structure for qpms_arc_cylinder().
438 typedef struct qpms_arc_cylinder_params_t {
439 double R; ///< Cylinder radius.
440 double h; ///< Cylinder height.
441 } qpms_arc_cylinder_params_t;
443 /// Arc parametrisation of cylindrical particle; for qpms_arc_function_t.
444 qpms_arc_function_retval_t qpms_arc_cylinder(double theta,
445 const void *params ///< Points to qpms_arc_cylinder_params_t
448 /// Arc parametrisation of spherical particle; for qpms_arc_function_t.
449 /** Useful mostly only for benchmarks or debugging, as one can use the Mie-Lorentz solution. */
450 qpms_arc_function_retval_t qpms_arc_sphere(double theta,
451 const void *R ///< Points to double containing particle's radius.
455 /// Replaces T-matrix contents with those of a particle with \f$ C_\infty \f$ symmetry.
456 /**
457 * N.B. currently, this might crash for higher values of lMax (lMax > 5).
458 * Also, it seems that I am doing something wrong, as the result is accurate for sphere
459 * with lMax = 1 and for higher l the accuracy decreases.
461 qpms_errno_t qpms_tmatrix_axialsym_fill(
462 qpms_tmatrix_t *t, ///< T-matrix whose contents are to be replaced. Not NULL.
463 complex double omega, ///< Angular frequency.
464 qpms_epsmu_t outside, ///< Optical properties of the outside medium.
465 qpms_epsmu_t inside, ///< Optical properties of the particle's material.
466 qpms_arc_function_t shape, ///< Particle surface parametrisation.
467 /** If `lMax_extend > t->bspec->lMax`, then the internal \a Q, \a R matrices will be
468 * trimmed at this larger lMax; the final T-matrix will then be trimmed
469 * according to bspec.
471 qpms_l_t lMax_extend
474 /// Creates a new T-matrix of a particle with \f$ C_\infty \f$ symmetry.
475 static inline qpms_tmatrix_t *qpms_tmatrix_axialsym(
476 const qpms_vswf_set_spec_t *bspec,
477 complex double omega, ///< Angular frequency.
478 qpms_epsmu_t outside, ///< Optical properties of the outside medium.
479 qpms_epsmu_t inside, ///< Optical properties of the particle's material.
480 qpms_arc_function_t shape, ///< Particle surface parametrisation.
481 /// Internal lMax to improve precision of the result.
482 /** If `lMax_extend > bspec->lMax`, then the internal \a Q, \a R matrices will be
483 * trimmed at this larger lMax; the final T-matrix will then be trimmed
484 * according to bspec.
486 qpms_l_t lMax_extend
488 qpms_tmatrix_t *t = qpms_tmatrix_init(bspec);
489 qpms_tmatrix_axialsym_fill(t, omega, outside, inside, shape, lMax_extend);
490 return t;
493 /// Parameter structure for qpms_tmatrix_generator_axialsym.
494 typedef struct qpms_tmatrix_generator_axialsym_param_t {
495 qpms_epsmu_generator_t outside;
496 qpms_epsmu_generator_t inside;
497 qpms_arc_function_t shape;
498 qpms_l_t lMax_extend;
499 } qpms_tmatrix_generator_axialsym_param_t;
502 /// qpms_tmatrix_axialsym for qpms_tmatrix_generator_t
503 qpms_errno_t qpms_tmatrix_generator_axialsym(qpms_tmatrix_t *t, ///< T-matrix to fill.
504 complex double omega, ///< Angular frequency.
505 const void *params ///< Parameters of type qpms_tmatrix_generator_axialsym_param_t.
509 /// Computes the (reduced) transposed R or Q matrix for axially symmetric particle (useful for debugging).
510 qpms_errno_t qpms_tmatrix_generator_axialsym_RQ_transposed_fill(complex double *target,
511 complex double omega,
512 const qpms_tmatrix_generator_axialsym_param_t *param,
513 qpms_normalisation_t norm,
514 qpms_bessel_t J
517 /// Computes the (reduced) transposed R or Q matrix for axially symmetric particle (useful mostly for debugging).
518 qpms_errno_t qpms_tmatrix_axialsym_RQ_transposed_fill(complex double *target,
519 complex double omega, qpms_epsmu_t outside, qpms_epsmu_t inside,
520 qpms_arc_function_t shape, qpms_l_t lMaxQR, qpms_normalisation_t norm,
521 qpms_bessel_t J ///< Use QPMS_BESSEL_REGULAR to calculate \f$ R^T\f$ or QPMS_HANKEL_PLUS to calculate \f$ Q^T\f$.
525 /// An "abstract" T-matrix, contains a T-matrix generator instead of actual data.
526 typedef struct qpms_tmatrix_function_t {
527 /** \brief VSWF basis specification, NOT owned by qpms_tmatrix_t by default.
529 * Usually not checked for meaningfulness by the functions (methods),
530 * so the caller should take care that \a spec->ilist does not
531 * contain any duplicities and that for each wave with order \a m
532 * there is also one with order \a −m.
534 const qpms_vswf_set_spec_t *spec;
535 const qpms_tmatrix_generator_t *gen; ///< A T-matrix generator function.
536 } qpms_tmatrix_function_t;
538 /// Convenience function to create a new T-matrix from qpms_tmatrix_function_t.
539 // FIXME the name is not very intuitive.
540 static inline qpms_tmatrix_t *qpms_tmatrix_init_from_function(qpms_tmatrix_function_t f, complex double omega) {
541 return qpms_tmatrix_init_from_generator(f.spec, *f.gen, omega);
544 /// Specifies different kinds of operations done on T-matrices
545 typedef enum {
546 QPMS_TMATRIX_OPERATION_NOOP, ///< Identity operation.
547 QPMS_TMATRIX_OPERATION_LRMATRIX, ///< General matrix transformation \f$ T' = MTM^\dagger \f$; see @ref qpms_tmatrix_operation_lrmatrix.
548 QPMS_TMATRIX_OPERATION_IROT3, ///< Single rotoreflection specified by a qpms_irot3_t.
549 QPMS_TMATRIX_OPERATION_IROT3ARR, ///< Symmetrise using an array of rotoreflection; see @ref qpms_tmatrix_operation_irot3arr.
550 QPMS_TMATRIX_OPERATION_COMPOSE_SUM, ///< Apply several transformations and sum the results, see @ref qpms_tmatrix_operation_compose_sum.
551 QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN, ///< Chain several different transformations; see @ref qpms_tmatrix_operation_compose_chain.
552 QPMS_TMATRIX_OPERATION_SCMULZ, ///< Elementwise scalar multiplication with a complex matrix; see @ref qpms_tmatrix_operation_scmulz.
553 //QPMS_TMATRIX_OPERATION_POINTGROUP, ///< TODO the new point group in pointgroup.h
554 QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE ///< Legacy qpms_finite_group_t with filled rep3d; see @ref qpms_tmatrix_operation_finite_group.
555 } qpms_tmatrix_operation_kind_t;
557 /// General matrix transformation \f[ T' = MTM^\dagger \f] spec for qpms_tmatrix_operation_t.
558 struct qpms_tmatrix_operation_lrmatrix {
559 /// Raw matrix data of \a M in row-major order.
560 /** The matrix must be taylored for the given bspec! */
561 complex double *m;
562 size_t m_size; ///< Total size of \a m matrix in terms of sizeof(complex double).
563 bool owns_m; ///< Whether \a m is owned by this;
566 struct qpms_tmatrix_operation_t; // Forward declaration for the composed operations.
568 /// Specifies a composed operation of type \f$ T' = c\sum_i f_i'(T) \f$ for qpms_tmatrix_operation_t.
569 struct qpms_tmatrix_operation_compose_sum {
570 size_t n; ///< Number of operations in ops;
571 const struct qpms_tmatrix_operation_t **ops; ///< Operations array. (Pointers array \a ops[] always owned by this.)
572 double factor; ///< Overall factor \a c.
573 /// (Optional) operations buffer into which elements of \a ops point.
574 /** Owned by this or NULL. If not NULL, all \a ops members are assumed to point into
575 * the memory held by \a opmem and to be properly initialised.
576 * Each \a ops member has to point to _different_ elements of \a opmem.
578 struct qpms_tmatrix_operation_t *opmem;
581 /// Specifies a composed operation of type \f$ T' = f_{n-1}(f_{n-2}(\dots f_0(T)\dots))) \f$ for qpms_tmatrix_operation_t.
582 struct qpms_tmatrix_operation_compose_chain {
583 size_t n; ///< Number of operations in ops;
584 const struct qpms_tmatrix_operation_t **ops; ///< Operations array. (Pointers owned by this.)
585 struct qpms_tmatrix_operation_t *opmem; ///< (Optional) operations buffer into which elements of \a ops point. (Owned by this or NULL.)
586 size_t opmem_size; ///< Length of the opmem array.
587 _Bool *ops_owned; ///< True for all sub operations owned by this and saved in opmem. NULL if opmem is NULL.
590 /// Specifies an elementwise complex multiplication of type \f$ T'_{ij} = M_{ij}T_{ij} \f$ for qpms_tmatrix_operation_t.
591 struct qpms_tmatrix_operation_scmulz {
592 /// Raw matrix data of \a M in row-major order.
593 /** The matrix must be taylored for the given bspec! */
594 complex double *m;
595 size_t m_size; ///< Total size of \a m matrix in terms of sizeof(complex double).
596 bool owns_m; ///< Whether \a m is owned by this.
599 /// Specifies a symmetrisation using a set of rotoreflections (with equal weights) for qpms_tmatrix_operation_t.
600 /** Internally, this is evaluated using a call to qpms_symmetrise_tmdata_irot3arr()
601 * or qpms_symmetrise_tmdata_irot3arr_inplace(). */
602 struct qpms_tmatrix_operation_irot3arr {
603 size_t n; ///< Number of rotoreflections;
604 qpms_irot3_t *ops; ///< Rotoreflection array of size \a n.
605 bool owns_ops; ///< Whether \a ops array is owned by this.
608 /// A generic T-matrix transformation operator.
609 typedef struct qpms_tmatrix_operation_t {
610 /// Specifies the operation kind to be performed and which type \op actually contains.
611 qpms_tmatrix_operation_kind_t typ;
612 union {
613 struct qpms_tmatrix_operation_lrmatrix lrmatrix;
614 struct qpms_tmatrix_operation_scmulz scmulz;
615 qpms_irot3_t irot3; ///< Single rotoreflection; \a typ = QPMS_TMATRIX_OPERATION_IROT3
616 struct qpms_tmatrix_operation_irot3arr irot3arr;
617 struct qpms_tmatrix_operation_compose_sum compose_sum;
618 struct qpms_tmatrix_operation_compose_chain compose_chain;
619 /// Finite group for QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE.
620 /** Not owned by this; \a rep3d must be filled. */
621 const qpms_finite_group_t *finite_group;
622 } op; ///< Operation data; actual type is determined by \a typ.
623 } qpms_tmatrix_operation_t;
625 static const qpms_tmatrix_operation_t qpms_tmatrix_operation_noop = {.typ = QPMS_TMATRIX_OPERATION_NOOP};
627 /// Apply an operation on a T-matrix, returning a newly allocated result.
628 qpms_tmatrix_t *qpms_tmatrix_apply_operation(const qpms_tmatrix_operation_t *op, const qpms_tmatrix_t *orig);
630 /// Apply an operation on a T-matrix and replace it with the result.
631 qpms_tmatrix_t *qpms_tmatrix_apply_operation_inplace(const qpms_tmatrix_operation_t *op, qpms_tmatrix_t *orig);
633 /// Apply an operation on a T-matrix and replace another one it with the result.
634 qpms_tmatrix_t *qpms_tmatrix_apply_operation_replace(qpms_tmatrix_t *dest,
635 const qpms_tmatrix_operation_t *op, const qpms_tmatrix_t *orig);
637 /// (Recursively) deallocates internal data of qpms_tmatrix_operation_t.
638 /** It does NOT clear the memory pointed to by op it self, but only
639 * heap-allocated data of certain operations, if labeled as owned by it.
640 * In case of compose operations, the recursion stops if the children are
641 * not owned by them (the opmem pointer is NULL).
643 void qpms_tmatrix_operation_clear(qpms_tmatrix_operation_t *f);
645 /// (Recursively) copies an qpms_tmatrix_operation_t.
646 /** Makes copies of all the internal data and takes ownership over them if needed */
647 void qpms_tmatrix_operation_copy(qpms_tmatrix_operation_t *target, const qpms_tmatrix_operation_t *src);
649 /// Inits a new "chain" of composed operations, some of which might be owned.
650 void qpms_tmatrix_operation_compose_chain_init(
651 qpms_tmatrix_operation_t *target, ///< The operation structure that will be set to the chain.
652 size_t nops, ///< Number of chained operations (length of the \a ops array)
653 size_t opmem_size ///< Size of the own operations buffer (length of the \a opmem array)
656 #if 0
657 // Abstract types that describe T-matrix/particle/scatsystem symmetries
658 // To be implemented later. See also the thoughts in the beginning of groups.h.
660 typedef *char qpms_tmatrix_id_t; ///< Maybe I want some usual integer type instead.
662 ///Abstract T-matrix type draft.
664 * TODO.
666 typedef struct qpms_abstract_tmatrix_t{
667 qpms_tmatrix_id_t id;
668 /// Generators of the discrete point group under which T-matrix is invariant.
669 qpms_irot3_t *invar_gens;
670 /// Length of invar_gens.
671 qpms_gmi_t invar_gens_size;
673 } qpms_abstract_tmatrix_t;
676 typedef struct qpms_abstract_particle_t{
677 } qpms_abstract_particle_t;
679 /// An abstract particle, defined by its position and abstract T-matrix.
680 typedef struct qpms_abstract_particle_t {
681 cart3_t pos; ///< Particle position in cartesian coordinates.
682 const qpms_abstract_tmatrix_t *tmatrix; ///< T-matrix; not owned by this.
683 } qpms_abstract_particle_t;
686 /** This is just an alias, as the same index can be used for
687 * abstract T-matrices as well.
689 typedef qpms_particle_tid_t qpms_abstract_particle_tid_t;
690 #endif // 0
692 #endif //TMATRICES_H