Molecule

Multi-level molecule (monomer)

The molecule is defined by the vector of energies of its states and by the transition dipole moments between allowed transitions.

>>> m = Molecule([0.0, 1.0])
>>> print(m.Nel)
2

Information about the molecule can be obtained simply by printing it

>>> print(m)
<BLANKLINE>
quantarhei.Molecule object
==========================
   name =
   position = None
   number of electronic states = 2
   # State properties
   State nr: 0 (ground state)
      electronic energy = 0.0 1/fs
      number of vibrational modes = 0
<BLANKLINE>
   State nr: 1
      electronic energy = 1.0 1/fs
      transition 0 -> 1
      transition dipole moment = [0.0, 0.0, 0.0]
      number of vibrational modes = 0
<BLANKLINE>
>>> import quantarhei as qr
>>> mol = qr.TestMolecule("two-levels-1-mode")
>>> print(mol)
<BLANKLINE>
quantarhei.Molecule object
==========================
   name = two-levels-1-mode
   position = None
   number of electronic states = 2
   # State properties
   State nr: 0 (ground state)
      electronic energy = 0.0 1/fs
      number of vibrational modes = 1
      # Mode properties
      mode no. = 0
         frequency = 1.0 1/fs
         shift = 0.0
         nmax = 2
   State nr: 1
      electronic energy = 1.0 1/fs
      transition 0 -> 1
      transition dipole moment = [0.0, 0.0, 0.0]
      number of vibrational modes = 1
      # Mode properties
      mode no. = 0
         frequency = 1.0 1/fs
         shift = 0.0
         nmax = 2

Class Details

class quantarhei.builders.molecules.Molecule(elenergies=[0.0, 1.0], name=None)[source]

Multi-level molecule (monomer)

Parameters:
  • name (str) – Monomer descriptor; a string identifying the monomer
  • elenergies (list of real numbers) – List of electronic energies, one per state. It includes ground state energy. It wise to chose the ground state energy as zero.
Attributes:
dmoments
elenergies
nal
nel
nmod
position

Methods

