Octopus
v_ks.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module v_ks_oct_m
22 use accel_oct_m
23 use comm_oct_m
25 use debug_oct_m
28 use energy_oct_m
32 use global_oct_m
33 use grid_oct_m
37 use ions_oct_m
38 use, intrinsic :: iso_fortran_env
42 use lda_u_oct_m
46 use mesh_oct_m
49 use mpi_oct_m
52 use parser_oct_m
57 use pseudo_oct_m
60 use sort_oct_m
61 use space_oct_m
69 use xc_oct_m
70 use xc_f03_lib_m
71 use xc_fbe_oct_m
72 use xc_functional_oct_m
75 use xc_oep_oct_m
76 use xc_sic_oct_m
78 use xc_vdw_oct_m
80
81 ! from the dftd3 library
82 use dftd3_api
83
84 implicit none
85
86 private
87 public :: &
88 v_ks_t, &
89 v_ks_init, &
90 v_ks_end, &
93 v_ks_calc, &
99
100 type v_ks_calc_t
101 private
102 logical :: calculating
103 logical :: time_present
104 real(real64) :: time
105 real(real64), allocatable :: density(:, :)
106 logical :: total_density_alloc
107 real(real64), pointer, contiguous :: total_density(:)
108 type(energy_t), allocatable :: energy
109 type(states_elec_t), pointer :: hf_st
110 real(real64), allocatable :: vxc(:, :)
111 real(real64), allocatable :: vtau(:, :)
112 real(real64), allocatable :: axc(:, :, :)
113 real(real64), allocatable :: a_ind(:, :)
114 real(real64), allocatable :: b_ind(:, :)
115 logical :: calc_energy
116 end type v_ks_calc_t
117
118 type v_ks_t
119 private
120 integer, public :: theory_level = -1
121
122 logical, public :: frozen_hxc = .false.
123
124 integer, public :: xc_family = 0
125 integer, public :: xc_flags = 0
126 integer, public :: xc_photon = 0
127 type(xc_t), public :: xc
128 type(xc_photons_t), public :: xc_photons
129 type(xc_oep_t), public :: oep
130 type(xc_oep_photon_t), public :: oep_photon
131 type(xc_ks_inversion_t), public :: ks_inversion
132 type(xc_sic_t), public :: sic
133 type(xc_vdw_t), public :: vdw
134 type(grid_t), pointer, public :: gr
135 type(v_ks_calc_t) :: calc
136 logical :: calculate_current = .false.
137 type(current_t) :: current_calculator
138 logical :: include_td_field = .false.
139 logical, public :: has_photons = .false.
140 logical :: xc_photon_include_hartree = .true.
141
142 real(real64), public :: stress_xc_gga(3, 3)
143 type(photon_mode_t), pointer, public :: pt => null()
144 type(mf_t), public :: pt_mx
145 end type v_ks_t
146
147contains
148
149 ! ---------------------------------------------------------
150 subroutine v_ks_init(ks, namespace, gr, st, ions, mc, space, kpoints)
151 type(v_ks_t), intent(inout) :: ks
152 type(namespace_t), intent(in) :: namespace
153 type(grid_t), target, intent(inout) :: gr
154 type(states_elec_t), intent(in) :: st
155 type(ions_t), intent(inout) :: ions
156 type(multicomm_t), intent(in) :: mc
157 class(space_t), intent(in) :: space
158 type(kpoints_t), intent(in) :: kpoints
159
160 integer :: x_id, c_id, xk_id, ck_id, default, val
161 logical :: parsed_theory_level, using_hartree_fock
162 integer :: pseudo_x_functional, pseudo_c_functional
163 integer :: oep_type
164
165 push_sub(v_ks_init)
166
167 ! We need to parse TheoryLevel and XCFunctional, this is
168 ! complicated because they are interdependent.
169
170 !%Variable TheoryLevel
171 !%Type integer
172 !%Section Hamiltonian
173 !%Description
174 !% The calculations can be run with different "theory levels" that
175 !% control how electrons are simulated. The default is
176 !% <tt>dft</tt>. When hybrid functionals are requested, through
177 !% the <tt>XCFunctional</tt> variable, the default is
178 !% <tt>hartree_fock</tt>.
179 !%Option independent_particles 2
180 !% Particles will be considered as independent, <i>i.e.</i> as non-interacting.
181 !% This mode is mainly used for testing purposes, as the code is usually
182 !% much faster with <tt>independent_particles</tt>.
183 !%Option hartree 1
184 !% Calculation within the Hartree method (experimental). Note that, contrary to popular
185 !% belief, the Hartree potential is self-interaction-free. Therefore, this run
186 !% mode will not yield the same result as <tt>kohn-sham</tt> without exchange-correlation.
187 !%Option hartree_fock 3
188 !% This is the traditional Hartree-Fock scheme. Like the Hartree scheme, it is fully
189 !% self-interaction-free.
190 !%Option kohn_sham 4
191 !% This is the default density-functional theory scheme. Note that you can also use
192 !% hybrid functionals in this scheme, but they will be handled the "DFT" way, <i>i.e.</i>,
193 !% solving the OEP equation.
194 !%Option generalized_kohn_sham 5
195 !% This is similar to the <tt>kohn-sham</tt> scheme, except that this allows for nonlocal operators.
196 !% This is the default mode to run hybrid functionals, meta-GGA functionals, or DFT+U.
197 !% It can be more convenient to use <tt>kohn-sham</tt> DFT within the OEP scheme to get similar (but not the same) results.
198 !% Note that within this scheme you can use a correlation functional, or a hybrid
199 !% functional (see <tt>XCFunctional</tt>). In the latter case, you will be following the
200 !% quantum-chemistry recipe to use hybrids.
201 !%Option rdmft 7
202 !% (Experimental) Reduced Density Matrix functional theory.
203 !%End
205 ks%xc_family = xc_family_none
206 ks%sic%level = sic_none
207 ks%oep%level = oep_level_none
208 ks%oep_photon%level = oep_level_none
209
210 ks%theory_level = kohn_sham_dft
211 parsed_theory_level = .false.
212
213 ! the user knows what he wants, give her that
214 if (parse_is_defined(namespace, 'TheoryLevel')) then
215 call parse_variable(namespace, 'TheoryLevel', kohn_sham_dft, ks%theory_level)
216 if (.not. varinfo_valid_option('TheoryLevel', ks%theory_level)) call messages_input_error(namespace, 'TheoryLevel')
218 parsed_theory_level = .true.
219 end if
221 ! parse the XC functional
223 call get_functional_from_pseudos(pseudo_x_functional, pseudo_c_functional)
225 default = 0
226 if (ks%theory_level == kohn_sham_dft .or. ks%theory_level == generalized_kohn_sham_dft) then
227 default = xc_get_default_functional(space%dim, pseudo_x_functional, pseudo_c_functional)
228 end if
230 if (.not. parse_is_defined(namespace, 'XCFunctional') &
231 .and. (pseudo_x_functional /= pseudo_exchange_any .or. pseudo_c_functional /= pseudo_correlation_any)) then
232 call messages_write('Info: the XCFunctional has been selected to match the pseudopotentials', new_line = .true.)
233 call messages_write(' used in the calculation.')
234 call messages_info(namespace=namespace)
235 end if
237 ! The description of this variable can be found in file src/xc/functionals_list.F90
238 call parse_variable(namespace, 'XCFunctional', default, val)
239
240 ! the first 3 digits of the number indicate the X functional and
241 ! the next 3 the C functional.
242 c_id = val / libxc_c_index
243 x_id = val - c_id * libxc_c_index
244
245 if ((x_id /= pseudo_x_functional .and. pseudo_x_functional /= pseudo_exchange_any) .or. &
246 (c_id /= pseudo_c_functional .and. pseudo_c_functional /= pseudo_exchange_any)) then
247 call messages_write('The XCFunctional that you selected does not match the one used', new_line = .true.)
248 call messages_write('to generate the pseudopotentials.')
249 call messages_warning(namespace=namespace)
250 end if
251
252 ! FIXME: we rarely need this. We should only parse when necessary.
253
254 !%Variable XCKernel
255 !%Type integer
256 !%Default -1
257 !%Section Hamiltonian::XC
258 !%Description
259 !% Defines the exchange-correlation kernel. Only LDA kernels are available currently.
260 !% The options are the same as <tt>XCFunctional</tt>.
261 !% Note: the kernel is only needed for Casida, Sternheimer, or optimal-control calculations.
262 !%Option xc_functional -1
263 !% The same functional defined by <tt>XCFunctional</tt>. By default, this is the case.
264 !%End
265 call parse_variable(namespace, 'XCKernel', -1, val)
266 if (-1 == val) then
267 ck_id = c_id
268 xk_id = x_id
269 else
270 ck_id = val / libxc_c_index
271 xk_id = val - ck_id * libxc_c_index
272 end if
273
274 call messages_obsolete_variable(namespace, 'XFunctional', 'XCFunctional')
275 call messages_obsolete_variable(namespace, 'CFunctional', 'XCFunctional')
276
277 !%Variable XCPhotonFunctional
278 !%Type integer
279 !%Default 0
280 !%Section Hamiltonian::XC
281 !%Description
282 !% Defines the exchange and correlation functionals to be used for the QEDFT
283 !% description of the electron-photon system.
284 !%Option none 0
285 !% No functional is used
286 !%Option photon_x_lda 10
287 !% Exchange-only local density approcimation
288 !%Option photon_xc_lda 11
289 !% Exchange-correlation local density approcimation
290 !%Option photon_x_wfn 20
291 !% Exchange-only based on wave functions
292 !%Option photon_xc_wfn 21
293 !% Exchange-correlation based on wave functions
294 !%End
295
296 call parse_variable(namespace, 'XCPhotonFunctional', option__xcphotonfunctional__none, ks%xc_photon)
297
298 !%Variable XCPhotonIncludeHartree
299 !%Type logical
300 !%Default yes
301 !%Section Hamiltonian::XC
302 !%Description
303 !% Use the Hartree potential and energy in calculations
304 !%End
305
306 call parse_variable(namespace, 'XCPhotonIncludeHartree', .true., ks%xc_photon_include_hartree)
307
308 if (.not. ks%xc_photon_include_hartree) then
309 call messages_write('turn off hartree potential and energy')
310 call messages_warning(namespace=namespace)
311 end if
312
313 ! initialize XC modules
314
315 ! This is a bit ugly, theory_level might not be generalized KS or HF now
316 ! but it might become generalized KS or HF later. This is safe because it
317 ! becomes generalized KS in the cases where the functional is hybrid
318 ! and the ifs inside check for both conditions.
319 using_hartree_fock = (ks%theory_level == hartree_fock) &
320 .or. (ks%theory_level == generalized_kohn_sham_dft .and. family_is_hybrid(ks%xc))
321 call xc_init(ks%xc, namespace, space%dim, space%periodic_dim, st%qtot, &
322 x_id, c_id, xk_id, ck_id, hartree_fock = using_hartree_fock, ispin=st%d%ispin)
323
324 ks%xc_family = ks%xc%family
325 ks%xc_flags = ks%xc%flags
326
327 if (.not. parsed_theory_level) then
328 default = kohn_sham_dft
329
330 ! the functional is a hybrid, use Hartree-Fock as theory level by default
331 if (family_is_hybrid(ks%xc) .or. family_is_mgga_with_exc(ks%xc)) then
333 end if
334
335 ! In principle we do not need to parse. However we do it for consistency
336 call parse_variable(namespace, 'TheoryLevel', default, ks%theory_level)
337 if (.not. varinfo_valid_option('TheoryLevel', ks%theory_level)) call messages_input_error(namespace, 'TheoryLevel')
338
339 end if
340
341 ! In case we need OEP, we need to find which type of OEP it is
342 oep_type = -1
343 if (family_is_mgga_with_exc(ks%xc)) then
344 call messages_experimental('MGGA energy functionals')
345
346 if (ks%theory_level == kohn_sham_dft) then
347 call messages_experimental("MGGA within the Kohn-Sham scheme")
348 ks%xc_family = ior(ks%xc_family, xc_family_oep)
349 oep_type = oep_type_mgga
350 end if
351 end if
352
353 call messages_obsolete_variable(namespace, 'NonInteractingElectrons', 'TheoryLevel')
354 call messages_obsolete_variable(namespace, 'HartreeFock', 'TheoryLevel')
355
356 ! Due to how the code is made, we need to set this to have theory level other than DFT
357 ! correct...
358 ks%sic%amaldi_factor = m_one
359
360 select case (ks%theory_level)
362
363 case (hartree)
364 call messages_experimental("Hartree theory level")
365 if (space%periodic_dim == space%dim) then
366 call messages_experimental("Hartree in fully periodic system")
367 end if
368 if (kpoints%full%npoints > 1) then
369 call messages_not_implemented("Hartree with k-points", namespace=namespace)
370 end if
371
372 case (hartree_fock)
373 if (kpoints%full%npoints > 1) then
374 call messages_experimental("Hartree-Fock with k-points")
375 end if
376
378 if (kpoints%full%npoints > 1 .and. family_is_hybrid(ks%xc)) then
379 call messages_experimental("Hybrid functionals with k-points")
380 end if
381
382 case (rdmft)
383 call messages_experimental('RDMFT theory level')
384
385 case (kohn_sham_dft)
386
387 ! check for SIC
388 if (bitand(ks%xc_family, xc_family_lda + xc_family_gga) /= 0) then
389 call xc_sic_init(ks%sic, namespace, gr, st, mc, space)
390 end if
391
392 if (bitand(ks%xc_family, xc_family_oep) /= 0) then
393 select case (ks%xc%functional(func_x,1)%id)
394 case (xc_oep_x_slater)
395 if (kpoints%reduced%npoints > 1) then
396 call messages_not_implemented("Slater with k-points", namespace=namespace)
397 end if
398 ks%oep%level = oep_level_none
399 case (xc_oep_x_fbe)
400 if (kpoints%reduced%npoints > 1) then
401 call messages_not_implemented("FBE functional with k-points", namespace=namespace)
402 end if
403 ks%oep%level = oep_level_none
404 case default
405 if((.not. ks%has_photons) .or. (ks%xc_photon /= 0)) then
406 if(oep_type == -1) then ! Else we have a MGGA
407 oep_type = oep_type_exx
408 end if
409 call xc_oep_init(ks%oep, namespace, gr, st, mc, space, oep_type)
410 end if
411 end select
412 else
413 ks%oep%level = oep_level_none
414 end if
415
416 if (bitand(ks%xc_family, xc_family_ks_inversion) /= 0) then
417 call xc_ks_inversion_init(ks%ks_inversion, namespace, gr, ions, st, ks%xc, mc, space, kpoints)
418 end if
419
420 end select
421
422 if (ks%theory_level /= kohn_sham_dft .and. parse_is_defined(namespace, "SICCorrection")) then
423 message(1) = "SICCorrection can only be used with Kohn-Sham DFT"
424 call messages_fatal(1, namespace=namespace)
425 end if
426
427 if (st%d%ispin == spinors) then
428 if (bitand(ks%xc_family, xc_family_mgga + xc_family_hyb_mgga) /= 0) then
429 call messages_not_implemented("MGGA with spinors", namespace=namespace)
430 end if
431 end if
432
433 ks%frozen_hxc = .false.
434
435 call v_ks_write_info(ks, namespace=namespace)
436
437 ks%gr => gr
438 ks%calc%calculating = .false.
439
440 !The value of ks%calculate_current is set to false or true by Output
441 call current_init(ks%current_calculator, namespace)
442
443 call ks%vdw%init(namespace, space, gr, ks%xc, ions, x_id, c_id)
444
445 if (ks%xc_photon /= 0) then
446 ! initilize the photon free variables
447 call ks%xc_photons%init(namespace, ks%xc_photon , space, gr, st)
448 ! remornalize the electron mass due to light-matter interaction; here we only deal with it in free space
449 ks%oep_photon%level = oep_level_none
450 end if
451
452
453 pop_sub(v_ks_init)
454
455 contains
456
458 subroutine get_functional_from_pseudos(x_functional, c_functional)
459 integer, intent(out) :: x_functional
460 integer, intent(out) :: c_functional
461
462 integer :: xf, cf, ispecies
463 logical :: warned_inconsistent
464
465 x_functional = pseudo_exchange_any
466 c_functional = pseudo_correlation_any
467
468 warned_inconsistent = .false.
469 do ispecies = 1, ions%nspecies
470 select type(spec=>ions%species(ispecies)%s)
471 class is(pseudopotential_t)
472 xf = spec%x_functional()
473 cf = spec%c_functional()
474
475 if (xf == pseudo_exchange_unknown .or. cf == pseudo_correlation_unknown) then
476 call messages_write("Unknown XC functional for species '"//trim(ions%species(ispecies)%s%get_label())//"'")
477 call messages_warning(namespace=namespace)
478 cycle
479 end if
480
481 if (x_functional == pseudo_exchange_any) then
482 x_functional = xf
483 else
484 if (xf /= x_functional .and. .not. warned_inconsistent) then
485 call messages_write('Inconsistent XC functional detected between species')
486 call messages_warning(namespace=namespace)
487 warned_inconsistent = .true.
488 end if
489 end if
490
491 if (c_functional == pseudo_correlation_any) then
492 c_functional = cf
493 else
494 if (cf /= c_functional .and. .not. warned_inconsistent) then
495 call messages_write('Inconsistent XC functional detected between species')
496 call messages_warning(namespace=namespace)
497 warned_inconsistent = .true.
498 end if
499 end if
500
501 class default
504 end select
505
506 end do
507
508 assert(x_functional /= pseudo_exchange_unknown)
509 assert(c_functional /= pseudo_correlation_unknown)
510
511 end subroutine get_functional_from_pseudos
512 end subroutine v_ks_init
513 ! ---------------------------------------------------------
514
515 ! ---------------------------------------------------------
516 subroutine v_ks_end(ks)
517 type(v_ks_t), intent(inout) :: ks
518
519 push_sub(v_ks_end)
520
521 call ks%vdw%end()
522
523 select case (ks%theory_level)
524 case (kohn_sham_dft)
525 if (bitand(ks%xc_family, xc_family_ks_inversion) /= 0) then
526 call xc_ks_inversion_end(ks%ks_inversion)
527 end if
528 if (bitand(ks%xc_family, xc_family_oep) /= 0) then
529 call xc_oep_end(ks%oep)
530 end if
531 call xc_end(ks%xc)
533 call xc_end(ks%xc)
534 end select
535
536 call xc_sic_end(ks%sic)
537
538 if (ks%xc_photon /= 0) then
539 call ks%xc_photons%end()
540 end if
541
542 pop_sub(v_ks_end)
543 end subroutine v_ks_end
544 ! ---------------------------------------------------------
545
546
547 ! ---------------------------------------------------------
548 subroutine v_ks_write_info(ks, iunit, namespace)
549 type(v_ks_t), intent(in) :: ks
550 integer, optional, intent(in) :: iunit
551 type(namespace_t), optional, intent(in) :: namespace
552
553 push_sub(v_ks_write_info)
554
555 call messages_print_with_emphasis(msg="Theory Level", iunit=iunit, namespace=namespace)
556 call messages_print_var_option("TheoryLevel", ks%theory_level, iunit=iunit, namespace=namespace)
557
558 select case (ks%theory_level)
560 call messages_info(iunit=iunit, namespace=namespace)
561 call xc_write_info(ks%xc, iunit, namespace)
562
563 case (kohn_sham_dft)
564 call messages_info(iunit=iunit, namespace=namespace)
565 call xc_write_info(ks%xc, iunit, namespace)
566
567 call messages_info(iunit=iunit, namespace=namespace)
568
569 call xc_sic_write_info(ks%sic, iunit, namespace)
570 call xc_oep_write_info(ks%oep, iunit, namespace)
571 call xc_ks_inversion_write_info(ks%ks_inversion, iunit, namespace)
572
573 end select
574
575 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
576
577 pop_sub(v_ks_write_info)
578 end subroutine v_ks_write_info
579 ! ---------------------------------------------------------
580
581
582 !----------------------------------------------------------
583 subroutine v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
584 type(namespace_t), intent(in) :: namespace
585 type(electron_space_t), intent(in) :: space
586 type(grid_t), intent(in) :: gr
587 type(ions_t), intent(in) :: ions
588 type(partner_list_t), intent(in) :: ext_partners
589 type(states_elec_t), intent(inout) :: st
590 type(v_ks_t), intent(inout) :: ks
591 type(hamiltonian_elec_t), intent(inout) :: hm
592 logical, optional, intent(in) :: calc_eigenval
593 logical, optional, intent(in) :: calc_current
594
595 integer, allocatable :: ind(:)
596 integer :: ist, ik
597 real(real64), allocatable :: copy_occ(:)
598 logical :: calc_eigenval_
599 logical :: calc_current_
600
601 push_sub(v_ks_h_setup)
602
603 calc_eigenval_ = optional_default(calc_eigenval, .true.)
604 calc_current_ = optional_default(calc_current, .true.)
605 call states_elec_fermi(st, namespace, gr)
606 call density_calc(st, gr, st%rho)
607 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval = calc_eigenval_, calc_current = calc_current_) ! get potentials
608
609 if (st%restart_reorder_occs .and. .not. st%fromScratch) then
610 message(1) = "Reordering occupations for restart."
611 call messages_info(1, namespace=namespace)
612
613 safe_allocate(ind(1:st%nst))
614 safe_allocate(copy_occ(1:st%nst))
615
616 do ik = 1, st%nik
617 call sort(st%eigenval(:, ik), ind)
618 copy_occ(1:st%nst) = st%occ(1:st%nst, ik)
619 do ist = 1, st%nst
620 st%occ(ist, ik) = copy_occ(ind(ist))
621 end do
622 end do
623
624 safe_deallocate_a(ind)
625 safe_deallocate_a(copy_occ)
626 end if
627
628 if (calc_eigenval_) call states_elec_fermi(st, namespace, gr) ! occupations
629 call energy_calc_total(namespace, space, hm, gr, st, ext_partners)
630
631 pop_sub(v_ks_h_setup)
632 end subroutine v_ks_h_setup
633
634 ! ---------------------------------------------------------
635 subroutine v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
636 calc_eigenval, time, calc_energy, calc_current, force_semilocal)
637 type(v_ks_t), intent(inout) :: ks
638 type(namespace_t), intent(in) :: namespace
639 type(electron_space_t), intent(in) :: space
640 type(hamiltonian_elec_t), intent(inout) :: hm
641 type(states_elec_t), intent(inout) :: st
642 type(ions_t), intent(in) :: ions
643 type(partner_list_t), intent(in) :: ext_partners
644 logical, optional, intent(in) :: calc_eigenval
645 real(real64), optional, intent(in) :: time
646 logical, optional, intent(in) :: calc_energy
647 logical, optional, intent(in) :: calc_current
648 logical, optional, intent(in) :: force_semilocal
649
650 logical :: calc_current_
651
652 push_sub(v_ks_calc)
653
654 calc_current_ = optional_default(calc_current, .true.)
655
656 call v_ks_calc_start(ks, namespace, space, hm, st, ions, hm%kpoints%latt, ext_partners, time, &
657 calc_energy, calc_current_, force_semilocal=force_semilocal)
658 call v_ks_calc_finish(ks, hm, namespace, space, hm%kpoints%latt, st, &
659 ext_partners, force_semilocal=force_semilocal)
660
661 if (optional_default(calc_eigenval, .false.)) then
662 call energy_calc_eigenvalues(namespace, hm, ks%gr%der, st)
663 end if
664
665 ! Update the magnetic constrain
666 call magnetic_constrain_update(hm%magnetic_constrain, ks%gr, st%d, space, hm%kpoints%latt, ions%pos, st%rho)
667 ! We add the potential to vxc, as this way the potential gets mixed together with vxc
668 ! While this is not ideal, this is a simple practical solution
669 if (hm%magnetic_constrain%level /= constrain_none) then
670 call lalg_axpy(ks%gr%np, st%d%nspin, m_one, hm%magnetic_constrain%pot, hm%ks_pot%vhxc)
671 end if
672
673 pop_sub(v_ks_calc)
674 end subroutine v_ks_calc
675
676 ! ---------------------------------------------------------
677
682 subroutine v_ks_calc_start(ks, namespace, space, hm, st, ions, latt, ext_partners, time, &
683 calc_energy, calc_current, force_semilocal)
684 type(v_ks_t), target, intent(inout) :: ks
685 type(namespace_t), intent(in) :: namespace
686 class(space_t), intent(in) :: space
687 type(hamiltonian_elec_t), target, intent(in) :: hm
688 type(states_elec_t), intent(inout) :: st
689 type(ions_t), intent(in) :: ions
690 type(lattice_vectors_t), intent(in) :: latt
691 type(partner_list_t), intent(in) :: ext_partners
692 real(real64), optional, intent(in) :: time
693 logical, optional, intent(in) :: calc_energy
694 logical, optional, intent(in) :: calc_current
695 logical, optional, intent(in) :: force_semilocal
696
697 logical :: calc_current_
698
699 push_sub(v_ks_calc_start)
700
701 calc_current_ = optional_default(calc_current, .true.) &
702 .and. ks%calculate_current &
703 .and. states_are_complex(st) &
705
706 call profiling_in("KOHN_SHAM_CALC")
707
708 assert(.not. ks%calc%calculating)
709 ks%calc%calculating = .true.
710
711 write(message(1), '(a)') 'Debug: Calculating Kohn-Sham potential.'
712 call messages_info(1, namespace=namespace, debug_only=.true.)
713
714 ks%calc%time_present = present(time)
715 ks%calc%time = optional_default(time, m_zero)
716
717 ks%calc%calc_energy = optional_default(calc_energy, .true.)
718
719 ! If the Hxc term is frozen, there is nothing more to do (WARNING: MISSING ks%calc%energy%intnvxc)
720 if (ks%frozen_hxc) then
721 if (calc_current_) then
722 call states_elec_allocate_current(st, space, ks%gr)
723 call current_calculate(ks%current_calculator, namespace, ks%gr, hm, space, st)
724 end if
725
726 call profiling_out("KOHN_SHAM_CALC")
727 pop_sub(v_ks_calc_start)
728 return
729 end if
730
731 allocate(ks%calc%energy)
732
733 call energy_copy(hm%energy, ks%calc%energy)
734
735 ks%calc%energy%intnvxc = m_zero
736
737 nullify(ks%calc%total_density)
738
739 if (ks%theory_level /= independent_particles .and. abs(ks%sic%amaldi_factor) > m_epsilon) then
740
741 call calculate_density()
742
743 if (poisson_is_async(hm%psolver)) then
744 call dpoisson_solve_start(hm%psolver, ks%calc%total_density)
745 end if
746
747 if (ks%theory_level /= hartree .and. ks%theory_level /= rdmft) call v_a_xc(hm, force_semilocal)
748 else
749 ks%calc%total_density_alloc = .false.
750 end if
751
752 if (calc_current_) then
753 call states_elec_allocate_current(st, space, ks%gr)
754 call current_calculate(ks%current_calculator, namespace, ks%gr, hm, space, st)
755 end if
756
757 nullify(ks%calc%hf_st)
758 if (ks%theory_level == hartree .or. ks%theory_level == hartree_fock &
759 .or. ks%theory_level == rdmft .or. (ks%theory_level == generalized_kohn_sham_dft &
760 .and. family_is_hybrid(ks%xc))) then
761 safe_allocate(ks%calc%hf_st)
762 call states_elec_copy(ks%calc%hf_st, st)
763
764 if (st%parallel_in_states) then
765 if (accel_is_enabled()) then
766 call messages_write('State parallelization of Hartree-Fock exchange is not supported')
767 call messages_new_line()
768 call messages_write('when running with OpenCL/CUDA. Please use domain parallelization')
769 call messages_new_line()
770 call messages_write("or disable acceleration using 'DisableAccel = yes'.")
771 call messages_fatal(namespace=namespace)
772 end if
774 end if
775 end if
776
777
778 ! Calculate the vector potential induced by the electronic current.
779 ! WARNING: calculating the self-induced magnetic field here only makes
780 ! sense if it is going to be used in the Hamiltonian, which does not happen
781 ! now. Otherwise one could just calculate it at the end of the calculation.
782 if (hm%self_induced_magnetic) then
783 safe_allocate(ks%calc%a_ind(1:ks%gr%np_part, 1:space%dim))
784 safe_allocate(ks%calc%b_ind(1:ks%gr%np_part, 1:space%dim))
785 call magnetic_induced(namespace, ks%gr, st, hm%psolver, hm%kpoints, ks%calc%a_ind, ks%calc%b_ind)
786 end if
787
788 if ((ks%has_photons) .and. (ks%calc%time_present) .and. (ks%xc_photon == 0) ) then
789 call mf_calc(ks%pt_mx, ks%gr, st, ions, ks%pt, time)
790 end if
791
792 ! if (ks%has_vibrations) then
793 ! call vibrations_eph_coup(ks%vib, ks%gr, hm, ions, st)
794 ! end if
795
796 call profiling_out("KOHN_SHAM_CALC")
797 pop_sub(v_ks_calc_start)
798
799 contains
800
801 subroutine calculate_density()
802 integer :: ip
803
805
806 ! get density taking into account non-linear core corrections
807 safe_allocate(ks%calc%density(1:ks%gr%np, 1:st%d%nspin))
808 call states_elec_total_density(st, ks%gr, ks%calc%density)
809
810 ! Amaldi correction
811 if (ks%sic%level == sic_amaldi) then
812 call lalg_scal(ks%gr%np, st%d%nspin, ks%sic%amaldi_factor, ks%calc%density)
813 end if
814
815 nullify(ks%calc%total_density)
816 if (allocated(st%rho_core) .or. hm%d%spin_channels > 1) then
817 ks%calc%total_density_alloc = .true.
818
819 safe_allocate(ks%calc%total_density(1:ks%gr%np))
820
821 do ip = 1, ks%gr%np
822 ks%calc%total_density(ip) = sum(ks%calc%density(ip, 1:hm%d%spin_channels))
823 end do
824
825 ! remove non-local core corrections
826 if (allocated(st%rho_core)) then
827 call lalg_axpy(ks%gr%np, -ks%sic%amaldi_factor, st%rho_core, ks%calc%total_density)
828 end if
829 else
830 ks%calc%total_density_alloc = .false.
831 ks%calc%total_density => ks%calc%density(:, 1)
832 end if
833
835 end subroutine calculate_density
836
837 ! ---------------------------------------------------------
838 subroutine v_a_xc(hm, force_semilocal)
839 type(hamiltonian_elec_t), intent(in) :: hm
840 logical, optional, intent(in) :: force_semilocal
841
842 integer :: ispin
843
844 push_sub(v_ks_calc_start.v_a_xc)
845 call profiling_in("XC")
846
847 ks%calc%energy%exchange = m_zero
848 ks%calc%energy%correlation = m_zero
849 ks%calc%energy%xc_j = m_zero
850 ks%calc%energy%vdw = m_zero
851
852 allocate(ks%calc%vxc(1:ks%gr%np, 1:st%d%nspin))
853 ks%calc%vxc = m_zero
854
855 if (family_is_mgga_with_exc(hm%xc)) then
856 safe_allocate(ks%calc%vtau(1:ks%gr%np, 1:st%d%nspin))
857 ks%calc%vtau = m_zero
858 end if
859
860 ! Get the *local* XC term
861 if (ks%calc%calc_energy) then
862 if (family_is_mgga_with_exc(hm%xc)) then
863 call xc_get_vxc(ks%gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, ks%calc%density, st%d%ispin, &
864 latt%rcell_volume, ks%calc%vxc, ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation, &
865 deltaxc = ks%calc%energy%delta_xc, vtau = ks%calc%vtau, force_orbitalfree=force_semilocal)
866 else
867 call xc_get_vxc(ks%gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, ks%calc%density, st%d%ispin, &
868 latt%rcell_volume, ks%calc%vxc, ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation, &
869 deltaxc = ks%calc%energy%delta_xc, stress_xc=ks%stress_xc_gga, force_orbitalfree=force_semilocal)
870 end if
871 else
872 if (family_is_mgga_with_exc(hm%xc)) then
873 call xc_get_vxc(ks%gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, ks%calc%density, &
874 st%d%ispin, latt%rcell_volume, ks%calc%vxc, vtau = ks%calc%vtau, force_orbitalfree=force_semilocal)
875 else
876 call xc_get_vxc(ks%gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, ks%calc%density, &
877 st%d%ispin, latt%rcell_volume, ks%calc%vxc, stress_xc=ks%stress_xc_gga, force_orbitalfree=force_semilocal)
878 end if
879 end if
880
881 !Noncollinear functionals
882 if (bitand(hm%xc%family, xc_family_nc_lda + xc_family_nc_mgga) /= 0) then
883 if (st%d%ispin /= spinors) then
884 message(1) = "Noncollinear functionals can only be used with spinor wavefunctions."
885 call messages_fatal(1)
886 end if
887
888 if (optional_default(force_semilocal, .false.)) then
889 message(1) = "Cannot perform LCAO for noncollinear MGGAs."
890 message(2) = "Please perform a LDA calculation first."
891 call messages_fatal(2)
892 end if
893
894 if (ks%calc%calc_energy) then
895 if (family_is_mgga_with_exc(hm%xc)) then
896 call xc_get_nc_vxc(ks%gr, ks%xc, st, hm%kpoints, space, namespace, ks%calc%density, ks%calc%vxc, &
897 vtau = ks%calc%vtau, ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation)
898 else
899 call xc_get_nc_vxc(ks%gr, ks%xc, st, hm%kpoints, space, namespace, ks%calc%density, ks%calc%vxc, &
900 ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation)
901 end if
902 else
903 if (family_is_mgga_with_exc(hm%xc)) then
904 call xc_get_nc_vxc(ks%gr, ks%xc, st, hm%kpoints, space, namespace, ks%calc%density, &
905 ks%calc%vxc, vtau = ks%calc%vtau)
906 else
907 call xc_get_nc_vxc(ks%gr, ks%xc, st, hm%kpoints, space, namespace, ks%calc%density, ks%calc%vxc)
908 end if
909 end if
910 end if
911
912 if (optional_default(force_semilocal, .false.)) then
913 call profiling_out("XC")
914 pop_sub(v_ks_calc_start.v_a_xc)
915 return
916 end if
917
918 ! ADSIC correction
919 if (ks%sic%level == sic_adsic) then
920 if (family_is_mgga(hm%xc%family)) then
921 call messages_not_implemented('ADSIC with MGGAs', namespace=namespace)
922 end if
923 if (ks%calc%calc_energy) then
924 call xc_sic_calc_adsic(ks%sic, namespace, space, ks%gr, st, hm, ks%xc, ks%calc%density, &
925 ks%calc%vxc, ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation)
926 else
927 call xc_sic_calc_adsic(ks%sic, namespace, space, ks%gr, st, hm, ks%xc, ks%calc%density, &
928 ks%calc%vxc)
929 end if
930 end if
931 !PZ SIC is done in the finish routine as OEP full needs to update the Hamiltonian
932
933 if (ks%theory_level == kohn_sham_dft) then
934 ! The OEP family has to be handled specially
935 ! Note that OEP is done in the finish state, as it requires updating the Hamiltonian and needs the new Hartre and vxc term
936 if (bitand(ks%xc_family, xc_family_oep) /= 0 .or. family_is_mgga_with_exc(ks%xc)) then
937
938 if (ks%xc%functional(func_x,1)%id == xc_oep_x_slater) then
939 call x_slater_calc(namespace, ks%gr, space, hm%exxop, st, hm%kpoints, ks%calc%energy%exchange, &
940 vxc = ks%calc%vxc)
941 else if (ks%xc%functional(func_x,1)%id == xc_oep_x_fbe .or. ks%xc%functional(func_x,1)%id == xc_oep_x_fbe_sl) then
942 call x_fbe_calc(ks%xc%functional(func_x,1)%id, namespace, hm%psolver, ks%gr, st, space, &
943 ks%calc%energy%exchange, vxc = ks%calc%vxc)
944
945 else if (ks%xc%functional(func_c,1)%id == xc_lda_c_fbe_sl) then
946
947 call fbe_c_lda_sl(namespace, hm%psolver, ks%gr, st, space, &
948 ks%calc%energy%correlation, vxc = ks%calc%vxc)
949
950 else
951 if (ks%oep_photon%level /= oep_level_none) then
952 if (states_are_real(st)) then
953 call dxc_oep_photon_calc(ks%oep_photon, namespace, ks%xc, ks%gr, &
954 hm, st, space, ks%calc%energy%exchange, ks%calc%energy%correlation, vxc = ks%calc%vxc)
955 else
956 call zxc_oep_photon_calc(ks%oep_photon, namespace, ks%xc, ks%gr, &
957 hm, st, space, ks%calc%energy%exchange, ks%calc%energy%correlation, vxc = ks%calc%vxc)
958 end if
959 ks%calc%energy%photon_exchange = ks%oep_photon%pt%ex
960 end if
961 end if
962
963 end if
964
965 if (bitand(ks%xc_family, xc_family_ks_inversion) /= 0) then
966 ! Also treat KS inversion separately (not part of libxc)
967 call xc_ks_inversion_calc(ks%ks_inversion, namespace, space, ks%gr, hm, ext_partners, st, vxc = ks%calc%vxc, &
968 time = ks%calc%time)
969 end if
970
971 ! compute the photon-free photon exchange potential and energy
972 if (ks%xc_photon /= 0) then
973
974 call ks%xc_photons%v_ks(namespace, ks%calc%total_density, ks%gr, space, hm%psolver, hm%ep, st)
975
976 ! add the photon-free px potential into the xc potential
977 do ispin = 1, hm%d%spin_channels
978 call lalg_axpy(ks%gr%np, m_one, ks%xc_photons%vpx(1:ks%gr%np), ks%calc%vxc(1:ks%gr%np, ispin) )
979 end do
980
981 ! photon-exchange energy
982 ks%calc%energy%photon_exchange = ks%xc_photons%ex
983 end if
984
985 end if
986
987 call ks%vdw%calc(namespace, space, latt, ions%atom, ions%natoms, ions%pos, &
988 ks%gr, st, ks%calc%energy%vdw, ks%calc%vxc)
989
990 if (ks%calc%calc_energy) then
991 ! MGGA vtau contribution is done after copying vtau to hm%vtau
992
993 if (hm%lda_u_level /= dft_u_none) then
994 if (states_are_real(st)) then
995 ks%calc%energy%int_dft_u = denergy_calc_electronic(namespace, hm, ks%gr%der, st, terms = term_dft_u)
996 else
997 ks%calc%energy%int_dft_u = zenergy_calc_electronic(namespace, hm, ks%gr%der, st, terms = term_dft_u)
998 end if
999 end if
1000 end if
1001
1002 call profiling_out("XC")
1003 pop_sub(v_ks_calc_start.v_a_xc)
1004 end subroutine v_a_xc
1005
1006 end subroutine v_ks_calc_start
1007 ! ---------------------------------------------------------
1008
1009 subroutine v_ks_calc_finish(ks, hm, namespace, space, latt, st, ext_partners, force_semilocal)
1010 type(v_ks_t), target, intent(inout) :: ks
1011 type(hamiltonian_elec_t), intent(inout) :: hm
1012 type(namespace_t), intent(in) :: namespace
1013 class(space_t), intent(in) :: space
1014 type(lattice_vectors_t), intent(in) :: latt
1015 type(states_elec_t), intent(inout) :: st
1016 type(partner_list_t), intent(in) :: ext_partners
1017 logical, optional, intent(in) :: force_semilocal
1018
1019 integer :: ip, ispin
1020 type(states_elec_t) :: xst
1022 real(real64) :: exx_energy
1023 real(real64) :: factor
1024
1025 push_sub(v_ks_calc_finish)
1026
1027 assert(ks%calc%calculating)
1028 ks%calc%calculating = .false.
1029
1030 if (ks%frozen_hxc) then
1031 pop_sub(v_ks_calc_finish)
1032 return
1033 end if
1034
1035 !change the pointer to the energy object
1036 safe_deallocate_a(hm%energy)
1037 call move_alloc(ks%calc%energy, hm%energy)
1038
1039 if (hm%self_induced_magnetic) then
1040 hm%a_ind(1:ks%gr%np, 1:space%dim) = ks%calc%a_ind(1:ks%gr%np, 1:space%dim)
1041 hm%b_ind(1:ks%gr%np, 1:space%dim) = ks%calc%b_ind(1:ks%gr%np, 1:space%dim)
1042
1043 safe_deallocate_a(ks%calc%a_ind)
1044 safe_deallocate_a(ks%calc%b_ind)
1045 end if
1046
1047 if (allocated(hm%v_static)) then
1048 hm%energy%intnvstatic = dmf_dotp(ks%gr, ks%calc%total_density, hm%v_static)
1049 else
1050 hm%energy%intnvstatic = m_zero
1051 end if
1052
1053 if (ks%theory_level == independent_particles .or. abs(ks%sic%amaldi_factor) <= m_epsilon) then
1054
1055 hm%ks_pot%vhxc = m_zero
1056 hm%energy%intnvxc = m_zero
1057 hm%energy%hartree = m_zero
1058 hm%energy%exchange = m_zero
1059 hm%energy%exchange_hf = m_zero
1060 hm%energy%correlation = m_zero
1061 else
1062
1063 hm%energy%hartree = m_zero
1064 call v_ks_hartree(namespace, ks, space, hm, ext_partners)
1065
1066 if (.not. optional_default(force_semilocal, .false.)) then
1067 !PZ-SIC
1068 if(ks%sic%level == sic_pz_oep) then
1069 if (states_are_real(st)) then
1070 call dxc_oep_calc(ks%sic%oep, namespace, ks%xc, ks%gr, hm, st, space, &
1071 latt%rcell_volume, hm%energy%exchange, hm%energy%correlation, vxc = ks%calc%vxc)
1072 else
1073 call zxc_oep_calc(ks%sic%oep, namespace, ks%xc, ks%gr, hm, st, space, &
1074 latt%rcell_volume, hm%energy%exchange, hm%energy%correlation, vxc = ks%calc%vxc)
1075 end if
1076 end if
1077
1078 ! OEP for exchange ad MGGAs (within Kohn-Sham DFT)
1079 if (ks%theory_level == kohn_sham_dft .and. ks%oep%level /= oep_level_none) then
1080 ! The OEP family has to be handled specially
1081 if (ks%xc%functional(func_x,1)%id == xc_oep_x .or. family_is_mgga_with_exc(ks%xc)) then
1082 if (states_are_real(st)) then
1083 call dxc_oep_calc(ks%oep, namespace, ks%xc, ks%gr, hm, st, space, &
1084 latt%rcell_volume, hm%energy%exchange, hm%energy%correlation, vxc = ks%calc%vxc)
1085 else
1086 call zxc_oep_calc(ks%oep, namespace, ks%xc, ks%gr, hm, st, space, &
1087 latt%rcell_volume, hm%energy%exchange, hm%energy%correlation, vxc = ks%calc%vxc)
1088 end if
1089 end if
1090 end if
1091 end if
1092
1093 if (ks%calc%calc_energy) then
1094 ! Now we calculate Int[n vxc] = energy%intnvxc
1095 hm%energy%intnvxc = m_zero
1096
1097 if (ks%theory_level /= independent_particles .and. ks%theory_level /= hartree .and. ks%theory_level /= rdmft) then
1098 do ispin = 1, hm%d%nspin
1099 if (ispin <= 2) then
1100 factor = m_one
1101 else
1102 factor = m_two
1103 end if
1104 hm%energy%intnvxc = hm%energy%intnvxc + &
1105 factor*dmf_dotp(ks%gr, st%rho(:, ispin), ks%calc%vxc(:, ispin), reduce = .false.)
1106 end do
1107 if (ks%gr%parallel_in_domains) call ks%gr%allreduce(hm%energy%intnvxc)
1108 end if
1109 end if
1110
1111
1112 if (ks%theory_level /= hartree .and. ks%theory_level /= rdmft) then
1113 ! move allocation of vxc from ks%calc to hm
1114 safe_deallocate_a(hm%ks_pot%vxc)
1115 call move_alloc(ks%calc%vxc, hm%ks_pot%vxc)
1116
1117 if (family_is_mgga_with_exc(hm%xc)) then
1118 call hm%ks_pot%set_vtau(ks%calc%vtau)
1119 safe_deallocate_a(ks%calc%vtau)
1120
1121 ! We need to evaluate the energy after copying vtau to hm%vtau
1122 if (ks%theory_level == generalized_kohn_sham_dft .and. ks%calc%calc_energy) then
1123 ! MGGA vtau contribution
1124 if (states_are_real(st)) then
1125 hm%energy%intnvxc = hm%energy%intnvxc &
1126 + denergy_calc_electronic(namespace, hm, ks%gr%der, st, terms = term_mgga)
1127 else
1128 hm%energy%intnvxc = hm%energy%intnvxc &
1129 + zenergy_calc_electronic(namespace, hm, ks%gr%der, st, terms = term_mgga)
1130 end if
1131 end if
1132 end if
1133
1134 else
1135 hm%ks_pot%vxc = m_zero
1136 end if
1137
1138 if (.not. ks%xc_photon_include_hartree) then
1139 hm%energy%hartree = m_zero
1140 hm%ks_pot%vhartree = m_zero
1141 end if
1142
1143 ! Build Hartree + XC potential
1144
1145 do ip = 1, ks%gr%np
1146 hm%ks_pot%vhxc(ip, 1) = hm%ks_pot%vxc(ip, 1) + hm%ks_pot%vhartree(ip)
1147 end do
1148 if (allocated(hm%vberry)) then
1149 do ip = 1, ks%gr%np
1150 hm%ks_pot%vhxc(ip, 1) = hm%ks_pot%vhxc(ip, 1) + hm%vberry(ip, 1)
1151 end do
1152 end if
1153
1154 if (hm%d%ispin > unpolarized) then
1155 do ip = 1, ks%gr%np
1156 hm%ks_pot%vhxc(ip, 2) = hm%ks_pot%vxc(ip, 2) + hm%ks_pot%vhartree(ip)
1157 end do
1158 if (allocated(hm%vberry)) then
1159 do ip = 1, ks%gr%np
1160 hm%ks_pot%vhxc(ip, 2) = hm%ks_pot%vhxc(ip, 2) + hm%vberry(ip, 2)
1161 end do
1162 end if
1163 end if
1164
1165 if (hm%d%ispin == spinors) then
1166 do ispin=3, 4
1167 do ip = 1, ks%gr%np
1168 hm%ks_pot%vhxc(ip, ispin) = hm%ks_pot%vxc(ip, ispin)
1169 end do
1170 end do
1171 end if
1172
1173 ! Note: this includes hybrids calculated with the Fock operator instead of OEP
1174 hm%energy%exchange_hf = m_zero
1175 if ((ks%theory_level == hartree .or. ks%theory_level == hartree_fock &
1176 .or. ks%theory_level == rdmft .or. ks%theory_level == generalized_kohn_sham_dft)) then
1177
1178 ! swap the states object
1179 if (associated(hm%exxop%st)) then
1180 if (hm%exxop%st%parallel_in_states) call states_elec_parallel_remote_access_stop(hm%exxop%st)
1181 call states_elec_end(hm%exxop%st)
1182 safe_deallocate_p(hm%exxop%st)
1183 end if
1184 if (associated(ks%calc%hf_st) .and. hm%exxop%useACE) then
1185 if (ks%calc%hf_st%parallel_in_states) call states_elec_parallel_remote_access_stop(ks%calc%hf_st)
1186 end if
1187
1188 !At the moment this block is called before the reinit call. This way the LCAO does not call the
1189 !exchange operator.
1190 !This should be changed and the CAM parameters should also be obtained from the restart information
1191 !Maybe the parameters should be mixed too.
1192 exx_energy = m_zero
1193 if (.not. optional_default(force_semilocal, .false.)) then
1194 if ((ks%theory_level == hartree_fock .or. ks%theory_level == rdmft &
1195 .or. (ks%theory_level == generalized_kohn_sham_dft &
1196 .and. family_is_hybrid(ks%xc))) .and. hm%exxop%useACE) then
1197 call xst%nullify()
1198 if (states_are_real(ks%calc%hf_st)) then
1199 call dexchange_operator_compute_potentials(hm%exxop, namespace, space, ks%gr, &
1200 ks%calc%hf_st, xst, hm%kpoints, exx_energy)
1201 call dexchange_operator_ace(hm%exxop, namespace, ks%gr, ks%calc%hf_st, xst)
1202 else
1203 call zexchange_operator_compute_potentials(hm%exxop, namespace, space, ks%gr, &
1204 ks%calc%hf_st, xst, hm%kpoints, exx_energy)
1205 if (hm%phase%is_allocated()) then
1206 call zexchange_operator_ace(hm%exxop, namespace, ks%gr, ks%calc%hf_st, xst, hm%phase)
1207 else
1208 call zexchange_operator_ace(hm%exxop, namespace, ks%gr, ks%calc%hf_st, xst)
1209 end if
1210 end if
1211 call states_elec_end(xst)
1212 exx_energy = exx_energy + hm%exxop%singul%energy
1213 end if
1214 end if
1215
1216 select case (ks%theory_level)
1218 if (family_is_hybrid(ks%xc)) then
1219 ! Add the energy only the ACE case. In the non-ACE case, the exchange energy is not correct
1220 hm%energy%exchange_hf = hm%energy%exchange_hf + exx_energy
1221 call exchange_operator_reinit(hm%exxop, ks%xc%cam_omega, ks%xc%cam_alpha, ks%xc%cam_beta, ks%calc%hf_st)
1222 end if
1223 case (hartree_fock)
1224 ! Add the energy only the ACE case. In the non-ACE case, the exchange energy is not correct
1225 hm%energy%exchange_hf = hm%energy%exchange_hf + exx_energy
1226 call exchange_operator_reinit(hm%exxop, ks%xc%cam_omega, ks%xc%cam_alpha, ks%xc%cam_beta, ks%calc%hf_st)
1227 case (hartree)
1228 call exchange_operator_reinit(hm%exxop, m_zero, m_one, m_zero, ks%calc%hf_st)
1229 case (rdmft)
1230 call exchange_operator_reinit(hm%exxop, m_zero, m_one, m_zero, ks%calc%hf_st)
1231 end select
1232 end if
1233
1234 end if
1235
1236 ! Because of the intent(in) in v_ks_calc_start, we need to update the parameters of hybrids for OEP
1237 ! here
1238 if (ks%theory_level == kohn_sham_dft .and. bitand(ks%xc_family, xc_family_oep) /= 0) then
1239 if (ks%xc%functional(func_x,1)%id /= xc_oep_x_slater .and. ks%xc%functional(func_x,1)%id /= xc_oep_x_fbe) then
1240 call exchange_operator_reinit(hm%exxop, ks%xc%cam_omega, ks%xc%cam_alpha, ks%xc%cam_beta)
1241 end if
1242 end if
1243
1244 if (ks%has_photons .and. (ks%xc_photon == 0)) then
1245 if (associated(ks%pt_mx%vmf)) then
1246 forall(ip = 1:ks%gr%np) hm%ks_pot%vhxc(ip, 1) = hm%ks_pot%vhxc(ip, 1) + ks%pt_mx%vmf(ip)
1247 if (hm%d%ispin > unpolarized) then
1248 forall(ip = 1:ks%gr%np) hm%ks_pot%vhxc(ip, 2) = hm%ks_pot%vhxc(ip, 2) + ks%pt_mx%vmf(ip)
1249 end if
1250 end if
1251 hm%ep%photon_forces(1:space%dim) = ks%pt_mx%fmf(1:space%dim)
1252 end if
1253
1254 if (ks%vdw%vdw_correction /= option__vdwcorrection__none) then
1255 hm%ep%vdw_forces = ks%vdw%forces
1256 hm%ep%vdw_stress = ks%vdw%stress
1257 safe_deallocate_a(ks%vdw%forces)
1258 else
1259 hm%ep%vdw_forces = 0.0_real64
1260 end if
1261
1262 if (ks%calc%time_present .or. hm%time_zero) then
1263 call hm%update(ks%gr, namespace, space, ext_partners, time = ks%calc%time)
1264 else
1265 call hamiltonian_elec_update_pot(hm, ks%gr)
1266 end if
1267
1268
1269 safe_deallocate_a(ks%calc%density)
1270 if (ks%calc%total_density_alloc) then
1271 safe_deallocate_p(ks%calc%total_density)
1272 end if
1273 nullify(ks%calc%total_density)
1274
1275 pop_sub(v_ks_calc_finish)
1276 end subroutine v_ks_calc_finish
1277
1278 ! ---------------------------------------------------------
1279 !
1283 !
1284 subroutine v_ks_hartree(namespace, ks, space, hm, ext_partners)
1285 type(namespace_t), intent(in) :: namespace
1286 type(v_ks_t), intent(inout) :: ks
1287 class(space_t), intent(in) :: space
1288 type(hamiltonian_elec_t), intent(inout) :: hm
1289 type(partner_list_t), intent(in) :: ext_partners
1290
1291 push_sub(v_ks_hartree)
1292
1293 if (.not. poisson_is_async(hm%psolver)) then
1294 ! solve the Poisson equation
1295 call dpoisson_solve(hm%psolver, namespace, hm%ks_pot%vhartree, ks%calc%total_density, reset=.false.)
1296 else
1297 ! The calculation was started by v_ks_calc_start.
1298 call dpoisson_solve_finish(hm%psolver, hm%ks_pot%vhartree)
1299 end if
1300
1301 if (ks%calc%calc_energy) then
1302 ! Get the Hartree energy
1303 hm%energy%hartree = m_half*dmf_dotp(ks%gr, ks%calc%total_density, hm%ks_pot%vhartree)
1304 end if
1305
1307 if(ks%calc%time_present) then
1308 if(hamiltonian_elec_has_kick(hm)) then
1309 call pcm_hartree_potential(hm%pcm, space, ks%gr, hm%psolver, ext_partners, hm%ks_pot%vhartree, &
1310 ks%calc%total_density, hm%energy%pcm_corr, kick=hm%kick, time=ks%calc%time)
1311 else
1312 call pcm_hartree_potential(hm%pcm, space, ks%gr, hm%psolver, ext_partners, hm%ks_pot%vhartree, &
1313 ks%calc%total_density, hm%energy%pcm_corr, time=ks%calc%time)
1314 end if
1315 else
1316 if(hamiltonian_elec_has_kick(hm)) then
1317 call pcm_hartree_potential(hm%pcm, space, ks%gr, hm%psolver, ext_partners, hm%ks_pot%vhartree, &
1318 ks%calc%total_density, hm%energy%pcm_corr, kick=hm%kick)
1319 else
1320 call pcm_hartree_potential(hm%pcm, space, ks%gr, hm%psolver, ext_partners, hm%ks_pot%vhartree, &
1321 ks%calc%total_density, hm%energy%pcm_corr)
1322 end if
1323 end if
1324
1325 pop_sub(v_ks_hartree)
1326 end subroutine v_ks_hartree
1327 ! ---------------------------------------------------------
1328
1329
1330 ! ---------------------------------------------------------
1331 subroutine v_ks_freeze_hxc(ks)
1332 type(v_ks_t), intent(inout) :: ks
1333
1334 push_sub(v_ks_freeze_hxc)
1335
1336 ks%frozen_hxc = .true.
1337
1338 pop_sub(v_ks_freeze_hxc)
1339 end subroutine v_ks_freeze_hxc
1340 ! ---------------------------------------------------------
1341
1342 subroutine v_ks_calculate_current(this, calc_cur)
1343 type(v_ks_t), intent(inout) :: this
1344 logical, intent(in) :: calc_cur
1345
1346 push_sub(v_ks_calculate_current)
1347
1348 this%calculate_current = calc_cur
1349
1350 pop_sub(v_ks_calculate_current)
1351 end subroutine v_ks_calculate_current
1352
1353end module v_ks_oct_m
1354
1355!! Local Variables:
1356!! mode: f90
1357!! coding: utf-8
1358!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:171
scales a vector by a constant
Definition: lalg_basic.F90:157
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:149
pure logical function, public accel_is_enabled()
Definition: accel.F90:400
subroutine, public current_calculate(this, namespace, gr, hm, space, st)
Compute total electronic current density.
Definition: current.F90:368
subroutine, public current_init(this, namespace)
Definition: current.F90:178
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
Definition: density.F90:850
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:609
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
subroutine, public energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit, full)
This subroutine calculates the total energy of the system. Basically, it adds up the KS eigenvalues,...
real(real64) function, public zenergy_calc_electronic(namespace, hm, der, st, terms)
real(real64) function, public denergy_calc_electronic(namespace, hm, der, st, terms)
subroutine, public energy_calc_eigenvalues(namespace, hm, der, st)
subroutine, public energy_copy(ein, eout)
Definition: energy.F90:167
subroutine, public dexchange_operator_ace(this, namespace, mesh, st, xst, phase)
subroutine, public zexchange_operator_compute_potentials(this, namespace, space, gr, st, xst, kpoints, ex, F_out)
subroutine, public dexchange_operator_compute_potentials(this, namespace, space, gr, st, xst, kpoints, ex, F_out)
subroutine, public zexchange_operator_ace(this, namespace, mesh, st, xst, phase)
subroutine, public exchange_operator_reinit(this, omega, alpha, beta, st)
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_epsilon
Definition: global.F90:204
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
integer, parameter, public term_mgga
integer, parameter, public term_dft_u
logical function, public hamiltonian_elec_has_kick(hm)
logical function, public hamiltonian_elec_needs_current(hm, states_are_real)
subroutine, public hamiltonian_elec_update_pot(this, mesh, accumulate)
Update the KS potential of the electronic Hamiltonian.
This module defines classes and functions for interaction partners.
A module to handle KS potential, without the external potential.
integer, parameter, public rdmft
integer, parameter, public hartree
integer, parameter, public hartree_fock
integer, parameter, public independent_particles
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public kohn_sham_dft
integer, parameter, public dft_u_none
Definition: lda_u.F90:201
This modules implements the routines for doing constrain DFT for noncollinear magnetism.
integer, parameter, public constrain_none
subroutine, public magnetic_constrain_update(this, mesh, std, space, latt, pos, rho)
Recomputes the magnetic contraining potential.
subroutine, public magnetic_induced(namespace, gr, st, psolver, kpoints, a_ind, b_ind)
This subroutine receives as input a current, and produces as an output the vector potential that it i...
Definition: magnetic.F90:383
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_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:920
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1113
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1045
subroutine, public messages_new_line()
Definition: messages.F90:1134
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1085
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
subroutine, public pcm_hartree_potential(pcm, space, mesh, psolver, ext_partners, vhartree, density, pcm_corr, kick, time)
PCM reaction field due to the electronic density.
subroutine, public mf_calc(this, gr, st, ions, pt_mode, time)
subroutine, public dpoisson_solve_start(this, rho)
Definition: poisson.F90:2115
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:892
subroutine, public dpoisson_solve_finish(this, pot)
Definition: poisson.F90:2138
logical pure function, public poisson_is_async(this)
Definition: poisson.F90:1247
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
integer, parameter, public pseudo_exchange_unknown
Definition: pseudo.F90:188
integer, parameter, public pseudo_correlation_unknown
Definition: pseudo.F90:192
integer, parameter, public pseudo_correlation_any
Definition: pseudo.F90:192
integer, parameter, public pseudo_exchange_any
Definition: pseudo.F90:188
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
integer, parameter, private libxc_c_index
Definition: species.F90:276
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_fermi(st, namespace, mesh, compute_spin)
calculate the Fermi level for the states in this object
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
subroutine, public states_elec_allocate_current(st, space, mesh)
This module provides routines for communicating states when using states parallelization.
subroutine, public states_elec_parallel_remote_access_stop(this)
stop remote memory access for states on other processors
subroutine, public states_elec_parallel_remote_access_start(this)
start remote memory access for states on other processors
subroutine v_ks_hartree(namespace, ks, space, hm, ext_partners)
Hartree contribution to the KS potential. This function is designed to be used by v_ks_calc_finish an...
Definition: v_ks.F90:1378
subroutine, public v_ks_calc_finish(ks, hm, namespace, space, latt, st, ext_partners, force_semilocal)
Definition: v_ks.F90:1103
subroutine, public v_ks_freeze_hxc(ks)
Definition: v_ks.F90:1425
subroutine, public v_ks_end(ks)
Definition: v_ks.F90:610
subroutine, public v_ks_calculate_current(this, calc_cur)
Definition: v_ks.F90:1436
subroutine, public v_ks_write_info(ks, iunit, namespace)
Definition: v_ks.F90:642
subroutine, public v_ks_calc_start(ks, namespace, space, hm, st, ions, latt, ext_partners, time, calc_energy, calc_current, force_semilocal)
This routine starts the calculation of the Kohn-Sham potential. The routine v_ks_calc_finish must be ...
Definition: v_ks.F90:777
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:730
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
Definition: v_ks.F90:677
subroutine, public v_ks_init(ks, namespace, gr, st, ions, mc, space, kpoints)
Definition: v_ks.F90:244
subroutine, public x_slater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
Interface to X(slater_calc)
Definition: x_slater.F90:145
subroutine, public x_fbe_calc(id, namespace, psolver, gr, st, space, ex, vxc)
Interface to X(x_fbe_calc) Two possible run modes possible: adiabatic and Sturm-Liouville....
Definition: xc_fbe.F90:170
subroutine, public fbe_c_lda_sl(namespace, psolver, gr, st, space, ec, vxc)
Sturm-Liouville version of the FBE local-density correlation functional.
Definition: xc_fbe.F90:436
integer, parameter, public func_c
integer, parameter, public func_x
subroutine, public xc_ks_inversion_end(ks_inv)
subroutine, public xc_ks_inversion_write_info(ks_inversion, iunit, namespace)
subroutine, public xc_ks_inversion_init(ks_inv, namespace, gr, ions, st, xc, mc, space, kpoints)
subroutine, public xc_ks_inversion_calc(ks_inversion, namespace, space, gr, hm, ext_partners, st, vxc, time)
integer, parameter, public oep_type_mgga
Definition: xc_oep.F90:184
integer, parameter, public oep_level_none
the OEP levels
Definition: xc_oep.F90:172
subroutine, public xc_oep_end(oep)
Definition: xc_oep.F90:356
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:2281
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:1444
subroutine, public xc_oep_write_info(oep, iunit, namespace)
Definition: xc_oep.F90:378
integer, parameter, public oep_type_exx
The different types of OEP that we can work with.
Definition: xc_oep.F90:184
subroutine, public xc_oep_init(oep, namespace, gr, st, mc, space, oep_type)
Definition: xc_oep.F90:217
subroutine, public zxc_oep_photon_calc(oep, namespace, xcs, gr, hm, st, space, ex, ec, vxc)
This file handles the evaluation of the OEP potential, in the KLI or full OEP as described in S....
subroutine, public dxc_oep_photon_calc(oep, namespace, xcs, gr, hm, st, space, ex, ec, vxc)
This file handles the evaluation of the OEP potential, in the KLI or full OEP as described in S....
This module implements the "photon-free" electron-photon exchange-correlation functional.
Definition: xc_photons.F90:121
integer, parameter, public sic_none
no self-interaction correction
Definition: xc_sic.F90:148
subroutine, public xc_sic_write_info(sic, iunit, namespace)
Definition: xc_sic.F90:251
integer, parameter, public sic_adsic
Averaged density SIC.
Definition: xc_sic.F90:148
subroutine, public xc_sic_init(sic, namespace, gr, st, mc, space)
initialize the SIC object
Definition: xc_sic.F90:168
subroutine, public xc_sic_end(sic)
finalize the SIC and, if needed, the included OEP
Definition: xc_sic.F90:237
integer, parameter, public sic_pz_oep
Perdew-Zunger SIC (OEP way)
Definition: xc_sic.F90:148
integer, parameter, public sic_amaldi
Amaldi correction term.
Definition: xc_sic.F90:148
subroutine, public xc_sic_calc_adsic(sic, namespace, space, gr, st, hm, xc, density, vxc, ex, ec)
Computes the ADSIC potential and energy.
Definition: xc_sic.F90:282
A module that takes care of xc contribution from vdW interactions.
Definition: xc_vdw.F90:116
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
The states_elec_t class contains all electronic wave functions.
int true(void)
subroutine get_functional_from_pseudos(x_functional, c_functional)
Tries to find out the functional from the pseudopotential.
Definition: v_ks.F90:552
subroutine v_a_xc(hm, force_semilocal)
Definition: v_ks.F90:932
subroutine calculate_density()
Definition: v_ks.F90:895