Density matrices

class rhodent.density_matrices.ConvolutionDensityMatricesFromDisk(ksd, pulserho_fmt, pulses, times, derivative_order_s=[0], real=True, imag=True, calc_size=1)[source]

Collection of density matrices in the Kohn-Sham basis after convolution with one or several pulses, for different times. Read from disk

Parameters:
  • ksd (KohnShamDecomposition | str) – KohnShamDecomposition object or filename

  • pulserho_fmt (str) –

    The pulserho_fmt is a formatting string for the density matrices saved to disk. Example:

    pulserho_fmt = ‘pulserho/t{time:09.1f}{tag}.npy’

  • pulses (Collection[GaussianPulse]) – Convolute the density matrices with these pulses

  • times (list[float] | NDArray[np.float64]) – Compute density matrices for these times (or as close to them as possible). In as

  • derivative_order_s (list[int]) – Compute density matrix derivatives of the following orders. 0 for plain density matrix and positive integers for derivatives

  • real (bool) – Calculate real part of density matrices

  • imag (bool) – Calculate imaginary part of density matrices

  • calc_size (int) – Size of the calculation communicator

class rhodent.density_matrices.ConvolutionDensityMatricesFromFrequency(ksd, frho_fmt, pulses, times, derivative_order_s=[0], real=True, imag=True, calc_size=1)[source]

Collection of density matrices in the Kohn-Sham basis after convolution with one or several pulses, for different times. Obtained from the the density matrices in frequency space, which are read from disk.

Parameters:
  • ksd (KohnShamDecomposition | str) – KohnShamDecomposition object or filename

  • frho_fmt (str) –

    The pulserho_fmt is a formatting string for the density matrices in frequency space saved to disk. Example:

    frho_fmt = ‘frho/w{freq:05.2f}-{reim}.npy’

    Should accept variables {freq} and {reim}

  • pulses (Collection[GaussianPulse]) – Convolute the density matrices with these pulses

  • times (list[float] | NDArray[np.float64]) – Compute density matrices for these times (or as close to them as possible). In as

  • derivative_order_s (list[int]) – Compute density matrix derivatives of the following orders. 0 for plain density matrix and positive integers for derivatives

  • real (bool) – Calculate real part of density matrices

  • imag (bool) – Calculate imaginary part of density matrices

  • calc_size (int) – Size of the calculation communicator

class rhodent.density_matrices.ConvolutionDensityMatricesFromWaveFunctions(ksd, wfs_fname, perturbation, pulses, times, derivative_order_s=[0], real=True, imag=True, calc_size=1, stride_opts=None, stridet=1)[source]

Collection of density matrices in the Kohn-Sham basis after convolution with one or several pulses, for different times. Obtained from the wave functions dump file, which is read from disk.

Parameters:
  • ksd (KohnShamDecomposition | str) – KohnShamDecomposition object or filename

  • wfs_fname (str) – Filename of the GPAW wave functions dump file

  • perturbation (Perturbation | dict) – Perturbation that was present during time propagation

  • pulses (Collection[GaussianPulse]) – Convolute the density matrices with these pulses

  • times (list[float] | NDArray[np.float64]) – Compute density matrices for these times (or as close to them as possible). In as

  • derivative_order_s (list[int]) – Compute density matrix derivatives of the following orders. 0 for plain density matrix and positive integers for derivatives

  • real (bool) – Calculate real part of density matrices

  • imag (bool) – Calculate imaginary part of density matrices

  • calc_size (int) – Size of the calculation communicator

  • stride_opts – Options passed to the GPAW wave functions reader

  • stridet (int) – Skip this many steps when reading wave function dump file

property myt: list[int]

List of indices corresponding to the time indices on held on this rank

parallel_prepare()[source]

Read everything necessary synchronously on all ranks

work_loop(rank)[source]

The work to be done by a certain rank of the loop communicator.

Parameters:

rank (int) – Rank of the loop communicator

Yields:
  • Objects representing the time, frequency or pulse to be computed by rank rank.

  • None is yielded when rank does not do any work while other ranks are doing work.

Return type:

Generator[ConvolutionDensityMatrixMetadata | None, None, None]

class rhodent.density_matrices.DensityMatrix(ksd, matrices, data_is_ravelled=True)[source]

Wrapper for the density matrix in the Kohn-Sham basis at one moment in time or at one frequency.

The plain density matrix and/or derivatives thereof may be stored.

Parameters:
  • ksd (gpaw.lcaotddft.ksddecomposition.KohnShamDecomposition) – KohnShamDecomposition object

  • data_is_ravelled (bool) – Whether the data is stored in a ravelled form (single index \(p\) for electron-hole pairs) or not (two indices \(ia\) for hole and electron).

  • matrices (dict[int, ndarray[Any, dtype[complex64]] | None]) – Dictionary mapping derivative orders (0, 1, 2) for zeroth, first, second derivative, .. to file names where density matrices are stored.

