(git:7a2a232)
Loading...
Searching...
No Matches
kpoint_methods.F
Go to the documentation of this file.
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! **************************************************************************************************
17 USE cell_types, ONLY: cell_type,&
21 USE cp_cfm_types, ONLY: cp_cfm_create,&
27 USE cp_dbcsr_api, ONLY: &
31 dbcsr_p_type, dbcsr_set, dbcsr_type, dbcsr_type_antisymmetric, dbcsr_type_no_symmetry, &
32 dbcsr_type_symmetric
40 USE cp_fm_types, ONLY: &
45 USE cryssym, ONLY: crys_sym_gen,&
46 csym_type,&
54 smear_mp,&
56 USE kinds, ONLY: dp
57 USE kpoint_types, ONLY: get_kpoint_info,&
65 USE mathconstants, ONLY: twopi
67 USE message_passing, ONLY: mp_cart_type,&
75 USE qs_mo_types, ONLY: allocate_mo_set,&
88 USE smearing_utils, ONLY: smearkp,&
90 USE util, ONLY: get_limit
91#include "./base/base_uses.f90"
92
93 IMPLICIT NONE
94
95 PRIVATE
96
97 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'kpoint_methods'
98
103
104! **************************************************************************************************
105
106CONTAINS
107
108! **************************************************************************************************
109!> \brief Generate the kpoints and initialize the kpoint environment
110!> \param kpoint The kpoint environment
111!> \param particle_set Particle types and coordinates
112!> \param cell Computational cell information
113! **************************************************************************************************
114 SUBROUTINE kpoint_initialize(kpoint, particle_set, cell)
115
116 TYPE(kpoint_type), POINTER :: kpoint
117 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
118 TYPE(cell_type), POINTER :: cell
119
120 CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_initialize'
121
122 INTEGER :: handle, i, ic, ik, iounit, ir, ira, is, &
123 j, natom, nkind, nr, ns
124 INTEGER, ALLOCATABLE, DIMENSION(:) :: atype
125 LOGICAL :: spez
126 REAL(kind=dp) :: wsum
127 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: coord, scoord
128 TYPE(csym_type) :: crys_sym
129 TYPE(kpoint_sym_type), POINTER :: kpsym
130
131 CALL timeset(routinen, handle)
132
133 cpassert(ASSOCIATED(kpoint))
134
135 SELECT CASE (kpoint%kp_scheme)
136 CASE ("NONE")
137 ! do nothing
138 CASE ("GAMMA")
139 kpoint%nkp = 1
140 ALLOCATE (kpoint%xkp(3, 1), kpoint%wkp(1))
141 kpoint%xkp(1:3, 1) = 0.0_dp
142 kpoint%wkp(1) = 1.0_dp
143 ALLOCATE (kpoint%kp_sym(1))
144 NULLIFY (kpoint%kp_sym(1)%kpoint_sym)
145 CALL kpoint_sym_create(kpoint%kp_sym(1)%kpoint_sym)
146 CASE ("MONKHORST-PACK", "MACDONALD")
147
148 IF (.NOT. kpoint%symmetry) THEN
149 ! we set up a random molecule to avoid any possible symmetry
150 natom = 10
151 ALLOCATE (coord(3, natom), scoord(3, natom), atype(natom))
152 DO i = 1, natom
153 atype(i) = i
154 coord(1, i) = sin(i*0.12345_dp)
155 coord(2, i) = cos(i*0.23456_dp)
156 coord(3, i) = sin(i*0.34567_dp)
157 CALL real_to_scaled(scoord(1:3, i), coord(1:3, i), cell)
158 END DO
159 ELSE
160 natom = SIZE(particle_set)
161 ALLOCATE (scoord(3, natom), atype(natom))
162 DO i = 1, natom
163 CALL get_atomic_kind(atomic_kind=particle_set(i)%atomic_kind, kind_number=atype(i))
164 CALL real_to_scaled(scoord(1:3, i), particle_set(i)%r(1:3), cell)
165 END DO
166 END IF
167 IF (kpoint%verbose) THEN
169 ELSE
170 iounit = -1
171 END IF
172 ! kind type list
173 ALLOCATE (kpoint%atype(natom))
174 kpoint%atype = atype
175
176 CALL crys_sym_gen(crys_sym, scoord, atype, cell%hmat, delta=kpoint%eps_geo, iounit=iounit)
177 CALL kpoint_gen(crys_sym, kpoint%nkp_grid, symm=kpoint%symmetry, shift=kpoint%kp_shift, &
178 full_grid=kpoint%full_grid)
179 kpoint%nkp = crys_sym%nkpoint
180 ALLOCATE (kpoint%xkp(3, kpoint%nkp), kpoint%wkp(kpoint%nkp))
181 wsum = sum(crys_sym%wkpoint)
182 DO ik = 1, kpoint%nkp
183 kpoint%xkp(1:3, ik) = crys_sym%xkpoint(1:3, ik)
184 kpoint%wkp(ik) = crys_sym%wkpoint(ik)/wsum
185 END DO
186
187 ! print output
188 IF (kpoint%symmetry) CALL print_crys_symmetry(crys_sym)
189 IF (kpoint%symmetry) CALL print_kp_symmetry(crys_sym)
190
191 ! transfer symmetry information
192 ALLOCATE (kpoint%kp_sym(kpoint%nkp))
193 DO ik = 1, kpoint%nkp
194 NULLIFY (kpoint%kp_sym(ik)%kpoint_sym)
195 CALL kpoint_sym_create(kpoint%kp_sym(ik)%kpoint_sym)
196 kpsym => kpoint%kp_sym(ik)%kpoint_sym
197 IF (crys_sym%symlib .AND. .NOT. crys_sym%fullgrid .AND. crys_sym%istriz == 1) THEN
198 ! set up the symmetrization information
199 kpsym%nwght = nint(crys_sym%wkpoint(ik))
200 ns = kpsym%nwght
201 !
202 IF (ns > 1) THEN
203 kpsym%apply_symmetry = .true.
204 natom = SIZE(particle_set)
205 ALLOCATE (kpsym%rot(3, 3, ns))
206 ALLOCATE (kpsym%xkp(3, ns))
207 ALLOCATE (kpsym%rotp(ns))
208 ALLOCATE (kpsym%f0(natom, ns))
209 nr = 0
210 DO is = 1, SIZE(crys_sym%kplink, 2)
211 IF (crys_sym%kplink(2, is) == ik) THEN
212 nr = nr + 1
213 ir = crys_sym%kpop(is)
214 ira = abs(ir)
215 kpsym%rotp(nr) = ir
216 IF (ir > 0) THEN
217 kpsym%rot(1:3, 1:3, nr) = crys_sym%rt(1:3, 1:3, ira)
218 ELSE
219 kpsym%rot(1:3, 1:3, nr) = -crys_sym%rt(1:3, 1:3, ira)
220 END IF
221 kpsym%xkp(1:3, nr) = crys_sym%kpmesh(1:3, is)
222 DO ic = 1, crys_sym%nrtot
223 IF (crys_sym%ibrot(ic) == ira) THEN
224 kpsym%f0(1:natom, nr) = crys_sym%f0(1:natom, ic)
225 EXIT
226 END IF
227 END DO
228 END IF
229 END DO
230 ! Reduce inversion?
231 kpsym%nwred = nr
232 END IF
233 END IF
234 END DO
235 IF (kpoint%symmetry) THEN
236 nkind = maxval(atype)
237 ns = crys_sym%nrtot
238 ALLOCATE (kpoint%kind_rotmat(ns, nkind))
239 DO i = 1, ns
240 DO j = 1, nkind
241 NULLIFY (kpoint%kind_rotmat(i, j)%rmat)
242 END DO
243 END DO
244 ALLOCATE (kpoint%ibrot(ns))
245 kpoint%ibrot(1:ns) = crys_sym%ibrot(1:ns)
246 END IF
247
248 CALL release_csym_type(crys_sym)
249 DEALLOCATE (scoord, atype)
250
251 CASE ("GENERAL")
252 ! default: no symmetry settings
253 ALLOCATE (kpoint%kp_sym(kpoint%nkp))
254 DO i = 1, kpoint%nkp
255 NULLIFY (kpoint%kp_sym(i)%kpoint_sym)
256 CALL kpoint_sym_create(kpoint%kp_sym(i)%kpoint_sym)
257 END DO
258 CASE DEFAULT
259 cpassert(.false.)
260 END SELECT
261
262 ! check for consistency of options
263 SELECT CASE (kpoint%kp_scheme)
264 CASE ("NONE")
265 ! don't use k-point code
266 CASE ("GAMMA")
267 cpassert(kpoint%nkp == 1)
268 cpassert(sum(abs(kpoint%xkp)) <= 1.e-12_dp)
269 cpassert(kpoint%wkp(1) == 1.0_dp)
270 cpassert(.NOT. kpoint%symmetry)
271 CASE ("GENERAL")
272 cpassert(.NOT. kpoint%symmetry)
273 cpassert(kpoint%nkp >= 1)
274 CASE ("MONKHORST-PACK", "MACDONALD")
275 cpassert(kpoint%nkp >= 1)
276 END SELECT
277 IF (kpoint%use_real_wfn) THEN
278 ! what about inversion symmetry?
279 ikloop: DO ik = 1, kpoint%nkp
280 DO i = 1, 3
281 spez = (kpoint%xkp(i, ik) == 0.0_dp .OR. kpoint%xkp(i, ik) == 0.5_dp)
282 IF (.NOT. spez) EXIT ikloop
283 END DO
284 END DO ikloop
285 IF (.NOT. spez) THEN
286 ! Warning: real wfn might be wrong for this system
287 CALL cp_warn(__location__, &
288 "A calculation using real wavefunctions is requested. "// &
289 "We could not determine if the symmetry of the system allows real wavefunctions. ")
290 END IF
291 END IF
292
293 CALL timestop(handle)
294
295 END SUBROUTINE kpoint_initialize
296
297! **************************************************************************************************
298!> \brief Initialize the kpoint environment
299!> \param kpoint Kpoint environment
300!> \param para_env ...
301!> \param blacs_env ...
302!> \param with_aux_fit ...
303! **************************************************************************************************
304 SUBROUTINE kpoint_env_initialize(kpoint, para_env, blacs_env, with_aux_fit)
305
306 TYPE(kpoint_type), INTENT(INOUT) :: kpoint
307 TYPE(mp_para_env_type), INTENT(IN), TARGET :: para_env
308 TYPE(cp_blacs_env_type), INTENT(IN), TARGET :: blacs_env
309 LOGICAL, INTENT(IN), OPTIONAL :: with_aux_fit
310
311 CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_env_initialize'
312
313 INTEGER :: handle, igr, ik, ikk, ngr, niogrp, nkp, &
314 nkp_grp, nkp_loc, npe, unit_nr
315 INTEGER, DIMENSION(2) :: dims, pos
316 LOGICAL :: aux_fit
317 TYPE(kpoint_env_p_type), DIMENSION(:), POINTER :: kp_aux_env, kp_env
318 TYPE(kpoint_env_type), POINTER :: kp
319 TYPE(mp_cart_type) :: comm_cart
320 TYPE(mp_para_env_type), POINTER :: para_env_inter_kp, para_env_kp
321
322 CALL timeset(routinen, handle)
323
324 IF (PRESENT(with_aux_fit)) THEN
325 aux_fit = with_aux_fit
326 ELSE
327 aux_fit = .false.
328 END IF
329
330 kpoint%para_env => para_env
331 CALL kpoint%para_env%retain()
332 kpoint%blacs_env_all => blacs_env
333 CALL kpoint%blacs_env_all%retain()
334
335 cpassert(.NOT. ASSOCIATED(kpoint%kp_env))
336 IF (aux_fit) THEN
337 cpassert(.NOT. ASSOCIATED(kpoint%kp_aux_env))
338 END IF
339
340 NULLIFY (kp_env, kp_aux_env)
341 nkp = kpoint%nkp
342 npe = para_env%num_pe
343 IF (npe == 1) THEN
344 ! only one process available -> owns all kpoints
345 ALLOCATE (kp_env(nkp))
346 DO ik = 1, nkp
347 NULLIFY (kp_env(ik)%kpoint_env)
348 CALL kpoint_env_create(kp_env(ik)%kpoint_env)
349 kp => kp_env(ik)%kpoint_env
350 kp%nkpoint = ik
351 kp%wkp = kpoint%wkp(ik)
352 kp%xkp(1:3) = kpoint%xkp(1:3, ik)
353 kp%is_local = .true.
354 END DO
355 kpoint%kp_env => kp_env
356
357 IF (aux_fit) THEN
358 ALLOCATE (kp_aux_env(nkp))
359 DO ik = 1, nkp
360 NULLIFY (kp_aux_env(ik)%kpoint_env)
361 CALL kpoint_env_create(kp_aux_env(ik)%kpoint_env)
362 kp => kp_aux_env(ik)%kpoint_env
363 kp%nkpoint = ik
364 kp%wkp = kpoint%wkp(ik)
365 kp%xkp(1:3) = kpoint%xkp(1:3, ik)
366 kp%is_local = .true.
367 END DO
368
369 kpoint%kp_aux_env => kp_aux_env
370 END IF
371
372 ALLOCATE (kpoint%kp_dist(2, 1))
373 kpoint%kp_dist(1, 1) = 1
374 kpoint%kp_dist(2, 1) = nkp
375 kpoint%kp_range(1) = 1
376 kpoint%kp_range(2) = nkp
377
378 ! parallel environments
379 kpoint%para_env_kp => para_env
380 CALL kpoint%para_env_kp%retain()
381 kpoint%para_env_inter_kp => para_env
382 CALL kpoint%para_env_inter_kp%retain()
383 kpoint%iogrp = .true.
384 kpoint%nkp_groups = 1
385 ELSE
386 IF (kpoint%parallel_group_size == -1) THEN
387 ! maximum parallelization over kpoints
388 ! making sure that the group size divides the npe and the nkp_grp the nkp
389 ! in the worst case, there will be no parallelism over kpoints.
390 DO igr = npe, 1, -1
391 IF (mod(npe, igr) /= 0) cycle
392 nkp_grp = npe/igr
393 IF (mod(nkp, nkp_grp) /= 0) cycle
394 ngr = igr
395 END DO
396 ELSE IF (kpoint%parallel_group_size == 0) THEN
397 ! no parallelization over kpoints
398 ngr = npe
399 ELSE IF (kpoint%parallel_group_size > 0) THEN
400 ngr = min(kpoint%parallel_group_size, npe)
401 ELSE
402 cpassert(.false.)
403 END IF
404 nkp_grp = npe/ngr
405 ! processor dimensions
406 dims(1) = ngr
407 dims(2) = nkp_grp
408 cpassert(mod(nkp, nkp_grp) == 0)
409 nkp_loc = nkp/nkp_grp
410
411 IF ((dims(1)*dims(2) /= npe)) THEN
412 cpabort("Number of processors is not divisible by the kpoint group size.")
413 END IF
414
415 ! Create the subgroups, one for each k-point group and one interconnecting group
416 CALL comm_cart%create(comm_old=para_env, ndims=2, dims=dims)
417 pos = comm_cart%mepos_cart
418 ALLOCATE (para_env_kp)
419 CALL para_env_kp%from_split(comm_cart, pos(2))
420 ALLOCATE (para_env_inter_kp)
421 CALL para_env_inter_kp%from_split(comm_cart, pos(1))
422 CALL comm_cart%free()
423
424 niogrp = 0
425 IF (para_env%is_source()) niogrp = 1
426 CALL para_env_kp%sum(niogrp)
427 kpoint%iogrp = (niogrp == 1)
428
429 ! parallel groups
430 kpoint%para_env_kp => para_env_kp
431 kpoint%para_env_inter_kp => para_env_inter_kp
432
433 ! distribution of kpoints
434 ALLOCATE (kpoint%kp_dist(2, nkp_grp))
435 DO igr = 1, nkp_grp
436 kpoint%kp_dist(1:2, igr) = get_limit(nkp, nkp_grp, igr - 1)
437 END DO
438 ! local kpoints
439 kpoint%kp_range(1:2) = kpoint%kp_dist(1:2, para_env_inter_kp%mepos + 1)
440
441 ALLOCATE (kp_env(nkp_loc))
442 DO ik = 1, nkp_loc
443 NULLIFY (kp_env(ik)%kpoint_env)
444 ikk = kpoint%kp_range(1) + ik - 1
445 CALL kpoint_env_create(kp_env(ik)%kpoint_env)
446 kp => kp_env(ik)%kpoint_env
447 kp%nkpoint = ikk
448 kp%wkp = kpoint%wkp(ikk)
449 kp%xkp(1:3) = kpoint%xkp(1:3, ikk)
450 kp%is_local = (ngr == 1)
451 END DO
452 kpoint%kp_env => kp_env
453
454 IF (aux_fit) THEN
455 ALLOCATE (kp_aux_env(nkp_loc))
456 DO ik = 1, nkp_loc
457 NULLIFY (kp_aux_env(ik)%kpoint_env)
458 ikk = kpoint%kp_range(1) + ik - 1
459 CALL kpoint_env_create(kp_aux_env(ik)%kpoint_env)
460 kp => kp_aux_env(ik)%kpoint_env
461 kp%nkpoint = ikk
462 kp%wkp = kpoint%wkp(ikk)
463 kp%xkp(1:3) = kpoint%xkp(1:3, ikk)
464 kp%is_local = (ngr == 1)
465 END DO
466 kpoint%kp_aux_env => kp_aux_env
467 END IF
468
470
471 IF (unit_nr > 0 .AND. kpoint%verbose) THEN
472 WRITE (unit_nr, *)
473 WRITE (unit_nr, fmt="(T2,A,T71,I10)") "KPOINTS| Number of kpoint groups ", nkp_grp
474 WRITE (unit_nr, fmt="(T2,A,T71,I10)") "KPOINTS| Size of each kpoint group", ngr
475 WRITE (unit_nr, fmt="(T2,A,T71,I10)") "KPOINTS| Number of kpoints per group", nkp_loc
476 END IF
477 kpoint%nkp_groups = nkp_grp
478
479 END IF
480
481 CALL timestop(handle)
482
483 END SUBROUTINE kpoint_env_initialize
484
485! **************************************************************************************************
486!> \brief Initialize a set of MOs and density matrix for each kpoint (kpoint group)
487!> \param kpoint Kpoint environment
488!> \param mos Reference MOs (global)
489!> \param added_mos ...
490!> \param for_aux_fit ...
491! **************************************************************************************************
492 SUBROUTINE kpoint_initialize_mos(kpoint, mos, added_mos, for_aux_fit)
493
494 TYPE(kpoint_type), POINTER :: kpoint
495 TYPE(mo_set_type), DIMENSION(:), INTENT(INOUT) :: mos
496 INTEGER, INTENT(IN), OPTIONAL :: added_mos
497 LOGICAL, OPTIONAL :: for_aux_fit
498
499 CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_initialize_mos'
500
501 INTEGER :: handle, ic, ik, is, nadd, nao, nc, &
502 nelectron, nkp_loc, nmo, nmorig(2), &
503 nspin
504 LOGICAL :: aux_fit
505 REAL(kind=dp) :: flexible_electron_count, maxocc, n_el_f
506 TYPE(cp_blacs_env_type), POINTER :: blacs_env
507 TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_ao_fm_pools
508 TYPE(cp_fm_struct_type), POINTER :: matrix_struct
509 TYPE(cp_fm_type), POINTER :: fmlocal
510 TYPE(kpoint_env_type), POINTER :: kp
511 TYPE(qs_matrix_pools_type), POINTER :: mpools
512
513 CALL timeset(routinen, handle)
514
515 IF (PRESENT(for_aux_fit)) THEN
516 aux_fit = for_aux_fit
517 ELSE
518 aux_fit = .false.
519 END IF
520
521 cpassert(ASSOCIATED(kpoint))
522
523 IF (.true. .OR. ASSOCIATED(mos(1)%mo_coeff)) THEN
524 IF (aux_fit) THEN
525 cpassert(ASSOCIATED(kpoint%kp_aux_env))
526 END IF
527
528 IF (PRESENT(added_mos)) THEN
529 nadd = added_mos
530 ELSE
531 nadd = 0
532 END IF
533
534 IF (kpoint%use_real_wfn) THEN
535 nc = 1
536 ELSE
537 nc = 2
538 END IF
539 nspin = SIZE(mos, 1)
540 nkp_loc = kpoint%kp_range(2) - kpoint%kp_range(1) + 1
541 IF (nkp_loc > 0) THEN
542 IF (aux_fit) THEN
543 cpassert(SIZE(kpoint%kp_aux_env) == nkp_loc)
544 ELSE
545 cpassert(SIZE(kpoint%kp_env) == nkp_loc)
546 END IF
547 ! allocate the mo sets, correct number of kpoints (local), real/complex, spin
548 DO ik = 1, nkp_loc
549 IF (aux_fit) THEN
550 kp => kpoint%kp_aux_env(ik)%kpoint_env
551 ELSE
552 kp => kpoint%kp_env(ik)%kpoint_env
553 END IF
554 ALLOCATE (kp%mos(nc, nspin))
555 DO is = 1, nspin
556 CALL get_mo_set(mos(is), nao=nao, nmo=nmo, nelectron=nelectron, &
557 n_el_f=n_el_f, maxocc=maxocc, flexible_electron_count=flexible_electron_count)
558 nmo = min(nao, nmo + nadd)
559 DO ic = 1, nc
560 CALL allocate_mo_set(kp%mos(ic, is), nao, nmo, nelectron, n_el_f, maxocc, &
561 flexible_electron_count)
562 END DO
563 END DO
564 END DO
565
566 ! generate the blacs environment for the kpoint group
567 ! we generate a blacs env for each kpoint group in parallel
568 ! we assume here that the group para_env_inter_kp will connect
569 ! equivalent parts of fm matrices, i.e. no reshuffeling of processors
570 NULLIFY (blacs_env)
571 IF (ASSOCIATED(kpoint%blacs_env)) THEN
572 blacs_env => kpoint%blacs_env
573 ELSE
574 CALL cp_blacs_env_create(blacs_env=blacs_env, para_env=kpoint%para_env_kp)
575 kpoint%blacs_env => blacs_env
576 END IF
577
578 ! set possible new number of MOs
579 DO is = 1, nspin
580 CALL get_mo_set(mos(is), nmo=nmorig(is))
581 nmo = min(nao, nmorig(is) + nadd)
582 CALL set_mo_set(mos(is), nmo=nmo)
583 END DO
584 ! matrix pools for the kpoint group, information on MOs is transferred using
585 ! generic mos structure
586 NULLIFY (mpools)
587 CALL mpools_create(mpools=mpools)
588 CALL mpools_rebuild_fm_pools(mpools=mpools, mos=mos, &
589 blacs_env=blacs_env, para_env=kpoint%para_env_kp)
590
591 IF (aux_fit) THEN
592 kpoint%mpools_aux_fit => mpools
593 ELSE
594 kpoint%mpools => mpools
595 END IF
596
597 ! reset old number of MOs
598 DO is = 1, nspin
599 CALL set_mo_set(mos(is), nmo=nmorig(is))
600 END DO
601
602 ! allocate density matrices
603 CALL mpools_get(mpools, ao_ao_fm_pools=ao_ao_fm_pools)
604 ALLOCATE (fmlocal)
605 CALL fm_pool_create_fm(ao_ao_fm_pools(1)%pool, fmlocal)
606 CALL cp_fm_get_info(fmlocal, matrix_struct=matrix_struct)
607 DO ik = 1, nkp_loc
608 IF (aux_fit) THEN
609 kp => kpoint%kp_aux_env(ik)%kpoint_env
610 ELSE
611 kp => kpoint%kp_env(ik)%kpoint_env
612 END IF
613 ! density matrix
614 CALL cp_fm_release(kp%pmat)
615 ALLOCATE (kp%pmat(nc, nspin))
616 DO is = 1, nspin
617 DO ic = 1, nc
618 CALL cp_fm_create(kp%pmat(ic, is), matrix_struct)
619 END DO
620 END DO
621 ! energy weighted density matrix
622 CALL cp_fm_release(kp%wmat)
623 ALLOCATE (kp%wmat(nc, nspin))
624 DO is = 1, nspin
625 DO ic = 1, nc
626 CALL cp_fm_create(kp%wmat(ic, is), matrix_struct)
627 END DO
628 END DO
629 END DO
630 CALL fm_pool_give_back_fm(ao_ao_fm_pools(1)%pool, fmlocal)
631 DEALLOCATE (fmlocal)
632
633 END IF
634
635 END IF
636
637 CALL timestop(handle)
638
639 END SUBROUTINE kpoint_initialize_mos
640
641! **************************************************************************************************
642!> \brief ...
643!> \param kpoint ...
644! **************************************************************************************************
645 SUBROUTINE kpoint_initialize_mo_set(kpoint)
646 TYPE(kpoint_type), POINTER :: kpoint
647
648 CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_initialize_mo_set'
649
650 INTEGER :: handle, ic, ik, ikk, ispin
651 TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_mo_fm_pools
652 TYPE(cp_fm_type), POINTER :: mo_coeff
653 TYPE(mo_set_type), DIMENSION(:, :), POINTER :: moskp
654
655 CALL timeset(routinen, handle)
656
657 DO ik = 1, SIZE(kpoint%kp_env)
658 CALL mpools_get(kpoint%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
659 moskp => kpoint%kp_env(ik)%kpoint_env%mos
660 ikk = kpoint%kp_range(1) + ik - 1
661 cpassert(ASSOCIATED(moskp))
662 DO ispin = 1, SIZE(moskp, 2)
663 DO ic = 1, SIZE(moskp, 1)
664 CALL get_mo_set(moskp(ic, ispin), mo_coeff=mo_coeff)
665 IF (.NOT. ASSOCIATED(mo_coeff)) THEN
666 CALL init_mo_set(moskp(ic, ispin), &
667 fm_pool=ao_mo_fm_pools(ispin)%pool, name="kpoints")
668 END IF
669 END DO
670 END DO
671 END DO
672
673 CALL timestop(handle)
674
675 END SUBROUTINE kpoint_initialize_mo_set
676
677! **************************************************************************************************
678!> \brief Generates the mapping of cell indices and linear RS index
679!> CELL (0,0,0) is always mapped to index 1
680!> \param kpoint Kpoint environment
681!> \param sab_nl Defining neighbour list
682!> \param para_env Parallel environment
683!> \param nimages [output]
684! **************************************************************************************************
685 SUBROUTINE kpoint_init_cell_index(kpoint, sab_nl, para_env, nimages)
686
687 TYPE(kpoint_type), POINTER :: kpoint
688 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
689 POINTER :: sab_nl
690 TYPE(mp_para_env_type), POINTER :: para_env
691 INTEGER, INTENT(OUT) :: nimages
692
693 CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_init_cell_index'
694
695 INTEGER :: handle, i1, i2, i3, ic, icount, it, &
696 ncount
697 INTEGER, DIMENSION(3) :: cell, itm
698 INTEGER, DIMENSION(:, :), POINTER :: index_to_cell, list
699 INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index, cti
700 LOGICAL :: new
702 DIMENSION(:), POINTER :: nl_iterator
703
704 NULLIFY (cell_to_index, index_to_cell)
705
706 CALL timeset(routinen, handle)
707
708 cpassert(ASSOCIATED(kpoint))
709
710 ALLOCATE (list(3, 125))
711 list = 0
712 icount = 1
713
714 CALL neighbor_list_iterator_create(nl_iterator, sab_nl)
715 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
716 CALL get_iterator_info(nl_iterator, cell=cell)
717
718 new = .true.
719 DO ic = 1, icount
720 IF (cell(1) == list(1, ic) .AND. cell(2) == list(2, ic) .AND. &
721 cell(3) == list(3, ic)) THEN
722 new = .false.
723 EXIT
724 END IF
725 END DO
726 IF (new) THEN
727 icount = icount + 1
728 IF (icount > SIZE(list, 2)) THEN
729 CALL reallocate(list, 1, 3, 1, 2*SIZE(list, 2))
730 END IF
731 list(1:3, icount) = cell(1:3)
732 END IF
733
734 END DO
735 CALL neighbor_list_iterator_release(nl_iterator)
736
737 itm(1) = maxval(abs(list(1, 1:icount)))
738 itm(2) = maxval(abs(list(2, 1:icount)))
739 itm(3) = maxval(abs(list(3, 1:icount)))
740 CALL para_env%max(itm)
741 it = maxval(itm(1:3))
742 IF (ASSOCIATED(kpoint%cell_to_index)) THEN
743 DEALLOCATE (kpoint%cell_to_index)
744 END IF
745 ALLOCATE (kpoint%cell_to_index(-itm(1):itm(1), -itm(2):itm(2), -itm(3):itm(3)))
746 cell_to_index => kpoint%cell_to_index
747 cti => cell_to_index
748 cti(:, :, :) = 0
749 DO ic = 1, icount
750 i1 = list(1, ic)
751 i2 = list(2, ic)
752 i3 = list(3, ic)
753 cti(i1, i2, i3) = ic
754 END DO
755 CALL para_env%sum(cti)
756 ncount = 0
757 DO i1 = -itm(1), itm(1)
758 DO i2 = -itm(2), itm(2)
759 DO i3 = -itm(3), itm(3)
760 IF (cti(i1, i2, i3) == 0) THEN
761 cti(i1, i2, i3) = 1000000
762 ELSE
763 ncount = ncount + 1
764 cti(i1, i2, i3) = (abs(i1) + abs(i2) + abs(i3))*1000 + abs(i3)*100 + abs(i2)*10 + abs(i1)
765 cti(i1, i2, i3) = cti(i1, i2, i3) + (i1 + i2 + i3)
766 END IF
767 END DO
768 END DO
769 END DO
770
771 IF (ASSOCIATED(kpoint%index_to_cell)) THEN
772 DEALLOCATE (kpoint%index_to_cell)
773 END IF
774 ALLOCATE (kpoint%index_to_cell(3, ncount))
775 index_to_cell => kpoint%index_to_cell
776 DO ic = 1, ncount
777 cell = minloc(cti)
778 i1 = cell(1) - 1 - itm(1)
779 i2 = cell(2) - 1 - itm(2)
780 i3 = cell(3) - 1 - itm(3)
781 cti(i1, i2, i3) = 1000000
782 index_to_cell(1, ic) = i1
783 index_to_cell(2, ic) = i2
784 index_to_cell(3, ic) = i3
785 END DO
786 cti(:, :, :) = 0
787 DO ic = 1, ncount
788 i1 = index_to_cell(1, ic)
789 i2 = index_to_cell(2, ic)
790 i3 = index_to_cell(3, ic)
791 cti(i1, i2, i3) = ic
792 END DO
793
794 ! keep pointer to this neighborlist
795 kpoint%sab_nl => sab_nl
796
797 ! set number of images
798 nimages = SIZE(index_to_cell, 2)
799
800 DEALLOCATE (list)
801
802 CALL timestop(handle)
803
804 END SUBROUTINE kpoint_init_cell_index
805
806! **************************************************************************************************
807!> \brief Transformation of real space matrices to a kpoint
808!> \param rmatrix Real part of kpoint matrix
809!> \param cmatrix Complex part of kpoint matrix (optional)
810!> \param rsmat Real space matrices
811!> \param ispin Spin index
812!> \param xkp Kpoint coordinates
813!> \param cell_to_index mapping of cell indices to RS index
814!> \param sab_nl Defining neighbor list
815!> \param is_complex Matrix to be transformed is imaginary
816!> \param rs_sign Matrix to be transformed is csaled by rs_sign
817! **************************************************************************************************
818 SUBROUTINE rskp_transform(rmatrix, cmatrix, rsmat, ispin, &
819 xkp, cell_to_index, sab_nl, is_complex, rs_sign)
820
821 TYPE(dbcsr_type) :: rmatrix
822 TYPE(dbcsr_type), OPTIONAL :: cmatrix
823 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: rsmat
824 INTEGER, INTENT(IN) :: ispin
825 REAL(kind=dp), DIMENSION(3), INTENT(IN) :: xkp
826 INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
827 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
828 POINTER :: sab_nl
829 LOGICAL, INTENT(IN), OPTIONAL :: is_complex
830 REAL(kind=dp), INTENT(IN), OPTIONAL :: rs_sign
831
832 CHARACTER(LEN=*), PARAMETER :: routinen = 'rskp_transform'
833
834 INTEGER :: handle, iatom, ic, icol, irow, jatom, &
835 nimg
836 INTEGER, DIMENSION(3) :: cell
837 LOGICAL :: do_symmetric, found, my_complex, &
838 wfn_real_only
839 REAL(kind=dp) :: arg, coskl, fsign, fsym, sinkl
840 REAL(kind=dp), DIMENSION(:, :), POINTER :: cblock, rblock, rsblock
842 DIMENSION(:), POINTER :: nl_iterator
843
844 CALL timeset(routinen, handle)
845
846 my_complex = .false.
847 IF (PRESENT(is_complex)) my_complex = is_complex
848
849 fsign = 1.0_dp
850 IF (PRESENT(rs_sign)) fsign = rs_sign
851
852 wfn_real_only = .true.
853 IF (PRESENT(cmatrix)) wfn_real_only = .false.
854
855 nimg = SIZE(rsmat, 2)
856
857 CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
858
859 CALL neighbor_list_iterator_create(nl_iterator, sab_nl)
860 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
861 CALL get_iterator_info(nl_iterator, iatom=iatom, jatom=jatom, cell=cell)
862
863 ! fsym = +- 1 is due to real space matrices being non-symmetric (although in a symmtric type)
864 ! with the link S_mu^0,nu^b = S_nu^0,mu^-b, and the KP matrices beeing Hermitian
865 fsym = 1.0_dp
866 irow = iatom
867 icol = jatom
868 IF (do_symmetric .AND. (iatom > jatom)) THEN
869 irow = jatom
870 icol = iatom
871 fsym = -1.0_dp
872 END IF
873
874 ic = cell_to_index(cell(1), cell(2), cell(3))
875 IF (ic < 1 .OR. ic > nimg) cycle
876
877 arg = real(cell(1), dp)*xkp(1) + real(cell(2), dp)*xkp(2) + real(cell(3), dp)*xkp(3)
878 IF (my_complex) THEN
879 coskl = fsign*fsym*cos(twopi*arg)
880 sinkl = fsign*sin(twopi*arg)
881 ELSE
882 coskl = fsign*cos(twopi*arg)
883 sinkl = fsign*fsym*sin(twopi*arg)
884 END IF
885
886 CALL dbcsr_get_block_p(matrix=rsmat(ispin, ic)%matrix, row=irow, col=icol, &
887 block=rsblock, found=found)
888 IF (.NOT. found) cycle
889
890 IF (wfn_real_only) THEN
891 CALL dbcsr_get_block_p(matrix=rmatrix, row=irow, col=icol, &
892 block=rblock, found=found)
893 IF (.NOT. found) cycle
894 rblock = rblock + coskl*rsblock
895 ELSE
896 CALL dbcsr_get_block_p(matrix=rmatrix, row=irow, col=icol, &
897 block=rblock, found=found)
898 IF (.NOT. found) cycle
899 CALL dbcsr_get_block_p(matrix=cmatrix, row=irow, col=icol, &
900 block=cblock, found=found)
901 IF (.NOT. found) cycle
902 rblock = rblock + coskl*rsblock
903 cblock = cblock + sinkl*rsblock
904 END IF
905
906 END DO
907 CALL neighbor_list_iterator_release(nl_iterator)
908
909 CALL timestop(handle)
910
911 END SUBROUTINE rskp_transform
912
913! **************************************************************************************************
914!> \brief Given the eigenvalues of all kpoints, calculates the occupation numbers
915!> \param kpoint Kpoint environment
916!> \param smear Smearing information
917!> \param probe ...
918! **************************************************************************************************
919 SUBROUTINE kpoint_set_mo_occupation(kpoint, smear, probe)
920
921 TYPE(kpoint_type), POINTER :: kpoint
922 TYPE(smear_type) :: smear
923 TYPE(hairy_probes_type), DIMENSION(:), OPTIONAL, &
924 POINTER :: probe
925
926 CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_set_mo_occupation'
927
928 INTEGER :: handle, ik, ikpgr, ispin, kplocal, nao, &
929 nb, ncol_global, ne_a, ne_b, &
930 nelectron, nkp, nmo, nrow_global, nspin
931 INTEGER, DIMENSION(2) :: kp_range
932 REAL(kind=dp) :: kts, mu, mus(2), nel
933 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: smatrix
934 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: weig, wocc
935 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :, :) :: icoeff, rcoeff
936 REAL(kind=dp), DIMENSION(:), POINTER :: eigenvalues, occupation, wkp
937 TYPE(cp_fm_type), POINTER :: mo_coeff
938 TYPE(kpoint_env_type), POINTER :: kp
939 TYPE(mo_set_type), POINTER :: mo_set
940 TYPE(mp_para_env_type), POINTER :: para_env_inter_kp
941
942 CALL timeset(routinen, handle)
943
944 ! first collect all the eigenvalues
945 CALL get_kpoint_info(kpoint, nkp=nkp)
946 kp => kpoint%kp_env(1)%kpoint_env
947 nspin = SIZE(kp%mos, 2)
948 mo_set => kp%mos(1, 1)
949 CALL get_mo_set(mo_set, nmo=nmo, nao=nao, nelectron=nelectron)
950 ne_a = nelectron
951 IF (nspin == 2) THEN
952 CALL get_mo_set(kp%mos(1, 2), nmo=nb, nelectron=ne_b)
953 cpassert(nmo == nb)
954 END IF
955 ALLOCATE (weig(nmo, nkp, nspin), wocc(nmo, nkp, nspin))
956 weig = 0.0_dp
957 wocc = 0.0_dp
958 IF (PRESENT(probe)) THEN
959 ALLOCATE (rcoeff(nao, nmo, nkp, nspin), icoeff(nao, nmo, nkp, nspin))
960 rcoeff = 0.0_dp !coeff, real part
961 icoeff = 0.0_dp !coeff, imaginary part
962 END IF
963 CALL get_kpoint_info(kpoint, kp_range=kp_range)
964 kplocal = kp_range(2) - kp_range(1) + 1
965 DO ikpgr = 1, kplocal
966 ik = kp_range(1) + ikpgr - 1
967 kp => kpoint%kp_env(ikpgr)%kpoint_env
968 DO ispin = 1, nspin
969 mo_set => kp%mos(1, ispin)
970 CALL get_mo_set(mo_set, eigenvalues=eigenvalues)
971 weig(1:nmo, ik, ispin) = eigenvalues(1:nmo)
972 IF (PRESENT(probe)) THEN
973 CALL get_mo_set(mo_set, mo_coeff=mo_coeff)
974 CALL cp_fm_get_info(mo_coeff, &
975 nrow_global=nrow_global, &
976 ncol_global=ncol_global)
977 ALLOCATE (smatrix(nrow_global, ncol_global))
978 CALL cp_fm_get_submatrix(mo_coeff, smatrix)
979
980 rcoeff(1:nao, 1:nmo, ik, ispin) = smatrix(1:nrow_global, 1:ncol_global)
981
982 DEALLOCATE (smatrix)
983
984 mo_set => kp%mos(2, ispin)
985
986 CALL get_mo_set(mo_set, mo_coeff=mo_coeff)
987 CALL cp_fm_get_info(mo_coeff, &
988 nrow_global=nrow_global, &
989 ncol_global=ncol_global)
990 ALLOCATE (smatrix(nrow_global, ncol_global))
991 CALL cp_fm_get_submatrix(mo_coeff, smatrix)
992
993 icoeff(1:nao, 1:nmo, ik, ispin) = smatrix(1:nrow_global, 1:ncol_global)
994
995 mo_set => kp%mos(1, ispin)
996
997 DEALLOCATE (smatrix)
998 END IF
999 END DO
1000 END DO
1001 CALL get_kpoint_info(kpoint, para_env_inter_kp=para_env_inter_kp)
1002 CALL para_env_inter_kp%sum(weig)
1003
1004 IF (PRESENT(probe)) THEN
1005 CALL para_env_inter_kp%sum(rcoeff)
1006 CALL para_env_inter_kp%sum(icoeff)
1007 END IF
1008
1009 CALL get_kpoint_info(kpoint, wkp=wkp)
1010
1011!calling of HP module HERE, before smear
1012 IF (PRESENT(probe)) THEN
1013 smear%do_smear = .false. !ensures smearing is switched off
1014
1015 IF (nspin == 1) THEN
1016 nel = real(nelectron, kind=dp)
1017 CALL probe_occupancy_kp(wocc(:, :, :), mus(1), kts, weig(:, :, :), rcoeff(:, :, :, :), icoeff(:, :, :, :), 2.0d0, &
1018 probe, nel, wkp)
1019 ELSE
1020 nel = real(ne_a, kind=dp) + real(ne_b, kind=dp)
1021 CALL probe_occupancy_kp(wocc(:, :, :), mu, kts, weig(:, :, :), rcoeff(:, :, :, :), icoeff(:, :, :, :), 1.0d0, &
1022 probe, nel, wkp)
1023 kts = kts/2._dp
1024 mus(1:2) = mu
1025 END IF
1026
1027 DO ikpgr = 1, kplocal
1028 ik = kp_range(1) + ikpgr - 1
1029 kp => kpoint%kp_env(ikpgr)%kpoint_env
1030 DO ispin = 1, nspin
1031 mo_set => kp%mos(1, ispin)
1032 CALL get_mo_set(mo_set, eigenvalues=eigenvalues, occupation_numbers=occupation)
1033 eigenvalues(1:nmo) = weig(1:nmo, ik, ispin)
1034 occupation(1:nmo) = wocc(1:nmo, ik, ispin)
1035 mo_set%kTS = kts
1036 mo_set%mu = mus(ispin)
1037
1038 CALL get_mo_set(mo_set, mo_coeff=mo_coeff)
1039 !get smatrix for kpoint_env ikp
1040 CALL cp_fm_get_info(mo_coeff, &
1041 nrow_global=nrow_global, &
1042 ncol_global=ncol_global)
1043 ALLOCATE (smatrix(nrow_global, ncol_global))
1044 CALL cp_fm_get_submatrix(mo_coeff, smatrix)
1045
1046 smatrix(1:nrow_global, 1:ncol_global) = rcoeff(1:nao, 1:nmo, ik, ispin)
1047 DEALLOCATE (smatrix)
1048
1049 mo_set => kp%mos(2, ispin)
1050
1051 CALL get_mo_set(mo_set, mo_coeff=mo_coeff)
1052 !get smatrix for kpoint_env ikp
1053 CALL cp_fm_get_info(mo_coeff, &
1054 nrow_global=nrow_global, &
1055 ncol_global=ncol_global)
1056 ALLOCATE (smatrix(nrow_global, ncol_global))
1057 CALL cp_fm_get_submatrix(mo_coeff, smatrix)
1058
1059 smatrix(1:nrow_global, 1:ncol_global) = icoeff(1:nao, 1:nmo, ik, ispin)
1060 DEALLOCATE (smatrix)
1061
1062 mo_set => kp%mos(1, ispin)
1063
1064 END DO
1065 END DO
1066
1067 DEALLOCATE (weig, wocc, rcoeff, icoeff)
1068
1069 END IF
1070
1071 IF (PRESENT(probe) .EQV. .false.) THEN
1072 IF (smear%do_smear) THEN
1073 SELECT CASE (smear%method)
1074 CASE (smear_fermi_dirac)
1075 ! finite electronic temperature
1076 IF (nspin == 1) THEN
1077 nel = real(nelectron, kind=dp)
1078 CALL smearkp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, &
1079 smear%electronic_temperature, 2.0_dp, smear_fermi_dirac)
1080 ELSE IF (smear%fixed_mag_mom > 0.0_dp) THEN
1081 nel = real(ne_a, kind=dp)
1082 CALL smearkp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, &
1083 smear%electronic_temperature, 1.0_dp, smear_fermi_dirac)
1084 nel = real(ne_b, kind=dp)
1085 CALL smearkp(wocc(:, :, 2), mus(2), kts, weig(:, :, 2), nel, wkp, &
1086 smear%electronic_temperature, 1.0_dp, smear_fermi_dirac)
1087 ELSE
1088 nel = real(ne_a, kind=dp) + real(ne_b, kind=dp)
1089 CALL smearkp2(wocc(:, :, :), mu, kts, weig(:, :, :), nel, wkp, &
1090 smear%electronic_temperature, smear_fermi_dirac)
1091 kts = kts/2._dp
1092 mus(1:2) = mu
1093 END IF
1095 IF (nspin == 1) THEN
1096 nel = real(nelectron, kind=dp)
1097 CALL smearkp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, &
1098 smear%smearing_width, 2.0_dp, smear%method)
1099 ELSE IF (smear%fixed_mag_mom > 0.0_dp) THEN
1100 nel = real(ne_a, kind=dp)
1101 CALL smearkp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, &
1102 smear%smearing_width, 1.0_dp, smear%method)
1103 nel = real(ne_b, kind=dp)
1104 CALL smearkp(wocc(:, :, 2), mus(2), kts, weig(:, :, 2), nel, wkp, &
1105 smear%smearing_width, 1.0_dp, smear%method)
1106 ELSE
1107 nel = real(ne_a, kind=dp) + real(ne_b, kind=dp)
1108 CALL smearkp2(wocc(:, :, :), mu, kts, weig(:, :, :), nel, wkp, &
1109 smear%smearing_width, smear%method)
1110 kts = kts/2._dp
1111 mus(1:2) = mu
1112 END IF
1113 CASE DEFAULT
1114 cpabort("kpoints: Selected smearing not (yet) supported")
1115 END SELECT
1116 ELSE
1117 ! fixed occupations (2/1)
1118 IF (nspin == 1) THEN
1119 nel = real(nelectron, kind=dp)
1120 CALL smearkp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, &
1121 0.0_dp, 2.0_dp, smear_gaussian)
1122 ELSE
1123 nel = real(ne_a, kind=dp)
1124 CALL smearkp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, &
1125 0.0_dp, 1.0_dp, smear_gaussian)
1126 nel = real(ne_b, kind=dp)
1127 CALL smearkp(wocc(:, :, 2), mus(2), kts, weig(:, :, 2), nel, wkp, &
1128 0.0_dp, 1.0_dp, smear_gaussian)
1129 END IF
1130 END IF
1131 DO ikpgr = 1, kplocal
1132 ik = kp_range(1) + ikpgr - 1
1133 kp => kpoint%kp_env(ikpgr)%kpoint_env
1134 DO ispin = 1, nspin
1135 mo_set => kp%mos(1, ispin)
1136 CALL get_mo_set(mo_set, eigenvalues=eigenvalues, occupation_numbers=occupation)
1137 eigenvalues(1:nmo) = weig(1:nmo, ik, ispin)
1138 occupation(1:nmo) = wocc(1:nmo, ik, ispin)
1139 mo_set%kTS = kts
1140 mo_set%mu = mus(ispin)
1141 END DO
1142 END DO
1143
1144 DEALLOCATE (weig, wocc)
1145
1146 END IF
1147
1148 CALL timestop(handle)
1149
1150 END SUBROUTINE kpoint_set_mo_occupation
1151
1152! **************************************************************************************************
1153!> \brief Calculate kpoint density matrices (rho(k), owned by kpoint groups)
1154!> \param kpoint kpoint environment
1155!> \param energy_weighted calculate energy weighted density matrix
1156!> \param for_aux_fit ...
1157! **************************************************************************************************
1158 SUBROUTINE kpoint_density_matrices(kpoint, energy_weighted, for_aux_fit)
1159
1160 TYPE(kpoint_type), POINTER :: kpoint
1161 LOGICAL, OPTIONAL :: energy_weighted, for_aux_fit
1162
1163 CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_density_matrices'
1164
1165 INTEGER :: handle, ikpgr, ispin, kplocal, nao, nmo, &
1166 nspin
1167 INTEGER, DIMENSION(2) :: kp_range
1168 LOGICAL :: aux_fit, wtype
1169 REAL(kind=dp), DIMENSION(:), POINTER :: eigenvalues, occupation
1170 TYPE(cp_fm_struct_type), POINTER :: matrix_struct
1171 TYPE(cp_fm_type) :: fwork
1172 TYPE(cp_fm_type), POINTER :: cpmat, pmat, rpmat
1173 TYPE(kpoint_env_type), POINTER :: kp
1174 TYPE(mo_set_type), POINTER :: mo_set
1175
1176 CALL timeset(routinen, handle)
1177
1178 IF (PRESENT(energy_weighted)) THEN
1179 wtype = energy_weighted
1180 ELSE
1181 ! default is normal density matrix
1182 wtype = .false.
1183 END IF
1184
1185 IF (PRESENT(for_aux_fit)) THEN
1186 aux_fit = for_aux_fit
1187 ELSE
1188 aux_fit = .false.
1189 END IF
1190
1191 IF (aux_fit) THEN
1192 cpassert(ASSOCIATED(kpoint%kp_aux_env))
1193 END IF
1194
1195 ! work matrix
1196 IF (aux_fit) THEN
1197 mo_set => kpoint%kp_aux_env(1)%kpoint_env%mos(1, 1)
1198 ELSE
1199 mo_set => kpoint%kp_env(1)%kpoint_env%mos(1, 1)
1200 END IF
1201 CALL get_mo_set(mo_set, nao=nao, nmo=nmo)
1202 CALL cp_fm_get_info(mo_set%mo_coeff, matrix_struct=matrix_struct)
1203 CALL cp_fm_create(fwork, matrix_struct)
1204
1205 CALL get_kpoint_info(kpoint, kp_range=kp_range)
1206 kplocal = kp_range(2) - kp_range(1) + 1
1207 DO ikpgr = 1, kplocal
1208 IF (aux_fit) THEN
1209 kp => kpoint%kp_aux_env(ikpgr)%kpoint_env
1210 ELSE
1211 kp => kpoint%kp_env(ikpgr)%kpoint_env
1212 END IF
1213 nspin = SIZE(kp%mos, 2)
1214 DO ispin = 1, nspin
1215 mo_set => kp%mos(1, ispin)
1216 IF (wtype) THEN
1217 CALL get_mo_set(mo_set, eigenvalues=eigenvalues)
1218 END IF
1219 IF (kpoint%use_real_wfn) THEN
1220 IF (wtype) THEN
1221 pmat => kp%wmat(1, ispin)
1222 ELSE
1223 pmat => kp%pmat(1, ispin)
1224 END IF
1225 CALL get_mo_set(mo_set, occupation_numbers=occupation)
1226 CALL cp_fm_to_fm(mo_set%mo_coeff, fwork)
1227 CALL cp_fm_column_scale(fwork, occupation)
1228 IF (wtype) THEN
1229 CALL cp_fm_column_scale(fwork, eigenvalues)
1230 END IF
1231 CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 0.0_dp, pmat)
1232 ELSE
1233 IF (wtype) THEN
1234 rpmat => kp%wmat(1, ispin)
1235 cpmat => kp%wmat(2, ispin)
1236 ELSE
1237 rpmat => kp%pmat(1, ispin)
1238 cpmat => kp%pmat(2, ispin)
1239 END IF
1240 CALL get_mo_set(mo_set, occupation_numbers=occupation)
1241 CALL cp_fm_to_fm(mo_set%mo_coeff, fwork)
1242 CALL cp_fm_column_scale(fwork, occupation)
1243 IF (wtype) THEN
1244 CALL cp_fm_column_scale(fwork, eigenvalues)
1245 END IF
1246 ! Re(c)*Re(c)
1247 CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 0.0_dp, rpmat)
1248 mo_set => kp%mos(2, ispin)
1249 ! Im(c)*Re(c)
1250 CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 0.0_dp, cpmat)
1251 ! Re(c)*Im(c)
1252 CALL parallel_gemm("N", "T", nao, nao, nmo, -1.0_dp, fwork, mo_set%mo_coeff, 1.0_dp, cpmat)
1253 CALL cp_fm_to_fm(mo_set%mo_coeff, fwork)
1254 CALL cp_fm_column_scale(fwork, occupation)
1255 IF (wtype) THEN
1256 CALL cp_fm_column_scale(fwork, eigenvalues)
1257 END IF
1258 ! Im(c)*Im(c)
1259 CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 1.0_dp, rpmat)
1260 END IF
1261 END DO
1262 END DO
1263
1264 CALL cp_fm_release(fwork)
1265
1266 CALL timestop(handle)
1267
1268 END SUBROUTINE kpoint_density_matrices
1269
1270! **************************************************************************************************
1271!> \brief Calculate Lowdin transformation of density matrix S^1/2 P S^1/2
1272!> Integrate diagonal elements over k-points to get Lowdin charges
1273!> \param kpoint kpoint environment
1274!> \param pmat_diag Sum over kpoints of diagonal elements
1275!> \par History
1276!> 04.2026 created [JGH]
1277! **************************************************************************************************
1278 SUBROUTINE lowdin_kp_trans(kpoint, pmat_diag)
1279
1280 TYPE(kpoint_type), POINTER :: kpoint
1281 REAL(kind=dp), DIMENSION(:, :), INTENT(INOUT) :: pmat_diag
1282
1283 CHARACTER(LEN=*), PARAMETER :: routinen = 'lowdin_kp_trans'
1284 COMPLEX(KIND=dp), PARAMETER :: cone = (1.0_dp, 0.0_dp), &
1285 czero = (0.0_dp, 0.0_dp)
1286
1287 INTEGER :: handle, ikpgr, ispin, kplocal, nao, nspin
1288 INTEGER, DIMENSION(2) :: kp_range
1289 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: dele
1290 TYPE(cp_cfm_type) :: cf1work, cf2work
1291 TYPE(cp_cfm_type), POINTER :: cshalf
1292 TYPE(cp_fm_struct_type), POINTER :: matrix_struct
1293 TYPE(cp_fm_type) :: f1work, f2work
1294 TYPE(cp_fm_type), POINTER :: cpmat, pmat, rpmat, shalf
1295 TYPE(kpoint_env_type), POINTER :: kp
1296 TYPE(mp_para_env_type), POINTER :: para_env_inter_kp
1297
1298 CALL timeset(routinen, handle)
1299
1300 nspin = SIZE(pmat_diag, 2)
1301 pmat_diag = 0.0_dp
1302
1303 ! work matrix
1304 CALL cp_fm_get_info(kpoint%kp_env(1)%kpoint_env%pmat(1, 1), &
1305 matrix_struct=matrix_struct, nrow_global=nao)
1306 IF (kpoint%use_real_wfn) THEN
1307 CALL cp_fm_create(f1work, matrix_struct, nrow=nao, ncol=nao)
1308 CALL cp_fm_create(f2work, matrix_struct, nrow=nao, ncol=nao)
1309 ELSE
1310 CALL cp_fm_create(f2work, matrix_struct, nrow=nao, ncol=nao)
1311 CALL cp_cfm_create(cf1work, matrix_struct, nrow=nao, ncol=nao)
1312 CALL cp_cfm_create(cf2work, matrix_struct, nrow=nao, ncol=nao)
1313 END IF
1314 ALLOCATE (dele(nao))
1315
1316 CALL get_kpoint_info(kpoint, kp_range=kp_range)
1317 kplocal = kp_range(2) - kp_range(1) + 1
1318 DO ikpgr = 1, kplocal
1319 kp => kpoint%kp_env(ikpgr)%kpoint_env
1320 DO ispin = 1, nspin
1321 IF (kpoint%use_real_wfn) THEN
1322 pmat => kp%pmat(1, ispin)
1323 shalf => kp%shalf
1324 CALL parallel_gemm("N", "N", nao, nao, nao, 1.0_dp, pmat, shalf, 0.0_dp, f1work)
1325 CALL parallel_gemm("N", "N", nao, nao, nao, 1.0_dp, shalf, f1work, 0.0_dp, f2work)
1326 ELSE
1327 rpmat => kp%pmat(1, ispin)
1328 cpmat => kp%pmat(2, ispin)
1329 cshalf => kp%cshalf
1330 CALL cp_fm_to_cfm(rpmat, cpmat, cf1work)
1331 CALL parallel_gemm("N", "N", nao, nao, nao, cone, cf1work, cshalf, czero, cf2work)
1332 CALL parallel_gemm("N", "N", nao, nao, nao, cone, cshalf, cf2work, czero, cf1work)
1333 CALL cp_cfm_to_fm(cf1work, mtargetr=f2work)
1334 END IF
1335 CALL cp_fm_get_diag(f2work, dele)
1336 pmat_diag(1:nao, ispin) = pmat_diag(1:nao, ispin) + kp%wkp*dele(1:nao)
1337 END DO
1338 END DO
1339
1340 CALL get_kpoint_info(kpoint, para_env_inter_kp=para_env_inter_kp)
1341 CALL para_env_inter_kp%sum(pmat_diag)
1342
1343 IF (kpoint%use_real_wfn) THEN
1344 CALL cp_fm_release(f1work)
1345 CALL cp_fm_release(f2work)
1346 ELSE
1347 CALL cp_fm_release(f2work)
1348 CALL cp_cfm_release(cf1work)
1349 CALL cp_cfm_release(cf2work)
1350 END IF
1351 DEALLOCATE (dele)
1352
1353 CALL timestop(handle)
1354
1355 END SUBROUTINE lowdin_kp_trans
1356
1357! **************************************************************************************************
1358!> \brief generate real space density matrices in DBCSR format
1359!> \param kpoint Kpoint environment
1360!> \param denmat Real space (DBCSR) density matrices
1361!> \param wtype True = energy weighted density matrix
1362!> False = normal density matrix
1363!> \param tempmat DBCSR matrix to be used as template
1364!> \param sab_nl ...
1365!> \param fmwork FM work matrices (kpoint group)
1366!> \param for_aux_fit ...
1367!> \param pmat_ext ...
1368! **************************************************************************************************
1369 SUBROUTINE kpoint_density_transform(kpoint, denmat, wtype, tempmat, sab_nl, fmwork, for_aux_fit, pmat_ext)
1370
1371 TYPE(kpoint_type), POINTER :: kpoint
1372 TYPE(dbcsr_p_type), DIMENSION(:, :) :: denmat
1373 LOGICAL, INTENT(IN) :: wtype
1374 TYPE(dbcsr_type), POINTER :: tempmat
1375 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1376 POINTER :: sab_nl
1377 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fmwork
1378 LOGICAL, OPTIONAL :: for_aux_fit
1379 TYPE(cp_fm_type), DIMENSION(:, :, :), INTENT(IN), &
1380 OPTIONAL :: pmat_ext
1381
1382 CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_density_transform'
1383
1384 INTEGER :: handle, ic, ik, ikk, indx, ir, ira, is, &
1385 ispin, jr, nc, nimg, nkp, nspin
1386 INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
1387 LOGICAL :: aux_fit, do_ext, do_symmetric, my_kpgrp, &
1388 real_only
1389 REAL(kind=dp) :: wkpx
1390 REAL(kind=dp), DIMENSION(:), POINTER :: wkp
1391 REAL(kind=dp), DIMENSION(:, :), POINTER :: xkp
1392 TYPE(copy_info_type), ALLOCATABLE, DIMENSION(:) :: info
1393 TYPE(cp_fm_type) :: fmdummy
1394 TYPE(dbcsr_type), POINTER :: cpmat, rpmat, scpmat, srpmat
1395 TYPE(kind_rotmat_type), DIMENSION(:), POINTER :: kind_rot
1396 TYPE(kpoint_env_type), POINTER :: kp
1397 TYPE(kpoint_sym_type), POINTER :: kpsym
1398 TYPE(mp_para_env_type), POINTER :: para_env
1399
1400 CALL timeset(routinen, handle)
1401
1402 CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
1403
1404 IF (PRESENT(for_aux_fit)) THEN
1405 aux_fit = for_aux_fit
1406 ELSE
1407 aux_fit = .false.
1408 END IF
1409
1410 do_ext = .false.
1411 IF (PRESENT(pmat_ext)) do_ext = .true.
1412
1413 IF (aux_fit) THEN
1414 cpassert(ASSOCIATED(kpoint%kp_aux_env))
1415 END IF
1416
1417 ! work storage
1418 ALLOCATE (rpmat)
1419 CALL dbcsr_create(rpmat, template=tempmat, &
1420 matrix_type=merge(dbcsr_type_symmetric, dbcsr_type_no_symmetry, do_symmetric))
1421 CALL cp_dbcsr_alloc_block_from_nbl(rpmat, sab_nl)
1422 CALL dbcsr_set(rpmat, 0.0_dp)
1423 ALLOCATE (cpmat)
1424 CALL dbcsr_create(cpmat, template=tempmat, &
1425 matrix_type=merge(dbcsr_type_antisymmetric, dbcsr_type_no_symmetry, do_symmetric))
1426 CALL cp_dbcsr_alloc_block_from_nbl(cpmat, sab_nl)
1427 CALL dbcsr_set(cpmat, 0.0_dp)
1428 IF (.NOT. kpoint%full_grid) THEN
1429 ALLOCATE (srpmat)
1430 CALL dbcsr_create(srpmat, template=rpmat)
1431 CALL cp_dbcsr_alloc_block_from_nbl(srpmat, sab_nl)
1432 CALL dbcsr_set(srpmat, 0.0_dp)
1433 ALLOCATE (scpmat)
1434 CALL dbcsr_create(scpmat, template=cpmat)
1435 CALL cp_dbcsr_alloc_block_from_nbl(scpmat, sab_nl)
1436 CALL dbcsr_set(scpmat, 0.0_dp)
1437 END IF
1438
1439 CALL get_kpoint_info(kpoint, nkp=nkp, xkp=xkp, wkp=wkp, &
1440 cell_to_index=cell_to_index)
1441 ! initialize real space density matrices
1442 IF (aux_fit) THEN
1443 kp => kpoint%kp_aux_env(1)%kpoint_env
1444 ELSE
1445 kp => kpoint%kp_env(1)%kpoint_env
1446 END IF
1447 nspin = SIZE(kp%mos, 2)
1448 nc = SIZE(kp%mos, 1)
1449 nimg = SIZE(denmat, 2)
1450 real_only = (nc == 1)
1451
1452 para_env => kpoint%blacs_env_all%para_env
1453 ALLOCATE (info(nspin*nkp*nc))
1454
1455 ! Start all the communication
1456 indx = 0
1457 DO ispin = 1, nspin
1458 DO ic = 1, nimg
1459 CALL dbcsr_set(denmat(ispin, ic)%matrix, 0.0_dp)
1460 END DO
1461 !
1462 DO ik = 1, nkp
1463 my_kpgrp = (ik >= kpoint%kp_range(1) .AND. ik <= kpoint%kp_range(2))
1464 IF (my_kpgrp) THEN
1465 ikk = ik - kpoint%kp_range(1) + 1
1466 IF (aux_fit) THEN
1467 kp => kpoint%kp_aux_env(ikk)%kpoint_env
1468 ELSE
1469 kp => kpoint%kp_env(ikk)%kpoint_env
1470 END IF
1471 ELSE
1472 NULLIFY (kp)
1473 END IF
1474 ! collect this density matrix on all processors
1475 cpassert(SIZE(fmwork) >= nc)
1476
1477 IF (my_kpgrp) THEN
1478 DO ic = 1, nc
1479 indx = indx + 1
1480 IF (do_ext) THEN
1481 CALL cp_fm_start_copy_general(pmat_ext(ikk, ic, ispin), fmwork(ic), para_env, info(indx))
1482 ELSE
1483 IF (wtype) THEN
1484 CALL cp_fm_start_copy_general(kp%wmat(ic, ispin), fmwork(ic), para_env, info(indx))
1485 ELSE
1486 CALL cp_fm_start_copy_general(kp%pmat(ic, ispin), fmwork(ic), para_env, info(indx))
1487 END IF
1488 END IF
1489 END DO
1490 ELSE
1491 DO ic = 1, nc
1492 indx = indx + 1
1493 CALL cp_fm_start_copy_general(fmdummy, fmwork(ic), para_env, info(indx))
1494 END DO
1495 END IF
1496 END DO
1497 END DO
1498
1499 ! Finish communication and transform the received matrices
1500 indx = 0
1501 DO ispin = 1, nspin
1502 DO ik = 1, nkp
1503 DO ic = 1, nc
1504 indx = indx + 1
1505 CALL cp_fm_finish_copy_general(fmwork(ic), info(indx))
1506 END DO
1507
1508 ! reduce to dbcsr storage
1509 IF (real_only) THEN
1510 CALL copy_fm_to_dbcsr(fmwork(1), rpmat, keep_sparsity=.true.)
1511 ELSE
1512 CALL copy_fm_to_dbcsr(fmwork(1), rpmat, keep_sparsity=.true.)
1513 CALL copy_fm_to_dbcsr(fmwork(2), cpmat, keep_sparsity=.true.)
1514 END IF
1515
1516 ! symmetrization
1517 kpsym => kpoint%kp_sym(ik)%kpoint_sym
1518 cpassert(ASSOCIATED(kpsym))
1519
1520 IF (kpsym%apply_symmetry) THEN
1521 wkpx = wkp(ik)/real(kpsym%nwght, kind=dp)
1522 DO is = 1, kpsym%nwght
1523 ir = abs(kpsym%rotp(is))
1524 ira = 0
1525 DO jr = 1, SIZE(kpoint%ibrot)
1526 IF (ir == kpoint%ibrot(jr)) ira = jr
1527 END DO
1528 cpassert(ira > 0)
1529 kind_rot => kpoint%kind_rotmat(ira, :)
1530 IF (real_only) THEN
1531 CALL symtrans(srpmat, rpmat, kind_rot, kpsym%rot(1:3, 1:3, is), &
1532 kpsym%f0(:, is), kpoint%atype, symmetric=.true.)
1533 ELSE
1534 CALL symtrans(srpmat, rpmat, kind_rot, kpsym%rot(1:3, 1:3, is), &
1535 kpsym%f0(:, is), kpoint%atype, symmetric=.true.)
1536 CALL symtrans(scpmat, cpmat, kind_rot, kpsym%rot(1:3, 1:3, is), &
1537 kpsym%f0(:, is), kpoint%atype, antisymmetric=.true.)
1538 END IF
1539 CALL transform_dmat(denmat, srpmat, scpmat, ispin, real_only, sab_nl, &
1540 cell_to_index, kpsym%xkp(1:3, is), wkpx)
1541 END DO
1542 ELSE
1543 ! transformation
1544 CALL transform_dmat(denmat, rpmat, cpmat, ispin, real_only, sab_nl, &
1545 cell_to_index, xkp(1:3, ik), wkp(ik))
1546 END IF
1547 END DO
1548 END DO
1549
1550 ! Clean up communication
1551 indx = 0
1552 DO ispin = 1, nspin
1553 DO ik = 1, nkp
1554 my_kpgrp = (ik >= kpoint%kp_range(1) .AND. ik <= kpoint%kp_range(2))
1555 IF (my_kpgrp) THEN
1556 ikk = ik - kpoint%kp_range(1) + 1
1557 IF (aux_fit) THEN
1558 kp => kpoint%kp_aux_env(ikk)%kpoint_env
1559 ELSE
1560 kp => kpoint%kp_env(ikk)%kpoint_env
1561 END IF
1562
1563 DO ic = 1, nc
1564 indx = indx + 1
1565 CALL cp_fm_cleanup_copy_general(info(indx))
1566 END DO
1567 ELSE
1568 ! calls with dummy arguments, so not included
1569 ! therefore just increment counter by trip count
1570 indx = indx + nc
1571 END IF
1572 END DO
1573 END DO
1574
1575 ! All done
1576 DEALLOCATE (info)
1577
1578 CALL dbcsr_deallocate_matrix(rpmat)
1579 CALL dbcsr_deallocate_matrix(cpmat)
1580 IF (.NOT. kpoint%full_grid) THEN
1581 CALL dbcsr_deallocate_matrix(srpmat)
1582 CALL dbcsr_deallocate_matrix(scpmat)
1583 END IF
1584
1585 CALL timestop(handle)
1586
1587 END SUBROUTINE kpoint_density_transform
1588
1589! **************************************************************************************************
1590!> \brief real space density matrices in DBCSR format
1591!> \param denmat Real space (DBCSR) density matrix
1592!> \param rpmat ...
1593!> \param cpmat ...
1594!> \param ispin ...
1595!> \param real_only ...
1596!> \param sab_nl ...
1597!> \param cell_to_index ...
1598!> \param xkp ...
1599!> \param wkp ...
1600! **************************************************************************************************
1601 SUBROUTINE transform_dmat(denmat, rpmat, cpmat, ispin, real_only, sab_nl, cell_to_index, xkp, wkp)
1602
1603 TYPE(dbcsr_p_type), DIMENSION(:, :) :: denmat
1604 TYPE(dbcsr_type), POINTER :: rpmat, cpmat
1605 INTEGER, INTENT(IN) :: ispin
1606 LOGICAL, INTENT(IN) :: real_only
1607 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1608 POINTER :: sab_nl
1609 INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
1610 REAL(kind=dp), DIMENSION(3), INTENT(IN) :: xkp
1611 REAL(kind=dp), INTENT(IN) :: wkp
1612
1613 CHARACTER(LEN=*), PARAMETER :: routinen = 'transform_dmat'
1614
1615 INTEGER :: handle, iatom, icell, icol, irow, jatom, &
1616 nimg
1617 INTEGER, DIMENSION(3) :: cell
1618 LOGICAL :: do_symmetric, found
1619 REAL(kind=dp) :: arg, coskl, fc, sinkl
1620 REAL(kind=dp), DIMENSION(:, :), POINTER :: cblock, dblock, rblock
1622 DIMENSION(:), POINTER :: nl_iterator
1623
1624 CALL timeset(routinen, handle)
1625
1626 nimg = SIZE(denmat, 2)
1627
1628 ! transformation
1629 CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
1630 CALL neighbor_list_iterator_create(nl_iterator, sab_nl)
1631 DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1632 CALL get_iterator_info(nl_iterator, iatom=iatom, jatom=jatom, cell=cell)
1633
1634 !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
1635 !Therefore, we have: S(R) = sum_k Re(S(k))*cos(k*R) -i^2*Im(S(k))*sin(k*R)
1636 ! = sum_k Re(S(k))*cos(k*R) + Im(S(k))*sin(k*R)
1637 !fc = +- 1 is due to the usual non-symmetric real-space matrices stored as symmetric ones
1638
1639 irow = iatom
1640 icol = jatom
1641 fc = 1.0_dp
1642 IF (do_symmetric .AND. iatom > jatom) THEN
1643 irow = jatom
1644 icol = iatom
1645 fc = -1.0_dp
1646 END IF
1647
1648 icell = cell_to_index(cell(1), cell(2), cell(3))
1649 IF (icell < 1 .OR. icell > nimg) cycle
1650
1651 arg = real(cell(1), dp)*xkp(1) + real(cell(2), dp)*xkp(2) + real(cell(3), dp)*xkp(3)
1652 coskl = wkp*cos(twopi*arg)
1653 sinkl = wkp*fc*sin(twopi*arg)
1654
1655 CALL dbcsr_get_block_p(matrix=denmat(ispin, icell)%matrix, row=irow, col=icol, &
1656 block=dblock, found=found)
1657 IF (.NOT. found) cycle
1658
1659 IF (real_only) THEN
1660 CALL dbcsr_get_block_p(matrix=rpmat, row=irow, col=icol, block=rblock, found=found)
1661 IF (.NOT. found) cycle
1662 dblock = dblock + coskl*rblock
1663 ELSE
1664 CALL dbcsr_get_block_p(matrix=rpmat, row=irow, col=icol, block=rblock, found=found)
1665 IF (.NOT. found) cycle
1666 CALL dbcsr_get_block_p(matrix=cpmat, row=irow, col=icol, block=cblock, found=found)
1667 IF (.NOT. found) cycle
1668 dblock = dblock + coskl*rblock
1669 dblock = dblock + sinkl*cblock
1670 END IF
1671 END DO
1672 CALL neighbor_list_iterator_release(nl_iterator)
1673
1674 CALL timestop(handle)
1675
1676 END SUBROUTINE transform_dmat
1677
1678! **************************************************************************************************
1679!> \brief Symmetrization of density matrix - transform to new k-point
1680!> \param smat density matrix at new kpoint
1681!> \param pmat reference density matrix
1682!> \param kmat Kind type rotation matrix
1683!> \param rot Rotation matrix
1684!> \param f0 Permutation of atoms under transformation
1685!> \param atype Atom to kind pointer
1686!> \param symmetric Symmetric matrix
1687!> \param antisymmetric Anti-Symmetric matrix
1688! **************************************************************************************************
1689 SUBROUTINE symtrans(smat, pmat, kmat, rot, f0, atype, symmetric, antisymmetric)
1690 TYPE(dbcsr_type), POINTER :: smat, pmat
1691 TYPE(kind_rotmat_type), DIMENSION(:), POINTER :: kmat
1692 REAL(kind=dp), DIMENSION(3, 3), INTENT(IN) :: rot
1693 INTEGER, DIMENSION(:), INTENT(IN) :: f0, atype
1694 LOGICAL, INTENT(IN), OPTIONAL :: symmetric, antisymmetric
1695
1696 CHARACTER(LEN=*), PARAMETER :: routinen = 'symtrans'
1697
1698 INTEGER :: handle, iatom, icol, ikind, ip, irow, &
1699 jcol, jkind, jp, jrow, natom, numnodes
1700 LOGICAL :: asym, dorot, found, perm, sym, trans
1701 REAL(kind=dp) :: dr, fsign
1702 REAL(kind=dp), DIMENSION(:, :), POINTER :: kroti, krotj, pblock, sblock
1703 TYPE(dbcsr_distribution_type) :: dist
1704 TYPE(dbcsr_iterator_type) :: iter
1705
1706 CALL timeset(routinen, handle)
1707
1708 ! check symmetry options
1709 sym = .false.
1710 IF (PRESENT(symmetric)) sym = symmetric
1711 asym = .false.
1712 IF (PRESENT(antisymmetric)) asym = antisymmetric
1713
1714 cpassert(.NOT. (sym .AND. asym))
1715 cpassert((sym .OR. asym))
1716
1717 ! do we have permutation of atoms
1718 natom = SIZE(f0)
1719 perm = .false.
1720 DO iatom = 1, natom
1721 IF (f0(iatom) == iatom) cycle
1722 perm = .true.
1723 EXIT
1724 END DO
1725
1726 ! do we have a real rotation
1727 dorot = .false.
1728 IF (abs(sum(abs(rot)) - 3.0_dp) > 1.e-12_dp) dorot = .true.
1729 dr = abs(rot(1, 1) - 1.0_dp) + abs(rot(2, 2) - 1.0_dp) + abs(rot(3, 3) - 1.0_dp)
1730 IF (abs(dr) > 1.e-12_dp) dorot = .true.
1731
1732 fsign = 1.0_dp
1733 IF (asym) fsign = -1.0_dp
1734
1735 IF (dorot .OR. perm) THEN
1736 CALL cp_abort(__location__, "k-points need FULL_GRID currently. "// &
1737 "Reduced grids not yet working correctly")
1738 CALL dbcsr_set(smat, 0.0_dp)
1739 IF (perm) THEN
1740 CALL dbcsr_get_info(pmat, distribution=dist)
1741 CALL dbcsr_distribution_get(dist, numnodes=numnodes)
1742 IF (numnodes == 1) THEN
1743 ! the matrices are local to this process
1744 CALL dbcsr_iterator_start(iter, pmat)
1745 DO WHILE (dbcsr_iterator_blocks_left(iter))
1746 CALL dbcsr_iterator_next_block(iter, irow, icol, pblock)
1747 ip = f0(irow)
1748 jp = f0(icol)
1749 IF (ip <= jp) THEN
1750 jrow = ip
1751 jcol = jp
1752 trans = .false.
1753 ELSE
1754 jrow = jp
1755 jcol = ip
1756 trans = .true.
1757 END IF
1758 CALL dbcsr_get_block_p(matrix=smat, row=jrow, col=jcol, block=sblock, found=found)
1759 cpassert(found)
1760 ikind = atype(irow)
1761 jkind = atype(icol)
1762 kroti => kmat(ikind)%rmat
1763 krotj => kmat(jkind)%rmat
1764 ! rotation
1765 IF (trans) THEN
1766 sblock = fsign*matmul(transpose(krotj), matmul(transpose(pblock), kroti))
1767 ELSE
1768 sblock = fsign*matmul(transpose(kroti), matmul(pblock, krotj))
1769 END IF
1770 END DO
1771 CALL dbcsr_iterator_stop(iter)
1772 !
1773 ELSE
1774 ! distributed matrices, most general code needed
1775 CALL cp_abort(__location__, "k-points need FULL_GRID currently. "// &
1776 "Reduced grids not yet working correctly")
1777 END IF
1778 ELSE
1779 ! no atom permutations, this is always local
1780 CALL dbcsr_copy(smat, pmat)
1781 CALL dbcsr_iterator_start(iter, smat)
1782 DO WHILE (dbcsr_iterator_blocks_left(iter))
1783 CALL dbcsr_iterator_next_block(iter, irow, icol, sblock)
1784 ip = f0(irow)
1785 jp = f0(icol)
1786 IF (ip <= jp) THEN
1787 jrow = ip
1788 jcol = jp
1789 trans = .false.
1790 ELSE
1791 jrow = jp
1792 jcol = ip
1793 trans = .true.
1794 END IF
1795 ikind = atype(irow)
1796 jkind = atype(icol)
1797 kroti => kmat(ikind)%rmat
1798 krotj => kmat(jkind)%rmat
1799 ! rotation
1800 IF (trans) THEN
1801 sblock = fsign*matmul(transpose(krotj), matmul(transpose(sblock), kroti))
1802 ELSE
1803 sblock = fsign*matmul(transpose(kroti), matmul(sblock, krotj))
1804 END IF
1805 END DO
1806 CALL dbcsr_iterator_stop(iter)
1807 !
1808 END IF
1809 ELSE
1810 ! this is the identity operation, just copy the matrix
1811 CALL dbcsr_copy(smat, pmat)
1812 END IF
1813
1814 CALL timestop(handle)
1815
1816 END SUBROUTINE symtrans
1817
1818! **************************************************************************************************
1819!> \brief ...
1820!> \param mat ...
1821! **************************************************************************************************
1822 SUBROUTINE matprint(mat)
1823 TYPE(dbcsr_type), POINTER :: mat
1824
1825 INTEGER :: i, icol, iounit, irow
1826 REAL(kind=dp), DIMENSION(:, :), POINTER :: mblock
1827 TYPE(dbcsr_iterator_type) :: iter
1828
1830 CALL dbcsr_iterator_start(iter, mat)
1831 DO WHILE (dbcsr_iterator_blocks_left(iter))
1832 CALL dbcsr_iterator_next_block(iter, irow, icol, mblock)
1833 !
1834 IF (iounit > 0) THEN
1835 WRITE (iounit, '(A,2I4)') 'BLOCK ', irow, icol
1836 DO i = 1, SIZE(mblock, 1)
1837 WRITE (iounit, '(8F12.6)') mblock(i, :)
1838 END DO
1839 END IF
1840 !
1841 END DO
1842 CALL dbcsr_iterator_stop(iter)
1843
1844 END SUBROUTINE matprint
1845! **************************************************************************************************
1846
1847END MODULE kpoint_methods
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
Definition cell_types.F:491
methods related to the blacs parallel environment
subroutine, public cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
allocates and initializes a type that represent a blacs context
Represents a complex full matrix distributed on many processors.
subroutine, public cp_cfm_release(matrix)
Releases a full matrix.
subroutine, public cp_fm_to_cfm(msourcer, msourcei, mtarget)
Construct a complex full matrix by taking its real and imaginary parts from two separate real-value f...
subroutine, public cp_cfm_create(matrix, matrix_struct, name, nrow, ncol, set_zero)
Creates a new full matrix with the given structure.
subroutine, public cp_cfm_to_fm(msource, mtargetr, mtargeti)
Copy real and imaginary parts of a complex full matrix into separate real-value full matrices.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_deallocate_matrix(matrix)
...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_distribution_get(dist, row_dist, col_dist, nrows, ncols, has_threads, group, mynode, numnodes, nprows, npcols, myprow, mypcol, pgrid, subgroups_defined, prow_group, pcol_group)
...
DBCSR operations in CP2K.
subroutine, public copy_fm_to_dbcsr(fm, matrix, keep_sparsity)
Copy a BLACS matrix to a dbcsr matrix.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_column_scale(matrixa, scaling)
scales column i of matrix a with scaling(i)
pool for for elements that are retained and released
subroutine, public fm_pool_create_fm(pool, element, name)
returns an element, allocating it if none is in the pool
subroutine, public fm_pool_give_back_fm(pool, element)
returns the element to the pool
represent the structure of a full matrix
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_get_diag(matrix, diag)
returns the diagonal elements of a fm
subroutine, public cp_fm_start_copy_general(source, destination, para_env, info)
Initiates the copy operation: get distribution data, post MPI isend and irecvs.
subroutine, public cp_fm_cleanup_copy_general(info)
Completes the copy operation: wait for comms clean up MPI state.
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
subroutine, public cp_fm_finish_copy_general(destination, info)
Completes the copy operation: wait for comms, unpack, clean up MPI state.
subroutine, public cp_fm_get_submatrix(fm, target_m, start_row, start_col, n_rows, n_cols, transpose)
gets a submatrix of a full matrix op(target_m)(1:n_rows,1:n_cols) =fm(start_row:start_row+n_rows,...
various routines to log and control the output. The idea is that decisions about where to log should ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
K-points and crystal symmetry routines.
Definition cryssym.F:12
subroutine, public print_crys_symmetry(csym)
...
Definition cryssym.F:579
subroutine, public crys_sym_gen(csym, scoor, types, hmat, delta, iounit)
...
Definition cryssym.F:129
subroutine, public release_csym_type(csym)
Release the CSYM type.
Definition cryssym.F:84
subroutine, public kpoint_gen(csym, nk, symm, shift, full_grid)
...
Definition cryssym.F:224
subroutine, public print_kp_symmetry(csym)
...
Definition cryssym.F:612
subroutine, public probe_occupancy_kp(occ, fermi, kts, energies, rcoeff, icoeff, maxocc, probe, n, wk)
subroutine to calculate occupation number and 'Fermi' level using the
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public smear_fermi_dirac
integer, parameter, public smear_gaussian
integer, parameter, public smear_mv
integer, parameter, public smear_mp
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Routines needed for kpoint calculation.
subroutine, public kpoint_initialize_mo_set(kpoint)
...
subroutine, public rskp_transform(rmatrix, cmatrix, rsmat, ispin, xkp, cell_to_index, sab_nl, is_complex, rs_sign)
Transformation of real space matrices to a kpoint.
subroutine, public kpoint_init_cell_index(kpoint, sab_nl, para_env, nimages)
Generates the mapping of cell indices and linear RS index CELL (0,0,0) is always mapped to index 1.
subroutine, public kpoint_initialize_mos(kpoint, mos, added_mos, for_aux_fit)
Initialize a set of MOs and density matrix for each kpoint (kpoint group)
subroutine, public kpoint_initialize(kpoint, particle_set, cell)
Generate the kpoints and initialize the kpoint environment.
subroutine, public kpoint_density_transform(kpoint, denmat, wtype, tempmat, sab_nl, fmwork, for_aux_fit, pmat_ext)
generate real space density matrices in DBCSR format
subroutine, public kpoint_density_matrices(kpoint, energy_weighted, for_aux_fit)
Calculate kpoint density matrices (rho(k), owned by kpoint groups)
subroutine, public lowdin_kp_trans(kpoint, pmat_diag)
Calculate Lowdin transformation of density matrix S^1/2 P S^1/2 Integrate diagonal elements over k-po...
subroutine, public kpoint_env_initialize(kpoint, para_env, blacs_env, with_aux_fit)
Initialize the kpoint environment.
subroutine, public kpoint_set_mo_occupation(kpoint, smear, probe)
Given the eigenvalues of all kpoints, calculates the occupation numbers.
Types and basic routines needed for a kpoint calculation.
subroutine, public kpoint_sym_create(kp_sym)
Create a single kpoint symmetry environment.
subroutine, public kpoint_env_create(kp_env)
Create a single kpoint environment.
subroutine, public get_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, sab_nl, sab_nl_nosym)
Retrieve information from a kpoint environment.
K-points and crystal symmetry routines based on.
Definition kpsym.F:28
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition list.F:24
Definition of mathematical constants and functions.
real(kind=dp), parameter, public twopi
Utility routines for the memory handling.
Interface to the message passing library MPI.
basic linear algebra operations for full matrixes
Define the data structure for the particle information.
wrapper for the pools of matrixes
subroutine, public mpools_create(mpools)
creates a mpools
subroutine, public mpools_rebuild_fm_pools(mpools, mos, blacs_env, para_env, nmosub)
rebuilds the pools of the (ao x mo, ao x ao , mo x mo) full matrixes
subroutine, public mpools_get(mpools, ao_mo_fm_pools, ao_ao_fm_pools, mo_mo_fm_pools, ao_mosub_fm_pools, mosub_mosub_fm_pools, maxao_maxmo_fm_pool, maxao_maxao_fm_pool, maxmo_maxmo_fm_pool)
returns various attributes of the mpools (notably the pools contained in it)
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public set_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, uniform_occupation, kts, mu, flexible_electron_count)
Set the components of a MO set data structure.
subroutine, public allocate_mo_set(mo_set, nao, nmo, nelectron, n_el_f, maxocc, flexible_electron_count)
Allocates a mo set and partially initializes it (nao,nmo,nelectron, and flexible_electron_count are v...
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
subroutine, public init_mo_set(mo_set, fm_pool, fm_ref, fm_struct, name)
initializes an allocated mo_set. eigenvalues, mo_coeff, occupation_numbers are valid only after this ...
Define the neighbor list data types and the corresponding functionality.
subroutine, public neighbor_list_iterator_create(iterator_set, nl, search, nthread)
Neighbor list iterator functions.
subroutine, public neighbor_list_iterator_release(iterator_set)
...
subroutine, public get_neighbor_list_set_p(neighbor_list_sets, nlist, symmetric)
Return the components of the first neighbor list set.
integer function, public neighbor_list_iterate(iterator_set, mepos)
...
subroutine, public get_iterator_info(iterator_set, mepos, ikind, jkind, nkind, ilist, nlist, inode, nnode, iatom, jatom, r, cell)
...
parameters that control an scf iteration
Unified smearing module supporting four methods: smear_fermi_dirac — Fermi-Dirac distribution smear_g...
subroutine, public smearkp(f, mu, kts, e, nel, wk, sigma, maxocc, method)
Bisection search for mu given a target electron count (k-point case, single spin channel or spin-dege...
subroutine, public smearkp2(f, mu, kts, e, nel, wk, sigma, method)
Bisection search for mu (k-point, spin-polarised with a shared chemical potential across both spin ch...
All kind of helpful little routines.
Definition util.F:14
pure integer function, dimension(2), public get_limit(m, n, me)
divide m entries into n parts, return size of part me
Definition util.F:333
Type defining parameters related to the simulation cell.
Definition cell_types.F:60
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
Represent a complex full matrix.
keeps the information about the structure of a full matrix
Stores the state of a copy between cp_fm_start_copy_general and cp_fm_finish_copy_general.
represent a full matrix
CSM type.
Definition cryssym.F:42
Rotation matrices for basis sets.
Keeps information about a specific k-point.
Keeps symmetry information about a specific k-point.
Contains information about kpoints.
stores all the informations relevant to an mpi environment
container for the pools of matrixes used by qs
contains the parameters needed by a scf run