예제 #1
0
def pwe_to_ff_conversion(vacuum_wavelength, plane_wave_expansion):
    """Compute the far field of a plane wave expansion object.

    Args:
        vacuum_wavelength (float):                 Vacuum wavelength in length units.
        plane_wave_expansion (PlaneWaveExpansion): Plane wave expansion to convert into far field object.

    Returns:
        A FarField object containing the far field intensity.
    """
    omega = flds.angular_frequency(vacuum_wavelength)
    k = plane_wave_expansion.k
    kp = plane_wave_expansion.k_parallel
    if plane_wave_expansion.kind == 'upgoing':
        polar_angles = np.arcsin(kp / k)
    elif plane_wave_expansion.kind == 'downgoing':
        polar_angles = np.pi - np.arcsin(kp / k)
    else:
        raise ValueError('PWE type not specified')
    if any(polar_angles.imag):
        raise ValueError('complex angles are not allowed')
    azimuthal_angles = plane_wave_expansion.azimuthal_angles
    kkz2 = flds.k_z(k_parallel=kp, k=k)**2 * k
    intens = (2 * np.pi**2 / omega * kkz2[np.newaxis, :, np.newaxis] *
              abs(plane_wave_expansion.coefficients)**2).real
    srt_idcs = np.argsort(polar_angles)  # reversing order in case of downgoing
    ff = FarField(polar_angles=polar_angles[srt_idcs],
                  azimuthal_angles=azimuthal_angles)
    ff.signal = intens[:, srt_idcs, :]
    return ff
예제 #2
0
    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] * flds.angular_frequency(
            vacuum_wavelength=vacuum_wavelength)
예제 #3
0
def internal_field_piecewise_expansion(vacuum_wavelength, particle_list, layer_system):
    """Compute a piecewise field expansion of the internal field of spheres.

    Args:
        vacuum_wavelength (float):                  vacuum wavelength
        particle_list (list):                       list of smuthi.particles.Particle objects
        layer_system (smuthi.layers.LayerSystem):   stratified medium

    Returns:
        internal field as smuthi.field_expansion.PiecewiseFieldExpansion object

    """
    intfld = fldex.PiecewiseFieldExpansion()

    for particle in particle_list:
        if type(particle).__name__ == 'Sphere':
            i_part = layer_system.layer_number(particle.position[2])
            k_medium = flds.angular_frequency(vacuum_wavelength) * layer_system.refractive_indices[i_part]
            k_particle = flds.angular_frequency(vacuum_wavelength) * particle.refractive_index

            internal_field = fldex.SphericalWaveExpansion(
                k_particle,
                l_max=particle.l_max,
                m_max=particle.m_max,
                kind='regular',
                reference_point=particle.position)

            internal_field.validity_conditions.append(particle.is_inside)

            for tau in range(2):
                for l in range(1, particle.l_max+1):
                    for m in range(-l, l+1):
                        n = flds.multi_to_single_index(tau, l, m, particle.l_max, particle.m_max)
                        b_to_c = (tmt.internal_mie_coefficient(tau, l, k_medium, k_particle, particle.radius)
                                  / tmt.mie_coefficient(tau, l, k_medium, k_particle, particle.radius))
                        internal_field.coefficients[n] = particle.scattered_field.coefficients[n] * b_to_c

            intfld.expansion_list.append(internal_field)
            particle.internal_field = internal_field

    return intfld
예제 #4
0
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.fields.default_polar_angles
        azimuthal_angles (numpy.ndarray or str):    azimuthal angle values (radian)
                                                    if 'default', use smuthi.fields.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 = flds.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 np.cos(beta_P) > 0:
        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
