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
23 use batch_oct_m
25 use blas_oct_m
26 use comm_oct_m
27 use debug_oct_m
29 use epot_oct_m
30 use global_oct_m
34 use math_oct_m
35 use mesh_oct_m
37 use mpi_oct_m
41 use ps_oct_m
43 use space_oct_m
47 use types_oct_m
49
50 implicit none
51
52 private
53
54 public :: &
57
61 private
62 type(projector_matrix_t), allocatable, public :: projector_matrices(:)
63 integer, public :: nprojector_matrices
64 logical, public :: apply_projector_matrices
65 logical, public :: has_non_local_potential
66 integer :: full_projection_size
67 integer, public :: max_npoints
68 integer, public :: total_points
69 integer :: max_nprojs
70 logical :: projector_mix
71 complex(real64), allocatable, public :: projector_phases(:, :, :, :)
72 integer, allocatable, public :: projector_to_atom(:)
73 integer :: nregions
74 integer, allocatable :: regions(:)
75 integer, public :: nphase
80 type(accel_mem_t) :: buff_offsets
81 type(accel_mem_t) :: buff_matrices
82 type(accel_mem_t) :: buff_maps
83 type(accel_mem_t) :: buff_scals
84 type(accel_mem_t) :: buff_position
85 type(accel_mem_t) :: buff_pos
86 type(accel_mem_t) :: buff_invmap
87 type(accel_mem_t), public :: buff_projector_phases
88 type(accel_mem_t) :: buff_mix
89 logical :: projector_self_overlap
90 real(real64), pointer, public :: spin(:,:,:) => null()
91 contains
92
93 procedure :: init => nonlocal_pseudopotential_init
94
95 procedure :: build => nonlocal_pseudopotential_build_proj
96
98
99 procedure :: has_self_overlap => nonlocal_pseudopotential_self_overlap
100
101 procedure :: dstart => dnonlocal_pseudopotential_start
102
103 procedure :: zstart => znonlocal_pseudopotential_start
104
105 procedure :: dfinish => dnonlocal_pseudopotential_finish
106
107 procedure :: zfinish => znonlocal_pseudopotential_finish
108
109 procedure :: dforce => dnonlocal_pseudopotential_force
110
111 procedure :: zforce => znonlocal_pseudopotential_force
112
113 procedure :: dposition_commutator => dnonlocal_pseudopotential_position_commutator
115 procedure :: zposition_commutator => znonlocal_pseudopotential_position_commutator
116
117 procedure :: dr_vn_local => dnonlocal_pseudopotential_r_vnlocal
118
119 procedure :: zr_vn_local => znonlocal_pseudopotential_r_vnlocal
120
122
124 !
125 type projection_t
126 private
127 real(real64), allocatable :: dprojection(:, :)
128 complex(real64), allocatable :: zprojection(:, :)
129 type(accel_mem_t) :: buff_projection
130 type(accel_mem_t) :: buff_spin_to_phase
131 end type projection_t
132
133contains
134
135 ! ---------------------------------------------------------
138 subroutine nonlocal_pseudopotential_init(this)
139 class(nonlocal_pseudopotential_t), intent(inout) :: this
140
142
143 this%apply_projector_matrices = .false.
144 this%has_non_local_potential = .false.
145 this%nprojector_matrices = 0
146
147 this%projector_self_overlap = .false.
148
150 end subroutine nonlocal_pseudopotential_init
151
152
153 !--------------------------------------------------------
156 class(nonlocal_pseudopotential_t), intent(inout) :: this
158 integer :: iproj
162 if (allocated(this%projector_matrices)) then
165 call accel_release_buffer(this%buff_offsets)
166 call accel_release_buffer(this%buff_matrices)
167 call accel_release_buffer(this%buff_maps)
168 call accel_release_buffer(this%buff_scals)
169 call accel_release_buffer(this%buff_position)
170 call accel_release_buffer(this%buff_pos)
171 call accel_release_buffer(this%buff_invmap)
172 if (this%projector_mix) call accel_release_buffer(this%buff_mix)
173 if (allocated(this%projector_phases)) call accel_release_buffer(this%buff_projector_phases)
174 end if
176 do iproj = 1, this%nprojector_matrices
177 call projector_matrix_deallocate(this%projector_matrices(iproj))
178 end do
179 safe_deallocate_a(this%regions)
180 safe_deallocate_a(this%projector_matrices)
181 safe_deallocate_a(this%projector_phases)
182 safe_deallocate_a(this%projector_to_atom)
183 end if
184
187
188 !-----------------------------------------------------------------
195 subroutine nonlocal_pseudopotential_build_proj(this, space, mesh, epot)
196 class(nonlocal_pseudopotential_t), target, intent(inout) :: this
197 class(space_t), intent(in) :: space
198 class(mesh_t), intent(in) :: mesh
199 type(epot_t), target, intent(in) :: epot
201 integer :: iatom, iproj, ll, lmax, lloc, mm, ic, jc
202 integer :: nmat, imat, ip, iorder
203 integer :: nregion, jatom, katom, iregion
204 integer, allocatable :: order(:), head(:), region_count(:)
205 logical, allocatable :: atom_counted(:)
206 logical :: overlap
207 type(projector_matrix_t), pointer :: pmat
208 type(kb_projector_t), pointer :: kb_p
209 type(rkb_projector_t), pointer :: rkb_p
210 type(hgh_projector_t), pointer :: hgh_p
211
213
214 call profiling_in("ATOM_COLORING")
215
216 ! this is most likely a very inefficient algorithm, O(natom**2) or
217 ! O(natom**3), probably it should be replaced by something better.
219 safe_allocate(order(1:epot%natoms)) ! order(iregion) = ?
220 safe_allocate(head(1:epot%natoms + 1)) ! head(iregion) points to the first atom in region iregion
221 safe_allocate(region_count(1:epot%natoms)) ! region_count(iregion): number of atoms in that region
222 safe_allocate(atom_counted(1:epot%natoms))
224 this%projector_self_overlap = .false.
225 atom_counted = .false.
226 order = -1
227
228 head(1) = 1
229 nregion = 0
230 do
231 nregion = nregion + 1
232 assert(nregion <= epot%natoms)
233
234 region_count(nregion) = 0
235
236 do iatom = 1, epot%natoms
237 if (atom_counted(iatom)) cycle
238
239 overlap = .false.
240
241 if (.not. projector_is(epot%proj(iatom), proj_none)) then
242 assert(associated(epot%proj(iatom)%sphere%mesh))
243 do jatom = 1, region_count(nregion)
244 katom = order(head(nregion) + jatom - 1)
245 if (projector_is(epot%proj(katom), proj_none)) cycle
246 overlap = submesh_overlap(epot%proj(iatom)%sphere, epot%proj(katom)%sphere, space)
247 if (overlap) exit
248 end do
249 end if
250
251 if (.not. overlap) then
252 ! iatom did not overlap with any previously counted atoms:
253 ! iatom will be added to the current region
254 region_count(nregion) = region_count(nregion) + 1
255 order(head(nregion) - 1 + region_count(nregion)) = iatom
256 atom_counted(iatom) = .true.
257 end if
258
259 end do
260
261 head(nregion + 1) = head(nregion) + region_count(nregion)
262
263 if (all(atom_counted)) exit
264 end do
265
266 safe_deallocate_a(atom_counted)
267 safe_deallocate_a(region_count)
268
269 if (debug%info) then
270 call messages_write('The atoms can be separated in ')
271 call messages_write(nregion)
272 call messages_write(' non-overlapping groups.')
273 call messages_info()
274 end if
275
276 do iregion = 1, nregion
277 do iatom = head(iregion), head(iregion + 1) - 1
278 if (.not. projector_is(epot%proj(order(iatom)), proj_kb)) cycle
279 do jatom = head(iregion), iatom - 1
280 if (.not. projector_is(epot%proj(order(jatom)), proj_kb)) cycle
281 assert(.not. submesh_overlap(epot%proj(order(iatom))%sphere, epot%proj(order(jatom))%sphere, space))
282 end do
283 end do
284 end do
285
286 call profiling_out("ATOM_COLORING")
287
288 ! deallocate previous projectors
289 call this%end()
290
291 ! count projectors
292 this%nprojector_matrices = 0
293 this%apply_projector_matrices = .false.
294 this%has_non_local_potential = .false.
295 this%nregions = nregion
296
297 !We determine if we have only local potential or not.
298 do iorder = 1, epot%natoms
299 iatom = order(iorder)
300
301 if (.not. projector_is_null(epot%proj(iatom))) then
302 this%has_non_local_potential = .true.
303 exit
304 end if
305 end do
306
307 do iorder = 1, epot%natoms
308 iatom = order(iorder)
309
310 if (.not. projector_is_null(epot%proj(iatom))) then
311 this%nprojector_matrices = this%nprojector_matrices + 1
312 this%apply_projector_matrices = .true.
313 !The HGH pseudopotentials are now supporting the SOC
314 if (epot%reltype /= norel .and. projector_is(epot%proj(iatom), proj_kb)) then
315 this%apply_projector_matrices = .false.
316 exit
317 end if
318 end if
319 end do
320
321 if (mesh%use_curvilinear) this%apply_projector_matrices = .false.
322
323 if (.not. this%apply_projector_matrices) then
324 safe_deallocate_a(order)
325 safe_deallocate_a(head)
326
328 return
329 end if
330
331
332 safe_allocate(this%projector_matrices(1:this%nprojector_matrices))
333 safe_allocate(this%regions(1:this%nprojector_matrices + 1))
334 safe_allocate(this%projector_to_atom(1:epot%natoms))
335
336 this%full_projection_size = 0
337 this%regions(this%nregions + 1) = this%nprojector_matrices + 1
338
339 this%projector_mix = .false.
340
341 iproj = 0
342 do iregion = 1, this%nregions
343 this%regions(iregion) = iproj + 1
344 do iorder = head(iregion), head(iregion + 1) - 1
345
346 iatom = order(iorder)
347
348 if (projector_is(epot%proj(iatom), proj_none)) cycle
349
350 iproj = iproj + 1
351
352 pmat => this%projector_matrices(iproj)
353
354 this%projector_to_atom(iproj) = iatom
355
356 lmax = epot%proj(iatom)%lmax
357 lloc = epot%proj(iatom)%lloc
358
359 if (projector_is(epot%proj(iatom), proj_kb)) then
360
361 ! count the number of projectors for this matrix
362 nmat = 0
363 do ll = 0, lmax
364 if (ll == lloc) cycle
365 do mm = -ll, ll
366 nmat = nmat + epot%proj(iatom)%kb_p(ll, mm)%n_c
367 end do
368 end do
369
370 call projector_matrix_allocate(pmat, nmat, epot%proj(iatom)%sphere, has_mix_matrix = .false.)
371
372 ! generate the matrix
373 pmat%dprojectors = m_zero
374 imat = 1
375 do ll = 0, lmax
376 if (ll == lloc) cycle
377 do mm = -ll, ll
378 kb_p => epot%proj(iatom)%kb_p(ll, mm)
379 do ic = 1, kb_p%n_c
380 do ip = 1, pmat%npoints
381 pmat%dprojectors(ip, imat) = kb_p%p(ip, ic)
382 end do
383 pmat%scal(imat) = kb_p%e(ic)*mesh%vol_pp(1)
384 imat = imat + 1
385 end do
386 end do
387 end do
388
389 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
390
391 else if (projector_is(epot%proj(iatom), proj_hgh)) then
392
393 this%projector_mix = .true.
394
395 ! count the number of projectors for this matrix
396 nmat = 0
397 do ll = 0, lmax
398 if (ll == lloc) cycle
399 do mm = -ll, ll
400 nmat = nmat + 3
401 end do
402 end do
403
404 call projector_matrix_allocate(pmat, nmat, epot%proj(iatom)%sphere, &
405 has_mix_matrix = .true., is_cmplx = (epot%reltype == spin_orbit))
406
407 ! generate the matrix
408 if (epot%reltype == spin_orbit) then
409 pmat%zprojectors = m_zero
410 pmat%zmix = m_zero
411 else
412 pmat%dprojectors = m_zero
413 pmat%dmix = m_zero
414 end if
415
416 imat = 1
417 do ll = 0, lmax
418 if (ll == lloc) cycle
419 do mm = -ll, ll
420 hgh_p => epot%proj(iatom)%hgh_p(ll, mm)
421
422 ! HGH pseudos mix different components, so we need to
423 ! generate a matrix that mixes the projections
424 if (epot%reltype == spin_orbit) then
425 do ic = 1, 3
426 do jc = 1, 3
427 pmat%zmix(imat - 1 + ic, imat - 1 + jc, 1) = hgh_p%h(ic, jc) + m_half*mm*hgh_p%k(ic, jc)
428 pmat%zmix(imat - 1 + ic, imat - 1 + jc, 2) = hgh_p%h(ic, jc) - m_half*mm*hgh_p%k(ic, jc)
429
430 if (mm < ll) then
431 pmat%zmix(imat - 1 + ic, imat + 3 - 1 + jc, 3) = m_half*hgh_p%k(ic, jc) * &
432 sqrt(real(ll*(ll+1)-mm*(mm+1), real64))
433 end if
434
435 if (-mm < ll) then
436 pmat%zmix(imat - 1 + ic, imat - 3 - 1 + jc, 4) = m_half*hgh_p%k(ic, jc) * &
437 sqrt(real(ll*(ll+1)-mm*(mm-1), real64))
438 end if
439 end do
440 end do
441 else
442 do ic = 1, 3
443 do jc = 1, 3
444 pmat%dmix(imat - 1 + ic, imat - 1 + jc) = hgh_p%h(ic, jc)
445 end do
446 end do
447 end if
448
449 do ic = 1, 3
450 if (epot%reltype == spin_orbit) then
451 do ip = 1, pmat%npoints
452 pmat%zprojectors(ip, imat) = hgh_p%zp(ip, ic)
453 end do
454 else
455 do ip = 1, pmat%npoints
456 pmat%dprojectors(ip, imat) = hgh_p%dp(ip, ic)
457 end do
458 end if
459 pmat%scal(imat) = mesh%volume_element
460 imat = imat + 1
461 end do
462
463 end do
464 end do
465
466 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
467
468 else if (projector_is(epot%proj(iatom), proj_rkb)) then
469 assert(epot%reltype == spin_orbit)
470
471 this%projector_mix = .true.
472
473 ! count the number of projectors for this matrix
474 nmat = 0
475 if (lloc /= 0) nmat = nmat + epot%proj(iatom)%kb_p(1, 1)%n_c
476
477 do ll = 1, lmax
478 if (ll == lloc) cycle
479 do mm = -ll, ll
480 nmat = nmat + epot%proj(iatom)%rkb_p(ll, mm)%n_c
481 end do
482 end do
483
484 call projector_matrix_allocate(pmat, nmat, epot%proj(iatom)%sphere, &
485 has_mix_matrix = .true., is_cmplx = .true.)
486
487 pmat%zprojectors = m_zero
488 pmat%zmix = m_zero
489
490 imat = 1
491 if (lloc /= 0) then
492 kb_p => epot%proj(iatom)%kb_p(1, 1)
493
494 do ic = 1, kb_p%n_c
495 pmat%zmix(ic, ic, 1:2) = kb_p%e(ic)
496 do ip = 1, pmat%npoints
497 pmat%zprojectors(ip, ic) = kb_p%p(ip, ic)
498 end do
499 pmat%scal(ic) = mesh%volume_element
500 end do
501 imat = kb_p%n_c + 1
502 nullify(kb_p)
503 end if
504
505 do ll = 1, lmax
506 if (ll == lloc) cycle
507 do mm = -ll, ll
508 rkb_p => epot%proj(iatom)%rkb_p(ll, mm)
509
510 ! See rkb_projector.F90 for understanding the indices
511 do ic = 0, rkb_p%n_c/2-1
512 pmat%zmix(imat + ic*2, imat + ic*2, 1) = rkb_p%f(ic*2+1, 1, 1)
513 pmat%zmix(imat + ic*2, imat + ic*2, 2) = rkb_p%f(ic*2+1, 2, 2)
514
515 pmat%zmix(imat + ic*2+1, imat + ic*2+1, 1) = rkb_p%f(ic*2+2, 1, 1)
516 pmat%zmix(imat + ic*2+1, imat + ic*2+1, 2) = rkb_p%f(ic*2+2, 2, 2)
517
518 if (mm < ll) then
519 pmat%zmix(imat + ic*2+rkb_p%n_c, imat + ic*2, 4) = rkb_p%f(ic*2+1, 2, 1)
520 pmat%zmix(imat + ic*2+1+rkb_p%n_c, imat + ic*2+1, 4) = rkb_p%f(ic*2+2, 2, 1)
521 end if
522
523 if (-mm < ll) then
524 pmat%zmix(imat + ic*2-rkb_p%n_c, imat + ic*2, 3) = rkb_p%f(ic*2+1, 1, 2)
525 pmat%zmix(imat + ic*2+1-rkb_p%n_c, imat + ic*2+1, 3) = rkb_p%f(ic*2+2, 1, 2)
526 end if
527 end do
528
529 do ic = 1, rkb_p%n_c
530 do ip = 1, pmat%npoints
531 pmat%zprojectors(ip, imat) = rkb_p%ket(ip, ic, 1, 1)
532 end do
533 pmat%scal(imat) = mesh%volume_element
534 imat = imat + 1
535 end do
536 end do
537
538 nullify(rkb_p)
539 end do
540
541 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
542
543 else
544 cycle
545 end if
546
547 pmat%map => epot%proj(iatom)%sphere%map
548 pmat%position => epot%proj(iatom)%sphere%rel_x
549
550 pmat%regions = epot%proj(iatom)%sphere%regions
551
552 this%full_projection_size = this%full_projection_size + pmat%nprojs
553
554 end do
555 end do
556
557 if (mesh%parallel_in_domains) then
558 call mesh%mpi_grp%allreduce_inplace(this%projector_self_overlap, 1, mpi_logical, mpi_lor)
559 end if
560
561 safe_deallocate_a(order)
562 safe_deallocate_a(head)
563
564 this%total_points = 0
565 this%max_npoints = 0
566 this%max_nprojs = 0
567 do imat = 1, this%nprojector_matrices
568 pmat => this%projector_matrices(imat)
569
570 this%max_npoints = max(this%max_npoints, pmat%npoints)
571 this%max_nprojs = max(this%max_nprojs, pmat%nprojs)
572 this%total_points = this%total_points + pmat%npoints
573 end do
574
575 if (accel_is_enabled()) call build_accel()
576
578
579 contains
580
581 subroutine build_accel()
582
583 integer :: matrix_size, scal_size
584 integer, allocatable :: cnt(:), invmap(:, :), invmap2(:), pos(:)
585 integer, allocatable :: offsets(:, :)
586 integer, parameter :: OFFSET_SIZE = 6 ! also defined in share/opencl/projectors.cl
587 integer, parameter :: POINTS = 1, projs = 2, matrix = 3, map = 4, scal = 5, mix = 6 ! update OFFSET_SIZE
588 integer :: ip, is, ii, ipos, mix_offset
589
591
592 safe_allocate(offsets(1:offset_size, 1:this%nprojector_matrices))
593 safe_allocate(cnt(1:mesh%np))
594
595 cnt = 0
596
597 ! Here we construct the offsets for accessing various arrays within the GPU kernels.
598 ! The offset(:,:) array contains a number of sizes and offsets, describing how to address the arrays.
599 ! This allows to transfer all these numbers to the GPU in one memory transfer.
600 !
601 ! For each projection matrix (addressed by imap), we have:
602 !
603 ! offset(POINTS, imap) : number of points of the sphere imap
604 ! offset(PROJS, imap) : number of projectors for imap
605 ! offset(MATRIX, imap) : address offset: cumulative of pmat%npoints * pmat%nprojs
606 ! offset(MAP, imap) : address offset: cumulative of pmat%npoints for each imap
607 ! offset(SCAL, imap) : address_offset: cumulative of pmat%nprojs
608 ! offset(MIX, imap) : address_offset: cumulative of pmat%nprojs**2 or 4*pmat%nprojs**2 for complex mixing
609
610 ! first we count
611 matrix_size = 0
612 this%total_points = 0
613 scal_size = 0
614 this%max_npoints = 0
615 this%max_nprojs = 0
616 mix_offset = 0
617 do imat = 1, this%nprojector_matrices
618 pmat => this%projector_matrices(imat)
619
620 this%max_npoints = max(this%max_npoints, pmat%npoints)
621 this%max_nprojs = max(this%max_nprojs, pmat%nprojs)
622
623 offsets(points, imat) = pmat%npoints
624 offsets(projs, imat) = pmat%nprojs
625
626 offsets(matrix, imat) = matrix_size
627 matrix_size = matrix_size + pmat%npoints*pmat%nprojs
628
629 offsets(map, imat) = this%total_points
630 this%total_points = this%total_points + pmat%npoints
631
632 offsets(scal, imat) = scal_size
633 scal_size = scal_size + pmat%nprojs
634
635 offsets(mix, imat) = mix_offset
636 if (allocated(pmat%dmix)) then
637 mix_offset = mix_offset + pmat%nprojs**2
638 else if (allocated(pmat%zmix)) then
639 mix_offset = mix_offset + 4 * pmat%nprojs**2
640 else
641 offsets(mix, imat) = -1
642 end if
643
644 do is = 1, pmat%npoints
645 ip = pmat%map(is)
646 cnt(ip) = cnt(ip) + 1
647 end do
648 end do
649
650 safe_allocate(invmap(1:max(maxval(cnt), 1), 1:mesh%np))
651 safe_allocate(invmap2(1:max(maxval(cnt)*mesh%np, 1)))
652 safe_allocate(pos(1:mesh%np + 1))
653
654 cnt = 0
655 ii = 0
656 do imat = 1, this%nprojector_matrices
657 pmat => this%projector_matrices(imat)
658 do is = 1, pmat%npoints
659 ip = pmat%map(is)
660 cnt(ip) = cnt(ip) + 1
661 invmap(cnt(ip), ip) = ii
662 ii = ii + 1
663 end do
664 end do
665
666 ipos = 0
667 pos(1) = 0
668 do ip = 1, mesh%np
669 do ii = 1, cnt(ip)
670 ipos = ipos + 1
671 invmap2(ipos) = invmap(ii, ip)
672 end do
673 pos(ip + 1) = ipos
674 end do
675
676 ! allocate
677 if (this%projector_matrices(1)%is_cmplx) then
678 call accel_create_buffer(this%buff_matrices, accel_mem_read_only, type_cmplx, matrix_size)
679 else
680 call accel_create_buffer(this%buff_matrices, accel_mem_read_only, type_float, matrix_size)
681 end if
682 call accel_create_buffer(this%buff_maps, accel_mem_read_only, type_integer, this%total_points)
683 call accel_create_buffer(this%buff_position, accel_mem_read_only, type_float, 3*this%total_points)
684 call accel_create_buffer(this%buff_scals, accel_mem_read_only, type_float, scal_size)
685
686 if (mix_offset > 0) then
687 if (allocated(this%projector_matrices(1)%zmix)) then
688 call accel_create_buffer(this%buff_mix, accel_mem_read_only, type_cmplx, mix_offset)
689 else
690 call accel_create_buffer(this%buff_mix, accel_mem_read_only, type_float, mix_offset)
691 end if
692 end if
693
694 ! now copy
695 do imat = 1, this%nprojector_matrices
696 pmat => this%projector_matrices(imat)
697 ! in parallel some spheres might not have points
698 if (pmat%npoints > 0) then
699 if (pmat%is_cmplx) then
700 call accel_write_buffer(this%buff_matrices, pmat%nprojs*pmat%npoints, pmat%zprojectors, offset = offsets(matrix, imat))
701 else
702 call accel_write_buffer(this%buff_matrices, pmat%nprojs*pmat%npoints, pmat%dprojectors, offset = offsets(matrix, imat))
703 end if
704 call accel_write_buffer(this%buff_maps, pmat%npoints, pmat%map, offset = offsets(map, imat))
705 call accel_write_buffer(this%buff_position, 3*pmat%npoints, pmat%position, offset = 3*offsets(map, imat))
706 end if
707 call accel_write_buffer(this%buff_scals, pmat%nprojs, pmat%scal, offset = offsets(scal, imat))
708 if (offsets(mix, imat) /= -1) then
709 if (allocated(pmat%zmix)) then
710 call accel_write_buffer(this%buff_mix, 4*pmat%nprojs**2, pmat%zmix, offset = offsets(mix, imat))
711 else
712 call accel_write_buffer(this%buff_mix, pmat%nprojs**2, pmat%dmix, offset = offsets(mix, imat))
713 end if
714 end if
715 end do
716
717 ! write the offsets
718 call accel_create_buffer(this%buff_offsets, accel_mem_read_only, type_integer, offset_size*this%nprojector_matrices)
719 call accel_write_buffer(this%buff_offsets, offset_size*this%nprojector_matrices, offsets)
720
721 ! the inverse map
722 call accel_create_buffer(this%buff_pos, accel_mem_read_only, type_integer, mesh%np + 1)
723 call accel_write_buffer(this%buff_pos, mesh%np + 1, pos)
724
725 call accel_create_buffer(this%buff_invmap, accel_mem_read_only, type_integer, ipos)
726 call accel_write_buffer(this%buff_invmap, ipos, invmap2)
727
728 safe_deallocate_a(offsets)
729 safe_deallocate_a(cnt)
730 safe_deallocate_a(invmap)
731 safe_deallocate_a(invmap2)
732 safe_deallocate_a(pos)
733
735 end subroutine build_accel
736
738
739 ! ----------------------------------------------------------------------------------
741 !
742 logical pure function nonlocal_pseudopotential_self_overlap(this) result(projector_self_overlap)
743 class(nonlocal_pseudopotential_t), intent(in) :: this
744
745 projector_self_overlap = this%projector_self_overlap
747
748#include "undef.F90"
749#include "real.F90"
750#include "nonlocal_pseudopotential_inc.F90"
751
752#include "undef.F90"
753#include "complex.F90"
754#include "nonlocal_pseudopotential_inc.F90"
755
757
758!! Local Variables:
759!! mode: f90
760!! coding: utf-8
761!! End:
double sqrt(double __x) __attribute__((__nothrow__
subroutine, public accel_release_buffer(this)
Definition: accel.F90:1246
pure logical function, public accel_is_enabled()
Definition: accel.F90:400
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
type(accel_kernel_t), pointer head
Definition: accel.F90:394
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:118
type(debug_t), save, public debug
Definition: debug.F90:154
integer, parameter, public norel
Definition: epot.F90:166
integer, parameter, public spin_orbit
Definition: epot.F90:166
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_half
Definition: global.F90:193
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
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
subroutine nonlocal_pseudopotential_destroy_proj(this)
Destroy the data of nonlocal_pseudopotential_t.
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 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 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_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:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
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:114
integer, parameter, public proj_hgh
Definition: ps.F90:165
integer, parameter, public proj_rkb
Definition: ps.F90:165
integer, parameter, public proj_none
Definition: ps.F90:165
integer, parameter, public proj_kb
Definition: ps.F90:165
This module handles spin dimensions of the states and the k-point distribution.
logical function, public submesh_overlap(sm1, sm2, space)
Definition: submesh.F90:778
type(type_t), public type_float
Definition: types.F90:133
type(type_t), public type_cmplx
Definition: types.F90:134
type(type_t), public type_integer
Definition: types.F90:135
subroutine build_accel()
Describes mesh distribution to nodes.
Definition: mesh.F90:186
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)