Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \brief Types and basic routines needed for a kpoint calculation
10 : !> \par History
11 : !> 2014.07 created [JGH]
12 : !> 2014.11 unified k-point and gamma-point code [Ole Schuett]
13 : !> \author JGH
14 : ! **************************************************************************************************
15 : MODULE kpoint_types
16 : USE cell_types, ONLY: cell_type
17 : USE cp_blacs_env, ONLY: cp_blacs_env_release,&
18 : cp_blacs_env_type
19 : USE cp_cfm_types, ONLY: cp_cfm_release,&
20 : cp_cfm_type
21 : USE cp_fm_types, ONLY: cp_fm_release,&
22 : cp_fm_type
23 : USE cp_log_handling, ONLY: cp_get_default_logger,&
24 : cp_logger_get_default_unit_nr,&
25 : cp_logger_type
26 : USE cp_output_handling, ONLY: cp_print_key_finished_output,&
27 : cp_print_key_unit_nr
28 : USE input_cp2k_kpoints, ONLY: use_complex_wfn,&
29 : use_k290_kpoint_backend,&
30 : use_k290_kpoint_symmetry,&
31 : use_real_wfn,&
32 : use_spglib_kpoint_backend,&
33 : use_spglib_kpoint_symmetry
34 : USE input_section_types, ONLY: section_vals_get,&
35 : section_vals_type,&
36 : section_vals_val_get
37 : USE kinds, ONLY: default_string_length,&
38 : dp
39 : USE mathconstants, ONLY: twopi
40 : USE message_passing, ONLY: mp_para_env_release,&
41 : mp_para_env_type
42 : USE physcon, ONLY: angstrom
43 : USE qs_diis_types, ONLY: qs_diis_b_release_kp,&
44 : qs_diis_buffer_type_kp
45 : USE qs_matrix_pools, ONLY: mpools_release,&
46 : qs_matrix_pools_type
47 : USE qs_mo_types, ONLY: deallocate_mo_set,&
48 : mo_set_type
49 : USE qs_neighbor_list_types, ONLY: neighbor_list_set_p_type
50 : USE string_utilities, ONLY: uppercase
51 : #include "./base/base_uses.f90"
52 :
53 : IMPLICIT NONE
54 :
55 : PRIVATE
56 :
57 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'kpoint_types'
58 :
59 : PUBLIC :: kpoint_type
60 : PUBLIC :: kpoint_create, kpoint_release, kpoint_reset_initialization, get_kpoint_info, set_kpoint_info
61 : PUBLIC :: read_kpoint_section, write_kpoint_info
62 : PUBLIC :: kpoint_env_type, kpoint_env_p_type
63 : PUBLIC :: kpoint_env_create, get_kpoint_env
64 : PUBLIC :: kind_rotmat_type
65 : PUBLIC :: kpoint_sym_type
66 : PUBLIC :: kpoint_sym_create
67 :
68 : ! **************************************************************************************************
69 : !> \brief Keeps information about a specific k-point
70 : !> \param nkpoint the kpoint index
71 : !> \param wkp weight of this kpoint
72 : !> \param xkp kpoint coordinates in units of b-vector
73 : !> \param is_local if this kpoint is calculated on a single thread
74 : !> \param mos associated MOs (r/i,spin)
75 : !> \param pmat associated density matrix (r/i,spin)
76 : !> \param wmat associated energy weighted density matrix (r/i,spin)
77 : !> \param smat associated overlap matrix (for ADMM) (r/i,spin)
78 : !> \param amat associated ADMM basis projection matrix (r/i,spin)
79 : !> \param shalf S(K)^(1/2) DFT+U Lowdin method (real wfn)
80 : !> \param cshalf S(K)^(1/2) DFT+U Lowdin method (complex wfn)
81 : !> \author JGH
82 : ! **************************************************************************************************
83 : TYPE kpoint_env_type
84 : INTEGER :: nkpoint = -1
85 : REAL(KIND=dp) :: wkp = 0.0_dp
86 : REAL(KIND=dp), DIMENSION(3) :: xkp = 0.0_dp
87 : LOGICAL :: is_local = .FALSE.
88 : TYPE(mo_set_type), DIMENSION(:, :), POINTER :: mos => NULL()
89 : TYPE(cp_fm_type), DIMENSION(:, :), POINTER :: pmat => NULL()
90 : TYPE(cp_fm_type), DIMENSION(:, :), POINTER :: wmat => NULL()
91 : TYPE(cp_fm_type), DIMENSION(:, :), POINTER :: smat => NULL()
92 : TYPE(cp_fm_type), DIMENSION(:, :), POINTER :: amat => NULL()
93 : TYPE(cp_fm_type) :: shalf
94 : TYPE(cp_cfm_type) :: cshalf
95 : END TYPE kpoint_env_type
96 :
97 : TYPE kpoint_env_p_type
98 : TYPE(kpoint_env_type), POINTER :: kpoint_env => NULL()
99 : END TYPE kpoint_env_p_type
100 :
101 : ! **************************************************************************************************
102 : !> \brief Rotation matrices for basis sets
103 : !> \param rmat atom basis function rotation matrix
104 : !> \author JGH
105 : ! **************************************************************************************************
106 : TYPE kind_rotmat_type
107 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: rmat => NULL()
108 : END TYPE kind_rotmat_type
109 :
110 : ! **************************************************************************************************
111 : !> \brief Keeps symmetry information about a specific k-point
112 : !> \param apply_symmetry ...
113 : !> \param nwght kpoint multiplicity
114 : !> \param xkp kpoint coordinates
115 : !> \param rot rotation matrices
116 : !> \param f0 atom permutation
117 : !> \param fcell atom cell shifts generated by the symmetry operation
118 : !> \param kgphase atom Bloch gauge from reciprocal-lattice folding of the mapped k-point
119 : !> \author JGH
120 : ! **************************************************************************************************
121 : TYPE kpoint_sym_type
122 : LOGICAL :: apply_symmetry = .FALSE.
123 : INTEGER :: nwght = -1
124 : INTEGER :: nwred = -1
125 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: xkp => NULL()
126 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: rot => NULL()
127 : INTEGER, DIMENSION(:), POINTER :: rotp => NULL()
128 : INTEGER, DIMENSION(:, :), POINTER :: f0 => NULL()
129 : INTEGER, DIMENSION(:, :, :), POINTER :: fcell => NULL()
130 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: kgphase => NULL()
131 : END TYPE kpoint_sym_type
132 :
133 : TYPE kpoint_sym_p_type
134 : TYPE(kpoint_sym_type), POINTER :: kpoint_sym => NULL()
135 : END TYPE kpoint_sym_p_type
136 :
137 : ! **************************************************************************************************
138 : !> \brief Contains information about kpoints
139 : !> \par History
140 : !> 2014.07 created [JGH]
141 : !> \param kp_scheme [input] Type of kpoint grid
142 : !> \param nkp_grid [input] Grid points
143 : !> \param kp_shift [input] Shift of the grid
144 : !> \param use_real_wfn [input] real/complex wfn
145 : !> \param symmetry [input] use symmetry (atoms) to reduce kpoints
146 : !> \param full_grid [input] don't reduce kpoints at all
147 : !> \param inversion_symmetry_only [input] reduce kpoints only by inversion symmetry
148 : !> \param symmetry_backend [input] k-point symmetry backend
149 : !> \param symmetry_reduction_method [input] k-point symmetry reduction method
150 : !> \param verbose [input] more output information
151 : !> \param eps_geo [input] accuracy of atom symmetry detection
152 : !> \param parallel_group_size [input] kpoint group size
153 : !> \param nkp number of kpoints
154 : !> \param xkp kpoint coordinates
155 : !> \param wkp kpoint weights
156 : !> \param xkp_input explicit GENERAL kpoint coordinates as read from the input
157 : !> \param wkp_input explicit GENERAL kpoint weights as read from the input
158 : !> \param para_env 'global' parallel environment
159 : !> \param para_env_kp parallel environment of the kpoint calculation
160 : !> \param para_env_inter_kp parallel environment between kpoints
161 : !> \param iogrp this kpoint group has the IO processor
162 : !> \param nkp_groups number of kpoint groups
163 : !> \param kp_dist kpoints distribution on groups
164 : !> \param kp_range kpoints distribution for local processor
165 : !> \param blacs_env BLACS env for the kpoint group
166 : !> \param opmats Operator matrices
167 : !> \param kp_env Information for each kpoint
168 : !> \param mpools FM matrix pools for kpoint groups
169 : !> \author JGH
170 : ! **************************************************************************************************
171 : TYPE kpoint_type
172 : CHARACTER(LEN=default_string_length) :: kp_scheme = ""
173 : INTEGER, DIMENSION(3) :: nkp_grid = -1
174 : REAL(KIND=dp), DIMENSION(3) :: kp_shift = 0.0_dp
175 : LOGICAL :: gamma_centered = .FALSE.
176 : LOGICAL :: use_real_wfn = .FALSE.
177 : LOGICAL :: symmetry = .FALSE.
178 : LOGICAL :: full_grid = .FALSE.
179 : LOGICAL :: inversion_symmetry_only = .FALSE.
180 : INTEGER :: symmetry_backend = use_k290_kpoint_backend
181 : INTEGER :: symmetry_reduction_method = use_k290_kpoint_symmetry
182 : LOGICAL :: verbose = .FALSE.
183 : REAL(KIND=dp) :: eps_geo = 0.0_dp
184 : INTEGER :: parallel_group_size = -1
185 : INTEGER :: nkp = -1
186 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: xkp => Null()
187 : REAL(KIND=dp), DIMENSION(:), POINTER :: wkp => Null()
188 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: xkp_input => Null()
189 : REAL(KIND=dp), DIMENSION(:), POINTER :: wkp_input => Null()
190 : ! parallel environment
191 : TYPE(mp_para_env_type), POINTER :: para_env => Null()
192 : TYPE(cp_blacs_env_type), POINTER :: blacs_env_all => Null()
193 : TYPE(mp_para_env_type), POINTER :: para_env_kp => Null(), &
194 : para_env_inter_kp => Null()
195 : LOGICAL :: iogrp = .FALSE.
196 : INTEGER :: nkp_groups = -1
197 : INTEGER, DIMENSION(:, :), POINTER :: kp_dist => Null()
198 : INTEGER, DIMENSION(2) :: kp_range = -1
199 : TYPE(cp_blacs_env_type), POINTER :: blacs_env => Null()
200 : INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index => Null()
201 : INTEGER, DIMENSION(:, :), POINTER :: index_to_cell => Null()
202 : TYPE(neighbor_list_set_p_type), &
203 : DIMENSION(:), POINTER :: sab_nl => Null(), &
204 : sab_nl_nosym => Null()
205 : ! environment
206 : TYPE(kpoint_env_p_type), DIMENSION(:), &
207 : POINTER :: kp_env => Null()
208 : TYPE(kpoint_env_p_type), DIMENSION(:), &
209 : POINTER :: kp_aux_env => Null()
210 : ! symmetry
211 : TYPE(kpoint_sym_p_type), DIMENSION(:), &
212 : POINTER :: kp_sym => Null()
213 : INTEGER, DIMENSION(:), POINTER :: atype => Null()
214 : INTEGER, DIMENSION(:), POINTER :: ibrot => Null()
215 : TYPE(kind_rotmat_type), DIMENSION(:, :), &
216 : POINTER :: kind_rotmat => Null()
217 : ! pools
218 : TYPE(qs_matrix_pools_type), POINTER :: mpools => Null()
219 : TYPE(qs_diis_buffer_type_kp), POINTER :: scf_diis_buffer => Null()
220 : TYPE(qs_matrix_pools_type), POINTER :: mpools_aux_fit => Null()
221 : END TYPE kpoint_type
222 :
223 : ! **************************************************************************************************
224 :
225 : CONTAINS
226 :
227 : ! **************************************************************************************************
228 : !> \brief Create a kpoint environment
229 : !> \param kpoint All the kpoint information
230 : !> \author JGH
231 : ! **************************************************************************************************
232 8656 : SUBROUTINE kpoint_create(kpoint)
233 : TYPE(kpoint_type), POINTER :: kpoint
234 :
235 8656 : CPASSERT(.NOT. ASSOCIATED(kpoint))
236 :
237 95216 : ALLOCATE (kpoint)
238 :
239 8656 : kpoint%kp_scheme = ""
240 34624 : kpoint%nkp_grid = 0
241 34624 : kpoint%kp_shift = 0.0_dp
242 8656 : kpoint%gamma_centered = .FALSE.
243 8656 : kpoint%symmetry = .FALSE.
244 8656 : kpoint%verbose = .FALSE.
245 8656 : kpoint%full_grid = .FALSE.
246 8656 : kpoint%inversion_symmetry_only = .FALSE.
247 8656 : kpoint%symmetry_backend = use_k290_kpoint_backend
248 8656 : kpoint%symmetry_reduction_method = use_k290_kpoint_symmetry
249 8656 : kpoint%use_real_wfn = .FALSE.
250 8656 : kpoint%eps_geo = 1.0e-6_dp
251 8656 : kpoint%parallel_group_size = -1
252 :
253 8656 : kpoint%nkp = 0
254 :
255 8656 : NULLIFY (kpoint%xkp, kpoint%wkp)
256 8656 : NULLIFY (kpoint%xkp_input, kpoint%wkp_input)
257 8656 : NULLIFY (kpoint%kp_dist)
258 :
259 8656 : NULLIFY (kpoint%para_env)
260 8656 : NULLIFY (kpoint%blacs_env_all)
261 8656 : NULLIFY (kpoint%para_env_kp, kpoint%para_env_inter_kp)
262 8656 : NULLIFY (kpoint%blacs_env)
263 8656 : kpoint%nkp_groups = 0
264 8656 : kpoint%iogrp = .FALSE.
265 25968 : kpoint%kp_range = 0
266 :
267 8656 : NULLIFY (kpoint%kp_env)
268 8656 : NULLIFY (kpoint%mpools)
269 :
270 8656 : ALLOCATE (kpoint%cell_to_index(0:0, 0:0, 0:0))
271 34624 : kpoint%cell_to_index(:, :, :) = 1
272 :
273 8656 : ALLOCATE (kpoint%index_to_cell(0:0, 0:0))
274 25968 : kpoint%index_to_cell(:, :) = 0
275 :
276 8656 : END SUBROUTINE kpoint_create
277 :
278 : ! **************************************************************************************************
279 : !> \brief Release a kpoint environment, deallocate all data
280 : !> \param kpoint The kpoint environment
281 : !> \author JGH
282 : ! **************************************************************************************************
283 17083 : SUBROUTINE kpoint_release(kpoint)
284 : TYPE(kpoint_type), POINTER :: kpoint
285 :
286 : INTEGER :: i, ik, j
287 :
288 17083 : IF (ASSOCIATED(kpoint)) THEN
289 :
290 8656 : IF (ASSOCIATED(kpoint%xkp)) THEN
291 686 : DEALLOCATE (kpoint%xkp)
292 : END IF
293 8656 : IF (ASSOCIATED(kpoint%wkp)) THEN
294 686 : DEALLOCATE (kpoint%wkp)
295 : END IF
296 8656 : IF (ASSOCIATED(kpoint%xkp_input)) THEN
297 18 : DEALLOCATE (kpoint%xkp_input)
298 : END IF
299 8656 : IF (ASSOCIATED(kpoint%wkp_input)) THEN
300 18 : DEALLOCATE (kpoint%wkp_input)
301 : END IF
302 8656 : IF (ASSOCIATED(kpoint%kp_dist)) THEN
303 562 : DEALLOCATE (kpoint%kp_dist)
304 : END IF
305 :
306 8656 : CALL mpools_release(kpoint%mpools)
307 8656 : CALL mpools_release(kpoint%mpools_aux_fit)
308 :
309 8656 : CALL cp_blacs_env_release(kpoint%blacs_env)
310 8656 : CALL cp_blacs_env_release(kpoint%blacs_env_all)
311 :
312 8656 : CALL mp_para_env_release(kpoint%para_env)
313 8656 : CALL mp_para_env_release(kpoint%para_env_kp)
314 8656 : CALL mp_para_env_release(kpoint%para_env_inter_kp)
315 :
316 8656 : IF (ASSOCIATED(kpoint%cell_to_index)) DEALLOCATE (kpoint%cell_to_index)
317 8656 : IF (ASSOCIATED(kpoint%index_to_cell)) DEALLOCATE (kpoint%index_to_cell)
318 :
319 8656 : IF (ASSOCIATED(kpoint%kp_env)) THEN
320 3462 : DO ik = 1, SIZE(kpoint%kp_env)
321 3462 : CALL kpoint_env_release(kpoint%kp_env(ik)%kpoint_env)
322 : END DO
323 562 : DEALLOCATE (kpoint%kp_env)
324 : END IF
325 :
326 8656 : IF (ASSOCIATED(kpoint%kp_aux_env)) THEN
327 250 : DO ik = 1, SIZE(kpoint%kp_aux_env)
328 250 : CALL kpoint_env_release(kpoint%kp_aux_env(ik)%kpoint_env)
329 : END DO
330 32 : DEALLOCATE (kpoint%kp_aux_env)
331 : END IF
332 :
333 8656 : IF (ASSOCIATED(kpoint%kp_sym)) THEN
334 15460 : DO ik = 1, SIZE(kpoint%kp_sym)
335 15460 : CALL kpoint_sym_release(kpoint%kp_sym(ik)%kpoint_sym)
336 : END DO
337 564 : DEALLOCATE (kpoint%kp_sym)
338 : END IF
339 :
340 8656 : IF (ASSOCIATED(kpoint%atype)) DEALLOCATE (kpoint%atype)
341 8656 : IF (ASSOCIATED(kpoint%ibrot)) DEALLOCATE (kpoint%ibrot)
342 :
343 8656 : IF (ASSOCIATED(kpoint%kind_rotmat)) THEN
344 16194 : DO i = 1, SIZE(kpoint%kind_rotmat, 1)
345 32072 : DO j = 1, SIZE(kpoint%kind_rotmat, 2)
346 31756 : IF (ASSOCIATED(kpoint%kind_rotmat(i, j)%rmat)) THEN
347 15462 : DEALLOCATE (kpoint%kind_rotmat(i, j)%rmat)
348 : END IF
349 : END DO
350 : END DO
351 316 : DEALLOCATE (kpoint%kind_rotmat)
352 : END IF
353 :
354 8656 : IF (ASSOCIATED(kpoint%scf_diis_buffer)) THEN
355 296 : CALL qs_diis_b_release_kp(kpoint%scf_diis_buffer)
356 296 : DEALLOCATE (kpoint%scf_diis_buffer)
357 : END IF
358 :
359 8656 : DEALLOCATE (kpoint)
360 :
361 : END IF
362 :
363 17083 : END SUBROUTINE kpoint_release
364 :
365 : ! **************************************************************************************************
366 : !> \brief Reset all data derived from a concrete k-point initialization.
367 : !> Input options such as the scheme, grid, shifts and symmetry settings are kept.
368 : !> \param kpoint The kpoint environment
369 : ! **************************************************************************************************
370 2330 : SUBROUTINE kpoint_reset_initialization(kpoint)
371 : TYPE(kpoint_type), INTENT(INOUT) :: kpoint
372 :
373 : INTEGER :: i, ik, j
374 :
375 2330 : IF (ASSOCIATED(kpoint%xkp)) THEN
376 2330 : DEALLOCATE (kpoint%xkp)
377 : NULLIFY (kpoint%xkp)
378 : END IF
379 2330 : IF (ASSOCIATED(kpoint%wkp)) THEN
380 2330 : DEALLOCATE (kpoint%wkp)
381 : NULLIFY (kpoint%wkp)
382 : END IF
383 2330 : IF (kpoint%kp_scheme == "GENERAL" .AND. ASSOCIATED(kpoint%xkp_input) .AND. &
384 : ASSOCIATED(kpoint%wkp_input)) THEN
385 18 : kpoint%nkp = SIZE(kpoint%wkp_input)
386 90 : ALLOCATE (kpoint%xkp(3, kpoint%nkp), kpoint%wkp(kpoint%nkp))
387 594 : kpoint%xkp(1:3, 1:kpoint%nkp) = kpoint%xkp_input(1:3, 1:kpoint%nkp)
388 162 : kpoint%wkp(1:kpoint%nkp) = kpoint%wkp_input(1:kpoint%nkp)
389 : END IF
390 2330 : IF (ASSOCIATED(kpoint%kp_dist)) THEN
391 2112 : DEALLOCATE (kpoint%kp_dist)
392 : NULLIFY (kpoint%kp_dist)
393 : END IF
394 :
395 2330 : CALL mpools_release(kpoint%mpools)
396 2330 : CALL mpools_release(kpoint%mpools_aux_fit)
397 :
398 2330 : CALL cp_blacs_env_release(kpoint%blacs_env)
399 2330 : CALL cp_blacs_env_release(kpoint%blacs_env_all)
400 :
401 2330 : CALL mp_para_env_release(kpoint%para_env)
402 2330 : CALL mp_para_env_release(kpoint%para_env_kp)
403 2330 : CALL mp_para_env_release(kpoint%para_env_inter_kp)
404 :
405 2330 : IF (ASSOCIATED(kpoint%cell_to_index)) THEN
406 2330 : DEALLOCATE (kpoint%cell_to_index)
407 : NULLIFY (kpoint%cell_to_index)
408 : END IF
409 2330 : IF (ASSOCIATED(kpoint%index_to_cell)) THEN
410 2330 : DEALLOCATE (kpoint%index_to_cell)
411 : NULLIFY (kpoint%index_to_cell)
412 : END IF
413 :
414 2330 : IF (ASSOCIATED(kpoint%kp_env)) THEN
415 5776 : DO ik = 1, SIZE(kpoint%kp_env)
416 5776 : CALL kpoint_env_release(kpoint%kp_env(ik)%kpoint_env)
417 : END DO
418 2112 : DEALLOCATE (kpoint%kp_env)
419 : NULLIFY (kpoint%kp_env)
420 : END IF
421 :
422 2330 : IF (ASSOCIATED(kpoint%kp_aux_env)) THEN
423 0 : DO ik = 1, SIZE(kpoint%kp_aux_env)
424 0 : CALL kpoint_env_release(kpoint%kp_aux_env(ik)%kpoint_env)
425 : END DO
426 0 : DEALLOCATE (kpoint%kp_aux_env)
427 : NULLIFY (kpoint%kp_aux_env)
428 : END IF
429 :
430 2330 : IF (ASSOCIATED(kpoint%kp_sym)) THEN
431 10066 : DO ik = 1, SIZE(kpoint%kp_sym)
432 10066 : CALL kpoint_sym_release(kpoint%kp_sym(ik)%kpoint_sym)
433 : END DO
434 2330 : DEALLOCATE (kpoint%kp_sym)
435 : NULLIFY (kpoint%kp_sym)
436 : END IF
437 :
438 2330 : IF (ASSOCIATED(kpoint%atype)) THEN
439 2328 : DEALLOCATE (kpoint%atype)
440 : NULLIFY (kpoint%atype)
441 : END IF
442 2330 : IF (ASSOCIATED(kpoint%ibrot)) THEN
443 2322 : DEALLOCATE (kpoint%ibrot)
444 : NULLIFY (kpoint%ibrot)
445 : END IF
446 :
447 2330 : IF (ASSOCIATED(kpoint%kind_rotmat)) THEN
448 24446 : DO i = 1, SIZE(kpoint%kind_rotmat, 1)
449 46722 : DO j = 1, SIZE(kpoint%kind_rotmat, 2)
450 44400 : IF (ASSOCIATED(kpoint%kind_rotmat(i, j)%rmat)) THEN
451 5794 : DEALLOCATE (kpoint%kind_rotmat(i, j)%rmat)
452 5794 : NULLIFY (kpoint%kind_rotmat(i, j)%rmat)
453 : END IF
454 : END DO
455 : END DO
456 2322 : DEALLOCATE (kpoint%kind_rotmat)
457 : NULLIFY (kpoint%kind_rotmat)
458 : END IF
459 :
460 2330 : IF (ASSOCIATED(kpoint%scf_diis_buffer)) THEN
461 1962 : CALL qs_diis_b_release_kp(kpoint%scf_diis_buffer)
462 1962 : DEALLOCATE (kpoint%scf_diis_buffer)
463 : NULLIFY (kpoint%scf_diis_buffer)
464 : END IF
465 :
466 2330 : NULLIFY (kpoint%sab_nl)
467 2330 : NULLIFY (kpoint%sab_nl_nosym)
468 :
469 2330 : ALLOCATE (kpoint%cell_to_index(0:0, 0:0, 0:0))
470 9320 : kpoint%cell_to_index(:, :, :) = 1
471 :
472 2330 : ALLOCATE (kpoint%index_to_cell(0:0, 0:0))
473 6990 : kpoint%index_to_cell(:, :) = 0
474 :
475 2330 : IF (.NOT. ASSOCIATED(kpoint%wkp)) kpoint%nkp = 0
476 2330 : kpoint%nkp_groups = 0
477 6990 : kpoint%kp_range = 0
478 2330 : kpoint%iogrp = .FALSE.
479 :
480 2330 : END SUBROUTINE kpoint_reset_initialization
481 :
482 : ! **************************************************************************************************
483 : !> \brief Retrieve information from a kpoint environment
484 : !> \param kpoint The kpoint environment
485 : !> \param kp_scheme Type of kpoint grid
486 : !> \param nkp_grid Grid points
487 : !> \param kp_shift Shift of the grid
488 : !> \param symmetry use symmetry (atoms) to reduce kpoints
489 : !> \param verbose more output information
490 : !> \param full_grid don't reduce kpoints at all
491 : !> \param use_real_wfn real/complex wfn
492 : !> \param eps_geo accuracy of atom symmetry detection
493 : !> \param parallel_group_size kpoint group size
494 : !> \param kp_range kpoints distribution for local processor
495 : !> \param nkp number of kpoints
496 : !> \param xkp kpoint coordinates in units of b-vector
497 : !> \param wkp kpoint weights
498 : !> \param para_env 'global' parallel environment
499 : !> \param blacs_env_all BLACS env for the total environment
500 : !> \param para_env_kp parallel environment of the kpoint calculation
501 : !> \param para_env_inter_kp parallel environment between kpoints
502 : !> \param blacs_env BLACS env for the kpoint group
503 : !> \param kp_env Information for each kpoint
504 : !> \param kp_aux_env ...
505 : !> \param mpools FM matrix pools for kpoint groups
506 : !> \param iogrp this kpoint group has the IO processor
507 : !> \param nkp_groups number of kpoint groups
508 : !> \param kp_dist kpoints distribution on groups
509 : !> \param cell_to_index given a cell triple, returns the real space index
510 : !> \param index_to_cell ...
511 : !> \param sab_nl neighbourlist that defines real space matrices
512 : !> \param sab_nl_nosym neighbourlist that defines real space matrices, non-symmetric
513 : !> \param inversion_symmetry_only reduce kpoints only by inversion symmetry
514 : !> \param symmetry_backend k-point symmetry backend
515 : !> \param symmetry_reduction_method k-point symmetry reduction method
516 : !> \param gamma_centered ...
517 : !> \author JGH
518 : ! **************************************************************************************************
519 2064542 : SUBROUTINE get_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, &
520 : full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, &
521 : para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, &
522 : kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, &
523 : sab_nl, sab_nl_nosym, inversion_symmetry_only, symmetry_backend, &
524 : symmetry_reduction_method, gamma_centered)
525 : TYPE(kpoint_type), INTENT(IN) :: kpoint
526 : CHARACTER(LEN=*), OPTIONAL :: kp_scheme
527 : INTEGER, DIMENSION(3), OPTIONAL :: nkp_grid
528 : REAL(KIND=dp), DIMENSION(3), OPTIONAL :: kp_shift
529 : LOGICAL, OPTIONAL :: symmetry, verbose, full_grid, &
530 : use_real_wfn
531 : REAL(KIND=dp), OPTIONAL :: eps_geo
532 : INTEGER, OPTIONAL :: parallel_group_size
533 : INTEGER, DIMENSION(2), OPTIONAL :: kp_range
534 : INTEGER, OPTIONAL :: nkp
535 : REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: xkp
536 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: wkp
537 : TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env
538 : TYPE(cp_blacs_env_type), OPTIONAL, POINTER :: blacs_env_all
539 : TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env_kp, para_env_inter_kp
540 : TYPE(cp_blacs_env_type), OPTIONAL, POINTER :: blacs_env
541 : TYPE(kpoint_env_p_type), DIMENSION(:), OPTIONAL, &
542 : POINTER :: kp_env, kp_aux_env
543 : TYPE(qs_matrix_pools_type), OPTIONAL, POINTER :: mpools
544 : LOGICAL, OPTIONAL :: iogrp
545 : INTEGER, OPTIONAL :: nkp_groups
546 : INTEGER, DIMENSION(:, :), OPTIONAL, POINTER :: kp_dist
547 : INTEGER, DIMENSION(:, :, :), OPTIONAL, POINTER :: cell_to_index
548 : INTEGER, DIMENSION(:, :), OPTIONAL, POINTER :: index_to_cell
549 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
550 : OPTIONAL, POINTER :: sab_nl, sab_nl_nosym
551 : LOGICAL, OPTIONAL :: inversion_symmetry_only
552 : INTEGER, OPTIONAL :: symmetry_backend, &
553 : symmetry_reduction_method
554 : LOGICAL, OPTIONAL :: gamma_centered
555 :
556 6696 : IF (PRESENT(kp_scheme)) kp_scheme = kpoint%kp_scheme
557 2065670 : IF (PRESENT(nkp_grid)) nkp_grid = kpoint%nkp_grid
558 2064686 : IF (PRESENT(kp_shift)) kp_shift = kpoint%kp_shift
559 2064542 : IF (PRESENT(gamma_centered)) gamma_centered = kpoint%gamma_centered
560 2064542 : IF (PRESENT(symmetry)) symmetry = kpoint%symmetry
561 2064542 : IF (PRESENT(verbose)) verbose = kpoint%verbose
562 2064542 : IF (PRESENT(full_grid)) full_grid = kpoint%full_grid
563 2064542 : IF (PRESENT(inversion_symmetry_only)) inversion_symmetry_only = kpoint%inversion_symmetry_only
564 2064542 : IF (PRESENT(symmetry_backend)) symmetry_backend = kpoint%symmetry_backend
565 2064542 : IF (PRESENT(symmetry_reduction_method)) &
566 0 : symmetry_reduction_method = kpoint%symmetry_reduction_method
567 2064542 : IF (PRESENT(use_real_wfn)) use_real_wfn = kpoint%use_real_wfn
568 2064542 : IF (PRESENT(eps_geo)) eps_geo = kpoint%eps_geo
569 2064542 : IF (PRESENT(parallel_group_size)) parallel_group_size = kpoint%parallel_group_size
570 :
571 2064542 : IF (PRESENT(nkp)) nkp = kpoint%nkp
572 2064542 : IF (PRESENT(wkp)) wkp => kpoint%wkp
573 2064542 : IF (PRESENT(xkp)) xkp => kpoint%xkp
574 :
575 2064542 : IF (PRESENT(para_env)) para_env => kpoint%para_env
576 2064542 : IF (PRESENT(para_env_kp)) para_env_kp => kpoint%para_env_kp
577 2064542 : IF (PRESENT(para_env_inter_kp)) para_env_inter_kp => kpoint%para_env_inter_kp
578 2064542 : IF (PRESENT(blacs_env_all)) blacs_env_all => kpoint%blacs_env_all
579 2064542 : IF (PRESENT(blacs_env)) blacs_env => kpoint%blacs_env
580 :
581 2064542 : IF (PRESENT(iogrp)) iogrp = kpoint%iogrp
582 2353976 : IF (PRESENT(kp_range)) kp_range = kpoint%kp_range
583 2064542 : IF (PRESENT(nkp_groups)) nkp_groups = kpoint%nkp_groups
584 2064542 : IF (PRESENT(kp_dist)) kp_dist => kpoint%kp_dist
585 :
586 2064542 : IF (PRESENT(kp_env)) kp_env => kpoint%kp_env
587 2064542 : IF (PRESENT(kp_aux_env)) kp_aux_env => kpoint%kp_aux_env
588 2064542 : IF (PRESENT(mpools)) mpools => kpoint%mpools
589 :
590 2064542 : IF (PRESENT(cell_to_index)) cell_to_index => kpoint%cell_to_index
591 2064542 : IF (PRESENT(index_to_cell)) index_to_cell => kpoint%index_to_cell
592 2064542 : IF (PRESENT(sab_nl)) sab_nl => kpoint%sab_nl
593 2064542 : IF (PRESENT(sab_nl_nosym)) sab_nl_nosym => kpoint%sab_nl_nosym
594 :
595 2064542 : END SUBROUTINE get_kpoint_info
596 :
597 : ! **************************************************************************************************
598 : !> \brief Set information in a kpoint environment
599 : !> \param kpoint The kpoint environment
600 : !> \param kp_scheme Type of kpoint grid
601 : !> \param nkp_grid Grid points
602 : !> \param kp_shift Shift of the grid
603 : !> \param symmetry use symmetry (atoms) to reduce kpoints
604 : !> \param verbose more output information
605 : !> \param full_grid don't reduce kpoints at all
606 : !> \param use_real_wfn real/complex wfn
607 : !> \param eps_geo accuracy of atom symmetry detection
608 : !> \param parallel_group_size kpoint group size
609 : !> \param kp_range kpoints distribution for local processor
610 : !> \param nkp number of kpoints
611 : !> \param xkp kpoint coordinates
612 : !> \param wkp kpoint weights
613 : !> \param para_env 'global' parallel environment
614 : !> \param blacs_env_all BLACS env for the total environment
615 : !> \param para_env_kp parallel environment of the kpoint calculation
616 : !> \param para_env_inter_kp parallel environment between kpoints
617 : !> \param blacs_env BLACS env for the kpoint group
618 : !> \param kp_env Information for each kpoint
619 : !> \param kp_aux_env ...
620 : !> \param mpools FM matrix pools for kpoint groups
621 : !> \param iogrp this kpoint group has the IO processor
622 : !> \param nkp_groups number of kpoint groups
623 : !> \param kp_dist kpoints distribution on groups
624 : !> \param cell_to_index given a cell triple, returns the real space index
625 : !> \param index_to_cell ...
626 : !> \param sab_nl neighbourlist that defines real space matrices
627 : !> \param sab_nl_nosym neighbourlist that defines real space matrices
628 : !> \param inversion_symmetry_only reduce kpoints only by inversion symmetry
629 : !> \param symmetry_backend k-point symmetry backend
630 : !> \param symmetry_reduction_method k-point symmetry reduction method
631 : !> \param gamma_centered ...
632 : !> \author JGH
633 : ! **************************************************************************************************
634 5706 : SUBROUTINE set_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, &
635 : full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, &
636 : para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, &
637 : kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, &
638 : sab_nl, sab_nl_nosym, inversion_symmetry_only, symmetry_backend, &
639 : symmetry_reduction_method, gamma_centered)
640 : TYPE(kpoint_type), INTENT(INOUT) :: kpoint
641 : CHARACTER(LEN=*), OPTIONAL :: kp_scheme
642 : INTEGER, DIMENSION(3), OPTIONAL :: nkp_grid
643 : REAL(KIND=dp), DIMENSION(3), OPTIONAL :: kp_shift
644 : LOGICAL, OPTIONAL :: symmetry, verbose, full_grid, &
645 : use_real_wfn
646 : REAL(KIND=dp), OPTIONAL :: eps_geo
647 : INTEGER, OPTIONAL :: parallel_group_size
648 : INTEGER, DIMENSION(2), OPTIONAL :: kp_range
649 : INTEGER, OPTIONAL :: nkp
650 : REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: xkp
651 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: wkp
652 : TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env
653 : TYPE(cp_blacs_env_type), OPTIONAL, POINTER :: blacs_env_all
654 : TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env_kp, para_env_inter_kp
655 : TYPE(cp_blacs_env_type), OPTIONAL, POINTER :: blacs_env
656 : TYPE(kpoint_env_p_type), DIMENSION(:), OPTIONAL, &
657 : POINTER :: kp_env, kp_aux_env
658 : TYPE(qs_matrix_pools_type), OPTIONAL, POINTER :: mpools
659 : LOGICAL, OPTIONAL :: iogrp
660 : INTEGER, OPTIONAL :: nkp_groups
661 : INTEGER, DIMENSION(:, :), OPTIONAL, POINTER :: kp_dist
662 : INTEGER, DIMENSION(:, :, :), OPTIONAL, POINTER :: cell_to_index
663 : INTEGER, DIMENSION(:, :), OPTIONAL, POINTER :: index_to_cell
664 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
665 : OPTIONAL, POINTER :: sab_nl, sab_nl_nosym
666 : LOGICAL, OPTIONAL :: inversion_symmetry_only
667 : INTEGER, OPTIONAL :: symmetry_backend, &
668 : symmetry_reduction_method
669 : LOGICAL, OPTIONAL :: gamma_centered
670 :
671 0 : IF (PRESENT(kp_scheme)) kpoint%kp_scheme = kp_scheme
672 5706 : IF (PRESENT(nkp_grid)) kpoint%nkp_grid = nkp_grid
673 5706 : IF (PRESENT(kp_shift)) kpoint%kp_shift = kp_shift
674 5706 : IF (PRESENT(gamma_centered)) kpoint%gamma_centered = gamma_centered
675 5706 : IF (PRESENT(symmetry)) kpoint%symmetry = symmetry
676 5706 : IF (PRESENT(verbose)) kpoint%verbose = verbose
677 5706 : IF (PRESENT(full_grid)) kpoint%full_grid = full_grid
678 5706 : IF (PRESENT(inversion_symmetry_only)) kpoint%inversion_symmetry_only = inversion_symmetry_only
679 5706 : IF (PRESENT(symmetry_backend)) kpoint%symmetry_backend = symmetry_backend
680 5706 : IF (PRESENT(symmetry_reduction_method)) &
681 0 : kpoint%symmetry_reduction_method = symmetry_reduction_method
682 5706 : IF (PRESENT(use_real_wfn)) kpoint%use_real_wfn = use_real_wfn
683 5706 : IF (PRESENT(eps_geo)) kpoint%eps_geo = eps_geo
684 5706 : IF (PRESENT(parallel_group_size)) kpoint%parallel_group_size = parallel_group_size
685 :
686 5706 : IF (PRESENT(nkp)) kpoint%nkp = nkp
687 5706 : IF (PRESENT(wkp)) kpoint%wkp => wkp
688 5706 : IF (PRESENT(xkp)) kpoint%xkp => xkp
689 :
690 5706 : IF (PRESENT(para_env)) kpoint%para_env => para_env
691 5706 : IF (PRESENT(para_env_kp)) kpoint%para_env_kp => para_env_kp
692 5706 : IF (PRESENT(para_env_inter_kp)) kpoint%para_env_inter_kp => para_env_inter_kp
693 5706 : IF (PRESENT(blacs_env_all)) kpoint%blacs_env_all => blacs_env_all
694 5706 : IF (PRESENT(blacs_env)) kpoint%blacs_env => blacs_env
695 :
696 5706 : IF (PRESENT(iogrp)) kpoint%iogrp = iogrp
697 5706 : IF (PRESENT(kp_range)) kpoint%kp_range = kp_range
698 5706 : IF (PRESENT(nkp_groups)) kpoint%nkp_groups = nkp_groups
699 5706 : IF (PRESENT(kp_dist)) kpoint%kp_dist => kp_dist
700 :
701 5706 : IF (PRESENT(kp_env)) kpoint%kp_env => kp_env
702 0 : IF (PRESENT(kp_env)) kpoint%kp_aux_env => kp_aux_env
703 5706 : IF (PRESENT(mpools)) kpoint%mpools => mpools
704 5706 : IF (PRESENT(sab_nl)) kpoint%sab_nl => sab_nl
705 5706 : IF (PRESENT(sab_nl_nosym)) kpoint%sab_nl_nosym => sab_nl_nosym
706 :
707 5706 : IF (PRESENT(cell_to_index)) THEN
708 0 : IF (ASSOCIATED(kpoint%cell_to_index)) DEALLOCATE (kpoint%cell_to_index)
709 0 : kpoint%cell_to_index => cell_to_index
710 : END IF
711 :
712 5706 : IF (PRESENT(index_to_cell)) THEN
713 0 : IF (ASSOCIATED(kpoint%index_to_cell)) DEALLOCATE (kpoint%index_to_cell)
714 0 : kpoint%index_to_cell => index_to_cell
715 : END IF
716 :
717 5706 : END SUBROUTINE set_kpoint_info
718 :
719 : ! **************************************************************************************************
720 : !> \brief Read the kpoint input section
721 : !> \param kpoint The kpoint environment
722 : !> \param kpoint_section The input section
723 : !> \param a_vec ...
724 : !> \param cell ...
725 : !> \author JGH
726 : ! **************************************************************************************************
727 8442 : SUBROUTINE read_kpoint_section(kpoint, kpoint_section, a_vec, cell)
728 : TYPE(kpoint_type), INTENT(INOUT) :: kpoint
729 : TYPE(section_vals_type), POINTER :: kpoint_section
730 : REAL(KIND=dp), DIMENSION(3, 3), INTENT(IN) :: a_vec
731 : TYPE(cell_type), OPTIONAL, POINTER :: cell
732 :
733 : REAL(KIND=dp), PARAMETER :: eps_cell = 1.0E-14_dp
734 :
735 : CHARACTER(LEN=default_string_length) :: ustr
736 : CHARACTER(LEN=default_string_length), &
737 8442 : DIMENSION(:), POINTER :: tmpstringlist
738 : INTEGER :: i, n_rep, nval, wfntype
739 : LOGICAL :: available, backend_explicit, &
740 : non_lower_triangular_cell, &
741 : non_orthogonal_cell, reduction_explicit
742 : REAL(KIND=dp) :: ff
743 : REAL(KIND=dp), DIMENSION(3, 3) :: cart_a_vec
744 8442 : REAL(KIND=dp), DIMENSION(:), POINTER :: reallist
745 :
746 8442 : CALL section_vals_get(kpoint_section, explicit=available)
747 8442 : cart_a_vec(:, :) = a_vec(:, :)
748 8442 : IF (PRESENT(cell)) THEN
749 8442 : IF (ASSOCIATED(cell)) THEN
750 9138 : IF (cell%input_cell_canonicalized) cart_a_vec(:, :) = cell%input_hmat(:, :)
751 : END IF
752 : END IF
753 :
754 8442 : IF (available) THEN
755 500 : CALL section_vals_val_get(kpoint_section, "SCHEME", c_vals=tmpstringlist)
756 500 : nval = SIZE(tmpstringlist)
757 500 : CPASSERT(nval >= 1)
758 500 : kpoint%kp_scheme = tmpstringlist(1)
759 500 : CALL uppercase(kpoint%kp_scheme)
760 :
761 : ! SCHEME [None, Gamma, Monkhorst-Pack, MacDonald, General]
762 426 : SELECT CASE (kpoint%kp_scheme)
763 : CASE ("NONE")
764 : ! do nothing
765 : CASE ("GAMMA")
766 : ! do nothing
767 : CASE ("MONKHORST-PACK")
768 426 : CPASSERT(nval >= 4)
769 1704 : DO i = 2, 4
770 1704 : READ (tmpstringlist(i), *) kpoint%nkp_grid(i - 1)
771 : END DO
772 : CASE ("MACDONALD")
773 8 : CPASSERT(nval >= 7)
774 32 : DO i = 2, 4
775 32 : READ (tmpstringlist(i), *) kpoint%nkp_grid(i - 1)
776 : END DO
777 32 : DO i = 5, 7
778 32 : READ (tmpstringlist(i), *) kpoint%kp_shift(i - 4)
779 : END DO
780 : CASE ("GENERAL")
781 18 : CALL section_vals_val_get(kpoint_section, "UNITS", c_val=ustr)
782 18 : CALL uppercase(ustr)
783 18 : CALL section_vals_val_get(kpoint_section, "KPOINT", n_rep_val=n_rep)
784 18 : kpoint%nkp = n_rep
785 18 : IF (ASSOCIATED(kpoint%xkp_input)) THEN
786 0 : DEALLOCATE (kpoint%xkp_input)
787 : NULLIFY (kpoint%xkp_input)
788 : END IF
789 18 : IF (ASSOCIATED(kpoint%wkp_input)) THEN
790 0 : DEALLOCATE (kpoint%wkp_input)
791 : NULLIFY (kpoint%wkp_input)
792 : END IF
793 90 : ALLOCATE (kpoint%xkp(3, n_rep), kpoint%wkp(n_rep))
794 126 : DO i = 1, n_rep
795 : CALL section_vals_val_get(kpoint_section, "KPOINT", i_rep_val=i, &
796 108 : r_vals=reallist)
797 108 : nval = SIZE(reallist)
798 108 : CPASSERT(nval >= 4)
799 108 : SELECT CASE (ustr)
800 : CASE ("B_VECTOR")
801 432 : kpoint%xkp(1:3, i) = reallist(1:3)
802 : CASE ("CART_ANGSTROM")
803 : kpoint%xkp(1:3, i) = (reallist(1)*cart_a_vec(1, 1:3) + &
804 : reallist(2)*cart_a_vec(2, 1:3) + &
805 0 : reallist(3)*cart_a_vec(3, 1:3))/twopi*angstrom
806 : CASE ("CART_BOHR")
807 : kpoint%xkp(1:3, i) = (reallist(1)*cart_a_vec(1, 1:3) + &
808 : reallist(2)*cart_a_vec(2, 1:3) + &
809 0 : reallist(3)*cart_a_vec(3, 1:3))/twopi
810 : CASE DEFAULT
811 108 : CPABORT("Unknown Unit for kpoint definition")
812 : END SELECT
813 126 : kpoint%wkp(i) = reallist(4)
814 : END DO
815 126 : ff = 1.0_dp/SUM(kpoint%wkp(:))
816 126 : kpoint%wkp(:) = ff*kpoint%wkp(:)
817 54 : ALLOCATE (kpoint%xkp_input(3, n_rep), kpoint%wkp_input(n_rep))
818 450 : kpoint%xkp_input(1:3, 1:n_rep) = kpoint%xkp(1:3, 1:n_rep)
819 144 : kpoint%wkp_input(1:n_rep) = kpoint%wkp(1:n_rep)
820 : CASE DEFAULT
821 500 : CPABORT("")
822 : END SELECT
823 :
824 500 : CALL section_vals_val_get(kpoint_section, "GAMMA_CENTERED", l_val=kpoint%gamma_centered)
825 500 : IF (kpoint%gamma_centered .AND. kpoint%kp_scheme /= "MONKHORST-PACK" .AND. &
826 : kpoint%kp_scheme /= "MACDONALD") THEN
827 : CALL cp_abort(__LOCATION__, &
828 0 : "KPOINTS%GAMMA_CENTERED is only supported with SCHEME MONKHORST-PACK or MACDONALD")
829 : END IF
830 :
831 500 : CALL section_vals_val_get(kpoint_section, "SYMMETRY", l_val=kpoint%symmetry)
832 500 : CALL section_vals_val_get(kpoint_section, "WAVEFUNCTIONS", i_val=wfntype)
833 500 : CALL section_vals_val_get(kpoint_section, "VERBOSE", l_val=kpoint%verbose)
834 500 : CALL section_vals_val_get(kpoint_section, "FULL_GRID", l_val=kpoint%full_grid)
835 : CALL section_vals_val_get(kpoint_section, "INVERSION_SYMMETRY_ONLY", &
836 500 : l_val=kpoint%inversion_symmetry_only)
837 : CALL section_vals_val_get(kpoint_section, "SYMMETRY_BACKEND", &
838 500 : i_val=kpoint%symmetry_backend, explicit=backend_explicit)
839 : CALL section_vals_val_get(kpoint_section, "SYMMETRY_REDUCTION_METHOD", &
840 500 : i_val=kpoint%symmetry_reduction_method, explicit=reduction_explicit)
841 500 : CALL resolve_kpoint_symmetry_settings(kpoint, backend_explicit, reduction_explicit)
842 500 : CALL section_vals_val_get(kpoint_section, "EPS_SYMMETRY", r_val=kpoint%eps_geo)
843 : IF ((kpoint%kp_scheme == "MONKHORST-PACK" .OR. kpoint%kp_scheme == "MACDONALD") .AND. &
844 500 : kpoint%symmetry .AND. .NOT. kpoint%full_grid .AND. &
845 : .NOT. kpoint%inversion_symmetry_only) THEN
846 : non_lower_triangular_cell = (ABS(a_vec(2, 1)) > eps_cell) .OR. &
847 : (ABS(a_vec(3, 1)) > eps_cell) .OR. &
848 424 : (ABS(a_vec(3, 2)) > eps_cell)
849 : non_orthogonal_cell = (ABS(DOT_PRODUCT(a_vec(:, 1), a_vec(:, 2))) > eps_cell) .OR. &
850 : (ABS(DOT_PRODUCT(a_vec(:, 1), a_vec(:, 3))) > eps_cell) .OR. &
851 2140 : (ABS(DOT_PRODUCT(a_vec(:, 2), a_vec(:, 3))) > eps_cell)
852 : ! The full point-group density-matrix transform is not charge-conserving for skew cells.
853 : IF (non_orthogonal_cell) THEN
854 20 : kpoint%inversion_symmetry_only = .TRUE.
855 : CALL cp_warn(__LOCATION__, &
856 : "Full atomic k-point symmetry was requested for a non-orthogonal cell. "// &
857 : "Falling back to KPOINTS%INVERSION_SYMMETRY_ONLY because the full "// &
858 : "point-group density-matrix symmetrization is only charge-conserving for "// &
859 20 : "orthogonal cell matrices.")
860 194 : ELSEIF (non_lower_triangular_cell) THEN
861 4 : kpoint%inversion_symmetry_only = .TRUE.
862 : CALL cp_warn(__LOCATION__, &
863 : "Full atomic k-point symmetry was requested for a cell matrix that does "// &
864 : "not follow the CP2K lower-triangular convention. Falling back to "// &
865 : "KPOINTS%INVERSION_SYMMETRY_ONLY. Use ABC/ALPHA_BETA_GAMMA or canonical "// &
866 4 : "A/B/C vectors to enable full point-group k-point reduction.")
867 : END IF
868 : END IF
869 : CALL section_vals_val_get(kpoint_section, "PARALLEL_GROUP_SIZE", &
870 500 : i_val=kpoint%parallel_group_size)
871 16 : SELECT CASE (wfntype)
872 : CASE (use_real_wfn)
873 16 : kpoint%use_real_wfn = .TRUE.
874 : CASE (use_complex_wfn)
875 484 : kpoint%use_real_wfn = .FALSE.
876 : CASE DEFAULT
877 500 : CPABORT("")
878 : END SELECT
879 :
880 : ELSE
881 7942 : kpoint%kp_scheme = "NONE"
882 : END IF
883 :
884 8442 : END SUBROUTINE read_kpoint_section
885 :
886 : ! **************************************************************************************************
887 : !> \brief Resolve legacy and backend k-point symmetry settings
888 : !> \param kpoint ...
889 : !> \param backend_explicit whether SYMMETRY_BACKEND was given
890 : !> \param reduction_explicit whether SYMMETRY_REDUCTION_METHOD was given
891 : ! **************************************************************************************************
892 500 : SUBROUTINE resolve_kpoint_symmetry_settings(kpoint, backend_explicit, reduction_explicit)
893 : TYPE(kpoint_type), INTENT(INOUT) :: kpoint
894 : LOGICAL, INTENT(IN) :: backend_explicit, reduction_explicit
895 :
896 500 : IF (backend_explicit .AND. .NOT. reduction_explicit) THEN
897 136 : SELECT CASE (kpoint%symmetry_backend)
898 : CASE (use_k290_kpoint_backend)
899 30 : kpoint%symmetry_reduction_method = use_k290_kpoint_symmetry
900 : CASE (use_spglib_kpoint_backend)
901 76 : kpoint%symmetry_reduction_method = use_spglib_kpoint_symmetry
902 : CASE DEFAULT
903 106 : CPABORT("Unknown k-point symmetry backend")
904 : END SELECT
905 : END IF
906 :
907 500 : IF (kpoint%symmetry_backend == use_spglib_kpoint_backend .AND. &
908 : kpoint%symmetry_reduction_method /= use_spglib_kpoint_symmetry) THEN
909 : CALL cp_abort(__LOCATION__, &
910 0 : "SYMMETRY_BACKEND SPGLIB requires SYMMETRY_REDUCTION_METHOD SPGLIB")
911 : END IF
912 :
913 500 : END SUBROUTINE resolve_kpoint_symmetry_settings
914 :
915 : ! **************************************************************************************************
916 : !> \brief Write information on the kpoints to output
917 : !> \param kpoint The kpoint environment
918 : !> \param iounit output unit
919 : !> \param dft_section DFT section information for output unit
920 : !> \author JGH
921 : ! **************************************************************************************************
922 8440 : SUBROUTINE write_kpoint_info(kpoint, iounit, dft_section)
923 : TYPE(kpoint_type), INTENT(IN) :: kpoint
924 : INTEGER, INTENT(IN), OPTIONAL :: iounit
925 : TYPE(section_vals_type), INTENT(IN), OPTIONAL :: dft_section
926 :
927 : INTEGER :: i, punit
928 : TYPE(cp_logger_type), POINTER :: logger
929 :
930 8440 : NULLIFY (logger)
931 8440 : logger => cp_get_default_logger()
932 :
933 8440 : IF (PRESENT(dft_section)) THEN
934 8428 : punit = cp_print_key_unit_nr(logger, dft_section, "PRINT%KPOINTS", extension=".Log")
935 12 : ELSE IF (PRESENT(iounit)) THEN
936 12 : punit = iounit
937 : ELSE
938 0 : punit = cp_logger_get_default_unit_nr(logger)
939 : END IF
940 :
941 8440 : IF (punit > 0) THEN
942 :
943 2055 : IF (kpoint%kp_scheme /= "NONE") THEN
944 91 : WRITE (punit, '(/," ",79("*"),/,T37,A,/," ",79("*"))') "Kpoints"
945 : END IF
946 1 : SELECT CASE (kpoint%kp_scheme)
947 : CASE ("NONE")
948 : ! be silent
949 : CASE ("GAMMA")
950 1 : WRITE (punit, '(A,T57,A)') ' BRILLOUIN|', ' Gamma-point calculation'
951 : CASE ("MONKHORST-PACK")
952 90 : WRITE (punit, '(A,T61,A20)') ' BRILLOUIN| K-point scheme ', ' Monkhorst-Pack'
953 90 : WRITE (punit, '(A,T66,3I5)') ' BRILLOUIN| K-Point grid', kpoint%nkp_grid
954 90 : IF (kpoint%gamma_centered) THEN
955 1 : WRITE (punit, '(A,T76,A)') ' BRILLOUIN| Gamma-centered k-point mesh', ' ON'
956 : END IF
957 : WRITE (punit, '(A,T66,G15.6)') &
958 90 : ' BRILLOUIN| K-point symmetry accuracy', kpoint%eps_geo
959 : CASE ("MACDONALD")
960 0 : WRITE (punit, '(A,T71,A10)') ' BRILLOUIN| K-point scheme ', ' MacDonald'
961 0 : WRITE (punit, '(A,T66,3I5)') ' BRILLOUIN| K-Point grid', kpoint%nkp_grid
962 0 : WRITE (punit, '(A,T51,3F10.4)') ' BRILLOUIN| K-Point shift', kpoint%kp_shift
963 0 : IF (kpoint%gamma_centered) THEN
964 0 : WRITE (punit, '(A,T76,A)') ' BRILLOUIN| Gamma-centered k-point mesh', ' ON'
965 : END IF
966 : WRITE (punit, '(A,T66,G15.6)') &
967 0 : ' BRILLOUIN| K-point symmetry accuracy', kpoint%eps_geo
968 : CASE ("GENERAL")
969 0 : WRITE (punit, '(A,T71,A10)') ' BRILLOUIN| K-point scheme ', ' General'
970 : CASE DEFAULT
971 2055 : CPABORT("")
972 : END SELECT
973 2055 : IF (kpoint%kp_scheme /= "NONE") THEN
974 91 : IF (kpoint%symmetry) THEN
975 59 : WRITE (punit, '(A,T76,A)') ' BRILLOUIN| K-Point point group symmetrization', ' ON'
976 : ELSE
977 32 : WRITE (punit, '(A,T76,A)') ' BRILLOUIN| K-Point point group symmetrization', ' OFF'
978 : END IF
979 91 : IF (kpoint%inversion_symmetry_only) THEN
980 5 : WRITE (punit, '(A,T76,A)') ' BRILLOUIN| Restrict symmetry to inversion', ' ON'
981 : END IF
982 : IF (kpoint%symmetry .AND. .NOT. kpoint%full_grid .AND. &
983 91 : .NOT. kpoint%inversion_symmetry_only .AND. &
984 : (kpoint%kp_scheme == "MONKHORST-PACK" .OR. kpoint%kp_scheme == "MACDONALD" .OR. &
985 : kpoint%kp_scheme == "GENERAL")) THEN
986 58 : SELECT CASE (kpoint%symmetry_backend)
987 : CASE (use_k290_kpoint_backend)
988 21 : WRITE (punit, '(A,T73,A)') ' BRILLOUIN| Symmetry backend', ' K290'
989 : CASE (use_spglib_kpoint_backend)
990 16 : WRITE (punit, '(A,T73,A)') ' BRILLOUIN| Symmetry backend', ' SPGLIB'
991 : CASE DEFAULT
992 37 : CPABORT("Unknown k-point symmetry backend")
993 : END SELECT
994 58 : SELECT CASE (kpoint%symmetry_reduction_method)
995 : CASE (use_k290_kpoint_symmetry)
996 21 : WRITE (punit, '(A,T73,A)') ' BRILLOUIN| Symmetry reduction method', ' K290'
997 : CASE (use_spglib_kpoint_symmetry)
998 16 : WRITE (punit, '(A,T73,A)') ' BRILLOUIN| Symmetry reduction method', ' SPGLIB'
999 : CASE DEFAULT
1000 91 : CPABORT("Unknown k-point symmetry reduction method")
1001 : END SELECT
1002 : END IF
1003 91 : IF (kpoint%use_real_wfn) THEN
1004 0 : WRITE (punit, '(A,T76,A)') ' BRILLOUIN| Wavefunction type', ' REAL'
1005 : ELSE
1006 91 : WRITE (punit, '(A,T73,A)') ' BRILLOUIN| Wavefunction type', ' COMPLEX'
1007 : END IF
1008 91 : IF (kpoint%full_grid) THEN
1009 17 : WRITE (punit, '(A,T76,A)') ' BRILLOUIN| Use full k-point grid '
1010 : END IF
1011 91 : IF (kpoint%kp_scheme /= "GAMMA") THEN
1012 90 : WRITE (punit, '(A,T71,I10)') ' BRILLOUIN| List of Kpoints [2 Pi/Bohr]', kpoint%nkp
1013 : WRITE (punit, '(A,T30,A,T48,A,T63,A,T78,A)') &
1014 90 : ' BRILLOUIN| Number ', 'Weight', 'X', 'Y', 'Z'
1015 822 : DO i = 1, kpoint%nkp
1016 732 : WRITE (punit, '(A,I5,3X,4F15.5)') ' BRILLOUIN| ', i, kpoint%wkp(i), &
1017 1554 : kpoint%xkp(1, i), kpoint%xkp(2, i), kpoint%xkp(3, i)
1018 : END DO
1019 : END IF
1020 91 : WRITE (punit, '(" ",79("*"))')
1021 : END IF
1022 :
1023 : END IF
1024 :
1025 8440 : IF (PRESENT(dft_section)) THEN
1026 8428 : CALL cp_print_key_finished_output(punit, logger, dft_section, "PRINT%KPOINTS")
1027 : END IF
1028 :
1029 8440 : END SUBROUTINE write_kpoint_info
1030 :
1031 : ! **************************************************************************************************
1032 : !> \brief Create a single kpoint environment
1033 : !> \param kp_env Single kpoint environment
1034 : !> \author JGH
1035 : ! **************************************************************************************************
1036 6782 : SUBROUTINE kpoint_env_create(kp_env)
1037 : TYPE(kpoint_env_type), POINTER :: kp_env
1038 :
1039 6782 : CPASSERT(.NOT. ASSOCIATED(kp_env))
1040 :
1041 33910 : ALLOCATE (kp_env)
1042 :
1043 6782 : kp_env%nkpoint = 0
1044 : kp_env%wkp = 0.0_dp
1045 27128 : kp_env%xkp = 0.0_dp
1046 6782 : kp_env%is_local = .FALSE.
1047 :
1048 6782 : NULLIFY (kp_env%mos)
1049 6782 : NULLIFY (kp_env%pmat)
1050 6782 : NULLIFY (kp_env%wmat)
1051 6782 : NULLIFY (kp_env%smat)
1052 6782 : NULLIFY (kp_env%amat)
1053 :
1054 6782 : END SUBROUTINE kpoint_env_create
1055 :
1056 : ! **************************************************************************************************
1057 : !> \brief Release a single kpoint environment
1058 : !> \param kp_env Single kpoint environment
1059 : !> \author JGH
1060 : ! **************************************************************************************************
1061 6782 : SUBROUTINE kpoint_env_release(kp_env)
1062 : TYPE(kpoint_env_type), POINTER :: kp_env
1063 :
1064 : INTEGER :: ic, is
1065 :
1066 6782 : IF (ASSOCIATED(kp_env)) THEN
1067 :
1068 6782 : IF (ASSOCIATED(kp_env%mos)) THEN
1069 13864 : DO is = 1, SIZE(kp_env%mos, 2)
1070 28008 : DO ic = 1, SIZE(kp_env%mos, 1)
1071 21226 : CALL deallocate_mo_set(kp_env%mos(ic, is))
1072 : END DO
1073 : END DO
1074 6782 : DEALLOCATE (kp_env%mos)
1075 : END IF
1076 :
1077 6782 : CALL cp_fm_release(kp_env%pmat)
1078 6782 : CALL cp_fm_release(kp_env%wmat)
1079 6782 : CALL cp_fm_release(kp_env%smat)
1080 6782 : CALL cp_fm_release(kp_env%amat)
1081 :
1082 6782 : CALL cp_fm_release(kp_env%shalf)
1083 6782 : CALL cp_cfm_release(kp_env%cshalf)
1084 :
1085 6782 : DEALLOCATE (kp_env)
1086 :
1087 : END IF
1088 :
1089 6782 : END SUBROUTINE kpoint_env_release
1090 :
1091 : ! **************************************************************************************************
1092 : !> \brief Get information from a single kpoint environment
1093 : !> \param kpoint_env Single kpoint environment
1094 : !> \param nkpoint Index of kpoint
1095 : !> \param wkp Weight of kpoint
1096 : !> \param xkp Coordinates of kpoint
1097 : !> \param is_local Is this kpoint local (single cpu group)
1098 : !> \param mos MOs of this kpoint
1099 : !> \author JGH
1100 : ! **************************************************************************************************
1101 7588 : SUBROUTINE get_kpoint_env(kpoint_env, nkpoint, wkp, xkp, is_local, mos)
1102 : TYPE(kpoint_env_type), INTENT(IN) :: kpoint_env
1103 : INTEGER, OPTIONAL :: nkpoint
1104 : REAL(KIND=dp), OPTIONAL :: wkp
1105 : REAL(KIND=dp), DIMENSION(3), OPTIONAL :: xkp
1106 : LOGICAL, OPTIONAL :: is_local
1107 : TYPE(mo_set_type), DIMENSION(:, :), OPTIONAL, &
1108 : POINTER :: mos
1109 :
1110 7588 : IF (PRESENT(nkpoint)) nkpoint = kpoint_env%nkpoint
1111 7588 : IF (PRESENT(wkp)) wkp = kpoint_env%wkp
1112 7588 : IF (PRESENT(xkp)) xkp = kpoint_env%xkp
1113 7588 : IF (PRESENT(is_local)) is_local = kpoint_env%is_local
1114 7588 : IF (PRESENT(mos)) mos => kpoint_env%mos
1115 :
1116 7588 : END SUBROUTINE get_kpoint_env
1117 :
1118 : ! **************************************************************************************************
1119 : !> \brief Create a single kpoint symmetry environment
1120 : !> \param kp_sym ...
1121 : !> \author JGH
1122 : ! **************************************************************************************************
1123 22632 : SUBROUTINE kpoint_sym_create(kp_sym)
1124 : TYPE(kpoint_sym_type), POINTER :: kp_sym
1125 :
1126 22632 : CPASSERT(.NOT. ASSOCIATED(kp_sym))
1127 :
1128 22632 : ALLOCATE (kp_sym)
1129 :
1130 22632 : kp_sym%nwght = 0
1131 22632 : kp_sym%nwred = 0
1132 : kp_sym%apply_symmetry = .FALSE.
1133 :
1134 : NULLIFY (kp_sym%rot)
1135 : NULLIFY (kp_sym%xkp)
1136 : NULLIFY (kp_sym%rotp)
1137 : NULLIFY (kp_sym%f0)
1138 : NULLIFY (kp_sym%fcell)
1139 : NULLIFY (kp_sym%kgphase)
1140 :
1141 22632 : END SUBROUTINE kpoint_sym_create
1142 :
1143 : ! **************************************************************************************************
1144 : !> \brief Release a single kpoint symmetry environment
1145 : !> \param kp_sym ...
1146 : !> \author JGH
1147 : ! **************************************************************************************************
1148 22632 : SUBROUTINE kpoint_sym_release(kp_sym)
1149 : TYPE(kpoint_sym_type), POINTER :: kp_sym
1150 :
1151 22632 : IF (ASSOCIATED(kp_sym)) THEN
1152 :
1153 22632 : IF (ASSOCIATED(kp_sym%rot)) THEN
1154 1090 : DEALLOCATE (kp_sym%rot)
1155 : END IF
1156 22632 : IF (ASSOCIATED(kp_sym%xkp)) THEN
1157 1090 : DEALLOCATE (kp_sym%xkp)
1158 : END IF
1159 22632 : IF (ASSOCIATED(kp_sym%f0)) THEN
1160 1090 : DEALLOCATE (kp_sym%f0)
1161 : END IF
1162 22632 : IF (ASSOCIATED(kp_sym%fcell)) THEN
1163 1090 : DEALLOCATE (kp_sym%fcell)
1164 : END IF
1165 22632 : IF (ASSOCIATED(kp_sym%kgphase)) THEN
1166 1090 : DEALLOCATE (kp_sym%kgphase)
1167 : END IF
1168 22632 : IF (ASSOCIATED(kp_sym%rotp)) THEN
1169 1090 : DEALLOCATE (kp_sym%rotp)
1170 : END IF
1171 :
1172 22632 : DEALLOCATE (kp_sym)
1173 :
1174 : END IF
1175 :
1176 22632 : END SUBROUTINE kpoint_sym_release
1177 :
1178 : ! **************************************************************************************************
1179 :
1180 0 : END MODULE kpoint_types
|