Loading [MathJax]/extensions/tex2jax.js
 (git:ed6f26b)
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros
qs_core_hamiltonian.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Calculation of the core Hamiltonian integral matrix <a|H|b> over
10!> Cartesian Gaussian-type functions.
11!>
12!> <a|H|b> = <a|T|b> + <a|V|b>
13!>
14!> Kinetic energy:
15!>
16!> <a|T|b> = <a|-nabla**2/2|b>
17!> \_______________/
18!> |
19!> kinetic
20!>
21!> Nuclear potential energy:
22!>
23!> a) Allelectron calculation:
24!>
25!> erfc(r)
26!> <a|V|b> = -Z*<a|---------|b>
27!> r
28!>
29!> 1 - erf(r)
30!> = -Z*<a|------------|b>
31!> r
32!>
33!> 1 erf(r)
34!> = -Z*(<a|---|b> - <a|--------|b>)
35!> r r
36!>
37!> 1
38!> = -Z*(<a|---|b> - N*<ab||c>)
39!> r
40!>
41!> -Z
42!> = <a|---|b> + Z*N*<ab||c>
43!> r
44!> \_______/ \_____/
45!> | |
46!> nuclear coulomb
47!>
48!> b) Pseudopotential calculation (Goedecker, Teter and Hutter; GTH):
49!>
50!> <a|V|b> = <a|(V(local) + V(non-local))|b>
51!>
52!> = <a|(V(local)|b> + <a|V(non-local))|b>
53!>
54!> <a|V(local)|b> = <a|-Z(eff)*erf(SQRT(2)*alpha*r)/r +
55!> (C1 + C2*(alpha*r)**2 + C3*(alpha*r)**4 +
56!> C4*(alpha*r)**6)*exp(-(alpha*r)**2/2))|b>
57!>
58!> <a|V(non-local)|b> = <a|p(l,i)>*h(i,j)*<p(l,j)|b>
59!> \par Literature
60!> S. Goedecker, M. Teter and J. Hutter, Phys. Rev. B 54, 1703 (1996)
61!> C. Hartwigsen, S. Goedecker and J. Hutter, Phys. Rev. B 58, 3641 (1998)
62!> M. Krack and M. Parrinello, Phys. Chem. Chem. Phys. 2, 2105 (2000)
63!> S. Obara and A. Saika, J. Chem. Phys. 84, 3963 (1986)
64!> \par History
65!> - Joost VandeVondele (April 2003) : added LSD forces
66!> - Non-redundant calculation of the non-local part of the GTH PP
67!> (22.05.2003,MK)
68!> - New parallelization scheme (27.06.2003,MK)
69!> - OpenMP version (07.12.2003,JGH)
70!> - Binary search loop for VPPNL operators (09.01.2004,JGH,MK)
71!> - Refactoring of pseudopotential and nuclear attraction integrals (25.02.2009,JGH)
72!> - General refactoring (01.10.2010,JGH)
73!> - Refactoring related to the new kinetic energy and overlap routines (07.2014,JGH)
74!> - k-point functionality (07.2015,JGH)
75!> \author Matthias Krack (14.09.2000,21.03.02)
76! **************************************************************************************************
80 USE core_ae, ONLY: build_core_ae
81 USE core_ppl, ONLY: build_core_ppl
82 USE core_ppnl, ONLY: build_core_ppnl
85 USE cp_dbcsr_api, ONLY: &
88 dbcsr_p_type, dbcsr_set, dbcsr_type, dbcsr_type_antisymmetric
96 USE cp_output_handling, ONLY: cp_p_file,&
103 rel_none,&
108 USE kinds, ONLY: default_string_length,&
109 dp
110 USE kpoint_types, ONLY: get_kpoint_info,&
115 USE qs_condnum, ONLY: overlap_condnum
120 USE qs_kind_types, ONLY: get_qs_kind,&
123 USE qs_ks_types, ONLY: get_ks_env,&
128 USE qs_oce_types, ONLY: allocate_oce_set,&
132 USE qs_rho_types, ONLY: qs_rho_get,&
134 USE virial_types, ONLY: virial_type
135#include "./base/base_uses.f90"
136
137 IMPLICIT NONE
138
139 PRIVATE
140
141 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_core_hamiltonian'
142
145
146CONTAINS
147
148! **************************************************************************************************
149!> \brief Cosntruction of the QS Core Hamiltonian Matrix
150!> \param qs_env ...
151!> \param calculate_forces ...
152!> \author Creation (11.03.2002,MK)
153!> Non-redundant calculation of the non-local part of the GTH PP (22.05.2003,MK)
154!> New parallelization scheme (27.06.2003,MK)
155! **************************************************************************************************
156 SUBROUTINE build_core_hamiltonian_matrix(qs_env, calculate_forces)
157
158 TYPE(qs_environment_type), POINTER :: qs_env
159 LOGICAL, INTENT(IN) :: calculate_forces
160
161 CHARACTER(LEN=*), PARAMETER :: routinen = 'build_core_hamiltonian_matrix'
162
163 INTEGER :: handle, ic, img, iw, nder, nders, &
164 nimages, nkind
165 INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
166 LOGICAL :: h_is_complex, norml1, norml2, ofdft, &
167 use_arnoldi, use_virial
168 REAL(kind=dp) :: eps_filter, eps_fit
169 REAL(kind=dp), DIMENSION(2) :: condnum
170 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
171 TYPE(cp_blacs_env_type), POINTER :: blacs_env
172 TYPE(cp_logger_type), POINTER :: logger
173 TYPE(dbcsr_distribution_type), POINTER :: dbcsr_dist
174 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_h, matrix_p, matrix_s, matrix_t, &
175 matrix_w
176 TYPE(dft_control_type), POINTER :: dft_control
177 TYPE(kg_environment_type), POINTER :: kg_env
178 TYPE(kpoint_type), POINTER :: kpoints
179 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
180 POINTER :: sab_orb, sap_oce
181 TYPE(oce_matrix_type), POINTER :: oce
182 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
183 TYPE(qs_force_type), DIMENSION(:), POINTER :: force
184 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
185 TYPE(qs_ks_env_type), POINTER :: ks_env
186 TYPE(qs_rho_type), POINTER :: rho
187 TYPE(virial_type), POINTER :: virial
188
189 IF (calculate_forces) THEN
190 CALL timeset(routinen//"_forces", handle)
191 ELSE
192 CALL timeset(routinen, handle)
193 END IF
194
195 NULLIFY (logger)
196 logger => cp_get_default_logger()
197
198 NULLIFY (dft_control)
199 CALL get_qs_env(qs_env=qs_env, dft_control=dft_control)
200
201 ! is this a orbital-free method calculation
202 ofdft = dft_control%qs_control%ofgpw
203
204 nimages = dft_control%nimages
205 IF (ofdft) THEN
206 cpassert(nimages == 1)
207 END IF
208
209 nders = 0
210 IF (calculate_forces) THEN
211 nder = 1
212 ELSE
213 IF (cp_print_key_should_output(logger%iter_info, qs_env%input, &
214 "DFT%PRINT%AO_MATRICES/DERIVATIVES") /= 0) THEN
215 nder = 1
216 ELSE
217 nder = 0
218 END IF
219 END IF
220
221 IF ((cp_print_key_should_output(logger%iter_info, qs_env%input, &
222 "DFT%PRINT%AO_MATRICES/OVERLAP") /= 0 .AND. &
223 btest(cp_print_key_should_output(logger%iter_info, qs_env%input, &
224 "DFT%PRINT%AO_MATRICES/DERIVATIVES"), cp_p_file))) THEN
225 nders = 1
226 END IF
227
228 ! the delta pulse in the periodic case needs the momentum operator,
229 ! which is equivalent to the derivative of the overlap matrix
230 IF (ASSOCIATED(dft_control%rtp_control)) THEN
231 IF (dft_control%rtp_control%apply_delta_pulse .AND. &
232 dft_control%rtp_control%periodic) THEN
233 nders = 1
234 END IF
235 END IF
236
237 IF (dft_control%tddfpt2_control%enabled) THEN
238 nders = 1
239 IF (dft_control%do_admm) THEN
240 IF (dft_control%admm_control%purification_method /= do_admm_purify_none) &
241 CALL cp_abort(__location__, &
242 "Only purification method NONE is possible with TDDFT at the moment")
243 END IF
244 END IF
245
246 ! filter for new matrices
247 eps_filter = dft_control%qs_control%eps_filter_matrix
248 !
249 NULLIFY (ks_env)
250 CALL get_qs_env(qs_env=qs_env, ks_env=ks_env)
251 NULLIFY (matrix_s, matrix_t)
252 CALL get_qs_env(qs_env=qs_env, kinetic_kp=matrix_t, matrix_s_kp=matrix_s)
253 NULLIFY (sab_orb)
254 CALL get_qs_env(qs_env=qs_env, sab_orb=sab_orb)
255 NULLIFY (rho, force, matrix_p, matrix_w)
256 IF (calculate_forces) THEN
257 CALL get_qs_env(qs_env=qs_env, force=force, matrix_w_kp=matrix_w)
258 CALL get_qs_env(qs_env=qs_env, rho=rho)
259 CALL qs_rho_get(rho, rho_ao_kp=matrix_p)
260 ! *** If LSD, then combine alpha density and beta density to
261 ! *** total density: alpha <- alpha + beta and
262 ! *** spin density: beta <- alpha - beta
263 ! (since all things can be computed based on the sum of these matrices anyway)
264 ! (matrix_p is restored at the end of the run, matrix_w is left in its modified state
265 ! (as it should not be needed afterwards)
266 IF (SIZE(matrix_p, 1) == 2) THEN
267 DO img = 1, nimages
268 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
269 alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
270 CALL dbcsr_add(matrix_p(2, img)%matrix, matrix_p(1, img)%matrix, &
271 alpha_scalar=-2.0_dp, beta_scalar=1.0_dp)
272 CALL dbcsr_add(matrix_w(1, img)%matrix, matrix_w(2, img)%matrix, &
273 alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
274 END DO
275 END IF
276 ! S matrix
277 CALL build_overlap_matrix(ks_env, nderivative=nders, matrixkp_s=matrix_s, &
278 matrix_name="OVERLAP MATRIX", &
279 basis_type_a="ORB", &
280 basis_type_b="ORB", &
281 sab_nl=sab_orb, calculate_forces=.true., &
282 matrixkp_p=matrix_w)
283 ! T matrix
284 IF (.NOT. ofdft) &
285 CALL build_kinetic_matrix(ks_env, matrixkp_t=matrix_t, &
286 matrix_name="KINETIC ENERGY MATRIX", &
287 basis_type="ORB", &
288 sab_nl=sab_orb, calculate_forces=.true., &
289 matrixkp_p=matrix_p, &
290 eps_filter=eps_filter)
291 ! *** If LSD, then recover alpha density and beta density ***
292 ! *** from the total density (1) and the spin density (2) ***
293 ! *** The W matrix is neglected, since it will be destroyed ***
294 ! *** in the calling force routine after leaving this routine ***
295 IF (SIZE(matrix_p, 1) == 2) THEN
296 DO img = 1, nimages
297 CALL dbcsr_add(matrix_p(1, img)%matrix, matrix_p(2, img)%matrix, &
298 alpha_scalar=0.5_dp, beta_scalar=0.5_dp)
299 CALL dbcsr_add(matrix_p(2, img)%matrix, matrix_p(1, img)%matrix, &
300 alpha_scalar=-1.0_dp, beta_scalar=1.0_dp)
301 END DO
302 END IF
303 ELSE
304 ! S matrix
305 CALL build_overlap_matrix(ks_env, nderivative=nders, matrixkp_s=matrix_s, &
306 matrix_name="OVERLAP MATRIX", &
307 basis_type_a="ORB", &
308 basis_type_b="ORB", &
309 sab_nl=sab_orb)
310 ! T matrix
311 IF (.NOT. ofdft) &
312 CALL build_kinetic_matrix(ks_env, matrixkp_t=matrix_t, &
313 matrix_name="KINETIC ENERGY MATRIX", &
314 basis_type="ORB", &
315 sab_nl=sab_orb, &
316 eps_filter=eps_filter)
317 END IF
318
319 ! (Re-)allocate H matrix based on overlap matrix
320 CALL get_ks_env(ks_env, complex_ks=h_is_complex)
321 CALL qs_matrix_h_allocate(qs_env, matrix_s(1, 1)%matrix, is_complex=h_is_complex)
322
323 NULLIFY (matrix_h)
324 CALL get_qs_env(qs_env, matrix_h_kp=matrix_h)
325
326 IF (.NOT. ofdft) THEN
327 DO img = 1, nimages
328 CALL dbcsr_copy(matrix_h(1, img)%matrix, matrix_t(1, img)%matrix, &
329 keep_sparsity=.true., name="CORE HAMILTONIAN MATRIX")
330 END DO
331 END IF
332
333 NULLIFY (qs_kind_set, atomic_kind_set, particle_set)
334 CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set, atomic_kind_set=atomic_kind_set, &
335 particle_set=particle_set)
336
337 IF (.NOT. ofdft) THEN
338 ! relativistic atomic correction to kinetic energy
339 IF (qs_env%rel_control%rel_method /= rel_none) THEN
340 IF (qs_env%rel_control%rel_transformation == rel_trans_atom) THEN
341 IF (nimages == 1) THEN
342 ic = 1
343 ELSE
344 CALL get_ks_env(ks_env=ks_env, kpoints=kpoints)
345 CALL get_kpoint_info(kpoint=kpoints, cell_to_index=cell_to_index)
346 ic = cell_to_index(0, 0, 0)
347 END IF
348 CALL build_atomic_relmat(matrix_h(1, 1)%matrix, atomic_kind_set, qs_kind_set)
349 ELSE
350 cpabort("Relativistic corrections of this type are currently not implemented")
351 END IF
352 END IF
353 END IF
354
355 ! *** core and pseudopotentials
356 CALL core_matrices(qs_env, matrix_h, matrix_p, calculate_forces, nder)
357
358 ! *** GAPW one-center-expansion (oce) matrices
359 NULLIFY (sap_oce)
360 CALL get_qs_env(qs_env=qs_env, sap_oce=sap_oce)
361 NULLIFY (oce)
362 IF (dft_control%qs_control%gapw .OR. dft_control%qs_control%gapw_xc) THEN
363 CALL get_qs_env(qs_env=qs_env, oce=oce)
364 CALL create_oce_set(oce)
365 nkind = SIZE(atomic_kind_set)
366 CALL allocate_oce_set(oce, nkind)
367 eps_fit = dft_control%qs_control%gapw_control%eps_fit
368 IF (ASSOCIATED(sap_oce)) &
369 CALL build_oce_matrices(oce%intac, calculate_forces, nder, qs_kind_set, particle_set, &
370 sap_oce, eps_fit)
371 END IF
372
373 ! *** KG atomic potentials for nonadditive kinetic energy
374 IF (dft_control%qs_control%do_kg) THEN
375 IF (qs_env%kg_env%tnadd_method == kg_tnadd_atomic) THEN
376 CALL get_qs_env(qs_env=qs_env, kg_env=kg_env, virial=virial, dbcsr_dist=dbcsr_dist)
377 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
378 CALL build_tnadd_mat(kg_env, matrix_p, force, virial, calculate_forces, use_virial, &
379 qs_kind_set, atomic_kind_set, particle_set, sab_orb, dbcsr_dist)
380 END IF
381 END IF
382
383 ! *** Put the core Hamiltonian matrix in the QS environment ***
384 CALL set_qs_env(qs_env, oce=oce)
385 CALL set_ks_env(ks_env, matrix_s_kp=matrix_s, kinetic_kp=matrix_t, matrix_h_kp=matrix_h)
386
387 ! *** Print matrices if requested
388 CALL dump_info_core_hamiltonian(qs_env, calculate_forces)
389
390 ! *** Overlap condition number
391 IF (.NOT. calculate_forces) THEN
392 IF (cp_print_key_should_output(logger%iter_info, qs_env%input, &
393 "DFT%PRINT%OVERLAP_CONDITION") .NE. 0) THEN
394 iw = cp_print_key_unit_nr(logger, qs_env%input, "DFT%PRINT%OVERLAP_CONDITION", &
395 extension=".Log")
396 CALL section_vals_val_get(qs_env%input, "DFT%PRINT%OVERLAP_CONDITION%1-NORM", l_val=norml1)
397 CALL section_vals_val_get(qs_env%input, "DFT%PRINT%OVERLAP_CONDITION%DIAGONALIZATION", l_val=norml2)
398 CALL section_vals_val_get(qs_env%input, "DFT%PRINT%OVERLAP_CONDITION%ARNOLDI", l_val=use_arnoldi)
399 CALL get_qs_env(qs_env=qs_env, blacs_env=blacs_env)
400 CALL overlap_condnum(matrix_s, condnum, iw, norml1, norml2, use_arnoldi, blacs_env)
401 END IF
402 END IF
403
404 CALL timestop(handle)
405
406 END SUBROUTINE build_core_hamiltonian_matrix
407
408! **************************************************************************************************
409!> \brief ...
410!> \param qs_env ...
411!> \param matrix_h ...
412!> \param matrix_p ...
413!> \param calculate_forces ...
414!> \param nder ...
415!> \param atcore ...
416! **************************************************************************************************
417 SUBROUTINE core_matrices(qs_env, matrix_h, matrix_p, calculate_forces, nder, atcore)
418
419 TYPE(qs_environment_type), POINTER :: qs_env
420 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_h, matrix_p
421 LOGICAL, INTENT(IN) :: calculate_forces
422 INTEGER, INTENT(IN) :: nder
423 REAL(kind=dp), DIMENSION(:), OPTIONAL :: atcore
424
425 INTEGER :: natom, nimages
426 INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
427 LOGICAL :: all_present, my_gt_nl, ppl_present, &
428 ppnl_present, use_virial
429 REAL(kind=dp) :: eps_ppnl
430 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
431 TYPE(dft_control_type), POINTER :: dft_control
432 TYPE(kpoint_type), POINTER :: kpoints
433 TYPE(lri_environment_type), POINTER :: lri_env
434 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
435 POINTER :: sab_orb, sac_ae, sac_ppl, sap_ppnl
436 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
437 TYPE(qs_force_type), DIMENSION(:), POINTER :: force
438 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
439 TYPE(qs_ks_env_type), POINTER :: ks_env
440 TYPE(virial_type), POINTER :: virial
441
442 NULLIFY (dft_control)
443 CALL get_qs_env(qs_env=qs_env, ks_env=ks_env, dft_control=dft_control, natom=natom)
444 nimages = dft_control%nimages
445 IF (PRESENT(atcore)) THEN
446 cpassert(SIZE(atcore) >= natom)
447 END IF
448
449 ! check whether a gauge transformed version of the non-local potential part has to be used
450 my_gt_nl = .false.
451 IF (qs_env%run_rtp) THEN
452 cpassert(ASSOCIATED(dft_control%rtp_control))
453 IF (dft_control%rtp_control%velocity_gauge) THEN
454 my_gt_nl = dft_control%rtp_control%nl_gauge_transform
455 END IF
456 END IF
457
458 ! prepare for k-points
459 NULLIFY (cell_to_index)
460 IF (nimages > 1) THEN
461 CALL get_ks_env(ks_env=ks_env, kpoints=kpoints)
462 CALL get_kpoint_info(kpoint=kpoints, cell_to_index=cell_to_index)
463 dft_control%qs_control%do_ppl_method = do_ppl_analytic
464 END IF
465 ! force analytic ppl calculation for GAPW methods
466 IF (dft_control%qs_control%gapw .OR. dft_control%qs_control%gapw_xc) THEN
467 dft_control%qs_control%do_ppl_method = do_ppl_analytic
468 END IF
469
470 ! force
471 NULLIFY (force)
472 IF (calculate_forces) CALL get_qs_env(qs_env=qs_env, force=force)
473 ! virial
474 CALL get_qs_env(qs_env=qs_env, virial=virial)
475 use_virial = virial%pv_availability .AND. (.NOT. virial%pv_numer)
476
477 NULLIFY (qs_kind_set, atomic_kind_set, particle_set)
478 CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set, atomic_kind_set=atomic_kind_set, &
479 particle_set=particle_set)
480
481 NULLIFY (sab_orb, sac_ae, sac_ppl, sap_ppnl)
482 CALL get_qs_env(qs_env=qs_env, &
483 sab_orb=sab_orb, &
484 sac_ae=sac_ae, &
485 sac_ppl=sac_ppl, &
486 sap_ppnl=sap_ppnl)
487
488 ! *** compute the nuclear attraction contribution to the core hamiltonian ***
489 all_present = ASSOCIATED(sac_ae)
490 IF (all_present) THEN
491 CALL build_core_ae(matrix_h, matrix_p, force, virial, calculate_forces, use_virial, nder, &
492 qs_kind_set, atomic_kind_set, particle_set, sab_orb, sac_ae, &
493 nimages, cell_to_index, atcore=atcore)
494 END IF
495
496 ! *** compute the ppl contribution to the core hamiltonian ***
497 ppl_present = ASSOCIATED(sac_ppl)
498 IF (ppl_present) THEN
499 IF (dft_control%qs_control%do_ppl_method == do_ppl_analytic) THEN
500 IF (dft_control%qs_control%lrigpw) THEN
501 CALL get_qs_env(qs_env, lri_env=lri_env)
502 IF (lri_env%ppl_ri) THEN
503 IF (lri_env%exact_1c_terms) THEN
504 cpabort("not implemented")
505 END IF
506 ELSE
507 CALL build_core_ppl(matrix_h, matrix_p, force, virial, calculate_forces, use_virial, nder, &
508 qs_kind_set, atomic_kind_set, particle_set, sab_orb, sac_ppl, &
509 nimages, cell_to_index, "ORB", atcore=atcore)
510 END IF
511 ELSE
512 CALL build_core_ppl(matrix_h, matrix_p, force, virial, calculate_forces, use_virial, nder, &
513 qs_kind_set, atomic_kind_set, particle_set, sab_orb, sac_ppl, &
514 nimages, cell_to_index, "ORB", atcore=atcore)
515 END IF
516 END IF
517 END IF
518
519 ! *** compute the ppnl contribution to the core hamiltonian ***
520 eps_ppnl = dft_control%qs_control%eps_ppnl
521 ppnl_present = ASSOCIATED(sap_ppnl)
522 IF (ppnl_present) THEN
523 IF (.NOT. my_gt_nl) THEN
524 CALL build_core_ppnl(matrix_h, matrix_p, force, virial, calculate_forces, use_virial, nder, &
525 qs_kind_set, atomic_kind_set, particle_set, sab_orb, sap_ppnl, eps_ppnl, &
526 nimages, cell_to_index, "ORB", atcore=atcore)
527 END IF
528 END IF
529
530 END SUBROUTINE core_matrices
531
532! **************************************************************************************************
533!> \brief Adds atomic blocks of relativistic correction for the kinetic energy
534!> \param matrix_h ...
535!> \param atomic_kind_set ...
536!> \param qs_kind_set ...
537! **************************************************************************************************
538 SUBROUTINE build_atomic_relmat(matrix_h, atomic_kind_set, qs_kind_set)
539 TYPE(dbcsr_type), POINTER :: matrix_h
540 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
541 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
542
543 INTEGER :: iatom, ikind, jatom
544 INTEGER, ALLOCATABLE, DIMENSION(:) :: kind_of
545 REAL(kind=dp), DIMENSION(:, :), POINTER :: hblock, reltmat
546 TYPE(dbcsr_iterator_type) :: iter
547
548 CALL get_atomic_kind_set(atomic_kind_set, kind_of=kind_of)
549
550 CALL dbcsr_iterator_start(iter, matrix_h)
551 DO WHILE (dbcsr_iterator_blocks_left(iter))
552 CALL dbcsr_iterator_next_block(iter, iatom, jatom, hblock)
553 IF (iatom == jatom) THEN
554 ikind = kind_of(iatom)
555 CALL get_qs_kind(qs_kind_set(ikind), reltmat=reltmat)
556 IF (ASSOCIATED(reltmat)) hblock = hblock + reltmat
557 END IF
558 END DO
559 CALL dbcsr_iterator_stop(iter)
560
561 END SUBROUTINE build_atomic_relmat
562
563! **************************************************************************************************
564!> \brief Possibly prints matrices after the construction of the Core
565!> Hamiltonian Matrix
566!> \param qs_env ...
567!> \param calculate_forces ...
568! **************************************************************************************************
569 SUBROUTINE dump_info_core_hamiltonian(qs_env, calculate_forces)
570 TYPE(qs_environment_type), POINTER :: qs_env
571 LOGICAL, INTENT(IN) :: calculate_forces
572
573 CHARACTER(LEN=*), PARAMETER :: routinen = 'dump_info_core_hamiltonian'
574
575 INTEGER :: after, handle, i, ic, iw, output_unit
576 LOGICAL :: omit_headers
577 TYPE(cp_logger_type), POINTER :: logger
578 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_v
579 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrixkp_h, matrixkp_s, matrixkp_t
580 TYPE(mp_para_env_type), POINTER :: para_env
581
582 CALL timeset(routinen, handle)
583
584 NULLIFY (logger, matrix_v, para_env)
585 logger => cp_get_default_logger()
586 CALL get_qs_env(qs_env, para_env=para_env)
587
588 ! Print the distribution of the overlap matrix blocks
589 ! this duplicates causes duplicate printing at the force calc
590 IF (.NOT. calculate_forces) THEN
591 IF (btest(cp_print_key_should_output(logger%iter_info, &
592 qs_env%input, "PRINT%DISTRIBUTION"), cp_p_file)) THEN
593 output_unit = cp_print_key_unit_nr(logger, qs_env%input, "PRINT%DISTRIBUTION", &
594 extension=".distribution")
595 CALL get_qs_env(qs_env, matrix_s_kp=matrixkp_s)
596 CALL cp_dbcsr_write_matrix_dist(matrixkp_s(1, 1)%matrix, output_unit, para_env)
597 CALL cp_print_key_finished_output(output_unit, logger, qs_env%input, "PRINT%DISTRIBUTION")
598 END IF
599 END IF
600
601 CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%OMIT_HEADERS", l_val=omit_headers)
602 ! Print the overlap integral matrix, if requested
603 IF (btest(cp_print_key_should_output(logger%iter_info, &
604 qs_env%input, "DFT%PRINT%AO_MATRICES/OVERLAP"), cp_p_file)) THEN
605 iw = cp_print_key_unit_nr(logger, qs_env%input, "DFT%PRINT%AO_MATRICES/OVERLAP", &
606 extension=".Log")
607 CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%NDIGITS", i_val=after)
608 after = min(max(after, 1), 16)
609 CALL get_qs_env(qs_env, matrix_s_kp=matrixkp_s)
610 IF (ASSOCIATED(matrixkp_s)) THEN
611 DO ic = 1, SIZE(matrixkp_s, 2)
612 CALL cp_dbcsr_write_sparse_matrix(matrixkp_s(1, ic)%matrix, 4, after, qs_env, para_env, &
613 output_unit=iw, omit_headers=omit_headers)
614 END DO
615 IF (btest(cp_print_key_should_output(logger%iter_info, qs_env%input, &
616 "DFT%PRINT%AO_MATRICES/DERIVATIVES"), cp_p_file)) THEN
617 DO ic = 1, SIZE(matrixkp_s, 2)
618 DO i = 2, SIZE(matrixkp_s, 1)
619 CALL cp_dbcsr_write_sparse_matrix(matrixkp_s(i, ic)%matrix, 4, after, qs_env, para_env, &
620 output_unit=iw, omit_headers=omit_headers)
621 END DO
622 END DO
623 END IF
624 END IF
625 CALL cp_print_key_finished_output(iw, logger, qs_env%input, &
626 "DFT%PRINT%AO_MATRICES/OVERLAP")
627 END IF
628
629 ! Print the kinetic energy integral matrix, if requested
630 IF (btest(cp_print_key_should_output(logger%iter_info, &
631 qs_env%input, "DFT%PRINT%AO_MATRICES/KINETIC_ENERGY"), cp_p_file)) THEN
632 iw = cp_print_key_unit_nr(logger, qs_env%input, "DFT%PRINT%AO_MATRICES/KINETIC_ENERGY", &
633 extension=".Log")
634 CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%NDIGITS", i_val=after)
635 after = min(max(after, 1), 16)
636 CALL get_qs_env(qs_env, kinetic_kp=matrixkp_t)
637 IF (ASSOCIATED(matrixkp_t)) THEN
638 DO ic = 1, SIZE(matrixkp_t, 2)
639 CALL cp_dbcsr_write_sparse_matrix(matrixkp_t(1, ic)%matrix, 4, after, qs_env, para_env, &
640 output_unit=iw, omit_headers=omit_headers)
641 END DO
642 END IF
643 CALL cp_print_key_finished_output(iw, logger, qs_env%input, &
644 "DFT%PRINT%AO_MATRICES/KINETIC_ENERGY")
645 END IF
646
647 ! Print the potential energy matrix, if requested
648 IF (btest(cp_print_key_should_output(logger%iter_info, &
649 qs_env%input, "DFT%PRINT%AO_MATRICES/POTENTIAL_ENERGY"), cp_p_file)) THEN
650 iw = cp_print_key_unit_nr(logger, qs_env%input, "DFT%PRINT%AO_MATRICES/POTENTIAL_ENERGY", &
651 extension=".Log")
652 CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%NDIGITS", i_val=after)
653 after = min(max(after, 1), 16)
654 CALL get_qs_env(qs_env, matrix_h_kp=matrixkp_h, kinetic_kp=matrixkp_t)
655 IF (ASSOCIATED(matrixkp_h)) THEN
656 IF (SIZE(matrixkp_h, 2) == 1) THEN
657 CALL dbcsr_allocate_matrix_set(matrix_v, 1)
658 ALLOCATE (matrix_v(1)%matrix)
659 CALL dbcsr_copy(matrix_v(1)%matrix, matrixkp_h(1, 1)%matrix, name="POTENTIAL ENERGY MATRIX")
660 CALL dbcsr_add(matrix_v(1)%matrix, matrixkp_t(1, 1)%matrix, &
661 alpha_scalar=1.0_dp, beta_scalar=-1.0_dp)
662 CALL cp_dbcsr_write_sparse_matrix(matrix_v(1)%matrix, 4, after, qs_env, &
663 para_env, output_unit=iw, omit_headers=omit_headers)
664 CALL dbcsr_deallocate_matrix_set(matrix_v)
665 ELSE
666 cpwarn("Printing of potential energy matrix not implemented for k-points")
667 END IF
668 END IF
669 CALL cp_print_key_finished_output(iw, logger, qs_env%input, &
670 "DFT%PRINT%AO_MATRICES/POTENTIAL_ENERGY")
671 END IF
672
673 ! Print the core Hamiltonian matrix, if requested
674 IF (btest(cp_print_key_should_output(logger%iter_info, &
675 qs_env%input, "DFT%PRINT%AO_MATRICES/CORE_HAMILTONIAN"), cp_p_file)) THEN
676 iw = cp_print_key_unit_nr(logger, qs_env%input, "DFT%PRINT%AO_MATRICES/CORE_HAMILTONIAN", &
677 extension=".Log")
678 CALL section_vals_val_get(qs_env%input, "DFT%PRINT%AO_MATRICES%NDIGITS", i_val=after)
679 after = min(max(after, 1), 16)
680 CALL get_qs_env(qs_env, matrix_h_kp=matrixkp_h)
681 IF (ASSOCIATED(matrixkp_h)) THEN
682 DO ic = 1, SIZE(matrixkp_h, 2)
683 CALL cp_dbcsr_write_sparse_matrix(matrixkp_h(1, ic)%matrix, 4, after, qs_env, para_env, &
684 output_unit=iw, omit_headers=omit_headers)
685 END DO
686 END IF
687 CALL cp_print_key_finished_output(iw, logger, qs_env%input, &
688 "DFT%PRINT%AO_MATRICES/CORE_HAMILTONIAN")
689 END IF
690
691 CALL timestop(handle)
692
693 END SUBROUTINE dump_info_core_hamiltonian
694
695! **************************************************************************************************
696!> \brief (Re-)allocate matrix_h based on the template (typically the overlap matrix)
697!> \param qs_env ...
698!> \param template ...
699!> \param is_complex ...
700! **************************************************************************************************
701 SUBROUTINE qs_matrix_h_allocate(qs_env, template, is_complex)
702 TYPE(qs_environment_type) :: qs_env
703 TYPE(dbcsr_type), INTENT(in) :: template
704 LOGICAL, INTENT(in) :: is_complex
705
706 CHARACTER(LEN=default_string_length) :: headline
707 INTEGER :: img, nimages
708 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_h, matrix_h_im
709 TYPE(dft_control_type), POINTER :: dft_control
710 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
711 POINTER :: sab_orb
712 TYPE(qs_ks_env_type), POINTER :: ks_env
713
714 NULLIFY (matrix_h, matrix_h_im, sab_orb, dft_control, ks_env)
715 CALL get_qs_env(qs_env=qs_env, &
716 matrix_h_kp=matrix_h, &
717 matrix_h_im_kp=matrix_h_im, &
718 sab_orb=sab_orb, &
719 dft_control=dft_control, &
720 ks_env=ks_env)
721
722 nimages = dft_control%nimages
723 CALL dbcsr_allocate_matrix_set(matrix_h, 1, nimages)
724 headline = "CORE HAMILTONIAN MATRIX"
725 DO img = 1, nimages
726 ALLOCATE (matrix_h(1, img)%matrix)
727 CALL dbcsr_create(matrix_h(1, img)%matrix, name=trim(headline), template=template)
728 CALL cp_dbcsr_alloc_block_from_nbl(matrix_h(1, img)%matrix, sab_orb)
729 CALL dbcsr_set(matrix_h(1, img)%matrix, 0.0_dp)
730 END DO
731 CALL set_ks_env(ks_env, matrix_h_kp=matrix_h)
732
733 IF (is_complex) THEN
734 headline = "IMAGINARY PART OF CORE HAMILTONIAN MATRIX"
735 CALL dbcsr_allocate_matrix_set(matrix_h_im, 1, nimages)
736 DO img = 1, nimages
737 ALLOCATE (matrix_h_im(1, img)%matrix)
738 CALL dbcsr_create(matrix_h_im(1, img)%matrix, name=trim(headline), template=template, &
739 matrix_type=dbcsr_type_antisymmetric)
740 CALL cp_dbcsr_alloc_block_from_nbl(matrix_h_im(1, img)%matrix, sab_orb)
741 CALL dbcsr_set(matrix_h_im(1, img)%matrix, 0.0_dp)
742 END DO
743 CALL set_ks_env(ks_env, matrix_h_im_kp=matrix_h_im)
744 END IF
745
746 END SUBROUTINE qs_matrix_h_allocate
747
748! **************************************************************************************************
749!> \brief (Re-)allocates matrix_h_im from matrix_h
750!> \param qs_env ...
751! **************************************************************************************************
753 TYPE(qs_environment_type) :: qs_env
754
755 CHARACTER(LEN=default_string_length) :: headline
756 INTEGER :: image, nimages
757 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_h, matrix_h_im
758 TYPE(dbcsr_type), POINTER :: template
759 TYPE(dft_control_type), POINTER :: dft_control
760 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
761 POINTER :: sab_orb
762 TYPE(qs_ks_env_type), POINTER :: ks_env
763
764 NULLIFY (matrix_h_im, matrix_h, dft_control, template, sab_orb, ks_env)
765
766 CALL get_qs_env(qs_env, &
767 matrix_h_im_kp=matrix_h_im, &
768 matrix_h_kp=matrix_h, &
769 dft_control=dft_control, &
770 sab_orb=sab_orb, &
771 ks_env=ks_env)
772
773 nimages = dft_control%nimages
774
775 cpassert(nimages .EQ. SIZE(matrix_h, 2))
776
777 CALL dbcsr_allocate_matrix_set(matrix_h_im, 1, nimages)
778
779 DO image = 1, nimages
780 headline = "IMAGINARY CORE HAMILTONIAN MATRIX"
781 ALLOCATE (matrix_h_im(1, image)%matrix)
782 template => matrix_h(1, image)%matrix ! base on real part, but anti-symmetric
783 CALL dbcsr_create(matrix=matrix_h_im(1, image)%matrix, template=template, &
784 name=trim(headline), matrix_type=dbcsr_type_antisymmetric)
785 CALL cp_dbcsr_alloc_block_from_nbl(matrix_h_im(1, image)%matrix, sab_orb)
786 CALL dbcsr_set(matrix_h_im(1, image)%matrix, 0.0_dp)
787 END DO
788 CALL set_ks_env(ks_env, matrix_h_im_kp=matrix_h_im)
789
791
792END MODULE qs_core_hamiltonian
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
Calculation of the nuclear attraction contribution to the core Hamiltonian <a|erfc|b> :we only calcul...
Definition core_ae.F:14
subroutine, public build_core_ae(matrix_h, matrix_p, force, virial, calculate_forces, use_virial, nder, qs_kind_set, atomic_kind_set, particle_set, sab_orb, sac_ae, nimages, cell_to_index, atcore)
...
Definition core_ae.F:86
Calculation of the local pseudopotential contribution to the core Hamiltonian <a|V(local)|b> = <a|Sum...
Definition core_ppl.F:18
subroutine, public build_core_ppl(matrix_h, matrix_p, force, virial, calculate_forces, use_virial, nder, qs_kind_set, atomic_kind_set, particle_set, sab_orb, sac_ppl, nimages, cell_to_index, basis_type, deltar, atcore)
...
Definition core_ppl.F:97
Calculation of the non-local pseudopotential contribution to the core Hamiltonian <a|V(non-local)|b> ...
Definition core_ppnl.F:15
subroutine, public build_core_ppnl(matrix_h, matrix_p, force, virial, calculate_forces, use_virial, nder, qs_kind_set, atomic_kind_set, particle_set, sab_orb, sap_ppnl, eps_ppnl, nimages, cell_to_index, basis_type, deltar, matrix_l, atcore)
...
Definition core_ppnl.F:90
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
DBCSR operations in CP2K.
DBCSR output in CP2K.
subroutine, public cp_dbcsr_write_matrix_dist(matrix, output_unit, para_env)
Print the distribution of a sparse matrix.
subroutine, public cp_dbcsr_write_sparse_matrix(sparse_matrix, before, after, qs_env, para_env, first_row, last_row, first_col, last_col, scale, output_unit, omit_headers)
...
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
integer, parameter, public cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_admm_purify_none
integer, parameter, public kg_tnadd_atomic
integer, parameter, public rel_trans_atom
integer, parameter, public rel_none
integer, parameter, public do_ppl_analytic
objects that represent the structure of input sections and the data contained in an input section
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
Types needed for a Kim-Gordon-like partitioning into molecular subunits.
Calculation of the local potential contribution of the nonadditive kinetic energy <a|V(local)|b> = <a...
subroutine, public build_tnadd_mat(kg_env, matrix_p, force, virial, calculate_forces, use_virial, qs_kind_set, atomic_kind_set, particle_set, sab_orb, dbcsr_dist)
...
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
Types and basic routines needed for a kpoint calculation.
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.
contains the types and subroutines for dealing with the lri_env lri : local resolution of the identit...
Interface to the message passing library MPI.
Define the data structure for the particle information.
Calculation of overlap matrix condition numbers.
Definition qs_condnum.F:13
subroutine, public overlap_condnum(matrixkp_s, condnum, iunit, norml1, norml2, use_arnoldi, blacs_env)
Calculation of the overlap matrix Condition Number.
Definition qs_condnum.F:66
Calculation of the core Hamiltonian integral matrix <a|H|b> over Cartesian Gaussian-type functions.
subroutine, public build_core_hamiltonian_matrix(qs_env, calculate_forces)
Cosntruction of the QS Core Hamiltonian Matrix.
subroutine, public dump_info_core_hamiltonian(qs_env, calculate_forces)
Possibly prints matrices after the construction of the Core Hamiltonian Matrix.
subroutine, public qs_matrix_h_allocate_imag_from_real(qs_env)
(Re-)allocates matrix_h_im from matrix_h
subroutine, public core_matrices(qs_env, matrix_h, matrix_p, calculate_forces, nder, atcore)
...
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_pp, 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, harris_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, eeq, rhs)
Get the QUICKSTEP environment.
subroutine, public set_qs_env(qs_env, super_cell, mos, qmmm, qmmm_periodic, ewald_env, ewald_pw, mpools, rho_external, external_vxc, mask, scf_control, rel_control, qs_charges, ks_env, ks_qmmm_env, wf_history, scf_env, active_space, input, oce, rho_atom_set, rho0_atom_set, rho0_mpole, run_rtp, rtp, rhoz_set, rhoz_tot, ecoul_1c, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, efield, linres_control, xas_env, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, ls_scf_env, do_transport, transport_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, harris_env, gcp_env, mp2_env, bs_env, kg_env, force, kpoints, wanniercentres, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs)
Set the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
Calculation of kinetic energy matrix and forces.
Definition qs_kinetic.F:15
subroutine, public build_kinetic_matrix(ks_env, matrix_t, matrixkp_t, matrix_name, basis_type, sab_nl, calculate_forces, matrix_p, matrixkp_p, eps_filter, nderivative)
Calculation of the kinetic energy matrix over Cartesian Gaussian functions.
Definition qs_kinetic.F:101
subroutine, public set_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, complex_ks, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, kinetic, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_ks_im_kp, vppl, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, kpoints, sab_orb, sab_all, sac_ae, sac_ppl, sac_lri, sap_ppnl, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_vdw, sab_scp, sab_almo, sab_kp, sab_kp_nosym, task_list, task_list_soft, subsys, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env)
...
subroutine, public get_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, complex_ks, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, kinetic, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_ks_im_kp, rho, rho_xc, vppl, rho_core, rho_nlcc, rho_nlcc_g, vee, neighbor_list_id, sab_orb, sab_all, sac_ae, sac_ppl, sac_lri, sap_ppnl, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_vdw, sab_scp, sab_almo, sab_kp, sab_kp_nosym, task_list, task_list_soft, kpoints, do_kpoints, atomic_kind_set, qs_kind_set, cell, cell_ref, use_ref_cell, particle_set, energy, force, local_particles, local_molecules, molecule_kind_set, molecule_set, subsys, cp_subsys, virial, results, atprop, nkind, natom, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env, nelectron_total, nelectron_spin)
...
Define the neighbor list data types and the corresponding functionality.
Routines for the construction of the coefficients for the expansion of the atomic densities rho1_hard...
subroutine, public build_oce_matrices(intac, calculate_forces, nder, qs_kind_set, particle_set, sap_oce, eps_fit)
Set up the sparse matrix for the coefficients of one center expansions This routine uses the same log...
subroutine, public allocate_oce_set(oce_set, nkind)
Allocate and initialize the matrix set of oce coefficients.
subroutine, public create_oce_set(oce_set)
...
Calculation of overlap matrix, its derivatives and forces.
Definition qs_overlap.F:19
subroutine, public build_overlap_matrix(ks_env, matrix_s, matrixkp_s, matrix_name, nderivative, basis_type_a, basis_type_b, sab_nl, calculate_forces, matrix_p, matrixkp_p)
Calculation of the overlap matrix over Cartesian Gaussian functions.
Definition qs_overlap.F:120
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
Provides all information about an atomic kind.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
type of a logger, at the moment it contains just a print level starting at which level it should be l...
Contains all the info needed for KG runs...
Contains information about kpoints.
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
keeps the density in various representations, keeping track of which ones are valid.