add_Mode(mod) Adds a vibrational mode to the monomer
build() Building routine for the molecule
check_temperature_consistent() Checks that the temperature is the same for all components
copy() Returns a shallow copy of the self
deepcopy() Returns a deep copy of the self
get_FoersterRateMatrix() Returns Förster rate matrix for the open system
get_Hamiltonian([multi, recalculate]) Returns the Hamiltonian of the Molecule object
get_KTHierarchy([depth]) Returns the Kubo-Tanimura hierarchy of an open system
get_KTHierarchyPropagator([depth]) Returns a propagator based on the Kubo-Tanimura hierarchy
get_Mode(N) Returns the Nth mode of the Molecule object
get_RedfieldRateMatrix() Returns Redfield rate matrix
get_ReducedDensityMatrixPropagator(timeaxis) Returns propagator of the density matrix
get_RelaxationTensor(timeaxis[, …]) Returns a relaxation tensor corresponding to the system
get_SystemBathInteraction() Returns a SystemBathInteraction object of the molecule
get_TransitionDipoleMoment([multi]) Returns the transition dipole moment operator
get_adiabatic_coupling(state1, state2) Returns adiabatic coupling between two states
get_diabatic_coupling(element) Returns list of coupling descriptors
get_dipole(N[, M]) Returns the dipole vector for a given electronic transition
get_egcf(transition)
get_electronic_natural_lifetime(N[, epsilon_r]) Returns natural lifetime of a given electronic state
get_electronic_natural_linewidth(N) Returns natural linewidth of a given electronic state
get_energy(N) Returns energy of the Nth state of the molecule
get_excited_density_matrix([condition, …]) Returns the density matrix corresponding to excitation condition
get_mode_environment(mode, elstate) Returns mode environment
get_number_of_modes() Retruns the number of modes in this molecule
get_potential_1D(mode, points[, other_modes]) Returns the one dimensional diabatic potentials
get_potential_2D(modes, points[, other_modes]) Returns the two dimensional diabatic potentials
get_temperature() Returns temperature of the molecule
get_thermal_ReducedDensityMatrix() Returns equilibrium density matrix for a give temparature
get_transition_dephasing(transition) Returns the dephasing rate of a given transition
get_transition_environment(transition) Returns energy gap correlation function of a monomer
get_transition_width(transition) Returns the transition width
liouville_pathways_1([eUt, ham, dtol, ptol, …]) Generator of the first order Liouville pathways
load(filename[, test]) Loads an object from a file and returns it
loaddir(dirname) Returns a directory of objects saved into a directory
plot_dressed_sticks([dfce, xlims, nsteps, …]) Plots a stick spectrum dessed by a supplied function
plot_potential_1D(mode, points[, …]) Plots the potentials
plot_stick_spectrum([xlims, ylims, …]) Plots the stick spectrum of the molecule
save(filename[, comment, test]) Saves the object with all its content into a file
savedir(dirname[, tag, comment, test]) Saves an object into directory containing a file with unique name
scopy() Creates a copy of the object by saving and loading it
set_adiabatic_coupling(state1, state2, coupl) Sets adiabatic coupling between two states
set_diabatic_coupling(element, factor[, shifts]) Sets off-diagobal elements of the diabatic potential matrix
set_dipole(N, M[, vec]) Sets transition dipole moment for an electronic transition
set_egcf(transition, egcf)
set_egcf_mapping(transition, …) Sets a correlation function mapping for a selected transition.
set_electronic_rwa(rwa_indices) Sets which electronic states belong to different blocks
set_energy(N, en) Sets the energy of the Nth state of the molecule
set_mode_environment([mode, elstate, corfunc]) Sets mode environment
set_transition_dephasing(transition, deph) Sets the dephasing rate of a given transition
set_transition_environment(transition, egcf) Sets a correlation function for a transition on this monomer
set_transition_width(transition, width) Sets the width of a given transition
unset_transition_environment(transition) Unsets correlation function from a transition on this monomer
convert_energy_2_current_u  
convert_energy_2_internal_u  
convert_length_2_current_u  
convert_length_2_internal_u  
get_diabatic_shifts  
set_electronic_natural_lifetime  
unit_repr  
unit_repr_latex  
add_Mode(mod)[source]

Adds a vibrational mode to the monomer

Parameters:mod (quantarhei.Mode) – Intramolecular vibrational mode

Examples

>>> mode = Mode()
>>> mol = Molecule([0.0, 2.0])
>>> print(mol.get_number_of_modes())
0
>>> mol.add_Mode(mode)
>>> print(mol.get_number_of_modes())
1
build()[source]

Building routine for the molecule

Unlike with the Aggregate, it is not compulsory to call build() before we start using the Molecule

check_temperature_consistent()[source]

Checks that the temperature is the same for all components

Examples

Isolated molecule has always consistent temperature

>>> m = Molecule([0.0, 1.0, 1.2])
>>> print(m.check_temperature_consistent())
True
>>> from ..qm.corfunctions import CorrelationFunction
>>> from .. import TimeAxis
>>> ta = TimeAxis(0.0,1000,1.0)
>>> params1 = dict(ftype="OverdampedBrownian",
...                reorg=20,cortime=100,T=300)
>>> params2 = dict(ftype="OverdampedBrownian",
...                reorg=20,cortime=100,T=200)
>>> with energy_units("1/cm"):
...     cf1 = CorrelationFunction(ta, params1)
...     cf2 = CorrelationFunction(ta, params2)
>>> m.set_transition_environment((0,1), cf1)
>>> m.set_transition_environment((0,2), cf2)
>>> print(m.get_temperature())
Traceback (most recent call last):
...
Exception: Molecular environment has an inconsisten temperature
get_Hamiltonian(multi=True, recalculate=False)[source]

Returns the Hamiltonian of the Molecule object

Examples

After the molecule is created, its Hamiltonian can be obtained

>>> import quantarhei as qr
>>> mol = qr.TestMolecule("two-levels-1-mode")
>>> H = mol.get_Hamiltonian()
>>> print(H.dim)
4
>>> print(H.data)
[[ 0.  0.  0.  0.]
 [ 0.  1.  0.  0.]
 [ 0.  0.  1.  0.]
 [ 0.  0.  0.  2.]]

