Voronoi

class rhodent.voronoi.VoronoiWeights(comm=None, log=None)[source]

Abstract base class for Voronoi weights in Kohn-Sham basis.

property root: bool

Whether this rank is the root rank.

property log: Logger

Logger object.

abstract property atom_projections: Sequence[list[int] | ndarray[tuple[int, ...], dtype[float64]]]

Atom projections.

abstract property nn: int

Number of bands.

property nproj: int

Number of projections.

property comm: gpaw.mpi.Communicator

MPI Communicator.

abstract property saved_fields: dict[str, Any]

Saved data fields associated with the object.

If this object is read from a ULM file there might be some extra data in that file. Return that data, or an empty dict if there is none.

class rhodent.voronoi.VoronoiWeightCalculator(voronoi_lcao)[source]

Calculate Voronoi weights in the basis of ground state Kohn-Sham orbitals.

For each atomic projection, calculates

\[W_{nn'} = \left<\psi_n | \hat{w} | \psi_{n'}\right> = \left<\tilde{\psi}_n | \hat{w} | \tilde{\psi}_{n'}\right> + \sum_{aij} \left<\tilde{\psi}_n | \tilde{p}_i^a\right> \Delta S_{ij}^a \left<\tilde{p}_i^a | \tilde{\psi}_{n'}\right>\]

where the operator \(\hat{w} = w(\vec{r})\) is 1 in the Voronoi region of atoms and 0 outside. The sum over atoms is restricted to atoms in the region \(w(\vec{r}) = 1\).

The PAW projectors

\[P_{ni}^a = \left<\tilde{p}_i^a | \tilde{\psi}_n \right>\]

and overlap matrix

\[\Delta S_{ij}^a = \left<\phi_i^a|\phi_j^a\right> - \left<\tilde{\phi}_i^a|\tilde{\phi}_j^a\right>\]

are obtained from the GPAW calculator and the smooth weight matrix is obtained by a basis change of the smooth LCAO weights

\[\left<\tilde{\psi}_n | \hat{w} | \tilde{\psi}_{n'}\right> = \sum_{\mu\nu} C_{n\mu} C_{n'\nu} \left<\tilde{\phi}_\mu | \hat{w} | \tilde{\phi}_{\nu}\right>.\]
Parameters:

voronoi_lcao (VoronoiLCAOWeights) – Object that calculates or loads the LCAO weights from file.

property atom_projections: Sequence[list[int] | ndarray[tuple[int, ...], dtype[float64]]]

Atom projections.

property nn: int

Number of bands.

property comm: gpaw.mpi.Communicator

MPI Communicator.

property saved_fields: dict[str, Any]

Saved data fields associated with the object.

If this object is read from a ULM file there might be some extra data in that file. Return that data, or an empty dict if there is none.

calculate_and_write(out_fname, write_extra={})[source]

Calculate the Voronoi weights in the Kohn-Sham basis and write to file.

The weights are saved in a numpy archive if the file extension is .npz or in a ULM file if the file extension is .ulm.

Parameters:
  • out_fname (str) – File name of the resulting data file.

  • write_extra (dict[str, Any]) – Dictionary of additional data written to the file.

classmethod from_gpw(atom_projections, gpw_file, voronoi_grid=None, comm=None)[source]

Set up calculation from GPAW file.

Parameters:
  • atom_projections (AtomProjectionsType) – List of atom groups. Each atom group is a list of integers (of any length).

  • gpw_file (str) – Filename of GPAW ground state file.

  • voronoi_grid (VoronoiGrid | None | str) – Voronoi grid, or None to calculate it, or a filename to read it from file.

  • comm (Communicator | None) – Communicator.

Return type:

VoronoiWeightCalculator

classmethod from_lcao_file(voronoi_lcao_file, comm=None)[source]

Set up calculation from file of already computed weights in LCAO basis.

Parameters:
  • voronoi_lcao_file (str) – Filename of file with Voronoi weights in LCAO basis.

  • comm (Communicator | None) – Communicator.

Return type:

VoronoiWeightCalculator

class rhodent.voronoi.VoronoiReader(ulm_fname, comm=None)[source]

Read Voronoi weights from ulm file.

Parameters:
  • ulm_fname (str) – Filename.

  • comm (Communicator | None) – GPAW MPI communicator object. Defaults to world.

property atom_projections: Sequence[list[int] | ndarray[tuple[int, ...], dtype[float64]]]

Atom projections.

property nn: int

Number of bands.

property saved_fields: dict[str, Any]

Saved data fields associated with the object

If this object is read from a ULM file there might be some extra data in that file. Return that data, or an empty dict if there is none.

class rhodent.voronoi.EmptyVoronoiWeights(comm=None, log=None)[source]

Object representing a lack of weights.

property atom_projections

Atom projections.

property nn

Number of bands.

property saved_fields

Saved data fields associated with the object.

If this object is read from a ULM file there might be some extra data in that file. Return that data, or an empty dict if there is none.

class rhodent.voronoi.VoronoiLCAOWeights[source]

Abstract base class for Voronoi weights in LCAO basis.

log_parallel(*args, **kwargs)[source]

Log message with communicator information.

Return type:

Logger

abstract property nn: int

Number of bands.

abstract property nM: int

Number of basis functions.