예제 #5
0
    def testGMRES(self):
        # Parameter input ----------------------------
        vacuum_wavelength = 550
        beam_polar_angle = np.pi * 7 / 8
        beam_azimuthal_angle = np.pi * 1 / 3
        beam_polarization = 0
        beam_amplitude = 1
        beam_neff_array = np.linspace(0, 2, 501, endpoint=False)
        beam_waist = 1000
        beam_focal_point = [200, 200, 200]
        #neff_waypoints = [0, 0.5, 0.8 - 0.01j, 2 - 0.01j, 2.5, 5]
        #neff_discr = 1e-2
        # --------------------------------------------

        #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 * flds.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', log_to_terminal='nose2' not in sys.modules.keys())  # suppress output if called by nose
        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,
                                            log_to_terminal=(
                                                not sys.argv[0].endswith('nose2')))  # suppress output if called by nose
        simulation_gmres.run()
        coefficients_gmres = particle_list[0].scattered_field.coefficients

        np.testing.assert_allclose(np.linalg.norm(coefficients_lu), np.linalg.norm(coefficients_gmres), rtol=1e-5)
예제 #6
0
    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 = flds.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)]
            particle.scattered_field.validity_conditions.append(
                particle.is_outside)
            self.initial_field.validity_conditions.append(particle.is_outside)
예제 #7
0
def extinction_cross_section(initial_field, particle_list, layer_system):
    """Evaluate the extinction cross section.

    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
            - 'top':      Extinction in the positinve z-direction (top layer)
            - 'bottom':     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 = flds.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 np.cos(beta_P) > 0:
        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, _ = sf.scattered_field_pwe(vacuum_wavelength, particle_list,
                                             layer_system, i_top, kappa_P,
                                             alpha_P)

    _, pwe_scat_bottom = sf.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 = flds.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 = flds.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
예제 #8
0
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.fields.default_polar_angles
        azimuthal_angles (numpy.ndarray or str):    azimuthal angle values (radian)
                                                    if 'default', use smuthi.fields.default_azimuthal_angles

    Returns:
        A smuthi.field_expansion.FarField object of the scattered field.
    """
    omega = flds.angular_frequency(vacuum_wavelength)
    if type(polar_angles) == str and polar_angles == 'default':
        polar_angles = flds.default_polar_angles
    if type(azimuthal_angles) == str and azimuthal_angles == 'default':
        azimuthal_angles = flds.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, _ = sf.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 = 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 = sf.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 = 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
