Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 : ! **************************************************************************************************
8 : MODULE qs_rho0_ggrid
9 : USE atomic_kind_types, ONLY: atomic_kind_type,&
10 : get_atomic_kind
11 : USE basis_set_types, ONLY: get_gto_basis_set,&
12 : gto_basis_set_type
13 : USE cell_types, ONLY: cell_type,&
14 : pbc
15 : USE cp_control_types, ONLY: dft_control_type
16 : USE gaussian_gridlevels, ONLY: gaussian_gridlevel
17 : USE grid_api, ONLY: GRID_FUNC_AB,&
18 : collocate_pgf_product
19 : USE kinds, ONLY: dp
20 : USE memory_utilities, ONLY: reallocate
21 : USE message_passing, ONLY: mp_para_env_type
22 : USE orbital_pointers, ONLY: indco,&
23 : nco,&
24 : ncoset,&
25 : nso,&
26 : nsoset
27 : USE orbital_transformation_matrices, ONLY: orbtramat
28 : USE particle_types, ONLY: particle_type
29 : USE pw_env_types, ONLY: pw_env_get,&
30 : pw_env_type
31 : USE pw_methods, ONLY: pw_axpy,&
32 : pw_copy,&
33 : pw_integrate_function,&
34 : pw_scale,&
35 : pw_transfer,&
36 : pw_zero
37 : USE pw_pool_types, ONLY: pw_pool_p_type,&
38 : pw_pool_type
39 : USE pw_types, ONLY: pw_c1d_gs_type,&
40 : pw_r3d_rs_type
41 : USE qs_cneo_ggrid, ONLY: integrate_vhgg_rspace,&
42 : rhoz_cneo_s_grid_create
43 : USE qs_cneo_types, ONLY: cneo_potential_type,&
44 : rhoz_cneo_type
45 : USE qs_environment_types, ONLY: get_qs_env,&
46 : qs_environment_type
47 : USE qs_force_types, ONLY: qs_force_type
48 : USE qs_harmonics_atom, ONLY: get_none0_cg_list,&
49 : harmonics_atom_type
50 : USE qs_integrate_potential, ONLY: integrate_pgf_product
51 : USE qs_kind_types, ONLY: get_qs_kind,&
52 : qs_kind_type
53 : USE qs_local_rho_types, ONLY: get_local_rho,&
54 : local_rho_type
55 : USE qs_rho0_types, ONLY: get_rho0_mpole,&
56 : rho0_mpole_type
57 : USE qs_rho_atom_types, ONLY: get_rho_atom,&
58 : rho_atom_coeff,&
59 : rho_atom_type
60 : USE realspace_grid_types, ONLY: realspace_grid_desc_p_type,&
61 : realspace_grid_desc_type,&
62 : realspace_grid_type,&
63 : rs_grid_create,&
64 : rs_grid_release,&
65 : rs_grid_zero,&
66 : transfer_pw2rs,&
67 : transfer_rs2pw
68 : USE util, ONLY: get_limit
69 : USE virial_types, ONLY: virial_type
70 : #include "./base/base_uses.f90"
71 :
72 : IMPLICIT NONE
73 :
74 : PRIVATE
75 :
76 : ! Global parameters (only in this module)
77 :
78 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_rho0_ggrid'
79 :
80 : ! Public subroutines
81 :
82 : PUBLIC :: put_rho0_on_grid, rho0_s_grid_create, integrate_vhg0_rspace
83 :
84 : CONTAINS
85 :
86 : ! **************************************************************************************************
87 : !> \brief ...
88 : !> \param qs_env ...
89 : !> \param rho0 ...
90 : !> \param tot_rs_int ...
91 : !> \param my_pools ...
92 : !> \param my_rs_grids ...
93 : !> \param my_rs_descs ...
94 : ! **************************************************************************************************
95 18394 : SUBROUTINE put_rho0_on_grid(qs_env, rho0, tot_rs_int, my_pools, my_rs_grids, my_rs_descs)
96 :
97 : TYPE(qs_environment_type), POINTER :: qs_env
98 : TYPE(rho0_mpole_type), POINTER :: rho0
99 : REAL(KIND=dp), INTENT(OUT) :: tot_rs_int
100 : TYPE(pw_pool_p_type), DIMENSION(:), OPTIONAL, &
101 : POINTER :: my_pools
102 : TYPE(realspace_grid_type), DIMENSION(:), &
103 : OPTIONAL, POINTER :: my_rs_grids
104 : TYPE(realspace_grid_desc_p_type), DIMENSION(:), &
105 : OPTIONAL, POINTER :: my_rs_descs
106 :
107 : CHARACTER(LEN=*), PARAMETER :: routineN = 'put_rho0_on_grid'
108 :
109 : INTEGER :: auxbas_grid, handle, iat, iatom, igrid, &
110 : ikind, ithread, j, l0_ikind, lmax0, &
111 : nat, nch_ik, nch_max, npme
112 18394 : INTEGER, DIMENSION(:), POINTER :: atom_list, cores
113 : LOGICAL :: paw_atom
114 : REAL(KIND=dp) :: eps_rho_rspace, rpgf0, zet0
115 : REAL(KIND=dp), DIMENSION(3) :: ra
116 18394 : REAL(KIND=dp), DIMENSION(:), POINTER :: Qlm_c
117 18394 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: pab
118 18394 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
119 : TYPE(cell_type), POINTER :: cell
120 : TYPE(dft_control_type), POINTER :: dft_control
121 : TYPE(mp_para_env_type), POINTER :: para_env
122 18394 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
123 : TYPE(pw_c1d_gs_type) :: coeff_gspace
124 : TYPE(pw_c1d_gs_type), POINTER :: rho0_s_gs
125 : TYPE(pw_env_type), POINTER :: pw_env
126 18394 : TYPE(pw_pool_p_type), DIMENSION(:), POINTER :: pw_pools
127 : TYPE(pw_pool_type), POINTER :: pw_pool
128 : TYPE(pw_r3d_rs_type) :: coeff_rspace, rho0_r_tmp
129 : TYPE(pw_r3d_rs_type), POINTER :: rho0_s_rs
130 18394 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
131 : TYPE(realspace_grid_desc_p_type), DIMENSION(:), &
132 18394 : POINTER :: descs
133 : TYPE(realspace_grid_desc_type), POINTER :: desc
134 18394 : TYPE(realspace_grid_type), DIMENSION(:), POINTER :: grids
135 : TYPE(realspace_grid_type), POINTER :: rs_grid
136 :
137 18394 : CALL timeset(routineN, handle)
138 :
139 18394 : NULLIFY (atomic_kind_set, qs_kind_set, cores, pab, Qlm_c)
140 :
141 18394 : NULLIFY (dft_control, pw_env, particle_set, para_env, cell, rho0_s_gs, rho0_s_rs)
142 : CALL get_qs_env(qs_env=qs_env, dft_control=dft_control, &
143 : particle_set=particle_set, &
144 : atomic_kind_set=atomic_kind_set, &
145 : qs_kind_set=qs_kind_set, &
146 : para_env=para_env, &
147 18394 : pw_env=pw_env, cell=cell)
148 18394 : eps_rho_rspace = dft_control%qs_control%eps_rho_rspace
149 :
150 18394 : NULLIFY (descs, pw_pools)
151 18394 : CALL pw_env_get(pw_env=pw_env, rs_descs=descs, rs_grids=grids, pw_pools=pw_pools)
152 18394 : auxbas_grid = pw_env%auxbas_grid
153 :
154 18394 : NULLIFY (rho0_s_gs, rho0_s_rs)
155 : CALL get_rho0_mpole(rho0_mpole=rho0, lmax_0=lmax0, &
156 : zet0_h=zet0, igrid_zet0_s=igrid, &
157 : rho0_s_gs=rho0_s_gs, &
158 18394 : rho0_s_rs=rho0_s_rs)
159 :
160 : ! *** set up the rs grid at level igrid
161 18394 : NULLIFY (rs_grid, desc, pw_pool)
162 : ! IF present, overwrite qs grid for new pool
163 18394 : IF (PRESENT(my_pools)) THEN
164 1304 : desc => my_rs_descs(igrid)%rs_desc
165 1304 : rs_grid => my_rs_grids(igrid)
166 1304 : pw_pool => my_pools(igrid)%pool
167 : ELSE
168 17090 : desc => descs(igrid)%rs_desc
169 17090 : rs_grid => grids(igrid)
170 17090 : pw_pool => pw_pools(igrid)%pool
171 : END IF
172 :
173 18394 : CPASSERT(ASSOCIATED(desc))
174 18394 : CPASSERT(ASSOCIATED(pw_pool))
175 :
176 18394 : IF (igrid /= auxbas_grid) THEN
177 12 : CALL pw_pool%create_pw(coeff_rspace)
178 12 : CALL pw_pool%create_pw(coeff_gspace)
179 : END IF
180 18394 : CALL rs_grid_zero(rs_grid)
181 :
182 18394 : nch_max = ncoset(lmax0)
183 :
184 55182 : ALLOCATE (pab(nch_max, 1))
185 :
186 54702 : DO ikind = 1, SIZE(atomic_kind_set)
187 36308 : CALL get_atomic_kind(atomic_kind_set(ikind), atom_list=atom_list, natom=nat)
188 36308 : CALL get_qs_kind(qs_kind_set(ikind), paw_atom=paw_atom)
189 :
190 36308 : IF (.NOT. paw_atom .AND. dft_control%qs_control%gapw_control%nopaw_as_gpw) CYCLE
191 :
192 : CALL get_rho0_mpole(rho0_mpole=rho0, ikind=ikind, l0_ikind=l0_ikind, &
193 34668 : rpgf0_s=rpgf0)
194 :
195 34668 : nch_ik = ncoset(l0_ikind)
196 430072 : pab = 0.0_dp
197 :
198 34668 : CALL reallocate(cores, 1, nat)
199 34668 : npme = 0
200 90532 : cores = 0
201 :
202 90532 : DO iat = 1, nat
203 55864 : iatom = atom_list(iat)
204 55864 : ra(:) = pbc(particle_set(iatom)%r, cell)
205 90532 : IF (rs_grid%desc%parallel .AND. .NOT. rs_grid%desc%distributed) THEN
206 : ! replicated realspace grid, split the atoms up between procs
207 55864 : IF (MODULO(nat, rs_grid%desc%group_size) == rs_grid%desc%my_pos) THEN
208 27932 : npme = npme + 1
209 27932 : cores(npme) = iat
210 : END IF
211 : ELSE
212 0 : npme = npme + 1
213 0 : cores(npme) = iat
214 : END IF
215 :
216 : END DO
217 :
218 : ithread = 0
219 151970 : DO j = 1, npme
220 :
221 27932 : iat = cores(j)
222 27932 : iatom = atom_list(iat)
223 :
224 27932 : CALL get_rho0_mpole(rho0_mpole=rho0, iat=iatom, Qlm_car=Qlm_c)
225 :
226 522822 : pab(1:nch_ik, 1) = Qlm_c(1:nch_ik)
227 :
228 27932 : ra(:) = pbc(particle_set(iatom)%r, cell)
229 :
230 : CALL collocate_pgf_product( &
231 : l0_ikind, zet0, 0, 0, 0.0_dp, 0, &
232 : ra, [0.0_dp, 0.0_dp, 0.0_dp], 1.0_dp, pab, 0, 0, &
233 : rs_grid, ga_gb_function=GRID_FUNC_AB, radius=rpgf0, &
234 64240 : use_subpatch=.TRUE., subpatch_pattern=0)
235 :
236 : END DO ! j
237 :
238 : END DO ! ikind
239 :
240 18394 : IF (ASSOCIATED(cores)) THEN
241 18250 : DEALLOCATE (cores)
242 : END IF
243 :
244 18394 : DEALLOCATE (pab)
245 :
246 18394 : IF (igrid /= auxbas_grid) THEN
247 12 : CALL transfer_rs2pw(rs_grid, coeff_rspace)
248 12 : CALL pw_zero(rho0_s_gs)
249 12 : CALL pw_transfer(coeff_rspace, coeff_gspace)
250 12 : CALL pw_axpy(coeff_gspace, rho0_s_gs)
251 :
252 12 : tot_rs_int = pw_integrate_function(coeff_rspace, isign=-1)
253 :
254 12 : CALL pw_pool%give_back_pw(coeff_rspace)
255 12 : CALL pw_pool%give_back_pw(coeff_gspace)
256 : ELSE
257 :
258 18382 : CALL pw_pool%create_pw(rho0_r_tmp)
259 :
260 18382 : CALL transfer_rs2pw(rs_grid, rho0_r_tmp)
261 :
262 18382 : tot_rs_int = pw_integrate_function(rho0_r_tmp, isign=-1)
263 :
264 18382 : CALL pw_transfer(rho0_r_tmp, rho0_s_rs)
265 18382 : CALL pw_pool%give_back_pw(rho0_r_tmp)
266 :
267 18382 : CALL pw_zero(rho0_s_gs)
268 18382 : CALL pw_transfer(rho0_s_rs, rho0_s_gs)
269 : END IF
270 18394 : CALL timestop(handle)
271 :
272 36788 : END SUBROUTINE put_rho0_on_grid
273 :
274 : ! **************************************************************************************************
275 : !> \brief ...
276 : !> \param pw_env ...
277 : !> \param rho0_mpole ...
278 : ! **************************************************************************************************
279 3052 : SUBROUTINE rho0_s_grid_create(pw_env, rho0_mpole)
280 :
281 : TYPE(pw_env_type), POINTER :: pw_env
282 : TYPE(rho0_mpole_type), POINTER :: rho0_mpole
283 :
284 : CHARACTER(len=*), PARAMETER :: routineN = 'rho0_s_grid_create'
285 :
286 : INTEGER :: handle
287 : TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
288 :
289 3052 : CALL timeset(routineN, handle)
290 :
291 3052 : CPASSERT(ASSOCIATED(pw_env))
292 :
293 3052 : NULLIFY (auxbas_pw_pool)
294 3052 : CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
295 3052 : CPASSERT(ASSOCIATED(auxbas_pw_pool))
296 :
297 : ! reallocate rho0 on the global grid in real and reciprocal space
298 3052 : CPASSERT(ASSOCIATED(rho0_mpole))
299 :
300 : ! rho0 density in real space
301 3052 : IF (ASSOCIATED(rho0_mpole%rho0_s_rs)) THEN
302 884 : CALL rho0_mpole%rho0_s_rs%release()
303 : ELSE
304 2168 : ALLOCATE (rho0_mpole%rho0_s_rs)
305 : END IF
306 3052 : CALL auxbas_pw_pool%create_pw(rho0_mpole%rho0_s_rs)
307 :
308 : ! rho0 density in reciprocal space
309 3052 : IF (ASSOCIATED(rho0_mpole%rho0_s_gs)) THEN
310 884 : CALL rho0_mpole%rho0_s_gs%release()
311 : ELSE
312 2168 : ALLOCATE (rho0_mpole%rho0_s_gs)
313 : END IF
314 3052 : CALL auxbas_pw_pool%create_pw(rho0_mpole%rho0_s_gs)
315 :
316 : ! Find the grid level suitable for rho0_soft
317 3052 : rho0_mpole%igrid_zet0_s = gaussian_gridlevel(pw_env%gridlevel_info, 2.0_dp*rho0_mpole%zet0_h)
318 :
319 3052 : CALL timestop(handle)
320 :
321 3052 : IF (rho0_mpole%do_cneo) THEN
322 16 : CALL rhoz_cneo_s_grid_create(pw_env, rho0_mpole)
323 : END IF
324 :
325 3052 : END SUBROUTINE rho0_s_grid_create
326 :
327 : ! **************************************************************************************************
328 : !> \brief ...
329 : !> \param qs_env ...
330 : !> \param v_rspace ...
331 : !> \param para_env ...
332 : !> \param calculate_forces ...
333 : !> \param local_rho_set ...
334 : !> \param local_rho_set_2nd ...
335 : !> \param atener ...
336 : !> \param kforce ...
337 : !> \param my_pools ...
338 : !> \param my_rs_descs ...
339 : ! **************************************************************************************************
340 17104 : SUBROUTINE integrate_vhg0_rspace(qs_env, v_rspace, para_env, calculate_forces, local_rho_set, &
341 17104 : local_rho_set_2nd, atener, kforce, my_pools, my_rs_descs)
342 :
343 : TYPE(qs_environment_type), POINTER :: qs_env
344 : TYPE(pw_r3d_rs_type), INTENT(IN) :: v_rspace
345 : TYPE(mp_para_env_type), POINTER :: para_env
346 : LOGICAL, INTENT(IN) :: calculate_forces
347 : TYPE(local_rho_type), OPTIONAL, POINTER :: local_rho_set, local_rho_set_2nd
348 : REAL(KIND=dp), DIMENSION(:), OPTIONAL :: atener
349 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: kforce
350 : TYPE(pw_pool_p_type), DIMENSION(:), OPTIONAL, &
351 : POINTER :: my_pools
352 : TYPE(realspace_grid_desc_p_type), DIMENSION(:), &
353 : OPTIONAL, POINTER :: my_rs_descs
354 :
355 : CHARACTER(LEN=*), PARAMETER :: routineN = 'integrate_vhg0_rspace'
356 :
357 : INTEGER :: auxbas_grid, bo(2), handle, i, iat, iatom, ic, icg, ico, ig1, ig2, igrid, ii, &
358 : ikind, ipgf1, ipgf2, is, iset1, iset2, iso, iso1, iso2, ispin, j, l0_ikind, llmax, &
359 : llmax_nuc, lmax0, lshell, lx, ly, lz, m1, m2, max_iso_not0_local, max_s_harm, &
360 : max_s_harm_nuc, maxl, maxl_nuc, maxso, maxso_nuc, mepos, n1, n2, nat, nch_ik, nch_max, &
361 : ncurr, nset, nset_nuc, nsotot, nsotot_nuc, nspins, num_pe
362 17104 : INTEGER, ALLOCATABLE, DIMENSION(:) :: cg_n_list
363 17104 : INTEGER, ALLOCATABLE, DIMENSION(:, :, :) :: cg_list
364 17104 : INTEGER, DIMENSION(:), POINTER :: atom_list, lmax, lmax_nuc, lmin, &
365 17104 : lmin_nuc, npgf, npgf_nuc
366 : LOGICAL :: grid_distributed, paw_atom, use_virial
367 : REAL(KIND=dp) :: eps_rho_rspace, force_tmp(3), fscale, &
368 : ra(3), rpgf0, zet0
369 : REAL(KIND=dp), DIMENSION(3, 3) :: my_virial_a, my_virial_b
370 17104 : REAL(KIND=dp), DIMENSION(:), POINTER :: hab_sph, norm_l, Qlm
371 17104 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: hab, hdab_sph, intloc, intloc_nuc, pab
372 17104 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: a_hdab_sph, hdab, Qlm_gg
373 17104 : REAL(KIND=dp), DIMENSION(:, :, :, :), POINTER :: a_hdab
374 17104 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
375 : TYPE(cell_type), POINTER :: cell
376 : TYPE(cneo_potential_type), POINTER :: cneo_potential
377 : TYPE(dft_control_type), POINTER :: dft_control
378 : TYPE(gto_basis_set_type), POINTER :: basis_1c_set, nuc_basis_set
379 : TYPE(harmonics_atom_type), POINTER :: harmonics
380 17104 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
381 : TYPE(pw_c1d_gs_type) :: coeff_gaux, coeff_gspace
382 : TYPE(pw_env_type), POINTER :: pw_env
383 17104 : TYPE(pw_pool_p_type), DIMENSION(:), POINTER :: pw_pools
384 : TYPE(pw_pool_type), POINTER :: pw_aux, pw_pool
385 : TYPE(pw_r3d_rs_type) :: coeff_raux, coeff_rspace
386 17104 : TYPE(qs_force_type), DIMENSION(:), POINTER :: force
387 17104 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
388 : TYPE(realspace_grid_desc_p_type), DIMENSION(:), &
389 17104 : POINTER :: rs_descs
390 : TYPE(realspace_grid_desc_type), POINTER :: rs_desc
391 342080 : TYPE(realspace_grid_type) :: rs_v
392 : TYPE(rho0_mpole_type), POINTER :: rho0_mpole
393 17104 : TYPE(rho_atom_coeff), DIMENSION(:), POINTER :: int_local_h, int_local_s
394 17104 : TYPE(rho_atom_type), DIMENSION(:), POINTER :: rho_atom_set
395 : TYPE(rho_atom_type), POINTER :: rho_atom
396 17104 : TYPE(rhoz_cneo_type), DIMENSION(:), POINTER :: rhoz_cneo_set
397 : TYPE(rhoz_cneo_type), POINTER :: rhoz_cneo
398 : TYPE(virial_type), POINTER :: virial
399 :
400 17104 : CALL timeset(routineN, handle)
401 :
402 : ! In case of the external density computed forces probably also
403 : ! need to be stored outside qs_env. We can then remove the
404 : ! attribute 'OPTIONAL' from the argument 'local_rho_set'.
405 :
406 : ! CPASSERT(.NOT. (calculate_forces .AND. PRESENT(local_rho_set)))
407 : ! IF (calculate_forces .AND. PRESENT(local_rho_set)) THEN
408 : ! CPWARN("Forces and External Density!")
409 : ! END IF
410 :
411 17104 : NULLIFY (atomic_kind_set, qs_kind_set, dft_control, particle_set)
412 17104 : NULLIFY (cell, force, pw_env, rho0_mpole, rho_atom_set, rhoz_cneo_set)
413 :
414 : CALL get_qs_env(qs_env=qs_env, &
415 : atomic_kind_set=atomic_kind_set, &
416 : qs_kind_set=qs_kind_set, &
417 : cell=cell, &
418 : dft_control=dft_control, &
419 : force=force, pw_env=pw_env, &
420 : rho0_mpole=rho0_mpole, &
421 : rho_atom_set=rho_atom_set, &
422 : rhoz_cneo_set=rhoz_cneo_set, &
423 : particle_set=particle_set, &
424 17104 : virial=virial)
425 :
426 17104 : use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
427 :
428 17104 : nspins = dft_control%nspins
429 :
430 : ! The aim of the following code was to return immediately if the subroutine
431 : ! was called for triplet excited states in spin-restricted case. This check
432 : ! is also performed before invocation of this subroutine. It should be save
433 : ! to remove the optional argument 'do_triplet' from the subroutine interface.
434 : !my_tddft = PRESENT(local_rho_set)
435 : !IF (my_tddft) THEN
436 : ! IF (PRESENT(do_triplet)) THEN
437 : ! IF (nspins == 1 .AND. do_triplet) RETURN
438 : ! ELSE
439 : ! IF (nspins == 1 .AND. dft_control%tddfpt_control%res_etype /= tddfpt_singlet) RETURN
440 : ! END IF
441 : !END IF
442 :
443 17104 : IF (PRESENT(local_rho_set)) &
444 : CALL get_local_rho(local_rho_set, rho0_mpole=rho0_mpole, rho_atom_set=rho_atom_set, &
445 2498 : rhoz_cneo_set=rhoz_cneo_set)
446 : ! Q from rho0_mpole of local_rho_set
447 : ! for TDDFT forces we need mixed potential / integral space
448 : ! potential stored on local_rho_set_2nd
449 17104 : IF (PRESENT(local_rho_set_2nd)) THEN
450 164 : CALL get_local_rho(local_rho_set_2nd, rho_atom_set=rho_atom_set)
451 : END IF
452 : CALL get_rho0_mpole(rho0_mpole=rho0_mpole, lmax_0=lmax0, &
453 : zet0_h=zet0, igrid_zet0_s=igrid, &
454 17104 : norm_g0l_h=norm_l)
455 :
456 : ! Setup of the potential on the multigrids
457 17104 : NULLIFY (rs_descs, pw_pools)
458 17104 : CPASSERT(ASSOCIATED(pw_env))
459 17104 : CALL pw_env_get(pw_env, rs_descs=rs_descs, pw_pools=pw_pools)
460 :
461 : ! Assign from pw_env
462 17104 : auxbas_grid = pw_env%auxbas_grid
463 :
464 : ! IF present, overwrite qs grid for new pool
465 : ! Get the potential on the right grid
466 17104 : IF (PRESENT(my_pools)) THEN
467 50 : rs_desc => my_rs_descs(igrid)%rs_desc
468 50 : pw_pool => my_pools(igrid)%pool
469 : ELSE
470 17054 : rs_desc => rs_descs(igrid)%rs_desc
471 17054 : pw_pool => pw_pools(igrid)%pool
472 : END IF
473 :
474 17104 : CALL pw_pool%create_pw(coeff_gspace)
475 17104 : CALL pw_pool%create_pw(coeff_rspace)
476 :
477 17104 : IF (igrid /= auxbas_grid) THEN
478 12 : pw_aux => pw_pools(auxbas_grid)%pool
479 12 : CALL pw_aux%create_pw(coeff_gaux)
480 12 : CALL pw_transfer(v_rspace, coeff_gaux)
481 12 : CALL pw_copy(coeff_gaux, coeff_gspace)
482 12 : CALL pw_transfer(coeff_gspace, coeff_rspace)
483 12 : CALL pw_aux%give_back_pw(coeff_gaux)
484 12 : CALL pw_aux%create_pw(coeff_raux)
485 12 : fscale = coeff_rspace%pw_grid%dvol/coeff_raux%pw_grid%dvol
486 12 : CALL pw_scale(coeff_rspace, fscale)
487 12 : CALL pw_aux%give_back_pw(coeff_raux)
488 : ELSE
489 :
490 17092 : IF (coeff_gspace%pw_grid%spherical) THEN
491 0 : CALL pw_transfer(v_rspace, coeff_gspace)
492 0 : CALL pw_transfer(coeff_gspace, coeff_rspace)
493 : ELSE
494 17092 : CALL pw_copy(v_rspace, coeff_rspace)
495 : END IF
496 : END IF
497 17104 : CALL pw_pool%give_back_pw(coeff_gspace)
498 :
499 : ! Setup the rs grid at level igrid
500 17104 : CALL rs_grid_create(rs_v, rs_desc)
501 17104 : CALL rs_grid_zero(rs_v)
502 17104 : CALL transfer_pw2rs(rs_v, coeff_rspace)
503 :
504 17104 : CALL pw_pool%give_back_pw(coeff_rspace)
505 :
506 : ! Now the potential is on the right grid => integration
507 :
508 17104 : eps_rho_rspace = dft_control%qs_control%eps_rho_rspace
509 :
510 : ! Allocate work storage
511 :
512 17104 : NULLIFY (hab, hab_sph, hdab, hdab_sph, pab, a_hdab, a_hdab_sph)
513 17104 : nch_max = ncoset(lmax0)
514 17104 : CALL reallocate(hab, 1, nch_max, 1, 1)
515 17104 : CALL reallocate(hab_sph, 1, nch_max)
516 17104 : CALL reallocate(hdab, 1, 3, 1, nch_max, 1, 1)
517 17104 : CALL reallocate(hdab_sph, 1, 3, 1, nch_max)
518 17104 : CALL reallocate(a_hdab, 1, 3, 1, 3, 1, nch_max, 1, 1)
519 17104 : CALL reallocate(a_hdab_sph, 1, 3, 1, 3, 1, nch_max)
520 17104 : CALL reallocate(pab, 1, nch_max, 1, 1)
521 :
522 17104 : ncurr = -1
523 :
524 17104 : grid_distributed = rs_v%desc%distributed
525 :
526 17104 : fscale = 1.0_dp
527 17104 : IF (PRESENT(kforce)) THEN
528 40 : fscale = kforce
529 : END IF
530 :
531 50540 : DO ikind = 1, SIZE(atomic_kind_set, 1)
532 33436 : NULLIFY (basis_1c_set, atom_list, harmonics, cneo_potential)
533 33436 : CALL get_atomic_kind(atomic_kind_set(ikind), atom_list=atom_list, natom=nat)
534 : CALL get_qs_kind(qs_kind_set(ikind), &
535 : basis_set=basis_1c_set, basis_type="GAPW_1C", &
536 : paw_atom=paw_atom, &
537 : harmonics=harmonics, &
538 33436 : cneo_potential=cneo_potential)
539 :
540 33436 : IF (.NOT. paw_atom) CYCLE
541 :
542 31846 : NULLIFY (Qlm_gg, lmax, npgf)
543 : CALL get_rho0_mpole(rho0_mpole=rho0_mpole, ikind=ikind, &
544 : l0_ikind=l0_ikind, Qlm_gg=Qlm_gg, & ! Qs different
545 31846 : rpgf0_s=rpgf0)
546 :
547 : CALL get_gto_basis_set(gto_basis_set=basis_1c_set, &
548 : lmax=lmax, lmin=lmin, &
549 : maxso=maxso, maxl=maxl, &
550 31846 : nset=nset, npgf=npgf)
551 :
552 31846 : nsotot = maxso*nset
553 127384 : ALLOCATE (intloc(nsotot, nsotot))
554 :
555 : ! Initialize the local KS integrals
556 :
557 31846 : nch_ik = ncoset(l0_ikind)
558 392008 : pab = 1.0_dp
559 31846 : max_s_harm = harmonics%max_s_harm
560 31846 : llmax = harmonics%llmax
561 :
562 31846 : NULLIFY (intloc_nuc)
563 31846 : maxl_nuc = -1
564 31846 : max_s_harm_nuc = 0
565 31846 : llmax_nuc = -1
566 31846 : IF (ASSOCIATED(cneo_potential)) THEN
567 48 : NULLIFY (nuc_basis_set)
568 : CALL get_qs_kind(qs_kind_set(ikind), &
569 : basis_set=nuc_basis_set, &
570 48 : basis_type="NUC")
571 :
572 : CALL get_gto_basis_set(gto_basis_set=nuc_basis_set, &
573 : lmax=lmax_nuc, lmin=lmin_nuc, &
574 : maxso=maxso_nuc, maxl=maxl_nuc, &
575 48 : nset=nset_nuc, npgf=npgf_nuc)
576 48 : nsotot_nuc = maxso_nuc*nset_nuc
577 192 : ALLOCATE (intloc_nuc(nsotot_nuc, nsotot_nuc))
578 48 : max_s_harm_nuc = cneo_potential%harmonics%max_s_harm
579 96 : llmax_nuc = cneo_potential%harmonics%llmax
580 : END IF
581 :
582 0 : ALLOCATE (cg_list(2, nsoset(MAX(maxl, maxl_nuc))**2, &
583 : MAX(max_s_harm, max_s_harm_nuc)), &
584 191076 : cg_n_list(MAX(max_s_harm, max_s_harm_nuc)))
585 :
586 31846 : num_pe = para_env%num_pe
587 31846 : mepos = para_env%mepos
588 95538 : DO j = 0, num_pe - 1
589 63692 : bo = get_limit(nat, num_pe, j)
590 63692 : IF (.NOT. grid_distributed .AND. j /= mepos) CYCLE
591 :
592 89580 : DO iat = bo(1), bo(2)
593 25888 : iatom = atom_list(iat)
594 25888 : ra(:) = pbc(particle_set(iatom)%r, cell)
595 :
596 25888 : NULLIFY (Qlm)
597 25888 : CALL get_rho0_mpole(rho0_mpole=rho0_mpole, iat=iatom, Qlm_tot=Qlm)
598 :
599 314948 : hab = 0.0_dp
600 1104464 : hdab = 0.0_dp
601 75190324 : intloc = 0._dp
602 25888 : IF (use_virial) THEN
603 294 : my_virial_a = 0.0_dp
604 294 : my_virial_b = 0.0_dp
605 38808 : a_hdab = 0.0_dp
606 : END IF
607 :
608 : CALL integrate_pgf_product( &
609 : l0_ikind, zet0, 0, 0, 0.0_dp, 0, &
610 : ra, [0.0_dp, 0.0_dp, 0.0_dp], rs_v, &
611 : hab, pab, o1=0, o2=0, &
612 : radius=rpgf0, &
613 : calculate_forces=calculate_forces, &
614 : use_virial=use_virial, my_virial_a=my_virial_a, my_virial_b=my_virial_b, &
615 25888 : hdab=hdab, a_hdab=a_hdab, use_subpatch=.TRUE., subpatch_pattern=0)
616 :
617 : ! Convert from cartesian to spherical
618 92978 : DO lshell = 0, l0_ikind
619 285466 : DO is = 1, nso(lshell)
620 192488 : iso = is + nsoset(lshell - 1)
621 192488 : hab_sph(iso) = 0.0_dp
622 769952 : hdab_sph(1:3, iso) = 0.0_dp
623 2502344 : a_hdab_sph(1:3, 1:3, iso) = 0.0_dp
624 1126089 : DO ic = 1, nco(lshell)
625 866511 : ico = ic + ncoset(lshell - 1)
626 866511 : lx = indco(1, ico)
627 866511 : ly = indco(2, ico)
628 866511 : lz = indco(3, ico)
629 : hab_sph(iso) = hab_sph(iso) + &
630 : norm_l(lshell)* &
631 : orbtramat(lshell)%slm(is, ic)* &
632 866511 : hab(ico, 1)
633 866511 : IF (calculate_forces) THEN
634 : hdab_sph(1:3, iso) = hdab_sph(1:3, iso) + &
635 : norm_l(lshell)* &
636 : orbtramat(lshell)%slm(is, ic)* &
637 313456 : hdab(1:3, ico, 1)
638 : END IF
639 1058999 : IF (use_virial) THEN
640 47040 : DO ii = 1, 3
641 152880 : DO i = 1, 3
642 : a_hdab_sph(i, ii, iso) = a_hdab_sph(i, ii, iso) + &
643 : norm_l(lshell)* &
644 : orbtramat(lshell)%slm(is, ic)* &
645 141120 : a_hdab(i, ii, ico, 1)
646 : END DO
647 : END DO
648 : END IF
649 :
650 : END DO ! ic
651 : END DO ! is
652 : END DO ! lshell
653 :
654 : m1 = 0
655 90953 : DO iset1 = 1, nset
656 :
657 : m2 = 0
658 282618 : DO iset2 = 1, nset
659 : CALL get_none0_cg_list(harmonics%my_CG, lmin(iset1), lmax(iset1), lmin(iset2), lmax(iset2), &
660 217553 : max_s_harm, llmax, cg_list, cg_n_list, max_iso_not0_local)
661 217553 : n1 = nsoset(lmax(iset1))
662 739586 : DO ipgf1 = 1, npgf(iset1)
663 522033 : n2 = nsoset(lmax(iset2))
664 2217342 : DO ipgf2 = 1, npgf(iset2)
665 :
666 9306370 : DO iso = 1, MIN(nsoset(l0_ikind), max_iso_not0_local)
667 21206813 : DO icg = 1, cg_n_list(iso)
668 12422476 : iso1 = cg_list(1, icg, iso)
669 12422476 : iso2 = cg_list(2, icg, iso)
670 :
671 12422476 : ig1 = iso1 + n1*(ipgf1 - 1) + m1
672 12422476 : ig2 = iso2 + n2*(ipgf2 - 1) + m2
673 :
674 19729057 : intloc(ig1, ig2) = intloc(ig1, ig2) + Qlm_gg(ig1, ig2, iso)*hab_sph(iso) ! potential times Q
675 :
676 : END DO ! icg
677 : END DO ! iso
678 :
679 : END DO ! ipgf2
680 : END DO ! ipgf1
681 500171 : m2 = m2 + maxso
682 : END DO ! iset2
683 90953 : m1 = m1 + maxso
684 : END DO ! iset1
685 :
686 25888 : IF (ASSOCIATED(cneo_potential)) THEN
687 305578 : intloc_nuc = 0.0_dp
688 46 : m1 = 0
689 460 : DO iset1 = 1, nset_nuc
690 414 : n1 = nsoset(lmax_nuc(iset1))
691 414 : m2 = 0
692 4140 : DO iset2 = 1, nset_nuc
693 3726 : n2 = nsoset(lmax_nuc(iset2))
694 : CALL get_none0_cg_list(cneo_potential%harmonics%my_CG, lmin_nuc(iset1), &
695 : lmax_nuc(iset1), lmin_nuc(iset2), lmax_nuc(iset2), &
696 : max_s_harm_nuc, llmax_nuc, cg_list, cg_n_list, &
697 3726 : max_iso_not0_local)
698 7452 : DO ipgf1 = 1, npgf_nuc(iset1)
699 11178 : DO ipgf2 = 1, npgf_nuc(iset2)
700 :
701 29578 : DO iso = 1, MIN(nsoset(l0_ikind), max_iso_not0_local)
702 50968 : DO icg = 1, cg_n_list(iso)
703 25116 : iso1 = cg_list(1, icg, iso)
704 25116 : iso2 = cg_list(2, icg, iso)
705 :
706 25116 : ig1 = iso1 + n1*(ipgf1 - 1) + m1
707 25116 : ig2 = iso2 + n2*(ipgf2 - 1) + m2
708 :
709 : intloc_nuc(ig1, ig2) = intloc_nuc(ig1, ig2) - cneo_potential%zeff* &
710 47242 : cneo_potential%Qlm_gg(ig1, ig2, iso)*hab_sph(iso)
711 :
712 : END DO ! icg
713 : END DO ! iso
714 :
715 : END DO ! ipgf2
716 : END DO ! ipgf1
717 7866 : m2 = m2 + maxso_nuc
718 : END DO ! iset2
719 460 : m1 = m1 + maxso_nuc
720 : END DO ! iset1
721 : END IF
722 :
723 25888 : IF (grid_distributed) THEN
724 : ! Sum result over all processors
725 0 : CALL para_env%sum(intloc)
726 0 : IF (ASSOCIATED(cneo_potential)) THEN
727 0 : CALL para_env%sum(intloc_nuc)
728 : END IF
729 : END IF
730 :
731 25888 : IF (j == mepos) THEN
732 25888 : rho_atom => rho_atom_set(iatom)
733 25888 : CALL get_rho_atom(rho_atom=rho_atom, ga_Vlocal_gb_h=int_local_h, ga_Vlocal_gb_s=int_local_s)
734 56373 : DO ispin = 1, nspins
735 169530986 : int_local_h(ispin)%r_coef = int_local_h(ispin)%r_coef + intloc
736 169556874 : int_local_s(ispin)%r_coef = int_local_s(ispin)%r_coef + intloc
737 : END DO
738 25888 : IF (ASSOCIATED(cneo_potential)) THEN
739 46 : rhoz_cneo => rhoz_cneo_set(iatom)
740 611156 : rhoz_cneo%ga_Vlocal_gb_h = rhoz_cneo%ga_Vlocal_gb_h + intloc_nuc
741 611156 : rhoz_cneo%ga_Vlocal_gb_s = rhoz_cneo%ga_Vlocal_gb_s + intloc_nuc
742 : END IF
743 : END IF
744 :
745 25888 : IF (PRESENT(atener)) THEN
746 178 : DO iso = 1, nsoset(l0_ikind)
747 178 : atener(iatom) = atener(iatom) + 0.5_dp*Qlm(iso)*hab_sph(iso)
748 : END DO
749 : END IF
750 :
751 25888 : IF (calculate_forces) THEN
752 997 : force_tmp(1:3) = 0.0_dp
753 9410 : DO iso = 1, nsoset(l0_ikind)
754 8413 : force_tmp(1) = force_tmp(1) + Qlm(iso)*hdab_sph(1, iso) ! Q here is from local_rho_set
755 8413 : force_tmp(2) = force_tmp(2) + Qlm(iso)*hdab_sph(2, iso)
756 9410 : force_tmp(3) = force_tmp(3) + Qlm(iso)*hdab_sph(3, iso)
757 : END DO
758 3988 : force(ikind)%g0s_Vh_elec(1:3, iat) = force(ikind)%g0s_Vh_elec(1:3, iat) + fscale*force_tmp(1:3)
759 : END IF
760 89580 : IF (use_virial) THEN
761 294 : my_virial_a = 0.0_dp
762 2940 : DO iso = 1, nsoset(l0_ikind)
763 10878 : DO ii = 1, 3
764 34398 : DO i = 1, 3
765 : ! Q from local_rho_set
766 23814 : virial%pv_gapw(i, ii) = virial%pv_gapw(i, ii) + fscale*Qlm(iso)*a_hdab_sph(i, ii, iso)
767 31752 : virial%pv_virial(i, ii) = virial%pv_virial(i, ii) + fscale*Qlm(iso)*a_hdab_sph(i, ii, iso)
768 : END DO
769 : END DO
770 : END DO
771 : END IF
772 :
773 : END DO
774 : END DO
775 :
776 31846 : DEALLOCATE (intloc)
777 31846 : IF (ASSOCIATED(intloc_nuc)) DEALLOCATE (intloc_nuc)
778 115822 : DEALLOCATE (cg_list, cg_n_list)
779 :
780 : END DO ! ikind
781 :
782 17104 : CALL rs_grid_release(rs_v)
783 :
784 17104 : DEALLOCATE (hab, hdab, hab_sph, hdab_sph, pab, a_hdab, a_hdab_sph)
785 :
786 17104 : CALL timestop(handle)
787 :
788 17104 : IF (rho0_mpole%do_cneo) THEN
789 48 : IF (PRESENT(kforce)) THEN
790 : CALL integrate_vhgg_rspace(qs_env, v_rspace, para_env, calculate_forces, &
791 0 : rhoz_cneo_set, kforce)
792 : ELSE
793 : CALL integrate_vhgg_rspace(qs_env, v_rspace, para_env, calculate_forces, &
794 48 : rhoz_cneo_set)
795 : END IF
796 : END IF
797 :
798 34208 : END SUBROUTINE integrate_vhg0_rspace
799 :
800 : END MODULE qs_rho0_ggrid
|