TABLE OF CONTENTS


ABINIT/m_wvl_denspot [ Modules ]

[ Top ] [ Modules ]

NAME

  m_wvl_denspot

FUNCTION

COPYRIGHT

  Copyright (C) 2008-2024 ABINIT group (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 .

SOURCE

16 #if defined HAVE_CONFIG_H
17 #include "config.h"
18 #endif
19 
20 #include "abi_common.h"
21 
22 module m_wvl_denspot
23 
24  use defs_basis
25  use m_errors
26  use m_abicore
27  use m_xmpi
28 
29  use defs_datatypes, only : pseudopotential_gth_type
30  use m_geometry,   only : xred2xcart
31 
32  implicit none
33 
34  private

ABINIT/wvl_denspot_free [ Functions ]

[ Top ] [ Functions ]

NAME

  wvl_denspot_free

FUNCTION

INPUTS

OUTPUT

SOURCE

170 subroutine wvl_denspot_free(den)
171 
172  use defs_wvltypes
173 #if defined HAVE_BIGDFT
174  use BigDFT_API, only: deallocate_rho_descriptors, &
175       & deallocate_denspot_distribution, denspot_free_history
176  use dynamic_memory
177 #endif
178  implicit none
179 
180 !Arguments ------------------------------------
181  type(wvl_denspot_type), intent(inout) :: den
182 
183 !Local variables-------------------------------
184 
185 ! *************************************************************************
186 
187 !DEBUG
188 !write (std_out,*) ' wvl_denspot_free : enter'
189 !ENDDEBUG
190 
191 #if defined HAVE_BIGDFT
192  if(associated(den%denspot%rhov)) then
193    call f_free_ptr(den%denspot%rhov)
194  end if
195  if(associated(den%denspot%rho_psi)) then
196    call f_free_ptr(den%denspot%rho_psi)
197  end if
198  if(associated(den%denspot%rho_C)) then
199    call f_free_ptr(den%denspot%rho_C)
200  end if
201  if(associated(den%denspot%V_ext)) then
202    call f_free_ptr(den%denspot%V_ext)
203  end if
204  if(associated(den%denspot%V_XC)) then
205    call f_free_ptr(den%denspot%V_XC)
206  end if
207  if(associated(den%denspot%Vloc_KS)) then
208    call f_free_ptr(den%denspot%Vloc_KS)
209  end if
210  if(associated(den%denspot%f_XC)) then
211    call f_free_ptr(den%denspot%f_XC)
212  end if
213  if(associated(den%denspot%rho_work)) then
214    call f_free_ptr(den%denspot%rho_work)
215  end if
216  if(associated(den%denspot%pot_work)) then
217    call f_free_ptr(den%denspot%pot_work)
218  end if
219  nullify(den%denspot%rhov)
220  nullify(den%denspot%rho_psi)
221  nullify(den%denspot%rho_C)
222  nullify(den%denspot%V_ext)
223  nullify(den%denspot%V_XC)
224  nullify(den%denspot%Vloc_KS)
225  nullify(den%denspot%f_XC)
226  nullify(den%denspot%rho_work)
227  nullify(den%denspot%pot_work)
228  !
229  call deallocate_rho_descriptors(den%denspot%rhod)
230  call deallocate_denspot_distribution(den%denspot%dpbox)
231  call denspot_free_history(den%denspot)
232 #else
233  BIGDFT_NOTENABLED_ERROR()
234  if (.false.) write(std_out,*) den%symObj
235 #endif
236 
237 end subroutine wvl_denspot_free

ABINIT/wvl_denspot_set [ Functions ]

[ Top ] [ Functions ]

NAME

  wvl_denspot_set

FUNCTION

  Fill in denspot datatype with information related
  to density and potential data.

INPUTS

  argin(sizein)=description

OUTPUT

  argout(sizeout)=description

SIDE EFFECTS

NOTES

SOURCE

 64 subroutine wvl_denspot_set(den,gth_params,ixc,natom,nsppol,rprimd,wvl,&
 65 &                          wvl_crmult,wvl_frmult,wvl_mpi_comm,xred)
 66 
 67  use defs_wvltypes
 68 
 69 #if defined HAVE_BIGDFT
 70  use BigDFT_API,only: initialize_DFT_local_fields,allocateRhoPot, &
 71 &                     input_variables,dpbox_set,density_descriptors
 72 #endif
 73  implicit none
 74 
 75 !Arguments ------------------------------------
 76   integer,intent(in):: ixc,natom,nsppol,wvl_mpi_comm
 77   real(dp), intent(in) :: rprimd(3, 3)
 78   real(dp), intent(in) :: wvl_frmult,wvl_crmult
 79   real(dp), intent(inout)  :: xred(3,natom)
 80   type(wvl_denspot_type), intent(out) :: den
 81   type(wvl_internal_type),intent(in)  :: wvl
 82   type(pseudopotential_gth_type),intent(in)::gth_params
 83 
 84 !Local variables-------------------------------
 85 #if defined HAVE_BIGDFT
 86   integer :: groupsize,me,nproc
 87   real(dp), allocatable :: xcart(:,:)
 88   character(len=3),parameter :: rho_commun='DBL'
 89   character(len=500) :: message
 90   character(len=4) :: SICapproach
 91   type(local_zone_descriptors) :: Lzd
 92 #endif
 93 
 94   ! *************************************************************************
 95 
 96 !DEBUG
 97 !write (std_out,*) ' wvl_denspot_set : enter'
 98 !ENDDEBUG
 99 
100 #if defined HAVE_BIGDFT
101 
102  write(message, '(a,a)' ) ch10,&
103 & ' wvl_denspot_set: Create wavelet type denspot.'
104  call wrtout(std_out,message,'COLL')
105 
106  nproc=xmpi_comm_size(wvl_mpi_comm)
107  me=xmpi_comm_rank(wvl_mpi_comm)
108  groupsize=0
109 
110 !Store xcart for each atom
111  ABI_MALLOC(xcart,(3, natom))
112  call xred2xcart(natom, rprimd, xcart, xred)
113 
114  call initialize_DFT_local_fields(den%denspot, ixc, nsppol)
115 
116 !number of planes for the density
117 !dpbox%nscatterarr(jproc, 1) = ngfft3_density
118 !number of planes for the potential
119 !dpbox%nscatterarr(jproc, 2) = ngfft3_potential
120 !starting offset for the potential
121 !dpbox%nscatterarr(jproc, 3) = density_start + potential_shift - 1
122 !GGA XC shift between density and potential
123 !dpbox%nscatterarr(jproc, 4) = potential_shift
124 
125  SICapproach="NONE"
126  Lzd%hgrids(1:3)=wvl%h(1:3)
127  Lzd%Glr%d%n1i=wvl%Glr%d%n1i
128  Lzd%Glr%d%n2i=wvl%Glr%d%n2i
129  Lzd%Glr%d%n3i=wvl%Glr%d%n3i
130  call dpbox_set(den%denspot%dpbox,Lzd,den%denspot%xc,me,nproc,wvl_mpi_comm,groupsize,&
131 & SICapproach,wvl%atoms%astruct%geocode,nsppol)
132 
133 !here dpbox can be put as input
134  call density_descriptors(me,nproc,den%denspot%xc,nsppol,wvl_crmult,wvl_frmult,wvl%atoms,&
135  den%denspot%dpbox,rho_commun,xcart,gth_params%radii_cf,den%denspot%rhod)
136 
137 !Note: change allocateRhoPot
138  call allocateRhoPot(wvl%Glr,nsppol,wvl%atoms,xcart,den%denspot)
139 
140 !Aditional information.
141  den%symObj = wvl%atoms%astruct%sym%symObj
142 
143  ABI_FREE(xcart)
144 
145 #else
146  BIGDFT_NOTENABLED_ERROR()
147  if (.false.) write(std_out,*) ixc,natom,nsppol,wvl_mpi_comm,rprimd(1,1),wvl_frmult,wvl_crmult,&
148 & xred(1,1),den%symObj,wvl%h(1),gth_params%psppar
149 #endif
150 
151 !DEBUG
152 !write (std_out,*) ' wvl_denspot_set : exit'
153 !ENDDEBUG
154 
155 end subroutine wvl_denspot_set