property nproj: int

Number of atomic projections

property atom_projections: Sequence[list[int] | ndarray[tuple[int, ...], dtype[float64]]]

List of atom projections.

property calc: GPAWCalculator | None

GPAW calculator instance.

property comm: gpaw.mpi.Communicator

MPI Communicator.

property root: bool

Whether this rank is the root rank.

property C_nM: ndarray[tuple[int, ...], dtype[float64]] | None

LCAO wave function coefficients on root rank, None on other ranks.

property P_ani: dict[int, ndarray[tuple[int, ...], dtype[float64]]] | None

PAW projectors \(P_{ni}^a\) on the root rank, None on other ranks.

\[P_{ni}^a = \left<\tilde{p}_i^a | \tilde{\psi}_n \right>\]
property dS_aii: dict[int, ndarray[tuple[int, ...], dtype[float64]]] | None

Overlap matrix PAW corrections \(\Delta S_{ij}^a\) on the root rank, None on other ranks.

\[\Delta S_{ij}^a = \left<\phi_i^a|\phi_j^a\right> - \left<\tilde{\phi}_i^a|\tilde{\phi}_j^a\right>\]
abstract property saved_fields: dict[str, Any]

Saved data fields associated with the object.

If this object is read from a ULM file there might be some extra data in that file. Return that data, or an empty dict if there is none.

class rhodent.voronoi.VoronoiLCAOReader(ulm_fname, comm=None)[source]

Read Voronoi weights in the LCAO basis from ULM file.

Parameters:
  • ulm_fname (str) – Filename.

  • comm – GPAW MPI communicator object. Defaults to world.

property nn: int

Number of bands.

property nM: int

Number of basis functions.

property saved_fields: dict[str, Any]

Saved data fields associated with the object.

If this object is read from a ULM file there might be some extra data in that file. Return that data, or an empty dict if there is none.

class rhodent.voronoi.VoronoiLCAOWeightCalculator(atom_projections, gpw_file, voronoi_grid=None, domain=-1, comm=None)[source]

Calculate Voronoi weights in LCAO basis.

For each atomic projection, calculates

\[\tilde{v}_{\mu\nu} = \left<\tilde{\phi}_{\mu} | \hat{w} | \tilde{\phi}_{\nu}\right> = \int w(\vec{r}) \tilde{\phi}_{\mu}^*(\vec{r}) \tilde{\phi}_{\nu}(\vec{r}) d\vec{r}\]

where the operator \(\hat{w} = w(\vec{r})\) is 1 in the Voronoi region of the atomic projections and 0 outside, and \(\tilde{\phi}_{\nu}\) are the smooth LCAO basis functions.

Parameters:
  • atom_projections (AtomProjectionsType) – List of atom groups. Each atom group is a list of integers (of any length).

  • gpw_file (str) – Filename of GPAW ground state file.

  • voronoi_grid (VoronoiGrid | None | str) – Voronoi grid, or None to calculate it, or a filename to read it from file.

  • domain (int) – Domain size.

  • comm (Communicator | None) – Communicator.

property nn: int

Number of bands.

property nM: int

Number of basis functions.

property calc: gpaw.GPAW

GPAW calculator instance.

property saved_fields: dict[str, Any]

Saved data fields associated with the object.

If this object is read from a ULM file there might be some extra data in that file. Return that data, or an empty dict if there is none.

calculate_and_write(out_fname, write_extra={})[source]

Calculate the Voronoi weights in the LCAO basis and write to file.

The weights are saved in a numpy archive if the file extension is .npz or in a ULM file if the file extension is .ulm.

Parameters:
  • out_fname (str) – File name of the resulting data file.

  • write_extra (dict[str, Any]) – Dictionary of additional data written to the file.

class rhodent.voronoi.VoronoiGrid[source]
abstract a_G(calc, log)[source]

Voronoi grid (on the coarse grid of the GPAW calculator) distributed on the domain communicator.

Each element in the grid is an integer corresponding to the closest atom.

Parameters:

calculator. (GPAW)

Return type:

ndarray[tuple[int, ...], dtype[int64]]

write(filename, calc, log)[source]

Write grid to file.

Parameters:
  • filename (str) – Filename.

  • calc (gpaw.GPAW) – GPAW calculator.

class rhodent.voronoi.VoronoiGridReader(filename)[source]

Read Voronoi grid from file.

a_G(calc, log)[source]

Voronoi grid (on the coarse grid of the GPAW calculator) distributed on the domain communicator.

Each element in the grid is an integer corresponding to the closest atom.

Parameters:

calculator. (GPAW)

Return type:

ndarray[tuple[int, ...], dtype[int64]]

class rhodent.voronoi.VoronoiGridCalculator[source]

Calculate the Voronoi grid.

a_G(calc, log)[source]

Voronoi grid (on the coarse grid of the GPAW calculator) distributed on the domain communicator.

Each element in the grid is an integer corresponding to the closest atom.

Parameters:

calculator. (GPAW)

Return type:

ndarray[tuple[int, ...], dtype[int64]]

rhodent.voronoi.proj_as_dict_on_master(proj, n1, n2)[source]

Collect the projectors to a dictionary on the root rank.

Return type:

dict[int, ndarray[tuple[int, ...], dtype[float64]]]