Octopus
xc_oep.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2012-2013 M. Gruning, P. Melo, M. Oliveira
3!! Copyright (C) 2022 N. Tancogne-Dejean
4!!
5!! This program is free software; you can redistribute it and/or modify
6!! it under the terms of the GNU General Public License as published by
7!! the Free Software Foundation; either version 2, or (at your option)
8!! any later version.
9!!
10!! This program is distributed in the hope that it will be useful,
11!! but WITHOUT ANY WARRANTY; without even the implied warranty of
12!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13!! GNU General Public License for more details.
14!!
15!! You should have received a copy of the GNU General Public License
16!! along with this program; if not, write to the Free Software
17!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
18!! 02110-1301, USA.
19!!
20
21#include "global.h"
22
23module xc_oep_oct_m
24 use comm_oct_m
25 use debug_oct_m
29 use global_oct_m
30 use grid_oct_m
33 use, intrinsic :: iso_fortran_env
40 use math_oct_m
42 use mesh_oct_m
44 use mpi_oct_m
47 use parser_oct_m
50 use scdm_oct_m
52 use space_oct_m
58 use xc_oct_m
59 use xc_f03_lib_m
61
62 implicit none
63
64 private
65 public :: &
66 xc_oep_t, &
68 xc_oep_end, &
76
78 integer, public, parameter :: &
79 OEP_LEVEL_NONE = 1, &
80 oep_level_kli = 3, &
82
84 integer, public, parameter :: &
85 OEP_MIXING_SCHEME_CONST = 1, &
88
90 integer, public, parameter :: &
91 OEP_TYPE_EXX = 1, &
92 oep_type_mgga = 2, &
93 oep_type_sic = 3, &
94 oep_type_photons = 4 ! one-photon OEP
95
96 type xc_oep_t
97 integer, public :: level
98 real(real64) :: mixing
99 integer :: mixing_scheme
100 type(lr_t) :: lr
101 type(linear_solver_t) :: solver
102 type(scf_tol_t) :: scftol
103 integer :: eigen_n
104 integer, allocatable :: eigen_type(:), eigen_index(:)
105 real(real64) :: socc, sfact
106 real(real64), allocatable, public :: vxc(:,:), uxc_bar(:,:,:)
107 real(real64), allocatable :: dlxc(:, :, :, :)
108 complex(real64), allocatable :: zlxc(:, :, :, :)
109 real(real64), public :: norm2ss
110 real(real64), allocatable :: vxc_old(:,:), ss_old(:,:)
111 integer :: noccst
112
113 type(scdm_t) :: scdm
114 logical :: scdm_for_pzsic
115
116 integer, public :: type = -1
117 end type xc_oep_t
118
119contains
120
121 ! ---------------------------------------------------------
122 subroutine xc_oep_init(oep, namespace, gr, st, mc, space, oep_type)
123 type(xc_oep_t), intent(inout) :: oep
124 type(namespace_t), intent(in) :: namespace
125 type(grid_t), intent(inout) :: gr
126 type(states_elec_t), intent(in) :: st
127 type(multicomm_t), intent(in) :: mc
128 class(space_t), intent(in) :: space
129 integer, intent(in) :: oep_type
130
131 push_sub(xc_oep_init)
132
133 !%Variable OEPLevel
134 !%Type integer
135 !%Default oep_kli
136 !%Section Hamiltonian::XC
137 !%Description
138 !% At what level shall <tt>Octopus</tt> handle the optimized effective potential (OEP) equation.
139 !%Option oep_none 1
140 !% Do not solve OEP equation.
141 !%Option oep_kli 3
142 !% Krieger-Li-Iafrate (KLI) approximation.
143 !% Ref: JB Krieger, Y Li, GJ Iafrate, <i>Phys. Lett. A</i> <b>146</b>, 256 (1990).
144 !%Option oep_full 5
145 !% (Experimental) Full solution of OEP equation using the Sternheimer approach.
146 !% The linear solver will be controlled by the variables in section <tt>Linear Response::Solver</tt>,
147 !% and the iterations for OEP by <tt>Linear Response::SCF in LR calculations</tt> and variable
148 !% <tt>OEPMixing</tt>. Note that default for <tt>LRMaximumIter</tt> is set to 10.
149 !% Ref: S. Kuemmel and J. Perdew, <i>Phys. Rev. Lett.</i> <b>90</b>, 043004 (2003).
150 !%End
151 call messages_obsolete_variable(namespace, 'OEP_Level', 'OEPLevel')
152 call parse_variable(namespace, 'OEPLevel', oep_level_kli, oep%level)
153 if (.not. varinfo_valid_option('OEPLevel', oep%level)) call messages_input_error(namespace, 'OEPLevel')
154
155 if (oep%level == oep_level_none) then
156 pop_sub(xc_oep_init)
157 return
158 end if
159
160 oep%type = oep_type
161
162 if (oep%level == oep_level_full) then
163
164 if (st%d%ispin == spinors) then
165 call messages_not_implemented("Full OEP with spinors", namespace=namespace)
166 end if
167
168 call messages_experimental("Full OEP")
169
170 !%Variable OEPMixing
171 !%Type float
172 !%Default 1.0
173 !%Section Hamiltonian::XC
174 !%Description
175 !% The linear mixing factor used to solve the Sternheimer
176 !% equation in the full OEP procedure.
177 !%End
178 call messages_obsolete_variable(namespace, 'OEP_Mixing', 'OEPMixing')
179 call parse_variable(namespace, 'OEPMixing', m_one, oep%mixing)
180
181 !%Variable OEPMixingScheme
182 !%Type integer
183 !%Default oep_mixing_scheme_const
184 !%Section Hamiltonian::XC
185 !%Description
186 !% Different Mixing Schemes are possible
187 !%Option oep_mixing_scheme_const 1
188 !% Use a constant
189 !% Reference: S. Kuemmel and J. Perdew, <i>Phys. Rev. Lett.</i> <b>90</b>, 4, 043004 (2003)
190 !%Option oep_mixing_scheme_bb 2
191 !% Use the Barzilai-Borwein (BB) Method
192 !% Reference: T. W. Hollins, S. J. Clark, K. Refson, and N. I. Gidopoulos,
193 !% <i>Phys. Rev. B</i> <b>85</b>, 235126 (2012)
194 !%Option oep_mixing_scheme_dens 3
195 !% Use the inverse of the electron density
196 !% Reference: S. Kuemmel and J. Perdew, <i>Phys. Rev. B</i> <b>68</b>, 035103 (2003)
197 !%End
198 call parse_variable(namespace, 'OEPMixingScheme', oep_mixing_scheme_const, oep%mixing_scheme)
200 if (oep%mixing_scheme == oep_mixing_scheme_bb) then
201 safe_allocate(oep%vxc_old(1:gr%np,st%d%ispin))
202 safe_allocate(oep%ss_old(1:gr%np,st%d%ispin))
203 oep%vxc_old = m_zero
204 oep%ss_old = m_zero
205 end if
207 oep%norm2ss = m_zero
208 end if
210 ! obtain the spin factors
211 call xc_oep_spinfactor(oep, st%d%nspin)
212
213 ! This variable will keep vxc across iterations
214 if ((st%d%ispin == spinors) .or. oep%level == oep_level_full) then
215 safe_allocate(oep%vxc(1:gr%np,st%d%nspin))
216 else
217 safe_allocate(oep%vxc(1:gr%np,1:min(st%d%nspin, 2)))
218 end if
219 oep%vxc = m_zero
220
221 ! when performing full OEP, we need to solve a linear equation
222 if (oep%level == oep_level_full) then
223 call scf_tol_init(oep%scftol, namespace, st%qtot, def_maximumiter=10)
224 call linear_solver_init(oep%solver, namespace, gr, states_are_real(st), mc, space)
225 call lr_init(oep%lr)
226 end if
227
228 if (oep%type == oep_type_sic) then
229 !%Variable SCDMforPZSIC
230 !%Type logical
231 !%Default .false.
232 !%Section Hamiltonian::XC
233 !%Description
234 !% If set to .true., the code will use the SCDM method to Wannierize the orbitals before
235 !% computing the orbital-dependent SIC correction
236 !%End
237 call parse_variable(namespace, 'SCDMforPZSIC', .false., oep%scdm_for_pzsic)
238 if (oep%scdm_for_pzsic) call scdm_init(oep%scdm, namespace)
239 end if
240
241 ! the linear equation has to be more converged if we are to attain the required precision
242 !oep%lr%conv_abs_dens = oep%lr%conv_abs_dens / (oep%mixing)
243
244 if (st%d%nspin == spinors) then
245 call messages_experimental("OEP with spinors")
246 end if
247
248 if (st%d%kpt%parallel .and. oep%level == oep_level_full) then
249 call messages_experimental("OEP Full parallel in spin/k-points", namespace=namespace)
250 end if
251
252 if (st%d%kpt%parallel .and. oep%type == oep_type_sic) then
253 call messages_not_implemented("PZ-SIC with k-point parallelization")
254 end if
255
256 pop_sub(xc_oep_init)
257 end subroutine xc_oep_init
258
259
260 ! ---------------------------------------------------------
261 subroutine xc_oep_end(oep)
262 class(xc_oep_t), intent(inout) :: oep
263
264 push_sub(xc_oep_end)
265
266 if (oep%level /= oep_level_none) then
267 safe_deallocate_a(oep%vxc)
268 if (oep%level == oep_level_full .or. oep%type == oep_type_photons) then
269 call lr_dealloc(oep%lr)
270 call linear_solver_end(oep%solver)
271 end if
272 if (oep%level == oep_level_full .and. oep%mixing_scheme == oep_mixing_scheme_bb) then
273 safe_deallocate_a(oep%vxc_old)
274 safe_deallocate_a(oep%ss_old)
275 end if
276 end if
277
278 pop_sub(xc_oep_end)
279 end subroutine xc_oep_end
280
281
282 ! ---------------------------------------------------------
283 subroutine xc_oep_write_info(oep, iunit, namespace)
284 type(xc_oep_t), intent(in) :: oep
285 integer, optional, intent(in) :: iunit
286 type(namespace_t), optional, intent(in) :: namespace
287
288 if (oep%level == oep_level_none) return
289
290 push_sub(xc_oep_write_info)
291
292 call messages_print_var_option('OEPLevel', oep%level, iunit=iunit, namespace=namespace)
293
294 pop_sub(xc_oep_write_info)
295 end subroutine xc_oep_write_info
296
297
298 ! ---------------------------------------------------------
300 ! ---------------------------------------------------------
301 subroutine xc_oep_spinfactor(oep, nspin)
302 class(xc_oep_t), intent(inout) :: oep
303 integer, intent(in) :: nspin
304
305 push_sub(xc_oep_spinfactor)
306
307 select case (nspin)
308 case (1) ! we need to correct or the spin occupancies
309 oep%socc = m_half
310 oep%sfact = m_two
311 case (2, 4)
312 oep%socc = m_one
313 oep%sfact = m_one
314 case default ! cannot handle any other case
315 assert(.false.)
316 end select
317
318 pop_sub(xc_oep_spinfactor)
319 end subroutine xc_oep_spinfactor
320
321
322 ! ---------------------------------------------------------
323 subroutine xc_oep_analyzeeigen(oep, st, is)
324 class(xc_oep_t), intent(inout) :: oep
325 type(states_elec_t), intent(in) :: st
326 integer, intent(in) :: is
327
328 integer :: ist
329 real(real64) :: max_eigen
330
331 push_sub(xc_oep_analyzeeigen)
332
333 ! The arrays are reuse for different spins, we reset them
334 oep%eigen_type = -1
335 oep%eigen_index = -1
336 oep%noccst = 0
337 oep%eigen_n = 0
338
339 ! find out the top occupied state, to correct for the asymptotics
340 ! of the potential
341 max_eigen = minval(st%eigenval(:, is))
342 do ist = 1, st%nst
343 if ((st%occ(ist, is) > m_min_occ).and.(st%eigenval(ist, is) > max_eigen)) then
344 max_eigen = st%eigenval(ist, is)
345 end if
346 end do
347
348 do ist = 1, st%nst
349 ! criterion for degeneracy
350 if (abs(st%eigenval(ist, is)-max_eigen) <= 1e-3_real64) then
351 oep%eigen_type(ist) = 2
352 else
353 if (st%occ(ist, is) > m_min_occ) then
354 oep%eigen_type(ist) = 1
355 oep%eigen_index(oep%eigen_n+1) = ist
356 oep%eigen_n = oep%eigen_n + 1
357 else
358 oep%eigen_type(ist) = 0
359 end if
360 end if
361 end do
362
363 ! find how many states are occupied.
364 do ist = 1, st%nst
365 if (st%occ(ist, is) > m_min_occ) oep%noccst = oep%noccst + 1
366 end do
367
368 ! We check that all states have been assigned
369 assert(all(oep%eigen_type >= 0))
370
371 pop_sub(xc_oep_analyzeeigen)
372 end subroutine xc_oep_analyzeeigen
373
374
375#include "xc_kli_pauli_inc.F90"
376#include "xc_oep_sic_pauli_inc.F90"
377
378#include "undef.F90"
379#include "real.F90"
380#include "xc_kli_inc.F90"
381#include "xc_oep_sic_inc.F90"
382#include "xc_oep_inc.F90"
383
384#include "undef.F90"
385#include "complex.F90"
386#include "xc_kli_inc.F90"
387#include "xc_oep_sic_inc.F90"
388#include "xc_oep_inc.F90"
389
390end module xc_oep_oct_m
391
392!! Local Variables:
393!! mode: f90
394!! coding: utf-8
395!! End:
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public spinors
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_min_occ
Minimal occupation that is considered to be non-zero.
Definition: global.F90:211
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public lr_init(lr)
subroutine, public lr_dealloc(lr)
subroutine, public linear_solver_end(this)
subroutine, public linear_solver_init(this, namespace, gr, states_are_real, mc, space)
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1046
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1086
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
This module provides routines for perform the Selected Column of Density Matrix (SCDM) method This fo...
Definition: scdm.F90:117
subroutine, public scdm_init(this, namespace)
Definition: scdm.F90:154
subroutine, public scf_tol_init(this, namespace, qtot, def_maximumiter, tol_scheme)
Definition: scf_tol.F90:160
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
Definition: xc.F90:114
integer, parameter, public oep_type_mgga
Definition: xc_oep.F90:183
subroutine, public xc_oep_spinfactor(oep, nspin)
A couple of auxiliary functions for oep.
Definition: xc_oep.F90:395
subroutine, public xc_oep_end(oep)
Definition: xc_oep.F90:355
subroutine, public zxc_oep_calc(oep, namespace, xcs, gr, hm, st, space, rcell_volume, ex, ec, vxc)
This file handles the evaluation of the OEP potential, in the KLI or full OEP as described in S....
Definition: xc_oep.F90:2311
integer, parameter, public oep_type_photons
Definition: xc_oep.F90:183
subroutine, public xc_oep_analyzeeigen(oep, st, is)
Definition: xc_oep.F90:417
integer, parameter, public oep_mixing_scheme_dens
Definition: xc_oep.F90:177
subroutine, public dxc_oep_calc(oep, namespace, xcs, gr, hm, st, space, rcell_volume, ex, ec, vxc)
This file handles the evaluation of the OEP potential, in the KLI or full OEP as described in S....
Definition: xc_oep.F90:1451
integer, parameter, public oep_mixing_scheme_bb
Definition: xc_oep.F90:177
subroutine, public xc_oep_write_info(oep, iunit, namespace)
Definition: xc_oep.F90:377
subroutine, public dxc_oep_mix(oep, mesh, ss, rho, is)
A routine that takes care of mixing the potential.
Definition: xc_oep.F90:1795
integer, parameter, public oep_level_full
Definition: xc_oep.F90:171
subroutine, public xc_oep_init(oep, namespace, gr, st, mc, space, oep_type)
Definition: xc_oep.F90:216
subroutine, public zxc_oep_mix(oep, mesh, ss, rho, is)
A routine that takes care of mixing the potential.
Definition: xc_oep.F90:2655
integer, parameter, public oep_type_sic
Definition: xc_oep.F90:183
integer, parameter, public oep_level_kli
Definition: xc_oep.F90:171
The states_elec_t class contains all electronic wave functions.