TABLE OF CONTENTS


ABINIT/gwls_sternheimer [ Functions ]

[ Top ] [ Functions ]

NAME

 gwls_sternheimer

FUNCTION

 .

INPUTS

 dtset <type(dataset_type)>=all input variables in this dataset

OUTPUT

SOURCE

 49 subroutine gwls_sternheimer(acell_img,amu_img,codvsn,cpui,dtfil,dtset,etotal_img,fcart_img,&
 50 &                           gred_img,iexit,intgres_img,mixalch_img,mpi_enreg,nimage,npwtot,occ_img,&
 51 &                           pawang,pawrad,pawtab,psps,rprim_img,strten_img,vel_cell_img,vel_img,xred_img,&
 52 &                           filnam,filstat,idtset,jdtset,ndtset) ! optional arguments
 53 
 54 use m_gwls_utility,                   only : master_debug, files_status_new, files_status_old
 55 ! use m_gwls_wf,                        only : norm_k
 56 use m_gwls_hamiltonian,               only : destroy_H, exchange, g_to_r, eig
 57 use m_gwls_TimingLog
 58 use m_gwls_valenceWavefunctions
 59 use m_gwls_ComputeCorrelationEnergy
 60 use m_gwls_GenerateEpsilon
 61 use m_dtset
 62 
 63 use defs_basis
 64 use defs_wvltypes
 65 use m_pawang
 66 use m_pawrad
 67 use m_pawtab
 68 use m_abicore
 69 use m_errors
 70 use m_dtfil
 71 
 72 use defs_datatypes, only : pseudopotential_type
 73 use defs_abitypes, only : MPI_type
 74 use m_time,      only : timab
 75 use m_gstateimg, only : gstateimg
 76 
 77 !Arguments ------------------------------------
 78 !scalars
 79 integer,intent(in) :: nimage
 80 integer,optional,intent(in) :: idtset,ndtset
 81 integer,intent(inout) :: iexit
 82 real(dp),intent(in) :: cpui
 83 character(len=8),intent(in) :: codvsn
 84 character(len=fnlen),optional,intent(in) :: filstat
 85 type(MPI_type),intent(inout) :: mpi_enreg
 86 type(datafiles_type),target,intent(inout) :: dtfil
 87 type(dataset_type),intent(inout) :: dtset
 88 type(pawang_type),intent(inout) :: pawang
 89 type(pseudopotential_type),intent(inout) :: psps
 90 !arrays
 91 integer,optional,intent(in) :: jdtset(:)
 92 integer,intent(out) :: npwtot(dtset%nkpt)
 93 character(len=fnlen),optional,intent(in) :: filnam(:)
 94 real(dp), intent(out) :: etotal_img(nimage),fcart_img(3,dtset%natom,nimage)
 95 real(dp), intent(out) :: gred_img(3,dtset%natom,nimage),strten_img(6,nimage)
 96 real(dp), intent(out) :: intgres_img(dtset%nspden,dtset%natom,nimage)
 97 real(dp),intent(inout) :: acell_img(3,nimage),amu_img(dtset%ntypat,nimage),occ_img(dtset%mband*dtset%nkpt*dtset%nsppol,nimage)
 98 real(dp),intent(inout) :: mixalch_img(dtset%npspalch,dtset%ntypalch,nimage)
 99 real(dp),intent(inout) :: rprim_img(3,3,nimage),vel_cell_img(3,3,nimage),vel_img(3,dtset%natom,nimage)
