TABLE OF CONTENTS


ABINIT/mkgrid_fft [ Functions ]

[ Top ] [ Functions ]

NAME

  mkgrid_fft

FUNCTION

  It sets the grid of fft (or real space) points to be treated.

COPYRIGHT

  Copyright (C) 2013-2018 ABINIT group (T.Rangel, DC)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

OUTPUT

SIDE EFFECTS

NOTES

PARENTS

      mkcore_paw,mklocl_realspace

CHILDREN

      xred2xcart

SOURCE

31 #if defined HAVE_CONFIG_H
32 #include "config.h"
33 #endif
34 
35 subroutine mkgrid_fft(ffti3_local,fftn3_distrib,gridcart,nfft,ngfft,rprimd)
36 
37  use defs_basis
38 
39  use m_geometry,       only : xred2xcart
40 
41 !This section has been created automatically by the script Abilint (TD).
42 !Do not modify the following lines by hand.
43 #undef ABI_FUNC
44 #define ABI_FUNC 'mkgrid_fft'
45 !End of the abilint section
46 
47  implicit none
48 
49 !Arguments ------------------------------------
50  integer, intent(in) :: nfft
51  integer,intent(in) :: ngfft(18)
52  integer, dimension(*), intent(in) :: ffti3_local,fftn3_distrib
53  real(dp), dimension(3,nfft), intent(out) :: gridcart
54  real(dp),intent(in) :: rprimd(3,3)
55 
56 !Local variables-------------------------------
57  integer :: ind,i1,i2,i3,i3loc,me,nproc
58  integer :: n1,n2,n3
59  real(dp), dimension(3) :: coord
60  real(dp), dimension(3,nfft) :: gridred
61 !character(len=500) :: message                   ! to be uncommented, if needed
62 
63 ! *************************************************************************
64 
65  n1    = ngfft(1)
66  n2    = ngfft(2)
67  n3    = ngfft(3)
68  nproc = ngfft(10)
69  me    = ngfft(11)
70 
71  do i3 = 1, n3, 1
72    if(fftn3_distrib(i3) == me) then !MPI
73      i3loc=ffti3_local(i3)
74      coord(3) = real(i3 - 1, dp) / real(n3, dp)
75      do i2 = 1, n2, 1
76        coord(2) = real(i2 - 1, dp) / real(n2, dp)
77        do i1 = 1, n1, 1
78          ind=i1+(i2-1)*n1+(i3loc-1)*n1*n2
79          coord(1) = real(i1 - 1, dp) / real(n1, dp)
80          gridred(:, ind) = coord(:)
81        end do
82      end do
83    end if
84  end do
85  call xred2xcart(nfft, rprimd, gridcart, gridred)
86 
87 end subroutine mkgrid_fft