The Hamiltonian is stored, and next time it is retrieved we obtain the same object.

>>> print(H.has_rwa)
False

We can manipulate with the Hamiltonian in between retrievals >>> H.set_rwa([0,1]) >>> print(H.has_rwa) True

The Hamiltonian the obtain the second time is affected by the changes performed outside the Molecule class

>>> H1 = mol.get_Hamiltonian()
>>> print(H1.has_rwa)
True

The newly obtained Hamiltonina IS the Hamiltonian obtained earlier.

>>> H1 is H
True
get_Mode(N)[source]

Returns the Nth mode of the Molecule object

Parameters:N (int) – Index of the mode to be returned

Examples

>>> import quantarhei as qr
>>> mol = qr.TestMolecule("two-levels-1-mode")
>>> mod = mol.get_Mode(1)
Traceback (most recent call last):
    ...
IndexError: list index out of range
>>> mod = mol.get_Mode(0)
>>> mod.get_energy(1)
1.0
get_SystemBathInteraction()[source]

Returns a SystemBathInteraction object of the molecule

get_TransitionDipoleMoment(multi=True)[source]

Returns the transition dipole moment operator

get_adiabatic_coupling(state1, state2)[source]

Returns adiabatic coupling between two states

get_diabatic_coupling(element)[source]

Returns list of coupling descriptors

get_dipole(N, M=None)[source]

Returns the dipole vector for a given electronic transition

There are two ways how to use this function:

  1. recommended

    get_dipole((0,1),[1.0, 0.0, 0.0])

    here N represents a transition by a tuple, M is the dipole

  2. deprecated (used in earlier versions of quantarhei)

    get_dipole(0,1,[1.0, 0.0, 0.0])

    here the transition is characterized by two integers and the last arguments is the vector

Examples

>>> m = Molecule([0.0, 1.0])
>>> m.set_dipole((0,1),[1.0, 0.0, 0.0])
>>> m.get_dipole((0,1))
array([ 1.,  0.,  0.])
get_electronic_natural_lifetime(N, epsilon_r=1.0)[source]

Returns natural lifetime of a given electronic state

get_electronic_natural_linewidth(N)[source]

Returns natural linewidth of a given electronic state

get_energy(N)[source]

Returns energy of the Nth state of the molecule

Parameters:N (int) – Index of the state

Examples

>>> import quantarhei as qr
>>> mol = qr.TestMolecule("two-levels-1-mode")
>>> mol.get_energy(1)
1.0

This methods reacts to the energy_units context manager

>>> with qr.energy_units("1/cm"):
...     print("{0:.8f}".format(mol.get_energy(1)))
5308.83745888
get_mode_environment(mode, elstate)[source]

Returns mode environment

Returns the environment (bath correlation function) interacting with a a given mode in a given electronic state.

Parameters:
  • mode (int) – index of the mode
  • elstate (int) – index of the electronic state
get_number_of_modes()[source]

Retruns the number of modes in this molecule

Examples

>>> m = Molecule([0.0, 1.0])
>>> m.get_number_of_modes()
0
>>> import quantarhei as qr
>>> m = qr.TestMolecule("two-levels-1-mode")
>>> m.get_number_of_modes()
1
get_potential_1D(mode, points, other_modes=None)[source]

Returns the one dimensional diabatic potentials

get_potential_2D(modes, points, other_modes=None)[source]

Returns the two dimensional diabatic potentials

get_temperature()[source]

Returns temperature of the molecule

Checks if the setting of environments is consistent and than takes the temperature from one of the energy gaps. If no environment (correlation function) is assigned to this molecule, we assume zero temperature.

Examples

Default temperature is 0 K

>>> import quantarhei as qr
>>> mol = qr.TestMolecule("two-levels-1-mode")
>>> mol.get_temperature()
0.0

Molecule gets its temperature from the environment

>>> from ..qm.corfunctions import CorrelationFunction
>>> from .. import TimeAxis
>>> ta = TimeAxis(0.0,1000,1.0)
>>> params = dict(ftype="OverdampedBrownian",reorg=20,cortime=100,T=300)
>>> with energy_units("1/cm"):
...     cf = CorrelationFunction(ta, params)
>>> m = Molecule([0.0, 1.0])
>>> m.set_transition_environment((0,1), cf)
>>> print(m.get_temperature())
300
get_transition_dephasing(transition)[source]

