Disable "hard" examples in CI
[qpms.git] / qpms / qpms_c.pyx
blob68ba71e231321f59f011207aaaf33d86c73dbd1e
1 """@package qpms_c
2 Cythonized parts of QPMS; mostly wrappers over the C data structures
3 to make them available in Python.
4 """
6 # Cythonized parts of QPMS here
7 # -----------------------------
9 import numpy as np
10 from .qpms_cdefs cimport *
11 from .cyquaternions cimport IRot3, CQuat
12 from .cybspec cimport BaseSpec
13 from .cycommon cimport make_c_string
14 from .cycommon import string_c2py, PointGroupClass, BesselType
15 from .cytmatrices cimport CTMatrix, TMatrixFunction, TMatrixGenerator, TMatrixInterpolator
16 from .cymaterials cimport EpsMuGenerator, EpsMu
17 from libc.stdlib cimport malloc, free, calloc
18 import warnings
20 def qpms_library_version():
21 return "git_" + qpms_version_git_sha1.decode('utf-8')
23 from cython.parallel cimport prange, parallel
24 from cython cimport boundscheck, wraparound
26 # Set custom GSL error handler. N.B. this is obviously not thread-safe.
27 cdef char *pgsl_err_reason
28 cdef char *pgsl_err_file
29 cdef int pgsl_err_line
30 cdef int pgsl_errno = 0
31 cdef int *pgsl_errno_ignorelist = NULL # list of ignored error codes, terminated by zero
33 # This error handler only sets the variables above
34 cdef void pgsl_error_handler(const char *reason, const char *_file, const int line, const int gsl_errno):
35 global pgsl_err_reason, pgsl_err_file, pgsl_err_line, pgsl_errno, pgsl_errno_ignorelist
36 cdef size_t i
37 if(pgsl_errno_ignorelist):
38 i = 0
39 while pgsl_errno_ignorelist[i] != 0:
40 if gsl_errno == pgsl_errno_ignorelist[i]:
41 return
42 i += 1
43 pgsl_err_file = _file
44 pgsl_err_reason = reason
45 pgsl_errno = gsl_errno
46 pgsl_err_line = line
47 return
49 cdef const int* pgsl_set_ignorelist(const int *new_ignorelist):
50 global pgsl_errno_ignorelist
51 cdef const int *oldlist = pgsl_errno_ignorelist
52 pgsl_errno_ignorelist = new_ignorelist
53 return oldlist
55 cdef class pgsl_ignore_error():
56 '''Context manager for setting a temporary list of errnos ignored by pgsl_error_handler.
57 Always sets pgsl_error_handler.
59 Performs pgsl_check_err() on exit unless
60 '''
61 cdef const int *ignorelist_old
62 cdef gsl_error_handler_t *old_handler
63 cdef bint final_check
64 cdef object ignorelist_python
66 cdef int *ignorelist
67 def __cinit__(self, *ignorelist, **kwargs):
68 self.ignorelist = <int*>calloc((len(ignorelist)+1), sizeof(int))
69 self.ignorelist_python = ignorelist
70 for i in range(len(ignorelist)):
71 self.ignorelist[i] = ignorelist[i]
72 if "final_check" in kwargs.keys() and not kwargs["final_check"]:
73 final_check = True
74 final_check = False
76 def __enter__(self):
77 global pgsl_error_handler
78 self.ignorelist_old = pgsl_set_ignorelist(self.ignorelist)
79 self.old_handler = gsl_set_error_handler(pgsl_error_handler)
80 return
82 def __exit__(self, type, value, traceback):
83 global pgsl_errno_ignorelist, pgsl_error_handler
84 pgsl_set_ignorelist(self.ignorelist_old)
85 gsl_set_error_handler(self.old_handler)
86 if self.final_check:
87 pgsl_check_err(retval = None, ignore = self.ignorelist_python)
89 def __dealloc__(self):
90 free(self.ignorelist)
92 def pgsl_check_err(retval = None, ignorelist = None):
93 global pgsl_err_reason, pgsl_err_file, pgsl_err_line, pgsl_errno
94 '''Check for possible errors encountered by pgsl_error_handler.
95 Takes return value of a function as an optional argument, which is now ignored.
96 '''
97 cdef int errno_was
98 if (pgsl_errno != 0):
99 errno_was = pgsl_errno
100 pgsl_errno = 0
101 raise RuntimeError("Error %d in GSL calculation in %s:%d: %s" % (errno_was,
102 string_c2py(pgsl_err_file), pgsl_err_line, string_c2py(pgsl_err_reason)))
103 if (retval is not None and retval != 0 and ignorelist is not None and retval not in ignorelist):
104 warnings.warn("Got non-zero return value %d" % retval)
105 if retval is not None:
106 return retval
107 else:
108 return 0
110 def set_gsl_pythonic_error_handling():
112 Sets pgsl_error_handler as the GSL error handler to avoid crashing.
114 gsl_set_error_handler(pgsl_error_handler)
116 cdef class PointGroup:
117 cdef readonly qpms_pointgroup_t G
119 def __init__(self, cls, qpms_gmi_t n = 0, IRot3 orientation = IRot3()):
120 cls = PointGroupClass(cls)
121 self.G.c = cls
122 if n <= 0 and qpms_pg_is_finite_axial(cls):
123 raise ValueError("For finite axial groups, n argument must be positive")
124 self.G.n = n
125 self.G.orientation = orientation.qd
127 def __len__(self):
128 return qpms_pg_order(self.G.c, self.G.n);
130 def __le__(PointGroup self, PointGroup other):
131 if qpms_pg_is_subgroup(self.G, other.G):
132 return True
133 else:
134 return False
135 def __ge__(PointGroup self, PointGroup other):
136 if qpms_pg_is_subgroup(other.G, self.G):
137 return True
138 else:
139 return False
140 def __lt__(PointGroup self, PointGroup other):
141 return qpms_pg_is_subgroup(self.G, other.G) and not qpms_pg_is_subgroup(other.G, self.G)
142 def __eq__(PointGroup self, PointGroup other):
143 return qpms_pg_is_subgroup(self.G, other.G) and qpms_pg_is_subgroup(other.G, self.G)
144 def __gt__(PointGroup self, PointGroup other):
145 return not qpms_pg_is_subgroup(self.G, other.G) and qpms_pg_is_subgroup(other.G, self.G)
147 def elems(self):
148 els = list()
149 cdef qpms_irot3_t *arr
150 arr = qpms_pg_elems(NULL, self.G)
151 cdef IRot3 q
152 for i in range(len(self)):
153 q = IRot3()
154 q.cset(arr[i])
155 els.append(q)
156 free(arr)
157 return els
159 cdef class FinitePointGroup:
161 Wrapper over the qpms_finite_group_t structure.
163 TODO more functionality to make it better usable in Python
164 (group element class at least)
167 def __cinit__(self, info):
168 '''Constructs a FinitePointGroup from PointGroupInfo'''
169 # TODO maybe I might use a try..finally statement to avoid leaks
170 # First, generate all basic data from info
171 permlist = info.deterministic_elemlist()
172 cdef int order = len(permlist)
173 permindices = {perm: i for i, perm in enumerate(permlist)} # 'invert' permlist
174 identity = info.permgroup.identity
175 # We use calloc to avoid calling free to unitialized pointers
176 self.G = <qpms_finite_group_t *>calloc(1,sizeof(qpms_finite_group_t))
177 if not self.G: raise MemoryError
178 self.G[0].name = make_c_string(info.name)
179 self.G[0].order = order
180 self.G[0].idi = permindices[identity]
181 self.G[0].mt = <qpms_gmi_t *>malloc(sizeof(qpms_gmi_t) * order * order)
182 if not self.G[0].mt: raise MemoryError
183 for i in range(order):
184 for j in range(order):
185 self.G[0].mt[i*order + j] = permindices[permlist[i] * permlist[j]]
186 self.G[0].invi = <qpms_gmi_t *>malloc(sizeof(qpms_gmi_t) * order)
187 if not self.G[0].invi: raise MemoryError
188 for i in range(order):
189 self.G[0].invi[i] = permindices[permlist[i]**-1]
190 self.G[0].ngens = len(info.permgroupgens)
191 self.G[0].gens = <qpms_gmi_t *>malloc(sizeof(qpms_gmi_t) * self.G[0].ngens)
192 if not self.G[0].gens: raise MemoryError
193 for i in range(self.G[0].ngens):
194 self.G[0].gens[i] = permindices[info.permgroupgens[i]]
195 self.G[0].permrep = <char **>calloc(order, sizeof(char *))
196 if not self.G[0].permrep: raise MemoryError
197 for i in range(order):
198 self.G[0].permrep[i] = make_c_string(str(permlist[i]))
199 if not self.G[0].permrep[i]: raise MemoryError
200 self.G[0].permrep_nelem = info.permgroup.degree
201 if info.rep3d is not None:
202 self.G[0].rep3d = <qpms_irot3_t *>malloc(order * sizeof(qpms_irot3_t))
203 for i in range(order):
204 self.G[0].rep3d[i] = info.rep3d[permlist[i]].qd
205 self.G[0].nirreps = len(info.irreps)
206 self.G[0].irreps = <qpms_finite_group_irrep_t *>calloc(self.G[0].nirreps, sizeof(qpms_finite_group_irrep_t))
207 if not self.G[0].irreps: raise MemoryError
208 cdef int dim
209 for iri, irname in enumerate(sorted(info.irreps.keys())):
210 irrep = info.irreps[irname]
211 is1d = isinstance(irrep[identity], (int, float, complex))
212 dim = 1 if is1d else irrep[identity].shape[0]
213 self.G[0].irreps[iri].dim = dim
214 self.G[0].irreps[iri].name = <char *>make_c_string(irname)
215 if not self.G[0].irreps[iri].name: raise MemoryError
216 self.G[0].irreps[iri].m = <cdouble *>malloc(dim*dim*sizeof(cdouble)*order)
217 if not self.G[0].irreps[iri].m: raise MemoryError
218 if is1d:
219 for i in range(order):
220 self.G[0].irreps[iri].m[i] = irrep[permlist[i]]
221 else:
222 for i in range(order):
223 for row in range(dim):
224 for col in range(dim):
225 self.G[0].irreps[iri].m[i*dim*dim + row*dim + col] = irrep[permlist[i]][row,col]
226 self.G[0].elemlabels = <char **> 0 # Elem labels not yet implemented
227 self.owns_data = True
229 def __dealloc__(self):
230 cdef qpms_gmi_t order
231 if self.owns_data:
232 if self.G:
233 order = self.G[0].order
234 free(self.G[0].name)
235 free(self.G[0].mt)
236 free(self.G[0].invi)
237 free(self.G[0].gens)
238 if self.G[0].permrep:
239 for i in range(order): free(self.G[0].permrep[i])
240 free(self.G[0].permrep)
241 if self.G[0].elemlabels: # this is not even contructed right now
242 for i in range(order): free(self.G[0].elemlabels[i])
243 if self.G[0].irreps:
244 for iri in range(self.G[0].nirreps):
245 free(self.G[0].irreps[iri].name)
246 free(self.G[0].irreps[iri].m)
247 free(self.G[0].irreps)
248 free(self.G)
249 self.G = <qpms_finite_group_t *>0
250 self.owns_data = False
252 property order: # LPTODO might instead be __len__() if iterable at some point
253 def __get__(self):
254 return self.G[0].order
256 @staticmethod
257 def TRIVIAL(): # quite ugly
258 from .symmetries import point_group_info
259 return FinitePointGroup(point_group_info['trivial_g'])
262 cdef class FinitePointGroupElement:
263 '''TODO'''
264 cdef readonly FinitePointGroup G
265 cdef readonly qpms_gmi_t gmi
266 def __cinit__(self, FinitePointGroup G, qpms_gmi_t gmi):
267 self.G = G
268 self.gmi = gmi
270 cdef class Particle:
272 Wrapper over the qpms_particle_t structure.
274 cdef qpms_particle_t p
275 cdef readonly TMatrixFunction f # Reference to ensure correct reference counting
278 def __cinit__(Particle self, pos, t, bspec = None):
279 """Particle object constructor.
281 Parameters
282 ----------
283 pos : (x, y, z)
284 Particle position in cartesian coordinates
286 t : TMatrixGenerator or CTMatrix or TMatrixInterpolator
287 T-matrix specification for the particle.
289 bspec : BaseSpec
290 WSWF basis specification for the particle. Might be omitted if
291 t is a CTMatrix.
293 cdef TMatrixGenerator tgen
294 cdef BaseSpec spec
295 if(len(pos)>=2 and len(pos) < 4):
296 self.p.pos.x = pos[0]
297 self.p.pos.y = pos[1]
298 self.p.pos.z = pos[2] if len(pos)==3 else 0
299 else:
300 raise ValueError("Position argument has to contain 3 or 2 cartesian coordinates")
301 if isinstance(t, CTMatrix):
302 tgen = TMatrixGenerator(t)
303 elif isinstance(t, TMatrixInterpolator):
304 tgen = TMatrixGenerator(t)
305 warnings.warn("Initialising a particle with interpolated T-matrix values. Imaginary frequencies will be discarded and mode search algorithm will yield nonsense (just saying).")
306 elif isinstance(t, TMatrixGenerator):
307 tgen = <TMatrixGenerator>t
308 else: raise TypeError('t must be either CTMatrix or TMatrixGenerator, was %s' % str(type(t)))
309 if bspec is not None:
310 spec = bspec
311 else:
312 if isinstance(tgen.holder, CTMatrix):
313 spec = (<CTMatrix>tgen.holder).spec
314 else:
315 raise ValueError("bspec argument must be specified separately for %s" % str(type(t)))
316 self.f = TMatrixFunction(tgen, spec)
317 self.p.tmg = self.f.rawpointer()
318 # TODO non-trivial transformations later; if modified, do not forget to update ScatteringSystem constructor
319 self.p.op = qpms_tmatrix_operation_noop
321 def __dealloc__(self):
322 qpms_tmatrix_operation_clear(&self.p.op)
324 cdef qpms_particle_t *rawpointer(Particle self):
325 '''Pointer to the qpms_particle_p structure.
327 return &(self.p)
328 property rawpointer:
329 def __get__(self):
330 return <uintptr_t> &(self.p)
332 cdef qpms_particle_t cval(Particle self):
333 '''Provides a copy for assigning in cython code'''
334 return self.p
336 property x:
337 def __get__(self):
338 return self.p.pos.x
339 def __set__(self,x):
340 self.p.pos.x = x
341 property y:
342 def __get__(self):
343 return self.p.pos.y
344 def __set__(self,y):
345 self.p.pos.y = y
346 property z:
347 def __get__(self):
348 return self.p.pos.z
349 def __set__(self,z):
350 self.p.pos.z = z
351 property pos:
352 def __get__(self):
353 return (self.p.pos.x, self.p.pos.y, self.p.pos.z)
354 def __set__(self, pos):
355 if(len(pos)>=2 and len(pos) < 4):
356 self.p.pos.x = pos[0]
357 self.p.pos.y = pos[1]
358 self.p.pos.z = pos[2] if len(pos)==3 else 0
359 else:
360 raise ValueError("Position argument has to contain 3 or 2 cartesian coordinates")
362 cpdef void scatsystem_set_nthreads(long n):
363 qpms_scatsystem_set_nthreads(n)
364 return
367 cdef class ScatteringSystem:
369 Wrapper over the C qpms_scatsys_t structure, representing a collection
370 of scatterers, finite or periodic.
372 Currently, it does not have a standard constructor. Use the
373 ScatteringSystem.create() method instead.
375 cdef list tmgobjs # here we keep the references to occuring TMatrixFunctions (and hence BaseSpecs and TMatrixGenerators)
376 #cdef list Tmatrices # Here we keep the references to occuring T-matrices
377 cdef EpsMuGenerator medium_holder # Here we keep the reference to medium generator
378 cdef qpms_scatsys_t *s
379 cdef FinitePointGroup sym
381 cdef qpms_iri_t iri_py2c(self, iri, allow_None = True):
382 if iri is None and allow_None:
383 return QPMS_NO_IRREP
384 cdef qpms_iri_t nir = self.nirreps
385 cdef qpms_iri_t ciri = iri
386 if ciri < 0 or ciri > nir:
387 raise ValueError("Invalid irrep index %s (of %d irreps)", str(iri), self.nirreps)
388 return ciri
390 def check_s(self): # cdef instead?
391 if self.s == <qpms_scatsys_t *>NULL:
392 raise ValueError("ScatteringSystem's s-pointer not set. You must not use the default constructor; use the create() method instead")
393 #TODO is there a way to disable the constructor outside this module?
395 @staticmethod # We don't have any "standard" constructor for this right now
396 def create(particles, medium, cdouble omega, FinitePointGroup sym = FinitePointGroup.TRIVIAL(),
397 latticebasis = None): # TODO tolerances
398 """(Non-standard) constructor of ScatteringSystem
400 Parameters
401 ----------
402 particles : list of Particles objects
403 Scatterers to be included in the system. These scatterers are then
404 copied around by the constructor using the point group symmetry defined in sym.
405 medium : EpsMu or EpsMuGenerator
406 Material properties of the background medium.
407 omega : complex
408 Any valid angular frequency for the initial T-matrix evaluation.
409 It must be a value which all the T-matrix generators and interpolators
410 referenced by particles can evaluate.
411 sym : FinitePointGroup
412 Symmetry group for a finite system.
413 Defaults to the trivial point group FinitePointGroup.TRIVIAL().
414 Currently, this must be left trivial for periodic systems.
415 latticebasis : None or array_like of float type.
416 Lattice base vectors of a periodic system in cartesian coordinates.
417 The shape must be [d][3], where d = 1, 2 or 3 determines the lattice dimension.
418 Defaults to None, i.e. a finite system.
420 Returns
421 -------
422 ss : ScatteringSystem
423 Object representing the system of compact scatterers.
424 ssw : _ScatteringSystemAtOmega
425 Object representing ss evaluated at omega (angular frequency used for the
426 initial evaluation).
428 Note
429 ----
430 Currently, the ScatterinSystem's T-matrices need to be evaluated for at least
431 one frequency in order to initialise the symmetry structure. For this reason,
432 this "constructor" creates an instance of ScatteringSystem and an instance
433 of _ScatteringSystemAtOmega at the same time, so that the evaluated T-matrices
434 can be used right away if needed. This is why this approach is used instead
435 of the usual __init__ method.
438 if latticebasis is not None and sym.order != 1:
439 raise NotImplementedError("Periodic systems don't currently support nontrivial point group symmetries")
441 # These we are going to construct
442 cdef ScatteringSystem self
443 cdef _ScatteringSystemAtOmega pyssw
445 cdef qpms_scatsys_t orig # This should be automatically init'd to 0 (CHECKME)
446 cdef qpms_ss_pi_t pi, p_count = len(particles)
447 cdef qpms_ss_tmi_t tmi, tm_count = 0
448 cdef qpms_ss_tmgi_t tmgi, tmg_count = 0
450 cdef qpms_scatsys_at_omega_t *ssw
451 cdef qpms_scatsys_t *ss
453 cdef Particle p
455 tmgindices = dict()
456 tmgobjs = list()
457 tmindices = dict()
458 tmlist = list()
459 for p in particles: # find and enumerate unique t-matrix generators
460 if p.p.op.typ != QPMS_TMATRIX_OPERATION_NOOP:
461 raise NotImplementedError("currently, only no-op T-matrix operations are allowed in ScatteringSystem constructor")
462 #tmg_key = id(p.f) # This causes a different generator for each particle -> SUPER SLOW
463 tmg_key = (id(p.f.generator), id(p.f.spec))
464 if tmg_key not in tmgindices:
465 tmgindices[tmg_key] = tmg_count
466 tmgobjs.append(p.f) # Save the references on BaseSpecs and TMatrixGenerators (via TMatrixFunctions)
467 tmg_count += 1
468 # Following lines have to be adjusted when nontrivial operations allowed:
469 tm_derived_key = (tmg_key, None) # TODO unique representation of p.p.op instead of None
470 if tm_derived_key not in tmindices:
471 tmindices[tm_derived_key] = tm_count
472 tmlist.append(tm_derived_key)
473 tm_count += 1
474 cdef EpsMuGenerator mediumgen = EpsMuGenerator(medium)
475 orig.medium = mediumgen.g
476 orig.tmg_count = tmg_count
477 orig.tm_count = tm_count
478 orig.p_count = p_count
479 try:
480 orig.tmg = <qpms_tmatrix_function_t *>malloc(orig.tmg_count * sizeof(orig.tmg[0]))
481 if not orig.tmg: raise MemoryError
482 orig.tm = <qpms_ss_derived_tmatrix_t *>malloc(orig.tm_count * sizeof(orig.tm[0]))
483 if not orig.tm: raise MemoryError
484 orig.p = <qpms_particle_tid_t *>malloc(orig.p_count * sizeof(orig.p[0]))
485 if not orig.p: raise MemoryError
486 for tmgi in range(orig.tmg_count):
487 orig.tmg[tmgi] = (<TMatrixFunction?>tmgobjs[tmgi]).raw()
488 for tmi in range(tm_count):
489 tm_derived_key = tmlist[tmi]
490 tmgi = tmgindices[tm_derived_key[0]]
491 orig.tm[tmi].tmgi = tmgi
492 orig.tm[tmi].op = qpms_tmatrix_operation_noop # TODO adjust when notrivial operations allowed
493 for pi in range(p_count):
494 p = particles[pi]
495 tmg_key = (id(p.f.generator), id(p.f.spec))
496 tm_derived_key = (tmg_key, None) # TODO unique representation of p.p.op instead of None
497 orig.p[pi].pos = p.cval().pos
498 orig.p[pi].tmatrix_id = tmindices[tm_derived_key]
499 if latticebasis is not None: # periodic system
500 assert(len(latticebasis) <= 3 and len(latticebasis) > 0)
501 orig.lattice_dimension = len(latticebasis)
502 for d in range(len(latticebasis)):
503 orig.per.lattice_basis[d] = {'x' : latticebasis[d][0], 'y' : latticebasis[d][1], 'z' : latticebasis[d][2] if len(latticebasis[d]) >= 3 else 0}
504 else: orig.lattice_dimension = 0
505 ssw = qpms_scatsys_apply_symmetry(&orig, sym.rawpointer(), omega, &QPMS_TOLERANCE_DEFAULT)
506 ss = ssw[0].ss
507 finally:
508 free(orig.tmg)
509 free(orig.tm)
510 free(orig.p)
511 self = ScatteringSystem()
512 self.medium_holder = mediumgen
513 self.s = ss
514 self.tmgobjs = tmgobjs
515 self.sym = sym
516 pyssw = _ScatteringSystemAtOmega()
517 pyssw.ssw = ssw
518 pyssw.ss_pyref = self
519 return self, pyssw
521 def __call__(self, cdouble omega):
522 self.check_s()
523 cdef _ScatteringSystemAtOmega pyssw = _ScatteringSystemAtOmega()
524 pyssw.ssw = qpms_scatsys_at_omega(self.s, omega)
525 pyssw.ss_pyref = self
526 return pyssw
528 def __dealloc__(self):
529 if(self.s):
530 qpms_scatsys_free(self.s)
532 property particles_tmi:
533 def __get__(self):
534 self.check_s()
535 r = list()
536 cdef qpms_ss_pi_t pi
537 for pi in range(self.s[0].p_count):
538 r.append(self.s[0].p[pi])
539 return r
541 property particle_count:
542 """Number of particles in the system (or unit cell)"""
543 def __get__(self):
544 self.check_s()
545 return self.s[0].p_count
546 property fecv_size:
547 """Length of the full excitation coefficient vector"""
548 def __get__(self):
549 self.check_s()
550 return self.s[0].fecv_size
551 property saecv_sizes:
552 """Lengths of the partial symmetry-adapted excitation coefficient vectors"""
553 def __get__(self):
554 self.check_s()
555 return [self.s[0].saecv_sizes[i]
556 for i in range(self.s[0].sym[0].nirreps)]
557 property irrep_names:
558 def __get__(self):
559 self.check_s()
560 return [string_c2py(self.s[0].sym[0].irreps[iri].name)
561 if (self.s[0].sym[0].irreps[iri].name) else None
562 for iri in range(self.s[0].sym[0].nirreps)]
563 property nirreps:
564 """Number of irreducible representations of the scattering system point symmetry group"""
565 def __get__(self):
566 self.check_s()
567 return self.s[0].sym[0].nirreps
568 property lattice_dimension:
569 def __get__(self):
570 return self.s[0].lattice_dimension
572 property lattice_basis:
573 """Direct lattice basis vectors for periodic system in 3D cartesian coordinates"""
574 def __get__(self):
575 cdef int d = self.lattice_dimension
576 if d == 0:
577 return None
578 else:
579 return np.array([[self.s[0].per.lattice_basis[i].x,
580 self.s[0].per.lattice_basis[i].y,
581 self.s[0].per.lattice_basis[i].z,
582 ] for i in range(d)])
584 property reciprocal_basis:
585 """Reciprocal lattice basis vectors for periodic system in 3D cartesian coordinates"""
586 def __get__(self):
587 cdef int d = self.lattice_dimension
588 if d == 0:
589 return None
590 else:
591 return np.array([[self.s[0].per.reciprocal_basis[i].x,
592 self.s[0].per.reciprocal_basis[i].y,
593 self.s[0].per.reciprocal_basis[i].z,
594 ] for i in range(d)])
596 def bspec_pi(self, qpms_ss_pi_t pi):
597 """Gives the BaseSpec of a particle defined by an index in the internal particles' order, in the form of an array that can be fed to BaseSpec constructor"""
598 self.check_s()
599 cdef qpms_ss_pi_t pcount = self.s[0].p_count
600 if pi >= pcount or pi < 0: raise IndexError("Invalid particle index")
601 cdef const qpms_vswf_set_spec_t *bspec = qpms_ss_bspec_pi(self.s, pi)
602 return np.array([ bspec[0].ilist[i] for i in range(bspec[0].n)])
603 property bspecs:
604 """Gives the BaseSpecs of all particles in the internal particles' order, in the form of arrays that can be fed to BaseSpec constructor"""
605 def __get__(self):
606 self.check_s()
607 cdef qpms_ss_pi_t pi, pcount = self.s[0].p_count
608 return [self.bspec_pi(pi) for pi in range(pcount)]
610 property unitcell_volume:
611 def __get__(self):
612 self.check_s()
613 if self.lattice_dimension:
614 return self.s[0].per.unitcell_volume
615 else:
616 return None
618 property eta:
619 """Ewald parameter η (only relevant for periodic systems)
621 This is the default value that will be copied into
622 _ScatteringSystemAtOmega by __call__"""
623 def __get__(self):
624 self.check_s()
625 if self.lattice_dimension:
626 return self.s[0].per.eta
627 else:
628 return None
630 def __set__(self, eta):
631 self.check_s()
632 if self.lattice_dimension:
633 self.s[0].per.eta = eta
634 else:
635 raise AttributeError("Cannot set Ewald parameter for finite system") # different exception?
638 def pack_vector(self, vect, iri):
639 """Converts (projects) a full excitation coefficient vector into an irrep subspace.
641 Parameters
642 ----------
643 vect : array_like of shape (self.fecv_size,)
644 The full excitation coefficient vector to be converted.
645 iri : int
646 Index of the irreducible representation.
648 Returns
649 -------
650 packed : ndarray of shape (self.saecv_sizes[iri],)
651 "irrep-packed" excitation vector: the part of vect belonging to the
652 irrep indexed by iri, in terms of given irrep's internal coordinates.
654 See Also
655 --------
656 unpack_vector : The beckward (not exactly inverse) operation.
657 pack_matrix : Corresponding conversion for matrices.
659 self.check_s()
660 if len(vect) != self.fecv_size:
661 raise ValueError("Length of a full vector has to be %d, not %d"
662 % (self.fecv_size, len(vect)))
663 vect = np.array(vect, dtype=complex, copy=False, order='C')
664 cdef cdouble[::1] vect_view = vect;
665 cdef np.ndarray[np.complex_t, ndim=1] target_np = np.empty(
666 (self.saecv_sizes[iri],), dtype=complex, order='C')
667 cdef cdouble[::1] target_view = target_np
668 qpms_scatsys_irrep_pack_vector(&target_view[0], &vect_view[0], self.s, iri)
669 return target_np
671 def unpack_vector(self, packed, iri):
672 """Unpacks an "irrep-packed" excitation coefficient vector to full coordinates.
674 Parameters
675 ----------
676 packed : array_like of shape (self.saecv_sizes[iri],)
677 Excitation coefficient vector part belonging to the irreducible representation
678 indexed by iri, in terms of given irrep's internal coordinates.
679 iri : int
680 Index of the irreducible representation.
682 Returns
683 -------
684 vect : ndarray of shape (self.fecv_size,)
685 The contribution to a full excitation coefficient vector from iri'th irrep.
687 See Also
688 --------
689 pack_vector : The inverse operation.
690 unpack_matrix : Corresponding conversion for matrices.
692 self.check_s()
693 if len(packed) != self.saecv_sizes[iri]:
694 raise ValueError("Length of %d. irrep-packed vector has to be %d, not %d"
695 % (iri, self.saecv_sizes, len(packed)))
696 packed = np.array(packed, dtype=complex, copy=False, order='C')
697 cdef cdouble[::1] packed_view = packed
698 cdef np.ndarray[np.complex_t, ndim=1] target_np = np.empty(
699 (self.fecv_size,), dtype=complex)
700 cdef cdouble[::1] target_view = target_np
701 qpms_scatsys_irrep_unpack_vector(&target_view[0], &packed_view[0],
702 self.s, iri, 0)
703 return target_np
705 def pack_matrix(self, fullmatrix, iri):
706 """Converts (projects) a matrix into an irrep subspace.
708 Parameters
709 ----------
710 fullmatrix : array_like of shape (self.fecv_size, self.fecv_size)
711 The full matrix (operating on the exctitation coefficient vectors) to be converted.
712 iri : int
713 Index of the irreducible representation.
715 Returns
716 -------
717 packedmatrix : ndarray of shape (self.saecv_sizes[iri], self.saecv_sizes[iri])
718 "irrep-packed" matrix: the part of fullmatrix belonging to the
719 irrep indexed by iri, in terms of given irrep's internal coordinates.
721 See Also
722 --------
723 unpack_matrix : The beckward (not exactly inverse) operation.
724 pack_vector : Corresponding conversion for excitation coefficient vectors.
726 self.check_s()
727 cdef size_t flen = self.s[0].fecv_size
728 cdef size_t rlen = self.saecv_sizes[iri]
729 fullmatrix = np.array(fullmatrix, dtype=complex, copy=False, order='C')
730 if fullmatrix.shape != (flen, flen):
731 raise ValueError("Full matrix shape should be (%d,%d), is %s."
732 % (flen, flen, repr(fullmatrix.shape)))
733 cdef cdouble[:,::1] fullmatrix_view = fullmatrix
734 cdef np.ndarray[np.complex_t, ndim=2] target_np = np.empty(
735 (rlen, rlen), dtype=complex, order='C')
736 cdef cdouble[:,::1] target_view = target_np
737 qpms_scatsys_irrep_pack_matrix(&target_view[0][0], &fullmatrix_view[0][0],
738 self.s, iri)
739 return target_np
740 def unpack_matrix(self, packedmatrix, iri):
741 """Unpacks an "irrep-packed" excitation coefficient vector to full coordinates.
743 Parameters
744 ----------
745 packedmatrix : array_like of shape (self.saecv_sizes[iri], self.saecv_sizes[iri])
746 A matrix (operating on the excitation coefficient vectors) part belonging to the
747 irreducible representation indexed by iri, in terms of given irrep's internal
748 coordinates.
749 iri : int
750 Index of the irreducible representation.
752 Returns
753 -------
754 fullmatrix : ndarray of shape (self.fecv_size, self.fecv_size)
755 The iri'th irrep contribution to a full matrix.
757 See Also
758 --------
759 pack_matrix : The inverse operation.
760 unpack_vector : Corresponding conversion for excitation coefficient vectors.
762 self.check_s()
763 cdef size_t flen = self.s[0].fecv_size
764 cdef size_t rlen = self.saecv_sizes[iri]
765 packedmatrix = np.array(packedmatrix, dtype=complex, copy=False, order='C')
766 if packedmatrix.shape != (rlen, rlen):
767 raise ValueError("Packed matrix shape should be (%d,%d), is %s."
768 % (rlen, rlen, repr(packedmatrix.shape)))
769 cdef cdouble[:,::1] packedmatrix_view = packedmatrix
770 cdef np.ndarray[np.complex_t, ndim=2] target_np = np.empty(
771 (flen, flen), dtype=complex, order='C')
772 cdef cdouble[:,::1] target_view = target_np
773 qpms_scatsys_irrep_unpack_matrix(&target_view[0][0], &packedmatrix_view[0][0],
774 self.s, iri, 0)
775 return target_np
777 def translation_matrix_full(self, cdouble wavenumber, blochvector = None, J = QPMS_HANKEL_PLUS, double eta = 0):
778 """Constructs full translation matrix of a scattering system.
780 This method enables to use any wave number for the background medium ignoring the
781 background EpsMuGenerator), using only system's particle positions (and lattice
782 basis for infinite system).
784 Parameters
785 ----------
786 wavenumber : complex
787 Wave number of the medium
789 blochvector : array_like of shape (3,)
790 Bloch vector (only for periodic systems)
792 J : BesselType
793 Optionally, one can replace Hankel functions of the first kind with different
794 Bessel functions.
796 eta : float
797 Ewald parameter η, only relevant for periodic lattices;
798 if set to 0 (default), the actual value is determined
799 automatically.
801 See Also
802 --------
803 ScatteringSystemAtOmega.translation_matrix_full : Translation matrix at a given frequency.
805 self.check_s()
806 cdef size_t flen = self.s[0].fecv_size
807 cdef np.ndarray[np.complex_t, ndim=2] target = np.empty(
808 (flen,flen),dtype=complex, order='C')
809 cdef cdouble[:,::1] target_view = target
810 cdef cart3_t blochvector_c
811 if self.lattice_dimension == 0:
812 if blochvector is None:
813 qpms_scatsys_build_translation_matrix_e_full(&target_view[0][0], self.s, wavenumber, J)
814 else: raise ValueError("Can't use blochvector with non-periodic system")
815 else:
816 if blochvector is None: raise ValueError("Valid blochvector must be specified for periodic system")
817 else:
818 if J != QPMS_HANKEL_PLUS:
819 raise NotImplementedError("Translation operators based on other than Hankel+ functions not supperted in periodic systems")
820 blochvector_c = {'x': blochvector[0], 'y': blochvector[1], 'z': blochvector[2]}
821 with pgsl_ignore_error(15):
822 qpms_scatsys_periodic_build_translation_matrix_full(&target_view[0][0], self.s, wavenumber, &blochvector_c, eta)
823 return target
825 def translation_matrix_packed(self, cdouble wavenumber, qpms_iri_t iri, J = QPMS_HANKEL_PLUS):
826 self.check_s()
827 cdef size_t rlen = self.saecv_sizes[iri]
828 cdef np.ndarray[np.complex_t, ndim=2] target = np.empty(
829 (rlen,rlen),dtype=complex, order='C')
830 cdef cdouble[:,::1] target_view = target
831 qpms_scatsys_build_translation_matrix_e_irrep_packed(&target_view[0][0],
832 self.s, iri, wavenumber, J)
833 return target
835 property fullvec_psizes:
836 def __get__(self):
837 self.check_s()
838 cdef np.ndarray[int32_t, ndim=1] ar = np.empty((self.s[0].p_count,), dtype=np.int32)
839 cdef int32_t[::1] ar_view = ar
840 for pi in range(self.s[0].p_count):
841 ar_view[pi] = qpms_ss_bspec_pi(self.s, pi)[0].n
842 return ar
845 property fullvec_poffsets:
846 def __get__(self):
847 self.check_s()
848 cdef np.ndarray[intptr_t, ndim=1] ar = np.empty((self.s[0].p_count,), dtype=np.intp)
849 cdef intptr_t[::1] ar_view = ar
850 cdef intptr_t offset = 0
851 for pi in range(self.s[0].p_count):
852 ar_view[pi] = offset
853 offset += qpms_ss_bspec_pi(self.s, pi)[0].n
854 return ar
856 property positions:
857 def __get__(self):
858 self.check_s()
859 cdef np.ndarray[np.double_t, ndim=2] ar = np.empty((self.s[0].p_count, 3), dtype=float)
860 cdef np.double_t[:,::1] ar_view = ar
861 for pi in range(self.s[0].p_count):
862 ar_view[pi,0] = self.s[0].p[pi].pos.x
863 ar_view[pi,1] = self.s[0].p[pi].pos.y
864 ar_view[pi,2] = self.s[0].p[pi].pos.z
865 return ar
867 def planewave_full(self, k_cart, E_cart):
868 self.check_s()
869 k_cart = np.array(k_cart)
870 E_cart = np.array(E_cart)
871 if k_cart.shape != (3,) or E_cart.shape != (3,):
872 raise ValueError("k_cart and E_cart must be ndarrays of shape (3,)")
873 cdef qpms_incfield_planewave_params_t p
874 p.use_cartesian = 1
875 p.k.cart.x = <cdouble>k_cart[0]
876 p.k.cart.y = <cdouble>k_cart[1]
877 p.k.cart.z = <cdouble>k_cart[2]
878 p.E.cart.x = <cdouble>E_cart[0]
879 p.E.cart.y = <cdouble>E_cart[1]
880 p.E.cart.z = <cdouble>E_cart[2]
881 cdef np.ndarray[np.complex_t, ndim=1] target_np = np.empty(
882 (self.fecv_size,), dtype=complex)
883 cdef cdouble[::1] target_view = target_np
884 qpms_scatsys_incident_field_vector_full(&target_view[0],
885 self.s, qpms_incfield_planewave, <void *>&p, 0)
886 return target_np
888 def find_modes(self, cdouble omega_centre, double omega_rr, double omega_ri, iri = None,
889 blochvector = None,
890 size_t contour_points = 20, double rank_tol = 1e-4, size_t rank_min_sel=1,
891 double res_tol = 0):
893 Attempts to find the eigenvalues and eigenvectors using Beyn's algorithm.
896 cdef beyn_result_t *res
897 cdef double blochvector_c[3]
898 if self.lattice_dimension == 0:
899 assert(blochvector is None)
900 res = qpms_scatsys_finite_find_eigenmodes(self.s, self.iri_py2c(iri),
901 omega_centre, omega_rr, omega_ri, contour_points,
902 rank_tol, rank_min_sel, res_tol)
903 else:
904 if iri is not None: raise NotImplementedError("Irreps decomposition not yet supported for periodic systems")
905 blochvector_c[0] = blochvector[0]
906 blochvector_c[1] = blochvector[1]
907 blochvector_c[2] = blochvector[2] if len(blochvector) > 2 else 0
908 res = qpms_scatsys_periodic_find_eigenmodes(self.s, blochvector_c,
909 omega_centre, omega_rr, omega_ri, contour_points, rank_tol, rank_min_sel, res_tol)
910 if res == NULL: raise RuntimeError
912 cdef size_t neig = res[0].neig, i, j
913 cdef size_t vlen = res[0].vlen # should be equal to self.s.fecv_size
915 cdef np.ndarray[complex, ndim=1] eigval = np.empty((neig,), dtype=complex)
916 cdef cdouble[::1] eigval_v = eigval
917 cdef np.ndarray[complex, ndim=1] eigval_err = np.empty((neig,), dtype=complex)
918 cdef cdouble[::1] eigval_err_v = eigval_err
919 cdef np.ndarray[double, ndim=1] residuals = np.empty((neig,), dtype=np.double)
920 cdef double[::1] residuals_v = residuals
921 cdef np.ndarray[complex, ndim=2] eigvec = np.empty((neig,vlen),dtype=complex)
922 cdef cdouble[:,::1] eigvec_v = eigvec
923 cdef np.ndarray[double, ndim=1] ranktest_SV = np.empty((vlen), dtype=np.double)
924 cdef double[::1] ranktest_SV_v = ranktest_SV
926 for i in range(neig):
927 eigval_v[i] = res[0].eigval[i]
928 eigval_err_v[i] = res[0].eigval_err[i]
929 residuals_v[i] = res[0].residuals[i]
930 for j in range(vlen):
931 eigvec_v[i,j] = res[0].eigvec[i*vlen + j]
932 for i in range(vlen):
933 ranktest_SV_v[i] = res[0].ranktest_SV[i]
935 zdist = eigval - omega_centre
936 eigval_inside_metric = np.hypot(zdist.real / omega_rr, zdist.imag / omega_ri)
938 beyn_result_free(res)
939 retdict = {
940 'eigval':eigval,
941 'eigval_inside_metric':eigval_inside_metric,
942 'eigvec':eigvec,
943 'residuals':residuals,
944 'eigval_err':eigval_err,
945 'ranktest_SV':ranktest_SV,
946 'iri': iri,
947 'blochvector': blochvector
950 return retdict
952 @boundscheck(False)
953 def scattered_E(self, cdouble wavenumber, scatcoeffvector_full, evalpos, btyp=BesselType.HANKEL_PLUS, bint alt=False):
954 # TODO DOC, periodic systems
955 if self.lattice_dimension != 0: # TODO enable also periodic systems
956 raise NotImplementedError("For periodic system, use _ScatteringSystemAtOmegaK.scattered_E")
957 cdef qpms_bessel_t btyp_c = BesselType(btyp)
958 evalpos = np.array(evalpos, dtype=float, copy=False)
959 if evalpos.shape[-1] != 3:
960 raise ValueError("Last dimension of evalpos has to be 3")
961 cdef np.ndarray[double,ndim=2] evalpos_a = evalpos.reshape(-1,3)
962 cdef np.ndarray[dtype=complex, ndim=1] scv = np.array(scatcoeffvector_full, copy=False, dtype=complex)
963 cdef cdouble[::1] scv_view = scv
964 cdef np.ndarray[complex, ndim=2] results = np.empty((evalpos_a.shape[0],3), dtype=complex)
965 cdef ccart3_t res
966 cdef cart3_t pos
967 cdef Py_ssize_t i
968 with nogil, parallel(), boundscheck(False), wraparound(False):
969 for i in prange(evalpos_a.shape[0]):
970 pos.x = evalpos_a[i,0]
971 pos.y = evalpos_a[i,1]
972 pos.z = evalpos_a[i,2]
973 if alt:
974 res = qpms_scatsys_scattered_E__alt(self.s, btyp_c, wavenumber, &scv_view[0], pos)
975 else:
976 res = qpms_scatsys_scattered_E(self.s, btyp_c, wavenumber, &scv_view[0], pos)
977 results[i,0] = res.x
978 results[i,1] = res.y
979 results[i,2] = res.z
980 return results.reshape(evalpos.shape)
982 def empty_lattice_modes_xy(EpsMu epsmu, reciprocal_basis, wavevector, double maxomega):
983 '''Empty (2D, xy-plane) lattice mode (diffraction order) frequencies of a non-dispersive medium.
985 reciprocal_basis is of the "mutliplied by 2Đż" type.
987 cdef double *omegas_c
988 cdef size_t n
989 cdef cart2_t k_, b1, b2
990 k_.x = wavevector[0]
991 k_.y = wavevector[1]
992 b1.x = reciprocal_basis[0][0]
993 b1.y = reciprocal_basis[0][1]
994 b2.x = reciprocal_basis[1][0]
995 b2.y = reciprocal_basis[1][1]
996 if(epsmu.n.imag != 0):
997 warnings.warn("Got complex refractive index", epsmu.n, "ignoring the imaginary part")
998 refindex = epsmu.n.real
999 n = qpms_emptylattice2_modes_maxfreq(&omegas_c, b1, b2, BASIS_RTOL,
1000 k_, GSL_CONST_MKSA_SPEED_OF_LIGHT / refindex, maxomega)
1001 cdef np.ndarray[double, ndim=1] omegas = np.empty((n,), dtype=np.double)
1002 cdef double[::1] omegas_v = omegas
1003 cdef size_t i
1004 for i in range(n):
1005 omegas_v[i] = omegas_c[i]
1006 free(omegas_c)
1007 return omegas
1009 cdef class _ScatteringSystemAtOmegaK:
1011 Wrapper over the C qpms_scatsys_at_omega_k_t structure
1013 This represents an infinite periodic system with a given (fixed) frequency and Bloch vector.
1015 cdef qpms_scatsys_at_omega_k_t sswk
1016 cdef _ScatteringSystemAtOmega ssw_pyref
1018 cdef qpms_scatsys_at_omega_k_t *rawpointer(self):
1019 return &self.sswk
1021 property eta:
1022 """Ewald parameter η"""
1023 def __get__(self):
1024 return self.sswk.eta
1025 def __set__(self, double eta):
1026 self.sswk.eta = eta
1028 @boundscheck(False)
1029 def scattered_E(self, scatcoeffvector_full, evalpos, btyp=QPMS_HANKEL_PLUS):
1030 """Evaluate electric field for a given excitation coefficient vector (periodic system)
1032 Parameters
1033 ----------
1034 scatcoeffvector_full: array_like of shape (self.fecv_size,)
1035 Onedimensional excitation coefficient vector, describing SVWFs leaving
1036 the scatterers inside the representative unit cell.
1037 evalpos: array_like of floats and shape (..., 3)
1038 Evaluation points in cartesian coordinates.
1040 Returns
1041 -------
1042 ndarray of complex, with the same shape as `evalpos`
1043 Electric field at the positions given in `evalpos`, in cartesian coordinates.
1045 if(btyp != QPMS_HANKEL_PLUS):
1046 raise NotImplementedError("Only first kind Bessel function-based fields are supported")
1047 cdef qpms_bessel_t btyp_c = BesselType(btyp)
1048 evalpos = np.array(evalpos, dtype=float, copy=False)
1049 if evalpos.shape[-1] != 3:
1050 raise ValueError("Last dimension of evalpos has to be 3")
1051 cdef np.ndarray[double,ndim=2] evalpos_a = evalpos.reshape(-1,3)
1052 cdef np.ndarray[dtype=complex, ndim=1] scv = np.array(scatcoeffvector_full, dtype=complex, copy=False)
1053 cdef cdouble[::1] scv_view = scv
1054 cdef np.ndarray[complex, ndim=2] results = np.empty((evalpos_a.shape[0],3), dtype=complex)
1055 cdef ccart3_t res
1056 cdef cart3_t pos
1057 cdef Py_ssize_t i
1058 with nogil, wraparound(False), parallel():
1059 for i in prange(evalpos_a.shape[0]):
1060 pos.x = evalpos_a[i,0]
1061 pos.y = evalpos_a[i,1]
1062 pos.z = evalpos_a[i,2]
1063 res = qpms_scatsyswk_scattered_E(&self.sswk, btyp_c, &scv_view[0], pos)
1064 results[i,0] = res.x
1065 results[i,1] = res.y
1066 results[i,2] = res.z
1067 return results.reshape(evalpos.shape)
1069 def scattered_field_basis(self, evalpos, btyp=QPMS_HANKEL_PLUS):
1070 # TODO examples
1071 """Evaluate scattered field "basis" (periodic system)
1073 This function enables the evaluation of "scattered" fields
1074 generated by the system for many different excitation
1075 coefficients vectors, without the expensive re-evaluation of the respective
1076 translation operators for each excitation coefficient vector.
1078 Parameters
1079 ----------
1080 evalpos: array_like of floats and shape (..., 3)
1081 Evaluation points in cartesian coordinates.
1083 Returns
1084 -------
1085 ndarray of complex, with the shape `evalpos.shape[:-1] + (self.fecv_size, 3)`
1086 "Basis" fields at the positions given in `evalpos`, in cartesian coordinates.
1088 if(btyp != QPMS_HANKEL_PLUS):
1089 raise NotImplementedError("Only first kind Bessel function-based fields are supported")
1090 cdef qpms_bessel_t btyp_c = BesselType(btyp)
1091 cdef Py_ssize_t fecv_size = self.fecv_size
1092 evalpos = np.array(evalpos, dtype=float, copy=False)
1093 if evalpos.shape[-1] != 3:
1094 raise ValueError("Last dimension of evalpos has to be 3")
1095 cdef np.ndarray[double,ndim=2] evalpos_a = evalpos.reshape(-1,3)
1096 cdef np.ndarray[complex, ndim=3] results = np.empty((evalpos_a.shape[0], fecv_size, 3), dtype=complex)
1097 cdef ccart3_t *res
1098 res = <ccart3_t *> malloc(fecv_size*sizeof(ccart3_t))
1099 cdef cart3_t pos
1100 cdef Py_ssize_t i, j
1101 with nogil, wraparound(False), parallel():
1102 for i in prange(evalpos_a.shape[0]):
1103 pos.x = evalpos_a[i,0]
1104 pos.y = evalpos_a[i,1]
1105 pos.z = evalpos_a[i,2]
1106 qpms_scatsyswk_scattered_field_basis(res, &self.sswk, btyp_c, pos)
1107 for j in range(fecv_size):
1108 results[i,j,0] = res[j].x
1109 results[i,j,1] = res[j].y
1110 results[i,j,2] = res[j].z
1111 free(res)
1112 return results.reshape(evalpos.shape[:-1] + (self.fecv_size, 3))
1114 property fecv_size:
1115 def __get__(self): return self.ssw_pyref.ss_pyref.fecv_size
1117 cdef class _ScatteringSystemAtOmega:
1119 Wrapper over the C qpms_scatsys_at_omega_t structure
1120 that keeps the T-matrix and background data evaluated
1121 at specific frequency.
1123 cdef qpms_scatsys_at_omega_t *ssw
1124 cdef ScatteringSystem ss_pyref
1126 def check(self): # cdef instead?
1127 if not self.ssw:
1128 raise ValueError("_ScatteringSystemAtOmega's ssw-pointer not set. You must not use the default constructor; ScatteringSystem.create() instead")
1129 self.ss_pyref.check_s()
1130 #TODO is there a way to disable the constructor outside this module?
1133 def ensure_finite(self):
1134 if self.ssw[0].ss[0].lattice_dimension != 0:
1135 raise NotImplementedError("Operation not supported for periodic systems")
1137 def __dealloc__(self):
1138 if (self.ssw):
1139 qpms_scatsys_at_omega_free(self.ssw)
1141 def apply_Tmatrices_full(self, a):
1142 self.check()
1143 if len(a) != self.fecv_size:
1144 raise ValueError("Length of a full vector has to be %d, not %d"
1145 % (self.fecv_size, len(a)))
1146 a = np.array(a, dtype=complex, copy=False, order='C')
1147 cdef cdouble[::1] a_view = a;
1148 cdef np.ndarray[np.complex_t, ndim=1] target_np = np.empty(
1149 (self.fecv_size,), dtype=complex, order='C')
1150 cdef cdouble[::1] target_view = target_np
1151 qpms_scatsysw_apply_Tmatrices_full(&target_view[0], &a_view[0], self.ssw)
1152 return target_np
1154 cdef qpms_scatsys_at_omega_t *rawpointer(self):
1155 return self.ssw
1157 def scatter_solver(self, iri=None, k=None):
1158 self.check()
1159 cdef _ScatteringSystemAtOmegaK sswk # used only for periodic systems
1160 if(self.ssw[0].ss[0].lattice_dimension == 0):
1161 return ScatteringMatrix(self, iri=iri)
1162 else:
1163 if iri is not None:
1164 raise NotImplementedError("Irrep decomposition not (yet) supported for periodic systems")
1165 sswk = self._sswk(k)
1166 return ScatteringMatrix(ssw=self, sswk=sswk, iri=None)
1168 def _sswk(self, k):
1169 cdef _ScatteringSystemAtOmegaK sswk = _ScatteringSystemAtOmegaK()
1170 sswk.sswk.ssw = self.ssw
1171 sswk.ssw_pyref=self
1172 sswk.sswk.k[0] = k[0]
1173 sswk.sswk.k[1] = k[1]
1174 sswk.sswk.k[2] = k[2]
1175 sswk.eta = qpms_ss_adjusted_eta(self.ssw[0].ss, self.ssw[0].wavenumber, sswk.sswk.k)
1176 return sswk
1178 property fecv_size:
1179 def __get__(self): return self.ss_pyref.fecv_size
1180 property saecv_sizes:
1181 def __get__(self): return self.ss_pyref.saecv_sizes
1182 property irrep_names:
1183 def __get__(self): return self.ss_pyref.irrep_names
1184 property nirreps:
1185 def __get__(self): return self.ss_pyref.nirreps
1186 property wavenumber:
1187 def __get__(self): return self.ssw[0].wavenumber
1190 def modeproblem_matrix_full(self, k=None):
1191 self.check()
1192 cdef size_t flen = self.ss_pyref.s[0].fecv_size
1193 cdef np.ndarray[np.complex_t, ndim=2] target = np.empty(
1194 (flen,flen),dtype=complex, order='C')
1195 cdef cdouble[:,::1] target_view = target
1196 cdef _ScatteringSystemAtOmegaK sswk
1197 if k is not None:
1198 sswk = self._sswk(k)
1199 qpms_scatsyswk_build_modeproblem_matrix_full(&target_view[0][0], &sswk.sswk)
1200 else:
1201 qpms_scatsysw_build_modeproblem_matrix_full(&target_view[0][0], self.ssw)
1202 return target
1204 def modeproblem_matrix_packed(self, qpms_iri_t iri, version='pR'):
1205 self.check()
1206 self.ensure_finite()
1207 cdef size_t rlen = self.saecv_sizes[iri]
1208 cdef np.ndarray[np.complex_t, ndim=2] target = np.empty(
1209 (rlen,rlen),dtype=complex, order='C')
1210 cdef cdouble[:,::1] target_view = target
1211 if (version == 'R'):
1212 qpms_scatsysw_build_modeproblem_matrix_irrep_packed_orbitorderR(&target_view[0][0], self.ssw, iri)
1213 elif (version == 'pR'):
1214 with nogil:
1215 qpms_scatsysw_build_modeproblem_matrix_irrep_packed(&target_view[0][0], self.ssw, iri)
1216 else:
1217 qpms_scatsysw_build_modeproblem_matrix_irrep_packed_serial(&target_view[0][0], self.ssw, iri)
1218 return target
1220 def translation_matrix_full(self, blochvector = None):
1221 """Constructs full translation matrix of a scattering system at given frequency.
1223 Parameters
1224 ----------
1225 blochvector : array_like of shape (3,)
1226 Bloch vector (only for periodic systems)
1228 See Also
1229 --------
1230 ScatteringSystem.translation_matrix_full: translation matrix for any wavenumber
1233 return self.ss_pyref.translation_matrix_full(wavenumber=self.wavenumber, blochvector=blochvector)
1235 def translation_matrix_packed(self, iri, J = QPMS_HANKEL_PLUS):
1236 return self.ss_pyref.translation_matrix_packed(wavenumber=self.wavenumber, iri=iri, J=J)
1238 @boundscheck(False)
1239 def scattered_E(self, scatcoeffvector_full, evalpos, blochvector=None, btyp=QPMS_HANKEL_PLUS, bint alt=False):
1240 # FIXME TODO this obviously does not work for periodic systems
1241 """Evaluate electric field for a given excitation coefficient vector
1243 Parameters
1244 ----------
1245 scatcoeffvector_full: array_like of shape (self.fecv_size,)
1246 Onedimensional excitation coefficient vector, describing SVWFs leaving
1247 the scatterers.
1248 evalpos: array_like of floats and shape (..., 3)
1249 Evaluation points in cartesian coordinates.
1250 blochvector: array_like or None
1251 Bloch vector, must be supplied (non-None) for periodic systems, else None.
1253 Returns
1254 -------
1255 array_like of complex, with the same shape as `evalpos`
1256 Electric field at the positions given in `evalpos`.
1258 See Also
1259 --------
1260 _ScatteringSystemAtOmegaK.scattered_E: Evaluate scattered field in an infinite periodic system.
1262 if(self.ssw[0].ss[0].lattice_dimension != 0):
1263 if blochvector is None:
1264 raise ValueError("For a periodic system, blochvector must be specified.")
1265 return self._sswk(blochvector).scattered_E(scatcoeffvector_full, evalpos=evalpos, btyp=btyp)
1266 cdef qpms_bessel_t btyp_c = BesselType(btyp)
1267 evalpos = np.array(evalpos, dtype=float, copy=False)
1268 if evalpos.shape[-1] != 3:
1269 raise ValueError("Last dimension of evalpos has to be 3")
1270 cdef np.ndarray[double,ndim=2] evalpos_a = evalpos.reshape(-1,3)
1271 cdef np.ndarray[dtype=complex, ndim=1] scv = np.array(scatcoeffvector_full, dtype=complex, copy=False)
1272 cdef cdouble[::1] scv_view = scv
1273 cdef np.ndarray[complex, ndim=2] results = np.empty((evalpos_a.shape[0],3), dtype=complex)
1274 cdef ccart3_t res
1275 cdef cart3_t pos
1276 cdef Py_ssize_t i
1277 with wraparound(False), nogil, parallel():
1278 for i in prange(evalpos_a.shape[0]):
1279 pos.x = evalpos_a[i,0]
1280 pos.y = evalpos_a[i,1]
1281 pos.z = evalpos_a[i,2]
1282 if alt:
1283 res = qpms_scatsysw_scattered_E__alt(self.ssw, btyp_c, &scv_view[0], pos)
1284 else:
1285 res = qpms_scatsysw_scattered_E(self.ssw, btyp_c, &scv_view[0], pos)
1286 results[i,0] = res.x
1287 results[i,1] = res.y
1288 results[i,2] = res.z
1289 return results.reshape(evalpos.shape)
1292 cdef class ScatteringMatrix:
1294 Wrapper over the C qpms_ss_LU structure that keeps the factorised mode problem matrix.
1296 cdef _ScatteringSystemAtOmega ssw # Here we keep the reference to the parent scattering system
1297 cdef _ScatteringSystemAtOmegaK sswk
1298 cdef qpms_ss_LU lu
1300 def __cinit__(self, _ScatteringSystemAtOmega ssw, sswk=None, iri=None):
1301 ssw.check()
1302 self.ssw = ssw
1303 if sswk is None:
1304 ssw.ensure_finite()
1305 # TODO? pre-allocate the matrix with numpy to make it transparent?
1306 if iri is None:
1307 self.lu = qpms_scatsysw_build_modeproblem_matrix_full_LU(
1308 NULL, NULL, ssw.rawpointer())
1309 else:
1310 self.lu = qpms_scatsysw_build_modeproblem_matrix_irrep_packed_LU(
1311 NULL, NULL, ssw.rawpointer(), iri)
1312 else:
1313 # TODO check sswk validity
1314 self.sswk = sswk
1315 self.lu = qpms_scatsyswk_build_modeproblem_matrix_full_LU(NULL, NULL, self.sswk.rawpointer())
1317 def __dealloc__(self):
1318 qpms_ss_LU_free(self.lu)
1320 property iri:
1321 def __get__(self):
1322 return None if self.lu.full else self.lu.iri
1324 def __call__(self, a_inc):
1325 cdef size_t vlen
1326 cdef qpms_iri_t iri = -1;
1327 if self.lu.full:
1328 vlen = self.lu.ssw[0].ss[0].fecv_size
1329 if len(a_inc) != vlen:
1330 raise ValueError("Length of a full coefficient vector has to be %d, not %d"
1331 % (vlen, len(a_inc)))
1332 else:
1333 iri = self.lu.iri
1334 vlen = self.lu.ssw[0].ss[0].saecv_sizes[iri]
1335 if len(a_inc) != vlen:
1336 raise ValueError("Length of a %d. irrep packed coefficient vector has to be %d, not %d"
1337 % (iri, vlen, len(a_inc)))
1338 a_inc = np.array(a_inc, dtype=complex, copy=False, order='C')
1339 cdef const cdouble[::1] a_view = a_inc;
1340 cdef np.ndarray f = np.empty((vlen,), dtype=complex, order='C')
1341 cdef cdouble[::1] f_view = f
1342 qpms_scatsys_scatter_solve(&f_view[0], &a_view[0], self.lu)
1343 return f
1345 def pitau(double theta, qpms_l_t lMax, double csphase = -1):
1346 if(abs(csphase) != 1):
1347 raise ValueError("csphase must be 1 or -1, is %g" % csphase)
1348 cdef size_t nelem = qpms_lMax2nelem(lMax)
1349 cdef np.ndarray[np.float_t, ndim=1] lega = np.empty((nelem,), dtype=float)
1350 cdef np.ndarray[np.float_t, ndim=1] pia = np.empty((nelem,), dtype=float)
1351 cdef np.ndarray[np.float_t, ndim=1] taua = np.empty((nelem,), dtype=float)
1352 cdef double[::1] leg = lega
1353 cdef double[::1] pi = pia
1354 cdef double[::1] tau = taua
1355 qpms_pitau_fill(&leg[0], &pi[0], &tau[0], theta, lMax, csphase)
1356 return (lega, pia, taua)
1358 def linton_gamma(cdouble x):
1359 return clilgamma(x)
1361 def linton_gamma_real(double x):
1362 return lilgamma(x)
1364 def gamma_inc(double a, cdouble x, int m = 0):
1365 cdef qpms_csf_result res
1366 with pgsl_ignore_error(15): #15 is underflow
1367 complex_gamma_inc_e(a, x, m, &res)
1368 return (res.val, res.err)
1370 def gamma_inc_series(double a, cdouble x):
1371 cdef qpms_csf_result res
1372 with pgsl_ignore_error(15): #15 is underflow
1373 cx_gamma_inc_series_e(a, x, &res)
1374 return (res.val, res.err)
1376 def gamma_inc_CF(double a, cdouble x):
1377 cdef qpms_csf_result res
1378 with pgsl_ignore_error(15): #15 is underflow
1379 cx_gamma_inc_CF_e(a, x, &res)
1380 return (res.val, res.err)
1382 def lll_reduce(basis, double delta=0.75):
1384 Lattice basis reduction with the Lenstra-Lenstra-Lovász algorithm.
1386 basis is array_like with dimensions (n, d), where
1387 n is the size of the basis (dimensionality of the lattice)
1388 and d is the dimensionality of the space into which the lattice
1389 is embedded.
1391 basis = np.array(basis, copy=True, order='C', dtype=np.double)
1392 if len(basis.shape) != 2:
1393 raise ValueError("Expected two-dimensional array (got %d-dimensional)"
1394 % len(basis.shape))
1395 cdef size_t n, d
1396 n, d = basis.shape
1397 if n > d:
1398 raise ValueError("Real space dimensionality (%d) cannot be smaller than"
1399 "the dimensionality of the lattice (%d) embedded into it."
1400 % (d, n))
1401 cdef double [:,:] basis_view = basis
1402 if 0 != qpms_reduce_lattice_basis(&basis_view[0,0], n, d, delta):
1403 raise RuntimeError("Something weird happened")
1404 return basis
1407 cdef PGen get_PGen_direct(direct_basis, bint include_origin=False, double layers=30):
1408 dba = np.array(direct_basis)
1409 if not (dba.shape == (2,2)):
1410 raise NotImplementedError
1411 cdef cart2_t b1, b2
1412 b1.x = dba[0,0]
1413 b1.y = dba[0,1]
1414 b2.x = dba[1,0]
1415 b2.y = dba[0,1]
1416 cdef double maxR = layers*max(cart2norm(b1), cart2norm(b2))
1417 return PGen_xyWeb_new(b1, b2, BASIS_RTOL, CART2_ZERO, 0, include_origin, maxR, False)
1419 cdef double get_unitcell_volume(direct_basis):
1420 dba = np.array(direct_basis)
1421 if not (dba.shape == (2,2)):
1422 raise NotImplementedError
1423 cdef cart2_t b1, b2
1424 b1.x = dba[0,0]
1425 b1.y = dba[0,1]
1426 b2.x = dba[1,0]
1427 b2.y = dba[0,1]
1428 return l2d_unitcell_area(b1, b2)
1430 cdef PGen get_PGen_reciprocal2pi(direct_basis, double layers = 30):
1431 dba = np.array(direct_basis)
1432 if not (dba.shape == (2,2)):
1433 raise NotImplementedError
1434 cdef cart2_t b1, b2, rb1, rb2
1435 b1.x = dba[0,0]
1436 b1.y = dba[0,1]
1437 b2.x = dba[1,0]
1438 b2.y = dba[0,1]
1439 if(l2d_reciprocalBasis2pi(b1, b2, &rb1, &rb2) != 0):
1440 raise RuntimeError
1441 cdef double maxK = layers*max(cart2norm(rb1), cart2norm(rb2))
1442 return PGen_xyWeb_new(rb1, rb2, BASIS_RTOL, CART2_ZERO,
1443 0, True, maxK, False)
1445 cdef class Ewald3Calculator:
1446 '''Wrapper class over qpms_ewald3_constants_t.
1448 Mainly for testing low-level scalar Ewald summation functionality.'''
1449 cdef qpms_ewald3_constants_t *c
1451 def __cinit__(self, qpms_l_t lMax, int csphase = -1):
1452 if (csphase != -1 and csphase != 1):
1453 raise ValueError("csphase must be +1 or -1, not %d" % csphase)
1454 self.c = qpms_ewald3_constants_init(lMax, csphase)
1456 def __dealloc__(self):
1457 qpms_ewald3_constants_free(self.c)
1459 def sigma0(self, double eta, cdouble wavenumber, do_err = False):
1460 cdef int retval
1461 cdef double err
1462 cdef cdouble result
1463 retval = ewald3_sigma0(&result, &err, self.c, eta, wavenumber)
1464 if retval:
1465 raise RuntimeError("ewald3_sigma0 returned non-zero value (%d)" % retval)
1466 if do_err:
1467 return (result, err)
1468 else:
1469 return result
1471 def sigma_short(self, double eta, cdouble wavenumber, direct_basis, wavevector, particle_shift, do_err=False):
1472 # FIXME now only 2d XY lattice in 3D is implemented here, we don't even do proper dimensionality checks.
1473 cdef cart3_t beta, pshift
1474 beta.x = wavevector[0]
1475 beta.y = wavevector[1]
1476 beta.z = 0
1477 pshift.x = particle_shift[0]
1478 pshift.y = particle_shift[1]
1479 pshift.z = 0
1480 cdef qpms_l_t n = self.c[0].nelem_sc
1481 cdef np.ndarray[complex, ndim=1] result = np.empty((n,), dtype=complex)
1482 cdef cdouble[::1] result_v = result
1483 cdef np.ndarray[double, ndim=1] err
1484 cdef double[::1] err_v
1485 if do_err:
1486 err = np.empty((n,), dtype=np.double)
1487 err_v = err
1488 cdef bint include_origin = not (particle_shift[0] == 0 and particle_shift[1] == 0)
1489 cdef PGen rgen = get_PGen_direct(direct_basis, include_origin)
1490 cdef int retval = ewald3_sigma_short(&result_v[0], &err_v[0] if do_err else NULL,
1491 self.c, eta, wavenumber, LAT_2D_IN_3D_XYONLY, &rgen, False, beta, pshift)
1492 if rgen.stateData: PGen_destroy(&rgen)
1493 if retval: raise RuntimeError("ewald3_sigma_short returned %d" % retval)
1494 if do_err:
1495 return (result, err)
1496 else:
1497 return result
1499 def sigma_long(self, double eta, cdouble wavenumber, direct_basis, wavevector, particle_shift, do_err=False):
1500 # FIXME now only 2d XY lattice in 3D is implemented here, we don't even do proper dimensionality checks.
1501 cdef cart3_t beta, pshift
1502 beta.x = wavevector[0]
1503 beta.y = wavevector[1]
1504 beta.z = 0
1505 pshift.x = particle_shift[0]
1506 pshift.y = particle_shift[1]
1507 pshift.z = 0
1508 cdef qpms_l_t n = self.c[0].nelem_sc
1509 cdef np.ndarray[complex, ndim=1] result = np.empty((n,), dtype=complex)
1510 cdef cdouble[::1] result_v = result
1511 cdef np.ndarray[double, ndim=1] err
1512 cdef double[::1] err_v
1513 if do_err:
1514 err = np.empty((n,), dtype=np.double)
1515 err_v = err
1516 cdef PGen kgen = get_PGen_reciprocal2pi(direct_basis)
1517 cdef double unitcell_volume = get_unitcell_volume(direct_basis)
1518 cdef int retval = ewald3_sigma_long(&result_v[0], &err_v[0] if do_err else NULL,
1519 self.c, eta, wavenumber, unitcell_volume, LAT_2D_IN_3D_XYONLY, &kgen, False, beta, pshift)
1521 if kgen.stateData: PGen_destroy(&kgen)
1522 if retval: raise RuntimeError("ewald3_sigma_long returned %d" % retval)
1523 if do_err:
1524 return (result, err)
1525 else:
1526 return result