[QE-developers] R: Parallelization question for G-wise projector-wavefunctions product computation
Stefano de Gironcoli
degironc at sissa.it
Sat Feb 8 19:41:14 CET 2025
not sure I understand what you need.
calbec computes the dot product <\beta_{R,l}|\psi_{ps}> and I don't
think in the paw formalism there is any use for the individual Fourier
components.
if you want the shape of a pw, say |k+G>, "dressed" with paw terms
(like one would have, conceptually, in the augmented plane waves method,
where pw in the interstitial region are matched to atomic solutions
inside the atomic sphere) you need to call calbec with a fake psi where
the |k+G> component is set to 1 and all the rest is set to 0. and then
apply the formula that you wrote... except that one does not know where
to represent that expression... certainly not on the fft grid because
the atomic partial waves (especially the AE ones) have Fourier
components that exceed the cutoff of the fft grid, and in fact are kept
on the radial logarithmic grid and never used anywhere else.
The matrix elements of a physical quantities of interest (if the
operator is local like the potential or the kinetic energy) are computed
exploiting the (assumed) completeness of the partial wave expansion in
the atomic sphere
I_\Omega_{R} x \sum_I |\phi_{R,l}^{ps}> <\beta_{R,l}| \approx
I_ \Omega_{R}
where I_\Omega_{R} is the identity inside the atomic sphere of atom R
and zero outside.
One can then write
<\psi_{ae}_i|OP|\psi_{ae}_j> = <\psi_{ps}_i|OP|\psi_{ps}_j> +
\sum_{R,l,J} [ <\phi_{R,l}^{ae}|OP|\phis_{R,J}^{ae}> -
<\phi_{R,l}^{ps}|OP|\phi_{R,J}^{ps}> ] x
<\psi_{ps}_i|\beta_{R,I}><\beta_{R,J}|\psi_{ps}_j>
where the first row is computed on the FFT grid, the second on the
atomic radial grids and the third row terms by calbec (and accumulated
if needed in sumbec)
electrostatics is a bit more complicated
stefano
On 08/02/25 18:22, Marin Luca wrote:
> Dear Stefano,
>
> Thanks a lot for the kind and fast reply! This is a good question, I
> also wondered that maybe the storage part was useless. I need this
> term because it enters the PAW equation for the "all-electron
> reconstruction" of the pseudopotential KS-states:
>
> \psi_{ae} = \psi_{ps} + \sum_{R,l} [\phi_{R,l}^{ae} - \phi_{R,l}^{ps}]
> x <\beta_{R,l}|\psi_{ps}>
>
> where R and l are the atomic species and projectors indices, \phi the
> "atomic-like" states from the PP calculation and \psi the KS state
> (there might be slight mistake in the formula).
>
> If I now want to apply this correction to each of the u_{k+G} pw
> coefficients that contributes to \psi_{ps}, I think I also need the
> <\beta_{R,l}|\psi_{ps}> in a "G-wise" fashion, as I was naming it in
> the question.
>
> I was therefore implementing this routine inside "gipaw", which
> already has precious stuff written for this PAW reconstruction.
>
> I hope this gives you more information, and please let me know if you
> think there is already a good place either in pw or gipaw where to
> perform this!
>
> Best,
> Luca
> ETH Zurich, PhD student
> ------------------------------------------------------------------------
> *Da:* developers <developers-bounces at lists.quantum-espresso.org> per
> conto di Stefano de Gironcoli <degironc at sissa.it>
> *Inviato:* sabato 8 febbraio 2025 15:55
> *A:* developers at lists.quantum-espresso.org
> <developers at lists.quantum-espresso.org>
> *Oggetto:* Re: [QE-developers] Parallelization question for G-wise
> projector-wavefunctions product computation
>
> Dear Luca,
>
> the array betapsi_g(1:npw,1:nkb,1:nbnd) looks very big to me,
> unless the system is very small
>
> what do you need it for ? maybe one can compute the target quantity
> on the fly rather than storing this intermediate result
>
>
> stefano
>
>
> On 08/02/25 15:31, Marin Luca wrote:
>> Dear QE developers,
>> I had a conceptual question of how I should correctly modify a
>> subroutine in order to avoid parallelization and memory’s
>> absurdities. Specifically, my idea is to slightly modify one of the
>> subroutines in the “calbec” (becmod.f90) interface that is computing
>> <\Beta|\psi> products, and I am using QE v.6.4 because I’m working on
>> a specific code repository using this version. By calling such or
>> similar line:
>> CALL ZGEMV( 'C', npw, nkb, (1.0_DP,0.0_DP), beta, npwx, psi, 1, &
>> (0.0_DP, 0.0_DP), betapsi, 1 )
>> and subsequently summing through bands:
>> CALL mp_sum( betapsi( :, 1:m ), intra_bgrp_comm )
>> “calbec_k” computes the betapsi product as a 2D matrix with dimension
>> (number of projectors, number of bands). My goal is to compute the
>> same product but accessing the elements at each G vectors, i.e.
>> without summing over the G vectors themselves. To this end, my (very
>> basic, probably terrible) idea was to define another subroutine with
>> an extra variable like this (without breaking the requirement of
>> having unambiguous interface procedures in “calbec”):
>> calbec_k_modified ( npw, beta, psi, betapsi, betapsi_g, nbnd )
>> COMPLEX (DP), INTENT (out) :: betapsi_g(:,:,:)
>> and explicitly computing betapsi_g with nested do loops:
>> DO k = 1, npw
>> DO i = 1, nkb
>> DO j = 1, m
>> betapsi_g(k, i, j) = CONJG(beta(k, i)) * psi(k, j)
>> ENDDO
>> ENDDO
>> ENDDO
>> I am very unsure of whether this is a feasible approach and how I
>> should handle the MPI parallelization in this case (I have little
>> experience in QE development). I read the parallelization section of
>> the developer’s manual, but I’m still quite confused. Any
>> suggestions, critics or explanations will be very appreciated.
>> Best regards,
>> Luca Marin
>> ETH Zurich, PhD student
>>
>>
>>
>> ________________________________________________
>> The Quantum ESPRESSO Foundation stands in solidarity with all civilians worldwide who are victims of terrorism, military aggression, and indiscriminate warfare.
>> _______________________________________________
>> developers mailing list
>> developers at lists.quantum-espresso.org <mailto:developers at lists.quantum-espresso.org>
>> https://lists.quantum-espresso.org/mailman/listinfo/developers <https://lists.quantum-espresso.org/mailman/listinfo/developers>
>
> ________________________________________________
> The Quantum ESPRESSO Foundation stands in solidarity with all civilians worldwide who are victims of terrorism, military aggression, and indiscriminate warfare.
> _______________________________________________
> developers mailing list
> developers at lists.quantum-espresso.org
> https://lists.quantum-espresso.org/mailman/listinfo/developers
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.quantum-espresso.org/pipermail/developers/attachments/20250208/659a7a09/attachment-0001.html>
More information about the developers
mailing list