TABLE OF CONTENTS


ABINIT/m_gwls_LanczosBasis [ Modules ]

[ Top ] [ Modules ]

NAME

 m_gwls_LanczosBasis

FUNCTION

  .

COPYRIGHT

 Copyright (C) 2009-2024 ABINIT group (JLJ, BR, MC)
 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 
23 module m_gwls_LanczosBasis
24 !----------------------------------------------------------------------------------------------------
25 ! This module contains the static Lanczos basis, which should be computed once and for all,
26 ! and then be made available to other modules.
27 !----------------------------------------------------------------------------------------------------
28 ! local modules
29 use m_gwls_utility
30 use m_gwls_wf
31 use m_gwls_hamiltonian
32 use m_gwls_GenerateEpsilon
33 use m_dtset
34 
35 ! Abinit modules
36 use defs_basis
37 use m_abicore
38 
39 implicit none
40 save
41 private

m_hamiltonian/cleanup_Lanczos_basis [ Functions ]

[ Top ] [ m_hamiltonian ] [ Functions ]

NAME

  cleanup_Lanczos_basis

FUNCTION

  .

INPUTS

OUTPUT

SOURCE

123 subroutine cleanup_Lanczos_basis()
124 
125 
126 ! *************************************************************************
127 
128 ! if(allocated()) ABI_FREE() can cause a segfault if used in this form, without the THEN.
129 ! This is because ABI_MALLOC is expanded at compilation time in many statements; and the if() can only prevent the execution of
130 ! the first.
131 ! So, the deallocation takes place even if allocated() returns .false. without the THEN.
132 if (allocated(Lbasis_lanczos))  then
133   ABI_FREE(Lbasis_lanczos)
134 end if
135 if (allocated(Lbasis_model_lanczos))  then
136   ABI_FREE(Lbasis_model_lanczos)
137 end if
138 
139 end subroutine cleanup_Lanczos_basis

m_hamiltonian/modify_Lbasis_Coulomb [ Functions ]

[ Top ] [ m_hamiltonian ] [ Functions ]

NAME

  modify_Lbasis_Coulomb

FUNCTION

  .

INPUTS

OUTPUT

SOURCE

