13 USE iso_c_binding,
ONLY: &
14 c_associated, c_bool, c_char, c_double, c_funloc, c_funptr, c_int, c_loc, c_null_ptr, c_ptr
19#include "../base/base_uses.f90"
25 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'grid_api'
80 TYPE(C_PTR) :: c_ptr = c_null_ptr
85 TYPE(C_PTR) :: c_ptr = c_null_ptr
112 lb_max, zetb, lb_min, &
113 ra, rab, scale, pab, o1, o2, &
115 ga_gb_function, radius, &
116 use_subpatch, subpatch_pattern)
118 INTEGER,
INTENT(IN) :: la_max
119 REAL(kind=
dp),
INTENT(IN) :: zeta
120 INTEGER,
INTENT(IN) :: la_min, lb_max
121 REAL(kind=
dp),
INTENT(IN) :: zetb
122 INTEGER,
INTENT(IN) :: lb_min
123 REAL(kind=
dp),
DIMENSION(3),
INTENT(IN),
TARGET :: ra, rab
124 REAL(kind=
dp),
INTENT(IN) :: scale
125 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: pab
126 INTEGER,
INTENT(IN) :: o1, o2
128 INTEGER,
INTENT(IN) :: ga_gb_function
129 REAL(kind=
dp),
INTENT(IN) :: radius
130 LOGICAL,
OPTIONAL :: use_subpatch
131 INTEGER,
INTENT(IN),
OPTIONAL :: subpatch_pattern
133 INTEGER :: border_mask
134 INTEGER,
DIMENSION(3),
TARGET :: border_width, npts_global, npts_local, &
136 LOGICAL(KIND=C_BOOL) :: orthorhombic
137 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: grid
139 SUBROUTINE grid_cpu_collocate_pgf_product_c(orthorhombic, &
141 la_max, la_min, lb_max, lb_min, &
142 zeta, zetb, rscale, dh, dh_inv, ra, rab, &
143 npts_global, npts_local, shift_local, border_width, &
144 radius, o1, o2, n1, n2, pab, &
146 BIND(C, name="grid_cpu_collocate_pgf_product")
147 IMPORT :: c_ptr, c_int, c_double, c_bool
148 LOGICAL(KIND=C_BOOL),
VALUE :: orthorhombic
149 INTEGER(KIND=C_INT),
VALUE :: border_mask
150 INTEGER(KIND=C_INT),
VALUE :: func
151 INTEGER(KIND=C_INT),
VALUE :: la_max
152 INTEGER(KIND=C_INT),
VALUE :: la_min
153 INTEGER(KIND=C_INT),
VALUE :: lb_max
154 INTEGER(KIND=C_INT),
VALUE :: lb_min
155 REAL(kind=c_double),
VALUE :: zeta
156 REAL(kind=c_double),
VALUE :: zetb
157 REAL(kind=c_double),
VALUE :: rscale
158 TYPE(c_ptr),
VALUE :: dh
159 TYPE(c_ptr),
VALUE :: dh_inv
160 TYPE(c_ptr),
VALUE :: ra
161 TYPE(c_ptr),
VALUE :: rab
162 TYPE(c_ptr),
VALUE :: npts_global
163 TYPE(c_ptr),
VALUE :: npts_local
164 TYPE(c_ptr),
VALUE :: shift_local
165 TYPE(c_ptr),
VALUE :: border_width
166 REAL(kind=c_double),
VALUE :: radius
167 INTEGER(KIND=C_INT),
VALUE :: o1
168 INTEGER(KIND=C_INT),
VALUE :: o2
169 INTEGER(KIND=C_INT),
VALUE :: n1
170 INTEGER(KIND=C_INT),
VALUE :: n2
171 TYPE(c_ptr),
VALUE :: pab
172 TYPE(c_ptr),
VALUE :: grid
173 END SUBROUTINE grid_cpu_collocate_pgf_product_c
177 IF (
PRESENT(use_subpatch))
THEN
178 IF (use_subpatch)
THEN
179 cpassert(
PRESENT(subpatch_pattern))
180 border_mask = iand(63, not(subpatch_pattern))
184 orthorhombic =
LOGICAL(rsgrid%desc%orthorhombic, c_bool)
186 cpassert(lbound(pab, 1) == 1)
187 cpassert(lbound(pab, 2) == 1)
189 CALL get_rsgrid_properties(rsgrid, npts_global=npts_global, &
190 npts_local=npts_local, &
191 shift_local=shift_local, &
192 border_width=border_width)
194 grid(1:, 1:, 1:) => rsgrid%r(:, :, :)
197 cpassert(is_contiguous(rsgrid%desc%dh))
198 cpassert(is_contiguous(rsgrid%desc%dh_inv))
199 cpassert(is_contiguous(ra))
200 cpassert(is_contiguous(rab))
201 cpassert(is_contiguous(npts_global))
202 cpassert(is_contiguous(npts_local))
203 cpassert(is_contiguous(shift_local))
204 cpassert(is_contiguous(border_width))
205 cpassert(is_contiguous(pab))
206 cpassert(is_contiguous(grid))
211 CALL grid_cpu_collocate_pgf_product_c(orthorhombic=orthorhombic, &
212 border_mask=border_mask, &
213 func=ga_gb_function, &
221 dh=c_loc(rsgrid%desc%dh(1, 1)), &
222 dh_inv=c_loc(rsgrid%desc%dh_inv(1, 1)), &
225 npts_global=c_loc(npts_global(1)), &
226 npts_local=c_loc(npts_local(1)), &
227 shift_local=c_loc(shift_local(1)), &
228 border_width=c_loc(border_width(1)), &
234 pab=c_loc(pab(1, 1)), &
235 grid=c_loc(grid(1, 1, 1)))
269 lb_max, zetb, lb_min, &
273 calculate_forces, force_a, force_b, &
275 use_virial, my_virial_a, &
276 my_virial_b, hdab, hadb, a_hdab, use_subpatch, subpatch_pattern)
278 INTEGER,
INTENT(IN) :: la_max
279 REAL(kind=
dp),
INTENT(IN) :: zeta
280 INTEGER,
INTENT(IN) :: la_min, lb_max
281 REAL(kind=
dp),
INTENT(IN) :: zetb
282 INTEGER,
INTENT(IN) :: lb_min
283 REAL(kind=
dp),
DIMENSION(3),
INTENT(IN),
TARGET :: ra, rab
285 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: hab
286 REAL(kind=
dp),
DIMENSION(:, :),
OPTIONAL,
POINTER :: pab
287 INTEGER,
INTENT(IN) :: o1, o2
288 REAL(kind=
dp),
INTENT(IN) :: radius
289 LOGICAL,
INTENT(IN) :: calculate_forces
290 REAL(kind=
dp),
DIMENSION(3),
INTENT(INOUT), &
291 OPTIONAL :: force_a, force_b
292 LOGICAL,
INTENT(IN),
OPTIONAL :: compute_tau, use_virial
293 REAL(kind=
dp),
DIMENSION(3, 3),
OPTIONAL :: my_virial_a, my_virial_b
294 REAL(kind=
dp),
DIMENSION(:, :, :),
OPTIONAL, &
295 POINTER :: hdab, hadb
296 REAL(kind=
dp),
DIMENSION(:, :, :, :),
OPTIONAL, &
298 LOGICAL,
OPTIONAL :: use_subpatch
299 INTEGER,
INTENT(IN),
OPTIONAL :: subpatch_pattern
301 INTEGER :: border_mask
302 INTEGER,
DIMENSION(3),
TARGET :: border_width, npts_global, npts_local, &
304 LOGICAL :: my_use_virial
305 LOGICAL(KIND=C_BOOL) :: my_compute_tau, orthorhombic
306 REAL(kind=
dp),
DIMENSION(3, 2),
TARGET :: forces
307 REAL(kind=
dp),
DIMENSION(3, 3, 2),
TARGET :: virials
308 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: grid
309 TYPE(c_ptr) :: a_hdab_cptr, forces_cptr, hadb_cptr, &
310 hdab_cptr, pab_cptr, virials_cptr
312 SUBROUTINE grid_cpu_integrate_pgf_product_c(orthorhombic, compute_tau, &
314 la_max, la_min, lb_max, lb_min, &
315 zeta, zetb, dh, dh_inv, ra, rab, &
316 npts_global, npts_local, shift_local, border_width, &
317 radius, o1, o2, n1, n2, grid, hab, pab, &
318 forces, virials, hdab, hadb, a_hdab) &
319 BIND(C, name="grid_cpu_integrate_pgf_product")
320 IMPORT :: c_ptr, c_int, c_double, c_bool
321 LOGICAL(KIND=C_BOOL),
VALUE :: orthorhombic
322 LOGICAL(KIND=C_BOOL),
VALUE :: compute_tau
323 INTEGER(KIND=C_INT),
VALUE :: border_mask
324 INTEGER(KIND=C_INT),
VALUE :: la_max
325 INTEGER(KIND=C_INT),
VALUE :: la_min
326 INTEGER(KIND=C_INT),
VALUE :: lb_max
327 INTEGER(KIND=C_INT),
VALUE :: lb_min
328 REAL(kind=c_double),
VALUE :: zeta
329 REAL(kind=c_double),
VALUE :: zetb
330 TYPE(c_ptr),
VALUE :: dh
331 TYPE(c_ptr),
VALUE :: dh_inv
332 TYPE(c_ptr),
VALUE :: ra
333 TYPE(c_ptr),
VALUE :: rab
334 TYPE(c_ptr),
VALUE :: npts_global
335 TYPE(c_ptr),
VALUE :: npts_local
336 TYPE(c_ptr),
VALUE :: shift_local
337 TYPE(c_ptr),
VALUE :: border_width
338 REAL(kind=c_double),
VALUE :: radius
339 INTEGER(KIND=C_INT),
VALUE :: o1
340 INTEGER(KIND=C_INT),
VALUE :: o2
341 INTEGER(KIND=C_INT),
VALUE :: n1
342 INTEGER(KIND=C_INT),
VALUE :: n2
343 TYPE(c_ptr),
VALUE :: grid
344 TYPE(c_ptr),
VALUE :: hab
345 TYPE(c_ptr),
VALUE :: pab
346 TYPE(c_ptr),
VALUE :: forces
347 TYPE(c_ptr),
VALUE :: virials
348 TYPE(c_ptr),
VALUE :: hdab
349 TYPE(c_ptr),
VALUE :: hadb
350 TYPE(c_ptr),
VALUE :: a_hdab
351 END SUBROUTINE grid_cpu_integrate_pgf_product_c
354 IF (radius == 0.0_dp)
THEN
359 IF (
PRESENT(use_subpatch))
THEN
360 IF (use_subpatch)
THEN
361 cpassert(
PRESENT(subpatch_pattern))
362 border_mask = iand(63, not(subpatch_pattern))
367 IF (
PRESENT(compute_tau))
THEN
368 my_compute_tau =
LOGICAL(compute_tau, c_bool)
370 my_compute_tau = .false.
373 IF (
PRESENT(use_virial))
THEN
374 my_use_virial = use_virial
376 my_use_virial = .false.
379 IF (calculate_forces)
THEN
380 cpassert(
PRESENT(pab))
381 pab_cptr = c_loc(pab(1, 1))
382 forces(:, :) = 0.0_dp
383 forces_cptr = c_loc(forces(1, 1))
385 pab_cptr = c_null_ptr
386 forces_cptr = c_null_ptr
389 IF (calculate_forces .AND. my_use_virial)
THEN
390 virials(:, :, :) = 0.0_dp
391 virials_cptr = c_loc(virials(1, 1, 1))
393 virials_cptr = c_null_ptr
396 IF (calculate_forces .AND.
PRESENT(hdab))
THEN
397 hdab_cptr = c_loc(hdab(1, 1, 1))
399 hdab_cptr = c_null_ptr
402 IF (calculate_forces .AND.
PRESENT(hadb))
THEN
403 hadb_cptr = c_loc(hadb(1, 1, 1))
405 hadb_cptr = c_null_ptr
408 IF (calculate_forces .AND. my_use_virial .AND.
PRESENT(a_hdab))
THEN
409 a_hdab_cptr = c_loc(a_hdab(1, 1, 1, 1))
411 a_hdab_cptr = c_null_ptr
414 orthorhombic =
LOGICAL(rsgrid%desc%orthorhombic, c_bool)
416 CALL get_rsgrid_properties(rsgrid, npts_global=npts_global, &
417 npts_local=npts_local, &
418 shift_local=shift_local, &
419 border_width=border_width)
421 grid(1:, 1:, 1:) => rsgrid%r(:, :, :)
424 cpassert(is_contiguous(rsgrid%desc%dh))
425 cpassert(is_contiguous(rsgrid%desc%dh_inv))
426 cpassert(is_contiguous(ra))
427 cpassert(is_contiguous(rab))
428 cpassert(is_contiguous(npts_global))
429 cpassert(is_contiguous(npts_local))
430 cpassert(is_contiguous(shift_local))
431 cpassert(is_contiguous(border_width))
432 cpassert(is_contiguous(grid))
433 cpassert(is_contiguous(hab))
434 cpassert(is_contiguous(forces))
435 cpassert(is_contiguous(virials))
436 IF (
PRESENT(pab))
THEN
437 cpassert(is_contiguous(pab))
439 IF (
PRESENT(hdab))
THEN
440 cpassert(is_contiguous(hdab))
442 IF (
PRESENT(a_hdab))
THEN
443 cpassert(is_contiguous(a_hdab))
447 CALL grid_cpu_integrate_pgf_product_c(orthorhombic=orthorhombic, &
448 compute_tau=my_compute_tau, &
449 border_mask=border_mask, &
456 dh=c_loc(rsgrid%desc%dh(1, 1)), &
457 dh_inv=c_loc(rsgrid%desc%dh_inv(1, 1)), &
460 npts_global=c_loc(npts_global(1)), &
461 npts_local=c_loc(npts_local(1)), &
462 shift_local=c_loc(shift_local(1)), &
463 border_width=c_loc(border_width(1)), &
469 grid=c_loc(grid(1, 1, 1)), &
470 hab=c_loc(hab(1, 1)), &
472 forces=forces_cptr, &
473 virials=virials_cptr, &
478 IF (
PRESENT(force_a) .AND. c_associated(forces_cptr)) &
479 force_a = force_a + forces(:, 1)
480 IF (
PRESENT(force_b) .AND. c_associated(forces_cptr)) &
481 force_b = force_b + forces(:, 2)
482 IF (
PRESENT(my_virial_a) .AND. c_associated(virials_cptr)) &
483 my_virial_a = my_virial_a + virials(:, :, 1)
484 IF (
PRESENT(my_virial_b) .AND. c_associated(virials_cptr)) &
485 my_virial_b = my_virial_b + virials(:, :, 2)
498 SUBROUTINE get_rsgrid_properties(rsgrid, npts_global, npts_local, shift_local, border_width)
500 INTEGER,
DIMENSION(:) :: npts_global, npts_local, shift_local, &
506 cpassert(lbound(rsgrid%r, 1) == rsgrid%lb_local(1))
507 cpassert(ubound(rsgrid%r, 1) == rsgrid%ub_local(1))
508 cpassert(lbound(rsgrid%r, 2) == rsgrid%lb_local(2))
509 cpassert(ubound(rsgrid%r, 2) == rsgrid%ub_local(2))
510 cpassert(lbound(rsgrid%r, 3) == rsgrid%lb_local(3))
511 cpassert(ubound(rsgrid%r, 3) == rsgrid%ub_local(3))
518 npts_global = rsgrid%desc%ub - rsgrid%desc%lb + 1
521 npts_local = rsgrid%ub_local - rsgrid%lb_local + 1
524 shift_local = rsgrid%lb_local - rsgrid%desc%lb
528 IF (rsgrid%desc%perd(i) == 1)
THEN
530 cpassert(npts_local(i) == npts_global(i))
531 cpassert(shift_local(i) == 0)
536 cpassert(npts_local(i) <= npts_global(i))
538 cpassert(rsgrid%lb_real(i) == rsgrid%lb_local(i) + rsgrid%desc%border)
539 cpassert(rsgrid%ub_real(i) == rsgrid%ub_local(i) - rsgrid%desc%border)
541 border_width(i) = rsgrid%desc%border
544 END SUBROUTINE get_rsgrid_properties
563 lmin, lmax, npgf, nsgf_set, first_sgf, sphi, zet, &
565 INTEGER,
INTENT(IN) :: nset, nsgf, maxco, maxpgf
566 INTEGER,
DIMENSION(:),
INTENT(IN),
TARGET :: lmin, lmax, npgf, nsgf_set
567 INTEGER,
DIMENSION(:, :),
INTENT(IN) :: first_sgf
568 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN),
TARGET :: sphi, zet
571 CHARACTER(LEN=*),
PARAMETER :: routinen =
'grid_create_basis_set'
574 INTEGER,
DIMENSION(nset),
TARGET :: my_first_sgf
575 TYPE(c_ptr) :: first_sgf_c, lmax_c, lmin_c, npgf_c, &
576 nsgf_set_c, sphi_c, zet_c
578 SUBROUTINE grid_create_basis_set_c(nset, nsgf, maxco, maxpgf, &
579 lmin, lmax, npgf, nsgf_set, first_sgf, sphi, zet, &
581 BIND(C, name="grid_create_basis_set")
582 IMPORT :: c_ptr, c_int
583 INTEGER(KIND=C_INT),
VALUE :: nset
584 INTEGER(KIND=C_INT),
VALUE :: nsgf
585 INTEGER(KIND=C_INT),
VALUE :: maxco
586 INTEGER(KIND=C_INT),
VALUE :: maxpgf
587 TYPE(c_ptr),
VALUE :: lmin
588 TYPE(c_ptr),
VALUE :: lmax
589 TYPE(c_ptr),
VALUE :: npgf
590 TYPE(c_ptr),
VALUE :: nsgf_set
591 TYPE(c_ptr),
VALUE :: first_sgf
592 TYPE(c_ptr),
VALUE :: sphi
593 TYPE(c_ptr),
VALUE :: zet
594 TYPE(c_ptr) :: basis_set
595 END SUBROUTINE grid_create_basis_set_c
598 CALL timeset(routinen, handle)
600 cpassert(
SIZE(lmin) == nset)
601 cpassert(
SIZE(lmin) == nset)
602 cpassert(
SIZE(lmax) == nset)
603 cpassert(
SIZE(npgf) == nset)
604 cpassert(
SIZE(nsgf_set) == nset)
605 cpassert(
SIZE(first_sgf, 2) == nset)
606 cpassert(
SIZE(sphi, 1) == maxco .AND.
SIZE(sphi, 2) == nsgf)
607 cpassert(
SIZE(zet, 1) == maxpgf .AND.
SIZE(zet, 2) == nset)
608 cpassert(.NOT. c_associated(basis_set%c_ptr))
611 cpassert(is_contiguous(lmin))
612 cpassert(is_contiguous(lmax))
613 cpassert(is_contiguous(npgf))
614 cpassert(is_contiguous(nsgf_set))
615 cpassert(is_contiguous(my_first_sgf))
616 cpassert(is_contiguous(sphi))
617 cpassert(is_contiguous(zet))
623 nsgf_set_c = c_null_ptr
624 first_sgf_c = c_null_ptr
630 lmin_c = c_loc(lmin(1))
631 lmax_c = c_loc(lmax(1))
632 npgf_c = c_loc(npgf(1))
633 nsgf_set_c = c_loc(nsgf_set(1))
635 IF (
SIZE(first_sgf) > 0)
THEN
636 my_first_sgf(:) = first_sgf(1, :)
637 first_sgf_c = c_loc(my_first_sgf(1))
639 IF (
SIZE(sphi) > 0)
THEN
640 sphi_c = c_loc(sphi(1, 1))
642 IF (
SIZE(zet) > 0)
THEN
643 zet_c = c_loc(zet(1, 1))
646 CALL grid_create_basis_set_c(nset=nset, &
653 nsgf_set=nsgf_set_c, &
654 first_sgf=first_sgf_c, &
657 basis_set=basis_set%c_ptr)
658 cpassert(c_associated(basis_set%c_ptr))
660 CALL timestop(handle)
671 CHARACTER(LEN=*),
PARAMETER :: routinen =
'grid_free_basis_set'
675 SUBROUTINE grid_free_basis_set_c(basis_set) &
676 BIND(C, name="grid_free_basis_set")
678 TYPE(c_ptr),
VALUE :: basis_set
679 END SUBROUTINE grid_free_basis_set_c
682 CALL timeset(routinen, handle)
684 cpassert(c_associated(basis_set%c_ptr))
686 CALL grid_free_basis_set_c(basis_set%c_ptr)
688 basis_set%c_ptr = c_null_ptr
690 CALL timestop(handle)
719 block_offsets, atom_positions, atom_kinds, basis_sets, &
720 level_list, iatom_list, jatom_list, &
721 iset_list, jset_list, ipgf_list, jpgf_list, &
722 border_mask_list, block_num_list, &
723 radius_list, rab_list, rs_grids, task_list)
725 INTEGER,
INTENT(IN) :: ntasks, natoms, nkinds, nblocks
726 INTEGER,
DIMENSION(:),
INTENT(IN),
TARGET :: block_offsets
727 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN),
TARGET :: atom_positions
728 INTEGER,
DIMENSION(:),
INTENT(IN),
TARGET :: atom_kinds
730 INTENT(IN),
TARGET :: basis_sets
731 INTEGER,
DIMENSION(:),
INTENT(IN),
TARGET :: level_list, iatom_list, jatom_list, &
732 iset_list, jset_list, ipgf_list, &
733 jpgf_list, border_mask_list, &
735 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN),
TARGET :: radius_list
736 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN),
TARGET :: rab_list
738 INTENT(IN) :: rs_grids
741 CHARACTER(LEN=*),
PARAMETER :: routinen =
'grid_create_task_list'
743 INTEGER :: handle, ikind, ilevel, nlevels
744 INTEGER,
ALLOCATABLE,
DIMENSION(:, :),
TARGET :: border_width, npts_global, npts_local, &
746 LOGICAL(KIND=C_BOOL) :: orthorhombic
747 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :), &
749 TYPE(c_ptr) :: block_num_list_c, block_offsets_c, border_mask_list_c, iatom_list_c, &
750 ipgf_list_c, iset_list_c, jatom_list_c, jpgf_list_c, jset_list_c, level_list_c, &
751 rab_list_c, radius_list_c
752 TYPE(c_ptr),
ALLOCATABLE,
DIMENSION(:),
TARGET :: basis_sets_c
754 SUBROUTINE grid_create_task_list_c(orthorhombic, &
755 ntasks, nlevels, natoms, nkinds, nblocks, &
756 block_offsets, atom_positions, atom_kinds, basis_sets, &
757 level_list, iatom_list, jatom_list, &
758 iset_list, jset_list, ipgf_list, jpgf_list, &
759 border_mask_list, block_num_list, &
760 radius_list, rab_list, &
761 npts_global, npts_local, shift_local, &
762 border_width, dh, dh_inv, task_list) &
763 BIND(C, name="grid_create_task_list")
764 IMPORT :: c_ptr, c_int, c_bool
765 LOGICAL(KIND=C_BOOL),
VALUE :: orthorhombic
766 INTEGER(KIND=C_INT),
VALUE :: ntasks
767 INTEGER(KIND=C_INT),
VALUE :: nlevels
768 INTEGER(KIND=C_INT),
VALUE :: natoms
769 INTEGER(KIND=C_INT),
VALUE :: nkinds
770 INTEGER(KIND=C_INT),
VALUE :: nblocks
771 TYPE(c_ptr),
VALUE :: block_offsets
772 TYPE(c_ptr),
VALUE :: atom_positions
773 TYPE(c_ptr),
VALUE :: atom_kinds
774 TYPE(c_ptr),
VALUE :: basis_sets
775 TYPE(c_ptr),
VALUE :: level_list
776 TYPE(c_ptr),
VALUE :: iatom_list
777 TYPE(c_ptr),
VALUE :: jatom_list
778 TYPE(c_ptr),
VALUE :: iset_list
779 TYPE(c_ptr),
VALUE :: jset_list
780 TYPE(c_ptr),
VALUE :: ipgf_list
781 TYPE(c_ptr),
VALUE :: jpgf_list
782 TYPE(c_ptr),
VALUE :: border_mask_list
783 TYPE(c_ptr),
VALUE :: block_num_list
784 TYPE(c_ptr),
VALUE :: radius_list
785 TYPE(c_ptr),
VALUE :: rab_list
786 TYPE(c_ptr),
VALUE :: npts_global
787 TYPE(c_ptr),
VALUE :: npts_local
788 TYPE(c_ptr),
VALUE :: shift_local
789 TYPE(c_ptr),
VALUE :: border_width
790 TYPE(c_ptr),
VALUE :: dh
791 TYPE(c_ptr),
VALUE :: dh_inv
792 TYPE(c_ptr) :: task_list
793 END SUBROUTINE grid_create_task_list_c
796 CALL timeset(routinen, handle)
798 cpassert(
SIZE(block_offsets) == nblocks)
799 cpassert(
SIZE(atom_positions, 1) == 3 .AND.
SIZE(atom_positions, 2) == natoms)
800 cpassert(
SIZE(atom_kinds) == natoms)
801 cpassert(
SIZE(basis_sets) == nkinds)
802 cpassert(
SIZE(level_list) == ntasks)
803 cpassert(
SIZE(iatom_list) == ntasks)
804 cpassert(
SIZE(jatom_list) == ntasks)
805 cpassert(
SIZE(iset_list) == ntasks)
806 cpassert(
SIZE(jset_list) == ntasks)
807 cpassert(
SIZE(ipgf_list) == ntasks)
808 cpassert(
SIZE(jpgf_list) == ntasks)
809 cpassert(
SIZE(border_mask_list) == ntasks)
810 cpassert(
SIZE(block_num_list) == ntasks)
811 cpassert(
SIZE(radius_list) == ntasks)
812 cpassert(
SIZE(rab_list, 1) == 3 .AND.
SIZE(rab_list, 2) == ntasks)
814 ALLOCATE (basis_sets_c(nkinds))
816 basis_sets_c(ikind) = basis_sets(ikind)%c_ptr
819 nlevels =
SIZE(rs_grids)
820 cpassert(nlevels > 0)
821 orthorhombic =
LOGICAL(rs_grids(1)%desc%orthorhombic, c_bool)
823 ALLOCATE (npts_global(3, nlevels), npts_local(3, nlevels))
824 ALLOCATE (shift_local(3, nlevels), border_width(3, nlevels))
825 ALLOCATE (dh(3, 3, nlevels), dh_inv(3, 3, nlevels))
826 DO ilevel = 1, nlevels
827 associate(rsgrid => rs_grids(ilevel))
828 CALL get_rsgrid_properties(rsgrid=rsgrid, &
829 npts_global=npts_global(:, ilevel), &
830 npts_local=npts_local(:, ilevel), &
831 shift_local=shift_local(:, ilevel), &
832 border_width=border_width(:, ilevel))
833 cpassert(rsgrid%desc%orthorhombic .EQV. orthorhombic)
834 dh(:, :, ilevel) = rsgrid%desc%dh(:, :)
835 dh_inv(:, :, ilevel) = rsgrid%desc%dh_inv(:, :)
840 cpassert(is_contiguous(block_offsets))
841 cpassert(is_contiguous(atom_positions))
842 cpassert(is_contiguous(atom_kinds))
843 cpassert(is_contiguous(basis_sets))
844 cpassert(is_contiguous(level_list))
845 cpassert(is_contiguous(iatom_list))
846 cpassert(is_contiguous(jatom_list))
847 cpassert(is_contiguous(iset_list))
848 cpassert(is_contiguous(jset_list))
849 cpassert(is_contiguous(ipgf_list))
850 cpassert(is_contiguous(jpgf_list))
851 cpassert(is_contiguous(border_mask_list))
852 cpassert(is_contiguous(block_num_list))
853 cpassert(is_contiguous(radius_list))
854 cpassert(is_contiguous(rab_list))
855 cpassert(is_contiguous(npts_global))
856 cpassert(is_contiguous(npts_local))
857 cpassert(is_contiguous(shift_local))
858 cpassert(is_contiguous(border_width))
859 cpassert(is_contiguous(dh))
860 cpassert(is_contiguous(dh_inv))
864 block_offsets_c = c_loc(block_offsets(1))
865 level_list_c = c_loc(level_list(1))
866 iatom_list_c = c_loc(iatom_list(1))
867 jatom_list_c = c_loc(jatom_list(1))
868 iset_list_c = c_loc(iset_list(1))
869 jset_list_c = c_loc(jset_list(1))
870 ipgf_list_c = c_loc(ipgf_list(1))
871 jpgf_list_c = c_loc(jpgf_list(1))
872 border_mask_list_c = c_loc(border_mask_list(1))
873 block_num_list_c = c_loc(block_num_list(1))
874 radius_list_c = c_loc(radius_list(1))
875 rab_list_c = c_loc(rab_list(1, 1))
878 block_offsets_c = c_null_ptr
879 level_list_c = c_null_ptr
880 iatom_list_c = c_null_ptr
881 jatom_list_c = c_null_ptr
882 iset_list_c = c_null_ptr
883 jset_list_c = c_null_ptr
884 ipgf_list_c = c_null_ptr
885 jpgf_list_c = c_null_ptr
886 border_mask_list_c = c_null_ptr
887 block_num_list_c = c_null_ptr
888 radius_list_c = c_null_ptr
889 rab_list_c = c_null_ptr
893 CALL grid_create_task_list_c(orthorhombic=orthorhombic, &
899 block_offsets=block_offsets_c, &
900 atom_positions=c_loc(atom_positions(1, 1)), &
901 atom_kinds=c_loc(atom_kinds(1)), &
902 basis_sets=c_loc(basis_sets_c(1)), &
903 level_list=level_list_c, &
904 iatom_list=iatom_list_c, &
905 jatom_list=jatom_list_c, &
906 iset_list=iset_list_c, &
907 jset_list=jset_list_c, &
908 ipgf_list=ipgf_list_c, &
909 jpgf_list=jpgf_list_c, &
910 border_mask_list=border_mask_list_c, &
911 block_num_list=block_num_list_c, &
912 radius_list=radius_list_c, &
913 rab_list=rab_list_c, &
914 npts_global=c_loc(npts_global(1, 1)), &
915 npts_local=c_loc(npts_local(1, 1)), &
916 shift_local=c_loc(shift_local(1, 1)), &
917 border_width=c_loc(border_width(1, 1)), &
918 dh=c_loc(dh(1, 1, 1)), &
919 dh_inv=c_loc(dh_inv(1, 1, 1)), &
920 task_list=task_list%c_ptr)
922 cpassert(c_associated(task_list%c_ptr))
924 CALL timestop(handle)
935 CHARACTER(LEN=*),
PARAMETER :: routinen =
'grid_free_task_list'
939 SUBROUTINE grid_free_task_list_c(task_list) &
940 BIND(C, name="grid_free_task_list")
942 TYPE(c_ptr),
VALUE :: task_list
943 END SUBROUTINE grid_free_task_list_c
946 CALL timeset(routinen, handle)
948 IF (c_associated(task_list%c_ptr))
THEN
949 CALL grid_free_task_list_c(task_list%c_ptr)
952 task_list%c_ptr = c_null_ptr
954 CALL timestop(handle)
967 INTEGER,
INTENT(IN) :: ga_gb_function
968 TYPE(offload_buffer_type),
INTENT(IN) :: pab_blocks
969 TYPE(realspace_grid_type),
DIMENSION(:), &
970 INTENT(IN) :: rs_grids
972 CHARACTER(LEN=*),
PARAMETER :: routinen =
'grid_collocate_task_list'
974 INTEGER :: handle, ilevel, nlevels
975 INTEGER,
ALLOCATABLE,
DIMENSION(:, :),
TARGET :: npts_local
976 TYPE(c_ptr),
ALLOCATABLE,
DIMENSION(:),
TARGET :: grids_c
978 SUBROUTINE grid_collocate_task_list_c(task_list, func, nlevels, &
979 npts_local, pab_blocks, grids) &
980 BIND(C, name="grid_collocate_task_list")
981 IMPORT :: c_ptr, c_int, c_bool
982 TYPE(c_ptr),
VALUE :: task_list
983 INTEGER(KIND=C_INT),
VALUE :: func
984 INTEGER(KIND=C_INT),
VALUE :: nlevels
985 TYPE(c_ptr),
VALUE :: npts_local
986 TYPE(c_ptr),
VALUE :: pab_blocks
987 TYPE(c_ptr),
VALUE :: grids
988 END SUBROUTINE grid_collocate_task_list_c
991 CALL timeset(routinen, handle)
993 nlevels =
SIZE(rs_grids)
994 cpassert(nlevels > 0)
996 ALLOCATE (grids_c(nlevels))
997 ALLOCATE (npts_local(3, nlevels))
998 DO ilevel = 1, nlevels
999 associate(rsgrid => rs_grids(ilevel))
1000 npts_local(:, ilevel) = rsgrid%ub_local - rsgrid%lb_local + 1
1001 grids_c(ilevel) = rsgrid%buffer%c_ptr
1006 cpassert(is_contiguous(npts_local))
1007 cpassert(is_contiguous(grids_c))
1010 cpassert(c_associated(task_list%c_ptr))
1011 cpassert(c_associated(pab_blocks%c_ptr))
1013 CALL grid_collocate_task_list_c(task_list=task_list%c_ptr, &
1014 func=ga_gb_function, &
1016 npts_local=c_loc(npts_local(1, 1)), &
1017 pab_blocks=pab_blocks%c_ptr, &
1018 grids=c_loc(grids_c(1)))
1020 CALL timestop(handle)
1037 pab_blocks, rs_grids, hab_blocks, forces, virial)
1039 LOGICAL,
INTENT(IN) :: compute_tau, calculate_forces, &
1041 TYPE(offload_buffer_type),
INTENT(IN) :: pab_blocks
1042 TYPE(realspace_grid_type),
DIMENSION(:), &
1043 INTENT(IN) :: rs_grids
1044 TYPE(offload_buffer_type),
INTENT(INOUT) :: hab_blocks
1045 REAL(kind=dp),
DIMENSION(:, :),
INTENT(INOUT), &
1047 REAL(kind=dp),
DIMENSION(3, 3),
INTENT(INOUT), &
1050 CHARACTER(LEN=*),
PARAMETER :: routinen =
'grid_integrate_task_list'
1052 INTEGER :: handle, ilevel, nlevels
1053 INTEGER,
ALLOCATABLE,
DIMENSION(:, :),
TARGET :: npts_local
1054 TYPE(c_ptr) :: forces_c, virial_c
1055 TYPE(c_ptr),
ALLOCATABLE,
DIMENSION(:),
TARGET :: grids_c
1057 SUBROUTINE grid_integrate_task_list_c(task_list, compute_tau, natoms, &
1058 nlevels, npts_local, &
1059 pab_blocks, grids, hab_blocks, forces, virial) &
1060 BIND(C, name="grid_integrate_task_list")
1061 IMPORT :: c_ptr, c_int, c_bool
1062 TYPE(c_ptr),
VALUE :: task_list
1063 LOGICAL(KIND=C_BOOL),
VALUE :: compute_tau
1064 INTEGER(KIND=C_INT),
VALUE :: natoms
1065 INTEGER(KIND=C_INT),
VALUE :: nlevels
1066 TYPE(c_ptr),
VALUE :: npts_local
1067 TYPE(c_ptr),
VALUE :: pab_blocks
1068 TYPE(c_ptr),
VALUE :: grids
1069 TYPE(c_ptr),
VALUE :: hab_blocks
1070 TYPE(c_ptr),
VALUE :: forces
1071 TYPE(c_ptr),
VALUE :: virial
1072 END SUBROUTINE grid_integrate_task_list_c
1075 CALL timeset(routinen, handle)
1077 nlevels =
SIZE(rs_grids)
1078 cpassert(nlevels > 0)
1080 ALLOCATE (grids_c(nlevels))
1081 ALLOCATE (npts_local(3, nlevels))
1082 DO ilevel = 1, nlevels
1083 associate(rsgrid => rs_grids(ilevel))
1084 npts_local(:, ilevel) = rsgrid%ub_local - rsgrid%lb_local + 1
1085 grids_c(ilevel) = rsgrid%buffer%c_ptr
1089 IF (calculate_forces)
THEN
1090 forces_c = c_loc(forces(1, 1))
1092 forces_c = c_null_ptr
1095 IF (calculate_virial)
THEN
1096 virial_c = c_loc(virial(1, 1))
1098 virial_c = c_null_ptr
1102 cpassert(is_contiguous(npts_local))
1103 cpassert(is_contiguous(grids_c))
1104 cpassert(is_contiguous(forces))
1105 cpassert(is_contiguous(virial))
1108 cpassert(
SIZE(forces, 1) == 3)
1109 cpassert(c_associated(task_list%c_ptr))
1110 cpassert(c_associated(hab_blocks%c_ptr))
1111 cpassert(c_associated(pab_blocks%c_ptr) .OR. .NOT. calculate_forces)
1112 cpassert(c_associated(pab_blocks%c_ptr) .OR. .NOT. calculate_virial)
1114 CALL grid_integrate_task_list_c(task_list=task_list%c_ptr, &
1115 compute_tau=
LOGICAL(compute_tau, C_BOOL), &
1116 natoms=size(forces, 2), &
1118 npts_local=c_loc(npts_local(1, 1)), &
1119 pab_blocks=pab_blocks%c_ptr, &
1120 grids=c_loc(grids_c(1)), &
1121 hab_blocks=hab_blocks%c_ptr, &
1125 CALL timestop(handle)
1134 SUBROUTINE grid_library_init_c()
BIND(C, name="grid_library_init")
1135 END SUBROUTINE grid_library_init_c
1138 CALL grid_library_init_c()
1148 SUBROUTINE grid_library_finalize_c()
BIND(C, name="grid_library_finalize")
1149 END SUBROUTINE grid_library_finalize_c
1152 CALL grid_library_finalize_c()
1164 INTEGER,
INTENT(IN) :: backend
1168 SUBROUTINE grid_library_set_config_c(backend, validate, apply_cutoff) &
1169 BIND(C, name="grid_library_set_config")
1170 IMPORT :: c_int, c_bool
1171 INTEGER(KIND=C_INT),
VALUE :: backend
1172 LOGICAL(KIND=C_BOOL),
VALUE :: validate
1174 END SUBROUTINE grid_library_set_config_c
1177 CALL grid_library_set_config_c(backend=backend, &
1178 validate=
LOGICAL(validate, C_BOOL), &
1190 TYPE(mp_comm_type) :: mpi_comm
1191 INTEGER,
INTENT(IN) :: output_unit
1194 SUBROUTINE grid_library_print_stats_c(mpi_comm, print_func, output_unit) &
1195 BIND(C, name="grid_library_print_stats")
1196 IMPORT :: c_funptr, c_int
1197 INTEGER(KIND=C_INT),
VALUE :: mpi_comm
1199 INTEGER(KIND=C_INT),
VALUE :: output_unit
1200 END SUBROUTINE grid_library_print_stats_c
1204 CALL grid_library_print_stats_c(mpi_comm=mpi_comm%get_handle(), &
1206 output_unit=output_unit)
1217 SUBROUTINE print_func(msg, msglen, output_unit)
BIND(C, name="grid_api_print_func")
1218 CHARACTER(KIND=C_CHAR),
INTENT(IN) :: msg(*)
1219 INTEGER(KIND=C_INT),
INTENT(IN),
VALUE :: msglen, output_unit
1221 IF (output_unit <= 0)
RETURN
1222 WRITE (output_unit, fmt=
"(100A)", advance=
"NO") msg(1:msglen)
static void print_func(const char *msg, int msglen, int output_unit)
Wrapper for printf, passed to dbm_library_print_stats.
void grid_create_basis_set(const int nset, const int nsgf, const int maxco, const int maxpgf, const int lmin[nset], const int lmax[nset], const int npgf[nset], const int nsgf_set[nset], const int first_sgf[nset], const double sphi[nsgf][maxco], const double zet[nset][maxpgf], grid_basis_set **basis_set_out)
Allocates a basis set which can be passed to grid_create_task_list. See grid_task_list....
void grid_free_basis_set(grid_basis_set *basis_set)
Deallocates given basis set.
void apply_cutoff(void *ptr)
void grid_library_finalize(void)
Finalizes the grid library.
void grid_library_init(void)
Initializes the grid library.
void grid_library_set_config(const enum grid_backend backend, const bool validate, const bool apply_cutoff)
Configures the grid library.
void grid_library_print_stats(const int fortran_comm, void(*print_func)(const char *, int, int), const int output_unit)
Prints statistics gathered by the grid library.
Fortran API for the grid package, which is written in C.
integer, parameter, public grid_func_adbmdab_z
integer, parameter, public grid_func_core_x
integer, parameter, public grid_func_adbmdab_y
integer, parameter, public grid_func_ardbmdarb_yx
integer, parameter, public grid_func_dab_z
subroutine, public grid_collocate_task_list(task_list, ga_gb_function, pab_blocks, rs_grids)
Collocate all tasks of in given list onto given grids.
integer, parameter, public grid_func_dzdx
integer, parameter, public grid_func_ardbmdarb_zz
integer, parameter, public grid_backend_auto
integer, parameter, public grid_backend_gpu
subroutine, public grid_free_task_list(task_list)
Deallocates given task list, basis_sets have to be freed separately.
integer, parameter, public grid_func_dzdz
integer, parameter, public grid_func_dydz
integer, parameter, public grid_func_adb_y
integer, parameter, public grid_func_dxdy
integer, parameter, public grid_func_dabpadb_y
integer, parameter, public grid_func_ardbmdarb_xy
integer, parameter, public grid_func_dab_y
integer, parameter, public grid_backend_hip
subroutine, public grid_create_task_list(ntasks, natoms, nkinds, nblocks, block_offsets, atom_positions, atom_kinds, basis_sets, level_list, iatom_list, jatom_list, iset_list, jset_list, ipgf_list, jpgf_list, border_mask_list, block_num_list, radius_list, rab_list, rs_grids, task_list)
Allocates a task list which can be passed to grid_collocate_task_list.
integer, parameter, public grid_func_adb_z
integer, parameter, public grid_func_ardbmdarb_zx
integer, parameter, public grid_func_adb_x
integer, parameter, public grid_func_dxdx
integer, parameter, public grid_func_ardbmdarb_xx
integer, parameter, public grid_func_dadb
integer, parameter, public grid_backend_dgemm
integer, parameter, public grid_func_dydy
integer, parameter, public grid_func_dabpadb_z
integer, parameter, public grid_backend_cpu
integer, parameter, public grid_func_dabpadb_x
integer, parameter, public grid_func_dx
integer, parameter, public grid_func_dz
integer, parameter, public grid_func_ardbmdarb_yz
integer, parameter, public grid_func_ab
subroutine, public integrate_pgf_product(la_max, zeta, la_min, lb_max, zetb, lb_min, ra, rab, rsgrid, hab, pab, o1, o2, radius, calculate_forces, force_a, force_b, compute_tau, use_virial, my_virial_a, my_virial_b, hdab, hadb, a_hdab, use_subpatch, subpatch_pattern)
low level function to compute matrix elements of primitive gaussian functions
integer, parameter, public grid_func_ardbmdarb_yy
subroutine, public grid_integrate_task_list(task_list, compute_tau, calculate_forces, calculate_virial, pab_blocks, rs_grids, hab_blocks, forces, virial)
Integrate all tasks of in given list from given grids.
integer, parameter, public grid_func_core_y
integer, parameter, public grid_backend_ref
integer, parameter, public grid_func_adbmdab_x
integer, parameter, public grid_func_dab_x
subroutine, public collocate_pgf_product(la_max, zeta, la_min, lb_max, zetb, lb_min, ra, rab, scale, pab, o1, o2, rsgrid, ga_gb_function, radius, use_subpatch, subpatch_pattern)
low level collocation of primitive gaussian functions
integer, parameter, public grid_func_ardbmdarb_zy
integer, parameter, public grid_func_core_z
integer, parameter, public grid_func_dy
integer, parameter, public grid_func_ardbmdarb_xz
Defines the basic variable types.
integer, parameter, public dp
Interface to the message passing library MPI.
Fortran API for the offload package, which is written in C.