Rewriting scatsystem: compiles without errors now.
[qpms.git] / qpms / tmatrices.h
blob8b4d2e6bd2cb77c160250f14a93af95beb311697
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 /// Destroys a T-matrix.
31 void qpms_tmatrix_free(qpms_tmatrix_t *t);
33 /// Check T-matrix equality/similarity.
34 /**
35 * This function actually checks for identical vswf specs.
36 * TODO define constants with "default" atol, rtol for this function.
38 bool qpms_tmatrix_isclose(const qpms_tmatrix_t *T1, const qpms_tmatrix_t *T2,
39 const double rtol, const double atol);
41 /// Creates a T-matrix from another matrix and a symmetry operation.
42 /** The symmetry operation is expected to be a unitary (square)
43 * matrix \a M with the same
44 * VSWF basis spec as the T-matrix (i.e. \a t->spec). The new T-matrix will then
45 * correspond to CHECKME \f[ T' = MTM^\dagger \f]
47 qpms_tmatrix_t *qpms_tmatrix_apply_symop(
48 const qpms_tmatrix_t *T, ///< the original T-matrix
49 const complex double *M ///< the symmetry op matrix in row-major format
52 /// Applies a symmetry operation onto a T-matrix, rewriting the original T-matrix data.
53 /** The symmetry operation is expected to be a unitary (square)
54 * matrix \a M with the same
55 * VSWF basis spec as the T-matrix (i.e. \a t->spec). The new T-matrix will then
56 * correspond to CHECKME \f[ T' = MTM^\dagger \f]
58 qpms_tmatrix_t *qpms_tmatrix_apply_symop_inplace(
59 qpms_tmatrix_t *T, ///< the original T-matrix
60 const complex double *M ///< the symmetry op matrix in row-major format
63 /// Symmetrizes a T-matrix with an involution symmetry operation.
64 /** The symmetry operation is expected to be a unitary (square)
65 * matrix \a M with the same
66 * VSWF basis spec as the T-matrix (i.e. \a t->spec). The new T-matrix will then
67 * correspond to CHECKME \f[ T' = \frac{T + MTM^\dagger}{2} \f]
69 qpms_tmatrix_t *qpms_tmatrix_symmetrise_involution(
70 const qpms_tmatrix_t *T, ///< the original T-matrix
71 const complex double *M ///< the symmetry op matrix in row-major format
74 /// Creates a \f$ C_\infty \f$ -symmetrized version of a T-matrix.
75 /**
76 * \f[ {T'}_{tlm}^{\tau\lambda\mu} = T_{tlm}^{\tau\lambda\mu} \delta_{m\mu} \f]
78 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_inf(
79 const qpms_tmatrix_t *T ///< the original T-matrix
81 /// Creates a \f$ C_N \f$ -symmetrized version of a T-matrix.
82 /**
83 * \f[ {T'}_{tlm}^{\tau\lambda\mu} = \begin{cases}
84 * T{}_{lm}^{\lambda\mu} & (m-\mu)\mod N=0\\
85 * 0 & (m-\mu)\mod N\ne0
86 * \end{cases} . \f]
88 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_N(
89 const qpms_tmatrix_t *T, ///< the original T-matrix
90 int N ///< number of z-axis rotations in the group
93 /// Symmetrizes a T-matrix with an involution symmetry operation, rewriting the original one.
94 /** The symmetry operation is expected to be a unitary (square)
95 * matrix \a M with the same
96 * VSWF basis spec as the T-matrix (i.e. \a t->spec). The new T-matrix will then
97 * correspond to CHECKME \f[ T' = \frac{T + MTM^\dagger}{2} \f]
99 qpms_tmatrix_t *qpms_tmatrix_symmetrise_involution_inplace(
100 qpms_tmatrix_t *T, ///< the original T-matrix
101 const complex double *M ///< the symmetry op matrix in row-major format
104 /// Creates a \f$ C_\infty \f$ -symmetrized version of a T-matrix, rewriting the original one.
106 * \f[ {T'}_{tlm}^{\tau\lambda\mu} = T_{tlm}^{\tau\lambda\mu} \delta_{m\mu} \f]
108 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_inf_inplace(
109 qpms_tmatrix_t *T ///< the original T-matrix
111 /// Creates a \f$ C_N \f$ -symmetrized version of a T-matrix, rewriting the original one.
113 * \f[ {T'}_{tlm}^{\tau\lambda\mu} = \begin{cases}
114 * T{}_{lm}^{\lambda\mu} & (m-\mu)\mod N=0\\
115 * 0 & (m-\mu)\mod N\ne0
116 * \end{cases} . \f]
118 qpms_tmatrix_t *qpms_tmatrix_symmetrise_C_N_inplace(
119 qpms_tmatrix_t *T, ///< the original T-matrix
120 int N ///< number of z-axis rotations in the group
123 /// Reads an open scuff-tmatrix generated file.
124 /**
125 * \a *freqs, \a *freqs_su, \a *tmatrices_array and \a *tmdata
126 * arrays are allocated by this function
127 * and have to be freed by the caller after use.
128 * \a freqs_su and \a tmatrices_array can be NULL, in that case
129 * the respective arrays are not filled nor allocated.
131 * The contents of tmatrices_array is NOT
132 * supposed to be freed element per element.
134 * TODO more checks and options regarding NANs etc.
137 qpms_errno_t qpms_load_scuff_tmatrix(
138 const char *path, ///< Path to the TMatrix file
139 const qpms_vswf_set_spec_t *bspec, ///< VSWF set spec
140 size_t *n, ///< Number of successfully loaded t-matrices
141 double **freqs, ///< Frequencies in SI units..
142 double **freqs_su, ///< Frequencies in SCUFF units (optional).
143 /// The resulting T-matrices (optional).
144 qpms_tmatrix_t **tmatrices_array,
145 complex double **tmdata ///< The T-matrices raw contents
148 /// Tells whether qpms_load_scuff_tmatrix should crash if fopen() fails.
149 /** If true (default), the function causes the program
150 * die e.g. when the tmatrix file
151 * does not exists.
153 * If false, it does nothing and returns an appropriate error value instead.
154 * This is desirable e.g. when used in Python (so that proper exception can
155 * be thrown).
157 extern bool qpms_load_scuff_tmatrix_crash_on_failure;
159 /// Loads a scuff-tmatrix generated file.
160 /** A simple wrapper over qpms_read_scuff_tmatrix() that needs a
161 * path instead of open FILE.
163 * The T-matrix is transformed from the VSWF basis defined by
164 * QPMS_NORMALISATION_CONVENTION_SCUFF into the basis defined
165 * by convention bspec->norm.
167 * Right now, bspec->norm with real or "reversed complex" spherical
168 * harmonics are not supported.
170 qpms_errno_t qpms_read_scuff_tmatrix(
171 FILE *f, ///< An open stream with the T-matrix data.
172 const qpms_vswf_set_spec_t *bspec, ///< VSWF set spec
173 size_t *n, ///< Number of successfully loaded t-matrices
174 double **freqs, ///< Frequencies in SI units.
175 double **freqs_su, ///< Frequencies in SCUFF units (optional).
176 /// The resulting T-matrices (optional).
177 qpms_tmatrix_t **tmatrices_array,
178 /// The T-matrices raw contents.
179 /** The coefficient of outgoing wave defined by
180 * \a bspec->ilist[desti] as a result of incoming wave
181 * \a bspec->ilist[srci] at frequency \a (*freqs)[fi]
182 * is accessed via
183 * (*tmdata)[bspec->n*bspec->n*fi + desti*bspec->n + srci].
185 complex double ** tmdata
188 /// In-place application of point group elements on raw T-matrix data.
189 /** \a tmdata can be e.g. obtained by qpms_load_scuff_tmatrix().
190 * The \a symops array should always contain all elements of a finite
191 * point (sub)group, including the identity operation.
193 * TODO more doc.
195 qpms_errno_t qpms_symmetrise_tmdata_irot3arr(
196 complex double *tmdata, const size_t tmcount,
197 const qpms_vswf_set_spec_t *bspec,
198 size_t n_symops,
199 const qpms_irot3_t *symops
202 /// In-place application of a point group on raw T-matrix data.
203 /** This does the same as qpms_symmetrise_tmdata_irot3arr(),
204 * but takes a valid finite point group as an argument.
206 * TODO more doc.
208 qpms_errno_t qpms_symmetrise_tmdata_finite_group(
209 complex double *tmdata, const size_t tmcount,
210 const qpms_vswf_set_spec_t *bspec,
211 const qpms_finite_group_t *pointgroup
214 /// In-place application of point group elements on a T-matrix.
215 /** The \a symops array should always contain all elements of a finite
216 * point (sub)group, including the identity operation.
218 * TODO more doc.
220 qpms_tmatrix_t *qpms_tmatrix_symmetrise_irot3arr_inplace(
221 qpms_tmatrix_t *T,
222 size_t n_symops,
223 const qpms_irot3_t *symops
226 /// In-place application of point group elements on a T-matrix.
227 /** This does the same as qpms_tmatrix_symmetrise_irot3arr(),
228 * but takes a valid finite point group as an argument.
230 * TODO more doc.
232 qpms_tmatrix_t *qpms_tmatrix_symmetrise_finite_group_inplace(
233 qpms_tmatrix_t *T,
234 const qpms_finite_group_t *pointgroup
237 /// Application of T-matrix on a vector of incident field coefficients, \f$ f = Ta \f$.
238 complex double *qpms_apply_tmatrix(
239 complex double *f_target, ///< Scattered field coefficient array of size T->spec->n; if NULL, a new one is allocated.
240 const complex double *a, ///< Incident field coefficient array of size T->spec->n.
241 const qpms_tmatrix_t *T ///< T-matrix \a T to apply.
244 /// Generic T-matrix generator function that fills a pre-initialised qpms_tmatrix_t given a frequency.
245 /** Implemented by:
246 * qpms_tmatrix_generator_axialsym()
247 * qpms_tmatrix_generator_interpolator()
248 * qpms_tmatrix_generator_sphere()
249 * qpms_tmatrix_generator_constant()
251 typedef struct qpms_tmatrix_generator_t {
252 qpms_errno_t (*function) (qpms_tmatrix_t *t, ///< T-matrix to fill.
253 complex double omega, ///< Angular frequency.
254 const void *params ///< Implementation dependent parameters.
256 const void *params; ///< Parameter pointer passed to the function.
257 } qpms_tmatrix_generator_t;
259 /// Initialises and evaluates a new T-matrix using a generator.
260 qpms_tmatrix_t *qpms_tmatrix_init_from_generator(
261 const qpms_vswf_set_spec_t *bspec,
262 qpms_tmatrix_generator_t gen,
263 complex double omega);
266 /// Implementation of qpms_matrix_generator_t that just copies a constant matrix.
267 /** N.B. this does almost no checks at all, so use it only for t-matrices with
268 * the same base spec.
270 qpms_errno_t qpms_tmatrix_generator_constant(qpms_tmatrix_t *t,
271 complex double omega,
272 /// Source T-matrix, real type is (const qpms_tmatrix_t*).
273 const void *tmatrix_orig
276 /* Fuck this, include the whole <gsl/gsl_spline.h>
277 typedef struct gsl_spline gsl_spline; // Forward declaration for the interpolator struct
278 typedef struct gsl_interp_type gsl_interp_type;
279 extern const gsl_interp_type * gsl_interp_linear;
280 extern const gsl_interp_type * gsl_interp_polynomial;
281 extern const gsl_interp_type * gsl_interp_cspline;
282 extern const gsl_interp_type * gsl_interp_cspline_periodic;
283 extern const gsl_interp_type * gsl_interp_akima;
284 extern const gsl_interp_type * gsl_interp_akima_periodic;
285 extern const gsl_interp_type * gsl_interp_steffen;
288 // struct gsl_interp_accel; // use if lookup proves to be too slow
289 typedef struct qpms_tmatrix_interpolator_t {
290 const qpms_vswf_set_spec_t *bspec;
291 //bool owns_bspec;
292 gsl_spline **splines_real; ///< There will be a spline object for each nonzero element
293 gsl_spline **splines_imag; ///< There will be a spline object for each nonzero element
294 // gsl_interp_accel **accel_real;
295 // gsl_interp_accel **accel_imag;
296 } qpms_tmatrix_interpolator_t;
298 /// Free a T-matrix interpolator.
299 void qpms_tmatrix_interpolator_free(qpms_tmatrix_interpolator_t *interp);
301 /// Fills an existing T-matrix with new interpolated values.
302 qpms_errno_t qpms_tmatrix_interpolator_eval_fill(qpms_tmatrix_t *target, ///< T-matrix to be updated, not NULL.
303 const qpms_tmatrix_interpolator_t *interp, double freq);
305 /// Evaluate a T-matrix interpolated value.
306 /** The result is to be freed using qpms_tmatrix_free().*/
307 qpms_tmatrix_t *qpms_tmatrix_interpolator_eval(const qpms_tmatrix_interpolator_t *interp, double freq);
309 /// Create a T-matrix interpolator from frequency and T-matrix arrays.
310 qpms_tmatrix_interpolator_t *qpms_tmatrix_interpolator_create(size_t n, ///< Number of freqs and T-matrices provided.
311 const double *freqs, const qpms_tmatrix_t *tmatrices_array, ///< N.B. array of qpms_tmatrix_t, not pointers!
312 const gsl_interp_type *iptype
313 //, bool copy_bspec ///< if true, copies its own copy of basis spec from the first T-matrix.
314 /*, ...? */);
317 /// qpms_tmatrix_interpolator for qpms_tmatrix_generator_t.
318 /** As in qpms_tmatrix_interpolator_eval(), the imaginary part of frequency is discarded!
320 qpms_errno_t qpms_tmatrix_generator_interpolator(qpms_tmatrix_t *t, ///< T-matrix to fill.
321 complex double omega, ///< Angular frequency.
322 const void *interpolator ///< Parameter of type qpms_tmatrix_interpolator_t *.
325 /// Calculates the reflection Mie-Lorentz coefficients for a spherical particle.
327 * This function is based on the previous python implementation mie_coefficients() from qpms_p.py,
328 * so any bugs therein should affect this function as well and perhaps vice versa.
330 * Most importantly, the case of magnetic material, \a mu_i != 0 or \a mu_e != 0 has never been tested
331 * and might give wrong results.
333 * \return Array with the Mie-Lorentz reflection coefficients in the order determined by bspec.
334 * If \a target was not NULL, this is target, otherwise a newly allocated array.
336 * TODO better doc.
338 complex double *qpms_mie_coefficients_reflection(
339 complex double *target, ///< Target array of length bspec->n. If NULL, a new one will be allocated.
340 const qpms_vswf_set_spec_t *bspec, ///< Defines which of the coefficients are calculated.
341 double a, ///< Radius of the sphere.
342 complex double k_i, ///< Wave number of the internal material of the sphere.
343 complex double k_e, ///< Wave number of the surrounding medium.
344 complex double mu_i, ///< Relative permeability of the sphere material.
345 complex double mu_e, ///< Relative permeability of the surrounding medium.
346 qpms_bessel_t J_ext, ///< Kind of the "incoming" waves. Most likely QPMS_BESSEL_REGULAR.
347 qpms_bessel_t J_scat ///< Kind of the "scattered" waves. Most likely QPMS_HANKEL_PLUS.
350 /// Replaces the contents of an existing T-matrix with that of a spherical nanoparticle calculated using the Lorentz-mie theory.
351 qpms_errno_t qpms_tmatrix_spherical_fill(qpms_tmatrix_t *t, ///< T-matrix whose contents are to be replaced. Not NULL.
352 double a, ///< Radius of the sphere.
353 complex double k_i, ///< Wave number of the internal material of the sphere.
354 complex double k_e, ///< Wave number of the surrounding medium.
355 complex double mu_i, ///< Relative permeability of the sphere material.
356 complex double mu_e ///< Relative permeability of the surrounding medium.
359 /// Parameter structure for qpms_tmatrix_generator_sphere().
360 typedef struct qpms_tmatrix_generator_sphere_param_t {
361 qpms_epsmu_generator_t outside;
362 qpms_epsmu_generator_t inside;
363 double radius;
364 } qpms_tmatrix_generator_sphere_param_t;
366 /// T-matrix generator for spherical particles using Lorentz-Mie solution.
367 qpms_errno_t qpms_tmatrix_generator_sphere(qpms_tmatrix_t *t,
368 complex double omega,
369 const void *params ///< Of type qpms_tmatrix_generator_sphere_param_t.
372 /// Creates a new T-matrix of a spherical particle using the Lorentz-Mie theory.
373 static inline qpms_tmatrix_t *qpms_tmatrix_spherical(
374 const qpms_vswf_set_spec_t *bspec,
375 double a, ///< Radius of the sphere.
376 complex double k_i, ///< Wave number of the internal material of the sphere.
377 complex double k_e, ///< Wave number of the surrounding medium.
378 complex double mu_i, ///< Relative permeability of the sphere material.
379 complex double mu_e ///< Relative permeability of the surrounding medium.
381 qpms_tmatrix_t *t = qpms_tmatrix_init(bspec);
382 qpms_tmatrix_spherical_fill(t, a, k_i, k_e, mu_i, mu_e);
383 return t;
386 /// Convenience function to calculate T-matrix of a non-magnetic spherical particle using the permittivity values, replacing existing T-matrix data.
387 qpms_errno_t qpms_tmatrix_spherical_mu0_fill(
388 qpms_tmatrix_t *t, ///< T-matrix whose contents are to be replaced. Not NULL.
389 double a, ///< Radius of the sphere.
390 double omega, ///< Angular frequency.
391 complex double epsilon_fg, ///< Relative permittivity of the sphere.
392 complex double epsilon_bg ///< Relative permittivity of the background medium.
395 /// Convenience function to calculate T-matrix of a non-magnetic spherical particle using the permittivity values.
396 static inline qpms_tmatrix_t *qpms_tmatrix_spherical_mu0(
397 const qpms_vswf_set_spec_t *bspec,
398 double a, ///< Radius of the sphere.
399 double omega, ///< Angular frequency.
400 complex double epsilon_fg, ///< Relative permittivity of the sphere.
401 complex double epsilon_bg ///< Relative permittivity of the background medium.
403 qpms_tmatrix_t *t = qpms_tmatrix_init(bspec);
404 qpms_tmatrix_spherical_mu0_fill(t, a, omega, epsilon_fg, epsilon_bg);
405 return t;
408 /// Return value type for qpms_arc_function_t.
409 typedef struct qpms_arc_function_retval_t {
410 double r; ///< Distance from the origin.
411 double beta; ///< Angle between surface normal and radial direction.
412 } qpms_arc_function_retval_t;
414 /// Prototype for general parametrisation of \f$ C_\infty \f$-symmetric particle's surface.
415 typedef struct qpms_arc_function_t {
416 /// Arc parametrisation function.
417 /** TODO link to notes.
419 * Implemented by:
420 * qpms_arc_cylinder(),
421 * qpms_arc_sphere().
423 qpms_arc_function_retval_t (*function) (
424 double theta, ///< Polar angle from interval \f$ [0, \pi] \f$.
425 const void *params ///< Pointer to implementation specific parameters.
427 const void *params;
428 } qpms_arc_function_t;
430 /// Parameter structure for qpms_arc_cylinder().
431 typedef struct qpms_arc_cylinder_params_t {
432 double R; ///< Cylinder radius.
433 double h; ///< Cylinder height.
434 } qpms_arc_cylinder_params_t;
436 /// Arc parametrisation of cylindrical particle; for qpms_arc_function_t.
437 qpms_arc_function_retval_t qpms_arc_cylinder(double theta,
438 const void *params; ///< Points to qpms_arc_cylinder_params_t
441 /// Arc parametrisation of spherical particle; for qpms_arc_function_t.
442 /** Useful mostly only for benchmarks or debugging, as one can use the Mie-Lorentz solution. */
443 qpms_arc_function_retval_t qpms_arc_sphere(double theta,
444 const void *R; ///< Points to double containing particle's radius.
448 /// Replaces T-matrix contents with those of a particle with \f$ C_\infty \f$ symmetry.
449 /**
450 * N.B. currently, this might crash for higher values of lMax (lMax > 5).
451 * Also, it seems that I am doing something wrong, as the result is accurate for sphere
452 * with lMax = 1 and for higher l the accuracy decreases.
454 qpms_errno_t qpms_tmatrix_axialsym_fill(
455 qpms_tmatrix_t *t, ///< T-matrix whose contents are to be replaced. Not NULL.
456 complex double omega, ///< Angular frequency.
457 qpms_epsmu_t outside, ///< Optical properties of the outside medium.
458 qpms_epsmu_t inside, ///< Optical properties of the particle's material.
459 qpms_arc_function_t shape, ///< Particle surface parametrisation.
460 /** If `lMax_extend > t->bspec->lMax`, then the internal \a Q, \a R matrices will be
461 * trimmed at this larger lMax; the final T-matrix will then be trimmed
462 * according to bspec.
464 qpms_l_t lMax_extend
467 /// Creates a new T-matrix of a particle with \f$ C_\infty \f$ symmetry.
468 static inline qpms_tmatrix_t *qpms_tmatrix_axialsym(
469 const qpms_vswf_set_spec_t *bspec,
470 complex double omega, ///< Angular frequency.
471 qpms_epsmu_t outside, ///< Optical properties of the outside medium.
472 qpms_epsmu_t inside, ///< Optical properties of the particle's material.
473 qpms_arc_function_t shape, ///< Particle surface parametrisation.
474 /// Internal lMax to improve precision of the result.
475 /** If `lMax_extend > bspec->lMax`, then the internal \a Q, \a R matrices will be
476 * trimmed at this larger lMax; the final T-matrix will then be trimmed
477 * according to bspec.
479 qpms_l_t lMax_extend
481 qpms_tmatrix_t *t = qpms_tmatrix_init(bspec);
482 qpms_tmatrix_axialsym_fill(t, omega, outside, inside, shape, lMax_extend);
483 return t;
486 /// Parameter structure for qpms_tmatrix_generator_axialsym.
487 typedef struct qpms_tmatrix_generator_axialsym_param_t {
488 qpms_epsmu_generator_t outside;
489 qpms_epsmu_generator_t inside;
490 qpms_arc_function_t shape;
491 qpms_l_t lMax_extend;
492 } qpms_tmatrix_generator_axialsym_param_t;
495 /// qpms_tmatrix_axialsym for qpms_tmatrix_generator_t
496 qpms_errno_t qpms_tmatrix_generator_axialsym(qpms_tmatrix_t *t, ///< T-matrix to fill.
497 complex double omega, ///< Angular frequency.
498 const void *params ///< Parameters of type qpms_tmatrix_generator_axialsym_param_t.
502 /// Computes the (reduced) transposed R or Q matrix for axially symmetric particle (useful for debugging).
503 qpms_errno_t qpms_tmatrix_generator_axialsym_RQ_transposed_fill(complex double *target,
504 complex double omega,
505 const qpms_tmatrix_generator_axialsym_param_t *param,
506 qpms_normalisation_t norm,
507 qpms_bessel_t J
510 /// Computes the (reduced) transposed R or Q matrix for axially symmetric particle (useful mostly for debugging).
511 qpms_errno_t qpms_tmatrix_axialsym_RQ_transposed_fill(complex double *target,
512 complex double omega, qpms_epsmu_t outside, qpms_epsmu_t inside,
513 qpms_arc_function_t shape, qpms_l_t lMaxQR, qpms_normalisation_t norm,
514 qpms_bessel_t J ///< Use QPMS_BESSEL_REGULAR to calculate \f$ R^T\f$ or QPMS_HANKEL_PLUS to calculate \f$ Q^T\f$.
518 /// An "abstract" T-matrix, contains a T-matrix generator instead of actual data.
519 typedef struct qpms_tmatrix_function_t {
520 /** \brief VSWF basis specification, NOT owned by qpms_tmatrix_t by default.
522 * Usually not checked for meaningfulness by the functions (methods),
523 * so the caller should take care that \a spec->ilist does not
524 * contain any duplicities and that for each wave with order \a m
525 * there is also one with order \a −m.
527 const qpms_vswf_set_spec_t *spec;
528 const qpms_tmatrix_generator_t *gen; ///< A T-matrix generator function.
529 } qpms_tmatrix_function_t;
531 /// Convenience function to create a new T-matrix from qpms_tmatrix_function_t.
532 // FIXME the name is not very intuitive.
533 static inline qpms_tmatrix_t *qpms_tmatrix_init_from_function(qpms_tmatrix_function_t f, complex double omega) {
534 return qpms_tmatrix_init_from_generator(f.spec, *f.gen, omega);
537 /// Specifies different kinds of operations done on T-matrices
538 typedef enum {
539 QPMS_TMATRIX_OPERATION_NOOP, ///< Identity operation.
540 QPMS_TMATRIX_OPERATION_LRMATRIX, ///< General matrix transformation \f$ T' = MTM^\dagger \f$; see @ref qpms_tmatrix_operation_lrmatrix.
541 QPMS_TMATRIX_OPERATION_IROT3, ///< Single rotoreflection specified by a qpms_irot3_t.
542 QPMS_TMATRIX_OPERATION_IROT3ARR, ///< Symmetrise using an array of rotoreflection; see @ref qpms_tmatrix_operation_irot3arr.
543 QPMS_TMATRIX_OPERATION_COMPOSE_SUM, ///< Apply several transformations and sum the results, see @ref qpms_tmatrix_operation_compose_sum.
544 QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN, ///< Chain several different transformations; see @ref qpms_tmatrix_operation_compose_chain.
545 QPMS_TMATRIX_OPERATION_SCMULZ, ///< Elementwise scalar multiplication with a complex matrix; see @ref qpms_tmatrix_operation_scmulz.
546 //QPMS_TMATRIX_OPERATION_POINTGROUP, ///< TODO the new point group in pointgroup.h
547 QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE ///< Legacy qpms_finite_group_t with filled rep3d; see @ref qpms_tmatrix_operation_finite_group.
548 } qpms_tmatrix_operation_kind_t;
550 /// General matrix transformation \f[ T' = MTM^\dagger \f] spec for qpms_tmatrix_operation_t.
551 struct qpms_tmatrix_operation_lrmatrix {
552 /// Raw matrix data of \a M in row-major order.
553 /** The matrix must be taylored for the given bspec! */
554 complex double *m;
555 size_t m_size; ///< Total size of \a m matrix in terms of sizeof(complex double).
556 bool owns_m; ///< Whether \a m is owned by this;
559 struct qpms_tmatrix_operation_t; // Forward declaration for the composed operations.
561 /// Specifies a composed operation of type \f$ T' = c\sum_i f_i'(T) \f$ for qpms_tmatrix_operation_t.
562 struct qpms_tmatrix_operation_compose_sum {
563 size_t n; ///< Number of operations in ops;
564 const struct qpms_tmatrix_operation_t **ops; ///< Operations array. (Pointers array \a ops[] always owned by this.)
565 double factor; ///< Overall factor \a c.
566 /// (Optional) operations buffer into which elements of \a ops point.
567 /** Owned by this or NULL. If not NULL, all \a ops members are assumed to point into
568 * the memory held by \a opmem and to be properly initialised.
569 * Each \a ops member has to point to _different_ elements of \a opmem.
571 struct qpms_tmatrix_operation_t *opmem;
574 /// 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.
575 struct qpms_tmatrix_operation_compose_chain {
576 size_t n; ///< Number of operations in ops;
577 const struct qpms_tmatrix_operation_t **ops; ///< Operations array. (Pointers owned by this.)
578 struct qpms_tmatrix_operation_t *opmem; ///< (Optional) operations buffer into which elements of \a ops point. (Owned by this or NULL.)
579 size_t opmem_size; ///< Length of the opmem array.
580 _Bool *ops_owned; ///< True for all sub operations owned by this and saved in opmem. NULL if opmem is NULL.
583 /// Specifies an elementwise complex multiplication of type \f$ T'_{ij} = M_{ij}T_{ij} \f$ for qpms_tmatrix_operation_t.
584 struct qpms_tmatrix_operation_scmulz {
585 /// Raw matrix data of \a M in row-major order.
586 /** The matrix must be taylored for the given bspec! */
587 complex double *m;
588 size_t m_size; ///< Total size of \a m matrix in terms of sizeof(complex double).
589 bool owns_m; ///< Whether \a m is owned by this.
592 /// Specifies a symmetrisation using a set of rotoreflections (with equal weights) for qpms_tmatrix_operation_t.
593 /** Internally, this is evaluated using a call to qpms_symmetrise_tmdata_irot3arr()
594 * or qpms_symmetrise_tmdata_irot3arr_inplace(). */
595 struct qpms_tmatrix_operation_irot3arr {
596 size_t n; ///< Number of rotoreflections;
597 qpms_irot3_t *ops; ///< Rotoreflection array of size \a n.
598 bool owns_ops; ///< Whether \a ops array is owned by this.
601 /// A generic T-matrix transformation operator.
602 typedef struct qpms_tmatrix_operation_t {
603 /// Specifies the operation kind to be performed and which type \op actually contains.
604 qpms_tmatrix_operation_kind_t typ;
605 union {
606 struct qpms_tmatrix_operation_lrmatrix lrmatrix;
607 struct qpms_tmatrix_operation_scmulz scmulz;
608 qpms_irot3_t irot3; ///< Single rotoreflection; \a typ = QPMS_TMATRIX_OPERATION_IROT3
609 struct qpms_tmatrix_operation_irot3arr irot3arr;
610 struct qpms_tmatrix_operation_compose_sum compose_sum;
611 struct qpms_tmatrix_operation_compose_chain compose_chain;
612 /// Finite group for QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE.
613 /** Not owned by this; \a rep3d must be filled. */
614 const qpms_finite_group_t *finite_group;
615 } op; ///< Operation data; actual type is determined by \a typ.
616 } qpms_tmatrix_operation_t;
618 static const qpms_tmatrix_operation_t qpms_tmatrix_operation_noop = {.typ = QPMS_TMATRIX_OPERATION_NOOP};
620 /// Apply an operation on a T-matrix, returning a newly allocated result.
621 qpms_tmatrix_t *qpms_tmatrix_apply_operation(const qpms_tmatrix_operation_t *op, const qpms_tmatrix_t *orig);
623 /// Apply an operation on a T-matrix and replace it with the result.
624 qpms_tmatrix_t *qpms_tmatrix_apply_operation_inplace(const qpms_tmatrix_operation_t *op, qpms_tmatrix_t *orig);
626 /// (Recursively) deallocates internal data of qpms_tmatrix_operation_t.
627 /** It does NOT clear the memory pointed to by op it self, but only
628 * heap-allocated data of certain operations, if labeled as owned by it.
629 * In case of compose operations, the recursion stops if the children are
630 * not owned by them (the opmem pointer is NULL).
632 void qpms_tmatrix_operation_clear(qpms_tmatrix_operation_t *f);
634 /// (Recursively) copies an qpms_tmatrix_operation_t.
635 /** Makes copies of all the internal data and takes ownership over them if needed */
636 void qpms_tmatrix_operation_copy(qpms_tmatrix_operation_t *target, const qpms_tmatrix_operation_t *src);
638 /// Inits a new "chain" of composed operations, some of which might be owned.
639 void qpms_tmatrix_operation_compose_chain_init(
640 qpms_tmatrix_operation_t *target, ///< The operation structure that will be set to the chain.
641 size_t nops, ///< Number of chained operations (length of the \a ops array)
642 size_t opmem_size ///< Size of the own operations buffer (length of the \a opmem array)
645 #if 0
646 // Abstract types that describe T-matrix/particle/scatsystem symmetries
647 // To be implemented later. See also the thoughts in the beginning of groups.h.
649 typedef *char qpms_tmatrix_id_t; ///< Maybe I want some usual integer type instead.
651 ///Abstract T-matrix type draft.
653 * TODO.
655 typedef struct qpms_abstract_tmatrix_t{
656 qpms_tmatrix_id_t id;
657 /// Generators of the discrete point group under which T-matrix is invariant.
658 qpms_irot3_t *invar_gens;
659 /// Length of invar_gens.
660 qpms_gmi_t invar_gens_size;
662 } qpms_abstract_tmatrix_t;
665 typedef struct qpms_abstract_particle_t{
666 } qpms_abstract_particle_t;
668 /// An abstract particle, defined by its position and abstract T-matrix.
669 typedef struct qpms_abstract_particle_t {
670 cart3_t pos; ///< Particle position in cartesian coordinates.
671 const qpms_abstract_tmatrix_t *tmatrix; ///< T-matrix; not owned by this.
672 } qpms_abstract_particle_t;
675 /** This is just an alias, as the same index can be used for
676 * abstract T-matrices as well.
678 typedef qpms_particle_tid_t qpms_abstract_particle_tid_t;
679 #endif // 0
681 #endif //TMATRICES_H