155 subroutine modify_Lbasis_Coulomb(psie_k, lmax, lmax_model)
156 !----------------------------------------------------------------------------------------------------
157 ! This subroutine computes, once and for all, the vectors (V^1/2.L)^*.psie
158 !----------------------------------------------------------------------------------------------------
159 real(dp),intent(in)  :: psie_k(2,npw_k)
160 integer ,intent(in)  :: lmax, lmax_model
161 
162 
163 real(dp), allocatable  :: psik_wrk(:,:), psikb_wrk(:,:), psikg_wrk(:,:)
164 real(dp), allocatable  :: psikg_e(:,:)
165 
166 integer  :: iblk_lanczos, nbdblock_lanczos
167 
168 integer   :: l, mb
169 
170 ! *************************************************************************
171 
172 ABI_MALLOC(psik_wrk  ,(2,npw_k))
173 ABI_MALLOC(psikb_wrk ,(2,npw_kb))
174 ABI_MALLOC(psikg_wrk ,(2,npw_g))
175 ABI_MALLOC(psikg_e ,(2,npw_g))
176 
177 ! copy the valence state on every row of FFT processors
178 do mb = 1, blocksize
179 psikb_wrk(:,(mb-1)*npw_k+1:mb*npw_k) = psie_k(:,:)
180 end do ! mb
181 ! Transform to FFT representation
182 call wf_block_distribute(psikb_wrk,  psikg_e,1) ! LA -> FFT
183 
184 
185 ! Number of blocks of lanczos vectors
186 nbdblock_lanczos = lmax/blocksize
187 if (modulo(lmax,blocksize) /= 0) nbdblock_lanczos = nbdblock_lanczos + 1
188 
189 ! loop on all blocks of lanczos vectors
190 do iblk_lanczos = 1, nbdblock_lanczos
191 
192 ! loop on all states within this block
193 do mb = 1, blocksize
194 ! Determine the index of the Lanczos vector
195 l = (iblk_lanczos-1)*blocksize + mb
196 
197 if ( l <= lmax) then
198   psik_wrk(1,:) = dble (Lbasis_lanczos(:,l))
199   psik_wrk(2,:) = dimag(Lbasis_lanczos(:,l))
200 else
201   psik_wrk(:,:) = zero
202 end if
203 ! Apply coulomb potential
204 call sqrt_vc_k(psik_wrk)
205 
206 ! Store in array of blocks of wavefunctions
207 psikb_wrk(:,(mb-1)*npw_k+1:mb*npw_k) = psik_wrk(:,:)
208 
209 end do ! mb
210 
211 ! Transform to FFT representation
212 call wf_block_distribute(psikb_wrk,  psikg_wrk,1) ! LA -> FFT
213 
214 ! fourier transform, conjugate
215 call g_to_r(psir1,psikg_wrk)
216 psir1(2,:,:,:) = -psir1(2,:,:,:)
217 
218 ! compute the product with the state "e"
219 call gr_to_g(psikg_wrk, psir1, psikg_e)
220 
221 ! return to LA configuration
222 
223 ! Transform to LA representation
224 call wf_block_distribute(psikb_wrk,  psikg_wrk, 2) ! FFT -> LA
225 
226 do mb = 1, blocksize
227 ! Determine the index of the Lanczos vector
228 l = (iblk_lanczos-1)*blocksize + mb
229 
230 if ( l <= lmax) then
231   psik_wrk(:,:) = psikb_wrk(:,(mb-1)*npw_k+1:mb*npw_k)
232   Lbasis_lanczos(:,l) = cmplx_1*psik_wrk(1,:)+cmplx_i*psik_wrk(2,:)
233 end if
234 
235 end do ! mb
236 
237 end do ! iblk_lanczos
238 
239 ! repeat, for model basis!
240 ! Number of blocks of lanczos vectors
241 nbdblock_lanczos = lmax_model/blocksize
242 if (modulo(lmax_model,blocksize) /= 0) nbdblock_lanczos = nbdblock_lanczos + 1
243 
244 ! loop on all blocks of lanczos vectors
245 do iblk_lanczos = 1, nbdblock_lanczos
246 
247 ! loop on all states within this block
248 do mb = 1, blocksize
249 ! Determine the index of the Lanczos vector
250 l = (iblk_lanczos-1)*blocksize + mb
251 
252 if ( l <= lmax_model) then
253   psik_wrk(1,:) = dble (Lbasis_model_lanczos(:,l))
254   psik_wrk(2,:) = dimag(Lbasis_model_lanczos(:,l))
255 else
256   psik_wrk(:,:) = zero
257 end if
258 ! Apply coulomb potential
259 call sqrt_vc_k(psik_wrk)
260 
261 ! Store in array of blocks of wavefunctions
262 psikb_wrk(:,(mb-1)*npw_k+1:mb*npw_k) = psik_wrk(:,:)
263 
264 end do ! mb
265 
266 ! Transform to FFT representation
267 call wf_block_distribute(psikb_wrk,  psikg_wrk,1) ! LA -> FFT
268 
269 ! fourier transform, conjugate
270 call g_to_r(psir1,psikg_wrk)
271 psir1(2,:,:,:) = -psir1(2,:,:,:)
272 
273 ! compute the product with the state "e"
274 call gr_to_g(psikg_wrk, psir1, psikg_e)
275 
276 ! return to LA configuration
277 ! Transform to LA representation
278 call wf_block_distribute(psikb_wrk,  psikg_wrk, 2) ! FFT -> LA
279 
280 do mb = 1, blocksize
281 ! Determine the index of the Lanczos vector
282 l = (iblk_lanczos-1)*blocksize + mb
283 
284 if ( l <= lmax_model) then
285   psik_wrk(:,:) = psikb_wrk(:,(mb-1)*npw_k+1:mb*npw_k)
286   Lbasis_model_lanczos(:,l) = cmplx_1*psik_wrk(1,:)+cmplx_i*psik_wrk(2,:)
287 end if
288 
289 end do ! mb
290 
291 end do
292 
293 ABI_FREE(psik_wrk  )
294 ABI_FREE(psikb_wrk )
295 ABI_FREE(psikg_wrk )
296 ABI_FREE(psikg_e )
297 
298 
299 end subroutine modify_Lbasis_Coulomb

m_hamiltonian/setup_Lanczos_basis [ Functions ]

[ Top ] [ m_hamiltonian ] [ Functions ]

NAME

  setup_Lanczos_basis

FUNCTION

  .

INPUTS

OUTPUT

SOURCE

 92 subroutine setup_Lanczos_basis(lmax,lmax_model)
 93 !----------------------------------------------------------------------------------------------------
 94 ! Set up the lanczos basis
 95 !----------------------------------------------------------------------------------------------------
 96 
 97 integer, intent(in) :: lmax, lmax_model
 98 
 99 ! *************************************************************************
100 
101 ABI_MALLOC(Lbasis_lanczos,(npw_k,lmax))
102 
103 if (lmax_model > 0) then
104   ABI_MALLOC(Lbasis_model_lanczos,(npw_k,lmax_model))
105 end if
106 
107 end subroutine setup_Lanczos_basis