Returns the dephasing rate of a given transition

Parameters:transition ({tuple, list}) – Quantum numbers of the states between which the transition occurs
get_transition_environment(transition)[source]

Returns energy gap correlation function of a monomer

Parameters:transition (tuple) – A tuple describing a transition in the molecule, e.g. (0,1) is a transition from the ground state to the first excited state.

Example

>>> m = Molecule([0.0, 1.0])

Environment of the transition has to be set first

>>> cc = m.get_transition_environment((0,1))
Traceback (most recent call last):
    ...
Exception: No environment set for the transition

Environment is characterized by the bath correlation function

>>> from ..qm.corfunctions import CorrelationFunction
>>> from .. import TimeAxis
>>> ta = TimeAxis(0.0,1000,1.0)
>>> params = dict(ftype="OverdampedBrownian",reorg=20,cortime=100,T=300)
>>> with energy_units("1/cm"):
...     cf = CorrelationFunction(ta, params)
>>> m.set_transition_environment((0,1), cf)
>>> cc = m.get_transition_environment((0,1))
>>> cc == cf
True

When non-existent transition is tried, exception is raised

>>> cc = m.get_transition_environment((0,2))
Traceback (most recent call last):
    ...
Exception: Index out of range
get_transition_width(transition)[source]

Returns the transition width

Returns the full width at half maximum (FWHM) of a Gaussian lineshape, or half width at half maximum (HWHM) of a Lorentzian lineshape

Parameters:transition ({tuple, list}) – Quantum numbers of the states between which the transition occurs
liouville_pathways_1(eUt=None, ham=None, dtol=0.01, ptol=0.001, etol=1e-06, verbose=0, lab=None)[source]

Generator of the first order Liouville pathways

Generator of the pathways for an absorption spectrum calculation.

Parameters:eUt (EvolutionSuperOperator) – Evolution superoperator representing the evolution of optical coherence in the system
dtol : float
Minimum acceptable strength of the transition from ground to excited state, relative to the maximum dipole strength available in the system
ptol : float
Minimum acceptable population of the ground state (e.g. states not thermally populated are excluded)
lab : LaboratorySetup
Object representing laboratory setup - number of pulses, polarization etc.
Returns:lst – List of LiouvillePathway objects
Return type:list
plot_dressed_sticks(dfce=None, xlims=[0, 1], nsteps=1000, show_zero_coupling=False, show=True)[source]

Plots a stick spectrum dessed by a supplied function

plot_potential_1D(mode, points, other_modes=None, nonint=True, states=None, energies=True, show=True, ylims=None)[source]

Plots the potentials

plot_stick_spectrum(xlims=[0.0, 1.0], ylims=None, show_zero_coupling=False, show=True)[source]

Plots the stick spectrum of the molecule

set_adiabatic_coupling(state1, state2, coupl)[source]

Sets adiabatic coupling between two states

set_diabatic_coupling(element, factor, shifts=None)[source]

Sets off-diagobal elements of the diabatic potential matrix

set_dipole(N, M, vec=None)[source]

Sets transition dipole moment for an electronic transition

There are two ways how to use this function:

  1. recommended

    set_dipole((0,1),[1.0, 0.0, 0.0])

    here N represents a transition by a tuple, M is the dipole

  2. deprecated (used in earlier versions of quantarhei)

    set_dipole(0,1,[1.0, 0.0, 0.0])

    here the transition is characterized by two integers and the last arguments is the vector

Examples

>>> m = Molecule([0.0, 1.0])
>>> m.set_dipole((0,1),[1.0, 0.0, 0.0])
>>> m.get_dipole((0,1))
array([ 1.,  0.,  0.])
set_egcf_mapping(transition, correlation_matrix, position)[source]

Sets a correlation function mapping for a selected transition.

The monomer can either have a correlation function assigned to it, or it can be a part of a correlation matrix. Here the mapping to the correlation matrix is specified.

