2 * \brief Common qpms types.
12 #define M_PI_2 (1.570796326794896619231321691639751442098584699687552910487L)
15 #define M_PI (3.14159265358979323846264338327950288419716939937510582097494L)
18 #define M_SQRT2 (1.41421356237309504880168872420969807856967187537694807317668L)
21 #define M_SQRTPI (1.77245385090551602729816748334114518279754945612238712821381L)
24 // integer index types
25 typedef int qpms_lm_t
;
26 /// Type for spherical harmonic degree \a l.
27 typedef int qpms_l_t
; /* can't be unsigned because of the behaviour under the - operator;
28 also -1 needed as an invalid value for scalar waves. */
30 /// Type for spherical harmonic order \a m.
31 typedef qpms_lm_t qpms_m_t
;
33 /// Type for the (\a l, \a m) multiindex of transversal (\a M or \a N -type) VSWFs.
34 /** This corresponds to the typical memory layout for various coefficient etc.
35 * Corresponds to the l-primary, m-secondary ordering, i.e.
36 * \f[ y = 0: l = 1, m = -1, \f]
37 * \f[ y = 1: l = 1, m = 0, \f]
38 * \f[ y = 2: l = 1, m = +1, \f]
39 * \f[ y = 3: l = 2, m = -2, \f]
42 * See also indexing.h.
44 typedef size_t qpms_y_t
;
46 /// Type for the (\a l, \a m) multiindex of spherical harmonics, including (0, 0).
47 /** This differs from qpms_y_t by being shifted by one and including
48 * the \a l = 0 option. Suitable also for scalar and longitudinal waves.
49 * Corresponds to the \a l -primary, \a m -secondary ordering, i.e.
50 * \f[ y = 0: l = 0, m = 0, \f]
51 * \f[ y = 1: l = 1, m = -1, \f]
52 * \f[ y = 2: l = 1, m = 0, \f]
53 * \f[ y = 3: l = 1, m = +1, \f]
54 * \f[ y = 4: l = 2, m = -2, \f]
57 * See also indexing.h.
59 typedef size_t qpms_y_sc_t
;
61 /// Codes of the VSWF types (electric/N, magnetic/M, longitudinal/L).
63 QPMS_VSWF_ELECTRIC
= 2, ///< "Electric" (\a N -type) transversal wave.
64 QPMS_VSWF_MAGNETIC
= 1, ///< "Magnetic" (\a M -type) transversal wave.
65 QPMS_VSWF_LONGITUDINAL
= 0 ///< Longitudinal (\a L -type) wave (not relevant for radiation).
68 // FIXME DOC The references to the functions do not work in doxygen:
69 /// Exhaustive index type for VSWF basis functions.
70 /** Carries information about the wave being of \a M/N/L (magnetic, electric,
71 * or longitudinal) type, as well as the wave's degree and order (\a l, \a m).
73 * The formula is 4 * (qpms_y_sc_t) y_sc + (qmps_vswf_type_t) type_code,
74 * but don't rely on this and use the functions
75 * \ref qpms_tmn2uvswfi() and \ref qpms_uvswfi2tmn()
76 * from qpms_types.h instead
77 * as the formula might change in future versions.
79 * See also indexing.h.
81 typedef unsigned long long qpms_uvswfi_t
;
83 /// Error codes / return values for certain numerical functions.
84 /** Those with values between -2 and 32 are subset of the GSL error codes. */
86 QPMS_SUCCESS
= 0, ///< Success.
87 QPMS_ERROR
= -1, ///< Unspecified error.
88 QPMS_FAILURE
= QPMS_ERROR
,
89 QPMS_ENOMEM
= 8, ///< Out of memory.
90 QPMS_ESING
= 21, ///< Apparent singularity detected.
91 QPMS_NAN_ENCOUNTERED
= 1024 ///< NaN value encountered in data processing.
96 /// Vector spherical wavefuction normalisation and phase convention codes.
98 * Throughout the literature, various conventions for VSWF bases are used.
99 * These bit flags are used by the functions declared in normalisation.h
100 * that return the appropriate convention-dependent factors.
102 * See @ref vswf_conventions for comparison of the various conventions used.
105 QPMS_NORMALISATION_UNDEF
= 0, ///< Convention undefined. This should not happen.
106 /// Flag indicating that qpms_normalisition_factor_* should actually return values inverse to the default.
107 QPMS_NORMALISATION_INVERSE
= 1,
108 /** Flag indicating inversion of the asimuthal phase for complex spherical harmonics (i.e. \f$ e^{-im\phi} \f$
109 * instead of \f$ e^{im\phi} \f$.
111 QPMS_NORMALISATION_REVERSE_AZIMUTHAL_PHASE
= 2,
112 /// Flag indicating use of the real spherical harmonics.
113 /** If QPMS_NORMALISATION_REVERSE_AZIMUTHAL_PHASE is unset, negative \a m
114 * correspond to sine in the asimuthal factor; if set, undefined behaviour.
116 QPMS_NORMALISATION_SPHARM_REAL
= 4,
117 /// Flag indicating usage of Condon-Shortley phase.
118 /** If set, the Ferrers functions and everything derived from them
119 * (spherical harmonics, VSWFs) will include a \f$ (-1)^m \f$ factor.
121 * On implementation level, this means that the relevant `gsl_sf_legendre_*_e()`
122 * functions will be called with argument `csphase = -1.` instead of `+1.`.
124 QPMS_NORMALISATION_CSPHASE
= 8,
125 QPMS_NORMALISATION_M_I
= 16, ///< Include an additional \a i -factor into the magnetic waves.
126 QPMS_NORMALISATION_M_MINUS
= 32, ///< Include an additional \f$-1\f$ -factor into the magnetic waves.
127 QPMS_NORMALISATION_N_I
= 64, ///< Include an additional \a i -factor into the electric waves.
128 QPMS_NORMALISATION_N_MINUS
= 128, ///< Include an additional \f$-1\f$ -factor into the magnetic waves.
129 QPMS_NORMALISATION_L_I
= 256, ///< Include an additional \a i -factor into the longitudinal waves.
130 QPMS_NORMALISATION_L_MINUS
= 512, ///< Include an additional \f$-1\f$ -factor into the longitudinal waves.
131 QPMS_NORMALISATION_NORM_BITSTART
= 65536,
132 /// The VSWFs shall be power-normalised. This is the "default".
134 * Power normalisation is used e.g. in \cite kristensson_spherical_2014 (complex spherical
135 * harmonics with Condon-Shortley phase) or \cite kristensson_scattering_2016 (real
136 * spherical harmonics). This is also the reference for all the other normalisation conventions,
137 * meaning that qpms_normalisation_factor_M() and qpms_normalisation_factor_N() shall
138 * always return `1. + 0.*I` if `norm == QPMS_NORMALISATION_NORM_POWER`.
140 QPMS_NORMALISATION_NORM_POWER
= QPMS_NORMALISATION_NORM_BITSTART
* 1,
141 /// The VSWFs shall be normalised as in \cite taylor_optical_2011 .
142 /** This includes a \f$ \sqrt{l(l+1)} \f$ factor compared to the power normalisation. */
143 QPMS_NORMALISATION_NORM_SPHARM
= QPMS_NORMALISATION_NORM_BITSTART
* 3,
144 /// The VSWFs shall be created using spherical harmonics without any normalisation. Do not use.
145 /** This includes a \f[
146 * \sqrt{l(l+1)} \left(\frac{(2l+1)}{4\pi}\frac{(l-m)!}{(l+m)!}\right)^{-\frac{1}{2}}
147 * \f] factor compared to the power normalisation.
149 * Note that this has no sense whatsoever for real spherical harmonics.
150 * Again, do not use this.
152 QPMS_NORMALISATION_NORM_NONE
= QPMS_NORMALISATION_NORM_BITSTART
* 2,
153 QPMS_NORMALISATION_NORM_BITS
= QPMS_NORMALISATION_NORM_POWER
154 | QPMS_NORMALISATION_NORM_NONE
| QPMS_NORMALISATION_NORM_SPHARM
,
156 /// VSWF convention used in \cite kristensson_scattering_2016
157 QPMS_NORMALISATION_CONVENTION_KRISTENSSON_REAL
= QPMS_NORMALISATION_NORM_POWER
158 | QPMS_NORMALISATION_SPHARM_REAL
,
159 /// VSWF convention used in \cite kristensson_spherical_2014
160 QPMS_NORMALISATION_CONVENTION_KRISTENSSON
= QPMS_NORMALISATION_NORM_POWER
161 | QPMS_NORMALISATION_CSPHASE
,
162 /// VSWF convention used in SCUFF-EM \cite reid_electromagnetism_2016
163 QPMS_NORMALISATION_CONVENTION_SCUFF
= QPMS_NORMALISATION_NORM_POWER
164 | QPMS_NORMALISATION_CSPHASE
| QPMS_NORMALISATION_M_I
165 | QPMS_NORMALISATION_N_MINUS
,
166 /// Default VSWF convention. We might encourage the compiler to expect this one.
167 QPMS_NORMALISATION_DEFAULT
= QPMS_NORMALISATION_CONVENTION_KRISTENSSON
168 } qpms_normalisation_t
;
170 /// Determine whether the convention includes Condon-Shortley phase (-1) or not (+1).
171 static inline int qpms_normalisation_t_csphase(qpms_normalisation_t norm
) {
172 return (norm
& QPMS_NORMALISATION_CSPHASE
)? -1 : 1;
175 /// Bessel function kinds.
177 QPMS_BESSEL_REGULAR
= 1, ///< regular (spherical) Bessel function \a j (Bessel function of the first kind)
178 QPMS_BESSEL_SINGULAR
= 2, ///< singular (spherical) Bessel function \a y (Bessel function of the second kind)
179 QPMS_HANKEL_PLUS
= 3, ///< (spherical) Hankel function \f$ h_1 = j + iy \f$
180 QPMS_HANKEL_MINUS
= 4, ///< (spherical) Hankel function \f$ h_2 = j - iy \f$
181 QPMS_BESSEL_UNDEF
= 0 ///< invalid / unspecified kind
184 // coordinate system types
185 /// 3D cartesian coordinates. See also vectors.h.
186 typedef struct cart3_t
{
190 /// 3D complex (actually 6D) coordinates. See also vectors.h.
191 typedef struct ccart3_t
{
192 complex double x
, y
, z
;
195 /// 3D complex vector pair (represents the E, H fields).
196 typedef struct ccart3_pair
{
200 /// 2D cartesian coordinates. See also vectors.h.
201 /** See also vectors.h */
202 typedef struct cart2_t
{
206 /// Spherical coordinates. See also vectors.h.
207 typedef struct sph_t
{
208 double r
, theta
, phi
;
211 /// Spherical coordinates with complex radial component. See also vectors.h.
212 typedef struct csph_t
{ // Do I really need this???
217 /// 3D complex vector components in local spherical basis. See also vectors.h.
218 typedef struct csphvec_t
{
219 complex double rc
, thetac
, phic
;
222 /// 2D polar coordinates. See also vectors.h.
223 typedef struct pol_t
{
227 /// Union type capable to contain various 1D, 2D and 3D coordinates.
228 /** Usually combined with qpms_coord_system_t. */
229 typedef union anycoord_point_t
{
230 double z
; ///< 1D cartesian coordinate.
237 /// Enum codes for common coordinate systems.
239 // IF EVER CHANGING THE CONSTANT VALUES HERE,
240 // CHECK THAT THEY DO NOT CLASH WITH THOSE IN PGenPointFlags!
241 QPMS_COORDS_CART1
= 64, ///< 1D cartesian (= double).
242 QPMS_COORDS_POL
= 128, ///< 2D polar.
243 QPMS_COORDS_SPH
= 256, ///< 3D spherical.
244 QPMS_COORDS_CART2
= 512, ///< 2D cartesian.
245 QPMS_COORDS_CART3
= 1024, ///< 3D cartesian.
246 /// Convenience bitmask (not a valid flag!).
247 QPMS_COORDS_BITRANGE
= QPMS_COORDS_CART1
252 } qpms_coord_system_t
;
256 * Internaly represented as a pair of complex numbers,
257 * \f$ Q_a = Q_1 + iQ_z, Q_b = Q_y + i Q_x\f$.
259 * See quaternions.h for "methods".
261 typedef struct qpms_quat_t
{
265 /// Quaternion type as four doubles.
266 /** See quaternions.h for "methods".
268 typedef struct qpms_quat4d_t
{
269 double c1
, ci
, cj
, ck
;
272 /// 3D improper rotations represented as a quaternion and a sign of the determinant.
273 /** See quaternions.h for "methods".
275 typedef struct qpms_irot3_t
{
276 qpms_quat_t rot
; ///< Quaternion representing the rotation part.
277 short det
; ///< Determinant of the transformation (valid values are 1 (rotation) or -1 (improper rotation)
280 /// Specifies a finite set of VSWFs.
282 * When for example not all the M and N -type waves up to a degree lMax
283 * need to be computed, this will specify the subset.
285 * A typical use case would be when only even/odd fields wrt. z-plane
286 * mirror symmetry are considered.
288 * See vswf.h for "methods".
290 typedef struct qpms_vswf_set_spec_t
{
291 size_t n
; ///< Actual number of VSWF indices included in ilist.
292 qpms_uvswfi_t
*ilist
; ///< List of wave indices.
293 /// Maximum degree of the waves specified in ilist overall.
294 /** `max(lMax_M, lMax_N, lMax_L)` */
296 /// Maximum degree of the magnetic (M-type) waves.
297 /** Set to 0 if no magnetic waves present. */
299 /// Maximum degree of the electric (N-type) waves.
300 /** Set to 0 if no electric waves present. */
302 /// Maximum degree of the longitudinal (L-type) waves.
303 /** Set to -1 if no longitudinal waves present. */
305 size_t capacity
; ///< Allocated capacity of ilist.
306 qpms_normalisation_t norm
; ///< Normalisation convention. To be set manually if needed.
307 } qpms_vswf_set_spec_t
;
309 /// T-matrix index used in qpms_scatsys_t and related structures.
310 typedef int32_t qpms_ss_tmi_t
;
312 /// T-matrix generator index used in qpms_scatsys_t and related structures.
313 typedef int32_t qpms_ss_tmgi_t
;
315 /// Particle index used in qpms_scatsys_t and related structures.
316 typedef int32_t qpms_ss_pi_t
;
318 // These types are mainly used in groups.h:
319 /// Finite group member index. See also groups.h.
320 typedef int qpms_gmi_t
;
322 /// Irreducible representation index. See also groups.h.
323 typedef int qpms_iri_t
;
324 /// Constant labeling that no irrep decomposition is done (full system solved instead).
325 #define QPMS_NO_IRREP ((qpms_iri_t) -1)
327 /// Permutation representation of a group element.
328 /** For now, it's just a string of the form "(0,1)(3,4,5)"
330 typedef const char * qpms_permutation_t
;
333 /** In the future, I might rather use a more abstract approach in which T-matrix
334 * is a mapping (function) of the field expansion coefficients.
335 * So the interface might change.
336 * For now, let me stick to the square dense matrix representation.
338 typedef struct qpms_tmatrix_t
{
339 /** \brief VSWF basis specification, NOT owned by qpms_tmatrix_t by default.
341 * Usually not checked for meaningfulness by the functions (methods),
342 * so the caller should take care that \a spec->ilist does not
343 * contain any duplicities and that for each wave with order \a m
344 * there is also one with order \a −m.
346 const qpms_vswf_set_spec_t
*spec
;
347 complex double *m
; ///< Matrix elements in row-major order.
348 bool owns_m
; ///< Information wheter m shall be deallocated with qpms_tmatrix_free()
352 /// Classification of a 3D point group.
356 QPMS_PGS_CN
, ///< Rotational symmetry \f$ \mathrm{C_{n}} \f$.
357 QPMS_PGS_S2N
, ///< Rotoreflectional symmetry \f$ \mathrm{S_{2n}} \f$.
358 QPMS_PGS_CNH
, ///< Rotational symmetry with horizontal reflection \f$ \mathrm{C_{nh}} \f$.
359 QPMS_PGS_CNV
, ///< Pyramidal symmetry \f$ \mathrm{C_{nv}} \f$.
360 QPMS_PGS_DN
, ///< Dihedral symmetry \f$ \mathrm{D_{n}} \f$.
361 QPMS_PGS_DND
, ///< Antiprismatic symmetry \f$ \mathrm{D_{nd}} \f$.
362 QPMS_PGS_DNH
, ///< Prismatic symmetry \f$ \mathrm{D_{nh}} \f$.
364 // Remaining polyhedral groups
365 QPMS_PGS_T
, ///< Chiral tetrahedral symmetry \f$ \mathrm{T} \f$.
366 QPMS_PGS_TD
, ///< Full tetrahedral symmetry \f$ \mathrm{T_d} \f$.
367 QPMS_PGS_TH
, ///< Pyritohedral symmetry \f$ \mathrm{T_h} \f$.
368 QPMS_PGS_O
, ///< Chiral octahedral symmetry \f$ \mathrm{O_h} \f$.
369 QPMS_PGS_OH
, ///< Full octahedral symmetry \f$ \mathrm{O_h} \f$.
370 QPMS_PGS_I
, ///< Chiral icosahedral symmetry \f$ \mathrm{I} \f$.
371 QPMS_PGS_IH
, ///< Full icosahedral symmetry \f$ \mathrm{I_h} \f$.
373 // Continuous axial groups
374 QPMS_PGS_CINF
, ///< \f$ \mathrm{C_\infty} \f$
375 QPMS_PGS_CINFH
, ///< \f$ \mathrm{C_{\infty h}} \f$
376 QPMS_PGS_CINFV
, ///< \f$ \mathrm{C_{\infty v}} \f$
377 QPMS_PGS_DINF
, ///< \f$ \mathrm{D_\infty} \f$
378 QPMS_PGS_DINFH
, ///< \f$ \mathrm{D_{\infty h}} \f$
380 // Remaining continuous groups
381 QPMS_PGS_SO3
, ///< Special orthogonal group \f$ \mathrm{SO(3)}.
382 QPMS_PGS_O3
, ///< Orthogonal group \f$ \mathrm{O(3)}.
383 } qpms_pointgroup_class
;
385 /// Full characterisation of a 3D point group.
386 typedef struct qpms_pointgroup_t
{
387 qpms_pointgroup_class c
; ///< Point group classification.
388 qpms_gmi_t n
; ///< Order of the rotational subgroup \f$ \mathrm{C_n} \f$ of finite axial groups.
389 /// Transformation between this point group and the "canonical" point group of the same type.
391 * Each 3D point group of a given type (specified by the \a c
392 * and \a n members) has its isomorphous "canonical instance",
393 * typically with the main rotation axis identical to the \a z
394 * cartesian axis and a mirror plane (if applicable) orthogonal
395 * to the \a x cartesian axis.
397 * If \f$ o \f$ is a transformation specified by \a orientation,
398 * then an element \f$ g \f$ of this group can be written
400 * g = o g_\mathrm{can.} o^{-1}
401 * \f] where \f$ g_\mathrm{can.} \f$ is a corresponding element
402 * from the canonical instance.
404 * CHECKME \f$ o \f$ transforms the cartesian \a z axis to the
405 * main rotation axis of this group (if applicable).
407 * TODO more detailed specification about the canonical instances
408 * and in which direction \a orientation goes.
410 qpms_irot3_t orientation
;
414 /// An abstract T-matrix without actual elements, but with info about particle symmetry.
415 typedef struct qpms_abstract_tmatrix_t
{
416 const qpms_vswf_set_spec_t
*spec
;
417 qpms_pointgroup_t sym
;
418 } qpms_abstract_tmatrix_t
;
421 /// A type holding electric permittivity and magnetic permeability of a material.
422 typedef struct qpms_epsmu_t
{
423 complex double eps
; ///< Relative permittivity.
424 complex double mu
; ///< Relative permeability.
427 struct qpms_tolerance_spec_t
; // See tolerances.h
429 #define lmcheck(l,m) assert((l) >= 1 && abs(m) <= (l))