Octopus
nonlocal_pseudopotential.F90
Go to the documentation of this file.
1!! Copyright (C) 2009 X. Andrade
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
22 use accel_oct_m
24 use batch_oct_m
26 use blas_oct_m
27 use comm_oct_m
28 use debug_oct_m
30 use epot_oct_m
31 use global_oct_m
36 use math_oct_m
37 use mesh_oct_m
39 use mpi_oct_m
43 use ps_oct_m
45 use space_oct_m
49 use types_oct_m
51
52 implicit none
53
54 private
55
56 public :: &
61
65 private
66 type(projector_matrix_t), allocatable, public :: projector_matrices(:)
67 integer, public :: nprojector_matrices
68 logical, public :: apply_projector_matrices
69 logical, public :: has_non_local_potential
70 integer :: full_projection_size
71 integer, public :: max_npoints
72 integer, public :: total_points
73 integer :: max_nprojs
74 logical :: projector_mix
75 complex(real64), allocatable, public :: projector_phases(:, :, :, :)
76 integer, allocatable, public :: projector_to_atom(:)
77 integer :: nregions
78 integer, allocatable :: regions(:)
79 integer, public :: nphase
84 type(accel_mem_t) :: buff_offsets
85 type(accel_mem_t) :: buff_matrices
86 type(accel_mem_t) :: buff_maps
87 type(accel_mem_t) :: buff_scals
88 type(accel_mem_t) :: buff_position
89 type(accel_mem_t) :: buff_pos
90 type(accel_mem_t) :: buff_invmap
91 type(accel_mem_t), public :: buff_projector_phases
92 type(accel_mem_t) :: buff_mix
93 logical :: projector_self_overlap
94 real(real64), pointer, public :: spin(:,:,:) => null()
95 contains
96
97 procedure :: init => nonlocal_pseudopotential_init
98
99 procedure :: build => nonlocal_pseudopotential_build_proj
100
101 procedure :: end => nonlocal_pseudopotential_destroy_proj
102
103 procedure :: has_self_overlap => nonlocal_pseudopotential_self_overlap
104
105 procedure :: dstart => dnonlocal_pseudopotential_start
106
107 procedure :: zstart => znonlocal_pseudopotential_start
108
109 procedure :: dfinish => dnonlocal_pseudopotential_finish
110
111 procedure :: zfinish => znonlocal_pseudopotential_finish
112
113 procedure :: dforce => dnonlocal_pseudopotential_force
114
115 procedure :: zforce => znonlocal_pseudopotential_force
117 procedure :: dposition_commutator => dnonlocal_pseudopotential_position_commutator
118
119 procedure :: zposition_commutator => znonlocal_pseudopotential_position_commutator
120
121 procedure :: dr_vn_local => dnonlocal_pseudopotential_r_vnlocal
122
123 procedure :: zr_vn_local => znonlocal_pseudopotential_r_vnlocal
124
126
128 !
129 type projection_t
130 private
131 real(real64), allocatable :: dprojection(:, :)
132 complex(real64), allocatable :: zprojection(:, :)
133 type(accel_mem_t) :: buff_projection
134 type(accel_mem_t) :: buff_spin_to_phase
135 type(accel_mem_t), allocatable :: buff_phasepsi(:)
136 type(accel_mem_t), allocatable :: buff_projection_temp(:)
137 integer :: cuda_stream_projection_DtH
138 end type projection_t
139
140contains
141
142 ! ---------------------------------------------------------
145 subroutine nonlocal_pseudopotential_init(this)
146 class(nonlocal_pseudopotential_t), intent(inout) :: this
147
149
150 this%apply_projector_matrices = .false.
151 this%has_non_local_potential = .false.
152 this%nprojector_matrices = 0
153
154 this%projector_self_overlap = .false.
155
157 end subroutine nonlocal_pseudopotential_init
158
159 ! ---------------------------------------------------------
165 class(nonlocal_pseudopotential_t), target, intent(inout) :: this
166 type(epot_t), target, intent(in) :: epot
168 integer :: imat, iatom
169 type(projector_matrix_t), pointer :: pmat
173 if (.not. allocated(this%projector_matrices) .or. this%nprojector_matrices == 0) then
175 return
176 end if
177
178 assert(allocated(this%projector_to_atom))
179 assert(size(this%projector_matrices) >= this%nprojector_matrices)
180 assert(size(this%projector_to_atom) >= this%nprojector_matrices)
182 do imat = 1, this%nprojector_matrices
183 pmat => this%projector_matrices(imat)
184 iatom = this%projector_to_atom(imat)
186 assert(iatom >= 1 .and. iatom <= epot%natoms)
187 assert(allocated(epot%proj))
188 assert(allocated(epot%proj(iatom)%sphere%map))
189 assert(allocated(epot%proj(iatom)%sphere%rel_x))
190
191 pmat%map => epot%proj(iatom)%sphere%map
192 pmat%position => epot%proj(iatom)%sphere%rel_x
193 end do
197
199 !--------------------------------------------------------
202 class(nonlocal_pseudopotential_t), intent(inout) :: this
203
204 integer :: iproj
205
207
208 if (allocated(this%projector_matrices)) then
209
211 call accel_free_buffer(this%buff_offsets)
212 call accel_free_buffer(this%buff_matrices)
213 call accel_free_buffer(this%buff_maps)
214 call accel_free_buffer(this%buff_scals)
215 call accel_free_buffer(this%buff_position)
216 call accel_free_buffer(this%buff_pos)
217 call accel_free_buffer(this%buff_invmap)
218 if (this%projector_mix) call accel_free_buffer(this%buff_mix)
219 if (allocated(this%projector_phases)) call accel_free_buffer(this%buff_projector_phases)
220 end if
221
222 do iproj = 1, this%nprojector_matrices
223 call projector_matrix_deallocate(this%projector_matrices(iproj))
224 end do
225 safe_deallocate_a(this%regions)
226 safe_deallocate_a(this%projector_matrices)
227 safe_deallocate_a(this%projector_phases)
228 safe_deallocate_a(this%projector_to_atom)
229 end if
233
234 !-----------------------------------------------------------------
241 subroutine nonlocal_pseudopotential_build_proj(this, space, mesh, epot)
242 class(nonlocal_pseudopotential_t), target, intent(inout) :: this
243 class(space_t), intent(in) :: space
244 class(mesh_t), intent(in) :: mesh
245 type(epot_t), target, intent(in) :: epot
246
247 integer :: iatom, iproj, ll, lmax, lloc, mm, ic, jc
248 integer :: nmat, imat, ip, iorder
249 integer :: nregion, jatom, katom, iregion
250 integer, allocatable :: order(:), head(:), region_count(:)
251 logical, allocatable :: atom_counted(:)
252 logical :: overlap
253 type(projector_matrix_t), pointer :: pmat
254 type(kb_projector_t), pointer :: kb_p
255 type(rkb_projector_t), pointer :: rkb_p
256 type(hgh_projector_t), pointer :: hgh_p
257
260 call profiling_in("ATOM_COLORING")
261
262 ! this is most likely a very inefficient algorithm, O(natom**2) or
263 ! O(natom**3), probably it should be replaced by something better.
264
265 safe_allocate(order(1:epot%natoms)) ! order(iregion) = ?
266 safe_allocate(head(1:epot%natoms + 1)) ! head(iregion) points to the first atom in region iregion
267 safe_allocate(region_count(1:epot%natoms)) ! region_count(iregion): number of atoms in that region
268 safe_allocate(atom_counted(1:epot%natoms))
269
270 this%projector_self_overlap = .false.
271 atom_counted = .false.
272 order = -1
273
274 head(1) = 1
275 nregion = 0
276 do
277 nregion = nregion + 1
278 assert(nregion <= epot%natoms)
279
280 region_count(nregion) = 0
281
282 do iatom = 1, epot%natoms
283 if (atom_counted(iatom)) cycle
284
285 overlap = .false.
286
287 if (.not. projector_is(epot%proj(iatom), proj_none)) then
288 assert(associated(epot%proj(iatom)%sphere%mesh))
289 do jatom = 1, region_count(nregion)
290 katom = order(head(nregion) + jatom - 1)
291 if (projector_is(epot%proj(katom), proj_none)) cycle
292 overlap = submesh_overlap(epot%proj(iatom)%sphere, epot%proj(katom)%sphere, space)
293 if (overlap) exit
294 end do
295 end if
297 if (.not. overlap) then
298 ! iatom did not overlap with any previously counted atoms:
299 ! iatom will be added to the current region
300 region_count(nregion) = region_count(nregion) + 1
301 order(head(nregion) - 1 + region_count(nregion)) = iatom
302 atom_counted(iatom) = .true.
303 end if
304
305 end do
306
307 head(nregion + 1) = head(nregion) + region_count(nregion)
308
309 if (all(atom_counted)) exit
310 end do
311
312 safe_deallocate_a(atom_counted)
313 safe_deallocate_a(region_count)
314
315 call messages_write('The atoms can be separated in ')
316 call messages_write(nregion)
317 call messages_write(' non-overlapping groups.')
318 call messages_info(debug_only=.true.)
319
320 do iregion = 1, nregion
321 do iatom = head(iregion), head(iregion + 1) - 1
322 if (.not. projector_is(epot%proj(order(iatom)), proj_kb)) cycle
323 do jatom = head(iregion), iatom - 1
324 if (.not. projector_is(epot%proj(order(jatom)), proj_kb)) cycle
325 assert(.not. submesh_overlap(epot%proj(order(iatom))%sphere, epot%proj(order(jatom))%sphere, space))
326 end do
327 end do
328 end do
329
330 call profiling_out("ATOM_COLORING")
331
332 ! deallocate previous projectors
333 call this%end()
334
335 ! count projectors
336 this%nprojector_matrices = 0
337 this%apply_projector_matrices = .false.
338 this%has_non_local_potential = .false.
339 this%nregions = nregion
340
341 !We determine if we have only local potential or not.
342 do iorder = 1, epot%natoms
343 iatom = order(iorder)
344
345 if (.not. projector_is_null(epot%proj(iatom))) then
346 this%has_non_local_potential = .true.
347 exit
348 end if
349 end do
350
351 do iorder = 1, epot%natoms
352 iatom = order(iorder)
353
354 if (.not. projector_is_null(epot%proj(iatom))) then
355 this%nprojector_matrices = this%nprojector_matrices + 1
356 this%apply_projector_matrices = .true.
357 end if
358 end do
359
360 ! This is currently the only not supported case
361 if (mesh%use_curvilinear) this%apply_projector_matrices = .false.
362
363 if (.not. this%apply_projector_matrices) then
364 safe_deallocate_a(order)
365 safe_deallocate_a(head)
366
368 return
369 end if
370
371
372 safe_allocate(this%projector_matrices(1:this%nprojector_matrices))
373 safe_allocate(this%regions(1:this%nprojector_matrices + 1))
374 safe_allocate(this%projector_to_atom(1:epot%natoms))
375
376 this%full_projection_size = 0
377 this%regions(this%nregions + 1) = this%nprojector_matrices + 1
378
379 this%projector_mix = .false.
380
381 iproj = 0
382 do iregion = 1, this%nregions
383 this%regions(iregion) = iproj + 1
384 do iorder = head(iregion), head(iregion + 1) - 1
385
386 iatom = order(iorder)
387
388 if (projector_is(epot%proj(iatom), proj_none)) cycle
389
390 iproj = iproj + 1
391
392 pmat => this%projector_matrices(iproj)
393
394 this%projector_to_atom(iproj) = iatom
395
396 lmax = epot%proj(iatom)%lmax
397 lloc = epot%proj(iatom)%lloc
398
399 if (projector_is(epot%proj(iatom), proj_kb)) then
400
401 ! count the number of projectors for this matrix
402 nmat = 0
403 do ll = 0, lmax
404 if (ll == lloc) cycle
405 do mm = -ll, ll
406 nmat = nmat + epot%proj(iatom)%kb_p(ll, mm)%n_c
407 end do
408 end do
409
410 call projector_matrix_allocate(pmat, nmat, epot%proj(iatom)%sphere, has_mix_matrix = .false.)
411
412 ! generate the matrix
413 pmat%dprojectors = m_zero
414 imat = 1
415 do ll = 0, lmax
416 if (ll == lloc) cycle
417 do mm = -ll, ll
418 kb_p => epot%proj(iatom)%kb_p(ll, mm)
419 do ic = 1, kb_p%n_c
420 call lalg_copy(pmat%npoints, kb_p%p(:, ic), pmat%dprojectors(:, imat))
421 pmat%scal(imat) = kb_p%e(ic)*mesh%vol_pp(1)
422 imat = imat + 1
423 end do
424 end do
425 end do
426
427 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
428
429 else if (projector_is(epot%proj(iatom), proj_hgh)) then
430
431 this%projector_mix = .true.
432
433 ! count the number of projectors for this matrix
434 nmat = 0
435 do ll = 0, lmax
436 if (ll == lloc) cycle
437 do mm = -ll, ll
438 nmat = nmat + 3
439 end do
440 end do
441
442 call projector_matrix_allocate(pmat, nmat, epot%proj(iatom)%sphere, &
443 has_mix_matrix = .true., is_cmplx = (epot%reltype == spin_orbit))
444
445 ! generate the matrix
446 if (epot%reltype == spin_orbit) then
447 pmat%zprojectors = m_zero
448 pmat%zmix = m_zero
449 else
450 pmat%dprojectors = m_zero
451 pmat%dmix = m_zero
452 end if
453
454 imat = 1
455 do ll = 0, lmax
456 if (ll == lloc) cycle
457 do mm = -ll, ll
458 hgh_p => epot%proj(iatom)%hgh_p(ll, mm)
459
460 ! HGH pseudos mix different components, so we need to
461 ! generate a matrix that mixes the projections
462 if (epot%reltype == spin_orbit .or. epot%reltype == fully_relativistic_zora) then
463 do ic = 1, 3
464 do jc = 1, 3
465 pmat%zmix(imat - 1 + ic, imat - 1 + jc, 1) = hgh_p%h(ic, jc) + m_half*mm*hgh_p%k(ic, jc)
466 pmat%zmix(imat - 1 + ic, imat - 1 + jc, 2) = hgh_p%h(ic, jc) - m_half*mm*hgh_p%k(ic, jc)
467
468 if (mm < ll) then
469 pmat%zmix(imat - 1 + ic, imat + 3 - 1 + jc, 3) = m_half*hgh_p%k(ic, jc) * &
470 sqrt(real(ll*(ll+1)-mm*(mm+1), real64))
471 end if
472
473 if (-mm < ll) then
474 pmat%zmix(imat - 1 + ic, imat - 3 - 1 + jc, 4) = m_half*hgh_p%k(ic, jc) * &
475 sqrt(real(ll*(ll+1)-mm*(mm-1), real64))
476 end if
477 end do
478 end do
479 else
480 do ic = 1, 3
481 do jc = 1, 3
482 pmat%dmix(imat - 1 + ic, imat - 1 + jc) = hgh_p%h(ic, jc)
483 end do
484 end do
485 end if
486
487 do ic = 1, 3
488 if (epot%reltype == spin_orbit .or. epot%reltype == fully_relativistic_zora) then
489 call lalg_copy(pmat%npoints, hgh_p%zp(:, ic), pmat%zprojectors(:, imat))
490 else
491 call lalg_copy(pmat%npoints, hgh_p%dp(:, ic), pmat%dprojectors(:, imat))
492 end if
493 pmat%scal(imat) = mesh%volume_element
494 imat = imat + 1
495 end do
496
497 end do
498 end do
499
500 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
501
502 else if (projector_is(epot%proj(iatom), proj_rkb)) then
503 assert(epot%reltype == spin_orbit)
504
505 this%projector_mix = .true.
506
507 ! count the number of projectors for this matrix
508 nmat = 0
509 if (lloc /= 0) nmat = nmat + epot%proj(iatom)%kb_p(1, 1)%n_c
510
511 do ll = 1, lmax
512 if (ll == lloc) cycle
513 do mm = -ll, ll
514 nmat = nmat + epot%proj(iatom)%rkb_p(ll, mm)%n_c
515 end do
516 end do
517
518 call projector_matrix_allocate(pmat, nmat, epot%proj(iatom)%sphere, &
519 has_mix_matrix = .true., is_cmplx = .true.)
520
521 pmat%zprojectors = m_zero
522 pmat%zmix = m_zero
523
524 imat = 1
525 if (lloc /= 0) then
526 kb_p => epot%proj(iatom)%kb_p(1, 1)
527
528 do ic = 1, kb_p%n_c
529 pmat%zmix(ic, ic, 1:2) = kb_p%e(ic)
530 do ip = 1, pmat%npoints
531 pmat%zprojectors(ip, ic) = kb_p%p(ip, ic)
532 end do
533 pmat%scal(ic) = mesh%volume_element
534 end do
535 imat = kb_p%n_c + 1
536 nullify(kb_p)
537 end if
538
539 do ll = 1, lmax
540 if (ll == lloc) cycle
541 do mm = -ll, ll
542 rkb_p => epot%proj(iatom)%rkb_p(ll, mm)
543
544 ! See rkb_projector.F90 for understanding the indices
545 do ic = 0, rkb_p%n_c/2-1
546 pmat%zmix(imat + ic*2, imat + ic*2, 1) = rkb_p%f(ic*2+1, 1, 1)
547 pmat%zmix(imat + ic*2, imat + ic*2, 2) = rkb_p%f(ic*2+1, 2, 2)
548
549 pmat%zmix(imat + ic*2+1, imat + ic*2+1, 1) = rkb_p%f(ic*2+2, 1, 1)
550 pmat%zmix(imat + ic*2+1, imat + ic*2+1, 2) = rkb_p%f(ic*2+2, 2, 2)
551
552 if (mm < ll) then
553 pmat%zmix(imat + ic*2+rkb_p%n_c, imat + ic*2, 4) = rkb_p%f(ic*2+1, 2, 1)
554 pmat%zmix(imat + ic*2+1+rkb_p%n_c, imat + ic*2+1, 4) = rkb_p%f(ic*2+2, 2, 1)
555 end if
556
557 if (-mm < ll) then
558 pmat%zmix(imat + ic*2-rkb_p%n_c, imat + ic*2, 3) = rkb_p%f(ic*2+1, 1, 2)
559 pmat%zmix(imat + ic*2+1-rkb_p%n_c, imat + ic*2+1, 3) = rkb_p%f(ic*2+2, 1, 2)
560 end if
561 end do
562
563 do ic = 1, rkb_p%n_c
564 call lalg_copy(pmat%npoints, rkb_p%ket(:, ic, 1, 1), pmat%zprojectors(:, imat))
565 pmat%scal(imat) = mesh%volume_element
566 imat = imat + 1
567 end do
568 end do
569
570 nullify(rkb_p)
571 end do
572
573 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
574
575 else
576 cycle
577 end if
578
579 pmat%map => epot%proj(iatom)%sphere%map
580 pmat%position => epot%proj(iatom)%sphere%rel_x
581
582 pmat%regions = epot%proj(iatom)%sphere%regions
583
584 this%full_projection_size = this%full_projection_size + pmat%nprojs
585
586 end do
587 end do
588
589 if (mesh%parallel_in_domains) then
590 call mesh%mpi_grp%allreduce_inplace(this%projector_self_overlap, 1, mpi_logical, mpi_lor)
591 end if
592
593 safe_deallocate_a(order)
594 safe_deallocate_a(head)
595
596 this%total_points = 0
597 this%max_npoints = 0
598 this%max_nprojs = 0
599 do imat = 1, this%nprojector_matrices
600 pmat => this%projector_matrices(imat)
601
602 this%max_npoints = max(this%max_npoints, pmat%npoints)
603 this%max_nprojs = max(this%max_nprojs, pmat%nprojs)
604 this%total_points = this%total_points + pmat%npoints
605 end do
606
607 if (accel_is_enabled()) then
610 end if
611
613
615
616 ! ----------------------------------------------------------------------------------
618 subroutine nonlocal_pseudopotential_accel_rebuild(this, space, mesh)
619 class(nonlocal_pseudopotential_t), target, intent(inout) :: this
620 class(space_t), intent(in) :: space
621 class(mesh_t), intent(in) :: mesh
622
624
625 if (.not. accel_is_enabled()) then
627 return
628 end if
629
631
632 if (.not. allocated(this%projector_matrices) .or. this%nprojector_matrices <= 0) then
634 return
635 end if
636
639
642
643 ! ----------------------------------------------------------------------------------
646 class(nonlocal_pseudopotential_t), intent(inout) :: this
647
649
650 call accel_detach_buffer(this%buff_offsets)
651 call accel_detach_buffer(this%buff_matrices)
652 call accel_detach_buffer(this%buff_maps)
653 call accel_detach_buffer(this%buff_scals)
654 call accel_detach_buffer(this%buff_position)
655 call accel_detach_buffer(this%buff_pos)
656 call accel_detach_buffer(this%buff_invmap)
657 call accel_detach_buffer(this%buff_projector_phases)
658 call accel_detach_buffer(this%buff_mix)
659
662
663 ! ----------------------------------------------------------------------------------
665 subroutine nonlocal_pseudopotential_build_accel_buffers(this, space, mesh)
666 class(nonlocal_pseudopotential_t), target, intent(inout) :: this
667 class(space_t), intent(in) :: space
668 class(mesh_t), intent(in) :: mesh
669
670 integer, parameter :: OFFSET_SIZE = 6
671 integer, parameter :: POINTS = 1, projs = 2, matrix = 3, map = 4, scal = 5, mix = 6
672 integer :: imat, matrix_size, scal_size
673 integer :: ip, is, ii, ipos, mix_offset
674 integer, allocatable :: cnt(:), invmap(:, :), invmap2(:), pos(:)
675 integer, allocatable :: offsets(:, :)
676 type(projector_matrix_t), pointer :: pmat
677
679
680 assert(allocated(this%projector_matrices))
681 assert(this%nprojector_matrices > 0)
682
683 safe_allocate(offsets(1:offset_size, 1:this%nprojector_matrices))
684 safe_allocate(cnt(1:mesh%np))
685
686 cnt = 0
687
688 matrix_size = 0
689 this%total_points = 0
690 scal_size = 0
691 this%max_npoints = 0
692 this%max_nprojs = 0
693 mix_offset = 0
694 do imat = 1, this%nprojector_matrices
695 pmat => this%projector_matrices(imat)
696
697 this%max_npoints = max(this%max_npoints, pmat%npoints)
698 this%max_nprojs = max(this%max_nprojs, pmat%nprojs)
699
700 offsets(points, imat) = pmat%npoints
701 offsets(projs, imat) = pmat%nprojs
702
703 offsets(matrix, imat) = matrix_size
704 matrix_size = matrix_size + pmat%npoints*pmat%nprojs
705
706 offsets(map, imat) = this%total_points
707 this%total_points = this%total_points + pmat%npoints
708
709 offsets(scal, imat) = scal_size
710 scal_size = scal_size + pmat%nprojs
711
712 offsets(mix, imat) = mix_offset
713 if (allocated(pmat%dmix)) then
714 mix_offset = mix_offset + pmat%nprojs**2
715 else if (allocated(pmat%zmix)) then
716 mix_offset = mix_offset + 4*pmat%nprojs**2
717 else
718 offsets(mix, imat) = -1
719 end if
720
721 do is = 1, pmat%npoints
722 ip = pmat%map(is)
723 cnt(ip) = cnt(ip) + 1
724 end do
725 end do
726
727 safe_allocate(invmap(1:max(maxval(cnt), 1), 1:mesh%np))
728 safe_allocate(invmap2(1:max(maxval(cnt)*mesh%np, 1)))
729 safe_allocate(pos(1:mesh%np + 1))
730
731 cnt = 0
732 ii = 0
733 do imat = 1, this%nprojector_matrices
734 pmat => this%projector_matrices(imat)
735 do is = 1, pmat%npoints
736 ip = pmat%map(is)
737 cnt(ip) = cnt(ip) + 1
738 invmap(cnt(ip), ip) = ii
739 ii = ii + 1
740 end do
741 end do
742
743 ipos = 0
744 pos(1) = 0
745 do ip = 1, mesh%np
746 do ii = 1, cnt(ip)
747 ipos = ipos + 1
748 invmap2(ipos) = invmap(ii, ip)
749 end do
750 pos(ip + 1) = ipos
751 end do
752
753 if (this%projector_matrices(1)%is_cmplx) then
754 call accel_create_buffer(this%buff_matrices, accel_mem_read_only, type_cmplx, matrix_size)
755 else
756 call accel_create_buffer(this%buff_matrices, accel_mem_read_only, type_float, matrix_size)
757 end if
758 call accel_create_buffer(this%buff_maps, accel_mem_read_only, type_integer, this%total_points)
759 call accel_create_buffer(this%buff_position, accel_mem_read_only, type_float, 3*this%total_points)
760 call accel_create_buffer(this%buff_scals, accel_mem_read_only, type_float, scal_size)
761
762 if (mix_offset > 0) then
763 if (allocated(this%projector_matrices(1)%zmix)) then
764 call accel_create_buffer(this%buff_mix, accel_mem_read_only, type_cmplx, mix_offset)
765 else
766 call accel_create_buffer(this%buff_mix, accel_mem_read_only, type_float, mix_offset)
767 end if
768 end if
769
770 do imat = 1, this%nprojector_matrices
771 pmat => this%projector_matrices(imat)
772 if (pmat%npoints > 0) then
773 if (pmat%is_cmplx) then
774 call accel_write_buffer(this%buff_matrices, pmat%npoints, pmat%nprojs, pmat%zprojectors, &
775 offset = offsets(matrix, imat))
776 else
777 call accel_write_buffer(this%buff_matrices, pmat%npoints, pmat%nprojs, pmat%dprojectors, &
778 offset = offsets(matrix, imat))
779 end if
780 call accel_write_buffer(this%buff_maps, pmat%npoints, pmat%map, offset = offsets(map, imat))
781 call accel_write_buffer(this%buff_position, space%dim, pmat%npoints, pmat%position, &
782 offset = 3*offsets(map, imat))
783 end if
784 call accel_write_buffer(this%buff_scals, pmat%nprojs, pmat%scal, offset = offsets(scal, imat))
785 if (offsets(mix, imat) /= -1) then
786 if (allocated(pmat%zmix)) then
787 call accel_write_buffer(this%buff_mix, pmat%nprojs, pmat%nprojs, 4, pmat%zmix, offset = offsets(mix, imat))
788 else
789 call accel_write_buffer(this%buff_mix, pmat%nprojs, pmat%nprojs, pmat%dmix, offset = offsets(mix, imat))
790 end if
791 end if
792 end do
793
794 call accel_create_buffer(this%buff_offsets, accel_mem_read_only, type_integer, offset_size*this%nprojector_matrices)
795 call accel_write_buffer(this%buff_offsets, offset_size, this%nprojector_matrices, offsets)
796
797 call accel_create_buffer(this%buff_pos, accel_mem_read_only, type_integer, mesh%np + 1)
798 call accel_write_buffer(this%buff_pos, mesh%np + 1, pos)
799
800 call accel_create_buffer(this%buff_invmap, accel_mem_read_only, type_integer, ipos)
801 call accel_write_buffer(this%buff_invmap, ipos, invmap2)
802
803 safe_deallocate_a(offsets)
804 safe_deallocate_a(cnt)
805 safe_deallocate_a(invmap)
806 safe_deallocate_a(invmap2)
807 safe_deallocate_a(pos)
808
811
812 ! ----------------------------------------------------------------------------------
815 class(nonlocal_pseudopotential_t), intent(inout) :: this
816
817 integer :: ik, imat, iphase, nphase, offset, npoints
818
820
821 if (.not. allocated(this%projector_phases)) then
823 return
824 end if
825
826 nphase = size(this%projector_phases, 2)
827 this%nphase = nphase
828
829 call accel_create_buffer(this%buff_projector_phases, accel_mem_read_only, type_cmplx, &
830 this%total_points*nphase*size(this%projector_phases, 4))
831
832 offset = 0
833 do ik = lbound(this%projector_phases, 4), ubound(this%projector_phases, 4)
834 do imat = 1, this%nprojector_matrices
835 npoints = this%projector_matrices(imat)%npoints
836 do iphase = 1, nphase
837 if (npoints > 0) then
838 call accel_write_buffer(this%buff_projector_phases, npoints, this%projector_phases(1:, iphase, imat, ik), &
839 offset = offset, async=.true.)
840 end if
841 offset = offset + npoints
842 end do
843 end do
844 end do
845 call accel_finish()
846
849
850 ! ----------------------------------------------------------------------------------
852 !
853 logical pure function nonlocal_pseudopotential_self_overlap(this) result(projector_self_overlap)
854 class(nonlocal_pseudopotential_t), intent(in) :: this
855
856 projector_self_overlap = this%projector_self_overlap
858
859#include "undef.F90"
860#include "real.F90"
861#include "nonlocal_pseudopotential_inc.F90"
862
863#include "undef.F90"
864#include "complex.F90"
865#include "nonlocal_pseudopotential_inc.F90"
866
868
869!! Local Variables:
870!! mode: f90
871!! coding: utf-8
872!! End:
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:188
double sqrt(double __x) __attribute__((__nothrow__
subroutine, public accel_free_buffer(this, async)
Definition: accel.F90:1005
subroutine, public accel_finish()
Definition: accel.F90:1098
subroutine, public accel_detach_buffer(this)
Clear a buffer handle without freeing device memory.
Definition: accel.F90:1049
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
integer, parameter, public accel_mem_read_only
Definition: accel.F90:185
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:120
integer, parameter, public spin_orbit
Definition: epot.F90:169
integer, parameter, public fully_relativistic_zora
Definition: epot.F90:169
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_half
Definition: global.F90:197
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
subroutine nonlocal_pseudopotential_destroy_proj(this)
Destroy the data of nonlocal_pseudopotential_t.
subroutine nonlocal_pseudopotential_detach_accel_buffers(this)
Clear copied accelerator handles so they do not alias source buffers.
subroutine dnonlocal_pseudopotential_force(this, mesh, st, spiral_bnd, iqn, ndim, psi1b, psi2b, force)
calculate contribution to forces, from non-local potentials
subroutine znonlocal_pseudopotential_position_commutator(this, mesh, std, spiral_bnd, psib, commpsib, async)
apply the commutator between the non-local potential and the position to the wave functions.
subroutine znonlocal_pseudopotential_force(this, mesh, st, spiral_bnd, iqn, ndim, psi1b, psi2b, force)
calculate contribution to forces, from non-local potentials
subroutine dnonlocal_pseudopotential_start(this, mesh, std, spiral_bnd, psib, projection, async)
Start application of non-local potentials (stored in the Hamiltonian) to the wave functions.
subroutine dnonlocal_pseudopotential_finish(this, mesh, spiral_bnd, std, projection, vpsib)
finish the application of non-local potentials.
subroutine, public nonlocal_pseudopotential_accel_rebuild(this, space, mesh)
Rebuild accelerator buffers after an intrinsic copy.
subroutine dnonlocal_pseudopotential_position_commutator(this, mesh, std, spiral_bnd, psib, commpsib, async)
apply the commutator between the non-local potential and the position to the wave functions.
subroutine, public nonlocal_pseudopotential_rebind_projectors(this, epot)
Rebind projector matrix pointers (map, position) to a target epot.
subroutine nonlocal_pseudopotential_build_accel_buffers(this, space, mesh)
Build accelerator buffers for projectors from host-side projector matrices.
subroutine znonlocal_pseudopotential_r_vnlocal(this, mesh, std, spiral_bnd, psib, commpsib)
Accumulates to commpsib the result of x V_{nl} | psib >
logical pure function nonlocal_pseudopotential_self_overlap(this)
Returns .true. if the Hamiltonian contains projectors, which overlap with themself.
subroutine nonlocal_pseudopotential_init(this)
initialize the nonlocal_pseudopotential_t object
subroutine dnonlocal_pseudopotential_r_vnlocal(this, mesh, std, spiral_bnd, psib, commpsib)
Accumulates to commpsib the result of x V_{nl} | psib >
subroutine znonlocal_pseudopotential_start(this, mesh, std, spiral_bnd, psib, projection, async)
Start application of non-local potentials (stored in the Hamiltonian) to the wave functions.
subroutine nonlocal_pseudopotential_build_projector_phase_accel_buffer(this)
Rebuild projector phase accelerator buffer from host-side projector phases.
subroutine nonlocal_pseudopotential_build_proj(this, space, mesh, epot)
build the projectors for the application of pseudo-potentials
subroutine znonlocal_pseudopotential_finish(this, mesh, spiral_bnd, std, projection, vpsib)
finish the application of non-local potentials.
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 projector_matrix_deallocate(this)
subroutine, public projector_matrix_allocate(this, nprojs, sphere, has_mix_matrix, is_cmplx)
logical elemental function, public projector_is(p, type)
Definition: projector.F90:210
logical elemental function, public projector_is_null(p)
Definition: projector.F90:203
Definition: ps.F90:116
integer, parameter, public proj_hgh
Definition: ps.F90:171
integer, parameter, public proj_rkb
Definition: ps.F90:171
integer, parameter, public proj_none
Definition: ps.F90:171
integer, parameter, public proj_kb
Definition: ps.F90:171
This module handles spin dimensions of the states and the k-point distribution.
logical function, public submesh_overlap(sm1, sm2, space)
Definition: submesh.F90:702
type(type_t), public type_float
Definition: types.F90:135
type(type_t), public type_cmplx
Definition: types.F90:136
type(type_t), public type_integer
Definition: types.F90:137
Describes mesh distribution to nodes.
Definition: mesh.F90:187
Class for projections of wave functions.
A set of projectors defined on a submesh.
The rkb_projector data type holds the KB projectors build with total angular momentum eigenfunctions....
int true(void)