12 #include <gsl/gsl_sf_legendre.h>
19 complex
double *target,
52 for(
size_t i = 0; i < bspec->n; ++i)
53 if (bspec->ilist[i] == index)
85 const complex
double *coeffs,
121 complex
double *target,
156 qpms_errno_t qpms_legendre_deriv_y_get(
double **result,
double **result_deriv,
double x,
qpms_l_t lMax,
157 gsl_sf_legendre_t lnorm,
double csphase);
158 qpms_errno_t qpms_legendre_deriv_y_fill(
double *where,
double *where_deriv,
double x,
159 qpms_l_t lMax, gsl_sf_legendre_t lnorm,
double csphase);
230 complex
double *targt_longcoeff, complex
double *target_mgcoeff, complex
double *target_elcoeff,
233 complex
double *targt_longcoeff, complex
double *target_mgcoeff, complex
double *target_elcoeff,
238 complex
double *longcoeffs, complex
double *mgcoeffs, complex
double *elcoeffs,
242 complex
double *longcoeffs, complex
double *mgcoeffs, complex
double *elcoeffs,
qpms_normalisation_t
Vector spherical wavefuction normalisation and phase convention codes.
Definition: qpms_types.h:104
int qpms_l_t
Type for spherical harmonic degree l.
Definition: qpms_types.h:27
unsigned long long qpms_uvswfi_t
Exhaustive index type for VSWF basis functions.
Definition: qpms_types.h:81
qpms_errno_t
Error codes / return values for certain numerical functions.
Definition: qpms_types.h:85
qpms_bessel_t
Bessel function kinds.
Definition: qpms_types.h:176
qpms_lm_t qpms_m_t
Type for spherical harmonic order m.
Definition: qpms_types.h:31
3D cartesian coordinates. See also vectors.h.
Definition: qpms_types.h:186
3D complex (actually 6D) coordinates. See also vectors.h.
Definition: qpms_types.h:191
Spherical coordinates with complex radial component. See also vectors.h.
Definition: qpms_types.h:212
3D complex vector components in local spherical basis. See also vectors.h.
Definition: qpms_types.h:218
Parameter structure for qpms_incfield_planewave()
Definition: vswf.h:93
union qpms_incfield_planewave_params_t::@1 k
Wave vector.
union qpms_incfield_planewave_params_t::@2 E
Electric field amplitude at origin.
bool use_cartesian
If true, wave direction k and amplitude E are specified in cartesian coordinates (via k....
Definition: vswf.h:94
Specifies a finite set of VSWFs.
Definition: qpms_types.h:290
Set of electric and magnetic VSWF values in spherical coordinate basis.
Definition: vswf.h:148
Spherical coordinates. See also vectors.h.
Definition: qpms_types.h:207
csphvec_t qpms_vswf_single_el_csph(qpms_m_t m, qpms_l_t n, csph_t kdlj, qpms_bessel_t btyp, qpms_normalisation_t norm)
Electric wave N, complex wave number version.
Definition: vswf.c:144
qpms_errno_t qpms_vswf_fill_csph(csphvec_t *resultL, csphvec_t *resultM, csphvec_t *resultN, qpms_l_t lMax, csph_t kdrj, qpms_bessel_t btyp, qpms_normalisation_t norm)
Evaluate VSWFs at a given point from l = 1 up to a given degree lMax (complex kr version).
Definition: vswf.c:230
bool qpms_vswf_set_spec_isidentical(const qpms_vswf_set_spec_t *a, const qpms_vswf_set_spec_t *b)
Compares two vswf basis specs.
Definition: vswf.c:55
qpms_errno_t qpms_vecspharm_fill(csphvec_t *const a1target, csphvec_t *const a2target, csphvec_t *const a3target, qpms_l_t lMax, sph_t dir, qpms_normalisation_t norm)
Evaluate vector spherical harmonics at dir.
Definition: vswf.c:348
csphvec_t qpms_eval_uvswf(const qpms_vswf_set_spec_t *setspec, const complex double *coeffs, csph_t kr, qpms_bessel_t btyp)
Evaluates field specified by SVWF coefficients at a given point.
Definition: vswf.c:611
static ssize_t qpms_vswf_set_spec_find_uvswfi(const qpms_vswf_set_spec_t *bspec, const qpms_uvswfi_t index)
Finds the position of a given index in the bspec's ilist.
Definition: vswf.h:50
qpms_errno_t qpms_vswf_set_spec_append(qpms_vswf_set_spec_t *self, qpms_uvswfi_t u)
Appends a VSWF index to a qpms_vswf_set_spec_t, also updating metadata.
Definition: vswf.c:24
struct qpms_vswfset_sph_t qpms_vswfset_sph_t
Set of electric and magnetic VSWF values in spherical coordinate basis.
void qpms_vswf_set_spec_free(qpms_vswf_set_spec_t *)
Destroys a qpms_vswf_set_spec_t.
Definition: vswf.c:94
qpms_errno_t qpms_vswf_fill(csphvec_t *resultL, csphvec_t *resultM, csphvec_t *resultN, qpms_l_t lMax, sph_t kdrj, qpms_bessel_t btyp, qpms_normalisation_t norm)
Evaluate VSWFs at a given point from l = 1 up to a given degree lMax.
Definition: vswf.c:290
qpms_vswf_set_spec_t * qpms_vswf_set_spec_copy(const qpms_vswf_set_spec_t *orig)
Copies an instance of qpms_vswf_set_spec_t.
Definition: vswf.c:66
csphvec_t qpms_vswf_single_el(qpms_m_t m, qpms_l_t n, sph_t kdlj, qpms_bessel_t btyp, qpms_normalisation_t norm)
Electric wave N.
Definition: vswf.c:192
csphvec_t qpms_vswf_L00(csph_t kdrj, qpms_bessel_t btyp, qpms_normalisation_t norm)
Evaluate the zeroth-degree longitudinal VSWF .
Definition: vswf.c:220
qpms_errno_t qpms_uvswf_fill(csphvec_t *const target, const qpms_vswf_set_spec_t *setspec, csph_t kr, qpms_bessel_t btyp)
Evaluates a set of VSWF basis functions at a given point.
Definition: vswf.c:569
csphvec_t qpms_vswf_single_mg_csph(qpms_m_t m, qpms_l_t n, csph_t kdlj, qpms_bessel_t btyp, qpms_normalisation_t norm)
Magnetic wave M, complex wave number version..
Definition: vswf.c:169
qpms_vswf_set_spec_t * qpms_vswf_set_spec_init(void)
Creates a qpms_vswf_set_spec_t structure with an empty list of wave indices.
Definition: vswf.c:14
qpms_errno_t(* qpms_incfield_t)(complex double *target, const qpms_vswf_set_spec_t *bspec, const cart3_t evalpoint, const void *args, bool add)
Calculates the (regular VSWF) expansion coefficients of an external incident field.
Definition: vswf.h:17
qpms_vswf_set_spec_t * qpms_vswf_set_spec_from_lMax(qpms_l_t lMax, qpms_normalisation_t norm)
Creates an instance of qpms_vswf_set_spec_t in the 'traditional' layout.
Definition: vswf.c:76
size_t * qpms_vswf_set_reindex(const qpms_vswf_set_spec_t *small, const qpms_vswf_set_spec_t *big)
Creates an index mapping between two bspecs.
Definition: vswf.c:111
struct qpms_incfield_planewave_params_t qpms_incfield_planewave_params_t
Parameter structure for qpms_incfield_planewave()
csphvec_t qpms_vswf_single_mg(qpms_m_t m, qpms_l_t n, sph_t kdlj, qpms_bessel_t btyp, qpms_normalisation_t norm)
Magnetic wave M.
Definition: vswf.c:197
qpms_errno_t qpms_incfield_planewave(complex double *target, const qpms_vswf_set_spec_t *bspec, const cart3_t evalpoint, const void *args, bool add)
Calculates the (regular VSWF) expansion coefficients of a plane wave.
Definition: vswf.c:480