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, &
76
77 integer, parameter, public :: &
78 CURRENT_GRADIENT = 1, &
81
82contains
83
84 subroutine current_init(this, namespace)
85 type(current_t), intent(out) :: this
86 type(namespace_t), intent(in) :: namespace
87
88 push_sub(current_init)
89
90 !%Variable CurrentDensity
91 !%Default gradient_corrected
92 !%Type integer
93 !%Section Hamiltonian
94 !%Description
95 !% This variable selects the method used to
96 !% calculate the current density. For the moment this variable is
97 !% for development purposes and users should not need to use
98 !% it.
99 !%Option gradient 1
100 !% The calculation of current is done using the gradient operator. (Experimental)
101 !%Option gradient_corrected 2
102 !% The calculation of current is done using the gradient operator
103 !% with additional corrections for the total current from non-local operators.
104 !%Option hamiltonian 3
105 !% The current density is obtained from the commutator of the
106 !% Hamiltonian with the position operator. (Experimental)
107 !%End
108
109 call parse_variable(namespace, 'CurrentDensity', current_gradient_corr, this%method)
110 if (.not. varinfo_valid_option('CurrentDensity', this%method)) call messages_input_error(namespace, 'CurrentDensity')
111 if (this%method /= current_gradient_corr) then
112 call messages_experimental("CurrentDensity /= gradient_corrected")
113 end if
114
115 !%Variable CalculateDiamagneticCurrent
116 !%Default no
117 !%Type logical
118 !%Section Hamiltonian
119 !%Description
120 !% This variable decides whether the current density arising from the non-uniform
121 !% vector potential, defined as:
122 !% <math> \vec{J}_{dmc}(\vec{r}, t)=-\frac{e^2}{m_e c_0} n(\vec{r}, t) \vec{A}(\vec{r},t)$ </math>
123 !% is included in the total current density.
124 !%End
125 call parse_variable(namespace, 'CalculateDiamagneticCurrent', .false., this%include_diamag)
126
127 pop_sub(current_init)
128 end subroutine current_init
129
130 ! ---------------------------------------------------------
131
132 subroutine current_batch_accumulate(st, der, ik, ib, psib, gpsib)
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(:)
139
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
144 real(real64) :: ww
145 real(real64), allocatable :: weight(:)
146 type(accel_mem_t) :: buff_weight, buff_current
147 type(accel_kernel_t), save :: kernel
148
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))
151
152 if (st%d%ispin == spinors .or. (psib%status() == batch_device_packed .and. der%dim /= 3)) then
153
154 do idir = 1, der%dim
155 do ist = states_elec_block_min(st, ib), states_elec_block_max(st, ib)
156
157 ww = st%kweights(ik)*st%occ(ist, ik)
158 if (abs(ww) <= m_epsilon) cycle
159
160 do idim = 1, st%d%dim
161 ii = st%group%psib(ib, ik)%inv_index((/ist, idim/))
162 call batch_get_state(psib, ii, der%mesh%np, psi(:, idim))
163 call batch_get_state(gpsib(idir), ii, der%mesh%np, gpsi(:, idim))
164 end do
165
166 if (st%d%ispin /= spinors) then
167 !$omp parallel do
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))
170 end do
171 !$omp end parallel do
172 else
173 !$omp parallel do private(c_tmp)
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)
180 end do
181 !$omp end parallel do
182 end if
183
184 end do
185 end do
186
187 else if (psib%status() == batch_device_packed) then
188
189 assert(der%dim == 3)
190
191 safe_allocate(weight(1:psib%nst))
192 do ist = 1, psib%nst
193 weight(ist) = st%kweights(ik)*st%occ(psib%ist(ist), ik)
194 end do
195
196
197 call accel_create_buffer(buff_weight, accel_mem_read_only, type_float, psib%nst)
198 call accel_write_buffer(buff_weight, psib%nst, weight, async=.true.)
199
200 call accel_create_buffer(buff_current, accel_mem_write_only, type_float, der%mesh%np*3)
201
202 call accel_kernel_start_call(kernel, 'density.cu', 'current_accumulate')
203
204 call accel_set_kernel_arg(kernel, 0, psib%nst)
205 call accel_set_kernel_arg(kernel, 1, der%mesh%np)
206 call accel_set_kernel_arg(kernel, 2, buff_weight)
207 call accel_set_kernel_arg(kernel, 3, psib%ff_device)
208 call accel_set_kernel_arg(kernel, 4, log2(int(psib%pack_size(1), int32)))
209 call accel_set_kernel_arg(kernel, 5, gpsib(1)%ff_device)
210 call accel_set_kernel_arg(kernel, 6, gpsib(2)%ff_device)
211 call accel_set_kernel_arg(kernel, 7, gpsib(3)%ff_device)
212 call accel_set_kernel_arg(kernel, 8, log2(int(gpsib(1)%pack_size(1), int32)))
213 call accel_set_kernel_arg(kernel, 9, buff_current)
214
215 ! Compute the grid size
216 bsize = accel_kernel_block_size(kernel)
217 call accel_grid_size(der%mesh%np, bsize, gsize)
218
219 call accel_kernel_run(kernel, (/gsize/), (/bsize/))
220
221 safe_allocate(current_tmp(1:der%mesh%np, 1:der%dim))
222
223 call accel_read_buffer(buff_current, der%mesh%np, der%dim, current_tmp)
224
225 call lalg_axpy(der%mesh%np, der%dim, m_one, current_tmp, st%current_kpt(:,:,ik))
226
227 safe_deallocate_a(current_tmp)
228
229 call accel_free_buffer(buff_weight)
230 call accel_free_buffer(buff_current)
231
232 safe_deallocate_a(weight)
233
234 else
235
236 assert(psib%is_packed() .eqv. gpsib(1)%is_packed())
237
238 !$omp parallel private(ip, ist, ww, idir)
239 do ii = 1, psib%nst
240 ist = states_elec_block_min(st, ib) + ii - 1
241 ww = st%kweights(ik)*st%occ(ist, ik)
242 if (abs(ww) <= m_epsilon) cycle
243
244 if (psib%is_packed()) then
245 do idir = 1, der%dim
246 !$omp do
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))
250 end do
251 !$omp end do nowait
252 end do
253 else
254 do idir = 1, der%dim
255 !$omp do
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))
259 end do
260 !$omp end do nowait
261 end do
262 end if
263 end do
264 !$omp end parallel
265
266 end if
267
268 safe_deallocate_a(psi)
269 safe_deallocate_a(gpsi)
270
271 end subroutine current_batch_accumulate
272
273
274 ! ---------------------------------------------------------
276 subroutine current_calculate(this, namespace, gr, hm, space, st)
277 type(current_t), intent(in) :: this
278 type(namespace_t), intent(in) :: namespace
279 type(grid_t), intent(inout) :: gr
280 type(hamiltonian_elec_t), intent(inout) :: hm
281 class(space_t), intent(in) :: space
282 type(states_elec_t), intent(inout) :: st
283
284
285 call profiling_in("CURRENT_TOTAL")
286 push_sub(current_calculate)
287
288 call current_calculate_para(this, namespace, gr, hm, space, st)
289 st%current = st%current_para
290
291 if (this%include_diamag) then
292 call current_calculate_dia_non_unif_vec_pot(gr%der, hm, st)
293 st%current = st%current + st%current_dia
294 end if
295
296 call profiling_out("CURRENT_TOTAL")
297 pop_sub(current_calculate)
298
299 end subroutine current_calculate
300
301 ! ---------------------------------------------------------
304 subroutine current_calculate_dia_non_unif_vec_pot(der, hm, st)
305 type(derivatives_t), intent(inout) :: der
306 type(hamiltonian_elec_t), intent(in) :: hm
307 type(states_elec_t), intent(inout) :: st
308
309 integer :: ispin, idir, ip
310
311 call profiling_in("CURRENT_DIA_NON_UNIF_A")
313
314 st%current_dia = m_zero
315
316 if(allocated(hm%hm_base%vector_potential)) then
317 do ispin = 1, st%d%nspin
318 do idir = 1, der%dim
319 !$omp parallel do
320 do ip = 1, der%mesh%np
321 ! the vector potential is assumed to be devided by c_0 already
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)
324 end do
325 !$omp end parallel do
326 end do
327 end do
328 end if
329
330 call profiling_out("CURRENT_DIA_NON_UNIF_A")
332
334
335 ! ---------------------------------------------------------
339 subroutine current_calculate_mag(der, st)
340 type(derivatives_t), intent(inout) :: der
341 type(states_elec_t), intent(inout) :: st
342
343 real(real64), allocatable :: magnetization_density(:, :), curl_mag(:, :)
344
345 call profiling_in("CURRENT_MAG")
346 push_sub(current_calculate_mag)
347
348 st%current_mag = m_zero
349
350 if (st%d%ispin /= unpolarized) then
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))
353
354 call magnetic_density(der%mesh, st%d, st%rho, magnetization_density)
355 call dderivatives_curl(der, magnetization_density, curl_mag)
356
357 call lalg_axpy(der%mesh%np, der%dim, m_half, curl_mag, st%current_mag(:, :, 1))
358
359 safe_deallocate_a(magnetization_density)
360 safe_deallocate_a(curl_mag)
361 end if
362
363 call profiling_out("CURRENT_MAG")
364 pop_sub(current_calculate_mag)
365
366 end subroutine current_calculate_mag
367
368 ! ---------------------------------------------------------
371 subroutine current_calculate_para(this, namespace, gr, hm, space, st)
372 type(current_t), intent(in) :: this
373 type(namespace_t), intent(in) :: namespace
374 type(grid_t), intent(inout) :: gr
375 type(hamiltonian_elec_t), intent(inout) :: hm
376 class(space_t), intent(in) :: space
377 type(states_elec_t), target, intent(inout) :: st
378
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
382 class(wfs_elec_t), allocatable :: commpsib(:)
383 real(real64) :: ww
384 complex(real64) :: c_tmp
385
386 call profiling_in("CURRENT_PARA")
387 push_sub(current_calculate_para)
388
389 ! spin not implemented or tested
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/)))
393
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))
401
402 !$omp parallel private(idir, ip, ispin)
403 do ik = st%d%kpt%start, st%d%kpt%end
404 do idir = 1, gr%der%dim
405 !$omp do
406 do ip = 1, gr%np
407 st%current_kpt(ip, idir, ik) = m_zero
408 end do
409 !$omp end do nowait
410 end do
411 end do
412 do ispin = 1, st%d%nspin
413 do idir = 1, gr%der%dim
414 !$omp do
415 do ip = 1, gr%np
416 st%current_para(ip, idir, ispin) = m_zero
417 end do
418 !$omp end do nowait
419 end do
420 end do
421 !$omp end parallel
422
423 select case (this%method)
424
426 if (family_is_hybrid(hm%xc)) then
427 if (.not. hm%exxop%useACE) then
428 hm%exxop%st => st
430 end if
431 end if
432
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
436
437 call st%group%psib(ib, ik)%do_pack(copy = .true.)
438
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)
443
444 call boundaries_set(gr%der%boundaries, gr, st%group%psib(ib, ik))
445 call zhamiltonian_elec_apply_batch(hm, namespace, gr, st%group%psib(ib, ik), hpsib, set_bc = .false.)
446
447 do idir = 1, gr%der%dim
448
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)
451
452 call zhamiltonian_elec_apply_batch(hm, namespace, gr, rpsib, hrpsib, set_bc = .false.)
453
454 do ist = states_elec_block_min(st, ib), states_elec_block_max(st, ib)
455 ww = st%kweights(ik)*st%occ(ist, ik)
456 if (ww <= m_epsilon) cycle
457
458 do idim = 1, st%d%dim
459 ii = st%group%psib(ib, ik)%inv_index((/ist, idim/))
460 call batch_get_state(st%group%psib(ib, ik), ii, gr%np, psi(:, idim))
461 call batch_get_state(hrpsib, ii, gr%np, hrpsi(:, idim))
462 call batch_get_state(rhpsib, ii, gr%np, rhpsi(:, idim))
463 end do
464
465 if (st%d%ispin /= spinors) then
466 !$omp parallel do
467 do ip = 1, gr%np
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))
470 end do
471 !$omp end parallel do
472 else
473 !$omp parallel do private(c_tmp)
474 do ip = 1, gr%np
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)
483 end do
484 !$omp end parallel do
485 end if
486
487 end do
488
489 end do
490
491 call st%group%psib(ib, ik)%do_unpack(copy = .false.)
492
493 call hpsib%end()
494 call rhpsib%end()
495 call rpsib%end()
496 call hrpsib%end()
497
498 end do
499 end do
500
501 if (family_is_hybrid(hm%xc)) then
502 if (.not. hm%exxop%useACE) then
504 nullify(hm%exxop%st)
505 end if
506 end if
507
508 case (current_gradient, current_gradient_corr)
509
510 if (this%method == current_gradient_corr .and. .not. family_is_mgga_with_exc(hm%xc) &
511 .and. hm%theory_level /= hartree_fock &
512 .and. hm%theory_level /= generalized_kohn_sham_dft &
513 .and. hm%theory_level /= rdmft &
514 .and. hm%vnl%apply_projector_matrices) then
515
516 ! we can use the packed version
517
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
521
522 ! copy st%group%psib(ib, ik) to epsib and set the phase
523 call hm%phase%copy_and_set_phase(gr, st%d%kpt, st%group%psib(ib, ik), epsib)
524
525 ! this now takes non-orthogonal axis into account
526 call zderivatives_batch_grad(gr%der, epsib, commpsib, set_bc=.false.)
527
528 call hm%vnl%zposition_commutator(gr, st%d, gr%der%boundaries%spiral, epsib, commpsib, async=.true.)
529
530 call zlda_u_commute_r(hm%lda_u, gr, space, st%d, namespace, epsib, commpsib)
531
532 call current_batch_accumulate(st, gr%der, ik, ib, epsib, commpsib)
533
534 do idir = 1, gr%der%dim
535 call commpsib(idir)%end()
536 end do
537
538 call epsib%end()
539
540 call accel_finish()
541 end do
542 end do
543
544 else
545
546 ! use the slow non-packed version
547
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
551
552 ww = st%kweights(ik)*st%occ(ist, ik)
553 if (abs(ww) <= m_epsilon) cycle
554
555 call states_elec_get_state(st, gr, ist, ik, psi)
556
557 if (hm%phase%is_allocated()) then
558 call hm%phase%apply_to_single(psi, gr%np, st%d%dim, ik, .false.)
559 ! apply phase correction while setting boundary -> memory needs to be
560 ! accessed only once
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)))
564 end do
565 else
566 do idim = 1, st%d%dim
567 call boundaries_set(gr%der%boundaries, gr, psi(:, idim))
568 end do
569 end if
570
571 do idim = 1, st%d%dim
572 call zderivatives_grad(gr%der, psi(:, idim), gpsi(:, :, idim), set_bc = .false.)
573 end do
574
575 if (this%method == current_gradient_corr) then
576 !A nonlocal contribution from the MGGA potential must be included
577 !This must be done first, as this is like a position-dependent mass
578 call hm%ks_pot%zcurrent_mass_renormalization(gpsi, gr%der%dim, st%d%dim, ispin)
579
580 !A nonlocal contribution from the pseudopotential must be included
581 call zprojector_commute_r_allatoms_alldir(hm%ep%proj, hm%ions, gr, st%d%dim, &
582 gr%der%boundaries, ik, psi, gpsi)
583 !A nonlocal contribution from the scissor must be included
584 if (hm%scissor%apply) then
585 call scissor_commute_r(hm%scissor, gr, ik, psi, gpsi)
586 end if
587
588 call zlda_u_commute_r_single(hm%lda_u, gr, space, st%d, namespace, ist, ik, &
589 psi, gpsi, hm%phase%is_allocated())
590
591 call zexchange_operator_commute_r(hm%exxop, namespace, gr, st%d, ik, psi, gpsi)
592
593 end if
594
595 if (st%d%ispin /= spinors) then
596 do idir = 1, gr%der%dim
597 !$omp parallel do
598 do ip = 1, gr%np
599 st%current_kpt(ip, idir, ik) = st%current_kpt(ip, idir, ik) + &
600 ww*aimag(conjg(psi(ip, 1))*gpsi(ip, idir, 1))
601 end do
602 !$omp end parallel do
603 end do
604 else
605 do idir = 1, gr%der%dim
606 !$omp parallel do private(c_tmp)
607 do ip = 1, gr%np
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)
615 end do
616 !$omp end parallel do
617 end do
618 end if
619
620 end do
621 end do
622
623 end if
624
625 case default
626
627 assert(.false.)
628
629 end select
630
631 if (st%d%ispin /= spinors) then
632 !We sum the current over k-points
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))
636 end do
637 end if
638
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/))
641 end if
642
643 if (st%symmetrize_density) then
644 do ispin = 1, st%d%nspin
645 call dgrid_symmetrize_vector_field(gr, st%current_para(:, :, ispin), suppress_warning = .true.)
646 end do
647 end if
648
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)
656
657 call profiling_out("CURRENT_PARA")
659
660 end subroutine current_calculate_para
661
662 ! ---------------------------------------------------------
663 subroutine current_heat_calculate(space, der, hm, st, current)
664 class(space_t), intent(in) :: space
665 type(derivatives_t), intent(in) :: der
666 type(hamiltonian_elec_t), intent(in) :: hm
667 type(states_elec_t), intent(in) :: st
668 real(real64), intent(out) :: current(:, :, :)
669
670 integer :: ik, ist, idir, idim, ip, ispin, ndim
671 complex(real64), allocatable :: gpsi(:, :, :), psi(:, :), g2psi(:, :, :, :)
672 complex(real64) :: tmp
673
674 push_sub(current_heat_calculate)
675
676 assert(space%is_periodic())
677 assert(st%d%dim == 1)
678
679 ndim = space%dim
680
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))
684
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)
687 end do
688
689
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
693
694 if (abs(st%kweights(ik)*st%occ(ist, ik)) <= m_epsilon) cycle
695
696 call states_elec_get_state(st, der%mesh, ist, ik, psi)
697 do idim = 1, st%d%dim
698 call boundaries_set(der%boundaries, der%mesh, psi(:, idim))
699 end do
700
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.)
703 end if
704
705 do idim = 1, st%d%dim
706 call zderivatives_grad(der, psi(:, idim), gpsi(:, :, idim), set_bc = .false.)
707 end do
708 do idir = 1, ndim
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.)
711 end if
712
713 do idim = 1, st%d%dim
714 call boundaries_set(der%boundaries, der%mesh, gpsi(:,idir, idim))
715 end do
716
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.)
719 end if
720
721 do idim = 1, st%d%dim
722 call zderivatives_grad(der, gpsi(:, idir, idim), g2psi(:, :, idir, idim), set_bc = .false.)
723 end do
724 end do
725 idim = 1
726 do ip = 1, der%mesh%np
727 do idir = 1, ndim
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
733 end do
734 end do
735 end do
736 end do
737
738
740
741 end subroutine current_heat_calculate
742
743end module current_oct_m
744
745!! Local Variables:
746!! mode: f90
747!! coding: utf-8
748!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:173
integer function, public accel_kernel_block_size(kernel)
Definition: accel.F90:1188
subroutine, public accel_free_buffer(this, async)
Definition: accel.F90:1005
subroutine, public accel_kernel_start_call(this, file_name, kernel_name, flags)
Definition: accel.F90:1413
subroutine, public accel_finish()
Definition: accel.F90:1098
integer, parameter, public accel_mem_write_only
Definition: accel.F90:185
integer, parameter, public accel_mem_read_only
Definition: accel.F90:185
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:228
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:400
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:435
integer, parameter, public current_hamiltonian
Definition: current.F90:172
integer, parameter, public current_gradient_corr
Definition: current.F90:172
subroutine, public current_heat_calculate(space, der, hm, st, current)
Definition: current.F90:759
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:467
subroutine, public current_calculate(this, namespace, gr, hm, space, st)
Compute total electronic current density.
Definition: current.F90:372
subroutine, public current_init(this, namespace)
Definition: current.F90:180
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:241
integer, parameter, public hartree_fock
Definition: global.F90:241
integer, parameter, public generalized_kohn_sham_dft
Definition: global.F90:241
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:5088
subroutine, public zlda_u_commute_r_single(this, mesh, space, d, namespace, ist, ik, psi, gpsi, has_phase)
Definition: lda_u.F90:5037
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:690
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
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)