34 use,
intrinsic :: iso_fortran_env
66 logical :: include_diamag
77 integer,
parameter,
public :: &
78 CURRENT_GRADIENT = 1, &
85 type(current_t),
intent(out) :: this
86 type(namespace_t),
intent(in) :: namespace
125 call parse_variable(namespace,
'CalculateDiamagneticCurrent', .false., this%include_diamag)
133 type(states_elec_t),
intent(inout) :: st
134 type(derivatives_t),
intent(inout) :: der
135 integer,
intent(in) :: ik
136 integer,
intent(in) :: ib
137 type(wfs_elec_t),
intent(in) :: psib
138 class(wfs_elec_t),
intent(in) :: gpsib(:)
140 integer :: ist, idir, ii, ip, idim, bsize, gsize
141 complex(real64),
allocatable :: psi(:, :), gpsi(:, :)
142 real(real64),
allocatable :: current_tmp(:, :)
143 complex(real64) :: c_tmp
145 real(real64),
allocatable :: weight(:)
146 type(accel_mem_t) :: buff_weight, buff_current
147 type(accel_kernel_t),
save :: kernel
149 safe_allocate(psi(1:der%mesh%np_part, 1:st%d%dim))
150 safe_allocate(gpsi(1:der%mesh%np_part, 1:st%d%dim))
157 ww = st%kweights(ik)*st%occ(ist, ik)
160 do idim = 1, st%d%dim
161 ii = st%group%psib(ib, ik)%inv_index((/ist, idim/))
166 if (st%d%ispin /=
spinors)
then
168 do ip = 1, der%mesh%np
169 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) + ww*aimag(conjg(psi(ip, 1))*gpsi(ip, 1))
174 do ip = 1, der%mesh%np
175 st%current_para(ip, idir, 1) = st%current_para(ip, idir, 1) + ww*aimag(conjg(psi(ip, 1))*gpsi(ip, 1))
176 st%current_para(ip, idir, 2) = st%current_para(ip, idir, 2) + ww*aimag(conjg(psi(ip, 2))*gpsi(ip, 2))
177 c_tmp = -
m_half*(conjg(psi(ip, 2))*gpsi(ip, 1) - psi(ip, 1)*conjg(gpsi(ip, 2)))
178 st%current_para(ip, idir, 3) = st%current_para(ip, idir, 3) - ww*aimag(c_tmp)
179 st%current_para(ip, idir, 4) = st%current_para(ip, idir, 4) + ww*real(c_tmp, real64)
191 safe_allocate(weight(1:psib%nst))
193 weight(ist) = st%kweights(ik)*st%occ(psib%ist(ist), ik)
221 safe_allocate(current_tmp(1:der%mesh%np, 1:der%dim))
225 call lalg_axpy(der%mesh%np, der%dim,
m_one, current_tmp, st%current_kpt(:,:,ik))
227 safe_deallocate_a(current_tmp)
232 safe_deallocate_a(weight)
236 assert(psib%is_packed() .eqv. gpsib(1)%is_packed())
241 ww = st%kweights(ik)*st%occ(ist, ik)
244 if (psib%is_packed())
then
247 do ip = 1, der%mesh%np
248 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) &
249 + ww*aimag(conjg(psib%zff_pack(ii, ip))*gpsib(idir)%zff_pack(ii, ip))
256 do ip = 1, der%mesh%np
257 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) &
258 + ww*aimag(conjg(psib%zff(ip, 1, ii))*gpsib(idir)%zff(ip, 1, ii))
268 safe_deallocate_a(psi)
269 safe_deallocate_a(gpsi)
279 type(
grid_t),
intent(inout) :: gr
281 class(
space_t),
intent(in) :: space
289 st%current = st%current_para
291 if (this%include_diamag)
then
293 st%current = st%current + st%current_dia
309 integer :: ispin, idir, ip
316 if(
allocated(hm%hm_base%vector_potential))
then
317 do ispin = 1, st%d%nspin
320 do ip = 1, der%mesh%np
322 st%current_dia(ip, idir, ispin) = st%current_dia(ip, idir, ispin) + &
323 st%rho(ip, ispin)*hm%hm_base%vector_potential(idir, ip)
343 real(real64),
allocatable :: magnetization_density(:, :), curl_mag(:, :)
351 safe_allocate(magnetization_density(1:der%mesh%np_part, 1:der%dim))
352 safe_allocate(curl_mag(1:der%mesh%np_part, 1:der%dim))
357 call lalg_axpy(der%mesh%np, der%dim,
m_half, curl_mag, st%current_mag(:, :, 1))
359 safe_deallocate_a(magnetization_density)
360 safe_deallocate_a(curl_mag)
374 type(
grid_t),
intent(inout) :: gr
376 class(
space_t),
intent(in) :: space
379 integer :: ik, ist, idir, idim, ip, ib, ii, ispin
380 complex(real64),
allocatable :: gpsi(:, :, :), psi(:, :), hpsi(:, :), rhpsi(:, :), rpsi(:, :), hrpsi(:, :)
381 type(
wfs_elec_t) :: hpsib, rhpsib, rpsib, hrpsib, epsib
384 complex(real64) :: c_tmp
390 assert(all(ubound(st%current_para) == (/gr%np_part, gr%der%dim, st%d%nspin/)))
391 assert(all(ubound(st%current_kpt) == (/gr%np, gr%der%dim, st%d%kpt%end/)))
392 assert(all(lbound(st%current_kpt) == (/1, 1, st%d%kpt%start/)))
394 safe_allocate(psi(1:gr%np_part, 1:st%d%dim))
395 safe_allocate(gpsi(1:gr%np, 1:gr%der%dim, 1:st%d%dim))
396 safe_allocate(hpsi(1:gr%np_part, 1:st%d%dim))
397 safe_allocate(rhpsi(1:gr%np_part, 1:st%d%dim))
398 safe_allocate(rpsi(1:gr%np_part, 1:st%d%dim))
399 safe_allocate(hrpsi(1:gr%np_part, 1:st%d%dim))
400 safe_allocate_type_array(
wfs_elec_t, commpsib, (1:gr%der%dim))
403 do ik = st%d%kpt%start, st%d%kpt%end
404 do idir = 1, gr%der%dim
407 st%current_kpt(ip, idir, ik) =
m_zero
412 do ispin = 1, st%d%nspin
413 do idir = 1, gr%der%dim
416 st%current_para(ip, idir, ispin) =
m_zero
423 select case (this%method)
427 if (.not. hm%exxop%useACE)
then
433 do ik = st%d%kpt%start, st%d%kpt%end
434 ispin = st%d%get_spin_index(ik)
435 do ib = st%group%block_start, st%group%block_end
437 call st%group%psib(ib, ik)%do_pack(copy = .
true.)
439 call st%group%psib(ib, ik)%copy_to(hpsib)
440 call st%group%psib(ib, ik)%copy_to(rhpsib)
441 call st%group%psib(ib, ik)%copy_to(rpsib)
442 call st%group%psib(ib, ik)%copy_to(hrpsib)
447 do idir = 1, gr%der%dim
449 call batch_mul(gr%np, gr%x_t(:, idir), hpsib, rhpsib)
450 call batch_mul(gr%np_part, gr%x_t(:, idir), st%group%psib(ib, ik), rpsib)
455 ww = st%kweights(ik)*st%occ(ist, ik)
458 do idim = 1, st%d%dim
459 ii = st%group%psib(ib, ik)%inv_index((/ist, idim/))
465 if (st%d%ispin /=
spinors)
then
468 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) &
469 - ww*aimag(conjg(psi(ip, 1))*hrpsi(ip, 1) - conjg(psi(ip, 1))*rhpsi(ip, 1))
475 st%current_para(ip, idir, 1) = st%current_para(ip, idir, 1) + &
476 ww*aimag(conjg(psi(ip, 1))*hrpsi(ip, 1) - conjg(psi(ip, 1))*rhpsi(ip, 1))
477 st%current_para(ip, idir, 2) = st%current_para(ip, idir, 2) + &
478 ww*aimag(conjg(psi(ip, 2))*hrpsi(ip, 2) - conjg(psi(ip, 2))*rhpsi(ip, 2))
479 c_tmp = -
m_half*
m_zi*(conjg(psi(ip, 2))*hrpsi(ip, 1) - conjg(psi(ip, 2))*rhpsi(ip, 1) &
480 -psi(ip, 1)*conjg(hrpsi(ip, 2)) - psi(ip, 1)*conjg(rhpsi(ip, 2)))
481 st%current_para(ip, idir, 3) = st%current_para(ip, idir, 3) + ww*real(c_tmp, real64)
482 st%current_para(ip, idir, 4) = st%current_para(ip, idir, 4) + ww*aimag(c_tmp)
491 call st%group%psib(ib, ik)%do_unpack(copy = .false.)
502 if (.not. hm%exxop%useACE)
then
513 .and. hm%theory_level /=
rdmft &
514 .and. hm%vnl%apply_projector_matrices)
then
518 do ik = st%d%kpt%start, st%d%kpt%end
519 ispin = st%d%get_spin_index(ik)
520 do ib = st%group%block_start, st%group%block_end
523 call hm%phase%copy_and_set_phase(gr, st%d%kpt, st%group%psib(ib, ik), epsib)
528 call hm%vnl%zposition_commutator(gr, st%d, gr%der%boundaries%spiral, epsib, commpsib, async=.
true.)
530 call zlda_u_commute_r(hm%lda_u, gr, space, st%d, namespace, epsib, commpsib)
534 do idir = 1, gr%der%dim
535 call commpsib(idir)%end()
548 do ik = st%d%kpt%start, st%d%kpt%end
549 ispin = st%d%get_spin_index(ik)
550 do ist = st%st_start, st%st_end
552 ww = st%kweights(ik)*st%occ(ist, ik)
557 if (hm%phase%is_allocated())
then
558 call hm%phase%apply_to_single(psi, gr%np, st%d%dim, ik, .false.)
561 do idim = 1, st%d%dim
562 call boundaries_set(gr%der%boundaries, gr, psi(:, idim), phase_correction = hm%phase%phase_corr(:, ik), &
563 buff_phase_corr = hm%phase%buff_phase_corr, offset=int((ik-st%d%kpt%start)*(gr%np_part-gr%np)))
566 do idim = 1, st%d%dim
571 do idim = 1, st%d%dim
578 call hm%ks_pot%zcurrent_mass_renormalization(gpsi, gr%der%dim, st%d%dim, ispin)
582 gr%der%boundaries, ik, psi, gpsi)
584 if (hm%scissor%apply)
then
589 psi, gpsi, hm%phase%is_allocated())
595 if (st%d%ispin /=
spinors)
then
596 do idir = 1, gr%der%dim
599 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) + &
600 ww*aimag(conjg(psi(ip, 1))*gpsi(ip, idir, 1))
605 do idir = 1, gr%der%dim
608 st%current_para(ip, idir, 1) = st%current_para(ip, idir, 1) + &
609 ww*aimag(conjg(psi(ip, 1))*gpsi(ip, idir, 1))
610 st%current_para(ip, idir, 2) = st%current_para(ip, idir, 2) + &
611 ww*aimag(conjg(psi(ip, 2))*gpsi(ip, idir, 2))
612 c_tmp = -
m_half*
m_zi*(conjg(psi(ip, 2))*gpsi(ip, idir, 1) - psi(ip, 1)*conjg(gpsi(ip, idir, 2)))
613 st%current_para(ip, idir, 3) = st%current_para(ip, idir, 3) + ww*real(c_tmp, real64)
614 st%current_para(ip, idir, 4) = st%current_para(ip, idir, 4) + ww*aimag(c_tmp)
631 if (st%d%ispin /=
spinors)
then
633 do ik = st%d%kpt%start, st%d%kpt%end
634 ispin = st%d%get_spin_index(ik)
635 call lalg_axpy(gr%np, gr%der%dim,
m_one, st%current_kpt(:, :, ik), st%current_para(:, :, ispin))
639 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
640 call comm_allreduce(st%st_kpt_mpi_grp, st%current_para, dim = (/gr%np, gr%der%dim, st%d%nspin/))
643 if (st%symmetrize_density)
then
644 do ispin = 1, st%d%nspin
649 safe_deallocate_a(gpsi)
650 safe_deallocate_a(psi)
651 safe_deallocate_a(hpsi)
652 safe_deallocate_a(rhpsi)
653 safe_deallocate_a(rpsi)
654 safe_deallocate_a(hrpsi)
655 safe_deallocate_a(commpsib)
664 class(
space_t),
intent(in) :: space
668 real(real64),
intent(out) :: current(:, :, :)
670 integer :: ik, ist, idir, idim, ip, ispin, ndim
671 complex(real64),
allocatable :: gpsi(:, :, :), psi(:, :), g2psi(:, :, :, :)
672 complex(real64) :: tmp
676 assert(space%is_periodic())
677 assert(st%d%dim == 1)
681 safe_allocate(psi(1:der%mesh%np_part, 1:st%d%dim))
682 safe_allocate(gpsi(1:der%mesh%np_part, 1:ndim, 1:st%d%dim))
683 safe_allocate(g2psi(1:der%mesh%np, 1:ndim, 1:ndim, 1:st%d%dim))
685 do ip = 1, der%mesh%np
686 current(ip, 1:ndim, 1:st%d%nspin) = st%current(ip, 1:ndim, 1:st%d%nspin)*hm%ep%vpsl(ip)
690 do ik = st%d%kpt%start, st%d%kpt%end
691 ispin = st%d%get_spin_index(ik)
692 do ist = st%st_start, st%st_end
694 if (abs(st%kweights(ik)*st%occ(ist, ik)) <=
m_epsilon) cycle
697 do idim = 1, st%d%dim
701 if (hm%phase%is_allocated())
then
702 call hm%phase%apply_to_single(psi, der%mesh%np_part, st%d%dim, ik, conjugate = .false.)
705 do idim = 1, st%d%dim
709 if (hm%phase%is_allocated())
then
710 call hm%phase%apply_to_single(gpsi(:, idir, :), der%mesh%np, st%d%dim, ik, conjugate = .
true.)
713 do idim = 1, st%d%dim
717 if (hm%phase%is_allocated())
then
718 call hm%phase%apply_to_single(gpsi(:, idir, :), der%mesh%np_part, st%d%dim, ik, conjugate = .false.)
721 do idim = 1, st%d%dim
722 call zderivatives_grad(der, gpsi(:, idir, idim), g2psi(:, :, idir, idim), set_bc = .false.)
726 do ip = 1, der%mesh%np
728 tmp = sum(conjg(g2psi(ip, 1:ndim, idir, idim))*gpsi(ip, 1:ndim, idim)) - &
729 sum(conjg(gpsi(ip, 1:ndim, idim))*g2psi(ip, 1:ndim, idir, idim))
730 tmp = tmp - conjg(gpsi(ip, idir, idim))*sum(g2psi(ip, 1:ndim, 1:ndim, idim)) + &
731 sum(conjg(g2psi(ip, 1:ndim, 1:ndim, idim)))*gpsi(ip, idir, idim)
732 current(ip, idir, ispin) = current(ip, idir, ispin) + st%kweights(ik)*st%occ(ist, ik)*aimag(tmp)/8.0_real64
constant times a vector plus a vector
integer function, public accel_kernel_block_size(kernel)
subroutine, public accel_free_buffer(this, async)
subroutine, public accel_kernel_start_call(this, file_name, kernel_name, flags)
subroutine, public accel_finish()
integer, parameter, public accel_mem_write_only
integer, parameter, public accel_mem_read_only
This module implements batches of mesh functions.
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
This module implements common operations on batches of mesh functions.
Module implementing boundary conditions in Octopus.
subroutine current_batch_accumulate(st, der, ik, ib, psib, gpsib)
subroutine current_calculate_dia_non_unif_vec_pot(der, hm, st)
Compute diamagnetic current density from non-uniform vector potential (the part coming from the unifo...
subroutine, public current_calculate_mag(der, st)
Compute magnetization current Note: due to the the numerical curl, the magnetization current could de...
integer, parameter, public current_hamiltonian
integer, parameter, public current_gradient_corr
subroutine, public current_heat_calculate(space, der, hm, st, current)
subroutine current_calculate_para(this, namespace, gr, hm, space, st)
Compute paramagnetic current density (including full diamagnetic term if method = Hamiltonian us used...
subroutine, public current_calculate(this, namespace, gr, hm, space, st)
Compute total electronic current density.
subroutine, public current_init(this, namespace)
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_curl(der, ff, op_ff, ghost_update, set_bc)
apply the curl operator to a vector of mesh functions
subroutine, public zderivatives_batch_grad(der, ffb, opffb, ghost_update, set_bc, to_cartesian, factor)
apply the gradient to a batch of mesh functions
subroutine, public zderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
subroutine, public zexchange_operator_commute_r(this, namespace, mesh, st_d, ik, psi, gpsi)
real(real64), parameter, public m_zero
integer, parameter, public rdmft
integer, parameter, public hartree_fock
integer, parameter, public generalized_kohn_sham_dft
complex(real64), parameter, public m_zi
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
subroutine, public dgrid_symmetrize_vector_field(gr, field, suppress_warning)
subroutine, public zhamiltonian_elec_apply_batch(hm, namespace, mesh, psib, hpsib, terms, set_bc)
A module to handle KS potential, without the external potential.
subroutine, public zlda_u_commute_r(this, mesh, space, d, namespace, psib, gpsib)
This routine computes [r,V_lda+u] .
subroutine, public zlda_u_commute_r_single(this, mesh, space, d, namespace, ist, ik, psi, gpsi, has_phase)
subroutine, public magnetic_density(mesh, std, rho, md)
This module is intended to contain "only mathematical" functions and procedures.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
subroutine, public zprojector_commute_r_allatoms_alldir(pj, ions, mesh, dim, bnd, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
integer pure function, public states_elec_block_max(st, ib)
return index of last state in block ib
integer pure function, public states_elec_block_min(st, ib)
return index of first state in block ib
This module provides routines for communicating states when using states parallelization.
subroutine, public states_elec_parallel_remote_access_stop(this)
stop remote memory access for states on other processors
subroutine, public states_elec_parallel_remote_access_start(this)
start remote memory access for states on other processors
type(type_t), public type_float
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
class representing derivatives
Description of the grid, containing information on derivatives, stencil, and symmetries.
The states_elec_t class contains all electronic wave functions.
batches of electronic states