Voronoi¶
- class rhodent.voronoi.EmptyVoronoiWeights[source]¶
- property atom_projections¶
Atom projections
- property dtype¶
dtype of wave functions
- property nn¶
Global 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.VoronoiLCAOReader(ulm_fname, comm_is_domain=False, comm=None)[source]¶
Read Voronoi weights in LCAO basis from ulm file.
- Parameters:
ulm_fname (
str
) – filenamebroadcast_weights – If true, the array of weights is broadcasted and yielded on all ranks
comm – GPAW MPI communicator object. Defaults to world
- property C_nM¶
LCAO wave function coefficients. Distributed over bands.
- property P_ani: ArrayDict¶
PAW projectors. Distributed over domains and bands.
\[P_{ni}^a = \left<\tilde{p}_i^a | \tilde{\psi}_n \right>\]
- property nM: int¶
Global number of bands
- property nn: int¶
Global 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.VoronoiLCAOWeightCalculator(atom_projections, gpw_file, voronoi_grid_file=None, recalculate_grid=False, domain=-1, use_pblas=False, comm=None)[source]¶
Loads Voronoi grid and calculates Voronoi weights from ground state file using LCAO basis function overlaps and PAW corrections.
\[W_{nn'} = \left<\psi_n|\hat{w}|\psi_{n'}\right> = \int w(\vec{r}) \psi_n^*(\vec{r}) \psi_{n'}(\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.
- Parameters:
- property C_nM¶
LCAO wave function coefficients. Distributed over bands.
- property P_ani¶
PAW projectors. Distributed over domains and bands.
\[P_{ni}^a = \left<\tilde{p}_i^a | \tilde{\psi}_n \right>\]
- property a_G¶
Returns the voronoi decomposition on the coarse grid, distributed over domain_comm
If this quantity is not stored, the voronoi decomposition on the master rank is distributed
- distribute_a_G()[source]¶
Distribute the voronoi decomposition on the coarse grid over domain_comm
Afterwards delete the non-distributed array
- property nM: int¶
Global number of bands
- property nn: int¶
Global 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.VoronoiLCAOWeights[source]¶
- abstract property C_nM¶
LCAO wave function coefficients. Distributed over bands.
- abstract property P_ani: ArrayDict¶
PAW projectors. Distributed over domains and bands.
\[P_{ni}^a = \left<\tilde{p}_i^a | \tilde{\psi}_n \right>\]
- property comm: gpaw.mpi.Communicator¶
Communicator
- property dS_aii¶
Overlap matrix PAW corrections. Same data on all ranks
\[\Delta S_{ij}^a= = \left<\phi_i^a|\phi_j^a\right> - \left<\tilde{\phi}_i^a|\tilde{\phi}_j^a\right>\]
- property dist_C_nM¶
LCAO wave function coefficients. Distributed in block cyclic form.
- property gather_P_ani: ArrayDict¶
PAW projectors. Gathered to band comm rank 0. None on all other band comm ranks.
- abstract property nM: int¶
Global number of bands
- abstract property nn: int¶
Global number of bands
- property nproj: int¶
Number of atomic projections
- 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.VoronoiReader(ulm_fname, broadcast_weights=False, comm=None)[source]¶
Read Voronoi weights from ulm file.
- Parameters:
ulm_fname (
str
) – filenamebroadcast_weights (
bool
) – If true, the array of weights is broadcasted and yielded on all rankscomm – GPAW MPI communicator object. Defaults to world
- property nn: int¶
Global 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.VoronoiWeightCalculator(voronoi_lcao_gen, use_pblas=False)[source]¶
Calculates KS wave function overlaps under operator \(\hat{w}\) using LCAO basis function overlaps and PAW corrections.
\[W_{nn'} = \left<\psi_n | \hat{w} | \psi_{n'}\right> = \left<\tilde{\psi}_n | \hat{T}^\dagger \hat{w} \hat{T} | \tilde{\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
proj_atoms
and 0 outside.The sum over atoms is restricted to atoms in the region \(w(\vec{r}) = 1\).
\(\hat{T}\) is the PAW transformation operator.
- Parameters:
voronoi_lcao_gen (
VoronoiLCAOWeights
) – Object that calculates or loads the LCAO weights from fileuse_pblas (
bool
) – Whether PBLAS/ScaLAPACK should be used (more efficient parallelization)
- calculate_weight_nn(proj_atoms, weight_MM, out=None)[source]¶
Calculates weights.
- Parameters:
proj_atoms (
list
[int
] |ndarray
[Any
,dtype
[float64
]]) – Do the Voronoi decomposition for these atomsweight_MM (
ndarray
[Any
,dtype
[float64
]] |None
) – The weights in LCAO without PAW correctionsout (
Optional
[ndarray
[Any
,dtype
[float64
]]]) – If this is specified the result is written to here on rank 0
- Return type:
- Returns:
The matrix elements \(W_{nn'}\) on root and
None
on other ranks
- calculate_weight_nn_pblas(proj_atoms, weight_MM, out=None)[source]¶
Calculates weights using PBLAS.
- Parameters:
proj_atoms (
list
[int
] |ndarray
[Any
,dtype
[float64
]]) – Do the Voronoi decomposition for these atomsweight_MM (
ndarray
[Any
,dtype
[float64
]] |None
) –The weights in LCAO without PAW corrections.
The weights should be summed to the rank 0 of the domain communicator The weights on other ranks than 0 on the band communicator are ignored (they should be identical on all ranks)
out (
Optional
[ndarray
[Any
,dtype
[float64
]]]) – If this is specified the result is written to here on rank 0
- Return type:
- Returns:
The matrix elements \(W_{nn'}\) on root and
None
on other ranks
- property comm: gpaw.mpi.Communicator¶
Communicator
- property nn: int¶
Global 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.VoronoiWeights[source]¶
- abstract property atom_projections: Sequence[list[int] | ndarray[Any, dtype[float64]]]¶
Atom projections
- property comm: gpaw.mpi.Communicator¶
Communicator
- abstract property nn: int¶
Global number of bands
- 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.
- rhodent.voronoi.calc_correction(correction_on_atoms, P_ani, dS_aii, Nn, out=None)[source]¶
Calculate the PAW correction to KS wave function overlaps on this domain.
\[c_{nn'} = \sum_a \left<\tilde{\psi}_n | \tilde{p}_i^a\right> \Delta S_{ij}^a \left<\tilde{p}_j^a | \tilde{\psi}_{n'}\right>,\]where the sum over atoms only includes a selection of atoms.
- Parameters:
correction_on_atoms (
list
[int
] |ndarray
[Any
,dtype
[float64
]]) – Which atoms to loop overP_ani (
dict
[int
,ndarray
[Any
,dtype
[float64
]]] |ArrayDict
) – The matrix of projectors \(\left<\tilde{p}_i^a | \tilde{\psi}_{n}\right>\)dS_aii (
dict
[int
,ndarray
[Any
,dtype
[float64
]]]) – \(\Delta S_{ij}^a\)Nn (
int
) – Global number of bandsout (
Optional
[ndarray
[Any
,dtype
[float64
]]]) – If this is specified the result is appended to here
- Return type:
- Returns:
\(c_{nn'}'\)
- rhodent.voronoi.calc_vt_MM(calc, w_G)[source]¶
Calculate LCAO basis function overlaps.
\[\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}\]- Parameters:
- Return type:
- Returns:
The contribution to \(\tilde{v}_{\mu\nu}\) from this domain
- rhodent.voronoi.calc_vt_nn(C_nM, vt_MM)[source]¶
Calculate pseudo KS wave function overlaps from LCAO overlaps.
\[\widetilde{W}_{nn'} = \left<\tilde{\psi}_n | \hat{w} | \tilde{\psi}_{n'}\right> = \sum_{\mu\nu} C_{\mu n}^* C_{\nu n'} \left<\tilde{\phi}_{\mu} | \hat{w} | \tilde{\phi}_{\nu}\right> = \sum_{\mu\nu} C_{\mu n}^* C_{\nu n'} \tilde{v}_{\mu\nu}\]
- rhodent.voronoi.calc_vt_nn_pblas(dist_C_nM, dist_vt_MM, dist_vt_nn)[source]¶
See
calc_vt_nn()
.- Parameters:
dist_C_nM – LCAO coefficients, in block cyclic form
dist_vt_MM – The contribution to :math:tilde{v}_{munu}` on domain, in block cyclic form
dist_vt_nn – Buffer in block cyclic form, into which the contribution to \(\widetilde{W}_{nn}\) on this domain will be written
- rhodent.voronoi.calculate_a_g(gpw_file=None, fine_grid=True, calc=None)[source]¶
Calculate Voronoi grid
- Parameters:
gpw_file (str | None) – Filename of ground state file. Ignored if calc is not None
calc (GPAWCalculator | None) – Initialized GPAW calculator
fine_grid (bool) – True if Voronoi grid on fine grid should be calculated, False if on coarse
- Return type:
Voronoi grid as nx,ny,nz array