34 use,
intrinsic :: iso_fortran_env
66 logical :: include_diamag
78 integer,
parameter,
public :: &
79 CURRENT_GRADIENT = 1, &
86 type(current_t),
intent(out) :: this
87 type(namespace_t),
intent(in) :: namespace
126 call parse_variable(namespace,
'CalculateDiamagneticCurrent', .false., this%include_diamag)
134 type(states_elec_t),
intent(inout) :: st
135 type(derivatives_t),
intent(inout) :: der
136 integer,
intent(in) :: ik
137 integer,
intent(in) :: ib
138 type(wfs_elec_t),
intent(in) :: psib
139 class(wfs_elec_t),
intent(in) :: gpsib(:)
141 integer :: ist, idir, ii, ip, idim, bsize, gsize
142 complex(real64),
allocatable :: psi(:, :), gpsi(:, :)
143 real(real64),
allocatable :: current_tmp(:, :)
144 complex(real64) :: c_tmp
146 real(real64),
allocatable :: weight(:)
147 type(accel_mem_t) :: buff_weight, buff_current
148 type(accel_kernel_t),
save :: kernel
150 safe_allocate(psi(1:der%mesh%np_part, 1:st%d%dim))
151 safe_allocate(gpsi(1:der%mesh%np_part, 1:st%d%dim))
158 ww = st%kweights(ik)*st%occ(ist, ik)
161 do idim = 1, st%d%dim
162 ii = st%group%psib(ib, ik)%inv_index((/ist, idim/))
167 if (st%d%ispin /=
spinors)
then
169 do ip = 1, der%mesh%np
170 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) + ww*aimag(conjg(psi(ip, 1))*gpsi(ip, 1))
175 do ip = 1, der%mesh%np
176 st%current_para(ip, idir, 1) = st%current_para(ip, idir, 1) + ww*aimag(conjg(psi(ip, 1))*gpsi(ip, 1))
177 st%current_para(ip, idir, 2) = st%current_para(ip, idir, 2) + ww*aimag(conjg(psi(ip, 2))*gpsi(ip, 2))
178 c_tmp = -
m_half*(conjg(psi(ip, 2))*gpsi(ip, 1) - psi(ip, 1)*conjg(gpsi(ip, 2)))
179 st%current_para(ip, idir, 3) = st%current_para(ip, idir, 3) - ww*aimag(c_tmp)
180 st%current_para(ip, idir, 4) = st%current_para(ip, idir, 4) + ww*real(c_tmp, real64)
192 safe_allocate(weight(1:psib%nst))
194 weight(ist) = st%kweights(ik)*st%occ(psib%ist(ist), ik)
222 safe_allocate(current_tmp(1:der%mesh%np, 1:der%dim))
226 call lalg_axpy(der%mesh%np, der%dim,
m_one, current_tmp, st%current_kpt(:,:,ik))
228 safe_deallocate_a(current_tmp)
233 safe_deallocate_a(weight)
237 assert(psib%is_packed() .eqv. gpsib(1)%is_packed())
242 ww = st%kweights(ik)*st%occ(ist, ik)
245 if (psib%is_packed())
then
248 do ip = 1, der%mesh%np
249 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) &
250 + ww*aimag(conjg(psib%zff_pack(ii, ip))*gpsib(idir)%zff_pack(ii, ip))
257 do ip = 1, der%mesh%np
258 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) &
259 + ww*aimag(conjg(psib%zff(ip, 1, ii))*gpsib(idir)%zff(ip, 1, ii))
269 safe_deallocate_a(psi)
270 safe_deallocate_a(gpsi)
280 type(
grid_t),
intent(inout) :: gr
282 class(
space_t),
intent(in) :: space
290 st%current = st%current_para
292 if (this%include_diamag)
then
294 st%current = st%current + st%current_dia
310 integer :: ispin, idir, ip
317 if(
allocated(hm%hm_base%vector_potential))
then
318 do ispin = 1, st%d%nspin
321 do ip = 1, der%mesh%np
323 st%current_dia(ip, idir, ispin) = st%current_dia(ip, idir, ispin) + &
324 st%rho(ip, ispin)*hm%hm_base%vector_potential(idir, ip)
344 real(real64),
allocatable :: magnetization_density(:, :), curl_mag(:, :)
352 safe_allocate(magnetization_density(1:der%mesh%np_part, 1:der%dim))
353 safe_allocate(curl_mag(1:der%mesh%np_part, 1:der%dim))
358 call lalg_axpy(der%mesh%np, der%dim,
m_half, curl_mag, st%current_mag(:, :, 1))
360 safe_deallocate_a(magnetization_density)
361 safe_deallocate_a(curl_mag)
375 type(
grid_t),
intent(inout) :: gr
377 class(
space_t),
intent(in) :: space
380 integer :: ik, ist, idir, idim, ip, ib, ii, ispin
381 complex(real64),
allocatable :: gpsi(:, :, :), psi(:, :), hpsi(:, :), rhpsi(:, :), rpsi(:, :), hrpsi(:, :)
382 type(
wfs_elec_t) :: hpsib, rhpsib, rpsib, hrpsib, epsib
385 complex(real64) :: c_tmp
391 assert(all(ubound(st%current_para) == (/gr%np_part, gr%der%dim, st%d%nspin/)))
392 assert(all(ubound(st%current_kpt) == (/gr%np, gr%der%dim, st%d%kpt%end/)))
393 assert(all(lbound(st%current_kpt) == (/1, 1, st%d%kpt%start/)))
395 safe_allocate(psi(1:gr%np_part, 1:st%d%dim))
396 safe_allocate(gpsi(1:gr%np, 1:gr%der%dim, 1:st%d%dim))
397 safe_allocate(hpsi(1:gr%np_part, 1:st%d%dim))
398 safe_allocate(rhpsi(1:gr%np_part, 1:st%d%dim))
399 safe_allocate(rpsi(1:gr%np_part, 1:st%d%dim))
400 safe_allocate(hrpsi(1:gr%np_part, 1:st%d%dim))
401 safe_allocate_type_array(
wfs_elec_t, commpsib, (1:gr%der%dim))
404 do ik = st%d%kpt%start, st%d%kpt%end
405 do idir = 1, gr%der%dim
408 st%current_kpt(ip, idir, ik) =
m_zero
413 do ispin = 1, st%d%nspin
414 do idir = 1, gr%der%dim
417 st%current_para(ip, idir, ispin) =
m_zero
424 select case (this%method)
428 if (.not. hm%exxop%useACE)
then
434 do ik = st%d%kpt%start, st%d%kpt%end
435 ispin = st%d%get_spin_index(ik)
436 do ib = st%group%block_start, st%group%block_end
438 call st%group%psib(ib, ik)%do_pack(copy = .
true.)
440 call st%group%psib(ib, ik)%copy_to(hpsib)
441 call st%group%psib(ib, ik)%copy_to(rhpsib)
442 call st%group%psib(ib, ik)%copy_to(rpsib)
443 call st%group%psib(ib, ik)%copy_to(hrpsib)
448 do idir = 1, gr%der%dim
450 call batch_mul(gr%np, gr%x_t(:, idir), hpsib, rhpsib)
451 call batch_mul(gr%np_part, gr%x_t(:, idir), st%group%psib(ib, ik), rpsib)
456 ww = st%kweights(ik)*st%occ(ist, ik)
459 do idim = 1, st%d%dim
460 ii = st%group%psib(ib, ik)%inv_index((/ist, idim/))
466 if (st%d%ispin /=
spinors)
then
469 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) &
470 - ww*aimag(conjg(psi(ip, 1))*hrpsi(ip, 1) - conjg(psi(ip, 1))*rhpsi(ip, 1))
476 st%current_para(ip, idir, 1) = st%current_para(ip, idir, 1) + &
477 ww*aimag(conjg(psi(ip, 1))*hrpsi(ip, 1) - conjg(psi(ip, 1))*rhpsi(ip, 1))
478 st%current_para(ip, idir, 2) = st%current_para(ip, idir, 2) + &
479 ww*aimag(conjg(psi(ip, 2))*hrpsi(ip, 2) - conjg(psi(ip, 2))*rhpsi(ip, 2))
480 c_tmp = -
m_half*
m_zi*(conjg(psi(ip, 2))*hrpsi(ip, 1) - conjg(psi(ip, 2))*rhpsi(ip, 1) &
481 -psi(ip, 1)*conjg(hrpsi(ip, 2)) - psi(ip, 1)*conjg(rhpsi(ip, 2)))
482 st%current_para(ip, idir, 3) = st%current_para(ip, idir, 3) + ww*real(c_tmp, real64)
483 st%current_para(ip, idir, 4) = st%current_para(ip, idir, 4) + ww*aimag(c_tmp)
492 call st%group%psib(ib, ik)%do_unpack(copy = .false.)
503 if (.not. hm%exxop%useACE)
then
514 .and. hm%theory_level /=
rdmft &
515 .and. hm%vnl%apply_projector_matrices)
then
519 do ik = st%d%kpt%start, st%d%kpt%end
520 ispin = st%d%get_spin_index(ik)
521 do ib = st%group%block_start, st%group%block_end
524 call hm%phase%copy_and_set_phase(gr, st%d%kpt, st%group%psib(ib, ik), epsib)
529 call hm%vnl%zposition_commutator(gr, st%d, gr%der%boundaries%spiral, epsib, commpsib, async=.
true.)
531 call zlda_u_commute_r(hm%lda_u, gr, space, st%d, namespace, epsib, commpsib)
535 do idir = 1, gr%der%dim
536 call commpsib(idir)%end()
549 do ik = st%d%kpt%start, st%d%kpt%end
550 ispin = st%d%get_spin_index(ik)
551 do ist = st%st_start, st%st_end
553 ww = st%kweights(ik)*st%occ(ist, ik)
558 if (hm%phase%is_allocated())
then
559 call hm%phase%apply_to_single(psi, gr%np, st%d%dim, ik, .false.)
562 do idim = 1, st%d%dim
563 call boundaries_set(gr%der%boundaries, gr, psi(:, idim), phase_correction = hm%phase%phase_corr(:, ik), &
564 buff_phase_corr = hm%phase%buff_phase_corr, offset=int((ik-st%d%kpt%start)*(gr%np_part-gr%np)))
567 do idim = 1, st%d%dim
572 do idim = 1, st%d%dim
579 call hm%ks_pot%zcurrent_mass_renormalization(gpsi, gr%der%dim, st%d%dim, ispin)
583 gr%der%boundaries, ik, psi, gpsi)
585 if (hm%scissor%apply)
then
590 psi, gpsi, hm%phase%is_allocated())
596 if (st%d%ispin /=
spinors)
then
597 do idir = 1, gr%der%dim
600 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) + &
601 ww*aimag(conjg(psi(ip, 1))*gpsi(ip, idir, 1))
606 do idir = 1, gr%der%dim
609 st%current_para(ip, idir, 1) = st%current_para(ip, idir, 1) + &
610 ww*aimag(conjg(psi(ip, 1))*gpsi(ip, idir, 1))
611 st%current_para(ip, idir, 2) = st%current_para(ip, idir, 2) + &
612 ww*aimag(conjg(psi(ip, 2))*gpsi(ip, idir, 2))
613 c_tmp = -
m_half*
m_zi*(conjg(psi(ip, 2))*gpsi(ip, idir, 1) - psi(ip, 1)*conjg(gpsi(ip, idir, 2)))
614 st%current_para(ip, idir, 3) = st%current_para(ip, idir, 3) + ww*real(c_tmp, real64)
615 st%current_para(ip, idir, 4) = st%current_para(ip, idir, 4) + ww*aimag(c_tmp)
632 if (st%d%ispin /=
spinors)
then
634 do ik = st%d%kpt%start, st%d%kpt%end
635 ispin = st%d%get_spin_index(ik)
636 call lalg_axpy(gr%np, gr%der%dim,
m_one, st%current_kpt(:, :, ik), st%current_para(:, :, ispin))
640 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
641 call comm_allreduce(st%st_kpt_mpi_grp, st%current_para, dim = (/gr%np, gr%der%dim, st%d%nspin/))
644 if (st%symmetrize_density)
then
645 do ispin = 1, st%d%nspin
650 safe_deallocate_a(gpsi)
651 safe_deallocate_a(psi)
652 safe_deallocate_a(hpsi)
653 safe_deallocate_a(rhpsi)
654 safe_deallocate_a(rpsi)
655 safe_deallocate_a(hrpsi)
656 safe_deallocate_a(commpsib)
673 complex(real64),
intent(in) :: psi_i(:,:)
674 complex(real64),
intent(in) :: psi_j(:,:)
675 integer,
intent(in) :: ik
676 complex(real64),
intent(out) :: cmel(:,:)
678 integer :: idir, idim, ispin
679 complex(real64),
allocatable :: gpsi_j(:, :, :), ppsi_j(:,:), gpsi_i(:, :, :), ppsi_i(:,:)
683 safe_allocate(gpsi_i(1:der%mesh%np, 1:der%dim, 1:hm%d%dim))
684 safe_allocate(ppsi_i(1:der%mesh%np_part,1:hm%d%dim))
685 safe_allocate(gpsi_j(1:der%mesh%np, 1:der%dim, 1:hm%d%dim))
686 safe_allocate(ppsi_j(1:der%mesh%np_part,1:hm%d%dim))
690 ispin = hm%d%get_spin_index(ik)
692 ppsi_i(1:der%mesh%np,:) = psi_i(1:der%mesh%np,:)
694 ppsi_j(1:der%mesh%np,:) = psi_j(1:der%mesh%np,:)
697 do idim = 1, hm%d%dim
702 if (hm%phase%is_allocated())
then
704 call hm%phase%apply_to_single(ppsi_i, der%mesh%np_part, hm%d%dim, ik, .false.)
705 call hm%phase%apply_to_single(ppsi_j, der%mesh%np_part, hm%d%dim, ik, .false.)
708 do idim = 1, hm%d%dim
709 call zderivatives_grad(der, ppsi_i(:, idim), gpsi_i(:, :, idim), set_bc = .false.)
710 call zderivatives_grad(der, ppsi_j(:, idim), gpsi_j(:, :, idim), set_bc = .false.)
715 call hm%ks_pot%zcurrent_mass_renormalization(gpsi_i, der%dim, hm%d%dim, ispin)
716 call hm%ks_pot%zcurrent_mass_renormalization(gpsi_j, der%dim, hm%d%dim, ispin)
720 der%boundaries, ik, ppsi_i, gpsi_i)
722 der%boundaries, ik, ppsi_j, gpsi_j)
724 if (hm%scissor%apply)
then
732 do idim = 1, hm%d%dim
734 cmel(idir,ispin) =
m_zi *
zmf_dotp(der%mesh, psi_i(:, idim), gpsi_j(:, idir,idim), reduce = .false.)
735 cmel(idir,ispin) = cmel(idir,ispin) -
m_zi *
zmf_dotp(der%mesh, gpsi_i(:, idir, idim), psi_j(:, idim), reduce = .false.)
740 call der%mesh%allreduce(cmel)
742 safe_deallocate_a(gpsi_i)
743 safe_deallocate_a(ppsi_i)
744 safe_deallocate_a(gpsi_j)
745 safe_deallocate_a(ppsi_j)
753 class(
space_t),
intent(in) :: space
757 real(real64),
intent(out) :: current(:, :, :)
759 integer :: ik, ist, idir, idim, ip, ispin, ndim
760 complex(real64),
allocatable :: gpsi(:, :, :), psi(:, :), g2psi(:, :, :, :)
761 complex(real64) :: tmp
765 assert(space%is_periodic())
766 assert(st%d%dim == 1)
770 safe_allocate(psi(1:der%mesh%np_part, 1:st%d%dim))
771 safe_allocate(gpsi(1:der%mesh%np_part, 1:ndim, 1:st%d%dim))
772 safe_allocate(g2psi(1:der%mesh%np, 1:ndim, 1:ndim, 1:st%d%dim))
774 do ip = 1, der%mesh%np
775 current(ip, 1:ndim, 1:st%d%nspin) = st%current(ip, 1:ndim, 1:st%d%nspin)*hm%ep%vpsl(ip)
779 do ik = st%d%kpt%start, st%d%kpt%end
780 ispin = st%d%get_spin_index(ik)
781 do ist = st%st_start, st%st_end
783 if (abs(st%kweights(ik)*st%occ(ist, ik)) <=
m_epsilon) cycle
786 do idim = 1, st%d%dim
790 if (hm%phase%is_allocated())
then
791 call hm%phase%apply_to_single(psi, der%mesh%np_part, st%d%dim, ik, conjugate = .false.)
794 do idim = 1, st%d%dim
798 if (hm%phase%is_allocated())
then
799 call hm%phase%apply_to_single(gpsi(:, idir, :), der%mesh%np, st%d%dim, ik, conjugate = .
true.)
802 do idim = 1, st%d%dim
806 if (hm%phase%is_allocated())
then
807 call hm%phase%apply_to_single(gpsi(:, idir, :), der%mesh%np_part, st%d%dim, ik, conjugate = .false.)
810 do idim = 1, st%d%dim
811 call zderivatives_grad(der, gpsi(:, idir, idim), g2psi(:, :, idir, idim), set_bc = .false.)
815 do ip = 1, der%mesh%np
817 tmp = sum(conjg(g2psi(ip, 1:ndim, idir, idim))*gpsi(ip, 1:ndim, idim)) - &
818 sum(conjg(gpsi(ip, 1:ndim, idim))*g2psi(ip, 1:ndim, idir, idim))
819 tmp = tmp - conjg(gpsi(ip, idir, idim))*sum(g2psi(ip, 1:ndim, 1:ndim, idim)) + &
820 sum(conjg(g2psi(ip, 1:ndim, 1:ndim, idim)))*gpsi(ip, idir, idim)
821 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_mel(der, hm, psi_i, psi_j, ik, cmel)
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
complex(real64), parameter, public m_z0
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