예제 #9
0
def direct_coupling_block_pvwf_mediated(vacuum_wavelength, receiving_particle,
                                        emitting_particle, layer_system,
                                        k_parallel):
    """Direct particle coupling matrix :math:`W` for two particles (via plane vector wave functions).
    For details, see: 
    Dominik Theobald et al., Phys. Rev. A 96, 033822, DOI: 10.1103/PhysRevA.96.033822 or arXiv:1708.04808 


    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
        k_parallel (numpy.array):                           In-plane wavenumber for plane wave expansion

    Returns:
        Direct coupling matrix block (numpy array).
    """
    if type(receiving_particle).__name__ != 'Spheroid' or type(
            emitting_particle).__name__ != 'Spheroid':
        raise NotImplementedError(
            'plane wave coupling currently implemented only for spheroids')

    lmax1 = receiving_particle.l_max
    mmax1 = receiving_particle.m_max
    assert lmax1 == mmax1, 'PVWF coupling requires lmax == mmax for each particle.'
    lmax2 = emitting_particle.l_max
    mmax2 = emitting_particle.m_max
    assert lmax2 == mmax2, 'PVWF coupling requires lmax == mmax for each particle.'
    lmax = max([lmax1, lmax2])
    m_max = max([mmax1, mmax2])
    blocksize1 = flds.blocksize(lmax1, mmax1)
    blocksize2 = flds.blocksize(lmax2, mmax2)

    n_medium = layer_system.refractive_indices[layer_system.layer_number(
        receiving_particle.position[2])]

    # finding the orientation of a plane separating the spheroids
    _, _, alpha, beta = spheroids_closest_points(
        emitting_particle.semi_axis_a, emitting_particle.semi_axis_c,
        emitting_particle.position, emitting_particle.euler_angles,
        receiving_particle.semi_axis_a, receiving_particle.semi_axis_c,
        receiving_particle.position, receiving_particle.euler_angles)

    # positions
    r1 = np.array(receiving_particle.position)
    r2 = np.array(emitting_particle.position)
    r21_lab = r1 - r2  # laboratory coordinate system

    # distance vector in rotated coordinate system
    r21_rot = np.dot(
        np.dot([[np.cos(beta), 0, np.sin(beta)], [0, 1, 0],
                [-np.sin(beta), 0, np.cos(beta)]],
               [[np.cos(alpha), -np.sin(alpha), 0],
                [np.sin(alpha), np.cos(alpha), 0], [0, 0, 1]]), r21_lab)
    rho21 = (r21_rot[0]**2 + r21_rot[1]**2)**0.5
    phi21 = np.arctan2(r21_rot[1], r21_rot[0])
    z21 = r21_rot[2]

    # wavenumbers
    omega = flds.angular_frequency(vacuum_wavelength)
    k = omega * n_medium
    kz = flds.k_z(k_parallel=k_parallel,
                  vacuum_wavelength=vacuum_wavelength,
                  refractive_index=n_medium)
    if z21 < 0:
        kz_var = -kz
    else:
        kz_var = kz

    # Bessel lookup
    bessel_list = []
    for dm in range(mmax1 + mmax2 + 1):
        bessel_list.append(scipy.special.jn(dm, k_parallel * rho21))

    # legendre function lookups
    ct = kz_var / k
    st = k_parallel / k
    _, pilm_list, taulm_list = sf.legendre_normalized(ct, st, lmax)

    # initialize result
    w = np.zeros((blocksize1, blocksize2), dtype=complex)

    # prefactor
    const_arr = k_parallel / (kz * k) * np.exp(1j * (kz_var * z21))

    for m1 in range(-mmax1, mmax1 + 1):
        for m2 in range(-mmax2, mmax2 + 1):
            jmm_eimphi_bessel = 4 * 1j**abs(m2 - m1) * np.exp(
                1j * phi21 * (m2 - m1)) * bessel_list[abs(m2 - m1)]
            prefactor = const_arr * jmm_eimphi_bessel
            for l1 in range(max(1, abs(m1)), lmax1 + 1):
                for l2 in range(max(1, abs(m2)), lmax2 + 1):
                    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)
                            for pol in range(2):
                                B_dag = trf.transformation_coefficients_vwf(
                                    tau1,
                                    l1,
                                    m1,
                                    pol,
                                    pilm_list=pilm_list,
                                    taulm_list=taulm_list,
                                    dagger=True)
                                B = trf.transformation_coefficients_vwf(
                                    tau2,
                                    l2,
                                    m2,
                                    pol,
                                    pilm_list=pilm_list,
                                    taulm_list=taulm_list,
                                    dagger=False)
                                integrand = prefactor * B * B_dag
                                w[n1, n2] += np.trapz(integrand, k_parallel)

    rot_mat_1 = trf.block_rotation_matrix_D_svwf(lmax1, mmax1, 0, beta, alpha)
    rot_mat_2 = trf.block_rotation_matrix_D_svwf(lmax2, mmax2, -alpha, -beta,
                                                 0)

    return np.dot(np.dot(np.transpose(rot_mat_1), w), np.transpose(rot_mat_2))
예제 #10
0
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
예제 #11
0
rho13 = np.sqrt(
    sum([(sphere1.position[i] - sphere3.position[i])**2 for i in range(3)]))
wr13 = layer_mediated_coupling_block(vacuum_wavelength, sphere1, sphere3,
                                     lay_sys)
w13 = 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 *
                             flds.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=(not sys.argv[0].endswith('nose2')
                     ))  # suppress output if called by nose
simulation_direct.run()
coefficients_direct = particle_list[0].scattered_field.coefficients