property P_ia: rhodent.utils.DistributedArray

\(P\) in non-ravelled form (indices \(ia\) for electron-hole pairs)

\[P_{ia} = \frac{2 \mathrm{Im}\:\delta\rho_{ia}}{\sqrt{2 f_{ia}}}\]

where \(f_{ia}\) is the occupation number difference of pair \(ia\).

property P_p: rhodent.utils.DistributedArray

\(P\) in ravelled form (common index \(p\) for electron-hole pairs)

\[P_p = \frac{2 \mathrm{Im}\:\delta\rho_p}{\sqrt{2 f_p}}\]

where \(f_p\) is the occupation number difference of pair \(p\).

property Q_ia: rhodent.utils.DistributedArray

\(Q\) in non-ravelled form (indices \(ia\) for electron-hole pairs)

\[Q_{ia} = \frac{2 \mathrm{Re}\:\delta\rho_{ia}}{\sqrt{2 f_{ia}}}\]

where \(f_{ia}\) is the occupation number difference of pair \(ia\).

property Q_p: rhodent.utils.DistributedArray

\(Q\) in ravelled form (common index \(p\) for electron-hole pairs)

\[Q_p = \frac{2 \mathrm{Re}\:\delta\rho_p}{\sqrt{2 f_p}}\]

where \(f_p\) is the occupation number difference of pair \(p\).

classmethod broadcast(density_matrix, ksd, root, comm)[source]

Broadcast a density matrix object which is on one rank to all other ranks

Parameters:
  • density_matrix (DensityMatrix | None) – The density matrix to be broadcast on the root rank, and None on other ranks

  • ksd (gpaw.lcaotddft.ksddecomposition.KohnShamDecomposition) – KohnShamDecomposition object

  • root (int) – Rank of the process that has the original data

  • comm – MPI communicator

Return type:

DensityMatrix

copy()[source]

Copy the density matrix

Return type:

DensityMatrix

property dP_ia: rhodent.utils.DistributedArray

First time derivative of \(P\) in non-ravelled form

property dP_p: rhodent.utils.DistributedArray

First time derivative of \(P\) in ravelled form

property dQ_ia: rhodent.utils.DistributedArray

First time derivative of \(Q\) in non-ravelled form

property dQ_p: rhodent.utils.DistributedArray

First time derivative of \(Q\) in ravelled form

property ddP_ia: rhodent.utils.DistributedArray

Second time derivative of \(P\) in non-ravelled form

property ddP_p: rhodent.utils.DistributedArray

Second time derivative of \(P\) in ravelled form

property ddQ_ia: rhodent.utils.DistributedArray

Second time derivative of \(Q\) in non-ravelled form

property ddQ_p: rhodent.utils.DistributedArray

Second time derivative of \(Q\) in ravelled form

property ddrho_ia: rhodent.utils.DistributedArray

Second time derivative of \(\delta rho_{ia}\)

In non-ravelled form (indices \(ia\) for electron-hole pairs).

property ddrho_p: rhodent.utils.DistributedArray

Second time derivative of \(\delta \rho_p\)

In ravelled form (common index \(p\) for electron-hole pairs).

property drho_ia: rhodent.utils.DistributedArray

First time derivative of \(\delta rho_{ia}\)

In non-ravelled form (indices \(ia\) for electron-hole pairs).

property drho_p: rhodent.utils.DistributedArray

First time derivative of \(\delta \rho_p\)

In ravelled form (common index \(p\) for electron-hole pairs).

property f_ia: rhodent.utils.DistributedArray

Occupation number difference \(f_{ia}\)

In non-ravelled form (indices \(ia\) for electron-hole pairs).

property f_p: rhodent.utils.DistributedArray

Occupation number difference \(f_p\)

In ravelled form (common index \(p\) for electron-hole pairs).

property ksd: gpaw.lcaotddft.ksddecomposition.KohnShamDecomposition

Kohn-Sham decomposition object

property rank: int

MPI rank of the ksd object

property rho_ia: rhodent.utils.DistributedArray

Electron-hole part of induced density matrix \(\delta rho_{ia}\)

In non-ravelled form (indices \(ia\) for electron-hole pairs).

property rho_p: rhodent.utils.DistributedArray

Electron-hole part of induced density matrix \(\delta \rho_p\)

In ravelled form (common index \(p\) for electron-hole pairs).

class rhodent.density_matrices.DensityMatrixBuffer(nnshape, xshape, dtype, re_buffers={}, im_buffers={})[source]

Buffer for the density matrix.

Objects of this class can hold buffers for real and imaginary parts and various derivatives at the same time.

Each buffer has two dimensions corresponding to (part of) the density matrix, and optionally additional dimensions for e.g. time.

