Dear Dyer,
you are right... the factor of two you find is the difference
between the radius and the diameter of a sphere.
the bit of code in FFTXlib/src/fft_type.f90 I pointed to earlier
defines the ranges nr1,nr2,nr3 of a parallelepiped box that contains the
G-sphere of sqrt(gcutm) radius.
these are then passed to grid_set( dfft, bg, gcutm, dfft%nr1,
dfft%nr2, dfft%nr3 ) where the actual size is determined in three nested
loops that range from - the estimated limit to + the estimated limit
DO k = -nr3, nr3
DO j = -nr2, nr2
DO i = -nr1, nr1
where the maximum ACTUAL absolute index values are set as
nb(1),nb(2),nb(3), which should coincide with nr1,nr2,nr3 if they were
computed properly, and finally set the grid dimensions as
nr1 = 2 * nb(1) + 1
nr2 = 2 * nb(2) + 1
nr3 = 2 * nb(3) + 1
stefano
On 02/03/25 06:04, Dyer, Brock wrote:
I've been tracing down all the variables that are required to generate
a good FFT grid, and I seem to be off by a factor of 2 in the end. I
can't quite figure out what the issue is, but my guess may be that
while I am storing 'ecutwfc' in Rydbergs there is a conversion to
Hartrees somewhere in the code that I haven't seen yet.
My current process looks like this (with some values from a run I had
recently so I can compare):
ecutwfc = 100 Ry
Ecutrho = 400 Ry
# Unit cell dimensions, given in Bohrs
v1 = [44.09733757, 0.0, 0.0]
v2 = [0.0, 44.09733757, 0.0]
v3 = [0.0, 0.0, 44.09733757]
alat = sqrt(v1[0]**2 + v1[1]**2 + v1[2]**2)
tpiba = (2.0 * pi) / alat
gcutm = ecutrho / (tpiba**2)
at = [v1/alat, v2/alat, v3/alat]
nr1 = int(sqrt(gcutm) * sqrt(at[0][0]**2 + at[0][1]**2 + at[0][2]**2)) + 1
nr2 = int(sqrt(gcutm) * sqrt(at[1][0]**2 + at[1][1]**2 + at[1][2]**2)) + 1
nr3 = int(sqrt(gcutm) * sqrt(at[2][0]**2 + at[2][1]**2 + at[2][2]**2)) + 1
These last lines are where I've noticed the problem. From looking at
the output of my run with the given cell sizes, I expect an FFT grid
of 288x288x288, however if I were to run this code (and the code that
checks if it's a good size) I'd get an FFT grid that is only half
that. I'd love some advice on this if it is at all possible. I also
can send some more formatted code if it would help (I decided to cut
down the python code so it looked a bit more like the original f90 code).
------------------------------------------------------------------------
*From:* users on behalf of Stefano de Gironcoli
*Sent:* Thursday, February 27, 2025 2:50 PM
*To:* users@lists.quantum-espresso.org
*Subject:* Re: [QE-users] Details of Automatic FFT Grid Size Calculation
it's in SUBROUTINE realspace_grid_init in FFTXlib/src/file fft_types.f90
!
! ... calculate the size of the real-space dense grid for FFT
! ... first, an estimate of nr1,nr2,nr3, based on the max values
! ... of n_i indices in: G = i*b_1 + j*b_2 + k*b_3
! ... We use G*a_i = n_i => n_i .le. |Gmax||a_i|
!
dfft%nr1 = int ( sqrt (gcutm) * sqrt (at(1, 1)**2 + at(2,
1)**2 + at(3, 1)**2) ) + 1
dfft%nr2 = int ( sqrt (gcutm) * sqrt (at(1, 2)**2 + at(2,
2)**2 + at(3, 2)**2) ) + 1
dfft%nr3 = int ( sqrt (gcutm) * sqrt (at(1, 3)**2 + at(2,
3)**2 + at(3, 3)**2) ) + 1
stefano
On 27/02/25 20:11, Dyer, Brock wrote:
Hello all, I have been working quite a bit lately on automating my
QE workflow, and as part of that workflow I check the
automatically calculated FFT grid sizes for the level of theory
that I have been using in order to improve my parallelization.
I have tried tracing down and reading the code that calculates the
FFT grid sizes, however I cannot find/understand the actual code
to calculate the grid sizes. My current understanding is that the
initial parameters to calculate the grid size are 'ecutwfc' and/or
'ecutrho', and the unit cell size, and then there seems to be some
more math, and perhaps at the end the final dimensions get
calculated in 'fft_ggen.f90'.
What I am looking for ideally is a mathematical formula that
includes all of the input parameters and operations required to
calculate the FFT grid sizes so that I can implement it into my
workflow and hopefully not have to run double calculations to
properly parallelize.
Thanks in advance for the help!
Brock Dyer, Ursinus College Class of 2025
_______________________________________________________________________________
The Quantum ESPRESSO Foundation stands in solidarity with all
civilians worldwide who are victims of terrorism, military
aggression, and indiscriminate warfare.
--------------------------------------------------------------------------------
Quantum ESPRESSO is supported by MaX (www.max-centre.eu
<https://linkprotect.cudasvc.com/url?a=http%3a%2f%2fwww.max-centre.eu&c=E,1,o6mRyDmIfQXWXGr_JikolRfWqjEMwbk8xt5Q6J2l_N7Fky0gnXgD_aQU9TdYIxh51SmHmmgN9S5AJxHAS6vpqDSxTijqfOUXWuQTbchy3TZnaqcHjGvJ&typo=1>)
users mailing list users@lists.quantum-espresso.org
<mailto:users@lists.quantum-espresso.org>
https://lists.quantum-espresso.org/mailman/listinfo/users
<https://linkprotect.cudasvc.com/url?a=https%3a%2f%2flists.quantum-espresso.org%2fmailman%2flistinfo%2fusers&c=E,1,MVtAKpGoN3VwtTvpYzhgehNBw3lc6Ccvz18zJQYi8KwO2z3w1lLip_vRNXU25kTl-Lr0eMQKhqjTUVWShMvELzUclPM0hjrp4RgEFfOWL7-pZnGUAGdhpqG_ZBdI&typo=1>
_______________________________________________________________________________
The Quantum ESPRESSO Foundation stands in solidarity with all civilians
worldwide who are victims of terrorism, military aggression, and indiscriminate
warfare.
--------------------------------------------------------------------------------
Quantum ESPRESSO is supported by MaX (www.max-centre.eu)
users mailing listus...@lists.quantum-espresso.org
https://lists.quantum-espresso.org/mailman/listinfo/users
_______________________________________________________________________________
The Quantum ESPRESSO Foundation stands in solidarity with all civilians
worldwide who are victims of terrorism, military aggression, and indiscriminate
warfare.
--------------------------------------------------------------------------------
Quantum ESPRESSO is supported by MaX (www.max-centre.eu)
users mailing list users@lists.quantum-espresso.org
https://lists.quantum-espresso.org/mailman/listinfo/users