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_f03_lib_m
97 use zora_oct_m
98
99 implicit none
100
101 private
102 public :: &
112 magnus, &
113 dvmask, &
114 zvmask, &
129
130
131 type, extends(hamiltonian_abst_t) :: hamiltonian_elec_t
132 ! Components are public by default
133
136 type(space_t), private :: space
137 type(states_elec_dim_t) :: d
138 type(hamiltonian_elec_base_t) :: hm_base
139 type(phase_t) :: phase
140 type(energy_t), allocatable :: energy
141 type(absorbing_boundaries_t) :: abs_boundaries
142 type(ks_potential_t) :: ks_pot
143 real(real64), allocatable :: vberry(:,:)
144
145 type(derivatives_t), pointer, private :: der
146
147 type(nonlocal_pseudopotential_t) :: vnl
148
149 type(ions_t), pointer :: ions
150 real(real64) :: exx_coef
151
152 type(poisson_t) :: psolver
153
155 logical :: self_induced_magnetic
156 real(real64), allocatable :: a_ind(:, :)
157 real(real64), allocatable :: b_ind(:, :)
158
159 integer :: theory_level
160 type(xc_t), pointer :: xc
161 type(xc_photons_t), pointer :: xc_photons
162
163 type(epot_t) :: ep
164 type(pcm_t) :: pcm
165
167 logical, private :: adjoint
168
170 real(real64), private :: mass
171
173 logical, private :: inh_term
174 type(states_elec_t) :: inh_st
175
178 type(oct_exchange_t) :: oct_exchange
179
180 type(scissor_t) :: scissor
181
182 real(real64) :: current_time
183 logical, private :: is_applied_packed
184
186 type(lda_u_t) :: lda_u
187 integer :: lda_u_level
188
189 logical, public :: time_zero
190
191 type(exchange_operator_t), public :: exxop
192
193 type(kpoints_t), pointer, public :: kpoints => null()
194
195 type(partner_list_t) :: external_potentials
196 real(real64), allocatable, public :: v_ext_pot(:)
197 real(real64), allocatable, public :: v_static(:)
198
199 type(ion_electron_local_potential_t) :: v_ie_loc
200 type(nlcc_t) :: nlcc
201
202 type(magnetic_constrain_t) :: magnetic_constrain
203
205 type(kick_t) :: kick
206
208 type(mxll_coupling_t) :: mxll
209 type(zora_t), pointer :: zora => null()
210
211 contains
212 procedure :: update => hamiltonian_elec_update
213 procedure :: apply_packed => hamiltonian_elec_apply_packed
214 procedure :: update_span => hamiltonian_elec_span
215 procedure :: dapply => dhamiltonian_elec_apply
216 procedure :: zapply => zhamiltonian_elec_apply
217 procedure :: dmagnus_apply => dhamiltonian_elec_magnus_apply
218 procedure :: zmagnus_apply => zhamiltonian_elec_magnus_apply
219 procedure :: is_hermitian => hamiltonian_elec_hermitian
220 procedure :: set_mass => hamiltonian_elec_set_mass
221 end type hamiltonian_elec_t
222
223 integer, public, parameter :: &
224 LENGTH = 1, &
225 velocity = 2
226
227
228contains
230 ! ---------------------------------------------------------
231 subroutine hamiltonian_elec_init(hm, namespace, space, gr, ions, ext_partners, st, theory_level, xc, &
232 mc, kpoints, need_exchange, xc_photons)
233 type(hamiltonian_elec_t), target, intent(inout) :: hm
234 type(namespace_t), intent(in) :: namespace
235 class(space_t), intent(in) :: space
236 type(grid_t), target, intent(inout) :: gr
237 type(ions_t), target, intent(inout) :: ions
238 type(partner_list_t), intent(inout) :: ext_partners
239 type(states_elec_t), target, intent(inout) :: st
240 integer, intent(in) :: theory_level
241 type(xc_t), target, intent(in) :: xc
242 type(multicomm_t), intent(in) :: mc
243 type(kpoints_t), target, intent(in) :: kpoints
244 logical, optional, intent(in) :: need_exchange
245 type(xc_photons_t), optional, target, intent(in) :: xc_photons
246
247
248 logical :: need_exchange_
249 real(real64) :: rashba_coupling
251 push_sub(hamiltonian_elec_init)
252 call profiling_in('HAMILTONIAN_ELEC_INIT')
254 ! make a couple of local copies
255 hm%space = space
256 hm%theory_level = theory_level
257 call states_elec_dim_copy(hm%d, st%d)
258
259 hm%kpoints => kpoints
261 !%Variable ParticleMass
262 !%Type float
263 !%Default 1.0
264 !%Section Hamiltonian
265 !%Description
266 !% It is possible to make calculations for a particle with a mass
267 !% different from one (atomic unit of mass, or mass of the electron).
268 !% This is useful to describe non-electronic systems, or for
269 !% esoteric purposes.
270 !%End
271 call parse_variable(namespace, 'ParticleMass', m_one, hm%mass)
272
273 !%Variable RashbaSpinOrbitCoupling
274 !%Type float
275 !%Default 0.0
276 !%Section Hamiltonian
277 !%Description
278 !% (Experimental.) For systems described in 2D (electrons confined to 2D in semiconductor structures), one
279 !% may add the Bychkov-Rashba spin-orbit coupling term [Bychkov and Rashba, <i>J. Phys. C: Solid
280 !% State Phys.</i> <b>17</b>, 6031 (1984)]. This variable determines the strength
281 !% of this perturbation, and has dimensions of energy times length.
282 !%End
283 call parse_variable(namespace, 'RashbaSpinOrbitCoupling', m_zero, rashba_coupling, units_inp%energy*units_inp%length)
284 if (parse_is_defined(namespace, 'RashbaSpinOrbitCoupling')) then
285 if (space%dim /= 2) then
286 write(message(1),'(a)') 'Rashba spin-orbit coupling can only be used for two-dimensional systems.'
287 call messages_fatal(1, namespace=namespace)
288 end if
289 call messages_experimental('RashbaSpinOrbitCoupling', namespace=namespace)
290 end if
291
292 call hm%hm_base%init(hm%d%nspin, hm%mass, rashba_coupling)
293 call hm%vnl%init()
294
295 assert(associated(gr%der%lapl))
296 hm%hm_base%kinetic => gr%der%lapl
297
298 safe_allocate(hm%energy)
299
300 !Keep pointers to derivatives, geometry and xc
301 hm%der => gr%der
302 hm%ions => ions
303 hm%xc => xc
304
305 if(present(xc_photons)) then
306 hm%xc_photons => xc_photons
307 else
308 hm%xc_photons => null()
309 end if
311 ! allocate potentials and density of the cores
312 ! In the case of spinors, vxc_11 = hm%vxc(:, 1), vxc_22 = hm%vxc(:, 2), Re(vxc_12) = hm%vxc(:. 3);
313 ! Im(vxc_12) = hm%vxc(:, 4)
314 call hm%ks_pot%init(gr%np, gr%np_part, hm%d%nspin, hm%theory_level, family_is_mgga_with_exc(hm%xc))
315
316 !Initialize Poisson solvers
317 call poisson_init(hm%psolver, namespace, space, gr%der, mc, gr%stencil, st%qtot)
318
319 ! Initialize external potential
320 call epot_init(hm%ep, namespace, gr, hm%ions, hm%psolver, hm%d%ispin, hm%xc%family, hm%kpoints)
321 call kick_init(hm%kick, namespace, space, hm%kpoints, hm%d%ispin)
322
323 hm%zora => zora_t(namespace, hm%der, hm%d, hm%ep, hm%mass)
325 !Temporary construction of the ion-electron interactions
326 call hm%v_ie_loc%init(gr, hm%psolver, hm%ions, namespace)
327 if (hm%ep%nlcc) then
328 call hm%nlcc%init(gr, hm%ions)
329 safe_allocate(st%rho_core(1:gr%np))
330 st%rho_core(:) = m_zero
331 end if
332
333 !Static magnetic field or rashba spin-orbit interaction requires complex wavefunctions
334 if (parse_is_defined(namespace, 'StaticMagneticField') .or. list_has_gauge_field(ext_partners) .or. &
335 parse_is_defined(namespace, 'RashbaSpinOrbitCoupling')) then
336 call states_set_complex(st)
337 end if
338
339 !%Variable CalculateSelfInducedMagneticField
340 !%Type logical
341 !%Default no
342 !%Section Hamiltonian
343 !%Description
344 !% The existence of an electronic current implies the creation of a self-induced magnetic
345 !% field, which may in turn back-react on the system. Of course, a fully consistent treatment
346 !% of this kind of effect should be done in QED theory, but we will attempt a first
347 !% approximation to the problem by considering the lowest-order relativistic terms
348 !% plugged into the normal Hamiltonian equations (spin-other-orbit coupling terms, etc.).
349 !% For the moment being, none of this is done, but a first step is taken by calculating
350 !% the induced magnetic field of a system that has a current, by considering the magnetostatic
351 !% approximation and Biot-Savart law:
352 !%
353 !% <math> \nabla^2 \vec{A} + 4\pi\alpha \vec{J} = 0</math>
354 !%
355 !% <math> \vec{B} = \vec{\nabla} \times \vec{A}</math>
356 !%
357 !% If <tt>CalculateSelfInducedMagneticField</tt> is set to yes, this <i>B</i> field is
358 !% calculated at the end of a <tt>gs</tt> calculation (nothing is done -- yet -- in the <tt>td</tt>case)
359 !% and printed out, if the <tt>Output</tt> variable contains the <tt>potential</tt> keyword (the prefix
360 !% of the output files is <tt>Bind</tt>).
361 !%End
362 call parse_variable(namespace, 'CalculateSelfInducedMagneticField', .false., hm%self_induced_magnetic)
363 if (hm%self_induced_magnetic) then
364 safe_allocate(hm%a_ind(1:gr%np_part, 1:space%dim))
365 safe_allocate(hm%b_ind(1:gr%np_part, 1:space%dim))
366
367 !(for dim = we could save some memory, but it is better to keep it simple)
368 end if
369
370 ! Absorbing boundaries
371 call absorbing_boundaries_init(hm%abs_boundaries, namespace, space, gr)
372
373 hm%inh_term = .false.
374 call oct_exchange_remove(hm%oct_exchange)
375
376 hm%adjoint = .false.
377
378 call hm%phase%init(gr, hm%d%kpt, hm%kpoints, st%d, space)
379
380 !%Variable DFTULevel
381 !%Type integer
382 !%Default no
383 !%Section Hamiltonian::XC
384 !%Description
385 !% This variable selects which DFT+U expression is added to the Hamiltonian.
386 !%Option dft_u_none 0
387 !% No +U term is not applied.
388 !%Option dft_u_empirical 1
389 !% An empiricial Hubbard U is added on the orbitals specified in the block species
390 !% with hubbard_l and hubbard_u
391 !%Option dft_u_acbn0 2
392 !% Octopus determines the effective U term using the
393 !% ACBN0 functional as defined in PRX 5, 011006 (2015)
394 !%End
395 call parse_variable(namespace, 'DFTULevel', dft_u_none, hm%lda_u_level)
396 call messages_print_var_option('DFTULevel', hm%lda_u_level, namespace=namespace)
397 if (hm%lda_u_level /= dft_u_none) then
398 call lda_u_init(hm%lda_u, namespace, space, hm%lda_u_level, gr, ions, st, mc, &
399 hm%kpoints, hm%phase%is_allocated())
400
401 !In the present implementation of DFT+U, in case of spinors, we have off-diagonal terms
402 !in spin space which break the assumption of the generalized Bloch theorem
403 if (kick_get_type(hm%kick) == kick_magnon_mode .and. gr%der%boundaries%spiral) then
404 call messages_not_implemented("DFT+U with generalized Bloch theorem and magnon kick", namespace=namespace)
405 end if
406
407 ! We rebuild the phase for the orbital projection, similarly to the one of the pseudopotentials
408 if(hm%lda_u_level /= dft_u_none .and. hm%phase%is_allocated()) then
409 call lda_u_build_phase_correction(hm%lda_u, space, hm%d, gr%der%boundaries, namespace, hm%kpoints)
410 end if
411 end if
412
413 !%Variable HamiltonianApplyPacked
414 !%Type logical
415 !%Default yes
416 !%Section Execution::Optimization
417 !%Description
418 !% If set to yes (the default), Octopus will 'pack' the
419 !% wave-functions when operating with them. This might involve some
420 !% additional copying but makes operations more efficient.
421 !% See also the related <tt>StatesPack</tt> variable.
422 !%End
423 call parse_variable(namespace, 'HamiltonianApplyPacked', .true., hm%is_applied_packed)
424
425 if (hm%theory_level == hartree_fock .and. st%parallel_in_states) then
426 call messages_experimental('Hartree-Fock parallel in states', namespace=namespace)
427 end if
428
429 if (hm%theory_level == generalized_kohn_sham_dft .and. family_is_hybrid(hm%xc) &
430 .and. st%parallel_in_states) then
431 call messages_experimental('Hybrid functionals parallel in states', namespace=namespace)
432 end if
433
434
435 !%Variable TimeZero
436 !%Type logical
437 !%Default no
438 !%Section Hamiltonian
439 !%Description
440 !% (Experimental) If set to yes, the ground state and other time
441 !% dependent calculation will assume that they are done at time
442 !% zero, so that all time depedent field at that time will be
443 !% included.
444 !%End
445 call parse_variable(namespace, 'TimeZero', .false., hm%time_zero)
446 if (hm%time_zero) call messages_experimental('TimeZero', namespace=namespace)
447
448 !Cam parameters are irrelevant here and are updated later
449 need_exchange_ = optional_default(need_exchange, .false.)
450 if (hm%theory_level == hartree_fock .or. hm%theory_level == hartree &
451 .or. hm%theory_level == rdmft .or. need_exchange_ .or. &
452 (hm%theory_level == generalized_kohn_sham_dft .and. family_is_hybrid(hm%xc)) &
453 .or. hm%xc%functional(func_x,1)%id == xc_oep_x_slater &
454 .or. bitand(hm%xc%family, xc_family_oep) /= 0) then
455 !We test Slater before OEP, as Slater is treated as OEP for the moment....
456 if (hm%xc%functional(func_x,1)%id == xc_oep_x_slater) then
457 call exchange_operator_init(hm%exxop, namespace, space, st, gr%der, mc, gr%stencil, &
458 hm%kpoints, m_zero, m_one, m_zero)
459 else if (bitand(hm%xc%family, xc_family_oep) /= 0 .or. hm%theory_level == rdmft) then
460 call exchange_operator_init(hm%exxop, namespace, space, st, gr%der, mc, gr%stencil, &
461 hm%kpoints, hm%xc%cam_omega, hm%xc%cam_alpha, hm%xc%cam_beta)
462 if (hm%theory_level == rdmft) hm%exxop%useACE = .false.
463 else
464 call exchange_operator_init(hm%exxop, namespace, space, st, gr%der, mc, gr%stencil, &
465 hm%kpoints, m_zero, m_one, m_zero)
466 end if
467 end if
468
469 if (hm%is_applied_packed .and. accel_is_enabled()) then
470 ! Check if we can actually apply the hamiltonian packed
471 if (gr%use_curvilinear) then
472 if (accel_allow_cpu_only()) then
473 hm%is_applied_packed = .false.
474 call messages_write('Cannot use CUDA or OpenCL as curvilinear coordinates are used.')
475 call messages_warning(namespace=namespace)
476 else
477 call messages_write('Cannot use CUDA or OpenCL as curvilinear coordinates are used.', new_line = .true.)
478 call messages_write('Calculation will not be continued. To force execution, set AllowCPUonly = yes.')
479 call messages_fatal(namespace=namespace)
480 end if
481 end if
482 end if
483
484 !We are building the list of external potentials
485 !This is done here at the moment, because we pass directly the mesh
486 !TODO: Once the abstract Hamiltonian knows about an abstract basis, we might move this to the
487 ! abstract Hamiltonian
488 call load_external_potentials(hm%external_potentials, namespace)
489
490 !Some checks which are electron specific, like k-points
492
493 !At the moment we do only have static external potential, so we never update them
495
496 !Build the resulting interactions
497 !TODO: This will be moved to the actual interactions
498 call build_interactions()
499
500 ! Constrained DFT for noncollinear magnetism
501 if (hm%theory_level /= independent_particles) then
502 call magnetic_constrain_init(hm%magnetic_constrain, namespace, gr, st%d, ions%natoms, ions%min_distance())
503 end if
504
505 ! init maxwell-electrons coupling
506 call mxll_coupling_init(hm%mxll, st%d, gr, namespace, hm%mass)
507
508 if (associated(hm%xc_photons)) then
509 if (hm%xc_photons%wants_to_renormalize_mass()) then
510 ! remornalize the electron mass due to light-matter interaction; here we only deal with it in free space
511 call hm%set_mass(namespace, hm%xc_photons%get_renormalized_mass())
512 end if
513 end if
514
515
516 call profiling_out('HAMILTONIAN_ELEC_INIT')
517 pop_sub(hamiltonian_elec_init)
518
519 contains
520
521 ! ---------------------------------------------------------
522 subroutine build_external_potentials()
523 type(list_iterator_t) :: iter
524 class(*), pointer :: potential
525 integer :: iop
526
528
529 safe_allocate(hm%v_ext_pot(1:gr%np))
530 hm%v_ext_pot(1:gr%np) = m_zero
531
532 call iter%start(hm%external_potentials)
533 do while (iter%has_next())
534 potential => iter%get_next()
535 select type (potential)
536 class is (external_potential_t)
537
538 call potential%allocate_memory(gr)
539 call potential%calculate(namespace, gr, hm%psolver)
540 !To preserve the old behavior, we are adding the various potentials
541 !to the corresponding arrays
542 select case (potential%type)
544 call lalg_axpy(gr%np, m_one, potential%pot, hm%v_ext_pot)
545
547 if (states_are_real(st)) then
548 message(1) = "Cannot use static magnetic field with real wavefunctions"
549 call messages_fatal(1, namespace=namespace)
550 end if
551
552 if (.not. allocated(hm%ep%b_field)) then
553 safe_allocate(hm%ep%b_field(1:3)) !Cannot be space%dim
554 hm%ep%b_field(1:3) = m_zero
555 end if
556 hm%ep%b_field(1:3) = hm%ep%b_field(1:3) + potential%b_field(1:3)
557
558 if (.not. allocated(hm%ep%a_static)) then
559 safe_allocate(hm%ep%a_static(1:gr%np, 1:space%dim))
560 hm%ep%a_static(1:gr%np, 1:space%dim) = m_zero
561 end if
562 call lalg_axpy(gr%np, space%dim, m_one, potential%a_static, hm%ep%a_static)
563
565 if (.not. allocated(hm%ep%e_field)) then
566 safe_allocate(hm%ep%e_field(1:space%dim))
567 hm%ep%e_field(1:space%dim) = m_zero
568 end if
569 hm%ep%e_field(1:space%dim) = hm%ep%e_field(1:space%dim) + potential%e_field(1:space%dim)
570
571 !In the fully periodic case, we use Berry phases
572 if (space%periodic_dim < space%dim) then
573 if (.not. allocated(hm%v_static)) then
574 safe_allocate(hm%v_static(1:gr%np))
575 hm%v_static(1:gr%np) = m_zero
576 end if
577 if (.not. allocated(hm%ep%v_ext)) then
578 safe_allocate(hm%ep%v_ext(1:gr%np_part))
579 hm%ep%v_ext(1:gr%np_part) = m_zero
580 end if
581 call lalg_axpy(gr%np, m_one, potential%pot, hm%v_static)
582 call lalg_axpy(gr%np, m_one, potential%v_ext, hm%ep%v_ext)
583 end if
584
585 if (hm%kpoints%use_symmetries) then
586 do iop = 1, symmetries_number(hm%kpoints%symm)
587 if (iop == symmetries_identity_index(hm%kpoints%symm)) cycle
588 if (.not. symm_op_invariant_cart(hm%kpoints%symm%ops(iop), hm%ep%e_field, 1e-5_real64)) then
589 message(1) = "The StaticElectricField breaks (at least) one of the symmetries used to reduce the k-points."
590 message(2) = "Set SymmetryBreakDir equal to StaticElectricField."
591 call messages_fatal(2, namespace=namespace)
592 end if
593 end do
594 end if
595
596 end select
597 call potential%deallocate_memory()
598
599 class default
600 assert(.false.)
601 end select
602 end do
603
605 end subroutine build_external_potentials
606
607 ! ---------------------------------------------------------
608 subroutine external_potentials_checks()
609 type(list_iterator_t) :: iter
610 class(*), pointer :: potential
611
613
614 call iter%start(hm%external_potentials)
615 do while (iter%has_next())
616 potential => iter%get_next()
617 select type (potential)
618 class is (external_potential_t)
619
620 if (potential%type == external_pot_static_efield .and. hm%kpoints%reduced%npoints > 1) then
621 message(1) = "Applying StaticElectricField in a periodic direction is only accurate for large supercells."
622 message(2) = "Single-point Berry phase is not appropriate when k-point sampling is needed."
623 call messages_warning(2, namespace=namespace)
624 end if
625
626 class default
627 assert(.false.)
628 end select
629 end do
630
632 end subroutine external_potentials_checks
633
634
635 !The code in this routines needs to know about the external potentials.
636 !This will be treated in the future by the interactions directly.
637 subroutine build_interactions()
638 logical :: external_potentials_present
639 logical :: kick_present
640
642
643 if (allocated(hm%ep%e_field) .and. space%is_periodic() .and. .not. list_has_gauge_field(ext_partners)) then
644 ! only need vberry if there is a field in a periodic direction
645 ! and we are not setting a gauge field
646 if (any(abs(hm%ep%e_field(1:space%periodic_dim)) > m_epsilon)) then
647 safe_allocate(hm%vberry(1:gr%np, 1:hm%d%nspin))
648 hm%vberry = m_zero
649 end if
650 end if
651
652 external_potentials_present = epot_have_external_potentials(hm%ep) .or. &
653 list_has_lasers(ext_partners) .or. allocated(hm%v_static)
654
655
656 kick_present = hamiltonian_elec_has_kick(hm)
657
658 call pcm_init(hm%pcm, namespace, space, ions, gr, st%qtot, st%val_charge, external_potentials_present, kick_present)
659 if (hm%pcm%run_pcm) then
660 if (hm%theory_level /= kohn_sham_dft) call messages_not_implemented("PCM for TheoryLevel /= kohn_sham", namespace=namespace)
661 end if
662
664
665 end subroutine build_interactions
666
667
668 end subroutine hamiltonian_elec_init
669
670
671 ! ---------------------------------------------------------
672 subroutine hamiltonian_elec_end(hm)
673 type(hamiltonian_elec_t), target, intent(inout) :: hm
674
675 type(partner_iterator_t) :: iter
676 class(interaction_partner_t), pointer :: potential
677
678 push_sub(hamiltonian_elec_end)
679
680 call hm%hm_base%end()
681 call hm%vnl%end()
682
683 call hm%phase%end()
684
685 call hm%ks_pot%end()
686 safe_deallocate_a(hm%vberry)
687 safe_deallocate_a(hm%a_ind)
688 safe_deallocate_a(hm%b_ind)
689 safe_deallocate_a(hm%v_ext_pot)
690
691 safe_deallocate_p(hm%zora)
692
693 call poisson_end(hm%psolver)
694
695 nullify(hm%xc)
696
697 call kick_end(hm%kick)
698 call epot_end(hm%ep)
699 nullify(hm%ions)
700
701 call absorbing_boundaries_end(hm%abs_boundaries)
702
703 call states_elec_dim_end(hm%d)
704
705 if (hm%scissor%apply) call scissor_end(hm%scissor)
706
707 call exchange_operator_end(hm%exxop)
708 call lda_u_end(hm%lda_u)
709
710 safe_deallocate_a(hm%energy)
711
712 if (hm%pcm%run_pcm) call pcm_end(hm%pcm)
713
714 call hm%v_ie_loc%end()
715 call hm%nlcc%end()
716
717 call iter%start(hm%external_potentials)
718 do while (iter%has_next())
719 potential => iter%get_next()
720 safe_deallocate_p(potential)
721 end do
722 call hm%external_potentials%empty()
723 safe_deallocate_a(hm%v_static)
724
725 call magnetic_constrain_end(hm%magnetic_constrain)
726
727 call mxll_coupling_end(hm%mxll)
728
729 pop_sub(hamiltonian_elec_end)
730 end subroutine hamiltonian_elec_end
731
732
733 ! ---------------------------------------------------------
734 ! True if the Hamiltonian is Hermitian, false otherwise
735 logical function hamiltonian_elec_hermitian(hm)
736 class(hamiltonian_elec_t), intent(in) :: hm
737
739 hamiltonian_elec_hermitian = .not.((hm%abs_boundaries%abtype == imaginary_absorbing) .or. &
740 oct_exchange_enabled(hm%oct_exchange))
741
743 end function hamiltonian_elec_hermitian
744
745
746 ! ---------------------------------------------------------
747 subroutine hamiltonian_elec_span(hm, delta, emin, namespace)
748 class(hamiltonian_elec_t), intent(inout) :: hm
749 real(real64), intent(in) :: delta(:)
750 real(real64), intent(in) :: emin
751 type(namespace_t), intent(in) :: namespace
752
753 real(real64) :: emax
754
755 push_sub(hamiltonian_elec_span)
756
757 ! estimate maximum energy of discrete kinetic operator
758 ! this neglects possible contributions from the non-local part of the pseudopotentials
760
761 hm%spectral_middle_point = (emax + emin) / m_two
762 hm%spectral_half_span = (emax - emin) / m_two
763
764 pop_sub(hamiltonian_elec_span)
766
767
768 ! ---------------------------------------------------------
769 pure logical function hamiltonian_elec_inh_term(hm) result(inh)
770 type(hamiltonian_elec_t), intent(in) :: hm
771
772 inh = hm%inh_term
773 end function hamiltonian_elec_inh_term
774
775
776 ! ---------------------------------------------------------
777 subroutine hamiltonian_elec_set_inh(hm, st)
778 type(hamiltonian_elec_t), intent(inout) :: hm
779 type(states_elec_t), intent(in) :: st
780
782
783 if (hm%inh_term) call states_elec_end(hm%inh_st)
784 call states_elec_copy(hm%inh_st, st)
785 hm%inh_term = .true.
786
788 end subroutine hamiltonian_elec_set_inh
789
790
791 ! ---------------------------------------------------------
792 subroutine hamiltonian_elec_remove_inh(hm)
793 type(hamiltonian_elec_t), intent(inout) :: hm
794
796
797 if (hm%inh_term) then
798 call states_elec_end(hm%inh_st)
799 hm%inh_term = .false.
800 end if
801
803 end subroutine hamiltonian_elec_remove_inh
804
805 ! ---------------------------------------------------------
806 subroutine hamiltonian_elec_adjoint(hm)
807 type(hamiltonian_elec_t), intent(inout) :: hm
808
810
811 if (.not. hm%adjoint) then
812 hm%adjoint = .true.
813 if (hm%abs_boundaries%abtype == imaginary_absorbing) then
814 hm%abs_boundaries%mf = -hm%abs_boundaries%mf
815 end if
816 end if
817
819 end subroutine hamiltonian_elec_adjoint
820
821
822 ! ---------------------------------------------------------
823 subroutine hamiltonian_elec_not_adjoint(hm)
824 type(hamiltonian_elec_t), intent(inout) :: hm
825
827
828 if (hm%adjoint) then
829 hm%adjoint = .false.
830 if (hm%abs_boundaries%abtype == imaginary_absorbing) then
831 hm%abs_boundaries%mf = -hm%abs_boundaries%mf
832 end if
833 end if
834
836 end subroutine hamiltonian_elec_not_adjoint
837
838
839 ! ---------------------------------------------------------
841 subroutine hamiltonian_elec_update(this, mesh, namespace, space, ext_partners, time)
842 class(hamiltonian_elec_t), intent(inout) :: this
843 class(mesh_t), intent(in) :: mesh
844 type(namespace_t), intent(in) :: namespace
845 class(space_t), intent(in) :: space
846 type(partner_list_t), intent(in) :: ext_partners
847 real(real64), optional, intent(in) :: time
848
849 integer :: ispin, ip, idir, iatom, ilaser
850 real(real64) :: aa(space%dim), time_
851 real(real64), allocatable :: vp(:,:)
852 type(lasers_t), pointer :: lasers
853 type(gauge_field_t), pointer :: gfield
854 real(real64) :: am(space%dim)
855
857 call profiling_in("HAMILTONIAN_ELEC_UPDATE")
858
859 this%current_time = m_zero
860 if (present(time)) this%current_time = time
861
862 time_ = optional_default(time, 0.0_real64)
863
864 ! set everything to zero
865 call this%hm_base%clear(mesh%np)
866
867 ! alllocate the scalar potential for the xc, hartree and external potentials
868 call this%hm_base%allocate_field(mesh, field_potential, &
869 complex_potential = this%abs_boundaries%abtype == imaginary_absorbing)
871 ! the lasers
872 if (present(time) .or. this%time_zero) then
873
874 lasers => list_get_lasers(ext_partners)
875 if(associated(lasers)) then
876 do ilaser = 1, lasers%no_lasers
877 select case (laser_kind(lasers%lasers(ilaser)))
879 do ispin = 1, this%d%spin_channels
880 call laser_potential(lasers%lasers(ilaser), mesh, &
881 this%hm_base%potential(:, ispin), time_)
882 end do
883 case (e_field_magnetic)
884 call this%hm_base%allocate_field(mesh, field_vector_potential + field_uniform_magnetic_field, &
885 .false.)
886 ! get the vector potential
887 safe_allocate(vp(1:mesh%np, 1:space%dim))
888 vp(1:mesh%np, 1:space%dim) = m_zero
889 call laser_vector_potential(lasers%lasers(ilaser), mesh, vp, time_)
890 !$omp parallel do private(idir) schedule(static)
891 do ip = 1, mesh%np
892 do idir = 1, space%dim
893 this%hm_base%vector_potential(idir, ip) = this%hm_base%vector_potential(idir, ip) + vp(ip, idir)/p_c
894 end do
895 end do
896 ! and the magnetic field
897 call laser_field(lasers%lasers(ilaser), this%hm_base%uniform_magnetic_field(1:space%dim), time_)
898 safe_deallocate_a(vp)
900 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
901 ! get the uniform vector potential associated with a magnetic field
902 aa = m_zero
903 call laser_field(lasers%lasers(ilaser), aa, time_)
904 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) - aa/p_c
905 end select
906 end do
907
908 if (lasers_with_nondipole_field(lasers)) then
909 assert( allocated(this%hm_base%uniform_vector_potential))
910 call lasers_nondipole_laser_field_step(lasers, am, time_)
911 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) - am/p_c
912 end if
913 end if
914
915 ! the gauge field
916 gfield => list_get_gauge_field(ext_partners)
917 if (associated(gfield)) then
918 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
919 call gauge_field_get_vec_pot(gfield, aa)
920 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) - aa/p_c
921 end if
922
923 ! the electric field for a periodic system through the gauge field
924 if (allocated(this%ep%e_field) .and. associated(gfield)) then
925 this%hm_base%uniform_vector_potential(1:space%periodic_dim) = &
926 this%hm_base%uniform_vector_potential(1:space%periodic_dim) - time_*this%ep%e_field(1:space%periodic_dim)
927 end if
928
929 ! add the photon-free mean-field vector potential
930 if (associated(this%xc_photons)) then
931 if(this%xc_photons%lpfmf .and. allocated(this%xc_photons%mf_vector_potential)) then
932 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
933 ! here we put a minus sign in front of the mean field term to get the right answer (need to check the formula)
934 this%hm_base%uniform_vector_potential(1:space%dim) = &
935 this%hm_base%uniform_vector_potential(1:space%dim) - this%xc_photons%mf_vector_potential(1:space%dim)/p_c
936 end if
937 end if
938
939 end if
940
941 ! the vector potential of a static magnetic field
942 if (allocated(this%ep%a_static)) then
943 call this%hm_base%allocate_field(mesh, field_vector_potential, .false.)
944 !ep%a_static contains 1/c A(r)
945 !$omp parallel do private(idir) schedule(static)
946 do ip = 1, mesh%np
947 do idir = 1, space%dim
948 this%hm_base%vector_potential(idir, ip) = this%hm_base%vector_potential(idir, ip) + this%ep%a_static(ip, idir)
949 end do
950 end do
951 end if
952
953 ! add Maxwell coupling to Hamiltonian and sets the magnetic field for the Zeeman term added below
954 call mxll_coupling_calc(this%mxll, this%hm_base, mesh, this%d, space)
955
956 !The electric field was added to the KS potential
957 call this%hm_base%accel_copy_pot(mesh)
958
959 ! and the static magnetic field
960 if (allocated(this%ep%b_field)) then
961 call this%hm_base%allocate_field(mesh, field_uniform_magnetic_field, .false.)
962 do idir = 1, 3
963 this%hm_base%uniform_magnetic_field(idir) = this%hm_base%uniform_magnetic_field(idir) + this%ep%b_field(idir)
964 end do
965 end if
966
967 ! Combine the uniform and non-uniform fields and compute the Zeeman term
968 call this%hm_base%update_magnetic_terms(mesh, this%ep%gyromagnetic_ratio, this%d%ispin)
969
970 ! This needs to be called at the end as the zeeman term enters the potential
971 call hamiltonian_elec_update_pot(this, mesh, accumulate = .true.)
972
973 if (this%mxll%test_equad) then
974 call set_electric_quadrupole_pot(this%mxll, mesh)
975 end if
976
977 call build_phase()
978
979 call profiling_out("HAMILTONIAN_ELEC_UPDATE")
981
982 contains
983
984 subroutine build_phase()
985 integer :: ik, imat, nmat, max_npoints, offset
986 integer :: ip
987 integer :: iphase, nphase
988
990
991 if ((.not. this%kpoints%gamma_only()) .or. allocated(this%hm_base%uniform_vector_potential)) then
992
993 call profiling_in('UPDATE_PHASES')
994 ! now regenerate the phases for the pseudopotentials
995 do iatom = 1, this%ep%natoms
996 call projector_init_phases(this%ep%proj(iatom), space%dim, this%d, this%der%boundaries, this%kpoints, &
997 vec_pot = this%hm_base%uniform_vector_potential, vec_pot_var = this%hm_base%vector_potential)
998 end do
999
1000 call profiling_out('UPDATE_PHASES')
1001 end if
1002
1003 if (allocated(this%hm_base%uniform_vector_potential)) then
1004
1005 call this%phase%update(mesh, this%d%kpt, this%kpoints, this%d, space, this%hm_base%uniform_vector_potential)
1006
1007 ! We rebuild the phase for the orbital projection, similarly to the one of the pseudopotentials
1008 if (this%lda_u_level /= dft_u_none) then
1009 call lda_u_build_phase_correction(this%lda_u, space, this%d, this%der%boundaries, namespace, this%kpoints, &
1010 vec_pot = this%hm_base%uniform_vector_potential, vec_pot_var = this%hm_base%vector_potential)
1011 end if
1012 end if
1013
1014 max_npoints = this%vnl%max_npoints
1015 nmat = this%vnl%nprojector_matrices
1016
1017
1018 if (this%phase%is_allocated() .and. allocated(this%vnl%projector_matrices)) then
1019
1020 nphase = 1
1021 if (this%der%boundaries%spiralBC) nphase = 3
1022
1023 if (.not. allocated(this%vnl%projector_phases)) then
1024 safe_allocate(this%vnl%projector_phases(1:max_npoints, 1:nphase, nmat, this%d%kpt%start:this%d%kpt%end))
1025 if (accel_is_enabled()) then
1026 call accel_create_buffer(this%vnl%buff_projector_phases, accel_mem_read_only, &
1027 type_cmplx, this%vnl%total_points*nphase*this%d%kpt%nlocal)
1028 ! We need to save nphase, with which the array has been build,
1029 ! as the number might change throughout the run
1030 this%vnl%nphase = nphase
1031 end if
1032 end if
1033
1034 offset = 0
1035 do ik = this%d%kpt%start, this%d%kpt%end
1036 do imat = 1, this%vnl%nprojector_matrices
1037 iatom = this%vnl%projector_to_atom(imat)
1038 do iphase = 1, nphase
1039 !$omp parallel do simd schedule(static)
1040 do ip = 1, this%vnl%projector_matrices(imat)%npoints
1041 this%vnl%projector_phases(ip, iphase, imat, ik) = this%ep%proj(iatom)%phase(ip, iphase, ik)
1042 end do
1043
1044 if (accel_is_enabled() .and. this%vnl%projector_matrices(imat)%npoints > 0) then
1045 call accel_write_buffer(this%vnl%buff_projector_phases, &
1046 this%vnl%projector_matrices(imat)%npoints, this%vnl%projector_phases(1:, iphase, imat, ik), &
1047 offset = offset, async=.true.)
1048 end if
1049 offset = offset + this%vnl%projector_matrices(imat)%npoints
1050 end do
1051 end do
1052 end do
1053
1054 end if
1055
1056 call accel_finish()
1057
1059 end subroutine build_phase
1060
1061 end subroutine hamiltonian_elec_update
1062
1063
1064 !----------------------------------------------------------------
1067 ! TODO: See Issue #1064
1068 subroutine hamiltonian_elec_update_pot(this, mesh, accumulate)
1069 type(hamiltonian_elec_t), intent(inout) :: this
1070 class(mesh_t), intent(in) :: mesh
1071 logical, optional, intent(in) :: accumulate
1072
1073 integer :: ispin, ip
1074
1076
1077 ! By default we nullify first the result
1078 if (.not. optional_default(accumulate, .false.)) then
1079 !$omp parallel private(ip, ispin)
1080 do ispin = 1, this%d%nspin
1081 !$omp do simd schedule(static)
1082 do ip = 1, mesh%np
1083 this%hm_base%potential(ip, ispin) = m_zero
1084 end do
1085 end do
1086 !$omp end parallel
1087 if (this%abs_boundaries%abtype == imaginary_absorbing) then
1088 this%hm_base%Impotential = m_zero
1089 end if
1090 end if
1091
1092 !$omp parallel private(ip, ispin)
1093 do ispin = 1, this%d%nspin
1094 if (ispin <= 2) then
1095 !$omp do simd schedule(static)
1096 ! this%vhxc(ip, ispin) is added after the calculation of ZORA potential
1097 do ip = 1, mesh%np
1098 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + this%ep%vpsl(ip) + this%v_ext_pot(ip)
1099 end do
1100
1102 if (this%pcm%run_pcm) then
1103 if (this%pcm%solute) then
1104 !$omp do simd schedule(static)
1105 do ip = 1, mesh%np
1106 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + &
1107 this%pcm%v_e_rs(ip) + this%pcm%v_n_rs(ip)
1108 end do
1109 !$omp end do simd nowait
1110 end if
1111 if (this%pcm%localf) then
1112 !$omp do simd schedule(static)
1113 do ip = 1, mesh%np
1114 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + &
1115 this%pcm%v_ext_rs(ip)
1116 end do
1117 !$omp end do simd nowait
1118 end if
1119 end if
1120
1122 if (this%abs_boundaries%abtype == imaginary_absorbing) then
1123 !$omp do simd schedule(static)
1124 do ip = 1, mesh%np
1125 this%hm_base%Impotential(ip, ispin) = this%hm_base%Impotential(ip, ispin) + this%abs_boundaries%mf(ip)
1126 end do
1127 !$omp end do simd nowait
1128 end if
1129 end if
1130 end do
1131 !$omp end parallel
1132
1133 ! scalar relativistic ZORA contribution
1134 ! \boldsymbol{p} \frac{c^2}{2c^2 - V} \boldsymbol{p} \Phi^\mathrm{ZORA}
1135 if (this%ep%reltype == scalar_relativistic_zora .or. this%ep%reltype == fully_relativistic_zora) then
1136 call this%zora%update(this%der, this%hm_base%potential)
1137 end if
1138
1139 !$omp parallel private(ip, ispin)
1140 do ispin = 1, this%d%nspin
1141 ! Adding Zeeman potential to hm_base%potential
1142 if (allocated(this%hm_base%zeeman_pot)) then
1143 !$omp do simd schedule(static)
1144 do ip = 1, mesh%np
1145 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + this%hm_base%zeeman_pot(ip, ispin)
1146 end do
1147 !$omp end do simd nowait
1148 end if
1149
1150 ! Adding Quadrupole potential from static E-field (test)
1151 if (this%mxll%test_equad) then
1152 !$omp do simd schedule(static)
1153 do ip = 1, mesh%np
1154 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + this%mxll%e_quadrupole_pot(ip)
1155 end do
1156 !$omp end do simd
1157 end if
1158 end do
1159 !$omp end parallel
1160
1162 ! Add the Hartree and KS potential
1163 call this%ks_pot%add_vhxc(this%hm_base%potential)
1164
1165 call this%hm_base%accel_copy_pot(mesh)
1166
1168 end subroutine hamiltonian_elec_update_pot
1169
1170 ! ---------------------------------------------------------
1171 subroutine hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
1172 type(hamiltonian_elec_t), intent(inout) :: this
1173 type(namespace_t), intent(in) :: namespace
1174 class(electron_space_t), intent(in) :: space
1175 type(grid_t), intent(in) :: gr
1176 type(ions_t), target, intent(inout) :: ions
1177 type(partner_list_t), intent(in) :: ext_partners
1178 type(states_elec_t), intent(inout) :: st
1179 real(real64), optional, intent(in) :: time
1180
1182
1183 this%ions => ions
1184 call epot_generate(this%ep, namespace, gr, this%ions, this%d)
1185
1186 ! Interation terms are treated below
1187
1188 ! First we add the static electric field
1189 if (allocated(this%ep%e_field) .and. space%periodic_dim < space%dim) then
1190 call lalg_axpy(gr%np, m_one, this%v_static, this%ep%vpsl)
1191 end if
1192
1193 ! Here we need to pass this again, else test are failing.
1194 ! This is not a real problem, as the multisystem framework will indeed to this anyway
1195 this%v_ie_loc%atoms_dist => ions%atoms_dist
1196 this%v_ie_loc%atom => ions%atom
1197 call this%v_ie_loc%calculate()
1198
1199 ! At the moment we need to add this to ep%vpsl, to keep the behavior of the code
1200 call lalg_axpy(gr%np, m_one, this%v_ie_loc%potential(:,1), this%ep%vpsl)
1201
1202 ! Here we need to reinit the NLCC object
1203 ! This is not a real problem, as the multisystem framework will indeed to this anyway
1204 if (this%ep%nlcc) then
1205 call this%nlcc%end()
1206 call this%nlcc%init(gr, ions)
1207 call this%nlcc%calculate()
1208 call lalg_copy(gr%np, this%nlcc%density(:,1), st%rho_core)
1209 end if
1210
1211 call this%vnl%build(space, gr, this%ep)
1212 call hamiltonian_elec_update(this, gr, namespace, space, ext_partners, time)
1213
1214 ! Check if projectors are still compatible with apply_packed on GPU
1215 if (this%is_applied_packed .and. accel_is_enabled()) then
1216 if (this%ep%non_local .and. .not. this%vnl%apply_projector_matrices) then
1217 if (accel_allow_cpu_only()) then
1218 call messages_write('Relativistic pseudopotentials have not been fully implemented for CUDA or OpenCL.')
1219 call messages_warning(namespace=namespace)
1220 else
1221 call messages_write('Relativistic pseudopotentials have not been fully implemented for CUDA or OpenCL.',&
1222 new_line = .true.)
1223 call messages_write('Calculation will not be continued. To force execution, set AllowCPUonly = yes.')
1224 call messages_fatal(namespace=namespace)
1225 end if
1226 end if
1227
1228 end if
1229
1230 if (this%pcm%run_pcm) then
1233 if (this%pcm%solute) then
1234 call pcm_calc_pot_rs(this%pcm, gr, this%psolver, ions = ions)
1235 end if
1236
1239 ! Interpolation is needed, hence gr%np_part -> 1:gr%np
1240 if (this%pcm%localf .and. allocated(this%v_static)) then
1241 call pcm_calc_pot_rs(this%pcm, gr, this%psolver, v_ext = this%ep%v_ext(1:gr%np_part))
1242 end if
1243
1244 end if
1245
1246 call lda_u_update_basis(this%lda_u, space, gr, ions, st, this%psolver, namespace, this%kpoints, &
1247 this%phase%is_allocated())
1248
1250 end subroutine hamiltonian_elec_epot_generate
1251
1252 ! -----------------------------------------------------------------
1253
1254 real(real64) function hamiltonian_elec_get_time(this) result(time)
1255 type(hamiltonian_elec_t), intent(inout) :: this
1256
1257 time = this%current_time
1258 end function hamiltonian_elec_get_time
1259
1260 ! -----------------------------------------------------------------
1261
1262 pure logical function hamiltonian_elec_apply_packed(this) result(apply)
1263 class(hamiltonian_elec_t), intent(in) :: this
1265 apply = this%is_applied_packed
1266
1268
1269
1270 ! -----------------------------------------------------------------
1271 subroutine zhamiltonian_elec_apply_atom (hm, namespace, space, latt, species, pos, ia, mesh, psi, vpsi)
1272 type(hamiltonian_elec_t), intent(in) :: hm
1273 type(namespace_t), intent(in) :: namespace
1274 class(space_t), intent(in) :: space
1275 type(lattice_vectors_t), intent(in) :: latt
1276 class(species_t), intent(in) :: species
1277 real(real64), intent(in) :: pos(1:space%dim)
1278 integer, intent(in) :: ia
1279 class(mesh_t), intent(in) :: mesh
1280 complex(real64), intent(in) :: psi(:,:)
1281 complex(real64), intent(out) :: vpsi(:,:)
1282
1283 integer :: idim
1284 real(real64), allocatable :: vlocal(:)
1286
1287 safe_allocate(vlocal(1:mesh%np_part))
1288 vlocal = m_zero
1289 call epot_local_potential(hm%ep, namespace, space, latt, mesh, species, pos, ia, vlocal)
1290
1291 do idim = 1, hm%d%dim
1292 vpsi(1:mesh%np, idim) = vlocal(1:mesh%np) * psi(1:mesh%np, idim)
1293 end do
1294
1295 safe_deallocate_a(vlocal)
1297 end subroutine zhamiltonian_elec_apply_atom
1298
1299 ! ---------------------------------------------------------
1304 subroutine hamiltonian_elec_update_with_ext_pot(this, mesh, space, ext_partners, time, mu)
1305 type(hamiltonian_elec_t), intent(inout) :: this
1306 class(space_t), intent(in) :: space
1307 class(mesh_t), intent(in) :: mesh
1308 type(partner_list_t), intent(in) :: ext_partners
1309 real(real64), intent(in) :: time(1:2)
1310 real(real64), intent(in) :: mu(1:2)
1311
1312 integer :: ispin, ip, idir, iatom, ilaser, itime
1313 real(real64) :: aa(space%dim), bb(space%dim), time_
1314 real(real64), allocatable :: vp(:,:)
1315 real(real64), allocatable :: velectric(:)
1316 type(lasers_t), pointer :: lasers
1317 type(gauge_field_t), pointer :: gfield
1318
1320 call profiling_in("HAMILTONIAN_ELEC_UPDATE_EXT_POT")
1321
1322 this%current_time = m_zero
1323 this%current_time = time(1)
1324
1325 ! set everything to zero
1326 call this%hm_base%clear(mesh%np)
1327
1328 ! the xc, hartree and external potentials
1329 call this%hm_base%allocate_field(mesh, field_potential, &
1330 complex_potential = this%abs_boundaries%abtype == imaginary_absorbing)
1331
1332 do itime = 1, 2
1333 time_ = time(itime)
1334
1335 lasers => list_get_lasers(ext_partners)
1336 if(associated(lasers)) then
1337 do ilaser = 1, lasers%no_lasers
1338 select case (laser_kind(lasers%lasers(ilaser)))
1339 case (e_field_scalar_potential, e_field_electric)
1340 safe_allocate(velectric(1:mesh%np))
1341 do ispin = 1, this%d%spin_channels
1342 velectric = m_zero
1343 call laser_potential(lasers%lasers(ilaser), mesh, velectric, time_)
1344 !$omp parallel do simd schedule(static)
1345 do ip = 1, mesh%np
1346 this%hm_base%potential(ip, ispin) = this%hm_base%potential(ip, ispin) + mu(itime) * velectric(ip)
1347 end do
1348 end do
1349 safe_deallocate_a(velectric)
1350 case (e_field_magnetic)
1351 call this%hm_base%allocate_field(mesh, field_vector_potential + field_uniform_magnetic_field, .false.)
1352 ! get the vector potential
1353 safe_allocate(vp(1:mesh%np, 1:space%dim))
1354 vp(1:mesh%np, 1:space%dim) = m_zero
1355 call laser_vector_potential(lasers%lasers(ilaser), mesh, vp, time_)
1356 do idir = 1, space%dim
1357 !$omp parallel do schedule(static)
1358 do ip = 1, mesh%np
1359 this%hm_base%vector_potential(idir, ip) = this%hm_base%vector_potential(idir, ip) &
1360 - mu(itime) * vp(ip, idir)/p_c
1361 end do
1362 end do
1363 ! and the magnetic field
1364 bb = m_zero
1365 call laser_field(lasers%lasers(ilaser), bb(1:space%dim), time_)
1366 this%hm_base%uniform_magnetic_field = this%hm_base%uniform_magnetic_field &
1367 - mu(itime) * bb
1368 safe_deallocate_a(vp)
1369 case (e_field_vector_potential)
1370 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
1371 ! get the uniform vector potential associated with a magnetic field
1372 aa = m_zero
1373 call laser_field(lasers%lasers(ilaser), aa, time_)
1374 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) &
1375 - mu(itime) * aa/p_c
1376 end select
1377 end do
1378 end if
1379
1380 ! the gauge field
1381 gfield => list_get_gauge_field(ext_partners)
1382 if (associated(gfield)) then
1383 call this%hm_base%allocate_field(mesh, field_uniform_vector_potential, .false.)
1384 call gauge_field_get_vec_pot(gfield, aa)
1385 this%hm_base%uniform_vector_potential(1:space%dim) = this%hm_base%uniform_vector_potential(1:space%dim) - aa/p_c
1386 end if
1387
1388 ! the electric field for a periodic system through the gauge field
1389 ! TODO: The condition is wrong here: the e_field should be in non-periodic dims as E field
1390 ! and as a gauge field in the periodic dim, unless we use a Bery phase, in which, we do not use it
1391 ! this way. But this is unrelated to the gauge field
1392 if (allocated(this%ep%e_field) .and. associated(gfield)) then
1393 this%hm_base%uniform_vector_potential(1:space%periodic_dim) = &
1394 this%hm_base%uniform_vector_potential(1:space%periodic_dim) - time_*this%ep%e_field(1:space%periodic_dim)
1395 end if
1396
1397 end do
1398
1399 ! the vector potential of a static magnetic field
1400 if (allocated(this%ep%a_static)) then
1401 call this%hm_base%allocate_field(mesh, field_vector_potential, .false.)
1402 !ep%a_static contains 1/c A(r)
1403 !$omp parallel do schedule(static) private(idir)
1404 do ip = 1, mesh%np
1405 do idir = 1, space%dim
1406 this%hm_base%vector_potential(idir, ip) = this%hm_base%vector_potential(idir, ip) + this%ep%a_static(ip, idir)
1407 end do
1408 end do
1409 end if
1410
1411 !The electric field is added to the KS potential
1412 call this%hm_base%accel_copy_pot(mesh)
1413
1414 ! and the static magnetic field
1415 if (allocated(this%ep%b_field)) then
1416 call this%hm_base%allocate_field(mesh, field_uniform_magnetic_field, .false.)
1417 do idir = 1, 3
1418 this%hm_base%uniform_magnetic_field(idir) = this%hm_base%uniform_magnetic_field(idir) + this%ep%b_field(idir)
1419 end do
1420 end if
1421
1422 call this%hm_base%update_magnetic_terms(mesh, this%ep%gyromagnetic_ratio, this%d%ispin)
1423
1424 call hamiltonian_elec_update_pot(this, mesh)
1425
1426 call build_phase()
1427
1428 call profiling_out("HAMILTONIAN_ELEC_UPDATE_EXT_POT")
1430
1431 contains
1432
1433 subroutine build_phase()
1434 integer :: ik, imat, nmat, max_npoints, offset, iphase, nphase
1435
1437
1438 if ((.not. this%kpoints%gamma_only()) .or. allocated(this%hm_base%uniform_vector_potential)) then
1439
1440 call profiling_in('UPDATE_PHASES')
1441 ! now regenerate the phases for the pseudopotentials
1442 do iatom = 1, this%ep%natoms
1443 call projector_init_phases(this%ep%proj(iatom), space%dim, this%d, this%der%boundaries, this%kpoints, &
1444 vec_pot = this%hm_base%uniform_vector_potential, vec_pot_var = this%hm_base%vector_potential)
1445 end do
1446
1447 call profiling_out('UPDATE_PHASES')
1448 end if
1449
1450 if (allocated(this%hm_base%uniform_vector_potential)) then
1451 call this%phase%update(mesh, this%d%kpt, this%kpoints, this%d, space, this%hm_base%uniform_vector_potential)
1452 end if
1453
1454 max_npoints = this%vnl%max_npoints
1455 nmat = this%vnl%nprojector_matrices
1456
1457
1458 if (this%phase%is_allocated() .and. allocated(this%vnl%projector_matrices)) then
1459
1460 nphase = 1
1461 if (this%der%boundaries%spiralBC) nphase = 3
1462
1463 if (.not. allocated(this%vnl%projector_phases)) then
1464 safe_allocate(this%vnl%projector_phases(1:max_npoints, nphase, nmat, this%d%kpt%start:this%d%kpt%end))
1465 if (accel_is_enabled()) then
1466 call accel_create_buffer(this%vnl%buff_projector_phases, accel_mem_read_only, &
1467 type_cmplx, this%vnl%total_points*nphase*this%d%kpt%nlocal)
1468 end if
1469 end if
1470
1471 offset = 0
1472 do ik = this%d%kpt%start, this%d%kpt%end
1473 do imat = 1, this%vnl%nprojector_matrices
1474 iatom = this%vnl%projector_to_atom(imat)
1475 do iphase = 1, nphase
1476 !$omp parallel do schedule(static)
1477 do ip = 1, this%vnl%projector_matrices(imat)%npoints
1478 this%vnl%projector_phases(ip, iphase, imat, ik) = this%ep%proj(iatom)%phase(ip, iphase, ik)
1479 end do
1480
1481 if (accel_is_enabled() .and. this%vnl%projector_matrices(imat)%npoints > 0) then
1482 call accel_write_buffer(this%vnl%buff_projector_phases, &
1483 this%vnl%projector_matrices(imat)%npoints, this%vnl%projector_phases(1:, iphase, imat, ik), &
1484 offset = offset)
1485 end if
1486 offset = offset + this%vnl%projector_matrices(imat)%npoints
1487 end do
1488 end do
1489 end do
1490
1491 end if
1492
1494 end subroutine build_phase
1495
1497
1498 logical function hamiltonian_elec_needs_current(hm, states_are_real)
1499 type(hamiltonian_elec_t), intent(in) :: hm
1500 logical, intent(in) :: states_are_real
1501
1503
1504 if (hm%self_induced_magnetic) then
1505 if (.not. states_are_real) then
1507 else
1508 message(1) = 'No current density for real states since it is identically zero.'
1509 call messages_warning(1)
1510 end if
1511 end if
1512
1514
1515 ! ---------------------------------------------------------
1516 subroutine zhamiltonian_elec_apply_all(hm, namespace, mesh, st, hst)
1517 type(hamiltonian_elec_t), intent(inout) :: hm
1518 type(namespace_t), intent(in) :: namespace
1519 class(mesh_t), intent(in) :: mesh
1520 type(states_elec_t), intent(inout) :: st
1521 type(states_elec_t), intent(inout) :: hst
1522
1523 integer :: ik, ib, ist
1524 complex(real64), allocatable :: psi(:, :)
1525 complex(real64), allocatable :: psiall(:, :, :, :)
1526
1528
1529 do ik = st%d%kpt%start, st%d%kpt%end
1530 do ib = st%group%block_start, st%group%block_end
1531 call zhamiltonian_elec_apply_batch(hm, namespace, mesh, st%group%psib(ib, ik), hst%group%psib(ib, ik))
1532 end do
1533 end do
1534
1535 if (oct_exchange_enabled(hm%oct_exchange)) then
1536
1537 safe_allocate(psiall(mesh%np_part, 1:hst%d%dim, st%st_start:st%st_end, st%d%kpt%start:st%d%kpt%end))
1538
1539 call states_elec_get_state(st, mesh, psiall)
1540
1541 call oct_exchange_prepare(hm%oct_exchange, mesh, psiall, hm%xc, hm%psolver, namespace)
1542
1543 safe_deallocate_a(psiall)
1544
1545 safe_allocate(psi(mesh%np_part, 1:hst%d%dim))
1546
1547 do ik = 1, st%nik
1548 do ist = 1, st%nst
1549 call states_elec_get_state(hst, mesh, ist, ik, psi)
1550 call oct_exchange_operator(hm%oct_exchange, namespace, mesh, psi, ist, ik)
1551 call states_elec_set_state(hst, mesh, ist, ik, psi)
1552 end do
1553 end do
1554
1555 safe_deallocate_a(psi)
1556
1557 end if
1558
1560 end subroutine zhamiltonian_elec_apply_all
1561
1562
1563 ! ---------------------------------------------------------
1564
1565 subroutine magnus(hm, namespace, mesh, psi, hpsi, ik, vmagnus, set_phase)
1566 type(hamiltonian_elec_t), intent(in) :: hm
1567 type(namespace_t), intent(in) :: namespace
1568 class(mesh_t), intent(in) :: mesh
1569 complex(real64), contiguous, intent(inout) :: psi(:,:)
1570 complex(real64), contiguous, intent(out) :: hpsi(:,:)
1571 integer, intent(in) :: ik
1572 real(real64), intent(in) :: vmagnus(:, :, :)
1573 logical, optional, intent(in) :: set_phase
1574
1575 complex(real64), allocatable :: auxpsi(:, :), aux2psi(:, :)
1576 integer :: idim, ispin
1577
1578 push_sub(magnus)
1579
1580 ! We will assume, for the moment, no spinors.
1581 if (hm%d%dim /= 1) then
1582 call messages_not_implemented("Magnus with spinors", namespace=namespace)
1583 end if
1584
1585 safe_allocate( auxpsi(1:mesh%np_part, 1:hm%d%dim))
1586 safe_allocate(aux2psi(1:mesh%np, 1:hm%d%dim))
1587
1588 ispin = hm%d%get_spin_index(ik)
1589
1590 ! Compute (T + Vnl)|psi> and store it
1591 call zhamiltonian_elec_apply_single(hm, namespace, mesh, psi, auxpsi, 1, ik, &
1592 terms = term_kinetic + term_non_local_potential, set_phase = set_phase)
1593
1594 ! H|psi> = (T + Vnl)|psi> + Vpsl|psi> + Vmagnus(t2)|psi> + Vborders
1595 do idim = 1, hm%d%dim
1596 call lalg_copy(mesh%np, auxpsi(:, idim), hpsi(:, idim))
1597 hpsi(1:mesh%np, idim) = hpsi(1:mesh%np, idim) + hm%ep%Vpsl(1:mesh%np)*psi(1:mesh%np,idim)
1598 call vborders(mesh, hm, psi(:, idim), hpsi(:, idim))
1599 end do
1600 hpsi(1:mesh%np, 1) = hpsi(1:mesh%np, 1) + vmagnus(1:mesh%np, ispin, 2)*psi(1:mesh%np, 1)
1601
1602 ! Add first term of the commutator: - i Vmagnus(t1) (T + Vnl) |psi>
1603 hpsi(1:mesh%np, 1) = hpsi(1:mesh%np, 1) - m_zi*vmagnus(1:mesh%np, ispin, 1)*auxpsi(1:mesh%np, 1)
1604
1605 ! Add second term of commutator: i (T + Vnl) Vmagnus(t1) |psi>
1606 auxpsi(1:mesh%np, 1) = vmagnus(1:mesh%np, ispin, 1)*psi(1:mesh%np, 1)
1607 call zhamiltonian_elec_apply_single(hm, namespace, mesh, auxpsi, aux2psi, 1, ik, &
1608 terms = term_kinetic + term_non_local_potential, set_phase = set_phase)
1609 hpsi(1:mesh%np, 1) = hpsi(1:mesh%np, 1) + m_zi*aux2psi(1:mesh%np, 1)
1610
1611 safe_deallocate_a(auxpsi)
1612 safe_deallocate_a(aux2psi)
1613 pop_sub(magnus)
1614 end subroutine magnus
1615
1616 ! ---------------------------------------------------------
1617 subroutine vborders (mesh, hm, psi, hpsi)
1618 class(mesh_t), intent(in) :: mesh
1619 type(hamiltonian_elec_t), intent(in) :: hm
1620 complex(real64), intent(in) :: psi(:)
1621 complex(real64), intent(inout) :: hpsi(:)
1622
1623 integer :: ip
1624
1625 push_sub(vborders)
1626
1627 if (hm%abs_boundaries%abtype == imaginary_absorbing) then
1628 do ip = 1, mesh%np
1629 hpsi(ip) = hpsi(ip) + m_zi*hm%abs_boundaries%mf(ip)*psi(ip)
1630 end do
1631 end if
1632
1633 pop_sub(vborders)
1634 end subroutine vborders
1635
1636 ! ---------------------------------------------------------
1637 logical function hamiltonian_elec_has_kick(hm)
1638 type(hamiltonian_elec_t), intent(in) :: hm
1639
1641
1642 hamiltonian_elec_has_kick = (abs(hm%kick%delta_strength) > m_epsilon)
1643
1645 end function hamiltonian_elec_has_kick
1646
1648 !
1649 subroutine hamiltonian_elec_set_mass(this, namespace, mass)
1650 class(hamiltonian_elec_t) , intent(inout) :: this
1651 type(namespace_t), intent(in) :: namespace
1652 real(real64), intent(in) :: mass
1653
1655
1656 if (parse_is_defined(namespace, 'ParticleMass')) then
1657 message(1) = 'Attempting to redefine a non-unit electron mass'
1658 call messages_fatal(1)
1659 else
1660 this%mass = mass
1661 end if
1662
1664 end subroutine hamiltonian_elec_set_mass
1665
1666 !----------------------------------------------------------
1673 subroutine hamiltonian_elec_copy_and_set_phase(hm, gr, kpt, psib, psib_with_phase)
1674 type(hamiltonian_elec_t), intent(in) :: hm
1675 type(grid_t), intent(in) :: gr
1676 type(distributed_t), intent(in) :: kpt
1677 type(wfs_elec_t), intent(in) :: psib
1678 type(wfs_elec_t), intent(out) :: psib_with_phase
1679
1680 integer :: k_offset, n_boundary_points
1681
1683
1684 call psib%copy_to(psib_with_phase)
1685 if (hm%phase%is_allocated()) then
1686 call hm%phase%apply_to(gr, gr%np, conjugate = .false., psib = psib_with_phase, src = psib, async=.true.)
1687 ! apply phase correction while setting boundary -> memory needs to be
1688 ! accessed only once
1689 k_offset = psib%ik - kpt%start
1690 n_boundary_points = int(gr%np_part - gr%np)
1691 call boundaries_set(gr%der%boundaries, gr, psib_with_phase, phase_correction = hm%phase%phase_corr(:, psib%ik), &
1692 buff_phase_corr = hm%phase%buff_phase_corr, offset=k_offset*n_boundary_points, async=.true.)
1693 else
1694 call psib%copy_data_to(gr%np, psib_with_phase)
1695 call boundaries_set(gr%der%boundaries, gr, psib_with_phase)
1696 end if
1697
1698 call psib_with_phase%do_pack(copy = .true.)
1699
1702
1703 ! ---------------------------------------------------------
1704 subroutine hamiltonian_elec_diagonal (hm, mesh, diag, ik)
1705 type(hamiltonian_elec_t), intent(in) :: hm
1706 class(mesh_t), intent(in) :: mesh
1707 real(real64), contiguous, intent(out) :: diag(:,:)
1708 integer, intent(in) :: ik
1709
1710 integer :: idim, ip, ispin
1711
1712 real(real64), allocatable :: ldiag(:)
1713
1715
1716 safe_allocate(ldiag(1:mesh%np))
1717
1718 diag = m_zero
1719
1720 call derivatives_lapl_diag(hm%der, ldiag)
1721
1722 do idim = 1, hm%d%dim
1723 diag(1:mesh%np, idim) = -m_half/hm%mass*ldiag(1:mesh%np)
1724 end do
1725
1726 select case (hm%d%ispin)
1727
1728 case (unpolarized, spin_polarized)
1729 ispin = hm%d%get_spin_index(ik)
1730 diag(1:mesh%np, 1) = diag(1:mesh%np, 1) + hm%ep%vpsl(1:mesh%np)
1731
1732 case (spinors)
1733 do ip = 1, mesh%np
1734 diag(ip, 1) = diag(ip, 1) + hm%ep%vpsl(ip)
1735 diag(ip, 2) = diag(ip, 2) + hm%ep%vpsl(ip)
1736 end do
1737
1738 end select
1739
1740 call hm%ks_pot%add_vhxc(diag)
1741
1743 end subroutine hamiltonian_elec_diagonal
1744
1745
1746
1747#include "undef.F90"
1748#include "real.F90"
1749#include "hamiltonian_elec_inc.F90"
1750
1751#include "undef.F90"
1752#include "complex.F90"
1753#include "hamiltonian_elec_inc.F90"
1754
1755end module hamiltonian_elec_oct_m
1756
1757!! Local Variables:
1758!! mode: f90
1759!! coding: utf-8
1760!! 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:171
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:186
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:411
subroutine, public accel_finish()
Definition: accel.F90:1298
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:122
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,...
logical function, public epot_have_external_potentials(ep)
Definition: epot.F90:602
integer, parameter, public scalar_relativistic_zora
Definition: epot.F90:166
subroutine, public epot_end(ep)
Definition: epot.F90:381
integer, parameter, public fully_relativistic_zora
Definition: epot.F90:166
subroutine, public epot_init(ep, namespace, gr, ions, psolver, ispin, xc_family, kpoints)
Definition: epot.F90:218
subroutine, public epot_generate(ep, namespace, mesh, ions, st_d)
Definition: epot.F90:417
subroutine, public exchange_operator_init(this, namespace, space, st, der, mc, stencil, kpoints, omega, alpha, beta)
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 gauge_field_get_vec_pot(this, vec_pot)
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:224
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
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 field_potential
subroutine vborders(mesh, hm, psi, hpsi)
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 hamiltonian_elec_copy_and_set_phase(hm, gr, kpt, psib, psib_with_phase)
Copy a batch to another batch and apply the Bloch phase to it.
subroutine, public zhamiltonian_elec_apply_single(hm, namespace, mesh, psi, hpsi, ist, ik, terms, set_bc, set_phase)
logical function, public hamiltonian_elec_has_kick(hm)
subroutine, public magnus(hm, namespace, mesh, psi, hpsi, ik, vmagnus, set_phase)
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 zhamiltonian_elec_apply_all(hm, namespace, mesh, st, hst)
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 zhamiltonian_elec_magnus_apply(hm, namespace, mesh, psib, hpsib, vmagnus)
subroutine, public zhamiltonian_elec_apply_atom(hm, namespace, space, latt, species, pos, ia, mesh, psi, vpsi)
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 dhamiltonian_elec_magnus_apply(hm, namespace, mesh, psib, hpsib, vmagnus)
subroutine hamiltonian_elec_update(this, mesh, namespace, space, ext_partners, time)
(re-)build the Hamiltonian for the next application:
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 zhamiltonian_elec_apply(hm, namespace, mesh, psib, hpsib, terms, set_bc)
This module defines classes and functions for interaction partners.
Definition: io.F90:114
integer, parameter, public kick_magnon_mode
Definition: kick.F90:163
subroutine, public kick_end(kick)
Definition: kick.F90:795
subroutine, public kick_init(kick, namespace, space, kpoints, nspin)
Definition: kick.F90:223
pure integer function, public kick_get_type(kick)
Definition: kick.F90:1364
A module to handle KS potential, without the external potential.
integer, parameter, public rdmft
integer, parameter, public hartree
integer, parameter, public hartree_fock
integer, parameter, public independent_particles
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public kohn_sham_dft
subroutine, public laser_vector_potential(laser, mesh, aa, time)
Definition: lasers.F90:1082
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:1154
logical function, public lasers_with_nondipole_field(lasers)
Check if a nondipole SFA correction should be computed for the given laser.
Definition: lasers.F90:741
integer, parameter, public e_field_electric
Definition: lasers.F90:177
integer, parameter, public e_field_vector_potential
Definition: lasers.F90:177
subroutine, public laser_potential(laser, mesh, pot, time)
Definition: lasers.F90:1047
integer, parameter, public e_field_scalar_potential
Definition: lasers.F90:177
integer pure elemental function, public laser_kind(laser)
Definition: lasers.F90:717
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:1119
integer, parameter, public e_field_magnetic
Definition: lasers.F90:177
integer, parameter, public dft_u_none
Definition: lda_u.F90:201
subroutine, public lda_u_update_basis(this, space, gr, ions, st, psolver, namespace, kpoints, has_phase)
Definition: lda_u.F90:696
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:824
subroutine, public lda_u_end(this)
Definition: lda_u.F90:654
subroutine, public lda_u_init(this, namespace, space, level, gr, ions, st, mc, kpoints, has_phase)
Definition: lda_u.F90:281
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_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:115
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_not_implemented(feature, namespace)
Definition: messages.F90:1114
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1086
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
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.
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:115
Some general things and nomenclature:
Definition: par_vec.F90:171
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
subroutine, public pcm_calc_pot_rs(pcm, mesh, psolver, ions, v_h, v_ext, kick, time_present, kick_time)
Definition: pcm.F90:1215
subroutine, public pcm_end(pcm)
Definition: pcm.F90:3071
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:294
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:709
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
subroutine, public projector_init_phases(this, dim, std, bnd, kpoints, vec_pot, vec_pot_var)
Definition: projector.F90:264
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:589
integer pure function, public symmetries_number(this)
Definition: symmetries.F90:543
type(type_t), parameter, public type_cmplx
Definition: types.F90:134
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
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
integer, parameter, public xc_oep_x_slater
Slater approximation to the exact exchange.
integer, parameter, public func_x
Definition: xc.F90:114
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
Definition: xc.F90:595
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
Definition: xc.F90:610
This module implements the "photon-free" electron-photon exchange-correlation functional.
Definition: xc_photons.F90:121
This module implements the ZORA terms for the Hamoiltonian.
Definition: zora.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 abstract Hamiltonian class defines a skeleton for specific implementations.
abstract class for general interaction partners
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Stores all communicators and groups.
Definition: multicomm.F90:206
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:156
This class is responsible for calculating and applying the ZORA.
Definition: zora.F90:145
int true(void)