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 end type projection_t
138
139contains
140
141 ! ---------------------------------------------------------
144 subroutine nonlocal_pseudopotential_init(this)
145 class(nonlocal_pseudopotential_t), intent(inout) :: this
146
148
149 this%apply_projector_matrices = .false.
150 this%has_non_local_potential = .false.
151 this%nprojector_matrices = 0
152
153 this%projector_self_overlap = .false.
154
156 end subroutine nonlocal_pseudopotential_init
157
158 ! ---------------------------------------------------------
164 class(nonlocal_pseudopotential_t), target, intent(inout) :: this
165 type(epot_t), target, intent(in) :: epot
167 integer :: imat, iatom
168 type(projector_matrix_t), pointer :: pmat
172 if (.not. allocated(this%projector_matrices) .or. this%nprojector_matrices == 0) then
174 return
175 end if
176
177 assert(allocated(this%projector_to_atom))
178 assert(size(this%projector_matrices) >= this%nprojector_matrices)
179 assert(size(this%projector_to_atom) >= this%nprojector_matrices)
181 do imat = 1, this%nprojector_matrices
182 pmat => this%projector_matrices(imat)
183 iatom = this%projector_to_atom(imat)
185 assert(iatom >= 1 .and. iatom <= epot%natoms)
186 assert(allocated(epot%proj))
187 assert(allocated(epot%proj(iatom)%sphere%map))
188 assert(allocated(epot%proj(iatom)%sphere%rel_x))
190 pmat%map => epot%proj(iatom)%sphere%map
191 pmat%position => epot%proj(iatom)%sphere%rel_x
192 end do
193
197
198 !--------------------------------------------------------
201 class(nonlocal_pseudopotential_t), intent(inout) :: this
203 integer :: iproj
207 if (allocated(this%projector_matrices)) then
209 if (accel_is_enabled()) then
210 call accel_free_buffer(this%buff_offsets)
211 call accel_free_buffer(this%buff_matrices)
212 call accel_free_buffer(this%buff_maps)
213 call accel_free_buffer(this%buff_scals)
214 call accel_free_buffer(this%buff_position)
215 call accel_free_buffer(this%buff_pos)
216 call accel_free_buffer(this%buff_invmap)
217 if (this%projector_mix) call accel_free_buffer(this%buff_mix)
218 if (allocated(this%projector_phases)) call accel_free_buffer(this%buff_projector_phases)
219 end if
220
221 do iproj = 1, this%nprojector_matrices
222 call projector_matrix_deallocate(this%projector_matrices(iproj))
223 end do
224 safe_deallocate_a(this%regions)
225 safe_deallocate_a(this%projector_matrices)
226 safe_deallocate_a(this%projector_phases)
227 safe_deallocate_a(this%projector_to_atom)
228 end if
232
233 !-----------------------------------------------------------------
240 subroutine nonlocal_pseudopotential_build_proj(this, space, mesh, epot)
241 class(nonlocal_pseudopotential_t), target, intent(inout) :: this
242 class(space_t), intent(in) :: space
243 class(mesh_t), intent(in) :: mesh
244 type(epot_t), target, intent(in) :: epot
245
246 integer :: iatom, iproj, ll, lmax, lloc, mm, ic, jc
247 integer :: nmat, imat, ip, iorder
248 integer :: nregion, jatom, katom, iregion
249 integer, allocatable :: order(:), head(:), region_count(:)
250 logical, allocatable :: atom_counted(:)
251 logical :: overlap
252 type(projector_matrix_t), pointer :: pmat
253 type(kb_projector_t), pointer :: kb_p
254 type(rkb_projector_t), pointer :: rkb_p
255 type(hgh_projector_t), pointer :: hgh_p
256
259 call profiling_in("ATOM_COLORING")
260
261 ! this is most likely a very inefficient algorithm, O(natom**2) or
262 ! O(natom**3), probably it should be replaced by something better.
263
264 safe_allocate(order(1:epot%natoms)) ! order(iregion) = ?
265 safe_allocate(head(1:epot%natoms + 1)) ! head(iregion) points to the first atom in region iregion
266 safe_allocate(region_count(1:epot%natoms)) ! region_count(iregion): number of atoms in that region
267 safe_allocate(atom_counted(1:epot%natoms))
268
269 this%projector_self_overlap = .false.
270 atom_counted = .false.
271 order = -1
272
273 head(1) = 1
274 nregion = 0
275 do
276 nregion = nregion + 1
277 assert(nregion <= epot%natoms)
278
279 region_count(nregion) = 0
280
281 do iatom = 1, epot%natoms
282 if (atom_counted(iatom)) cycle
283
284 overlap = .false.
285
286 if (.not. projector_is(epot%proj(iatom), proj_none)) then
287 assert(associated(epot%proj(iatom)%sphere%mesh))
288 do jatom = 1, region_count(nregion)
289 katom = order(head(nregion) + jatom - 1)
290 if (projector_is(epot%proj(katom), proj_none)) cycle
291 overlap = submesh_overlap(epot%proj(iatom)%sphere, epot%proj(katom)%sphere, space)
292 if (overlap) exit
293 end do
294 end if
296 if (.not. overlap) then
297 ! iatom did not overlap with any previously counted atoms:
298 ! iatom will be added to the current region
299 region_count(nregion) = region_count(nregion) + 1
300 order(head(nregion) - 1 + region_count(nregion)) = iatom
301 atom_counted(iatom) = .true.
302 end if
303
304 end do
305
306 head(nregion + 1) = head(nregion) + region_count(nregion)
307
308 if (all(atom_counted)) exit
309 end do
310
311 safe_deallocate_a(atom_counted)
312 safe_deallocate_a(region_count)
313
314 call messages_write('The atoms can be separated in ')
315 call messages_write(nregion)
316 call messages_write(' non-overlapping groups.')
317 call messages_info(debug_only=.true.)
318
319 do iregion = 1, nregion
320 do iatom = head(iregion), head(iregion + 1) - 1
321 if (.not. projector_is(epot%proj(order(iatom)), proj_kb)) cycle
322 do jatom = head(iregion), iatom - 1
323 if (.not. projector_is(epot%proj(order(jatom)), proj_kb)) cycle
324 assert(.not. submesh_overlap(epot%proj(order(iatom))%sphere, epot%proj(order(jatom))%sphere, space))
325 end do
326 end do
327 end do
328
329 call profiling_out("ATOM_COLORING")
330
331 ! deallocate previous projectors
332 call this%end()
333
334 ! count projectors
335 this%nprojector_matrices = 0
336 this%apply_projector_matrices = .false.
337 this%has_non_local_potential = .false.
338 this%nregions = nregion
339
340 !We determine if we have only local potential or not.
341 do iorder = 1, epot%natoms
342 iatom = order(iorder)
343
344 if (.not. projector_is_null(epot%proj(iatom))) then
345 this%has_non_local_potential = .true.
346 exit
347 end if
348 end do
349
350 do iorder = 1, epot%natoms
351 iatom = order(iorder)
352
353 if (.not. projector_is_null(epot%proj(iatom))) then
354 this%nprojector_matrices = this%nprojector_matrices + 1
355 this%apply_projector_matrices = .true.
356 end if
357 end do
358
359 ! This is currently the only not supported case
360 if (mesh%use_curvilinear) this%apply_projector_matrices = .false.
361
362 if (.not. this%apply_projector_matrices) then
363 safe_deallocate_a(order)
364 safe_deallocate_a(head)
365
367 return
368 end if
369
370
371 safe_allocate(this%projector_matrices(1:this%nprojector_matrices))
372 safe_allocate(this%regions(1:this%nprojector_matrices + 1))
373 safe_allocate(this%projector_to_atom(1:epot%natoms))
374
375 this%full_projection_size = 0
376 this%regions(this%nregions + 1) = this%nprojector_matrices + 1
377
378 this%projector_mix = .false.
379
380 iproj = 0
381 do iregion = 1, this%nregions
382 this%regions(iregion) = iproj + 1
383 do iorder = head(iregion), head(iregion + 1) - 1
384
385 iatom = order(iorder)
386
387 if (projector_is(epot%proj(iatom), proj_none)) cycle
388
389 iproj = iproj + 1
390
391 pmat => this%projector_matrices(iproj)
392
393 this%projector_to_atom(iproj) = iatom
394
395 lmax = epot%proj(iatom)%lmax
396 lloc = epot%proj(iatom)%lloc
397
398 if (projector_is(epot%proj(iatom), proj_kb)) then
399
400 ! count the number of projectors for this matrix
401 nmat = 0
402 do ll = 0, lmax
403 if (ll == lloc) cycle
404 do mm = -ll, ll
405 nmat = nmat + epot%proj(iatom)%kb_p(ll, mm)%n_c
406 end do
407 end do
408
409 call projector_matrix_allocate(pmat, nmat, epot%proj(iatom)%sphere, has_mix_matrix = .false.)
410
411 ! generate the matrix
412 pmat%dprojectors = m_zero
413 imat = 1
414 do ll = 0, lmax
415 if (ll == lloc) cycle
416 do mm = -ll, ll
417 kb_p => epot%proj(iatom)%kb_p(ll, mm)
418 do ic = 1, kb_p%n_c
419 call lalg_copy(pmat%npoints, kb_p%p(:, ic), pmat%dprojectors(:, imat))
420 pmat%scal(imat) = kb_p%e(ic)*mesh%vol_pp(1)
421 imat = imat + 1
422 end do
423 end do
424 end do
425
426 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
427
428 else if (projector_is(epot%proj(iatom), proj_hgh)) then
429
430 this%projector_mix = .true.
431
432 ! count the number of projectors for this matrix
433 nmat = 0
434 do ll = 0, lmax
435 if (ll == lloc) cycle
436 do mm = -ll, ll
437 nmat = nmat + 3
438 end do
439 end do
440
441 call projector_matrix_allocate(pmat, nmat, epot%proj(iatom)%sphere, &
442 has_mix_matrix = .true., is_cmplx = (epot%reltype == spin_orbit))
443
444 ! generate the matrix
445 if (epot%reltype == spin_orbit) then
446 pmat%zprojectors = m_zero
447 pmat%zmix = m_zero
448 else
449 pmat%dprojectors = m_zero
450 pmat%dmix = m_zero
451 end if
452
453 imat = 1
454 do ll = 0, lmax
455 if (ll == lloc) cycle
456 do mm = -ll, ll
457 hgh_p => epot%proj(iatom)%hgh_p(ll, mm)
458
459 ! HGH pseudos mix different components, so we need to
460 ! generate a matrix that mixes the projections
461 if (epot%reltype == spin_orbit .or. epot%reltype == fully_relativistic_zora) then
462 do ic = 1, 3
463 do jc = 1, 3
464 pmat%zmix(imat - 1 + ic, imat - 1 + jc, 1) = hgh_p%h(ic, jc) + m_half*mm*hgh_p%k(ic, jc)
465 pmat%zmix(imat - 1 + ic, imat - 1 + jc, 2) = hgh_p%h(ic, jc) - m_half*mm*hgh_p%k(ic, jc)
466
467 if (mm < ll) then
468 pmat%zmix(imat - 1 + ic, imat + 3 - 1 + jc, 3) = m_half*hgh_p%k(ic, jc) * &
469 sqrt(real(ll*(ll+1)-mm*(mm+1), real64))
470 end if
471
472 if (-mm < ll) then
473 pmat%zmix(imat - 1 + ic, imat - 3 - 1 + jc, 4) = m_half*hgh_p%k(ic, jc) * &
474 sqrt(real(ll*(ll+1)-mm*(mm-1), real64))
475 end if
476 end do
477 end do
478 else
479 do ic = 1, 3
480 do jc = 1, 3
481 pmat%dmix(imat - 1 + ic, imat - 1 + jc) = hgh_p%h(ic, jc)
482 end do
483 end do
484 end if
485
486 do ic = 1, 3
487 if (epot%reltype == spin_orbit .or. epot%reltype == fully_relativistic_zora) then
488 call lalg_copy(pmat%npoints, hgh_p%zp(:, ic), pmat%zprojectors(:, imat))
489 else
490 call lalg_copy(pmat%npoints, hgh_p%dp(:, ic), pmat%dprojectors(:, imat))
491 end if
492 pmat%scal(imat) = mesh%volume_element
493 imat = imat + 1
494 end do
495
496 end do
497 end do
498
499 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
500
501 else if (projector_is(epot%proj(iatom), proj_rkb)) then
502 assert(epot%reltype == spin_orbit)
503
504 this%projector_mix = .true.
505
506 ! count the number of projectors for this matrix
507 nmat = 0
508 if (lloc /= 0) nmat = nmat + epot%proj(iatom)%kb_p(1, 1)%n_c
509
510 do ll = 1, lmax
511 if (ll == lloc) cycle
512 do mm = -ll, ll
513 nmat = nmat + epot%proj(iatom)%rkb_p(ll, mm)%n_c
514 end do
515 end do
516
517 call projector_matrix_allocate(pmat, nmat, epot%proj(iatom)%sphere, &
518 has_mix_matrix = .true., is_cmplx = .true.)
519
520 pmat%zprojectors = m_zero
521 pmat%zmix = m_zero
522
523 imat = 1
524 if (lloc /= 0) then
525 kb_p => epot%proj(iatom)%kb_p(1, 1)
526
527 do ic = 1, kb_p%n_c
528 pmat%zmix(ic, ic, 1:2) = kb_p%e(ic)
529 do ip = 1, pmat%npoints
530 pmat%zprojectors(ip, ic) = kb_p%p(ip, ic)
531 end do
532 pmat%scal(ic) = mesh%volume_element
533 end do
534 imat = kb_p%n_c + 1
535 nullify(kb_p)
536 end if
537
538 do ll = 1, lmax
539 if (ll == lloc) cycle
540 do mm = -ll, ll
541 rkb_p => epot%proj(iatom)%rkb_p(ll, mm)
542
543 ! See rkb_projector.F90 for understanding the indices
544 do ic = 0, rkb_p%n_c/2-1
545 pmat%zmix(imat + ic*2, imat + ic*2, 1) = rkb_p%f(ic*2+1, 1, 1)
546 pmat%zmix(imat + ic*2, imat + ic*2, 2) = rkb_p%f(ic*2+1, 2, 2)
547
548 pmat%zmix(imat + ic*2+1, imat + ic*2+1, 1) = rkb_p%f(ic*2+2, 1, 1)
549 pmat%zmix(imat + ic*2+1, imat + ic*2+1, 2) = rkb_p%f(ic*2+2, 2, 2)
550
551 if (mm < ll) then
552 pmat%zmix(imat + ic*2+rkb_p%n_c, imat + ic*2, 4) = rkb_p%f(ic*2+1, 2, 1)
553 pmat%zmix(imat + ic*2+1+rkb_p%n_c, imat + ic*2+1, 4) = rkb_p%f(ic*2+2, 2, 1)
554 end if
555
556 if (-mm < ll) then
557 pmat%zmix(imat + ic*2-rkb_p%n_c, imat + ic*2, 3) = rkb_p%f(ic*2+1, 1, 2)
558 pmat%zmix(imat + ic*2+1-rkb_p%n_c, imat + ic*2+1, 3) = rkb_p%f(ic*2+2, 1, 2)
559 end if
560 end do
561
562 do ic = 1, rkb_p%n_c
563 call lalg_copy(pmat%npoints, rkb_p%ket(:, ic, 1, 1), pmat%zprojectors(:, imat))
564 pmat%scal(imat) = mesh%volume_element
565 imat = imat + 1
566 end do
567 end do
568
569 nullify(rkb_p)
570 end do
571
572 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
573
574 else
575 cycle
576 end if
577
578 pmat%map => epot%proj(iatom)%sphere%map
579 pmat%position => epot%proj(iatom)%sphere%rel_x
580
581 pmat%regions = epot%proj(iatom)%sphere%regions
582
583 this%full_projection_size = this%full_projection_size + pmat%nprojs
584
585 end do
586 end do
587
588 if (mesh%parallel_in_domains) then
589 call mesh%mpi_grp%allreduce_inplace(this%projector_self_overlap, 1, mpi_logical, mpi_lor)
590 end if
591
592 safe_deallocate_a(order)
593 safe_deallocate_a(head)
594
595 this%total_points = 0
596 this%max_npoints = 0
597 this%max_nprojs = 0
598 do imat = 1, this%nprojector_matrices
599 pmat => this%projector_matrices(imat)
600
601 this%max_npoints = max(this%max_npoints, pmat%npoints)
602 this%max_nprojs = max(this%max_nprojs, pmat%nprojs)
603 this%total_points = this%total_points + pmat%npoints
604 end do
605
606 if (accel_is_enabled()) then
609 end if
610
612
614
615 ! ----------------------------------------------------------------------------------
617 subroutine nonlocal_pseudopotential_accel_rebuild(this, space, mesh)
618 class(nonlocal_pseudopotential_t), target, intent(inout) :: this
619 class(space_t), intent(in) :: space
620 class(mesh_t), intent(in) :: mesh
621
623
624 if (.not. accel_is_enabled()) then
626 return
627 end if
628
630
631 if (.not. allocated(this%projector_matrices) .or. this%nprojector_matrices <= 0) then
633 return
634 end if
635
638
641
642 ! ----------------------------------------------------------------------------------
645 class(nonlocal_pseudopotential_t), intent(inout) :: this
646
648
649 call accel_detach_buffer(this%buff_offsets)
650 call accel_detach_buffer(this%buff_matrices)
651 call accel_detach_buffer(this%buff_maps)
652 call accel_detach_buffer(this%buff_scals)
653 call accel_detach_buffer(this%buff_position)
654 call accel_detach_buffer(this%buff_pos)
655 call accel_detach_buffer(this%buff_invmap)
656 call accel_detach_buffer(this%buff_projector_phases)
657 call accel_detach_buffer(this%buff_mix)
658
661
662 ! ----------------------------------------------------------------------------------
664 subroutine nonlocal_pseudopotential_build_accel_buffers(this, space, mesh)
665 class(nonlocal_pseudopotential_t), target, intent(inout) :: this
666 class(space_t), intent(in) :: space
667 class(mesh_t), intent(in) :: mesh
668
669 integer, parameter :: OFFSET_SIZE = 6
670 integer, parameter :: POINTS = 1, projs = 2, matrix = 3, map = 4, scal = 5, mix = 6
671 integer :: imat, matrix_size, scal_size
672 integer :: ip, is, ii, ipos, mix_offset
673 integer, allocatable :: cnt(:), invmap(:, :), invmap2(:), pos(:)
674 integer, allocatable :: offsets(:, :)
675 type(projector_matrix_t), pointer :: pmat
676
678
679 assert(allocated(this%projector_matrices))
680 assert(this%nprojector_matrices > 0)
681
682 safe_allocate(offsets(1:offset_size, 1:this%nprojector_matrices))
683 safe_allocate(cnt(1:mesh%np))
684
685 cnt = 0
686
687 matrix_size = 0
688 this%total_points = 0
689 scal_size = 0
690 this%max_npoints = 0
691 this%max_nprojs = 0
692 mix_offset = 0
693 do imat = 1, this%nprojector_matrices
694 pmat => this%projector_matrices(imat)
695
696 this%max_npoints = max(this%max_npoints, pmat%npoints)
697 this%max_nprojs = max(this%max_nprojs, pmat%nprojs)
698
699 offsets(points, imat) = pmat%npoints
700 offsets(projs, imat) = pmat%nprojs
701
702 offsets(matrix, imat) = matrix_size
703 matrix_size = matrix_size + pmat%npoints*pmat%nprojs
704
705 offsets(map, imat) = this%total_points
706 this%total_points = this%total_points + pmat%npoints
707
708 offsets(scal, imat) = scal_size
709 scal_size = scal_size + pmat%nprojs
710
711 offsets(mix, imat) = mix_offset
712 if (allocated(pmat%dmix)) then
713 mix_offset = mix_offset + pmat%nprojs**2
714 else if (allocated(pmat%zmix)) then
715 mix_offset = mix_offset + 4*pmat%nprojs**2
716 else
717 offsets(mix, imat) = -1
718 end if
719
720 do is = 1, pmat%npoints
721 ip = pmat%map(is)
722 cnt(ip) = cnt(ip) + 1
723 end do
724 end do
725
726 safe_allocate(invmap(1:max(maxval(cnt), 1), 1:mesh%np))
727 safe_allocate(invmap2(1:max(maxval(cnt)*mesh%np, 1)))
728 safe_allocate(pos(1:mesh%np + 1))
729
730 cnt = 0
731 ii = 0
732 do imat = 1, this%nprojector_matrices
733 pmat => this%projector_matrices(imat)
734 do is = 1, pmat%npoints
735 ip = pmat%map(is)
736 cnt(ip) = cnt(ip) + 1
737 invmap(cnt(ip), ip) = ii
738 ii = ii + 1
739 end do
740 end do
741
742 ipos = 0
743 pos(1) = 0
744 do ip = 1, mesh%np
745 do ii = 1, cnt(ip)
746 ipos = ipos + 1
747 invmap2(ipos) = invmap(ii, ip)
748 end do
749 pos(ip + 1) = ipos
750 end do
751
752 if (this%projector_matrices(1)%is_cmplx) then
753 call accel_create_buffer(this%buff_matrices, accel_mem_read_only, type_cmplx, matrix_size)
754 else
755 call accel_create_buffer(this%buff_matrices, accel_mem_read_only, type_float, matrix_size)
756 end if
757 call accel_create_buffer(this%buff_maps, accel_mem_read_only, type_integer, this%total_points)
758 call accel_create_buffer(this%buff_position, accel_mem_read_only, type_float, 3*this%total_points)
759 call accel_create_buffer(this%buff_scals, accel_mem_read_only, type_float, scal_size)
760
761 if (mix_offset > 0) then
762 if (allocated(this%projector_matrices(1)%zmix)) then
763 call accel_create_buffer(this%buff_mix, accel_mem_read_only, type_cmplx, mix_offset)
764 else
765 call accel_create_buffer(this%buff_mix, accel_mem_read_only, type_float, mix_offset)
766 end if
767 end if
768
769 do imat = 1, this%nprojector_matrices
770 pmat => this%projector_matrices(imat)
771 if (pmat%npoints > 0) then
772 if (pmat%is_cmplx) then
773 call accel_write_buffer(this%buff_matrices, pmat%npoints, pmat%nprojs, pmat%zprojectors, &
774 offset = offsets(matrix, imat))
775 else
776 call accel_write_buffer(this%buff_matrices, pmat%npoints, pmat%nprojs, pmat%dprojectors, &
777 offset = offsets(matrix, imat))
778 end if
779 call accel_write_buffer(this%buff_maps, pmat%npoints, pmat%map, offset = offsets(map, imat))
780 call accel_write_buffer(this%buff_position, space%dim, pmat%npoints, pmat%position, &
781 offset = 3*offsets(map, imat))
782 end if
783 call accel_write_buffer(this%buff_scals, pmat%nprojs, pmat%scal, offset = offsets(scal, imat))
784 if (offsets(mix, imat) /= -1) then
785 if (allocated(pmat%zmix)) then
786 call accel_write_buffer(this%buff_mix, pmat%nprojs, pmat%nprojs, 4, pmat%zmix, offset = offsets(mix, imat))
787 else
788 call accel_write_buffer(this%buff_mix, pmat%nprojs, pmat%nprojs, pmat%dmix, offset = offsets(mix, imat))
789 end if
790 end if
791 end do
792
793 call accel_create_buffer(this%buff_offsets, accel_mem_read_only, type_integer, offset_size*this%nprojector_matrices)
794 call accel_write_buffer(this%buff_offsets, offset_size, this%nprojector_matrices, offsets)
795
796 call accel_create_buffer(this%buff_pos, accel_mem_read_only, type_integer, mesh%np + 1)
797 call accel_write_buffer(this%buff_pos, mesh%np + 1, pos)
798
799 call accel_create_buffer(this%buff_invmap, accel_mem_read_only, type_integer, ipos)
800 call accel_write_buffer(this%buff_invmap, ipos, invmap2)
801
802 safe_deallocate_a(offsets)
803 safe_deallocate_a(cnt)
804 safe_deallocate_a(invmap)
805 safe_deallocate_a(invmap2)
806 safe_deallocate_a(pos)
807
810
811 ! ----------------------------------------------------------------------------------
814 class(nonlocal_pseudopotential_t), intent(inout) :: this
815
816 integer :: ik, imat, iphase, nphase, offset, npoints
817
819
820 if (.not. allocated(this%projector_phases)) then
822 return
823 end if
824
825 nphase = size(this%projector_phases, 2)
826 this%nphase = nphase
827
828 call accel_create_buffer(this%buff_projector_phases, accel_mem_read_only, type_cmplx, &
829 this%total_points*nphase*size(this%projector_phases, 4))
830
831 offset = 0
832 do ik = lbound(this%projector_phases, 4), ubound(this%projector_phases, 4)
833 do imat = 1, this%nprojector_matrices
834 npoints = this%projector_matrices(imat)%npoints
835 do iphase = 1, nphase
836 if (npoints > 0) then
837 call accel_write_buffer(this%buff_projector_phases, npoints, this%projector_phases(1:, iphase, imat, ik), &
838 offset = offset, async=.true.)
839 end if
840 offset = offset + npoints
841 end do
842 end do
843 end do
844 call accel_finish()
845
848
849 ! ----------------------------------------------------------------------------------
851 !
852 logical pure function nonlocal_pseudopotential_self_overlap(this) result(projector_self_overlap)
853 class(nonlocal_pseudopotential_t), intent(in) :: this
854
855 projector_self_overlap = this%projector_self_overlap
857
858#include "undef.F90"
859#include "real.F90"
860#include "nonlocal_pseudopotential_inc.F90"
861
862#include "undef.F90"
863#include "complex.F90"
864#include "nonlocal_pseudopotential_inc.F90"
865
867
868!! Local Variables:
869!! mode: f90
870!! coding: utf-8
871!! 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:594
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:774
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)