Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \brief Routines needed for kpoint calculation
10 : !> \par History
11 : !> 2014.07 created [JGH]
12 : !> 2014.11 unified k-point and gamma-point code [Ole Schuett]
13 : !> \author JGH
14 : ! **************************************************************************************************
15 : MODULE kpoint_methods
16 : USE atomic_kind_types, ONLY: get_atomic_kind
17 : USE cell_types, ONLY: cell_type,&
18 : real_to_scaled
19 : USE cp_blacs_env, ONLY: cp_blacs_env_create,&
20 : cp_blacs_env_type
21 : USE cp_cfm_types, ONLY: cp_cfm_create,&
22 : cp_cfm_get_info,&
23 : cp_cfm_release,&
24 : cp_cfm_to_fm,&
25 : cp_cfm_type,&
26 : cp_fm_to_cfm
27 : USE cp_control_types, ONLY: hairy_probes_type
28 : USE cp_dbcsr_api, ONLY: &
29 : dbcsr_copy, dbcsr_create, dbcsr_deallocate_matrix, dbcsr_distribute, &
30 : dbcsr_distribution_get, dbcsr_distribution_type, dbcsr_get_block_p, dbcsr_get_info, &
31 : dbcsr_get_stored_coordinates, dbcsr_iterator_blocks_left, dbcsr_iterator_next_block, &
32 : dbcsr_iterator_start, dbcsr_iterator_stop, dbcsr_iterator_type, dbcsr_p_type, &
33 : dbcsr_replicate_all, dbcsr_set, dbcsr_type, dbcsr_type_antisymmetric, &
34 : dbcsr_type_no_symmetry, dbcsr_type_symmetric
35 : USE cp_dbcsr_cp2k_link, ONLY: cp_dbcsr_alloc_block_from_nbl
36 : USE cp_dbcsr_operations, ONLY: copy_fm_to_dbcsr
37 : USE cp_fm_basic_linalg, ONLY: cp_fm_column_scale
38 : USE cp_fm_pool_types, ONLY: cp_fm_pool_p_type,&
39 : fm_pool_create_fm,&
40 : fm_pool_give_back_fm
41 : USE cp_fm_struct, ONLY: cp_fm_struct_type
42 : USE cp_fm_types, ONLY: &
43 : copy_info_type, cp_fm_cleanup_copy_general, cp_fm_create, cp_fm_finish_copy_general, &
44 : cp_fm_get_diag, cp_fm_get_info, cp_fm_get_submatrix, cp_fm_release, &
45 : cp_fm_start_copy_general, cp_fm_to_fm, cp_fm_type
46 : USE cp_log_handling, ONLY: cp_logger_get_default_io_unit
47 : USE cryssym, ONLY: crys_sym_gen,&
48 : csym_type,&
49 : kpoint_gen,&
50 : kpoint_gen_general,&
51 : print_crys_symmetry,&
52 : print_kp_symmetry,&
53 : release_csym_type
54 : USE hairy_probes, ONLY: probe_occupancy_kp
55 : USE input_constants, ONLY: smear_fermi_dirac,&
56 : smear_gaussian,&
57 : smear_mp,&
58 : smear_mv
59 : USE input_cp2k_kpoints, ONLY: use_spglib_kpoint_backend,&
60 : use_spglib_kpoint_symmetry
61 : USE kinds, ONLY: dp
62 : USE kpoint_types, ONLY: get_kpoint_info,&
63 : kind_rotmat_type,&
64 : kpoint_env_create,&
65 : kpoint_env_p_type,&
66 : kpoint_env_type,&
67 : kpoint_sym_create,&
68 : kpoint_sym_type,&
69 : kpoint_type
70 : USE mathconstants, ONLY: twopi
71 : USE mathlib, ONLY: inv_3x3
72 : USE memory_utilities, ONLY: reallocate
73 : USE message_passing, ONLY: mp_cart_type,&
74 : mp_para_env_type
75 : USE parallel_gemm_api, ONLY: parallel_gemm
76 : USE particle_types, ONLY: particle_type
77 : USE qs_matrix_pools, ONLY: mpools_create,&
78 : mpools_get,&
79 : mpools_rebuild_fm_pools,&
80 : qs_matrix_pools_type
81 : USE qs_mo_types, ONLY: allocate_mo_set,&
82 : get_mo_set,&
83 : init_mo_set,&
84 : mo_set_type,&
85 : set_mo_set
86 : USE qs_neighbor_list_types, ONLY: get_iterator_info,&
87 : get_neighbor_list_set_p,&
88 : neighbor_list_iterate,&
89 : neighbor_list_iterator_create,&
90 : neighbor_list_iterator_p_type,&
91 : neighbor_list_iterator_release,&
92 : neighbor_list_set_p_type
93 : USE scf_control_types, ONLY: smear_type
94 : USE smearing_utils, ONLY: Smearkp,&
95 : Smearkp2
96 : USE util, ONLY: get_limit
97 : #include "./base/base_uses.f90"
98 :
99 : IMPLICIT NONE
100 :
101 : PRIVATE
102 :
103 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'kpoint_methods'
104 :
105 : PUBLIC :: kpoint_initialize, kpoint_env_initialize, kpoint_initialize_mos, kpoint_initialize_mo_set
106 : PUBLIC :: kpoint_init_cell_index, kpoint_set_mo_occupation
107 : PUBLIC :: kpoint_density_matrices, kpoint_density_transform
108 : PUBLIC :: rskp_transform, lowdin_kp_trans, lowdin_kp_mo_coeff
109 :
110 : ! **************************************************************************************************
111 :
112 : CONTAINS
113 :
114 : ! **************************************************************************************************
115 : !> \brief Generate the kpoints and initialize the kpoint environment
116 : !> \param kpoint The kpoint environment
117 : !> \param particle_set Particle types and coordinates
118 : !> \param cell Computational cell information
119 : ! **************************************************************************************************
120 10804 : SUBROUTINE kpoint_initialize(kpoint, particle_set, cell)
121 :
122 : TYPE(kpoint_type), POINTER :: kpoint
123 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
124 : TYPE(cell_type), POINTER :: cell
125 :
126 : CHARACTER(LEN=*), PARAMETER :: routineN = 'kpoint_initialize'
127 :
128 : INTEGER :: handle, i, ic, ik, iounit, ir, ira, is, &
129 : isign, j, natom, nkind, nr, ns
130 10804 : INTEGER, ALLOCATABLE, DIMENSION(:) :: atype
131 10804 : INTEGER, ALLOCATABLE, DIMENSION(:, :) :: agauge
132 : INTEGER, DIMENSION(3, 3) :: frot, krot
133 : LOGICAL :: spez
134 : REAL(KIND=dp) :: eps_kpoint, wsum
135 10804 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: coord, scoord
136 : REAL(KIND=dp), DIMENSION(3) :: diff, kgvec, srot
137 : REAL(KIND=dp), DIMENSION(3, 3) :: srotmat
138 10804 : REAL(KIND=dp), DIMENSION(:), POINTER :: wkp_full
139 10804 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: xkp_full
140 205276 : TYPE(csym_type) :: crys_sym
141 : TYPE(kpoint_sym_type), POINTER :: kpsym
142 :
143 10804 : CALL timeset(routineN, handle)
144 :
145 10804 : CPASSERT(ASSOCIATED(kpoint))
146 :
147 10822 : SELECT CASE (kpoint%kp_scheme)
148 : CASE ("NONE")
149 : ! do nothing
150 : CASE ("GAMMA")
151 18 : kpoint%nkp = 1
152 18 : ALLOCATE (kpoint%xkp(3, 1), kpoint%wkp(1))
153 72 : kpoint%xkp(1:3, 1) = 0.0_dp
154 18 : kpoint%wkp(1) = 1.0_dp
155 36 : ALLOCATE (kpoint%kp_sym(1))
156 18 : NULLIFY (kpoint%kp_sym(1)%kpoint_sym)
157 18 : CALL kpoint_sym_create(kpoint%kp_sym(1)%kpoint_sym)
158 : CASE ("MONKHORST-PACK", "MACDONALD")
159 :
160 2776 : IF (.NOT. kpoint%symmetry) THEN
161 : ! we set up a random molecule to avoid any possible symmetry
162 164 : natom = 10
163 164 : ALLOCATE (coord(3, natom), scoord(3, natom), atype(natom))
164 1804 : DO i = 1, natom
165 1640 : atype(i) = i
166 1640 : coord(1, i) = SIN(i*0.12345_dp)
167 1640 : coord(2, i) = COS(i*0.23456_dp)
168 1640 : coord(3, i) = SIN(i*0.34567_dp)
169 1804 : CALL real_to_scaled(scoord(1:3, i), coord(1:3, i), cell)
170 : END DO
171 : ELSE
172 2612 : natom = SIZE(particle_set)
173 13060 : ALLOCATE (scoord(3, natom), atype(natom))
174 14894 : DO i = 1, natom
175 12282 : CALL get_atomic_kind(atomic_kind=particle_set(i)%atomic_kind, kind_number=atype(i))
176 14894 : CALL real_to_scaled(scoord(1:3, i), particle_set(i)%r(1:3), cell)
177 : END DO
178 : END IF
179 2776 : IF (kpoint%verbose) THEN
180 2102 : iounit = cp_logger_get_default_io_unit()
181 : ELSE
182 674 : iounit = -1
183 : END IF
184 : ! kind type list
185 8328 : ALLOCATE (kpoint%atype(natom))
186 16698 : kpoint%atype = atype
187 : ! Match the periodic gauge used by the k-space matrices.
188 8328 : ALLOCATE (agauge(3, natom))
189 16698 : DO i = 1, natom
190 58464 : agauge(1:3, i) = -FLOOR(scoord(1:3, i) + 0.5_dp - kpoint%eps_geo)
191 : END DO
192 :
193 : CALL crys_sym_gen(crys_sym, scoord, atype, cell%hmat, delta=kpoint%eps_geo, iounit=iounit, &
194 2776 : use_spglib=kpoint%symmetry)
195 : CALL kpoint_gen(crys_sym, kpoint%nkp_grid, symm=kpoint%symmetry, shift=kpoint%kp_shift, &
196 : full_grid=kpoint%full_grid, gamma_centered=kpoint%gamma_centered, &
197 : inversion_symmetry_only=kpoint%inversion_symmetry_only, &
198 : use_spglib_reduction= &
199 : kpoint%symmetry_reduction_method == use_spglib_kpoint_symmetry, &
200 2776 : use_spglib_backend=kpoint%symmetry_backend == use_spglib_kpoint_backend)
201 2776 : kpoint%nkp = crys_sym%nkpoint
202 13880 : ALLOCATE (kpoint%xkp(3, kpoint%nkp), kpoint%wkp(kpoint%nkp))
203 13814 : wsum = SUM(crys_sym%wkpoint)
204 13814 : DO ik = 1, kpoint%nkp
205 44152 : kpoint%xkp(1:3, ik) = crys_sym%xkpoint(1:3, ik)
206 13814 : kpoint%wkp(ik) = crys_sym%wkpoint(ik)/wsum
207 : END DO
208 :
209 2776 : eps_kpoint = MAX(1.e-12_dp, 10.0_dp*kpoint%eps_geo)
210 : ! print output
211 2776 : IF (kpoint%symmetry) CALL print_crys_symmetry(crys_sym)
212 2776 : IF (kpoint%symmetry) CALL print_kp_symmetry(crys_sym)
213 :
214 : ! transfer symmetry information
215 19366 : ALLOCATE (kpoint%kp_sym(kpoint%nkp))
216 13814 : DO ik = 1, kpoint%nkp
217 11038 : NULLIFY (kpoint%kp_sym(ik)%kpoint_sym)
218 11038 : CALL kpoint_sym_create(kpoint%kp_sym(ik)%kpoint_sym)
219 11038 : kpsym => kpoint%kp_sym(ik)%kpoint_sym
220 : IF (crys_sym%nrtot > 0 .AND. .NOT. crys_sym%fullgrid .AND. &
221 13814 : crys_sym%istriz == 1 .AND. .NOT. crys_sym%inversion_only) THEN
222 : ! set up the symmetrization information
223 1074 : kpsym%nwght = NINT(crys_sym%wkpoint(ik))
224 1074 : ns = kpsym%nwght
225 : !
226 1074 : IF (ns > 1) THEN
227 42338 : DO is = 1, SIZE(crys_sym%kplink, 2)
228 42338 : IF (crys_sym%kplink(2, is) == ik) THEN
229 1232836 : DO ic = 1, crys_sym%nrtot
230 96715592 : srotmat = MATMUL(cell%h_inv, MATMUL(crys_sym%rt(1:3, 1:3, ic), cell%hmat))
231 15915224 : frot(1:3, 1:3) = NINT(srotmat(1:3, 1:3))
232 31830448 : krot(1:3, 1:3) = NINT(TRANSPOSE(inv_3x3(REAL(frot(1:3, 1:3), KIND=dp))))
233 3681332 : DO isign = 1, 2
234 2448496 : ir = MERGE(crys_sym%ibrot(ic), -crys_sym%ibrot(ic), isign == 1)
235 2448496 : IF (ir == crys_sym%kpop(is)) CYCLE
236 : kgvec(1:3) = crys_sym%kpmesh(1:3, is) - &
237 : MATMUL(REAL(MERGE(krot(1:3, 1:3), -krot(1:3, 1:3), &
238 : isign == 1), KIND=dp), &
239 68317424 : kpoint%xkp(1:3, ik))
240 9759632 : diff(1:3) = kgvec(1:3) - ANINT(kgvec(1:3))
241 4984792 : IF (ALL(ABS(diff(1:3)) < eps_kpoint)) ns = ns + 1
242 : END DO
243 : END DO
244 : END IF
245 : END DO
246 1064 : kpsym%apply_symmetry = .TRUE.
247 1064 : natom = SIZE(particle_set)
248 3192 : ALLOCATE (kpsym%rot(3, 3, ns))
249 3192 : ALLOCATE (kpsym%xkp(3, ns))
250 3192 : ALLOCATE (kpsym%rotp(ns))
251 4256 : ALLOCATE (kpsym%f0(natom, ns))
252 4256 : ALLOCATE (kpsym%fcell(3, natom, ns))
253 4256 : ALLOCATE (kpsym%kgphase(natom, ns))
254 1064 : nr = 0
255 42338 : DO is = 1, SIZE(crys_sym%kplink, 2)
256 42338 : IF (crys_sym%kplink(2, is) == ik) THEN
257 8588 : nr = nr + 1
258 8588 : ir = crys_sym%kpop(is)
259 8588 : ira = ABS(ir)
260 61200 : DO ic = 1, crys_sym%nrtot
261 61200 : IF (crys_sym%ibrot(ic) == ira) THEN
262 8588 : kpsym%rotp(nr) = ir
263 111644 : kpsym%rot(1:3, 1:3, nr) = crys_sym%rt(1:3, 1:3, ic)
264 678452 : srotmat = MATMUL(cell%h_inv, MATMUL(crys_sym%rt(1:3, 1:3, ic), cell%hmat))
265 111644 : frot(1:3, 1:3) = NINT(srotmat(1:3, 1:3))
266 34352 : kpsym%xkp(1:3, nr) = crys_sym%kpmesh(1:3, is)
267 223288 : krot(1:3, 1:3) = NINT(TRANSPOSE(inv_3x3(REAL(frot(1:3, 1:3), KIND=dp))))
268 60116 : IF (ir < 0) krot(1:3, 1:3) = -krot(1:3, 1:3)
269 : kgvec(1:3) = kpsym%xkp(1:3, nr) - &
270 : MATMUL(REAL(krot(1:3, 1:3), KIND=dp), &
271 240464 : kpoint%xkp(1:3, ik))
272 34352 : kgvec(1:3) = ANINT(kgvec(1:3))
273 72412 : kpsym%f0(1:natom, nr) = crys_sym%f0(1:natom, ic)
274 72412 : DO j = 1, natom
275 1021184 : srot(1:3) = MATMUL(srotmat, scoord(1:3, j)) + crys_sym%vt(1:3, ic)
276 : kpsym%fcell(1:3, j, nr) = &
277 : NINT(srot(1:3) - scoord(1:3, kpsym%f0(j, nr))) + &
278 : MATMUL(frot(1:3, 1:3), agauge(1:3, j)) - &
279 1212656 : agauge(1:3, kpsym%f0(j, nr))
280 : kpsym%kgphase(j, nr) = DOT_PRODUCT(kgvec(1:3), &
281 : scoord(1:3, j) + &
282 263884 : REAL(agauge(1:3, j), KIND=dp))
283 : END DO
284 : EXIT
285 : END IF
286 : END DO
287 8588 : CPASSERT(ic <= crys_sym%nrtot)
288 : END IF
289 : END DO
290 42338 : DO is = 1, SIZE(crys_sym%kplink, 2)
291 42338 : IF (crys_sym%kplink(2, is) == ik) THEN
292 1232836 : DO ic = 1, crys_sym%nrtot
293 96715592 : srotmat = MATMUL(cell%h_inv, MATMUL(crys_sym%rt(1:3, 1:3, ic), cell%hmat))
294 15915224 : frot(1:3, 1:3) = NINT(srotmat(1:3, 1:3))
295 31830448 : krot(1:3, 1:3) = NINT(TRANSPOSE(inv_3x3(REAL(frot(1:3, 1:3), KIND=dp))))
296 3681332 : DO isign = 1, 2
297 2448496 : ir = MERGE(crys_sym%ibrot(ic), -crys_sym%ibrot(ic), isign == 1)
298 2448496 : IF (ir == crys_sym%kpop(is)) CYCLE
299 : kgvec(1:3) = crys_sym%kpmesh(1:3, is) - &
300 : MATMUL(REAL(MERGE(krot(1:3, 1:3), -krot(1:3, 1:3), &
301 : isign == 1), KIND=dp), &
302 68317424 : kpoint%xkp(1:3, ik))
303 9759632 : diff(1:3) = kgvec(1:3) - ANINT(kgvec(1:3))
304 4984792 : IF (ALL(ABS(diff(1:3)) < eps_kpoint)) THEN
305 159068 : nr = nr + 1
306 159068 : kpsym%rotp(nr) = ir
307 2067884 : kpsym%rot(1:3, 1:3, nr) = crys_sym%rt(1:3, 1:3, ic)
308 636272 : kpsym%xkp(1:3, nr) = crys_sym%kpmesh(1:3, is)
309 636272 : kgvec(1:3) = ANINT(kgvec(1:3))
310 1412268 : kpsym%f0(1:natom, nr) = crys_sym%f0(1:natom, ic)
311 1412268 : DO j = 1, natom
312 20051200 : srot(1:3) = MATMUL(srotmat, scoord(1:3, j)) + crys_sym%vt(1:3, ic)
313 : kpsym%fcell(1:3, j, nr) = &
314 : NINT(srot(1:3) - scoord(1:3, kpsym%f0(j, nr))) + &
315 : MATMUL(frot(1:3, 1:3), agauge(1:3, j)) - &
316 23810800 : agauge(1:3, kpsym%f0(j, nr))
317 : kpsym%kgphase(j, nr) = DOT_PRODUCT(kgvec(1:3), &
318 : scoord(1:3, j) + &
319 5171868 : REAL(agauge(1:3, j), KIND=dp))
320 : END DO
321 : END IF
322 : END DO
323 : END DO
324 : END IF
325 : END DO
326 1064 : kpsym%nwred = nr
327 : END IF
328 : END IF
329 : END DO
330 2776 : IF (kpoint%symmetry) THEN
331 14894 : nkind = MAXVAL(atype)
332 2612 : ns = crys_sym%nrtot
333 46088 : ALLOCATE (kpoint%kind_rotmat(ns, nkind))
334 36774 : DO i = 1, ns
335 71088 : DO j = 1, nkind
336 68476 : NULLIFY (kpoint%kind_rotmat(i, j)%rmat)
337 : END DO
338 : END DO
339 5654 : ALLOCATE (kpoint%ibrot(ns))
340 36774 : kpoint%ibrot(1:ns) = crys_sym%ibrot(1:ns)
341 : END IF
342 :
343 2776 : CALL release_csym_type(crys_sym)
344 2776 : DEALLOCATE (scoord, atype)
345 2776 : DEALLOCATE (agauge)
346 :
347 : CASE ("GENERAL")
348 : NULLIFY (xkp_full, wkp_full)
349 36 : IF (ASSOCIATED(kpoint%xkp_input)) THEN
350 36 : xkp_full => kpoint%xkp_input
351 36 : wkp_full => kpoint%wkp_input
352 : ELSE
353 0 : xkp_full => kpoint%xkp
354 0 : wkp_full => kpoint%wkp
355 : END IF
356 36 : CPASSERT(ASSOCIATED(xkp_full))
357 36 : CPASSERT(ASSOCIATED(wkp_full))
358 36 : IF (.NOT. ASSOCIATED(kpoint%xkp_input)) THEN
359 0 : ALLOCATE (kpoint%xkp_input(3, SIZE(wkp_full)), kpoint%wkp_input(SIZE(wkp_full)))
360 0 : kpoint%xkp_input(1:3, 1:SIZE(wkp_full)) = xkp_full(1:3, 1:SIZE(wkp_full))
361 0 : kpoint%wkp_input(1:SIZE(wkp_full)) = wkp_full(1:SIZE(wkp_full))
362 0 : xkp_full => kpoint%xkp_input
363 0 : wkp_full => kpoint%wkp_input
364 : END IF
365 36 : IF (.NOT. kpoint%symmetry) THEN
366 10 : IF (.NOT. ASSOCIATED(kpoint%xkp)) THEN
367 0 : kpoint%nkp = SIZE(wkp_full)
368 0 : ALLOCATE (kpoint%xkp(3, kpoint%nkp), kpoint%wkp(kpoint%nkp))
369 0 : kpoint%xkp(1:3, 1:kpoint%nkp) = xkp_full(1:3, 1:kpoint%nkp)
370 0 : kpoint%wkp(1:kpoint%nkp) = wkp_full(1:kpoint%nkp)
371 : END IF
372 : ! default: no symmetry settings
373 74 : ALLOCATE (kpoint%kp_sym(kpoint%nkp))
374 54 : DO i = 1, kpoint%nkp
375 44 : NULLIFY (kpoint%kp_sym(i)%kpoint_sym)
376 54 : CALL kpoint_sym_create(kpoint%kp_sym(i)%kpoint_sym)
377 : END DO
378 : ELSE
379 26 : IF (kpoint%verbose) THEN
380 16 : iounit = cp_logger_get_default_io_unit()
381 : ELSE
382 10 : iounit = -1
383 : END IF
384 26 : natom = SIZE(particle_set)
385 130 : ALLOCATE (scoord(3, natom), atype(natom))
386 234 : DO i = 1, natom
387 208 : CALL get_atomic_kind(atomic_kind=particle_set(i)%atomic_kind, kind_number=atype(i))
388 234 : CALL real_to_scaled(scoord(1:3, i), particle_set(i)%r(1:3), cell)
389 : END DO
390 52 : ALLOCATE (kpoint%atype(natom))
391 234 : kpoint%atype = atype
392 78 : ALLOCATE (agauge(3, natom))
393 234 : DO i = 1, natom
394 858 : agauge(1:3, i) = -FLOOR(scoord(1:3, i) + 0.5_dp - kpoint%eps_geo)
395 : END DO
396 :
397 : CALL crys_sym_gen(crys_sym, scoord, atype, cell%hmat, delta=kpoint%eps_geo, iounit=iounit, &
398 : use_spglib=(kpoint%symmetry_backend == use_spglib_kpoint_backend .OR. &
399 30 : kpoint%symmetry_reduction_method == use_spglib_kpoint_symmetry))
400 : CALL kpoint_gen_general(crys_sym, xkp_full, wkp_full, symm=kpoint%symmetry, &
401 : full_grid=kpoint%full_grid, &
402 : inversion_symmetry_only=kpoint%inversion_symmetry_only, &
403 : use_spglib_reduction= &
404 : kpoint%symmetry_reduction_method == use_spglib_kpoint_symmetry, &
405 26 : use_spglib_backend=kpoint%symmetry_backend == use_spglib_kpoint_backend)
406 26 : IF (ASSOCIATED(kpoint%xkp)) THEN
407 26 : DEALLOCATE (kpoint%xkp)
408 26 : NULLIFY (kpoint%xkp)
409 : END IF
410 26 : IF (ASSOCIATED(kpoint%wkp)) THEN
411 26 : DEALLOCATE (kpoint%wkp)
412 26 : NULLIFY (kpoint%wkp)
413 : END IF
414 26 : kpoint%nkp = crys_sym%nkpoint
415 130 : ALLOCATE (kpoint%xkp(3, kpoint%nkp), kpoint%wkp(kpoint%nkp))
416 52 : wsum = SUM(crys_sym%wkpoint)
417 52 : DO ik = 1, kpoint%nkp
418 104 : kpoint%xkp(1:3, ik) = crys_sym%xkpoint(1:3, ik)
419 52 : kpoint%wkp(ik) = crys_sym%wkpoint(ik)/wsum
420 : END DO
421 :
422 26 : eps_kpoint = MAX(1.e-12_dp, 10.0_dp*kpoint%eps_geo)
423 26 : CALL print_crys_symmetry(crys_sym)
424 26 : CALL print_kp_symmetry(crys_sym)
425 :
426 104 : ALLOCATE (kpoint%kp_sym(kpoint%nkp))
427 52 : DO ik = 1, kpoint%nkp
428 26 : NULLIFY (kpoint%kp_sym(ik)%kpoint_sym)
429 26 : CALL kpoint_sym_create(kpoint%kp_sym(ik)%kpoint_sym)
430 26 : kpsym => kpoint%kp_sym(ik)%kpoint_sym
431 : IF (crys_sym%nrtot > 0 .AND. .NOT. crys_sym%fullgrid .AND. &
432 52 : crys_sym%istriz == 1 .AND. .NOT. crys_sym%inversion_only) THEN
433 26 : kpsym%nwght = NINT(crys_sym%wkpoint(ik))
434 26 : ns = kpsym%nwght
435 26 : IF (ns > 1) THEN
436 234 : DO is = 1, SIZE(crys_sym%kplink, 2)
437 234 : IF (crys_sym%kplink(2, is) == ik) THEN
438 30928 : DO ic = 1, crys_sym%nrtot
439 2426880 : srotmat = MATMUL(cell%h_inv, MATMUL(crys_sym%rt(1:3, 1:3, ic), cell%hmat))
440 399360 : frot(1:3, 1:3) = NINT(srotmat(1:3, 1:3))
441 798720 : krot(1:3, 1:3) = NINT(TRANSPOSE(inv_3x3(REAL(frot(1:3, 1:3), KIND=dp))))
442 92368 : DO isign = 1, 2
443 61440 : ir = MERGE(crys_sym%ibrot(ic), -crys_sym%ibrot(ic), isign == 1)
444 61440 : IF (ir == crys_sym%kpop(is)) CYCLE
445 : kgvec(1:3) = crys_sym%kpmesh(1:3, is) - &
446 : MATMUL(REAL(MERGE(krot(1:3, 1:3), -krot(1:3, 1:3), &
447 : isign == 1), KIND=dp), &
448 1714496 : kpoint%xkp(1:3, ik))
449 244928 : diff(1:3) = kgvec(1:3) - ANINT(kgvec(1:3))
450 145088 : IF (ALL(ABS(diff(1:3)) < eps_kpoint)) ns = ns + 1
451 : END DO
452 : END DO
453 : END IF
454 : END DO
455 26 : kpsym%apply_symmetry = .TRUE.
456 78 : ALLOCATE (kpsym%rot(3, 3, ns))
457 78 : ALLOCATE (kpsym%xkp(3, ns))
458 78 : ALLOCATE (kpsym%rotp(ns))
459 104 : ALLOCATE (kpsym%f0(natom, ns))
460 104 : ALLOCATE (kpsym%fcell(3, natom, ns))
461 104 : ALLOCATE (kpsym%kgphase(natom, ns))
462 26 : nr = 0
463 234 : DO is = 1, SIZE(crys_sym%kplink, 2)
464 234 : IF (crys_sym%kplink(2, is) == ik) THEN
465 208 : nr = nr + 1
466 208 : ir = crys_sym%kpop(is)
467 208 : ira = ABS(ir)
468 628 : DO ic = 1, crys_sym%nrtot
469 628 : IF (crys_sym%ibrot(ic) == ira) THEN
470 208 : kpsym%rotp(nr) = ir
471 2704 : kpsym%rot(1:3, 1:3, nr) = crys_sym%rt(1:3, 1:3, ic)
472 16432 : srotmat = MATMUL(cell%h_inv, MATMUL(crys_sym%rt(1:3, 1:3, ic), cell%hmat))
473 2704 : frot(1:3, 1:3) = NINT(srotmat(1:3, 1:3))
474 832 : kpsym%xkp(1:3, nr) = crys_sym%kpmesh(1:3, is)
475 5408 : krot(1:3, 1:3) = NINT(TRANSPOSE(inv_3x3(REAL(frot(1:3, 1:3), KIND=dp))))
476 1456 : IF (ir < 0) krot(1:3, 1:3) = -krot(1:3, 1:3)
477 : kgvec(1:3) = kpsym%xkp(1:3, nr) - &
478 : MATMUL(REAL(krot(1:3, 1:3), KIND=dp), &
479 5824 : kpoint%xkp(1:3, ik))
480 832 : kgvec(1:3) = ANINT(kgvec(1:3))
481 1872 : kpsym%f0(1:natom, nr) = crys_sym%f0(1:natom, ic)
482 1872 : DO j = 1, natom
483 26624 : srot(1:3) = MATMUL(srotmat, scoord(1:3, j)) + crys_sym%vt(1:3, ic)
484 : kpsym%fcell(1:3, j, nr) = &
485 : NINT(srot(1:3) - scoord(1:3, kpsym%f0(j, nr))) + &
486 : MATMUL(frot(1:3, 1:3), agauge(1:3, j)) - &
487 31616 : agauge(1:3, kpsym%f0(j, nr))
488 : kpsym%kgphase(j, nr) = DOT_PRODUCT(kgvec(1:3), &
489 : scoord(1:3, j) + &
490 6864 : REAL(agauge(1:3, j), KIND=dp))
491 : END DO
492 : EXIT
493 : END IF
494 : END DO
495 208 : CPASSERT(ic <= crys_sym%nrtot)
496 : END IF
497 : END DO
498 234 : DO is = 1, SIZE(crys_sym%kplink, 2)
499 234 : IF (crys_sym%kplink(2, is) == ik) THEN
500 30928 : DO ic = 1, crys_sym%nrtot
501 2426880 : srotmat = MATMUL(cell%h_inv, MATMUL(crys_sym%rt(1:3, 1:3, ic), cell%hmat))
502 399360 : frot(1:3, 1:3) = NINT(srotmat(1:3, 1:3))
503 798720 : krot(1:3, 1:3) = NINT(TRANSPOSE(inv_3x3(REAL(frot(1:3, 1:3), KIND=dp))))
504 92368 : DO isign = 1, 2
505 61440 : ir = MERGE(crys_sym%ibrot(ic), -crys_sym%ibrot(ic), isign == 1)
506 61440 : IF (ir == crys_sym%kpop(is)) CYCLE
507 : kgvec(1:3) = crys_sym%kpmesh(1:3, is) - &
508 : MATMUL(REAL(MERGE(krot(1:3, 1:3), -krot(1:3, 1:3), &
509 : isign == 1), KIND=dp), &
510 1714496 : kpoint%xkp(1:3, ik))
511 244928 : diff(1:3) = kgvec(1:3) - ANINT(kgvec(1:3))
512 145088 : IF (ALL(ABS(diff(1:3)) < eps_kpoint)) THEN
513 7472 : nr = nr + 1
514 7472 : kpsym%rotp(nr) = ir
515 97136 : kpsym%rot(1:3, 1:3, nr) = crys_sym%rt(1:3, 1:3, ic)
516 29888 : kpsym%xkp(1:3, nr) = crys_sym%kpmesh(1:3, is)
517 29888 : kgvec(1:3) = ANINT(kgvec(1:3))
518 67248 : kpsym%f0(1:natom, nr) = crys_sym%f0(1:natom, ic)
519 67248 : DO j = 1, natom
520 956416 : srot(1:3) = MATMUL(srotmat, scoord(1:3, j)) + crys_sym%vt(1:3, ic)
521 : kpsym%fcell(1:3, j, nr) = &
522 : NINT(srot(1:3) - scoord(1:3, kpsym%f0(j, nr))) + &
523 : MATMUL(frot(1:3, 1:3), agauge(1:3, j)) - &
524 1135744 : agauge(1:3, kpsym%f0(j, nr))
525 : kpsym%kgphase(j, nr) = DOT_PRODUCT(kgvec(1:3), &
526 : scoord(1:3, j) + &
527 246576 : REAL(agauge(1:3, j), KIND=dp))
528 : END DO
529 : END IF
530 : END DO
531 : END DO
532 : END IF
533 : END DO
534 26 : kpsym%nwred = nr
535 : END IF
536 : END IF
537 : END DO
538 234 : nkind = MAXVAL(atype)
539 26 : ns = crys_sym%nrtot
540 3970 : ALLOCATE (kpoint%kind_rotmat(ns, nkind))
541 3866 : DO i = 1, ns
542 7706 : DO j = 1, nkind
543 7680 : NULLIFY (kpoint%kind_rotmat(i, j)%rmat)
544 : END DO
545 : END DO
546 78 : ALLOCATE (kpoint%ibrot(ns))
547 3866 : kpoint%ibrot(1:ns) = crys_sym%ibrot(1:ns)
548 :
549 26 : CALL release_csym_type(crys_sym)
550 26 : DEALLOCATE (scoord, atype)
551 26 : DEALLOCATE (agauge)
552 : END IF
553 : CASE DEFAULT
554 10804 : CPASSERT(.FALSE.)
555 : END SELECT
556 :
557 : ! check for consistency of options
558 10822 : SELECT CASE (kpoint%kp_scheme)
559 : CASE ("NONE")
560 : ! don't use k-point code
561 : CASE ("GAMMA")
562 18 : CPASSERT(kpoint%nkp == 1)
563 90 : CPASSERT(SUM(ABS(kpoint%xkp)) <= 1.e-12_dp)
564 18 : CPASSERT(kpoint%wkp(1) == 1.0_dp)
565 18 : CPASSERT(.NOT. kpoint%symmetry)
566 : CASE ("GENERAL")
567 36 : CPASSERT(kpoint%nkp >= 1)
568 : CASE ("MONKHORST-PACK", "MACDONALD")
569 10804 : CPASSERT(kpoint%nkp >= 1)
570 : END SELECT
571 10804 : IF (kpoint%use_real_wfn) THEN
572 : ! what about inversion symmetry?
573 40 : ikloop: DO ik = 1, kpoint%nkp
574 100 : DO i = 1, 3
575 60 : spez = (kpoint%xkp(i, ik) == 0.0_dp .OR. kpoint%xkp(i, ik) == 0.5_dp)
576 20 : IF (.NOT. spez) EXIT ikloop
577 : END DO
578 : END DO ikloop
579 20 : IF (.NOT. spez) THEN
580 : ! Warning: real wfn might be wrong for this system
581 : CALL cp_warn(__LOCATION__, &
582 : "A calculation using real wavefunctions is requested. "// &
583 0 : "We could not determine if the symmetry of the system allows real wavefunctions. ")
584 : END IF
585 : END IF
586 :
587 10804 : CALL timestop(handle)
588 :
589 21608 : END SUBROUTINE kpoint_initialize
590 :
591 : ! **************************************************************************************************
592 : !> \brief Initialize the kpoint environment
593 : !> \param kpoint Kpoint environment
594 : !> \param para_env ...
595 : !> \param blacs_env ...
596 : !> \param with_aux_fit ...
597 : ! **************************************************************************************************
598 2674 : SUBROUTINE kpoint_env_initialize(kpoint, para_env, blacs_env, with_aux_fit)
599 :
600 : TYPE(kpoint_type), INTENT(INOUT) :: kpoint
601 : TYPE(mp_para_env_type), INTENT(IN), TARGET :: para_env
602 : TYPE(cp_blacs_env_type), INTENT(IN), TARGET :: blacs_env
603 : LOGICAL, INTENT(IN), OPTIONAL :: with_aux_fit
604 :
605 : CHARACTER(LEN=*), PARAMETER :: routineN = 'kpoint_env_initialize'
606 :
607 : INTEGER :: handle, igr, ik, ikk, ngr, niogrp, nkp, &
608 : nkp_grp, nkp_loc, npe, unit_nr
609 : INTEGER, DIMENSION(2) :: dims, pos
610 : LOGICAL :: aux_fit
611 2674 : TYPE(kpoint_env_p_type), DIMENSION(:), POINTER :: kp_aux_env, kp_env
612 : TYPE(kpoint_env_type), POINTER :: kp
613 2674 : TYPE(mp_cart_type) :: comm_cart
614 : TYPE(mp_para_env_type), POINTER :: para_env_inter_kp, para_env_kp
615 :
616 2674 : CALL timeset(routineN, handle)
617 :
618 2674 : IF (PRESENT(with_aux_fit)) THEN
619 2568 : aux_fit = with_aux_fit
620 : ELSE
621 : aux_fit = .FALSE.
622 : END IF
623 :
624 2674 : kpoint%para_env => para_env
625 2674 : CALL kpoint%para_env%retain()
626 2674 : kpoint%blacs_env_all => blacs_env
627 2674 : CALL kpoint%blacs_env_all%retain()
628 :
629 2674 : CPASSERT(.NOT. ASSOCIATED(kpoint%kp_env))
630 2674 : IF (aux_fit) THEN
631 32 : CPASSERT(.NOT. ASSOCIATED(kpoint%kp_aux_env))
632 : END IF
633 :
634 2674 : NULLIFY (kp_env, kp_aux_env)
635 2674 : nkp = kpoint%nkp
636 2674 : npe = para_env%num_pe
637 2674 : IF (npe == 1) THEN
638 : ! only one process available -> owns all kpoints
639 0 : ALLOCATE (kp_env(nkp))
640 0 : DO ik = 1, nkp
641 0 : NULLIFY (kp_env(ik)%kpoint_env)
642 0 : CALL kpoint_env_create(kp_env(ik)%kpoint_env)
643 0 : kp => kp_env(ik)%kpoint_env
644 0 : kp%nkpoint = ik
645 0 : kp%wkp = kpoint%wkp(ik)
646 0 : kp%xkp(1:3) = kpoint%xkp(1:3, ik)
647 0 : kp%is_local = .TRUE.
648 : END DO
649 0 : kpoint%kp_env => kp_env
650 :
651 0 : IF (aux_fit) THEN
652 0 : ALLOCATE (kp_aux_env(nkp))
653 0 : DO ik = 1, nkp
654 0 : NULLIFY (kp_aux_env(ik)%kpoint_env)
655 0 : CALL kpoint_env_create(kp_aux_env(ik)%kpoint_env)
656 0 : kp => kp_aux_env(ik)%kpoint_env
657 0 : kp%nkpoint = ik
658 0 : kp%wkp = kpoint%wkp(ik)
659 0 : kp%xkp(1:3) = kpoint%xkp(1:3, ik)
660 0 : kp%is_local = .TRUE.
661 : END DO
662 :
663 0 : kpoint%kp_aux_env => kp_aux_env
664 : END IF
665 :
666 0 : ALLOCATE (kpoint%kp_dist(2, 1))
667 0 : kpoint%kp_dist(1, 1) = 1
668 0 : kpoint%kp_dist(2, 1) = nkp
669 0 : kpoint%kp_range(1) = 1
670 0 : kpoint%kp_range(2) = nkp
671 :
672 : ! parallel environments
673 0 : kpoint%para_env_kp => para_env
674 0 : CALL kpoint%para_env_kp%retain()
675 0 : kpoint%para_env_inter_kp => para_env
676 0 : CALL kpoint%para_env_inter_kp%retain()
677 0 : kpoint%iogrp = .TRUE.
678 0 : kpoint%nkp_groups = 1
679 : ELSE
680 2674 : IF (kpoint%parallel_group_size == -1) THEN
681 : ! maximum parallelization over kpoints
682 : ! making sure that the group size divides the npe and the nkp_grp the nkp
683 : ! in the worst case, there will be no parallelism over kpoints.
684 7200 : DO igr = npe, 1, -1
685 4800 : IF (MOD(npe, igr) /= 0) CYCLE
686 4800 : nkp_grp = npe/igr
687 4800 : IF (MOD(nkp, nkp_grp) /= 0) CYCLE
688 7200 : ngr = igr
689 : END DO
690 274 : ELSE IF (kpoint%parallel_group_size == 0) THEN
691 : ! no parallelization over kpoints
692 186 : ngr = npe
693 88 : ELSE IF (kpoint%parallel_group_size > 0) THEN
694 88 : ngr = MIN(kpoint%parallel_group_size, npe)
695 : ELSE
696 0 : CPASSERT(.FALSE.)
697 : END IF
698 2674 : nkp_grp = npe/ngr
699 : ! processor dimensions
700 2674 : dims(1) = ngr
701 2674 : dims(2) = nkp_grp
702 2674 : CPASSERT(MOD(nkp, nkp_grp) == 0)
703 2674 : nkp_loc = nkp/nkp_grp
704 :
705 2674 : IF ((dims(1)*dims(2) /= npe)) THEN
706 0 : CPABORT("Number of processors is not divisible by the kpoint group size.")
707 : END IF
708 :
709 : ! Create the subgroups, one for each k-point group and one interconnecting group
710 2674 : CALL comm_cart%create(comm_old=para_env, ndims=2, dims=dims)
711 8022 : pos = comm_cart%mepos_cart
712 2674 : ALLOCATE (para_env_kp)
713 2674 : CALL para_env_kp%from_split(comm_cart, pos(2))
714 2674 : ALLOCATE (para_env_inter_kp)
715 2674 : CALL para_env_inter_kp%from_split(comm_cart, pos(1))
716 2674 : CALL comm_cart%free()
717 :
718 2674 : niogrp = 0
719 2674 : IF (para_env%is_source()) niogrp = 1
720 2674 : CALL para_env_kp%sum(niogrp)
721 2674 : kpoint%iogrp = (niogrp == 1)
722 :
723 : ! parallel groups
724 2674 : kpoint%para_env_kp => para_env_kp
725 2674 : kpoint%para_env_inter_kp => para_env_inter_kp
726 :
727 : ! distribution of kpoints
728 8022 : ALLOCATE (kpoint%kp_dist(2, nkp_grp))
729 7084 : DO igr = 1, nkp_grp
730 15904 : kpoint%kp_dist(1:2, igr) = get_limit(nkp, nkp_grp, igr - 1)
731 : END DO
732 : ! local kpoints
733 8022 : kpoint%kp_range(1:2) = kpoint%kp_dist(1:2, para_env_inter_kp%mepos + 1)
734 :
735 14586 : ALLOCATE (kp_env(nkp_loc))
736 9238 : DO ik = 1, nkp_loc
737 6564 : NULLIFY (kp_env(ik)%kpoint_env)
738 6564 : ikk = kpoint%kp_range(1) + ik - 1
739 6564 : CALL kpoint_env_create(kp_env(ik)%kpoint_env)
740 6564 : kp => kp_env(ik)%kpoint_env
741 6564 : kp%nkpoint = ikk
742 6564 : kp%wkp = kpoint%wkp(ikk)
743 26256 : kp%xkp(1:3) = kpoint%xkp(1:3, ikk)
744 9238 : kp%is_local = (ngr == 1)
745 : END DO
746 2674 : kpoint%kp_env => kp_env
747 :
748 2674 : IF (aux_fit) THEN
749 282 : ALLOCATE (kp_aux_env(nkp_loc))
750 250 : DO ik = 1, nkp_loc
751 218 : NULLIFY (kp_aux_env(ik)%kpoint_env)
752 218 : ikk = kpoint%kp_range(1) + ik - 1
753 218 : CALL kpoint_env_create(kp_aux_env(ik)%kpoint_env)
754 218 : kp => kp_aux_env(ik)%kpoint_env
755 218 : kp%nkpoint = ikk
756 218 : kp%wkp = kpoint%wkp(ikk)
757 872 : kp%xkp(1:3) = kpoint%xkp(1:3, ikk)
758 250 : kp%is_local = (ngr == 1)
759 : END DO
760 32 : kpoint%kp_aux_env => kp_aux_env
761 : END IF
762 :
763 2674 : unit_nr = cp_logger_get_default_io_unit()
764 :
765 2674 : IF (unit_nr > 0 .AND. kpoint%verbose) THEN
766 1060 : WRITE (unit_nr, *)
767 1060 : WRITE (unit_nr, FMT="(T2,A,T71,I10)") "KPOINTS| Number of kpoint groups ", nkp_grp
768 1060 : WRITE (unit_nr, FMT="(T2,A,T71,I10)") "KPOINTS| Size of each kpoint group", ngr
769 1060 : WRITE (unit_nr, FMT="(T2,A,T71,I10)") "KPOINTS| Number of kpoints per group", nkp_loc
770 : END IF
771 2674 : kpoint%nkp_groups = nkp_grp
772 :
773 : END IF
774 :
775 2674 : CALL timestop(handle)
776 :
777 5348 : END SUBROUTINE kpoint_env_initialize
778 :
779 : ! **************************************************************************************************
780 : !> \brief Initialize a set of MOs and density matrix for each kpoint (kpoint group)
781 : !> \param kpoint Kpoint environment
782 : !> \param mos Reference MOs (global)
783 : !> \param added_mos ...
784 : !> \param for_aux_fit ...
785 : ! **************************************************************************************************
786 2706 : SUBROUTINE kpoint_initialize_mos(kpoint, mos, added_mos, for_aux_fit)
787 :
788 : TYPE(kpoint_type), POINTER :: kpoint
789 : TYPE(mo_set_type), DIMENSION(:), INTENT(INOUT) :: mos
790 : INTEGER, INTENT(IN), OPTIONAL :: added_mos
791 : LOGICAL, OPTIONAL :: for_aux_fit
792 :
793 : CHARACTER(LEN=*), PARAMETER :: routineN = 'kpoint_initialize_mos'
794 :
795 : INTEGER :: handle, ic, ik, is, nadd, nao, nc, &
796 : nelectron, nkp_loc, nmo, nmorig(2), &
797 : nspin
798 : LOGICAL :: aux_fit
799 : REAL(KIND=dp) :: flexible_electron_count, maxocc, n_el_f
800 : TYPE(cp_blacs_env_type), POINTER :: blacs_env
801 2706 : TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_ao_fm_pools
802 : TYPE(cp_fm_struct_type), POINTER :: matrix_struct
803 : TYPE(cp_fm_type), POINTER :: fmlocal
804 : TYPE(kpoint_env_type), POINTER :: kp
805 : TYPE(qs_matrix_pools_type), POINTER :: mpools
806 :
807 2706 : CALL timeset(routineN, handle)
808 :
809 2706 : IF (PRESENT(for_aux_fit)) THEN
810 32 : aux_fit = for_aux_fit
811 : ELSE
812 : aux_fit = .FALSE.
813 : END IF
814 :
815 2706 : CPASSERT(ASSOCIATED(kpoint))
816 :
817 : IF (.TRUE. .OR. ASSOCIATED(mos(1)%mo_coeff)) THEN
818 2706 : IF (aux_fit) THEN
819 32 : CPASSERT(ASSOCIATED(kpoint%kp_aux_env))
820 : END IF
821 :
822 2706 : IF (PRESENT(added_mos)) THEN
823 90 : nadd = added_mos
824 : ELSE
825 : nadd = 0
826 : END IF
827 :
828 2706 : IF (kpoint%use_real_wfn) THEN
829 : nc = 1
830 : ELSE
831 2688 : nc = 2
832 : END IF
833 2706 : nspin = SIZE(mos, 1)
834 2706 : nkp_loc = kpoint%kp_range(2) - kpoint%kp_range(1) + 1
835 2706 : IF (nkp_loc > 0) THEN
836 2706 : IF (aux_fit) THEN
837 32 : CPASSERT(SIZE(kpoint%kp_aux_env) == nkp_loc)
838 : ELSE
839 2674 : CPASSERT(SIZE(kpoint%kp_env) == nkp_loc)
840 : END IF
841 : ! allocate the mo sets, correct number of kpoints (local), real/complex, spin
842 9488 : DO ik = 1, nkp_loc
843 6782 : IF (aux_fit) THEN
844 218 : kp => kpoint%kp_aux_env(ik)%kpoint_env
845 : ELSE
846 6564 : kp => kpoint%kp_env(ik)%kpoint_env
847 : END IF
848 48354 : ALLOCATE (kp%mos(nc, nspin))
849 16570 : DO is = 1, nspin
850 : CALL get_mo_set(mos(is), nao=nao, nmo=nmo, nelectron=nelectron, &
851 7082 : n_el_f=n_el_f, maxocc=maxocc, flexible_electron_count=flexible_electron_count)
852 7082 : nmo = MIN(nao, nmo + nadd)
853 28008 : DO ic = 1, nc
854 : CALL allocate_mo_set(kp%mos(ic, is), nao, nmo, nelectron, n_el_f, maxocc, &
855 21226 : flexible_electron_count)
856 : END DO
857 : END DO
858 : END DO
859 :
860 : ! generate the blacs environment for the kpoint group
861 : ! we generate a blacs env for each kpoint group in parallel
862 : ! we assume here that the group para_env_inter_kp will connect
863 : ! equivalent parts of fm matrices, i.e. no reshuffeling of processors
864 2706 : NULLIFY (blacs_env)
865 2706 : IF (ASSOCIATED(kpoint%blacs_env)) THEN
866 32 : blacs_env => kpoint%blacs_env
867 : ELSE
868 2674 : CALL cp_blacs_env_create(blacs_env=blacs_env, para_env=kpoint%para_env_kp)
869 2674 : kpoint%blacs_env => blacs_env
870 : END IF
871 :
872 : ! set possible new number of MOs
873 5460 : DO is = 1, nspin
874 2754 : CALL get_mo_set(mos(is), nmo=nmorig(is))
875 2754 : nmo = MIN(nao, nmorig(is) + nadd)
876 5460 : CALL set_mo_set(mos(is), nmo=nmo)
877 : END DO
878 : ! matrix pools for the kpoint group, information on MOs is transferred using
879 : ! generic mos structure
880 2706 : NULLIFY (mpools)
881 2706 : CALL mpools_create(mpools=mpools)
882 : CALL mpools_rebuild_fm_pools(mpools=mpools, mos=mos, &
883 2706 : blacs_env=blacs_env, para_env=kpoint%para_env_kp)
884 :
885 2706 : IF (aux_fit) THEN
886 32 : kpoint%mpools_aux_fit => mpools
887 : ELSE
888 2674 : kpoint%mpools => mpools
889 : END IF
890 :
891 : ! reset old number of MOs
892 5460 : DO is = 1, nspin
893 5460 : CALL set_mo_set(mos(is), nmo=nmorig(is))
894 : END DO
895 :
896 : ! allocate density matrices
897 2706 : CALL mpools_get(mpools, ao_ao_fm_pools=ao_ao_fm_pools)
898 2706 : ALLOCATE (fmlocal)
899 2706 : CALL fm_pool_create_fm(ao_ao_fm_pools(1)%pool, fmlocal)
900 2706 : CALL cp_fm_get_info(fmlocal, matrix_struct=matrix_struct)
901 9488 : DO ik = 1, nkp_loc
902 6782 : IF (aux_fit) THEN
903 218 : kp => kpoint%kp_aux_env(ik)%kpoint_env
904 : ELSE
905 6564 : kp => kpoint%kp_env(ik)%kpoint_env
906 : END IF
907 : ! density matrix
908 6782 : CALL cp_fm_release(kp%pmat)
909 48354 : ALLOCATE (kp%pmat(nc, nspin))
910 13864 : DO is = 1, nspin
911 28008 : DO ic = 1, nc
912 21226 : CALL cp_fm_create(kp%pmat(ic, is), matrix_struct)
913 : END DO
914 : END DO
915 : ! energy weighted density matrix
916 6782 : CALL cp_fm_release(kp%wmat)
917 41572 : ALLOCATE (kp%wmat(nc, nspin))
918 16570 : DO is = 1, nspin
919 28008 : DO ic = 1, nc
920 21226 : CALL cp_fm_create(kp%wmat(ic, is), matrix_struct)
921 : END DO
922 : END DO
923 : END DO
924 2706 : CALL fm_pool_give_back_fm(ao_ao_fm_pools(1)%pool, fmlocal)
925 2706 : DEALLOCATE (fmlocal)
926 :
927 : END IF
928 :
929 : END IF
930 :
931 2706 : CALL timestop(handle)
932 :
933 2706 : END SUBROUTINE kpoint_initialize_mos
934 :
935 : ! **************************************************************************************************
936 : !> \brief ...
937 : !> \param kpoint ...
938 : ! **************************************************************************************************
939 106 : SUBROUTINE kpoint_initialize_mo_set(kpoint)
940 : TYPE(kpoint_type), POINTER :: kpoint
941 :
942 : CHARACTER(LEN=*), PARAMETER :: routineN = 'kpoint_initialize_mo_set'
943 :
944 : INTEGER :: handle, ic, ik, ikk, ispin
945 106 : TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_mo_fm_pools
946 : TYPE(cp_fm_type), POINTER :: mo_coeff
947 106 : TYPE(mo_set_type), DIMENSION(:, :), POINTER :: moskp
948 :
949 106 : CALL timeset(routineN, handle)
950 :
951 988 : DO ik = 1, SIZE(kpoint%kp_env)
952 882 : CALL mpools_get(kpoint%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
953 882 : moskp => kpoint%kp_env(ik)%kpoint_env%mos
954 882 : ikk = kpoint%kp_range(1) + ik - 1
955 882 : CPASSERT(ASSOCIATED(moskp))
956 1904 : DO ispin = 1, SIZE(moskp, 2)
957 3630 : DO ic = 1, SIZE(moskp, 1)
958 1832 : CALL get_mo_set(moskp(ic, ispin), mo_coeff=mo_coeff)
959 2748 : IF (.NOT. ASSOCIATED(mo_coeff)) THEN
960 : CALL init_mo_set(moskp(ic, ispin), &
961 1832 : fm_pool=ao_mo_fm_pools(ispin)%pool, name="kpoints")
962 : END IF
963 : END DO
964 : END DO
965 : END DO
966 :
967 106 : CALL timestop(handle)
968 :
969 106 : END SUBROUTINE kpoint_initialize_mo_set
970 :
971 : ! **************************************************************************************************
972 : !> \brief Generates the mapping of cell indices and linear RS index
973 : !> CELL (0,0,0) is always mapped to index 1
974 : !> \param kpoint Kpoint environment
975 : !> \param sab_nl Defining neighbour list
976 : !> \param para_env Parallel environment
977 : !> \param nimages [output]
978 : ! **************************************************************************************************
979 3314 : SUBROUTINE kpoint_init_cell_index(kpoint, sab_nl, para_env, nimages)
980 :
981 : TYPE(kpoint_type), POINTER :: kpoint
982 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
983 : POINTER :: sab_nl
984 : TYPE(mp_para_env_type), POINTER :: para_env
985 : INTEGER, INTENT(OUT) :: nimages
986 :
987 : CHARACTER(LEN=*), PARAMETER :: routineN = 'kpoint_init_cell_index'
988 :
989 : INTEGER :: handle, i1, i2, i3, ic, icount, it, &
990 : ncount
991 : INTEGER, DIMENSION(3) :: cell, itm
992 3314 : INTEGER, DIMENSION(:, :), POINTER :: index_to_cell, list
993 3314 : INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index, cti
994 : LOGICAL :: new
995 : TYPE(neighbor_list_iterator_p_type), &
996 3314 : DIMENSION(:), POINTER :: nl_iterator
997 :
998 3314 : NULLIFY (cell_to_index, index_to_cell)
999 :
1000 3314 : CALL timeset(routineN, handle)
1001 :
1002 3314 : CPASSERT(ASSOCIATED(kpoint))
1003 :
1004 3314 : ALLOCATE (list(3, 125))
1005 1660314 : list = 0
1006 3314 : icount = 1
1007 :
1008 3314 : CALL neighbor_list_iterator_create(nl_iterator, sab_nl)
1009 1221277 : DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1010 1217963 : CALL get_iterator_info(nl_iterator, cell=cell)
1011 :
1012 1217963 : new = .TRUE.
1013 69008333 : DO ic = 1, icount
1014 68838728 : IF (cell(1) == list(1, ic) .AND. cell(2) == list(2, ic) .AND. &
1015 169605 : cell(3) == list(3, ic)) THEN
1016 : new = .FALSE.
1017 : EXIT
1018 : END IF
1019 : END DO
1020 1221277 : IF (new) THEN
1021 169605 : icount = icount + 1
1022 169605 : IF (icount > SIZE(list, 2)) THEN
1023 520 : CALL reallocate(list, 1, 3, 1, 2*SIZE(list, 2))
1024 : END IF
1025 678420 : list(1:3, icount) = cell(1:3)
1026 : END IF
1027 :
1028 : END DO
1029 3314 : CALL neighbor_list_iterator_release(nl_iterator)
1030 :
1031 176233 : itm(1) = MAXVAL(ABS(list(1, 1:icount)))
1032 176233 : itm(2) = MAXVAL(ABS(list(2, 1:icount)))
1033 176233 : itm(3) = MAXVAL(ABS(list(3, 1:icount)))
1034 3314 : CALL para_env%max(itm)
1035 13256 : it = MAXVAL(itm(1:3))
1036 3314 : IF (ASSOCIATED(kpoint%cell_to_index)) THEN
1037 3310 : DEALLOCATE (kpoint%cell_to_index)
1038 : END IF
1039 16570 : ALLOCATE (kpoint%cell_to_index(-itm(1):itm(1), -itm(2):itm(2), -itm(3):itm(3)))
1040 3314 : cell_to_index => kpoint%cell_to_index
1041 3314 : cti => cell_to_index
1042 455812 : cti(:, :, :) = 0
1043 176233 : DO ic = 1, icount
1044 172919 : i1 = list(1, ic)
1045 172919 : i2 = list(2, ic)
1046 172919 : i3 = list(3, ic)
1047 176233 : cti(i1, i2, i3) = ic
1048 : END DO
1049 908310 : CALL para_env%sum(cti)
1050 3314 : ncount = 0
1051 17200 : DO i1 = -itm(1), itm(1)
1052 85742 : DO i2 = -itm(2), itm(2)
1053 453098 : DO i3 = -itm(3), itm(3)
1054 439212 : IF (cti(i1, i2, i3) == 0) THEN
1055 171270 : cti(i1, i2, i3) = 1000000
1056 : ELSE
1057 199400 : ncount = ncount + 1
1058 199400 : cti(i1, i2, i3) = (ABS(i1) + ABS(i2) + ABS(i3))*1000 + ABS(i3)*100 + ABS(i2)*10 + ABS(i1)
1059 199400 : cti(i1, i2, i3) = cti(i1, i2, i3) + (i1 + i2 + i3)
1060 : END IF
1061 : END DO
1062 : END DO
1063 : END DO
1064 :
1065 3314 : IF (ASSOCIATED(kpoint%index_to_cell)) THEN
1066 3314 : DEALLOCATE (kpoint%index_to_cell)
1067 : END IF
1068 9942 : ALLOCATE (kpoint%index_to_cell(3, ncount))
1069 3314 : index_to_cell => kpoint%index_to_cell
1070 202714 : DO ic = 1, ncount
1071 73489452 : cell = MINLOC(cti)
1072 199400 : i1 = cell(1) - 1 - itm(1)
1073 199400 : i2 = cell(2) - 1 - itm(2)
1074 199400 : i3 = cell(3) - 1 - itm(3)
1075 199400 : cti(i1, i2, i3) = 1000000
1076 199400 : index_to_cell(1, ic) = i1
1077 199400 : index_to_cell(2, ic) = i2
1078 202714 : index_to_cell(3, ic) = i3
1079 : END DO
1080 455812 : cti(:, :, :) = 0
1081 202714 : DO ic = 1, ncount
1082 199400 : i1 = index_to_cell(1, ic)
1083 199400 : i2 = index_to_cell(2, ic)
1084 199400 : i3 = index_to_cell(3, ic)
1085 202714 : cti(i1, i2, i3) = ic
1086 : END DO
1087 :
1088 : ! keep pointer to this neighborlist
1089 3314 : kpoint%sab_nl => sab_nl
1090 :
1091 : ! set number of images
1092 3314 : nimages = SIZE(index_to_cell, 2)
1093 :
1094 3314 : DEALLOCATE (list)
1095 :
1096 3314 : CALL timestop(handle)
1097 :
1098 3314 : END SUBROUTINE kpoint_init_cell_index
1099 :
1100 : ! **************************************************************************************************
1101 : !> \brief Transformation of real space matrices to a kpoint
1102 : !> \param rmatrix Real part of kpoint matrix
1103 : !> \param cmatrix Complex part of kpoint matrix (optional)
1104 : !> \param rsmat Real space matrices
1105 : !> \param ispin Spin index
1106 : !> \param xkp Kpoint coordinates
1107 : !> \param cell_to_index mapping of cell indices to RS index
1108 : !> \param sab_nl Defining neighbor list
1109 : !> \param is_complex Matrix to be transformed is imaginary
1110 : !> \param rs_sign Matrix to be transformed is csaled by rs_sign
1111 : ! **************************************************************************************************
1112 480508 : SUBROUTINE rskp_transform(rmatrix, cmatrix, rsmat, ispin, &
1113 : xkp, cell_to_index, sab_nl, is_complex, rs_sign)
1114 :
1115 : TYPE(dbcsr_type) :: rmatrix
1116 : TYPE(dbcsr_type), OPTIONAL :: cmatrix
1117 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: rsmat
1118 : INTEGER, INTENT(IN) :: ispin
1119 : REAL(KIND=dp), DIMENSION(3), INTENT(IN) :: xkp
1120 : INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
1121 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1122 : POINTER :: sab_nl
1123 : LOGICAL, INTENT(IN), OPTIONAL :: is_complex
1124 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: rs_sign
1125 :
1126 : CHARACTER(LEN=*), PARAMETER :: routineN = 'rskp_transform'
1127 :
1128 : INTEGER :: handle, iatom, ic, icol, irow, jatom, &
1129 : nimg
1130 : INTEGER, DIMENSION(3) :: cell
1131 : LOGICAL :: do_symmetric, found, my_complex, &
1132 : wfn_real_only
1133 : REAL(KIND=dp) :: arg, coskl, fsign, fsym, sinkl
1134 240254 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: cblock, rblock, rsblock
1135 : TYPE(neighbor_list_iterator_p_type), &
1136 240254 : DIMENSION(:), POINTER :: nl_iterator
1137 :
1138 240254 : CALL timeset(routineN, handle)
1139 :
1140 240254 : my_complex = .FALSE.
1141 240254 : IF (PRESENT(is_complex)) my_complex = is_complex
1142 :
1143 240254 : fsign = 1.0_dp
1144 240254 : IF (PRESENT(rs_sign)) fsign = rs_sign
1145 :
1146 240254 : wfn_real_only = .TRUE.
1147 240254 : IF (PRESENT(cmatrix)) wfn_real_only = .FALSE.
1148 :
1149 240254 : nimg = SIZE(rsmat, 2)
1150 :
1151 240254 : CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
1152 :
1153 240254 : CALL neighbor_list_iterator_create(nl_iterator, sab_nl)
1154 104304598 : DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1155 104064344 : CALL get_iterator_info(nl_iterator, iatom=iatom, jatom=jatom, cell=cell)
1156 :
1157 : ! fsym = +- 1 is due to real space matrices being non-symmetric (although in a symmtric type)
1158 : ! with the link S_mu^0,nu^b = S_nu^0,mu^-b, and the KP matrices beeing Hermitian
1159 104064344 : fsym = 1.0_dp
1160 104064344 : irow = iatom
1161 104064344 : icol = jatom
1162 104064344 : IF (do_symmetric .AND. (iatom > jatom)) THEN
1163 45411900 : irow = jatom
1164 45411900 : icol = iatom
1165 45411900 : fsym = -1.0_dp
1166 : END IF
1167 :
1168 104064344 : ic = cell_to_index(cell(1), cell(2), cell(3))
1169 104064344 : IF (ic < 1 .OR. ic > nimg) CYCLE
1170 :
1171 104063600 : arg = REAL(cell(1), dp)*xkp(1) + REAL(cell(2), dp)*xkp(2) + REAL(cell(3), dp)*xkp(3)
1172 104063600 : IF (my_complex) THEN
1173 3451728 : coskl = fsign*fsym*COS(twopi*arg)
1174 3451728 : sinkl = fsign*SIN(twopi*arg)
1175 : ELSE
1176 100611872 : coskl = fsign*COS(twopi*arg)
1177 100611872 : sinkl = fsign*fsym*SIN(twopi*arg)
1178 : END IF
1179 :
1180 : CALL dbcsr_get_block_p(matrix=rsmat(ispin, ic)%matrix, row=irow, col=icol, &
1181 104063600 : block=rsblock, found=found)
1182 104063600 : IF (.NOT. found) CYCLE
1183 :
1184 104303854 : IF (wfn_real_only) THEN
1185 : CALL dbcsr_get_block_p(matrix=rmatrix, row=irow, col=icol, &
1186 529850 : block=rblock, found=found)
1187 529850 : IF (.NOT. found) CYCLE
1188 249444630 : rblock = rblock + coskl*rsblock
1189 : ELSE
1190 : CALL dbcsr_get_block_p(matrix=rmatrix, row=irow, col=icol, &
1191 103533750 : block=rblock, found=found)
1192 103533750 : IF (.NOT. found) CYCLE
1193 : CALL dbcsr_get_block_p(matrix=cmatrix, row=irow, col=icol, &
1194 103533750 : block=cblock, found=found)
1195 103533750 : IF (.NOT. found) CYCLE
1196 12372069650 : rblock = rblock + coskl*rsblock
1197 12372069650 : cblock = cblock + sinkl*rsblock
1198 : END IF
1199 :
1200 : END DO
1201 240254 : CALL neighbor_list_iterator_release(nl_iterator)
1202 :
1203 240254 : CALL timestop(handle)
1204 :
1205 240254 : END SUBROUTINE rskp_transform
1206 :
1207 : ! **************************************************************************************************
1208 : !> \brief Given the eigenvalues of all kpoints, calculates the occupation numbers
1209 : !> \param kpoint Kpoint environment
1210 : !> \param smear Smearing information
1211 : !> \param probe ...
1212 : ! **************************************************************************************************
1213 29810 : SUBROUTINE kpoint_set_mo_occupation(kpoint, smear, probe)
1214 :
1215 : TYPE(kpoint_type), POINTER :: kpoint
1216 : TYPE(smear_type) :: smear
1217 : TYPE(hairy_probes_type), DIMENSION(:), OPTIONAL, &
1218 : POINTER :: probe
1219 :
1220 : CHARACTER(LEN=*), PARAMETER :: routineN = 'kpoint_set_mo_occupation'
1221 :
1222 : INTEGER :: handle, ik, ikpgr, ispin, kplocal, nao, &
1223 : nb, ncol_global, ne_a, ne_b, &
1224 : nelectron, nkp, nmo, nrow_global, nspin
1225 : INTEGER, DIMENSION(2) :: kp_range
1226 : REAL(KIND=dp) :: kTS, mu, mus(2), nel
1227 29810 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: smatrix
1228 29810 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: weig, wocc
1229 29810 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :, :) :: icoeff, rcoeff
1230 29810 : REAL(KIND=dp), DIMENSION(:), POINTER :: eigenvalues, occupation, wkp
1231 : TYPE(cp_fm_type), POINTER :: mo_coeff
1232 : TYPE(kpoint_env_type), POINTER :: kp
1233 : TYPE(mo_set_type), POINTER :: mo_set
1234 : TYPE(mp_para_env_type), POINTER :: para_env_inter_kp
1235 :
1236 29810 : CALL timeset(routineN, handle)
1237 :
1238 : ! first collect all the eigenvalues
1239 29810 : CALL get_kpoint_info(kpoint, nkp=nkp)
1240 29810 : kp => kpoint%kp_env(1)%kpoint_env
1241 29810 : nspin = SIZE(kp%mos, 2)
1242 29810 : mo_set => kp%mos(1, 1)
1243 29810 : CALL get_mo_set(mo_set, nmo=nmo, nao=nao, nelectron=nelectron)
1244 29810 : ne_a = nelectron
1245 29810 : IF (nspin == 2) THEN
1246 748 : CALL get_mo_set(kp%mos(1, 2), nmo=nb, nelectron=ne_b)
1247 748 : CPASSERT(nmo == nb)
1248 : END IF
1249 238480 : ALLOCATE (weig(nmo, nkp, nspin), wocc(nmo, nkp, nspin))
1250 29810 : weig = 0.0_dp
1251 29810 : wocc = 0.0_dp
1252 29810 : IF (PRESENT(probe)) THEN
1253 0 : ALLOCATE (rcoeff(nao, nmo, nkp, nspin), icoeff(nao, nmo, nkp, nspin))
1254 0 : rcoeff = 0.0_dp !coeff, real part
1255 0 : icoeff = 0.0_dp !coeff, imaginary part
1256 : END IF
1257 29810 : CALL get_kpoint_info(kpoint, kp_range=kp_range)
1258 29810 : kplocal = kp_range(2) - kp_range(1) + 1
1259 88246 : DO ikpgr = 1, kplocal
1260 58436 : ik = kp_range(1) + ikpgr - 1
1261 58436 : kp => kpoint%kp_env(ikpgr)%kpoint_env
1262 149114 : DO ispin = 1, nspin
1263 60868 : mo_set => kp%mos(1, ispin)
1264 60868 : CALL get_mo_set(mo_set, eigenvalues=eigenvalues)
1265 1265158 : weig(1:nmo, ik, ispin) = eigenvalues(1:nmo)
1266 119304 : IF (PRESENT(probe)) THEN
1267 0 : CALL get_mo_set(mo_set, mo_coeff=mo_coeff)
1268 : CALL cp_fm_get_info(mo_coeff, &
1269 : nrow_global=nrow_global, &
1270 0 : ncol_global=ncol_global)
1271 0 : ALLOCATE (smatrix(nrow_global, ncol_global))
1272 0 : CALL cp_fm_get_submatrix(mo_coeff, smatrix)
1273 :
1274 0 : rcoeff(1:nao, 1:nmo, ik, ispin) = smatrix(1:nrow_global, 1:ncol_global)
1275 :
1276 0 : DEALLOCATE (smatrix)
1277 :
1278 0 : mo_set => kp%mos(2, ispin)
1279 :
1280 0 : CALL get_mo_set(mo_set, mo_coeff=mo_coeff)
1281 : CALL cp_fm_get_info(mo_coeff, &
1282 : nrow_global=nrow_global, &
1283 0 : ncol_global=ncol_global)
1284 0 : ALLOCATE (smatrix(nrow_global, ncol_global))
1285 0 : CALL cp_fm_get_submatrix(mo_coeff, smatrix)
1286 :
1287 0 : icoeff(1:nao, 1:nmo, ik, ispin) = smatrix(1:nrow_global, 1:ncol_global)
1288 :
1289 0 : mo_set => kp%mos(1, ispin)
1290 :
1291 0 : DEALLOCATE (smatrix)
1292 : END IF
1293 : END DO
1294 : END DO
1295 29810 : CALL get_kpoint_info(kpoint, para_env_inter_kp=para_env_inter_kp)
1296 29810 : CALL para_env_inter_kp%sum(weig)
1297 :
1298 29810 : IF (PRESENT(probe)) THEN
1299 0 : CALL para_env_inter_kp%sum(rcoeff)
1300 0 : CALL para_env_inter_kp%sum(icoeff)
1301 : END IF
1302 :
1303 29810 : CALL get_kpoint_info(kpoint, wkp=wkp)
1304 :
1305 : !calling of HP module HERE, before smear
1306 29810 : IF (PRESENT(probe)) THEN
1307 0 : smear%do_smear = .FALSE. !ensures smearing is switched off
1308 :
1309 0 : IF (nspin == 1) THEN
1310 0 : nel = REAL(nelectron, KIND=dp)
1311 : CALL probe_occupancy_kp(wocc(:, :, :), mus(1), kTS, weig(:, :, :), rcoeff(:, :, :, :), icoeff(:, :, :, :), 2.0d0, &
1312 0 : probe, nel, wkp)
1313 : ELSE
1314 0 : nel = REAL(ne_a, KIND=dp) + REAL(ne_b, KIND=dp)
1315 : CALL probe_occupancy_kp(wocc(:, :, :), mu, kTS, weig(:, :, :), rcoeff(:, :, :, :), icoeff(:, :, :, :), 1.0d0, &
1316 0 : probe, nel, wkp)
1317 0 : kTS = kTS/2._dp
1318 0 : mus(1:2) = mu
1319 : END IF
1320 :
1321 0 : DO ikpgr = 1, kplocal
1322 0 : ik = kp_range(1) + ikpgr - 1
1323 0 : kp => kpoint%kp_env(ikpgr)%kpoint_env
1324 0 : DO ispin = 1, nspin
1325 0 : mo_set => kp%mos(1, ispin)
1326 0 : CALL get_mo_set(mo_set, eigenvalues=eigenvalues, occupation_numbers=occupation)
1327 0 : eigenvalues(1:nmo) = weig(1:nmo, ik, ispin)
1328 0 : occupation(1:nmo) = wocc(1:nmo, ik, ispin)
1329 0 : mo_set%kTS = kTS
1330 0 : mo_set%mu = mus(ispin)
1331 :
1332 0 : CALL get_mo_set(mo_set, mo_coeff=mo_coeff)
1333 : !get smatrix for kpoint_env ikp
1334 : CALL cp_fm_get_info(mo_coeff, &
1335 : nrow_global=nrow_global, &
1336 0 : ncol_global=ncol_global)
1337 0 : ALLOCATE (smatrix(nrow_global, ncol_global))
1338 0 : CALL cp_fm_get_submatrix(mo_coeff, smatrix)
1339 :
1340 0 : smatrix(1:nrow_global, 1:ncol_global) = rcoeff(1:nao, 1:nmo, ik, ispin)
1341 0 : DEALLOCATE (smatrix)
1342 :
1343 0 : mo_set => kp%mos(2, ispin)
1344 :
1345 0 : CALL get_mo_set(mo_set, mo_coeff=mo_coeff)
1346 : !get smatrix for kpoint_env ikp
1347 : CALL cp_fm_get_info(mo_coeff, &
1348 : nrow_global=nrow_global, &
1349 0 : ncol_global=ncol_global)
1350 0 : ALLOCATE (smatrix(nrow_global, ncol_global))
1351 0 : CALL cp_fm_get_submatrix(mo_coeff, smatrix)
1352 :
1353 0 : smatrix(1:nrow_global, 1:ncol_global) = icoeff(1:nao, 1:nmo, ik, ispin)
1354 0 : DEALLOCATE (smatrix)
1355 :
1356 0 : mo_set => kp%mos(1, ispin)
1357 :
1358 : END DO
1359 : END DO
1360 :
1361 0 : DEALLOCATE (weig, wocc, rcoeff, icoeff)
1362 :
1363 : END IF
1364 :
1365 : IF (PRESENT(probe) .EQV. .FALSE.) THEN
1366 29810 : IF (smear%do_smear) THEN
1367 19988 : SELECT CASE (smear%method)
1368 : CASE (smear_fermi_dirac)
1369 : ! finite electronic temperature
1370 9930 : IF (nspin == 1) THEN
1371 9852 : nel = REAL(nelectron, KIND=dp)
1372 : CALL Smearkp(wocc(:, :, 1), mus(1), kTS, weig(:, :, 1), nel, wkp, &
1373 9852 : smear%electronic_temperature, 2.0_dp, smear_fermi_dirac)
1374 78 : ELSE IF (smear%fixed_mag_mom > 0.0_dp) THEN
1375 0 : nel = REAL(ne_a, KIND=dp)
1376 : CALL Smearkp(wocc(:, :, 1), mus(1), kTS, weig(:, :, 1), nel, wkp, &
1377 0 : smear%electronic_temperature, 1.0_dp, smear_fermi_dirac)
1378 0 : nel = REAL(ne_b, KIND=dp)
1379 : CALL Smearkp(wocc(:, :, 2), mus(2), kTS, weig(:, :, 2), nel, wkp, &
1380 0 : smear%electronic_temperature, 1.0_dp, smear_fermi_dirac)
1381 : ELSE
1382 78 : nel = REAL(ne_a, KIND=dp) + REAL(ne_b, KIND=dp)
1383 : CALL Smearkp2(wocc(:, :, :), mu, kTS, weig(:, :, :), nel, wkp, &
1384 78 : smear%electronic_temperature, smear_fermi_dirac)
1385 78 : kTS = kTS/2._dp
1386 234 : mus(1:2) = mu
1387 : END IF
1388 : CASE (smear_gaussian, smear_mp, smear_mv)
1389 128 : IF (nspin == 1) THEN
1390 96 : nel = REAL(nelectron, KIND=dp)
1391 : CALL Smearkp(wocc(:, :, 1), mus(1), kTS, weig(:, :, 1), nel, wkp, &
1392 96 : smear%smearing_width, 2.0_dp, smear%method)
1393 32 : ELSE IF (smear%fixed_mag_mom > 0.0_dp) THEN
1394 0 : nel = REAL(ne_a, KIND=dp)
1395 : CALL Smearkp(wocc(:, :, 1), mus(1), kTS, weig(:, :, 1), nel, wkp, &
1396 0 : smear%smearing_width, 1.0_dp, smear%method)
1397 0 : nel = REAL(ne_b, KIND=dp)
1398 : CALL Smearkp(wocc(:, :, 2), mus(2), kTS, weig(:, :, 2), nel, wkp, &
1399 0 : smear%smearing_width, 1.0_dp, smear%method)
1400 : ELSE
1401 32 : nel = REAL(ne_a, KIND=dp) + REAL(ne_b, KIND=dp)
1402 : CALL Smearkp2(wocc(:, :, :), mu, kTS, weig(:, :, :), nel, wkp, &
1403 32 : smear%smearing_width, smear%method)
1404 32 : kTS = kTS/2._dp
1405 96 : mus(1:2) = mu
1406 : END IF
1407 : CASE DEFAULT
1408 10058 : CPABORT("kpoints: Selected smearing not (yet) supported")
1409 : END SELECT
1410 : ELSE
1411 : ! fixed occupations (2/1)
1412 19752 : IF (nspin == 1) THEN
1413 19114 : nel = REAL(nelectron, KIND=dp)
1414 : CALL Smearkp(wocc(:, :, 1), mus(1), kTS, weig(:, :, 1), nel, wkp, &
1415 19114 : 0.0_dp, 2.0_dp, smear_gaussian)
1416 : ELSE
1417 638 : nel = REAL(ne_a, KIND=dp)
1418 : CALL Smearkp(wocc(:, :, 1), mus(1), kTS, weig(:, :, 1), nel, wkp, &
1419 638 : 0.0_dp, 1.0_dp, smear_gaussian)
1420 638 : nel = REAL(ne_b, KIND=dp)
1421 : CALL Smearkp(wocc(:, :, 2), mus(2), kTS, weig(:, :, 2), nel, wkp, &
1422 638 : 0.0_dp, 1.0_dp, smear_gaussian)
1423 : END IF
1424 : END IF
1425 88246 : DO ikpgr = 1, kplocal
1426 58436 : ik = kp_range(1) + ikpgr - 1
1427 58436 : kp => kpoint%kp_env(ikpgr)%kpoint_env
1428 149114 : DO ispin = 1, nspin
1429 60868 : mo_set => kp%mos(1, ispin)
1430 60868 : CALL get_mo_set(mo_set, eigenvalues=eigenvalues, occupation_numbers=occupation)
1431 1265158 : eigenvalues(1:nmo) = weig(1:nmo, ik, ispin)
1432 1265158 : occupation(1:nmo) = wocc(1:nmo, ik, ispin)
1433 60868 : mo_set%kTS = kTS
1434 119304 : mo_set%mu = mus(ispin)
1435 : END DO
1436 : END DO
1437 :
1438 29810 : DEALLOCATE (weig, wocc)
1439 :
1440 : END IF
1441 :
1442 29810 : CALL timestop(handle)
1443 :
1444 89430 : END SUBROUTINE kpoint_set_mo_occupation
1445 :
1446 : ! **************************************************************************************************
1447 : !> \brief Calculate kpoint density matrices (rho(k), owned by kpoint groups)
1448 : !> \param kpoint kpoint environment
1449 : !> \param energy_weighted calculate energy weighted density matrix
1450 : !> \param for_aux_fit ...
1451 : ! **************************************************************************************************
1452 90984 : SUBROUTINE kpoint_density_matrices(kpoint, energy_weighted, for_aux_fit)
1453 :
1454 : TYPE(kpoint_type), POINTER :: kpoint
1455 : LOGICAL, OPTIONAL :: energy_weighted, for_aux_fit
1456 :
1457 : CHARACTER(LEN=*), PARAMETER :: routineN = 'kpoint_density_matrices'
1458 :
1459 : INTEGER :: handle, ikpgr, ispin, kplocal, nao, nmo, &
1460 : nspin
1461 : INTEGER, DIMENSION(2) :: kp_range
1462 : LOGICAL :: aux_fit, wtype
1463 30328 : REAL(KIND=dp), DIMENSION(:), POINTER :: eigenvalues, occupation
1464 : TYPE(cp_fm_struct_type), POINTER :: matrix_struct
1465 : TYPE(cp_fm_type) :: fwork
1466 : TYPE(cp_fm_type), POINTER :: cpmat, pmat, rpmat
1467 : TYPE(kpoint_env_type), POINTER :: kp
1468 : TYPE(mo_set_type), POINTER :: mo_set
1469 :
1470 30328 : CALL timeset(routineN, handle)
1471 :
1472 30328 : IF (PRESENT(energy_weighted)) THEN
1473 378 : wtype = energy_weighted
1474 : ELSE
1475 : ! default is normal density matrix
1476 : wtype = .FALSE.
1477 : END IF
1478 :
1479 30328 : IF (PRESENT(for_aux_fit)) THEN
1480 124 : aux_fit = for_aux_fit
1481 : ELSE
1482 : aux_fit = .FALSE.
1483 : END IF
1484 :
1485 124 : IF (aux_fit) THEN
1486 124 : CPASSERT(ASSOCIATED(kpoint%kp_aux_env))
1487 : END IF
1488 :
1489 : ! work matrix
1490 30328 : IF (aux_fit) THEN
1491 124 : mo_set => kpoint%kp_aux_env(1)%kpoint_env%mos(1, 1)
1492 : ELSE
1493 30204 : mo_set => kpoint%kp_env(1)%kpoint_env%mos(1, 1)
1494 : END IF
1495 30328 : CALL get_mo_set(mo_set, nao=nao, nmo=nmo)
1496 30328 : CALL cp_fm_get_info(mo_set%mo_coeff, matrix_struct=matrix_struct)
1497 30328 : CALL cp_fm_create(fwork, matrix_struct)
1498 :
1499 30328 : CALL get_kpoint_info(kpoint, kp_range=kp_range)
1500 30328 : kplocal = kp_range(2) - kp_range(1) + 1
1501 92022 : DO ikpgr = 1, kplocal
1502 61694 : IF (aux_fit) THEN
1503 1876 : kp => kpoint%kp_aux_env(ikpgr)%kpoint_env
1504 : ELSE
1505 59818 : kp => kpoint%kp_env(ikpgr)%kpoint_env
1506 : END IF
1507 61694 : nspin = SIZE(kp%mos, 2)
1508 156400 : DO ispin = 1, nspin
1509 64378 : mo_set => kp%mos(1, ispin)
1510 64378 : IF (wtype) THEN
1511 1406 : CALL get_mo_set(mo_set, eigenvalues=eigenvalues)
1512 : END IF
1513 126072 : IF (kpoint%use_real_wfn) THEN
1514 168 : IF (wtype) THEN
1515 12 : pmat => kp%wmat(1, ispin)
1516 : ELSE
1517 156 : pmat => kp%pmat(1, ispin)
1518 : END IF
1519 168 : CALL get_mo_set(mo_set, occupation_numbers=occupation)
1520 168 : CALL cp_fm_to_fm(mo_set%mo_coeff, fwork)
1521 168 : CALL cp_fm_column_scale(fwork, occupation)
1522 168 : IF (wtype) THEN
1523 12 : CALL cp_fm_column_scale(fwork, eigenvalues)
1524 : END IF
1525 168 : CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 0.0_dp, pmat)
1526 : ELSE
1527 64210 : IF (wtype) THEN
1528 1394 : rpmat => kp%wmat(1, ispin)
1529 1394 : cpmat => kp%wmat(2, ispin)
1530 : ELSE
1531 62816 : rpmat => kp%pmat(1, ispin)
1532 62816 : cpmat => kp%pmat(2, ispin)
1533 : END IF
1534 64210 : CALL get_mo_set(mo_set, occupation_numbers=occupation)
1535 64210 : CALL cp_fm_to_fm(mo_set%mo_coeff, fwork)
1536 64210 : CALL cp_fm_column_scale(fwork, occupation)
1537 64210 : IF (wtype) THEN
1538 1394 : CALL cp_fm_column_scale(fwork, eigenvalues)
1539 : END IF
1540 : ! Re(c)*Re(c)
1541 64210 : CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 0.0_dp, rpmat)
1542 64210 : mo_set => kp%mos(2, ispin)
1543 : ! Im(c)*Re(c)
1544 64210 : CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 0.0_dp, cpmat)
1545 : ! Re(c)*Im(c)
1546 64210 : CALL parallel_gemm("N", "T", nao, nao, nmo, -1.0_dp, fwork, mo_set%mo_coeff, 1.0_dp, cpmat)
1547 64210 : CALL cp_fm_to_fm(mo_set%mo_coeff, fwork)
1548 64210 : CALL cp_fm_column_scale(fwork, occupation)
1549 64210 : IF (wtype) THEN
1550 1394 : CALL cp_fm_column_scale(fwork, eigenvalues)
1551 : END IF
1552 : ! Im(c)*Im(c)
1553 64210 : CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 1.0_dp, rpmat)
1554 : END IF
1555 : END DO
1556 : END DO
1557 :
1558 30328 : CALL cp_fm_release(fwork)
1559 :
1560 30328 : CALL timestop(handle)
1561 :
1562 30328 : END SUBROUTINE kpoint_density_matrices
1563 :
1564 : ! **************************************************************************************************
1565 : !> \brief Calculate Lowdin transformation of density matrix S^1/2 P S^1/2
1566 : !> Integrate diagonal elements over k-points to get Lowdin charges
1567 : !> \param kpoint kpoint environment
1568 : !> \param pmat_diag Sum over kpoints of diagonal elements
1569 : !> \par History
1570 : !> 04.2026 created [JGH]
1571 : ! **************************************************************************************************
1572 6 : SUBROUTINE lowdin_kp_trans(kpoint, pmat_diag)
1573 :
1574 : TYPE(kpoint_type), POINTER :: kpoint
1575 : REAL(KIND=dp), DIMENSION(:, :), INTENT(INOUT) :: pmat_diag
1576 :
1577 : CHARACTER(LEN=*), PARAMETER :: routineN = 'lowdin_kp_trans'
1578 : COMPLEX(KIND=dp), PARAMETER :: cone = (1.0_dp, 0.0_dp), &
1579 : czero = (0.0_dp, 0.0_dp)
1580 :
1581 : INTEGER :: handle, ikpgr, ispin, kplocal, nao, nspin
1582 : INTEGER, DIMENSION(2) :: kp_range
1583 6 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: dele
1584 : TYPE(cp_cfm_type) :: cf1work, cf2work
1585 : TYPE(cp_cfm_type), POINTER :: cshalf
1586 : TYPE(cp_fm_struct_type), POINTER :: matrix_struct
1587 : TYPE(cp_fm_type) :: f1work, f2work
1588 : TYPE(cp_fm_type), POINTER :: cpmat, pmat, rpmat, shalf
1589 : TYPE(kpoint_env_type), POINTER :: kp
1590 : TYPE(mp_para_env_type), POINTER :: para_env_inter_kp
1591 :
1592 6 : CALL timeset(routineN, handle)
1593 :
1594 6 : nspin = SIZE(pmat_diag, 2)
1595 336 : pmat_diag = 0.0_dp
1596 :
1597 : ! work matrix
1598 : CALL cp_fm_get_info(kpoint%kp_env(1)%kpoint_env%pmat(1, 1), &
1599 6 : matrix_struct=matrix_struct, nrow_global=nao)
1600 6 : IF (kpoint%use_real_wfn) THEN
1601 0 : CALL cp_fm_create(f1work, matrix_struct, nrow=nao, ncol=nao)
1602 0 : CALL cp_fm_create(f2work, matrix_struct, nrow=nao, ncol=nao)
1603 : ELSE
1604 6 : CALL cp_fm_create(f2work, matrix_struct, nrow=nao, ncol=nao)
1605 6 : CALL cp_cfm_create(cf1work, matrix_struct, nrow=nao, ncol=nao)
1606 6 : CALL cp_cfm_create(cf2work, matrix_struct, nrow=nao, ncol=nao)
1607 : END IF
1608 18 : ALLOCATE (dele(nao))
1609 :
1610 6 : CALL get_kpoint_info(kpoint, kp_range=kp_range)
1611 6 : kplocal = kp_range(2) - kp_range(1) + 1
1612 238 : DO ikpgr = 1, kplocal
1613 232 : kp => kpoint%kp_env(ikpgr)%kpoint_env
1614 486 : DO ispin = 1, nspin
1615 248 : IF (kpoint%use_real_wfn) THEN
1616 0 : pmat => kp%pmat(1, ispin)
1617 0 : shalf => kp%shalf
1618 0 : CALL parallel_gemm("N", "N", nao, nao, nao, 1.0_dp, pmat, shalf, 0.0_dp, f1work)
1619 0 : CALL parallel_gemm("N", "N", nao, nao, nao, 1.0_dp, shalf, f1work, 0.0_dp, f2work)
1620 : ELSE
1621 248 : rpmat => kp%pmat(1, ispin)
1622 248 : cpmat => kp%pmat(2, ispin)
1623 248 : cshalf => kp%cshalf
1624 248 : CALL cp_fm_to_cfm(rpmat, cpmat, cf1work)
1625 248 : CALL parallel_gemm("N", "N", nao, nao, nao, cone, cf1work, cshalf, czero, cf2work)
1626 248 : CALL parallel_gemm("N", "N", nao, nao, nao, cone, cshalf, cf2work, czero, cf1work)
1627 248 : CALL cp_cfm_to_fm(cf1work, mtargetr=f2work)
1628 : END IF
1629 248 : CALL cp_fm_get_diag(f2work, dele)
1630 2592 : pmat_diag(1:nao, ispin) = pmat_diag(1:nao, ispin) + kp%wkp*dele(1:nao)
1631 : END DO
1632 : END DO
1633 :
1634 6 : CALL get_kpoint_info(kpoint, para_env_inter_kp=para_env_inter_kp)
1635 666 : CALL para_env_inter_kp%sum(pmat_diag)
1636 :
1637 6 : IF (kpoint%use_real_wfn) THEN
1638 0 : CALL cp_fm_release(f1work)
1639 0 : CALL cp_fm_release(f2work)
1640 : ELSE
1641 6 : CALL cp_fm_release(f2work)
1642 6 : CALL cp_cfm_release(cf1work)
1643 6 : CALL cp_cfm_release(cf2work)
1644 : END IF
1645 6 : DEALLOCATE (dele)
1646 :
1647 6 : CALL timestop(handle)
1648 :
1649 12 : END SUBROUTINE lowdin_kp_trans
1650 :
1651 : ! **************************************************************************************************
1652 : !> \brief Calculate S(k)^1/2 C(k) for real or complex k-point wavefunctions
1653 : !> \param kp K-point environment for one local k point
1654 : !> \param ispin Spin index
1655 : !> \param use_real_wfn Use real k-point wavefunctions
1656 : !> \param shalfc Output matrix containing S(k)^1/2 C(k) for real wavefunctions
1657 : !> \param cshalfc Output matrix containing S(k)^1/2 C(k) for complex wavefunctions
1658 : ! **************************************************************************************************
1659 100 : SUBROUTINE lowdin_kp_mo_coeff(kp, ispin, use_real_wfn, shalfc, cshalfc)
1660 :
1661 : TYPE(kpoint_env_type), POINTER :: kp
1662 : INTEGER, INTENT(IN) :: ispin
1663 : LOGICAL, INTENT(IN) :: use_real_wfn
1664 : TYPE(cp_fm_type), INTENT(INOUT), OPTIONAL :: shalfc
1665 : TYPE(cp_cfm_type), INTENT(INOUT), OPTIONAL :: cshalfc
1666 :
1667 : INTEGER :: nao, nmo
1668 : TYPE(cp_fm_struct_type), POINTER :: matrix_struct_mo, matrix_struct_shalf
1669 : TYPE(cp_fm_type) :: cshalf_im, cshalf_re, shalf_im, shalf_re
1670 : TYPE(mo_set_type), POINTER :: mo_set, mo_set_im, mo_set_re
1671 :
1672 600 : IF (use_real_wfn) THEN
1673 0 : CPASSERT(PRESENT(shalfc))
1674 0 : mo_set => kp%mos(1, ispin)
1675 0 : CALL get_mo_set(mo_set, nao=nao, nmo=nmo)
1676 :
1677 : CALL parallel_gemm("N", "N", nao, nmo, nao, 1.0_dp, kp%shalf, &
1678 0 : mo_set%mo_coeff, 0.0_dp, shalfc)
1679 : ELSE
1680 100 : CPASSERT(PRESENT(cshalfc))
1681 100 : mo_set_re => kp%mos(1, ispin)
1682 100 : mo_set_im => kp%mos(2, ispin)
1683 100 : CALL get_mo_set(mo_set_re, nao=nao, nmo=nmo)
1684 100 : CALL cp_fm_get_info(mo_set_re%mo_coeff, matrix_struct=matrix_struct_mo)
1685 100 : CALL cp_cfm_get_info(kp%cshalf, matrix_struct=matrix_struct_shalf)
1686 :
1687 100 : CALL cp_fm_create(shalf_re, matrix_struct_shalf, nrow=nao, ncol=nao)
1688 100 : CALL cp_fm_create(shalf_im, matrix_struct_shalf, nrow=nao, ncol=nao)
1689 100 : CALL cp_fm_create(cshalf_re, matrix_struct_mo, nrow=nao, ncol=nmo)
1690 100 : CALL cp_fm_create(cshalf_im, matrix_struct_mo, nrow=nao, ncol=nmo)
1691 :
1692 100 : CALL cp_cfm_to_fm(kp%cshalf, mtargetr=shalf_re, mtargeti=shalf_im)
1693 :
1694 : ! Re[S(k)^1/2 C(k)] = Re[S(k)^1/2] C_re(k) - Im[S(k)^1/2] C_im(k)
1695 : CALL parallel_gemm("N", "N", nao, nmo, nao, 1.0_dp, shalf_re, &
1696 100 : mo_set_re%mo_coeff, 0.0_dp, cshalf_re)
1697 : CALL parallel_gemm("N", "N", nao, nmo, nao, -1.0_dp, shalf_im, &
1698 100 : mo_set_im%mo_coeff, 1.0_dp, cshalf_re)
1699 :
1700 : ! Im[S(k)^1/2 C(k)] = Re[S(k)^1/2] C_im(k) + Im[S(k)^1/2] C_re(k)
1701 : CALL parallel_gemm("N", "N", nao, nmo, nao, 1.0_dp, shalf_re, &
1702 100 : mo_set_im%mo_coeff, 0.0_dp, cshalf_im)
1703 : CALL parallel_gemm("N", "N", nao, nmo, nao, 1.0_dp, shalf_im, &
1704 100 : mo_set_re%mo_coeff, 1.0_dp, cshalf_im)
1705 :
1706 100 : CALL cp_fm_to_cfm(cshalf_re, cshalf_im, cshalfc)
1707 :
1708 100 : CALL cp_fm_release(shalf_re)
1709 100 : CALL cp_fm_release(shalf_im)
1710 100 : CALL cp_fm_release(cshalf_re)
1711 100 : CALL cp_fm_release(cshalf_im)
1712 : END IF
1713 :
1714 100 : END SUBROUTINE lowdin_kp_mo_coeff
1715 :
1716 : ! **************************************************************************************************
1717 : !> \brief generate real space density matrices in DBCSR format
1718 : !> \param kpoint Kpoint environment
1719 : !> \param denmat Real space (DBCSR) density matrices
1720 : !> \param wtype True = energy weighted density matrix
1721 : !> False = normal density matrix
1722 : !> \param tempmat DBCSR matrix to be used as template
1723 : !> \param sab_nl ...
1724 : !> \param fmwork FM work matrices (kpoint group)
1725 : !> \param for_aux_fit ...
1726 : !> \param pmat_ext ...
1727 : ! **************************************************************************************************
1728 30576 : SUBROUTINE kpoint_density_transform(kpoint, denmat, wtype, tempmat, sab_nl, fmwork, for_aux_fit, pmat_ext)
1729 :
1730 : TYPE(kpoint_type), POINTER :: kpoint
1731 : TYPE(dbcsr_p_type), DIMENSION(:, :) :: denmat
1732 : LOGICAL, INTENT(IN) :: wtype
1733 : TYPE(dbcsr_type), POINTER :: tempmat
1734 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1735 : POINTER :: sab_nl
1736 : TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fmwork
1737 : LOGICAL, OPTIONAL :: for_aux_fit
1738 : TYPE(cp_fm_type), DIMENSION(:, :, :), INTENT(IN), &
1739 : OPTIONAL :: pmat_ext
1740 :
1741 : CHARACTER(LEN=*), PARAMETER :: routineN = 'kpoint_density_transform'
1742 :
1743 : INTEGER :: handle, ic, ik, ikk, indx, ir, ira, is, &
1744 : ispin, jr, nc, nimg, nkp, nspin
1745 30576 : INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
1746 : LOGICAL :: aux_fit, do_ext, do_symmetric, my_kpgrp, &
1747 : real_only, reverse_phase
1748 : REAL(KIND=dp) :: wkpx
1749 30576 : REAL(KIND=dp), DIMENSION(:), POINTER :: wkp
1750 30576 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: xkp
1751 30576 : TYPE(copy_info_type), ALLOCATABLE, DIMENSION(:) :: info
1752 : TYPE(cp_fm_type) :: fmdummy
1753 : TYPE(dbcsr_type), POINTER :: cpmat, rpmat, scpmat, srpmat
1754 30576 : TYPE(kind_rotmat_type), DIMENSION(:), POINTER :: kind_rot
1755 : TYPE(kpoint_env_type), POINTER :: kp
1756 : TYPE(kpoint_sym_type), POINTER :: kpsym
1757 : TYPE(mp_para_env_type), POINTER :: para_env
1758 :
1759 30576 : CALL timeset(routineN, handle)
1760 :
1761 30576 : CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
1762 :
1763 30576 : IF (PRESENT(for_aux_fit)) THEN
1764 372 : aux_fit = for_aux_fit
1765 : ELSE
1766 : aux_fit = .FALSE.
1767 : END IF
1768 :
1769 30576 : do_ext = .FALSE.
1770 30576 : IF (PRESENT(pmat_ext)) do_ext = .TRUE.
1771 :
1772 30576 : IF (aux_fit) THEN
1773 216 : CPASSERT(ASSOCIATED(kpoint%kp_aux_env))
1774 : END IF
1775 :
1776 : ! work storage
1777 30576 : ALLOCATE (rpmat)
1778 : CALL dbcsr_create(rpmat, template=tempmat, &
1779 30636 : matrix_type=MERGE(dbcsr_type_symmetric, dbcsr_type_no_symmetry, do_symmetric))
1780 30576 : CALL cp_dbcsr_alloc_block_from_nbl(rpmat, sab_nl)
1781 30576 : CALL dbcsr_set(rpmat, 0.0_dp)
1782 30576 : ALLOCATE (cpmat)
1783 : CALL dbcsr_create(cpmat, template=tempmat, &
1784 30636 : matrix_type=MERGE(dbcsr_type_antisymmetric, dbcsr_type_no_symmetry, do_symmetric))
1785 30576 : CALL cp_dbcsr_alloc_block_from_nbl(cpmat, sab_nl)
1786 30576 : CALL dbcsr_set(cpmat, 0.0_dp)
1787 30576 : IF (.NOT. kpoint%full_grid) THEN
1788 28214 : ALLOCATE (srpmat)
1789 28214 : CALL dbcsr_create(srpmat, template=rpmat)
1790 28214 : CALL cp_dbcsr_alloc_block_from_nbl(srpmat, sab_nl)
1791 28214 : CALL dbcsr_set(srpmat, 0.0_dp)
1792 28214 : ALLOCATE (scpmat)
1793 28214 : CALL dbcsr_create(scpmat, template=cpmat)
1794 28214 : CALL cp_dbcsr_alloc_block_from_nbl(scpmat, sab_nl)
1795 28214 : CALL dbcsr_set(scpmat, 0.0_dp)
1796 : END IF
1797 :
1798 : CALL get_kpoint_info(kpoint, nkp=nkp, xkp=xkp, wkp=wkp, &
1799 30576 : cell_to_index=cell_to_index)
1800 : ! initialize real space density matrices
1801 30576 : IF (aux_fit) THEN
1802 216 : kp => kpoint%kp_aux_env(1)%kpoint_env
1803 : ELSE
1804 30360 : kp => kpoint%kp_env(1)%kpoint_env
1805 : END IF
1806 30576 : nspin = SIZE(kp%mos, 2)
1807 30576 : nc = SIZE(kp%mos, 1)
1808 30576 : nimg = SIZE(denmat, 2)
1809 30576 : real_only = (nc == 1)
1810 : ! Gamma-centered even grids store atom-cell shifts in the opposite Bloch-phase convention.
1811 41022 : reverse_phase = kpoint%gamma_centered .AND. ANY(MOD(kpoint%nkp_grid(1:3), 2) == 0)
1812 :
1813 30576 : para_env => kpoint%blacs_env_all%para_env
1814 555612 : ALLOCATE (info(nspin*nkp*nc))
1815 :
1816 : ! Start all the communication
1817 30576 : indx = 0
1818 61982 : DO ispin = 1, nspin
1819 1448164 : DO ic = 1, nimg
1820 1448164 : CALL dbcsr_set(denmat(ispin, ic)%matrix, 0.0_dp)
1821 : END DO
1822 : !
1823 171704 : DO ik = 1, nkp
1824 109722 : my_kpgrp = (ik >= kpoint%kp_range(1) .AND. ik <= kpoint%kp_range(2))
1825 : IF (my_kpgrp) THEN
1826 67358 : ikk = ik - kpoint%kp_range(1) + 1
1827 67358 : IF (aux_fit) THEN
1828 2714 : kp => kpoint%kp_aux_env(ikk)%kpoint_env
1829 : ELSE
1830 64644 : kp => kpoint%kp_env(ikk)%kpoint_env
1831 : END IF
1832 : ELSE
1833 : NULLIFY (kp)
1834 : END IF
1835 : ! collect this density matrix on all processors
1836 109722 : CPASSERT(SIZE(fmwork) >= nc)
1837 :
1838 141128 : IF (my_kpgrp) THEN
1839 201906 : DO ic = 1, nc
1840 134548 : indx = indx + 1
1841 201906 : IF (do_ext) THEN
1842 5428 : CALL cp_fm_start_copy_general(pmat_ext(ikk, ic, ispin), fmwork(ic), para_env, info(indx))
1843 : ELSE
1844 129120 : IF (wtype) THEN
1845 2800 : CALL cp_fm_start_copy_general(kp%wmat(ic, ispin), fmwork(ic), para_env, info(indx))
1846 : ELSE
1847 126320 : CALL cp_fm_start_copy_general(kp%pmat(ic, ispin), fmwork(ic), para_env, info(indx))
1848 : END IF
1849 : END IF
1850 : END DO
1851 : ELSE
1852 127092 : DO ic = 1, nc
1853 84728 : indx = indx + 1
1854 127092 : CALL cp_fm_start_copy_general(fmdummy, fmwork(ic), para_env, info(indx))
1855 : END DO
1856 : END IF
1857 : END DO
1858 : END DO
1859 :
1860 : ! Finish communication and transform the received matrices
1861 30576 : indx = 0
1862 61982 : DO ispin = 1, nspin
1863 171704 : DO ik = 1, nkp
1864 328998 : DO ic = 1, nc
1865 219276 : indx = indx + 1
1866 328998 : CALL cp_fm_finish_copy_general(fmwork(ic), info(indx))
1867 : END DO
1868 :
1869 : ! reduce to dbcsr storage
1870 109722 : IF (real_only) THEN
1871 168 : CALL copy_fm_to_dbcsr(fmwork(1), rpmat, keep_sparsity=.TRUE.)
1872 : ELSE
1873 109554 : CALL copy_fm_to_dbcsr(fmwork(1), rpmat, keep_sparsity=.TRUE.)
1874 109554 : CALL copy_fm_to_dbcsr(fmwork(2), cpmat, keep_sparsity=.TRUE.)
1875 : END IF
1876 :
1877 : ! symmetrization
1878 109722 : kpsym => kpoint%kp_sym(ik)%kpoint_sym
1879 109722 : CPASSERT(ASSOCIATED(kpsym))
1880 :
1881 141128 : IF (kpsym%apply_symmetry) THEN
1882 4154 : wkpx = wkp(ik)/REAL(kpsym%nwght, KIND=dp)
1883 33622 : DO is = 1, kpsym%nwght
1884 29468 : ir = ABS(kpsym%rotp(is))
1885 29468 : ira = 0
1886 4083332 : DO jr = 1, SIZE(kpoint%ibrot)
1887 4083332 : IF (ir == kpoint%ibrot(jr)) ira = jr
1888 : END DO
1889 29468 : CPASSERT(ira > 0)
1890 29468 : kind_rot => kpoint%kind_rotmat(ira, :)
1891 : CALL symtrans_phase(srpmat, scpmat, rpmat, cpmat, real_only, kind_rot, &
1892 : kpsym%rot(1:3, 1:3, is), kpsym%f0(:, is), &
1893 : kpsym%fcell(:, :, is), kpoint%atype, kpsym%xkp(1:3, is), &
1894 29468 : kpsym%rotp(is) < 0, reverse_phase)
1895 : CALL transform_dmat(denmat, srpmat, scpmat, ispin, real_only, sab_nl, &
1896 33622 : cell_to_index, kpsym%xkp(1:3, is), wkpx)
1897 : END DO
1898 : ELSE
1899 : ! transformation
1900 : CALL transform_dmat(denmat, rpmat, cpmat, ispin, real_only, sab_nl, &
1901 105568 : cell_to_index, xkp(1:3, ik), wkp(ik))
1902 : END IF
1903 : END DO
1904 : END DO
1905 :
1906 : ! Clean up communication
1907 30576 : indx = 0
1908 61982 : DO ispin = 1, nspin
1909 171704 : DO ik = 1, nkp
1910 109722 : my_kpgrp = (ik >= kpoint%kp_range(1) .AND. ik <= kpoint%kp_range(2))
1911 31406 : IF (my_kpgrp) THEN
1912 201906 : ikk = ik - kpoint%kp_range(1) + 1
1913 : IF (aux_fit) THEN
1914 201906 : kp => kpoint%kp_aux_env(ikk)%kpoint_env
1915 : ELSE
1916 201906 : kp => kpoint%kp_env(ikk)%kpoint_env
1917 : END IF
1918 :
1919 201906 : DO ic = 1, nc
1920 134548 : indx = indx + 1
1921 201906 : CALL cp_fm_cleanup_copy_general(info(indx))
1922 : END DO
1923 : ELSE
1924 : ! calls with dummy arguments, so not included
1925 : ! therefore just increment counter by trip count
1926 42364 : indx = indx + nc
1927 : END IF
1928 : END DO
1929 : END DO
1930 :
1931 : ! All done
1932 249852 : DEALLOCATE (info)
1933 :
1934 30576 : CALL dbcsr_deallocate_matrix(rpmat)
1935 30576 : CALL dbcsr_deallocate_matrix(cpmat)
1936 30576 : IF (.NOT. kpoint%full_grid) THEN
1937 28214 : CALL dbcsr_deallocate_matrix(srpmat)
1938 28214 : CALL dbcsr_deallocate_matrix(scpmat)
1939 : END IF
1940 :
1941 30576 : CALL timestop(handle)
1942 :
1943 30576 : END SUBROUTINE kpoint_density_transform
1944 :
1945 : ! **************************************************************************************************
1946 : !> \brief real space density matrices in DBCSR format
1947 : !> \param denmat Real space (DBCSR) density matrix
1948 : !> \param rpmat ...
1949 : !> \param cpmat ...
1950 : !> \param ispin ...
1951 : !> \param real_only ...
1952 : !> \param sab_nl ...
1953 : !> \param cell_to_index ...
1954 : !> \param xkp ...
1955 : !> \param wkp ...
1956 : ! **************************************************************************************************
1957 135036 : SUBROUTINE transform_dmat(denmat, rpmat, cpmat, ispin, real_only, sab_nl, cell_to_index, xkp, wkp)
1958 :
1959 : TYPE(dbcsr_p_type), DIMENSION(:, :) :: denmat
1960 : TYPE(dbcsr_type), POINTER :: rpmat, cpmat
1961 : INTEGER, INTENT(IN) :: ispin
1962 : LOGICAL, INTENT(IN) :: real_only
1963 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1964 : POINTER :: sab_nl
1965 : INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
1966 : REAL(KIND=dp), DIMENSION(3), INTENT(IN) :: xkp
1967 : REAL(KIND=dp), INTENT(IN) :: wkp
1968 :
1969 : CHARACTER(LEN=*), PARAMETER :: routineN = 'transform_dmat'
1970 :
1971 : INTEGER :: handle, iatom, icell, icol, irow, jatom, &
1972 : nimg
1973 : INTEGER, DIMENSION(3) :: cell
1974 : LOGICAL :: do_symmetric, found
1975 : REAL(KIND=dp) :: arg, coskl, fc, sinkl
1976 135036 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: cblock, dblock, rblock
1977 : TYPE(neighbor_list_iterator_p_type), &
1978 135036 : DIMENSION(:), POINTER :: nl_iterator
1979 :
1980 135036 : CALL timeset(routineN, handle)
1981 :
1982 135036 : nimg = SIZE(denmat, 2)
1983 :
1984 : ! transformation
1985 135036 : CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
1986 135036 : CALL neighbor_list_iterator_create(nl_iterator, sab_nl)
1987 63645835 : DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1988 63510799 : CALL get_iterator_info(nl_iterator, iatom=iatom, jatom=jatom, cell=cell)
1989 :
1990 : !We have a FT from KP to real-space: S(R) = sum_k S(k)*exp(-i*k*R), with S(k) a complex number
1991 : !Therefore, we have: S(R) = sum_k Re(S(k))*cos(k*R) -i^2*Im(S(k))*sin(k*R)
1992 : ! = sum_k Re(S(k))*cos(k*R) + Im(S(k))*sin(k*R)
1993 : !fc = +- 1 is due to the usual non-symmetric real-space matrices stored as symmetric ones
1994 :
1995 63510799 : irow = iatom
1996 63510799 : icol = jatom
1997 63510799 : fc = 1.0_dp
1998 63510799 : IF (do_symmetric .AND. iatom > jatom) THEN
1999 27697112 : irow = jatom
2000 27697112 : icol = iatom
2001 27697112 : fc = -1.0_dp
2002 : END IF
2003 :
2004 63510799 : icell = cell_to_index(cell(1), cell(2), cell(3))
2005 63510799 : IF (icell < 1 .OR. icell > nimg) CYCLE
2006 :
2007 63509521 : arg = REAL(cell(1), dp)*xkp(1) + REAL(cell(2), dp)*xkp(2) + REAL(cell(3), dp)*xkp(3)
2008 63509521 : coskl = wkp*COS(twopi*arg)
2009 63509521 : sinkl = wkp*fc*SIN(twopi*arg)
2010 :
2011 : CALL dbcsr_get_block_p(matrix=denmat(ispin, icell)%matrix, row=irow, col=icol, &
2012 63509521 : block=dblock, found=found)
2013 63509521 : IF (.NOT. found) CYCLE
2014 :
2015 63644557 : IF (real_only) THEN
2016 294113 : CALL dbcsr_get_block_p(matrix=rpmat, row=irow, col=icol, block=rblock, found=found)
2017 294113 : IF (.NOT. found) CYCLE
2018 142452095 : dblock = dblock + coskl*rblock
2019 : ELSE
2020 63215408 : CALL dbcsr_get_block_p(matrix=rpmat, row=irow, col=icol, block=rblock, found=found)
2021 63215408 : IF (.NOT. found) CYCLE
2022 63215408 : CALL dbcsr_get_block_p(matrix=cpmat, row=irow, col=icol, block=cblock, found=found)
2023 63215408 : IF (.NOT. found) CYCLE
2024 10087499886 : dblock = dblock + coskl*rblock
2025 10087499886 : dblock = dblock + sinkl*cblock
2026 : END IF
2027 : END DO
2028 135036 : CALL neighbor_list_iterator_release(nl_iterator)
2029 :
2030 135036 : CALL timestop(handle)
2031 :
2032 135036 : END SUBROUTINE transform_dmat
2033 :
2034 : ! **************************************************************************************************
2035 : !> \brief Allocate a dense work matrix with the requested shape
2036 : !> \param work dense work matrix
2037 : !> \param nrow number of rows
2038 : !> \param ncol number of columns
2039 : ! **************************************************************************************************
2040 838562 : SUBROUTINE ensure_work_matrix(work, nrow, ncol)
2041 :
2042 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :), &
2043 : INTENT(INOUT) :: work
2044 : INTEGER, INTENT(IN) :: nrow, ncol
2045 :
2046 838562 : IF (ALLOCATED(work)) THEN
2047 813437 : IF (SIZE(work, 1) == nrow .AND. SIZE(work, 2) == ncol) RETURN
2048 304 : DEALLOCATE (work)
2049 : END IF
2050 101716 : ALLOCATE (work(nrow, ncol))
2051 :
2052 : END SUBROUTINE ensure_work_matrix
2053 :
2054 : ! **************************************************************************************************
2055 : !> \brief Symmetrize a complex k-point matrix including Bloch phase shifts
2056 : !> \param srpmat real part of transformed matrix
2057 : !> \param scpmat imaginary part of transformed matrix
2058 : !> \param rpmat real part of reference matrix
2059 : !> \param cpmat imaginary part of reference matrix
2060 : !> \param real_only ...
2061 : !> \param kmat kind type rotation matrix
2062 : !> \param rot rotation matrix
2063 : !> \param f0 atom permutation
2064 : !> \param fcell atom cell shifts generated by the symmetry operation
2065 : !> \param atype atom to kind pointer
2066 : !> \param xkp target k-point coordinates
2067 : !> \param time_reversal ...
2068 : !> \param reverse_phase ...
2069 : ! **************************************************************************************************
2070 29468 : SUBROUTINE symtrans_phase(srpmat, scpmat, rpmat, cpmat, real_only, kmat, rot, f0, fcell, atype, &
2071 : xkp, time_reversal, reverse_phase)
2072 :
2073 : TYPE(dbcsr_type), POINTER :: srpmat, scpmat, rpmat, cpmat
2074 : LOGICAL, INTENT(IN) :: real_only
2075 : TYPE(kind_rotmat_type), DIMENSION(:), POINTER :: kmat
2076 : REAL(KIND=dp), DIMENSION(3, 3), INTENT(IN) :: rot
2077 : INTEGER, DIMENSION(:), INTENT(IN) :: f0
2078 : INTEGER, DIMENSION(:, :), INTENT(IN) :: fcell
2079 : INTEGER, DIMENSION(:), INTENT(IN) :: atype
2080 : REAL(KIND=dp), DIMENSION(3), INTENT(IN) :: xkp
2081 : LOGICAL, INTENT(IN) :: time_reversal, reverse_phase
2082 :
2083 : CHARACTER(LEN=*), PARAMETER :: routineN = 'symtrans_phase'
2084 :
2085 : INTEGER :: handle, iatom, icol, ikind, ip, irow, &
2086 : jcol, jkind, jp, jrow, mynode, natom, &
2087 : numnodes, owner
2088 : INTEGER, DIMENSION(3) :: shift
2089 : LOGICAL :: dorot, found, has_phase, perm, trans
2090 : REAL(KIND=dp) :: arg, coskl, dr, sinkl
2091 29468 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: cwork, rwork, twork
2092 29468 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: cblock, kroti, krotj, rblock, scblock, &
2093 29468 : srblock
2094 : TYPE(dbcsr_distribution_type) :: dist
2095 : TYPE(dbcsr_iterator_type) :: iter
2096 :
2097 29468 : CALL timeset(routineN, handle)
2098 :
2099 29468 : natom = SIZE(f0)
2100 29468 : perm = .FALSE.
2101 127288 : DO iatom = 1, natom
2102 118304 : IF (f0(iatom) == iatom) CYCLE
2103 : perm = .TRUE.
2104 106804 : EXIT
2105 : END DO
2106 :
2107 29468 : dorot = .FALSE.
2108 383084 : IF (ABS(SUM(ABS(rot)) - 3.0_dp) > 1.e-12_dp) dorot = .TRUE.
2109 29468 : dr = ABS(rot(1, 1) - 1.0_dp) + ABS(rot(2, 2) - 1.0_dp) + ABS(rot(3, 3) - 1.0_dp)
2110 29468 : IF (ABS(dr) > 1.e-12_dp) dorot = .TRUE.
2111 439848 : has_phase = ANY(fcell /= 0) .OR. time_reversal
2112 :
2113 29468 : IF (.NOT. (dorot .OR. perm .OR. has_phase)) THEN
2114 4154 : CALL dbcsr_copy(srpmat, rpmat)
2115 4154 : IF (.NOT. real_only) CALL dbcsr_copy(scpmat, cpmat)
2116 4154 : CALL timestop(handle)
2117 4154 : RETURN
2118 : END IF
2119 :
2120 25314 : CALL dbcsr_get_info(rpmat, distribution=dist)
2121 25314 : CALL dbcsr_distribution_get(dist, mynode=mynode, numnodes=numnodes)
2122 25314 : IF (numnodes /= 1 .AND. (perm .OR. has_phase)) THEN
2123 25270 : CALL dbcsr_replicate_all(rpmat)
2124 25270 : IF (.NOT. real_only) CALL dbcsr_replicate_all(cpmat)
2125 : END IF
2126 :
2127 25314 : CALL dbcsr_set(srpmat, 0.0_dp)
2128 25314 : IF (.NOT. real_only) CALL dbcsr_set(scpmat, 0.0_dp)
2129 :
2130 25314 : CALL dbcsr_iterator_start(iter, rpmat)
2131 863810 : DO WHILE (dbcsr_iterator_blocks_left(iter))
2132 838496 : CALL dbcsr_iterator_next_block(iter, irow, icol, rblock)
2133 838496 : IF (.NOT. ALLOCATED(rwork)) THEN
2134 101256 : ALLOCATE (rwork(SIZE(rblock, 1), SIZE(rblock, 2)))
2135 813182 : ELSEIF (SIZE(rwork, 1) /= SIZE(rblock, 1) .OR. SIZE(rwork, 2) /= SIZE(rblock, 2)) THEN
2136 608 : DEALLOCATE (rwork)
2137 2432 : ALLOCATE (rwork(SIZE(rblock, 1), SIZE(rblock, 2)))
2138 : END IF
2139 838496 : IF (.NOT. real_only) THEN
2140 838496 : IF (.NOT. ALLOCATED(cwork)) THEN
2141 101256 : ALLOCATE (cwork(SIZE(rblock, 1), SIZE(rblock, 2)))
2142 813182 : ELSEIF (SIZE(cwork, 1) /= SIZE(rblock, 1) .OR. SIZE(cwork, 2) /= SIZE(rblock, 2)) THEN
2143 608 : DEALLOCATE (cwork)
2144 2432 : ALLOCATE (cwork(SIZE(rblock, 1), SIZE(rblock, 2)))
2145 : END IF
2146 : END IF
2147 :
2148 838496 : ikind = atype(irow)
2149 838496 : jkind = atype(icol)
2150 838496 : kroti => kmat(ikind)%rmat
2151 838496 : krotj => kmat(jkind)%rmat
2152 :
2153 838496 : IF (reverse_phase) THEN
2154 0 : shift = fcell(1:3, irow) - fcell(1:3, icol)
2155 : ELSE
2156 3353984 : shift = fcell(1:3, icol) - fcell(1:3, irow)
2157 : END IF
2158 838496 : arg = REAL(shift(1), dp)*xkp(1) + REAL(shift(2), dp)*xkp(2) + REAL(shift(3), dp)*xkp(3)
2159 : ! The Bloch phase depends only on the mapped atom-cell shifts and the target k-point.
2160 838496 : coskl = COS(twopi*arg)
2161 838496 : sinkl = SIN(twopi*arg)
2162 838496 : IF (real_only) THEN
2163 0 : IF (ABS(sinkl) > 1.e-12_dp) THEN
2164 0 : CALL cp_abort(__LOCATION__, "Real k-point wavefunctions cannot represent symmetry phases")
2165 : END IF
2166 0 : rwork(:, :) = coskl*rblock
2167 : ELSE
2168 838496 : CALL dbcsr_get_block_p(matrix=cpmat, row=irow, col=icol, block=cblock, found=found)
2169 67348608 : rwork(:, :) = coskl*rblock
2170 838496 : IF (time_reversal) THEN
2171 37421772 : cwork(:, :) = -sinkl*rblock
2172 470660 : IF (found) THEN
2173 37421772 : rwork(:, :) = rwork - sinkl*cblock
2174 37421772 : cwork(:, :) = cwork - coskl*cblock
2175 : END IF
2176 : ELSE
2177 29926836 : cwork(:, :) = -sinkl*rblock
2178 367836 : IF (found) THEN
2179 29926836 : rwork(:, :) = rwork + sinkl*cblock
2180 29926836 : cwork(:, :) = cwork + coskl*cblock
2181 : END IF
2182 : END IF
2183 : END IF
2184 :
2185 838496 : ip = f0(irow)
2186 838496 : jp = f0(icol)
2187 838496 : IF (ip <= jp) THEN
2188 745920 : jrow = ip
2189 745920 : jcol = jp
2190 745920 : trans = .FALSE.
2191 : ELSE
2192 92576 : jrow = jp
2193 92576 : jcol = ip
2194 92576 : trans = .TRUE.
2195 : END IF
2196 :
2197 838496 : CALL dbcsr_get_block_p(matrix=srpmat, row=jrow, col=jcol, block=srblock, found=found)
2198 838496 : IF (.NOT. found) THEN
2199 419215 : CALL dbcsr_get_stored_coordinates(srpmat, jrow, jcol, owner)
2200 419215 : CPASSERT(owner /= mynode)
2201 : CYCLE
2202 : END IF
2203 419281 : IF (trans) THEN
2204 46288 : CALL ensure_work_matrix(twork, SIZE(krotj, 1), SIZE(rwork, 1))
2205 : CALL dgemm('N', 'T', SIZE(krotj, 1), SIZE(rwork, 1), SIZE(krotj, 2), &
2206 : 1.0_dp, krotj, SIZE(krotj, 1), rwork, SIZE(rwork, 1), &
2207 46288 : 0.0_dp, twork, SIZE(twork, 1))
2208 : CALL dgemm('N', 'T', SIZE(twork, 1), SIZE(kroti, 1), SIZE(twork, 2), &
2209 : 1.0_dp, twork, SIZE(twork, 1), kroti, SIZE(kroti, 1), &
2210 46288 : 1.0_dp, srblock, SIZE(srblock, 1))
2211 : ELSE
2212 372993 : CALL ensure_work_matrix(twork, SIZE(kroti, 1), SIZE(rwork, 2))
2213 : CALL dgemm('N', 'N', SIZE(kroti, 1), SIZE(rwork, 2), SIZE(kroti, 2), &
2214 : 1.0_dp, kroti, SIZE(kroti, 1), rwork, SIZE(rwork, 1), &
2215 372993 : 0.0_dp, twork, SIZE(twork, 1))
2216 : CALL dgemm('N', 'T', SIZE(twork, 1), SIZE(krotj, 1), SIZE(twork, 2), &
2217 : 1.0_dp, twork, SIZE(twork, 1), krotj, SIZE(krotj, 1), &
2218 372993 : 1.0_dp, srblock, SIZE(srblock, 1))
2219 : END IF
2220 :
2221 444595 : IF (.NOT. real_only) THEN
2222 419281 : CALL dbcsr_get_block_p(matrix=scpmat, row=jrow, col=jcol, block=scblock, found=found)
2223 419281 : CPASSERT(found)
2224 419281 : IF (trans) THEN
2225 46288 : CALL ensure_work_matrix(twork, SIZE(krotj, 1), SIZE(cwork, 1))
2226 : CALL dgemm('N', 'T', SIZE(krotj, 1), SIZE(cwork, 1), SIZE(krotj, 2), &
2227 : 1.0_dp, krotj, SIZE(krotj, 1), cwork, SIZE(cwork, 1), &
2228 46288 : 0.0_dp, twork, SIZE(twork, 1))
2229 : CALL dgemm('N', 'T', SIZE(twork, 1), SIZE(kroti, 1), SIZE(twork, 2), &
2230 : -1.0_dp, twork, SIZE(twork, 1), kroti, SIZE(kroti, 1), &
2231 46288 : 1.0_dp, scblock, SIZE(scblock, 1))
2232 : ELSE
2233 372993 : CALL ensure_work_matrix(twork, SIZE(kroti, 1), SIZE(cwork, 2))
2234 : CALL dgemm('N', 'N', SIZE(kroti, 1), SIZE(cwork, 2), SIZE(kroti, 2), &
2235 : 1.0_dp, kroti, SIZE(kroti, 1), cwork, SIZE(cwork, 1), &
2236 372993 : 0.0_dp, twork, SIZE(twork, 1))
2237 : CALL dgemm('N', 'T', SIZE(twork, 1), SIZE(krotj, 1), SIZE(twork, 2), &
2238 : 1.0_dp, twork, SIZE(twork, 1), krotj, SIZE(krotj, 1), &
2239 372993 : 1.0_dp, scblock, SIZE(scblock, 1))
2240 : END IF
2241 : END IF
2242 : END DO
2243 25314 : CALL dbcsr_iterator_stop(iter)
2244 25314 : IF (numnodes /= 1 .AND. (perm .OR. has_phase)) THEN
2245 25270 : CALL dbcsr_distribute(rpmat)
2246 25270 : IF (.NOT. real_only) CALL dbcsr_distribute(cpmat)
2247 : END IF
2248 :
2249 25314 : CALL timestop(handle)
2250 :
2251 58936 : END SUBROUTINE symtrans_phase
2252 :
2253 : ! **************************************************************************************************
2254 : !> \brief Symmetrization of density matrix - transform to new k-point
2255 : !> \param smat density matrix at new kpoint
2256 : !> \param pmat reference density matrix
2257 : !> \param kmat Kind type rotation matrix
2258 : !> \param rot Rotation matrix
2259 : !> \param f0 Permutation of atoms under transformation
2260 : !> \param atype Atom to kind pointer
2261 : !> \param symmetric Symmetric matrix
2262 : !> \param antisymmetric Anti-Symmetric matrix
2263 : ! **************************************************************************************************
2264 0 : SUBROUTINE symtrans(smat, pmat, kmat, rot, f0, atype, symmetric, antisymmetric)
2265 : TYPE(dbcsr_type), POINTER :: smat, pmat
2266 : TYPE(kind_rotmat_type), DIMENSION(:), POINTER :: kmat
2267 : REAL(KIND=dp), DIMENSION(3, 3), INTENT(IN) :: rot
2268 : INTEGER, DIMENSION(:), INTENT(IN) :: f0, atype
2269 : LOGICAL, INTENT(IN), OPTIONAL :: symmetric, antisymmetric
2270 :
2271 : CHARACTER(LEN=*), PARAMETER :: routineN = 'symtrans'
2272 :
2273 : INTEGER :: handle, iatom, icol, ikind, ip, irow, &
2274 : jcol, jkind, jp, jrow, natom, numnodes
2275 : LOGICAL :: asym, dorot, found, perm, sym, trans
2276 : REAL(KIND=dp) :: dr, fsign
2277 0 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: work
2278 0 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: kroti, krotj, pblock, sblock
2279 : TYPE(dbcsr_distribution_type) :: dist
2280 : TYPE(dbcsr_iterator_type) :: iter
2281 :
2282 0 : CALL timeset(routineN, handle)
2283 :
2284 : ! check symmetry options
2285 0 : sym = .FALSE.
2286 0 : IF (PRESENT(symmetric)) sym = symmetric
2287 0 : asym = .FALSE.
2288 0 : IF (PRESENT(antisymmetric)) asym = antisymmetric
2289 :
2290 0 : CPASSERT(.NOT. (sym .AND. asym))
2291 0 : CPASSERT((sym .OR. asym))
2292 :
2293 : ! do we have permutation of atoms
2294 0 : natom = SIZE(f0)
2295 0 : perm = .FALSE.
2296 0 : DO iatom = 1, natom
2297 0 : IF (f0(iatom) == iatom) CYCLE
2298 : perm = .TRUE.
2299 0 : EXIT
2300 : END DO
2301 :
2302 : ! do we have a real rotation
2303 0 : dorot = .FALSE.
2304 0 : IF (ABS(SUM(ABS(rot)) - 3.0_dp) > 1.e-12_dp) dorot = .TRUE.
2305 0 : dr = ABS(rot(1, 1) - 1.0_dp) + ABS(rot(2, 2) - 1.0_dp) + ABS(rot(3, 3) - 1.0_dp)
2306 0 : IF (ABS(dr) > 1.e-12_dp) dorot = .TRUE.
2307 :
2308 0 : fsign = 1.0_dp
2309 0 : IF (asym) fsign = -1.0_dp
2310 :
2311 0 : IF (dorot .OR. perm) THEN
2312 : CALL cp_abort(__LOCATION__, "k-points need FULL_GRID currently. "// &
2313 0 : "Reduced grids not yet working correctly")
2314 0 : CALL dbcsr_set(smat, 0.0_dp)
2315 0 : IF (perm) THEN
2316 0 : CALL dbcsr_get_info(pmat, distribution=dist)
2317 0 : CALL dbcsr_distribution_get(dist, numnodes=numnodes)
2318 0 : IF (numnodes == 1) THEN
2319 : ! the matrices are local to this process
2320 0 : CALL dbcsr_iterator_start(iter, pmat)
2321 0 : DO WHILE (dbcsr_iterator_blocks_left(iter))
2322 0 : CALL dbcsr_iterator_next_block(iter, irow, icol, pblock)
2323 0 : ip = f0(irow)
2324 0 : jp = f0(icol)
2325 0 : IF (ip <= jp) THEN
2326 0 : jrow = ip
2327 0 : jcol = jp
2328 0 : trans = .FALSE.
2329 : ELSE
2330 0 : jrow = jp
2331 0 : jcol = ip
2332 0 : trans = .TRUE.
2333 : END IF
2334 0 : CALL dbcsr_get_block_p(matrix=smat, row=jrow, col=jcol, BLOCK=sblock, found=found)
2335 0 : CPASSERT(found)
2336 0 : ikind = atype(irow)
2337 0 : jkind = atype(icol)
2338 0 : kroti => kmat(ikind)%rmat
2339 0 : krotj => kmat(jkind)%rmat
2340 : ! rotation
2341 0 : IF (trans) THEN
2342 0 : CALL ensure_work_matrix(work, SIZE(krotj, 2), SIZE(pblock, 1))
2343 : CALL dgemm('T', 'T', SIZE(krotj, 2), SIZE(pblock, 1), SIZE(krotj, 1), &
2344 : 1.0_dp, krotj, SIZE(krotj, 1), pblock, SIZE(pblock, 1), &
2345 0 : 0.0_dp, work, SIZE(work, 1))
2346 : CALL dgemm('N', 'N', SIZE(work, 1), SIZE(kroti, 2), SIZE(work, 2), &
2347 : fsign, work, SIZE(work, 1), kroti, SIZE(kroti, 1), &
2348 0 : 0.0_dp, sblock, SIZE(sblock, 1))
2349 : ELSE
2350 0 : CALL ensure_work_matrix(work, SIZE(kroti, 2), SIZE(pblock, 2))
2351 : CALL dgemm('T', 'N', SIZE(kroti, 2), SIZE(pblock, 2), SIZE(kroti, 1), &
2352 : 1.0_dp, kroti, SIZE(kroti, 1), pblock, SIZE(pblock, 1), &
2353 0 : 0.0_dp, work, SIZE(work, 1))
2354 : CALL dgemm('N', 'N', SIZE(work, 1), SIZE(krotj, 2), SIZE(work, 2), &
2355 : fsign, work, SIZE(work, 1), krotj, SIZE(krotj, 1), &
2356 0 : 0.0_dp, sblock, SIZE(sblock, 1))
2357 : END IF
2358 : END DO
2359 0 : CALL dbcsr_iterator_stop(iter)
2360 : !
2361 : ELSE
2362 : ! distributed matrices, most general code needed
2363 : CALL cp_abort(__LOCATION__, "k-points need FULL_GRID currently. "// &
2364 0 : "Reduced grids not yet working correctly")
2365 : END IF
2366 : ELSE
2367 : ! no atom permutations, this is always local
2368 0 : CALL dbcsr_copy(smat, pmat)
2369 0 : CALL dbcsr_iterator_start(iter, smat)
2370 0 : DO WHILE (dbcsr_iterator_blocks_left(iter))
2371 0 : CALL dbcsr_iterator_next_block(iter, irow, icol, sblock)
2372 0 : ip = f0(irow)
2373 0 : jp = f0(icol)
2374 0 : IF (ip <= jp) THEN
2375 : jrow = ip
2376 : jcol = jp
2377 0 : trans = .FALSE.
2378 : ELSE
2379 : jrow = jp
2380 : jcol = ip
2381 0 : trans = .TRUE.
2382 : END IF
2383 0 : ikind = atype(irow)
2384 0 : jkind = atype(icol)
2385 0 : kroti => kmat(ikind)%rmat
2386 0 : krotj => kmat(jkind)%rmat
2387 : ! rotation
2388 0 : IF (trans) THEN
2389 0 : CALL ensure_work_matrix(work, SIZE(krotj, 2), SIZE(sblock, 1))
2390 : CALL dgemm('T', 'T', SIZE(krotj, 2), SIZE(sblock, 1), SIZE(krotj, 1), &
2391 : 1.0_dp, krotj, SIZE(krotj, 1), sblock, SIZE(sblock, 1), &
2392 0 : 0.0_dp, work, SIZE(work, 1))
2393 : CALL dgemm('N', 'N', SIZE(work, 1), SIZE(kroti, 2), SIZE(work, 2), &
2394 : fsign, work, SIZE(work, 1), kroti, SIZE(kroti, 1), &
2395 0 : 0.0_dp, sblock, SIZE(sblock, 1))
2396 : ELSE
2397 0 : CALL ensure_work_matrix(work, SIZE(kroti, 2), SIZE(sblock, 2))
2398 : CALL dgemm('T', 'N', SIZE(kroti, 2), SIZE(sblock, 2), SIZE(kroti, 1), &
2399 : 1.0_dp, kroti, SIZE(kroti, 1), sblock, SIZE(sblock, 1), &
2400 0 : 0.0_dp, work, SIZE(work, 1))
2401 : CALL dgemm('N', 'N', SIZE(work, 1), SIZE(krotj, 2), SIZE(work, 2), &
2402 : fsign, work, SIZE(work, 1), krotj, SIZE(krotj, 1), &
2403 0 : 0.0_dp, sblock, SIZE(sblock, 1))
2404 : END IF
2405 : END DO
2406 0 : CALL dbcsr_iterator_stop(iter)
2407 : !
2408 : END IF
2409 : ELSE
2410 : ! this is the identity operation, just copy the matrix
2411 0 : CALL dbcsr_copy(smat, pmat)
2412 : END IF
2413 :
2414 0 : CALL timestop(handle)
2415 :
2416 0 : END SUBROUTINE symtrans
2417 :
2418 : ! **************************************************************************************************
2419 : !> \brief ...
2420 : !> \param mat ...
2421 : ! **************************************************************************************************
2422 0 : SUBROUTINE matprint(mat)
2423 : TYPE(dbcsr_type), POINTER :: mat
2424 :
2425 : INTEGER :: i, icol, iounit, irow
2426 0 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: mblock
2427 : TYPE(dbcsr_iterator_type) :: iter
2428 :
2429 0 : iounit = cp_logger_get_default_io_unit()
2430 0 : CALL dbcsr_iterator_start(iter, mat)
2431 0 : DO WHILE (dbcsr_iterator_blocks_left(iter))
2432 0 : CALL dbcsr_iterator_next_block(iter, irow, icol, mblock)
2433 : !
2434 0 : IF (iounit > 0) THEN
2435 0 : WRITE (iounit, '(A,2I4)') 'BLOCK ', irow, icol
2436 0 : DO i = 1, SIZE(mblock, 1)
2437 0 : WRITE (iounit, '(8F12.6)') mblock(i, :)
2438 : END DO
2439 : END IF
2440 : !
2441 : END DO
2442 0 : CALL dbcsr_iterator_stop(iter)
2443 :
2444 0 : END SUBROUTINE matprint
2445 : ! **************************************************************************************************
2446 :
2447 : END MODULE kpoint_methods
|