Parameters:
  • transition (tuple) – A tuple describing a transition in the molecule, e.g. (0,1) is a transition from the ground state to the first excited state.
  • correlation_matrix (CorrelationFunctionMatrix) – An instance of CorrelationFunctionMatrix
  • position (int) – Position in the CorrelationFunctionMatrix corresponding to the monomer.

Examples

A set of three monomers

>>> en1 = [0.0,12100, 13000] #*cm2int]
>>> en2 = [0.0,12000] #*cm2int]
>>> with energy_units("1/cm"):
...     m1 = Molecule(en1)
...     m2 = Molecule(en1)
...     m3 = Molecule(en2)

Bath correlation functions to describe molecular environment

>>> from .. import TimeAxis
>>> from .. import CorrelationFunction
>>> time = TimeAxis(0.0, 2000, 1.0) # in fs
>>> temperature = 300.0 # in Kelvins
>>> cfce_params1 = dict(ftype="OverdampedBrownian",
...           reorg=30.0,
...           cortime=60.0,
...           T=temperature)
>>> cfce_params2 = dict(ftype="OverdampedBrownian",
...           reorg=30.0,
...           cortime=60.0,
...           T=temperature)
>>> with energy_units("1/cm"):
...    cf1 = CorrelationFunction(time,cfce_params1)
...    cf2 = CorrelationFunction(time,cfce_params2)

Environment of the molecules is collected to a matrix. A smaller number of correlation functions can be assigned to a large number of “sites”.

>>> cm = CorrelationFunctionMatrix(time,3)
>>> ic1 = cm.set_correlation_function(cf1,[(0,0),(2,2)])
>>> ic2 = cm.set_correlation_function(cf2,[(1,1)])

The sites of the in the matrix are assigned to molecules.

>>> m1.set_egcf_mapping((0,1), cm, 0)
>>> m2.set_egcf_mapping((0,1), cm, 1)
>>> m3.set_egcf_mapping((0,1), cm, 2)

The environment cannot be set twice

>>> m1.set_egcf_mapping((0,1), cm, 2)
Traceback (most recent call last):
    ...
Exception: Monomer has a correlation function already
>>> m1.set_egcf_mapping((0,2), cm, 2)
>>> m1.set_egcf_mapping((0,2), cm, 2)
Traceback (most recent call last):
    ...
Exception: Monomer has a correlation function already
set_electronic_rwa(rwa_indices)[source]

Sets which electronic states belong to different blocks

Setting the indices of blocks should allow the construction of Rotating Wave Approximation for the Hamiltonian. This in turn is used for the calculations of optical spectra.

Examples

>>> with energy_units("1/cm"):
...     mol1 = Molecule([0.0, 10000.0])
>>> mol1.set_electronic_rwa([0, 1])
>>> H = mol1.get_Hamiltonian()
>>> H.has_rwa
True
>>> with energy_units("1/cm"):
...     mol2 = Molecule([0.0, 10000.0, 10500.0, 20000.0])
>>> with energy_units("1/cm"):
...     mod2 = Mode(frequency=200.0)
...     mol2.add_Mode(mod2)
>>> mod2.set_nmax(0, 5)
>>> mod2.set_nmax(1, 4)
>>> mod2.set_nmax(2, 5)
>>> mod2.set_nmax(3, 3)
>>> mol2.set_electronic_rwa([0, 1, 3])
>>> H = mol2.get_Hamiltonian()
>>> H.has_rwa
True
>>> print(H.rwa_indices)
[ 0  5 14]
set_energy(N, en)[source]

Sets the energy of the Nth state of the molecule

Parameters:
  • N (int) – Index of the state
  • en (float) – Energy to be assigned to the Nth state.

Examples

>>> import quantarhei as qr
>>> mol = qr.TestMolecule("two-levels-1-mode")
>>> mol.set_energy(1, 1.5)
>>> mol.get_energy(1)
1.5

This method reacts to the energy_units context manager

>>> import quantarhei as qr
>>> mol = qr.TestMolecule("two-levels-1-mode")
>>> with qr.energy_units("1/cm"):
...     mol.set_energy(1, 5308.8374588761453)
>>> mol.get_energy(1)
1.0
set_mode_environment(mode=0, elstate=0, corfunc=None)[source]