cu.enable_gpu()
simulation_lookup_linear_gpu = simul.Simulation(
def select_numerical_parameters(simulation,
                                detector="extinction cross section",
                                tolerance=1e-3,
                                max_iter=20,
                                neff_imag=1e-2,
                                neff_step=1e-2,
                                select_neff_max=True,
                                neff_max_increment=0.5,
                                neff_max=None,
                                select_neff_step=True,
                                select_multipole_cutoff=True,
                                relative_convergence=True,
                                suppress_simulation_output=True):
    """Trigger automatic selection routines for various numerical parameters.

    Args:
        simulation (smuthi.simulation.Simulation):    Simulation object from which parameters are read and into which
                                                      results are stored.
        detector (function or string):                Function that accepts a simulation object and returns a detector
                                                      value the change of which is used to define convergence.
                                                      Alternatively, use "extinction cross section" (default) to have
                                                      the extinction cross section as the detector value.
        tolerance (float):                            Relative tolerance for the detector value change.
        max_iter (int):                               Break convergence loops after that number of iterations, even if
                                                      no convergence has been achieved.
        neff_imag (float):                            Extent of the contour into the negative imaginary direction
                                                      (in terms of effective refractive index, n_eff=kappa/omega).
        neff_step (float):                            Discretization of the contour (in terms of eff. refractive
                                                      index) - if `select_neff_step` is true, this value will be
                                                      eventually overwritten. However, it is required in any case.
                                                      Default: 1e-2
        select_neff_max (logical):                    If set to true (default), the Sommerfeld integral truncation parameter
                                                      `neff_max` is determined automatically with the help of a
                                                      Cauchy convergence criterion.
        neff_max_increment (float):                   Only needed if `select_neff_max` is true.
                                                      Step size with which `neff_max` is incremented.
        neff_max (float):                             Only needed if `select_neff_max` is false.
                                                      Truncation value of contour (in terms of effective refractive
                                                      index).
        select_neff_step (logical):                   If set to true (default), the Sommerfeld integral discretization
                                                      parameter `neff_step` is determined automatically with the help of
                                                      a Cauchy convergence criterion.
        select_multipole_cutoff (logical):            If set to true (default), the multipole expansion cutoff
                                                      parameters `l_max` and `m_max` are determined automatically with
                                                      the help of a Cauchy convergence criterion.
        relative_convergence (logical):               If set to true (default), the `neff_max` convergence and the
                                                      `l_max` and `m_max` convergence routine are performed in the
                                                      spirit of relative convergence, i.e., the multipole expansion
                                                      convergence is checked again for each value of the Sommerfeld
                                                      integral truncation. This takes more time, but is required at
                                                      least in the case of flat particles near interfaces.
        suppress_simulation_output (logical):         If set to false (default), suppress any output from simulations
                                                      and only display convergence progress
    """
    print("")
    print("----------------------------------------")
    print("Starting automatic paramateter selection")
    print("----------------------------------------")

    if select_neff_max:
        converge_neff_max(simulation=simulation,
                          detector=detector,
                          tolerance=tolerance,
                          max_iter=max_iter,
                          neff_imag=neff_imag,
                          neff_step=neff_step,
                          neff_max_increment=neff_max_increment,
                          converge_lm=relative_convergence)
        neff_max = simulation.k_parallel[-1] / angular_frequency(
            simulation.initial_field.vacuum_wavelength)

    if select_multipole_cutoff:
        converge_multipole_cutoff(simulation=simulation,
                                  detector=detector,
                                  tolerance=tolerance,
                                  max_iter=max_iter)

    if select_neff_step:
        converge_neff_step(simulation=simulation,
                           detector=detector,
                           tolerance=tolerance,
                           max_iter=max_iter,
                           neff_imag=neff_imag,
                           neff_max=neff_max)
def converge_neff_max(simulation,
                      detector="extinction cross section",
                      tolerance=1e-3,
                      max_iter=20,
                      neff_imag=1e-2,
                      neff_step=2e-3,
                      neff_max_increment=0.5,
                      converge_lm=True):
    """Find a suitable truncation value for the multiple scattering Sommerfeld integral contour and update the
    simulation object accordingly.

    Args:
        simulation (smuthi.simulation.Simulation):    Simulation object
        detector (function or string):                Function that accepts a simulation object and returns a detector
                                                      value the change of which is used to define convergence.
                                                      Alternatively, use "extinction cross section" (default) to have
                                                      the extinction cross section as the detector value.
        tolerance (float):                            Relative tolerance for the detector value change.
        max_iter (int):                               Break convergence loops after that number of iterations, even if
                                                      no convergence has been achieved.
        neff_imag (float):                              Extent of the contour into the negative imaginary direction
                                                        (in terms of effective refractive index, n_eff=kappa/omega).
        neff_step (float):                              Discretization of the contour (in terms of eff. refractive
                                                        index).
        neff_max_increment (float):                     Increment the neff_max parameter with that step size
        converge_lm (logical):                          If set to true, update multipole truncation during each step
                                                        (this takes longer time, but is necessary for critical use cases
                                                        like flat particles on a substrate)

    Returns:
        Detector value for converged settings.
    """

    print("")
    print("---------------------------")
    log.write_blue("Searching suitable neff_max")

    update_contour(simulation=simulation,
                   neff_imag=neff_imag,
                   neff_max_offset=0,
                   neff_step=neff_step)

    simulation.k_parallel = flds.reasonable_Sommerfeld_kpar_contour(
        vacuum_wavelength=simulation.initial_field.vacuum_wavelength,
        layer_refractive_indices=simulation.layer_system.refractive_indices,
        neff_imag=neff_imag,
        neff_max_offset=0,
        neff_resolution=neff_step)

    neff_max = simulation.k_parallel[-1] / angular_frequency(
        simulation.initial_field.vacuum_wavelength)
    print("Starting value: neff_max=%f" % neff_max.real)

    if converge_lm:
        with log.LoggerIndented():
            current_value = converge_multipole_cutoff(
                simulation=simulation,
                detector=detector,
                tolerance=tolerance /
                2,  # otherwise, results flucutate by tolerance and convergence check is compromised
                max_iter=max_iter)
    else:
        current_value = evaluate(simulation, detector)

    for _ in range(max_iter):
        old_neff_max = neff_max
        neff_max = neff_max + neff_max_increment
        update_contour(simulation=simulation,
                       neff_imag=neff_imag,
                       neff_max=neff_max,
                       neff_step=neff_step)

        print("---------------------------------------")
        print("Try neff_max = %f" % neff_max.real)

        if converge_lm:
            with log.LoggerIndented():
                new_value = converge_multipole_cutoff(
                    simulation=simulation,
                    detector=detector,
                    tolerance=tolerance,
                    max_iter=max_iter,
                    current_value=current_value)
        else:
            new_value = evaluate(simulation, detector)

        rel_diff = abs(new_value - current_value) / abs(current_value)
        print("Old detector value:", current_value)
        print("New detector value:", new_value)
        print("Relative difference:", rel_diff)

        if rel_diff < tolerance:  # in this case: discard l_max increment
            neff_max = old_neff_max
            update_contour(simulation=simulation,
                           neff_imag=neff_imag,
                           neff_max=neff_max,
                           neff_step=neff_step)
            log.write_green(
                "Relative difference smaller than tolerance. Keep neff_max = %f"
                % neff_max.real)
            return current_value
        else:
            current_value = new_value

    log.write_red("No convergence achieved. Keep neff_max = %i" % neff_max)
    return None
예제 #14
0
lmax = 10

rx = 100
ry = -100
rz = 100

dx = 400
dy = 800
dz = -500

tau = 0
l = 2
m = -1
# --------------------------------------------

k = flds.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 = trf.translation_coefficients_svwf(tau, l, m, tau2, l2, m2, k, [dx, dy, dz])
            Ex2 += Mx * A
            Ey2 += My * A
            Ez2 += Mz * A