Octopus
current.F90
Go to the documentation of this file.
1!! Copyright (C) 2008-2019 X. Andrade, F. Bonafe, R. Jestaedt, H. Appel
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
21module current_oct_m
22 use accel_oct_m
23 use batch_oct_m
26 use comm_oct_m
27 use debug_oct_m
31 use global_oct_m
32 use grid_oct_m
34 use, intrinsic :: iso_fortran_env
37 use lda_u_oct_m
38 use math_oct_m
41 use mesh_oct_m
45 use parser_oct_m
49 use space_oct_m
52 use types_oct_m
53 use unit_oct_m
57 use xc_oct_m
58
59 implicit none
60
61 private
62
63 type current_t
64 private
65 integer :: method
66 logical :: include_diamag
67 end type current_t
68
69
70 public :: &
71 current_t, &
77
78 integer, parameter, public :: &
79 CURRENT_GRADIENT = 1, &
82
83contains
84
85 subroutine current_init(this, namespace)
86 type(current_t), intent(out) :: this
87 type(namespace_t), intent(in) :: namespace
88
89 push_sub(current_init)
90
91 !%Variable CurrentDensity
92 !%Default gradient_corrected
93 !%Type integer
94 !%Section Hamiltonian
95 !%Description
96 !% This variable selects the method used to
97 !% calculate the current density. For the moment this variable is
98 !% for development purposes and users should not need to use
99 !% it.
100 !%Option gradient 1
101 !% The calculation of current is done using the gradient operator. (Experimental)
102 !%Option gradient_corrected 2
103 !% The calculation of current is done using the gradient operator
104 !% with additional corrections for the total current from non-local operators.
105 !%Option hamiltonian 3
106 !% The current density is obtained from the commutator of the
107 !% Hamiltonian with the position operator. (Experimental)
108 !%End
109
110 call parse_variable(namespace, 'CurrentDensity', current_gradient_corr, this%method)
111 if (.not. varinfo_valid_option('CurrentDensity', this%method)) call messages_input_error(namespace, 'CurrentDensity')
112 if (this%method /= current_gradient_corr) then
113 call messages_experimental("CurrentDensity /= gradient_corrected")
114 end if
115
116 !%Variable CalculateDiamagneticCurrent
117 !%Default no
118 !%Type logical
119 !%Section Hamiltonian
120 !%Description
121 !% This variable decides whether the current density arising from the non-uniform
122 !% vector potential, defined as:
123 !% <math> \vec{J}_{dmc}(\vec{r}, t)=-\frac{e^2}{m_e c_0} n(\vec{r}, t) \vec{A}(\vec{r},t)$ </math>
124 !% is included in the total current density.
125 !%End
126 call parse_variable(namespace, 'CalculateDiamagneticCurrent', .false., this%include_diamag)
127
128 pop_sub(current_init)
129 end subroutine current_init
130
131 ! ---------------------------------------------------------
132
133 subroutine current_batch_accumulate(st, der, ik, ib, psib, gpsib)
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(:)
140
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
145 real(real64) :: ww
146 real(real64), allocatable :: weight(:)
147 type(accel_mem_t) :: buff_weight, buff_current
148 type(accel_kernel_t), save :: kernel
149
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))
152
153 if (st%d%ispin == spinors .or. (psib%status() == batch_device_packed .and. der%dim /= 3)) then
154
155 do idir = 1, der%dim
156 do ist = states_elec_block_min(st, ib), states_elec_block_max(st, ib)
157
158 ww = st%kweights(ik)*st%occ(ist, ik)
159 if (abs(ww) <= m_epsilon) cycle
161 do idim = 1, st%d%dim
162 ii = st%group%psib(ib, ik)%inv_index((/ist, idim/))
163 call batch_get_state(psib, ii, der%mesh%np, psi(:, idim))
164 call batch_get_state(gpsib(idir), ii, der%mesh%np, gpsi(:, idim))
165 end do
166
167 if (st%d%ispin /= spinors) then
168 !$omp parallel do
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))
171 end do
172 !$omp end parallel do
173 else
174 !$omp parallel do private(c_tmp)
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)
181 end do
182 !$omp end parallel do
183 end if
184
185 end do
186 end do
187
188 else if (psib%status() == batch_device_packed) then
189
190 assert(der%dim == 3)
191
192 safe_allocate(weight(1:psib%nst))
193 do ist = 1, psib%nst
194 weight(ist) = st%kweights(ik)*st%occ(psib%ist(ist), ik)
195 end do
196
197
198 call accel_create_buffer(buff_weight, accel_mem_read_only, type_float, psib%nst)
199 call accel_write_buffer(buff_weight, psib%nst, weight, async=.true.)
200
201 call accel_create_buffer(buff_current, accel_mem_write_only, type_float, der%mesh%np*3)
202
203 call accel_kernel_start_call(kernel, 'density.cu', 'current_accumulate')
204
205 call accel_set_kernel_arg(kernel, 0, psib%nst)
206 call accel_set_kernel_arg(kernel, 1, der%mesh%np)
207 call accel_set_kernel_arg(kernel, 2, buff_weight)
208 call accel_set_kernel_arg(kernel, 3, psib%ff_device)
209 call accel_set_kernel_arg(kernel, 4, log2(int(psib%pack_size(1), int32)))
210 call accel_set_kernel_arg(kernel, 5, gpsib(1)%ff_device)
211 call accel_set_kernel_arg(kernel, 6, gpsib(2)%ff_device)
212 call accel_set_kernel_arg(kernel, 7, gpsib(3)%ff_device)
213 call accel_set_kernel_arg(kernel, 8, log2(int(gpsib(1)%pack_size(1), int32)))
214 call accel_set_kernel_arg(kernel, 9, buff_current)
215
216 ! Compute the grid size
217 bsize = accel_kernel_block_size(kernel)
218 call accel_grid_size(der%mesh%np, bsize, gsize)
219
220 call accel_kernel_run(kernel, (/gsize/), (/bsize/))
221
222 safe_allocate(current_tmp(1:der%mesh%np, 1:der%dim))
223
224 call accel_read_buffer(buff_current, der%mesh%np, der%dim, current_tmp)
225
226 call lalg_axpy(der%mesh%np, der%dim, m_one, current_tmp, st%current_kpt(:,:,ik))
227
228 safe_deallocate_a(current_tmp)
229
230 call accel_free_buffer(buff_weight)
231 call accel_free_buffer(buff_current)
232
233 safe_deallocate_a(weight)
234
235 else
236
237 assert(psib%is_packed() .eqv. gpsib(1)%is_packed())
238
239 !$omp parallel private(ip, ist, ww, idir)
240 do ii = 1, psib%nst
241 ist = states_elec_block_min(st, ib) + ii - 1
242 ww = st%kweights(ik)*st%occ(ist, ik)
243 if (abs(ww) <= m_epsilon) cycle
244
245 if (psib%is_packed()) then
246 do idir = 1, der%dim
247 !$omp do
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))
251 end do
252 !$omp end do nowait
253 end do
254 else
255 do idir = 1, der%dim
256 !$omp do
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))
260 end do
261 !$omp end do nowait
262 end do
263 end if
264 end do
265 !$omp end parallel
266
267 end if
268
269 safe_deallocate_a(psi)
270 safe_deallocate_a(gpsi)
271
272 end subroutine current_batch_accumulate
273
274
275 ! ---------------------------------------------------------
277 subroutine current_calculate(this, namespace, gr, hm, space, st)
278 type(current_t), intent(in) :: this
279 type(namespace_t), intent(in) :: namespace
280 type(grid_t), intent(inout) :: gr
281 type(hamiltonian_elec_t), intent(inout) :: hm
282 class(space_t), intent(in) :: space
283 type(states_elec_t), intent(inout) :: st
284
285
286 call profiling_in("CURRENT_TOTAL")
287 push_sub(current_calculate)
288
289 call current_calculate_para(this, namespace, gr, hm, space, st)
290 st%current = st%current_para
291
292 if (this%include_diamag) then
293 call current_calculate_dia_non_unif_vec_pot(gr%der, hm, st)
294 st%current = st%current + st%current_dia
295 end if
296
297 call profiling_out("CURRENT_TOTAL")
298 pop_sub(current_calculate)
299
300 end subroutine current_calculate
301
302 ! ---------------------------------------------------------
305 subroutine current_calculate_dia_non_unif_vec_pot(der, hm, st)
306 type(derivatives_t), intent(inout) :: der
307 type(hamiltonian_elec_t), intent(in) :: hm
308 type(states_elec_t), intent(inout) :: st
309
310 integer :: ispin, idir, ip
311
312 call profiling_in("CURRENT_DIA_NON_UNIF_A")
314
315 st%current_dia = m_zero
316
317 if(allocated(hm%hm_base%vector_potential)) then
318 do ispin = 1, st%d%nspin
319 do idir = 1, der%dim
320 !$omp parallel do
321 do ip = 1, der%mesh%np
322 ! the vector potential is assumed to be devided by c_0 already
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)
325 end do
326 !$omp end parallel do
327 end do
328 end do
329 end if
330
331 call profiling_out("CURRENT_DIA_NON_UNIF_A")
333
335
336 ! ---------------------------------------------------------
340 subroutine current_calculate_mag(der, st)
341 type(derivatives_t), intent(inout) :: der
342 type(states_elec_t), intent(inout) :: st
343
344 real(real64), allocatable :: magnetization_density(:, :), curl_mag(:, :)
345
346 call profiling_in("CURRENT_MAG")
347 push_sub(current_calculate_mag)
348
349 st%current_mag = m_zero
350
351 if (st%d%ispin /= unpolarized) then
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))
354
355 call magnetic_density(der%mesh, st%d, st%rho, magnetization_density)
356 call dderivatives_curl(der, magnetization_density, curl_mag)
357
358 call lalg_axpy(der%mesh%np, der%dim, m_half, curl_mag, st%current_mag(:, :, 1))
359
360 safe_deallocate_a(magnetization_density)
361 safe_deallocate_a(curl_mag)
362 end if
363
364 call profiling_out("CURRENT_MAG")
365 pop_sub(current_calculate_mag)
366
367 end subroutine current_calculate_mag
368
369 ! ---------------------------------------------------------
372 subroutine current_calculate_para(this, namespace, gr, hm, space, st)
373 type(current_t), intent(in) :: this
374 type(namespace_t), intent(in) :: namespace
375 type(grid_t), intent(inout) :: gr
376 type(hamiltonian_elec_t), intent(inout) :: hm
377 class(space_t), intent(in) :: space
378 type(states_elec_t), target, intent(inout) :: st
379
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
383 class(wfs_elec_t), allocatable :: commpsib(:)
384 real(real64) :: ww
385 complex(real64) :: c_tmp
386
387 call profiling_in("CURRENT_PARA")
388 push_sub(current_calculate_para)
389
390 ! spin not implemented or tested
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/)))
394
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))
402
403 !$omp parallel private(idir, ip, ispin)
404 do ik = st%d%kpt%start, st%d%kpt%end
405 do idir = 1, gr%der%dim
406 !$omp do
407 do ip = 1, gr%np
408 st%current_kpt(ip, idir, ik) = m_zero
409 end do
410 !$omp end do nowait
411 end do
412 end do
413 do ispin = 1, st%d%nspin
414 do idir = 1, gr%der%dim
415 !$omp do
416 do ip = 1, gr%np
417 st%current_para(ip, idir, ispin) = m_zero
418 end do
419 !$omp end do nowait
420 end do
421 end do
422 !$omp end parallel
423
424 select case (this%method)
425
427 if (family_is_hybrid(hm%xc)) then
428 if (.not. hm%exxop%useACE) then
429 hm%exxop%st => st
431 end if
432 end if
433
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
437
438 call st%group%psib(ib, ik)%do_pack(copy = .true.)
439
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)
444
445 call boundaries_set(gr%der%boundaries, gr, st%group%psib(ib, ik))
446 call zhamiltonian_elec_apply_batch(hm, namespace, gr, st%group%psib(ib, ik), hpsib, set_bc = .false.)
447
448 do idir = 1, gr%der%dim
449
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)
452
453 call zhamiltonian_elec_apply_batch(hm, namespace, gr, rpsib, hrpsib, set_bc = .false.)
454
455 do ist = states_elec_block_min(st, ib), states_elec_block_max(st, ib)
456 ww = st%kweights(ik)*st%occ(ist, ik)
457 if (ww <= m_epsilon) cycle
458
459 do idim = 1, st%d%dim
460 ii = st%group%psib(ib, ik)%inv_index((/ist, idim/))
461 call batch_get_state(st%group%psib(ib, ik), ii, gr%np, psi(:, idim))
462 call batch_get_state(hrpsib, ii, gr%np, hrpsi(:, idim))
463 call batch_get_state(rhpsib, ii, gr%np, rhpsi(:, idim))
464 end do
465
466 if (st%d%ispin /= spinors) then
467 !$omp parallel do
468 do ip = 1, gr%np
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))
471 end do
472 !$omp end parallel do
473 else
474 !$omp parallel do private(c_tmp)
475 do ip = 1, gr%np
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)
484 end do
485 !$omp end parallel do
486 end if
487
488 end do
489
490 end do
491
492 call st%group%psib(ib, ik)%do_unpack(copy = .false.)
493
494 call hpsib%end()
495 call rhpsib%end()
496 call rpsib%end()
497 call hrpsib%end()
498
499 end do
500 end do
501
502 if (family_is_hybrid(hm%xc)) then
503 if (.not. hm%exxop%useACE) then
505 nullify(hm%exxop%st)
506 end if
507 end if
508
509 case (current_gradient, current_gradient_corr)
510
511 if (this%method == current_gradient_corr .and. .not. family_is_mgga_with_exc(hm%xc) &
512 .and. hm%theory_level /= hartree_fock &
513 .and. hm%theory_level /= generalized_kohn_sham_dft &
514 .and. hm%theory_level /= rdmft &
515 .and. hm%vnl%apply_projector_matrices) then
516
517 ! we can use the packed version
518
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
522
523 ! copy st%group%psib(ib, ik) to epsib and set the phase
524 call hm%phase%copy_and_set_phase(gr, st%d%kpt, st%group%psib(ib, ik), epsib)
525
526 ! this now takes non-orthogonal axis into account
527 call zderivatives_batch_grad(gr%der, epsib, commpsib, set_bc=.false.)
528
529 call hm%vnl%zposition_commutator(gr, st%d, gr%der%boundaries%spiral, epsib, commpsib, async=.true.)
530
531 call zlda_u_commute_r(hm%lda_u, gr, space, st%d, namespace, epsib, commpsib)
532
533 call current_batch_accumulate(st, gr%der, ik, ib, epsib, commpsib)
534
535 do idir = 1, gr%der%dim
536 call commpsib(idir)%end()
537 end do
538
539 call epsib%end()
540
541 call accel_finish()
542 end do
543 end do
544
545 else
546
547 ! use the slow non-packed version
548
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
552
553 ww = st%kweights(ik)*st%occ(ist, ik)
554 if (abs(ww) <= m_epsilon) cycle
555
556 call states_elec_get_state(st, gr, ist, ik, psi)
557
558 if (hm%phase%is_allocated()) then
559 call hm%phase%apply_to_single(psi, gr%np, st%d%dim, ik, .false.)
560 ! apply phase correction while setting boundary -> memory needs to be
561 ! accessed only once
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)))
565 end do
566 else
567 do idim = 1, st%d%dim
568 call boundaries_set(gr%der%boundaries, gr, psi(:, idim))
569 end do
570 end if
571
572 do idim = 1, st%d%dim
573 call zderivatives_grad(gr%der, psi(:, idim), gpsi(:, :, idim), set_bc = .false.)
574 end do
575
576 if (this%method == current_gradient_corr) then
577 !A nonlocal contribution from the MGGA potential must be included
578 !This must be done first, as this is like a position-dependent mass
579 call hm%ks_pot%zcurrent_mass_renormalization(gpsi, gr%der%dim, st%d%dim, ispin)
580
581 !A nonlocal contribution from the pseudopotential must be included
582 call zprojector_commute_r_allatoms_alldir(hm%ep%proj, hm%ions, gr, st%d%dim, &
583 gr%der%boundaries, ik, psi, gpsi)
584 !A nonlocal contribution from the scissor must be included
585 if (hm%scissor%apply) then
586 call scissor_commute_r(hm%scissor, gr, ik, psi, gpsi)
587 end if
588
589 call zlda_u_commute_r_single(hm%lda_u, gr, space, st%d, namespace, ist, ik, &
590 psi, gpsi, hm%phase%is_allocated())
591
592 call zexchange_operator_commute_r(hm%exxop, namespace, gr, st%d, ik, psi, gpsi)
593
594 end if
595
596 if (st%d%ispin /= spinors) then
597 do idir = 1, gr%der%dim
598 !$omp parallel do
599 do ip = 1, gr%np
600 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) + &
601 ww*aimag(conjg(psi(ip, 1))*gpsi(ip, idir, 1))
602 end do
603 !$omp end parallel do
604 end do
605 else
606 do idir = 1, gr%der%dim
607 !$omp parallel do private(c_tmp)
608 do ip = 1, gr%np
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)
616 end do
617 !$omp end parallel do
618 end do
619 end if
620
621 end do
622 end do
623
624 end if
625
626 case default
627
628 assert(.false.)
629
630 end select
631
632 if (st%d%ispin /= spinors) then
633 !We sum the current over k-points
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))
637 end do
638 end if
639
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/))
642 end if
643
644 if (st%symmetrize_density) then
645 do ispin = 1, st%d%nspin
646 call dgrid_symmetrize_vector_field(gr, st%current_para(:, :, ispin), suppress_warning = .true.)
647 end do
648 end if
649
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)
657
658 call profiling_out("CURRENT_PARA")
660
661 end subroutine current_calculate_para
662
663
664 ! ---------------------------------------------------------
665 ! Calculate the current matrix element between two states
666 ! I_{ij}(t) = <i| J(t) |j>
667 ! This is used only in the floquet_observables utility and
668 ! is highly experimental
669
670 subroutine current_calculate_mel(der, hm, psi_i, psi_j, ik, cmel)
671 type(derivatives_t), intent(inout) :: der
672 type(hamiltonian_elec_t), intent(in) :: hm
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(:,:) ! the current vector cmel(1:der%dim, 1:st%d%nspin)
677
678 integer :: idir, idim, ispin
679 complex(real64), allocatable :: gpsi_j(:, :, :), ppsi_j(:,:), gpsi_i(:, :, :), ppsi_i(:,:)
680
681 push_sub(current_calculate_mel)
682
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))
687
688 cmel = m_z0
689
690 ispin = hm%d%get_spin_index(ik)
691 ppsi_i(:,:) = m_z0
692 ppsi_i(1:der%mesh%np,:) = psi_i(1:der%mesh%np,:)
693 ppsi_j(:,:) = m_z0
694 ppsi_j(1:der%mesh%np,:) = psi_j(1:der%mesh%np,:)
695
696
697 do idim = 1, hm%d%dim
698 call boundaries_set(der%boundaries, der%mesh, ppsi_i(:, idim))
699 call boundaries_set(der%boundaries, der%mesh, ppsi_j(:, idim))
700 end do
701
702 if (hm%phase%is_allocated()) then
703 ! Apply the phase that contains both the k-point and vector-potential terms.
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.)
706 end if
707
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.)
711 end do
712
713 !A nonlocal contribution from the MGGA potential must be included
714 !This must be done first, as this is like a position-dependent mass
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)
717
718 !A nonlocal contribution from the pseudopotential must be included
719 call zprojector_commute_r_allatoms_alldir(hm%ep%proj, hm%ions, der%mesh, hm%d%dim, &
720 der%boundaries, ik, ppsi_i, gpsi_i)
721 call zprojector_commute_r_allatoms_alldir(hm%ep%proj, hm%ions, der%mesh, hm%d%dim, &
722 der%boundaries, ik, ppsi_j, gpsi_j)
723 !A nonlocal contribution from the scissor must be included
724 if (hm%scissor%apply) then
725 call scissor_commute_r(hm%scissor, der%mesh, ik, ppsi_i, gpsi_i)
726 call scissor_commute_r(hm%scissor, der%mesh, ik, ppsi_j, gpsi_j)
727 end if
728
729
730 do idir = 1, der%dim
731
732 do idim = 1, hm%d%dim
733
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.)
736
737 end do
738 end do
739
740 call der%mesh%allreduce(cmel)
741
742 safe_deallocate_a(gpsi_i)
743 safe_deallocate_a(ppsi_i)
744 safe_deallocate_a(gpsi_j)
745 safe_deallocate_a(ppsi_j)
746
747 pop_sub(current_calculate_mel)
748
749 end subroutine current_calculate_mel
750
751 ! ---------------------------------------------------------
752 subroutine current_heat_calculate(space, der, hm, st, current)
753 class(space_t), intent(in) :: space
754 type(derivatives_t), intent(in) :: der
755 type(hamiltonian_elec_t), intent(in) :: hm
756 type(states_elec_t), intent(in) :: st
757 real(real64), intent(out) :: current(:, :, :)
758
759 integer :: ik, ist, idir, idim, ip, ispin, ndim
760 complex(real64), allocatable :: gpsi(:, :, :), psi(:, :), g2psi(:, :, :, :)
761 complex(real64) :: tmp
762
763 push_sub(current_heat_calculate)
764
765 assert(space%is_periodic())
766 assert(st%d%dim == 1)
767
768 ndim = space%dim
769
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))
773
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)
776 end do
777
778
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
782
783 if (abs(st%kweights(ik)*st%occ(ist, ik)) <= m_epsilon) cycle
784
785 call states_elec_get_state(st, der%mesh, ist, ik, psi)
786 do idim = 1, st%d%dim
787 call boundaries_set(der%boundaries, der%mesh, psi(:, idim))
788 end do
789
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.)
792 end if
793
794 do idim = 1, st%d%dim
795 call zderivatives_grad(der, psi(:, idim), gpsi(:, :, idim), set_bc = .false.)
796 end do
797 do idir = 1, ndim
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.)
800 end if
801
802 do idim = 1, st%d%dim
803 call boundaries_set(der%boundaries, der%mesh, gpsi(:,idir, idim))
804 end do
805
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.)
808 end if
809
810 do idim = 1, st%d%dim
811 call zderivatives_grad(der, gpsi(:, idir, idim), g2psi(:, :, idir, idim), set_bc = .false.)
812 end do
813 end do
814 idim = 1
815 do ip = 1, der%mesh%np
816 do idir = 1, ndim
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
822 end do
823 end do
824 end do
825 end do
826
827
829
830 end subroutine current_heat_calculate
831
832end module current_oct_m
833
834!! Local Variables:
835!! mode: f90
836!! coding: utf-8
837!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:173
integer function, public accel_kernel_block_size(kernel)
Definition: accel.F90:1167
subroutine, public accel_free_buffer(this, async)
Definition: accel.F90:1004
subroutine, public accel_kernel_start_call(this, file_name, kernel_name, flags)
Definition: accel.F90:1392
subroutine, public accel_finish()
Definition: accel.F90:1077
integer, parameter, public accel_mem_write_only
Definition: accel.F90:184
integer, parameter, public accel_mem_read_only
Definition: accel.F90:184
This module implements batches of mesh functions.
Definition: batch.F90:135
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
Definition: batch.F90:286
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:124
subroutine current_batch_accumulate(st, der, ik, ib, psib, gpsib)
Definition: current.F90:229
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...
Definition: current.F90:401
subroutine, public current_calculate_mag(der, st)
Compute magnetization current Note: due to the the numerical curl, the magnetization current could de...
Definition: current.F90:436
integer, parameter, public current_hamiltonian
Definition: current.F90:173
integer, parameter, public current_gradient_corr
Definition: current.F90:173
subroutine, public current_heat_calculate(space, der, hm, st, current)
Definition: current.F90:848
subroutine current_calculate_para(this, namespace, gr, hm, space, st)
Compute paramagnetic current density (including full diamagnetic term if method = Hamiltonian us used...
Definition: current.F90:468
subroutine, public current_calculate_mel(der, hm, psi_i, psi_j, ik, cmel)
Definition: current.F90:766
subroutine, public current_calculate(this, namespace, gr, hm, space, st)
Compute total electronic current density.
Definition: current.F90:373
subroutine, public current_init(this, namespace)
Definition: current.F90:181
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
Definition: global.F90:191
integer, parameter, public rdmft
Definition: global.F90:237
integer, parameter, public hartree_fock
Definition: global.F90:237
complex(real64), parameter, public m_z0
Definition: global.F90:201
integer, parameter, public generalized_kohn_sham_dft
Definition: global.F90:237
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public dgrid_symmetrize_vector_field(gr, field, suppress_warning)
Definition: grid.F90:698
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] .
Definition: lda_u.F90:4994
subroutine, public zlda_u_commute_r_single(this, mesh, space, d, namespace, ist, ik, psi, gpsi, has_phase)
Definition: lda_u.F90:4943
subroutine, public magnetic_density(mesh, std, rho, md)
Definition: magnetic.F90:164
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public zprojector_commute_r_allatoms_alldir(pj, ions, mesh, dim, bnd, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
Definition: projector.F90:1758
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
Definition: types.F90:135
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
Definition: xc.F90:116
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
Definition: xc.F90:593
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
Definition: xc.F90:608
class representing derivatives
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:141
int true(void)