QPMS
Electromagnetic multiple scattering library and toolkit.
|
T-matrices for scattering systems. More...
Go to the source code of this file.
Data Structures | |
struct | qpms_tmatrix_generator_t |
Generic T-matrix generator function that fills a pre-initialised qpms_tmatrix_t given a frequency. More... | |
struct | qpms_tmatrix_interpolator_t |
struct | qpms_tmatrix_generator_sphere_param_t |
Parameter structure for qpms_tmatrix_generator_sphere(). More... | |
struct | qpms_arc_function_retval_t |
Return value type for qpms_arc_function_t. More... | |
struct | qpms_arc_function_t |
Prototype for general parametrisation of \( C_\infty \)-symmetric particle's surface. More... | |
struct | qpms_arc_cylinder_params_t |
Parameter structure for qpms_arc_cylinder(). More... | |
struct | qpms_tmatrix_generator_axialsym_param_t |
Parameter structure for qpms_tmatrix_generator_axialsym. More... | |
struct | qpms_tmatrix_function_t |
An "abstract" T-matrix, contains a T-matrix generator instead of actual data. More... | |
struct | qpms_tmatrix_operation_lrmatrix |
General matrix transformation. More... | |
struct | qpms_tmatrix_operation_compose_sum |
Specifies a composed operation of type \( T' = c\sum_i f_i'(T) \) for qpms_tmatrix_operation_t. More... | |
struct | qpms_tmatrix_operation_compose_chain |
Specifies a composed operation of type \( T' = f_{n-1}(f_{n-2}(\dots f_0(T)\dots))) \) for qpms_tmatrix_operation_t. More... | |
struct | qpms_tmatrix_operation_scmulz |
Specifies an elementwise complex multiplication of type \( T'_{ij} = M_{ij}T_{ij} \) for qpms_tmatrix_operation_t. More... | |
struct | qpms_tmatrix_operation_irot3arr |
Specifies a symmetrisation using a set of rotoreflections (with equal weights) for qpms_tmatrix_operation_t. More... | |
struct | qpms_tmatrix_operation_t |
A generic T-matrix transformation operator. More... | |
Typedefs | |
typedef struct qpms_finite_group_t | qpms_finite_group_t |
typedef struct qpms_tmatrix_generator_t | qpms_tmatrix_generator_t |
Generic T-matrix generator function that fills a pre-initialised qpms_tmatrix_t given a frequency. More... | |
typedef struct qpms_tmatrix_interpolator_t | qpms_tmatrix_interpolator_t |
typedef struct qpms_tmatrix_generator_sphere_param_t | qpms_tmatrix_generator_sphere_param_t |
Parameter structure for qpms_tmatrix_generator_sphere(). | |
typedef struct qpms_arc_function_retval_t | qpms_arc_function_retval_t |
Return value type for qpms_arc_function_t. | |
typedef struct qpms_arc_function_t | qpms_arc_function_t |
Prototype for general parametrisation of \( C_\infty \)-symmetric particle's surface. | |
typedef struct qpms_arc_cylinder_params_t | qpms_arc_cylinder_params_t |
Parameter structure for qpms_arc_cylinder(). | |
typedef struct qpms_tmatrix_generator_axialsym_param_t | qpms_tmatrix_generator_axialsym_param_t |
Parameter structure for qpms_tmatrix_generator_axialsym. | |
typedef struct qpms_tmatrix_function_t | qpms_tmatrix_function_t |
An "abstract" T-matrix, contains a T-matrix generator instead of actual data. | |
typedef struct qpms_tmatrix_operation_t | qpms_tmatrix_operation_t |
A generic T-matrix transformation operator. | |
Enumerations | |
enum | qpms_tmatrix_operation_kind_t { QPMS_TMATRIX_OPERATION_NOOP , QPMS_TMATRIX_OPERATION_LRMATRIX , QPMS_TMATRIX_OPERATION_IROT3 , QPMS_TMATRIX_OPERATION_IROT3ARR , QPMS_TMATRIX_OPERATION_COMPOSE_SUM , QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN , QPMS_TMATRIX_OPERATION_SCMULZ , QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE } |
Specifies different kinds of operations done on T-matrices. More... | |
Functions | |
static complex double * | qpms_tmatrix_row (qpms_tmatrix_t *t, size_t rowno) |
Returns a pointer to the beginning of the T-matrix row rowno. | |
qpms_tmatrix_t * | qpms_tmatrix_init (const qpms_vswf_set_spec_t *bspec) |
Initialises a zero T-matrix. More... | |
qpms_tmatrix_t * | qpms_tmatrix_copy (const qpms_tmatrix_t *T) |
Copies a T-matrix, allocating new array for the T-matrix data. | |
qpms_tmatrix_t * | qpms_tmatrix_mv (qpms_tmatrix_t *dest, const qpms_tmatrix_t *orig) |
Copies a T-matrix contents between two pre-allocated T-matrices. More... | |
void | qpms_tmatrix_free (qpms_tmatrix_t *t) |
Destroys a T-matrix. | |
bool | qpms_tmatrix_isclose (const qpms_tmatrix_t *T1, const qpms_tmatrix_t *T2, const double rtol, const double atol) |
Check T-matrix equality/similarity. More... | |
qpms_tmatrix_t * | qpms_tmatrix_apply_symop (const qpms_tmatrix_t *T, const complex double *M) |
Creates a T-matrix from another matrix and a symmetry operation. More... | |
qpms_tmatrix_t * | qpms_tmatrix_apply_symop_inplace (qpms_tmatrix_t *T, const complex double *M) |
Applies a symmetry operation onto a T-matrix, rewriting the original T-matrix data. More... | |
qpms_tmatrix_t * | qpms_tmatrix_symmetrise_involution (const qpms_tmatrix_t *T, const complex double *M) |
Symmetrizes a T-matrix with an involution symmetry operation. More... | |
qpms_tmatrix_t * | qpms_tmatrix_symmetrise_C_inf (const qpms_tmatrix_t *T) |
Creates a \( C_\infty \) -symmetrized version of a T-matrix. More... | |
qpms_tmatrix_t * | qpms_tmatrix_symmetrise_C_N (const qpms_tmatrix_t *T, int N) |
Creates a \( C_N \) -symmetrized version of a T-matrix. More... | |
qpms_tmatrix_t * | qpms_tmatrix_symmetrise_involution_inplace (qpms_tmatrix_t *T, const complex double *M) |
Symmetrizes a T-matrix with an involution symmetry operation, rewriting the original one. More... | |
qpms_tmatrix_t * | qpms_tmatrix_symmetrise_C_inf_inplace (qpms_tmatrix_t *T) |
Creates a \( C_\infty \) -symmetrized version of a T-matrix, rewriting the original one. More... | |
qpms_tmatrix_t * | qpms_tmatrix_symmetrise_C_N_inplace (qpms_tmatrix_t *T, int N) |
Creates a \( C_N \) -symmetrized version of a T-matrix, rewriting the original one. More... | |
qpms_errno_t | qpms_load_scuff_tmatrix (const char *path, const qpms_vswf_set_spec_t *bspec, size_t *n, double **freqs, double **freqs_su, qpms_tmatrix_t **tmatrices_array, complex double **tmdata) |
Reads an open scuff-tmatrix generated file. More... | |
qpms_errno_t | qpms_read_scuff_tmatrix (FILE *f, const qpms_vswf_set_spec_t *bspec, size_t *n, double **freqs, double **freqs_su, qpms_tmatrix_t **tmatrices_array, complex double **tmdata) |
Loads a scuff-tmatrix generated file. More... | |
qpms_errno_t | qpms_symmetrise_tmdata_irot3arr (complex double *tmdata, const size_t tmcount, const qpms_vswf_set_spec_t *bspec, size_t n_symops, const qpms_irot3_t *symops) |
In-place application of point group elements on raw T-matrix data. More... | |
qpms_errno_t | qpms_symmetrise_tmdata_finite_group (complex double *tmdata, const size_t tmcount, const qpms_vswf_set_spec_t *bspec, const qpms_finite_group_t *pointgroup) |
In-place application of a point group on raw T-matrix data. More... | |
qpms_tmatrix_t * | qpms_tmatrix_symmetrise_irot3arr_inplace (qpms_tmatrix_t *T, size_t n_symops, const qpms_irot3_t *symops) |
In-place application of point group elements on a T-matrix. More... | |
qpms_tmatrix_t * | qpms_tmatrix_symmetrise_finite_group_inplace (qpms_tmatrix_t *T, const qpms_finite_group_t *pointgroup) |
In-place application of point group elements on a T-matrix. More... | |
complex double * | qpms_apply_tmatrix (complex double *f_target, const complex double *a, const qpms_tmatrix_t *T) |
Application of T-matrix on a vector of incident field coefficients, \( f = Ta \). More... | |
qpms_tmatrix_t * | qpms_tmatrix_init_from_generator (const qpms_vswf_set_spec_t *bspec, qpms_tmatrix_generator_t gen, complex double omega) |
Initialises and evaluates a new T-matrix using a generator. | |
qpms_errno_t | qpms_tmatrix_generator_constant (qpms_tmatrix_t *t, complex double omega, const void *tmatrix_orig) |
Implementation of qpms_matrix_generator_t that just copies a constant matrix. More... | |
void | qpms_tmatrix_interpolator_free (qpms_tmatrix_interpolator_t *interp) |
Free a T-matrix interpolator. | |
qpms_errno_t | qpms_tmatrix_interpolator_eval_fill (qpms_tmatrix_t *target, const qpms_tmatrix_interpolator_t *interp, double freq) |
Fills an existing T-matrix with new interpolated values. More... | |
qpms_tmatrix_t * | qpms_tmatrix_interpolator_eval (const qpms_tmatrix_interpolator_t *interp, double freq) |
Evaluate a T-matrix interpolated value. More... | |
qpms_tmatrix_interpolator_t * | qpms_tmatrix_interpolator_create (size_t n, const double *freqs, const qpms_tmatrix_t *tmatrices_array, const gsl_interp_type *iptype) |
Create a T-matrix interpolator from frequency and T-matrix arrays. More... | |
qpms_errno_t | qpms_tmatrix_generator_interpolator (qpms_tmatrix_t *t, complex double omega, const void *interpolator) |
qpms_tmatrix_interpolator for qpms_tmatrix_generator_t. More... | |
complex double * | qpms_mie_coefficients_reflection (complex double *target, const qpms_vswf_set_spec_t *bspec, double a, complex double k_i, complex double k_e, complex double mu_i, complex double mu_e, qpms_bessel_t J_ext, qpms_bessel_t J_scat) |
Calculates the reflection Mie-Lorentz coefficients for a spherical particle. More... | |
qpms_errno_t | qpms_tmatrix_spherical_fill (qpms_tmatrix_t *t, double a, complex double k_i, complex double k_e, complex double mu_i, complex double mu_e) |
Replaces the contents of an existing T-matrix with that of a spherical nanoparticle calculated using the Lorentz-mie theory. More... | |
qpms_errno_t | qpms_tmatrix_generator_sphere (qpms_tmatrix_t *t, complex double omega, const void *params) |
T-matrix generator for spherical particles using Lorentz-Mie solution. More... | |
static qpms_tmatrix_t * | qpms_tmatrix_spherical (const qpms_vswf_set_spec_t *bspec, double a, complex double k_i, complex double k_e, complex double mu_i, complex double mu_e) |
Creates a new T-matrix of a spherical particle using the Lorentz-Mie theory. More... | |
qpms_errno_t | qpms_tmatrix_spherical_mu0_fill (qpms_tmatrix_t *t, double a, double omega, complex double epsilon_fg, complex double epsilon_bg) |
Convenience function to calculate T-matrix of a non-magnetic spherical particle using the permittivity values, replacing existing T-matrix data. More... | |
static qpms_tmatrix_t * | qpms_tmatrix_spherical_mu0 (const qpms_vswf_set_spec_t *bspec, double a, double omega, complex double epsilon_fg, complex double epsilon_bg) |
Convenience function to calculate T-matrix of a non-magnetic spherical particle using the permittivity values. More... | |
qpms_arc_function_retval_t | qpms_arc_cylinder (double theta, const void *params) |
Arc parametrisation of cylindrical particle; for qpms_arc_function_t. More... | |
qpms_arc_function_retval_t | qpms_arc_sphere (double theta, const void *R) |
Arc parametrisation of spherical particle; for qpms_arc_function_t. More... | |
qpms_errno_t | qpms_tmatrix_axialsym_fill (qpms_tmatrix_t *t, complex double omega, qpms_epsmu_t outside, qpms_epsmu_t inside, qpms_arc_function_t shape, qpms_l_t lMax_extend) |
Replaces T-matrix contents with those of a particle with \( C_\infty \) symmetry. More... | |
static qpms_tmatrix_t * | qpms_tmatrix_axialsym (const qpms_vswf_set_spec_t *bspec, complex double omega, qpms_epsmu_t outside, qpms_epsmu_t inside, qpms_arc_function_t shape, qpms_l_t lMax_extend) |
Creates a new T-matrix of a particle with \( C_\infty \) symmetry. More... | |
qpms_errno_t | qpms_tmatrix_generator_axialsym (qpms_tmatrix_t *t, complex double omega, const void *params) |
qpms_tmatrix_axialsym for qpms_tmatrix_generator_t More... | |
qpms_errno_t | qpms_tmatrix_generator_axialsym_RQ_transposed_fill (complex double *target, complex double omega, const qpms_tmatrix_generator_axialsym_param_t *param, qpms_normalisation_t norm, qpms_bessel_t J) |
Computes the (reduced) transposed R or Q matrix for axially symmetric particle (useful for debugging). | |
qpms_errno_t | qpms_tmatrix_axialsym_RQ_transposed_fill (complex double *target, complex double omega, qpms_epsmu_t outside, qpms_epsmu_t inside, qpms_arc_function_t shape, qpms_l_t lMaxQR, qpms_normalisation_t norm, qpms_bessel_t J) |
Computes the (reduced) transposed R or Q matrix for axially symmetric particle (useful mostly for debugging). More... | |
static qpms_tmatrix_t * | qpms_tmatrix_init_from_function (qpms_tmatrix_function_t f, complex double omega) |
Convenience function to create a new T-matrix from qpms_tmatrix_function_t. | |
qpms_tmatrix_t * | qpms_tmatrix_apply_operation (const qpms_tmatrix_operation_t *op, const qpms_tmatrix_t *orig) |
Apply an operation on a T-matrix, returning a newly allocated result. | |
qpms_tmatrix_t * | qpms_tmatrix_apply_operation_inplace (const qpms_tmatrix_operation_t *op, qpms_tmatrix_t *orig) |
Apply an operation on a T-matrix and replace it with the result. | |
qpms_tmatrix_t * | qpms_tmatrix_apply_operation_replace (qpms_tmatrix_t *dest, const qpms_tmatrix_operation_t *op, const qpms_tmatrix_t *orig) |
Apply an operation on a T-matrix and replace another one it with the result. | |
void | qpms_tmatrix_operation_clear (qpms_tmatrix_operation_t *f) |
(Recursively) deallocates internal data of qpms_tmatrix_operation_t. More... | |
void | qpms_tmatrix_operation_copy (qpms_tmatrix_operation_t *target, const qpms_tmatrix_operation_t *src) |
(Recursively) copies an qpms_tmatrix_operation_t. More... | |
void | qpms_tmatrix_operation_compose_chain_init (qpms_tmatrix_operation_t *target, size_t nops, size_t opmem_size) |
Inits a new "chain" of composed operations, some of which might be owned. More... | |
Variables | |
bool | qpms_load_scuff_tmatrix_crash_on_failure |
Tells whether qpms_load_scuff_tmatrix should crash if fopen() fails. More... | |
static const qpms_tmatrix_operation_t | qpms_tmatrix_operation_noop = {.typ = QPMS_TMATRIX_OPERATION_NOOP} |
T-matrices for scattering systems.
typedef struct qpms_tmatrix_generator_t qpms_tmatrix_generator_t |
Generic T-matrix generator function that fills a pre-initialised qpms_tmatrix_t given a frequency.
Implemented by: qpms_tmatrix_generator_axialsym() qpms_tmatrix_generator_interpolator() qpms_tmatrix_generator_sphere() qpms_tmatrix_generator_constant()
Specifies different kinds of operations done on T-matrices.
Enumerator | |
---|---|
QPMS_TMATRIX_OPERATION_NOOP | Identity operation. |
QPMS_TMATRIX_OPERATION_LRMATRIX | General matrix transformation \( T' = MTM^\dagger \); see qpms_tmatrix_operation_lrmatrix. |
QPMS_TMATRIX_OPERATION_IROT3 | Single rotoreflection specified by a qpms_irot3_t. |
QPMS_TMATRIX_OPERATION_IROT3ARR | Symmetrise using an array of rotoreflection; see qpms_tmatrix_operation_irot3arr. |
QPMS_TMATRIX_OPERATION_COMPOSE_SUM | Apply several transformations and sum the results, see qpms_tmatrix_operation_compose_sum. |
QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN | Chain several different transformations; see qpms_tmatrix_operation_compose_chain. |
QPMS_TMATRIX_OPERATION_SCMULZ | Elementwise scalar multiplication with a complex matrix; see qpms_tmatrix_operation_scmulz. |
QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE | Legacy qpms_finite_group_t with filled rep3d; see qpms_tmatrix_operation_finite_group. |
complex double* qpms_apply_tmatrix | ( | complex double * | f_target, |
const complex double * | a, | ||
const qpms_tmatrix_t * | T | ||
) |
Application of T-matrix on a vector of incident field coefficients, \( f = Ta \).
f_target | Scattered field coefficient array of size T->spec->n; if NULL, a new one is allocated. |
a | Incident field coefficient array of size T->spec->n. |
T | T-matrix T to apply. |
qpms_arc_function_retval_t qpms_arc_cylinder | ( | double | theta, |
const void * | params | ||
) |
Arc parametrisation of cylindrical particle; for qpms_arc_function_t.
params | Points to qpms_arc_cylinder_params_t |
qpms_arc_function_retval_t qpms_arc_sphere | ( | double | theta, |
const void * | R | ||
) |
Arc parametrisation of spherical particle; for qpms_arc_function_t.
Useful mostly only for benchmarks or debugging, as one can use the Mie-Lorentz solution.
R | Points to double containing particle's radius. |
qpms_errno_t qpms_load_scuff_tmatrix | ( | const char * | path, |
const qpms_vswf_set_spec_t * | bspec, | ||
size_t * | n, | ||
double ** | freqs, | ||
double ** | freqs_su, | ||
qpms_tmatrix_t ** | tmatrices_array, | ||
complex double ** | tmdata | ||
) |
Reads an open scuff-tmatrix generated file.
*freqs, *freqs_su, *tmatrices_array and *tmdata arrays are allocated by this function and have to be freed by the caller after use. freqs_su and tmatrices_array can be NULL, in that case the respective arrays are not filled nor allocated.
The contents of tmatrices_array is NOT supposed to be freed element per element.
TODO more checks and options regarding NANs etc.
path | Path to the TMatrix file |
bspec | VSWF set spec |
n | Number of successfully loaded t-matrices |
freqs | Frequencies in SI units.. |
freqs_su | Frequencies in SCUFF units (optional). |
tmatrices_array | The resulting T-matrices (optional). |
tmdata | The T-matrices raw contents |
complex double* qpms_mie_coefficients_reflection | ( | complex double * | target, |
const qpms_vswf_set_spec_t * | bspec, | ||
double | a, | ||
complex double | k_i, | ||
complex double | k_e, | ||
complex double | mu_i, | ||
complex double | mu_e, | ||
qpms_bessel_t | J_ext, | ||
qpms_bessel_t | J_scat | ||
) |
Calculates the reflection Mie-Lorentz coefficients for a spherical particle.
This function is based on the previous python implementation mie_coefficients() from qpms_p.py, so any bugs therein should affect this function as well and perhaps vice versa.
Most importantly, the case of magnetic material, mu_i != 0 or mu_e != 0 has never been tested and might give wrong results.
TODO better doc.
target | Target array of length bspec->n. If NULL, a new one will be allocated. |
bspec | Defines which of the coefficients are calculated. |
a | Radius of the sphere. |
k_i | Wave number of the internal material of the sphere. |
k_e | Wave number of the surrounding medium. |
mu_i | Relative permeability of the sphere material. |
mu_e | Relative permeability of the surrounding medium. |
J_ext | Kind of the "incoming" waves. Most likely QPMS_BESSEL_REGULAR. |
J_scat | Kind of the "scattered" waves. Most likely QPMS_HANKEL_PLUS. |
qpms_errno_t qpms_read_scuff_tmatrix | ( | FILE * | f, |
const qpms_vswf_set_spec_t * | bs, | ||
size_t *const | n, | ||
double **const | freqs, | ||
double **const | freqs_su, | ||
qpms_tmatrix_t **const | tmatrices_array, | ||
complex double **const | tmdata | ||
) |
Loads a scuff-tmatrix generated file.
A simple wrapper over qpms_read_scuff_tmatrix() that needs a path instead of open FILE.
The T-matrix is transformed from the VSWF basis defined by QPMS_NORMALISATION_CONVENTION_SCUFF into the basis defined by convention bspec->norm.
Right now, bspec->norm with real or "reversed complex" spherical harmonics are not supported.
Loads a scuff-tmatrix generated file.
f | file handle |
bs | VSWF set spec |
n | Number of successfully loaded t-matrices |
freqs | Frequencies in SI units |
freqs_su | Frequencies in SCUFF units (optional) |
tmatrices_array | The resulting T-matrices (optional). |
tmdata | The T-matrices raw contents. The coefficient of outgoing wave defined by bspec->ilist[desti] as a result of incoming wave bspec->ilist[srci] at frequency (*freqs)[fi] is accessed via (*tmdata)[bspec->n*bspec->n*fi + desti*bspec->n + srci]. |
qpms_errno_t qpms_symmetrise_tmdata_finite_group | ( | complex double * | tmdata, |
const size_t | tmcount, | ||
const qpms_vswf_set_spec_t * | bspec, | ||
const qpms_finite_group_t * | pointgroup | ||
) |
In-place application of a point group on raw T-matrix data.
This does the same as qpms_symmetrise_tmdata_irot3arr(), but takes a valid finite point group as an argument.
TODO more doc.
qpms_errno_t qpms_symmetrise_tmdata_irot3arr | ( | complex double * | tmdata, |
const size_t | tmcount, | ||
const qpms_vswf_set_spec_t * | bspec, | ||
size_t | n_symops, | ||
const qpms_irot3_t * | symops | ||
) |
In-place application of point group elements on raw T-matrix data.
tmdata can be e.g. obtained by qpms_load_scuff_tmatrix(). The symops array should always contain all elements of a finite point (sub)group, including the identity operation.
TODO more doc.
qpms_tmatrix_t* qpms_tmatrix_apply_symop | ( | const qpms_tmatrix_t * | T, |
const complex double * | M | ||
) |
Creates a T-matrix from another matrix and a symmetry operation.
The symmetry operation is expected to be a unitary (square) matrix M with the same VSWF basis spec as the T-matrix (i.e. t->spec). The new T-matrix will then correspond to CHECKME
\[ T' = MTM^\dagger \]
T | the original T-matrix |
M | the symmetry op matrix in row-major format |
qpms_tmatrix_t* qpms_tmatrix_apply_symop_inplace | ( | qpms_tmatrix_t * | T, |
const complex double * | M | ||
) |
Applies a symmetry operation onto a T-matrix, rewriting the original T-matrix data.
The symmetry operation is expected to be a unitary (square) matrix M with the same VSWF basis spec as the T-matrix (i.e. t->spec). The new T-matrix will then correspond to CHECKME
\[ T' = MTM^\dagger \]
T | the original T-matrix |
M | the symmetry op matrix in row-major format |
|
inlinestatic |
Creates a new T-matrix of a particle with \( C_\infty \) symmetry.
omega | Angular frequency. |
outside | Optical properties of the outside medium. |
inside | Optical properties of the particle's material. |
shape | Particle surface parametrisation. |
lMax_extend | Internal lMax to improve precision of the result. If `lMax_extend > bspec->lMax`, then the internal Q, R matrices will be trimmed at this larger lMax; the final T-matrix will then be trimmed according to bspec. |
qpms_errno_t qpms_tmatrix_axialsym_fill | ( | qpms_tmatrix_t * | t, |
complex double | omega, | ||
qpms_epsmu_t | outside, | ||
qpms_epsmu_t | inside, | ||
qpms_arc_function_t | shape, | ||
qpms_l_t | lMax_extend | ||
) |
Replaces T-matrix contents with those of a particle with \( C_\infty \) symmetry.
N.B. currently, this might crash for higher values of lMax (lMax > 5). Also, it seems that I am doing something wrong, as the result is accurate for sphere with lMax = 1 and for higher l the accuracy decreases.
t | T-matrix whose contents are to be replaced. Not NULL. |
omega | Angular frequency. |
outside | Optical properties of the outside medium. |
inside | Optical properties of the particle's material. |
shape | Particle surface parametrisation. |
lMax_extend | If `lMax_extend > t->bspec->lMax`, then the internal Q, R matrices will be trimmed at this larger lMax; the final T-matrix will then be trimmed according to bspec. |
qpms_errno_t qpms_tmatrix_axialsym_RQ_transposed_fill | ( | complex double * | target, |
complex double | omega, | ||
qpms_epsmu_t | outside, | ||
qpms_epsmu_t | inside, | ||
qpms_arc_function_t | shape, | ||
qpms_l_t | lMaxQR, | ||
qpms_normalisation_t | norm, | ||
qpms_bessel_t | J | ||
) |
Computes the (reduced) transposed R or Q matrix for axially symmetric particle (useful mostly for debugging).
J | Use QPMS_BESSEL_REGULAR to calculate \( R^T\) or QPMS_HANKEL_PLUS to calculate \( Q^T\). |
qpms_errno_t qpms_tmatrix_generator_axialsym | ( | qpms_tmatrix_t * | t, |
complex double | omega, | ||
const void * | params | ||
) |
qpms_tmatrix_axialsym for qpms_tmatrix_generator_t
t | T-matrix to fill. |
omega | Angular frequency. |
params | Parameters of type qpms_tmatrix_generator_axialsym_param_t. |
qpms_errno_t qpms_tmatrix_generator_constant | ( | qpms_tmatrix_t * | t, |
complex double | omega, | ||
const void * | tmatrix_orig | ||
) |
Implementation of qpms_matrix_generator_t that just copies a constant matrix.
N.B. this does almost no checks at all, so use it only for t-matrices with the same base spec.
tmatrix_orig | Source T-matrix, real type is (const qpms_tmatrix_t*). |
qpms_errno_t qpms_tmatrix_generator_interpolator | ( | qpms_tmatrix_t * | t, |
complex double | omega, | ||
const void * | interpolator | ||
) |
qpms_tmatrix_interpolator for qpms_tmatrix_generator_t.
As in qpms_tmatrix_interpolator_eval(), the imaginary part of frequency is discarded!
t | T-matrix to fill. |
omega | Angular frequency. |
interpolator | Parameter of type qpms_tmatrix_interpolator_t *. |
qpms_errno_t qpms_tmatrix_generator_sphere | ( | qpms_tmatrix_t * | t, |
complex double | omega, | ||
const void * | params | ||
) |
T-matrix generator for spherical particles using Lorentz-Mie solution.
params | Of type qpms_tmatrix_generator_sphere_param_t. |
qpms_tmatrix_t* qpms_tmatrix_init | ( | const qpms_vswf_set_spec_t * | bspec | ) |
Initialises a zero T-matrix.
qpms_tmatrix_interpolator_t* qpms_tmatrix_interpolator_create | ( | size_t | n, |
const double * | freqs, | ||
const qpms_tmatrix_t * | tmatrices_array, | ||
const gsl_interp_type * | iptype | ||
) |
Create a T-matrix interpolator from frequency and T-matrix arrays.
n | Number of freqs and T-matrices provided. |
tmatrices_array | N.B. array of qpms_tmatrix_t, not pointers! |
qpms_tmatrix_t* qpms_tmatrix_interpolator_eval | ( | const qpms_tmatrix_interpolator_t * | interp, |
double | freq | ||
) |
Evaluate a T-matrix interpolated value.
The result is to be freed using qpms_tmatrix_free().
qpms_errno_t qpms_tmatrix_interpolator_eval_fill | ( | qpms_tmatrix_t * | target, |
const qpms_tmatrix_interpolator_t * | interp, | ||
double | freq | ||
) |
Fills an existing T-matrix with new interpolated values.
target | T-matrix to be updated, not NULL. |
bool qpms_tmatrix_isclose | ( | const qpms_tmatrix_t * | T1, |
const qpms_tmatrix_t * | T2, | ||
const double | rtol, | ||
const double | atol | ||
) |
Check T-matrix equality/similarity.
This function actually checks for identical vswf specs. TODO define constants with "default" atol, rtol for this function.
qpms_tmatrix_t* qpms_tmatrix_mv | ( | qpms_tmatrix_t * | dest, |
const qpms_tmatrix_t * | orig | ||
) |
Copies a T-matrix contents between two pre-allocated T-matrices.
orig->spec and dest->spec must be identical.
void qpms_tmatrix_operation_clear | ( | qpms_tmatrix_operation_t * | f | ) |
(Recursively) deallocates internal data of qpms_tmatrix_operation_t.
It does NOT clear the memory pointed to by op it self, but only heap-allocated data of certain operations, if labeled as owned by it. In case of compose operations, the recursion stops if the children are not owned by them (the opmem pointer is NULL).
void qpms_tmatrix_operation_compose_chain_init | ( | qpms_tmatrix_operation_t * | target, |
size_t | nops, | ||
size_t | opmem_size | ||
) |
Inits a new "chain" of composed operations, some of which might be owned.
target | The operation structure that will be set to the chain. |
nops | Number of chained operations (length of the ops array) |
opmem_size | Size of the own operations buffer (length of the opmem array) |
void qpms_tmatrix_operation_copy | ( | qpms_tmatrix_operation_t * | target, |
const qpms_tmatrix_operation_t * | src | ||
) |
(Recursively) copies an qpms_tmatrix_operation_t.
Makes copies of all the internal data and takes ownership over them if needed
|
inlinestatic |
Creates a new T-matrix of a spherical particle using the Lorentz-Mie theory.
a | Radius of the sphere. |
k_i | Wave number of the internal material of the sphere. |
k_e | Wave number of the surrounding medium. |
mu_i | Relative permeability of the sphere material. |
mu_e | Relative permeability of the surrounding medium. |
qpms_errno_t qpms_tmatrix_spherical_fill | ( | qpms_tmatrix_t * | t, |
double | a, | ||
complex double | k_i, | ||
complex double | k_e, | ||
complex double | mu_i, | ||
complex double | mu_e | ||
) |
Replaces the contents of an existing T-matrix with that of a spherical nanoparticle calculated using the Lorentz-mie theory.
t | T-matrix whose contents are to be replaced. Not NULL. |
a | Radius of the sphere. |
k_i | Wave number of the internal material of the sphere. |
k_e | Wave number of the surrounding medium. |
mu_i | Relative permeability of the sphere material. |
mu_e | Relative permeability of the surrounding medium. |
|
inlinestatic |
Convenience function to calculate T-matrix of a non-magnetic spherical particle using the permittivity values.
a | Radius of the sphere. |
omega | Angular frequency. |
epsilon_fg | Relative permittivity of the sphere. |
epsilon_bg | Relative permittivity of the background medium. |
qpms_errno_t qpms_tmatrix_spherical_mu0_fill | ( | qpms_tmatrix_t * | t, |
double | a, | ||
double | omega, | ||
complex double | epsilon_fg, | ||
complex double | epsilon_bg | ||
) |
Convenience function to calculate T-matrix of a non-magnetic spherical particle using the permittivity values, replacing existing T-matrix data.
t | T-matrix whose contents are to be replaced. Not NULL. |
a | Radius of the sphere. |
omega | Angular frequency. |
epsilon_fg | Relative permittivity of the sphere. |
epsilon_bg | Relative permittivity of the background medium. |
qpms_tmatrix_t* qpms_tmatrix_symmetrise_C_inf | ( | const qpms_tmatrix_t * | T | ) |
Creates a \( C_\infty \) -symmetrized version of a T-matrix.
\[ {T'}_{tlm}^{\tau\lambda\mu} = T_{tlm}^{\tau\lambda\mu} \delta_{m\mu} \]
T | the original T-matrix |
qpms_tmatrix_t* qpms_tmatrix_symmetrise_C_inf_inplace | ( | qpms_tmatrix_t * | T | ) |
Creates a \( C_\infty \) -symmetrized version of a T-matrix, rewriting the original one.
\[ {T'}_{tlm}^{\tau\lambda\mu} = T_{tlm}^{\tau\lambda\mu} \delta_{m\mu} \]
T | the original T-matrix |
qpms_tmatrix_t* qpms_tmatrix_symmetrise_C_N | ( | const qpms_tmatrix_t * | T, |
int | N | ||
) |
Creates a \( C_N \) -symmetrized version of a T-matrix.
\[ {T'}_{tlm}^{\tau\lambda\mu} = \begin{cases} T{}_{lm}^{\lambda\mu} & (m-\mu)\mod N=0\\ 0 & (m-\mu)\mod N\ne0 \end{cases} . \]
T | the original T-matrix |
N | number of z-axis rotations in the group |
qpms_tmatrix_t* qpms_tmatrix_symmetrise_C_N_inplace | ( | qpms_tmatrix_t * | T, |
int | N | ||
) |
Creates a \( C_N \) -symmetrized version of a T-matrix, rewriting the original one.
\[ {T'}_{tlm}^{\tau\lambda\mu} = \begin{cases} T{}_{lm}^{\lambda\mu} & (m-\mu)\mod N=0\\ 0 & (m-\mu)\mod N\ne0 \end{cases} . \]
T | the original T-matrix |
N | number of z-axis rotations in the group |
qpms_tmatrix_t* qpms_tmatrix_symmetrise_finite_group_inplace | ( | qpms_tmatrix_t * | T, |
const qpms_finite_group_t * | pointgroup | ||
) |
In-place application of point group elements on a T-matrix.
This does the same as qpms_tmatrix_symmetrise_irot3arr(), but takes a valid finite point group as an argument.
TODO more doc.
qpms_tmatrix_t* qpms_tmatrix_symmetrise_involution | ( | const qpms_tmatrix_t * | T, |
const complex double * | M | ||
) |
Symmetrizes a T-matrix with an involution symmetry operation.
The symmetry operation is expected to be a unitary (square) matrix M with the same VSWF basis spec as the T-matrix (i.e. t->spec). The new T-matrix will then correspond to CHECKME
\[ T' = \frac{T + MTM^\dagger}{2} \]
T | the original T-matrix |
M | the symmetry op matrix in row-major format |
qpms_tmatrix_t* qpms_tmatrix_symmetrise_involution_inplace | ( | qpms_tmatrix_t * | T, |
const complex double * | M | ||
) |
Symmetrizes a T-matrix with an involution symmetry operation, rewriting the original one.
The symmetry operation is expected to be a unitary (square) matrix M with the same VSWF basis spec as the T-matrix (i.e. t->spec). The new T-matrix will then correspond to CHECKME
\[ T' = \frac{T + MTM^\dagger}{2} \]
T | the original T-matrix |
M | the symmetry op matrix in row-major format |
qpms_tmatrix_t* qpms_tmatrix_symmetrise_irot3arr_inplace | ( | qpms_tmatrix_t * | T, |
size_t | n_symops, | ||
const qpms_irot3_t * | symops | ||
) |
In-place application of point group elements on a T-matrix.
The symops array should always contain all elements of a finite point (sub)group, including the identity operation.
TODO more doc.
|
extern |
Tells whether qpms_load_scuff_tmatrix should crash if fopen() fails.
If true (default), the function causes the program die e.g. when the tmatrix file does not exists.
If false, it does nothing and returns an appropriate error value instead. This is desirable e.g. when used in Python (so that proper exception can be thrown).