Sets mode environment

Sets the environment (bath correlation function) interacting with a a given mode in a given electronic state.

Parameters:
  • mode (int) – index of the mode
  • elstate (int) – index of the electronic state
  • corfunc (quantarhei.CorrelationFunction) – CorrelationFunction object
set_transition_dephasing(transition, deph)[source]

Sets the dephasing rate of a given transition

Parameters:
  • transition ({tuple, list}) – Quantum numbers of the states between which the transition occurs
  • deph (float) – Dephasing rate of the transition
set_transition_environment(transition, egcf)[source]

Sets a correlation function for a transition on this monomer

Parameters:
  • transition (tuple) – A tuple describing a transition in the molecule, e.g. (0,1) is a transition from the ground state to the first excited state.
  • egcf (CorrelationFunction) – CorrelationFunction object

Example

>>> from ..qm.corfunctions import CorrelationFunction
>>> from .. import TimeAxis
>>> ta = TimeAxis(0.0,1000,1.0)
>>> params = dict(ftype="OverdampedBrownian",reorg=20,cortime=100,T=300)
>>> with energy_units("1/cm"):
...     cf = CorrelationFunction(ta, params)
>>> m = Molecule([0.0, 1.0])
>>> m.set_transition_environment((0,1), cf)
>>> print(m._has_system_bath_coupling)
True

When the environment is already set, the next attempt is refused

>>> m.set_transition_environment((0,1), cf)
Traceback (most recent call last):
...
Exception: Correlation function already speficied for this monomer

The environment cannot be set when the molecule is mapped on a correlation function matrix

>>> with energy_units("1/cm"):
...     cf1 = CorrelationFunction(ta, params)
>>> cm = CorrelationFunctionMatrix(ta,1)
>>> ic1 = cm.set_correlation_function(cf1,[(0,0)])
>>> m1 = Molecule([0.0, 1.0])
>>> m1.set_egcf_mapping((0,1), cm, 0)
>>> m1.set_transition_environment((0,1), cf1)
Traceback (most recent call last):
    ...
Exception: This monomer is mapped on a CorrelationFunctionMatrix
set_transition_width(transition, width)[source]

Sets the width of a given transition

Parameters:
  • transition ({tuple, list}) – Quantum numbers of the states between which the transition occurs
  • width (float) – The full width at half maximum (FWHM) of a Gaussian lineshape, or half width at half maximum (HWHM) of a Lorentzian lineshape
unset_transition_environment(transition)[source]

Unsets correlation function from a transition on this monomer

This is needed if the environment is to be replaced

Parameters:transition (tuple) – A tuple describing a transition in the molecule, e.g. (0,1) is a transition from the ground state to the first excited state.

Example

>>> from ..qm.corfunctions import CorrelationFunction
>>> from .. import TimeAxis
>>> ta = TimeAxis(0.0,1000,1.0)
>>> params = dict(ftype="OverdampedBrownian",reorg=20,cortime=100,T=300)
>>> with energy_units("1/cm"):
...     cf = CorrelationFunction(ta, params)
>>> m = Molecule([0.0, 1.0])
>>> m.set_transition_environment((0,1), cf)
>>> print(m._has_system_bath_coupling)
True
>>> m.unset_transition_environment((0,1))
>>> print(m._has_system_bath_coupling)
False

When the environment is unset, the next attempt to set is succesful

>>> m.set_transition_environment((0,1), cf)
>>> print(m._has_system_bath_coupling)
True

The environment cannot be unset when the molecule is mapped on a correlation function matrix

>>> with energy_units("1/cm"):
...     cf1 = CorrelationFunction(ta, params)
>>> cm = CorrelationFunctionMatrix(ta,1)
>>> ic1 = cm.set_correlation_function(cf1,[(0,0)])
>>> m1 = Molecule([0.0, 1.0])
>>> m1.set_egcf_mapping((0,1), cm, 0)
>>> m1.unset_transition_environment((0,1))
Traceback (most recent call last):
    ...
Exception: This monomer is mapped on a CorrelationFunctionMatrix