Octopus
hamiltonian_elec_base.F90
Go to the documentation of this file.
1!! Copyright (C) 2009 X. Andrade
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use accel_oct_m
23 use batch_oct_m
25 use blas_oct_m
26 use comm_oct_m
27 use debug_oct_m
30 use epot_oct_m
31 use global_oct_m
33 use math_oct_m
34 use mesh_oct_m
36 use mpi_oct_m
39 use ps_oct_m
40 use space_oct_m
43 use types_oct_m
45
46 implicit none
47
48 private
49
50 public :: &
55
65 private
66 integer :: nspin
67 real(real64) :: mass
68 real(real64) :: rashba_coupling
69 type(nl_operator_t), pointer, public :: kinetic
70 real(real64), allocatable, public :: potential(:, :)
71 real(real64), allocatable, public :: Impotential(:, :)
72 real(real64), allocatable, public :: uniform_magnetic_field(:)
74 real(real64), allocatable, public :: magnetic_field(:, :)
76 real(real64), allocatable, public :: zeeman_pot(:, :)
77 real(real64), allocatable, public :: uniform_vector_potential(:)
81 real(real64), allocatable, public :: vector_potential(:, :)
82 type(accel_mem_t) :: potential_accel
83 type(accel_mem_t) :: impotential_accel
84 contains
85
86 procedure :: init => hamiltonian_elec_base_init
87
88 procedure :: end => hamiltonian_elec_base_end
89
90 procedure :: clear => hamiltonian_elec_base_clear
91
92 procedure :: update_magnetic_terms => hamiltonian_elec_base_update_magnetic_terms
93
94 procedure :: allocate_field => hamiltonian_elec_base_allocate
95
96 procedure :: accel_copy_pot => hamiltonian_elec_base_accel_copy_pot
97
98 procedure :: has_magnetic => hamiltonian_elec_base_has_magnetic
99
100 procedure :: has_zeeman => hamiltonian_elec_base_has_zeeman
101
102 procedure :: has_vector_potential => hamiltonian_elec_base_has_vector_potential
103
104 procedure :: calc_rashba => hamiltonian_elec_base_rashba
105
106 procedure :: dcalc_magnetic => dhamiltonian_elec_base_magnetic
107
108 procedure :: zcalc_magnetic => zhamiltonian_elec_base_magnetic
109
110 procedure :: dcalc_local => dhamiltonian_elec_base_local
111
112 procedure :: zcalc_local => zhamiltonian_elec_base_local
113
115
116 integer, parameter, public :: &
117 TERM_ALL = huge(1), &
118 term_kinetic = 1, &
121 term_others = 8, &
122 term_local_external = 16, &
123 term_mgga = 32, &
124 term_dft_u = 64, &
125 term_rdmft_occ = 128
126
127 integer, parameter, public :: &
128 FIELD_POTENTIAL = 1, &
133
134contains
135
136 ! ---------------------------------------------------------
139 subroutine hamiltonian_elec_base_init(this, nspin, mass, rashba_coupling)
140 class(hamiltonian_elec_base_t), intent(inout) :: this
141 integer, intent(in) :: nspin
142 real(real64), intent(in) :: mass
143 real(real64), intent(in) :: rashba_coupling
144
146
147 this%nspin = nspin
148 this%mass = mass
149 this%rashba_coupling = rashba_coupling
150
152 end subroutine hamiltonian_elec_base_init
153
154 ! ---------------------------------------------------------
157 subroutine hamiltonian_elec_base_end(this)
158 class(hamiltonian_elec_base_t), intent(inout) :: this
162 if (allocated(this%potential) .and. accel_is_enabled()) then
163 call accel_free_buffer(this%potential_accel)
164 if (allocated(this%Impotential)) then
165 call accel_free_buffer(this%impotential_accel)
166 end if
167 end if
168
169 safe_deallocate_a(this%potential)
170 safe_deallocate_a(this%Impotential)
171 safe_deallocate_a(this%vector_potential)
172 safe_deallocate_a(this%uniform_vector_potential)
173 safe_deallocate_a(this%uniform_magnetic_field)
174 safe_deallocate_a(this%magnetic_field)
175 safe_deallocate_a(this%zeeman_pot)
179
180 ! ----------------------------------------------------------
184 !
185 subroutine hamiltonian_elec_base_clear(this, np)
186 class(hamiltonian_elec_base_t), intent(inout) :: this
187 integer, intent(in) :: np
188
189 integer :: ip, ispin
190
191 push_sub(hamiltonian_elec_clear)
192
193 if (allocated(this%potential)) then
194 do ispin = 1, this%nspin
195 !$omp parallel do schedule(static)
196 do ip = 1, np
197 this%potential(ip, ispin) = m_zero
198 end do
199 end do
200 end if
202 if (allocated(this%Impotential)) then
203 do ispin = 1, this%nspin
204 !$omp parallel do schedule(static)
205 do ip = 1, np
206 this%Impotential(ip, ispin) = m_zero
207 end do
208 end do
209 end if
210
211 if (allocated(this%uniform_vector_potential)) then
212 this%uniform_vector_potential = m_zero
213 end if
214
215 if (allocated(this%vector_potential)) then
216 this%vector_potential = m_zero
217 end if
218
219 if (allocated(this%uniform_magnetic_field)) then
220 this%uniform_magnetic_field = m_zero
221 end if
223 if (allocated(this%magnetic_field)) then
224 this%magnetic_field = m_zero
225 end if
226
227 if (allocated(this%zeeman_pot)) then
228 this%zeeman_pot = m_zero
229 end if
230
231 pop_sub(hamiltonian_elec_clear)
232 end subroutine hamiltonian_elec_base_clear
233
235 ! ---------------------------------------------------------------
237 !
238 subroutine hamiltonian_elec_base_allocate(this, mesh, field, complex_potential)
239 class(hamiltonian_elec_base_t), intent(inout) :: this
240 class(mesh_t), intent(in) :: mesh
241 integer, intent(in) :: field
242 logical, intent(in) :: complex_potential
243
244 integer :: ip, ispin
245
247
248 if (bitand(field_potential, field) /= 0) then
249 if (.not. allocated(this%potential)) then
250 safe_allocate(this%potential(1:mesh%np, 1:this%nspin))
251 do ispin = 1, this%nspin
252 !$omp parallel do schedule(static)
253 do ip = 1, mesh%np
254 this%potential(ip, ispin) = m_zero
255 end do
256 end do
257 if (complex_potential) then
258 safe_allocate(this%Impotential(1:mesh%np, 1:this%nspin))
259 do ispin = 1, this%nspin
260 !$omp parallel do schedule(static)
261 do ip = 1, mesh%np
262 this%Impotential(ip, ispin) = m_zero
263 end do
264 end do
265 end if
266 if (accel_is_enabled()) then
267 call accel_create_buffer(this%potential_accel, accel_mem_read_only, type_float, accel_padded_size(mesh%np)*this%nspin)
268 if (complex_potential) then
269 call accel_create_buffer(this%impotential_accel, accel_mem_read_only, type_float, &
270 accel_padded_size(mesh%np)*this%nspin)
271 end if
272 end if
273 end if
274 end if
275
276 if (bitand(field_uniform_vector_potential, field) /= 0) then
277 if (.not. allocated(this%uniform_vector_potential)) then
278 safe_allocate(this%uniform_vector_potential(1:mesh%box%dim))
279 this%uniform_vector_potential = m_zero
280 end if
281 end if
282
283 if (bitand(field_vector_potential, field) /= 0) then
284 if (.not. allocated(this%vector_potential)) then
285 safe_allocate(this%vector_potential(1:mesh%box%dim, 1:mesh%np))
286 this%vector_potential = m_zero
287 end if
288 end if
289
290 if (bitand(field_uniform_magnetic_field, field) /= 0) then
291 if (.not. allocated(this%uniform_magnetic_field)) then
292 safe_allocate(this%uniform_magnetic_field(1:max(mesh%box%dim, 3)))
293 this%uniform_magnetic_field = m_zero
294 end if
295 end if
296
297 if ((bitand(field_magnetic_field, field) /= 0) .or. &
298 bitand(field_uniform_magnetic_field, field) /= 0) then
299 if (.not. allocated(this%magnetic_field)) then
300 safe_allocate(this%magnetic_field(1:mesh%np, 1:max(mesh%box%dim, 3)))
301 this%magnetic_field = m_zero
302 safe_allocate(this%zeeman_pot(1:mesh%np, 1:this%nspin))
303 this%zeeman_pot = m_zero
304 end if
305 end if
306
308 end subroutine hamiltonian_elec_base_allocate
309
310 ! ----------------------------------------------------------
311 !
317 !
318 subroutine hamiltonian_elec_base_update_magnetic_terms(this, mesh, gyromagnetic_ratio, ispin)
319 class(hamiltonian_elec_base_t), intent(inout) :: this
320 class(mesh_t), intent(in) :: mesh
321 real(real64), intent(in) :: gyromagnetic_ratio
322 integer, intent(in) :: ispin
323
324 real(real64) :: b_norm2, cc
325 integer :: idir, ip
326
328
329 if (allocated(this%uniform_vector_potential) .and. allocated(this%vector_potential)) then
330 ! copy the uniform vector potential onto the non-uniform one
331 do idir = 1, mesh%box%dim
332 !$omp parallel do schedule(static)
333 do ip = 1, mesh%np
334 this%vector_potential(idir, ip) = &
335 this%vector_potential(idir, ip) + this%uniform_vector_potential(idir)
336 end do
337 end do
338 safe_deallocate_a(this%uniform_vector_potential)
339 end if
340
341
342 if (.not. allocated(this%magnetic_field)) then
344 return
345 end if
346
347 ! add the Zeeman term
348 ! It reads in the Pauli equation as $\frac{|e|\hbar}{2m} \sigma\cdot\vec{B}$
349 ! The gyromagnetic ratio is twice the magnetic moment of the electron in units of Bohr magneton
350 ! This is why we multiply by $\frac{|g_e|}{2}$.
351 ! The factor 1/c comes from the factor 1/c in front of the vector potential, which leads to
352 ! the Zeeman term starting from the Dirac Hamiltonian.
353 cc = m_half/p_c*gyromagnetic_ratio*m_half
354
355 if (allocated(this%uniform_magnetic_field) ) then
356 !$omp parallel private(ip)
357 do idir = 1, 3
358 !$omp do
359 do ip = 1, mesh%np
360 this%magnetic_field(ip, idir) = this%magnetic_field(ip, idir) + this%uniform_magnetic_field(idir)
361 end do
362 !$omp end do nowait
363 end do
364 !$omp end parallel
365 end if
366
367 select case (ispin)
368 case (spin_polarized)
369 !$omp parallel do private(b_norm2)
370 do ip = 1, mesh%np
371 b_norm2 = norm2(this%magnetic_field(ip, 1:max(mesh%box%dim, 3)))
372 this%zeeman_pot(ip, 1) = cc*b_norm2
373 this%zeeman_pot(ip, 2) = - cc*b_norm2
374 end do
375
376 case (spinors)
377 !$omp parallel do
378 do ip = 1, mesh%np
379 this%zeeman_pot(ip, 1) = cc*this%magnetic_field(ip, 3)
380 this%zeeman_pot(ip, 2) = - cc*this%magnetic_field(ip, 3)
381 this%zeeman_pot(ip, 3) = cc*this%magnetic_field(ip, 1)
382 this%zeeman_pot(ip, 4) = - cc*this%magnetic_field(ip, 2)
383 end do
384 end select
385 safe_deallocate_a(this%uniform_magnetic_field)
386
389
390
391 !--------------------------------------------------------
393 !
394 subroutine hamiltonian_elec_base_accel_copy_pot(this, mesh)
395 class(hamiltonian_elec_base_t), intent(inout) :: this
396 class(mesh_t), intent(in) :: mesh
397
398 integer :: offset, ispin
399
401
402 if (allocated(this%potential) .and. accel_is_enabled()) then
403 offset = 0
404 do ispin = 1, this%nspin
405 call accel_write_buffer(this%potential_accel, mesh%np, this%potential(:, ispin), offset = offset)
406 if(allocated(this%Impotential)) then
407 call accel_write_buffer(this%impotential_accel, mesh%np, this%Impotential(:, ispin), offset = offset)
408 end if
409 offset = offset + accel_padded_size(mesh%np)
410 end do
411 end if
412
415
416 ! ---------------------------------------------------------
418 subroutine hamiltonian_elec_base_accel_rebuild(this, mesh)
419 class(hamiltonian_elec_base_t), intent(inout) :: this
420 class(mesh_t), intent(in) :: mesh
421
422 logical :: complex_potential
423
425
426 if (.not. accel_is_enabled()) then
428 return
429 end if
430
431 call accel_detach_buffer(this%potential_accel)
432 call accel_detach_buffer(this%impotential_accel)
433
434 if (.not. allocated(this%potential)) then
436 return
437 end if
438
439 complex_potential = allocated(this%Impotential)
440
441 call accel_create_buffer(this%potential_accel, accel_mem_read_only, type_float, &
442 accel_padded_size(mesh%np)*this%nspin)
443 if (complex_potential) then
444 call accel_create_buffer(this%impotential_accel, accel_mem_read_only, type_float, &
445 accel_padded_size(mesh%np)*this%nspin)
446 end if
447
448 call this%accel_copy_pot(mesh)
449
452
453 ! ----------------------------------------------------------------------------------
455 !
456 logical pure function hamiltonian_elec_base_has_magnetic(this) result(has_magnetic)
457 class(hamiltonian_elec_base_t), intent(in) :: this
458
459 has_magnetic = allocated(this%vector_potential) &
460 .or. allocated(this%uniform_magnetic_field)
461
463
464 ! ----------------------------------------------------------------------------------
466 !
467 logical pure function hamiltonian_elec_base_has_zeeman(this) result(has_zeeman)
468 class(hamiltonian_elec_base_t), intent(in) :: this
469
470 has_zeeman = allocated(this%zeeman_pot)
471
473
474 ! ----------------------------------------------------------------------------------
476 !
477 logical pure function hamiltonian_elec_base_has_vector_potential(this) result(has_vector_potential)
478 class(hamiltonian_elec_base_t), intent(in) :: this
479
480 has_vector_potential = allocated(this%vector_potential) &
481 .or. allocated(this%uniform_vector_potential)
482
484
485 ! ---------------------------------------------------------------------------------------
486 subroutine hamiltonian_elec_base_rashba(this, mesh, der, std, psib, vpsib)
487 class(hamiltonian_elec_base_t), intent(in) :: this
488 class(mesh_t), intent(in) :: mesh
489 type(derivatives_t), intent(in) :: der
490 type(states_elec_dim_t), intent(in) :: std
491 type(wfs_elec_t), target, intent(in) :: psib
492 type(wfs_elec_t), target, intent(inout) :: vpsib
493
494 integer :: ist, idim, ip
495 complex(real64), allocatable :: psi(:, :), vpsi(:, :), grad(:, :, :)
496
498
499 if (abs(this%rashba_coupling) < m_epsilon) then
501 return
502 end if
503 assert(std%ispin == spinors)
504 assert(mesh%box%dim == 2)
505 assert(psib%type() == type_cmplx)
506 assert(vpsib%type() == type_cmplx)
507
508 safe_allocate(psi(1:mesh%np_part, 1:std%dim))
509 safe_allocate(vpsi(1:mesh%np, 1:std%dim))
510 safe_allocate(grad(1:mesh%np, 1:mesh%box%dim, 1:std%dim))
511
512 do ist = 1, psib%nst
513 call batch_get_state(psib, ist, mesh%np_part, psi)
514 call batch_get_state(vpsib, ist, mesh%np, vpsi)
515
516 do idim = 1, std%dim
517 call zderivatives_grad(der, psi(:, idim), grad(:, :, idim), ghost_update = .false., set_bc = .false.)
518 end do
519
520 if (allocated(this%vector_potential)) then
521 do ip = 1, mesh%np
522 vpsi(ip, 1) = vpsi(ip, 1) - &
523 (this%rashba_coupling) * cmplx(this%vector_potential(2, ip), this%vector_potential(1, ip), real64) * psi(ip, 2)
524 vpsi(ip, 2) = vpsi(ip, 2) - &
525 (this%rashba_coupling) * cmplx(this%vector_potential(2, ip), -this%vector_potential(1, ip), real64) * psi(ip, 1)
526 end do
527 end if
528
529 do ip = 1, mesh%np
530 vpsi(ip, 1) = vpsi(ip, 1) - &
531 this%rashba_coupling*( grad(ip, 1, 2) - m_zi*grad(ip, 2, 2))
532 vpsi(ip, 2) = vpsi(ip, 2) + &
533 this%rashba_coupling*( grad(ip, 1, 1) + m_zi*grad(ip, 2, 1))
534 end do
535
536 call batch_set_state(vpsib, ist, mesh%np, vpsi)
537 end do
538
539 safe_deallocate_a(grad)
540 safe_deallocate_a(vpsi)
541 safe_deallocate_a(psi)
542
544 end subroutine hamiltonian_elec_base_rashba
545
546#include "undef.F90"
547#include "real.F90"
548#include "hamiltonian_elec_base_inc.F90"
549
550#include "undef.F90"
551#include "complex.F90"
552#include "hamiltonian_elec_base_inc.F90"
553
555
556!! Local Variables:
557!! mode: f90
558!! coding: utf-8
559!! End:
subroutine, public accel_free_buffer(this, async)
Definition: accel.F90:1005
subroutine, public accel_detach_buffer(this)
Clear a buffer handle without freeing device memory.
Definition: accel.F90:1049
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
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:120
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:233
logical pure function hamiltonian_elec_base_has_vector_potential(this)
return .true. of the Hamiltonian contains any vector potential
subroutine hamiltonian_elec_base_accel_copy_pot(this, mesh)
copy the potential to the acceleration device buffer
subroutine, public dhamiltonian_elec_base_local_sub(potential, mesh, std, ispin, psib, vpsib, Impotential, potential_accel, impotential_accel, async)
apply a local potential to a set of states
subroutine hamiltonian_elec_base_update_magnetic_terms(this, mesh, gyromagnetic_ratio, ispin)
update the magnetic terms of the hamiltonian_elec_base_t.
integer, parameter, public term_local_external
subroutine hamiltonian_elec_base_rashba(this, mesh, der, std, psib, vpsib)
subroutine dhamiltonian_elec_base_local(this, mesh, std, ispin, psib, vpsib, async)
apply the local potential (stored in the hamiltonian) to the states
integer, parameter, public field_uniform_magnetic_field
subroutine zhamiltonian_elec_base_magnetic(this, mesh, der, std, ispin, psib, vpsib)
apply magnetic terms form the Hamiltonian to the wave functions
subroutine hamiltonian_elec_base_allocate(this, mesh, field, complex_potential)
This function ensures that the corresponding field is allocated.
integer, parameter, public field_uniform_vector_potential
integer, parameter, public field_vector_potential
subroutine hamiltonian_elec_base_clear(this, np)
This functions sets to zero all fields that are currently allocated.
logical pure function hamiltonian_elec_base_has_magnetic(this)
return .true. if the Hamiltonian contains any magnetic field
subroutine hamiltonian_elec_base_end(this)
Finalizer for hamiltonian_elec_base_t.
integer, parameter, public term_mgga
integer, parameter, public term_others
integer, parameter, public term_non_local_potential
integer, parameter, public term_rdmft_occ
integer, parameter, public term_kinetic
subroutine zhamiltonian_elec_base_local(this, mesh, std, ispin, psib, vpsib, async)
apply the local potential (stored in the hamiltonian) to the states
subroutine dhamiltonian_elec_base_magnetic(this, mesh, der, std, ispin, psib, vpsib)
apply magnetic terms form the Hamiltonian to the wave functions
logical pure function hamiltonian_elec_base_has_zeeman(this)
return .true. of the Hamiltonian contains a zeeman term
integer, parameter, public field_magnetic_field
subroutine, public zhamiltonian_elec_base_local_sub(potential, mesh, std, ispin, psib, vpsib, Impotential, potential_accel, impotential_accel, async)
apply a local potential to a set of states
subroutine hamiltonian_elec_base_init(this, nspin, mass, rashba_coupling)
initialize the hamiltonian_elec_base_t object
subroutine, public hamiltonian_elec_base_accel_rebuild(this, mesh)
Rebuild accelerator buffers after an intrinsic copy.
integer, parameter, public term_dft_u
integer, parameter, public term_local_potential
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
This module defines non-local operators.
Definition: ps.F90:116
This module handles spin dimensions of the states and the k-point distribution.
type(type_t), public type_float
Definition: types.F90:135
The basic Hamiltonian for electronic system.
Describes mesh distribution to nodes.
Definition: mesh.F90:187