Copyright 2016, Vinothan N. Manoharan, Victoria Hwang, Annie Stephenson
This file is part of the structural-color python package.
This package is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
This package is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with this package. If not, see http://www.gnu.org/licenses/.
One of the advantages of the bulk montecarlo model is that we can sample phase functions and scattering lengths for spheres that contain different particle assemblies. This means we can simulate the reflectance of bulk films made of mixtures of spheres with different colors, allowing us to simulate color mixing using the bulk Monte Carlo model.
Below is an example that calculates a reflectance spectrum from a bulk film made of a mixture of two types of spheres of different colors.
In [1]:
%matplotlib inline
import numpy as np
import time
import structcol as sc
import structcol.refractive_index as ri
from structcol import montecarlo as mc
from structcol import detector as det
from structcol import phase_func_sphere as pfs
import matplotlib.pyplot as plt
import seaborn as sns
from scipy.misc import factorial
import os
This is essentially the same as running MC for a sphere as described in montecarlo_tutorial.ipynb, only we return a few extra parameters from calc_refl_trans() and use them to calculate the phase function, scattering coefficient, and absorption coefficient for the bulk Monte Carlo simulation.
We have to set a few extra parameters for the bulk simulation
In [2]:
# Properties of the source
wavelengths = sc.Quantity(np.arange(400., 801.,10),'nm') # wavelengths at which to calculate reflectance
# Geometric properties of the sample
particle_radii = sc.Quantity([130, 160],'nm') # radii of the two species of particles
volume_fraction_bulk = sc.Quantity(0.63,'') # volume fraction of the spheres in the bulk film
volume_fraction_particles = sc.Quantity(0.55, '') # volume fraction of the particles in the sphere boundary
sphere_boundary_diameter = sc.Quantity('10 um') # diameter of sphere boundary in bulk film
bulk_thickness = sc.Quantity('50 um') # thickness of the bulk film
boundary = 'sphere' # geometry of sample
boundary_bulk = 'film' # geometry of the bulk sample
# Refractive indices
n_particle = ri.n('vacuum', wavelengths) # refractive index of particle
n_matrix = ri.n('polystyrene', wavelengths) + 2e-5*1j # refractive index of matrix
n_matrix_bulk = ri.n('vacuum', wavelengths) # refractive index of the bulk matrix
n_medium = ri.n('vacuum', wavelengths) # refractive index of medium outside the bulk sample.
# Monte Carlo parameters
ntrajectories = 2000 # number of trajectories to run with a spherical boundary
nevents = 300 # number of scattering events for each trajectory in a spherical boundary
ntrajectories_bulk = 2000 # number of trajectories to run in the bulk film
nevents_bulk = 300 # number of events to run in the bulk film
# Properties that should not need to be changed
z_low = sc.Quantity('0.0 um') # sets trajectories starting point
sns.set_style('white') # sets white plotting background
In [3]:
p_bulk = np.zeros((particle_radii.size, wavelengths.size, 200))
reflectance_sphere = np.zeros(wavelengths.size)
mu_scat_bulk = sc.Quantity(np.zeros((particle_radii.size, wavelengths.size)),'1/um')
mu_abs_bulk = sc.Quantity(np.zeros((particle_radii.size, wavelengths.size)),'1/um')
for j in range(particle_radii.size):
# print radius to keep track of where we are in calculation
print('particle radius: ' + str(particle_radii[j]))
for i in range(wavelengths.size):
# caculate the effective index of the sample
n_sample = ri.n_eff(n_particle[i], n_matrix[i], volume_fraction_particles)
# Calculate the phase function and scattering and absorption coefficients from the single scattering model
# (this absorption coefficient is of the scatterer, not of an absorber added to the system)
p, mu_scat, mu_abs = mc.calc_scat(particle_radii[j], n_particle[i], n_sample,
volume_fraction_particles, wavelengths[i])
# Initialize the trajectories
r0, k0, W0 = mc.initialize(nevents, ntrajectories, n_matrix_bulk[i], n_sample,
boundary, sample_diameter = sphere_boundary_diameter)
r0 = sc.Quantity(r0, 'um')
k0 = sc.Quantity(k0, '')
W0 = sc.Quantity(W0, '')
# Create trajectories object
trajectories = mc.Trajectory(r0, k0, W0)
# Generate a matrix of all the randomly sampled angles first
sintheta, costheta, sinphi, cosphi, _, _ = mc.sample_angles(nevents, ntrajectories, p)
# Create step size distribution
step = mc.sample_step(nevents, ntrajectories, mu_scat)
# Run photons
trajectories.absorb(mu_abs, step)
trajectories.scatter(sintheta, costheta, sinphi, cosphi)
trajectories.move(step)
# Calculate reflection and transmition
(refl_indices,
trans_indices,
_, _, _,
refl_per_traj, trans_per_traj,
_,_,_,_,
reflectance_sphere[i],
_,_, norm_refl, norm_trans) = det.calc_refl_trans(trajectories, sphere_boundary_diameter,
n_matrix_bulk[i], n_sample, boundary,
run_fresnel_traj = False,
return_extra = True)
### Calculate phase function and lscat ###
# use output of calc_refl_trans to calculate phase function, mu_scat, and mu_abs for the bulk
p_bulk[j,i,:], mu_scat_bulk[j,i], mu_abs_bulk[j,i] = pfs.calc_scat_bulk(refl_per_traj, trans_per_traj,
trans_indices,
norm_refl, norm_trans,
volume_fraction_bulk,
sphere_boundary_diameter,
n_matrix_bulk[i],
wavelengths[i],
plot=False, phi_dependent=False)
In [4]:
# sample
prob = np.array([0.5,0.5]) # fraction of each sphere color type
sphere_type_sampled = pfs.sample_concentration(prob, ntrajectories_bulk, nevents_bulk)
# plot
sns.distplot(np.ndarray.flatten(sphere_type_sampled), kde = False)
plt.xlim([1,2])
plt.ylabel('number sampled')
plt.xlabel('sphere type number')
Out[4]:
The only difference from a normal bulk reflectance calculation (see bulk_montecarlo_tutorial.ipynb) is that we use the function pfs.sample_angles_step_poly() instead of sample_angles() and sample_step()
Note that for mixtures of different sphere types, absorption only works in the bulk matrix, not in the spheres themselves. This is because sampling the different absorption lengths for different sphere types has not yet been implemented.
In [5]:
reflectance_bulk_mix = np.zeros(wavelengths.size)
for i in range(wavelengths.size):
# print the wavelength keep track of where we are in calculation
print('wavelength: ' + str(wavelengths[i]))
# Initialize the trajectories
r0, k0, W0 = mc.initialize(nevents_bulk, ntrajectories_bulk, n_medium[i], n_matrix_bulk[i], boundary_bulk)
r0 = sc.Quantity(r0, 'um')
W0 = sc.Quantity(W0, '')
k0 = sc.Quantity(k0, '')
# Sample angles and calculate step size based on sampled radii
(sintheta, costheta, sinphi, cosphi,
step, _, _) = pfs.sample_angles_step_poly(nevents_bulk, ntrajectories_bulk,
p_bulk[:,i,:],
sphere_type_sampled,
mu_scat_bulk[:,i])
# Create trajectories object
trajectories = mc.Trajectory(r0, k0, W0)
# Run photons
trajectories.absorb(mu_abs_bulk[0,i], step) # Note: we assume that all scattering events
# have the same amount of absorption
trajectories.scatter(sintheta, costheta, sinphi, cosphi)
trajectories.move(step)
# calculate reflectance
reflectance_bulk_mix[i], transmittance = det.calc_refl_trans(trajectories, bulk_thickness,
n_medium[i], n_matrix_bulk[i], boundary_bulk)
In [6]:
plt.figure()
plt.plot(wavelengths, reflectance_bulk_mix, linewidth = 3)
plt.ylim([0,1])
plt.xlim([400,800])
plt.xlabel('Wavelength (nm)')
plt.ylabel('Reflectance')
plt.title('Bulk Reflectance')
Out[6]:
In [ ]: