DielectricTensor

class DielectricTensor(configuration, dynamical_matrix, optical_spectrum, born_effective_charge, phonon_modes=None, q_point=None, charge_sum_rule_method=None, polarity_tolerance=None, broadening=None)

Constructor for the static dielectric tensor object.

Parameters:
  • configuration (BulkConfiguration) – The bulk configuration with an attached calculator for which to calculate the dielectric tensor.

  • dynamical_matrix (DynamicalMatrix) – The dynamical matrix to calculate the phonon modes from (incl. polar phonon splitting).

  • optical_spectrum (OpticalSpectrum) – The optical spectrum object from which to get the electronic contribution to the dielectric tensor.

  • born_effective_charge (BornEffectiveCharge) – The born effective charge object calculated by the user. This is used for the ionic contribution to the dielectric tensor.

  • phonon_modes (list with positive ints | All) – The phonon modes to include.
    Default: All (All phonon modes are included).

  • q_point (list of 3 floats) – The fractional q-point to use for phonon modes and polar correction.
    Default: [0.0, 0.0, 0.0]

  • charge_sum_rule_method (None | DistributeEvenly | DistributeRelative.) – Keyword specifying if a charge sum-rule should be applied to the Born effective charges. Can be None for no sum-rule application. Alternatively the charge error is either distributed evenly (DistributeEvenly) or relative (DistributeRelative) to the Born effective charges on the atoms.
    Default: None

  • polarity_tolerance (float) – Unit-less tolerance used to control when modes should be considered non-polar. This is helpful to avoid numerical errors for acoustic modes with vanishing frequency.
    Default: 1e-6

  • broadening (PhysicalQuantity of type energy) – The broadening parameter used for the spectra.
    Default: 0.25 * meV  (~2 cm⁻1)

absorption(contribution=<class 'NL.ComputerScienceUtilities.NLFlag._NLFlag.Sum'>)
Parameters:

contribution – NLFlag to determine the contributions to include. If Sum, the sum of the electronic and ionic contribution is returned. If Electronic or Ionic, only the electronic or ionic contributions are returned.

Returns:

The absorption tensor for the contribution. If contribution is Sum or Electronic, the dimension is (3, 3, N), where N is the number of energies. If contribution is Ionic, the dimension is (M, 3, 3, N), where N is the number of energies and M is the number of phonon modes.

Return type:

ndarray.

broadening()
Returns:

The broadening parameter used for the spectra.

Return type:

PhysicalQuantity of type energy.

dielectricTensor(contribution=<class 'NL.ComputerScienceUtilities.NLFlag._NLFlag.Sum'>)
Parameters:

contribution (Sum | Electronic | Ionic) – NLFlag to determine the contributions to include. If Sum, the sum of the electronic and ionic contribution is returned. If Electronic or Ionic, only the electronic or ionic contributions are returned.

Returns:

The real part of the dielectric constant tensor. If contribution is Sum or Electronic, the dimension is (3, 3, N), where N is the number of energies. If contribution is Ionic, the dimension is (M, 3, 3, N), where N is the number of energies and M is the number of phonon modes.

Return type:

Dimensionless ndarray.

electronicDielectricTensor()
Returns:

The electronic dielectric constant tensor - i.e. excluding the ionic contribution. Dimension is (3, 3, N), where N is the number of energies.

Return type:

Dimensionless ndarray.

energies()
Returns:

The energies used for the spectrum.

Return type:

PhysicalQuantity of type energy.

evaluate()
Returns:

The static dielectric constant tensor - i.e. including the ionic contribution. Dimension is (3, 3, N), where N is the number of energies.

Return type:

Dimensionless ndarray.

evaluateImaginaryPolarizability()

Function for evaluating the imaginary part of the polarizability.

Returns:

A list of imaginary parts of the polarizability tensors [electronic, static, ionic].

Return type:

PhysicalQuantity with the unit Hartree**3 * vacuum_permitivity

evaluatePolarizability()

Function for evaluating the real part of the polarizability.

Returns:

A list of real parts of the polarizability tensors [electronic, static, ionic].

Return type:

PhysicalQuantity with the unit Hartree**3 * vacuum_permitivity

extinctionCoefficient(contribution=<class 'NL.ComputerScienceUtilities.NLFlag._NLFlag.Sum'>)
Parameters:

contribution (Sum | Electronic | Ionic) – NLFlag to determine the contributions to include. If Sum, the sum of the electronic and ionic contribution is returned. If Electronic or Ionic, only the electronic or ionic contributions are returned.

Returns:

The extinction coefficient including the ionic contribution. If contribution is Sum or Electronic, the dimension is (3, 3, N), where N is the number of energies. If contribution is Ionic, the dimension is (M, 3, 3, N), where N is the number of energies and M is the number of phonon modes.

Return type:

Dimensionless ndarray.

imaginaryDielectricTensor(contribution=<class 'NL.ComputerScienceUtilities.NLFlag._NLFlag.Sum'>)
Parameters:

contribution – NLFlag to determine the contributions to include. If Sum, the sum of the electronic and ionic contribution is returned. If Electronic or Ionic, only the electronic or ionic contributions are returned.

Returns:

The imaginary part of the dielectric constant tensor. If contribution is Sum or Electronic, the dimension is (3, 3, N), where N is the number of energies. If contribution is Ionic, the dimension is (M, 3, 3, N), where N is the number of energies and M is the number of phonon modes.

Return type:

Dimensionless ndarray.

imaginaryElectronicDielectricTensor()
Returns:

The electronic dielectric constant tensor - i.e. excluding the ionic contribution. Dimension is (3, 3, N), where N is the number of energies.

Return type:

Dimensionless ndarray.

imaginaryIonicDielectricTensor()
Returns:

The ionic part of the static dielectric constant tensor. Dimension is (M, 3, 3, N), where N is the number of energies and M is the number of phonon modes.

Return type:

Dimensionless ndarray.

imaginaryStaticDielectricTensor()
Returns:

The static dielectric constant tensor - i.e. including the ionic contribution. Dimension is (3, 3, N), where N is the number of energies.

Return type:

Dimensionless ndarray.

imaginarySusceptibility(contribution=<class 'NL.ComputerScienceUtilities.NLFlag._NLFlag.Sum'>)
Parameters:

contribution – NLFlag to determine the contributions to include. If Sum, the sum of the electronic and ionic contribution is returned. If Electronic or Ionic, only the electronic or ionic contributions are returned.

Returns:

Imaginary part of the susceptibility for the contribution.

Return type:

numpy.ndarray.

infraredOscillatorStrength()
Returns:

The infrared oscillator strength tensor. Dimension is (M, 3, 3) where M is the number of phonon modes.

Return type:

Dimensionless complex ndarray.

ionicDielectricTensor()
Returns:

The ionic part of the static dielectric constant tensor. Dimension is (M, 3, 3, N), where N is the number of energies and M is the number of phonon modes.

Return type:

Dimensionless ndarray.

metatext()
Returns:

The metatext of the object or None if no metatext is present.

Return type:

str | None

nlprint(stream=None)

Print a string containing an ASCII table useful for plotting the AnalysisSpin object.

Parameters:

stream (python stream) – The stream the table should be written to.
Default: NLPrintLogger()

opticalConductivity(contribution=<class 'NL.ComputerScienceUtilities.NLFlag._NLFlag.Sum'>)
Parameters:

contribution – NLFlag to determine the contributions to include. If Sum, the sum of the electronic and ionic contribution is returned. If Electronic or Ionic, only the electronic or ionic contributions are returned.

Returns:

The optical conductivity tensor for the contribution. If contribution is Sum or Electronic, the dimension is (3, 3, N), where N is the number of energies. If contribution is Ionic, the dimension is (M, 3, 3, N), where N is the number of energies and M is the number of phonon modes.

Return type:

ndarray.

phononEigenvectors()
Returns:

The phonon eigenvectors.

Return type:

numpy.ndarray.

phononEnergies()
Returns:

The phonon eigenenergies.

Return type:

PhysicalQuantity of type energy.

phononModes()

Query method for the phonon modes.

Returns:

The list of phonon mode indices used in the calculation.

Return type:

list of int.

phononWaveNumbers()
Returns:

The phonon wave numbers calculated from the phonon frequencies.

Return type:

PhysicalQuantity of type inverse length.

polarity()
Returns:

The infrared polarity of the phonon modes. Dimension is (M, 3) where M is the number of phonon modes.

Return type:

PhysicalQuantity with the unit charge/(mass^0.5); specifically elementary_charge / electron_mass^(0.5).

qpoint()
Returns:

The q-point for which the phonon modes and Raman signals are calculated.

Return type:

ndarray with 3 floats.

realSusceptibility(contribution=<class 'NL.ComputerScienceUtilities.NLFlag._NLFlag.Sum'>)
Parameters:

contribution – NLFlag to determine the contributions to include. If Sum, the sum of the electronic and ionic contribution is returned. If Electronic or Ionic, only the electronic or ionic contributions are returned.

Returns:

Real part of the susceptibility for the contribution.

Return type:

numpy.ndarray.

reflectivity(contribution=<class 'NL.ComputerScienceUtilities.NLFlag._NLFlag.Sum'>)
Parameters:

contribution – NLFlag to determine the contributions to include. If Sum, the sum of the electronic and ionic contribution is returned. If Electronic or Ionic, only the electronic or ionic contributions are returned.

Returns:

The reflectivity tensor for the contribution. If contribution is Sum or Electronic, the dimension is (3, 3, N), where N is the number of energies. If contribution is Ionic, the dimension is (M, 3, 3, N), where N is the number of energies and M is the number of phonon modes.

Return type:

Dimensionless ndarray.

refractiveIndex(contribution=<class 'NL.ComputerScienceUtilities.NLFlag._NLFlag.Sum'>)
Parameters:

contribution (Sum | Electronic | Ionic) – NLFlag to determine the contributions to include. If Sum, the sum of the electronic and ionic contribution is returned. If Electronic or Ionic, only the electronic or ionic contributions are returned.

Returns:

The refractive index including the ionic contribution. If contribution is Sum or Electronic, the dimension is (3, 3, N), where N is the number of energies. If contribution is Ionic, the dimension is (M, 3, 3, N), where N is the number of energies and M is the number of phonon modes.

Return type:

Dimensionless ndarray.

setMetatext(metatext)

Set a given metatext string on the object.

Parameters:

metatext (str | None) – The metatext string that should be set. A value of “None” can be given to remove the current metatext.

staticDielectricTensor()
Returns:

The static dielectric constant tensor - i.e. including the ionic contribution. Dimension is (3, 3, N), where N is the number of energies.

Return type:

Dimensionless ndarray.

uniqueString()

Return a unique string representing the state of the object.

vibrationalMode()
Returns:

A VibrationalMode object for the calculated modes.

Return type:

VibrationalMode

Usage Examples

Calculate the static DielectricTensor of a configuration including both ionic and electronic contributions.

dielectric_tensor = DielectricTensor(
  configuration,
  dynamical_matrix,
  optical_spectrum,
  born_effective_charge)

The dielectric tensor will be evaluated at the same energies as defined in the input optical spectrum object.

Notes

The static (low frequency) dielectric response \(\epsilon(\omega)=\epsilon_0\epsilon_r(\omega)\) of polar materials consists of an electronic and a ionic part:

\[\epsilon_{r, \alpha \beta}(\omega) = \epsilon_{r, \alpha \beta}^{elec}(\omega) + \epsilon_{r, \alpha \beta}^{ionic}(\omega)\]

where \(\epsilon_r\) is the unit-less relative permittivity and \(\alpha, \beta\) are Cartesian directions. The electronic part can be calculated from the OpticalSpectrum while the ionic part is related to the BornEffectiveCharge and phonon modes. At high frequency only the electrons can follow the optical field and the dynamic (high frequency) dielectric constant \(\epsilon_{r, \alpha \beta}^{\infty} = \epsilon_{r, \alpha \beta}^{elec}(0)\) is obtained by evaluating the electronic dielectric constant at zero frequency.

