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