TABLE OF CONTENTS


ABINIT/gwls_sternheimer [ Functions ]

[ Top ] [ Functions ]

NAME

 gwls_sternheimer

FUNCTION

 .

COPYRIGHT

 Copyright (C) 1999-2018 ABINIT group (GMR, VO, LR, RWG, MG, RShaltaf)
 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

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

OUTPUT

PARENTS

      driver

CHILDREN

      cleanupvalencewavefunctions,close_timing_log
      compute_correlations_no_model_shift_lanczos
      compute_correlations_shift_lanczos
      compute_exchange_and_correlation_energies,destroy_h
      driver_generateprintdielectriceigenvalues,gstateimg
      preparevalencewavefunctions,setup_timing_log,timab

SOURCE

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