(git:1f285aa)
kpoint_methods.F
Go to the documentation of this file.
1 !--------------------------------------------------------------------------------------------------!
2 ! CP2K: A general program to perform molecular dynamics simulations !
3 ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 ! !
5 ! SPDX-License-Identifier: GPL-2.0-or-later !
6 !--------------------------------------------------------------------------------------------------!
7 
8 ! **************************************************************************************************
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,&
20  cp_blacs_env_type
21  USE cp_control_types, ONLY: dft_control_type
25  USE cp_fm_pool_types, ONLY: cp_fm_pool_p_type,&
28  USE cp_fm_struct, ONLY: cp_fm_struct_type
29  USE cp_fm_types, ONLY: &
31  cp_fm_get_info, cp_fm_release, cp_fm_start_copy_general, cp_fm_to_fm, cp_fm_type
33  USE cryssym, ONLY: apply_rotation_coord,&
34  crys_sym_gen,&
35  csym_type,&
36  kpoint_gen,&
40  USE dbcsr_api, ONLY: &
41  dbcsr_copy, dbcsr_create, dbcsr_deallocate_matrix, dbcsr_distribution_get, &
42  dbcsr_distribution_type, dbcsr_get_block_p, dbcsr_get_info, dbcsr_iterator_blocks_left, &
43  dbcsr_iterator_next_block, dbcsr_iterator_start, dbcsr_iterator_stop, dbcsr_iterator_type, &
44  dbcsr_p_type, dbcsr_set, dbcsr_type, dbcsr_type_antisymmetric, dbcsr_type_no_symmetry, &
45  dbcsr_type_symmetric
46  USE fermi_utils, ONLY: fermikp,&
47  fermikp2
49  USE kinds, ONLY: dp
50  USE kpoint_types, ONLY: get_kpoint_info,&
52  kpoint_env_p_type,&
53  kpoint_env_type,&
55  kpoint_sym_type,&
56  kpoint_type
57  USE mathconstants, ONLY: twopi
58  USE memory_utilities, ONLY: reallocate
59  USE message_passing, ONLY: mp_cart_type,&
60  mp_para_env_type
61  USE parallel_gemm_api, ONLY: parallel_gemm
62  USE particle_types, ONLY: particle_type
63  USE qs_matrix_pools, ONLY: mpools_create,&
64  mpools_get,&
66  qs_matrix_pools_type
67  USE qs_mo_types, ONLY: allocate_mo_set,&
68  get_mo_set,&
69  init_mo_set,&
70  mo_set_type,&
76  neighbor_list_iterator_p_type,&
78  neighbor_list_set_p_type
79  USE scf_control_types, ONLY: smear_type
80  USE util, ONLY: get_limit
81 #include "./base/base_uses.f90"
82 
83  IMPLICIT NONE
84 
85  PRIVATE
86 
87  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'kpoint_methods'
88 
92  PUBLIC :: rskp_transform
93 
94 ! **************************************************************************************************
95 
96 CONTAINS
97 
98 ! **************************************************************************************************
99 !> \brief Generate the kpoints and initialize the kpoint environment
100 !> \param kpoint The kpoint environment
101 !> \param particle_set Particle types and coordinates
102 !> \param cell Computational cell information
103 ! **************************************************************************************************
104  SUBROUTINE kpoint_initialize(kpoint, particle_set, cell)
105 
106  TYPE(kpoint_type), POINTER :: kpoint
107  TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
108  TYPE(cell_type), POINTER :: cell
109 
110  CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_initialize'
111 
112  INTEGER :: handle, i, ik, iounit, ir, is, natom, &
113  nr, ns
114  INTEGER, ALLOCATABLE, DIMENSION(:) :: atype
115  LOGICAL :: spez
116  REAL(kind=dp) :: wsum
117  REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: coord, scoord
118  TYPE(csym_type) :: crys_sym
119  TYPE(kpoint_sym_type), POINTER :: kpsym
120 
121  CALL timeset(routinen, handle)
122 
123  cpassert(ASSOCIATED(kpoint))
124 
125  SELECT CASE (kpoint%kp_scheme)
126  CASE ("NONE")
127  ! do nothing
128  CASE ("GAMMA")
129  kpoint%nkp = 1
130  ALLOCATE (kpoint%xkp(3, 1), kpoint%wkp(1))
131  kpoint%xkp(1:3, 1) = 0.0_dp
132  kpoint%wkp(1) = 1.0_dp
133  ALLOCATE (kpoint%kp_sym(1))
134  NULLIFY (kpoint%kp_sym(1)%kpoint_sym)
135  CALL kpoint_sym_create(kpoint%kp_sym(1)%kpoint_sym)
136  CASE ("MONKHORST-PACK", "MACDONALD")
137 
138  IF (.NOT. kpoint%symmetry) THEN
139  ! we set up a random molecule to avoid any possible symmetry
140  natom = 10
141  ALLOCATE (coord(3, natom), scoord(3, natom), atype(natom))
142  DO i = 1, natom
143  atype(i) = i
144  coord(1, i) = sin(i*0.12345_dp)
145  coord(2, i) = cos(i*0.23456_dp)
146  coord(3, i) = sin(i*0.34567_dp)
147  CALL real_to_scaled(scoord(1:3, i), coord(1:3, i), cell)
148  END DO
149  ELSE
150  natom = SIZE(particle_set)
151  ALLOCATE (scoord(3, natom), atype(natom))
152  DO i = 1, natom
153  CALL get_atomic_kind(atomic_kind=particle_set(i)%atomic_kind, kind_number=atype(i))
154  CALL real_to_scaled(scoord(1:3, i), particle_set(i)%r(1:3), cell)
155  END DO
156  END IF
157  IF (kpoint%verbose) THEN
159  ELSE
160  iounit = -1
161  END IF
162 
163  CALL crys_sym_gen(crys_sym, scoord, atype, cell%hmat, delta=kpoint%eps_geo, iounit=iounit)
164  CALL kpoint_gen(crys_sym, kpoint%nkp_grid, symm=kpoint%symmetry, shift=kpoint%kp_shift, &
165  full_grid=kpoint%full_grid)
166  kpoint%nkp = crys_sym%nkpoint
167  ALLOCATE (kpoint%xkp(3, kpoint%nkp), kpoint%wkp(kpoint%nkp))
168  wsum = sum(crys_sym%wkpoint)
169  DO ik = 1, kpoint%nkp
170  kpoint%xkp(1:3, ik) = crys_sym%xkpoint(1:3, ik)
171  kpoint%wkp(ik) = crys_sym%wkpoint(ik)/wsum
172  END DO
173 
174  ! print output
175  IF (kpoint%symmetry) CALL print_crys_symmetry(crys_sym)
176  IF (kpoint%symmetry) CALL print_kp_symmetry(crys_sym)
177 
178  ! transfer symmetry information
179  ALLOCATE (kpoint%kp_sym(kpoint%nkp))
180  DO ik = 1, kpoint%nkp
181  NULLIFY (kpoint%kp_sym(ik)%kpoint_sym)
182  CALL kpoint_sym_create(kpoint%kp_sym(ik)%kpoint_sym)
183  kpsym => kpoint%kp_sym(ik)%kpoint_sym
184  IF (crys_sym%symlib .AND. .NOT. crys_sym%fullgrid .AND. crys_sym%istriz == 1) THEN
185  ! set up the symmetrization information
186  kpsym%nwght = nint(crys_sym%wkpoint(ik))
187  ns = kpsym%nwght
188  ! to be done correctly
189  IF (ns > 1) THEN
190  kpsym%apply_symmetry = .true.
191  natom = SIZE(particle_set)
192  ALLOCATE (kpsym%rot(3, 3, ns))
193  ALLOCATE (kpsym%xkp(3, ns))
194  ALLOCATE (kpsym%f0(natom, ns))
195  nr = 0
196  DO is = 1, SIZE(crys_sym%kplink, 2)
197  IF (crys_sym%kplink(2, is) == ik) THEN
198  nr = nr + 1
199  ir = crys_sym%kpop(is)
200  kpsym%rot(1:3, 1:3, nr) = crys_sym%rotations(1:3, 1:3, ir)
201  kpsym%xkp(1:3, nr) = crys_sym%kpmesh(1:3, is)
202  CALL apply_rotation_coord(kpsym%f0(1:natom, nr), crys_sym, ir)
203  END IF
204  END DO
205  END IF
206  END IF
207  END DO
208 
209  CALL release_csym_type(crys_sym)
210  DEALLOCATE (scoord, atype)
211 
212  CASE ("GENERAL")
213  ! default: no symmetry settings
214  ALLOCATE (kpoint%kp_sym(kpoint%nkp))
215  DO i = 1, kpoint%nkp
216  NULLIFY (kpoint%kp_sym(i)%kpoint_sym)
217  CALL kpoint_sym_create(kpoint%kp_sym(i)%kpoint_sym)
218  END DO
219  CASE DEFAULT
220  cpassert(.false.)
221  END SELECT
222 
223  ! check for consistency of options
224  SELECT CASE (kpoint%kp_scheme)
225  CASE ("NONE")
226  ! don't use k-point code
227  CASE ("GAMMA")
228  cpassert(kpoint%nkp == 1)
229  cpassert(sum(abs(kpoint%xkp)) <= 1.e-12_dp)
230  cpassert(kpoint%wkp(1) == 1.0_dp)
231  cpassert(.NOT. kpoint%symmetry)
232  CASE ("GENERAL")
233  cpassert(.NOT. kpoint%symmetry)
234  cpassert(kpoint%nkp >= 1)
235  CASE ("MONKHORST-PACK", "MACDONALD")
236  cpassert(kpoint%nkp >= 1)
237  END SELECT
238  IF (kpoint%use_real_wfn) THEN
239  ! what about inversion symmetry?
240  ikloop: DO ik = 1, kpoint%nkp
241  DO i = 1, 3
242  spez = (kpoint%xkp(i, ik) == 0.0_dp .OR. kpoint%xkp(i, ik) == 0.5_dp)
243  IF (.NOT. spez) EXIT ikloop
244  END DO
245  END DO ikloop
246  IF (.NOT. spez) THEN
247  ! Warning: real wfn might be wrong for this system
248  CALL cp_warn(__location__, &
249  "A calculation using real wavefunctions is requested. "// &
250  "We could not determine if the symmetry of the system allows real wavefunctions. ")
251  END IF
252  END IF
253 
254  CALL timestop(handle)
255 
256  END SUBROUTINE kpoint_initialize
257 
258 ! **************************************************************************************************
259 !> \brief Initialize the kpoint environment
260 !> \param kpoint Kpoint environment
261 !> \param para_env ...
262 !> \param blacs_env ...
263 !> \param with_aux_fit ...
264 ! **************************************************************************************************
265  SUBROUTINE kpoint_env_initialize(kpoint, para_env, blacs_env, with_aux_fit)
266 
267  TYPE(kpoint_type), INTENT(INOUT) :: kpoint
268  TYPE(mp_para_env_type), INTENT(IN), TARGET :: para_env
269  TYPE(cp_blacs_env_type), INTENT(IN), TARGET :: blacs_env
270  LOGICAL, INTENT(IN), OPTIONAL :: with_aux_fit
271 
272  CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_env_initialize'
273 
274  INTEGER :: handle, igr, ik, ikk, ngr, niogrp, nkp, &
275  nkp_grp, nkp_loc, npe, unit_nr
276  INTEGER, DIMENSION(2) :: dims, pos
277  LOGICAL :: aux_fit
278  TYPE(kpoint_env_p_type), DIMENSION(:), POINTER :: kp_aux_env, kp_env
279  TYPE(kpoint_env_type), POINTER :: kp
280  TYPE(mp_cart_type) :: comm_cart
281  TYPE(mp_para_env_type), POINTER :: para_env_inter_kp, para_env_kp
282 
283  CALL timeset(routinen, handle)
284 
285  IF (PRESENT(with_aux_fit)) THEN
286  aux_fit = with_aux_fit
287  ELSE
288  aux_fit = .false.
289  END IF
290 
291  kpoint%para_env => para_env
292  CALL kpoint%para_env%retain()
293  kpoint%blacs_env_all => blacs_env
294  CALL kpoint%blacs_env_all%retain()
295 
296  cpassert(.NOT. ASSOCIATED(kpoint%kp_env))
297  IF (aux_fit) THEN
298  cpassert(.NOT. ASSOCIATED(kpoint%kp_aux_env))
299  END IF
300 
301  NULLIFY (kp_env, kp_aux_env)
302  nkp = kpoint%nkp
303  npe = para_env%num_pe
304  IF (npe == 1) THEN
305  ! only one process available -> owns all kpoints
306  ALLOCATE (kp_env(nkp))
307  DO ik = 1, nkp
308  NULLIFY (kp_env(ik)%kpoint_env)
309  CALL kpoint_env_create(kp_env(ik)%kpoint_env)
310  kp => kp_env(ik)%kpoint_env
311  kp%nkpoint = ik
312  kp%wkp = kpoint%wkp(ik)
313  kp%xkp(1:3) = kpoint%xkp(1:3, ik)
314  kp%is_local = .true.
315  END DO
316  kpoint%kp_env => kp_env
317 
318  IF (aux_fit) THEN
319  ALLOCATE (kp_aux_env(nkp))
320  DO ik = 1, nkp
321  NULLIFY (kp_aux_env(ik)%kpoint_env)
322  CALL kpoint_env_create(kp_aux_env(ik)%kpoint_env)
323  kp => kp_aux_env(ik)%kpoint_env
324  kp%nkpoint = ik
325  kp%wkp = kpoint%wkp(ik)
326  kp%xkp(1:3) = kpoint%xkp(1:3, ik)
327  kp%is_local = .true.
328  END DO
329 
330  kpoint%kp_aux_env => kp_aux_env
331  END IF
332 
333  ALLOCATE (kpoint%kp_dist(2, 1))
334  kpoint%kp_dist(1, 1) = 1
335  kpoint%kp_dist(2, 1) = nkp
336  kpoint%kp_range(1) = 1
337  kpoint%kp_range(2) = nkp
338 
339  ! parallel environments
340  kpoint%para_env_kp => para_env
341  CALL kpoint%para_env_kp%retain()
342  kpoint%para_env_inter_kp => para_env
343  CALL kpoint%para_env_inter_kp%retain()
344  kpoint%iogrp = .true.
345  kpoint%nkp_groups = 1
346  ELSE
347  IF (kpoint%parallel_group_size == -1) THEN
348  ! maximum parallelization over kpoints
349  ! making sure that the group size divides the npe and the nkp_grp the nkp
350  ! in the worst case, there will be no parallelism over kpoints.
351  DO igr = npe, 1, -1
352  IF (mod(npe, igr) .NE. 0) cycle
353  nkp_grp = npe/igr
354  IF (mod(nkp, nkp_grp) .NE. 0) cycle
355  ngr = igr
356  END DO
357  ELSE IF (kpoint%parallel_group_size == 0) THEN
358  ! no parallelization over kpoints
359  ngr = npe
360  ELSE IF (kpoint%parallel_group_size > 0) THEN
361  ngr = min(kpoint%parallel_group_size, npe)
362  ELSE
363  cpassert(.false.)
364  END IF
365  nkp_grp = npe/ngr
366  ! processor dimensions
367  dims(1) = ngr
368  dims(2) = nkp_grp
369  cpassert(mod(nkp, nkp_grp) == 0)
370  nkp_loc = nkp/nkp_grp
371 
372  IF ((dims(1)*dims(2) /= npe)) THEN
373  cpabort("Number of processors is not divisible by the kpoint group size.")
374  END IF
375 
376  ! Create the subgroups, one for each k-point group and one interconnecting group
377  CALL comm_cart%create(comm_old=para_env, ndims=2, dims=dims)
378  pos = comm_cart%mepos_cart
379  ALLOCATE (para_env_kp)
380  CALL para_env_kp%from_split(comm_cart, pos(2))
381  ALLOCATE (para_env_inter_kp)
382  CALL para_env_inter_kp%from_split(comm_cart, pos(1))
383  CALL comm_cart%free()
384 
385  niogrp = 0
386  IF (para_env%is_source()) niogrp = 1
387  CALL para_env_kp%sum(niogrp)
388  kpoint%iogrp = (niogrp == 1)
389 
390  ! parallel groups
391  kpoint%para_env_kp => para_env_kp
392  kpoint%para_env_inter_kp => para_env_inter_kp
393 
394  ! distribution of kpoints
395  ALLOCATE (kpoint%kp_dist(2, nkp_grp))
396  DO igr = 1, nkp_grp
397  kpoint%kp_dist(1:2, igr) = get_limit(nkp, nkp_grp, igr - 1)
398  END DO
399  ! local kpoints
400  kpoint%kp_range(1:2) = kpoint%kp_dist(1:2, para_env_inter_kp%mepos + 1)
401 
402  ALLOCATE (kp_env(nkp_loc))
403  DO ik = 1, nkp_loc
404  NULLIFY (kp_env(ik)%kpoint_env)
405  ikk = kpoint%kp_range(1) + ik - 1
406  CALL kpoint_env_create(kp_env(ik)%kpoint_env)
407  kp => kp_env(ik)%kpoint_env
408  kp%nkpoint = ikk
409  kp%wkp = kpoint%wkp(ikk)
410  kp%xkp(1:3) = kpoint%xkp(1:3, ikk)
411  kp%is_local = (ngr == 1)
412  END DO
413  kpoint%kp_env => kp_env
414 
415  IF (aux_fit) THEN
416  ALLOCATE (kp_aux_env(nkp_loc))
417  DO ik = 1, nkp_loc
418  NULLIFY (kp_aux_env(ik)%kpoint_env)
419  ikk = kpoint%kp_range(1) + ik - 1
420  CALL kpoint_env_create(kp_aux_env(ik)%kpoint_env)
421  kp => kp_aux_env(ik)%kpoint_env
422  kp%nkpoint = ikk
423  kp%wkp = kpoint%wkp(ikk)
424  kp%xkp(1:3) = kpoint%xkp(1:3, ikk)
425  kp%is_local = (ngr == 1)
426  END DO
427  kpoint%kp_aux_env => kp_aux_env
428  END IF
429 
431 
432  IF (unit_nr > 0 .AND. kpoint%verbose) THEN
433  WRITE (unit_nr, *)
434  WRITE (unit_nr, fmt="(T2,A,T71,I10)") "KPOINTS| Number of kpoint groups ", nkp_grp
435  WRITE (unit_nr, fmt="(T2,A,T71,I10)") "KPOINTS| Size of each kpoint group", ngr
436  WRITE (unit_nr, fmt="(T2,A,T71,I10)") "KPOINTS| Number of kpoints per group", nkp_loc
437  END IF
438  kpoint%nkp_groups = nkp_grp
439 
440  END IF
441 
442  CALL timestop(handle)
443 
444  END SUBROUTINE kpoint_env_initialize
445 
446 ! **************************************************************************************************
447 !> \brief Initialize a set of MOs and density matrix for each kpoint (kpoint group)
448 !> \param kpoint Kpoint environment
449 !> \param mos Reference MOs (global)
450 !> \param added_mos ...
451 !> \param for_aux_fit ...
452 ! **************************************************************************************************
453  SUBROUTINE kpoint_initialize_mos(kpoint, mos, added_mos, for_aux_fit)
454 
455  TYPE(kpoint_type), POINTER :: kpoint
456  TYPE(mo_set_type), DIMENSION(:), INTENT(INOUT) :: mos
457  INTEGER, INTENT(IN), OPTIONAL :: added_mos
458  LOGICAL, OPTIONAL :: for_aux_fit
459 
460  CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_initialize_mos'
461 
462  INTEGER :: handle, ic, ik, is, nadd, nao, nc, &
463  nelectron, nkp_loc, nmo, nmorig(2), &
464  nspin
465  LOGICAL :: aux_fit
466  REAL(kind=dp) :: flexible_electron_count, maxocc, n_el_f
467  TYPE(cp_blacs_env_type), POINTER :: blacs_env
468  TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_ao_fm_pools
469  TYPE(cp_fm_struct_type), POINTER :: matrix_struct
470  TYPE(cp_fm_type), POINTER :: fmlocal
471  TYPE(kpoint_env_type), POINTER :: kp
472  TYPE(qs_matrix_pools_type), POINTER :: mpools
473 
474  CALL timeset(routinen, handle)
475 
476  IF (PRESENT(for_aux_fit)) THEN
477  aux_fit = for_aux_fit
478  ELSE
479  aux_fit = .false.
480  END IF
481 
482  cpassert(ASSOCIATED(kpoint))
483 
484  IF (.true. .OR. ASSOCIATED(mos(1)%mo_coeff)) THEN
485  IF (aux_fit) THEN
486  cpassert(ASSOCIATED(kpoint%kp_aux_env))
487  END IF
488 
489  IF (PRESENT(added_mos)) THEN
490  nadd = added_mos
491  ELSE
492  nadd = 0
493  END IF
494 
495  IF (kpoint%use_real_wfn) THEN
496  nc = 1
497  ELSE
498  nc = 2
499  END IF
500  nspin = SIZE(mos, 1)
501  nkp_loc = kpoint%kp_range(2) - kpoint%kp_range(1) + 1
502  IF (nkp_loc > 0) THEN
503  IF (aux_fit) THEN
504  cpassert(SIZE(kpoint%kp_aux_env) == nkp_loc)
505  ELSE
506  cpassert(SIZE(kpoint%kp_env) == nkp_loc)
507  END IF
508  ! allocate the mo sets, correct number of kpoints (local), real/complex, spin
509  DO ik = 1, nkp_loc
510  IF (aux_fit) THEN
511  kp => kpoint%kp_aux_env(ik)%kpoint_env
512  ELSE
513  kp => kpoint%kp_env(ik)%kpoint_env
514  END IF
515  ALLOCATE (kp%mos(nc, nspin))
516  DO is = 1, nspin
517  CALL get_mo_set(mos(is), nao=nao, nmo=nmo, nelectron=nelectron, &
518  n_el_f=n_el_f, maxocc=maxocc, flexible_electron_count=flexible_electron_count)
519  nmo = min(nao, nmo + nadd)
520  DO ic = 1, nc
521  CALL allocate_mo_set(kp%mos(ic, is), nao, nmo, nelectron, n_el_f, maxocc, &
522  flexible_electron_count)
523  END DO
524  END DO
525  END DO
526 
527  ! generate the blacs environment for the kpoint group
528  ! we generate a blacs env for each kpoint group in parallel
529  ! we assume here that the group para_env_inter_kp will connect
530  ! equivalent parts of fm matrices, i.e. no reshuffeling of processors
531  NULLIFY (blacs_env)
532  IF (ASSOCIATED(kpoint%blacs_env)) THEN
533  blacs_env => kpoint%blacs_env
534  ELSE
535  CALL cp_blacs_env_create(blacs_env=blacs_env, para_env=kpoint%para_env_kp)
536  kpoint%blacs_env => blacs_env
537  END IF
538 
539  ! set possible new number of MOs
540  DO is = 1, nspin
541  CALL get_mo_set(mos(is), nmo=nmorig(is))
542  nmo = min(nao, nmorig(is) + nadd)
543  CALL set_mo_set(mos(is), nmo=nmo)
544  END DO
545  ! matrix pools for the kpoint group, information on MOs is transferred using
546  ! generic mos structure
547  NULLIFY (mpools)
548  CALL mpools_create(mpools=mpools)
549  CALL mpools_rebuild_fm_pools(mpools=mpools, mos=mos, &
550  blacs_env=blacs_env, para_env=kpoint%para_env_kp)
551 
552  IF (aux_fit) THEN
553  kpoint%mpools_aux_fit => mpools
554  ELSE
555  kpoint%mpools => mpools
556  END IF
557 
558  ! reset old number of MOs
559  DO is = 1, nspin
560  CALL set_mo_set(mos(is), nmo=nmorig(is))
561  END DO
562 
563  ! allocate density matrices
564  CALL mpools_get(mpools, ao_ao_fm_pools=ao_ao_fm_pools)
565  ALLOCATE (fmlocal)
566  CALL fm_pool_create_fm(ao_ao_fm_pools(1)%pool, fmlocal)
567  CALL cp_fm_get_info(fmlocal, matrix_struct=matrix_struct)
568  DO ik = 1, nkp_loc
569  IF (aux_fit) THEN
570  kp => kpoint%kp_aux_env(ik)%kpoint_env
571  ELSE
572  kp => kpoint%kp_env(ik)%kpoint_env
573  END IF
574  ! density matrix
575  CALL cp_fm_release(kp%pmat)
576  ALLOCATE (kp%pmat(nc, nspin))
577  DO is = 1, nspin
578  DO ic = 1, nc
579  CALL cp_fm_create(kp%pmat(ic, is), matrix_struct)
580  END DO
581  END DO
582  ! energy weighted density matrix
583  CALL cp_fm_release(kp%wmat)
584  ALLOCATE (kp%wmat(nc, nspin))
585  DO is = 1, nspin
586  DO ic = 1, nc
587  CALL cp_fm_create(kp%wmat(ic, is), matrix_struct)
588  END DO
589  END DO
590  END DO
591  CALL fm_pool_give_back_fm(ao_ao_fm_pools(1)%pool, fmlocal)
592  DEALLOCATE (fmlocal)
593 
594  END IF
595 
596  END IF
597 
598  CALL timestop(handle)
599 
600  END SUBROUTINE kpoint_initialize_mos
601 
602 ! **************************************************************************************************
603 !> \brief ...
604 !> \param kpoint ...
605 ! **************************************************************************************************
606  SUBROUTINE kpoint_initialize_mo_set(kpoint)
607  TYPE(kpoint_type), POINTER :: kpoint
608 
609  CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_initialize_mo_set'
610 
611  INTEGER :: handle, ic, ik, ikk, ispin
612  TYPE(cp_fm_pool_p_type), DIMENSION(:), POINTER :: ao_mo_fm_pools
613  TYPE(cp_fm_type), POINTER :: mo_coeff
614  TYPE(mo_set_type), DIMENSION(:, :), POINTER :: moskp
615 
616  CALL timeset(routinen, handle)
617 
618  DO ik = 1, SIZE(kpoint%kp_env)
619  CALL mpools_get(kpoint%mpools, ao_mo_fm_pools=ao_mo_fm_pools)
620  moskp => kpoint%kp_env(ik)%kpoint_env%mos
621  ikk = kpoint%kp_range(1) + ik - 1
622  cpassert(ASSOCIATED(moskp))
623  DO ispin = 1, SIZE(moskp, 2)
624  DO ic = 1, SIZE(moskp, 1)
625  CALL get_mo_set(moskp(ic, ispin), mo_coeff=mo_coeff)
626  IF (.NOT. ASSOCIATED(mo_coeff)) THEN
627  CALL init_mo_set(moskp(ic, ispin), &
628  fm_pool=ao_mo_fm_pools(ispin)%pool, name="kpoints")
629  END IF
630  END DO
631  END DO
632  END DO
633 
634  CALL timestop(handle)
635 
636  END SUBROUTINE kpoint_initialize_mo_set
637 
638 ! **************************************************************************************************
639 !> \brief Generates the mapping of cell indices and linear RS index
640 !> CELL (0,0,0) is always mapped to index 1
641 !> \param kpoint Kpoint environment
642 !> \param sab_nl Defining neighbour list
643 !> \param para_env Parallel environment
644 !> \param dft_control ...
645 ! **************************************************************************************************
646  SUBROUTINE kpoint_init_cell_index(kpoint, sab_nl, para_env, dft_control)
647 
648  TYPE(kpoint_type), POINTER :: kpoint
649  TYPE(neighbor_list_set_p_type), DIMENSION(:), &
650  POINTER :: sab_nl
651  TYPE(mp_para_env_type), POINTER :: para_env
652  TYPE(dft_control_type), POINTER :: dft_control
653 
654  CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_init_cell_index'
655 
656  INTEGER :: handle, i1, i2, i3, ic, icount, it, &
657  ncount
658  INTEGER, DIMENSION(3) :: cell, itm
659  INTEGER, DIMENSION(:, :), POINTER :: index_to_cell, list
660  INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index, cti
661  LOGICAL :: new
662  TYPE(neighbor_list_iterator_p_type), &
663  DIMENSION(:), POINTER :: nl_iterator
664 
665  NULLIFY (cell_to_index, index_to_cell)
666 
667  CALL timeset(routinen, handle)
668 
669  cpassert(ASSOCIATED(kpoint))
670 
671  ALLOCATE (list(3, 125))
672  list = 0
673  icount = 1
674 
675  CALL neighbor_list_iterator_create(nl_iterator, sab_nl)
676  DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
677  CALL get_iterator_info(nl_iterator, cell=cell)
678 
679  new = .true.
680  DO ic = 1, icount
681  IF (cell(1) == list(1, ic) .AND. cell(2) == list(2, ic) .AND. &
682  cell(3) == list(3, ic)) THEN
683  new = .false.
684  EXIT
685  END IF
686  END DO
687  IF (new) THEN
688  icount = icount + 1
689  IF (icount > SIZE(list, 2)) THEN
690  CALL reallocate(list, 1, 3, 1, 2*SIZE(list, 2))
691  END IF
692  list(1:3, icount) = cell(1:3)
693  END IF
694 
695  END DO
696  CALL neighbor_list_iterator_release(nl_iterator)
697 
698  itm(1) = maxval(abs(list(1, 1:icount)))
699  itm(2) = maxval(abs(list(2, 1:icount)))
700  itm(3) = maxval(abs(list(3, 1:icount)))
701  CALL para_env%max(itm)
702  it = maxval(itm(1:3))
703  IF (ASSOCIATED(kpoint%cell_to_index)) THEN
704  DEALLOCATE (kpoint%cell_to_index)
705  END IF
706  ALLOCATE (kpoint%cell_to_index(-itm(1):itm(1), -itm(2):itm(2), -itm(3):itm(3)))
707  cell_to_index => kpoint%cell_to_index
708  cti => cell_to_index
709  cti(:, :, :) = 0
710  DO ic = 1, icount
711  i1 = list(1, ic)
712  i2 = list(2, ic)
713  i3 = list(3, ic)
714  cti(i1, i2, i3) = ic
715  END DO
716  CALL para_env%sum(cti)
717  ncount = 0
718  DO i1 = -itm(1), itm(1)
719  DO i2 = -itm(2), itm(2)
720  DO i3 = -itm(3), itm(3)
721  IF (cti(i1, i2, i3) == 0) THEN
722  cti(i1, i2, i3) = 1000000
723  ELSE
724  ncount = ncount + 1
725  cti(i1, i2, i3) = (abs(i1) + abs(i2) + abs(i3))*1000 + abs(i3)*100 + abs(i2)*10 + abs(i1)
726  cti(i1, i2, i3) = cti(i1, i2, i3) + (i1 + i2 + i3)
727  END IF
728  END DO
729  END DO
730  END DO
731 
732  IF (ASSOCIATED(kpoint%index_to_cell)) THEN
733  DEALLOCATE (kpoint%index_to_cell)
734  END IF
735  ALLOCATE (kpoint%index_to_cell(3, ncount))
736  index_to_cell => kpoint%index_to_cell
737  DO ic = 1, ncount
738  cell = minloc(cti)
739  i1 = cell(1) - 1 - itm(1)
740  i2 = cell(2) - 1 - itm(2)
741  i3 = cell(3) - 1 - itm(3)
742  cti(i1, i2, i3) = 1000000
743  index_to_cell(1, ic) = i1
744  index_to_cell(2, ic) = i2
745  index_to_cell(3, ic) = i3
746  END DO
747  cti(:, :, :) = 0
748  DO ic = 1, ncount
749  i1 = index_to_cell(1, ic)
750  i2 = index_to_cell(2, ic)
751  i3 = index_to_cell(3, ic)
752  cti(i1, i2, i3) = ic
753  END DO
754 
755  ! keep pointer to this neighborlist
756  kpoint%sab_nl => sab_nl
757 
758  ! set number of images
759  dft_control%nimages = SIZE(index_to_cell, 2)
760 
761  DEALLOCATE (list)
762 
763  CALL timestop(handle)
764  END SUBROUTINE kpoint_init_cell_index
765 
766 ! **************************************************************************************************
767 !> \brief Transformation of real space matrices to a kpoint
768 !> \param rmatrix Real part of kpoint matrix
769 !> \param cmatrix Complex part of kpoint matrix (optional)
770 !> \param rsmat Real space matrices
771 !> \param ispin Spin index
772 !> \param xkp Kpoint coordinates
773 !> \param cell_to_index mapping of cell indices to RS index
774 !> \param sab_nl Defining neighbor list
775 !> \param is_complex Matrix to be transformed is imaginary
776 !> \param rs_sign Matrix to be transformed is csaled by rs_sign
777 ! **************************************************************************************************
778  SUBROUTINE rskp_transform(rmatrix, cmatrix, rsmat, ispin, &
779  xkp, cell_to_index, sab_nl, is_complex, rs_sign)
780 
781  TYPE(dbcsr_type) :: rmatrix
782  TYPE(dbcsr_type), OPTIONAL :: cmatrix
783  TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: rsmat
784  INTEGER, INTENT(IN) :: ispin
785  REAL(kind=dp), DIMENSION(3), INTENT(IN) :: xkp
786  INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
787  TYPE(neighbor_list_set_p_type), DIMENSION(:), &
788  POINTER :: sab_nl
789  LOGICAL, INTENT(IN), OPTIONAL :: is_complex
790  REAL(kind=dp), INTENT(IN), OPTIONAL :: rs_sign
791 
792  CHARACTER(LEN=*), PARAMETER :: routinen = 'rskp_transform'
793 
794  INTEGER :: handle, iatom, ic, icol, irow, jatom, &
795  nimg
796  INTEGER, DIMENSION(3) :: cell
797  LOGICAL :: do_symmetric, found, my_complex, &
798  wfn_real_only
799  REAL(kind=dp) :: arg, coskl, fsign, fsym, sinkl
800  REAL(kind=dp), DIMENSION(:, :), POINTER :: cblock, rblock, rsblock
801  TYPE(neighbor_list_iterator_p_type), &
802  DIMENSION(:), POINTER :: nl_iterator
803 
804  CALL timeset(routinen, handle)
805 
806  my_complex = .false.
807  IF (PRESENT(is_complex)) my_complex = is_complex
808 
809  fsign = 1.0_dp
810  IF (PRESENT(rs_sign)) fsign = rs_sign
811 
812  wfn_real_only = .true.
813  IF (PRESENT(cmatrix)) wfn_real_only = .false.
814 
815  nimg = SIZE(rsmat, 2)
816 
817  CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
818 
819  CALL neighbor_list_iterator_create(nl_iterator, sab_nl)
820  DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
821  CALL get_iterator_info(nl_iterator, iatom=iatom, jatom=jatom, cell=cell)
822 
823  ! fsym = +- 1 is due to real space matrices being non-symmetric (although in a symmtric type)
824  ! with the link S_mu^0,nu^b = S_nu^0,mu^-b, and the KP matrices beeing Hermitian
825  fsym = 1.0_dp
826  irow = iatom
827  icol = jatom
828  IF (do_symmetric .AND. (iatom > jatom)) THEN
829  irow = jatom
830  icol = iatom
831  fsym = -1.0_dp
832  END IF
833 
834  ic = cell_to_index(cell(1), cell(2), cell(3))
835  IF (ic < 1 .OR. ic > nimg) cycle
836 
837  arg = real(cell(1), dp)*xkp(1) + real(cell(2), dp)*xkp(2) + real(cell(3), dp)*xkp(3)
838  IF (my_complex) THEN
839  coskl = fsign*fsym*cos(twopi*arg)
840  sinkl = fsign*sin(twopi*arg)
841  ELSE
842  coskl = fsign*cos(twopi*arg)
843  sinkl = fsign*fsym*sin(twopi*arg)
844  END IF
845 
846  CALL dbcsr_get_block_p(matrix=rsmat(ispin, ic)%matrix, row=irow, col=icol, &
847  block=rsblock, found=found)
848  IF (.NOT. found) cycle
849 
850  IF (wfn_real_only) THEN
851  CALL dbcsr_get_block_p(matrix=rmatrix, row=irow, col=icol, &
852  block=rblock, found=found)
853  IF (.NOT. found) cycle
854  rblock = rblock + coskl*rsblock
855  ELSE
856  CALL dbcsr_get_block_p(matrix=rmatrix, row=irow, col=icol, &
857  block=rblock, found=found)
858  IF (.NOT. found) cycle
859  CALL dbcsr_get_block_p(matrix=cmatrix, row=irow, col=icol, &
860  block=cblock, found=found)
861  IF (.NOT. found) cycle
862  rblock = rblock + coskl*rsblock
863  cblock = cblock + sinkl*rsblock
864  END IF
865 
866  END DO
867  CALL neighbor_list_iterator_release(nl_iterator)
868 
869  CALL timestop(handle)
870 
871  END SUBROUTINE rskp_transform
872 
873 ! **************************************************************************************************
874 !> \brief Given the eigenvalues of all kpoints, calculates the occupation numbers
875 !> \param kpoint Kpoint environment
876 !> \param smear Smearing information
877 ! **************************************************************************************************
878  SUBROUTINE kpoint_set_mo_occupation(kpoint, smear)
879 
880  TYPE(kpoint_type), POINTER :: kpoint
881  TYPE(smear_type), POINTER :: smear
882 
883  CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_set_mo_occupation'
884 
885  INTEGER :: handle, ik, ikpgr, ispin, kplocal, nb, &
886  ne_a, ne_b, nelectron, nkp, nmo, nspin
887  INTEGER, DIMENSION(2) :: kp_range
888  REAL(kind=dp) :: kts, mu, mus(2), nel
889  REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: weig, wocc
890  REAL(kind=dp), DIMENSION(:), POINTER :: eigenvalues, occupation, wkp
891  TYPE(kpoint_env_type), POINTER :: kp
892  TYPE(mo_set_type), POINTER :: mo_set
893  TYPE(mp_para_env_type), POINTER :: para_env_inter_kp
894 
895  CALL timeset(routinen, handle)
896 
897  ! first collect all the eigenvalues
898  CALL get_kpoint_info(kpoint, nkp=nkp)
899  kp => kpoint%kp_env(1)%kpoint_env
900  nspin = SIZE(kp%mos, 2)
901  mo_set => kp%mos(1, 1)
902  CALL get_mo_set(mo_set, nmo=nmo, nelectron=nelectron)
903  ne_a = nelectron
904  IF (nspin == 2) THEN
905  CALL get_mo_set(kp%mos(1, 2), nmo=nb, nelectron=ne_b)
906  cpassert(nmo == nb)
907  END IF
908  ALLOCATE (weig(nmo, nkp, nspin), wocc(nmo, nkp, nspin))
909  weig = 0.0_dp
910  wocc = 0.0_dp
911  CALL get_kpoint_info(kpoint, kp_range=kp_range)
912  kplocal = kp_range(2) - kp_range(1) + 1
913  DO ikpgr = 1, kplocal
914  ik = kp_range(1) + ikpgr - 1
915  kp => kpoint%kp_env(ikpgr)%kpoint_env
916  DO ispin = 1, nspin
917  mo_set => kp%mos(1, ispin)
918  CALL get_mo_set(mo_set, eigenvalues=eigenvalues)
919  weig(1:nmo, ik, ispin) = eigenvalues(1:nmo)
920  END DO
921  END DO
922  CALL get_kpoint_info(kpoint, para_env_inter_kp=para_env_inter_kp)
923  CALL para_env_inter_kp%sum(weig)
924 
925  CALL get_kpoint_info(kpoint, wkp=wkp)
926  IF (smear%do_smear) THEN
927  ! finite electronic temperature
928  SELECT CASE (smear%method)
929  CASE (smear_fermi_dirac)
930  IF (nspin == 1) THEN
931  nel = real(nelectron, kind=dp)
932  CALL fermikp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, &
933  smear%electronic_temperature, 2.0_dp)
934  ELSE IF (smear%fixed_mag_mom > 0.0_dp) THEN
935  cpabort("kpoints: Smearing with fixed magnetic moments not (yet) supported")
936  nel = real(ne_a, kind=dp)
937  CALL fermikp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, &
938  smear%electronic_temperature, 1.0_dp)
939  nel = real(ne_b, kind=dp)
940  CALL fermikp(wocc(:, :, 2), mus(2), kts, weig(:, :, 2), nel, wkp, &
941  smear%electronic_temperature, 1.0_dp)
942  ELSE
943  nel = real(ne_a, kind=dp) + real(ne_b, kind=dp)
944  CALL fermikp2(wocc(:, :, :), mu, kts, weig(:, :, :), nel, wkp, &
945  smear%electronic_temperature)
946  kts = kts/2._dp
947  mus(1:2) = mu
948  END IF
949  CASE DEFAULT
950  cpabort("kpoints: Selected smearing not (yet) supported")
951  END SELECT
952  ELSE
953  ! fixed occupations (2/1)
954  IF (nspin == 1) THEN
955  nel = real(nelectron, kind=dp)
956  CALL fermikp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, 0.0_dp, 2.0_dp)
957  ELSE
958  nel = real(ne_a, kind=dp)
959  CALL fermikp(wocc(:, :, 1), mus(1), kts, weig(:, :, 1), nel, wkp, 0.0_dp, 1.0_dp)
960  nel = real(ne_b, kind=dp)
961  CALL fermikp(wocc(:, :, 2), mus(2), kts, weig(:, :, 2), nel, wkp, 0.0_dp, 1.0_dp)
962  END IF
963  END IF
964  DO ikpgr = 1, kplocal
965  ik = kp_range(1) + ikpgr - 1
966  kp => kpoint%kp_env(ikpgr)%kpoint_env
967  DO ispin = 1, nspin
968  mo_set => kp%mos(1, ispin)
969  CALL get_mo_set(mo_set, eigenvalues=eigenvalues, occupation_numbers=occupation)
970  eigenvalues(1:nmo) = weig(1:nmo, ik, ispin)
971  occupation(1:nmo) = wocc(1:nmo, ik, ispin)
972  mo_set%kTS = kts
973  mo_set%mu = mus(ispin)
974  END DO
975  END DO
976 
977  DEALLOCATE (weig, wocc)
978 
979  CALL timestop(handle)
980 
981  END SUBROUTINE kpoint_set_mo_occupation
982 
983 ! **************************************************************************************************
984 !> \brief Calculate kpoint density matrices (rho(k), owned by kpoint groups)
985 !> \param kpoint kpoint environment
986 !> \param energy_weighted calculate energy weighted density matrix
987 !> \param for_aux_fit ...
988 ! **************************************************************************************************
989  SUBROUTINE kpoint_density_matrices(kpoint, energy_weighted, for_aux_fit)
990 
991  TYPE(kpoint_type), POINTER :: kpoint
992  LOGICAL, OPTIONAL :: energy_weighted, for_aux_fit
993 
994  CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_density_matrices'
995 
996  INTEGER :: handle, ikpgr, ispin, kplocal, nao, nmo, &
997  nspin
998  INTEGER, DIMENSION(2) :: kp_range
999  LOGICAL :: aux_fit, wtype
1000  REAL(kind=dp), DIMENSION(:), POINTER :: eigenvalues, occupation
1001  TYPE(cp_fm_struct_type), POINTER :: matrix_struct
1002  TYPE(cp_fm_type) :: fwork
1003  TYPE(cp_fm_type), POINTER :: cpmat, pmat, rpmat
1004  TYPE(kpoint_env_type), POINTER :: kp
1005  TYPE(mo_set_type), POINTER :: mo_set
1006 
1007  CALL timeset(routinen, handle)
1008 
1009  IF (PRESENT(energy_weighted)) THEN
1010  wtype = energy_weighted
1011  ELSE
1012  ! default is normal density matrix
1013  wtype = .false.
1014  END IF
1015 
1016  IF (PRESENT(for_aux_fit)) THEN
1017  aux_fit = for_aux_fit
1018  ELSE
1019  aux_fit = .false.
1020  END IF
1021 
1022  IF (aux_fit) THEN
1023  cpassert(ASSOCIATED(kpoint%kp_aux_env))
1024  END IF
1025 
1026  ! work matrix
1027  IF (aux_fit) THEN
1028  mo_set => kpoint%kp_aux_env(1)%kpoint_env%mos(1, 1)
1029  ELSE
1030  mo_set => kpoint%kp_env(1)%kpoint_env%mos(1, 1)
1031  END IF
1032  CALL get_mo_set(mo_set, nao=nao, nmo=nmo)
1033  CALL cp_fm_get_info(mo_set%mo_coeff, matrix_struct=matrix_struct)
1034  CALL cp_fm_create(fwork, matrix_struct)
1035 
1036  CALL get_kpoint_info(kpoint, kp_range=kp_range)
1037  kplocal = kp_range(2) - kp_range(1) + 1
1038  DO ikpgr = 1, kplocal
1039  IF (aux_fit) THEN
1040  kp => kpoint%kp_aux_env(ikpgr)%kpoint_env
1041  ELSE
1042  kp => kpoint%kp_env(ikpgr)%kpoint_env
1043  END IF
1044  nspin = SIZE(kp%mos, 2)
1045  DO ispin = 1, nspin
1046  mo_set => kp%mos(1, ispin)
1047  IF (wtype) THEN
1048  CALL get_mo_set(mo_set, eigenvalues=eigenvalues)
1049  END IF
1050  IF (kpoint%use_real_wfn) THEN
1051  IF (wtype) THEN
1052  pmat => kp%wmat(1, ispin)
1053  ELSE
1054  pmat => kp%pmat(1, ispin)
1055  END IF
1056  CALL get_mo_set(mo_set, occupation_numbers=occupation)
1057  CALL cp_fm_to_fm(mo_set%mo_coeff, fwork)
1058  CALL cp_fm_column_scale(fwork, occupation)
1059  IF (wtype) THEN
1060  CALL cp_fm_column_scale(fwork, eigenvalues)
1061  END IF
1062  CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 0.0_dp, pmat)
1063  ELSE
1064  IF (wtype) THEN
1065  rpmat => kp%wmat(1, ispin)
1066  cpmat => kp%wmat(2, ispin)
1067  ELSE
1068  rpmat => kp%pmat(1, ispin)
1069  cpmat => kp%pmat(2, ispin)
1070  END IF
1071  CALL get_mo_set(mo_set, occupation_numbers=occupation)
1072  CALL cp_fm_to_fm(mo_set%mo_coeff, fwork)
1073  CALL cp_fm_column_scale(fwork, occupation)
1074  IF (wtype) THEN
1075  CALL cp_fm_column_scale(fwork, eigenvalues)
1076  END IF
1077  ! Re(c)*Re(c)
1078  CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 0.0_dp, rpmat)
1079  mo_set => kp%mos(2, ispin)
1080  ! Im(c)*Re(c)
1081  CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 0.0_dp, cpmat)
1082  ! Re(c)*Im(c)
1083  CALL parallel_gemm("N", "T", nao, nao, nmo, -1.0_dp, fwork, mo_set%mo_coeff, 1.0_dp, cpmat)
1084  CALL cp_fm_to_fm(mo_set%mo_coeff, fwork)
1085  CALL cp_fm_column_scale(fwork, occupation)
1086  IF (wtype) THEN
1087  CALL cp_fm_column_scale(fwork, eigenvalues)
1088  END IF
1089  ! Im(c)*Im(c)
1090  CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, mo_set%mo_coeff, fwork, 1.0_dp, rpmat)
1091  END IF
1092  END DO
1093  END DO
1094 
1095  CALL cp_fm_release(fwork)
1096 
1097  CALL timestop(handle)
1098 
1099  END SUBROUTINE kpoint_density_matrices
1100 
1101 ! **************************************************************************************************
1102 !> \brief generate real space density matrices in DBCSR format
1103 !> \param kpoint Kpoint environment
1104 !> \param denmat Real space (DBCSR) density matrices
1105 !> \param wtype True = energy weighted density matrix
1106 !> False = normal density matrix
1107 !> \param tempmat DBCSR matrix to be used as template
1108 !> \param sab_nl ...
1109 !> \param fmwork FM work matrices (kpoint group)
1110 !> \param for_aux_fit ...
1111 !> \param pmat_ext ...
1112 ! **************************************************************************************************
1113  SUBROUTINE kpoint_density_transform(kpoint, denmat, wtype, tempmat, sab_nl, fmwork, for_aux_fit, pmat_ext)
1114 
1115  TYPE(kpoint_type), POINTER :: kpoint
1116  TYPE(dbcsr_p_type), DIMENSION(:, :) :: denmat
1117  LOGICAL, INTENT(IN) :: wtype
1118  TYPE(dbcsr_type), POINTER :: tempmat
1119  TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1120  POINTER :: sab_nl
1121  TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fmwork
1122  LOGICAL, OPTIONAL :: for_aux_fit
1123  TYPE(cp_fm_type), DIMENSION(:, :, :), INTENT(IN), &
1124  OPTIONAL :: pmat_ext
1125 
1126  CHARACTER(LEN=*), PARAMETER :: routinen = 'kpoint_density_transform'
1127 
1128  INTEGER :: handle, ic, ik, ikk, indx, is, ispin, &
1129  nc, nimg, nkp, nspin
1130  INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
1131  LOGICAL :: aux_fit, do_ext, do_symmetric, my_kpgrp, &
1132  real_only
1133  REAL(kind=dp) :: wkpx
1134  REAL(kind=dp), DIMENSION(:), POINTER :: wkp
1135  REAL(kind=dp), DIMENSION(:, :), POINTER :: xkp
1136  TYPE(copy_info_type), ALLOCATABLE, DIMENSION(:) :: info
1137  TYPE(cp_fm_type) :: fmdummy
1138  TYPE(dbcsr_type), POINTER :: cpmat, rpmat, scpmat, srpmat
1139  TYPE(kpoint_env_type), POINTER :: kp
1140  TYPE(kpoint_sym_type), POINTER :: kpsym
1141  TYPE(mp_para_env_type), POINTER :: para_env
1142 
1143  CALL timeset(routinen, handle)
1144 
1145  CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
1146 
1147  IF (PRESENT(for_aux_fit)) THEN
1148  aux_fit = for_aux_fit
1149  ELSE
1150  aux_fit = .false.
1151  END IF
1152 
1153  do_ext = .false.
1154  IF (PRESENT(pmat_ext)) do_ext = .true.
1155 
1156  IF (aux_fit) THEN
1157  cpassert(ASSOCIATED(kpoint%kp_aux_env))
1158  END IF
1159 
1160  ! work storage
1161  ALLOCATE (rpmat)
1162  CALL dbcsr_create(rpmat, template=tempmat, &
1163  matrix_type=merge(dbcsr_type_symmetric, dbcsr_type_no_symmetry, do_symmetric))
1164  CALL cp_dbcsr_alloc_block_from_nbl(rpmat, sab_nl)
1165  ALLOCATE (cpmat)
1166  CALL dbcsr_create(cpmat, template=tempmat, &
1167  matrix_type=merge(dbcsr_type_antisymmetric, dbcsr_type_no_symmetry, do_symmetric))
1168  CALL cp_dbcsr_alloc_block_from_nbl(cpmat, sab_nl)
1169  IF (.NOT. kpoint%full_grid) THEN
1170  ALLOCATE (srpmat)
1171  CALL dbcsr_create(srpmat, template=rpmat)
1172  CALL cp_dbcsr_alloc_block_from_nbl(srpmat, sab_nl)
1173  ALLOCATE (scpmat)
1174  CALL dbcsr_create(scpmat, template=cpmat)
1175  CALL cp_dbcsr_alloc_block_from_nbl(scpmat, sab_nl)
1176  END IF
1177 
1178  CALL get_kpoint_info(kpoint, nkp=nkp, xkp=xkp, wkp=wkp, &
1179  cell_to_index=cell_to_index)
1180  ! initialize real space density matrices
1181  IF (aux_fit) THEN
1182  kp => kpoint%kp_aux_env(1)%kpoint_env
1183  ELSE
1184  kp => kpoint%kp_env(1)%kpoint_env
1185  END IF
1186  nspin = SIZE(kp%mos, 2)
1187  nc = SIZE(kp%mos, 1)
1188  nimg = SIZE(denmat, 2)
1189  real_only = (nc == 1)
1190 
1191  para_env => kpoint%blacs_env_all%para_env
1192  ALLOCATE (info(nspin*nkp*nc))
1193 
1194  ! Start all the communication
1195  indx = 0
1196  DO ispin = 1, nspin
1197  DO ic = 1, nimg
1198  CALL dbcsr_set(denmat(ispin, ic)%matrix, 0.0_dp)
1199  END DO
1200  !
1201  DO ik = 1, nkp
1202  my_kpgrp = (ik >= kpoint%kp_range(1) .AND. ik <= kpoint%kp_range(2))
1203  IF (my_kpgrp) THEN
1204  ikk = ik - kpoint%kp_range(1) + 1
1205  IF (aux_fit) THEN
1206  kp => kpoint%kp_aux_env(ikk)%kpoint_env
1207  ELSE
1208  kp => kpoint%kp_env(ikk)%kpoint_env
1209  END IF
1210  ELSE
1211  NULLIFY (kp)
1212  END IF
1213  ! collect this density matrix on all processors
1214  cpassert(SIZE(fmwork) >= nc)
1215 
1216  IF (my_kpgrp) THEN
1217  DO ic = 1, nc
1218  indx = indx + 1
1219  IF (do_ext) THEN
1220  CALL cp_fm_start_copy_general(pmat_ext(ikk, ic, ispin), fmwork(ic), para_env, info(indx))
1221  ELSE
1222  IF (wtype) THEN
1223  CALL cp_fm_start_copy_general(kp%wmat(ic, ispin), fmwork(ic), para_env, info(indx))
1224  ELSE
1225  CALL cp_fm_start_copy_general(kp%pmat(ic, ispin), fmwork(ic), para_env, info(indx))
1226  END IF
1227  END IF
1228  END DO
1229  ELSE
1230  DO ic = 1, nc
1231  indx = indx + 1
1232  CALL cp_fm_start_copy_general(fmdummy, fmwork(ic), para_env, info(indx))
1233  END DO
1234  END IF
1235  END DO
1236  END DO
1237 
1238  ! Finish communication and transform the received matrices
1239  indx = 0
1240  DO ispin = 1, nspin
1241  DO ik = 1, nkp
1242  DO ic = 1, nc
1243  indx = indx + 1
1244  CALL cp_fm_finish_copy_general(fmwork(ic), info(indx))
1245  END DO
1246 
1247  ! reduce to dbcsr storage
1248  IF (real_only) THEN
1249  CALL copy_fm_to_dbcsr(fmwork(1), rpmat, keep_sparsity=.true.)
1250  ELSE
1251  CALL copy_fm_to_dbcsr(fmwork(1), rpmat, keep_sparsity=.true.)
1252  CALL copy_fm_to_dbcsr(fmwork(2), cpmat, keep_sparsity=.true.)
1253  END IF
1254 
1255  ! symmetrization
1256  kpsym => kpoint%kp_sym(ik)%kpoint_sym
1257  cpassert(ASSOCIATED(kpsym))
1258 
1259  IF (kpsym%apply_symmetry) THEN
1260  wkpx = wkp(ik)/real(kpsym%nwght, kind=dp)
1261  DO is = 1, kpsym%nwght
1262  IF (real_only) THEN
1263  CALL symtrans(srpmat, rpmat, kpsym%rot(1:3, 1:3, is), kpsym%f0(:, is), symmetric=.true.)
1264  ELSE
1265  CALL symtrans(srpmat, rpmat, kpsym%rot(1:3, 1:3, is), kpsym%f0(:, is), symmetric=.true.)
1266  CALL symtrans(scpmat, cpmat, kpsym%rot(1:3, 1:3, is), kpsym%f0(:, is), antisymmetric=.true.)
1267  END IF
1268  CALL transform_dmat(denmat, srpmat, scpmat, ispin, real_only, sab_nl, &
1269  cell_to_index, kpsym%xkp(1:3, is), wkpx)
1270  END DO
1271  ELSE
1272  ! transformation
1273  CALL transform_dmat(denmat, rpmat, cpmat, ispin, real_only, sab_nl, &
1274  cell_to_index, xkp(1:3, ik), wkp(ik))
1275  END IF
1276  END DO
1277  END DO
1278 
1279  ! Clean up communication
1280  indx = 0
1281  DO ispin = 1, nspin
1282  DO ik = 1, nkp
1283  my_kpgrp = (ik >= kpoint%kp_range(1) .AND. ik <= kpoint%kp_range(2))
1284  IF (my_kpgrp) THEN
1285  ikk = ik - kpoint%kp_range(1) + 1
1286  IF (aux_fit) THEN
1287  kp => kpoint%kp_aux_env(ikk)%kpoint_env
1288  ELSE
1289  kp => kpoint%kp_env(ikk)%kpoint_env
1290  END IF
1291 
1292  DO ic = 1, nc
1293  indx = indx + 1
1294  CALL cp_fm_cleanup_copy_general(info(indx))
1295  END DO
1296  ELSE
1297  ! calls with dummy arguments, so not included
1298  ! therefore just increment counter by trip count
1299  indx = indx + nc
1300  END IF
1301  END DO
1302  END DO
1303 
1304  ! All done
1305  DEALLOCATE (info)
1306 
1307  CALL dbcsr_deallocate_matrix(rpmat)
1308  CALL dbcsr_deallocate_matrix(cpmat)
1309  IF (.NOT. kpoint%full_grid) THEN
1310  CALL dbcsr_deallocate_matrix(srpmat)
1311  CALL dbcsr_deallocate_matrix(scpmat)
1312  END IF
1313 
1314  CALL timestop(handle)
1315 
1316  END SUBROUTINE kpoint_density_transform
1317 
1318 ! **************************************************************************************************
1319 !> \brief real space density matrices in DBCSR format
1320 !> \param denmat Real space (DBCSR) density matrix
1321 !> \param rpmat ...
1322 !> \param cpmat ...
1323 !> \param ispin ...
1324 !> \param real_only ...
1325 !> \param sab_nl ...
1326 !> \param cell_to_index ...
1327 !> \param xkp ...
1328 !> \param wkp ...
1329 ! **************************************************************************************************
1330  SUBROUTINE transform_dmat(denmat, rpmat, cpmat, ispin, real_only, sab_nl, cell_to_index, xkp, wkp)
1331 
1332  TYPE(dbcsr_p_type), DIMENSION(:, :) :: denmat
1333  TYPE(dbcsr_type), POINTER :: rpmat, cpmat
1334  INTEGER, INTENT(IN) :: ispin
1335  LOGICAL, INTENT(IN) :: real_only
1336  TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1337  POINTER :: sab_nl
1338  INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
1339  REAL(kind=dp), DIMENSION(3), INTENT(IN) :: xkp
1340  REAL(kind=dp), INTENT(IN) :: wkp
1341 
1342  CHARACTER(LEN=*), PARAMETER :: routinen = 'transform_dmat'
1343 
1344  INTEGER :: handle, iatom, icell, icol, irow, jatom, &
1345  nimg
1346  INTEGER, DIMENSION(3) :: cell
1347  LOGICAL :: do_symmetric, found
1348  REAL(kind=dp) :: arg, coskl, fc, sinkl
1349  REAL(kind=dp), DIMENSION(:, :), POINTER :: cblock, dblock, rblock
1350  TYPE(neighbor_list_iterator_p_type), &
1351  DIMENSION(:), POINTER :: nl_iterator
1352 
1353  CALL timeset(routinen, handle)
1354 
1355  nimg = SIZE(denmat, 2)
1356 
1357  ! transformation
1358  CALL get_neighbor_list_set_p(neighbor_list_sets=sab_nl, symmetric=do_symmetric)
1359  CALL neighbor_list_iterator_create(nl_iterator, sab_nl)
1360  DO WHILE (neighbor_list_iterate(nl_iterator) == 0)
1361  CALL get_iterator_info(nl_iterator, iatom=iatom, jatom=jatom, cell=cell)
1362 
1363  !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
1364  !Therefore, we have: S(R) = sum_k Re(S(k))*cos(k*R) -i^2*Im(S(k))*sin(k*R)
1365  ! = sum_k Re(S(k))*cos(k*R) + Im(S(k))*sin(k*R)
1366  !fc = +- 1 is due to the usual non-symmetric real-sapce matrices stored as symmetric ones
1367 
1368  irow = iatom
1369  icol = jatom
1370  fc = 1.0_dp
1371  IF (do_symmetric .AND. iatom > jatom) THEN
1372  irow = jatom
1373  icol = iatom
1374  fc = -1.0_dp
1375  END IF
1376 
1377  icell = cell_to_index(cell(1), cell(2), cell(3))
1378  IF (icell < 1 .OR. icell > nimg) cycle
1379 
1380  arg = real(cell(1), dp)*xkp(1) + real(cell(2), dp)*xkp(2) + real(cell(3), dp)*xkp(3)
1381  coskl = wkp*cos(twopi*arg)
1382  sinkl = wkp*fc*sin(twopi*arg)
1383 
1384  CALL dbcsr_get_block_p(matrix=denmat(ispin, icell)%matrix, row=irow, col=icol, &
1385  block=dblock, found=found)
1386  IF (.NOT. found) cycle
1387 
1388  IF (real_only) THEN
1389  CALL dbcsr_get_block_p(matrix=rpmat, row=irow, col=icol, block=rblock, found=found)
1390  IF (.NOT. found) cycle
1391  dblock = dblock + coskl*rblock
1392  ELSE
1393  CALL dbcsr_get_block_p(matrix=rpmat, row=irow, col=icol, block=rblock, found=found)
1394  IF (.NOT. found) cycle
1395  CALL dbcsr_get_block_p(matrix=cpmat, row=irow, col=icol, block=cblock, found=found)
1396  IF (.NOT. found) cycle
1397  dblock = dblock + coskl*rblock
1398  dblock = dblock + sinkl*cblock
1399  END IF
1400  END DO
1401  CALL neighbor_list_iterator_release(nl_iterator)
1402 
1403  CALL timestop(handle)
1404 
1405  END SUBROUTINE transform_dmat
1406 
1407 ! **************************************************************************************************
1408 !> \brief Symmetrization of density matrix - transform to new k-point
1409 !> \param smat density matrix at new kpoint
1410 !> \param pmat reference density matrix
1411 !> \param rot Rotation matrix
1412 !> \param f0 Permutation of atoms under transformation
1413 !> \param symmetric Symmetric matrix
1414 !> \param antisymmetric Anti-Symmetric matrix
1415 ! **************************************************************************************************
1416  SUBROUTINE symtrans(smat, pmat, rot, f0, symmetric, antisymmetric)
1417  TYPE(dbcsr_type), POINTER :: smat, pmat
1418  REAL(kind=dp), DIMENSION(3, 3), INTENT(IN) :: rot
1419  INTEGER, DIMENSION(:), INTENT(IN) :: f0
1420  LOGICAL, INTENT(IN), OPTIONAL :: symmetric, antisymmetric
1421 
1422  CHARACTER(LEN=*), PARAMETER :: routinen = 'symtrans'
1423 
1424  INTEGER :: handle, iatom, icol, ip, irow, jcol, jp, &
1425  jrow, natom, numnodes
1426  LOGICAL :: asym, dorot, found, perm, sym, trans
1427  REAL(kind=dp) :: dr, fsign
1428  REAL(kind=dp), DIMENSION(:, :), POINTER :: pblock, sblock
1429  TYPE(dbcsr_distribution_type) :: dist
1430  TYPE(dbcsr_iterator_type) :: iter
1431 
1432  CALL timeset(routinen, handle)
1433 
1434  ! check symmetry options
1435  sym = .false.
1436  IF (PRESENT(symmetric)) sym = symmetric
1437  asym = .false.
1438  IF (PRESENT(antisymmetric)) asym = antisymmetric
1439 
1440  cpassert(.NOT. (sym .AND. asym))
1441  cpassert((sym .OR. asym))
1442 
1443  ! do we have permutation of atoms
1444  natom = SIZE(f0)
1445  perm = .false.
1446  DO iatom = 1, natom
1447  IF (f0(iatom) == iatom) cycle
1448  perm = .true.
1449  EXIT
1450  END DO
1451 
1452  ! do we have a real rotation
1453  dorot = .false.
1454  IF (abs(sum(abs(rot)) - 3.0_dp) > 1.e-12_dp) dorot = .true.
1455  dr = abs(rot(1, 1) - 1.0_dp) + abs(rot(2, 2) - 1.0_dp) + abs(rot(3, 3) - 1.0_dp)
1456  IF (abs(dr) > 1.e-12_dp) dorot = .true.
1457 
1458  fsign = 1.0_dp
1459  IF (asym) fsign = -1.0_dp
1460 
1461  IF (dorot .OR. perm) THEN
1462  CALL dbcsr_set(smat, 0.0_dp)
1463  IF (perm) THEN
1464  CALL dbcsr_get_info(pmat, distribution=dist)
1465  CALL dbcsr_distribution_get(dist, numnodes=numnodes)
1466  IF (numnodes == 1) THEN
1467  ! the matrices are local to this process
1468  CALL dbcsr_iterator_start(iter, pmat)
1469  DO WHILE (dbcsr_iterator_blocks_left(iter))
1470  CALL dbcsr_iterator_next_block(iter, irow, icol, pblock)
1471  ip = f0(irow)
1472  jp = f0(icol)
1473  IF (ip <= jp) THEN
1474  jrow = ip
1475  jcol = jp
1476  trans = .false.
1477  ELSE
1478  jrow = jp
1479  jcol = ip
1480  trans = .true.
1481  END IF
1482  CALL dbcsr_get_block_p(matrix=smat, row=jrow, col=jcol, block=sblock, found=found)
1483  IF (.NOT. found) cycle
1484  IF (trans) THEN
1485  sblock = fsign*transpose(pblock)
1486  ELSE
1487  sblock = pblock
1488  END IF
1489  END DO
1490  CALL dbcsr_iterator_stop(iter)
1491  !
1492  ELSE
1493  ! distributed matrices, most general code needed
1494  CALL cp_abort(__location__, "k-points need FULL_GRID currently. "// &
1495  "Reduced grids not yet working correctly")
1496  END IF
1497  ELSE
1498  ! no atom permutations, this is always local
1499  ! ignore rotations for now
1500  CALL dbcsr_copy(smat, pmat)
1501  CALL cp_abort(__location__, "k-points need FULL_GRID currently. "// &
1502  "Reduced grids not yet working correctly")
1503  END IF
1504  ELSE
1505  ! this is the identity operation, just copy the matrix
1506  CALL dbcsr_copy(smat, pmat)
1507  END IF
1508 
1509  CALL timestop(handle)
1510 
1511  END SUBROUTINE symtrans
1512 
1513 ! **************************************************************************************************
1514 
1515 END 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:486
methods related to the blacs parallel environment
Definition: cp_blacs_env.F:15
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
Definition: cp_blacs_env.F:123
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
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
Definition: cp_fm_struct.F:14
represent a full matrix distributed on many processors
Definition: cp_fm_types.F:15
subroutine, public cp_fm_start_copy_general(source, destination, para_env, info)
Initiates the copy operation: get distribution data, post MPI isend and irecvs.
Definition: cp_fm_types.F:1568
subroutine, public cp_fm_cleanup_copy_general(info)
Completes the copy operation: wait for comms clean up MPI state.
Definition: cp_fm_types.F:1946
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
Definition: cp_fm_types.F:1016
subroutine, public cp_fm_finish_copy_general(destination, info)
Completes the copy operation: wait for comms, unpack, clean up MPI state.
Definition: cp_fm_types.F:1882
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
Definition: cp_fm_types.F:167
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:510
subroutine, public crys_sym_gen(csym, scoor, types, hmat, delta, iounit)
...
Definition: cryssym.F:114
subroutine, public release_csym_type(csym)
Release the CSYM type.
Definition: cryssym.F:72
subroutine, public kpoint_gen(csym, nk, symm, shift, full_grid)
...
Definition: cryssym.F:209
subroutine, public apply_rotation_coord(f0, csym, ir)
...
Definition: cryssym.F:473
subroutine, public print_kp_symmetry(csym)
...
Definition: cryssym.F:543
deal with the Fermi distribution, compute it, fix mu, get derivs
Definition: fermi_utils.F:13
subroutine, public fermikp2(f, mu, kTS, e, nel, wk, t)
Bisection search to find mu for a given nel (kpoint case)
Definition: fermi_utils.F:336
subroutine, public fermikp(f, mu, kTS, e, nel, wk, t, maxocc)
Bisection search to find mu for a given nel (kpoint case)
Definition: fermi_utils.F:283
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public smear_fermi_dirac
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_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_init_cell_index(kpoint, sab_nl, para_env, dft_control)
Generates the mapping of cell indices and linear RS index CELL (0,0,0) is always mapped to index 1.
subroutine, public kpoint_density_matrices(kpoint, energy_weighted, for_aux_fit)
Calculate kpoint density matrices (rho(k), owned by kpoint groups)
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)
Given the eigenvalues of all kpoints, calculates the occupation numbers.
Types and basic routines needed for a kpoint calculation.
Definition: kpoint_types.F:15
subroutine, public kpoint_sym_create(kp_sym)
Create a single kpoint symmetry environment.
Definition: kpoint_types.F:769
subroutine, public kpoint_env_create(kp_env)
Create a single kpoint environment.
Definition: kpoint_types.F:685
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.
Definition: kpoint_types.F:333
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.
Definition: mathconstants.F:16
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.
Definition: qs_mo_types.F:452
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...
Definition: qs_mo_types.F:206
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.
Definition: qs_mo_types.F:397
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 ...
Definition: qs_mo_types.F:245
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
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