Fix saving lists of arrays with recent versions of numpy
[qpms.git] / qpms / qpms_cdefs.pxd
blob1ca0476acd63add640f49d61258521d0ff3f9c2f
1 cimport numpy as np
3 ctypedef double complex cdouble
5 from libc.stdint cimport *
7 cdef extern from "version.h":
8     const char *qpms_version_git_sha1
10 cdef extern from "gsl/gsl_errno.h":
11     ctypedef void gsl_error_handler_t (const char *reason, const char *file,
12             int line, int gsl_errno)
13     gsl_error_handler_t *gsl_set_error_handler(gsl_error_handler_t *new_handler)
14     gsl_error_handler_t *gsl_set_error_handler_off();
16 cdef extern from "gsl/gsl_const_mksa.h":
17     const double GSL_CONST_MKSA_SPEED_OF_LIGHT
19 cdef extern from "qpms_types.h":
20     cdef struct cart3_t:
21         double x
22         double y
23         double z
24     cdef struct ccart3_t:
25         cdouble x
26         cdouble y
27         cdouble z
28     cdef struct cart2_t:
29         double x
30         double y
31     cdef struct sph_t:
32         double r
33         double theta
34         double phi
35     cdef struct csph_t:
36         cdouble r
37         double theta
38         double phi
39     cdef struct csphvec_t:
40         cdouble rc
41         cdouble thetac
42         cdouble phic
43     cdef struct pol_t:
44         double r
45         double phi
46     cdef union anycoord_point_t:
47         double z
48         cart3_t cart3
49         cart2_t cart2
50         pol_t pol
51     ctypedef enum qpms_normalisation_t:
52         QPMS_NORMALISATION_UNDEF
53         QPMS_NORMALISATION_INVERSE
54         QPMS_NORMALISATION_REVERSE_AZIMUTHAL_PHASE
55         QPMS_NORMALISATION_SPHARM_REAL
56         QPMS_NORMALISATION_CSPHASE
57         QPMS_NORMALISATION_M_I
58         QPMS_NORMALISATION_M_MINUS
59         QPMS_NORMALISATION_N_I
60         QPMS_NORMALISATION_N_MINUS
61         QPMS_NORMALISATION_L_I
62         QPMS_NORMALISATION_L_MINUS
63         QPMS_NORMALISATION_NORM_BITSTART
64         QPMS_NORMALISATION_NORM_POWER
65         QPMS_NORMALISATION_NORM_SPHARM
66         QPMS_NORMALISATION_NORM_NONE
67         QPMS_NORMALISATION_NORM_BITS
68         QPMS_NORMALISATION_CONVENTION_KRISTENSSON_REAL
69         QPMS_NORMALISATION_CONVENTION_KRISTENSSON
70         QPMS_NORMALISATION_CONVENTION_SCUFF
71         QPMS_NORMALISATION_DEFAULT
72     ctypedef enum qpms_bessel_t:
73         QPMS_BESSEL_REGULAR
74         QPMS_BESSEL_SINGULAR
75         QPMS_HANKEL_PLUS
76         QPMS_HANKEL_MINUS
77         QPMS_BESSEL_UNDEF
78     ctypedef int qpms_lm_t
79     ctypedef int qpms_l_t
80     ctypedef int qpms_m_t
81     ctypedef size_t qpms_y_t
82     struct qpms_quat_t:
83         cdouble a
84         cdouble b
85     struct qpms_quat4d_t:
86         double c1
87         double ci
88         double cj
89         double ck
90     struct qpms_irot3_t:
91         qpms_quat_t rot
92         short det
93     ctypedef np.ulonglong_t qpms_uvswfi_t
94     struct qpms_vswf_set_spec_t:
95         size_t n
96         qpms_uvswfi_t *ilist
97         qpms_l_t lMax
98         qpms_l_t lMax_M
99         qpms_l_t lMax_N
100         qpms_l_t lMax_L
101         size_t capacity
102         qpms_normalisation_t norm
103     ctypedef enum qpms_errno_t:
104         QPMS_SUCCESS
105         QPMS_ERROR
106         # more if needed
107     ctypedef enum qpms_vswf_type_t:
108         QPMS_VSWF_ELECTRIC
109         QPMS_VSWF_MAGNETIC
110         QPMS_VSWF_LONGITUDINAL
111     ctypedef int32_t qpms_ss_tmgi_t
112     ctypedef int32_t qpms_ss_tmi_t
113     ctypedef int32_t qpms_ss_pi_t
114     ctypedef int qpms_gmi_t
115     ctypedef int qpms_iri_t
116     qpms_iri_t QPMS_NO_IRREP
117     ctypedef const char * qpms_permutation_t
118     struct qpms_tmatrix_t:
119         qpms_vswf_set_spec_t *spec
120         cdouble *m
121         bint owns_m # FIXME in fact bool
122     ctypedef enum qpms_pointgroup_class:
123         QPMS_PGS_CN
124         QPMS_PGS_S2N
125         QPMS_PGS_CNH
126         QPMS_PGS_CNV
127         QPMS_PGS_DN
128         QPMS_PGS_DND
129         QPMS_PGS_DNH
130         QPMS_PGS_T
131         QPMS_PGS_TD
132         QPMS_PGS_TH
133         QPMS_PGS_O
134         QPMS_PGS_OH
135         QPMS_PGS_I
136         QPMS_PGS_IH
137         QPMS_PGS_CINF
138         QPMS_PGS_CINFH
139         QPMS_PGS_CINFV
140         QPMS_PGS_DINF
141         QPMS_PGS_DINFH
142         QPMS_PGS_SO3
143         QPMS_PGS_O3
144     struct qpms_pointgroup_t:
145         qpms_pointgroup_class c
146         qpms_gmi_t n
147         qpms_irot3_t orientation
148     struct qpms_epsmu_t:
149         cdouble eps
150         cdouble mu
151     ctypedef enum qpms_coord_system_t:
152         pass
153     # maybe more if needed
155 cdef extern from "qpms_error.h":
156     ctypedef enum qpms_dbgmsg_flags:
157         QPMS_DBGMSG_MISC
158         QPMS_DBGMSG_THREADS
159         QPMS_DBGMSG_INTEGRATION
160     qpms_dbgmsg_flags qpms_dbgmsg_enable(qpms_dbgmsg_flags types)
161     qpms_dbgmsg_flags qpms_dbgmsg_disable(qpms_dbgmsg_flags types)
163 cdef extern from "tolerances.h":
164     struct qpms_tolerance_spec_t:
165         pass # TODO
166     const qpms_tolerance_spec_t QPMS_TOLERANCE_DEFAULT
169 # This is due to the fact that cython apparently cannot nest the unnamed struct/unions in an obvious way
170 ctypedef union qpms_incfield_planewave_params_k:
171     ccart3_t cart
172     csph_t sph
173 ctypedef union qpms_incfield_planewave_params_E:
174     ccart3_t cart
175     csphvec_t sph
177 cdef extern from "vswf.h":
178     bint qpms_vswf_set_spec_isidentical(const qpms_vswf_set_spec_t *a, const qpms_vswf_set_spec_t *b)
179     ctypedef qpms_errno_t (*qpms_incfield_t)(cdouble *target, const qpms_vswf_set_spec_t *bspec,
180             const cart3_t evalpoint, const void *args, bint add)
181     ctypedef struct qpms_incfield_planewave_params_t:
182         bint use_cartesian
183         qpms_incfield_planewave_params_k k
184         qpms_incfield_planewave_params_E E
185     qpms_errno_t qpms_incfield_planewave(cdouble *target, const qpms_vswf_set_spec_t *bspec,
186             const cart3_t evalpoint, const void *args, bint add)
187     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)
188     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)
190 cdef extern from "indexing.h":
191     qpms_y_t qpms_lMax2nelem(qpms_l_t lMax)
192     qpms_uvswfi_t qpms_tmn2uvswfi(qpms_vswf_type_t t, qpms_m_t m, qpms_l_t n)
193     qpms_errno_t qpms_uvswfi2tmn(qpms_uvswfi_t u, qpms_vswf_type_t* t, qpms_m_t* m, qpms_l_t* n)
194     qpms_m_t qpms_uvswfi2m(qpms_uvswfi_t u)
195     # maybe more if needed
197 # Point generators from lattices.h
198 cdef extern from "lattices.h":
199     ctypedef enum PGenPointFlags:
200         pass
201     struct PGenReturnData:
202         pass
203     struct PGenZReturnData:
204         pass
205     struct PGenPolReturnData:
206         pass
207     struct PGenSphReturnData:
208         pass
209     struct PGenCart2ReturnData:
210         pass
211     struct PGenCart3ReturnData:
212         pass
213     struct PGenClassInfo: # maybe important
214         pass
215     struct PGen: # probably important
216         PGenClassInfo* c
217         void *stateData
218     void PGen_destroy(PGen *g)
219     
220     int l2d_reciprocalBasis2pi(cart2_t b1, cart2_t b2, cart2_t *rb1, cart2_t *rb2);
221     double l2d_unitcell_area(cart2_t b1, cart2_t b2)
222     void l2d_reduceBasis(cart2_t in1, cart2_t in2, cart2_t *out1, cart2_t *out2)
224     const double BASIS_RTOL
225     size_t qpms_emptylattice2_modes_maxfreq(double **target_freqs, cart2_t b1_rec, cart2_t b2_rec,
226             double rtol, cart2_t k, double wave_speed, double maxomega)
227     void qpms_emptylattice2_modes_nearest(double *target_freqs, cart2_t b1_rec, cart2_t b2_rec,
228             double rtol, cart2_t k, double wave_speed, double omega)
230     # now the individual PGen implementations:
231     # FIXME Is bint always guaranteed to be equivalent to _Bool? (I dont't think so.)
232     PGen PGen_xyWeb_new(cart2_t b1, cart2_t b2, double rtol, cart2_t offset,
233             double minR, bint inc_minR, double maxR, bint inc_maxR)
234     ctypedef enum PGen_1D_incrementDirection:
235         PGEN_1D_INC_FROM_ORIGIN
236         PGEN_1D_INC_TOWARDS_ORIGIN
237     PGen PGen_1D_new_minMaxR(double period, double offset, double minR, bint inc_minR,
238             double maxR, bint inc_maxR, PGen_1D_incrementDirection incdir)
239     int qpms_reduce_lattice_basis(double *b, size_t bsize, size_t ndim, double delta)
240     ctypedef enum LatticeDimensionality:
241         LAT1D
242         LAT2D
243         LAT3D
244         SPACE1D
245         SPACE2D
246         SPACE3D
247         LAT_1D_IN_3D
248         LAT_2D_IN_3D
249         LAT_3D_IN_3D
250         LAT_ZONLY
251         LAT_XYONLY
252         LAT_1D_IN_3D_ZONLY
253         LAT_2D_IN_3D_XYONLY
256 cdef extern from "vectors.h":
257     cart2_t cart2_add(const cart2_t a, const cart2_t b)
258     cart2_t cart2_substract(const cart2_t a, const cart2_t b)
259     cart2_t cart2_scale(const double c, const cart2_t v)
260     double cart2_dot(const cart2_t a, const cart2_t b)
261     double cart2_normsq(const cart2_t a)
262     double cart2norm(const cart2_t v)
263     pol_t cart2pol(const cart2_t cart)
264     sph_t pol2sph_equator(const pol_t pol)
265     sph_t cart22sph(const cart2_t cart)
266     sph_t cart12sph_zaxis(double z)
267     cart3_t cart12cart3z(double z)
268     cart3_t cart22cart3xy(const cart2_t a)
269     cart2_t cart3xy2cart2(const cart3_t a)
270     double cart3_dot(const cart3_t a, const cart3_t b)
271     cart3_t cart3_vectorprod(const cart3_t a, const cart3_t b)
272     double cart3_tripleprod(const cart3_t a, const cart3_t b, const cart3_t c)
273     double cart3_normsq(const cart3_t a)
274     double cart3norm(const cart3_t v)
275     sph_t cart2sph(const cart3_t cart)
276     cart3_t sph2cart(const sph_t sph)
277     cart2_t pol2cart(const pol_t pol)
278     cart3_t pol2cart3_equator(const pol_t pol)
279     cart3_t cart3_add(const cart3_t a, const cart3_t b)
280     cart3_t cart3_substract(const cart3_t a, const cart3_t b)
281     cart3_t cart3_scale(const double c, const cart3_t v)
282     cart3_t cart3_divscale( const cart3_t v, const double c)
283     double cart3_dist(const cart3_t a, const cart3_t b)
284     bint cart3_isclose(const cart3_t a, const cart3_t b, double rtol, double atol)
285     ccart3_t ccart3_scale(const cdouble c, const ccart3_t v)
286     ccart3_t ccart3_add(const ccart3_t a, const ccart3_t b)
287     ccart3_t ccart3_substract(const ccart3_t a, const ccart3_t b)
288     cdouble ccart3_dotnc(const ccart3_t a, const ccart3_t b)
289     ccart3_t cart32ccart3(cart3_t c)
290     csphvec_t csphvec_add(const csphvec_t a, const csphvec_t b)
291     csphvec_t csphvec_substract(const csphvec_t a, const csphvec_t b)
292     csphvec_t csphvec_scale(cdouble c, const csphvec_t v)
293     cdouble csphvec_dotnc(const csphvec_t a, const csphvec_t b)
294     sph_t sph_scale(double c, const sph_t s)
295     csph_t sph_cscale(cdouble c, const sph_t s)
296     ccart3_t csphvec2ccart(const csphvec_t sphvec, const sph_t at)
297     ccart3_t csphvec2ccart_csph(const csphvec_t sphvec, const csph_t at)
298     csphvec_t ccart2csphvec(const ccart3_t cartvec, const sph_t at)
299     csph_t sph2csph(sph_t s)
300     sph_t csph2sph(csph_t s)
301     csph_t ccart2csph(const ccart3_t cart)
302     csph_t cart2csph(const cart3_t cart)
303     ccart3_t csph2ccart(const csph_t sph)
304     void ccart3_kahaninit(ccart3_t *sum, ccart3_t *compensation)
305     void csphvec_kahaninit(csphvec_t *sum, csphvec_t *compensation)
306     void ccart3_kahanadd(ccart3_t *sum, ccart3_t *compensation, const ccart3_t input)
307     void csphvec_kahanadd(csphvec_t *sum, csphvec_t *compensation, const csphvec_t input)
308     double csphvec_norm(const csphvec_t a)
309     double csphvec_reldiff_abstol(const csphvec_t a, const csphvec_t b, double tolerance)
310     double csphvec_reldiff(const csphvec_t a, const csphvec_t b)
311     sph_t anycoord2sph(anycoord_point_t p, qpms_coord_system_t t)
312     cart3_t anycoord2cart3(anycoord_point_t p, qpms_coord_system_t t)
313     double anycoord_norm(anycoord_point_t p, qpms_coord_system_t t)
314     pol_t anycoord2pol(anycoord_point_t p, qpms_coord_system_t t)
315     cart2_t anycoord2cart2(anycoord_point_t p, qpms_coord_system_t t)
316     double anycoord2cart1(anycoord_point_t p, qpms_coord_system_t t)
317     void qpms_array_coord_transform(void *dest, qpms_coord_system_t tdest,
318             const void *src, qpms_coord_system_t tsrc, size_t nmemb)
319     void anycoord_arr2something(void *dest, qpms_coord_system_t tdest,
320             const void *src, qpms_coord_system_t tsrc, size_t nmemb)
321     void cart3_to_double_array(double *a, cart3_t b)
322     cart3_t cart3_from_double_array(const double *a)
323     void cart2_to_double_array(double *a, cart2_t b)
324     cart2_t cart2_from_double_array(const double *a)
325     const cart2_t CART2_ZERO
328 cdef extern from "quaternions.h":
329     qpms_quat_t qpms_quat_2c_from_4d(qpms_quat4d_t q)
330     qpms_quat4d_t qpms_quat_4d_from_2c(qpms_quat_t q)
331     qpms_quat_t qpms_quat_mult(qpms_quat_t p, qpms_quat_t q)
332     qpms_quat_t qpms_quat_add(qpms_quat_t p, qpms_quat_t q)
333     qpms_quat_t qpms_quat_rscale(double s, qpms_quat_t q) 
334     qpms_quat_t qpms_quat_conj(qpms_quat_t q) 
335     double qpms_quat_norm(qpms_quat_t q) 
336     double qpms_quat_imnorm(qpms_quat_t q)
337     qpms_quat_t qpms_quat_normalise(qpms_quat_t q) 
338     qpms_quat_t qpms_quat_log(qpms_quat_t q)
339     qpms_quat_t qpms_quat_exp(qpms_quat_t q)
340     qpms_quat_t qpms_quat_pow(qpms_quat_t q, double exponent)
341     cdouble qpms_wignerD_elem(qpms_quat_t q, qpms_l_t l,
342                            qpms_m_t mp, qpms_m_t m)
343     qpms_irot3_t qpms_irot3_mult(qpms_irot3_t p, qpms_irot3_t q)
344     qpms_irot3_t qpms_irot3_pow(qpms_irot3_t p, int n)
345     qpms_quat_t qpms_quat_from_rotvector(cart3_t v)
347 cdef extern from "groups.h":
348     struct qpms_finite_group_irrep_t:
349         int dim
350         char *name
351         cdouble *m
352     struct qpms_finite_group_t:
353         char *name
354         qpms_gmi_t order
355         qpms_gmi_t idi
356         qpms_gmi_t *mt
357         qpms_gmi_t *invi
358         qpms_gmi_t *gens
359         int ngens
360         qpms_permutation_t *permrep
361         char **elemlabels
362         int permrep_nelem
363         qpms_irot3_t *rep3d
364         qpms_iri_t nirreps
365         qpms_finite_group_irrep_t *irreps
366     qpms_finite_group_t QPMS_FINITE_GROUP_TRIVIAL
367     qpms_finite_group_t QPMS_FINITE_GROUP_TRIVIAL_G
369 cdef extern from "symmetries.h":
370     cdouble *qpms_zflip_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec)
371     cdouble *qpms_yflip_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec)
372     cdouble *qpms_xflip_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec)
373     cdouble *qpms_zrot_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec, double phi)
374     cdouble *qpms_zrot_rational_uvswi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec, int N, int w)
375     cdouble *qpms_irot3_uvswfi_dense(cdouble *target, const qpms_vswf_set_spec_t *bspec, qpms_irot3_t transf)
376     size_t qpms_zero_roundoff_clean(double *arr, size_t nmemb, double atol)
377     size_t qpms_czero_roundoff_clean(cdouble *arr, size_t nmemb, double atol)
379 #cdef extern from "numpy/arrayobject.h":
380 #    cdef enum NPY_TYPES:
381 #        NPY_DOUBLE
382 #        NPY_CDOUBLE # complex double
383 #        NPY_LONG # int
384 #    ctypedef int npy_intp
387 cdef extern from "translations.h":
388     struct qpms_trans_calculator:
389         int lMax
390         size_t nelem
391         cdouble** A_multipliers
392         cdouble** B_multipliers
393     enum qpms_normalization_t:
394         pass
395     qpms_trans_calculator* qpms_trans_calculator_init(int lMax, int nt) # should be qpms_normalization_t
396     void qpms_trans_calculator_free(qpms_trans_calculator* c)
397     cdouble qpms_trans_calculator_get_A_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     cdouble qpms_trans_calculator_get_B_ext(const qpms_trans_calculator* c,
401             int m, int n, int mu, int nu, cdouble kdlj_r, double kdlj_th, double kdlj_phi,
402             int r_ge_d, int J) nogil
403     int qpms_trans_calculator_get_AB_p_ext(const qpms_trans_calculator* c,
404             cdouble *Adest, cdouble *Bdest,
405             int m, int n, int mu, int nu, cdouble kdlj_r, double kdlj_th, double kdlj_phi,
406             int r_ge_d, int J) nogil
407     int qpms_trans_calculator_get_AB_arrays_ext(const qpms_trans_calculator *c,
408             cdouble *Adest, cdouble *Bdest,
409             size_t deststride, size_t srcstride,
410             cdouble kdlj_r, double kdlj_theta, double kdlj_phi,
411             int r_ge_d, int J) nogil
412     int qpms_cython_trans_calculator_get_AB_arrays_loop(qpms_trans_calculator *c,
413             int J, int resnd,
414             int daxis, int saxis,
415             char *A_data, np.npy_intp *A_shape, np.npy_intp *A_strides,
416             char *B_data, np.npy_intp *B_shape, np.npy_intp *B_strides,
417             char *r_data, np.npy_intp *r_shape, np.npy_intp *r_strides,
418             char *theta_data, np.npy_intp *theta_shape, np.npy_intp *theta_strides,
419             char *phi_data, np.npy_intp *phi_shape, np.npy_intp *phi_strides,
420             char *r_ge_d_data, np.npy_intp *phi_shape, np.npy_intp *phi_strides) nogil
422     int qpms_trans_calculator_get_trans_array(const qpms_trans_calculator *c,
423                 cdouble *target,
424                 const qpms_vswf_set_spec_t *destspec, size_t deststride,
425                 const qpms_vswf_set_spec_t *srcspec, size_t srcstride,
426                 csph_t kdlj, bint r_ge_d, qpms_bessel_t J);
428     int qpms_trans_calculator_get_trans_array_lc3p(
429                 const qpms_trans_calculator *c,
430                 cdouble *target,
431                 const qpms_vswf_set_spec_t *destspec, size_t deststride,
432                 const qpms_vswf_set_spec_t *srcspec, size_t srcstride,
433                 cdouble k, cart3_t destpos, cart3_t srcpos,
434                 qpms_bessel_t J
435                 );
437 cdef extern from "qpms_specfunc.h":
438     struct qpms_pitau_t:
439         qpms_l_t lMax
440         double *leg
441         double *pi
442         double *tau
443     qpms_pitau_t qpms_pitau_get(double theta, qpms_l_t lMax, double csphase)
444     qpms_errno_t qpms_pitau_fill(double *leg, double *pi, double *tau, double theta, qpms_l_t lMax, double csphase)
445     void qpms_pitau_free(qpms_pitau_t pitau)
448 cdef extern from "gsl/gsl_interp.h":
449     struct gsl_interp_type:
450         pass
451     const gsl_interp_type *gsl_interp_linear
452     const gsl_interp_type *gsl_interp_cspline
453     # ^^^ These are probably the only relevant ones.
455 cdef extern from "materials.h":
456     struct qpms_epsmu_generator_t:
457         qpms_epsmu_t (*function)(cdouble omega, const void *params)
458         const void *params
459     qpms_epsmu_t qpms_epsmu_const_g(cdouble omega, const void *params)
460     qpms_epsmu_t qpms_permittivity_interpolator_epsmu_g(cdouble omega, const void *epsmu)
461     qpms_epsmu_t qpms_lorentzdrude_epsmu_g(cdouble omega, const void *ldparams)
463     struct qpms_permittivity_interpolator_t:
464         pass
465     qpms_permittivity_interpolator_t *qpms_permittivity_interpolator_create(const size_t incount,
466             cdouble *wavelength_m, cdouble *n, cdouble *k, const gsl_interp_type *iptype)
467     qpms_permittivity_interpolator_t *qpms_permittivity_interpolator_from_yml(const char *path,
468             const gsl_interp_type *iptype)
469     cdouble qpms_permittivity_interpolator_eps_at_omega(const qpms_permittivity_interpolator_t *interp, double omega_SI)
470     double qpms_permittivity_interpolator_omega_max(const qpms_permittivity_interpolator_t *interp)
471     double qpms_permittivity_interpolator_omega_min(const qpms_permittivity_interpolator_t *interp)
472     void qpms_permittivity_interpolator_free(qpms_permittivity_interpolator_t *interp)
473     struct qpms_ldparams_triple_t:
474         double f
475         double omega
476         double gamma
477     struct qpms_ldparams_t:
478         cdouble eps_inf
479         double omega_p
480         size_t n
481         qpms_ldparams_triple_t data[0]
482     cdouble qpms_lorentzdrude_eps(cdouble, const qpms_ldparams_t *)
484 cdef extern from "tmatrices.h":
485     bint qpms_load_scuff_tmatrix_crash_on_failure
486     struct qpms_tmatrix_generator_t:
487         qpms_errno_t (*function)(qpms_tmatrix_t *t, cdouble omega, const void *params)
488         const void *params
489     qpms_errno_t qpms_tmatrix_generator_axialsym(qpms_tmatrix_t *t, cdouble omega, const void *params)
490     qpms_errno_t qpms_tmatrix_generator_interpolator(qpms_tmatrix_t *t, cdouble omega, const void *params)
491     qpms_errno_t qpms_tmatrix_generator_sphere(qpms_tmatrix_t *t, cdouble omega, const void *params)
492     qpms_errno_t qpms_tmatrix_generator_constant(qpms_tmatrix_t *t, cdouble omega, const void *params)
493     struct qpms_tmatrix_generator_sphere_param_t:
494         qpms_epsmu_generator_t outside
495         qpms_epsmu_generator_t inside
496         double radius
497     struct qpms_arc_function_retval_t:
498         double r
499         double beta
500     struct qpms_arc_function_t:
501         qpms_arc_function_retval_t (*function)(double theta, const void *params)
502         const void *params
503     struct qpms_tmatrix_generator_axialsym_param_t:
504         qpms_epsmu_generator_t outside
505         qpms_epsmu_generator_t inside
506         qpms_arc_function_t shape
507         qpms_l_t lMax_extend
508     struct qpms_arc_cylinder_params_t:
509         double R
510         double h
511     qpms_arc_function_retval_t qpms_arc_cylinder(double theta, const void *params)
512     qpms_arc_function_retval_t qpms_arc_sphere(double theta, const void *params)
513     struct qpms_tmatrix_interpolator_t:
514         const qpms_vswf_set_spec_t *bspec
515     void qpms_tmatrix_interpolator_free(qpms_tmatrix_interpolator_t *interp)
516     qpms_tmatrix_t *qpms_tmatrix_interpolator_eval(const qpms_tmatrix_interpolator_t *interp, double freq)
517     qpms_tmatrix_interpolator_t *qpms_tmatrix_interpolator_create(size_t n, double *freqs, 
518             const qpms_tmatrix_t *tmatrices_array, const gsl_interp_type *iptype)
519     void qpms_tmatrix_free(qpms_tmatrix_t *tmatrix)
520     qpms_tmatrix_isclose(const qpms_tmatrix_t *A, const qpms_tmatrix_t *B,
521                 const double rtol, const double atol)
522     qpms_errno_t qpms_symmetrise_tmdata_irot3arr(
523             cdouble *tmdata, const size_t tmcount,
524             const qpms_vswf_set_spec_t *bspec,
525             size_t n_symops,
526             const qpms_irot3_t *symops
527             )
528     qpms_errno_t qpms_symmetrise_tmdata_finite_group(
529             cdouble *tmdata, const size_t tmcount,
530             const qpms_vswf_set_spec_t *bspec,
531             const qpms_finite_group_t *pointgroup
532             )
533     qpms_tmatrix_t *qpms_tmatrix_symmetrise_irot3arr_inplace(
534             qpms_tmatrix_t *T,
535             size_t n_symops,
536             const qpms_irot3_t *symops
537             )
538     qpms_tmatrix_t *qpms_tmatrix_symmetrise_finite_group_inplace(
539             qpms_tmatrix_t *T,
540             const qpms_finite_group_t *pointgroup
541             )
542     qpms_errno_t qpms_load_scuff_tmatrix(const char *path, const qpms_vswf_set_spec_t *bspec,
543             size_t *n, double **freqs, double **freqs_su, qpms_tmatrix_t **tmatrices_array,
544             cdouble **tmdata)
545     cdouble *qpms_mie_coefficients_reflection(cdouble *target, const qpms_vswf_set_spec_t *bspec,
546             double a, cdouble k_i, cdouble k_e, cdouble mu_i, cdouble mu_e, qpms_bessel_t J_ext, qpms_bessel_t J_scat)
547     qpms_tmatrix_t *qpms_tmatrix_spherical(const qpms_vswf_set_spec_t *bspec, double a, 
548             cdouble k_i, cdouble k_e, cdouble mu_i, cdouble mu_e)
549     qpms_errno_t qpms_tmatrix_spherical_fill(qpms_tmatrix_t *t, double a, 
550             cdouble k_i, cdouble k_e, cdouble mu_i, cdouble mu_e)
551     qpms_tmatrix_t *qpms_tmatrix_spherical(const qpms_vswf_set_spec_t *bspec,
552             double a, cdouble k_i, cdouble k_e, cdouble mu_i, cdouble mu_e)
553     cdouble qpms_drude_epsilon(cdouble eps_inf, cdouble omega_p, cdouble gamma_p, cdouble omega)
554     qpms_errno_t qpms_tmatrix_spherical_mu0_fill(qpms_tmatrix_t *t, double a, double omega,
555             cdouble epsilon_fg, cdouble epsilon_bg)
556     qpms_tmatrix_t *qpms_tmatrix_spherical_mu0(const qpms_vswf_set_spec_t *bspec, double a,
557             double omega, cdouble epsilon_fg, cdouble epsilon_bg)
558     qpms_errno_t qpms_tmatrix_generator_axialsym_RQ_transposed_fill(cdouble *target, cdouble omega,
559             const qpms_tmatrix_generator_axialsym_param_t *param, qpms_normalisation_t norm, qpms_bessel_t J)
560     struct qpms_tmatrix_function_t:
561         const qpms_vswf_set_spec_t *spec
562         const qpms_tmatrix_generator_t *gen
563     ctypedef enum qpms_tmatrix_operation_kind_t:
564         QPMS_TMATRIX_OPERATION_NOOP
565         QPMS_TMATRIX_OPERATION_LRMATRIX
566         QPMS_TMATRIX_OPERATION_IROT3
567         QPMS_TMATRIX_OPERATION_IROT3ARR
568         QPMS_TMATRIX_OPERATION_COMPOSE_SUM
569         QPMS_TMATRIX_OPERATION_COMPOSE_CHAIN
570         QPMS_TMATRIX_OPERATION_SCMULZ
571         QPMS_TMATRIX_OPERATION_FINITE_GROUP_SYMMETRISE
572     
573     struct qpms_tmatrix_operation_t:
574         qpms_tmatrix_operation_kind_t typ
575         pass # TODO add the op union later if needed
576     const qpms_tmatrix_operation_t qpms_tmatrix_operation_noop
577     void qpms_tmatrix_operation_clear(qpms_tmatrix_operation_t *)
579 cdef extern from "pointgroups.h":
580     bint qpms_pg_is_finite_axial(qpms_pointgroup_class cls)
581     double qpms_pg_quat_cmp_atol
582     int qpms_pg_irot3_cmp(const qpms_irot3_t *, const qpms_irot3_t *);
583     int qpms_pg_irot3_cmp_v(const void *, const void *);
584     int qpms_pg_irot3_approx_cmp(const qpms_irot3_t *a, const qpms_irot3_t *b, double atol)
585     int qpms_pg_irot3_approx_cmp_v(const void *a, const void *b)
587     qpms_gmi_t qpms_pg_order(qpms_pointgroup_class cls, qpms_gmi_t n)
588     qpms_irot3_t *qpms_pg_canonical_elems( qpms_irot3_t *target, qpms_pointgroup_class cls, qpms_gmi_t n)
589     qpms_gmi_t qpms_pg_genset_size(qpms_pointgroup_class cls, qpms_gmi_t n)
590     qpms_gmi_t qpms_pg_genset(qpms_pointgroup_class cls, qpms_gmi_t n, qpms_irot3_t *gen)
591     qpms_irot3_t *qpms_pg_elems(qpms_irot3_t *target, qpms_pointgroup_t g)
592     bint qpms_pg_is_subgroup(qpms_pointgroup_t a, qpms_pointgroup_t b);
594 cdef extern from "scatsystem.h":
595     void qpms_scatsystem_set_nthreads(long n)
596     struct qpms_particle_t:
597         cart3_t pos
598         const qpms_tmatrix_function_t *tmg
599         qpms_tmatrix_operation_t op
600     struct qpms_particle_tid_t:
601         cart3_t pos
602         qpms_ss_tmi_t tmatrix_id
603     struct qpms_ss_derived_tmatrix_t:
604         qpms_ss_tmgi_t tmgi
605         qpms_tmatrix_operation_t op
606     struct qpms_scatsys_periodic_info_t:
607         cart3_t lattice_basis[3]
608         cart3_t reciprocal_basis[3]
609         double unitcell_volume
610         double eta
611         #etc.
612     struct qpms_scatsys_t:
613         int lattice_dimension
614         qpms_epsmu_generator_t medium
615         qpms_tmatrix_function_t *tmg
616         qpms_ss_tmgi_t tmg_count
617         qpms_ss_derived_tmatrix_t *tm
618         qpms_ss_tmi_t tm_count
619         qpms_particle_tid_t *p
620         qpms_ss_pi_t p_count
621         # We shouldn't need more to construct a symmetric scatsystem ^^^
622         size_t fecv_size
623         size_t *saecv_sizes
624         const qpms_finite_group_t *sym
625         qpms_scatsys_periodic_info_t per
627         # per[] and other stuff not currently needed in cython
628     void qpms_scatsys_free(qpms_scatsys_t *s)
629     qpms_errno_t qpms_scatsys_dump(qpms_scatsys_t *ss, char *path) #NI
630     qpms_scatsys_t *qpms_scatsys_load(char *path) #NI
631     struct qpms_scatsys_at_omega_t:
632         const qpms_scatsys_t *ss
633         qpms_tmatrix_t **tm,
634         cdouble omega
635         qpms_epsmu_t medium
636         cdouble wavenumber
637     qpms_scatsys_at_omega_t *qpms_scatsys_apply_symmetry(const qpms_scatsys_t *orig, const qpms_finite_group_t *sym,
638             cdouble omega, const qpms_tolerance_spec_t *tol)
639     qpms_scatsys_at_omega_t *qpms_scatsys_at_omega(const qpms_scatsys_t *ss, cdouble omega)
640     void qpms_scatsys_at_omega_free(qpms_scatsys_at_omega_t *ssw)
641     cdouble *qpms_scatsys_irrep_pack_matrix(cdouble *target_packed,
642             const cdouble *orig_full, const qpms_scatsys_t *ss, qpms_iri_t iri)
643     cdouble *qpms_scatsys_irrep_unpack_matrix(cdouble *target_full, 
644             const cdouble *orig_packed, const qpms_scatsys_t *ss, qpms_iri_t iri, bint add)
645     cdouble *qpms_scatsys_irrep_pack_vector(cdouble *target_packed,
646             const cdouble *orig_full, const qpms_scatsys_t *ss, qpms_iri_t iri)
647     cdouble *qpms_scatsys_irrep_unpack_vector(cdouble *target_full,
648             const cdouble *orig_packed, const qpms_scatsys_t *ss, qpms_iri_t iri, bint add)
649     cdouble *qpms_scatsysw_build_modeproblem_matrix_full(cdouble *target,
650             const qpms_scatsys_at_omega_t *ssw)
651     cdouble *qpms_scatsys_build_translation_matrix_full(cdouble *target,
652             const qpms_scatsys_t *ss, cdouble k)
653     cdouble *qpms_scatsys_build_translation_matrix_e_full(cdouble *target,
654             const qpms_scatsys_t *ss, cdouble k, qpms_bessel_t J)
655     cdouble *qpms_scatsysw_build_modeproblem_matrix_irrep_packed(cdouble *target,
656             const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri) nogil
657     cdouble *qpms_scatsys_build_translation_matrix_e_irrep_packed(cdouble *target,
658             const qpms_scatsys_t *ss, qpms_iri_t iri, cdouble k, qpms_bessel_t J) nogil
659     cdouble *qpms_scatsysw_build_modeproblem_matrix_irrep_packed_orbitorderR(
660             cdouble *target, const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri) nogil
661     cdouble *qpms_scatsysw_build_modeproblem_matrix_irrep_packed_serial(
662             cdouble *target, const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri) nogil
663     cdouble *qpms_scatsys_incident_field_vector_full(cdouble *target_full,
664             const qpms_scatsys_t *ss, qpms_incfield_t field_at_point, 
665             const void *args, bint add) nogil
666     cdouble *qpms_scatsysw_apply_Tmatrices_full(cdouble *target_full, const cdouble *inc_full,
667             const qpms_scatsys_at_omega_t *ssw) nogil
668     struct qpms_ss_LU:
669         const qpms_scatsys_at_omega_t *ssw
670         const qpms_scatsys_at_omega_k_t *sswk
671         bint full
672         qpms_iri_t iri
673         cdouble *a
674         int *ipiv
675     void qpms_ss_LU_free(qpms_ss_LU lu)
676     qpms_ss_LU qpms_scatsysw_build_modeproblem_matrix_full_LU(cdouble *target,
677             int *target_piv, const qpms_scatsys_at_omega_t *ssw) nogil
678     qpms_ss_LU qpms_scatsysw_build_modeproblem_matrix_irrep_packed_LU(cdouble *target,
679             int *target_piv, const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri) nogil
680     qpms_ss_LU qpms_scatsysw_modeproblem_matrix_full_factorise(cdouble *modeproblem_matrix_full,
681             int *target_piv, const qpms_scatsys_at_omega_t *ssw) nogil
682     qpms_ss_LU qpms_scatsysw_modeproblem_matrix_irrep_packed_factorise(cdouble *modeproblem_matrix_full,
683             int *target_piv, const qpms_scatsys_at_omega_t *ssw, qpms_iri_t iri) nogil
684     cdouble *qpms_scatsys_scatter_solve(cdouble *target_f, const cdouble *a_inc, qpms_ss_LU ludata) nogil
685     const qpms_vswf_set_spec_t *qpms_ss_bspec_tmi(const qpms_scatsys_t *ss, qpms_ss_tmi_t tmi) nogil
686     const qpms_vswf_set_spec_t *qpms_ss_bspec_pi(const qpms_scatsys_t *ss, qpms_ss_pi_t pi) nogil
687     beyn_result_t *qpms_scatsys_finite_find_eigenmodes(const qpms_scatsys_t *ss, qpms_iri_t iri,
688             cdouble omega_centre, double omega_rr, double omega_ri, size_t contour_npoints, 
689             double rank_tol, size_t rank_sel_min, double res_tol) nogil
690     # periodic-related funs
691     struct qpms_scatsys_at_omega_k_t:
692         const qpms_scatsys_at_omega_t *ssw
693         double k[3]
694         double eta
695     cdouble *qpms_scatsyswk_build_modeproblem_matrix_full(cdouble *target, const qpms_scatsys_at_omega_k_t *sswk)
696     cdouble *qpms_scatsys_periodic_build_translation_matrix_full(cdouble *target, const qpms_scatsys_t *ss, cdouble wavenumber, const cart3_t *wavevector, double eta)
697     qpms_ss_LU qpms_scatsyswk_build_modeproblem_matrix_full_LU(cdouble *target, int *target_piv, const qpms_scatsys_at_omega_k_t *sswk)
698     beyn_result_t *qpms_scatsys_periodic_find_eigenmodes(const qpms_scatsys_t *ss, const double *k,
699             cdouble omega_centre, double omega_rr, double omega_ri, size_t contour_npoints, 
700             double rank_tol, size_t rank_sel_min, double res_tol) nogil
701     ccart3_t qpms_scatsys_scattered_E(const qpms_scatsys_t *ss, qpms_bessel_t btyp, cdouble wavenumber,
702             const cdouble *f_excitation_vector_full, cart3_t where) nogil
703     ccart3_t qpms_scatsysw_scattered_E(const qpms_scatsys_at_omega_t *ssw, qpms_bessel_t btyp,
704             const cdouble *f_excitation_vector_full, cart3_t where) nogil
705     ccart3_t qpms_scatsys_scattered_E__alt(const qpms_scatsys_t *ss, qpms_bessel_t btyp, cdouble wavenumber,
706             const cdouble *f_excitation_vector_full, cart3_t where) nogil
707     ccart3_t qpms_scatsysw_scattered_E__alt(const qpms_scatsys_at_omega_t *ssw, qpms_bessel_t btyp,
708             const cdouble *f_excitation_vector_full, cart3_t where) nogil
709     ccart3_t qpms_scatsyswk_scattered_E(const qpms_scatsys_at_omega_k_t *sswk, qpms_bessel_t btyp,
710             const cdouble *f_excitation_vector_full, cart3_t where) nogil
711     qpms_errno_t qpms_scatsyswk_scattered_field_basis(ccart3_t *target, const qpms_scatsys_at_omega_k_t *sswk,
712             qpms_bessel_t btyp, cart3_t where) nogil
713     double qpms_ss_adjusted_eta(const qpms_scatsys_t *ss, cdouble wavenumber, const double *wavevector) nogil
716 cdef extern from "ewald.h":
717     struct qpms_csf_result:
718         cdouble val
719         double err
721     ctypedef enum qpms_ewald_part:
722         QPMS_EWALD_LONG_RANGE
723         QPMS_EWALD_SHORT_RANGE
724         QPMS_EWALD_FULL
725         QPMS_EWALD_0TERM
727     struct qpms_ewald3_constants_t:
728         qpms_l_t lMax
729         qpms_y_t nelem_sc
730     qpms_ewald3_constants_t *qpms_ewald3_constants_init(qpms_l_t lMax, int csphase)
731     void qpms_ewald3_constants_free(qpms_ewald3_constants_t *c)
732     
733     cdouble lilgamma(double t)
734     cdouble clilgamma(cdouble z)
735     int cx_gamma_inc_series_e(double a, cdouble x, qpms_csf_result *result)
736     int cx_gamma_inc_CF_e(double a, cdouble x, qpms_csf_result *result)
737     int complex_gamma_inc_e(double a, cdouble x, int m, qpms_csf_result *result)
739     int ewald3_sigma0(cdouble *target, double *err, const qpms_ewald3_constants_t *c, double eta, cdouble wavenumber)
740     int ewald3_sigma_long(cdouble *target_sigmalr_y, double *target_sigmalr_y_err, const qpms_ewald3_constants_t *c, 
741             double eta, cdouble wavenumber, double unitcell_volume, LatticeDimensionality latdim, PGen *pgen_K,
742             bint pgen_generates_shifted_points, cart3_t k, cart3_t particle_shift)
743     int ewald3_sigma_short(cdouble *target_sigmasr_y, double *target_sigmasr_y_err, const qpms_ewald3_constants_t *c, 
744             double eta, cdouble wavenumber, LatticeDimensionality latdim, PGen *pgen_R, 
745             bint pgen_generates_shifted_points, cart3_t k, cart3_t particle_shift)
747 cdef extern from "beyn.h":
748     ctypedef struct beyn_contour_t:
749         bint (*inside_test)(beyn_contour_t *, cdouble z)
750         pass
751     ctypedef struct beyn_result_t:
752         size_t neig
753         size_t vlen
754         cdouble *eigval
755         cdouble *eigval_err
756         double *residuals
757         cdouble *eigvec
758         double *ranktest_SV
759     ctypedef enum beyn_contour_halfellipse_orientation:
760         BEYN_CONTOUR_HALFELLIPSE_RE_PLUS
761         BEYN_CONTOUR_HALFELLIPSE_IM_PLUS
762         BEYN_CONTOUR_HALFELLIPSE_RE_MINUS
763         BEYN_CONTOUR_HALFELLIPSE_IM_MINUS
765     ctypedef int (*beyn_function_M_t)(cdouble *target_M, size_t m, cdouble z, void *params)
766     ctypedef int (*beyn_function_M_inv_Vhat_t)(cdouble *target, size_t m, size_t l, const cdouble *Vhat, cdouble z, void *params)
768     void beyn_result_free(beyn_result_t *result)
770     beyn_result_t *beyn_solve(size_t m, size_t l, beyn_function_M_t M,
771             beyn_function_M_inv_Vhat_t M_inv_Vhat, void *params, const beyn_contour_t *contour,
772             double rank_tol, size_t rank_min_sel, double res_tol)
774     beyn_contour_t *beyn_contour_ellipse(cdouble centre, double halfax_re, double halfax_im, size_t npoints)
775     beyn_contour_t *beyn_contour_halfellipse(cdouble centre, double halfax_re, double halfax_im, size_t npoints,
776             beyn_contour_halfellipse_orientation ori)
777     beyn_contour_t *beyn_contour_kidney(cdouble centre, double halfax_re, double halfax_im, size_t npoints,
778             double rounding, beyn_contour_halfellipse_orientation ori)