def mie_coefficient(tau, l, k_medium, k_particle, radius): """Return the Mie coefficients of a sphere. Args: tau integer: spherical polarization, 0 for spherical TE and 1 for spherical TM l integer: l=1,... multipole degree (polar quantum number) k_medium float or complex: wavenumber in surrounding medium (inverse length unit) k_particle float or complex: wavenumber inside sphere (inverse length unit) radius float: radius of sphere (length unit) Returns: Mie coefficients as complex """ jlkr_medium = sf.spherical_bessel(l, k_medium * radius) jlkr_particle = sf.spherical_bessel(l, k_particle * radius) dxxj_medium = sf.dx_xj(l, k_medium * radius) dxxj_particle = sf.dx_xj(l, k_particle * radius) hlkr_medium = sf.spherical_hankel(l, k_medium * radius) dxxh_medium = sf.dx_xh(l, k_medium * radius) if tau == 0: q = (jlkr_medium * dxxj_particle - jlkr_particle * dxxj_medium) / ( jlkr_particle * dxxh_medium - hlkr_medium * dxxj_particle) elif tau == 1: q = ((k_medium**2 * jlkr_medium * dxxj_particle - k_particle**2 * jlkr_particle * dxxj_medium) / (k_particle**2 * jlkr_particle * dxxh_medium - k_medium**2 * hlkr_medium * dxxj_particle)) else: raise ValueError('tau must be 0 (spherical TE) or 1 (spherical TM)') return q
def sph_bessel(nu, l, Z): if nu == 1: return sf.spherical_bessel(l, Z) elif nu == 3: return sf.spherical_hankel(l, Z) else: return 0
def translation_coefficients_svwf(tau1, l1, m1, tau2, l2, m2, k, d, sph_hankel=None, legendre=None, exp_immphi=None): r"""Coefficients of the translation operator for the expansion of an outgoing spherical wave in terms of regular spherical waves with respect to a different origin: .. math:: \mathbf{\Psi}_{\tau l m}^{(3)}(\mathbf{r} + \mathbf{d} = \sum_{\tau'} \sum_{l'} \sum_{m'} A_{\tau l m, \tau' l' m'} (\mathbf{d}) \mathbf{\Psi}_{\tau' l' m'}^{(1)}(\mathbf{r}) for :math:`|\mathbf{r}|<|\mathbf{d}|`. See also section 2.3.3 and appendix B of [Egel 2018 diss]. Args: tau1 (int): tau1=0,1: Original wave's spherical polarization l1 (int): l=1,...: Original wave's SVWF multipole degree m1 (int): m=-l,...,l: Original wave's SVWF multipole order tau2 (int): tau2=0,1: Partial wave's spherical polarization l2 (int): l=1,...: Partial wave's SVWF multipole degree m2 (int): m=-l,...,l: Partial wave's SVWF multipole order k (float or complex): wavenumber (inverse length unit) d (list): translation vectors in format [dx, dy, dz] (length unit) dx, dy, dz can be scalars or ndarrays sph_hankel (list): Optional. sph_hankel[i] contains the spherical hankel funciton of degree i, evaluated at k*d where d is the norm of the distance vector(s) legendre (list): Optional. legendre[l][m] contains the legendre function of order l and degree m, evaluated at cos(theta) where theta is the polar angle(s) of the distance vector(s) Returns: translation coefficient A (complex) """ # spherical coordinates of d: dd = np.sqrt(d[0] ** 2 + d[1] ** 2 + d[2] ** 2) if exp_immphi is None: phid = np.arctan2(d[1], d[0]) eimph = np.exp(1j * (m1 - m2) * phid) else: eimph = exp_immphi[m1][m2] if sph_hankel is None: sph_hankel = [mathfunc.spherical_hankel(n, k * dd) for n in range(l1 + l2 + 1)] if legendre is None: costthetd = d[2] / dd sinthetd = np.sqrt(d[0] ** 2 + d[1] ** 2) / dd legendre, _, _ = mathfunc.legendre_normalized(costthetd, sinthetd, l1 + l2) A = complex(0) for ld in range(abs(l1 - l2), l1 + l2 + 1): a5, b5 = ab5_coefficients(l1, m1, l2, m2, ld) if tau1 == tau2: A += a5 * sph_hankel[ld] * legendre[ld][abs(m1 - m2)] else: A += b5 * sph_hankel[ld] * legendre[ld][abs(m1 - m2)] A = eimph * A return A
def spherical_vector_wave_function(x, y, z, k, nu, tau, l, m): """Electric field components of spherical vector wave function (SVWF). The conventions are chosen according to `A. Doicu, T. Wriedt, and Y. A. Eremin: "Light Scattering by Systems of Particles", Springer-Verlag, 2006 <https://doi.org/10.1007/978-3-540-33697-6>`_ See also section 2.3.2 of [Egel 2018 diss]. Args: x (numpy.ndarray): x-coordinate of position where to test the field (length unit) y (numpy.ndarray): y-coordinate of position where to test the field z (numpy.ndarray): z-coordinate of position where to test the field k (float or complex): wavenumber (inverse length unit) nu (int): 1 for regular waves, 3 for outgoing waves tau (int): spherical polarization, 0 for spherical TE and 1 for spherical TM l (int): l=1,... multipole degree (polar quantum number) m (int): m=-l,...,l multipole order (azimuthal quantum number) Returns: - x-coordinate of SVWF electric field (numpy.ndarray) - y-coordinate of SVWF electric field (numpy.ndarray) - z-coordinate of SVWF electric field (numpy.ndarray) """ x = np.array(x) y = np.array(y) z = np.array(z) r = np.sqrt(x**2 + y**2 + z**2) non0 = np.logical_not(r == 0) theta = np.zeros(x.shape) theta[non0] = np.arccos(z[non0] / r[non0]) phi = np.arctan2(y, x) # unit vector in r-direction er_x = np.zeros(x.shape) er_y = np.zeros(x.shape) er_z = np.ones(x.shape) er_x[non0] = x[non0] / r[non0] er_y[non0] = y[non0] / r[non0] er_z[non0] = z[non0] / r[non0] # unit vector in theta-direction eth_x = np.cos(theta) * np.cos(phi) eth_y = np.cos(theta) * np.sin(phi) eth_z = -np.sin(theta) # unit vector in phi-direction eph_x = -np.sin(phi) eph_y = np.cos(phi) eph_z = x - x cos_thet = np.cos(theta) sin_thet = np.sin(theta) plm_list, pilm_list, taulm_list = sf.legendre_normalized( cos_thet, sin_thet, l) plm = plm_list[l][abs(m)] pilm = pilm_list[l][abs(m)] taulm = taulm_list[l][abs(m)] kr = k * r if nu == 1: bes = sf.spherical_bessel(l, kr) dxxz = sf.dx_xj(l, kr) bes_kr = np.zeros(kr.shape, dtype=np.complex) dxxz_kr = np.zeros(kr.shape, dtype=np.complex) zero = (r == 0) nonzero = np.logical_not(zero) bes_kr[zero] = (l == 1) / 3 dxxz_kr[zero] = (l == 1) * 2 / 3 bes_kr[nonzero] = (bes[nonzero] / kr[nonzero]) dxxz_kr[nonzero] = (dxxz[nonzero] / kr[nonzero]) elif nu == 3: bes = sf.spherical_hankel(l, kr) dxxz = sf.dx_xh(l, kr) bes_kr = bes / kr dxxz_kr = dxxz / kr else: raise ValueError('nu must be 1 (regular SVWF) or 3 (outgoing SVWF)') eimphi = np.exp(1j * m * phi) prefac = 1 / np.sqrt(2 * l * (l + 1)) if tau == 0: Ex = prefac * bes * (1j * m * pilm * eth_x - taulm * eph_x) * eimphi Ey = prefac * bes * (1j * m * pilm * eth_y - taulm * eph_y) * eimphi Ez = prefac * bes * (1j * m * pilm * eth_z - taulm * eph_z) * eimphi elif tau == 1: Ex = prefac * (l * (l + 1) * bes_kr * plm * er_x + dxxz_kr * (taulm * eth_x + 1j * m * pilm * eph_x)) * eimphi Ey = prefac * (l * (l + 1) * bes_kr * plm * er_y + dxxz_kr * (taulm * eth_y + 1j * m * pilm * eph_y)) * eimphi Ez = prefac * (l * (l + 1) * bes_kr * plm * er_z + dxxz_kr * (taulm * eth_z + 1j * m * pilm * eph_z)) * eimphi else: raise ValueError('tau must be 0 (spherical TE) or 1 (spherical TM)') return Ex, Ey, Ez
def radial_coupling_lookup_table(vacuum_wavelength, particle_list, layer_system, k_parallel='default', resolution=None): """Prepare Sommerfeld integral lookup table to allow for a fast calculation of the coupling matrix by interpolation. This function is called when all particles are on the same z-position. Args: vacuum_wavelength (float): Vacuum wavelength in length units particle_list (list): List of particle objects layer_system (smuthi.layers.LayerSystem): Stratified medium k_parallel (numpy.ndarray or str): In-plane wavenumber for Sommerfeld integrals. If 'default', smuthi.fields.default_Sommerfeld_k_parallel_array resolution (float): Spatial resolution of lookup table in length units. (default: vacuum_wavelength / 100) Smaller means more accurate but higher memory footprint Returns: (tuple) tuple containing: lookup_table (ndarray): Coupling lookup, indices are [rho, n1, n2]. rho_array (ndarray): Values for the radial distance considered for the lookup (starting from negative numbers to allow for simpler cubic interpolation without distinction of cases at rho=0) """ sys.stdout.write('Prepare radial particle coupling lookup:\n') sys.stdout.flush() if resolution is None: resolution = vacuum_wavelength / 100 sys.stdout.write('Setting lookup resolution to %f\n' % resolution) sys.stdout.flush() l_max = max([particle.l_max for particle in particle_list]) m_max = max([particle.m_max for particle in particle_list]) blocksize = smuthi.fields.blocksize(l_max, m_max) x_array = np.array([particle.position[0] for particle in particle_list]) y_array = np.array([particle.position[1] for particle in particle_list]) rho_array = np.sqrt((x_array[:, None] - x_array[None, :])**2 + (y_array[:, None] - y_array[None, :])**2) radial_distance_array = np.arange(-3 * resolution, rho_array.max() + 3 * resolution, resolution) z = particle_list[0].position[2] i_s = layer_system.layer_number(z) k_is = layer_system.wavenumber(i_s, vacuum_wavelength) dz = z - layer_system.reference_z(i_s) len_rho = len(radial_distance_array) # direct ----------------------------------------------------------------------------------------------------------- w = np.zeros((len_rho, blocksize, blocksize), dtype=np.complex64) sys.stdout.write('Memory footprint: ' + size_format(w.nbytes) + '\n') sys.stdout.flush() ct = np.array([0.0]) st = np.array([1.0]) bessel_h = [] for n in range(2 * l_max + 1): bessel_h.append(sf.spherical_hankel(n, k_is * radial_distance_array)) bessel_h[-1][radial_distance_array <= 0] = np.nan legendre, _, _ = sf.legendre_normalized(ct, st, 2 * l_max) pbar = tqdm( total=blocksize**2, desc='Direct coupling ', file=sys.stdout, bar_format='{l_bar}{bar}| elapsed: {elapsed} remaining: {remaining}') for m1 in range(-m_max, m_max + 1): for m2 in range(-m_max, m_max + 1): for l1 in range(max(1, abs(m1)), l_max + 1): for l2 in range(max(1, abs(m2)), l_max + 1): A = np.zeros(len_rho, dtype=complex) B = np.zeros(len_rho, dtype=complex) for ld in range(max(abs(l1 - l2), abs(m1 - m2)), l1 + l2 + 1): # if ld<abs(m1-m2) then P=0 a5, b5 = trf.ab5_coefficients(l2, m2, l1, m1, ld) A = A + a5 * bessel_h[ld] * legendre[ld][abs(m1 - m2)] B = B + b5 * bessel_h[ld] * legendre[ld][abs(m1 - m2)] for tau1 in range(2): n1 = smuthi.fields.multi_to_single_index( tau1, l1, m1, l_max, m_max) for tau2 in range(2): n2 = smuthi.fields.multi_to_single_index( tau2, l2, m2, l_max, m_max) if tau1 == tau2: w[:, n1, n2] = A # remember that w = A.T else: w[:, n1, n2] = B pbar.update() pbar.close() close_to_zero = radial_distance_array < rho_array[ ~np.eye(rho_array.shape[0], dtype=bool)].min() / 2 w[close_to_zero, :, :] = 0 # switch off direct coupling contribution near rho=0 # layer mediated --------------------------------------------------------------------------------------------------- sys.stdout.write('Layer mediated coupling : ...') sys.stdout.flush() if type(k_parallel) == str and k_parallel == 'default': k_parallel = smuthi.fields.default_Sommerfeld_k_parallel_array kz_is = smuthi.fields.k_z(k_parallel=k_parallel, k=k_is) len_kp = len(k_parallel) # phase factors epl2jkz = np.exp(2j * kz_is * dz) emn2jkz = np.exp(-2j * kz_is * dz) # layer response L = np.zeros((2, 2, 2, len_kp), dtype=complex) # pol, pl/mn1, pl/mn2, kp for pol in range(2): L[pol, :, :, :] = lay.layersystem_response_matrix( pol, layer_system.thicknesses, layer_system.refractive_indices, k_parallel, smuthi.fields.angular_frequency(vacuum_wavelength), i_s, i_s) # transformation coefficients B_dag = np.zeros((2, 2, blocksize, len_kp), dtype=complex) # pol, pl/mn, n, kp B = np.zeros((2, 2, blocksize, len_kp), dtype=complex) # pol, pl/mn, n, kp ct = kz_is / k_is st = k_parallel / k_is _, pilm_pl, taulm_pl = sf.legendre_normalized(ct, st, l_max) _, pilm_mn, taulm_mn = sf.legendre_normalized(-ct, st, l_max) m_list = [None for n in range(blocksize)] for tau in range(2): for m in range(-m_max, m_max + 1): for l in range(max(1, abs(m)), l_max + 1): n = smuthi.fields.multi_to_single_index( tau, l, m, l_max, m_max) m_list[n] = m for pol in range(2): B_dag[pol, 0, n, :] = trf.transformation_coefficients_vwf( tau, l, m, pol, pilm_list=pilm_pl, taulm_list=taulm_pl, dagger=True) B_dag[pol, 1, n, :] = trf.transformation_coefficients_vwf( tau, l, m, pol, pilm_list=pilm_mn, taulm_list=taulm_mn, dagger=True) B[pol, 0, n, :] = trf.transformation_coefficients_vwf( tau, l, m, pol, pilm_list=pilm_pl, taulm_list=taulm_pl, dagger=False) B[pol, 1, n, :] = trf.transformation_coefficients_vwf( tau, l, m, pol, pilm_list=pilm_mn, taulm_list=taulm_mn, dagger=False) # pairs of (n1, n2), listed by abs(m1-m2) n1n2_combinations = [[] for dm in range(2 * m_max + 1)] for n1 in range(blocksize): m1 = m_list[n1] for n2 in range(blocksize): m2 = m_list[n2] n1n2_combinations[abs(m1 - m2)].append((n1, n2)) wr = np.zeros((len_rho, blocksize, blocksize), dtype=complex) dkp = np.diff(k_parallel) if cu.use_gpu: re_dkp_d = cu.gpuarray.to_gpu(np.float32(dkp.real)) im_dkp_d = cu.gpuarray.to_gpu(np.float32(dkp.imag)) kernel_source_code = cusrc.radial_lookup_assembly_code % ( blocksize, len_rho, len_kp) helper_function = cu.SourceModule(kernel_source_code).get_function( "helper") cuda_blocksize = 128 cuda_gridsize = (len_rho + cuda_blocksize - 1) // cuda_blocksize re_dwr_d = cu.gpuarray.to_gpu(np.zeros(len_rho, dtype=np.float32)) im_dwr_d = cu.gpuarray.to_gpu(np.zeros(len_rho, dtype=np.float32)) n1n2_combinations = [[] for dm in range(2 * m_max + 1)] for n1 in range(blocksize): m1 = m_list[n1] for n2 in range(blocksize): m2 = m_list[n2] n1n2_combinations[abs(m1 - m2)].append((n1, n2)) pbar = tqdm( total=blocksize**2, desc='Layer mediated coupling ', file=sys.stdout, bar_format='{l_bar}{bar}| elapsed: {elapsed} remaining: {remaining}') for dm in range(2 * m_max + 1): bessel = scipy.special.jv( dm, (k_parallel[None, :] * radial_distance_array[:, None])) besjac = bessel * (k_parallel / (kz_is * k_is))[None, :] for n1n2 in n1n2_combinations[dm]: n1 = n1n2[0] m1 = m_list[n1] n2 = n1n2[1] m2 = m_list[n2] belbe = np.zeros(len_kp, dtype=complex) # n1, n2, kp for pol in range(2): belbe += L[pol, 0, 0, :] * B_dag[pol, 0, n1, :] * B[pol, 0, n2, :] belbe += L[pol, 1, 0, :] * B_dag[pol, 1, n1, :] * B[pol, 0, n2, :] * emn2jkz belbe += L[pol, 0, 1, :] * B_dag[pol, 0, n1, :] * B[pol, 1, n2, :] * epl2jkz belbe += L[pol, 1, 1, :] * B_dag[pol, 1, n1, :] * B[pol, 1, n2, :] if cu.use_gpu: re_belbe_d = cu.gpuarray.to_gpu(np.float32( belbe[None, :].real)) im_belbe_d = cu.gpuarray.to_gpu(np.float32( belbe[None, :].imag)) re_besjac_d = cu.gpuarray.to_gpu(np.float32(besjac.real)) im_besjac_d = cu.gpuarray.to_gpu(np.float32(besjac.imag)) helper_function(re_besjac_d.gpudata, im_besjac_d.gpudata, re_belbe_d.gpudata, im_belbe_d.gpudata, re_dkp_d.gpudata, im_dkp_d.gpudata, re_dwr_d.gpudata, im_dwr_d.gpudata, block=(cuda_blocksize, 1, 1), grid=(cuda_gridsize, 1)) wr[:, n1, n2] = 4 * (1j)**abs(m2 - m1) * (re_dwr_d.get() + 1j * im_dwr_d.get()) else: integrand = besjac * belbe[None, :] # rho, kp wr[:, n1, n2] = 2 * (1j)**abs(m2 - m1) * ( (integrand[:, :-1] + integrand[:, 1:]) * dkp[None, :]).sum( axis=-1) # trapezoidal rule pbar.update() pbar.close() return w + wr, radial_distance_array
def volumetric_coupling_lookup_table(vacuum_wavelength, particle_list, layer_system, k_parallel='default', resolution=None): """Prepare Sommerfeld integral lookup table to allow for a fast calculation of the coupling matrix by interpolation. This function is called when not all particles are on the same z-position. Args: vacuum_wavelength (float): Vacuum wavelength in length units particle_list (list): List of particle objects layer_system (smuthi.layers.LayerSystem): Stratified medium k_parallel (numpy.ndarray or str): In-plane wavenumber for Sommerfeld integrals. If 'default', smuthi.fields.default_Sommerfeld_k_parallel_array resolution (float): Spatial resolution of lookup table in length units. (default: vacuum_wavelength / 100) Smaller means more accurate but higher memory footprint Returns: (tuple): tuple containing: w_pl (ndarray): Coupling lookup for z1 + z2, indices are [rho, z, n1, n2]. Includes layer mediated coupling. w_mn (ndarray): Coupling lookup for z1 + z2, indices are [rho, z, n1, n2]. Includes layer mediated and direct coupling. rho_array (ndarray): Values for the radial distance considered for the lookup (starting from negative numbers to allow for simpler cubic interpolation without distinction of cases for lookup edges sz_array (ndarray): Values for the sum of z-coordinates (z1 + z2) considered for the lookup dz_array (ndarray): Values for the difference of z-coordinates (z1 - z2) considered for the lookup """ sys.stdout.write('Prepare 3D particle coupling lookup:\n') sys.stdout.flush() if resolution is None: resolution = vacuum_wavelength / 100 sys.stdout.write('Setting lookup resolution to %f\n' % resolution) sys.stdout.flush() l_max = max([particle.l_max for particle in particle_list]) m_max = max([particle.m_max for particle in particle_list]) blocksize = smuthi.fields.blocksize(l_max, m_max) particle_x_array = np.array( [particle.position[0] for particle in particle_list]) particle_y_array = np.array( [particle.position[1] for particle in particle_list]) particle_z_array = np.array( [particle.position[2] for particle in particle_list]) particle_rho_array = np.sqrt( (particle_x_array[:, None] - particle_x_array[None, :])**2 + (particle_y_array[:, None] - particle_y_array[None, :])**2) dz_min = particle_z_array.min() - particle_z_array.max() dz_max = particle_z_array.max() - particle_z_array.min() sz_min = 2 * particle_z_array.min() sz_max = 2 * particle_z_array.max() rho_array = np.arange(-3 * resolution, particle_rho_array.max() + 3 * resolution, resolution) sz_array = np.arange(sz_min - 3 * resolution, sz_max + 3 * resolution, resolution) dz_array = np.arange(dz_min - 3 * resolution, dz_max + 3 * resolution, resolution) len_rho = len(rho_array) len_sz = len(sz_array) len_dz = len(dz_array) assert len_sz == len_dz i_s = layer_system.layer_number(particle_list[0].position[2]) k_is = layer_system.wavenumber(i_s, vacuum_wavelength) z_is = layer_system.reference_z(i_s) # direct ----------------------------------------------------------------------------------------------------------- w = np.zeros((len_rho, len_dz, blocksize, blocksize), dtype=np.complex64) sys.stdout.write('Lookup table memory footprint: ' + size_format(2 * w.nbytes) + '\n') sys.stdout.flush() r_array = np.sqrt(dz_array[None, :]**2 + rho_array[:, None]**2) r_array[r_array == 0] = 1e-20 ct = dz_array[None, :] / r_array st = rho_array[:, None] / r_array legendre, _, _ = sf.legendre_normalized(ct, st, 2 * l_max) bessel_h = [] for dm in tqdm( range(2 * l_max + 1), desc='Spherical Hankel lookup ', file=sys.stdout, bar_format='{l_bar}{bar}| elapsed: {elapsed} remaining: {remaining}' ): bessel_h.append(sf.spherical_hankel(dm, k_is * r_array)) pbar = tqdm( total=blocksize**2, desc='Direct coupling ', file=sys.stdout, bar_format='{l_bar}{bar}| elapsed: {elapsed} remaining: {remaining}') for m1 in range(-m_max, m_max + 1): for m2 in range(-m_max, m_max + 1): for l1 in range(max(1, abs(m1)), l_max + 1): for l2 in range(max(1, abs(m2)), l_max + 1): A = np.zeros((len_rho, len_dz), dtype=complex) B = np.zeros((len_rho, len_dz), dtype=complex) for ld in range(max(abs(l1 - l2), abs(m1 - m2)), l1 + l2 + 1): # if ld<abs(m1-m2) then P=0 a5, b5 = trf.ab5_coefficients( l2, m2, l1, m1, ld) # remember that w = A.T A += a5 * bessel_h[ld] * legendre[ld][abs( m1 - m2)] # remember that w = A.T B += b5 * bessel_h[ld] * legendre[ld][abs( m1 - m2)] # remember that w = A.T for tau1 in range(2): n1 = smuthi.fields.multi_to_single_index( tau1, l1, m1, l_max, m_max) for tau2 in range(2): n2 = smuthi.fields.multi_to_single_index( tau2, l2, m2, l_max, m_max) if tau1 == tau2: w[:, :, n1, n2] = A else: w[:, :, n1, n2] = B pbar.update() pbar.close() # switch off direct coupling contribution near rho=0: w[rho_array < particle_rho_array[~np.eye(particle_rho_array.shape[0], dtype=bool)].min( ) / 2, :, :, :] = 0 # layer mediated --------------------------------------------------------------------------------------------------- sys.stdout.write('Layer mediated coupling : ...') sys.stdout.flush() if type(k_parallel) == str and k_parallel == 'default': k_parallel = smuthi.fields.default_Sommerfeld_k_parallel_array kz_is = smuthi.fields.k_z(k_parallel=k_parallel, k=k_is) len_kp = len(k_parallel) # phase factors epljksz = np.exp(1j * kz_is[None, :] * (sz_array[:, None] - 2 * z_is)) # z, k emnjksz = np.exp(-1j * kz_is[None, :] * (sz_array[:, None] - 2 * z_is)) epljkdz = np.exp(1j * kz_is[None, :] * dz_array[:, None]) emnjkdz = np.exp(-1j * kz_is[None, :] * dz_array[:, None]) # layer response L = np.zeros((2, 2, 2, len_kp), dtype=complex) # pol, pl/mn1, pl/mn2, kp for pol in range(2): L[pol, :, :, :] = lay.layersystem_response_matrix( pol, layer_system.thicknesses, layer_system.refractive_indices, k_parallel, smuthi.fields.angular_frequency(vacuum_wavelength), i_s, i_s) # transformation coefficients B_dag = np.zeros((2, 2, blocksize, len_kp), dtype=complex) # pol, pl/mn, n, kp B = np.zeros((2, 2, blocksize, len_kp), dtype=complex) # pol, pl/mn, n, kp ct_k = kz_is / k_is st_k = k_parallel / k_is _, pilm_pl, taulm_pl = sf.legendre_normalized(ct_k, st_k, l_max) _, pilm_mn, taulm_mn = sf.legendre_normalized(-ct_k, st_k, l_max) m_list = [None for i in range(blocksize)] for tau in range(2): for m in range(-m_max, m_max + 1): for l in range(max(1, abs(m)), l_max + 1): n = smuthi.fields.multi_to_single_index( tau, l, m, l_max, m_max) m_list[n] = m for pol in range(2): B_dag[pol, 0, n, :] = trf.transformation_coefficients_vwf( tau, l, m, pol, pilm_list=pilm_pl, taulm_list=taulm_pl, dagger=True) B_dag[pol, 1, n, :] = trf.transformation_coefficients_vwf( tau, l, m, pol, pilm_list=pilm_mn, taulm_list=taulm_mn, dagger=True) B[pol, 0, n, :] = trf.transformation_coefficients_vwf( tau, l, m, pol, pilm_list=pilm_pl, taulm_list=taulm_pl, dagger=False) B[pol, 1, n, :] = trf.transformation_coefficients_vwf( tau, l, m, pol, pilm_list=pilm_mn, taulm_list=taulm_mn, dagger=False) # pairs of (n1, n2), listed by abs(m1-m2) n1n2_combinations = [[] for dm in range(2 * m_max + 1)] for n1 in range(blocksize): m1 = m_list[n1] for n2 in range(blocksize): m2 = m_list[n2] n1n2_combinations[abs(m1 - m2)].append((n1, n2)) wr_pl = np.zeros((len_rho, len_dz, blocksize, blocksize), dtype=np.complex64) wr_mn = np.zeros((len_rho, len_dz, blocksize, blocksize), dtype=np.complex64) dkp = np.diff(k_parallel) if cu.use_gpu: re_dkp_d = cu.gpuarray.to_gpu(np.float32(dkp.real)) im_dkp_d = cu.gpuarray.to_gpu(np.float32(dkp.imag)) kernel_source_code = cusrc.volume_lookup_assembly_code % ( blocksize, len_rho, len_sz, len_kp) helper_function = cu.SourceModule(kernel_source_code).get_function( "helper") cuda_blocksize = 128 cuda_gridsize = (len_rho * len_sz + cuda_blocksize - 1) // cuda_blocksize re_dwr_d = cu.gpuarray.to_gpu( np.zeros((len_rho, len_sz), dtype=np.float32)) im_dwr_d = cu.gpuarray.to_gpu( np.zeros((len_rho, len_sz), dtype=np.float32)) pbar = tqdm( total=blocksize**2, desc='Layer mediated coupling ', file=sys.stdout, bar_format='{l_bar}{bar}| elapsed: {elapsed} remaining: {remaining}') for dm in range(2 * m_max + 1): bessel = scipy.special.jv(dm, (k_parallel[None, :] * rho_array[:, None])) besjac = bessel * (k_parallel / (kz_is * k_is))[None, :] for n1n2 in n1n2_combinations[dm]: n1 = n1n2[0] m1 = m_list[n1] n2 = n1n2[1] m2 = m_list[n2] belbee_pl = np.zeros((len_dz, len_kp), dtype=complex) belbee_mn = np.zeros((len_dz, len_kp), dtype=complex) for pol in range(2): belbee_pl += ((L[pol, 0, 1, :] * B_dag[pol, 0, n1, :] * B[pol, 1, n2, :])[None, :] * epljksz + (L[pol, 1, 0, :] * B_dag[pol, 1, n1, :] * B[pol, 0, n2, :])[None, :] * emnjksz) belbee_mn += ((L[pol, 0, 0, :] * B_dag[pol, 0, n1, :] * B[pol, 0, n2, :])[None, :] * epljkdz + (L[pol, 1, 1, :] * B_dag[pol, 1, n1, :] * B[pol, 1, n2, :])[None, :] * emnjkdz) if cu.use_gpu: re_belbee_pl_d = cu.gpuarray.to_gpu( np.float32(belbee_pl[None, :, :].real)) im_belbee_pl_d = cu.gpuarray.to_gpu( np.float32(belbee_pl[None, :, :].imag)) re_belbee_mn_d = cu.gpuarray.to_gpu( np.float32(belbee_mn[None, :, :].real)) im_belbee_mn_d = cu.gpuarray.to_gpu( np.float32(belbee_mn[None, :, :].imag)) re_besjac_d = cu.gpuarray.to_gpu( np.float32(besjac[:, None, :].real)) im_besjac_d = cu.gpuarray.to_gpu( np.float32(besjac[:, None, :].imag)) helper_function(re_besjac_d.gpudata, im_besjac_d.gpudata, re_belbee_pl_d.gpudata, im_belbee_pl_d.gpudata, re_dkp_d.gpudata, im_dkp_d.gpudata, re_dwr_d.gpudata, im_dwr_d.gpudata, block=(cuda_blocksize, 1, 1), grid=(cuda_gridsize, 1)) wr_pl[:, :, n1, n2] = 4 * (1j)**abs(m2 - m1) * (re_dwr_d.get() + 1j * im_dwr_d.get()) helper_function(re_besjac_d.gpudata, im_besjac_d.gpudata, re_belbee_mn_d.gpudata, im_belbee_mn_d.gpudata, re_dkp_d.gpudata, im_dkp_d.gpudata, re_dwr_d.gpudata, im_dwr_d.gpudata, block=(cuda_blocksize, 1, 1), grid=(cuda_gridsize, 1)) wr_mn[:, :, n1, n2] = 4 * (1j)**abs(m2 - m1) * (re_dwr_d.get() + 1j * im_dwr_d.get()) else: integrand = besjac[:, None, :] * belbee_pl[None, :, :] wr_pl[:, :, n1, n2] = 2 * (1j)**abs(m2 - m1) * ( (integrand[:, :, :-1] + integrand[:, :, 1:]) * dkp[None, None, :]).sum(axis=-1) # trapezoidal rule integrand = besjac[:, None, :] * belbee_mn[None, :, :] wr_mn[:, :, n1, n2] = 2 * (1j)**abs(m2 - m1) * ( (integrand[:, :, :-1] + integrand[:, :, 1:]) * dkp[None, None, :]).sum(axis=-1) pbar.update() pbar.close() return wr_pl, w + wr_mn, rho_array, sz_array, dz_array
def direct_coupling_block(vacuum_wavelength, receiving_particle, emitting_particle, layer_system): """Direct particle coupling matrix :math:`W` for two particles. This routine is explicit, but slow. Args: vacuum_wavelength (float): Vacuum wavelength :math:`\lambda` (length unit) receiving_particle (smuthi.particles.Particle): Particle that receives the scattered field emitting_particle (smuthi.particles.Particle): Particle that emits the scattered field layer_system (smuthi.layers.LayerSystem): Stratified medium in which the coupling takes place Returns: Direct coupling matrix block as numpy array. """ omega = flds.angular_frequency(vacuum_wavelength) # index specs lmax1 = receiving_particle.l_max mmax1 = receiving_particle.m_max lmax2 = emitting_particle.l_max mmax2 = emitting_particle.m_max blocksize1 = flds.blocksize(lmax1, mmax1) blocksize2 = flds.blocksize(lmax2, mmax2) # initialize result w = np.zeros((blocksize1, blocksize2), dtype=complex) # check if particles are in same layer rS1 = receiving_particle.position rS2 = emitting_particle.position iS1 = layer_system.layer_number(rS1[2]) iS2 = layer_system.layer_number(rS2[2]) if iS1 == iS2 and not emitting_particle == receiving_particle: k = omega * layer_system.refractive_indices[iS1] dx = rS1[0] - rS2[0] dy = rS1[1] - rS2[1] dz = rS1[2] - rS2[2] d = np.sqrt(dx**2 + dy**2 + dz**2) cos_theta = dz / d sin_theta = np.sqrt(dx**2 + dy**2) / d phi = np.arctan2(dy, dx) # spherical functions bessel_h = [ sf.spherical_hankel(n, k * d) for n in range(lmax1 + lmax2 + 1) ] legendre, _, _ = sf.legendre_normalized(cos_theta, sin_theta, lmax1 + lmax2) # the particle coupling operator is the transpose of the SVWF translation operator # therefore, (l1,m1) and (l2,m2) are interchanged: for m1 in range(-mmax1, mmax1 + 1): for m2 in range(-mmax2, mmax2 + 1): eimph = np.exp(1j * (m2 - m1) * phi) for l1 in range(max(1, abs(m1)), lmax1 + 1): for l2 in range(max(1, abs(m2)), lmax2 + 1): A, B = complex(0), complex(0) for ld in range(max(abs(l1 - l2), abs(m1 - m2)), l1 + l2 + 1): # if ld<abs(m1-m2) then P=0 a5, b5 = trf.ab5_coefficients(l2, m2, l1, m1, ld) A += a5 * bessel_h[ld] * legendre[ld][abs(m1 - m2)] B += b5 * bessel_h[ld] * legendre[ld][abs(m1 - m2)] A, B = eimph * A, eimph * B for tau1 in range(2): n1 = flds.multi_to_single_index( tau1, l1, m1, lmax1, mmax1) for tau2 in range(2): n2 = flds.multi_to_single_index( tau2, l2, m2, lmax2, mmax2) if tau1 == tau2: w[n1, n2] = A else: w[n1, n2] = B return w