Density matrices
Buffers
These classes store one or many density matrices.
- 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
[tuple
[int
,...
],dtype
[complex128
]] |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.typing.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.typing.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.typing.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.typing.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, andNone
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:
- property dP_ia: rhodent.typing.DistributedArray
First time derivative of \(P\) in non-ravelled form.
- property dP_p: rhodent.typing.DistributedArray
First time derivative of \(P\) in ravelled form.
- property dQ_ia: rhodent.typing.DistributedArray
First time derivative of \(Q\) in non-ravelled form.
- property dQ_p: rhodent.typing.DistributedArray
First time derivative of \(Q\) in ravelled form.
- property ddP_ia: rhodent.typing.DistributedArray
Second time derivative of \(P\) in non-ravelled form.
- property ddP_p: rhodent.typing.DistributedArray
Second time derivative of \(P\) in ravelled form.
- property ddQ_ia: rhodent.typing.DistributedArray
Second time derivative of \(Q\) in non-ravelled form.
- property ddQ_p: rhodent.typing.DistributedArray
Second time derivative of \(Q\) in ravelled form.
- property ddrho_ia: rhodent.typing.DistributedArray
Second time derivative of \(\delta rho_{ia}\).
In non-ravelled form (indices \(ia\) for electron-hole pairs).
- property ddrho_p: rhodent.typing.DistributedArray
Second time derivative of \(\delta \rho_p\).
In ravelled form (common index \(p\) for electron-hole pairs).
- property drho_ia: rhodent.typing.DistributedArray
First time derivative of \(\delta rho_{ia}\).
In non-ravelled form (indices \(ia\) for electron-hole pairs).
- property drho_p: rhodent.typing.DistributedArray
First time derivative of \(\delta \rho_p\).
In ravelled form (common index \(p\) for electron-hole pairs).
- property f_ia: rhodent.typing.DistributedArray
Occupation number difference \(f_{ia}\).
In non-ravelled form (indices \(ia\) for electron-hole pairs).
- property f_p: rhodent.typing.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.typing.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.typing.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
[tuple
[int
,...
],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
[tuple
[int
,...
],dtype
[TypeVar
(DTypeT
, bound=generic
, covariant=True)]]]) – Same asre_buffers
but for imaginary part.
- 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 imag: ndarray[tuple[int, ...], dtype[DTypeT]]
Buffer of shape nnshape + xshape corresponding to imaginary data.
- property imag1: ndarray[tuple[int, ...], dtype[DTypeT]]
Buffer of shape nnshape + xshape corresponding to imaginary part of first derivative.
- property imag2: ndarray[tuple[int, ...], 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[tuple[int, ...], dtype[DTypeT]]
Buffer of shape nnshape + xshape corresponding to real data.
- property real1: ndarray[tuple[int, ...], dtype[DTypeT]]
Buffer of shape nnshape + xshape corresponding to real part of first derivative.
- property real2: ndarray[tuple[int, ...], 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
) – TargetDensityMatrixBuffer
.comm – MPI communicator.
source_indices_r (
Sequence
[Union
[tuple
[Union
[ndarray
[tuple
[int
,...
],dtype
[int64
]],list
[int
],int
,slice
],...
],ndarray
[tuple
[int
,...
],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 ranksrecv_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.
- 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.
- class rhodent.density_matrices.frequency.FrequencyDensityMatrixMetadata(density_matrices: FrequencyDensityMatrices, globalw: int, globalr: int, localw: int, localr: int)[source]
Metadata to the density matrices.
- Parameters:
density_matrices – Parent of this object.
globalw – Frequency index.
localw – Frequency index on this rank.
globalr – Real/imaginary index.
localr – Real/imaginary index on this rank.
- property freq: float
Frequency in units of eV.
- property global_indices
Unique index for this work.
- property reim: str
Returns real/imaginary tag
'Re'
or'Im'
.The tag corresponds to the Fourier transform of the real or imaginary part of the density matrix.
- class rhodent.density_matrices.time.ConvolutionDensityMatrixMetadata(density_matrices: ConvolutionDensityMatrices, globalt: int, globalp: int, localt: int, localp: int)[source]
Metadata to the density matrices.
Properties
- density_matrices
Parent of this object.
- globalt
Time index.
- localt
Time index on this rank.
- globalp
Pulse index.
- localp
Pulse index on this rank.
- property global_indices
Unique index for this work.
- property pulse: Perturbation
Pulse.
- property time: float
Simulation time in units of as.
Generators
These classes are used internally to iteratively obtain density matrices in
the form of DensityMatrix
object.
- class rhodent.density_matrices.ConvolutionDensityMatricesFromDisk(ksd, pulserho_fmt, pulses, times, derivative_order_s=[0], log=None, 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) –
Formatting string for the density matrices saved to disk.
Example:
pulserho_fmt =
pulserho/t{time:09.1f}{tag}.npy
.
Accepts variables
{time}
- Time in units of as.{tag}
- Derivative tag,''
,'-Iomega'
, or'-omega2'
.{pulsefreq}
- Pulse frequency in units of eV.{pulsefwhm}
- Pulse FWHM in units of fs.
pulses (Collection[PerturbationLike]) – 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 units of 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, perturbation, pulses, times, derivative_order_s=[0], log=None, 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) –
Formatting string for the density matrices in frequency space saved to disk.
Example:
frho_fmt =
frho/w{freq:05.2f}-{reim}.npy'
.
Accepts variables:
{freq}
- Frequency in units of eV.{reim}
- ‘Re’ or ‘Im’ for Fourier transform of real/imaginary part of density matrix.
perturbation (PerturbationLike) – Perturbation that was present during time propagation.
pulses (Collection[PerturbationLike]) – 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 units of 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, log=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 time-dependent wave functions file, which is read from disk.
- Parameters:
ksd (KohnShamDecomposition | str) – KohnShamDecomposition object or filename.
wfs_fname (str) – Filename of the time-dependent wave functions file.
perturbation (PerturbationLike) – Perturbation that was present during time propagation.
pulses (Collection[PerturbationLike]) – 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 units of 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.
stridet (int) – Skip this many steps when reading the time-dependent wave functions file.
- property myt: list[int]
List of indices corresponding to the time indices on held on this rank.
- 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.FrequencyDensityMatricesFromDisk(ksd, frho_fmt, perturbation, frequencies, real=True, imag=True, derivative_order_s=[0], calc_size=1, log=None)[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) –
Formatting string for the density matrices in frequency space saved to disk.
Example:
frho_fmt =
frho/w{freq:05.2f}-{reim}.npy
.
Accepts variables:
{freq}
- Frequency in units of eV.{reim}
-'Re'
or'Im'
for Fourier transform of real/imaginary part of density matrix.
perturbation (PerturbationLike) – The perturbation which the density matrices are a response to.
frequencies (list[float] | NDArray[np.float64]) – Compute density matrices for these frequencies (or as close to them as possible). In units of 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.
- class rhodent.density_matrices.FrequencyDensityMatricesFromWaveFunctions(ksd, wfs_fname, perturbation, frequencies, frequency_broadening=0, real=True, imag=True, derivative_order_s=[0], log=None, calc_size=1, stridet=1)[source]
Collection of density matrices in the Kohn-Sham basis in the frequency domain, for different frequencies. Obtained from the time-dependent wave functions 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 time-dependent wave functions file.
perturbation (PerturbationLike) – The perturbation which the density matrices are a response to.
frequencies (list[float] | NDArray[np.float64]) – Compute density matrices for these frequencies (or as close to them as possible). In units of eV.
frequency_broadening (float) – Gaussian broadening width in units of eV. Default (0) is no broadening.
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.
stridet (int) – Skip this many steps when reading the time-dependent wave functions file.
- property myw: list[int]
List of indices corresponding to the frequency indices on held on this rank.
- class rhodent.density_matrices.TimeDensityMatricesFromWaveFunctions(ksd, wfs_fname, times, real=True, imag=True, calc_size=1, log=None, stridet=1)[source]
Collection of density matrices in the Kohn-Sham basis for different times, obtained by reading the time-dependent wave functions file.
- Parameters:
ksd (KohnShamDecomposition | str) – KohnShamDecomposition object or filename.
wfs_fname (str) – Filename of the time-dependent wave functions file.
times (list[float] | NDArray[np.float64]) – Compute density matrices for these times (or as close to them as possible). In units of 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.
stridet (int) – Skip this many steps when reading the time-dependent wave functions file.
Readers
These classes are used internally to read and process density matrices in the form of
DensityMatrixBuffer
objects.
- class rhodent.density_matrices.readers.gpaw.BaseWfsReader(wfs_fname, comm=<MPI object>, yield_re=True, yield_im=True, stridet=1, tmax=-1, filter_times=None, log=None)[source]
Read wave functions or density matrices from the time-dependent wave functions file.
- Parameters:
wfs_fname (
str
) – Filename of the time-dependent wave functions file written by GPAW.comm – MPI communicator.
yield_re (
bool
) – Whether to read and yield the real part of wave functions/density matrices.yield_im (
bool
) – Whether to read and yield the imaginary part of wave functions/density matrices.stridet (
int
) – Skip this many steps when reading.tmax (
int
) – Last time index to read.filter_times (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – A list of times to read in atomic units. The closest times in the time-dependent wave functions file will be read. Applied after skipping the stridet steps to tmax.log (
Optional
[Logger
]) – Logger object.
- property comm: gpaw.mpi.Communicator
MPI communicator.
- property dt: float
Time step in atomic units.
- abstract iread(*args, **kwargs)[source]
Iteratively read wave functions or density matrices time by time.
- Return type:
Generator
[DensityMatrixBuffer
,None
,None
]
- abstract nnshape(*args, **kwargs)[source]
Shape of the density matrices or wave functions.
- Return type:
tuple
[int
,int
]
- property nt: int
Number of times to read.
- work_loop(rank)[source]
Yield the time indices that this rank will read.
New indices are yielded until the end of self.reader_t is reached (across all ranks).
- Yields:
Time index between 0 and len(self.reader_t) - 1 corresponding to
the time being read by this rank. Or None if this rank has nothing
to read.
- Return type:
Generator
[int
|None
,None
,None
]
- property yield_im: bool
Whether this object should read imaginary parts.
- property yield_re: bool
Whether this object should read real parts.
- class rhodent.density_matrices.readers.gpaw.KohnShamRhoWfsReader(wfs_fname, ksd, comm=<MPI object>, yield_re=True, yield_im=True, stridet=1, tmax=-1, filter_times=None, log=None, striden=0)[source]
Read density matrices from the time-dependent wave functions file.
Yield density matrices time by time.
- Parameters:
wfs_fname (str) – Filename of the time-dependent wave functions file.
ksd (str | KohnShamDecomposition) – KohnShamDecomposition object or filename to the ksd file.
comm – MPI communicator.
yield_re (bool) – Whether to read and yield the real part of wave functions/density matrices.
yield_im (bool) – Whether to read and yield the imaginary part of wave functions/density matrices.
stridet (int) – Skip this many steps when reading.
tmax (int) – Last time index to read.
filter_times (list[float] | NDArray[np.float64] | None) – A list of times to read in atomic units. The closest times in the time-dependent wave functions file will be read.
striden (int) – Option passed through to the LCAORhoWfsReader.
log (Logger | None) – Logger object.
- iread(s, k, n1, n2)[source]
Read the density matrices time by time.
- Parameters:
s (
int
) – Read these indices.k (
int
) – Read these indices.n1 (
slice
) – Read these indices.n2 (
slice
) – Read these indices.
- Return type:
Generator
[DensityMatrixBuffer
,None
,None
]
- property ksd: gpaw.lcaotddft.ksddecomposition.KohnShamDecomposition
Kohn-Sham decomposition object.
- class rhodent.density_matrices.readers.gpaw.LCAORhoWfsReader(wfs_fname, comm=<MPI object>, yield_re=True, yield_im=True, stridet=1, tmax=-1, filter_times=None, log=None, striden=4)[source]
Read density matrices in the LCAO basis from the time-dependent wave functions file.
Yield density matrices time by time.
- iread(s, k, M1, M2)[source]
Read the density matrices time by time.
- Parameters:
s (
int
) – Read these indices.k (
int
) – Read these indices.M1 (
slice
) – Read these indices.M2 (
slice
) – Read these indices.
- Return type:
Generator
[DensityMatrixBuffer
,None
,None
]
- class rhodent.density_matrices.readers.gpaw.WfsReader(wfs_fname, comm=<MPI object>, yield_re=True, yield_im=True, stridet=1, tmax=-1, filter_times=None, log=None)[source]
Read wave function LCAO coefficients from the time-dependent wave functions file.
Yield wave functions time by time.
- iread(s, k, n, M)[source]
Read the density matrices time by time.
- Parameters:
s (
int
) – Read these indices.k (
int
) – Read these indices.n (
slice
) – Read these indices.M (
slice
) – Read these indices.
- Return type:
Generator
[DensityMatrixBuffer
,None
,None
]
- rhodent.density_matrices.readers.gpaw.proxy_C_nM(reader, *indices)[source]
Proxy the wave function coefficients, with the correct units.
Make sure that the proxied array has at least two dimensions
- rhodent.density_matrices.readers.gpaw.proxy_coefficients(reader, *indices)[source]
Proxy the wave function coefficients, with the correct units.
- rhodent.density_matrices.readers.gpaw.proxy_occupations(reader, *indices)[source]
Proxy the wave function occupations with the correct scale.
- class rhodent.density_matrices.distributed.base.BaseDistributor(rho_reader, parameters=None, comm=None)[source]
Distribute density matrices over time, frequency or other dimensions across MPI ranks
- property comm: gpaw.mpi.Communicator
MPI communicator.
- classmethod from_parameters(wfs_fname, ksd, comm=<MPI object>, yield_re=True, yield_im=True, stridet=1, log=None, verbose=False, **kwargs)[source]
Set up this class, trying to enforce memory limit.
- Parameters:
wfs_fname (str) – File name of the time-dependent wave functions file.
ksd (KohnShamDecomposition | str) – KohnShamDecomposition object or filename to the ksd file.
comm – MPI communicator.
yield_re (bool) – Whether to read and yield the real part of wave functions/density matrices.
yield_im (bool) – Whether to read and yield the imaginary part of wave functions/density matrices.
stridet (int) – Skip this many steps when reading the time-dependent wave functions file.
log (Logger | None) – Logger object.
verbose (bool) – Be verbose in the attempts to satisfy memory requirement.
kwargs – Options passed through to
from_reader()
.
- abstract classmethod from_reader(rho_nn_reader, parameters, **kwargs)[source]
Set up this class from a density matrix reader and parameters object
- Return type:
- property ksd: gpaw.lcaotddft.ksddecomposition.KohnShamDecomposition
Kohn-Sham decomposition object.
- property maxnchunks: int
Maximum number of ranks participating in reading of chunks.
- property maxntimes: int
Maximum number of ranks participating in reading of times.
- property niters: int
Number of iterations needed to read all chunks.
- work_loop(rank)[source]
Like work_loop_by_rank but for one particular rank
- Return type:
Generator
[RhoIndices
|None
,None
,None
]
- work_loop_by_ranks()[source]
Yield slice objects corresponding to the chunk of the density matrix that is gathered on each rank.
New indices are yielded until the entire density matrix is processed (across all ranks).
- Yields:
List of slice objects corresponding to part of the density matrix
yielded on each ranks. None in place of the slice object if there is
nothing yielded for that rank.
- Return type:
Generator
[list
[RhoIndices
|None
],None
,None
]
- abstract property xshape: tuple[int, ...]
Shape of x-dimension in buffers.
- property yield_im: bool
Whether imaginary part of density matrices is calculated.
- property yield_re: bool
Whether real part of density matrices is calculated.
- class rhodent.density_matrices.distributed.base.RhoIndices(s, k, n1, n2)[source]
-
k:
int
Alias for field number 1
-
n1:
slice
Alias for field number 2
-
n2:
slice
Alias for field number 3
-
s:
int
Alias for field number 0
-
k:
- class rhodent.density_matrices.distributed.base.RhoParameters(ns: int, nk: int, n1min: int, n1max: int, n2min: int, n2max: int, striden1: int = 4, striden2: int = 4)[source]
Utility class to describe density matrix indices.
- Parameters:
ns – Number of spins.
nk – Number of kpoints.
n1min – Smallest index of row to read.
n1max – Largest index of row to read.
n2min – Smallest index of column to read.
n2max – Largest index of column to read.
striden1 – Stride for reading rows. Each chunk will be this size in the first dimension.
striden2 – Stride for reading columns. Each chunk will be this size in the second dimension.
- classmethod from_ksd(ksd, comm=None, only_ia=True, chunk_iterations=1, **kwargs)[source]
Initialize from KohnShamDecomposition.
- Parameters:
ksd (KohnShamDecomposition) – KohnShamDecomposition.
comm (Communicator | None) – MPI Communicator.
only_ia (bool) –
True
if the parameters should be set up such that the electron-hole part of the density matrix is read, otherwise full density matrix.chunk_iterations (int) – Attempt to set up the strides so that the total number of chunks is as close as possible but not more than the number of MPI ranks times
chunk_iterations
.kwargs – Options passed through to the constructor.
- Return type:
RhoParameters
- property full_nnshape: tuple[int, int]
Shape of the full density matrix to be read.
- iterate_indices()[source]
Iteratively yield indices slicing chunks of the density matrix.
- Return type:
Generator
[RhoIndices
,None
,None
]
-
n1max:
int
Alias for field number 3
-
n1min:
int
Alias for field number 2
- property n1size: int
Size of full density matrix in the first dimension.
-
n2max:
int
Alias for field number 5
-
n2min:
int
Alias for field number 4
- property n2size: int
Size of full density matrix in the first dimension.
-
nk:
int
Alias for field number 1
- property nnshape: tuple[int, int]
Shape of each density matrix chunk.
-
ns:
int
Alias for field number 0
-
striden1:
int
Alias for field number 6
-
striden2:
int
Alias for field number 7
- class rhodent.density_matrices.distributed.time.AlltoallvTimeDistributor(rho_reader, parameters=None)[source]
Iteratively read density matrices in the time domain
This class uses the KohnShamRhoWfsReader to iteratively read wave functions (one time per rank) and construct density matrices in the ground state Kohn-Sham basis for each time. When all ranks have read one time each, this class performs a redistribution of data, such that each rank only gets one chunk of the density matrices, but for all times. The density matrices are accumulated in a buffer and yielded when all times have been read.
- Parameters:
rho_reader (
KohnShamRhoWfsReader
) – Density matrices reader
- class rhodent.density_matrices.distributed.time.TimeDistributor(rho_reader, parameters=None, comm=None)[source]
Iteratively read density matrices in the time domain
This class uses the KohnShamRhoWfsReader to iteratively read wave functions (each rank reading the same times) and construct density matrices in the ground state Kohn-Sham basis for each time. The different ranks are reading different chunks of the density matrices. The density matrices are accumulated in a buffer and yielded when all times have been read.
- Parameters:
rho_reader (KohnShamRhoWfsReader) – Density matrices reader
comm (Communicator | None) – Communicator
log – Logger
- property dtype
Dtype of buffers.
- classmethod from_reader(rho_nn_reader, parameters, **kwargs)[source]
Set up this class from a density matrix reader and parameters object
- Return type:
- property xshape
Shape of x-dimension in buffers.
- class rhodent.density_matrices.distributed.frequency.FourierTransformer(rho_nn_reader, perturbation, filter_frequencies=None, frequency_broadening=0, result_on_ranks=[])[source]
Iteratively take the Fourier transform of density matrices.
- Parameters:
rho_nn_reader (
TimeDistributor
) – Object that can iteratively read density matrices in the time domain, distributed such that different ranks have different chunks of the density matrix, but each ranks has all times for the same chunk.perturbation (
Union
[Perturbation
,Laser
,dict
,None
]) – The perturbation which the density matrices are a response to.filter_frequencies (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – After Fourier transformation keep only these frequencies (or the frequencies closest to them). In atomic units.frequency_broadening (
float
) – Gaussian broadening width in atomic units. Default (0) is no broadening.result_on_ranks (
list
[int
]) – List of ranks among which the resulting arrays will be distributed. Empty list (default) to distribute among all ranks.
- create_out_buffer()[source]
Create the DensityMatrixBuffer to hold the temporary density matrix after each redistribution
- Return type:
- create_result_buffer()[source]
Create the DensityMatrixBuffer to hold the resulting density matrix
- Return type:
- property dist_buffer: DensityMatrixBuffer
Buffer of density matrices on this rank after redistribution
- property dtype
Dtype of buffers.
- classmethod from_reader(rho_nn_reader, parameters, *, perturbation, filter_frequencies=None, frequency_broadening=0, result_on_ranks=[])[source]
Set up this class from a density matrix reader and parameters object
- Return type:
- property nranks_result: int
Number of ranks that the resulting arrays will be distributed among
- redistribute()[source]
Perform the Fourier transform and redistribute the data
When the Fourier transform is performed, the data is distributed such that each rank stores the entire time/frequency series for one chunk of the density matrices, i.e. indices n1, n2.
This function then performs a redistribution of the data such that each rank stores full density matrices, for certain frequencies.
If the density matrices are split into more chunks than there are ranks, then the chunks are read, Fourier transformed and distributed in a loop several times until all data has been processed.
- Return type:
- Returns:
Density matrix buffer with x-dimensions (Number of local frequencies, )
where the Number of local frequencies variers between the ranks.
- property result_on_ranks: list[int]
Set of ranks among which the result will be distributed
- property xshape
Shape of x-dimension in buffers.
- class rhodent.density_matrices.distributed.pulse.PulseConvolver(rho_nn_reader, perturbation, pulses, derivative_order_s=[0], filter_times=None, result_on_ranks=[])[source]
Class performing pulse convolution of density matrices.
The procedure of the pulse convolution is the following:
The entire time series of (real and/or imaginary parts of) density matrices is read, in several chunks of indices (n1, n2). Each MPI rank works on different chunks.
Each chunk is Fourier transformed, divided by the Fourier transform of the original perturbation and multiplied by the Fourier transform of the new pulse(s).
Optionally, derivatives are computed by multiplying the density matrices in the frequency domain by factors of \(i \omega\).
Each chunk is inverse Fourier tranformed, and only selected times are kept.
Additionally, this class can redistribute the resulting convoluted density matrices so that the each rank holds the entire density matrix, for a few times.
- Parameters:
rho_nn_reader (
TimeDistributor
) – Object begin able to iteratively read density matrices in the time domain. Density matrices are split in chunks and distributed among ranks.perturbation (
Union
[Perturbation
,Laser
,dict
,None
]) – The perturbation which the density matrices read byrho_nn_reader
are a response to.pulses (
Collection
[Union
[Perturbation
,Laser
,dict
,None
]]) – List of pulses to perform to convolution with.derivative_order_s (
list
[int
]) – List of derivative orders to compute.filter_times (
Union
[list
[float
],ndarray
[tuple
[int
,...
],dtype
[float64
]],None
]) – After convolution keep only these times (or the times closest to them). In atomic units.result_on_ranks (
list
[int
]) – List of ranks among which the resulting arrays will be distributed. Empty list (default) to distribute among all ranks.
- create_out_buffer()[source]
Create the DensityMatrixBuffer to hold the temporary density matrix after each redistribution
- Return type:
- create_result_buffer()[source]
Create the DensityMatrixBuffer to hold the resulting density matrix
- Return type:
- property dist_buffer: DensityMatrixBuffer
Buffer of denisty matrices on this rank after redistribution.
- property dtype
Dtype of buffers.
- classmethod from_reader(rho_nn_reader, parameters, perturbation, pulses, derivative_order_s=[0], filter_times=None, result_on_ranks=[])[source]
Set up this class from a density matrix reader and parameters object
- Return type:
- property nlocalt: int
Number of times stored on this rank after redistribution of the result.
- property nranks_result: int
Number of ranks storing part of the result after redistribution.
- property nt: int
Number of times for which convoluted density matrices are calculated.
- redistribute()[source]
Perform the pulse convolution and redistribute the resulting density matrices.
During the pulse convolution step, the data is distributed such that each rank stores the entire time series for one chunk of the density matrices, i.e. indices n1, n2.
This function then performs a redistribution of the data such that each rank stores full density matrices, for certain times.
If the density matrices are split into more chunks than there are ranks, then the chunks are read, convoluted with pulses and distributed in a loop several times until all data has been processed.
- Return type:
- Returns:
Density matrix buffer with x-dimensions (number of pulses, number of local times)
where the number of local times variers between the ranks.
- property result_on_ranks: list[int]
Set of ranks among which the result will be distributed.
- property time_t: ndarray[tuple[int, ...], dtype[float64]]
Array of times corresponding to convoluted density matrices; in atomic units.
- property upscaling: int
Upscaling factor.
Data is upscaled in time by this factor during the Fourier transformation step, in order to calculate convoluted density matrices at a denser grid of times than what is present in the time-dependent wave functions file.
- property xshape
Shape of x-dimension in buffers.