At low frequency the ions can also respond to an optical field. The ionic part is calculated as [1]

\[\epsilon_{r, \alpha \beta}^{ionic}(\omega) = \frac{1}{\epsilon_0 \Omega_0} \sum_{\lambda} \frac{S^{\lambda}_{\alpha \beta}}{\omega_{\lambda}^2 - \omega^2 - i \Gamma \omega}\,,\]

where we introduced a broadening \(\Gamma\) for each phonon mode Lorentz oscillator, \(\Omega_0\) is the unit cell volume, \(\epsilon_0\) is the vacuum permittivity and \(S^{\lambda}_{\alpha \beta} = p^{*\lambda}_{\alpha} p^{\lambda}_{\beta}\) is the mode oscillator strength which is related to the infrared mode polarity:

\[p^{\lambda}_{\alpha} = \sum_{i, \beta} Z_{\alpha \beta}(i) e_{\beta}^{\lambda}/\sqrt{m_{i}}\]

\(Z_{\alpha \beta}(i)\) is the Born effective charge tensor for atom \(i\) with mass \(m_{i}\) and \(e_{\beta}^{\lambda}\) is the eigendisplacement of phonon mode \(\lambda\) with frequency \(\omega_{\lambda}\).

Infrared spectrum

The refractive index, \(n\) and extinction coefficient, \(\kappa\), are obtained from the real and imaginary parts of the square root of the dielectric tensor (see also OpticalSpectrum):

\[n + i\kappa = \sqrt{\epsilon_{r}}\,.\]

The reflectivity at the interface between a vacuum region and the material is given by

\[R = \frac{(1-n)^2 + \kappa^2}{(1+n)^2 + \kappa^2}\]

The infrared spectrum is obtained by performing a dense energy sampling in the range of the vibrational frequency range (typically up to 0.25 eV or correspondingly 2000 \(\textrm{cm}^{-1}\), depending on the material of interest).

The term infrared spectrum is used frequently for either the imaginary part of the dielectric tensor, \(\textrm{Im}(\epsilon_{r})\), the absorption coefficient \(\alpha = 2 \frac{\omega}{c} \kappa\) where \(c\) is the speed of light, and the transmittance, \(T = 1-r\) throughout literature.

When \(\alpha(\omega)=0\), with corresponding \(T(\omega)=1\), the material is transparent at that frequency.

Usage Example

This example script illustrates the workflow of the DielectricTensor calculation, using a classical force field for the phonons and DFT+1/2 for the electronic structure, to obtain the infrared spectrum of \(\textrm{SiO}_2\) (Quatz):

filename = 'Quartz.hdf5'

# -------------------------------------------------------------
# Bulk Configuration
# -------------------------------------------------------------
# Set up lattice
lattice = Hexagonal(4.916*Angstrom, 5.4054*Angstrom)

# Define elements
elements = [Silicon, Silicon, Silicon, Oxygen, Oxygen, Oxygen, Oxygen, Oxygen,
            Oxygen]

# Define coordinates
fractional_coordinates = [[ 0.4697        ,  0.            ,  0.            ],
                          [ 0.            ,  0.4697        ,  0.666666666667],
                          [ 0.5303        ,  0.5303        ,  0.333333333333],
                          [ 0.4135        ,  0.2669        ,  0.1191        ],
                          [ 0.2669        ,  0.4135        ,  0.547567      ],
                          [ 0.7331        ,  0.1466        ,  0.785767      ],
                          [ 0.5865        ,  0.8534        ,  0.214233      ],
                          [ 0.8534        ,  0.5865        ,  0.452433      ],
                          [ 0.1466        ,  0.7331        ,  0.8809        ]]

# Set up configuration
bulk_configuration = BulkConfiguration(
    bravais_lattice=lattice,
    elements=elements,
    fractional_coordinates=fractional_coordinates
    )

