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
51 use space_oct_m
56 use scdm_oct_m
59 use xc_oct_m
60 use xc_f03_lib_m
62 use xc_vxc_oct_m
63
64 implicit none
65
66 private
67 public :: &
68 xc_oep_t, &
70 xc_oep_end, &
78
80 integer, public, parameter :: &
81 OEP_LEVEL_NONE = 1, &
82 oep_level_kli = 3, &
84
86 integer, public, parameter :: &
87 OEP_MIXING_SCHEME_CONST = 1, &
90
92 integer, public, parameter :: &
93 OEP_TYPE_EXX = 1, &
94 oep_type_mgga = 2, &
95 oep_type_sic = 3, &
96 oep_type_photons = 4 ! one-photon OEP
97
98 type xc_oep_t
99 integer, public :: level
100 real(real64) :: mixing
101 integer :: mixing_scheme
102 type(lr_t) :: lr
103 type(linear_solver_t) :: solver
104 type(scf_tol_t) :: scftol
105 integer :: eigen_n
106 integer, allocatable :: eigen_type(:), eigen_index(:)
107 real(real64) :: socc, sfact
108 real(real64), allocatable, public :: vxc(:,:), uxc_bar(:,:,:)
109 real(real64), allocatable :: dlxc(:, :, :, :)
110 complex(real64), allocatable :: zlxc(:, :, :, :)
111 real(real64), public :: norm2ss
112 real(real64), allocatable :: vxc_old(:,:), ss_old(:,:)
113 integer :: noccst
114
115 type(scdm_t) :: scdm
116 logical :: scdm_for_pzsic
117
118 integer, public :: type = -1
119 end type xc_oep_t
120
121contains
122
123 ! ---------------------------------------------------------
124 subroutine xc_oep_init(oep, namespace, gr, st, mc, space, oep_type)
125 type(xc_oep_t), intent(inout) :: oep
126 type(namespace_t), intent(in) :: namespace
127 type(grid_t), intent(inout) :: gr
128 type(states_elec_t), intent(in) :: st
129 type(multicomm_t), intent(in) :: mc
130 class(space_t), intent(in) :: space
131 integer, intent(in) :: oep_type
132
133 push_sub(xc_oep_init)
134
135 !%Variable OEPLevel
136 !%Type integer
137 !%Default oep_kli
138 !%Section Hamiltonian::XC
139 !%Description
140 !% At what level shall <tt>Octopus</tt> handle the optimized effective potential (OEP) equation.
141 !%Option oep_none 1
142 !% Do not solve OEP equation.
143 !%Option oep_kli 3
144 !% Krieger-Li-Iafrate (KLI) approximation.
145 !% Ref: JB Krieger, Y Li, GJ Iafrate, <i>Phys. Lett. A</i> <b>146</b>, 256 (1990).
146 !%Option oep_full 5
147 !% (Experimental) Full solution of OEP equation using the Sternheimer approach.
148 !% The linear solver will be controlled by the variables in section <tt>Linear Response::Solver</tt>,
149 !% and the iterations for OEP by <tt>Linear Response::SCF in LR calculations</tt> and variable
150 !% <tt>OEPMixing</tt>. Note that default for <tt>LRMaximumIter</tt> is set to 10.
151 !% Ref: S. Kuemmel and J. Perdew, <i>Phys. Rev. Lett.</i> <b>90</b>, 043004 (2003).
152 !%End
153 call messages_obsolete_variable(namespace, 'OEP_Level', 'OEPLevel')
154 call parse_variable(namespace, 'OEPLevel', oep_level_kli, oep%level)
155 if (.not. varinfo_valid_option('OEPLevel', oep%level)) call messages_input_error(namespace, 'OEPLevel')
156
157 if (oep%level == oep_level_none) then
158 pop_sub(xc_oep_init)
159 return
160 end if
161
162 oep%type = oep_type
163
164 if (oep%level == oep_level_full) then
165
166 if (st%d%nspin == spinors) then
167 call messages_not_implemented("Full OEP with spinors", namespace=namespace)
168 end if
169
170 call messages_experimental("Full OEP")
171
172 !%Variable OEPMixing
173 !%Type float
174 !%Default 1.0
175 !%Section Hamiltonian::XC
176 !%Description
177 !% The linear mixing factor used to solve the Sternheimer
178 !% equation in the full OEP procedure.
179 !%End
180 call messages_obsolete_variable(namespace, 'OEP_Mixing', 'OEPMixing')
181 call parse_variable(namespace, 'OEPMixing', m_one, oep%mixing)
182
183 !%Variable OEPMixingScheme
184 !%Type integer
185 !%Default oep_mixing_scheme_const
186 !%Section Hamiltonian::XC
187 !%Description
188 !% Different Mixing Schemes are possible
189 !%Option oep_mixing_scheme_const 1
190 !% Use a constant
191 !% Reference: S. Kuemmel and J. Perdew, <i>Phys. Rev. Lett.</i> <b>90</b>, 4, 043004 (2003)
192 !%Option oep_mixing_scheme_bb 2
193 !% Use the Barzilai-Borwein (BB) Method
194 !% Reference: T. W. Hollins, S. J. Clark, K. Refson, and N. I. Gidopoulos,
195 !% <i>Phys. Rev. B</i> <b>85</b>, 235126 (2012)
196 !%Option oep_mixing_scheme_dens 3
197 !% Use the inverse of the electron density
198 !% Reference: S. Kuemmel and J. Perdew, <i>Phys. Rev. B</i> <b>68</b>, 035103 (2003)
199 !%End
200 call parse_variable(namespace, 'OEPMixingScheme', oep_mixing_scheme_const, oep%mixing_scheme)
202 if (oep%mixing_scheme == oep_mixing_scheme_bb) then
203 safe_allocate(oep%vxc_old(1:gr%np,st%d%ispin))
204 safe_allocate(oep%ss_old(1:gr%np,st%d%ispin))
205 oep%vxc_old = m_zero
206 oep%ss_old = m_zero
207 end if
209 oep%norm2ss = m_zero
210 end if
212 ! obtain the spin factors
213 call xc_oep_spinfactor(oep, st%d%nspin)
214
215 ! This variable will keep vxc across iterations
216 if ((st%d%ispin == spinors) .or. oep%level == oep_level_full) then
217 safe_allocate(oep%vxc(1:gr%np,st%d%nspin))
218 else
219 safe_allocate(oep%vxc(1:gr%np,1:min(st%d%nspin, 2)))
220 end if
221 oep%vxc = m_zero
222
223 ! when performing full OEP, we need to solve a linear equation
224 if (oep%level == oep_level_full) then
225 call scf_tol_init(oep%scftol, namespace, st%qtot, def_maximumiter=10)
226 call linear_solver_init(oep%solver, namespace, gr, states_are_real(st), mc, space)
227 call lr_init(oep%lr)
228 end if
229
230 if (oep%type == oep_type_sic) then
231 !%Variable SCDMforPZSIC
232 !%Type logical
233 !%Default .false.
234 !%Section Hamiltonian::XC
235 !%Description
236 !% If set to .true., the code will use the SCDM method to Wannierize the orbitals before
237 !% computing the orbital-dependent SIC correction
238 !%End
239 call parse_variable(namespace, 'SCDMforPZSIC', .false., oep%scdm_for_pzsic)
240 if (oep%scdm_for_pzsic) call scdm_init(oep%scdm, namespace)
241 end if
242
243 ! the linear equation has to be more converged if we are to attain the required precision
244 !oep%lr%conv_abs_dens = oep%lr%conv_abs_dens / (oep%mixing)
245
246 if (st%d%nspin == spinors) then
247 call messages_experimental("OEP with spinors")
248 end if
249
250 if (st%d%kpt%parallel .and. oep%level == oep_level_full) then
251 call messages_experimental("OEP Full parallel in spin/k-points", namespace=namespace)
252 end if
253
254 if (st%d%kpt%parallel .and. oep%type == oep_type_sic) then
255 call messages_not_implemented("PZ-SIC with k-point parallelization")
256 end if
257
258 pop_sub(xc_oep_init)
259 end subroutine xc_oep_init
260
261
262 ! ---------------------------------------------------------
263 subroutine xc_oep_end(oep)
264 class(xc_oep_t), intent(inout) :: oep
265
266 push_sub(xc_oep_end)
267
268 if (oep%level /= oep_level_none) then
269 safe_deallocate_a(oep%vxc)
270 if (oep%level == oep_level_full .or. oep%type == oep_type_photons) then
271 call lr_dealloc(oep%lr)
272 call linear_solver_end(oep%solver)
273 end if
274 if (oep%level == oep_level_full .and. oep%mixing_scheme == oep_mixing_scheme_bb) then
275 safe_deallocate_a(oep%vxc_old)
276 safe_deallocate_a(oep%ss_old)
277 end if
278 end if
279
280 pop_sub(xc_oep_end)
281 end subroutine xc_oep_end
282
283
284 ! ---------------------------------------------------------
285 subroutine xc_oep_write_info(oep, iunit, namespace)
286 type(xc_oep_t), intent(in) :: oep
287 integer, optional, intent(in) :: iunit
288 type(namespace_t), optional, intent(in) :: namespace
289
290 if (oep%level == oep_level_none) return
291
292 push_sub(xc_oep_write_info)
293
294 call messages_print_var_option('OEPLevel', oep%level, iunit=iunit, namespace=namespace)
295
296 pop_sub(xc_oep_write_info)
297 end subroutine xc_oep_write_info
298
299
300 ! ---------------------------------------------------------
302 ! ---------------------------------------------------------
303 subroutine xc_oep_spinfactor(oep, nspin)
304 class(xc_oep_t), intent(inout) :: oep
305 integer, intent(in) :: nspin
306
307 push_sub(xc_oep_spinfactor)
308
309 select case (nspin)
310 case (1) ! we need to correct or the spin occupancies
311 oep%socc = m_half
312 oep%sfact = m_two
313 case (2, 4)
314 oep%socc = m_one
315 oep%sfact = m_one
316 case default ! cannot handle any other case
317 assert(.false.)
318 end select
319
320 pop_sub(xc_oep_spinfactor)
321 end subroutine xc_oep_spinfactor
322
323
324 ! ---------------------------------------------------------
325 subroutine xc_oep_analyzeeigen(oep, st, is)
326 class(xc_oep_t), intent(inout) :: oep
327 type(states_elec_t), intent(in) :: st
328 integer, intent(in) :: is
329
330 integer :: ist
331 real(real64) :: max_eigen
332
333 push_sub(xc_oep_analyzeeigen)
334
335 ! The arrays are reuse for different spins, we reset them
336 oep%eigen_type = -1
337 oep%eigen_index = -1
338 oep%noccst = 0
339 oep%eigen_n = 0
340
341 ! find out the top occupied state, to correct for the asymptotics
342 ! of the potential
343 max_eigen = minval(st%eigenval(:, is))
344 do ist = 1, st%nst
345 if ((st%occ(ist, is) > m_min_occ).and.(st%eigenval(ist, is) > max_eigen)) then
346 max_eigen = st%eigenval(ist, is)
347 end if
348 end do
349
350 do ist = 1, st%nst
351 ! criterion for degeneracy
352 if (abs(st%eigenval(ist, is)-max_eigen) <= 1e-3_real64) then
353 oep%eigen_type(ist) = 2
354 else
355 if (st%occ(ist, is) > m_min_occ) then
356 oep%eigen_type(ist) = 1
357 oep%eigen_index(oep%eigen_n+1) = ist
358 oep%eigen_n = oep%eigen_n + 1
359 else
360 oep%eigen_type(ist) = 0
361 end if
362 end if
363 end do
364
365 ! find how many states are occupied.
366 do ist = 1, st%nst
367 if (st%occ(ist, is) > m_min_occ) oep%noccst = ist
368 end do
369
370 ! We check that all states have been assigned
371 assert(all(oep%eigen_type >= 0))
372
373 pop_sub(xc_oep_analyzeeigen)
374 end subroutine xc_oep_analyzeeigen
375
376
377#include "xc_kli_pauli_inc.F90"
378#include "xc_oep_sic_pauli_inc.F90"
379
380#include "undef.F90"
381#include "real.F90"
382#include "xc_kli_inc.F90"
383#include "xc_oep_sic_inc.F90"
384#include "xc_oep_inc.F90"
385
386#include "undef.F90"
387#include "complex.F90"
388#include "xc_kli_inc.F90"
389#include "xc_oep_sic_inc.F90"
390#include "xc_oep_inc.F90"
391
392end module xc_oep_oct_m
393
394!! Local Variables:
395!! mode: f90
396!! coding: utf-8
397!! 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:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_min_occ
Minimal occupation that is considered to be non-zero.
Definition: global.F90:218
This module implements the underlying real-space grid.
Definition: grid.F90:119
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:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:999
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:690
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
This module provides routines for perform the Selected Column of Density Matrix (SCDM) method This fo...
Definition: scdm.F90:119
subroutine, public scdm_init(this, namespace)
Definition: scdm.F90:156
subroutine, public scf_tol_init(this, namespace, qtot, def_maximumiter, tol_scheme)
Definition: scf_tol.F90:162
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
Definition: xc.F90:116
integer, parameter, public oep_type_mgga
Definition: xc_oep.F90:187
subroutine, public xc_oep_spinfactor(oep, nspin)
A couple of auxiliary functions for oep.
Definition: xc_oep.F90:399
subroutine, public xc_oep_end(oep)
Definition: xc_oep.F90:359
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:2340
integer, parameter, public oep_type_photons
Definition: xc_oep.F90:187
subroutine, public xc_oep_analyzeeigen(oep, st, is)
Definition: xc_oep.F90:421
integer, parameter, public oep_mixing_scheme_dens
Definition: xc_oep.F90:181
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:1448
integer, parameter, public oep_mixing_scheme_bb
Definition: xc_oep.F90:181
subroutine, public xc_oep_write_info(oep, iunit, namespace)
Definition: xc_oep.F90:381
subroutine, public dxc_oep_mix(oep, mesh, ss, rho, is)
A routine that takes care of mixing the potential.
Definition: xc_oep.F90:1826
integer, parameter, public oep_level_full
Definition: xc_oep.F90:175
subroutine, public xc_oep_init(oep, namespace, gr, st, mc, space, oep_type)
Definition: xc_oep.F90:220
subroutine, public zxc_oep_mix(oep, mesh, ss, rho, is)
A routine that takes care of mixing the potential.
Definition: xc_oep.F90:2718
integer, parameter, public oep_type_sic
Definition: xc_oep.F90:187
integer, parameter, public oep_level_kli
Definition: xc_oep.F90:175
The states_elec_t class contains all electronic wave functions.