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