Octopus
poisson_fft.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2011 M. Marques, A. Castro, A. Rubio, G. Bertsch, M. Oliveira
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
23 use cube_oct_m
24 use debug_oct_m
25 use fft_oct_m
27 use global_oct_m
28 use, intrinsic :: iso_fortran_env
31 use math_oct_m
33 use mesh_oct_m
36 use parser_oct_m
39 use space_oct_m
42 use unit_oct_m
44
45 implicit none
46 private
47 public :: &
54
55 integer, public, parameter :: &
56 POISSON_FFT_KERNEL_NONE = -1, &
63
64 type poisson_fft_t
65 ! Components are public by default
66 type(fourier_space_op_t) :: coulb
67 integer :: kernel
68 real(real64) :: soft_coulb_param
69 end type poisson_fft_t
70
71 real(real64), parameter :: TOL_VANISHING_Q = 1e-6_real64
72contains
73
74 subroutine poisson_fft_init(this, namespace, space, cube, kernel, soft_coulb_param, fullcube)
75 type(poisson_fft_t), intent(out) :: this
76 type(namespace_t), intent(in) :: namespace
77 class(space_t), intent(in) :: space
78 type(cube_t), intent(inout) :: cube
79 integer, intent(in) :: kernel
80 real(real64), optional, intent(in) :: soft_coulb_param
81 type(cube_t), optional, intent(in) :: fullcube
82
83 push_sub(poisson_fft_init)
84
85 this%kernel = kernel
86 this%soft_coulb_param = optional_default(soft_coulb_param, m_zero)
87
88 safe_allocate(this%coulb%qq(1:space%dim))
89 this%coulb%qq(1:space%periodic_dim) = m_zero
90 this%coulb%qq(space%periodic_dim+1:space%dim) = 1e-5_real64
91 this%coulb%singularity = m_zero
92 this%coulb%mu = m_zero
93 this%coulb%alpha = m_zero
94 this%coulb%beta = m_zero
95
96 call poisson_fft_get_kernel(namespace, space, cube, this%coulb, kernel, soft_coulb_param, fullcube)
97
98 pop_sub(poisson_fft_init)
99 end subroutine poisson_fft_init
100
101 subroutine poisson_fft_get_kernel(namespace, space, cube, coulb, kernel, soft_coulb_param, fullcube)
102 type(namespace_t), intent(in) :: namespace
103 class(space_t), intent(in) :: space
104 type(cube_t), intent(in) :: cube
105 type(fourier_space_op_t), intent(inout) :: coulb
106 integer, intent(in) :: kernel
107 real(real64), optional, intent(in) :: soft_coulb_param
108 type(cube_t), optional, intent(in) :: fullcube
109
110 push_sub(poisson_fft_get_kernel)
111
112 if (coulb%mu > m_epsilon) then
113 if (space%dim /= 3 .or. kernel /= poisson_fft_kernel_nocut) then
114 message(1) = "The screened Coulomb potential is only implemented in 3D for PoissonFFTKernel=fft_nocut."
115 call messages_fatal(1, namespace=namespace)
116 end if
117 end if
118
119
120 if (kernel == poisson_fft_kernel_hockney) then
121 if (.not. present(fullcube)) then
122 message(1) = "Hockney's FFT-kernel needs cube of full unit cell "
123 call messages_fatal(1, namespace=namespace)
124 else
125 if (.not. allocated(fullcube%fft)) then
126 message(1) = "Hockney's FFT-kernel needs PoissonSolver=fft"
127 call messages_fatal(1, namespace=namespace)
128 end if
129 end if
130 end if
131
132
133 select case (space%dim)
134 case (1)
135 assert(present(soft_coulb_param))
136 select case (kernel)
138 call poisson_fft_build_1d_0d(namespace, cube, coulb, soft_coulb_param)
140 call poisson_fft_build_1d_1d(cube, coulb, soft_coulb_param)
141 case default
142 message(1) = "Invalid Poisson FFT kernel for 1D."
143 call messages_fatal(1, namespace=namespace)
144 end select
145
146 case (2)
147 select case (kernel)
149 call poisson_fft_build_2d_0d(namespace, cube, coulb)
151 call poisson_fft_build_2d_1d(namespace, cube, coulb)
153 call poisson_fft_build_2d_2d(cube, coulb)
154 case default
155 message(1) = "Invalid Poisson FFT kernel for 2D."
156 call messages_fatal(1, namespace=namespace)
157 end select
158
159 case (3)
160 select case (kernel)
162 call poisson_fft_build_3d_0d(namespace, cube, kernel, coulb, space%is_periodic())
165 call poisson_fft_build_3d_1d(namespace, space, cube, coulb)
168 call poisson_fft_build_3d_2d(namespace, cube, coulb)
171 call poisson_fft_build_3d_3d(cube, coulb)
172
174 call poisson_fft_build_3d_3d_hockney(cube, coulb, fullcube)
175
176 case default
177 message(1) = "Invalid Poisson FFT kernel for 3D."
178 call messages_fatal(1, namespace=namespace)
179 end select
180 end select
181
183 end subroutine poisson_fft_get_kernel
184
185 !-----------------------------------------------------------------
186
187 subroutine get_cutoff(namespace, default_r_c, r_c)
188 type(namespace_t), intent(in) :: namespace
189 real(real64), intent(in) :: default_r_c
190 real(real64), intent(out) :: r_c
191
192 push_sub(get_cutoff)
193
194 call parse_variable(namespace, 'PoissonCutoffRadius', default_r_c, r_c, units_inp%length)
195
196 call messages_write('Info: Poisson Cutoff Radius =')
197 call messages_write(r_c, units = units_out%length, fmt = '(f6.1)')
198 call messages_info()
199
200 if (r_c > default_r_c + m_epsilon) then
201 call messages_write('Poisson cutoff radius is larger than cell size.', new_line = .true.)
202 call messages_write('You can see electrons in neighboring cell(s).')
203 call messages_warning()
204 end if
205
206 pop_sub(get_cutoff)
207 end subroutine get_cutoff
208
209
211 subroutine poisson_fft_build_3d_3d(cube, coulb)
212 type(cube_t), intent(in) :: cube
213 type(fourier_space_op_t), intent(inout) :: coulb
214
215 integer :: ix, iy, iz, ixx(3), db(3), n1, n2, n3, lx, ly, lz
216 real(real64) :: temp(3), modg2, inv_four_mu2, ecut
217 real(real64) :: gg(3)
218 real(real64) :: singularity_term
219 real(real64), allocatable :: fft_coulb_fs(:,:,:)
220
222
223 db(1:3) = cube%rs_n_global(1:3)
224
225 if (coulb%mu > m_epsilon) then
226 inv_four_mu2 = m_one/((m_two*coulb%mu)**2)
227 else
228 inv_four_mu2 = m_zero
229 end if
230
231 n1 = max(1, cube%fs_n(1))
232 n2 = max(1, cube%fs_n(2))
233 n3 = max(1, cube%fs_n(3))
234
235 ! Define q+G = 0 term
236 ! Screened short-range coulomb potential (erfc function)
237 if(coulb%mu > m_epsilon .and. abs(coulb%alpha) < m_epsilon) then
238 ! Analytical limit of 4pi/q^2 * (1-beta*exp(-|q|^2/4mu^2))
239 singularity_term = m_four * m_pi * inv_four_mu2 * coulb%beta
240 else ! We use the user-defined value of the singularity
241 ! Long-range screened singularity
242 if(abs(coulb%alpha) > m_epsilon) then
243 singularity_term = coulb%singularity*coulb%alpha + m_four*m_pi*inv_four_mu2 * coulb%beta
244 else
245 singularity_term = coulb%singularity
246 end if
247 ! 4pi/q^2 * alpha
248 end if
249
250 ! store the Fourier transform of the Coulomb interaction
251 safe_allocate(fft_coulb_fs(1:n1, 1:n2, 1:n3))
252 fft_coulb_fs = m_zero
253
254 temp(1:3) = m_two*m_pi/(db(1:3)*cube%spacing(1:3))
255
256 ecut = fft_get_ecut_from_box(cube%rs_n_global, cube%fs_istart, cube%latt, temp, 3, coulb%qq)
257
258 !$omp parallel do collapse(2) default(none) schedule(static) &
259 !$omp shared(cube, coulb, n3, n2, n1, db, temp, ecut, inv_four_mu2, singularity_term, fft_Coulb_FS) &
260 !$omp private(lz, ly, lx, iz, iy, ix, ixx, gg, modg2)
261 do lz = 1, n3
262 do ly = 1, n2
263 iz = cube%fs_istart(3) + lz - 1
264 ixx(3) = pad_feq(iz, db(3), .true.)
265 iy = cube%fs_istart(2) + ly - 1
266 ixx(2) = pad_feq(iy, db(2), .true.)
267 do lx = 1, n1
268 ix = cube%fs_istart(1) + lx - 1
269 ixx(1) = pad_feq(ix, db(1), .true.)
270
271 ! Integer G -> Cartesian (G + q)
272 call fft_gg_transform(ixx, temp, 3, cube%latt, coulb%qq, gg, modg2)
273
274 if(modg2 > m_two * ecut * 1.001_real64) cycle
275
276 !HH not very elegant
277 if (cube%fft%library.eq.fftlib_nfft) modg2=cube%Lfs(ix,1)**2 + cube%Lfs(iy,2)**2 + cube%Lfs(iz,3)**2
278
279 if (modg2 > tol_vanishing_q) then
280 !Screened coulomb potential (erfc function)
281 if (coulb%mu > m_epsilon) then
282 if(abs(coulb%alpha) > m_epsilon) then ! CAM
283 fft_coulb_fs(lx, ly, lz) = m_four * m_pi / modg2 * (coulb%alpha + coulb%beta * exp(-modg2*inv_four_mu2))
284 else ! purely short-range screened
285 fft_coulb_fs(lx, ly, lz) = m_four * m_pi / modg2 * coulb%beta * (m_one-exp(-modg2 * inv_four_mu2))
286 end if
287 else
288 if (abs(coulb%alpha) > m_epsilon) then ! global screened hybrids
289 fft_coulb_fs(lx, ly, lz) = m_four * m_pi / modg2 * coulb%alpha
290 else ! Bare interaction
291 fft_coulb_fs(lx, ly, lz) = m_four * m_pi / modg2
292 end if
293 end if
294 else ! This is the term q+G = 0
295 fft_coulb_fs(lx, ly, lz) = singularity_term
296 end if
297
298 end do
299 end do
300 end do
301
302 call dfourier_space_op_init(coulb, cube, fft_coulb_fs)
303
304 safe_deallocate_a(fft_coulb_fs)
307 !-----------------------------------------------------------------
308
309 !-----------------------------------------------------------------
314 subroutine poisson_fft_build_3d_3d_hockney(cube, coulb, fullcube)
315 type(cube_t), intent(in) :: cube
316 type(fourier_space_op_t), intent(inout) :: coulb
317 type(cube_t), intent(in) :: fullcube
318
319 integer :: ix, iy, iz, ixx(3), db(3), nfs(3), nrs(3), nfs_s(3), nrs_s(3), dnrs(3)
320 real(real64) :: temp(3), modg2, weight
321 real(real64) :: gg(3)
322 real(real64), allocatable :: fft_Coulb_small_RS(:,:,:)
323 real(real64), allocatable :: fft_Coulb_RS(:,:,:)
324 complex(real64), allocatable :: fft_Coulb_small_FS(:,:,:), fft_Coulb_FS(:,:,:)
325
327
328 assert(abs(coulb%mu) < m_epsilon)
329
330 ! dimensions of large boxes
331 nfs(1:3) = fullcube%fs_n_global(1:3)
332 nrs(1:3) = fullcube%rs_n_global(1:3)
333
334 safe_allocate(fft_coulb_fs(1:nfs(1),1:nfs(2),1:nfs(3)))
335 safe_allocate(fft_coulb_rs(1:nrs(1),1:nrs(2),1:nrs(3)))
336
337 ! dimensions of small boxes x_s
338 nfs_s(1:3) = cube%fs_n_global(1:3)
339 nrs_s(1:3) = cube%rs_n_global(1:3)
340
341 safe_allocate(fft_coulb_small_fs(1:nfs_s(1),1:nfs_s(2),1:nfs_s(3)))
342 safe_allocate(fft_coulb_small_rs(1:nrs_s(1),1:nrs_s(2),1:nrs_s(3)))
343
344 ! build full periodic Coulomb potenital in Fourier space
345 fft_coulb_fs = m_zero
346
347 db(1:3) = fullcube%rs_n_global(1:3)
348 temp(1:3) = m_two*m_pi/(db(1:3)*cube%spacing(1:3))
349
350 do iz = 1, nfs(3)
351 ixx(3) = pad_feq(iz, db(3), .true.)
352 do iy = 1, nfs(2)
353 ixx(2) = pad_feq(iy, db(2), .true.)
354 do ix = 1, nfs(1)
355 ixx(1) = pad_feq(ix, db(1), .true.)
356
357 call fft_gg_transform(ixx, temp, 3, cube%latt, coulb%qq, gg, modg2)
358
359 if (abs(modg2) > tol_vanishing_q) then
360 fft_coulb_fs(ix, iy, iz) = m_one/modg2
361 else
362 fft_coulb_fs(ix, iy, iz) = m_zero
363 end if
364 end do
365 end do
366 end do
367
368 ! Full range hybrids weight
369 weight = m_four*m_pi
370 if(coulb%alpha > m_epsilon) weight = weight * coulb%alpha
371
372 do iz = 1, nfs(3)
373 do iy = 1, nfs(2)
374 do ix = 1, nfs(1)
375 fft_coulb_fs(ix, iy, iz) = weight*fft_coulb_fs(ix, iy, iz)
376 end do
377 end do
378 end do
379
380 ! get periodic Coulomb potential in real space
381 call dfft_backward(fullcube%fft, fft_coulb_fs, fft_coulb_rs)
382
383 ! copy to small box by respecting this pattern
384 ! full periodic coulomb: |abc--------------------------xyz|
385 ! Hockney: |abcxyz|
386 dnrs = nrs - nrs_s
387
388 do iz = 1, nrs_s(3)
389 ixx(3) = iz
390 if (iz > nrs_s(3)/2+1) ixx(3) = ixx(3) + dnrs(3)
391 do iy = 1, nrs_s(2)
392 ixx(2) = iy
393 if (iy > nrs_s(2)/2+1) ixx(2) = ixx(2) + dnrs(2)
394 do ix = 1, nrs_s(1)
395 ixx(1) = ix
396 if (ix > nrs_s(1)/2+1) ixx(1) = ixx(1) + dnrs(1)
397 fft_coulb_small_rs(ix, iy, iz) = fft_coulb_rs(ixx(1),ixx(2),ixx(3))
398 end do
399 end do
400 end do
401 ! make Hockney kernel in Fourier space
402 call dfft_forward(cube%fft, fft_coulb_small_rs, fft_coulb_small_fs)
403 !dummy copy for type conversion
404 fft_coulb_small_rs(1:nfs_s(1),1:nfs_s(2),1:nfs_s(3)) = &
405 real( fft_Coulb_small_FS(1:nfs_s(1),1:nfs_s(2),1:nfs_s(3)), real64)
406
407
408 ! Restrict array to local part to support pfft
409 ! For FFTW this reduces simply to the full array
410 call dfourier_space_op_init(coulb, cube, &
411 fft_coulb_small_rs(cube%fs_istart(1):cube%fs_istart(1)+cube%fs_n(1), &
412 cube%fs_istart(2):cube%fs_istart(2)+cube%fs_n(2), &
413 cube%fs_istart(3):cube%fs_istart(3)+cube%fs_n(3)))
414
415 safe_deallocate_a(fft_coulb_fs)
416 safe_deallocate_a(fft_coulb_rs)
417 safe_deallocate_a(fft_coulb_small_fs)
418 safe_deallocate_a(fft_coulb_small_rs)
419
421
423
424 !-----------------------------------------------------------------
426 subroutine poisson_fft_build_3d_2d(namespace, cube, coulb)
427 type(namespace_t), intent(in) :: namespace
428 type(cube_t), intent(in) :: cube
429 type(fourier_space_op_t), intent(inout) :: coulb
430
431 integer :: ix, iy, iz, ixx(3), db(3)
432 integer :: lx, ly, lz, n1, n2, n3
433 real(real64) :: temp(3), modg2, ecut, weight
434 real(real64) :: gpar, gz, r_c, gg(3), default_r_c
435 real(real64), allocatable :: fft_coulb_FS(:,:,:)
436
438
439 db(1:3) = cube%rs_n_global(1:3)
440
441 assert(abs(coulb%mu) < m_epsilon)
442 ! Full range hybrids weight
443 weight = m_four*m_pi
444 if(coulb%alpha > m_epsilon) weight = weight * coulb%alpha
445
446
447 !%Variable PoissonCutoffRadius
448 !%Type float
449 !%Section Hamiltonian::Poisson
450 !%Description
451 !% When <tt>PoissonSolver = fft</tt> and <tt>PoissonFFTKernel</tt> is neither <tt>multipole_corrections</tt>
452 !% nor <tt>fft_nocut</tt>,
453 !% this variable controls the distance after which the electron-electron interaction goes to zero.
454 !% A warning will be written if the value is too large and will cause spurious interactions between images.
455 !% The default is half of the FFT box max dimension in a finite direction.
456 !%End
457
458 default_r_c = db(3)*cube%spacing(3)/m_two
459 call get_cutoff(namespace, default_r_c, r_c)
460
461 n1 = max(1, cube%fs_n(1))
462 n2 = max(1, cube%fs_n(2))
463 n3 = max(1, cube%fs_n(3))
464 ! store the Fourier transform of the Coulomb interaction
465 safe_allocate(fft_coulb_fs(1:n1, 1:n2, 1:n3))
466 fft_coulb_fs = m_zero
467
468 temp(1:3) = m_two*m_pi/(db(1:3)*cube%spacing(1:3))
469
470 ecut = fft_get_ecut_from_box(cube%rs_n_global, cube%fs_istart, cube%latt, temp, 2, coulb%qq)
471
472 do lz = 1, n3
473 iz = cube%fs_istart(3) + lz - 1
474 ixx(3) = pad_feq(iz, db(3), .true.)
475 do ly = 1, n2
476 iy = cube%fs_istart(2) + ly - 1
477 ixx(2) = pad_feq(iy, db(2), .true.)
478 do lx = 1, n1
479 ix = cube%fs_istart(1) + lx - 1
480 ixx(1) = pad_feq(ix, db(1), .true.)
481
482 call fft_gg_transform(ixx, temp, 2, cube%latt, coulb%qq, gg, modg2)
483
484 if(sum(gg(1:2)**2) > m_two*ecut*1.001_real64) cycle
485
486 if (abs(modg2) > tol_vanishing_q) then
487 gz = abs(gg(3))
488 gpar = hypot(gg(1), gg(2))
489 ! note: if gpar = 0, then modg2 = gz**2
490 fft_coulb_fs(lx, ly, lz) = poisson_cutoff_3d_2d(gpar,gz,r_c)/modg2
491 else
492 fft_coulb_fs(lx, ly, lz) = -m_half*r_c**2
493 end if
494 fft_coulb_fs(lx, ly, lz) = weight*fft_coulb_fs(lx, ly, lz)
495 end do
496 end do
497
498 end do
499
500 call dfourier_space_op_init(coulb, cube, fft_coulb_fs)
501
502 safe_deallocate_a(fft_coulb_fs)
504 end subroutine poisson_fft_build_3d_2d
505 !-----------------------------------------------------------------
506
507
508 !-----------------------------------------------------------------
510 subroutine poisson_fft_build_3d_1d(namespace, space, cube, coulb)
511 type(namespace_t), intent(in) :: namespace
512 class(space_t), intent(in) :: space
513 type(cube_t), intent(in) :: cube
514 type(fourier_space_op_t), intent(inout) :: coulb
515
516 type(spline_t) :: cylinder_cutoff_f
517 real(real64), allocatable :: x(:), y(:)
518 integer :: ix, iy, iz, ixx(3), db(3), k, ngp
519 integer :: lx, ly, lz, n1, n2, n3, lxx(3)
520 real(real64) :: temp(3), modg2, xmax, weight
521 real(real64) :: gperp, gx, gy, gz, r_c, gg(3), default_r_c
522 real(real64), allocatable :: fft_coulb_FS(:,:,:)
523
525
526 assert(abs(coulb%mu) < m_epsilon)
527 ! Full range hybrids weight
528 weight = m_four*m_pi
529 if(coulb%alpha > m_epsilon) weight = weight * coulb%alpha
530
531
532 db(1:3) = cube%rs_n_global(1:3)
533
534 default_r_c = maxval(db(2:3)*cube%spacing(2:3)/m_two)
535 call get_cutoff(namespace, default_r_c, r_c)
536
537 n1 = max(1, cube%fs_n(1))
538 n2 = max(1, cube%fs_n(2))
539 n3 = max(1, cube%fs_n(3))
540 ! store the Fourier transform of the Coulomb interaction
541 safe_allocate(fft_coulb_fs(1:n1, 1:n2, 1:n3))
542 fft_coulb_fs = m_zero
543
544 temp(1:3) = m_two*m_pi/(db(1:3)*cube%spacing(1:3))
545
546 if (.not. space%is_periodic()) then
547 ngp = 8*db(2)
548 safe_allocate(x(1:ngp))
549 safe_allocate(y(1:ngp))
550 end if
551
552 ! Note(Alex) This loop ordering results in bad memory access for fft_Coulb_FS
553 ! It should be refactored
554 do lx = 1, n1
555 ix = cube%fs_istart(1) + lx - 1
556 ixx(1) = pad_feq(ix, db(1), .true.)
557 lxx(1) = ixx(1) - cube%fs_istart(1) + 1
558 gx = temp(1)*ixx(1)
559
560 if (.not. space%is_periodic()) then
561 call spline_init(cylinder_cutoff_f)
562 xmax = norm2(temp(2:3)*db(2:3))/2
563 do k = 1, ngp
564 x(k) = (k-1)*(xmax/(ngp-1))
565 y(k) = poisson_cutoff_3d_1d_finite(gx, x(k), norm2(cube%latt%rlattice_primitive(:, 1)), &
566 maxval(norm2(cube%latt%rlattice_primitive(:, 2:3), dim=1)))
567 end do
568 call spline_fit(ngp, x, y, cylinder_cutoff_f, m_zero)
569 end if
570
571 do ly = 1, n2
572 iy = cube%fs_istart(2) + ly - 1
573 ixx(2) = pad_feq(iy, db(2), .true.)
574 lxx(2) = ixx(2) - cube%fs_istart(2) + 1
575 do lz = 1, n3
576 iz = cube%fs_istart(3) + lz - 1
577 ixx(3) = pad_feq(iz, db(3), .true.)
578 lxx(3) = ixx(3) - cube%fs_istart(3) + 1
579
580 call fft_gg_transform(ixx, temp, 1, cube%latt, coulb%qq, gg, modg2)
581
582 if (abs(modg2) > tol_vanishing_q) then
583 gperp = hypot(gg(2), gg(3))
584 if (space%periodic_dim == 1) then
585 if (gperp > r_c) then
586 fft_coulb_fs(lx, ly, lz) = m_zero
587 else
588 fft_coulb_fs(lx, ly, lz) = poisson_cutoff_3d_1d(abs(gx), gperp, r_c)/modg2
589 end if
590 else if (.not. space%is_periodic()) then
591 gy = gg(2)
592 gz = gg(3)
593 if ((gz >= m_zero) .and. (gy >= m_zero)) then
594 fft_coulb_fs(lx, ly, lz) = spline_eval(cylinder_cutoff_f, gperp)
595 end if
596 if ((gz >= m_zero) .and. (gy < m_zero)) then
597 fft_coulb_fs(lx, ly, lz) = fft_coulb_fs(lx, -lxx(2) + 1, lz)
598 end if
599 if ((gz < m_zero) .and. (gy >= m_zero)) then
600 fft_coulb_fs(lx, ly, lz) = fft_coulb_fs(lx, ly, -lxx(3) + 1)
601 end if
602 if ((gz < m_zero) .and. (gy < m_zero)) then
603 fft_coulb_fs(lx, ly, lz) = fft_coulb_fs(lx, -lxx(2) + 1, -lxx(3) + 1)
604 end if
605 end if
606
607 else
608 if (space%periodic_dim == 1) then
609 fft_coulb_fs(lx, ly, lz) = -(m_half*log(r_c) - m_fourth)*r_c**2
610 else if (.not. space%is_periodic()) then
611 fft_coulb_fs(lx, ly, lz) = poisson_cutoff_3d_1d_finite(m_zero, m_zero, &
612 norm2(cube%latt%rlattice_primitive(:, 1)), maxval(norm2(cube%latt%rlattice_primitive(:, 2:3), dim=1)))
613 end if
614
615 end if
616 fft_coulb_fs(lx, ly, lz) = weight*fft_coulb_fs(lx, ly, lz)
617 end do
618 end do
619
620 if (.not. space%is_periodic()) then
621 call spline_end(cylinder_cutoff_f)
622 end if
623 end do
624
625 call dfourier_space_op_init(coulb, cube, fft_coulb_fs)
626
627 safe_deallocate_a(fft_coulb_fs)
628 safe_deallocate_a(x)
629 safe_deallocate_a(y)
631 end subroutine poisson_fft_build_3d_1d
632 !-----------------------------------------------------------------
633
634
635 !-----------------------------------------------------------------
637 subroutine poisson_fft_build_3d_0d(namespace, cube, kernel, coulb, is_periodic)
638 type(namespace_t), intent(in) :: namespace
639 type(cube_t), intent(in) :: cube
640 integer, intent(in) :: kernel
641 type(fourier_space_op_t), intent(inout) :: coulb
642 logical, intent(in) :: is_periodic
643
644 integer :: ix, iy, iz, ixx(3), db(3), lx, ly, lz, n1, n2, n3
645 real(real64) :: temp(3), modg2, ecut, weight
646 real(real64) :: r_c, gg(3), default_r_c
647 real(real64), allocatable :: fft_coulb_FS(:,:,:)
648 real(real64) :: axis(3,3)
649
651
652 assert(abs(coulb%mu) < m_epsilon)
653 ! Full range hybrids weight
654 weight = m_four*m_pi
655 if(coulb%alpha > m_epsilon) weight = weight * coulb%alpha
656
657
658 db(1:3) = cube%rs_n_global(1:3)
659
660 if (kernel /= poisson_fft_kernel_corrected) then
661
662 ! This is the real-space cutoff
663 do ix = 1, 3
664 axis(:,ix) = cube%latt%rlattice_primitive(:, ix) * cube%spacing(ix) * db(ix) / m_two
665 end do
666
667 default_r_c = m_huge
668 do ix = 1, 3
669 iy = mod(ix, 3)+1
670 iz = mod(ix+1, 3)+1
671
672 ! For orthogonal cells, this is determined by the size of the cube
673 if (.not. cube%latt%nonorthogonal) then
674 temp(1:3) = axis(:, ix)
675 default_r_c = min(default_r_c, norm2(temp(1:3)))
676 else
677 ! At the moment, this codepath is only called for DFT+U submesh Poisson solver
678 !
679 ! For non-orthogonal cells, the relevant length is the distance between two planes of
680 ! the parallelepiped. This ensures that we draw a sphere that touches the borders of the box,
681 ! thus avoiding contribution from periodic replicas
682 ! This distance is given by the usual formula of the distance from a point to a plan
683 temp = dcross_product(axis(:, iy), axis(:, iz))
684 temp = temp / norm2(temp)
685 default_r_c = min(default_r_c, dot_product(temp, axis(:, ix)-axis(:, iy)))
686 end if
687 end do
688 call get_cutoff(namespace, default_r_c, r_c)
689 end if
690
691 n1 = max(1, cube%fs_n(1))
692 n2 = max(1, cube%fs_n(2))
693 n3 = max(1, cube%fs_n(3))
694
695 ! store the fourier transform of the Coulomb interaction
696 ! store only the relevant part if PFFT is used
697 safe_allocate(fft_coulb_fs(1:n1,1:n2,1:n3))
698 fft_coulb_fs = m_zero
699
700 temp(1:3) = m_two*m_pi/(db(1:3)*cube%spacing(1:3))
701
702 ecut = fft_get_ecut_from_box(cube%rs_n_global, cube%fs_istart, cube%latt, temp, 3, coulb%qq)
703
704 do lz = 1, n3
705 iz = cube%fs_istart(3) + lz - 1
706 ixx(3) = pad_feq(iz, db(3), .true.)
707 do ly = 1, n2
708 iy = cube%fs_istart(2) + ly - 1
709 ixx(2) = pad_feq(iy, db(2), .true.)
710 do lx = 1, n1
711 ix = cube%fs_istart(1) + lx - 1
712 ixx(1) = pad_feq(ix, db(1), .true.)
713
714 call fft_gg_transform(ixx, temp, 0, cube%latt, coulb%qq, gg, modg2)
715
716 !HH
717 if (cube%fft%library.eq.fftlib_nfft) then
718 modg2=cube%Lfs(ix,1)**2+cube%Lfs(iy,2)**2+cube%Lfs(iz,3)**2
719 r_c = default_r_c*m_two
720 end if
721
722 ! At the moment this is only done for periodic space, so for DFT+U Coulomb integrals
723 if(modg2 > m_two*ecut*1.001_real64 .and. is_periodic) cycle
724
725 if (abs(modg2) > tol_vanishing_q) then
726 select case (kernel)
728 fft_coulb_fs(lx, ly, lz) = weight*poisson_cutoff_3d_0d(sqrt(modg2),r_c)/modg2
730 fft_coulb_fs(lx, ly, lz) = weight/modg2
731 end select
732 else
733 select case (kernel)
735 fft_coulb_fs(lx, ly, lz) = weight*r_c**2/m_two
737 fft_coulb_fs(lx, ly, lz) = m_zero
738 end select
739 end if
740 end do
741 end do
742 end do
743
744 call dfourier_space_op_init(coulb, cube, fft_coulb_fs, in_device = (kernel /= poisson_fft_kernel_corrected))
745
746 safe_deallocate_a(fft_coulb_fs)
748 end subroutine poisson_fft_build_3d_0d
749 !-----------------------------------------------------------------
750
751
752 !-----------------------------------------------------------------
754 subroutine poisson_fft_build_2d_0d(namespace, cube, coulb)
755 type(namespace_t), intent(in) :: namespace
756 type(cube_t), intent(in) :: cube
757 type(fourier_space_op_t), intent(inout) :: coulb
758
759 type(spline_t) :: besselintf
760 integer :: i, ix, iy, ixx(2), db(2), npoints
761 real(real64) :: temp(2), vec, r_c, maxf, dk, default_r_c, weight
762 real(real64), allocatable :: x(:), y(:)
763 real(real64), allocatable :: fft_coulb_FS(:,:,:)
764
766
767 assert(abs(coulb%mu) < m_epsilon)
768 ! Full range hybrids weight
769 weight = m_one
770 if(coulb%alpha > m_epsilon) weight = weight * coulb%alpha
771
772 db(1:2) = cube%rs_n_global(1:2)
773
774 default_r_c = maxval(db(1:2)*cube%spacing(1:2)/m_two)
775 call get_cutoff(namespace, default_r_c, r_c)
776
777 call spline_init(besselintf)
778
779 ! store the fourier transform of the Coulomb interaction
780 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
781 fft_coulb_fs = m_zero
782 temp(1:2) = m_two*m_pi/(db(1:2)*cube%spacing(1:2))
783
784 maxf = r_c * norm2(temp(1:2)*db(1:2))/2
785 dk = 0.25_real64 ! This seems to be reasonable.
786 npoints = nint(maxf/dk)
787 safe_allocate(x(1:npoints))
788 safe_allocate(y(1:npoints))
789 x(1) = m_zero
790 y(1) = m_zero
791 do i = 2, npoints
792 x(i) = (i-1) * maxf / (npoints-1)
793 y(i) = y(i-1) + poisson_cutoff_2d_0d(x(i-1), x(i))
794 end do
795 call spline_fit(npoints, x, y, besselintf, m_zero)
796
797 do iy = 1, cube%fs_n_global(2)
798 ixx(2) = pad_feq(iy, db(2), .true.)
799 do ix = 1, cube%fs_n_global(1)
800 ixx(1) = pad_feq(ix, db(1), .true.)
801 vec = norm2(temp(1:2)*ixx(1:2))
802 ! extra check to avoid extrapolation which leads to an error in gsl
803 if (vec*r_c >= x(npoints)) then
804 fft_coulb_fs(ix, iy, 1) = weight * y(npoints)
805 else if (vec > m_zero) then
806 fft_coulb_fs(ix, iy, 1) = weight * (m_two * m_pi / vec) * spline_eval(besselintf, vec*r_c)
807 else
808 fft_coulb_fs(ix, iy, 1) = weight * m_two * m_pi * r_c
809 end if
810 end do
811 end do
812
813 call dfourier_space_op_init(coulb, cube, fft_coulb_fs)
814
815 safe_deallocate_a(fft_coulb_fs)
816 safe_deallocate_a(x)
817 safe_deallocate_a(y)
818 call spline_end(besselintf)
820 end subroutine poisson_fft_build_2d_0d
821 !-----------------------------------------------------------------
822
823
824 !-----------------------------------------------------------------
826 subroutine poisson_fft_build_2d_1d(namespace, cube, coulb)
827 type(namespace_t), intent(in) :: namespace
828 type(cube_t), intent(in) :: cube
829 type(fourier_space_op_t), intent(inout) :: coulb
830
831 integer :: ix, iy, ixx(2), db(2)
832 real(real64) :: temp(2), r_c, gx, gy, default_r_c, weight
833 real(real64), allocatable :: fft_coulb_FS(:,:,:)
834
836
837 assert(abs(coulb%mu) < m_epsilon)
838 ! Full range hybrids weight
839 weight = m_one
840 if(coulb%alpha > m_epsilon) weight = weight * coulb%alpha
841
842
843 db(1:2) = cube%rs_n_global(1:2)
844
845 default_r_c = db(2)*cube%spacing(2)/m_two
846 call get_cutoff(namespace, default_r_c, r_c)
847
848 ! store the fourier transform of the Coulomb interaction
849 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
850 fft_coulb_fs = m_zero
851 temp(1:2) = m_two*m_pi/(db(1:2)*cube%spacing(1:2))
852
853 ! First, the term ix = 0 => gx = 0.
854 fft_coulb_fs(1, 1, 1) = -m_four * r_c * (log(r_c)-m_one)
855 do iy = 2, cube%fs_n_global(2)
856 ixx(2) = pad_feq(iy, db(2), .true.)
857 gy = temp(2)*ixx(2)
858 fft_coulb_fs(1, iy, 1) = -m_four * poisson_cutoff_intcoslog(r_c, gy, m_one) * weight
859 end do
860
861 do iy = 1, db(2)
862 ixx(2) = pad_feq(iy, db(2), .true.)
863 gy = temp(2)*ixx(2)
864 do ix = 2, cube%fs_n_global(1)
865 ixx(1) = pad_feq(ix, db(1), .true.)
866 gx = temp(1)*ixx(1)
867 fft_coulb_fs(ix, iy, 1) = poisson_cutoff_2d_1d(gy, gx, r_c) * weight
868 end do
869 end do
870
871 call dfourier_space_op_init(coulb, cube, fft_coulb_fs)
872
873 safe_deallocate_a(fft_coulb_fs)
874
876 end subroutine poisson_fft_build_2d_1d
877 !-----------------------------------------------------------------
878
879
880 !-----------------------------------------------------------------
882 subroutine poisson_fft_build_2d_2d(cube, coulb)
883 type(cube_t), intent(in) :: cube
884 type(fourier_space_op_t), intent(inout) :: coulb
885
886 integer :: ix, iy, ixx(2), db(2)
887 real(real64) :: temp(2), vec, weight
888 real(real64), allocatable :: fft_coulb_FS(:,:,:)
889
891
892 assert(abs(coulb%mu) < m_epsilon)
893 ! Full range hybrids weight
894 weight = m_one
895 if(coulb%alpha > m_epsilon) weight = weight * coulb%alpha
896
897 db(1:2) = cube%rs_n_global(1:2)
898
899 ! store the fourier transform of the Coulomb interaction
900 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
901 fft_coulb_fs = m_zero
902 temp(1:2) = m_two*m_pi/(db(1:2)*cube%spacing(1:2))
903
904 do iy = 1, cube%fs_n_global(2)
905 ixx(2) = pad_feq(iy, db(2), .true.)
906 do ix = 1, cube%fs_n_global(1)
907 ixx(1) = pad_feq(ix, db(1), .true.)
908 vec = sqrt((temp(1) * ixx(1))**2 + (temp(2) * ixx(2))**2)
909 if (vec > m_zero) fft_coulb_fs(ix, iy, 1) = m_two * m_pi / vec * weight
910 end do
911 end do
912
913 call dfourier_space_op_init(coulb, cube, fft_coulb_fs)
914
915 safe_deallocate_a(fft_coulb_fs)
917 end subroutine poisson_fft_build_2d_2d
918 !-----------------------------------------------------------------
919
920
921 !-----------------------------------------------------------------
922 subroutine poisson_fft_build_1d_1d(cube, coulb, poisson_soft_coulomb_param)
923 type(cube_t), intent(in) :: cube
924 type(fourier_space_op_t), intent(inout) :: coulb
925 real(real64), intent(in) :: poisson_soft_coulomb_param
926
927 integer :: ix, ixx
928 real(real64) :: g
929 real(real64), allocatable :: fft_coulb_fs(:, :, :), weight
930
932
933 assert(abs(coulb%mu) < m_epsilon)
934 ! Full range hybrids weight
935 weight = m_one
936 if(coulb%alpha > m_epsilon) weight = weight * coulb%alpha
937
938 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
939 fft_coulb_fs = m_zero
940
941 ! Fourier transform of Soft Coulomb interaction.
942 do ix = 1, cube%fs_n_global(1)
943 ixx = pad_feq(ix, cube%rs_n_global(1), .true.)
944 g = (ixx + coulb%qq(1))*m_two*m_pi/abs(cube%latt%rlattice(1,1))
945 if (abs(g) > tol_vanishing_q) then
946 fft_coulb_fs(ix, 1, 1) = m_two * loct_bessel_k0(poisson_soft_coulomb_param*abs(g)) * weight
947 else
948 fft_coulb_fs(ix, 1, 1) = coulb%singularity * m_two * weight
949 end if
950 end do
951
952 call dfourier_space_op_init(coulb, cube, fft_coulb_fs)
953 safe_deallocate_a(fft_coulb_fs)
954
956 end subroutine poisson_fft_build_1d_1d
957 !-----------------------------------------------------------------
958
959
960 !-----------------------------------------------------------------
961 subroutine poisson_fft_build_1d_0d(namespace, cube, coulb, poisson_soft_coulomb_param)
962 type(namespace_t), intent(in) :: namespace
963 type(cube_t), intent(in) :: cube
964 type(fourier_space_op_t), intent(inout) :: coulb
965 real(real64), intent(in) :: poisson_soft_coulomb_param
966
967 integer :: box(1), ixx(1), ix
968 real(real64) :: temp(1), g, r_c, default_r_c, weight
969 real(real64), allocatable :: fft_coulb_fs(:, :, :)
970
972
973 assert(abs(coulb%mu) < m_epsilon)
974 ! Full range hybrids weight
975 weight = m_one
976 if(coulb%alpha > m_epsilon) weight = weight * coulb%alpha
978 box(1:1) = cube%rs_n_global(1:1)
979
980 default_r_c = box(1)*cube%spacing(1)/m_two
981 call get_cutoff(namespace, default_r_c, r_c)
982
983 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
984 fft_coulb_fs = m_zero
985 temp(1:1) = m_two*m_pi/(box(1:1)*cube%spacing(1:1))
986
987 ! Fourier transform of Soft Coulomb interaction.
988 do ix = 1, cube%fs_n_global(1)
989 ixx(1) = pad_feq(ix, box(1), .true.)
990 g = temp(1)*ixx(1)
991 fft_coulb_fs(ix, 1, 1) = poisson_cutoff_1d_0d(g, poisson_soft_coulomb_param, r_c) * weight
992 end do
993
994 call dfourier_space_op_init(coulb, cube, fft_coulb_fs)
995 safe_deallocate_a(fft_coulb_fs)
996
998 end subroutine poisson_fft_build_1d_0d
999 !-----------------------------------------------------------------
1000
1001
1002 !-----------------------------------------------------------------
1003 subroutine poisson_fft_end(this)
1004 type(poisson_fft_t), intent(inout) :: this
1005
1006 push_sub(poisson_fft_end)
1007
1008 call fourier_space_op_end(this%coulb)
1009
1010 pop_sub(poisson_fft_end)
1011 end subroutine poisson_fft_end
1012
1013#include "undef.F90"
1014#include "real.F90"
1015#include "poisson_fft_inc.F90"
1016#include "undef.F90"
1017#include "complex.F90"
1018#include "poisson_fft_inc.F90"
1019
1020
1021end module poisson_fft_oct_m
1022
1023!! Local Variables:
1024!! mode: f90
1025!! coding: utf-8
1026!! End:
Some operations may be done for one spline-function, or for an array of them.
Definition: splines.F90:181
double hypot(double __x, double __y) __attribute__((__nothrow__
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:120
real(real64) function, public fft_get_ecut_from_box(box_dim, fs_istart, latt, gspacing, periodic_dim, qq)
Given an fft box (fixed by the real-space grid), it returns the cutoff energy of the sphere that fits...
Definition: fft.F90:1017
pure integer function, public pad_feq(ii, nn, mode)
convert between array index and G-vector
Definition: fft.F90:876
integer, parameter, public fftlib_nfft
Definition: fft.F90:179
pure subroutine, public fft_gg_transform(gg_in, temp, periodic_dim, latt, qq, gg, modg2)
Convert FFT grid index into the Cartesian reciprocal-space vector .
Definition: fft.F90:972
subroutine, public fourier_space_op_end(this)
subroutine, public dfourier_space_op_init(this, cube, op, in_device)
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_huge
Definition: global.F90:209
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_fourth
Definition: global.F90:200
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
Definition: math.F90:1941
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
real(real64) function, public poisson_cutoff_3d_2d(p, z, r)
real(real64) function, public poisson_cutoff_3d_1d(x, p, rmax)
real(real64) function, public poisson_cutoff_3d_0d(x, r)
integer, parameter, public poisson_fft_kernel_hockney
subroutine poisson_fft_build_2d_0d(namespace, cube, coulb)
A. Castro et al., Phys. Rev. B 80, 033102 (2009)
subroutine poisson_fft_build_3d_1d(namespace, space, cube, coulb)
C. A. Rozzi et al., Phys. Rev. B 73, 205119 (2006), Table I.
subroutine poisson_fft_build_3d_2d(namespace, cube, coulb)
C. A. Rozzi et al., Phys. Rev. B 73, 205119 (2006), Table I.
integer, parameter, public poisson_fft_kernel_nocut
integer, parameter, public poisson_fft_kernel_cyl
subroutine poisson_fft_build_2d_1d(namespace, cube, coulb)
A. Castro et al., Phys. Rev. B 80, 033102 (2009)
subroutine poisson_fft_build_1d_0d(namespace, cube, coulb, poisson_soft_coulomb_param)
subroutine, public zpoisson_fft_solve(this, mesh, cube, pot, rho, mesh_cube_map, average_to_zero, kernel, sm)
subroutine poisson_fft_build_3d_0d(namespace, cube, kernel, coulb, is_periodic)
C. A. Rozzi et al., Phys. Rev. B 73, 205119 (2006), Table I.
subroutine get_cutoff(namespace, default_r_c, r_c)
subroutine poisson_fft_build_3d_3d(cube, coulb)
Compute the Coulomb kernel in reciprocal space, for a 3D FFT grid.
subroutine poisson_fft_build_1d_1d(cube, coulb, poisson_soft_coulomb_param)
subroutine, public poisson_fft_get_kernel(namespace, space, cube, coulb, kernel, soft_coulb_param, fullcube)
subroutine, public poisson_fft_end(this)
subroutine poisson_fft_build_2d_2d(cube, coulb)
A. Castro et al., Phys. Rev. B 80, 033102 (2009)
subroutine, public poisson_fft_init(this, namespace, space, cube, kernel, soft_coulb_param, fullcube)
integer, parameter, public poisson_fft_kernel_pla
integer, parameter, public poisson_fft_kernel_corrected
integer, parameter, public poisson_fft_kernel_sph
subroutine, public dpoisson_fft_solve(this, mesh, cube, pot, rho, mesh_cube_map, average_to_zero, kernel, sm)
subroutine poisson_fft_build_3d_3d_hockney(cube, coulb, fullcube)
Kernel for Hockneys algorithm that solves the poisson equation in a small box while respecting the pe...
subroutine, public spline_fit(nrc, rofi, ffit, spl, threshold)
Definition: splines.F90:415
real(real64) function, public spline_eval(spl, x)
Definition: splines.F90:443
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
the basic spline datatype
Definition: splines.F90:158
int true(void)