3 ctypedef double complex cdouble
5 from libc.stdint cimport *
7 cdef extern from "gsl/gsl_errno.h":
8 ctypedef void gsl_error_handler_t (const char *reason, const char *file,
9 int line, int gsl_errno)
10 gsl_error_handler_t *gsl_set_error_handler(gsl_error_handler_t *new_handler)
11 gsl_error_handler_t *gsl_set_error_handler_off();
13 cdef extern from "gsl/gsl_const_mksa.h":
14 const double GSL_CONST_MKSA_SPEED_OF_LIGHT
16 cdef extern from "qpms_types.h":
36 cdef struct csphvec_t:
43 cdef union anycoord_point_t:
48 ctypedef enum qpms_normalisation_t:
49 QPMS_NORMALISATION_UNDEF
50 QPMS_NORMALISATION_INVERSE
51 QPMS_NORMALISATION_REVERSE_AZIMUTHAL_PHASE
52 QPMS_NORMALISATION_SPHARM_REAL
53 QPMS_NORMALISATION_CSPHASE
54 QPMS_NORMALISATION_M_I
55 QPMS_NORMALISATION_M_MINUS
56 QPMS_NORMALISATION_N_I
57 QPMS_NORMALISATION_N_MINUS
58 QPMS_NORMALISATION_L_I
59 QPMS_NORMALISATION_L_MINUS
60 QPMS_NORMALISATION_NORM_BITSTART
61 QPMS_NORMALISATION_NORM_POWER
62 QPMS_NORMALISATION_NORM_SPHARM
63 QPMS_NORMALISATION_NORM_NONE
64 QPMS_NORMALISATION_NORM_BITS
65 QPMS_NORMALISATION_CONVENTION_KRISTENSSON_REAL
66 QPMS_NORMALISATION_CONVENTION_KRISTENSSON
67 QPMS_NORMALISATION_CONVENTION_SCUFF
68 QPMS_NORMALISATION_DEFAULT
69 ctypedef enum qpms_bessel_t:
75 ctypedef int qpms_lm_t
78 ctypedef size_t qpms_y_t
90 ctypedef np.ulonglong_t qpms_uvswfi_t
91 struct qpms_vswf_set_spec_t:
99 qpms_normalisation_t norm
100 ctypedef enum qpms_errno_t:
104 ctypedef enum qpms_vswf_type_t:
107 QPMS_VSWF_LONGITUDINAL
108 ctypedef int32_t qpms_ss_tmgi_t
109 ctypedef int32_t qpms_ss_tmi_t
110 ctypedef int32_t qpms_ss_pi_t
111 ctypedef int qpms_gmi_t
112 ctypedef int qpms_iri_t
113 qpms_iri_t QPMS_NO_IRREP
114 ctypedef const char * qpms_permutation_t
115 struct qpms_tmatrix_t:
116 qpms_vswf_set_spec_t *spec
118 bint owns_m # FIXME in fact bool
119 ctypedef enum qpms_pointgroup_class:
141 struct qpms_pointgroup_t:
142 qpms_pointgroup_class c
144 qpms_irot3_t orientation
148 ctypedef enum qpms_coord_system_t:
150 # maybe more if needed
152 cdef extern from "qpms_error.h":
153 ctypedef enum qpms_dbgmsg_flags:
156 QPMS_DBGMSG_INTEGRATION
157 qpms_dbgmsg_flags qpms_dbgmsg_enable(qpms_dbgmsg_flags types)
158 qpms_dbgmsg_flags qpms_dbgmsg_disable(qpms_dbgmsg_flags types)
160 cdef extern from "tolerances.h":
161 struct qpms_tolerance_spec_t:
163 const qpms_tolerance_spec_t QPMS_TOLERANCE_DEFAULT
166 # This is due to the fact that cython apparently cannot nest the unnamed struct/unions in an obvious way
167 ctypedef union qpms_incfield_planewave_params_k:
170 ctypedef union qpms_incfield_planewave_params_E:
174 cdef extern from "vswf.h":
175 bint qpms_vswf_set_spec_isidentical(const qpms_vswf_set_spec_t *a, const qpms_vswf_set_spec_t *b)
176 ctypedef qpms_errno_t (*qpms_incfield_t)(cdouble *target, const qpms_vswf_set_spec_t *bspec,
177 const cart3_t evalpoint, const void *args, bint add)
178 ctypedef struct qpms_incfield_planewave_params_t:
180 qpms_incfield_planewave_params_k k
181 qpms_incfield_planewave_params_E E
182 qpms_errno_t qpms_incfield_planewave(cdouble *target, const qpms_vswf_set_spec_t *bspec,
183 const cart3_t evalpoint, const void *args, bint add)
184 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)
185 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)
187 cdef extern from "indexing.h":
188 qpms_y_t qpms_lMax2nelem(qpms_l_t lMax)
189 qpms_uvswfi_t qpms_tmn2uvswfi(qpms_vswf_type_t t, qpms_m_t m, qpms_l_t n)
190 qpms_errno_t qpms_uvswfi2tmn(qpms_uvswfi_t u, qpms_vswf_type_t* t, qpms_m_t* m, qpms_l_t* n)
191 qpms_m_t qpms_uvswfi2m(qpms_uvswfi_t u)
192 # maybe more if needed
194 # Point generators from lattices.h
195 cdef extern from "lattices.h":
196 ctypedef enum PGenPointFlags:
198 struct PGenReturnData:
200 struct PGenZReturnData:
202 struct PGenPolReturnData:
204 struct PGenSphReturnData:
206 struct PGenCart2ReturnData:
208 struct PGenCart3ReturnData:
210 struct PGenClassInfo: # maybe important
212 struct PGen: # probably important
215 void PGen_destroy(PGen *g)
217 int l2d_reciprocalBasis2pi(cart2_t b1, cart2_t b2, cart2_t *rb1, cart2_t *rb2);
218 double l2d_unitcell_area(cart2_t b1, cart2_t b2)
219 void l2d_reduceBasis(cart2_t in1, cart2_t in2, cart2_t *out1, cart2_t *out2)
221 const double BASIS_RTOL
222 size_t qpms_emptylattice2_modes_maxfreq(double **target_freqs, cart2_t b1_rec, cart2_t b2_rec,
223 double rtol, cart2_t k, double wave_speed, double maxomega)
224 void qpms_emptylattice2_modes_nearest(double *target_freqs, cart2_t b1_rec, cart2_t b2_rec,
225 double rtol, cart2_t k, double wave_speed, double omega)
227 # now the individual PGen implementations:
228 # FIXME Is bint always guaranteed to be equivalent to _Bool? (I dont't think so.)
229 PGen PGen_xyWeb_new(cart2_t b1, cart2_t b2, double rtol, cart2_t offset,
230 double minR, bint inc_minR, double maxR, bint inc_maxR)
231 ctypedef enum PGen_1D_incrementDirection:
232 PGEN_1D_INC_FROM_ORIGIN
233 PGEN_1D_INC_TOWARDS_ORIGIN
234 PGen PGen_1D_new_minMaxR(double period, double offset, double minR, bint inc_minR,
235 double maxR, bint inc_maxR, PGen_1D_incrementDirection incdir)
236 int qpms_reduce_lattice_basis(double *b, size_t bsize, size_t ndim, double delta)
237 ctypedef enum LatticeDimensionality:
253 cdef extern from "vectors.h":
254 cart2_t cart2_add(const cart2_t a, const cart2_t b)
255 cart2_t cart2_substract(const cart2_t a, const cart2_t b)
256 cart2_t cart2_scale(const double c, const cart2_t v)
257 double cart2_dot(const cart2_t a, const cart2_t b)
258 double cart2_normsq(const cart2_t a)
259 double cart2norm(const cart2_t v)
260 pol_t cart2pol(const cart2_t cart)
261 sph_t pol2sph_equator(const pol_t pol)
262 sph_t cart22sph(const cart2_t cart)
263 sph_t cart12sph_zaxis(double z)
264 cart3_t cart12cart3z(double z)
265 cart3_t cart22cart3xy(const cart2_t a)
266 cart2_t cart3xy2cart2(const cart3_t a)
267 double cart3_dot(const cart3_t a, const cart3_t b)
268 cart3_t cart3_vectorprod(const cart3_t a, const cart3_t b)
269 double cart3_tripleprod(const cart3_t a, const cart3_t b, const cart3_t c)
270 double cart3_normsq(const cart3_t a)
271 double cart3norm(const cart3_t v)
272 sph_t cart2sph(const cart3_t cart)
273 cart3_t sph2cart(const sph_t sph)
274 cart2_t pol2cart(const pol_t pol)
275 cart3_t pol2cart3_equator(const pol_t pol)
276 cart3_t cart3_add(const cart3_t a, const cart3_t b)
277 cart3_t cart3_substract(const cart3_t a, const cart3_t b)
278 cart3_t cart3_scale(const double c, const cart3_t v)
279 cart3_t cart3_divscale( const cart3_t v, const double c)
280 double cart3_dist(const cart3_t a, const cart3_t b)
281 bint cart3_isclose(const cart3_t a, const cart3_t b, double rtol, double atol)
282 ccart3_t ccart3_scale(const cdouble c, const ccart3_t v)
283 ccart3_t ccart3_add(const ccart3_t a, const ccart3_t b)
284 ccart3_t ccart3_substract(const ccart3_t a, const ccart3_t b)
285 cdouble ccart3_dotnc(const ccart3_t a, const ccart3_t b)
286 ccart3_t cart32ccart3(cart3_t c)
287 csphvec_t csphvec_add(const csphvec_t a, const csphvec_t b)
288 csphvec_t csphvec_substract(const csphvec_t a, const csphvec_t b)
289 csphvec_t csphvec_scale(cdouble c, const csphvec_t v)
290 cdouble csphvec_dotnc(const csphvec_t a, const csphvec_t b)
291 sph_t sph_scale(double c, const sph_t s)
292 csph_t sph_cscale(cdouble c, const sph_t s)
293 ccart3_t csphvec2ccart(const csphvec_t sphvec, const sph_t at)
294 ccart3_t csphvec2ccart_csph(const csphvec_t sphvec, const csph_t at)
295 csphvec_t ccart2csphvec(const ccart3_t cartvec, const sph_t at)
296 csph_t sph2csph(sph_t s)
297 sph_t csph2sph(csph_t s)
298 csph_t ccart2csph(const ccart3_t cart)
299 csph_t cart2csph(const cart3_t cart)
300 ccart3_t csph2ccart(const csph_t sph)
301 void ccart3_kahaninit(ccart3_t *sum, ccart3_t *compensation)
302 void csphvec_kahaninit(csphvec_t *sum, csphvec_t *compensation)
303 void ccart3_kahanadd(ccart3_t *sum, ccart3_t *compensation, const ccart3_t input)
304 void csphvec_kahanadd(csphvec_t *sum, csphvec_t *compensation, const csphvec_t input)
305 double csphvec_norm(const csphvec_t a)
306 double csphvec_reldiff_abstol(const csphvec_t a, const csphvec_t b, double tolerance)
307 double csphvec_reldiff(const csphvec_t a, const csphvec_t b)
308 sph_t anycoord2sph(anycoord_point_t p, qpms_coord_system_t t)
309 cart3_t anycoord2cart3(anycoord_point_t p, qpms_coord_system_t t)
310 double anycoord_norm(anycoord_point_t p, qpms_coord_system_t t)
311 pol_t anycoord2pol(anycoord_point_t p, qpms_coord_system_t t)
312 cart2_t anycoord2cart2(anycoord_point_t p, qpms_coord_system_t t)
313 double anycoord2cart1(anycoord_point_t p, qpms_coord_system_t t)
314 void qpms_array_coord_transform(void *dest, qpms_coord_system_t tdest,
315 const void *src, qpms_coord_system_t tsrc, size_t nmemb)
316 void anycoord_arr2something(void *dest, qpms_coord_system_t tdest,
317 const void *src, qpms_coord_system_t tsrc, size_t nmemb)
318 void cart3_to_double_array(double *a, cart3_t b)
319 cart3_t cart3_from_double_array(const double *a)
320 void cart2_to_double_array(double *a, cart2_t b)
321 cart2_t cart2_from_double_array(const double *a)
322 const cart2_t CART2_ZERO
325 cdef extern from "quaternions.h":
326 qpms_quat_t qpms_quat_2c_from_4d(qpms_quat4d_t q)
327 qpms_quat4d_t qpms_quat_4d_from_2c(qpms_quat_t q)
328 qpms_quat_t qpms_quat_mult(qpms_quat_t p, qpms_quat_t q)
329 qpms_quat_t qpms_quat_add(qpms_quat_t p, qpms_quat_t q)
330 qpms_quat_t qpms_quat_rscale(double s, qpms_quat_t q)
331 qpms_quat_t qpms_quat_conj(qpms_quat_t q)
332 double qpms_quat_norm(qpms_quat_t q)
333 double qpms_quat_imnorm(qpms_quat_t q)
334 qpms_quat_t qpms_quat_normalise(qpms_quat_t q)
335 qpms_quat_t qpms_quat_log(qpms_quat_t q)
336 qpms_quat_t qpms_quat_exp(qpms_quat_t q)
337 qpms_quat_t qpms_quat_pow(qpms_quat_t q, double exponent)
338 cdouble qpms_wignerD_elem(qpms_quat_t q, qpms_l_t l,
339 qpms_m_t mp, qpms_m_t m)
340 qpms_irot3_t qpms_irot3_mult(qpms_irot3_t p, qpms_irot3_t q)
341 qpms_irot3_t qpms_irot3_pow(qpms_irot3_t p, int n)
342 qpms_quat_t qpms_quat_from_rotvector(cart3_t v)
344 cdef extern from "groups.h":
345 struct qpms_finite_group_irrep_t:
349 struct qpms_finite_group_t:
357 qpms_permutation_t *permrep
362 qpms_finite_group_irrep_t *irreps
363 qpms_finite_group_t QPMS_FINITE_GROUP_TRIVIAL
364 qpms_finite_group_t QPMS_FINITE_GROUP_TRIVIAL_G
366 cdef extern from "symmetries.h":
367 cdouble *qpms_zflip_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec)
368 cdouble *qpms_yflip_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec)
369 cdouble *qpms_xflip_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec)
370 cdouble *qpms_zrot_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec, double phi)
371 cdouble *qpms_zrot_rational_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec, int N, int w)
372 cdouble *qpms_irot3_uvswfi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec, qpms_irot3_t transf)
373 size_t qpms_zero_roundoff_clean(double *arr, size_t nmemb, double atol)
374 size_t qpms_czero_roundoff_clean(cdouble *arr, size_t nmemb, double atol)
376 #cdef extern from "numpy/arrayobject.h":
377 # cdef enum NPY_TYPES:
379 # NPY_CDOUBLE # complex double
381 # ctypedef int npy_intp
384 cdef extern from "translations.h":
385 struct qpms_trans_calculator:
388 cdouble** A_multipliers
389 cdouble** B_multipliers
390 enum qpms_normalization_t:
392 qpms_trans_calculator* qpms_trans_calculator_init(int lMax, int nt) # should be qpms_normalization_t
393 void qpms_trans_calculator_free(qpms_trans_calculator* c)
394 cdouble qpms_trans_calculator_get_A_ext(const qpms_trans_calculator* c,
395 int m, int n, int mu, int nu, cdouble kdlj_r, double kdlj_th, double kdlj_phi,
396 int r_ge_d, int J) nogil
397 cdouble qpms_trans_calculator_get_B_ext(const qpms_trans_calculator* c,
398 int m, int n, int mu, int nu, cdouble kdlj_r, double kdlj_th, double kdlj_phi,
399 int r_ge_d, int J) nogil
400 int qpms_trans_calculator_get_AB_p_ext(const qpms_trans_calculator* c,
401 cdouble *Adest, cdouble *Bdest,
402 int m, int n, int mu, int nu, cdouble kdlj_r, double kdlj_th, double kdlj_phi,
403 int r_ge_d, int J) nogil
404 int qpms_trans_calculator_get_AB_arrays_ext(const qpms_trans_calculator *c,
405 cdouble *Adest, cdouble *Bdest,
406 size_t deststride, size_t srcstride,
407 cdouble kdlj_r, double kdlj_theta, double kdlj_phi,
408 int r_ge_d, int J) nogil
409 int qpms_cython_trans_calculator_get_AB_arrays_loop(qpms_trans_calculator *c,
411 int daxis, int saxis,
412 char *A_data, np.npy_intp *A_shape, np.npy_intp *A_strides,
413 char *B_data, np.npy_intp *B_shape, np.npy_intp *B_strides,
414 char *r_data, np.npy_intp *r_shape, np.npy_intp *r_strides,
415 char *theta_data, np.npy_intp *theta_shape, np.npy_intp *theta_strides,
416 char *phi_data, np.npy_intp *phi_shape, np.npy_intp *phi_strides,
417 char *r_ge_d_data, np.npy_intp *phi_shape, np.npy_intp *phi_strides) nogil
419 int qpms_trans_calculator_get_trans_array(const qpms_trans_calculator *c,
421 const qpms_vswf_set_spec_t *destspec, size_t deststride,
422 const qpms_vswf_set_spec_t *srcspec, size_t srcstride,
423 csph_t kdlj, bint r_ge_d, qpms_bessel_t J);
425 int qpms_trans_calculator_get_trans_array_lc3p(
426 const qpms_trans_calculator *c,
428 const qpms_vswf_set_spec_t *destspec, size_t deststride,
429 const qpms_vswf_set_spec_t *srcspec, size_t srcstride,
430 cdouble k, cart3_t destpos, cart3_t srcpos,
434 cdef extern from "qpms_specfunc.h":
440 qpms_pitau_t qpms_pitau_get(double theta, qpms_l_t lMax, double csphase)
441 qpms_errno_t qpms_pitau_fill(double *leg, double *pi, double *tau, double theta, qpms_l_t lMax, double csphase)
442 void qpms_pitau_free(qpms_pitau_t pitau)
445 cdef extern from "gsl/gsl_interp.h":
446 struct gsl_interp_type:
448 const gsl_interp_type *gsl_interp_linear
449 const gsl_interp_type *gsl_interp_cspline
450 # ^^^ These are probably the only relevant ones.
452 cdef extern from "materials.h":
453 struct qpms_epsmu_generator_t:
454 qpms_epsmu_t (*function)(cdouble omega, const void *params)
456 qpms_epsmu_t qpms_epsmu_const_g(cdouble omega, const void *params)
457 qpms_epsmu_t qpms_permittivity_interpolator_epsmu_g(cdouble omega, const void *epsmu)
458 qpms_epsmu_t qpms_lorentzdrude_epsmu_g(cdouble omega, const void *ldparams)
460 struct qpms_permittivity_interpolator_t:
462 qpms_permittivity_interpolator_t *qpms_permittivity_interpolator_create(const size_t incount,
463 cdouble *wavelength_m, cdouble *n, cdouble *k, const gsl_interp_type *iptype)
464 qpms_permittivity_interpolator_t *qpms_permittivity_interpolator_from_yml(const char *path,
465 const gsl_interp_type *iptype)
466 cdouble qpms_permittivity_interpolator_eps_at_omega(const qpms_permittivity_interpolator_t *interp, double omega_SI)
467 double qpms_permittivity_interpolator_omega_max(const qpms_permittivity_interpolator_t *interp)
468 double qpms_permittivity_interpolator_omega_min(const qpms_permittivity_interpolator_t *interp)
469 void qpms_permittivity_interpolator_free(qpms_permittivity_interpolator_t *interp)
470 struct qpms_ldparams_triple_t:
474 struct qpms_ldparams_t:
478 qpms_ldparams_triple_t data[0]
479 cdouble qpms_lorentzdrude_eps(cdouble, const qpms_ldparams_t *)
481 cdef extern from "tmatrices.h":
482 bint qpms_load_scuff_tmatrix_crash_on_failure
483 struct qpms_tmatrix_generator_t:
484 qpms_errno_t (*function)(qpms_tmatrix_t *t, cdouble omega, const void *params)
486 qpms_errno_t qpms_tmatrix_generator_axialsym(qpms_tmatrix_t *t, cdouble omega, const void *params)
487 qpms_errno_t qpms_tmatrix_generator_interpolator(qpms_tmatrix_t *t, cdouble omega, const void *params)
488 qpms_errno_t qpms_tmatrix_generator_sphere(qpms_tmatrix_t *t, cdouble omega, const void *params)
489 qpms_errno_t qpms_tmatrix_generator_constant(qpms_tmatrix_t *t, cdouble omega, const void *params)
490 struct qpms_tmatrix_generator_sphere_param_t:
491 qpms_epsmu_generator_t outside
492 qpms_epsmu_generator_t inside
494 struct qpms_arc_function_retval_t:
497 struct qpms_arc_function_t:
498 qpms_arc_function_retval_t (*function)(double theta, const void *params)
500 struct qpms_tmatrix_generator_axialsym_param_t:
501 qpms_epsmu_generator_t outside
502 qpms_epsmu_generator_t inside
503 qpms_arc_function_t shape
505 struct qpms_arc_cylinder_params_t:
508 qpms_arc_function_retval_t qpms_arc_cylinder(double theta, const void *params)
509 qpms_arc_function_retval_t qpms_arc_sphere(double theta, const void *params)
510 struct qpms_tmatrix_interpolator_t:
511 const qpms_vswf_set_spec_t *bspec
512 void qpms_tmatrix_interpolator_free(qpms_tmatrix_interpolator_t *interp)
513 qpms_tmatrix_t *qpms_tmatrix_interpolator_eval(const qpms_tmatrix_interpolator_t *interp, double freq)
514 qpms_tmatrix_interpolator_t *qpms_tmatrix_interpolator_create(size_t n, double *freqs,
515 const qpms_tmatrix_t *tmatrices_array, const gsl_interp_type *iptype)
516 void qpms_tmatrix_free(qpms_tmatrix_t *tmatrix)
517 qpms_tmatrix_isclose(const qpms_tmatrix_t *A, const qpms_tmatrix_t *B,
518 const double rtol, const double atol)
519 qpms_errno_t qpms_symmetrise_tmdata_irot3arr(
520 cdouble *tmdata, const size_t tmcount,
521 const qpms_vswf_set_spec_t *bspec,
523 const qpms_irot3_t *symops
525 qpms_errno_t qpms_symmetrise_tmdata_finite_group(
526 cdouble *tmdata, const size_t tmcount,
527 const qpms_vswf_set_spec_t *bspec,
528 const qpms_finite_group_t *pointgroup
530 qpms_tmatrix_t *qpms_tmatrix_symmetrise_irot3arr_inplace(
533 const qpms_irot3_t *symops
535 qpms_tmatrix_t *qpms_tmatrix_symmetrise_finite_group_inplace(
537 const qpms_finite_group_t *pointgroup
539 qpms_errno_t qpms_load_scuff_tmatrix(const char *path, const qpms_vswf_set_spec_t *bspec,
540 size_t *n, double **freqs, double **freqs_su, qpms_tmatrix_t **tmatrices_array,
542 cdouble *qpms_mie_coefficients_reflection(cdouble *target, const qpms_vswf_set_spec_t *bspec,
543 double a, cdouble k_i, cdouble k_e, cdouble mu_i, cdouble mu_e, qpms_bessel_t J_ext, qpms_bessel_t J_scat)
544 qpms_tmatrix_t *qpms_tmatrix_spherical(const qpms_vswf_set_spec_t *bspec, double a,
545 cdouble k_i, cdouble k_e, cdouble mu_i, cdouble mu_e)
546 qpms_errno_t qpms_tmatrix_spherical_fill(qpms_tmatrix_t *t, double a,
547 cdouble k_i, cdouble k_e, cdouble mu_i, cdouble mu_e)
548 qpms_tmatrix_t *qpms_tmatrix_spherical(const qpms_vswf_set_spec_t *bspec,
549 double a, cdouble k_i, cdouble k_e, cdouble mu_i, cdouble mu_e)
550 cdouble qpms_drude_epsilon(cdouble eps_inf, cdouble omega_p, cdouble gamma_p, cdouble omega)
551 qpms_errno_t qpms_tmatrix_spherical_mu0_fill(qpms_tmatrix_t *t, double a, double omega,
552 cdouble epsilon_fg, cdouble epsilon_bg)
553 qpms_tmatrix_t *qpms_tmatrix_spherical_mu0(const qpms_vswf_set_spec_t *bspec, double a,
554 double omega, cdouble epsilon_fg, cdouble epsilon_bg)
555 qpms_errno_t qpms_tmatrix_generator_axialsym_RQ_transposed_fill(cdouble *target, cdouble omega,
556 const qpms_tmatrix_generator_axialsym_param_t *param, qpms_normalisation_t norm, qpms_bessel_t J)
557 struct qpms_tmatrix_function_t:
558 const qpms_vswf_set_spec_t *spec
559 const qpms_tmatrix_generator_t *gen
560 ctypedef enum qpms_tmatrix_operation_kind_t:
561 QPMS_TMATRIX_OPERATION_NOOP
562 QPMS_TMATRIX_OPERATION_LRMATRIX
563 QPMS_TMATRIX_OPERATION_IROT3
564 QPMS_TMATRIX_OPERATION_IROT3ARR
565 QPMS_TMATRIX_OPERATION_COMPOSE_SUM
566 QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN
567 QPMS_TMATRIX_OPERATION_SCMULZ
568 QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE
570 struct qpms_tmatrix_operation_t:
571 qpms_tmatrix_operation_kind_t typ
572 pass # TODO add the op union later if needed
573 const qpms_tmatrix_operation_t qpms_tmatrix_operation_noop
574 void qpms_tmatrix_operation_clear(qpms_tmatrix_operation_t *)
576 cdef extern from "pointgroups.h":
577 bint qpms_pg_is_finite_axial(qpms_pointgroup_class cls)
578 double qpms_pg_quat_cmp_atol
579 int qpms_pg_irot3_cmp(const qpms_irot3_t *, const qpms_irot3_t *);
580 int qpms_pg_irot3_cmp_v(const void *, const void *);
581 int qpms_pg_irot3_approx_cmp(const qpms_irot3_t *a, const qpms_irot3_t *b, double atol)
582 int qpms_pg_irot3_approx_cmp_v(const void *a, const void *b)
584 qpms_gmi_t qpms_pg_order(qpms_pointgroup_class cls, qpms_gmi_t n)
585 qpms_irot3_t *qpms_pg_canonical_elems( qpms_irot3_t *target, qpms_pointgroup_class cls, qpms_gmi_t n)
586 qpms_gmi_t qpms_pg_genset_size(qpms_pointgroup_class cls, qpms_gmi_t n)
587 qpms_gmi_t qpms_pg_genset(qpms_pointgroup_class cls, qpms_gmi_t n, qpms_irot3_t *gen)
588 qpms_irot3_t *qpms_pg_elems(qpms_irot3_t *target, qpms_pointgroup_t g)
589 bint qpms_pg_is_subgroup(qpms_pointgroup_t a, qpms_pointgroup_t b);
591 cdef extern from "scatsystem.h":
592 void qpms_scatsystem_set_nthreads(long n)
593 struct qpms_particle_t:
595 const qpms_tmatrix_function_t *tmg
596 qpms_tmatrix_operation_t op
597 struct qpms_particle_tid_t:
599 qpms_ss_tmi_t tmatrix_id
600 struct qpms_ss_derived_tmatrix_t:
602 qpms_tmatrix_operation_t op
603 struct qpms_scatsys_periodic_info_t:
604 cart3_t lattice_basis[3]
605 double unitcell_volume
608 struct qpms_scatsys_t:
609 int lattice_dimension
610 qpms_epsmu_generator_t medium
611 qpms_tmatrix_function_t *tmg
612 qpms_ss_tmgi_t tmg_count
613 qpms_ss_derived_tmatrix_t *tm
614 qpms_ss_tmi_t tm_count
615 qpms_particle_tid_t *p
617 # We shouldn't need more to construct a symmetric scatsystem ^^^
620 const qpms_finite_group_t *sym
621 qpms_scatsys_periodic_info_t per
623 # per[] and other stuff not currently needed in cython
624 void qpms_scatsys_free(qpms_scatsys_t *s)
625 qpms_errno_t qpms_scatsys_dump(qpms_scatsys_t *ss, char *path) #NI
626 qpms_scatsys_t *qpms_scatsys_load(char *path) #NI
627 struct qpms_scatsys_at_omega_t:
628 const qpms_scatsys_t *ss
633 qpms_scatsys_at_omega_t *qpms_scatsys_apply_symmetry(const qpms_scatsys_t *orig, const qpms_finite_group_t *sym,
634 cdouble omega, const qpms_tolerance_spec_t *tol)
635 qpms_scatsys_at_omega_t *qpms_scatsys_at_omega(const qpms_scatsys_t *ss, cdouble omega)
636 void qpms_scatsys_at_omega_free(qpms_scatsys_at_omega_t *ssw)
637 cdouble *qpms_scatsys_irrep_pack_matrix(cdouble *target_packed,
638 const cdouble *orig_full, const qpms_scatsys_t *ss, qpms_iri_t iri)
639 cdouble *qpms_scatsys_irrep_unpack_matrix(cdouble *target_full,
640 const cdouble *orig_packed, const qpms_scatsys_t *ss, qpms_iri_t iri, bint add)
641 cdouble *qpms_scatsys_irrep_pack_vector(cdouble *target_packed,
642 const cdouble *orig_full, const qpms_scatsys_t *ss, qpms_iri_t iri)
643 cdouble *qpms_scatsys_irrep_unpack_vector(cdouble *target_full,
644 const cdouble *orig_packed, const qpms_scatsys_t *ss, qpms_iri_t iri, bint add)
645 cdouble *qpms_scatsysw_build_modeproblem_matrix_full(cdouble *target,
646 const qpms_scatsys_at_omega_t *ssw)
647 cdouble *qpms_scatsys_build_translation_matrix_full(cdouble *target,
648 const qpms_scatsys_t *ss, cdouble k)
649 cdouble *qpms_scatsys_build_translation_matrix_e_full(cdouble *target,
650 const qpms_scatsys_t *ss, cdouble k, qpms_bessel_t J)
651 cdouble *qpms_scatsysw_build_modeproblem_matrix_irrep_packed(cdouble *target,
652 const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri) nogil
653 cdouble *qpms_scatsys_build_translation_matrix_e_irrep_packed(cdouble *target,
654 const qpms_scatsys_t *ss, qpms_iri_t iri, cdouble k, qpms_bessel_t J) nogil
655 cdouble *qpms_scatsysw_build_modeproblem_matrix_irrep_packed_orbitorderR(
656 cdouble *target, const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri) nogil
657 cdouble *qpms_scatsysw_build_modeproblem_matrix_irrep_packed_serial(
658 cdouble *target, const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri) nogil
659 cdouble *qpms_scatsys_incident_field_vector_full(cdouble *target_full,
660 const qpms_scatsys_t *ss, qpms_incfield_t field_at_point,
661 const void *args, bint add)
662 cdouble *qpms_scatsysw_apply_Tmatrices_full(cdouble *target_full, const cdouble *inc_full,
663 const qpms_scatsys_at_omega_t *ssw)
665 const qpms_scatsys_at_omega_t *ssw
666 const qpms_scatsys_at_omega_k_t *sswk
671 void qpms_ss_LU_free(qpms_ss_LU lu)
672 qpms_ss_LU qpms_scatsysw_build_modeproblem_matrix_full_LU(cdouble *target,
673 int *target_piv, const qpms_scatsys_at_omega_t *ssw)
674 qpms_ss_LU qpms_scatsysw_build_modeproblem_matrix_irrep_packed_LU(cdouble *target,
675 int *target_piv, const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri)
676 qpms_ss_LU qpms_scatsysw_modeproblem_matrix_full_factorise(cdouble *modeproblem_matrix_full,
677 int *target_piv, const qpms_scatsys_at_omega_t *ssw)
678 qpms_ss_LU qpms_scatsysw_modeproblem_matrix_irrep_packed_factorise(cdouble *modeproblem_matrix_full,
679 int *target_piv, const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri)
680 cdouble *qpms_scatsys_scatter_solve(cdouble *target_f, const cdouble *a_inc, qpms_ss_LU ludata)
681 const qpms_vswf_set_spec_t *qpms_ss_bspec_tmi(const qpms_scatsys_t *ss, qpms_ss_tmi_t tmi)
682 const qpms_vswf_set_spec_t *qpms_ss_bspec_pi(const qpms_scatsys_t *ss, qpms_ss_pi_t pi)
683 beyn_result_t *qpms_scatsys_finite_find_eigenmodes(const qpms_scatsys_t *ss, qpms_iri_t iri,
684 cdouble omega_centre, double omega_rr, double omega_ri, size_t contour_npoints,
685 double rank_tol, size_t rank_sel_min, double res_tol)
686 # periodic-related funs
687 struct qpms_scatsys_at_omega_k_t:
688 const qpms_scatsys_at_omega_t *ssw
691 cdouble *qpms_scatsyswk_build_modeproblem_matrix_full(cdouble *target, const qpms_scatsys_at_omega_k_t *sswk)
692 cdouble *qpms_scatsys_periodic_build_translation_matrix_full(cdouble *target, const qpms_scatsys_t *ss, cdouble wavenumber, const cart3_t *wavevector, double eta)
693 qpms_ss_LU qpms_scatsyswk_build_modeproblem_matrix_full_LU(cdouble *target, int *target_piv, const qpms_scatsys_at_omega_k_t *sswk)
694 beyn_result_t *qpms_scatsys_periodic_find_eigenmodes(const qpms_scatsys_t *ss, const double *k,
695 cdouble omega_centre, double omega_rr, double omega_ri, size_t contour_npoints,
696 double rank_tol, size_t rank_sel_min, double res_tol)
697 const qpms_vswf_set_spec_t *qpms_ss_bspec_pi(const qpms_scatsys_t *ss, qpms_ss_pi_t pi)
698 ccart3_t qpms_scatsys_scattered_E(const qpms_scatsys_t *ss, qpms_bessel_t btyp, cdouble wavenumber,
699 const cdouble *f_excitation_vector_full, cart3_t where)
700 ccart3_t qpms_scatsysw_scattered_E(const qpms_scatsys_at_omega_t *ssw, qpms_bessel_t btyp,
701 const cdouble *f_excitation_vector_full, cart3_t where)
702 ccart3_t qpms_scatsys_scattered_E__alt(const qpms_scatsys_t *ss, qpms_bessel_t btyp, cdouble wavenumber,
703 const cdouble *f_excitation_vector_full, cart3_t where)
704 ccart3_t qpms_scatsysw_scattered_E__alt(const qpms_scatsys_at_omega_t *ssw, qpms_bessel_t btyp,
705 const cdouble *f_excitation_vector_full, cart3_t where)
706 double qpms_ss_adjusted_eta(const qpms_scatsys_t *ss, cdouble wavenumber, const double *wavevector);
709 cdef extern from "ewald.h":
710 struct qpms_csf_result:
714 ctypedef enum qpms_ewald_part:
715 QPMS_EWALD_LONG_RANGE
716 QPMS_EWALD_SHORT_RANGE
720 struct qpms_ewald3_constants_t:
723 qpms_ewald3_constants_t *qpms_ewald3_constants_init(qpms_l_t lMax, int csphase)
724 void qpms_ewald3_constants_free(qpms_ewald3_constants_t *c)
726 cdouble lilgamma(double t)
727 cdouble clilgamma(cdouble z)
728 int cx_gamma_inc_series_e(double a, cdouble x, qpms_csf_result *result)
729 int cx_gamma_inc_CF_e(double a, cdouble x, qpms_csf_result *result)
730 int complex_gamma_inc_e(double a, cdouble x, int m, qpms_csf_result *result)
732 int ewald3_sigma0(cdouble *target, double *err, const qpms_ewald3_constants_t *c, double eta, cdouble wavenumber)
733 int ewald3_sigma_long(cdouble *target_sigmalr_y, double *target_sigmalr_y_err, const qpms_ewald3_constants_t *c,
734 double eta, cdouble wavenumber, double unitcell_volume, LatticeDimensionality latdim, PGen *pgen_K,
735 bint pgen_generates_shifted_points, cart3_t k, cart3_t particle_shift)
736 int ewald3_sigma_short(cdouble *target_sigmasr_y, double *target_sigmasr_y_err, const qpms_ewald3_constants_t *c,
737 double eta, cdouble wavenumber, LatticeDimensionality latdim, PGen *pgen_R,
738 bint pgen_generates_shifted_points, cart3_t k, cart3_t particle_shift)
741 cdef extern from "gsl/gsl_complex.h":
742 ctypedef struct gsl_complex:
745 cdef extern from "gsl/gsl_matrix.h":
746 ctypedef struct gsl_matrix_complex:
748 ctypedef struct gsl_vector:
750 ctypedef struct gsl_vector_complex:
753 cdef extern from "beyn.h":
754 ctypedef struct beyn_contour_t:
755 bint (*inside_test)(beyn_contour_t *, cdouble z)
757 ctypedef struct beyn_result_gsl_t:
759 ctypedef struct beyn_result_t:
767 beyn_result_gsl_t *gsl
768 ctypedef enum beyn_contour_halfellipse_orientation:
769 BEYN_CONTOUR_HALFELLIPSE_RE_PLUS
770 BEYN_CONTOUR_HALFELLIPSE_IM_PLUS
771 BEYN_CONTOUR_HALFELLIPSE_RE_MINUS
772 BEYN_CONTOUR_HALFELLIPSE_IM_MINUS
774 ctypedef int (*beyn_function_M_gsl_t)(gsl_matrix_complex *target_M, cdouble z, void *params)
775 ctypedef int (*beyn_function_M_inv_Vhat_gsl_t)(gsl_matrix_complex *target, const gsl_matrix_complex *Vhat, cdouble z, void *params)
776 ctypedef int (*beyn_function_M_t)(cdouble *target_M, size_t m, cdouble z, void *params)
777 ctypedef int (*beyn_function_M_inv_Vhat_t)(cdouble *target, size_t m, size_t l, const cdouble *Vhat, cdouble z, void *params)
779 void beyn_result_gsl_free(beyn_result_gsl_t *result)
780 void beyn_result_free(beyn_result_t *result)
782 beyn_result_gsl_t *beyn_solve_gsl(size_t m, size_t l, beyn_function_M_gsl_t M,
783 beyn_function_M_inv_Vhat_gsl_t M_inv_Vhat, void *params, const beyn_contour_t *contour,
784 double rank_tol, size_t rank_min_sel, double res_tol)
786 beyn_result_t *beyn_solve(size_t m, size_t l, beyn_function_M_t M,
787 beyn_function_M_inv_Vhat_t M_inv_Vhat, void *params, const beyn_contour_t *contour,
788 double rank_tol, size_t rank_min_sel, double res_tol)
790 beyn_contour_t *beyn_contour_ellipse(cdouble centre, double halfax_re, double halfax_im, size_t npoints)
791 beyn_contour_t *beyn_contour_halfellipse(cdouble centre, double halfax_re, double halfax_im, size_t npoints,
792 beyn_contour_halfellipse_orientation ori)
793 beyn_contour_t *beyn_contour_kidney(cdouble centre, double halfax_re, double halfax_im, size_t npoints,
794 double rounding, beyn_contour_halfellipse_orientation ori)
797 cdouble gsl_comlpex_tostd(gsl_complex z)
798 gsl_complex gsl_complex_fromstd(cdouble z)