Octopus
phase.F90
Go to the documentation of this file.
1!! Copyright (C) 2009 X. Andrade
2!! Copyright (C) 2024 N. Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22module phase_oct_m
23 use accel_oct_m
24 use batch_oct_m
26 use comm_oct_m
27 use debug_oct_m
31 use global_oct_m
32 use grid_oct_m
35 use math_oct_m
36 use mesh_oct_m
38 use mpi_oct_m
40 use space_oct_m
44 use types_oct_m
46
47 implicit none
48
49 private
50
51 public :: &
53
84 type phase_t
85 private
86 complex(real64), allocatable :: phase(:, :)
89 complex(real64), public, allocatable :: phase_corr(:,:)
92 complex(real64), allocatable :: phase_spiral(:,:)
95 type(accel_mem_t) :: buff_phase
96 type(accel_mem_t) :: buff_phase_spiral
97 type(accel_mem_t), public :: buff_phase_corr
98 integer :: buff_phase_qn_start
99 real(real64), public, pointer :: spin(:,:,:) => null()
100 contains
101 procedure :: init => phase_init_phases
102
103 procedure :: update => phase_update_phases
104
105 procedure :: end => phase_end
106
107 procedure :: set_phase_corr => phase_set_phase_corr
108
109 procedure :: unset_phase_corr => phase_unset_phase_corr
110
111 procedure :: apply_to => phase_apply_batch
112
113 procedure :: apply_to_single => phase_apply_mf
114
115 procedure :: apply_phase_spiral => phase_phase_spiral
116
117 procedure :: is_allocated => phase_is_allocated
118 end type phase_t
119
120contains
121
122 ! ---------------------------------------------------------
124 subroutine phase_init_phases(phase, gr, kpt, kpoints, d, space)
125 class(phase_t), intent(inout) :: phase
126 type(grid_t), intent(in) :: gr
127 type(distributed_t), intent(in) :: kpt
128 type(kpoints_t), intent(in) :: kpoints
129 type(states_elec_dim_t), intent(in) :: d
130 type(space_t), intent(in) :: space
131
132 integer :: ip, ik, sp
133 integer(int64) :: ip_inner_global
134 real(real64) :: kpoint(space%dim), x_global(space%dim)
135
136 push_sub(phase_init_phases)
137
138 ! no e^ik phase needed for Gamma-point-only periodic calculations
139 ! unless for velocity-gauge for lasers
140 if (accel_is_enabled()) then
141 phase%buff_phase_qn_start = kpt%start
142 end if
143 if(kpoints%gamma_only()) then
144 pop_sub(phase_init_phases)
145 return
146 end if
147
148 safe_allocate(phase%phase(1:gr%np_part, kpt%start:kpt%end))
149 safe_allocate(phase%phase_corr(gr%np+1:gr%np_part, kpt%start:kpt%end))
150 !$omp parallel private(ip)
151 do ik = kpt%start, kpt%end
152 !$omp do
153 do ip = gr%np + 1, gr%np_part
154 phase%phase_corr(ip, ik) = m_one
155 end do
156 !$omp end do nowait
157 end do
158 !$omp end parallel
159
160 if (gr%der%boundaries%spiralBC) then
161 sp = gr%np
162 if (gr%parallel_in_domains) sp = gr%np + gr%pv%np_ghost
163
164 ! We decided to allocate the array from 1:np_part-sp as this is less error prone when passing
165 ! the array to other routines, or in particular creating a C-style pointer from phase_spiral(1,1).
166 ! We will also update phase_corr and possible other similar arrays.
167
168 safe_allocate(phase%phase_spiral(1:gr%np_part-sp, 1:2))
169
170 ! loop over boundary points
171 do ip = sp + 1, gr%np_part
172 ! get corresponding inner point
173 ip_inner_global = mesh_periodic_point(gr, space, ip)
174 x_global = mesh_x_global(gr, ip_inner_global)
175 phase%phase_spiral(ip-sp, 1) = &
176 exp(m_zi * sum((gr%x(ip, 1:space%dim)-x_global(1:space%dim)) * gr%der%boundaries%spiral_q(1:space%dim)))
177 phase%phase_spiral(ip-sp, 2) = &
178 exp(-m_zi * sum((gr%x(ip, 1:space%dim)-x_global(1:space%dim)) * gr%der%boundaries%spiral_q(1:space%dim)))
179 end do
180
181 if (accel_is_enabled()) then
182 call accel_create_buffer(phase%buff_phase_spiral, accel_mem_read_only, type_cmplx, (gr%np_part-sp)*2)
183 call accel_write_buffer(phase%buff_phase_spiral, (gr%np_part-sp)*2, phase%phase_spiral)
184 end if
185 end if
186
187
188 kpoint(1:space%dim) = m_zero
190 sp = gr%np
191 if (gr%parallel_in_domains) sp = gr%np + gr%pv%np_ghost
193 do ik = kpt%start, kpt%end
194 kpoint(1:space%dim) = kpoints%get_point(d%get_kpoint_index(ik))
195 !$omp parallel private(ip, ip_inner_global, x_global)
196 !$omp do
197 do ip = 1, gr%np_part
198 phase%phase(ip, ik) = exp(-m_zi * sum(gr%x(ip, 1:space%dim) * kpoint(1:space%dim)))
199 end do
200 !$omp end do
201
202 ! loop over boundary points
203 !$omp do
204 do ip = sp + 1, gr%np_part
205 ! get corresponding inner point
206 ip_inner_global = mesh_periodic_point(gr, space, ip)
207
208 ! compute phase correction from global coordinate (opposite sign!)
209 x_global = mesh_x_global(gr, ip_inner_global)
210 phase%phase_corr(ip, ik) = phase%phase(ip, ik)* &
211 exp(m_zi * sum(x_global(1:space%dim) * kpoint(1:space%dim)))
212 end do
213 !$omp end do nowait
214 !$omp end parallel
215 end do
216
218 call accel_create_buffer(phase%buff_phase, accel_mem_read_only, type_cmplx, gr%np_part*kpt%nlocal)
219 call accel_write_buffer(phase%buff_phase, gr%np_part*kpt%nlocal, phase%phase)
220 call accel_create_buffer(phase%buff_phase_corr, accel_mem_read_only, type_cmplx, (gr%np_part - gr%np)*kpt%nlocal)
221 call accel_write_buffer(phase%buff_phase_corr, (gr%np_part - gr%np)*kpt%nlocal, phase%phase_corr)
222 end if
223
224 pop_sub(phase_init_phases)
225 end subroutine phase_init_phases
226
227 ! ----------------------------------------------------------------------------------
229 subroutine phase_update_phases(phase, mesh, kpt, kpoints, d, space, uniform_vector_potential)
230 class(phase_t), intent(inout) :: phase
231 class(mesh_t), intent(in) :: mesh
232 type(distributed_t), intent(in) :: kpt
233 type(kpoints_t), intent(in) :: kpoints
234 type(states_elec_dim_t), intent(in) :: d
235 type(space_t), intent(in) :: space
236 real(real64), allocatable, intent(in) :: uniform_vector_potential(:)
237
238 integer :: ik, ip, sp
239 integer(int64) :: ip_inner_global
240 real(real64) :: kpoint(space%dim)
241 real(real64), allocatable :: x_global(:,:)
242
243 push_sub_with_profile(phase_update_phases)
244
245 if (allocated(uniform_vector_potential)) then
246
247 if (.not. allocated(phase%phase)) then
248 safe_allocate(phase%phase(1:mesh%np_part, kpt%start:kpt%end))
249 if (accel_is_enabled()) then
250 call accel_create_buffer(phase%buff_phase, accel_mem_read_only, type_cmplx, &
251 mesh%np_part*kpt%nlocal)
252 end if
253 end if
254
255 if (.not. allocated(phase%phase_corr)) then
256 safe_allocate(phase%phase_corr(mesh%np+1:mesh%np_part, kpt%start:kpt%end))
257 if (accel_is_enabled()) then
258 call accel_create_buffer(phase%buff_phase_corr, accel_mem_read_only, type_cmplx, &
259 (mesh%np_part - mesh%np)*kpt%nlocal)
260 end if
261 end if
262
263 kpoint(1:space%dim) = m_zero
264 ! loop over boundary points
265 sp = mesh%np
266 ! skip ghost points
267 if (mesh%parallel_in_domains) sp = mesh%np + mesh%pv%np_ghost
268
269 safe_allocate(x_global(1:space%dim,(sp + 1):mesh%np_part))
270 !$omp parallel private(ik, ip, kpoint, ip_inner_global)
271 !$omp do schedule(static)
272 do ip = sp + 1, mesh%np_part
273 ! get corresponding inner point
274 ip_inner_global = mesh_periodic_point(mesh, space, ip)
275
276 ! compute the difference between the global coordinate and the local coordinate
277 x_global(:,ip) = mesh_x_global(mesh, ip_inner_global) - mesh%x(ip, 1:space%dim)
278 end do
279
280 do ik = kpt%start, kpt%end
281 kpoint(1:space%dim) = kpoints%get_point(d%get_kpoint_index(ik))
282 !We add the vector potential
283 kpoint(1:space%dim) = kpoint(1:space%dim) + uniform_vector_potential(1:space%dim)
284
285 !$omp do schedule(static)
286 do ip = 1, mesh%np_part
287 phase%phase(ip, ik) = exp(cmplx(m_zero, -sum(mesh%x(ip, 1:space%dim)*kpoint(1:space%dim)), real64))
288 end do
289 !$omp end do
290
291 !$omp do schedule(static)
292 do ip = sp + 1, mesh%np_part
293 phase%phase_corr(ip, ik) = exp(cmplx(m_zero, sum(x_global(1:space%dim, ip) * kpoint(1:space%dim)), real64))
294 end do
295 !$omp end do nowait
296 end do
297 !$omp end parallel
298 safe_deallocate_a(x_global)
299
300 if (accel_is_enabled()) then
301 call accel_write_buffer(phase%buff_phase, mesh%np_part*kpt%nlocal, phase%phase, async=.true.)
302 call accel_write_buffer(phase%buff_phase_corr, (mesh%np_part - mesh%np)*kpt%nlocal, &
303 phase%phase_corr, async=.true.)
304 end if
305 end if
306
307 pop_sub_with_profile(phase_update_phases)
308 end subroutine phase_update_phases
309
310 ! ----------------------------------------------------------------------------------
312 subroutine phase_end(phase)
313 class(phase_t), intent(inout) :: phase
314
315 push_sub(phase_end)
316
317 if (phase%is_allocated() .and. accel_is_enabled()) then
318 call accel_release_buffer(phase%buff_phase)
319 call accel_release_buffer(phase%buff_phase_corr)
320 end if
321
322 if (allocated(phase%phase_spiral) .and. accel_is_enabled()) then
323 call accel_release_buffer(phase%buff_phase_spiral)
324 end if
325
326 safe_deallocate_a(phase%phase)
327 safe_deallocate_a(phase%phase_corr)
328 safe_deallocate_a(phase%phase_spiral)
329
330 pop_sub(phase_end)
331 end subroutine phase_end
332
333 ! ----------------------------------------------------------------------------------
335 !
336 subroutine phase_set_phase_corr(phase, mesh, psib, async)
337 class(phase_t), intent(in) :: phase
338 class(mesh_t), intent(in) :: mesh
339 type(wfs_elec_t), intent(inout) :: psib
340 logical, optional, intent(in) :: async
341
342
343 logical :: phase_correction
344
345 push_sub(phase_set_phase_corr)
346
347 ! check if we only want a phase correction for the boundary points
348 phase_correction = phase%is_allocated()
349
350 !We apply the phase only to np points, and the phase for the np+1 to np_part points
351 !will be treated as a phase correction in the Hamiltonian
352 if (phase_correction) then
353 call phase%apply_to(mesh, mesh%np, .false., psib, async=async)
354 end if
355
356 pop_sub(phase_set_phase_corr)
357 end subroutine phase_set_phase_corr
358
359 ! ----------------------------------------------------------------------------------
361 !
362 subroutine phase_unset_phase_corr(phase, mesh, psib, async)
363 class(phase_t), intent(in) :: phase
364 class(mesh_t), intent(in) :: mesh
365 type(wfs_elec_t), intent(inout) :: psib
366 logical, optional, intent(in) :: async
367
368 logical :: phase_correction
369
370 push_sub(phase_unset_phase_corr)
371
372 ! check if we only want a phase correction for the boundary points
373 phase_correction = phase%is_allocated()
374
375 !We apply the phase only to np points, and the phase for the np+1 to np_part points
376 !will be treated as a phase correction in the Hamiltonian
377 if (phase_correction) then
378 call phase%apply_to(mesh, mesh%np, .true., psib, async=async)
379 end if
380
382 end subroutine phase_unset_phase_corr
383
384 ! ---------------------------------------------------------------------------------------
386 !
387 subroutine phase_apply_batch(this, mesh, np, conjugate, psib, src, async)
388 class(phase_t), intent(in) :: this
389 class(mesh_t), intent(in) :: mesh
390 integer, intent(in) :: np
391 logical, intent(in) :: conjugate
392 type(wfs_elec_t), target, intent(inout) :: psib
393 type(wfs_elec_t), optional, target, intent(in) :: src
394 logical, optional, intent(in) :: async
395
396 integer :: ip, ii, sp
397 type(wfs_elec_t), pointer :: src_
398 complex(real64) :: phase
399 integer(int64) :: wgsize, dim2, dim3
400 type(accel_kernel_t), save :: ker_phase
401
402 push_sub(phase_apply_batch)
403 call profiling_in("PHASE_APPLY_BATCH")
404
405 call profiling_count_operations(6*np*psib%nst_linear)
406
407 assert(np <= mesh%np_part)
408 assert(psib%type() == type_cmplx)
409 assert(psib%ik >= lbound(this%phase, dim=2))
410 assert(psib%ik <= ubound(this%phase, dim=2))
411
412 src_ => psib
413 if (present(src)) src_ => src
414
415 assert(src_%has_phase .eqv. conjugate)
416 assert(src_%ik == psib%ik)
417 assert(src_%type() == type_cmplx)
418
419 ! We want to skip the ghost points for setting the phase
420 sp = min(np, mesh%np)
421 if (np > mesh%np .and. mesh%parallel_in_domains) sp = mesh%np + mesh%pv%np_ghost
422
423 select case (psib%status())
424 case (batch_packed)
425
426 if (conjugate) then
427
428 !$omp parallel private(ii, phase)
429 !$omp do
430 do ip = 1, min(mesh%np, np)
431 phase = conjg(this%phase(ip, psib%ik))
432 !$omp simd
433 do ii = 1, psib%nst_linear
434 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
435 end do
436 end do
437 !$omp end do nowait
438
439 ! Boundary points, if requested
440 !$omp do
441 do ip = sp+1, np
442 phase = conjg(this%phase(ip, psib%ik))
443 !$omp simd
444 do ii = 1, psib%nst_linear
445 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
446 end do
447 end do
448 !$omp end parallel
449
450 else
451
452 !$omp parallel private(ii, phase)
453 !$omp do
454 do ip = 1, min(mesh%np, np)
455 phase = this%phase(ip, psib%ik)
456 !$omp simd
457 do ii = 1, psib%nst_linear
458 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
459 end do
460 end do
461 !$omp end do nowait
462
463 ! Boundary points, if requested
464 !$omp do
465 do ip = sp+1, np
466 phase = this%phase(ip, psib%ik)
467 !$omp simd
468 do ii = 1, psib%nst_linear
469 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
470 end do
471 end do
472 !$omp end parallel
473
474 end if
475
476 case (batch_not_packed)
477
478 if (conjugate) then
479
480 !$omp parallel private(ii, ip)
481 do ii = 1, psib%nst_linear
482 !$omp do simd
483 do ip = 1, min(mesh%np, np)
484 psib%zff_linear(ip, ii) = conjg(this%phase(ip, psib%ik))*src_%zff_linear(ip, ii)
485 end do
486 !$omp end do simd nowait
487
488 ! Boundary points, if requested
489 !$omp do simd
490 do ip = sp+1, np
491 psib%zff_linear(ip, ii) = conjg(this%phase(ip, psib%ik))*src_%zff_linear(ip, ii)
492 end do
493 !$omp end do simd nowait
494 end do
495 !$omp end parallel
496
497 else
498 !$omp parallel private(ii, ip)
499 do ii = 1, psib%nst_linear
500 !$omp do simd
501 do ip = 1, min(mesh%np, np)
502 psib%zff_linear(ip, ii) = this%phase(ip, psib%ik)*src_%zff_linear(ip, ii)
503 end do
504 !$omp end do simd nowait
505
506 ! Boundary points, if requested
507 !$omp do simd
508 do ip = sp+1, np
509 psib%zff_linear(ip, ii) = this%phase(ip, psib%ik)*src_%zff_linear(ip, ii)
510 end do
511 !$omp end do simd nowait
512
513 end do
514 !$omp end parallel
515
516 end if
517
519 call accel_kernel_start_call(ker_phase, 'phase.cl', 'phase_hamiltonian')
520
521 if (conjugate) then
522 call accel_set_kernel_arg(ker_phase, 0, 1_4)
523 else
524 call accel_set_kernel_arg(ker_phase, 0, 0_4)
525 end if
526
527 call accel_set_kernel_arg(ker_phase, 1, (psib%ik - this%buff_phase_qn_start)*mesh%np_part)
528 call accel_set_kernel_arg(ker_phase, 2, np)
529 call accel_set_kernel_arg(ker_phase, 3, this%buff_phase)
530 call accel_set_kernel_arg(ker_phase, 4, src_%ff_device)
531 call accel_set_kernel_arg(ker_phase, 5, log2(int(src_%pack_size(1), int32)))
532 call accel_set_kernel_arg(ker_phase, 6, psib%ff_device)
533 call accel_set_kernel_arg(ker_phase, 7, log2(int(psib%pack_size(1), int32)))
534
535 wgsize = accel_kernel_workgroup_size(ker_phase)/psib%pack_size(1)
536
537 dim3 = np/(accel_max_size_per_dim(2)*wgsize) + 1
538 dim2 = min(accel_max_size_per_dim(2)*wgsize, pad(np, wgsize))
539
540 call accel_kernel_run(ker_phase, (/psib%pack_size(1), dim2, dim3/), (/psib%pack_size(1), wgsize, 1_int64/))
541
542 if(.not. optional_default(async, .false.)) call accel_finish()
543 end select
544
545 psib%has_phase = .not. conjugate
546
547 call profiling_out("PHASE_APPLY_BATCH")
548 pop_sub(phase_apply_batch)
549 end subroutine phase_apply_batch
550
556 !
557 subroutine phase_apply_mf(this, psi, np, dim, ik, conjugate)
558 class(phase_t), intent(in) :: this
559 complex(real64), intent(inout) :: psi(:, :)
560 integer, intent(in) :: np
561 integer, intent(in) :: dim
562 integer, intent(in) :: ik
563 logical, intent(in) :: conjugate
564
565 integer :: idim, ip
566
567 push_sub(phase_apply_mf)
568
569 assert(ik >= lbound(this%phase, dim=2))
570 assert(ik <= ubound(this%phase, dim=2))
571
572 call profiling_in("PHASE_APPLY_SINGLE")
573
574 if (conjugate) then
575 ! Apply the phase that contains both the k-point and vector-potential terms.
576 do idim = 1, dim
577 !$omp parallel do
578 do ip = 1, np
579 psi(ip, idim) = conjg(this%phase(ip, ik))*psi(ip, idim)
580 end do
581 !$omp end parallel do
582 end do
583 else
584 ! Apply the conjugate of (i.e. remove) the phase that contains both the k-point and vector-potential terms.
585 do idim = 1, dim
586 !$omp parallel do
587 do ip = 1, np
588 psi(ip, idim) = this%phase(ip, ik)*psi(ip, idim)
589 end do
590 !$omp end parallel do
591 end do
592 end if
593
594 call profiling_out("PHASE_APPLY_SINGLE")
595
596 pop_sub(phase_apply_mf)
597 end subroutine phase_apply_mf
598
599
600 ! ---------------------------------------------------------------------------------------
602 !
603 subroutine phase_phase_spiral(this, der, psib)
604 class(phase_t), intent(in) :: this
605 type(derivatives_t), intent(in) :: der
606 class(wfs_elec_t), intent(inout) :: psib
607
608 integer :: ip, ii, sp
609 integer, allocatable :: spin_label(:)
610 type(accel_mem_t) :: spin_label_buffer
611 integer(int64) :: wgsize
612
613 push_sub(phase_phase_spiral)
614 call profiling_in("PBC_PHASE_SPIRAL")
615
616 call profiling_count_operations(6*(der%mesh%np_part-der%mesh%np)*psib%nst_linear)
617
618 assert(der%boundaries%spiral)
619 assert(psib%type() == type_cmplx)
620
621 sp = der%mesh%np
622 if (der%mesh%parallel_in_domains) sp = der%mesh%np + der%mesh%pv%np_ghost
623
624
625 select case (psib%status())
626 case (batch_packed)
627
628 !$omp parallel do private(ip, ii)
629 do ip = sp + 1, der%mesh%np_part
630 do ii = 1, psib%nst_linear, 2
631 if (this%spin(3,psib%linear_to_ist(ii), psib%ik)>0) then
632 psib%zff_pack(ii+1, ip) = psib%zff_pack(ii+1, ip)*this%phase_spiral(ip-sp, 1)
633 else
634 psib%zff_pack(ii, ip) = psib%zff_pack(ii, ip)*this%phase_spiral(ip-sp, 2)
635 end if
636 end do
637 end do
638 !$omp end parallel do
639
640 case (batch_not_packed)
641
642 !$omp parallel private(ii, ip)
643 do ii = 1, psib%nst_linear, 2
644 if (this%spin(3,psib%linear_to_ist(ii), psib%ik)>0) then
645 !$omp do
646 do ip = sp + 1, der%mesh%np_part
647 psib%zff_linear(ip, ii+1) = psib%zff_linear(ip, ii+1)*this%phase_spiral(ip-sp, 1)
648 end do
649 !$omp end do nowait
650 else
651 !$omp do
652 do ip = sp + 1, der%mesh%np_part
653 psib%zff_linear(ip, ii) = psib%zff_linear(ip, ii)*this%phase_spiral(ip-sp, 2)
654 end do
655 !$omp end do nowait
656 end if
657 end do
658 !$omp end parallel
659
661
662 assert(accel_is_enabled())
663
664 ! generate array of offsets for access of psib and phase_spiral:
665 ! TODO: Move this to the routine where spin(:,:,:) is generated
666 ! and also move the buffer to the GPU at this point to
667 ! avoid unecessary latency here!
668
669 safe_allocate(spin_label(1:psib%nst_linear))
670 spin_label = 0
671 do ii = 1, psib%nst_linear, 2
672 if (this%spin(3, psib%linear_to_ist(ii), psib%ik) > 0) spin_label(ii)=1
673 end do
674
675 call accel_create_buffer(spin_label_buffer, accel_mem_read_only, type_integer, psib%nst_linear)
676 call accel_write_buffer(spin_label_buffer, psib%nst_linear, spin_label)
677
678 call accel_kernel_start_call(kernel_phase_spiral, 'phase_spiral.cl', 'phase_spiral_apply')
679
682 call accel_set_kernel_arg(kernel_phase_spiral, 2, der%mesh%np_part)
683 call accel_set_kernel_arg(kernel_phase_spiral, 3, psib%ff_device)
684 call accel_set_kernel_arg(kernel_phase_spiral, 4, log2(psib%pack_size(1)))
685 call accel_set_kernel_arg(kernel_phase_spiral, 5, this%buff_phase_spiral)
686 call accel_set_kernel_arg(kernel_phase_spiral, 6, spin_label_buffer)
687
688 wgsize = accel_kernel_workgroup_size(kernel_phase_spiral)/psib%pack_size(1)
689
691 (/psib%pack_size(1)/2, pad(der%mesh%np_part - sp, 2*wgsize)/), &
692 (/psib%pack_size(1)/2, 2*wgsize/))
693
694 call accel_finish()
695
696 call accel_release_buffer(spin_label_buffer)
697
698 safe_deallocate_a(spin_label)
699
700 end select
701
702 call profiling_out("PBC_PHASE_SPIRAL")
703 pop_sub(phase_phase_spiral)
704 end subroutine phase_phase_spiral
705
706
707 ! ---------------------------------------------------------------------------------------
708 logical pure function phase_is_allocated(this)
709 class(phase_t), intent(in) :: this
710
711 phase_is_allocated = allocated(this%phase)
712 end function phase_is_allocated
713
714end module phase_oct_m
715
716!! Local Variables:
717!! mode: f90
718!! coding: utf-8
719!! End:
double exp(double __x) __attribute__((__nothrow__
subroutine, public accel_kernel_start_call(this, file_name, kernel_name, flags)
Definition: accel.F90:2144
subroutine, public accel_finish()
Definition: accel.F90:1300
integer pure function, public accel_max_size_per_dim(dim)
Definition: accel.F90:2179
subroutine, public accel_release_buffer(this)
Definition: accel.F90:1250
type(accel_kernel_t), target, save, public kernel_phase_spiral
Definition: accel.F90:291
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
integer function, public accel_kernel_workgroup_size(kernel)
Definition: accel.F90:1482
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
This module implements batches of mesh functions.
Definition: batch.F90:133
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
Definition: batch.F90:282
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
Definition: batch.F90:282
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
Definition: batch.F90:282
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64), parameter, public m_zero
Definition: global.F90:188
complex(real64), parameter, public m_zi
Definition: global.F90:202
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
integer(int64) function, public mesh_periodic_point(mesh, space, ip)
This function returns the point inside the grid corresponding to a boundary point when PBCs are used....
Definition: mesh.F90:716
real(real64) function, dimension(1:mesh%box%dim), public mesh_x_global(mesh, ipg)
Definition: mesh.F90:804
subroutine phase_phase_spiral(this, der, psib)
apply spiral phase
Definition: phase.F90:697
subroutine phase_unset_phase_corr(phase, mesh, psib, async)
unset the phase correction (if necessary)
Definition: phase.F90:456
subroutine phase_init_phases(phase, gr, kpt, kpoints, d, space)
Initiliaze the phase arrays and copy to GPU the data.
Definition: phase.F90:218
subroutine phase_end(phase)
Releases the memory of the phase object.
Definition: phase.F90:406
subroutine phase_update_phases(phase, mesh, kpt, kpoints, d, space, uniform_vector_potential)
Update the phases.
Definition: phase.F90:323
logical pure function phase_is_allocated(this)
Definition: phase.F90:802
subroutine phase_apply_batch(this, mesh, np, conjugate, psib, src, async)
apply (remove) the phase to the wave functions before (after) applying the Hamiltonian
Definition: phase.F90:481
subroutine phase_set_phase_corr(phase, mesh, psib, async)
set the phase correction (if necessary)
Definition: phase.F90:430
subroutine phase_apply_mf(this, psi, np, dim, ik, conjugate)
apply (or remove) the phase to a wave function psi
Definition: phase.F90:651
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
This module handles spin dimensions of the states and the k-point distribution.
type(type_t), public type_cmplx
Definition: types.F90:134
type(type_t), public type_integer
Definition: types.F90:135
class representing derivatives
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....
Describes mesh distribution to nodes.
Definition: mesh.F90:186
A container for the phase.
Definition: phase.F90:177
class for organizing spins and k-points
batches of electronic states
Definition: wfs_elec.F90:139
int true(void)