# -------------------------------------------------------------
# Optimize Geometry
# -------------------------------------------------------------
# Define a classical calculator
potentialSet = Tersoff_SiO_2007()
calculator_classical = TremoloXCalculator(parameters=potentialSet)
calculator_classical.setVerletListsDelta(0.25 * Angstrom)

bulk_configuration.setCalculator(calculator_classical)

constraints = [BravaisLatticeConstraint()]
bulk_configuration = OptimizeGeometry(
    bulk_configuration,
    max_forces=0.001 * eV / Angstrom,
    max_stress=0.01 * GPa,
    max_steps=500,
    max_step_length=0.2 * Angstrom,
    constraints=constraints,
    trajectory_filename=None,
)
nlsave(filename, bulk_configuration, object_id='bulk_configuration_classical')

# -------------------------------------------------------------
# Dynamical Matrix
# -------------------------------------------------------------
nrep = 5
dynamical_matrix = DynamicalMatrix(
    bulk_configuration,
    filename=filename,
    object_id='dynamical_matrix',
    repetitions=(nrep, nrep, nrep),
    atomic_displacement=0.01 * Angstrom,
    acoustic_sum_rule=True,
    force_tolerance=1e-08 * Hartree / Bohr**2,
    processes_per_displacement=1,
    log_filename_prefix=None,
    )
dynamical_matrix.update()

# -------------------------------------------------------------
# Set relaxed configuration with DFT calculator
# -------------------------------------------------------------
exchange_correlation = GGAHalf.PBE

k_point_sampling = MonkhorstPackGrid(
    na=6,
    nb=6,
    nc=5,
    )
numerical_accuracy_parameters = NumericalAccuracyParameters(
    density_mesh_cutoff=125.0*Hartree,
    k_point_sampling=k_point_sampling,
    occupation_method=FermiDirac(100.0*Kelvin*boltzmann_constant),
    )

basis_set = [
    BasisGGAPseudoDojo.Oxygen_Medium,
    BasisGGAPseudoDojo.Silicon_Medium,
    ]

calculator = LCAOCalculator(
    exchange_correlation=exchange_correlation,
    numerical_accuracy_parameters=numerical_accuracy_parameters,
    basis_set=basis_set,
    )

bulk_configuration.setCalculator(calculator)
nlprint(bulk_configuration)
bulk_configuration.update()
nlsave(filename, bulk_configuration, object_id='bulk_configuration_GGAHalf')

# -------------------------------------------------------------
# Optical spectrum
# -------------------------------------------------------------
energies = numpy.linspace(0, 0.17358, 1001) * eV
optical_spectrum = OpticalSpectrum(
    configuration=bulk_configuration,
    kpoints=MonkhorstPackGrid(26, 26, 20),
    energies=energies,
    broadening=0.1 * eV,
    bands_below_fermi_level=1000,
    bands_above_fermi_level=1000,
)
nlsave(filename, optical_spectrum, object_id='optical_spectrum')

# -------------------------------------------------------------
# Born effective charges
# -------------------------------------------------------------
# Setup Born effective class object.
nk, nk_off = 7, 3
# Perform calculation.
born_effective_charge = BornEffectiveCharge(
   configuration=bulk_configuration,
   kpoints_a=MonkhorstPackGrid(nk, nk_off, nk_off),
   kpoints_b=MonkhorstPackGrid(nk_off, nk, nk_off),
   kpoints_c=MonkhorstPackGrid(nk_off, nk_off, nk),
   atomic_displacement=0.01 * Angstrom)
nlsave(filename, born_effective_charge, object_id='born_effective_charge')

# -------------------------------------------------------------
# Dielectric tensor
# -------------------------------------------------------------
dielectric_tensor = DielectricTensor(
        configuration=bulk_configuration,
        optical_spectrum=optical_spectrum,
        dynamical_matrix=dynamical_matrix,
        born_effective_charge=born_effective_charge,
        )
nlsave(filename, dielectric_tensor, object_id='dielectric_tensor')

dielectric_tensor_infrared_spectrum.py

Notice that the energies at interest are defined in the input to the OpticalSpectrum object.

It is important to note that this example is only intended to give an overview of the capabilities of DielectricTensor.