Parameters:
  • nnshape (tuple[int, int]) – Shape of the dimension corresponding to the density matrix. Must have dimension 2

  • xshape (tuple[int, ...]) – Shape of the additional dimension corresponding to, e.g., time

  • dtype (Union[dtype[TypeVar(DTypeT, bound= generic, covariant=True)], type[TypeVar(DTypeT, bound= generic, covariant=True)], _SupportsDType[dtype[TypeVar(DTypeT, bound= generic, covariant=True)]]]) – Data type of the density matrices

  • re_buffers (dict[int, ndarray[Any, dtype[TypeVar(DTypeT, bound= generic, covariant=True)]]]) – Buffers for the different derivatives of the real part of the density matrix as dictionaries, where the keys is the derivative order (0, 1, 2, etc.) and the value is a numpy array of shape (nnshape, xshape)

  • im_buffers (dict[int, ndarray[Any, dtype[TypeVar(DTypeT, bound= generic, covariant=True)]]]) – Same as re_buffers but for imaginary part

broadcast_x(data_nnx)[source]

Broadcast the x dimensions of data_nnx

Return type:

ndarray[Any, dtype[TypeVar(DTypeT, bound= generic, covariant=True)]]

copy()[source]

Return a deep copy of this object (buffers are copied)

Return type:

DensityMatrixBuffer

property derivatives_imag: list[int]

Stored derivative order of real density matrices in sorted order

property derivatives_real: list[int]

Stored derivative order of real density matrices in sorted order

property dtype: dtype[DTypeT]

Dtype of the buffers

ensure_contiguous_buffers()[source]

Make the buffers contiguous if they are not already

property imag: ndarray[Any, dtype[DTypeT]]

Buffer of shape nnshape + xshape corresponding to imaginary data

property imag1: ndarray[Any, dtype[DTypeT]]

Buffer of shape nnshape + xshape corresponding to imaginary part of first derivative

property imag2: ndarray[Any, dtype[DTypeT]]

Buffer of shape nnshape + xshape corresponding to imaginary part of second derivative

property narrays: int

Number of arrays stored in this object

property nnshape: tuple[int, int]

Shape of the part of the density matrix

property real: ndarray[Any, dtype[DTypeT]]

Buffer of shape nnshape + xshape corresponding to real data

property real1: ndarray[Any, dtype[DTypeT]]

Buffer of shape nnshape + xshape corresponding to real part of first derivative

property real2: ndarray[Any, dtype[DTypeT]]

Buffer of shape nnshape + xshape corresponding to real part of second derivative

recv_arrays(comm, rank, log=None)[source]

Receive data to another MPI rank

Parameters:
  • comm – Communicator

  • rank (int) – Send to this rank

  • log (Optional[Logger]) – Optional logger

redistribute(target, comm, source_indices_r, target_indices_r, log=None)[source]

Redistribute this DensityMatrixBuffer to another according to user specified options

The nn dimensions of the self and target buffers should be the same, but the x dimensions can be different. However, self need to have the same shape on all ranks and target needs to have the same shape on all ranks.

Parameters:
  • target (DensityMatrixBuffer) – Target DensityMatrixBuffer

  • comm – MPI communicator

  • source_indices_r (Sequence[Union[tuple[Union[ndarray[Any, dtype[int64]], list[int], int, slice], ...], ndarray[Any, dtype[int64]], list[int], int, slice, None]]) – List of x-indices. The length of the list must equal to the communicator size. The x-index that is element r of the list corresponds to the data from the source that will be sent to rank r. If the x-index is None, then the rank corresponding to that element will not be receiving data. This argument must be the same on all ranks

  • recv_indices_r – List of x-indices. The length of the list must equal to the communicator size. The x-index that is element r of the list corresponds to the data in the target that will be received from rank r. If the x-index is None, then the rank corresponding to that element will not be sending data. This argument must be the same on all ranks

  • log (Optional[Logger]) – Optional logger

safe_add(real, derivative, data_nnx)[source]

Add data_nnx to the corrsponding buffer, if the dimensions of data_nnx are equal to or smaller than the buffer

If the dimensions of data_nnx are smaller than or equal to the dimensions of the buffer, add to the first elements of the buffer. Otherwise raise and error.

Parameters:
  • real (bool) – True if real, False if imaginary

  • derivative (int) – Derivative order

  • buffer_nnx – Data of shape (nnshape, xshape)

safe_fill(real, derivative, data_nnx)[source]

Write data_nnx to the corrsponding buffer, if the dimensions of data_nnx are equal to or smaller than the buffer

If the dimensions of data_nnx are smaller than or equal to the dimensions of the buffer, write to the first elements of the buffer. Otherwise raise and error.

Parameters:
  • real (bool) – True if real, False if imaginary

  • derivative (int) – Derivative order

  • buffer_nnx – Data of shape (nnshape, xshape)

