TABLE OF CONTENTS


ABINIT/psolver_kernel [ Functions ]

[ Top ] [ Functions ]

NAME

 psolver_kernel

FUNCTION

 Build, get or free the kernel matrix used by the Poisson solver to compute the
 the convolution between 1/r and rho. The kernel is a saved variable. If
 this routine is called for building while a kernel already exists, it is not
 recomputed if all parameters (grid step and data size) are unchanged. Otherwise
 the kernel is freed and recompute again. The build action has a returned variable
 which is a pointer on the kernel. The get action also returns the kernel, or
 NULL if none has been associated.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR, TRangel).
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  iaction=0 to free all kernel allocated array,
          1 to compute the kernel (parallel case),
          2 to get it (parallel case),
          3 to compute the kernel (sequential),
          4 to get the sequential kernel.
  rprimd(3,3)=dimensional primitive translations in real space (bohr)

OUTPUT

  kernel= associated kernel on build (iaction = 1) and get action (iaction = 2).

PARENTS

      gstate,mklocl_realspace,psolver_hartree,psolver_rhohxc
      wvl_wfsinp_reformat

CHILDREN

      deallocate_coulomb_operator,nullify_coulomb_operator,pkernel_set,wrtout

SOURCE

 42 #if defined HAVE_CONFIG_H
 43 #include "config.h"
 44 #endif
 45 
 46 #include "abi_common.h"
 47 
 48 subroutine psolver_kernel(hgrid, iaction,  icoulomb, &
 49 & iproc, kernel, mpi_comm, ngfft, nproc, nscforder)
 50   
 51  use defs_basis
 52  use m_profiling_abi
 53  use m_errors
 54 
 55 #if defined HAVE_BIGDFT
 56  use BigDFT_API,     only  : coulomb_operator,nullify_coulomb_operator, &
 57 &                            deallocate_coulomb_operator,mpi_environment
 58  use poisson_solver, only : pkernel_init,pkernel_set
 59 #else
 60  use defs_wvltypes,  only : coulomb_operator
 61 #endif
 62 
 63 !This section has been created automatically by the script Abilint (TD).
 64 !Do not modify the following lines by hand.
 65 #undef ABI_FUNC
 66 #define ABI_FUNC 'psolver_kernel'
 67  use interfaces_14_hidewrite
 68 !End of the abilint section
 69 
 70   implicit none
 71 
 72 !Arguments ------------------------------------
 73   !scalars
 74   integer,intent(in) :: iaction, icoulomb, mpi_comm, nscforder, iproc, nproc
 75   !arrays
 76   integer, intent(in) :: ngfft(3)
 77   type(coulomb_operator),intent(inout)::kernel
 78   real(dp),intent(in) :: hgrid(3)
 79 
 80 !Local variables-------------------------
 81 #if defined HAVE_BIGDFT
 82   !scalars
 83   integer,parameter :: igpu=0 !no GPUs
 84   !arrays
 85   integer, save :: kernel_scfOrder 
 86   integer, save :: kernel_icoulomb
 87   integer, save :: data_size(3) = (/ -2, -2, -2 /) 
 88   real(dp), save :: kernel_hgrid(3)   ! Grid step used when creating the kernel.
 89   character(len = 1) :: geocode
 90   character(len=500) :: message
 91   integer :: current_size(3)
 92   type(coulomb_operator),save :: pkernel, kernelseq 
 93   type(mpi_environment) :: mpi_env
 94 #endif
 95 
 96 ! *************************************************************************
 97 
 98 #if defined HAVE_BIGDFT
 99 
