def angular_frequency(self): """Angular frequency. Returns: Angular frequency (float) according to the vacuum wavelength in units of c=1. """ return coord.angular_frequency(self.vacuum_wavelength)
def plane_wave_expansion(self, layer_system, i): """Plane wave expansion for the plane wave including its layer system response. As it already is a plane wave, the plane wave expansion is somehow trivial (containing only one partial wave, i.e., a discrete plane wave expansion). Args: layer_system (smuthi.layers.LayerSystem): Layer system object i (int): layer number in which the plane wave expansion is valid Returns: Tuple of smuthi.field_expansion.PlaneWaveExpansion objects. The first element is an upgoing PWE, whereas the second element is a downgoing PWE. """ if np.cos(self.polar_angle) > 0: iP = 0 kind = 'upgoing' else: iP = layer_system.number_of_layers() - 1 kind = 'downgoing' niP = layer_system.refractive_indices[iP] neff = np.sin([self.polar_angle]) * niP alpha = np.array([self.azimuthal_angle]) angular_frequency = coord.angular_frequency(self.vacuum_wavelength) k_iP = niP * angular_frequency k_Px = k_iP * np.sin(self.polar_angle) * np.cos(self.azimuthal_angle) k_Py = k_iP * np.sin(self.polar_angle) * np.sin(self.azimuthal_angle) k_Pz = k_iP * np.cos(self.polar_angle) z_iP = layer_system.reference_z(iP) amplitude_iP = self.amplitude * np.exp( -1j * (k_Px * self.reference_point[0] + k_Py * self.reference_point[1] + k_Pz * (self.reference_point[2] - z_iP))) loz = layer_system.lower_zlimit(iP) upz = layer_system.upper_zlimit(iP) pwe_exc = fldex.PlaneWaveExpansion(k=k_iP, k_parallel=neff * angular_frequency, azimuthal_angles=alpha, kind=kind, reference_point=[0, 0, z_iP], lower_z=loz, upper_z=upz) pwe_exc.coefficients[self.polarization, 0, 0] = amplitude_iP pwe_up, pwe_down = layer_system.response(pwe_exc, from_layer=iP, to_layer=i) if iP == i: if kind == 'upgoing': pwe_up = pwe_up + pwe_exc elif kind == 'downgoing': pwe_down = pwe_down + pwe_exc return pwe_up, pwe_down
def wavenumber(self, layer_number, vacuum_wavelength): """ Args: layer_number (int): number of layer in question vacuum_wavelength (float): vacuum wavelength Returns: wavenumber in that layer as float """ return self.refractive_indices[layer_number] * coord.angular_frequency( vacuum_wavelength=vacuum_wavelength)
def scattering_cross_section(initial_field, particle_list, layer_system, polar_angles='default', azimuthal_angles='default'): """Evaluate and display the differential scattering cross section as a function of solid angle. Args: initial_field (smuthi.initial.PlaneWave): Initial Plane wave particle_list (list): scattering particles layer_system (smuthi.layers.LayerSystem): stratified medium polar_angles (numpy.ndarray or str): polar angles values (radian). if 'default', use smuthi.coordinates.default_polar_angles azimuthal_angles (numpy.ndarray or str): azimuthal angle values (radian) if 'default', use smuthi.coordinates.default_azimuthal_angles Returns: A tuple of smuthi.field_expansion.FarField objects, one for forward scattering (i.e., into the top hemisphere) and one for backward scattering (bottom hemisphere). """ if not type(initial_field).__name__ == 'PlaneWave': raise ValueError( 'Cross section only defined for plane wave excitation.') i_top = layer_system.number_of_layers() - 1 vacuum_wavelength = initial_field.vacuum_wavelength omega = coord.angular_frequency(vacuum_wavelength) k_bot = omega * layer_system.refractive_indices[0] k_top = omega * layer_system.refractive_indices[-1] # read plane wave parameters A_P = initial_field.amplitude beta_P = initial_field.polar_angle if beta_P < np.pi / 2: i_P = 0 n_P = layer_system.refractive_indices[i_P] else: i_P = i_top n_P = layer_system.refractive_indices[i_P] if n_P.imag: raise ValueError( 'plane wave from absorbing layer: cross section undefined') else: n_P = n_P.real initial_intensity = abs(A_P)**2 * abs(np.cos(beta_P)) * n_P / 2 dscs = scattered_far_field(vacuum_wavelength, particle_list, layer_system, polar_angles, azimuthal_angles) dscs.signal_type = 'differential scattering cross section' dscs.signal = dscs.signal / initial_intensity return dscs
rho13 = np.sqrt( sum([(sphere1.position[i] - sphere3.position[i])**2 for i in range(3)])) wr13 = coup.layer_mediated_coupling_block(vacuum_wavelength, sphere1, sphere3, lay_sys) w13 = coup.direct_coupling_block(vacuum_wavelength, sphere1, sphere3, lay_sys) # initialize initial field object init_fld = init.GaussianBeam(vacuum_wavelength=vacuum_wavelength, polar_angle=beam_polar_angle, azimuthal_angle=beam_azimuthal_angle, polarization=beam_polarization, amplitude=beam_amplitude, reference_point=beam_focal_point, beam_waist=beam_waist, k_parallel_array=beam_neff_array * coord.angular_frequency(vacuum_wavelength)) # initialize simulation object simulation_direct = simul.Simulation(layer_system=lay_sys, particle_list=particle_list, initial_field=init_fld, solver_type='LU', store_coupling_matrix=True, log_to_terminal=False) simulation_direct.run() coefficients_direct = particle_list[0].scattered_field.coefficients cu.enable_gpu() simulation_lookup_linear_gpu = simul.Simulation( layer_system=lay_sys, particle_list=particle_list,
def scattered_field_pwe(vacuum_wavelength, particle_list, layer_system, layer_number, k_parallel='default', azimuthal_angles='default', include_direct=True, include_layer_response=True): """Calculate the plane wave expansion of the scattered field of a set of particles. Args: vacuum_wavelength (float): Vacuum wavelength (length unit) particle_list (list): List of Particle objects layer_system (smuthi.layers.LayerSystem): Stratified medium layer_number (int): Layer number in which the plane wave expansion should be valid k_parallel (numpy.ndarray or str): in-plane wavenumbers array. if 'default', use smuthi.coordinates.default_k_parallel azimuthal_angles (numpy.ndarray or str): azimuthal angles array if 'default', use smuthi.coordinates.default_azimuthal_angles include_direct (bool): If True, include the direct scattered field include_layer_response (bool): If True, include the layer system response Returns: A tuple of PlaneWaveExpansion objects for upgoing and downgoing waves. """ sys.stdout.write( 'Evaluating scattered field plane wave expansion in layer number %i ...\n' % layer_number) sys.stdout.flush() omega = coord.angular_frequency(vacuum_wavelength) k = omega * layer_system.refractive_indices[layer_number] z = layer_system.reference_z(layer_number) vb = (layer_system.lower_zlimit(layer_number), layer_system.upper_zlimit(layer_number)) pwe_up = fldex.PlaneWaveExpansion(k=k, k_parallel=k_parallel, azimuthal_angles=azimuthal_angles, kind='upgoing', reference_point=[0, 0, z], lower_z=vb[0], upper_z=vb[1]) pwe_down = fldex.PlaneWaveExpansion(k=k, k_parallel=k_parallel, azimuthal_angles=azimuthal_angles, kind='downgoing', reference_point=[0, 0, z], lower_z=vb[0], upper_z=vb[1]) for iS, particle in enumerate( tqdm(particle_list, desc='Scatt. field pwe ', file=sys.stdout, bar_format='{l_bar}{bar}| elapsed: {elapsed} ' 'remaining: {remaining}')): i_iS = layer_system.layer_number(particle.position[2]) # direct contribution if i_iS == layer_number and include_direct: pu, pd = fldex.swe_to_pwe_conversion( swe=particle.scattered_field, k_parallel=k_parallel, azimuthal_angles=azimuthal_angles, layer_system=layer_system) pwe_up = pwe_up + pu pwe_down = pwe_down + pd # layer mediated contribution if include_layer_response: pu, pd = fldex.swe_to_pwe_conversion( swe=particle.scattered_field, k_parallel=k_parallel, azimuthal_angles=azimuthal_angles, layer_system=layer_system, layer_number=layer_number, layer_system_mediated=True) pwe_up = pwe_up + pu pwe_down = pwe_down + pd return pwe_up, pwe_down
def scattered_field_piecewise_expansion(vacuum_wavelength, particle_list, layer_system, k_parallel='default', azimuthal_angles='default', layer_numbers=None): """Compute a piecewise field expansion of the scattered field. Args: vacuum_wavelength (float): vacuum wavelength particle_list (list): list of smuthi.particles.Particle objects layer_system (smuthi.layers.LayerSystem): stratified medium k_parallel (numpy.ndarray or str): in-plane wavenumbers array. if 'default', use smuthi.coordinates.default_k_parallel azimuthal_angles (numpy.ndarray or str): azimuthal angles array if 'default', use smuthi.coordinates.default_azimuthal_angles layer_numbers (list): if specified, append only plane wave expansions for these layers Returns: scattered field as smuthi.field_expansion.PiecewiseFieldExpansion object """ if layer_numbers is None: layer_numbers = range(layer_system.number_of_layers()) sfld = fldex.PiecewiseFieldExpansion() for i in tqdm(layer_numbers, desc='Scatt. field expansion ', file=sys.stdout, bar_format='{l_bar}{bar}| elapsed: {elapsed} ' 'remaining: {remaining}'): # layer mediated scattered field --------------------------------------------------------------------------- k = coord.angular_frequency( vacuum_wavelength) * layer_system.refractive_indices[i] ref = [0, 0, layer_system.reference_z(i)] vb = (layer_system.lower_zlimit(i), layer_system.upper_zlimit(i)) pwe_up = fldex.PlaneWaveExpansion(k=k, k_parallel=k_parallel, azimuthal_angles=azimuthal_angles, kind='upgoing', reference_point=ref, lower_z=vb[0], upper_z=vb[1]) pwe_down = fldex.PlaneWaveExpansion(k=k, k_parallel=k_parallel, azimuthal_angles=azimuthal_angles, kind='downgoing', reference_point=ref, lower_z=vb[0], upper_z=vb[1]) for particle in particle_list: add_up, add_down = fldex.swe_to_pwe_conversion( particle.scattered_field, k_parallel, azimuthal_angles, layer_system, i, True) pwe_up = pwe_up + add_up pwe_down = pwe_down + add_down # in bottom_layer, suppress upgoing waves, and in top layer, suppress downgoing waves if i > 0: sfld.expansion_list.append(pwe_up) if i < layer_system.number_of_layers() - 1: sfld.expansion_list.append(pwe_down) # direct field --------------------------------------------------------------------------------------------- for particle in particle_list: sfld.expansion_list.append(particle.scattered_field) return sfld
def extinction_cross_section(initial_field, particle_list, layer_system): """Evaluate and display the differential scattering cross section as a function of solid angle. Args: initial_field (smuthi.initial_field.PlaneWave): Plane wave object particle_list (list): List of smuthi.particles.Particle objects layer_system (smuthi.layers.LayerSystem): Representing the stratified medium Returns: Dictionary with following entries - 'forward': Extinction in the positinve z-direction (top layer) - 'backward': Extinction in the negative z-direction (bottom layer) """ if not type(initial_field).__name__ == 'PlaneWave': raise ValueError( 'Cross section only defined for plane wave excitation.') i_top = layer_system.number_of_layers() - 1 vacuum_wavelength = initial_field.vacuum_wavelength omega = coord.angular_frequency(vacuum_wavelength) k_bot = omega * layer_system.refractive_indices[0] k_top = omega * layer_system.refractive_indices[-1] # read plane wave parameters pol_P = initial_field.polarization beta_P = initial_field.polar_angle alpha_P = initial_field.azimuthal_angle if beta_P < np.pi / 2: i_P = 0 n_P = layer_system.refractive_indices[i_P] k_P = k_bot else: i_P = i_top n_P = layer_system.refractive_indices[i_P] k_P = k_top if n_P.imag: raise ValueError( 'plane wave from absorbing layer: cross section undefined') else: n_P = n_P.real # complex amplitude of initial wave (including phase factor for reference point) kappa_P = np.sin(beta_P) * k_P kx = np.cos(alpha_P) * kappa_P ky = np.sin(alpha_P) * kappa_P pm_kz_P = k_P * np.cos(beta_P) kvec_P = np.array([kx, ky, pm_kz_P]) rvec_iP = np.array([0, 0, layer_system.reference_z(i_P)]) rvec_0 = np.array(initial_field.reference_point) ejkriP = np.exp(1j * np.dot(kvec_P, rvec_iP - rvec_0)) A_P = initial_field.amplitude * ejkriP initial_intensity = abs(A_P)**2 * abs(np.cos(beta_P)) * n_P / 2 pwe_scat_top, _ = scattered_field_pwe(vacuum_wavelength, particle_list, layer_system, i_top, kappa_P, alpha_P) _, pwe_scat_bottom = scattered_field_pwe(vacuum_wavelength, particle_list, layer_system, 0, kappa_P, alpha_P) # bottom extinction _, pwe_init_bottom = initial_field.plane_wave_expansion(layer_system, 0) kz_bot = coord.k_z(k_parallel=kappa_P, k=k_bot) gRPbot = np.squeeze(pwe_init_bottom.coefficients[pol_P]) g_scat_bottom = np.squeeze(pwe_scat_bottom.coefficients[pol_P]) P_bot_ext = 4 * np.pi**2 * kz_bot / omega * (gRPbot * np.conj(g_scat_bottom)).real bottom_extinction_cs = -P_bot_ext / initial_intensity # top extinction pwe_init_top, _ = initial_field.plane_wave_expansion(layer_system, i_top) gRPtop = np.squeeze(pwe_init_top.coefficients[pol_P]) kz_top = coord.k_z(k_parallel=kappa_P, k=k_top) g_scat_top = np.squeeze(pwe_scat_top.coefficients[pol_P]) P_top_ext = 4 * np.pi**2 * kz_top / omega * (gRPtop * np.conj(g_scat_top)).real top_extinction_cs = -P_top_ext / initial_intensity extinction_cs = {'top': top_extinction_cs, 'bottom': bottom_extinction_cs} return extinction_cs
def scattered_far_field(vacuum_wavelength, particle_list, layer_system, polar_angles='default', azimuthal_angles='default'): """ Evaluate the scattered far field. Args: vacuum_wavelength (float): in length units particle_list (list): list of smuthi.Particle objects layer_system (smuthi.layers.LayerSystem): represents the stratified medium polar_angles (numpy.ndarray or str): polar angles values (radian). if 'default', use smuthi.coordinates.default_polar_angles azimuthal_angles (numpy.ndarray or str): azimuthal angle values (radian) if 'default', use smuthi.coordinates.default_azimuthal_angles Returns: A smuthi.field_expansion.FarField object of the scattered field. """ omega = coord.angular_frequency(vacuum_wavelength) if type(polar_angles) == str and polar_angles == 'default': polar_angles = coord.default_polar_angles if type(azimuthal_angles) == str and azimuthal_angles == 'default': azimuthal_angles = coord.default_azimuthal_angles if any(polar_angles.imag): raise ValueError("complex angles not allowed in far field") i_top = layer_system.number_of_layers() - 1 top_polar_angles = polar_angles[polar_angles < (np.pi / 2)] bottom_polar_angles = polar_angles[polar_angles > (np.pi / 2)] neff_top = np.sort( np.sin(top_polar_angles) * layer_system.refractive_indices[i_top]) neff_bottom = np.sort( np.sin(bottom_polar_angles) * layer_system.refractive_indices[0]) if len(top_polar_angles ) > 1 and layer_system.refractive_indices[i_top].imag == 0: pwe_top, _ = scattered_field_pwe(vacuum_wavelength, particle_list, layer_system, i_top, k_parallel=neff_top * omega, azimuthal_angles=azimuthal_angles, include_direct=True, include_layer_response=True) top_far_field = fldex.pwe_to_ff_conversion( vacuum_wavelength=vacuum_wavelength, plane_wave_expansion=pwe_top) else: top_far_field = None if len(bottom_polar_angles ) > 1 and layer_system.refractive_indices[0].imag == 0: _, pwe_bottom = scattered_field_pwe(vacuum_wavelength, particle_list, layer_system, 0, k_parallel=neff_bottom * omega, azimuthal_angles=azimuthal_angles, include_direct=True, include_layer_response=True) bottom_far_field = fldex.pwe_to_ff_conversion( vacuum_wavelength=vacuum_wavelength, plane_wave_expansion=pwe_bottom) else: bottom_far_field = None if top_far_field is not None: far_field = top_far_field if bottom_far_field is not None: far_field.append(bottom_far_field) else: far_field = bottom_far_field far_field.polar_angles = far_field.polar_angles.real return far_field
import smuthi.layers import smuthi.particles import smuthi.coordinates as coord ld = 532 A = 1 beta = 0 alpha = 0.2 * np.pi pol = 1 rS = [100, -200, 300] laysys = smuthi.layers.LayerSystem(thicknesses=[0, 0], refractive_indices=[1, 1]) particle = smuthi.particles.Sphere(position=rS, l_max=3, m_max=3) particle_list = [particle] al_ar = np.linspace(0, 2 * np.pi, 1000) kp_ar = np.linspace(0, 0.99999, 1000) * coord.angular_frequency(ld) bw = 4000 ref = [-100, 100, 200] gauss_beam = init.GaussianBeam(vacuum_wavelength=ld, polar_angle=beta, azimuthal_angle=alpha, polarization=pol, amplitude=A, reference_point=ref, k_parallel_array=kp_ar, beam_waist=bw) particle.initial_field = gauss_beam.spherical_wave_expansion(particle, laysys) def test_SWE_coefficients_against_prototype(): aI = particle.initial_field.coefficients
coord.set_default_k_parallel(vacuum_wavelength, neff_waypoints, neff_discr) # initialize particle object sphere1 = part.Sphere(position=[100, 100, 150], refractive_index=2.4 + 0.0j, radius=110, l_max=4, m_max=4) sphere2 = part.Sphere(position=[-100, -100, 250], refractive_index=1.9 + 0.0j, radius=120, l_max=3, m_max=3) sphere3 = part.Sphere(position=[-200, 100, 300], refractive_index=1.7 + 0.0j, radius=90, l_max=3, m_max=3) particle_list = [sphere1, sphere2, sphere3] # initialize layer system object lay_sys = lay.LayerSystem([0, 400, 0], [2, 1.4, 2]) # initialize initial field object init_fld = init.GaussianBeam(vacuum_wavelength=vacuum_wavelength, polar_angle=beam_polar_angle, azimuthal_angle=beam_azimuthal_angle, polarization=beam_polarization, amplitude=beam_amplitude, reference_point=beam_focal_point, beam_waist=beam_waist, k_parallel_array=beam_neff_array*coord.angular_frequency(vacuum_wavelength)) # initialize simulation object simulation_lu = simul.Simulation(layer_system=lay_sys, particle_list=particle_list, initial_field=init_fld, solver_type='LU') simulation_lu.run() coefficients_lu = particle_list[0].scattered_field.coefficients simulation_gmres = simul.Simulation(layer_system=lay_sys, particle_list=particle_list, initial_field=init_fld, solver_type='gmres', solver_tolerance=1e-5) simulation_gmres.run() coefficients_gmres = particle_list[0].scattered_field.coefficients def test_result(): relerr = np.linalg.norm(coefficients_lu - coefficients_gmres) / np.linalg.norm(coefficients_lu) print('relative error: ', relerr)
lmax = 10 rx = 100 ry = -100 rz = 100 dx = 400 dy = 800 dz = -500 tau = 0 l = 2 m = -1 # -------------------------------------------- k = coord.angular_frequency( vacuum_wavelength) * surrounding_medium_refractive_index # outgoing wave Ex, Ey, Ez = vwf.spherical_vector_wave_function(rx + dx, ry + dy, rz + dz, k, 3, tau, l, m) # series of incoming waves Ex2, Ey2, Ez2 = complex(0), complex(0), complex(0) for tau2 in range(2): for l2 in range(1, lmax + 1): for m2 in range(-l2, l2 + 1): Mx, My, Mz = vwf.spherical_vector_wave_function( rx, ry, rz, k, 1, tau2, l2, m2) A = vwf.translation_coefficients_svwf(tau, l, m, tau2, l2, m2, k, [dx, dy, dz]) Ex2 += Mx * A
def solve(self): """Compute scattered field coefficients and store them in the particles' spherical wave expansion objects.""" if len(self.particle_list) > 0: if self.solver_type == 'LU': sys.stdout.write('Solve (LU decomposition) : ...') if not hasattr(self.master_matrix.linear_operator, 'A'): raise ValueError( 'LU factorization only possible ' 'with the option "store coupling matrix".') if not hasattr(self.master_matrix, 'LU_piv'): lu, piv = scipy.linalg.lu_factor( self.master_matrix.linear_operator.A, overwrite_a=False) self.master_matrix.LU_piv = (lu, piv) b = scipy.linalg.lu_solve(self.master_matrix.LU_piv, self.t_matrix.right_hand_side()) sys.stdout.write(' done\n') sys.stdout.flush() elif self.solver_type == 'gmres': rhs = self.t_matrix.right_hand_side() start_time = time.time() def status_msg(rk): global iter_num iter_msg = ('Solve (GMRES) : Iter ' + str(iter_num) + ' | Rel. residual: ' + "{:.2e}".format(np.linalg.norm(rk)) + ' | elapsed: ' + str(int(time.time() - start_time)) + 's') sys.stdout.write('\r' + iter_msg) iter_num += 1 global iter_num iter_num = 0 b, info = scipy.sparse.linalg.gmres( self.master_matrix.linear_operator, rhs, rhs, tol=self.solver_tolerance, callback=status_msg) # sys.stdout.write('\n') else: raise ValueError( 'This solver type is currently not implemented.') for iS, particle in enumerate(self.particle_list): i_iS = self.layer_system.layer_number(particle.position[2]) n_iS = self.layer_system.refractive_indices[i_iS] k = coord.angular_frequency( self.initial_field.vacuum_wavelength) * n_iS loz, upz = self.layer_system.lower_zlimit( i_iS), self.layer_system.upper_zlimit(i_iS) particle.scattered_field = fldex.SphericalWaveExpansion( k=k, l_max=particle.l_max, m_max=particle.m_max, kind='outgoing', reference_point=particle.position, lower_z=loz, upper_z=upz) particle.scattered_field.coefficients = b[ self.master_matrix.index_block(iS)]