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
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, &
100
101 type v_ks_calc_t
102 private
103 logical :: calculating
104 logical :: time_present
105 real(real64) :: time
106 real(real64), allocatable :: density(:, :)
107 logical :: total_density_alloc
108 real(real64), pointer, contiguous :: total_density(:)
109 type(energy_t), allocatable :: energy
110
111 type(states_elec_t), pointer :: hf_st
115
116 real(real64), allocatable :: vxc(:, :)
117 real(real64), allocatable :: vtau(:, :)
118 real(real64), allocatable :: axc(:, :, :)
119 real(real64), allocatable :: a_ind(:, :)
120 real(real64), allocatable :: b_ind(:, :)
121 logical :: calc_energy
122 end type v_ks_calc_t
123
124 type v_ks_t
125 private
126 integer, public :: theory_level = -1
127
128 logical, public :: frozen_hxc = .false.
129
130 integer, public :: xc_family = 0
131 integer, public :: xc_flags = 0
132 integer, public :: xc_photon = 0
133 type(xc_t), public :: xc
134 type(xc_photons_t), public :: xc_photons
135 type(xc_oep_t), public :: oep
136 type(xc_oep_photon_t), public :: oep_photon
137 type(xc_ks_inversion_t), public :: ks_inversion
138 type(xc_sic_t), public :: sic
139 type(xc_vdw_t), public :: vdw
140 type(grid_t), pointer, public :: gr
141 type(v_ks_calc_t) :: calc
142 logical :: calculate_current = .false.
143 type(current_t) :: current_calculator
144 logical :: include_td_field = .false.
145 logical, public :: has_photons = .false.
146 logical :: xc_photon_include_hartree = .true.
147
148 real(real64), public :: stress_xc_gga(3, 3)
149 type(photon_mode_t), pointer, public :: pt => null()
150 type(mf_t), public :: pt_mx
151 end type v_ks_t
152
153contains
154
155 ! ---------------------------------------------------------
156 subroutine v_ks_init(ks, namespace, gr, st, ions, mc, space, kpoints)
157 type(v_ks_t), intent(inout) :: ks
158 type(namespace_t), intent(in) :: namespace
159 type(grid_t), target, intent(inout) :: gr
160 type(states_elec_t), intent(in) :: st
161 type(ions_t), intent(inout) :: ions
162 type(multicomm_t), intent(in) :: mc
163 class(space_t), intent(in) :: space
164 type(kpoints_t), intent(in) :: kpoints
165
166 integer :: x_id, c_id, xk_id, ck_id, default, val
167 logical :: parsed_theory_level, using_hartree_fock
168 integer :: pseudo_x_functional, pseudo_c_functional
169 integer :: oep_type
170
171 push_sub(v_ks_init)
172
173 ! We need to parse TheoryLevel and XCFunctional, this is
174 ! complicated because they are interdependent.
175
176 !%Variable TheoryLevel
177 !%Type integer
178 !%Section Hamiltonian
179 !%Description
180 !% The calculations can be run with different "theory levels" that
181 !% control how electrons are simulated. The default is
182 !% <tt>dft</tt>. When hybrid functionals are requested, through
183 !% the <tt>XCFunctional</tt> variable, the default is
184 !% <tt>hartree_fock</tt>.
185 !%Option independent_particles 2
186 !% Particles will be considered as independent, <i>i.e.</i> as non-interacting.
187 !% This mode is mainly used for testing purposes, as the code is usually
188 !% much faster with <tt>independent_particles</tt>.
189 !%Option hartree 1
190 !% Calculation within the Hartree method (experimental). Note that, contrary to popular
191 !% belief, the Hartree potential is self-interaction-free. Therefore, this run
192 !% mode will not yield the same result as <tt>kohn-sham</tt> without exchange-correlation.
193 !%Option hartree_fock 3
194 !% This is the traditional Hartree-Fock scheme. Like the Hartree scheme, it is fully
195 !% self-interaction-free.
196 !%Option kohn_sham 4
197 !% This is the default density-functional theory scheme. Note that you can also use
198 !% hybrid functionals in this scheme, but they will be handled the "DFT" way, <i>i.e.</i>,
199 !% solving the OEP equation.
200 !%Option generalized_kohn_sham 5
201 !% This is similar to the <tt>kohn-sham</tt> scheme, except that this allows for nonlocal operators.
202 !% This is the default mode to run hybrid functionals, meta-GGA functionals, or DFT+U.
203 !% It can be more convenient to use <tt>kohn-sham</tt> DFT within the OEP scheme to get similar (but not the same) results.
204 !% Note that within this scheme you can use a correlation functional, or a hybrid
205 !% functional (see <tt>XCFunctional</tt>). In the latter case, you will be following the
206 !% quantum-chemistry recipe to use hybrids.
207 !%Option rdmft 7
208 !% (Experimental) Reduced Density Matrix functional theory.
209 !%End
211 ks%xc_family = xc_family_none
212 ks%sic%level = sic_none
213 ks%oep%level = oep_level_none
214 ks%oep_photon%level = oep_level_none
215
216 ks%theory_level = kohn_sham_dft
217 parsed_theory_level = .false.
218
219 ! the user knows what he wants, give her that
220 if (parse_is_defined(namespace, 'TheoryLevel')) then
221 call parse_variable(namespace, 'TheoryLevel', kohn_sham_dft, ks%theory_level)
222 if (.not. varinfo_valid_option('TheoryLevel', ks%theory_level)) call messages_input_error(namespace, 'TheoryLevel')
224 parsed_theory_level = .true.
225 end if
227 ! parse the XC functional
229 call get_functional_from_pseudos(pseudo_x_functional, pseudo_c_functional)
231 default = 0
232 if (ks%theory_level == kohn_sham_dft .or. ks%theory_level == generalized_kohn_sham_dft) then
233 default = xc_get_default_functional(space%dim, pseudo_x_functional, pseudo_c_functional)
234 end if
236 if (.not. parse_is_defined(namespace, 'XCFunctional') &
237 .and. (pseudo_x_functional /= pseudo_exchange_any .or. pseudo_c_functional /= pseudo_correlation_any)) then
238 call messages_write('Info: the XCFunctional has been selected to match the pseudopotentials', new_line = .true.)
239 call messages_write(' used in the calculation.')
240 call messages_info(namespace=namespace)
241 end if
243 ! The description of this variable can be found in file src/xc/functionals_list.F90
244 call parse_variable(namespace, 'XCFunctional', default, val)
245
246 ! the first 3 digits of the number indicate the X functional and
247 ! the next 3 the C functional.
248 c_id = val / libxc_c_index
249 x_id = val - c_id * libxc_c_index
250
251 if ((x_id /= pseudo_x_functional .and. pseudo_x_functional /= pseudo_exchange_any) .or. &
252 (c_id /= pseudo_c_functional .and. pseudo_c_functional /= pseudo_exchange_any)) then
253 call messages_write('The XCFunctional that you selected does not match the one used', new_line = .true.)
254 call messages_write('to generate the pseudopotentials.')
255 call messages_warning(namespace=namespace)
256 end if
257
258 ! FIXME: we rarely need this. We should only parse when necessary.
259
260 !%Variable XCKernel
261 !%Type integer
262 !%Default -1
263 !%Section Hamiltonian::XC
264 !%Description
265 !% Defines the exchange-correlation kernel. Only LDA kernels are available currently.
266 !% The options are the same as <tt>XCFunctional</tt>.
267 !% Note: the kernel is only needed for Casida, Sternheimer, or optimal-control calculations.
268 !%Option xc_functional -1
269 !% The same functional defined by <tt>XCFunctional</tt>. By default, this is the case.
270 !%End
271 call parse_variable(namespace, 'XCKernel', -1, val)
272 if (-1 == val) then
273 ck_id = c_id
274 xk_id = x_id
275 else
276 ck_id = val / libxc_c_index
277 xk_id = val - ck_id * libxc_c_index
278 end if
279
280 call messages_obsolete_variable(namespace, 'XFunctional', 'XCFunctional')
281 call messages_obsolete_variable(namespace, 'CFunctional', 'XCFunctional')
282
283 !%Variable XCPhotonFunctional
284 !%Type integer
285 !%Default 0
286 !%Section Hamiltonian::XC
287 !%Description
288 !% Defines the exchange and correlation functionals to be used for the QEDFT
289 !% description of the electron-photon system.
290 !%Option none 0
291 !% No functional is used
292 !%Option photon_x_lda 10
293 !% Exchange-only local density approcimation
294 !%Option photon_xc_lda 11
295 !% Exchange-correlation local density approcimation
296 !%Option photon_x_wfn 20
297 !% Exchange-only based on wave functions
298 !%Option photon_xc_wfn 21
299 !% Exchange-correlation based on wave functions
300 !%End
301
302 call parse_variable(namespace, 'XCPhotonFunctional', option__xcphotonfunctional__none, ks%xc_photon)
303
304 !%Variable XCPhotonIncludeHartree
305 !%Type logical
306 !%Default yes
307 !%Section Hamiltonian::XC
308 !%Description
309 !% Use the Hartree potential and energy in calculations
310 !%End
311
312 call parse_variable(namespace, 'XCPhotonIncludeHartree', .true., ks%xc_photon_include_hartree)
313
314 if (.not. ks%xc_photon_include_hartree) then
315 call messages_write('turn off hartree potential and energy')
316 call messages_warning(namespace=namespace)
317 end if
318
319 ! initialize XC modules
320
321 ! This is a bit ugly, theory_level might not be generalized KS or HF now
322 ! but it might become generalized KS or HF later. This is safe because it
323 ! becomes generalized KS in the cases where the functional is hybrid
324 ! and the ifs inside check for both conditions.
325 using_hartree_fock = (ks%theory_level == hartree_fock) &
326 .or. (ks%theory_level == generalized_kohn_sham_dft .and. family_is_hybrid(ks%xc))
327 call xc_init(ks%xc, namespace, space%dim, space%periodic_dim, st%qtot, &
328 x_id, c_id, xk_id, ck_id, hartree_fock = using_hartree_fock, ispin=st%d%ispin)
329
330 ks%xc_family = ks%xc%family
331 ks%xc_flags = ks%xc%flags
332
333 if (.not. parsed_theory_level) then
334 default = kohn_sham_dft
335
336 ! the functional is a hybrid, use Hartree-Fock as theory level by default
337 if (family_is_hybrid(ks%xc) .or. family_is_mgga_with_exc(ks%xc)) then
339 end if
340
341 ! In principle we do not need to parse. However we do it for consistency
342 call parse_variable(namespace, 'TheoryLevel', default, ks%theory_level)
343 if (.not. varinfo_valid_option('TheoryLevel', ks%theory_level)) call messages_input_error(namespace, 'TheoryLevel')
344
345 end if
346
347 ! In case we need OEP, we need to find which type of OEP it is
348 oep_type = -1
349 if (family_is_mgga_with_exc(ks%xc)) then
350 call messages_experimental('MGGA energy functionals')
351
352 if (ks%theory_level == kohn_sham_dft) then
353 call messages_experimental("MGGA within the Kohn-Sham scheme")
354 ks%xc_family = ior(ks%xc_family, xc_family_oep)
355 oep_type = oep_type_mgga
356 end if
357 end if
358
359 call messages_obsolete_variable(namespace, 'NonInteractingElectrons', 'TheoryLevel')
360 call messages_obsolete_variable(namespace, 'HartreeFock', 'TheoryLevel')
361
362 ! Due to how the code is made, we need to set this to have theory level other than DFT
363 ! correct...
364 ks%sic%amaldi_factor = m_one
365
366 select case (ks%theory_level)
368
369 case (hartree)
370 call messages_experimental("Hartree theory level")
371 if (space%periodic_dim == space%dim) then
372 call messages_experimental("Hartree in fully periodic system")
373 end if
374 if (kpoints%full%npoints > 1) then
375 call messages_not_implemented("Hartree with k-points", namespace=namespace)
376 end if
377
378 case (hartree_fock)
379 if (kpoints%full%npoints > 1) then
380 call messages_experimental("Hartree-Fock with k-points")
381 end if
382
384 if (kpoints%full%npoints > 1 .and. family_is_hybrid(ks%xc)) then
385 call messages_experimental("Hybrid functionals with k-points")
386 end if
387
388 case (rdmft)
389 call messages_experimental('RDMFT theory level')
390
391 case (kohn_sham_dft)
392
393 ! check for SIC
394 if (bitand(ks%xc_family, xc_family_lda + xc_family_gga) /= 0) then
395 call xc_sic_init(ks%sic, namespace, gr, st, mc, space)
396 end if
397
398 if (bitand(ks%xc_family, xc_family_oep) /= 0) then
399 select case (ks%xc%functional(func_x,1)%id)
400 case (xc_oep_x_slater)
401 if (kpoints%reduced%npoints > 1) then
402 call messages_not_implemented("Slater with k-points", namespace=namespace)
403 end if
404 ks%oep%level = oep_level_none
405 case (xc_oep_x_fbe)
406 if (kpoints%reduced%npoints > 1) then
407 call messages_not_implemented("FBE functional with k-points", namespace=namespace)
408 end if
409 ks%oep%level = oep_level_none
410 case default
411 if((.not. ks%has_photons) .or. (ks%xc_photon /= 0)) then
412 if(oep_type == -1) then ! Else we have a MGGA
413 oep_type = oep_type_exx
414 end if
415 call xc_oep_init(ks%oep, namespace, gr, st, mc, space, oep_type)
416 end if
417 end select
418 else
419 ks%oep%level = oep_level_none
420 end if
421
422 if (bitand(ks%xc_family, xc_family_ks_inversion) /= 0) then
423 call xc_ks_inversion_init(ks%ks_inversion, namespace, gr, ions, st, ks%xc, mc, space, kpoints)
424 end if
425
426 end select
427
428 if (ks%theory_level /= kohn_sham_dft .and. parse_is_defined(namespace, "SICCorrection")) then
429 message(1) = "SICCorrection can only be used with Kohn-Sham DFT"
430 call messages_fatal(1, namespace=namespace)
431 end if
432
433 if (st%d%ispin == spinors) then
434 if (bitand(ks%xc_family, xc_family_mgga + xc_family_hyb_mgga) /= 0) then
435 call messages_not_implemented("MGGA with spinors", namespace=namespace)
436 end if
437 end if
438
439 ks%frozen_hxc = .false.
440
441 call v_ks_write_info(ks, namespace=namespace)
442
443 ks%gr => gr
444 ks%calc%calculating = .false.
445
446 !The value of ks%calculate_current is set to false or true by Output
447 call current_init(ks%current_calculator, namespace)
448
449 call ks%vdw%init(namespace, space, gr, ks%xc, ions, x_id, c_id)
450
451 if (ks%xc_photon /= 0) then
452 ! initilize the photon free variables
453 call ks%xc_photons%init(namespace, ks%xc_photon , space, gr, st)
454 ! remornalize the electron mass due to light-matter interaction; here we only deal with it in free space
455 ks%oep_photon%level = oep_level_none
456 end if
457
458
459 pop_sub(v_ks_init)
460
461 contains
462
464 subroutine get_functional_from_pseudos(x_functional, c_functional)
465 integer, intent(out) :: x_functional
466 integer, intent(out) :: c_functional
467
468 integer :: xf, cf, ispecies
469 logical :: warned_inconsistent
470
471 x_functional = pseudo_exchange_any
472 c_functional = pseudo_correlation_any
473
474 warned_inconsistent = .false.
475 do ispecies = 1, ions%nspecies
476 select type(spec=>ions%species(ispecies)%s)
477 class is(pseudopotential_t)
478 xf = spec%x_functional()
479 cf = spec%c_functional()
480
481 if (xf == pseudo_exchange_unknown .or. cf == pseudo_correlation_unknown) then
482 call messages_write("Unknown XC functional for species '"//trim(ions%species(ispecies)%s%get_label())//"'")
483 call messages_warning(namespace=namespace)
484 cycle
485 end if
486
487 if (x_functional == pseudo_exchange_any) then
488 x_functional = xf
489 else
490 if (xf /= x_functional .and. .not. warned_inconsistent) then
491 call messages_write('Inconsistent XC functional detected between species')
492 call messages_warning(namespace=namespace)
493 warned_inconsistent = .true.
494 end if
495 end if
496
497 if (c_functional == pseudo_correlation_any) then
498 c_functional = cf
499 else
500 if (cf /= c_functional .and. .not. warned_inconsistent) then
501 call messages_write('Inconsistent XC functional detected between species')
502 call messages_warning(namespace=namespace)
503 warned_inconsistent = .true.
504 end if
505 end if
506
507 class default
510 end select
511
512 end do
513
514 assert(x_functional /= pseudo_exchange_unknown)
515 assert(c_functional /= pseudo_correlation_unknown)
516
517 end subroutine get_functional_from_pseudos
518 end subroutine v_ks_init
519 ! ---------------------------------------------------------
520
521 ! ---------------------------------------------------------
522 subroutine v_ks_end(ks)
523 type(v_ks_t), intent(inout) :: ks
524
525 push_sub(v_ks_end)
526
527 call ks%vdw%end()
528
529 select case (ks%theory_level)
530 case (kohn_sham_dft)
531 if (bitand(ks%xc_family, xc_family_ks_inversion) /= 0) then
532 call xc_ks_inversion_end(ks%ks_inversion)
533 end if
534 if (bitand(ks%xc_family, xc_family_oep) /= 0) then
535 call xc_oep_end(ks%oep)
536 end if
537 call xc_end(ks%xc)
539 call xc_end(ks%xc)
540 end select
541
542 call xc_sic_end(ks%sic)
543
544 if (ks%xc_photon /= 0) then
545 call ks%xc_photons%end()
546 end if
547
548 pop_sub(v_ks_end)
549 end subroutine v_ks_end
550 ! ---------------------------------------------------------
551
552
553 ! ---------------------------------------------------------
554 subroutine v_ks_write_info(ks, iunit, namespace)
555 type(v_ks_t), intent(in) :: ks
556 integer, optional, intent(in) :: iunit
557 type(namespace_t), optional, intent(in) :: namespace
558
559 push_sub(v_ks_write_info)
560
561 call messages_print_with_emphasis(msg="Theory Level", iunit=iunit, namespace=namespace)
562 call messages_print_var_option("TheoryLevel", ks%theory_level, iunit=iunit, namespace=namespace)
563
564 select case (ks%theory_level)
566 call messages_info(iunit=iunit, namespace=namespace)
567 call xc_write_info(ks%xc, iunit, namespace)
568
569 case (kohn_sham_dft)
570 call messages_info(iunit=iunit, namespace=namespace)
571 call xc_write_info(ks%xc, iunit, namespace)
572
573 call messages_info(iunit=iunit, namespace=namespace)
574
575 call xc_sic_write_info(ks%sic, iunit, namespace)
576 call xc_oep_write_info(ks%oep, iunit, namespace)
577 call xc_ks_inversion_write_info(ks%ks_inversion, iunit, namespace)
578
579 end select
580
581 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
582
583 pop_sub(v_ks_write_info)
584 end subroutine v_ks_write_info
585 ! ---------------------------------------------------------
586
587
588 !----------------------------------------------------------
589 subroutine v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
590 type(namespace_t), intent(in) :: namespace
591 type(electron_space_t), intent(in) :: space
592 type(grid_t), intent(in) :: gr
593 type(ions_t), intent(in) :: ions
594 type(partner_list_t), intent(in) :: ext_partners
595 type(states_elec_t), intent(inout) :: st
596 type(v_ks_t), intent(inout) :: ks
597 type(hamiltonian_elec_t), intent(inout) :: hm
598 logical, optional, intent(in) :: calc_eigenval
599 logical, optional, intent(in) :: calc_current
600
601 integer, allocatable :: ind(:)
602 integer :: ist, ik
603 real(real64), allocatable :: copy_occ(:)
604 logical :: calc_eigenval_
605 logical :: calc_current_
606
607 push_sub(v_ks_h_setup)
608
609 calc_eigenval_ = optional_default(calc_eigenval, .true.)
610 calc_current_ = optional_default(calc_current, .true.)
611 call states_elec_fermi(st, namespace, gr)
612 call density_calc(st, gr, st%rho)
613 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval = calc_eigenval_, calc_current = calc_current_) ! get potentials
614
615 if (st%restart_reorder_occs .and. .not. st%fromScratch) then
616 message(1) = "Reordering occupations for restart."
617 call messages_info(1, namespace=namespace)
618
619 safe_allocate(ind(1:st%nst))
620 safe_allocate(copy_occ(1:st%nst))
621
622 do ik = 1, st%nik
623 call sort(st%eigenval(:, ik), ind)
624 copy_occ(1:st%nst) = st%occ(1:st%nst, ik)
625 do ist = 1, st%nst
626 st%occ(ist, ik) = copy_occ(ind(ist))
627 end do
628 end do
629
630 safe_deallocate_a(ind)
631 safe_deallocate_a(copy_occ)
632 end if
633
634 if (calc_eigenval_) call states_elec_fermi(st, namespace, gr) ! occupations
635 call energy_calc_total(namespace, space, hm, gr, st, ext_partners)
636
637 pop_sub(v_ks_h_setup)
638 end subroutine v_ks_h_setup
639
640 ! ---------------------------------------------------------
641 subroutine v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
642 calc_eigenval, time, calc_energy, calc_current, force_semilocal)
643 type(v_ks_t), intent(inout) :: ks
644 type(namespace_t), intent(in) :: namespace
645 type(electron_space_t), intent(in) :: space
646 type(hamiltonian_elec_t), intent(inout) :: hm
647 type(states_elec_t), intent(inout) :: st
648 type(ions_t), intent(in) :: ions
649 type(partner_list_t), intent(in) :: ext_partners
650 logical, optional, intent(in) :: calc_eigenval
651 real(real64), optional, intent(in) :: time
652 logical, optional, intent(in) :: calc_energy
653 logical, optional, intent(in) :: calc_current
654 logical, optional, intent(in) :: force_semilocal
655
656 logical :: calc_current_
657
658 push_sub(v_ks_calc)
659
660 calc_current_ = optional_default(calc_current, .true.)
661
662 call v_ks_calc_start(ks, namespace, space, hm, st, ions, hm%kpoints%latt, ext_partners, time, &
663 calc_energy, calc_current_, force_semilocal=force_semilocal)
664 call v_ks_calc_finish(ks, hm, namespace, space, hm%kpoints%latt, st, &
665 ext_partners, force_semilocal=force_semilocal)
666
667 if (optional_default(calc_eigenval, .false.)) then
668 call energy_calc_eigenvalues(namespace, hm, ks%gr%der, st)
669 end if
670
671 ! Update the magnetic constrain
672 call magnetic_constrain_update(hm%magnetic_constrain, ks%gr, st%d, space, hm%kpoints%latt, ions%pos, st%rho)
673 ! We add the potential to vxc, as this way the potential gets mixed together with vxc
674 ! While this is not ideal, this is a simple practical solution
675 if (hm%magnetic_constrain%level /= constrain_none) then
676 call lalg_axpy(ks%gr%np, st%d%nspin, m_one, hm%magnetic_constrain%pot, hm%ks_pot%vhxc)
677 end if
678
679 pop_sub(v_ks_calc)
680 end subroutine v_ks_calc
681
682 ! ---------------------------------------------------------
683
688 subroutine v_ks_calc_start(ks, namespace, space, hm, st, ions, latt, ext_partners, time, &
689 calc_energy, calc_current, force_semilocal)
690 type(v_ks_t), target, intent(inout) :: ks
691 type(namespace_t), intent(in) :: namespace
692 class(space_t), intent(in) :: space
693 type(hamiltonian_elec_t), target, intent(in) :: hm
694 type(states_elec_t), target, intent(inout) :: st
695 type(ions_t), intent(in) :: ions
696 type(lattice_vectors_t), intent(in) :: latt
697 type(partner_list_t), intent(in) :: ext_partners
698 real(real64), optional, intent(in) :: time
699 logical, optional, intent(in) :: calc_energy
700 logical, optional, intent(in) :: calc_current
701 logical, optional, intent(in) :: force_semilocal
702
703 logical :: calc_current_
704
705 push_sub(v_ks_calc_start)
706
707 calc_current_ = optional_default(calc_current, .true.) &
708 .and. ks%calculate_current &
709 .and. states_are_complex(st) &
711
712 call profiling_in("KOHN_SHAM_CALC")
713
714 assert(.not. ks%calc%calculating)
715 ks%calc%calculating = .true.
716
717 write(message(1), '(a)') 'Debug: Calculating Kohn-Sham potential.'
718 call messages_info(1, namespace=namespace, debug_only=.true.)
719
720 ks%calc%time_present = present(time)
721 ks%calc%time = optional_default(time, m_zero)
722
723 ks%calc%calc_energy = optional_default(calc_energy, .true.)
724
725 ! If the Hxc term is frozen, there is nothing more to do (WARNING: MISSING ks%calc%energy%intnvxc)
726 if (ks%frozen_hxc) then
727 if (calc_current_) then
728 call states_elec_allocate_current(st, space, ks%gr)
729 call current_calculate(ks%current_calculator, namespace, ks%gr, hm, space, st)
730 end if
731
732 call profiling_out("KOHN_SHAM_CALC")
733 pop_sub(v_ks_calc_start)
734 return
735 end if
736
737 allocate(ks%calc%energy)
738
739 call energy_copy(hm%energy, ks%calc%energy)
740
741 ks%calc%energy%intnvxc = m_zero
742
743 nullify(ks%calc%total_density)
744
745 if (ks%theory_level /= independent_particles .and. abs(ks%sic%amaldi_factor) > m_epsilon) then
746
747 call calculate_density()
748
749 if (poisson_is_async(hm%psolver)) then
750 call dpoisson_solve_start(hm%psolver, ks%calc%total_density)
751 end if
752
753 if (ks%theory_level /= hartree .and. ks%theory_level /= rdmft) call v_a_xc(hm, force_semilocal)
754 else
755 ks%calc%total_density_alloc = .false.
756 end if
757
758 if (calc_current_) then
759 call states_elec_allocate_current(st, space, ks%gr)
760 call current_calculate(ks%current_calculator, namespace, ks%gr, hm, space, st)
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
801
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 if (optional_default(force_semilocal, .false.)) then
927 call profiling_out("XC")
928 pop_sub(v_ks_calc_start.v_a_xc)
929 return
930 end if
931
932 ! ADSIC correction
933 if (ks%sic%level == sic_adsic) then
934 if (family_is_mgga(hm%xc%family)) then
935 call messages_not_implemented('ADSIC with MGGAs', namespace=namespace)
936 end if
937 if (ks%calc%calc_energy) then
938 call xc_sic_calc_adsic(ks%sic, namespace, space, ks%gr, st, hm, ks%xc, ks%calc%density, &
939 ks%calc%vxc, ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation)
940 else
941 call xc_sic_calc_adsic(ks%sic, namespace, space, ks%gr, st, hm, ks%xc, ks%calc%density, &
942 ks%calc%vxc)
943 end if
944 end if
945 !PZ SIC is done in the finish routine as OEP full needs to update the Hamiltonian
946
947 if (ks%theory_level == kohn_sham_dft) then
948 ! The OEP family has to be handled specially
949 ! Note that OEP is done in the finish state, as it requires updating the Hamiltonian and needs the new Hartre and vxc term
950 if (bitand(ks%xc_family, xc_family_oep) /= 0 .or. family_is_mgga_with_exc(ks%xc)) then
951
952 if (ks%xc%functional(func_x,1)%id == xc_oep_x_slater) then
953 call x_slater_calc(namespace, ks%gr, space, hm%exxop, st, hm%kpoints, ks%calc%energy%exchange, &
954 vxc = ks%calc%vxc)
955 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
956 call x_fbe_calc(ks%xc%functional(func_x,1)%id, namespace, hm%psolver, ks%gr, st, space, &
957 ks%calc%energy%exchange, vxc = ks%calc%vxc)
958
959 else if (ks%xc%functional(func_c,1)%id == xc_lda_c_fbe_sl) then
960
961 call fbe_c_lda_sl(namespace, hm%psolver, ks%gr, st, space, &
962 ks%calc%energy%correlation, vxc = ks%calc%vxc)
963
964 end if
965
966 end if
967
968 if (bitand(ks%xc_family, xc_family_ks_inversion) /= 0) then
969 ! Also treat KS inversion separately (not part of libxc)
970 call xc_ks_inversion_calc(ks%ks_inversion, namespace, space, ks%gr, hm, ext_partners, st, vxc = ks%calc%vxc, &
971 time = ks%calc%time)
972 end if
973
974 ! compute the photon-free photon exchange potential and energy
975 if (ks%xc_photon /= 0) then
976
977 call ks%xc_photons%v_ks(namespace, ks%calc%total_density, ks%gr, space, hm%psolver, hm%ep, st)
978
979 ! add the photon-free px potential into the xc potential
980 do ispin = 1, hm%d%spin_channels
981 call lalg_axpy(ks%gr%np, m_one, ks%xc_photons%vpx(1:ks%gr%np), ks%calc%vxc(1:ks%gr%np, ispin) )
982 end do
983
984 ! photon-exchange energy
985 ks%calc%energy%photon_exchange = ks%xc_photons%ex
986 end if
987
988 end if
989
990 call ks%vdw%calc(namespace, space, latt, ions%atom, ions%natoms, ions%pos, &
991 ks%gr, st, ks%calc%energy%vdw, ks%calc%vxc)
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
1101
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_omega, ks%xc%cam_alpha, ks%xc%cam_beta, ks%calc%hf_st)
1209 end if
1210 case (hartree_fock)
1211 call exchange_operator_reinit(hm%exxop, ks%xc%cam_omega, ks%xc%cam_alpha, ks%xc%cam_beta, ks%calc%hf_st)
1212 case (hartree, rdmft)
1213 call exchange_operator_reinit(hm%exxop, m_zero, m_one, m_zero, 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, exx_energy)
1224 call dexchange_operator_ace(hm%exxop, namespace, ks%gr, ks%calc%hf_st, xst)
1225 else
1226 call zexchange_operator_compute_potentials(hm%exxop, namespace, space, ks%gr, &
1227 ks%calc%hf_st, xst, hm%kpoints, exx_energy)
1228 if (hm%phase%is_allocated()) then
1229 call zexchange_operator_ace(hm%exxop, namespace, ks%gr, ks%calc%hf_st, xst, hm%phase)
1230 else
1231 call zexchange_operator_ace(hm%exxop, namespace, ks%gr, ks%calc%hf_st, xst)
1232 end if
1233 end if
1234 call states_elec_end(xst)
1235 exx_energy = exx_energy + hm%exxop%singul%energy
1236 end if
1237
1238 ! Add the energy only the ACE case. In the non-ACE case, the singularity energy is added in energy_calc.F90
1239 select case (ks%theory_level)
1241 if (family_is_hybrid(ks%xc)) then
1242 hm%energy%exchange_hf = hm%energy%exchange_hf + exx_energy
1243 end if
1244 case (hartree_fock)
1245 hm%energy%exchange_hf = hm%energy%exchange_hf + exx_energy
1246 end select
1247 else
1248 ! If we ask for semilocal, we deactivate the exchange operator entirely
1249 call exchange_operator_reinit(hm%exxop, m_zero, m_zero, m_zero, ks%calc%hf_st)
1250 end if
1251 end if
1252
1253 end if
1254
1255 ! Because of the intent(in) in v_ks_calc_start, we need to update the parameters of hybrids for OEP
1256 ! here
1257 if (ks%theory_level == kohn_sham_dft .and. bitand(ks%xc_family, xc_family_oep) /= 0) then
1258 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
1259 call exchange_operator_reinit(hm%exxop, ks%xc%cam_omega, ks%xc%cam_alpha, ks%xc%cam_beta)
1260 end if
1261 end if
1262
1263 if (ks%has_photons .and. (ks%xc_photon == 0)) then
1264 if (associated(ks%pt_mx%vmf)) then
1265 forall(ip = 1:ks%gr%np) hm%ks_pot%vhxc(ip, 1) = hm%ks_pot%vhxc(ip, 1) + ks%pt_mx%vmf(ip)
1266 if (hm%d%ispin > unpolarized) then
1267 forall(ip = 1:ks%gr%np) hm%ks_pot%vhxc(ip, 2) = hm%ks_pot%vhxc(ip, 2) + ks%pt_mx%vmf(ip)
1268 end if
1269 end if
1270 hm%ep%photon_forces(1:space%dim) = ks%pt_mx%fmf(1:space%dim)
1271 end if
1272
1273 if (ks%vdw%vdw_correction /= option__vdwcorrection__none) then
1274 hm%ep%vdw_forces(:, :) = ks%vdw%forces(:, :)
1275 hm%ep%vdw_stress = ks%vdw%stress
1276 safe_deallocate_a(ks%vdw%forces)
1277 else
1278 hm%ep%vdw_forces = 0.0_real64
1279 end if
1280
1281 if (ks%calc%time_present .or. hm%time_zero) then
1282 call hm%update(ks%gr, namespace, space, ext_partners, time = ks%calc%time)
1283 else
1284 call hamiltonian_elec_update_pot(hm, ks%gr)
1285 end if
1286
1287
1288 safe_deallocate_a(ks%calc%density)
1289 if (ks%calc%total_density_alloc) then
1290 safe_deallocate_p(ks%calc%total_density)
1291 end if
1292 nullify(ks%calc%total_density)
1293
1294 pop_sub(v_ks_calc_finish)
1295 end subroutine v_ks_calc_finish
1296
1297 ! ---------------------------------------------------------
1298 !
1302 !
1303 subroutine v_ks_hartree(namespace, ks, space, hm, ext_partners)
1304 type(namespace_t), intent(in) :: namespace
1305 type(v_ks_t), intent(inout) :: ks
1306 class(space_t), intent(in) :: space
1307 type(hamiltonian_elec_t), intent(inout) :: hm
1308 type(partner_list_t), intent(in) :: ext_partners
1309
1310 push_sub(v_ks_hartree)
1311
1312 if (.not. poisson_is_async(hm%psolver)) then
1313 ! solve the Poisson equation
1314 call dpoisson_solve(hm%psolver, namespace, hm%ks_pot%vhartree, ks%calc%total_density, reset=.false.)
1315 else
1316 ! The calculation was started by v_ks_calc_start.
1317 call dpoisson_solve_finish(hm%psolver, hm%ks_pot%vhartree)
1318 end if
1319
1320 if (ks%calc%calc_energy) then
1321 ! Get the Hartree energy
1322 hm%energy%hartree = m_half*dmf_dotp(ks%gr, ks%calc%total_density, hm%ks_pot%vhartree)
1323 end if
1324
1326 if(ks%calc%time_present) then
1327 if(hamiltonian_elec_has_kick(hm)) then
1328 call pcm_hartree_potential(hm%pcm, space, ks%gr, hm%psolver, ext_partners, hm%ks_pot%vhartree, &
1329 ks%calc%total_density, hm%energy%pcm_corr, kick=hm%kick, time=ks%calc%time)
1330 else
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, time=ks%calc%time)
1333 end if
1334 else
1335 if(hamiltonian_elec_has_kick(hm)) then
1336 call pcm_hartree_potential(hm%pcm, space, ks%gr, hm%psolver, ext_partners, hm%ks_pot%vhartree, &
1337 ks%calc%total_density, hm%energy%pcm_corr, kick=hm%kick)
1338 else
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)
1341 end if
1342 end if
1343
1344 pop_sub(v_ks_hartree)
1345 end subroutine v_ks_hartree
1346 ! ---------------------------------------------------------
1347
1348
1349 ! ---------------------------------------------------------
1350 subroutine v_ks_freeze_hxc(ks)
1351 type(v_ks_t), intent(inout) :: ks
1352
1353 push_sub(v_ks_freeze_hxc)
1354
1355 ks%frozen_hxc = .true.
1356
1357 pop_sub(v_ks_freeze_hxc)
1358 end subroutine v_ks_freeze_hxc
1359 ! ---------------------------------------------------------
1360
1361 subroutine v_ks_calculate_current(this, calc_cur)
1362 type(v_ks_t), intent(inout) :: this
1363 logical, intent(in) :: calc_cur
1364
1365 push_sub(v_ks_calculate_current)
1366
1367 this%calculate_current = calc_cur
1368
1369 pop_sub(v_ks_calculate_current)
1370 end subroutine v_ks_calculate_current
1371
1373 subroutine v_ks_update_dftu_energy(ks, namespace, hm, st, int_dft_u)
1374 type(v_ks_t), intent(inout) :: ks
1375 type(hamiltonian_elec_t), intent(in) :: hm
1376 type(namespace_t), intent(in) :: namespace
1377 type(states_elec_t), intent(inout) :: st
1378 real(real64), intent(out) :: int_dft_u
1379
1380 if (hm%lda_u_level == dft_u_none) return
1381
1382 push_sub(v_ks_update_dftu_energy)
1383
1384 if (states_are_real(st)) then
1385 int_dft_u = denergy_calc_electronic(namespace, hm, ks%gr%der, st, terms = term_dft_u)
1386 else
1387 int_dft_u = zenergy_calc_electronic(namespace, hm, ks%gr%der, st, terms = term_dft_u)
1388 end if
1389
1391 end subroutine v_ks_update_dftu_energy
1392end module v_ks_oct_m
1393
1394!! Local Variables:
1395!! mode: f90
1396!! coding: utf-8
1397!! 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:401
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:369
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:1397
subroutine, public v_ks_calc_finish(ks, hm, namespace, space, latt, st, ext_partners, force_semilocal)
Definition: v_ks.F90:1100
subroutine, public v_ks_freeze_hxc(ks)
Definition: v_ks.F90:1444
subroutine, public v_ks_end(ks)
Definition: v_ks.F90:616
subroutine, public v_ks_calculate_current(this, calc_cur)
Definition: v_ks.F90:1455
subroutine, public v_ks_write_info(ks, iunit, namespace)
Definition: v_ks.F90:648
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:1467
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:783
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:736
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
Definition: v_ks.F90:683
subroutine, public v_ks_init(ks, namespace, gr, st, ions, mc, space, kpoints)
Definition: v_ks.F90:250
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 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:114
subroutine, public xc_write_info(xcs, iunit, namespace)
Definition: xc.F90:221
subroutine, public xc_init(xcs, namespace, ndim, periodic_dim, nel, x_id, c_id, xk_id, ck_id, hartree_fock, ispin)
Definition: xc.F90:251
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
Definition: xc.F90:582
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:595
subroutine, public xc_end(xcs)
Definition: xc.F90:520
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:757
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
Definition: xc.F90:610
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:2417
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:2301
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:558
subroutine v_a_xc(hm, force_semilocal)
Definition: v_ks.F90:946
subroutine calculate_density()
Definition: v_ks.F90:909