send_arrays(comm, rank, log=None)[source]

Send data to another MPI rank

Parameters:
  • comm – Communicator

  • rank (int) – Send to this rank

  • log (Optional[Logger]) – Optional logger

property shape: tuple[int, ...]

Shape of the buffers

store(real, derivative, buffer_nnx)[source]

Store buffer for part of density matrix

Parameters:
  • real (bool) – True if real, False if imaginary

  • derivative (int) – Derivative order

  • buffer_nnx (ndarray[Any, dtype[TypeVar(DTypeT, bound= generic, covariant=True)]]) – Buffer of shape (nnshape, xshape)

property xshape: tuple[int, ...]

Shape of the additional dimension of the buffers

zero_buffers(real, imag, derivative_order_s)[source]

Initialize many buffers at once to and write zeros

Parameters:
  • real (bool) – Initialize buffers for real parts

  • imag (bool) – Initialize buffers for imaginary parts

  • derivative_order_s (list[int]) – Initialize buffers for these derivative orders

zeros(real, derivative)[source]

Initialize buffer with zeros for part of density matrix

Parameters:
  • real (bool) – True if real, False if imaginary

  • derivative (int) – Derivative order

class rhodent.density_matrices.FrequencyDensityMatricesFromDisk(ksd, frho_fmt, frequencies, real=True, imag=True, derivative_order_s=[0], calc_size=1, kickstr=1e-05)[source]

Collection of density matrices in the Kohn-Sham basis in the frequency domain, for different frequencies. Read from disk.

Plain density matrices and/or derivatives thereof may be represented.

Parameters:
  • ksd (KohnShamDecomposition | str) – KohnShamDecomposition object or filename

  • frho_fmt (str) –

    The pulserho_fmt is a formatting string for the density matrices in frequency space saved to disk. Example:

    frho_fmt = ‘frho/w{freq:05.2f}-{reim}.npy’

    Should accept variables {freq} and {reim}

  • frequencies (list[float] | NDArray[np.float64]) – Compute density matrices for these frequencies (or as close to them as possible). In eV

  • derivative_order_s (list[int]) – Compute density matrix derivatives of the following orders. 0 for plain density matrix and positive integers for derivatives

  • real (bool) – Calculate the Fourier transform of the real part of the density matrix

  • imag (bool) – Calculate the Fourier transform of the imaginary part of the density matrix

  • calc_size (int) – Size of the calculation communicator

  • kickstr (float) – Strength of the delta kick used during time propagation

class rhodent.density_matrices.FrequencyDensityMatricesFromWaveFunctions(ksd, wfs_fname, frequencies, real=True, imag=True, derivative_order_s=[0], calc_size=1, stride_opts=None, stridet=1)[source]

Collection of density matrices in the Kohn-Sham basis in the frequency domain, for different frequencies. Obtained from the wave functions dump file, which is read from disk.

Plain density matrices and/or derivatives thereof may be represented.

Parameters:
  • ksd (KohnShamDecomposition | str) – KohnShamDecomposition object or filename

  • wfs_fname (str) – Filename of the GPAW wave functions dump file

  • frequencies (list[float] | NDArray[np.float64]) – Compute density matrices for these frequencies (or as close to them as possible). In eV

  • derivative_order_s (list[int]) – Compute density matrix derivatives of the following orders. 0 for plain density matrix and positive integers for derivatives

  • real (bool) – Calculate the real part of density matrices

  • imag (bool) – Calculate the imaginary part of density matrices

  • calc_size (int) – Size of the calculation communicator

  • stride_opts – Options passed to the GPAW wave functions reader

  • stridet (int) – Skip this many steps when reading wave function dump file

property myw: list[int]

List of indices corresponding to the frequency indices on held on this rank

parallel_prepare()[source]

Read everything necessary synchronously on all ranks

class rhodent.density_matrices.TimeDensityMatricesFromWaveFunctions(ksd, wfs_fname, times, real=True, imag=True, calc_size=1, stride_opts=None, stridet=1)[source]

Collection of density matrices in the Kohn-Sham basis for different times, obtained by reading the wave functions dump file.

Parameters:
  • ksd (KohnShamDecomposition | str) – KohnShamDecomposition object or filename

  • wfs_fname (str) – Filename of the GPAW wave functions dump file

  • times (list[float] | NDArray[np.float64]) – Compute density matrices for these times (or as close to them as possible). In as

  • real (bool) – Calculate the real part of density matrices

  • imag (bool) – Calculate the imaginary part of density matrices

  • calc_size (int) – Size of the calculation communicator

  • stride_opts – Options passed to the GPAW wave functions reader

  • stridet (int) – Skip this many steps when reading wave function dump file

parallel_prepare()[source]

Read everything necessary synchronously on all ranks