100 real(dp),intent(inout) :: xred_img(3,dtset%natom,nimage)
101 type(pawrad_type),intent(inout) :: pawrad(psps%ntypat*psps%usepaw)
102 type(pawtab_type),intent(inout) :: pawtab(psps%ntypat*psps%usepaw)
103 !Local variables ------------------------------
104 
105 integer  :: e_index
106 real(dp) :: exchange_energy
107 real(dp) :: vxc_energy
108 real(dp) :: DFT_band_energy
109  type(wvl_data) :: wvl
110 
111 ! timing
112 real(dp) :: tsec(2)
113 integer :: GWLS_TIMAB, OPTION_TIMAB
114 !END DEBUG
115 ! *************************************************************************
116 
117 ! initialize to zero
118 
119 ! start new increment of the GWLS timer
120  GWLS_TIMAB   = 1501
121  OPTION_TIMAB = 1
122  call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
123 
124 
125  master_debug = .false.
126 !master_debug = .true.
127 
128 !Governs the status of all opened files formerly with the status 'new' harcoded.
129 !It is initialized in the m_gwls_utility module as 'new'; this default is overriden here if desired.
130  files_status_new = 'unknown'
131 !Governs the status of all opened files formerly with the status 'old' harcoded.
132 !It is initialized in the m_gwls_utility module as 'old'; this default is overriden here if desired.
133  files_status_old = 'unknown'
134 
135 ! Test the input to make sure it is consistent with what
136 ! the code can do
137 !call test_input(dtset2,psps2)
138 
139 !Initializing variables and build the Hamiltonian using build_H
140 
141  GWLS_TIMAB   = 1521
142  OPTION_TIMAB = 1
143  call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
144 
145 
146  call gstateimg(acell_img,amu_img,codvsn,cpui,dtfil,dtset,etotal_img,fcart_img,&
147 & gred_img,iexit,intgres_img,mixalch_img,mpi_enreg,nimage,npwtot,occ_img,&
148 & pawang,pawrad,pawtab,psps,rprim_img,strten_img,vel_cell_img,vel_img,wvl,xred_img,&
149 & filnam,filstat,idtset,jdtset,ndtset) ! optional arguments
150 
151  OPTION_TIMAB = 2
152  call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
153 
154 
155 ! tabulate the valence wavefunctions, to be used throughout the code!
156 ! TODO : put in build_H and destroy_H
157  GWLS_TIMAB   = 1522
158  OPTION_TIMAB = 1
159  call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
160 
161  call prepareValenceWavefunctions()
162 
163  OPTION_TIMAB = 2
164  call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
165 
166 
167 
168  e_index         = dtset%gwls_band_index
169  DFT_band_energy = eig(e_index)
170 
171  write(std_out,10) '                               '
172  write(std_out,10) ' GWLS RESULTS                  '
173  write(std_out,10) ' ----------------------------- '
174  write(std_out,12) ' For orbital |psi_e>, with e : ',e_index
175  write(std_out,14) ' DFT eigenenergy             : ',DFT_band_energy,' Ha = ',DFT_band_energy*Ha_eV,' eV'
176  flush(std_out)
177 
178  write(ab_out,10) '                               '
179  write(ab_out,10) ' GWLS RESULTS                  '
180  write(ab_out,10) ' ----------------------------- '
181  write(ab_out,12) ' For orbital |psi_e>, with e : ',e_index
182  write(ab_out,14) ' DFT eigenenergy             : ',DFT_band_energy,' Ha = ',DFT_band_energy*Ha_eV,' eV'
183 
184  if (dtset%gwls_exchange /= 0) then
185 
186    GWLS_TIMAB   = 1502
187    OPTION_TIMAB = 1
188    call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
189 
190    call compute_Exchange_and_Correlation_energies(e_index, exchange_energy, Vxc_energy)
191 
192    OPTION_TIMAB = 2
193    call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
194 
195 
196 
197 
198    write(std_out,14) ' <psi_e |   V_xc  | psi_e>   : ',vxc_energy     ,' Ha = ',vxc_energy     *Ha_eV,' eV'
199    write(std_out,14) ' <psi_e | Sigma_x | psi_e>   : ',exchange_energy,' Ha = ',exchange_energy*Ha_eV,' eV'
200    flush(std_out)
201 
202    write(ab_out,14) ' <psi_e |   V_xc  | psi_e>   : ',vxc_energy     ,' Ha = ',vxc_energy     *Ha_eV,' eV'
203    write(ab_out,14) ' <psi_e | Sigma_x | psi_e>   : ',exchange_energy,' Ha = ',exchange_energy*Ha_eV,' eV'
204  else
205    exchange_energy  = zero
206    Vxc_energy       = zero
207 
208  end if
209 
210  if (dtset%gwls_correlation == 3) then
211 
212    call setup_timing_log()
213 
214    GWLS_TIMAB   = 1503
215    OPTION_TIMAB = 1
216 
217    call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
218 
219    call compute_correlations_shift_lanczos(dtset, exchange_energy, Vxc_energy,master_debug)
220 
221    OPTION_TIMAB = 2
222    call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
223 
224 
225    call close_timing_log()
226 
227  else if (dtset%gwls_correlation == 4) then
228    call setup_timing_log()
229    call compute_correlations_no_model_shift_lanczos(dtset, exchange_energy, Vxc_energy,master_debug)
230    call close_timing_log()
231 
232  else if (dtset%gwls_correlation == 5) then
233    call setup_timing_log()
234    call Driver_GeneratePrintDielectricEigenvalues(dtset)
235    call close_timing_log()
236 
237 
238  end if
239 
240  call cleanupValenceWavefunctions()
241 
242  call destroy_H()
243 
244  GWLS_TIMAB   = 1501
245  OPTION_TIMAB = 2
246  call timab(GWLS_TIMAB,OPTION_TIMAB,tsec)
247 
248 
249  10 format(A)
250  12 format(A,I6)
251  14 format(A,ES24.16,A,F16.8,A)
252 
253 end subroutine gwls_sternheimer

ABINIT/m_gwls_sternheimer [ Modules ]

[ Top ] [ Modules ]

NAME

   m_gwls_sternheimer

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 module m_gwls_sternheimer
23 
24  implicit none
25 
26  private