100  if (icoulomb == 0) then
101 !  The kernel is built with 'P'eriodic boundary counditions.
102    geocode = 'P'
103  else if (icoulomb == 1) then
104 !  The kernel is built with 'F'ree boundary counditions.
105    geocode = 'F'
106  else if (icoulomb == 2) then
107 !  The kernel is built with 'S'urface boundary counditions.
108    geocode = 'S'
109  end if
110  current_size(:) = ngfft(1:3)
111 
112 !Initialise kernel_array pointer.
113  if (maxval(data_size) == -2) then
114    call nullify_coulomb_operator(pkernel)
115    call nullify_coulomb_operator(kernelseq)
116  end if
117 
118 !If iaction == 0, we free the kernel.
119  if (iaction == 0) then
120    if (associated(pkernel%kernel)) then
121      write(message, "(A)") "Psolver_kernel() : deallocating pkernel..."
122      call wrtout(std_out, message,'COLL')
123      
124      call deallocate_coulomb_operator(pkernel)
125    end if
126    if (associated(kernelseq%kernel)) then
127      write(message, "(A)") "Psolver_kernel() : deallocating kernelseq..."
128      call wrtout(std_out, message,'COLL')
129 
130      call deallocate_coulomb_operator(kernelseq)
131    end if
132    data_size = (/ -1, -1, -1 /)
133    return
134  end if
135 
136 
137 !Action is build or get. We check the sizes before doing anything else.
138 
139 !!$!Get the size depending on wavelets calculations or not
140 !!$ if (dtset%usewvl == 0) then
141 !!$   hgrid(1) = rprimd(1, 1) / ngfft(1)
142 !!$   hgrid(2) = rprimd(2, 2) / ngfft(2)
143 !!$   hgrid(3) = rprimd(3, 3) / ngfft(3)
144 !!$
145 !!$ else
146 !!$   hgrid(:) = 0.5d0 * wvl%h(:)
147 !!$   current_size(1:3) = (/ wvl%Glr%d%n1i, wvl%Glr%d%n2i, wvl%Glr%d%n3i /)
148 !!$ end if
149 
150 !Compute a new kernel if grid size has changed or if the kernel
151 !has never been computed.
152  if ((iaction == 1 .and. .not. associated(pkernel%kernel))    .or. &
153 & (iaction == 3 .and. .not. associated(kernelseq%kernel)) .or. &
154 & kernel_icoulomb /= icoulomb  .or. &
155 & data_size(1)    /= current_size(1) .or. &
156 & data_size(2)    /= current_size(2) .or. &
157 & data_size(3)    /= current_size(3) .or. &
158 & kernel_hgrid(1) /= hgrid(1)        .or. &
159 & kernel_hgrid(2) /= hgrid(2)        .or. &
160 & kernel_hgrid(3) /= hgrid(3)        .or. &
161 & kernel_scfOrder /= nscforder) then
162    write(message, "(A,A,A,3I6)") "Psolver_kernel() : building kernel...", ch10, &
163 &   " | data dimensions:", current_size
164    call wrtout(std_out, message, 'COLL')
165 
166    if (iaction == 1 .or. iaction == 2) then
167      if (associated(pkernel%kernel)) then
168        call deallocate_coulomb_operator(pkernel)
169      end if
170      mpi_env%mpi_comm = mpi_comm
171      mpi_env%iproc    = iproc
172      mpi_env%nproc    = nproc
173      mpi_env%igroup   = 0 ! no task groups
174      mpi_env%ngroup   = 1 ! no task groups
175      pkernel= pkernel_init(.True.,iproc,nproc,igpu,geocode,&
176 &     current_size,hgrid,nscforder, mpi_env = mpi_env)
177      call pkernel_set(pkernel,.True.)
178    end if
179 
180    if (iaction == 3 .or. iaction == 4) then
181      if (associated(kernelseq%kernel)) then
182        call deallocate_coulomb_operator(kernelseq)
183      end if
184      mpi_env%mpi_comm = mpi_comm
185      mpi_env%iproc    = 0
186      mpi_env%nproc    = 1
187      mpi_env%igroup   = 0 ! no task groups
188      mpi_env%ngroup   = 1 ! no task groups
189      kernelseq= pkernel_init(.True.,iproc,nproc,igpu,geocode,&
190 &     current_size,hgrid,nscforder, mpi_env = mpi_env)
191      call pkernel_set(kernelseq,.True.)
192    end if
193 
194 !  Storing variables which were used to make the kernel
195    kernel_icoulomb = icoulomb
196    data_size(:)    = current_size(:)
197    kernel_hgrid(:) = hgrid(:)
198    kernel_scfOrder = nscforder
199  end if
200  
201  ! Shallow copy if kernel has been associated.
202  if (iaction == 1 .or. iaction == 2) then
203    kernel = pkernel
204  end if
205  if (iaction == 3 .or. iaction == 4) then
206    kernel = kernelseq
207  end if
208 
209 #else
210  BIGDFT_NOTENABLED_ERROR()
211  if (.false.) write(std_out,*)  iaction,icoulomb,mpi_comm,nscforder,iproc,nproc,ngfft(1),kernel%co,hgrid(1)
212 #endif
213 
214 end subroutine psolver_kernel