Octopus
hamiltonian_elec.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2020 M. Marques, A. Castro, A. Rubio, G. Bertsch,
2!! N. Tancogne-Dejean, M. Lueders
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
24 use accel_oct_m
26 use batch_oct_m
29 use comm_oct_m
30 use debug_oct_m
33 use energy_oct_m
38 use epot_oct_m
41 use global_oct_m
42 use grid_oct_m
46 use io_oct_m
47 use ions_oct_m
48 use kick_oct_m
49 use, intrinsic :: iso_fortran_env
53 use lasers_oct_m
55 use lda_u_oct_m
58 use math_oct_m
59 use mesh_oct_m
62 use mpi_oct_m
66 use nlcc_oct_m
70 use parser_oct_m
75 use pcm_oct_m
76 use phase_oct_m
79 use space_oct_m
87 use types_oct_m
88 use unit_oct_m
92 use xc_oct_m
93 use xc_cam_oct_m
94 use xc_f03_lib_m
98 use zora_oct_m
99
100 implicit none
101
102 private
103 public :: &
114 dvmask, &
115 zvmask, &
130
131
132 type, extends(hamiltonian_abst_t) :: hamiltonian_elec_t
133 ! Components are public by default
134
137 type(space_t), private :: space
138 type(states_elec_dim_t) :: d
139 type(hamiltonian_elec_base_t) :: hm_base
140 type(phase_t) :: phase
141 type(energy_t), allocatable :: energy
142 type(absorbing_boundaries_t) :: abs_boundaries
143 type(ks_potential_t) :: ks_pot
144 real(real64), allocatable :: vberry(:,:)
145
146 type(derivatives_t), pointer, private :: der
147
148 type(nonlocal_pseudopotential_t) :: vnl
149
150 type(ions_t), pointer :: ions
151 logical, private :: owns_ions = .false.
152 logical, private :: is_copy_snapshot = .false.
153 real(real64) :: exx_coef
154
155 type(poisson_t) :: psolver
156
158 logical :: self_induced_magnetic
159 real(real64), allocatable :: a_ind(:, :)
160 real(real64), allocatable :: b_ind(:, :)
161
162 integer :: theory_level
163 type(xc_t), pointer :: xc
164 type(xc_photons_t), pointer :: xc_photons
165
166 type(epot_t) :: ep
167 type(pcm_t) :: pcm
168
170 logical, private :: adjoint
171
173 real(real64), private :: mass
174
176 logical, private :: inh_term
177 type(states_elec_t) :: inh_st
178
181 type(oct_exchange_t) :: oct_exchange
182
183 type(scissor_t) :: scissor
184
185 real(real64) :: current_time
186 logical, private :: is_applied_packed
187
189 type(lda_u_t) :: lda_u
190 integer :: lda_u_level
191
192 logical, public :: time_zero
193
194 type(exchange_operator_t), public :: exxop
195
196 type(kpoints_t), pointer, public :: kpoints => null()
197
198 type(partner_list_t) :: external_potentials
199 real(real64), allocatable, public :: v_ext_pot(:)
200 real(real64), allocatable, public :: v_static(:)
201
202 type(ion_electron_local_potential_t) :: v_ie_loc
203 type(nlcc_t) :: nlcc
204
205 type(magnetic_constrain_t) :: magnetic_constrain
206
208 type(kick_t) :: kick
209
211 type(mxll_coupling_t) :: mxll
212 type(zora_t), pointer :: zora => null()
213
214 contains
215 procedure :: update => hamiltonian_elec_update
216 procedure :: apply_packed => hamiltonian_elec_apply_packed
217 procedure :: update_span => hamiltonian_elec_span
218 procedure :: dapply => dhamiltonian_elec_apply
219 procedure :: zapply => zhamiltonian_elec_apply
220 procedure :: is_hermitian => hamiltonian_elec_hermitian
221 procedure :: needs_mgga_term => hamiltonian_elec_needs_mgga_term
222 procedure :: set_mass => hamiltonian_elec_set_mass
223 end type hamiltonian_elec_t
224
225 integer, public, parameter :: &
226 LENGTH = 1, &
228
229
230contains
231
232 ! ---------------------------------------------------------
233 subroutine hamiltonian_elec_init(hm, namespace, space, gr, ions, ext_partners, st, theory_level, xc, &
234 mc, kpoints, need_exchange, xc_photons)
235 type(hamiltonian_elec_t), target, intent(inout) :: hm
236 type(namespace_t), intent(in) :: namespace
237 class(space_t), intent(in) :: space
238 type(grid_t), target, intent(inout) :: gr
239 type(ions_t), target, intent(inout) :: ions
240 type(partner_list_t), intent(inout) :: ext_partners
241 type(states_elec_t), target, intent(inout) :: st
242 integer, intent(in) :: theory_level
243 type(xc_t), target, intent(in) :: xc
244 type(multicomm_t), intent(in) :: mc
245 type(kpoints_t), target, intent(in) :: kpoints
246 logical, optional, intent(in) :: need_exchange
247 type(xc_photons_t), optional, target, intent(in) :: xc_photons
249
250 logical :: need_exchange_
251 real(real64) :: rashba_coupling
252
254 call profiling_in('HAMILTONIAN_ELEC_INIT')
256 ! make a couple of local copies
257 hm%space = space
258 hm%theory_level = theory_level
259 call states_elec_dim_copy(hm%d, st%d)
260
261 hm%kpoints => kpoints
263 !%Variable ParticleMass
264 !%Type float
265 !%Default 1.0
266 !%Section Hamiltonian
267 !%Description
268 !% It is possible to make calculations for a particle with a mass
269 !% different from one (atomic unit of mass, or mass of the electron).
270 !% This is useful to describe non-electronic systems, or for
271 !% esoteric purposes.
272 !%End
273 call parse_variable(namespace, 'ParticleMass', m_one, hm%mass)
274
275 !%Variable RashbaSpinOrbitCoupling
276 !%Type float
277 !%Default 0.0
278 !%Section Hamiltonian
279 !%Description
280 !% (Experimental.) For systems described in 2D (electrons confined to 2D in semiconductor structures), one
281 !% may add the Bychkov-Rashba spin-orbit coupling term [Bychkov and Rashba, <i>J. Phys. C: Solid
282 !% State Phys.</i> <b>17</b>, 6031 (1984)]. This variable determines the strength
283 !% of this perturbation, and has dimensions of energy times length.
284 !%End
285 call parse_variable(namespace, 'RashbaSpinOrbitCoupling', m_zero, rashba_coupling, units_inp%energy*units_inp%length)
286 if (parse_is_defined(namespace, 'RashbaSpinOrbitCoupling')) then
287 if (space%dim /= 2) then
288 write(message(1),'(a)') 'Rashba spin-orbit coupling can only be used for two-dimensional systems.'
289 call messages_fatal(1, namespace=namespace)
290 end if
291 call messages_experimental('RashbaSpinOrbitCoupling', namespace=namespace)
292 end if
294 call hm%hm_base%init(hm%d%nspin, hm%mass, rashba_coupling)
295 call hm%vnl%init()
296
297 assert(associated(gr%der%lapl))
298 hm%hm_base%kinetic => gr%der%lapl
299
300 safe_allocate(hm%energy)
301
302 !Keep pointers to derivatives, geometry and xc
303 hm%der => gr%der
304 hm%ions => ions
305 hm%owns_ions = .false.
306 hm%is_copy_snapshot = .false.
307 hm%xc => xc
308
309 if(present(xc_photons)) then
310 hm%xc_photons => xc_photons
311 else
312 hm%xc_photons => null()
313 end if
315 ! allocate potentials and density of the cores
316 ! In the case of spinors, vxc_11 = hm%vxc(:, 1), vxc_22 = hm%vxc(:, 2), Re(vxc_12) = hm%vxc(:. 3);
317 ! Im(vxc_12) = hm%vxc(:, 4)
318 call hm%ks_pot%init(gr%der, gr%np, gr%np_part, hm%d%nspin, hm%theory_level, family_is_mgga_with_exc(hm%xc))
319
320 !Initialize Poisson solvers
321 call poisson_init(hm%psolver, namespace, space, gr%der, mc, gr%stencil, st%qtot)
322
323 ! Initialize external potential
324 call epot_init(hm%ep, namespace, gr, hm%ions, hm%psolver, hm%d%ispin, hm%xc%family, hm%kpoints)
325 call kick_init(hm%kick, namespace, space, hm%kpoints, hm%d%ispin)
326
327 hm%zora => zora_t(namespace, hm%der, hm%d, hm%ep, hm%mass)
329 !Temporary construction of the ion-electron interactions
330 call hm%v_ie_loc%init(gr, hm%psolver, hm%ions, namespace)
331 if (hm%ep%nlcc) then
332 call hm%nlcc%init(gr, hm%ions)
333 safe_allocate(st%rho_core(1:gr%np))
334 st%rho_core(:) = m_zero
335 end if
336
337 !Static magnetic field or rashba spin-orbit interaction requires complex wavefunctions
338 if (parse_is_defined(namespace, 'StaticMagneticField') .or. list_has_gauge_field(ext_partners) .or. &
339 parse_is_defined(namespace, 'RashbaSpinOrbitCoupling')) then
340 call states_set_complex(st)
341 end if
342
343 !%Variable CalculateSelfInducedMagneticField
344 !%Type logical
345 !%Default no
346 !%Section Hamiltonian
347 !%Description
348 !% The existence of an electronic current implies the creation of a self-induced magnetic
349 !% field, which may in turn back-react on the system. Of course, a fully consistent treatment
350 !% of this kind of effect should be done in QED theory, but we will attempt a first
351 !% approximation to the problem by considering the lowest-order relativistic terms
352 !% plugged into the normal Hamiltonian equations (spin-other-orbit coupling terms, etc.).
353 !% For the moment being, none of this is done, but a first step is taken by calculating
354 !% the induced magnetic field of a system that has a current, by considering the magnetostatic
355 !% approximation and Biot-Savart law:
356 !%
357 !% <math> \nabla^2 \vec{A} + 4\pi\alpha \vec{J} = 0</math>
358 !%
359 !% <math> \vec{B} = \vec{\nabla} \times \vec{A}</math>
360 !%
361 !% If <tt>CalculateSelfInducedMagneticField</tt> is set to yes, this <i>B</i> field is
362 !% calculated at the end of a <tt>gs</tt> calculation (nothing is done -- yet -- in the <tt>td</tt>case)
363 !% and printed out, if the <tt>Output</tt> variable contains the <tt>potential</tt> keyword (the prefix
364 !% of the output files is <tt>Bind</tt>).
365 !%End
366 call parse_variable(namespace, 'CalculateSelfInducedMagneticField', .false., hm%self_induced_magnetic)
367 if (hm%self_induced_magnetic) then
368 safe_allocate(hm%a_ind(1:gr%np_part, 1:space%dim))
369 safe_allocate(hm%b_ind(1:gr%np_part, 1:space%dim))
370
371 !(for dim = we could save some memory, but it is better to keep it simple)
372 end if
373
374 ! Absorbing boundaries
375 call absorbing_boundaries_init(hm%abs_boundaries, namespace, space, gr)
376
377 hm%inh_term = .false.
378 call oct_exchange_remove(hm%oct_exchange)
379
380 hm%adjoint = .false.
381
382 call hm%phase%init(gr, hm%d%kpt, hm%kpoints, st%d, space)
383
384 !%Variable DFTULevel
385 !%Type integer
386 !%Default no
387 !%Section Hamiltonian::XC
388 !%Description
389 !% This variable selects which DFT+U expression is added to the Hamiltonian.
390 !%Option dft_u_none 0
391 !% No +U term is not applied.
392 !%Option dft_u_empirical 1
393 !% An empiricial Hubbard U is added on the orbitals specified in the block species
394 !% with hubbard_l and hubbard_u
395 !%Option dft_u_acbn0 2
396 !% Octopus determines the effective U term using the
397 !% ACBN0 functional as defined in PRX 5, 011006 (2015)
398 !%End
399 call parse_variable(namespace, 'DFTULevel', dft_u_none, hm%lda_u_level)
400 call messages_print_var_option('DFTULevel', hm%lda_u_level, namespace=namespace)
401 if (hm%lda_u_level /= dft_u_none) then
402 call lda_u_init(hm%lda_u, namespace, space, hm%lda_u_level, gr, ions, st, mc, hm%kpoints)
403
404 !In the present implementation of DFT+U, in case of spinors, we have off-diagonal terms
405 !in spin space which break the assumption of the generalized Bloch theorem
406 if (kick_get_type(hm%kick) == kick_magnon_mode .and. gr%der%boundaries%spiral) then
407 call messages_not_implemented("DFT+U with generalized Bloch theorem and magnon kick", namespace=namespace)
408 end if
409
410 ! We rebuild the phase for the orbital projection, similarly to the one of the pseudopotentials
411 if(hm%lda_u_level /= dft_u_none .and. hm%phase%is_allocated()) then
412 call lda_u_build_phase_correction(hm%lda_u, space, hm%d, gr%der%boundaries, namespace, hm%kpoints)
413 end if
414 end if
415
416 !%Variable HamiltonianApplyPacked
417 !%Type logical
418 !%Default yes
419 !%Section Execution::Optimization
420 !%Description
421 !% If set to yes (the default), Octopus will 'pack' the
422 !% wave-functions when operating with them. This might involve some
423 !% additional copying but makes operations more efficient.
424 !% See also the related <tt>StatesPack</tt> variable.
425 !%End
426 call parse_variable(namespace, 'HamiltonianApplyPacked', .true., hm%is_applied_packed)
427
428 if (hm%theory_level == hartree_fock .and. st%parallel_in_states) then
429 call messages_experimental('Hartree-Fock parallel in states', namespace=namespace)
430 end if
431
432 if (hm%theory_level == generalized_kohn_sham_dft .and. family_is_hybrid(hm%xc) &
433 .and. st%parallel_in_states) then
434 call messages_experimental('Hybrid functionals parallel in states', namespace=namespace)
435 end if
436
437 !%Variable TimeZero
438 !%Type logical
439 !%Default no
440 !%Section Hamiltonian
441 !%Description
442 !% (Experimental) If set to yes, the ground state and other time
443 !% dependent calculation will assume that they are done at time
444 !% zero, so that all time depedent field at that time will be
445 !% included.
446 !%End
447 call parse_variable(namespace, 'TimeZero', .false., hm%time_zero)
448 if (hm%time_zero) call messages_experimental('TimeZero', namespace=namespace)
449
450 !Cam parameters are irrelevant here and are updated later
451 need_exchange_ = optional_default(need_exchange, .false.)
452 if (hm%xc%compute_exchange(hm%theory_level) .or. need_exchange_) then
453 !We test Slater before OEP, as Slater is treated as OEP for the moment....
454 if (hm%xc%functional(func_x,1)%id == xc_oep_x_slater) then
455 call exchange_operator_init(hm%exxop, namespace, space, st, gr%der, mc, gr%stencil, &
456 hm%kpoints, cam_exact_exchange)
457 else if (bitand(hm%xc%family, xc_family_oep) /= 0 .or. hm%theory_level == rdmft) then
458 call exchange_operator_init(hm%exxop, namespace, space, st, gr%der, mc, gr%stencil, &
459 hm%kpoints, hm%xc%cam)
460 if (hm%theory_level == rdmft) hm%exxop%useACE = .false.
461 else
462 call exchange_operator_init(hm%exxop, namespace, space, st, gr%der, mc, gr%stencil, &
463 hm%kpoints, cam_exact_exchange)
464 end if
465 end if
466
467 if (hm%is_applied_packed .and. accel_is_enabled()) then
468 ! Check if we can actually apply the hamiltonian packed
469 if (gr%use_curvilinear) then
470 if (accel_allow_cpu_only()) then
471 hm%is_applied_packed = .false.
472 call messages_write('Cannot use GPUs as curvilinear coordinates are used.')
473 call messages_warning(namespace=namespace)
474 else
475 call messages_write('Cannot use GPUs as curvilinear coordinates are used.', new_line = .true.)
476 call messages_write('Calculation will not be continued. To force execution, set AllowCPUonly = yes.')
477 call messages_fatal(namespace=namespace)
478 end if
479 end if
480 end if
481
482 !We are building the list of external potentials
483 !This is done here at the moment, because we pass directly the mesh
484 !TODO: Once the abstract Hamiltonian knows about an abstract basis, we might move this to the
485 ! abstract Hamiltonian
486 call load_external_potentials(hm%external_potentials, namespace)
487
488 !Some checks which are electron specific, like k-points
490
491 !At the moment we do only have static external potential, so we never update them
493
494 !Build the resulting interactions
495 !TODO: This will be moved to the actual interactions
496 call build_interactions()
497
498 ! Constrained DFT for noncollinear magnetism
499 if (hm%theory_level /= independent_particles) then
500 call magnetic_constrain_init(hm%magnetic_constrain, namespace, gr, st%d, ions%natoms, ions%min_distance())
501 end if
502
503 ! init maxwell-electrons coupling
504 call mxll_coupling_init(hm%mxll, st%d, gr, namespace, hm%mass)
505
506 if (associated(hm%xc_photons)) then
507 if (hm%xc_photons%wants_to_renormalize_mass()) then
508 ! remornalize the electron mass due to light-matter interaction; here we only deal with it in free space
509 call hm%set_mass(namespace, hm%xc_photons%get_renormalized_mass())
510 end if
511 end if
512
513 if (hm%xc%compute_exchange(hm%theory_level) .or. need_exchange_) call hm%exxop%write_info(namespace)
514
515 call profiling_out('HAMILTONIAN_ELEC_INIT')
516 pop_sub(hamiltonian_elec_init)
517
518 contains
519
520 ! ---------------------------------------------------------
521 subroutine build_external_potentials()
522 type(list_iterator_t) :: iter
523 class(*), pointer :: potential
524 integer :: iop
525
527
528 safe_allocate(hm%v_ext_pot(1:gr%np))
529 hm%v_ext_pot(1:gr%np) = m_zero
530
531 call iter%start(hm%external_potentials)
532 do while (iter%has_next())
533 potential => iter%get_next()
534 select type (potential)
535 class is (external_potential_t)
536
537 call potential%allocate_memory(gr)
538 call potential%calculate(namespace, gr, hm%psolver)
539 !To preserve the old behavior, we are adding the various potentials
540 !to the corresponding arrays
541 select case (potential%type)
543 call lalg_axpy(gr%np, m_one, potential%pot, hm%v_ext_pot)
544
546 if (states_are_real(st)) then
547 message(1) = "Cannot use static magnetic field with real wavefunctions"
548 call messages_fatal(1, namespace=namespace)
549 end if
550
551 if (.not. allocated(hm%ep%b_field)) then
552 safe_allocate(hm%ep%b_field(1:3)) !Cannot be space%dim
553 hm%ep%b_field(1:3) = m_zero
554 end if
555 hm%ep%b_field(1:3) = hm%ep%b_field(1:3) + potential%b_field(1:3)
556
557 if (.not. allocated(hm%ep%a_static)) then
558 safe_allocate(hm%ep%a_static(1:gr%np, 1:space%dim))
559 hm%ep%a_static(1:gr%np, 1:space%dim) = m_zero
560 end if
561 call lalg_axpy(gr%np, space%dim, m_one, potential%a_static, hm%ep%a_static)
562
564 if (.not. allocated(hm%ep%e_field)) then
565 safe_allocate(hm%ep%e_field(1:space%dim))
566 hm%ep%e_field(1:space%dim) = m_zero
567 end if
568 hm%ep%e_field(1:space%dim) = hm%ep%e_field(1:space%dim) + potential%e_field(1:space%dim)
569
570 !In the fully periodic case, we use Berry phases
571 if (space%periodic_dim < space%dim) then
572 if (.not. allocated(hm%v_static)) then
573 safe_allocate(hm%v_static(1:gr%np))
574 hm%v_static(1:gr%np) = m_zero
575 end if
576 if (.not. allocated(hm%ep%v_ext)) then
577 safe_allocate(hm%ep%v_ext(1:gr%np_part))
578 hm%ep%v_ext(1:gr%np_part) = m_zero
579 end if
580 call lalg_axpy(gr%np, m_one, potential%pot, hm%v_static)
581 call lalg_axpy(gr%np, m_one, potential%v_ext, hm%ep%v_ext)
582 end if
583
584 if (hm%kpoints%use_symmetries) then
585 do iop = 1, symmetries_number(hm%kpoints%symm)
586 if (iop == symmetries_identity_index(hm%kpoints%symm)) cycle
587 if (.not. symm_op_invariant_cart(hm%kpoints%symm%ops(iop), hm%ep%e_field, 1e-5_real64)) then
588 message(1) = "The StaticElectricField breaks (at least) one of the symmetries used to reduce the k-points."
589 message(2) = "Set SymmetryBreakDir equal to StaticElectricField."
590 call messages_fatal(2, namespace=namespace)
591 end if
592 end do
593 end if
594
595 end select
596 call potential%deallocate_memory()
597
598 class default
599 assert(.false.)
600 end select
601 end do
602
604 end subroutine build_external_potentials
605
606 ! ---------------------------------------------------------
607 subroutine external_potentials_checks()
608 type(list_iterator_t) :: iter
609 class(*), pointer :: potential
610
612
613 call iter%start(hm%external_potentials)
614 do while (iter%has_next())
615 potential => iter%get_next()
616 select type (potential)
617 class is (external_potential_t)
618
619 if (potential%type == external_pot_static_efield .and. hm%kpoints%reduced%npoints > 1) then
620 message(1) = "Applying StaticElectricField in a periodic direction is only accurate for large supercells."
621 message(2) = "Single-point Berry phase is not appropriate when k-point sampling is needed."
622 call messages_warning(2, namespace=namespace)
623 end if
624
625 class default
626 assert(.false.)
627 end select
628 end do
629
631 end subroutine external_potentials_checks
632
633
634 !The code in this routines needs to know about the external potentials.
635 !This will be treated in the future by the interactions directly.
636 subroutine build_interactions()
637 logical :: external_potentials_present
638 logical :: kick_present
639
641
642 if (allocated(hm%ep%e_field) .and. space%is_periodic() .and. .not. list_has_gauge_field(ext_partners)) then
643 ! only need vberry if there is a field in a periodic direction
644 ! and we are not setting a gauge field
645 if (any(abs(hm%ep%e_field(1:space%periodic_dim)) > m_epsilon)) then
646 safe_allocate(hm%vberry(1:gr%np, 1:hm%d%nspin))
647 hm%vberry = m_zero
648 end if
649 end if
650
651 external_potentials_present = epot_have_external_potentials(hm%ep) .or. &
652 list_has_lasers(ext_partners) .or. allocated(hm%v_static)
653
654
655 kick_present = hamiltonian_elec_has_kick(hm)
656
657 call pcm_init(hm%pcm, namespace, space, ions, gr, st%qtot, st%val_charge, external_potentials_present, kick_present)
658 if (hm%pcm%run_pcm) then
659 if (hm%theory_level /= kohn_sham_dft) call messages_not_implemented("PCM for TheoryLevel /= kohn_sham", namespace=namespace)
660 end if
661
663
664 end subroutine build_interactions
665
666
667 end subroutine hamiltonian_elec_init
668
669
670
671
672 ! ---------------------------------------------------------
673 subroutine hamiltonian_elec_end(hm)
674 type(hamiltonian_elec_t), target, intent(inout) :: hm
675
676 type(partner_iterator_t) :: iter
677 class(interaction_partner_t), pointer :: potential
678
679 push_sub(hamiltonian_elec_end)
680
681 if (hm%is_copy_snapshot) then
682 ! release only pointer targets created by hamiltonian_elec_copy
684 pop_sub(hamiltonian_elec_end)
685 return
686 end if
687
688 call hm%hm_base%end()
689 call hm%vnl%end()
690
691 call hm%phase%end()
692
693 call hm%ks_pot%end()
694 safe_deallocate_a(hm%vberry)
695 safe_deallocate_a(hm%a_ind)
696 safe_deallocate_a(hm%b_ind)
697 safe_deallocate_a(hm%v_ext_pot)
698
699 safe_deallocate_p(hm%zora)
700
701 call poisson_end(hm%psolver)
703 nullify(hm%xc)
704
705 call kick_end(hm%kick)
706 call epot_end(hm%ep)
707 if (hm%owns_ions .and. associated(hm%ions)) then
708 ! ions_copy currently aliases species objects through species wrappers.
709 ! Detach wrapper storage to avoid double-finalizing shared species.
710 if (allocated(hm%ions%species)) then
711 safe_deallocate_a(hm%ions%species)
712 hm%ions%nspecies = 0
713 end if
714 safe_deallocate_p(hm%ions)
715 hm%owns_ions = .false.
716 else
717 nullify(hm%ions)
718 hm%owns_ions = .false.
719 end if
720
721 call absorbing_boundaries_end(hm%abs_boundaries)
722
723 call states_elec_dim_end(hm%d)
724
725 if (hm%scissor%apply) call scissor_end(hm%scissor)
726
727 call exchange_operator_end(hm%exxop)
728 call lda_u_end(hm%lda_u)
729
730 safe_deallocate_a(hm%energy)
732 if (hm%pcm%run_pcm) call pcm_end(hm%pcm)
733
734 call hm%v_ie_loc%end()
735 call hm%nlcc%end()
736
737 call iter%start(hm%external_potentials)
738 do while (iter%has_next())
739 potential => iter%get_next()
740 safe_deallocate_p(potential)
741 end do
742 call hm%external_potentials%empty()
743 safe_deallocate_a(hm%v_static)
744
745 call magnetic_constrain_end(hm%magnetic_constrain)
746
747 call mxll_coupling_end(hm%mxll)
748
749 hm%is_copy_snapshot = .false.
750
751 pop_sub(hamiltonian_elec_end)
752 end subroutine hamiltonian_elec_end
753
754 ! ---------------------------------------------------------
759 type(hamiltonian_elec_t), target, intent(inout) :: hm
760
761 type(partner_iterator_t) :: iter
762 class(interaction_partner_t), pointer :: potential
763
765
766 if (hm%owns_ions .and. associated(hm%ions)) then
767 if (allocated(hm%ions%species)) then
768 safe_deallocate_a(hm%ions%species)
769 hm%ions%nspecies = 0
770 end if
771 safe_deallocate_p(hm%ions)
772 hm%owns_ions = .false.
773 else
774 nullify(hm%ions)
775 hm%owns_ions = .false.
776 end if
777
778 if (associated(hm%zora)) then
779 safe_deallocate_p(hm%zora)
780 end if
781
782 call iter%start(hm%external_potentials)
783 do while (iter%has_next())
784 potential => iter%get_next()
785 safe_deallocate_p(potential)
786 end do
787 call hm%external_potentials%empty()
788
789 call hm%hm_base%end()
790 call hm%ks_pot%end()
791 call hm%phase%end()
792 call hm%vnl%end()
793 if (hm%lda_u_level /= dft_u_none) call lda_u_end(hm%lda_u)
794
795 hm%is_copy_snapshot = .false.
796
799
800 ! ---------------------------------------------------------
807 subroutine hamiltonian_elec_copy(hm_out, hm_in)
808 type(hamiltonian_elec_t), target, intent(out) :: hm_out
809 type(hamiltonian_elec_t), intent(inout) :: hm_in
810
811 type(ions_t), pointer :: ions_out
812 type(zora_t), pointer :: zora_out
813 type(partner_list_t) :: external_potentials_in
814 type(partner_list_t) :: external_potentials_out
815 type(list_iterator_t) :: iter
816 class(*), pointer :: ptr
817 class(external_potential_t), pointer :: ext_pot
818 logical :: runtime_initialized
819
820 push_sub(hamiltonian_elec_copy)
821
823 runtime_initialized = associated(hm_in%kpoints)
824
825 call iter%start(hm_in%external_potentials)
826 do while (iter%has_next())
827 ptr => iter%get_next()
828 select type (ptr)
829 class is (external_potential_t)
830 call external_potentials_in%add(ptr)
831 call external_potential_clone(ext_pot, ptr)
832 call external_potentials_out%add(ext_pot)
833 class default
834 call messages_not_implemented("hamiltonian_elec_copy with unsupported external potential partner")
835 end select
836 end do
837 call hm_in%external_potentials%empty()
838
839 hm_out = hm_in
840 hm_out%is_copy_snapshot = .true.
841
842 call epot_bind_poisson_solver(hm_out%ep, hm_out%psolver)
843 call nonlocal_pseudopotential_rebind_projectors(hm_out%vnl, hm_out%ep)
844
845 if (accel_is_enabled() .and. runtime_initialized .and. associated(hm_out%der)) then
846 call hamiltonian_elec_base_accel_rebuild(hm_out%hm_base, hm_out%der%mesh)
847 call ks_potential_accel_rebuild(hm_out%ks_pot)
848 call phase_accel_rebuild(hm_out%phase, hm_out%der%mesh, hm_out%d%kpt)
849 call nonlocal_pseudopotential_accel_rebuild(hm_out%vnl, hm_out%space, hm_out%der%mesh)
850 end if
851
852 call hm_out%external_potentials%empty()
853 call iter%start(external_potentials_out)
854 do while (iter%has_next())
855 ptr => iter%get_next()
856 select type (ptr)
857 class is (external_potential_t)
858 call hm_out%external_potentials%add(ptr)
859 class default
860 call messages_not_implemented("hamiltonian_elec_copy with unsupported external potential partner")
861 end select
862 end do
863 call iter%start(external_potentials_in)
864 do while (iter%has_next())
865 ptr => iter%get_next()
866 select type (ptr)
867 class is (external_potential_t)
868 call hm_in%external_potentials%add(ptr)
869 class default
870 call messages_not_implemented("hamiltonian_elec_copy with unsupported external potential partner")
871 end select
872 end do
873 call external_potentials_out%empty()
874 call external_potentials_in%empty()
875
876 if (runtime_initialized) then
877 call kick_end(hm_out%kick)
878 call kick_copy(hm_out%kick, hm_in%kick)
879
880 if (allocated(hm_in%energy)) then
881 if (.not. allocated(hm_out%energy)) safe_allocate(hm_out%energy)
882 call energy_copy(hm_in%energy, hm_out%energy)
883 else
884 safe_deallocate_a(hm_out%energy)
885 end if
886
887 call magnetic_constrain_copy(hm_out%magnetic_constrain, hm_in%magnetic_constrain)
888 call mxll_coupling_copy(hm_out%mxll, hm_in%mxll, hm_out%der)
889 end if
890
891 if (associated(hm_in%ions)) then
892 safe_allocate(ions_out)
893 ions_out = hm_in%ions
894 hm_out%ions => ions_out
895 hm_out%owns_ions = .true.
896 call hm_out%v_ie_loc%bind(hm_out%psolver, hm_out%ions)
897 call hm_out%nlcc%bind(hm_out%ions)
898 else
899 nullify(hm_out%ions)
900 hm_out%owns_ions = .false.
901 nullify(hm_out%v_ie_loc%atoms_dist)
902 nullify(hm_out%v_ie_loc%atom)
903 nullify(hm_out%v_ie_loc%pos)
904 nullify(hm_out%nlcc%atoms_dist)
905 nullify(hm_out%nlcc%atom)
906 nullify(hm_out%nlcc%pos)
907 end if
908
909 if (hm_out%lda_u_level /= dft_u_none .and. associated(hm_out%ions)) then
910 call lda_u_rebind_after_copy(hm_out%lda_u, hm_out%ions)
911 if (accel_is_enabled() .and. runtime_initialized) then
912 call lda_u_accel_rebuild(hm_out%lda_u, hm_out%d%kpt)
913 end if
914 end if
915
916 if (associated(hm_in%zora)) then
917 allocate(zora_out, source=hm_in%zora)
918 hm_out%zora => zora_out
919 else
920 nullify(hm_out%zora)
921 end if
922
923 pop_sub(hamiltonian_elec_copy)
924
925 end subroutine hamiltonian_elec_copy
926
927 ! ---------------------------------------------------------
930 type(hamiltonian_elec_t), intent(in) :: hm_in
931
932 logical :: runtime_initialized
933
935
936 runtime_initialized = associated(hm_in%kpoints)
937 if (.not. runtime_initialized) then
939 return
940 end if
941
942 if (hm_in%theory_level == hartree_fock .or. hm_in%theory_level == rdmft) then
943 call messages_not_implemented("hamiltonian_elec_copy with HF/RDMFT exchange operator paths")
944 end if
945
946 if (associated(hm_in%xc)) then
947 if (hm_in%xc%compute_exchange(hm_in%theory_level)) then
948 call messages_not_implemented("hamiltonian_elec_copy with hybrid/OEP exchange paths")
949 end if
950 end if
951
952 if (hamiltonian_elec_inh_term(hm_in)) then
953 call messages_not_implemented("hamiltonian_elec_copy with OCT inhomogeneous/source term")
954 end if
955
956 if (oct_exchange_enabled(hm_in%oct_exchange)) then
957 call messages_not_implemented("hamiltonian_elec_copy with OCT exchange term")
958 end if
959
960 if (hm_in%pcm%run_pcm) then
961 call messages_not_implemented("hamiltonian_elec_copy with PCM enabled")
962 end if
963
966
967
968 ! ---------------------------------------------------------
969 ! True if the Hamiltonian is Hermitian, false otherwise
970 logical function hamiltonian_elec_hermitian(hm)
971 class(hamiltonian_elec_t), intent(in) :: hm
972
974 hamiltonian_elec_hermitian = .not.((hm%abs_boundaries%abtype == imaginary_absorbing) .or. &
975 oct_exchange_enabled(hm%oct_exchange))
976
978 end function hamiltonian_elec_hermitian
979
980 ! ---------------------------------------------------------
981 ! True if the Hamiltonian needs to apply the mGGA term
982 pure logical function hamiltonian_elec_needs_mgga_term(hm, terms)
983 class(hamiltonian_elec_t), intent(in) :: hm
984 integer, intent(in) :: terms
985
987
988 if (bitand(term_mgga, terms) /= 0 .and. family_is_mgga_with_exc(hm%xc) &
989 .and. hm%theory_level == generalized_kohn_sham_dft) then
991 end if
992
993 !For OEP
994 if(terms == term_mgga .and. family_is_mgga_with_exc(hm%xc) .and. hm%theory_level == kohn_sham_dft) then
996 end if
997
999
1000
1001
1002 ! ---------------------------------------------------------
1003 subroutine hamiltonian_elec_span(hm, delta, emin, namespace)
1004 class(hamiltonian_elec_t), intent(inout) :: hm
1005 real(real64), intent(in) :: delta(:)
1006 real(real64), intent(in) :: emin
1007 type(namespace_t), intent(in) :: namespace
1008
1009 real(real64) :: emax
1010
1011 push_sub(hamiltonian_elec_span)
1012
1013 ! estimate maximum energy of discrete kinetic operator
1014 ! this neglects possible contributions from the non-local part of the pseudopotentials
1016
1017 hm%spectral_middle_point = (emax + emin) / m_two
1018 hm%spectral_half_span = (emax - emin) / m_two
1019
1020 pop_sub(hamiltonian_elec_span)
1021 end subroutine hamiltonian_elec_span
1022
1023
1024 ! ---------------------------------------------------------
1025 pure logical function hamiltonian_elec_inh_term(hm) result(inh)
1026 type(hamiltonian_elec_t), intent(in) :: hm
1027
1028 inh = hm%inh_term
1029 end function hamiltonian_elec_inh_term
1030
1031
1032 ! ---------------------------------------------------------
1033 subroutine hamiltonian_elec_set_inh(hm, st)
1034 type(hamiltonian_elec_t), intent(inout) :: hm
1035 type(states_elec_t), intent(in) :: st
1036
1037 push_sub(hamiltonian_elec_set_inh)
1038
1039 if (hm%inh_term) call states_elec_end(hm%inh_st)
1040 call states_elec_copy(hm%inh_st, st)
1041 hm%inh_term = .true.
1042
1044 end subroutine hamiltonian_elec_set_inh
1045
1046
1047 ! ---------------------------------------------------------
1048 subroutine hamiltonian_elec_remove_inh(hm)
1049 type(hamiltonian_elec_t), intent(inout) :: hm
1050
1052
1053 if (hm%inh_term) then
1054 call states_elec_end(hm%inh_st)
1055 hm%inh_term = .false.
1056 end if
1057
1059 end subroutine hamiltonian_elec_remove_inh
1060
1061 ! ---------------------------------------------------------
1062 subroutine hamiltonian_elec_adjoint(hm)
1063 type(hamiltonian_elec_t), intent(inout) :: hm
1064
1066
1067 if (.not. hm%adjoint) then
1068 hm%adjoint = .true.
1069 if (hm%abs_boundaries%abtype == imaginary_absorbing) then
1070 hm%abs_boundaries%mf = -hm%abs_boundaries%mf
1071 end if
1072 end if
1073
1075 end subroutine hamiltonian_elec_adjoint
1076
1078 ! ---------------------------------------------------------
1079 subroutine hamiltonian_elec_not_adjoint(hm)
1080 type(hamiltonian_elec_t), intent(inout) :: hm
1081
1083
1084 if (hm%adjoint) then
1085 hm%adjoint = .false.
1086 if (hm%abs_boundaries%abtype == imaginary_absorbing) then
1087 hm%abs_boundaries%mf = -hm%abs_boundaries%mf
1088 end if
1089 end if
1090
1092 end subroutine hamiltonian_elec_not_adjoint
1093
1094
1095 ! ---------------------------------------------------------
1097 subroutine hamiltonian_elec_update(this, mesh, namespace, space, ext_partners, time)
1098 class(hamiltonian_elec_t), intent(inout) :: this
1099 class(mesh_t), intent(in) :: mesh
1100 type(namespace_t), intent(in) :: namespace
1101 class(space_t), intent(in) :: space
1102 type(partner_list_t), intent(in) :: ext_partners
1103 real(real64), optional, intent(in) :: time
1104
1105 integer :: ispin, ip, idir, iatom, ilaser
1106 real(real64) :: aa(space%dim), time_
1107 real(real64), allocatable :: vp(:,:)
1108 type(lasers_t), pointer :: lasers
1109 type(gauge_field_t), pointer :: gfield
1110 real(real64) :: am(space%dim)
1111
1112 push_sub(hamiltonian_elec_update)
1113 call profiling_in("HAMILTONIAN_ELEC_UPDATE")
1114
1115 this%current_time = m_zero
1116 if (present(time)) this%current_time = time
1117
1118 time_ = optional_default(time, 0.0_real64)
1119
1120 ! set everything to zero
1121 call this%hm_base%clear(mesh%np)
1122
1123 ! alllocate the scalar potential for the xc, hartree and external potentials
1124 call this%hm_base%allocate_field(mesh, field_potential, &
1125 complex_potential = this%abs_boundaries%abtype == imaginary_absorbing)
1126
1127 ! the lasers
1128 if (present(time) .or. this%time_zero) then
1129
1130 lasers => list_get_lasers(ext_partners)
1131 if(associated(lasers)) then
1132 do ilaser = 1, lasers%no_lasers
1133 select case (laser_kind(lasers%lasers(ilaser)))
1135 do ispin = 1, this%d%spin_channels
1136 call laser_potential(lasers%lasers(ilaser), mesh, &
1137 this%hm_base%potential(:, ispin), time_)
1138 end do
1139 case (e_field_magnetic)
1140 call this%hm_base%allocate_field(mesh, field_vector_potential + field_uniform_magnetic_field, &
1141 .false.)
1142 ! get the vector potential
1143 safe_allocate(vp(1:mesh%np, 1:space%dim))
1144 vp(1:mesh%np, 1:space%dim) = m_zero
1145 call laser_vector_potential(lasers%lasers(ilaser), mesh, vp, time_)
1146 !$omp parallel do private(idir) schedule(static)
1147 do ip = 1, mesh%np
1148 do idir = 1, space%dim
1149 this%hm_base%vector_potential(idir, ip) = this%hm_base%vector_potential(idir, ip) + vp(ip, idir)/p_c
1150 end do
1151 end do
1152 ! and the magnetic field
1153 call laser_field(lasers%lasers(ilaser), this%hm_base%uniform_magnetic_field(1:space%dim), time_)
1154 safe_deallocate_a(vp)
1156 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
1157 ! get the uniform vector potential associated with a magnetic field
1158 aa = m_zero
1159 call laser_field(lasers%lasers(ilaser), aa, time_)
1160 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) - aa/p_c
1161 end select
1162 end do
1163
1164 if (lasers_with_nondipole_field(lasers)) then
1165 assert( allocated(this%hm_base%uniform_vector_potential))
1166 call lasers_nondipole_laser_field_step(lasers, am, time_)
1167 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) - am/p_c
1168 end if
1169 end if
1170
1171 ! the gauge field
1172 gfield => list_get_gauge_field(ext_partners)
1173 if (associated(gfield)) then
1174 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
1175 call gauge_field_get_vec_pot(gfield, aa)
1176 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) - aa/p_c
1177 end if
1178
1179 ! the electric field for a periodic system through the gauge field
1180 if (allocated(this%ep%e_field) .and. associated(gfield)) then
1181 this%hm_base%uniform_vector_potential(1:space%periodic_dim) = &
1182 this%hm_base%uniform_vector_potential(1:space%periodic_dim) - time_*this%ep%e_field(1:space%periodic_dim)
1183 end if
1184
1185 ! add the photon-free mean-field vector potential
1186 if (associated(this%xc_photons)) then
1187 if(this%xc_photons%lpfmf .and. allocated(this%xc_photons%mf_vector_potential)) then
1188 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
1189 ! here we put a minus sign in front of the mean field term to get the right answer (need to check the formula)
1190 this%hm_base%uniform_vector_potential(1:space%dim) = &
1191 this%hm_base%uniform_vector_potential(1:space%dim) - this%xc_photons%mf_vector_potential(1:space%dim)/p_c
1192 end if
1193 end if
1194
1195 end if
1196
1197 ! the vector potential of a static magnetic field
1198 if (allocated(this%ep%a_static)) then
1199 call this%hm_base%allocate_field(mesh, field_vector_potential, .false.)
1200 !ep%a_static contains 1/c A(r)
1201 !$omp parallel do private(idir) schedule(static)
1202 do ip = 1, mesh%np
1203 do idir = 1, space%dim
1204 this%hm_base%vector_potential(idir, ip) = this%hm_base%vector_potential(idir, ip) + this%ep%a_static(ip, idir)
1205 end do
1206 end do
1207 end if
1208
1209 ! add Maxwell coupling to Hamiltonian and sets the magnetic field for the Zeeman term added below
1210 call mxll_coupling_calc(this%mxll, this%hm_base, mesh, this%d, space)
1211
1212 !The electric field was added to the KS potential
1213 call this%hm_base%accel_copy_pot(mesh)
1214
1215 ! and the static magnetic field
1216 if (allocated(this%ep%b_field)) then
1217 call this%hm_base%allocate_field(mesh, field_uniform_magnetic_field, .false.)
1218 do idir = 1, 3
1219 this%hm_base%uniform_magnetic_field(idir) = this%hm_base%uniform_magnetic_field(idir) + this%ep%b_field(idir)
1220 end do
1221 end if
1222
1223 ! Combine the uniform and non-uniform fields and compute the Zeeman term
1224 call this%hm_base%update_magnetic_terms(mesh, this%ep%gyromagnetic_ratio, this%d%ispin)
1225
1226 ! This needs to be called at the end as the zeeman term enters the potential
1227 call hamiltonian_elec_update_pot(this, mesh, accumulate = .true.)
1228
1229 if (this%mxll%test_equad) then
1230 call set_electric_quadrupole_pot(this%mxll, mesh)
1231 end if
1232
1233 call build_phase()
1234
1235 call profiling_out("HAMILTONIAN_ELEC_UPDATE")
1237
1238 contains
1239
1240 subroutine build_phase()
1241 integer :: ik, imat, nmat, max_npoints, offset
1242 integer :: ip
1243 integer :: iphase, nphase
1244
1246
1247 if ((.not. this%kpoints%gamma_only()) .or. allocated(this%hm_base%uniform_vector_potential)) then
1248
1249 call profiling_in('UPDATE_PHASES')
1250 ! now regenerate the phases for the pseudopotentials
1251 do iatom = 1, this%ep%natoms
1252 call projector_init_phases(this%ep%proj(iatom), space%dim, this%d, this%der%boundaries, this%kpoints, &
1253 vec_pot = this%hm_base%uniform_vector_potential, vec_pot_var = this%hm_base%vector_potential)
1254 end do
1255
1256 call profiling_out('UPDATE_PHASES')
1257 end if
1258
1259 if (allocated(this%hm_base%uniform_vector_potential)) then
1260
1261 call this%phase%update(mesh, this%d%kpt, this%kpoints, this%d, space, this%hm_base%uniform_vector_potential)
1262
1263 ! We rebuild the phase for the orbital projection, similarly to the one of the pseudopotentials
1264 if (this%lda_u_level /= dft_u_none) then
1265 call lda_u_build_phase_correction(this%lda_u, space, this%d, this%der%boundaries, namespace, this%kpoints, &
1266 vec_pot = this%hm_base%uniform_vector_potential, vec_pot_var = this%hm_base%vector_potential)
1267 end if
1268 end if
1269
1270 max_npoints = this%vnl%max_npoints
1271 nmat = this%vnl%nprojector_matrices
1272
1273
1274 if (this%phase%is_allocated() .and. allocated(this%vnl%projector_matrices)) then
1275
1276 nphase = 1
1277 if (this%der%boundaries%spiralBC) nphase = 3
1278
1279 if (.not. allocated(this%vnl%projector_phases)) then
1280 safe_allocate(this%vnl%projector_phases(1:max_npoints, 1:nphase, nmat, this%d%kpt%start:this%d%kpt%end))
1281 if (accel_is_enabled()) then
1282 call accel_create_buffer(this%vnl%buff_projector_phases, accel_mem_read_only, &
1283 type_cmplx, this%vnl%total_points*nphase*this%d%kpt%nlocal)
1284 ! We need to save nphase, with which the array has been build,
1285 ! as the number might change throughout the run
1286 this%vnl%nphase = nphase
1287 end if
1288 end if
1289
1290 offset = 0
1291 do ik = this%d%kpt%start, this%d%kpt%end
1292 do imat = 1, this%vnl%nprojector_matrices
1293 iatom = this%vnl%projector_to_atom(imat)
1294 do iphase = 1, nphase
1295 !$omp parallel do simd schedule(static)
1296 do ip = 1, this%vnl%projector_matrices(imat)%npoints
1297 this%vnl%projector_phases(ip, iphase, imat, ik) = this%ep%proj(iatom)%phase(ip, iphase, ik)
1298 end do
1299
1300 if (accel_is_enabled() .and. this%vnl%projector_matrices(imat)%npoints > 0) then
1301 call accel_write_buffer(this%vnl%buff_projector_phases, &
1302 this%vnl%projector_matrices(imat)%npoints, this%vnl%projector_phases(1:, iphase, imat, ik), &
1303 offset = offset, async=.true.)
1304 end if
1305 offset = offset + this%vnl%projector_matrices(imat)%npoints
1306 end do
1307 end do
1308 end do
1309
1310 end if
1311
1312 call accel_finish()
1313
1315 end subroutine build_phase
1316
1317 end subroutine hamiltonian_elec_update
1318
1319
1320 !----------------------------------------------------------------
1323 ! TODO: See Issue #1064
1324 subroutine hamiltonian_elec_update_pot(this, mesh, accumulate)
1325 type(hamiltonian_elec_t), intent(inout) :: this
1326 class(mesh_t), intent(in) :: mesh
1327 logical, optional, intent(in) :: accumulate
1328
1329 integer :: ispin, ip
1330
1332
1333 ! By default we nullify first the result
1334 if (.not. optional_default(accumulate, .false.)) then
1335 !$omp parallel private(ip, ispin)
1336 do ispin = 1, this%d%nspin
1337 !$omp do simd schedule(static)
1338 do ip = 1, mesh%np
1339 this%hm_base%potential(ip, ispin) = m_zero
1340 end do
1341 end do
1342 !$omp end parallel
1343 if (this%abs_boundaries%abtype == imaginary_absorbing) then
1344 this%hm_base%Impotential = m_zero
1345 end if
1346 end if
1347
1348 !$omp parallel private(ip, ispin)
1349 do ispin = 1, this%d%nspin
1350 if (ispin <= 2) then
1351 !$omp do simd schedule(static)
1352 ! this%vhxc(ip, ispin) is added after the calculation of ZORA potential
1353 do ip = 1, mesh%np
1354 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + this%ep%vpsl(ip) + this%v_ext_pot(ip)
1355 end do
1356
1358 if (this%pcm%run_pcm) then
1359 if (this%pcm%solute) then
1360 !$omp do simd schedule(static)
1361 do ip = 1, mesh%np
1362 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + &
1363 this%pcm%v_e_rs(ip) + this%pcm%v_n_rs(ip)
1364 end do
1365 !$omp end do simd nowait
1366 end if
1367 if (this%pcm%localf) then
1368 !$omp do simd schedule(static)
1369 do ip = 1, mesh%np
1370 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + &
1371 this%pcm%v_ext_rs(ip)
1372 end do
1373 !$omp end do simd nowait
1374 end if
1375 end if
1376
1378 if (this%abs_boundaries%abtype == imaginary_absorbing) then
1379 !$omp do simd schedule(static)
1380 do ip = 1, mesh%np
1381 this%hm_base%Impotential(ip, ispin) = this%hm_base%Impotential(ip, ispin) + this%abs_boundaries%mf(ip)
1382 end do
1383 !$omp end do simd nowait
1384 end if
1385 end if
1386 end do
1387 !$omp end parallel
1388
1389 ! scalar relativistic ZORA contribution
1390 ! \boldsymbol{p} \frac{c^2}{2c^2 - V} \boldsymbol{p} \Phi^\mathrm{ZORA}
1391 if (this%ep%reltype == scalar_relativistic_zora .or. this%ep%reltype == fully_relativistic_zora) then
1392 call this%zora%update(this%der, this%hm_base%potential)
1393 end if
1394
1395 !$omp parallel private(ip, ispin)
1396 do ispin = 1, this%d%nspin
1397 ! Adding Zeeman potential to hm_base%potential
1398 if (allocated(this%hm_base%zeeman_pot)) then
1399 !$omp do simd schedule(static)
1400 do ip = 1, mesh%np
1401 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + this%hm_base%zeeman_pot(ip, ispin)
1402 end do
1403 !$omp end do simd nowait
1404 end if
1405
1406 ! Adding Quadrupole potential from static E-field (test)
1407 if (this%mxll%test_equad) then
1408 !$omp do simd schedule(static)
1409 do ip = 1, mesh%np
1410 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + this%mxll%e_quadrupole_pot(ip)
1411 end do
1412 !$omp end do simd
1413 end if
1414 end do
1415 !$omp end parallel
1416
1417
1418 ! Add the Hartree and KS potential
1419 call this%ks_pot%add_vhxc(this%hm_base%potential)
1420
1421 call this%hm_base%accel_copy_pot(mesh)
1422
1424 end subroutine hamiltonian_elec_update_pot
1425
1426 ! ---------------------------------------------------------
1427 subroutine hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
1428 type(hamiltonian_elec_t), intent(inout) :: this
1429 type(namespace_t), intent(in) :: namespace
1430 class(electron_space_t), intent(in) :: space
1431 type(grid_t), intent(in) :: gr
1432 type(ions_t), target, intent(inout) :: ions
1433 type(partner_list_t), intent(in) :: ext_partners
1434 type(states_elec_t), intent(inout) :: st
1435 real(real64), optional, intent(in) :: time
1436
1438
1439 this%ions => ions
1440 call epot_generate(this%ep, namespace, gr, this%ions, this%d)
1441
1442 ! Interation terms are treated below
1443
1444 ! First we add the static electric field
1445 if (allocated(this%ep%e_field) .and. space%periodic_dim < space%dim) then
1446 call lalg_axpy(gr%np, m_one, this%v_static, this%ep%vpsl)
1447 end if
1448
1449 ! Here we need to pass this again, else test are failing.
1450 ! This is not a real problem, as the multisystem framework will indeed to this anyway
1451 this%v_ie_loc%atoms_dist => ions%atoms_dist
1452 this%v_ie_loc%atom => ions%atom
1453 call this%v_ie_loc%calculate()
1454
1455 ! At the moment we need to add this to ep%vpsl, to keep the behavior of the code
1456 call lalg_axpy(gr%np, m_one, this%v_ie_loc%potential(:,1), this%ep%vpsl)
1457
1458 ! Here we need to reinit the NLCC object
1459 ! This is not a real problem, as the multisystem framework will indeed to this anyway
1460 if (this%ep%nlcc) then
1461 call this%nlcc%end()
1462 call this%nlcc%init(gr, ions)
1463 call this%nlcc%calculate()
1464 call lalg_copy(gr%np, this%nlcc%density(:,1), st%rho_core)
1465 end if
1466
1467 call this%vnl%build(space, gr, this%ep)
1468 call hamiltonian_elec_update(this, gr, namespace, space, ext_partners, time)
1469
1470 ! Check if projectors are still compatible with apply_packed on GPU
1471 if (this%is_applied_packed .and. accel_is_enabled()) then
1472 if (this%ep%non_local .and. .not. this%vnl%apply_projector_matrices) then
1473 if (accel_allow_cpu_only()) then
1474 call messages_write('Relativistic pseudopotentials have not been fully implemented for GPUs.')
1475 call messages_warning(namespace=namespace)
1476 else
1477 call messages_write('Relativistic pseudopotentials have not been fully implemented for GPUs.',&
1478 new_line = .true.)
1479 call messages_write('Calculation will not be continued. To force execution, set AllowCPUonly = yes.')
1480 call messages_fatal(namespace=namespace)
1481 end if
1482 end if
1483
1484 end if
1485
1486 if (this%pcm%run_pcm) then
1489 if (this%pcm%solute) then
1490 call pcm_calc_pot_rs(this%pcm, gr, this%psolver, ions = ions)
1491 end if
1492
1495 ! Interpolation is needed, hence gr%np_part -> 1:gr%np
1496 if (this%pcm%localf .and. allocated(this%v_static)) then
1497 call pcm_calc_pot_rs(this%pcm, gr, this%psolver, v_ext = this%ep%v_ext(1:gr%np_part))
1498 end if
1499
1500 end if
1501
1502 call lda_u_update_basis(this%lda_u, space, gr, ions, st, this%psolver, namespace, this%kpoints, &
1503 this%phase%is_allocated())
1504
1506 end subroutine hamiltonian_elec_epot_generate
1507
1508 ! -----------------------------------------------------------------
1509
1510 real(real64) function hamiltonian_elec_get_time(this) result(time)
1511 type(hamiltonian_elec_t), intent(inout) :: this
1512
1513 time = this%current_time
1514 end function hamiltonian_elec_get_time
1515
1516 ! -----------------------------------------------------------------
1517
1518 pure logical function hamiltonian_elec_apply_packed(this) result(apply)
1519 class(hamiltonian_elec_t), intent(in) :: this
1520
1521 apply = this%is_applied_packed
1524
1525
1526 ! -----------------------------------------------------------------
1527 subroutine zhamiltonian_elec_apply_atom (hm, namespace, space, latt, species, pos, ia, mesh, psi, vpsi)
1528 type(hamiltonian_elec_t), intent(in) :: hm
1529 type(namespace_t), intent(in) :: namespace
1530 class(space_t), intent(in) :: space
1531 type(lattice_vectors_t), intent(in) :: latt
1532 class(species_t), intent(in) :: species
1533 real(real64), intent(in) :: pos(1:space%dim)
1534 integer, intent(in) :: ia
1535 class(mesh_t), intent(in) :: mesh
1536 complex(real64), intent(in) :: psi(:,:)
1537 complex(real64), intent(out) :: vpsi(:,:)
1538
1539 integer :: idim
1540 real(real64), allocatable :: vlocal(:)
1542
1543 safe_allocate(vlocal(1:mesh%np_part))
1544 vlocal = m_zero
1545 call epot_local_potential(hm%ep, namespace, space, latt, mesh, species, pos, ia, vlocal)
1546
1547 do idim = 1, hm%d%dim
1548 vpsi(1:mesh%np, idim) = vlocal(1:mesh%np) * psi(1:mesh%np, idim)
1549 end do
1550
1551 safe_deallocate_a(vlocal)
1553 end subroutine zhamiltonian_elec_apply_atom
1554
1555 ! ---------------------------------------------------------
1560 subroutine hamiltonian_elec_update_with_ext_pot(this, mesh, space, ext_partners, time, mu)
1561 type(hamiltonian_elec_t), intent(inout) :: this
1562 class(space_t), intent(in) :: space
1563 class(mesh_t), intent(in) :: mesh
1564 type(partner_list_t), intent(in) :: ext_partners
1565 real(real64), intent(in) :: time(1:2)
1566 real(real64), intent(in) :: mu(1:2)
1567
1568 integer :: ispin, ip, idir, iatom, ilaser, itime
1569 real(real64) :: aa(space%dim), bb(space%dim), time_
1570 real(real64), allocatable :: vp(:,:)
1571 real(real64), allocatable :: velectric(:)
1572 type(lasers_t), pointer :: lasers
1573 type(gauge_field_t), pointer :: gfield
1574
1576 call profiling_in("HAMILTONIAN_ELEC_UPDATE_EXT_POT")
1577
1578 this%current_time = m_zero
1579 this%current_time = time(1)
1580
1581 ! set everything to zero
1582 call this%hm_base%clear(mesh%np)
1583
1584 ! the xc, hartree and external potentials
1585 call this%hm_base%allocate_field(mesh, field_potential, &
1586 complex_potential = this%abs_boundaries%abtype == imaginary_absorbing)
1587
1588 do itime = 1, 2
1589 time_ = time(itime)
1590
1591 lasers => list_get_lasers(ext_partners)
1592 if(associated(lasers)) then
1593 do ilaser = 1, lasers%no_lasers
1594 select case (laser_kind(lasers%lasers(ilaser)))
1595 case (e_field_scalar_potential, e_field_electric)
1596 safe_allocate(velectric(1:mesh%np))
1597 do ispin = 1, this%d%spin_channels
1598 velectric = m_zero
1599 call laser_potential(lasers%lasers(ilaser), mesh, velectric, time_)
1600 !$omp parallel do simd schedule(static)
1601 do ip = 1, mesh%np
1602 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + mu(itime) * velectric(ip)
1603 end do
1604 end do
1605 safe_deallocate_a(velectric)
1606 case (e_field_magnetic)
1607 call this%hm_base%allocate_field(mesh, field_vector_potential + field_uniform_magnetic_field, .false.)
1608 ! get the vector potential
1609 safe_allocate(vp(1:mesh%np, 1:space%dim))
1610 vp(1:mesh%np, 1:space%dim) = m_zero
1611 call laser_vector_potential(lasers%lasers(ilaser), mesh, vp, time_)
1612 do idir = 1, space%dim
1613 !$omp parallel do schedule(static)
1614 do ip = 1, mesh%np
1615 this%hm_base%vector_potential(idir, ip) = this%hm_base%vector_potential(idir, ip) &
1616 - mu(itime) * vp(ip, idir)/p_c
1617 end do
1618 end do
1619 ! and the magnetic field
1620 bb = m_zero
1621 call laser_field(lasers%lasers(ilaser), bb(1:space%dim), time_)
1622 this%hm_base%uniform_magnetic_field = this%hm_base%uniform_magnetic_field &
1623 - mu(itime) * bb
1624 safe_deallocate_a(vp)
1625 case (e_field_vector_potential)
1626 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
1627 ! get the uniform vector potential associated with a magnetic field
1628 aa = m_zero
1629 call laser_field(lasers%lasers(ilaser), aa, time_)
1630 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) &
1631 - mu(itime) * aa/p_c
1632 end select
1633 end do
1634 end if
1635
1636 ! the gauge field
1637 gfield => list_get_gauge_field(ext_partners)
1638 if (associated(gfield)) then
1639 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
1640 call gauge_field_get_vec_pot(gfield, aa)
1641 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) - aa/p_c
1642 end if
1643
1644 ! the electric field for a periodic system through the gauge field
1645 ! TODO: The condition is wrong here: the e_field should be in non-periodic dims as E field
1646 ! and as a gauge field in the periodic dim, unless we use a Bery phase, in which, we do not use it
1647 ! this way. But this is unrelated to the gauge field
1648 if (allocated(this%ep%e_field) .and. associated(gfield)) then
1649 this%hm_base%uniform_vector_potential(1:space%periodic_dim) = &
1650 this%hm_base%uniform_vector_potential(1:space%periodic_dim) - time_*this%ep%e_field(1:space%periodic_dim)
1651 end if
1652
1653 end do
1654
1655 ! the vector potential of a static magnetic field
1656 if (allocated(this%ep%a_static)) then
1657 call this%hm_base%allocate_field(mesh, field_vector_potential, .false.)
1658 !ep%a_static contains 1/c A(r)
1659 !$omp parallel do schedule(static) private(idir)
1660 do ip = 1, mesh%np
1661 do idir = 1, space%dim
1662 this%hm_base%vector_potential(idir, ip) = this%hm_base%vector_potential(idir, ip) + this%ep%a_static(ip, idir)
1663 end do
1664 end do
1665 end if
1666
1667 !The electric field is added to the KS potential
1668 call this%hm_base%accel_copy_pot(mesh)
1669
1670 ! and the static magnetic field
1671 if (allocated(this%ep%b_field)) then
1672 call this%hm_base%allocate_field(mesh, field_uniform_magnetic_field, .false.)
1673 do idir = 1, 3
1674 this%hm_base%uniform_magnetic_field(idir) = this%hm_base%uniform_magnetic_field(idir) + this%ep%b_field(idir)
1675 end do
1676 end if
1677
1678 call this%hm_base%update_magnetic_terms(mesh, this%ep%gyromagnetic_ratio, this%d%ispin)
1679
1680 call hamiltonian_elec_update_pot(this, mesh)
1681
1682 call build_phase()
1683
1684 call profiling_out("HAMILTONIAN_ELEC_UPDATE_EXT_POT")
1686
1687 contains
1688
1689 subroutine build_phase()
1690 integer :: ik, imat, nmat, max_npoints, offset, iphase, nphase
1691
1693
1694 if ((.not. this%kpoints%gamma_only()) .or. allocated(this%hm_base%uniform_vector_potential)) then
1695
1696 call profiling_in('UPDATE_PHASES')
1697 ! now regenerate the phases for the pseudopotentials
1698 do iatom = 1, this%ep%natoms
1699 call projector_init_phases(this%ep%proj(iatom), space%dim, this%d, this%der%boundaries, this%kpoints, &
1700 vec_pot = this%hm_base%uniform_vector_potential, vec_pot_var = this%hm_base%vector_potential)
1701 end do
1702
1703 call profiling_out('UPDATE_PHASES')
1704 end if
1705
1706 if (allocated(this%hm_base%uniform_vector_potential)) then
1707 call this%phase%update(mesh, this%d%kpt, this%kpoints, this%d, space, this%hm_base%uniform_vector_potential)
1708 end if
1709
1710 max_npoints = this%vnl%max_npoints
1711 nmat = this%vnl%nprojector_matrices
1712
1713
1714 if (this%phase%is_allocated() .and. allocated(this%vnl%projector_matrices)) then
1715
1716 nphase = 1
1717 if (this%der%boundaries%spiralBC) nphase = 3
1718
1719 if (.not. allocated(this%vnl%projector_phases)) then
1720 safe_allocate(this%vnl%projector_phases(1:max_npoints, nphase, nmat, this%d%kpt%start:this%d%kpt%end))
1721 if (accel_is_enabled()) then
1722 call accel_create_buffer(this%vnl%buff_projector_phases, accel_mem_read_only, &
1723 type_cmplx, this%vnl%total_points*nphase*this%d%kpt%nlocal)
1724 end if
1725 end if
1726
1727 offset = 0
1728 do ik = this%d%kpt%start, this%d%kpt%end
1729 do imat = 1, this%vnl%nprojector_matrices
1730 iatom = this%vnl%projector_to_atom(imat)
1731 do iphase = 1, nphase
1732 !$omp parallel do schedule(static)
1733 do ip = 1, this%vnl%projector_matrices(imat)%npoints
1734 this%vnl%projector_phases(ip, iphase, imat, ik) = this%ep%proj(iatom)%phase(ip, iphase, ik)
1735 end do
1736
1737 if (accel_is_enabled() .and. this%vnl%projector_matrices(imat)%npoints > 0) then
1738 call accel_write_buffer(this%vnl%buff_projector_phases, &
1739 this%vnl%projector_matrices(imat)%npoints, this%vnl%projector_phases(1:, iphase, imat, ik), &
1740 offset = offset)
1741 end if
1742 offset = offset + this%vnl%projector_matrices(imat)%npoints
1743 end do
1744 end do
1745 end do
1746
1747 end if
1748
1750 end subroutine build_phase
1751
1753
1754 logical function hamiltonian_elec_needs_current(hm, states_are_real)
1755 type(hamiltonian_elec_t), intent(in) :: hm
1756 logical, intent(in) :: states_are_real
1757
1759
1760 if (hm%self_induced_magnetic) then
1761 if (.not. states_are_real) then
1763 else
1764 message(1) = 'No current density for real states since it is identically zero.'
1765 call messages_warning(1)
1766 end if
1767 end if
1768
1770
1771 ! ---------------------------------------------------------
1772 subroutine zhamiltonian_elec_apply_all(hm, namespace, gr, st, hst)
1773 type(hamiltonian_elec_t), intent(inout) :: hm
1774 type(namespace_t), intent(in) :: namespace
1775 type(grid_t), intent(in) :: gr
1776 type(states_elec_t), intent(inout) :: st
1777 type(states_elec_t), intent(inout) :: hst
1778
1779 integer :: ik, ib, ist
1780 complex(real64), allocatable :: psi(:, :)
1781 complex(real64), allocatable :: psiall(:, :, :, :)
1782
1784
1785 do ik = st%d%kpt%start, st%d%kpt%end
1786 do ib = st%group%block_start, st%group%block_end
1787 call zhamiltonian_elec_apply_batch(hm, namespace, gr, st%group%psib(ib, ik), hst%group%psib(ib, ik))
1788 end do
1789 end do
1790
1791 if (oct_exchange_enabled(hm%oct_exchange)) then
1792
1793 safe_allocate(psiall(gr%np_part, 1:hst%d%dim, st%st_start:st%st_end, st%d%kpt%start:st%d%kpt%end))
1794
1795 call states_elec_get_state(st, gr, psiall)
1796
1797 call oct_exchange_prepare(hm%oct_exchange, gr, psiall, hm%xc, hm%psolver, namespace)
1798
1799 safe_deallocate_a(psiall)
1800
1801 safe_allocate(psi(gr%np_part, 1:hst%d%dim))
1802
1803 do ik = 1, st%nik
1804 do ist = 1, st%nst
1805 call states_elec_get_state(hst, gr, ist, ik, psi)
1806 call oct_exchange_operator(hm%oct_exchange, namespace, gr, psi, ist, ik)
1807 call states_elec_set_state(hst, gr, ist, ik, psi)
1808 end do
1809 end do
1810
1811 safe_deallocate_a(psi)
1812
1813 end if
1814
1816 end subroutine zhamiltonian_elec_apply_all
1817
1818 ! ---------------------------------------------------------
1819 logical function hamiltonian_elec_has_kick(hm)
1820 type(hamiltonian_elec_t), intent(in) :: hm
1821
1823
1824 hamiltonian_elec_has_kick = (abs(hm%kick%delta_strength) > m_epsilon)
1825
1827 end function hamiltonian_elec_has_kick
1828
1830 !
1831 subroutine hamiltonian_elec_set_mass(this, namespace, mass)
1832 class(hamiltonian_elec_t) , intent(inout) :: this
1833 type(namespace_t), intent(in) :: namespace
1834 real(real64), intent(in) :: mass
1835
1837
1838 if (parse_is_defined(namespace, 'ParticleMass')) then
1839 message(1) = 'Attempting to redefine a non-unit electron mass'
1840 call messages_fatal(1)
1841 else
1842 this%mass = mass
1843 end if
1844
1846 end subroutine hamiltonian_elec_set_mass
1847
1848 ! ---------------------------------------------------------
1849 subroutine hamiltonian_elec_diagonal (hm, mesh, diag, ik)
1850 type(hamiltonian_elec_t), intent(in) :: hm
1851 class(mesh_t), intent(in) :: mesh
1852 real(real64), contiguous, intent(out) :: diag(:,:)
1853 integer, intent(in) :: ik
1854
1855 integer :: idim, ip, ispin
1856
1857 real(real64), allocatable :: ldiag(:)
1858
1860
1861 safe_allocate(ldiag(1:mesh%np))
1862
1863 diag = m_zero
1864
1865 call derivatives_lapl_diag(hm%der, ldiag)
1866
1867 do idim = 1, hm%d%dim
1868 diag(1:mesh%np, idim) = -m_half/hm%mass*ldiag(1:mesh%np)
1869 end do
1870
1871 select case (hm%d%ispin)
1872
1873 case (unpolarized, spin_polarized)
1874 ispin = hm%d%get_spin_index(ik)
1875 diag(1:mesh%np, 1) = diag(1:mesh%np, 1) + hm%ep%vpsl(1:mesh%np)
1876
1877 case (spinors)
1878 do ip = 1, mesh%np
1879 diag(ip, 1) = diag(ip, 1) + hm%ep%vpsl(ip)
1880 diag(ip, 2) = diag(ip, 2) + hm%ep%vpsl(ip)
1881 end do
1882
1883 end select
1884
1885 call hm%ks_pot%add_vhxc(diag)
1886
1888 end subroutine hamiltonian_elec_diagonal
1889
1890
1891
1892#include "undef.F90"
1893#include "real.F90"
1894#include "hamiltonian_elec_inc.F90"
1895
1896#include "undef.F90"
1897#include "complex.F90"
1898#include "hamiltonian_elec_inc.F90"
1899
1900end module hamiltonian_elec_oct_m
1901
1902!! Local Variables:
1903!! mode: f90
1904!! coding: utf-8
1905!! End:
subroutine build_external_potentials()
subroutine build_phase()
subroutine external_potentials_checks()
subroutine build_interactions()
constant times a vector plus a vector
Definition: lalg_basic.F90:173
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:188
integer, parameter, public imaginary_absorbing
subroutine, public absorbing_boundaries_end(this)
subroutine, public absorbing_boundaries_init(this, namespace, space, gr)
pure logical function, public accel_allow_cpu_only()
Definition: accel.F90:412
subroutine, public accel_finish()
Definition: accel.F90:1098
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
integer, parameter, public accel_mem_read_only
Definition: accel.F90:185
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:124
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64) function, public derivatives_lapl_get_max_eigenvalue(this)
Get maximum eigenvalue of discrete Laplacian. For the star and star_general stencils,...
subroutine, public energy_copy(ein, eout)
Definition: energy.F90:170
subroutine, public epot_bind_poisson_solver(ep, psolver)
Bind the Poisson solver if the potential manages a density. The Poisson solver pointer is aliased whe...
Definition: epot.F90:423
logical function, public epot_have_external_potentials(ep)
Definition: epot.F90:624
integer, parameter, public scalar_relativistic_zora
Definition: epot.F90:169
subroutine, public epot_end(ep)
Definition: epot.F90:384
integer, parameter, public fully_relativistic_zora
Definition: epot.F90:169
subroutine, public epot_init(ep, namespace, gr, ions, psolver, ispin, xc_family, kpoints)
Definition: epot.F90:221
subroutine, public epot_generate(ep, namespace, mesh, ions, st_d)
Definition: epot.F90:439
subroutine, public exchange_operator_init(this, namespace, space, st, der, mc, stencil, kpoints, cam)
subroutine, public exchange_operator_end(this)
logical function, public list_has_gauge_field(partners)
type(gauge_field_t) function, pointer, public list_get_gauge_field(partners)
logical function, public list_has_lasers(partners)
type(lasers_t) function, pointer, public list_get_lasers(partners)
integer, parameter, public external_pot_from_file
potential, defined in a file
subroutine, public load_external_potentials(external_potentials, namespace)
integer, parameter, public external_pot_charge_density
user-defined function for charge density
integer, parameter, public external_pot_usdef
user-defined function for local potential
integer, parameter, public external_pot_static_efield
Static electric field.
integer, parameter, public external_pot_static_bfield
Static magnetic field.
subroutine, public external_potential_clone(pot_out, pot_in)
Deep-clone an external potential instance.
subroutine, public gauge_field_get_vec_pot(this, vec_pot)
real(real64), parameter, public m_two
Definition: global.F90:202
real(real64), parameter, public m_zero
Definition: global.F90:200
integer, parameter, public rdmft
Definition: global.F90:250
integer, parameter, public hartree_fock
Definition: global.F90:250
integer, parameter, public independent_particles
Theory level.
Definition: global.F90:250
integer, parameter, public generalized_kohn_sham_dft
Definition: global.F90:250
integer, parameter, public kohn_sham_dft
Definition: global.F90:250
real(real64), parameter, public m_epsilon
Definition: global.F90:216
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:242
real(real64), parameter, public m_one
Definition: global.F90:201
This module implements the underlying real-space grid.
Definition: grid.F90:119
This module defines an abstract class for Hamiltonians.
integer, parameter, public field_uniform_magnetic_field
integer, parameter, public field_uniform_vector_potential
integer, parameter, public field_vector_potential
integer, parameter, public term_mgga
subroutine, public hamiltonian_elec_base_accel_rebuild(this, mesh)
Rebuild accelerator buffers after an intrinsic copy.
integer, parameter, public field_potential
pure logical function, public hamiltonian_elec_apply_packed(this)
subroutine, public hamiltonian_elec_set_inh(hm, st)
subroutine, public zvmask(mesh, hm, st)
subroutine, public zhamiltonian_elec_apply_batch(hm, namespace, mesh, psib, hpsib, terms, set_bc)
subroutine, public hamiltonian_elec_adjoint(hm)
subroutine, public hamiltonian_elec_end(hm)
subroutine, public zhamiltonian_elec_apply_single(hm, namespace, mesh, psi, hpsi, ist, ik, terms, set_bc, set_phase)
pure logical function hamiltonian_elec_needs_mgga_term(hm, terms)
logical function, public hamiltonian_elec_has_kick(hm)
logical function hamiltonian_elec_hermitian(hm)
subroutine, public dhamiltonian_elec_apply_single(hm, namespace, mesh, psi, hpsi, ist, ik, terms, set_bc, set_phase)
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
real(real64) function, public hamiltonian_elec_get_time(this)
subroutine, public hamiltonian_elec_copy(hm_out, hm_in)
Deep-copy a hamiltonian_elec_t snapshot.
subroutine, public dvmask(mesh, hm, st)
logical function, public hamiltonian_elec_needs_current(hm, states_are_real)
subroutine, public hamiltonian_elec_remove_inh(hm)
subroutine, public zhamiltonian_elec_apply_atom(hm, namespace, space, latt, species, pos, ia, mesh, psi, vpsi)
subroutine hamiltonian_elec_release_copy_owned(hm)
Release owned pointer targets created by hamiltonian_elec_copy.
subroutine, public zhamiltonian_elec_apply_all(hm, namespace, gr, st, hst)
subroutine, public hamiltonian_elec_diagonal(hm, mesh, diag, ik)
subroutine, public hamiltonian_elec_update_pot(this, mesh, accumulate)
Update the KS potential of the electronic Hamiltonian.
integer, parameter, public velocity
subroutine hamiltonian_elec_update(this, mesh, namespace, space, ext_partners, time)
(re-)build the Hamiltonian for the next application:
subroutine hamiltonian_elec_copy_guard_supported(hm_in)
Hard-fail guard matrix for unsupported hamiltonian snapshots.
subroutine hamiltonian_elec_span(hm, delta, emin, namespace)
subroutine, public hamiltonian_elec_init(hm, namespace, space, gr, ions, ext_partners, st, theory_level, xc, mc, kpoints, need_exchange, xc_photons)
subroutine dhamiltonian_elec_apply(hm, namespace, mesh, psib, hpsib, terms, set_bc)
pure logical function, public hamiltonian_elec_inh_term(hm)
subroutine hamiltonian_elec_set_mass(this, namespace, mass)
set the effective electron mass, checking whether it was previously redefined.
subroutine, public dhamiltonian_elec_apply_batch(hm, namespace, mesh, psib, hpsib, terms, set_bc)
subroutine, public hamiltonian_elec_update_with_ext_pot(this, mesh, space, ext_partners, time, mu)
This is an extension of "hamiltonian_elec_update_pot" to be used by the CFM4 propagator....
subroutine, public hamiltonian_elec_not_adjoint(hm)
subroutine, public zhamiltonian_elec_external(this, mesh, psib, vpsib)
subroutine zhamiltonian_elec_apply(hm, namespace, mesh, psib, hpsib, terms, set_bc)
This module defines classes and functions for interaction partners.
Definition: io.F90:116
subroutine, public kick_copy(kick_out, kick_in)
Definition: kick.F90:753
integer, parameter, public kick_magnon_mode
Definition: kick.F90:165
subroutine, public kick_end(kick)
Definition: kick.F90:796
subroutine, public kick_init(kick, namespace, space, kpoints, nspin)
Definition: kick.F90:225
pure integer function, public kick_get_type(kick)
Definition: kick.F90:1365
A module to handle KS potential, without the external potential.
subroutine, public ks_potential_accel_rebuild(this)
Rebuild accelerator buffers after an intrinsic copy.
subroutine, public laser_vector_potential(laser, mesh, aa, time)
Definition: lasers.F90:1084
subroutine, public lasers_nondipole_laser_field_step(this, field, time)
Retrieves the NDSFA vector_potential correction. The nondipole field is obtained for consecutive time...
Definition: lasers.F90:1156
logical function, public lasers_with_nondipole_field(lasers)
Check if a nondipole SFA correction should be computed for the given laser.
Definition: lasers.F90:743
integer, parameter, public e_field_electric
Definition: lasers.F90:179
integer, parameter, public e_field_vector_potential
Definition: lasers.F90:179
subroutine, public laser_potential(laser, mesh, pot, time)
Definition: lasers.F90:1049
integer, parameter, public e_field_scalar_potential
Definition: lasers.F90:179
integer pure elemental function, public laser_kind(laser)
Definition: lasers.F90:719
subroutine, public laser_field(laser, field, time)
Retrieves the value of either the electric or the magnetic field. If the laser is given by a scalar p...
Definition: lasers.F90:1121
integer, parameter, public e_field_magnetic
Definition: lasers.F90:179
integer, parameter, public dft_u_none
Definition: lda_u.F90:205
subroutine, public lda_u_init(this, namespace, space, level, gr, ions, st, mc, kpoints)
Definition: lda_u.F90:286
subroutine, public lda_u_update_basis(this, space, gr, ions, st, psolver, namespace, kpoints, has_phase)
Definition: lda_u.F90:791
subroutine, public lda_u_accel_rebuild(this, kpt)
Rebuild DFT+U accelerator buffers after intrinsic assignment.
Definition: lda_u.F90:748
subroutine, public lda_u_rebind_after_copy(this, ions)
Rebind non-owning pointers after intrinsic assignment of lda_u_t.
Definition: lda_u.F90:703
subroutine, public lda_u_build_phase_correction(this, space, std, boundaries, namespace, kpoints, vec_pot, vec_pot_var)
Build the phase correction to the global phase for all orbitals.
Definition: lda_u.F90:918
subroutine, public lda_u_end(this)
Definition: lda_u.F90:657
This module implements fully polymorphic linked lists, and some specializations thereof.
This modules implements the routines for doing constrain DFT for noncollinear magnetism.
subroutine, public magnetic_constrain_copy(this_out, this_in)
Deep-copy magnetic constrain data.
subroutine, public magnetic_constrain_end(this)
Releases memory of the magnetic constrain.
subroutine, public magnetic_constrain_init(this, namespace, mesh, std, natoms, min_dist)
Initilializes the magnetic_constrain_t object.
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1068
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1040
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public mxll_coupling_init(this, d, gr, namespace, mass)
Parse variables and initialize Maxwell coupling.
subroutine, public set_electric_quadrupole_pot(this, mesh)
Computes the electric quadrupole potential where .
subroutine, public mxll_coupling_end(this)
Finalize and deallocate Maxwell coupling arrays.
subroutine, public mxll_coupling_calc(this, hm_base, mesh, d, space)
Add the Maxwell coupling to the electronic Hamiltonian.
subroutine, public mxll_coupling_copy(this_out, this_in, der_target)
Deep-copy Maxwell-electron coupling data and rebind derivatives.
subroutine, public nonlocal_pseudopotential_accel_rebuild(this, space, mesh)
Rebuild accelerator buffers after an intrinsic copy.
subroutine, public nonlocal_pseudopotential_rebind_projectors(this, epot)
Rebind projector matrix pointers (map, position) to a target epot.
logical function, public oct_exchange_enabled(this)
subroutine, public oct_exchange_remove(this)
this module contains the low-level part of the output system
Definition: output_low.F90:117
Some general things and nomenclature:
Definition: par_vec.F90:173
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:455
subroutine, public pcm_calc_pot_rs(pcm, mesh, psolver, ions, v_h, v_ext, kick, time_present, kick_time)
Definition: pcm.F90:1217
subroutine, public pcm_end(pcm)
Definition: pcm.F90:3072
real(real64), dimension(:,:), allocatable delta
D_E matrix in JCP 139, 024105 (2013).
Definition: pcm.F90:270
subroutine, public pcm_init(pcm, namespace, space, ions, grid, qtot, val_charge, external_potentials_present, kick_present)
Initializes the PCM calculation: reads the VdW molecular cavity and generates the PCM response matrix...
Definition: pcm.F90:296
subroutine, public phase_accel_rebuild(phase, mesh, kpt)
Rebuild phase accelerator buffers after an intrinsic copy.
Definition: phase.F90:492
subroutine, public poisson_init(this, namespace, space, der, mc, stencil, qtot, label, solver, verbose, force_serial, force_cmplx)
Definition: poisson.F90:233
subroutine, public poisson_end(this)
Definition: poisson.F90:689
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public projector_init_phases(this, dim, std, bnd, kpoints, vec_pot, vec_pot_var)
Definition: projector.F90:266
subroutine, public scissor_end(this)
Definition: scissor.F90:240
subroutine, public states_set_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_dim_copy(dout, din)
subroutine, public states_elec_dim_end(dim)
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
This module provides routines for communicating states when using states parallelization.
integer pure function, public symmetries_identity_index(this)
Definition: symmetries.F90:599
integer pure function, public symmetries_number(this)
Definition: symmetries.F90:553
type(type_t), public type_cmplx
Definition: types.F90:136
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
type(xc_cam_t), parameter, public cam_exact_exchange
Use only Hartree Fock exact exchange.
Definition: xc_cam.F90:155
integer, parameter, public xc_oep_x_slater
Slater approximation to the exact exchange.
integer, parameter, public func_x
Definition: xc.F90:116
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:599
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
Definition: xc.F90:614
This module implements the "photon-free" electron-photon exchange-correlation functional.
Definition: xc_photons.F90:123
This module implements the ZORA terms for the Hamoiltonian.
Definition: zora.F90:118
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
The abstract Hamiltonian class defines a skeleton for specific implementations.
abstract class for general interaction partners
Describes mesh distribution to nodes.
Definition: mesh.F90:187
Stores all communicators and groups.
Definition: multicomm.F90:208
The states_elec_t class contains all electronic wave functions.
This class described the 'photon-exchange' electron-photon xc functionals, based on QEDFT.
Definition: xc_photons.F90:160
This class is responsible for calculating and applying the ZORA.
Definition: zora.F90:147
int true(void)