PDielec.Calculator
#
Calculator module.
Module Contents#
Functions#
|
Routine to calculate the error in the Bruggeman method. |
|
Bruggeman method using scalar quantities (suitable for powell minimisation method). |
|
Bruggeman method using tensor quantities. |
|
Calculate the absorption from the frequencies and intensities using a Lorentzian. |
|
Calculate the effective constant permittivity using a Mie scattering approach. |
Return the averaged tensor. |
|
|
Calculate the effective constant permittivity using the averaged permittivity method. |
|
Calculate the effective constant permittivity using the method of Balan. |
|
Calculate an iteration of the Bruggeman method.. |
|
Calculate the effective constant permittivity using the method of Bruggeman (minimisation). |
|
Calculate the angle between a, b and c in degrees. |
|
Calculate the scattering from bubbles embedded in a possibly complex dielectric at a given frequency. |
|
Calculate centre of mass. |
|
Calculate the distance between a and b. |
|
Calculate energy distribution in the phonon modes. |
|
Calculate the permittivity from the refractive index. |
|
Calculate the refractive index from the dielectric constant. |
|
Calculate the refractive index from the dielectric constant. |
Calculate a size effect using Equations 10.38 and 10.39 in Sihvola. |
|
|
Calculate the torsion angle between a, b, c and d in degrees. |
Return a true element from the symbol. |
|
|
Calculate the effective constant permittivity using the Coherent method. |
|
Construct the projection operator for the molecule defined by atoms, xyz, masses. |
Determine the euler angles of a rotation matrix. |
|
|
Determine the unique direction of the shape from the shape data. |
|
Apply a passive Euler rotation to a vector. |
|
Return a pool of processors given the number of cpus and whether threading is requested. |
|
Apply a Hodrick Prescott filter to the spectrum in x, y. |
|
Calculate the IR intensities from the trace of the oscillator strengths. |
Initialise a complex 3x3 tensor with the given diagonal components. |
|
|
Initialise a real 3x3 tensor with the given diagonal components. |
|
Initialise a 3x3 tensor with the ellipsoid depolarisation matrix. |
Initialise a 3x3 tensor with the needle depolarisation matrix. |
|
Initialise a 3x3 tensor with the plate depolarisation matrix, returns a tensor. |
|
Initialise a 3x3 tensor with the sphere depolarisation matrix, returns a tensor. |
|
Initialise a 3x3 tensor to a unit tensor. |
|
|
Calculate the low frequency permittivity or zero frequency permittivity. |
|
Apply the nonanalytic correction to the dynamical matrix and calculate the LO frequencies. |
|
Calculate the effective constant permittivity using the Maxwell Garnett method. |
|
Calculate the effective constant permittivity using the Maxwell Garnett method. |
|
Calculate the effective constant permittivity using a Mie scattering approach. |
|
Transform from mass weighted coordinates to xyz. |
Orthogonalise the list of projection operators ps. |
|
|
Calculate oscillator strengths from the normal modes and the born charges. |
|
Calculate the atr s and p reflectance. |
|
Take the field directions in efield and use each direction to calculate a random rotation about that axis. |
When a new worker process is created, the affinity is set to all CPUs. |
|
|
Set default number of threads. |
|
Solve the effective medium equations. |
|
Calculate the effective wavenumber based on Waterman-Truell scattering model. |
- PDielec.Calculator._brug_iter_error(epsbr, eps1, eps2, shape, L, f1, size)[source]#
Routine to calculate the error in the Bruggeman method.
Parameters#
- epsbrlist of 2 floats
The real and imaginary components of the scalar permittivity
- eps1tensor (3x3)
Permittivity of phase 1
- eps2tensor (3x3)
Permittivity of phase 2
- shapestring
The shape descriptor
- L3x3 tensor
Depolarisation tensor
- f1float
The volume fraction of phase 1
- sizefloat
The size of the particle
Returns#
float : The error associated with the current values of the permittivities
- PDielec.Calculator._brug_minimise_scalar(variables, eps1, eps2, shape, L, f1, size)[source]#
Bruggeman method using scalar quantities (suitable for powell minimisation method).
Parameters#
- variableslist of 2 floats
The real and imaginary components of the scalar permittivity
- eps1tensor (3x3)
Permittivity of phase 1
- eps2tensor (3x3)
Permittivity of phase 2
- shapestring
The particle shape
- L3x3 tensor
Depolarisation tensor
- f1float
The volume fraction of phase 1
- sizefloat
The size of the particle
Returns#
float : The error associated with the current values of the permittivities
- PDielec.Calculator._brug_minimise_tensor(variables, eps1, eps2, shape, L, f1, size)[source]#
Bruggeman method using tensor quantities.
Parameters#
- variableslist of 2 floats
The real and imaginary components of the scalar permittivity
- eps1tensor (3x3)
Permittivity of phase 1
- eps2tensor (3x3)
Permittivity of phase 2
- shapestring
The particle shape
- L3x3 tensor
Depolarisation tensor
- f1float
The volume fraction of phase 1
- sizefloat
The size of the particle
Returns#
float : The error associated with the current values of the permittivities
- PDielec.Calculator.absorption_from_mode_intensities(f, modes, frequencies, sigmas, intensities)[source]#
Calculate the absorption from the frequencies and intensities using a Lorentzian.
Parameters#
- ffloat
The frequency of the absorption in cm-1.
- modeslist of ints
A list of the modes
- frequencieslist of reals
Mode frequencies (cm-1)
- sigmaslist of reals
Mode widths in cm-1
- intensitieslist of reals
Mode intensities (D2/A2/amu).
Returns#
- float
The molar absorption coefficient at f, in L/mol/cm.
Notes#
The number 4225.6 converts the units of D2/A2/amu to L mole-1 cm-1 cm-1.
- PDielec.Calculator.anisotropic_mie_scattering(dielectric_medium, crystal_permittivity, shape, L, vf, size, size_mu, size_distribution_sigma)[source]#
Calculate the effective constant permittivity using a Mie scattering approach.
Parameters#
- dielectric_mediumarray_like (3x3)
The dielectric constant tensor of the medium.
- crystal_permittivityarray_like (3x3)
The total frequency dielectric constant tensor at the current frequency.
- shapestr
The name of the current shape (NOT USED).
- Larray_like
The shape’s depolarisation matrix (NOT USED).
- sizefloat
The dimensionless size parameter for the frequency under consideration.
- size_mufloat
The particle size in microns
- size_distribution_sigmafloat
The log normal value of sigma.
- vffloat
The volume fraction of filler.
Notes#
Mie only works for spherical particles, so the shape, and L parameters are ignored.
Returns#
- effective_dielectric_constantfloat
The effective dielectric constant.
- PDielec.Calculator.average_tensor(t)[source]#
Return the averaged tensor.
Parameters#
- ttensor (3x3)
The tensor
Returns#
- outputtensor
The averaged tensor.
- PDielec.Calculator.averaged_permittivity(dielectric_medium, crystal_permittivity, shape, L, vf, size)[source]#
Calculate the effective constant permittivity using the averaged permittivity method.
Parameters#
- dielectric_mediumarray_like (3x3)
The dielectric constant tensor of the medium.
- crystal_permittivityarray_like (3x3)
The total frequency dielectric constant tensor at the current frequency.
- shapestr
The name of the current shape.
- Larray_like (3x3)
The shape’s depolarisation matrix.
- vffloat
Volume fraction
- sizefloat
The dimensionless size parameter for the frequency under consideration (not used).
Returns#
- array_like
The effective dielectric constant.
- PDielec.Calculator.balan(dielectric_medium, crystal_permittivity, shape, L, vf, size)[source]#
Calculate the effective constant permittivity using the method of Balan.
Parameters#
- dielectric_mediumarray_like (3x3)
The dielectric constant tensor of the medium.
- crystal_permittivityarray_like (3x3)
The total frequency dielectric constant tensor at the current frequency.
- shapestr
The name of the current shape.
- Larray_like (3x3)
The shape’s depolarisation matrix.
- vffloat
Volume fraction
- sizefloat
The dimensionless size parameter for the frequency under consideration (not used).
Returns#
- array_like
The effective dielectric constant.
- PDielec.Calculator.bruggeman_iter(eps1, eps2, shape, L, f2, size, epsbr)[source]#
Calculate an iteration of the Bruggeman method..
Parameters#
- eps1array_like (3x3)
The dielectric constant tensor of medium 1.
- eps2array_like (3x3)
The dielectric constant tensor of medium 2.
- shapestr
The name of the current shape.
- Larray_like (3x3)
The shape’s depolarisation matrix.
- f2float
The volume fraction of component 2.
- sizefloat
The dimensionless size parameter for the frequency under consideration.
- epsbrfloat (3x3)
An initial guess at the solution.
Returns#
- tensor (3x3)
The effective dielectric constant.
Notes#
This function applies homogenization formalisms to active dielectric composite materials as discussed in the work of Tom G. Mackay and Akhlesh Lakhtakia.
- PDielec.Calculator.bruggeman_minimise(eps1, eps2, shape, L, f2, size, epsbr)[source]#
Calculate the effective constant permittivity using the method of Bruggeman (minimisation).
Parameters#
- eps1array_like (3x3)
The dielectric constant tensor of medium 1.
- eps2array_like (3x3)
The dielectric constant tensor of medium 2.
- shapestr
The name of the current shape.
- Larray_like (3x3)
The shape’s depolarisation matrix.
- f2float
The volume fraction of component 2.
- sizefloat
The dimensionless size parameter for the frequency under consideration.
- epsbrfloat (3x3)
An initial guess at the solution.
Returns#
- tensor (3x3)
The effective dielectric constant.
Notes#
This function applies homogenization formalisms to active dielectric composite materials as discussed in the work of Tom G. Mackay and Akhlesh Lakhtakia.
- PDielec.Calculator.calculate_angle(a, b, c)[source]#
Calculate the angle between a, b and c in degrees.
The bond is a-b-c, b is the central atom
Parameters#
- alist of reals
Coordinates of a
- blist of reals
Coordinates of b
- clist of reals
Coordinates of c
Returns#
- float
The angle between a-b-c in degrees
- PDielec.Calculator.calculate_bubble_refractive_index(v_cm1, ri_medium, vf, radius_mu)[source]#
Calculate the scattering from bubbles embedded in a possibly complex dielectric at a given frequency.
Parameters#
- v_cm1float
The frequency in cm-1.
- ri_mediumfloat
The refractive index of the medium.
- vffloat
The volume fraction of bubbles.
- radius_mufloat
The radius of the bubbles in microns.
Returns#
- effective_dielectric_constantfloat
The effective dielectric constant.
- ri_mediumfloat
The refractive index of the medium
Notes#
This function calculates the scattering from bubbles embedded in a dielectric medium, which can have a complex refractive index. It considers the frequency of interest, the refractive index of the medium, the volume fraction of bubbles, and the radius of the bubbles to calculate the effective dielectric constant and its associated refractive index.
- PDielec.Calculator.calculate_centre_of_mass(xyzs, masses)[source]#
Calculate centre of mass.
Parameters#
- xyzslist of xyz coordinates of the atoms
The xyz coordinates
- masseslist of the atomic masses
The list of atomic masses in amu
Returns#
mass : float - The total mass cm : vector (3) - the coordinates of the centre of mass
- PDielec.Calculator.calculate_distance(a, b)[source]#
Calculate the distance between a and b.
Parameters#
- alist of reals
Coordinates of a
- blist of reals
Coordinates of b
Returns#
- float
The distance between a and b
- PDielec.Calculator.calculate_energy_distribution(cell, frequencies, normal_modes, debug=False)[source]#
Calculate energy distribution in the phonon modes.
Parameters#
- cellunit cell object
The unit cell object
- frequenciesarray_like
The frequencies in cm-1.
- normal_modesarray_like
The mass weighted normal modes.
- debugboolean
True for debugging
- PDielec.Calculator.calculate_permittivity(refractive_index, debug=False)[source]#
Calculate the permittivity from the refractive index.
Parameters#
- refractive_indexcomplex
The refractive index from which the permittivity is calculated.
- debugboolean
True for debugging information
Returns#
- complex
The calculated permittivity.
- PDielec.Calculator.calculate_refractive_index(dielectric, debug=False)[source]#
Calculate the refractive index from the dielectric constant.
Calculate the trace of the dielectric and calculate both square roots. Then choose the root with the largest imaginary component. This obeys the Kramers-Konig requirements.
Parameters#
- dielectriccomplex
The permittivity
- debugboolean
True for debugging information
Returns#
- complex
The refractive index calculated from the dielectric constant.
Notes#
The calculation of the refractive index from the dielectric constant involves the trace of the dielectric tensor and the selection of the square root with the largest imaginary component.
- PDielec.Calculator.calculate_refractive_index_scalar(dielectric_scalar, debug=False)[source]#
Calculate the refractive index from the dielectric constant.
Calculate the trace of the dielectric and calculate both square roots. Then choose the root with the largest imaginary component. This obeys the Konig-Kramer requirements.
Parameters#
- dielectric_scalarcomplex
The permittivity
- debugboolean
True for debugging information
Returns#
- complex
The refractive index calculated from the dielectric constant.
Notes#
The calculation of the refractive index from the dielectric constant involves the trace of the dielectric tensor and the selection of the square root with the largest imaginary component to satisfy the Konig-Kramer conditions.
- PDielec.Calculator.calculate_size_factor(x)[source]#
Calculate a size effect using Equations 10.38 and 10.39 in Sihvola.
If x is small the result is close to 1
Parameters#
- xfloat
The size of the particle
Returns#
- float
Size effect
- PDielec.Calculator.calculate_torsion(a, b, c, d)[source]#
Calculate the torsion angle between a, b, c and d in degrees.
Parameters#
- alist of reals
Coordinates of a
- blist of reals
Coordinates of b
- clist of reals
Coordinates of c
- dlist of reals
Coordinates of d
Returns#
- float
The torsion angle between a-b-c-d in degrees
- PDielec.Calculator.cleanup_symbol(s)[source]#
Return a true element from the symbol.
Parameters#
- sstr
The element symbol to be cleaned up
Returns#
- str
The cleaned symbol
- PDielec.Calculator.coherent(dielectric_medium, crystal_permittivity, shape, L, vf, size, dielectric_apparent)[source]#
Calculate the effective constant permittivity using the Coherent method.
Parameters#
- dielectric_mediumtensor (3x3)
The dielectric constant tensor of the medium.
- crystal_permittivitytensor (3x3)
The total frequency dielectric constant tensor at the current frequency.
- shapestr
The name of the current shape.
- Lmatrix (3x3)
The shape’s depolarisation matrix.
- sizefloat
The dimensionless size parameter for the frequency under consideration.
- vffloat
The volume fraction of filler.
- dielectric_apparent3x3 array of floats
The current estimate of the dielectric
Returns#
- tensor
The effective dielectric constant.
- PDielec.Calculator.construct_projection_operator(atoms, xyzs, masses, nats)[source]#
Construct the projection operator for the molecule defined by atoms, xyz, masses.
Parameters#
- atomslist strings
The atom types
- xyzslist of atom coordinates vector (3)
The list of coordinates
- masseslist of atom masses
The list of atomic masses in amu
- natsint
The number of atoms
Returns#
- type
Description of the returned object.
- PDielec.Calculator.determineEulerAngles(R)[source]#
Determine the euler angles of a rotation matrix.
Parameters#
- Rtensor (3,3)
The rotation matrix
Returns#
- theta, phi, psifloat
The Euler angles
- PDielec.Calculator.direction_from_shape(data, reader)[source]#
Determine the unique direction of the shape from the shape data.
Parameters#
- datalist of strings
Data may contain a miller indices which defines a surface, e.g., (1,1,-1), or a direction as a miller direction vector, e.g., [1,0,-1].
- readera reader objecy
The reader is used to get the unit-cell
Returns#
- outvector (3)
Description of the unique direction determined from the data.
- PDielec.Calculator.euler_rotation(vector, theta, phi, psi)[source]#
Apply a passive Euler rotation to a vector.
Parameters#
- vectorvector (3)
The vector to be rotated
- thetafloat
The angle theta
- phifloat
The angle phi
- psifloat
The angle psi
Returns#
vector (3)
Notes#
A passive Euler rotation refers to the rotation of the coordinate system while the vector remains fixed. This operation is often used in physics and engineering to describe the orientation of an object with respect to a reference coordinate system.
- PDielec.Calculator.get_pool(ncpus, threading, initializer=None, initargs=None, debugger=None)[source]#
Return a pool of processors given the number of cpus and whether threading is requested.
Parameters#
- ncpusint
the number of processors
- threadingbool
true if threading is to be used
- initializerfunction
Function to be called before getting the pool
- initargsfunction arguments
Any other parameters
- debuggera debugger object
A debugger object
Returns#
pool : the pool of processors
- PDielec.Calculator.hodrick_prescott_filter(y, damping, lambda_value, niters)[source]#
Apply a Hodrick Prescott filter to the spectrum in x, y.
Parameters#
- yarray_like
The experimental absorption data.
- dampingfloat
The damping factor used to damp the iterations.
- lambda_valuefloat
The chosen smoothing factor.
- nitersint
The number of iterations
Returns#
list of floats : The new spectrum
Notes#
Based on ideas in the thesis of Mayank Kaushik (University Adelaide).
- PDielec.Calculator.infrared_intensities(oscillator_strengths)[source]#
Calculate the IR intensities from the trace of the oscillator strengths.
Returns#
- np.array
An array of the calculated IR intensities in units of (D/A)^2/amu.
- PDielec.Calculator.initialise_complex_diagonal_tensor(reals)[source]#
Initialise a complex 3x3 tensor with the given diagonal components.
Parameters#
- realslist
A list of 3 real numbers for the diagonals.
Returns#
- array
The returned tensor is a complex 3x3 array.
- PDielec.Calculator.initialise_diagonal_tensor(reals)[source]#
Initialise a real 3x3 tensor with the given diagonal components.
Parameters#
- realslist
A list of 3 real numbers for the diagonals.
Returns#
- array
The returned tensor is a real 3x3 array.
- PDielec.Calculator.initialise_ellipsoid_depolarisation_matrix(unique, aoverb)[source]#
Initialise a 3x3 tensor with the ellipsoid depolarisation matrix.
Parameters#
- uniquelist of 3 floats
Unique direction for ellipsoid
- aoverbfloat
The ratio of a / b (the ratio of the principle axis lengths of the ellipsoid)
Returns#
- np.array
The Ellipsoid dpolarisation tensor
- PDielec.Calculator.initialise_needle_depolarisation_matrix(unique)[source]#
Initialise a 3x3 tensor with the needle depolarisation matrix.
Parameters#
- unique: np.array
The unique direction of the needle
Returns#
- ndarray
A 3x3 tensor initialized with the plate depolarisation matrix.
- PDielec.Calculator.initialise_plate_depolarisation_matrix(normal)[source]#
Initialise a 3x3 tensor with the plate depolarisation matrix, returns a tensor.
Parameters#
- normalnp.array
A 3 vector giving the normal direction of the plate
Returns#
- ndarray
A 3x3 tensor initialized with the plate depolarisation matrix.
- PDielec.Calculator.initialise_sphere_depolarisation_matrix()[source]#
Initialise a 3x3 tensor with the sphere depolarisation matrix, returns a tensor.
Parmeters#
None
Returns#
- ndarray
A 3x3 tensor representing the sphere depolarisation matrix.
- PDielec.Calculator.initialise_unit_tensor()[source]#
Initialise a 3x3 tensor to a unit tensor.
Parameters#
None
Returns#
- array
The returned tensor is a 3x3 array.
- PDielec.Calculator.ionic_permittivity(mode_list, oscillator_strengths, frequencies, volume)[source]#
Calculate the low frequency permittivity or zero frequency permittivity.
Parameters#
- mode_listlist
List of integers giving the active modes
- oscillator_strengthsarray_like
Oscillator strengths, in atomic units.
- frequenciesarray_like
Frequencies, in atomic units.
- volumefloat
Volume, in atomic units.
Returns#
- 3x3 np.array
The calculated low frequency permittivity or zero frequency permittivity.
Notes#
The calculation of low frequency permittivity or zero frequency permittivity requires oscillator strengths, frequencies, and volume all to be specified in atomic units.
- PDielec.Calculator.longitudinal_modes(frequencies, normal_modes, born_charges, masses, epsilon_inf, volume, qlist, reader)[source]#
Apply the nonanalytic correction to the dynamical matrix and calculate the LO frequencies.
Parameters#
- frequenciesarray_like
The frequencies (f) in atomic units.
- normal_modesarray_like
The mass weighted normal modes (U).
- born_chargesarray_like
The born charges (Z) stored as [[Z1x, Z1y, Z1z], [Z2x, Z2y, Z2z], [Z3x, Z3y, Z3z]], where 1, 2, 3 are the directions of the field and x, y, z are the coordinates of the atom.
- massesarray_like
The atomic masses in atomic units
- epsilon_inf3x3 array
Epsilon infinity
- volumefloat
volume in atomic units
- qlistlist
A list of direction vectors.
- readera reader object
a reader object
Returns#
- list
A list of (real) frequencies in atomic units. Any imaginary frequencies are set to 0.
Notes#
If projection was requested in the reader, the correction is modified to ensure translational invariance.
- PDielec.Calculator.maxwell(dielectric_medium, crystal_permittivity, shape, L, vf, size)[source]#
Calculate the effective constant permittivity using the Maxwell Garnett method.
Parameters#
- dielectric_mediumtensor (3x3)
The dielectric constant tensor of the medium.
- crystal_permittivitytensor (3x3)
The total frequency dielectric constant tensor at the current frequency.
- shapestr
The name of the current shape.
- Lmatrix (3x3)
The shape’s depolarisation matrix.
- sizefloat
The dimensionless size parameter for the frequency under consideration.
- vffloat
The volume fraction of filler.
Returns#
- tensor (3x3)
The effective dielectric constant.
- PDielec.Calculator.maxwell_sihvola(dielectric_medium, crystal_permittivity, shape, L, vf, size)[source]#
Calculate the effective constant permittivity using the Maxwell Garnett method.
Parameters#
- dielectric_mediumtensor (3x3)
The dielectric constant tensor of the medium.
- crystal_permittivitytensor (3x3)
The total frequency dielectric constant tensor at the current frequency.
- shapestr
The name of the current shape.
- Lmatrix (3x3)
The shape’s depolarisation matrix.
- sizefloat
The dimensionless size parameter for the frequency under consideration.
- vffloat
The volume fraction of filler.
Returns#
- tensor
The effective dielectric constant.
- PDielec.Calculator.mie_scattering(dielectric_medium, crystal_permittivity, shape, L, vf, size, size_mu, size_distribution_sigma)[source]#
Calculate the effective constant permittivity using a Mie scattering approach.
Parameters#
- dielectric_mediumarray_like (3x3)
Dielectric constant tensor of the medium.
- crystal_permittivityarray_like (3x3)
Total frequency dielectric constant tensor at the current frequency.
- shapestr, optional
The name of the current shape (NOT USED).
- Larray_like, optional
The shape’s depolarisation matrix (NOT USED).
- sizefloat
The dimensionless size parameter for the frequency under consideration.
- size_mufloat
The particle size in microns
- size_distribution_sigmafloat
The log normal value of sigma.
- vffloat
The volume fraction of filler.
Returns#
- float or array_like
The effective dielectric constant.
Notes#
In this method, the MG method is used to calculate the averaged effective permittivity. Then, the permittivity of the isotropic sphere that would give the same average permittivity is calculated. Finally, the Mie scattering of that sphere is calculated. The routine returns the effective dielectric constant.
- PDielec.Calculator.normal_modes(masses, mass_weighted_normal_modes)[source]#
Transform from mass weighted coordinates to xyz.
Note this returns an array object. The returned normal modes have NOT been renormalized. The input masses are in atomic units. The output normal modes are in atomic units.
Parameters#
- massesarray of floats
The masses in atomic units
- mass_weighted_normal_modesarray of floats
The mass weighted normal modes
Returns#
- np.array
An array of xyz coordinates derived from mass-weighted coordinates.
Notes#
The transformation to xyz coordinates does not involve renormalization of the normal modes. Both the input masses and the output normal modes are in atomic units.
- PDielec.Calculator.orthogonalise_projection_operator(ps)[source]#
Orthogonalise the list of projection operators ps.
Use Gramm Schmidt orthogonalisation to perform the operation
Parameters#
- pslist of operators
The projection operators to orthogonalise.
Returns#
- outlist of operators
The orthogonalised projection operators.
Notes#
This function orthogonalises a set of projection operators.
- PDielec.Calculator.oscillator_strengths(normal_modes, born_charges)[source]#
Calculate oscillator strengths from the normal modes and the born charges.
Parameters#
- normal_modesarray_like
Normal modes are in the mass weighted coordinate system and normalized.
- born_chargesarray_like
Born charges are in electrons, so in atomic units.
Returns#
- array_like
Oscillator strengths calculated from the given normal modes and born charges.
- PDielec.Calculator.reflectance_atr(ns, n0, theta, atrSPolFraction)[source]#
Calculate the atr s and p reflectance.
Parameters#
- nscomplex
The complex permittivity of the effective medium.
- n0complex
The permittivity of atr material.
- thetafloat
The angle of incidence in degrees.
- atrSPolFractionfloat
The fraction of S wave to be considered. The amount of P wave is 1 - atrSPolFraction.
Returns#
- rsfloat
The s-wave Fresnel amplitude.
- rpfloat
The p-wave Fresnel amplitude.
- PDielec.Calculator.rodridgues_rotations(efield)[source]#
Take the field directions in efield and use each direction to calculate a random rotation about that axis.
Use the field (which is a random unit vector in xyz space) to generate an orthogonal rotation matrix using the Rodrigues rotation formula A = I3.cos(theta) + (1-cos(theta)) e . eT + ex sin(theta), where I3 is a unit matrix, e is the direction, and ex is the cross product matrix:
\[\begin{split}\begin{align} &\text{ex} = \begin{pmatrix} 0 & -e3 & e2 \\ e3 & 0 & -e1 \\ -e2 & e1 & 0 \end{pmatrix} \end{align}\end{split}\]Parameters#
- efieldarray_like
The field directions, where each direction is used to calculate a random rotation about that axis. Assumes the field is real.
Returns#
- list of ndarray
A list of rotation matrices for each direction in efield.
Notes#
The input field is assumed to be real.
Output is a list of 3x3 rotation matrices.
- PDielec.Calculator.set_affinity_on_worker()[source]#
When a new worker process is created, the affinity is set to all CPUs.
- PDielec.Calculator.set_no_of_threads(nthreads)[source]#
Set default number of threads.
Parameters#
- nthreadsint the number of threads to be used
The number of threads to be applied
The environment is modified to set the most common environment variables for the number of threads a BLAS implementation will use. BLAS implementations include: MKL, OPENBLAS, OMP, NUMEXPR, BLIS and VECLIB
- PDielec.Calculator.solve_effective_medium_equations(method, vf, size_mu, size_distribution_sigma, matrixPermittivityFunction, shape, L, concentration, atrPermittivity, atrTheta, atrSPol, bubble_vf, bubble_radius, previous_solution_shared, atuple)[source]#
Solve the effective medium equations.
Parameters#
- methodstr
The method to be used, options include bruggeman, balan, maxwell, maxwell-garnet, averagedpermittivity, maxwell-sihvola, coherent, bruggeman-minimise, mie, anisotropic-mie.
- vffloat
The volume fraction of dielectric.
- size_mufloat
The particle size in micron
- size_distribution_sigmafloat
The width of the size distribution.
- matrixPermittivityFunctionfunction
Function returning the matrix permittivity at a frequency.
- shapestr
The shape of the particles.
- Larray
The depolarisation matrix.
- concentrationfloat
The concentration of particles.
- atrPermittivityfloat
The permittivity of the ATR substrate.
- atrThetafloat
The ATR angle of incidence.
- atrSPolstr
The ATR polarisation.
- bubble_vffloat
Volume fraction of bubbles.
- bubble_radiusfloat
The radius of bubbles.
- previous_solution_sharedbool
Use the previous solution to speed up iterations in the case of Bruggeman and coherent methods.
- atupletuple
A tuple containing frequency in cm-1 (v_cm1) and a rank 3 tensor of the permittivity of the crystal at a given frequency (crystalPermittivity).
Returns#
- Tuple of results
v_cm1 method size_mu size_distribution_sigma shape data trace absorption_coefficient molar_absorption_coefficient spatr
- PDielec.Calculator.waterman_truell_scattering(lambda_vacuum_nm, N_nm, radius_nm, ri_medium)[source]#
Calculate the effective wavenumber based on Waterman-Truell scattering model.
Parameters#
- lambda_vacuum_nmfloat
Wavelength of the incident light in vacuum, in nanometers.
- N_nmfloat
Number density of scatterers, in inverse cubic nanometers.
- radius_nmfloat
Radius of a scatterer, in nanometers.
- ri_mediumfloat
Refractive index of the medium.
Returns#
- new_kcomplex
The complex effective wavenumber calculated based on the Waterman-Truell scattering model, which accounts for the multiple scattering effects among the particles in the medium.
Notes#
The Waterman-Truell scattering model is a method to calculate the effective wavenumber of a medium containing scatterers. This calculation takes into account the size of the scatterers, their number density, and the refractive index of the medium.
The model utilizes the Mie scattering solutions to evaluate the forward scattering amplitudes, which are then used to approximate the effective wavenumber for wave propagation in the medium. This model is particularly useful in the study of wave scattering in composite materials and biological tissues.