(git:374b731)
Loading...
Searching...
No Matches
qs_loc_main.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 Driver for the localization that should be general
10!> for all the methods available and all the definition of the
11!> spread functional
12!> Write centers, spread and cubes only if required and for the
13!> selected states
14!> The localized functions are copied in the standard mos array
15!> for the next use
16!> \par History
17!> 01.2008 Teodoro Laino [tlaino] - University of Zurich
18!> - Merging the two localization codes and updating to new structures
19!> 04.2023 JGH Code isolation and refactoring
20!> \author MI (04.2005)
21! **************************************************************************************************
24 USE cell_types, ONLY: cell_type
33 USE cp_fm_types, ONLY: &
36 USE dbcsr_api, ONLY: dbcsr_create,&
37 dbcsr_p_type,&
38 dbcsr_set,&
39 dbcsr_type,&
40 dbcsr_type_symmetric
41 USE input_constants, ONLY: &
48 USE kinds, ONLY: default_string_length,&
49 dp
60 USE qs_loc_types, ONLY: get_qs_loc_env,&
65 USE qs_mo_types, ONLY: get_mo_set,&
68#include "./base/base_uses.f90"
69
70 IMPLICIT NONE
71
72 PRIVATE
73
74! *** Global parameters ***
75
76 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_loc_main'
77
78! *** Public ***
79 PUBLIC :: qs_loc_driver
80
81CONTAINS
82
83! **************************************************************************************************
84!> \brief set up the calculation of localized orbitals
85!> \param qs_env ...
86!> \param qs_loc_env ...
87!> \param print_loc_section ...
88!> \param myspin ...
89!> \param ext_mo_coeff ...
90!> \par History
91!> 04.2005 created [MI]
92!> 04.2023 refactored [JGH]
93!> \author MI
94! **************************************************************************************************
95 SUBROUTINE qs_loc_driver(qs_env, qs_loc_env, print_loc_section, myspin, ext_mo_coeff)
96
97 TYPE(qs_environment_type), POINTER :: qs_env
98 TYPE(qs_loc_env_type), POINTER :: qs_loc_env
99 TYPE(section_vals_type), POINTER :: print_loc_section
100 INTEGER, INTENT(IN) :: myspin
101 TYPE(cp_fm_type), INTENT(IN), OPTIONAL, TARGET :: ext_mo_coeff
102
103 CHARACTER(len=*), PARAMETER :: routinen = 'qs_loc_driver'
104
105 INTEGER :: dim_op, handle, i, imo, imoloc, j, lb, &
106 loc_method, nao, nmosub, restricted, ub
107 INTEGER, DIMENSION(:), POINTER :: ivec
108 LOGICAL, SAVE :: first_time = .true.
109 REAL(dp), DIMENSION(6) :: weights
110 REAL(kind=dp), DIMENSION(:, :), POINTER :: vecbuffer
111 TYPE(cell_type), POINTER :: cell
112 TYPE(cp_fm_struct_type), POINTER :: tmp_fm_struct
113 TYPE(cp_fm_type), DIMENSION(:), POINTER :: moloc_coeff
114 TYPE(cp_fm_type), DIMENSION(:, :), POINTER :: op_fm_set
115 TYPE(cp_fm_type), POINTER :: locorb
116 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: op_sm_set
117 TYPE(dft_control_type), POINTER :: dft_control
118 TYPE(localized_wfn_control_type), POINTER :: localized_wfn_control
119 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
120 TYPE(mp_para_env_type), POINTER :: para_env
121 TYPE(section_vals_type), POINTER :: input, low_spin_roks_section
122
123 CALL timeset(routinen, handle)
124 NULLIFY (para_env, mos, dft_control)
125 NULLIFY (cell, localized_wfn_control, moloc_coeff, op_sm_set, op_fm_set)
126 qs_loc_env%first_time = first_time
127 qs_loc_env%target_time = qs_env%target_time
128 qs_loc_env%start_time = qs_env%start_time
129
130 CALL get_qs_loc_env(qs_loc_env=qs_loc_env, &
131 localized_wfn_control=localized_wfn_control, &
132 moloc_coeff=moloc_coeff, op_sm_set=op_sm_set, op_fm_set=op_fm_set, cell=cell, &
133 weights=weights, dim_op=dim_op)
134
135 CALL get_qs_env(qs_env=qs_env, dft_control=dft_control, &
136 para_env=para_env, mos=mos, input=input)
137
138 !calculation of single occupied states to which unitary transformations should not be applied in LOW SPIN ROKS
139 IF (dft_control%restricted) THEN
140 low_spin_roks_section => section_vals_get_subs_vals(input, "DFT%LOW_SPIN_ROKS")
141 CALL section_vals_val_get(low_spin_roks_section, "SPIN_CONFIGURATION", i_rep_val=1, i_vals=ivec)
142 restricted = SIZE(ivec)
143 ELSE
144 restricted = 0
145 END IF
146
147 NULLIFY (locorb)
148 IF (PRESENT(ext_mo_coeff)) THEN
149 locorb => ext_mo_coeff
150 ELSE
151 CALL get_mo_set(mo_set=mos(myspin), mo_coeff=locorb)
152 END IF
153
154 loc_method = localized_wfn_control%localization_method
155
156 nmosub = localized_wfn_control%nloc_states(myspin)
157 IF (localized_wfn_control%operator_type == op_loc_berry) THEN
158 ! Here we allocate op_fm_set with the RIGHT size for uks
159 NULLIFY (tmp_fm_struct)
160 CALL cp_fm_struct_create(tmp_fm_struct, nrow_global=nmosub, &
161 ncol_global=nmosub, para_env=para_env, &
162 context=locorb%matrix_struct%context)
163 !
164 ALLOCATE (op_fm_set(2, dim_op))
165 DO i = 1, dim_op
166 DO j = 1, SIZE(op_fm_set, 1)
167 CALL cp_fm_create(op_fm_set(j, i), tmp_fm_struct)
168 CALL cp_fm_get_info(op_fm_set(j, i), nrow_global=nmosub)
169 CALL cp_fm_set_all(op_fm_set(j, i), 0.0_dp)
170 END DO
171 END DO
172 CALL cp_fm_struct_release(tmp_fm_struct)
173 END IF
174
175 IF (localized_wfn_control%do_mixed) THEN
176 CALL loc_mixed_method(qs_env, qs_loc_env, print_loc_section, myspin, op_fm_set)
177 ELSE
178 SELECT CASE (localized_wfn_control%operator_type)
179 CASE (op_loc_berry)
180 CALL optimize_loc_berry(loc_method, qs_loc_env, moloc_coeff(myspin), op_sm_set, &
181 op_fm_set, para_env, cell, weights, myspin, print_loc_section, &
182 restricted=restricted)
183 CASE (op_loc_boys)
184 cpabort("Boys localization not implemented")
185 CASE (op_loc_pipek)
186 CALL optimize_loc_pipek(qs_env, loc_method, qs_loc_env, moloc_coeff(myspin), &
187 op_fm_set, myspin, print_loc_section)
188 END SELECT
189 END IF
190
191 ! Here we dealloctate op_fm_set
192 IF (localized_wfn_control%operator_type == op_loc_berry) THEN
193 IF (ASSOCIATED(op_fm_set)) THEN
194 DO i = 1, dim_op
195 DO j = 1, SIZE(op_fm_set, 1)
196 CALL cp_fm_release(op_fm_set(j, i))
197 END DO
198 END DO
199 DEALLOCATE (op_fm_set)
200 END IF
201 END IF
202
203 ! give back the localized orbitals
204 CALL get_mo_set(mo_set=mos(myspin), nao=nao)
205 lb = localized_wfn_control%lu_bound_states(1, myspin)
206 ub = localized_wfn_control%lu_bound_states(2, myspin)
207
208 IF (localized_wfn_control%set_of_states == state_loc_list) THEN
209 ALLOCATE (vecbuffer(1, nao))
210 nmosub = SIZE(localized_wfn_control%loc_states, 1)
211 imoloc = 0
212 DO i = lb, ub
213 ! Get the index in the subset
214 imoloc = imoloc + 1
215 ! Get the index in the full set
216 imo = localized_wfn_control%loc_states(i, myspin)
217
218 CALL cp_fm_get_submatrix(moloc_coeff(myspin), vecbuffer, 1, imoloc, &
219 nao, 1, transpose=.true.)
220 CALL cp_fm_set_submatrix(locorb, vecbuffer, 1, imo, nao, 1, transpose=.true.)
221 END DO
222 DEALLOCATE (vecbuffer)
223 ELSE
224 nmosub = localized_wfn_control%nloc_states(myspin)
225 CALL cp_fm_to_fm(moloc_coeff(myspin), locorb, nmosub, 1, lb)
226 END IF
227
228 ! Write cube files if required
229 IF (localized_wfn_control%print_cubes) THEN
230 CALL loc_print(qs_env, qs_loc_env, moloc_coeff, myspin, print_loc_section)
231 END IF
232 first_time = .false.
233
234 CALL timestop(handle)
235
236 END SUBROUTINE qs_loc_driver
237
238! **************************************************************************************************
239!> \brief set up the calculation of localized orbitals
240!> \param qs_env ...
241!> \param qs_loc_env ...
242!> \param print_loc_section ...
243!> \param myspin ...
244!> \param op_fm_set ...
245!> \par History
246!> 04.2023 refactored [JGH]
247!> \author MI
248! **************************************************************************************************
249 SUBROUTINE loc_mixed_method(qs_env, qs_loc_env, print_loc_section, myspin, op_fm_set)
250
251 TYPE(qs_environment_type), POINTER :: qs_env
252 TYPE(qs_loc_env_type), POINTER :: qs_loc_env
253 TYPE(section_vals_type), POINTER :: print_loc_section
254 INTEGER, INTENT(IN) :: myspin
255 TYPE(cp_fm_type), DIMENSION(:, :), POINTER :: op_fm_set
256
257 CHARACTER(len=*), PARAMETER :: routinen = 'loc_mixed_method'
258
259 INTEGER :: dim_op, handle, jspin, loc_method, nao, &
260 ndummy, nextra, ngextra, nguess, nmo, &
261 nmosub, norextra, restricted
262 INTEGER, DIMENSION(2) :: nelectron_spin
263 INTEGER, DIMENSION(:), POINTER :: ivec
264 LOGICAL :: do_ortho, has_unit_metric, &
265 my_guess_atomic, my_guess_wan
266 REAL(dp), DIMENSION(6) :: weights
267 REAL(kind=dp), DIMENSION(:, :), POINTER :: tmp_mat
268 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
269 TYPE(cell_type), POINTER :: cell
270 TYPE(cp_fm_struct_type), POINTER :: tmp_fm_struct
271 TYPE(cp_fm_type) :: mos_guess, tmp_fm, tmp_fm_1, vectors_2
272 TYPE(cp_fm_type), DIMENSION(:), POINTER :: moloc_coeff
273 TYPE(cp_fm_type), POINTER :: mo_coeff
274 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: p_rmpv
275 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s_kp, op_sm_set
276 TYPE(dbcsr_type), POINTER :: refmatrix, tmatrix
277 TYPE(dft_control_type), POINTER :: dft_control
278 TYPE(localized_wfn_control_type), POINTER :: localized_wfn_control
279 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
280 TYPE(mp_para_env_type), POINTER :: para_env
281 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
282 POINTER :: sab_orb
283 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
284 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
285 TYPE(section_vals_type), POINTER :: input, low_spin_roks_section
286
287 CALL timeset(routinen, handle)
288
289 NULLIFY (moloc_coeff, op_sm_set)
290 CALL get_qs_env(qs_env=qs_env, dft_control=dft_control, para_env=para_env, mos=mos, input=input)
291
292 !calculation of single occupied states to which unitary transformations should not be applied in LOW SPIN ROKS
293 IF (dft_control%restricted) THEN
294 low_spin_roks_section => section_vals_get_subs_vals(input, "DFT%LOW_SPIN_ROKS")
295 CALL section_vals_val_get(low_spin_roks_section, "SPIN_CONFIGURATION", i_rep_val=1, i_vals=ivec)
296 restricted = SIZE(ivec)
297 ELSE
298 restricted = 0
299 END IF
300
301 CALL get_qs_loc_env(qs_loc_env=qs_loc_env, &
302 localized_wfn_control=localized_wfn_control, &
303 moloc_coeff=moloc_coeff, op_sm_set=op_sm_set, cell=cell, &
304 weights=weights, dim_op=dim_op)
305
306 CALL get_mo_set(mo_set=mos(myspin), nao=nao, nmo=nmo)
307 loc_method = localized_wfn_control%localization_method
308 nmosub = localized_wfn_control%nloc_states(myspin)
309
310 cpassert(localized_wfn_control%operator_type == op_loc_berry)
311 cpassert(localized_wfn_control%do_mixed)
312
313 my_guess_atomic = .false.
314 ! SGh-wan: if atomic guess and do_mixed and nextra > 0
315 ! read CPO_GUESS; CASE ATOMIC / RESTART / RANDOM (0/1/2)
316 ! read CPO_GUESS_SPACE if CASE ATOMIC; CASE ALL / WAN
317 nextra = localized_wfn_control%nextra
318 IF (nextra > 0) THEN
319 my_guess_atomic = .true.
320 my_guess_wan = .false.
321 do_ortho = .true.
322 SELECT CASE (localized_wfn_control%coeff_po_guess)
323
324 CASE (do_loc_cpo_atomic)
325 my_guess_atomic = .true.
326 NULLIFY (atomic_kind_set, qs_kind_set, particle_set, matrix_s_kp, sab_orb, p_rmpv, &
327 refmatrix, tmatrix)
328 CALL get_qs_env(qs_env=qs_env, &
329 atomic_kind_set=atomic_kind_set, &
330 qs_kind_set=qs_kind_set, &
331 particle_set=particle_set, &
332 matrix_s_kp=matrix_s_kp, &
333 has_unit_metric=has_unit_metric, &
334 nelectron_spin=nelectron_spin, &
335 sab_orb=sab_orb)
336
337 refmatrix => matrix_s_kp(1, 1)%matrix
338 ! create p_rmpv
339 CALL dbcsr_allocate_matrix_set(p_rmpv, dft_control%nspins)
340 DO jspin = 1, dft_control%nspins
341 ALLOCATE (p_rmpv(jspin)%matrix)
342 tmatrix => p_rmpv(jspin)%matrix
343 CALL dbcsr_create(matrix=tmatrix, template=refmatrix, &
344 matrix_type=dbcsr_type_symmetric, nze=0)
345 CALL cp_dbcsr_alloc_block_from_nbl(tmatrix, sab_orb)
346 CALL dbcsr_set(tmatrix, 0.0_dp)
347 END DO
348 CALL calculate_atomic_block_dm(p_rmpv, refmatrix, atomic_kind_set, qs_kind_set, &
349 dft_control%nspins, nelectron_spin, 0, para_env)
350 CASE (do_loc_cpo_restart)
351 my_guess_atomic = .false.
352 my_guess_wan = .true.
353 CASE (do_loc_cpo_random)
354 my_guess_atomic = .false.
355 END SELECT
356
357 norextra = nmo - nmosub
358 CALL get_mo_set(mo_set=mos(myspin), mo_coeff=mo_coeff)
359 CALL cp_fm_struct_create(tmp_fm_struct, nrow_global=nao, &
360 ncol_global=norextra, para_env=para_env, context=mo_coeff%matrix_struct%context)
361 CALL cp_fm_create(vectors_2, tmp_fm_struct)
362 CALL cp_fm_struct_release(tmp_fm_struct)
363 ALLOCATE (tmp_mat(nao, norextra))
364 CALL cp_fm_get_submatrix(mo_coeff, tmp_mat, 1, nmosub + 1)
365 CALL cp_fm_set_submatrix(vectors_2, tmp_mat)
366 DEALLOCATE (tmp_mat)
367
368 ! if guess "atomic" generate MOs based on atomic densities and
369 ! pass on to optimize_loc_berry
370 IF (my_guess_atomic .OR. my_guess_wan) THEN
371
372 SELECT CASE (localized_wfn_control%coeff_po_guess_mo_space)
373
375 ndummy = nmosub
377 ndummy = nmo
378 do_ortho = .false.
379
380 END SELECT
381
382 CALL cp_fm_struct_create(tmp_fm_struct, nrow_global=nao, &
383 ncol_global=ndummy, para_env=para_env, &
384 context=mo_coeff%matrix_struct%context)
385 CALL cp_fm_create(mos_guess, tmp_fm_struct)
386 CALL cp_fm_set_all(mos_guess, 0.0_dp)
387
388 IF (my_guess_atomic) THEN
389 CALL cp_fm_create(tmp_fm, tmp_fm_struct)
390 CALL cp_fm_create(tmp_fm_1, tmp_fm_struct)
391 CALL cp_fm_set_all(tmp_fm, 0.0_dp)
392 CALL cp_fm_set_all(tmp_fm_1, 0.0_dp)
393 CALL cp_fm_init_random(tmp_fm, ndummy)
394 IF (has_unit_metric) THEN
395 CALL cp_fm_to_fm(tmp_fm, tmp_fm_1)
396 ELSE
397 ! PS*C(:,1:nomo)+C(:,nomo+1:nmo) (nomo=NINT(nelectron/maxocc))
398 CALL cp_dbcsr_sm_fm_multiply(refmatrix, tmp_fm, tmp_fm_1, ndummy)
399 END IF
400 CALL cp_dbcsr_sm_fm_multiply(p_rmpv(myspin)%matrix, tmp_fm_1, mos_guess, ndummy)
401 CALL cp_fm_release(tmp_fm)
402 CALL cp_fm_release(tmp_fm_1)
403 CALL cp_fm_struct_release(tmp_fm_struct)
404 ELSEIF (my_guess_wan) THEN
405 nguess = localized_wfn_control%nguess(myspin)
406 ALLOCATE (tmp_mat(nao, nguess))
407 CALL cp_fm_get_submatrix(moloc_coeff(myspin), tmp_mat, 1, 1, nao, nguess)
408 CALL cp_fm_set_submatrix(mos_guess, tmp_mat, 1, 1, nao, nguess)
409 DEALLOCATE (tmp_mat)
410 ngextra = nmosub - nguess
411 !WRITE(*,*) 'nguess, ngextra = ', nguess, ngextra
412 CALL cp_fm_struct_release(tmp_fm_struct)
413 IF (ngextra > 0) THEN
414 CALL cp_fm_struct_create(tmp_fm_struct, nrow_global=nao, &
415 ncol_global=ngextra, para_env=para_env, &
416 context=mo_coeff%matrix_struct%context)
417 CALL cp_fm_create(tmp_fm, tmp_fm_struct)
418 CALL cp_fm_init_random(tmp_fm, ngextra)
419 ALLOCATE (tmp_mat(nao, ngextra))
420 CALL cp_fm_get_submatrix(tmp_fm, tmp_mat, 1, 1, nao, ngextra)
421 CALL cp_fm_set_submatrix(mos_guess, tmp_mat, 1, nguess + 1, nao, ngextra)
422 DEALLOCATE (tmp_mat)
423 CALL cp_fm_release(tmp_fm)
424 CALL cp_fm_struct_release(tmp_fm_struct)
425 ELSE
426 do_ortho = .false.
427 END IF
428 ALLOCATE (tmp_mat(nao, nmosub))
429 CALL cp_fm_get_submatrix(mo_coeff, tmp_mat, 1, 1, nao, nmosub)
430 CALL cp_fm_set_submatrix(moloc_coeff(myspin), tmp_mat)
431 DEALLOCATE (tmp_mat)
432 END IF
433
434 IF (do_ortho) THEN
435 IF ((my_guess_atomic) .OR. (my_guess_wan)) THEN
436 !! and ortho the result
437 IF (has_unit_metric) THEN
438 CALL make_basis_simple(mos_guess, ndummy)
439 ELSE
440 CALL make_basis_sm(mos_guess, ndummy, refmatrix)
441 END IF
442 END IF
443 END IF
444
445 CALL optimize_loc_berry(loc_method, qs_loc_env, moloc_coeff(myspin), op_sm_set, &
446 op_fm_set, para_env, cell, weights, myspin, print_loc_section, &
447 restricted=restricted, &
448 nextra=nextra, nmo=nmo, vectors_2=vectors_2, guess_mos=mos_guess)
449 CALL cp_fm_release(mos_guess)
450 ELSE
451 CALL optimize_loc_berry(loc_method, qs_loc_env, moloc_coeff(myspin), op_sm_set, &
452 op_fm_set, para_env, cell, weights, myspin, print_loc_section, &
453 restricted=restricted, &
454 nextra=nextra, nmo=nmo, vectors_2=vectors_2)
455 END IF
456 CALL cp_fm_release(vectors_2)
457 IF (my_guess_atomic) CALL dbcsr_deallocate_matrix_set(p_rmpv)
458 ELSE
459 CALL optimize_loc_berry(loc_method, qs_loc_env, moloc_coeff(myspin), op_sm_set, &
460 op_fm_set, para_env, cell, weights, myspin, print_loc_section, &
461 restricted=restricted, nextra=0)
462 END IF
463
464 CALL timestop(handle)
465
466 END SUBROUTINE loc_mixed_method
467
468! **************************************************************************************************
469!> \brief printing of Cube files of localized orbitals
470!> \param qs_env ...
471!> \param qs_loc_env ...
472!> \param moloc_coeff ...
473!> \param ispin ...
474!> \param print_loc_section ...
475! **************************************************************************************************
476 SUBROUTINE loc_print(qs_env, qs_loc_env, moloc_coeff, ispin, print_loc_section)
477
478 TYPE(qs_environment_type), POINTER :: qs_env
479 TYPE(qs_loc_env_type), POINTER :: qs_loc_env
480 TYPE(cp_fm_type), DIMENSION(:), POINTER :: moloc_coeff
481 INTEGER, INTENT(IN), OPTIONAL :: ispin
482 TYPE(section_vals_type), POINTER :: print_loc_section
483
484 CHARACTER(LEN=default_string_length) :: my_pos
485 INTEGER :: i, ir, istate, j, jstate, n_rep, ncubes, &
486 nmo
487 INTEGER, DIMENSION(:), POINTER :: bounds, list, list_cubes
488 LOGICAL :: append_cube, list_cubes_setup
489 REAL(kind=dp), DIMENSION(:, :), POINTER :: centers
490 TYPE(localized_wfn_control_type), POINTER :: localized_wfn_control
491 TYPE(section_vals_type), POINTER :: print_key
492
493 list_cubes_setup = .false.
494 NULLIFY (bounds, list, list_cubes)
495
496 CALL get_qs_loc_env(qs_loc_env=qs_loc_env, &
497 localized_wfn_control=localized_wfn_control)
498
499 ! Provides boundaries of MOs
500 CALL section_vals_val_get(print_loc_section, "WANNIER_CUBES%CUBES_LU_BOUNDS", &
501 i_vals=bounds)
502 ncubes = bounds(2) - bounds(1) + 1
503 IF (ncubes > 0) THEN
504 list_cubes_setup = .true.
505 ALLOCATE (list_cubes(ncubes))
506 DO ir = 1, ncubes
507 list_cubes(ir) = bounds(1) + (ir - 1)
508 END DO
509 END IF
510
511 ! Provides the list of MOs
512 CALL section_vals_val_get(print_loc_section, "WANNIER_CUBES%CUBES_LIST", &
513 n_rep_val=n_rep)
514 IF (.NOT. list_cubes_setup) THEN
515 ncubes = 0
516 DO ir = 1, n_rep
517 CALL section_vals_val_get(print_loc_section, "WANNIER_CUBES%CUBES_LIST", &
518 i_rep_val=ir, i_vals=list)
519 IF (ASSOCIATED(list)) THEN
520 CALL reallocate(list_cubes, 1, ncubes + SIZE(list))
521 DO i = 1, SIZE(list)
522 list_cubes(i + ncubes) = list(i)
523 END DO
524 ncubes = ncubes + SIZE(list)
525 END IF
526 END DO
527 IF (ncubes > 0) list_cubes_setup = .true.
528 END IF
529
530 ! Full list of Mos
531 IF (.NOT. list_cubes_setup) THEN
532 list_cubes_setup = .true.
533 ncubes = localized_wfn_control%nloc_states(1)
534 IF (ncubes > 0) THEN
535 ALLOCATE (list_cubes(ncubes))
536 END IF
537 DO i = 1, ncubes
538 list_cubes(i) = i
539 END DO
540 END IF
541
542 ncubes = SIZE(list_cubes)
543 ncubes = min(ncubes, nmo)
544 ALLOCATE (centers(6, ncubes))
545 DO i = 1, ncubes
546 istate = list_cubes(i)
547 DO j = 1, localized_wfn_control%nloc_states(ispin)
548 jstate = localized_wfn_control%loc_states(j, ispin)
549 IF (istate == jstate) THEN
550 centers(1:6, i) = localized_wfn_control%centers_set(ispin)%array(1:6, j)
551 EXIT
552 END IF
553 END DO
554 END DO ! ncubes
555
556 ! Real call for dumping the cube files
557 print_key => section_vals_get_subs_vals(print_loc_section, "WANNIER_CUBES")
558 append_cube = section_get_lval(print_loc_section, "WANNIER_CUBES%APPEND")
559 my_pos = "REWIND"
560 IF (append_cube) THEN
561 my_pos = "APPEND"
562 END IF
563
564 CALL qs_print_cubes(qs_env, moloc_coeff(ispin), ncubes, list_cubes, centers, &
565 print_key, "loc"//trim(adjustl(qs_loc_env%tag_mo)), &
566 ispin=ispin, file_position=my_pos)
567
568 DEALLOCATE (centers)
569 DEALLOCATE (list_cubes)
570
571 END SUBROUTINE loc_print
572
573END MODULE qs_loc_main
Define the atomic kind types and their sub types.
Handles all functions related to the CELL.
Definition cell_types.F:15
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_sm_fm_multiply(matrix, fm_in, fm_out, ncol, alpha, beta)
multiply a dbcsr with a fm matrix
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_set_submatrix(fm, new_values, start_row, start_col, n_rows, n_cols, alpha, beta, transpose)
sets a submatrix of a full matrix fm(start_row:start_row+n_rows,start_col:start_col+n_cols) = alpha*o...
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
subroutine, public cp_fm_get_submatrix(fm, target_m, start_row, start_col, n_rows, n_cols, transpose)
gets a submatrix of a full matrix op(target_m)(1:n_rows,1:n_cols) =fm(start_row:start_row+n_rows,...
subroutine, public cp_fm_init_random(matrix, ncol, start_col)
fills a matrix with random numbers
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public op_loc_pipek
integer, parameter, public do_loc_cpo_space_nmo
integer, parameter, public op_loc_berry
integer, parameter, public op_loc_boys
integer, parameter, public do_loc_cpo_random
integer, parameter, public state_loc_list
integer, parameter, public do_loc_cpo_restart
integer, parameter, public do_loc_cpo_space_wan
integer, parameter, public do_loc_cpo_atomic
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
logical function, public section_get_lval(section_vals, keyword_name)
...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition list.F:24
Utility routines for the memory handling.
Interface to the message passing library MPI.
Define the data structure for the particle information.
Routine to return block diagonal density matrix. Blocks correspond to the atomic densities.
subroutine, public calculate_atomic_block_dm(pmatrix, matrix_s, atomic_kind_set, qs_kind_set, nspin, nelectron_spin, ounit, para_env)
returns a block diagonal density matrix. Blocks correspond to the atomic densities.
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
Driver for the localization that should be general for all the methods available and all the definiti...
Definition qs_loc_main.F:22
subroutine, public qs_loc_driver(qs_env, qs_loc_env, print_loc_section, myspin, ext_mo_coeff)
set up the calculation of localized orbitals
Definition qs_loc_main.F:96
Driver for the localization that should be general for all the methods available and all the definiti...
subroutine, public optimize_loc_pipek(qs_env, method, qs_loc_env, vectors, zij_fm_set, ispin, print_loc_section)
...
subroutine, public optimize_loc_berry(method, qs_loc_env, vectors, op_sm_set, zij_fm_set, para_env, cell, weights, ispin, print_loc_section, restricted, nextra, nmo, vectors_2, guess_mos)
Calculate and optimize the spread functional as calculated from the selected mos and by the definitio...
subroutine, public qs_print_cubes(qs_env, mo_coeff, nstates, state_list, centers, print_key, root, ispin, idir, state0, file_position)
write the cube files for a set of selected states
New version of the module for the localization of the molecular orbitals This should be able to use d...
subroutine, public get_qs_loc_env(qs_loc_env, cell, local_molecules, localized_wfn_control, moloc_coeff, op_sm_set, op_fm_set, para_env, particle_set, weights, dim_op)
...
collects routines that perform operations directly related to MOs
subroutine, public make_basis_simple(vmatrix, ncol)
given a set of vectors, return an orthogonal (C^T C == 1) set spanning the same space (notice,...
subroutine, public make_basis_sm(vmatrix, ncol, matrix_s)
returns an S-orthonormal basis v (v^T S v ==1)
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
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.
Define the neighbor list data types and the corresponding functionality.
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
keeps the information about the structure of a full matrix
represent a full matrix
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
A type that holds controlling information for the calculation of the spread of wfn and the optimizati...
contains all the info needed by quickstep to calculate the spread of a selected set of orbitals and i...