Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 : ! **************************************************************************************************
8 : !> \brief Calculation of the non-local pseudopotential contribution to the core Hamiltonian
9 : !> <a|V(non-local)|b> = <a|p(l,i)>*h(i,j)*<p(l,j)|b>
10 : !> \par History
11 : !> - refactered from qs_core_hamiltian [Joost VandeVondele, 2008-11-01]
12 : !> - full rewrite [jhu, 2009-01-23]
13 : !> - Extended by the derivatives for DFPT [Sandra Luber, Edward Ditler, 2021]
14 : ! **************************************************************************************************
15 : MODULE core_ppnl
16 : USE ai_angmom, ONLY: angmom
17 : USE ai_overlap, ONLY: overlap
18 : USE atomic_kind_types, ONLY: atomic_kind_type,&
19 : get_atomic_kind_set
20 : USE basis_set_types, ONLY: gto_basis_set_p_type,&
21 : gto_basis_set_type
22 : USE cp_dbcsr_api, ONLY: dbcsr_add,&
23 : dbcsr_get_block_p,&
24 : dbcsr_p_type
25 : USE external_potential_types, ONLY: gth_potential_p_type,&
26 : gth_potential_type,&
27 : sgp_potential_p_type,&
28 : sgp_potential_type
29 : USE kinds, ONLY: dp,&
30 : int_8
31 : USE orbital_pointers, ONLY: init_orbital_pointers,&
32 : nco,&
33 : ncoset
34 : USE particle_types, ONLY: particle_type
35 : USE qs_force_types, ONLY: qs_force_type
36 : USE qs_kind_types, ONLY: get_qs_kind,&
37 : get_qs_kind_set,&
38 : qs_kind_type
39 : USE qs_neighbor_list_types, ONLY: neighbor_list_set_p_type
40 : USE sap_kind_types, ONLY: alist_type,&
41 : clist_type,&
42 : get_alist,&
43 : release_sap_int,&
44 : sap_int_type,&
45 : sap_sort
46 : USE virial_methods, ONLY: virial_pair_force
47 : USE virial_types, ONLY: virial_type
48 :
49 : !$ USE OMP_LIB, ONLY: omp_lock_kind, &
50 : !$ omp_init_lock, omp_set_lock, &
51 : !$ omp_unset_lock, omp_destroy_lock
52 :
53 : #include "./base/base_uses.f90"
54 :
55 : IMPLICIT NONE
56 :
57 : PRIVATE
58 :
59 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'core_ppnl'
60 :
61 : PUBLIC :: build_core_ppnl
62 :
63 : CONTAINS
64 :
65 : ! **************************************************************************************************
66 : !> \brief ...
67 : !> \param matrix_h ...
68 : !> \param matrix_p ...
69 : !> \param force ...
70 : !> \param virial ...
71 : !> \param calculate_forces ...
72 : !> \param use_virial ...
73 : !> \param nder ...
74 : !> \param qs_kind_set ...
75 : !> \param atomic_kind_set ...
76 : !> \param particle_set ...
77 : !> \param sab_orb ...
78 : !> \param sap_ppnl ...
79 : !> \param eps_ppnl ...
80 : !> \param nimages ...
81 : !> \param cell_to_index ...
82 : !> \param basis_type ...
83 : !> \param deltaR Weighting factors of the derivatives wrt. nuclear positions
84 : !> \param matrix_l ...
85 : !> \param atcore ...
86 : ! **************************************************************************************************
87 14476 : SUBROUTINE build_core_ppnl(matrix_h, matrix_p, force, virial, calculate_forces, use_virial, nder, &
88 : qs_kind_set, atomic_kind_set, particle_set, sab_orb, sap_ppnl, eps_ppnl, &
89 14476 : nimages, cell_to_index, basis_type, deltaR, matrix_l, atcore)
90 :
91 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_h, matrix_p
92 : TYPE(qs_force_type), DIMENSION(:), POINTER :: force
93 : TYPE(virial_type), POINTER :: virial
94 : LOGICAL, INTENT(IN) :: calculate_forces
95 : LOGICAL :: use_virial
96 : INTEGER :: nder
97 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
98 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
99 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
100 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
101 : POINTER :: sab_orb, sap_ppnl
102 : REAL(KIND=dp), INTENT(IN) :: eps_ppnl
103 : INTEGER, INTENT(IN) :: nimages
104 : INTEGER, DIMENSION(:, :, :), OPTIONAL, POINTER :: cell_to_index
105 : CHARACTER(LEN=*), INTENT(IN) :: basis_type
106 : REAL(KIND=dp), DIMENSION(:, :), INTENT(IN), &
107 : OPTIONAL :: deltaR
108 : TYPE(dbcsr_p_type), DIMENSION(:, :), OPTIONAL, &
109 : POINTER :: matrix_l
110 : REAL(KIND=dp), DIMENSION(:), INTENT(INOUT), &
111 : OPTIONAL :: atcore
112 :
113 : CHARACTER(LEN=*), PARAMETER :: routineN = 'build_core_ppnl'
114 :
115 : INTEGER :: atom_a, first_col, handle, i, i_dim, iab, iac, iatom, ib, ibc, icol, ikind, &
116 : ilist, img, irow, iset, j, jatom, jb, jkind, jneighbor, kac, katom, kbc, kkind, l, &
117 : lc_max, lc_min, ldai, ldsab, lppnl, maxco, maxder, maxl, maxlgto, maxlppnl, maxppnl, &
118 : maxsgf, na, natom, nb, ncoa, ncoc, nkind, nlist, nneighbor, nnl, np, nppnl, nprjc, nseta, &
119 : nsgfa, prjc, sgfa, slot
120 14476 : INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_of_kind, kind_of
121 : INTEGER, DIMENSION(3) :: cell_b, cell_c
122 14476 : INTEGER, DIMENSION(:), POINTER :: la_max, la_min, npgfa, nprj_ppnl, &
123 14476 : nsgf_seta
124 14476 : INTEGER, DIMENSION(:, :), POINTER :: first_sgfa
125 : LOGICAL :: do_dR, do_gth, do_kp, do_soc, doat, &
126 : found, ppnl_present
127 : REAL(KIND=dp) :: atk, dac, f0, ppnl_radius
128 14476 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: radp
129 14476 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: sab, work
130 14476 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: ai_work, lab, work_l
131 : REAL(KIND=dp), DIMENSION(1) :: rprjc, zetc
132 : REAL(KIND=dp), DIMENSION(3) :: fa, fb, rab, rac, rbc
133 : REAL(KIND=dp), DIMENSION(3, 3) :: pv_thread
134 : TYPE(gto_basis_set_type), POINTER :: orb_basis_set
135 14476 : TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: basis_set
136 : TYPE(gth_potential_type), POINTER :: gth_potential
137 14476 : TYPE(gth_potential_p_type), DIMENSION(:), POINTER :: gpotential
138 : TYPE(clist_type), POINTER :: clist
139 : TYPE(alist_type), POINTER :: alist_ac, alist_bc
140 28952 : REAL(KIND=dp), DIMENSION(SIZE(particle_set)) :: at_thread
141 14476 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: achint, acint, alkint, bchint, bcint, &
142 14476 : blkint
143 14476 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: cprj, h_block, l_block_x, l_block_y, &
144 14476 : l_block_z, p_block, r_2block, &
145 14476 : r_3block, rpgfa, sphi_a, vprj_ppnl, &
146 14476 : wprj_ppnl, zeta
147 14476 : REAL(KIND=dp), DIMENSION(:), POINTER :: a_nl, alpha_ppnl, hprj, set_radius_a
148 28952 : REAL(KIND=dp), DIMENSION(3, SIZE(particle_set)) :: force_thread
149 14476 : TYPE(sap_int_type), DIMENSION(:), POINTER :: sap_int
150 14476 : TYPE(sgp_potential_p_type), DIMENSION(:), POINTER :: spotential
151 : TYPE(sgp_potential_type), POINTER :: sgp_potential
152 :
153 : !$ INTEGER(kind=omp_lock_kind), &
154 14476 : !$ ALLOCATABLE, DIMENSION(:) :: locks
155 : !$ INTEGER(KIND=int_8) :: iatom8
156 : !$ INTEGER :: lock_num, hash
157 : !$ INTEGER, PARAMETER :: nlock = 501
158 :
159 : MARK_USED(int_8)
160 :
161 14476 : do_dR = .FALSE.
162 72 : IF (PRESENT(deltaR)) do_dR = .TRUE.
163 14476 : doat = .FALSE.
164 14476 : IF (PRESENT(atcore)) doat = .TRUE.
165 :
166 14476 : IF (calculate_forces) THEN
167 6171 : CALL timeset(routineN//"_forces", handle)
168 : ELSE
169 8305 : CALL timeset(routineN, handle)
170 : END IF
171 :
172 14476 : do_soc = PRESENT(matrix_l)
173 :
174 14476 : ppnl_present = ASSOCIATED(sap_ppnl)
175 :
176 14476 : IF (ppnl_present) THEN
177 :
178 14476 : nkind = SIZE(atomic_kind_set)
179 14476 : natom = SIZE(particle_set)
180 :
181 14476 : do_kp = (nimages > 1)
182 :
183 14476 : IF (do_kp) THEN
184 220 : CPASSERT(PRESENT(cell_to_index) .AND. ASSOCIATED(cell_to_index))
185 : END IF
186 :
187 14476 : IF (calculate_forces .OR. doat) THEN
188 6233 : IF (SIZE(matrix_p, 1) == 2) THEN
189 1894 : DO img = 1, nimages
190 : CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
191 1230 : alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
192 : CALL dbcsr_add(matrix_p(2, img)%matrix, matrix_p(1, img)%matrix, &
193 1894 : alpha_scalar=-2.0_dp, beta_scalar=1.0_dp)
194 : END DO
195 : END IF
196 : END IF
197 :
198 14476 : maxder = ncoset(nder)
199 :
200 : CALL get_qs_kind_set(qs_kind_set, &
201 : maxco=maxco, &
202 : maxlgto=maxlgto, &
203 : maxsgf=maxsgf, &
204 : maxlppnl=maxlppnl, &
205 : maxppnl=maxppnl, &
206 14476 : basis_type=basis_type)
207 :
208 14476 : maxl = MAX(maxlgto, maxlppnl)
209 14476 : CALL init_orbital_pointers(maxl + nder + 1)
210 :
211 14476 : ldsab = MAX(maxco, ncoset(maxlppnl), maxsgf, maxppnl)
212 14476 : ldai = ncoset(maxl + nder + 1)
213 :
214 : ! sap_int needs to be shared as multiple threads need to access this
215 99784 : ALLOCATE (sap_int(nkind*nkind))
216 70832 : DO i = 1, nkind*nkind
217 56356 : NULLIFY (sap_int(i)%alist, sap_int(i)%asort, sap_int(i)%aindex)
218 70832 : sap_int(i)%nalist = 0
219 : END DO
220 :
221 : ! Set up direct access to basis and potential
222 154532 : ALLOCATE (basis_set(nkind), gpotential(nkind), spotential(nkind))
223 41860 : DO ikind = 1, nkind
224 27384 : CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set, basis_type=basis_type)
225 27384 : IF (ASSOCIATED(orb_basis_set)) THEN
226 27384 : basis_set(ikind)%gto_basis_set => orb_basis_set
227 : ELSE
228 0 : NULLIFY (basis_set(ikind)%gto_basis_set)
229 : END IF
230 27384 : CALL get_qs_kind(qs_kind_set(ikind), gth_potential=gth_potential, sgp_potential=sgp_potential)
231 27384 : NULLIFY (gpotential(ikind)%gth_potential)
232 27384 : NULLIFY (spotential(ikind)%sgp_potential)
233 41860 : IF (ASSOCIATED(gth_potential)) THEN
234 27152 : gpotential(ikind)%gth_potential => gth_potential
235 27152 : IF (do_soc .AND. (.NOT. gth_potential%soc)) THEN
236 0 : CPABORT("Spin-orbit coupling selected, but GTH potential without SOC parameters provided")
237 : END IF
238 232 : ELSE IF (ASSOCIATED(sgp_potential)) THEN
239 10 : spotential(ikind)%sgp_potential => sgp_potential
240 : END IF
241 : END DO
242 :
243 : ! Allocate sap int
244 782634 : DO slot = 1, sap_ppnl(1)%nl_size
245 :
246 768158 : ikind = sap_ppnl(1)%nlist_task(slot)%ikind
247 768158 : kkind = sap_ppnl(1)%nlist_task(slot)%jkind
248 768158 : iatom = sap_ppnl(1)%nlist_task(slot)%iatom
249 768158 : katom = sap_ppnl(1)%nlist_task(slot)%jatom
250 768158 : nlist = sap_ppnl(1)%nlist_task(slot)%nlist
251 768158 : ilist = sap_ppnl(1)%nlist_task(slot)%ilist
252 768158 : nneighbor = sap_ppnl(1)%nlist_task(slot)%nnode
253 :
254 768158 : iac = ikind + nkind*(kkind - 1)
255 768158 : IF (.NOT. ASSOCIATED(basis_set(ikind)%gto_basis_set)) CYCLE
256 768158 : IF (.NOT. ASSOCIATED(gpotential(kkind)%gth_potential) .AND. &
257 : .NOT. ASSOCIATED(spotential(kkind)%sgp_potential)) CYCLE
258 768158 : IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) THEN
259 30826 : sap_int(iac)%a_kind = ikind
260 30826 : sap_int(iac)%p_kind = kkind
261 30826 : sap_int(iac)%nalist = nlist
262 154058 : ALLOCATE (sap_int(iac)%alist(nlist))
263 92406 : DO i = 1, nlist
264 61580 : NULLIFY (sap_int(iac)%alist(i)%clist)
265 61580 : sap_int(iac)%alist(i)%aatom = 0
266 92406 : sap_int(iac)%alist(i)%nclist = 0
267 : END DO
268 : END IF
269 782634 : IF (.NOT. ASSOCIATED(sap_int(iac)%alist(ilist)%clist)) THEN
270 61526 : sap_int(iac)%alist(ilist)%aatom = iatom
271 61526 : sap_int(iac)%alist(ilist)%nclist = nneighbor
272 1321892 : ALLOCATE (sap_int(iac)%alist(ilist)%clist(nneighbor))
273 829684 : DO i = 1, nneighbor
274 829684 : sap_int(iac)%alist(ilist)%clist(i)%catom = 0
275 : END DO
276 : END IF
277 : END DO
278 :
279 : ! Calculate the overlap integrals <a|p>
280 : !$OMP PARALLEL &
281 : !$OMP DEFAULT (NONE) &
282 : !$OMP SHARED (basis_set, gpotential, spotential, maxder, ncoset, &
283 : !$OMP sap_ppnl, sap_int, nkind, ldsab, ldai, nder, nco, do_soc ) &
284 : !$OMP PRIVATE (ikind, kkind, iatom, katom, nlist, ilist, nneighbor, jneighbor, &
285 : !$OMP cell_c, rac, iac, first_sgfa, la_max, la_min, npgfa, nseta, nsgfa, nsgf_seta, &
286 : !$OMP slot, sphi_a, zeta, cprj, hprj, lppnl, nppnl, nprj_ppnl, &
287 : !$OMP clist, iset, ncoa, sgfa, prjc, work, work_l, sab, lab, ai_work, nprjc, &
288 : !$OMP ppnl_radius, ncoc, rpgfa, first_col, vprj_ppnl, wprj_ppnl, i, j, l, do_gth, &
289 : !$OMP set_radius_a, rprjc, dac, lc_max, lc_min, zetc, alpha_ppnl, &
290 14476 : !$OMP na, nb, np, nnl, a_nl, radp, i_dim, ib, jb)
291 :
292 : ALLOCATE (sab(ldsab, ldsab*maxder), work(ldsab, ldsab*maxder))
293 : sab = 0.0_dp
294 : ALLOCATE (ai_work(ldai, ldai, ncoset(nder + 1)))
295 : ai_work = 0.0_dp
296 : IF (do_soc) THEN
297 : ALLOCATE (lab(ldsab, ldsab, 3), work_l(ldsab, ldsab, 3))
298 : lab = 0.0_dp
299 : END IF
300 :
301 : !$OMP DO SCHEDULE(GUIDED)
302 : DO slot = 1, sap_ppnl(1)%nl_size
303 :
304 : ikind = sap_ppnl(1)%nlist_task(slot)%ikind
305 : kkind = sap_ppnl(1)%nlist_task(slot)%jkind
306 : iatom = sap_ppnl(1)%nlist_task(slot)%iatom
307 : katom = sap_ppnl(1)%nlist_task(slot)%jatom
308 : nlist = sap_ppnl(1)%nlist_task(slot)%nlist
309 : ilist = sap_ppnl(1)%nlist_task(slot)%ilist
310 : nneighbor = sap_ppnl(1)%nlist_task(slot)%nnode
311 : jneighbor = sap_ppnl(1)%nlist_task(slot)%inode
312 : cell_c(:) = sap_ppnl(1)%nlist_task(slot)%cell(:)
313 : rac(1:3) = sap_ppnl(1)%nlist_task(slot)%r(1:3)
314 :
315 : iac = ikind + nkind*(kkind - 1)
316 : IF (.NOT. ASSOCIATED(basis_set(ikind)%gto_basis_set)) CYCLE
317 : ! Get definition of basis set
318 : first_sgfa => basis_set(ikind)%gto_basis_set%first_sgf
319 : la_max => basis_set(ikind)%gto_basis_set%lmax
320 : la_min => basis_set(ikind)%gto_basis_set%lmin
321 : npgfa => basis_set(ikind)%gto_basis_set%npgf
322 : nseta = basis_set(ikind)%gto_basis_set%nset
323 : nsgfa = basis_set(ikind)%gto_basis_set%nsgf
324 : nsgf_seta => basis_set(ikind)%gto_basis_set%nsgf_set
325 : rpgfa => basis_set(ikind)%gto_basis_set%pgf_radius
326 : set_radius_a => basis_set(ikind)%gto_basis_set%set_radius
327 : sphi_a => basis_set(ikind)%gto_basis_set%sphi
328 : zeta => basis_set(ikind)%gto_basis_set%zet
329 : ! Get definition of PP projectors
330 : IF (ASSOCIATED(gpotential(kkind)%gth_potential)) THEN
331 : ! GTH potential
332 : do_gth = .TRUE.
333 : alpha_ppnl => gpotential(kkind)%gth_potential%alpha_ppnl
334 : cprj => gpotential(kkind)%gth_potential%cprj
335 : lppnl = gpotential(kkind)%gth_potential%lppnl
336 : nppnl = gpotential(kkind)%gth_potential%nppnl
337 : nprj_ppnl => gpotential(kkind)%gth_potential%nprj_ppnl
338 : ppnl_radius = gpotential(kkind)%gth_potential%ppnl_radius
339 : vprj_ppnl => gpotential(kkind)%gth_potential%vprj_ppnl
340 : wprj_ppnl => gpotential(kkind)%gth_potential%wprj_ppnl
341 : ELSE IF (ASSOCIATED(spotential(kkind)%sgp_potential)) THEN
342 : ! SGP potential
343 : do_gth = .FALSE.
344 : nprjc = spotential(kkind)%sgp_potential%nppnl
345 : IF (nprjc == 0) CYCLE
346 : nnl = spotential(kkind)%sgp_potential%n_nonlocal
347 : lppnl = spotential(kkind)%sgp_potential%lmax
348 : a_nl => spotential(kkind)%sgp_potential%a_nonlocal
349 : ppnl_radius = spotential(kkind)%sgp_potential%ppnl_radius
350 : ALLOCATE (radp(nnl))
351 : radp(:) = ppnl_radius
352 : cprj => spotential(kkind)%sgp_potential%cprj_ppnl
353 : hprj => spotential(kkind)%sgp_potential%vprj_ppnl
354 : nppnl = SIZE(cprj, 2)
355 : ELSE
356 : CYCLE
357 : END IF
358 :
359 : dac = SQRT(SUM(rac*rac))
360 : clist => sap_int(iac)%alist(ilist)%clist(jneighbor)
361 : clist%catom = katom
362 : clist%cell = cell_c
363 : clist%rac = rac
364 : ALLOCATE (clist%acint(nsgfa, nppnl, maxder), &
365 : clist%achint(nsgfa, nppnl, maxder), &
366 : clist%alint(nsgfa, nppnl, 3), &
367 : clist%alkint(nsgfa, nppnl, 3))
368 : clist%acint = 0.0_dp
369 : clist%achint = 0.0_dp
370 : clist%alint = 0.0_dp
371 : clist%alkint = 0.0_dp
372 :
373 : clist%nsgf_cnt = 0
374 : NULLIFY (clist%sgf_list)
375 : DO iset = 1, nseta
376 : ncoa = npgfa(iset)*ncoset(la_max(iset))
377 : sgfa = first_sgfa(1, iset)
378 : IF (do_gth) THEN
379 : ! GTH potential
380 : prjc = 1
381 : work = 0.0_dp
382 : DO l = 0, lppnl
383 : nprjc = nprj_ppnl(l)*nco(l)
384 : IF (nprjc == 0) CYCLE
385 : rprjc(1) = ppnl_radius
386 : IF (set_radius_a(iset) + rprjc(1) < dac) CYCLE
387 : lc_max = l + 2*(nprj_ppnl(l) - 1)
388 : lc_min = l
389 : zetc(1) = alpha_ppnl(l)
390 : ncoc = ncoset(lc_max)
391 :
392 : ! Calculate the primitive overlap integrals
393 : CALL overlap(la_max(iset), la_min(iset), npgfa(iset), rpgfa(:, iset), zeta(:, iset), &
394 : lc_max, lc_min, 1, rprjc, zetc, rac, dac, sab, nder, .TRUE., ai_work, ldai)
395 : ! Transformation step projector functions (Cartesian -> spherical)
396 : na = ncoa
397 : nb = nprjc
398 : np = ncoc
399 : DO i = 1, maxder
400 : first_col = (i - 1)*ldsab
401 : ! CALL dgemm("N", "N", ncoa, nprjc, ncoc, 1.0_dp, sab(1, first_col + 1), SIZE(sab, 1), &
402 : ! cprj(1, prjc), SIZE(cprj, 1), 0.0_dp, work(1, first_col + prjc), ldsab)
403 : work(1:na, first_col + prjc:first_col + prjc + nb - 1) = &
404 : MATMUL(sab(1:na, first_col + 1:first_col + np), cprj(1:np, prjc:prjc + nb - 1))
405 : END DO
406 :
407 : IF (do_soc) THEN
408 : ! Calculate the primitive angular momentum integrals needed for spin-orbit coupling
409 : lab = 0.0_dp
410 : CALL angmom(la_max(iset), npgfa(iset), zeta(:, iset), rpgfa(:, iset), la_min(iset), &
411 : lc_max, 1, zetc, rprjc, -rac, (/0._dp, 0._dp, 0._dp/), lab)
412 : DO i_dim = 1, 3
413 : work_l(1:na, prjc:prjc + nb - 1, i_dim) = &
414 : MATMUL(lab(1:na, 1:np, i_dim), cprj(1:np, prjc:prjc + nb - 1))
415 : END DO
416 : END IF
417 :
418 : prjc = prjc + nprjc
419 :
420 : END DO
421 : na = nsgf_seta(iset)
422 : nb = nppnl
423 : np = ncoa
424 : DO i = 1, maxder
425 : first_col = (i - 1)*ldsab + 1
426 : ! Contraction step (basis functions)
427 : ! CALL dgemm("T", "N", nsgf_seta(iset), nppnl, ncoa, 1.0_dp, sphi_a(1, sgfa), SIZE(sphi_a, 1), &
428 : ! work(1, first_col), ldsab, 0.0_dp, clist%acint(sgfa, 1, i), nsgfa)
429 : clist%acint(sgfa:sgfa + na - 1, 1:nb, i) = &
430 : MATMUL(TRANSPOSE(sphi_a(1:np, sgfa:sgfa + na - 1)), work(1:np, first_col:first_col + nb - 1))
431 : ! Multiply with interaction matrix(h)
432 : ! CALL dgemm("N", "N", nsgf_seta(iset), nppnl, nppnl, 1.0_dp, clist%acint(sgfa, 1, i), nsgfa, &
433 : ! vprj_ppnl(1, 1), SIZE(vprj_ppnl, 1), 0.0_dp, clist%achint(sgfa, 1, i), nsgfa)
434 : clist%achint(sgfa:sgfa + na - 1, 1:nb, i) = &
435 : MATMUL(clist%acint(sgfa:sgfa + na - 1, 1:nb, i), vprj_ppnl(1:nb, 1:nb))
436 : END DO
437 : IF (do_soc) THEN
438 : DO i_dim = 1, 3
439 : clist%alint(sgfa:sgfa + na - 1, 1:nb, i_dim) = &
440 : MATMUL(TRANSPOSE(sphi_a(1:np, sgfa:sgfa + na - 1)), work_l(1:np, 1:nb, i_dim))
441 : clist%alkint(sgfa:sgfa + na - 1, 1:nb, i_dim) = &
442 : MATMUL(clist%alint(sgfa:sgfa + na - 1, 1:nb, i_dim), wprj_ppnl(1:nb, 1:nb))
443 : END DO
444 : END IF
445 : ELSE
446 : ! SGP potential
447 : ! Calculate the primitive overlap integrals
448 : CALL overlap(la_max(iset), la_min(iset), npgfa(iset), rpgfa(:, iset), zeta(:, iset), &
449 : lppnl, 0, nnl, radp, a_nl, rac, dac, sab, nder, .TRUE., ai_work, ldai)
450 : na = nsgf_seta(iset)
451 : nb = nppnl
452 : np = ncoa
453 : DO i = 1, maxder
454 : first_col = (i - 1)*ldsab + 1
455 : ! Transformation step projector functions (cartesian->spherical)
456 : ! CALL dgemm("N", "N", ncoa, nppnl, nprjc, 1.0_dp, sab(1, first_col), ldsab, &
457 : ! cprj(1, 1), SIZE(cprj, 1), 0.0_dp, work(1, 1), ldsab)
458 : work(1:np, 1:nb) = MATMUL(sab(1:np, first_col:first_col + nprjc - 1), cprj(1:nprjc, 1:nb))
459 : ! Contraction step (basis functions)
460 : ! CALL dgemm("T", "N", nsgf_seta(iset), nppnl, ncoa, 1.0_dp, sphi_a(1, sgfa), SIZE(sphi_a, 1), &
461 : ! work(1, 1), ldsab, 0.0_dp, clist%acint(sgfa, 1, i), nsgfa)
462 : clist%acint(sgfa:sgfa + na - 1, 1:nb, i) = &
463 : MATMUL(TRANSPOSE(sphi_a(1:np, sgfa:sgfa + na - 1)), work(1:np, 1:nb))
464 : ! *** Multiply with interaction matrix(h) ***
465 : ncoc = sgfa + nsgf_seta(iset) - 1
466 : DO j = 1, nppnl
467 : clist%achint(sgfa:ncoc, j, i) = clist%acint(sgfa:ncoc, j, i)*hprj(j)
468 : END DO
469 : END DO
470 : END IF
471 : END DO
472 : clist%maxac = MAXVAL(ABS(clist%acint(:, :, 1)))
473 : clist%maxach = MAXVAL(ABS(clist%achint(:, :, 1)))
474 : IF (.NOT. do_gth) DEALLOCATE (radp)
475 : END DO
476 :
477 : DEALLOCATE (sab, ai_work, work)
478 : IF (do_soc) DEALLOCATE (lab, work_l)
479 : !$OMP END PARALLEL
480 :
481 : ! Set up a sorting index
482 14476 : CALL sap_sort(sap_int)
483 : ! All integrals needed have been calculated and stored in sap_int
484 : ! We now calculate the Hamiltonian matrix elements
485 :
486 232484 : force_thread = 0.0_dp
487 68978 : at_thread = 0.0_dp
488 14476 : pv_thread = 0.0_dp
489 :
490 : !$OMP PARALLEL &
491 : !$OMP DEFAULT (NONE) &
492 : !$OMP SHARED (do_kp, basis_set, matrix_h, matrix_l, cell_to_index,&
493 : !$OMP sab_orb, matrix_p, sap_int, nkind, eps_ppnl, force, &
494 : !$OMP doat, do_dR, deltaR, maxder, nder, &
495 : !$OMP locks, virial, use_virial, calculate_forces, do_soc, natom) &
496 : !$OMP PRIVATE (ikind, jkind, iatom, jatom, cell_b, rab, &
497 : !$OMP slot, iab, atom_a, f0, irow, icol, h_block, &
498 : !$OMP l_block_x, l_block_y, l_block_z, &
499 : !$OMP r_2block, r_3block, atk, &
500 : !$OMP found,p_block, iac, ibc, alist_ac, alist_bc, acint, bcint, &
501 : !$OMP achint, bchint, alkint, blkint, &
502 : !$OMP na, np, nb, katom, j, fa, fb, rbc, rac, &
503 : !$OMP kkind, kac, kbc, i, img, hash, iatom8) &
504 14476 : !$OMP REDUCTION (+ : at_thread, pv_thread, force_thread )
505 :
506 : !$OMP SINGLE
507 : !$ ALLOCATE (locks(nlock))
508 : !$OMP END SINGLE
509 :
510 : !$OMP DO
511 : !$ DO lock_num = 1, nlock
512 : !$ call omp_init_lock(locks(lock_num))
513 : !$ END DO
514 : !$OMP END DO
515 :
516 : !$OMP DO SCHEDULE(GUIDED)
517 : DO slot = 1, sab_orb(1)%nl_size
518 :
519 : ikind = sab_orb(1)%nlist_task(slot)%ikind
520 : jkind = sab_orb(1)%nlist_task(slot)%jkind
521 : iatom = sab_orb(1)%nlist_task(slot)%iatom
522 : jatom = sab_orb(1)%nlist_task(slot)%jatom
523 : cell_b(:) = sab_orb(1)%nlist_task(slot)%cell(:)
524 : rab(1:3) = sab_orb(1)%nlist_task(slot)%r(1:3)
525 :
526 : IF (.NOT. ASSOCIATED(basis_set(ikind)%gto_basis_set)) CYCLE
527 : IF (.NOT. ASSOCIATED(basis_set(jkind)%gto_basis_set)) CYCLE
528 :
529 : iab = ikind + nkind*(jkind - 1)
530 :
531 : ! Use the symmetry of the first derivatives
532 : IF (iatom == jatom) THEN
533 : f0 = 1.0_dp
534 : ELSE
535 : f0 = 2.0_dp
536 : END IF
537 :
538 : IF (do_kp) THEN
539 : img = cell_to_index(cell_b(1), cell_b(2), cell_b(3))
540 : ELSE
541 : img = 1
542 : END IF
543 :
544 : ! Create matrix blocks for a new matrix block column
545 : IF (iatom <= jatom) THEN
546 : irow = iatom
547 : icol = jatom
548 : ELSE
549 : irow = jatom
550 : icol = iatom
551 : END IF
552 : NULLIFY (h_block)
553 : CALL dbcsr_get_block_p(matrix_h(1, img)%matrix, irow, icol, h_block, found)
554 : IF (do_soc) THEN
555 : NULLIFY (l_block_x, l_block_y, l_block_z)
556 : CALL dbcsr_get_block_p(matrix_l(1, img)%matrix, irow, icol, l_block_x, found)
557 : CALL dbcsr_get_block_p(matrix_l(2, img)%matrix, irow, icol, l_block_y, found)
558 : CALL dbcsr_get_block_p(matrix_l(3, img)%matrix, irow, icol, l_block_z, found)
559 : END IF
560 :
561 : IF (do_dR) THEN
562 : NULLIFY (r_2block, r_3block)
563 : CALL dbcsr_get_block_p(matrix_h(2, img)%matrix, irow, icol, r_2block, found)
564 : CALL dbcsr_get_block_p(matrix_h(3, img)%matrix, irow, icol, r_3block, found)
565 : END IF
566 :
567 : IF (calculate_forces .OR. doat) THEN
568 : NULLIFY (p_block)
569 : CALL dbcsr_get_block_p(matrix_p(1, img)%matrix, irow, icol, p_block, found)
570 : END IF
571 :
572 : ! loop over all kinds for projector atom
573 : IF (ASSOCIATED(h_block)) THEN
574 : !$ iatom8 = INT(iatom - 1, int_8)*INT(natom, int_8) + INT(jatom, int_8)
575 : !$ hash = INT(MOD(iatom8, INT(nlock, int_8)) + 1)
576 :
577 : DO kkind = 1, nkind
578 : iac = ikind + nkind*(kkind - 1)
579 : ibc = jkind + nkind*(kkind - 1)
580 : IF (.NOT. ASSOCIATED(sap_int(iac)%alist)) CYCLE
581 : IF (.NOT. ASSOCIATED(sap_int(ibc)%alist)) CYCLE
582 : CALL get_alist(sap_int(iac), alist_ac, iatom)
583 : CALL get_alist(sap_int(ibc), alist_bc, jatom)
584 : IF (.NOT. ASSOCIATED(alist_ac)) CYCLE
585 : IF (.NOT. ASSOCIATED(alist_bc)) CYCLE
586 : DO kac = 1, alist_ac%nclist
587 : DO kbc = 1, alist_bc%nclist
588 : IF (alist_ac%clist(kac)%catom /= alist_bc%clist(kbc)%catom) CYCLE
589 : IF (ALL(cell_b + alist_bc%clist(kbc)%cell - alist_ac%clist(kac)%cell == 0)) THEN
590 : IF (alist_ac%clist(kac)%maxac*alist_bc%clist(kbc)%maxach < eps_ppnl) CYCLE
591 : acint => alist_ac%clist(kac)%acint
592 : bcint => alist_bc%clist(kbc)%acint
593 : achint => alist_ac%clist(kac)%achint
594 : bchint => alist_bc%clist(kbc)%achint
595 : IF (do_soc) THEN
596 : alkint => alist_ac%clist(kac)%alkint
597 : blkint => alist_bc%clist(kbc)%alkint
598 : END IF
599 : na = SIZE(acint, 1)
600 : np = SIZE(acint, 2)
601 : nb = SIZE(bcint, 1)
602 : !$ CALL omp_set_lock(locks(hash))
603 : IF (.NOT. do_dR) THEN
604 : IF (iatom <= jatom) THEN
605 : h_block(1:na, 1:nb) = h_block(1:na, 1:nb) + &
606 : MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, 1)))
607 : ELSE
608 : h_block(1:nb, 1:na) = h_block(1:nb, 1:na) + &
609 : MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, 1)))
610 : END IF
611 : END IF
612 : IF (do_soc) THEN
613 : IF (iatom <= jatom) THEN
614 : l_block_x(1:na, 1:nb) = l_block_x(1:na, 1:nb) + &
615 : MATMUL(alkint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, 1)))
616 : l_block_y(1:na, 1:nb) = l_block_y(1:na, 1:nb) + &
617 : MATMUL(alkint(1:na, 1:np, 2), TRANSPOSE(bcint(1:nb, 1:np, 1)))
618 : l_block_z(1:na, 1:nb) = l_block_z(1:na, 1:nb) + &
619 : MATMUL(alkint(1:na, 1:np, 3), TRANSPOSE(bcint(1:nb, 1:np, 1)))
620 :
621 : ELSE
622 : l_block_x(1:nb, 1:na) = l_block_x(1:nb, 1:na) - &
623 : MATMUL(blkint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, 1)))
624 : l_block_y(1:nb, 1:na) = l_block_y(1:nb, 1:na) - &
625 : MATMUL(blkint(1:nb, 1:np, 2), TRANSPOSE(acint(1:na, 1:np, 1)))
626 : l_block_z(1:nb, 1:na) = l_block_z(1:nb, 1:na) - &
627 : MATMUL(blkint(1:nb, 1:np, 3), TRANSPOSE(acint(1:na, 1:np, 1)))
628 : END IF
629 : END IF
630 : !$ CALL omp_unset_lock(locks(hash))
631 : IF (calculate_forces) THEN
632 : IF (ASSOCIATED(p_block)) THEN
633 : katom = alist_ac%clist(kac)%catom
634 : DO i = 1, 3
635 : j = i + 1
636 : IF (iatom <= jatom) THEN
637 : fa(i) = SUM(p_block(1:na, 1:nb)* &
638 : MATMUL(acint(1:na, 1:np, j), TRANSPOSE(bchint(1:nb, 1:np, 1))))
639 : fb(i) = SUM(p_block(1:na, 1:nb)* &
640 : MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, j))))
641 : ELSE
642 : fa(i) = SUM(p_block(1:nb, 1:na)* &
643 : MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, j))))
644 : fb(i) = SUM(p_block(1:nb, 1:na)* &
645 : MATMUL(bcint(1:nb, 1:np, j), TRANSPOSE(achint(1:na, 1:np, 1))))
646 : END IF
647 : force_thread(i, iatom) = force_thread(i, iatom) + f0*fa(i)
648 : force_thread(i, katom) = force_thread(i, katom) - f0*fa(i)
649 : force_thread(i, jatom) = force_thread(i, jatom) + f0*fb(i)
650 : force_thread(i, katom) = force_thread(i, katom) - f0*fb(i)
651 : END DO
652 :
653 : IF (use_virial) THEN
654 : rac = alist_ac%clist(kac)%rac
655 : rbc = alist_bc%clist(kbc)%rac
656 : CALL virial_pair_force(pv_thread, f0, fa, rac)
657 : CALL virial_pair_force(pv_thread, f0, fb, rbc)
658 : END IF
659 : END IF
660 : END IF
661 :
662 : IF (do_dR) THEN
663 : i = 1; j = 2;
664 : katom = alist_ac%clist(kac)%catom
665 : IF (iatom <= jatom) THEN
666 : h_block(1:na, 1:nb) = h_block(1:na, 1:nb) + &
667 : (deltaR(i, iatom) - deltaR(i, katom))* &
668 : MATMUL(acint(1:na, 1:np, j), TRANSPOSE(bchint(1:nb, 1:np, 1)))
669 :
670 : h_block(1:na, 1:nb) = h_block(1:na, 1:nb) + &
671 : (deltaR(i, jatom) - deltaR(i, katom))* &
672 : MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, j)))
673 : ELSE
674 : h_block(1:nb, 1:na) = h_block(1:nb, 1:na) + &
675 : (deltaR(i, iatom) - deltaR(i, katom))* &
676 : MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, j)))
677 : h_block(1:nb, 1:na) = h_block(1:nb, 1:na) + &
678 : (deltaR(i, jatom) - deltaR(i, katom))* &
679 : MATMUL(bcint(1:nb, 1:np, j), TRANSPOSE(achint(1:na, 1:np, 1)))
680 : END IF
681 :
682 : i = 2; j = 3;
683 : katom = alist_ac%clist(kac)%catom
684 : IF (iatom <= jatom) THEN
685 : r_2block(1:na, 1:nb) = r_2block(1:na, 1:nb) + &
686 : (deltaR(i, iatom) - deltaR(i, katom))* &
687 : MATMUL(acint(1:na, 1:np, j), TRANSPOSE(bchint(1:nb, 1:np, 1)))
688 :
689 : r_2block(1:na, 1:nb) = r_2block(1:na, 1:nb) + &
690 : (deltaR(i, jatom) - deltaR(i, katom))* &
691 : MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, j)))
692 : ELSE
693 : r_2block(1:nb, 1:na) = r_2block(1:nb, 1:na) + &
694 : (deltaR(i, iatom) - deltaR(i, katom))* &
695 : MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, j)))
696 : r_2block(1:nb, 1:na) = r_2block(1:nb, 1:na) + &
697 : (deltaR(i, jatom) - deltaR(i, katom))* &
698 : MATMUL(bcint(1:nb, 1:np, j), TRANSPOSE(achint(1:na, 1:np, 1)))
699 : END IF
700 :
701 : i = 3; j = 4;
702 : katom = alist_ac%clist(kac)%catom
703 : IF (iatom <= jatom) THEN
704 : r_3block(1:na, 1:nb) = r_3block(1:na, 1:nb) + &
705 : (deltaR(i, iatom) - deltaR(i, katom))* &
706 : MATMUL(acint(1:na, 1:np, j), TRANSPOSE(bchint(1:nb, 1:np, 1)))
707 :
708 : r_3block(1:na, 1:nb) = r_3block(1:na, 1:nb) + &
709 : (deltaR(i, jatom) - deltaR(i, katom))* &
710 : MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, j)))
711 : ELSE
712 : r_3block(1:nb, 1:na) = r_3block(1:nb, 1:na) + &
713 : (deltaR(i, iatom) - deltaR(i, katom))* &
714 : MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, j)))
715 : r_3block(1:nb, 1:na) = r_3block(1:nb, 1:na) + &
716 : (deltaR(i, jatom) - deltaR(i, katom))* &
717 : MATMUL(bcint(1:nb, 1:np, j), TRANSPOSE(achint(1:na, 1:np, 1)))
718 : END IF
719 :
720 : END IF
721 : IF (doat) THEN
722 : IF (ASSOCIATED(p_block)) THEN
723 : katom = alist_ac%clist(kac)%catom
724 : IF (iatom <= jatom) THEN
725 : atk = SUM(p_block(1:na, 1:nb)* &
726 : MATMUL(achint(1:na, 1:np, 1), TRANSPOSE(bcint(1:nb, 1:np, 1))))
727 : ELSE
728 : atk = SUM(p_block(1:nb, 1:na)* &
729 : MATMUL(bchint(1:nb, 1:np, 1), TRANSPOSE(acint(1:na, 1:np, 1))))
730 : END IF
731 : at_thread(katom) = at_thread(katom) + f0*atk
732 : END IF
733 : END IF
734 : EXIT ! We have found a match and there can be only one single match
735 : END IF
736 : END DO
737 : END DO
738 : END DO
739 : END IF
740 : END DO
741 :
742 : !$OMP DO
743 : !$ DO lock_num = 1, nlock
744 : !$ call omp_destroy_lock(locks(lock_num))
745 : !$ END DO
746 : !$OMP END DO
747 :
748 : !$OMP SINGLE
749 : !$ DEALLOCATE (locks)
750 : !$OMP END SINGLE NOWAIT
751 :
752 : !$OMP END PARALLEL
753 :
754 14476 : CALL release_sap_int(sap_int)
755 :
756 14476 : DEALLOCATE (basis_set, gpotential, spotential)
757 14476 : IF (calculate_forces) THEN
758 6171 : CALL get_atomic_kind_set(atomic_kind_set, atom_of_kind=atom_of_kind, kind_of=kind_of)
759 : !$OMP DO
760 : DO iatom = 1, natom
761 22257 : atom_a = atom_of_kind(iatom)
762 22257 : ikind = kind_of(iatom)
763 89028 : force(ikind)%gth_ppnl(:, atom_a) = force(ikind)%gth_ppnl(:, atom_a) + force_thread(:, iatom)
764 : END DO
765 : !$OMP END DO
766 6171 : DEALLOCATE (atom_of_kind, kind_of)
767 : END IF
768 :
769 14476 : IF (calculate_forces .AND. use_virial) THEN
770 10296 : virial%pv_ppnl = virial%pv_ppnl + pv_thread
771 10296 : virial%pv_virial = virial%pv_virial + pv_thread
772 : END IF
773 :
774 14476 : IF (doat) THEN
775 280 : atcore(1:natom) = atcore(1:natom) + at_thread
776 : END IF
777 :
778 28952 : IF (calculate_forces .OR. doat) THEN
779 : ! If LSD, then recover alpha density and beta density
780 : ! from the total density (1) and the spin density (2)
781 6233 : IF (SIZE(matrix_p, 1) == 2) THEN
782 1894 : DO img = 1, nimages
783 : CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
784 1230 : alpha_scalar=0.5_dp, beta_scalar=0.5_dp)
785 : CALL dbcsr_add(matrix_p(2, img)%matrix, matrix_p(1, img)%matrix, &
786 1894 : alpha_scalar=-1.0_dp, beta_scalar=1.0_dp)
787 : END DO
788 : END IF
789 : END IF
790 :
791 : END IF !ppnl_present
792 :
793 14476 : CALL timestop(handle)
794 :
795 28952 : END SUBROUTINE build_core_ppnl
796 :
797 : END MODULE core_ppnl
|