(git:d5c4d39)
Loading...
Searching...
No Matches
xc_gauxc_functional.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
11 USE cell_types, ONLY: cell_type
13 USE cp_dbcsr_api, ONLY: dbcsr_add,&
29 USE iso_c_binding, ONLY: c_char,&
30 c_double,&
31 c_int,&
32 c_null_char
33 USE kinds, ONLY: default_path_length,&
35 dp
36 USE message_passing, ONLY: mp_comm_self,&
43 USE qs_kind_types, ONLY: get_qs_kind,&
44 has_nlcc,&
46 USE qs_ks_types, ONLY: qs_ks_env_type,&
48 USE qs_rho_types, ONLY: qs_rho_get,&
52 USE xc_gauxc_interface, ONLY: &
53 cp_gauxc_basisset_type, cp_gauxc_grid_type, cp_gauxc_integrator_type, &
54 cp_gauxc_molecule_type, cp_gauxc_status_type, cp_gauxc_xc_gradient_type, cp_gauxc_xc_type, &
55 gauxc_check_status, gauxc_compute_xc, gauxc_compute_xc_gradient, gauxc_create_basisset, &
56 gauxc_create_grid, gauxc_create_integrator, gauxc_create_molecule, gauxc_destroy_basisset, &
57 gauxc_destroy_grid, gauxc_destroy_integrator, gauxc_destroy_molecule, &
58 gauxc_write_basisset_hdf5, gauxc_write_molecule_hdf5
60#include "../base/base_uses.f90"
61
62 IMPLICIT NONE
63
64 PRIVATE
65
66 LOGICAL, PARAMETER :: debug_this_module = .true.
67 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'xc_gauxc_functional'
68
70
71 INTERFACE
72 INTEGER(c_int) FUNCTION c_setenv(name, value, overwrite) BIND(C, name="setenv")
73 IMPORT :: c_char, c_int
74 CHARACTER(KIND=c_char), DIMENSION(*), INTENT(IN) :: name, value
75 INTEGER(c_int), VALUE :: overwrite
76 END FUNCTION c_setenv
77
78 INTEGER(c_int) FUNCTION c_unsetenv(name) BIND(C, name="unsetenv")
79 IMPORT :: c_char, c_int
80 CHARACTER(KIND=c_char), DIMENSION(*), INTENT(IN) :: name
81 END FUNCTION c_unsetenv
82 END INTERFACE
83
84CONTAINS
85
86! **************************************************************************************************
87!> \brief Set the GauXC OneDFT atom chunk environment knob when the CP2K keyword is explicit.
88!> \param atom_chunk_size ...
89!> \param is_explicit ...
90! **************************************************************************************************
91 SUBROUTINE set_gauxc_onedft_atom_chunk_env(atom_chunk_size, is_explicit)
92 INTEGER, INTENT(IN) :: atom_chunk_size
93 LOGICAL, INTENT(IN) :: is_explicit
94
95 CHARACTER(LEN=32) :: chunk_value
96 INTEGER(c_int) :: ierr
97
98 IF (.NOT. is_explicit) RETURN
99
100 IF (atom_chunk_size < 0) THEN
101 ierr = c_unsetenv("GAUXC_ONEDFT_ATOM_CHUNK_SIZE"//c_null_char)
102 ELSE
103 WRITE (chunk_value, '(I0)') atom_chunk_size
104 ierr = c_setenv( &
105 "GAUXC_ONEDFT_ATOM_CHUNK_SIZE"//c_null_char, &
106 trim(chunk_value)//c_null_char, &
107 1_c_int)
108 END IF
109 IF (ierr /= 0_c_int) THEN
110 CALL cp_abort(__location__, &
111 "Could not set GAUXC_ONEDFT_ATOM_CHUNK_SIZE for GauXC OneDFT/SKALA.")
112 END IF
113 END SUBROUTINE set_gauxc_onedft_atom_chunk_env
114
115! **************************************************************************************************
116!> \brief ...
117!> \param dbcsr_mat ...
118!> \param dense_mat ...
119! **************************************************************************************************
120 SUBROUTINE dbcsr_to_dense(dbcsr_mat, dense_mat)
123 dbcsr_iterator_stop, dbcsr_type_antisymmetric, dbcsr_type_symmetric
124 TYPE(dbcsr_p_type), INTENT(IN) :: dbcsr_mat
125 REAL(c_double), ALLOCATABLE, DIMENSION(:, :), &
126 INTENT(INOUT) :: dense_mat
127
128 CHARACTER :: matrix_type
129 INTEGER :: col, col_end, col_start, icol, irow, &
130 nblkcols_total, nblkrows_total, ncols, &
131 nrows, row, row_end, row_start
132 INTEGER, ALLOCATABLE, DIMENSION(:) :: c_offset, r_offset
133 INTEGER, DIMENSION(:), POINTER :: col_blk_size, row_blk_size
134 LOGICAL :: transposed
135 REAL(c_double), POINTER :: block(:, :)
136 TYPE(dbcsr_iterator_type) :: iter
137
138 CALL dbcsr_get_info(dbcsr_mat%matrix, &
139 row_blk_size=row_blk_size, &
140 col_blk_size=col_blk_size, &
141 nblkrows_total=nblkrows_total, &
142 nblkcols_total=nblkcols_total, &
143 nfullrows_total=nrows, &
144 nfullcols_total=ncols)
145 matrix_type = dbcsr_get_matrix_type(dbcsr_mat%matrix)
146
147 IF (.NOT. ALLOCATED(dense_mat)) THEN
148 ALLOCATE (dense_mat(nrows, ncols))
149 ELSE IF (.NOT. all(shape(dense_mat) == [nrows, ncols])) THEN
150 DEALLOCATE (dense_mat)
151 ALLOCATE (dense_mat(nrows, ncols))
152 ELSE
153 cpassert(all(shape(dense_mat) == [nrows, ncols]))
154 END IF
155 dense_mat = 0._dp
156
157 ALLOCATE (r_offset(nblkrows_total), c_offset(nblkcols_total))
158
159 r_offset(1) = 1
160 DO row = 2, nblkrows_total
161 r_offset(row) = r_offset(row - 1) + row_blk_size(row - 1)
162 END DO
163 c_offset(1) = 1
164 DO col = 2, nblkcols_total
165 c_offset(col) = c_offset(col - 1) + col_blk_size(col - 1)
166 END DO
167
168 CALL dbcsr_iterator_start(iter, dbcsr_mat%matrix)
169 DO WHILE (dbcsr_iterator_blocks_left(iter))
170 CALL dbcsr_iterator_next_block(iter, irow, icol, block, transposed=transposed)
171 row_start = r_offset(irow)
172 row_end = row_start + row_blk_size(irow) - 1
173 col_start = c_offset(icol)
174 col_end = col_start + col_blk_size(icol) - 1
175 IF (transposed) THEN
176 dense_mat(row_start:row_end, col_start:col_end) = transpose(block)
177 IF (irow /= icol) THEN
178 IF (matrix_type == dbcsr_type_symmetric) THEN
179 dense_mat(col_start:col_end, row_start:row_end) = block
180 ELSE IF (matrix_type == dbcsr_type_antisymmetric) THEN
181 dense_mat(col_start:col_end, row_start:row_end) = -block
182 END IF
183 END IF
184 ELSE
185 dense_mat(row_start:row_end, col_start:col_end) = block
186 IF (irow /= icol) THEN
187 IF (matrix_type == dbcsr_type_symmetric) THEN
188 dense_mat(col_start:col_end, row_start:row_end) = transpose(block)
189 ELSE IF (matrix_type == dbcsr_type_antisymmetric) THEN
190 dense_mat(col_start:col_end, row_start:row_end) = -transpose(block)
191 END IF
192 END IF
193 END IF
194 END DO
195 CALL dbcsr_iterator_stop(iter)
196
197 DEALLOCATE (r_offset, c_offset)
198
199 END SUBROUTINE dbcsr_to_dense
200
201! ******, ***********************************************************************************
202!> \brief Convert a dense symmetric matrix to a DBCSR matrix with full upper block structure.
203!> This creates all upper-triangular blocks, not just those present in a template.
204!> This is needed because GauXC computes VXC for the full dense density matrix.
205!> \param dense_mat Input dense matrix
206!> \param template_dbcsr Template DBCSR matrix for distribution and block sizes
207!> \return dbcsr_mat Output DBCSR matrix with full upper block structure
208! **************************************************************************************************
209 FUNCTION dense_to_dbcsr(dense_mat, template_dbcsr) RESULT(dbcsr_mat)
210 USE cp_dbcsr_api, ONLY: &
211 dbcsr_create, &
217 dbcsr_init_p, &
220 dbcsr_type_symmetric, &
222 REAL(c_double), DIMENSION(:, :), INTENT(IN) :: dense_mat
223 TYPE(dbcsr_p_type), INTENT(IN) :: template_dbcsr
224 TYPE(dbcsr_p_type) :: dbcsr_mat
225
226 INTEGER :: col, icol, irow, mynode, nblkcols_total, &
227 nblkrows_total, ncols, nrows, owner, &
228 row
229 INTEGER, ALLOCATABLE, DIMENSION(:) :: c_offset, r_offset
230 INTEGER, DIMENSION(:), POINTER :: col_blk_size, row_blk_size
231 TYPE(dbcsr_distribution_type) :: dist
232
233 CALL dbcsr_get_info(template_dbcsr%matrix, &
234 row_blk_size=row_blk_size, &
235 col_blk_size=col_blk_size, &
236 nblkrows_total=nblkrows_total, &
237 nblkcols_total=nblkcols_total, &
238 nfullrows_total=nrows, &
239 nfullcols_total=ncols, &
240 distribution=dist)
241 CALL dbcsr_distribution_get(dist, mynode=mynode)
242
243 cpassert(nrows == SIZE(dense_mat, 1))
244 cpassert(ncols == SIZE(dense_mat, 2))
245
246 CALL dbcsr_init_p(dbcsr_mat%matrix)
247 CALL dbcsr_create(dbcsr_mat%matrix, &
248 template=template_dbcsr%matrix, &
249 name="VXC from GauXC (dense)", &
250 matrix_type=dbcsr_type_symmetric)
251 CALL dbcsr_work_create(dbcsr_mat%matrix, work_mutable=.true.)
252
253 ALLOCATE (r_offset(nblkrows_total), c_offset(nblkcols_total))
254
255 r_offset(1) = 1
256 DO row = 2, nblkrows_total
257 r_offset(row) = r_offset(row - 1) + row_blk_size(row - 1)
258 END DO
259 c_offset(1) = 1
260 DO col = 2, nblkcols_total
261 c_offset(col) = c_offset(col - 1) + col_blk_size(col - 1)
262 END DO
263
264 DO irow = 1, nblkrows_total
265 DO icol = 1, nblkcols_total
266 IF (irow > icol) cycle
267 CALL dbcsr_get_stored_coordinates(dbcsr_mat%matrix, irow, icol, owner)
268 IF (owner /= mynode) cycle
269 CALL dbcsr_put_block(dbcsr_mat%matrix, irow, icol, &
270 0.5_dp*( &
271 dense_mat(r_offset(irow):r_offset(irow) + row_blk_size(irow) - 1, &
272 c_offset(icol):c_offset(icol) + col_blk_size(icol) - 1) + &
273 transpose(dense_mat(r_offset(icol):r_offset(icol) + row_blk_size(icol) - 1, &
274 c_offset(irow):c_offset(irow) + col_blk_size(irow) - 1))))
275 END DO
276 END DO
277
278 CALL dbcsr_finalize(dbcsr_mat%matrix)
279
280 DEALLOCATE (r_offset, c_offset)
281
282 END FUNCTION dense_to_dbcsr
283
284! **************************************************************************************************
285!> \brief ...
286!> \param xc_section ...
287!> \return ...
288! **************************************************************************************************
289 FUNCTION get_gauxc_functional(xc_section) RESULT(gauxc_functional_section)
290 TYPE(section_vals_type), INTENT(in), POINTER :: xc_section
291 TYPE(section_vals_type), POINTER :: gauxc_functional_section
292
293 INTEGER :: ifun
294 TYPE(section_vals_type), POINTER :: functionals, xc_fun
295
296 NULLIFY (gauxc_functional_section)
297
298 functionals => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
299 IF (.NOT. ASSOCIATED(functionals)) THEN
300 cpabort("XC_FUNCTIONAL section not found")
301 END IF
302
303 ifun = 0
304 DO
305 ifun = ifun + 1
306 xc_fun => section_vals_get_subs_vals2(functionals, i_section=ifun)
307 IF (.NOT. ASSOCIATED(xc_fun)) EXIT
308 IF (xc_fun%section%name /= "GAUXC" .OR. ifun > 1) THEN
309 cpabort("GauXC functionals are mutually exclusive with any other functional.")
310 END IF
311 gauxc_functional_section => xc_fun
312 END DO
313
314 IF (.NOT. ASSOCIATED(gauxc_functional_section)) THEN
315 cpabort("No XC functional found in XC_FUNCTIONAL section")
316 END IF
317 END FUNCTION get_gauxc_functional
318
319! **************************************************************************************************
320!> \brief ...
321!> \param xc_section ...
322!> \return ...
323! **************************************************************************************************
324 FUNCTION xc_section_uses_gauxc(xc_section) RESULT(uses_gauxc)
325 TYPE(section_vals_type), INTENT(in), POINTER :: xc_section
326 LOGICAL :: uses_gauxc
327
328 INTEGER :: ifun
329 TYPE(section_vals_type), POINTER :: functionals, xc_fun
330
331 uses_gauxc = .false.
332 IF (.NOT. ASSOCIATED(xc_section)) RETURN
333
334 functionals => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
335 IF (.NOT. ASSOCIATED(functionals)) RETURN
336
337 ifun = 0
338 DO
339 ifun = ifun + 1
340 xc_fun => section_vals_get_subs_vals2(functionals, i_section=ifun)
341 IF (.NOT. ASSOCIATED(xc_fun)) EXIT
342 IF (xc_fun%section%name == "GAUXC") THEN
343 uses_gauxc = .true.
344 EXIT
345 END IF
346 END DO
347
348 END FUNCTION xc_section_uses_gauxc
349
350! **************************************************************************************************
351!> \brief Reject unsupported pseudopotential variants in GauXC GAPW mode
352!> \param qs_kind_set ...
353! **************************************************************************************************
354 SUBROUTINE ensure_gauxc_gapw_all_electron(qs_kind_set)
355 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
356
357 INTEGER :: ikind
358 TYPE(gth_potential_type), POINTER :: gth_potential
359 TYPE(sgp_potential_type), POINTER :: sgp_potential
360
361 cpassert(ASSOCIATED(qs_kind_set))
362
363 DO ikind = 1, SIZE(qs_kind_set)
364 NULLIFY (gth_potential, sgp_potential)
365 CALL get_qs_kind(qs_kind_set(ikind), &
366 gth_potential=gth_potential, &
367 sgp_potential=sgp_potential)
368 IF (ASSOCIATED(gth_potential) .OR. ASSOCIATED(sgp_potential)) THEN
369 CALL cp_abort(__location__, &
370 "GauXC with METHOD GAPW currently supports all-electron potentials only. "// &
371 "Use POTENTIAL ALL for GAPW validation or METHOD GPW with pseudopotentials.")
372 END IF
373 END DO
374
375 END SUBROUTINE ensure_gauxc_gapw_all_electron
376
377! **************************************************************************************************
378!> \brief Check the current periodic scope of the CP2K-GauXC bridge
379!> \param dft_control ...
380!> \param cell ...
381!> \param qs_kind_set ...
382!> \param do_kpoints ...
383!> \param periodic_reference ...
384!> \note This path keeps isolated validation cells usable under PERIODIC XYZ.
385!> It intentionally does not implement compact periodic GauXC quadrature.
386! **************************************************************************************************
387 SUBROUTINE ensure_gauxc_periodic_reference_scope( &
388 dft_control, cell, qs_kind_set, do_kpoints, periodic_reference)
389 TYPE(dft_control_type), POINTER :: dft_control
390 TYPE(cell_type), POINTER :: cell
391 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
392 LOGICAL, INTENT(IN) :: do_kpoints, periodic_reference
393
394 INTEGER :: ikind
395 LOGICAL :: is_periodic
396 TYPE(gth_potential_type), POINTER :: gth_potential
397 TYPE(sgp_potential_type), POINTER :: sgp_potential
398
399 cpassert(ASSOCIATED(dft_control))
400 cpassert(ASSOCIATED(qs_kind_set))
401
402 is_periodic = .false.
403 IF (ASSOCIATED(cell)) is_periodic = any(cell%perd /= 0)
404
405 IF (do_kpoints) THEN
406 CALL cp_abort(__location__, &
407 "GauXC currently supports only Gamma-only density matrices in CP2K. "// &
408 "Periodic k-point density matrices require a dedicated GauXC periodic interface.")
409 END IF
410 IF (dft_control%nimages /= 1) THEN
411 CALL cp_abort(__location__, &
412 "GauXC currently supports only a single AO image in CP2K. "// &
413 "Periodic neighbour-cell AO blocks require a dedicated GauXC periodic interface.")
414 END IF
415 IF (.NOT. is_periodic) RETURN
416
417 IF (.NOT. periodic_reference) THEN
418 CALL cp_abort(__location__, &
419 "Periodic GauXC calculations in CP2K require GAUXC%PERIODIC_REFERENCE T. "// &
420 "This opt-in documents that the current path is only an isolated-cell, "// &
421 "Gamma-only, single-image METHOD GPW reference path using GauXC molecular "// &
422 "quadrature, not a dedicated periodic GauXC interface.")
423 END IF
424
425 IF (.NOT. all(cell%perd == 1)) THEN
426 CALL cp_abort(__location__, &
427 "The current GauXC isolated-cell reference path supports only PERIODIC XYZ. "// &
428 "Partial periodicity requires a dedicated GauXC periodic interface.")
429 END IF
430 IF (.NOT. dft_control%qs_control%gpw) THEN
431 CALL cp_abort(__location__, &
432 "The current GauXC isolated-cell reference path is limited to METHOD GPW with GTH "// &
433 "pseudopotentials. GAPW, GAPW_XC, and other QS methods are not supported here.")
434 END IF
435
436 DO ikind = 1, SIZE(qs_kind_set)
437 NULLIFY (gth_potential, sgp_potential)
438 CALL get_qs_kind(qs_kind_set(ikind), &
439 gth_potential=gth_potential, &
440 sgp_potential=sgp_potential)
441 IF (.NOT. ASSOCIATED(gth_potential) .OR. ASSOCIATED(sgp_potential)) THEN
442 CALL cp_abort(__location__, &
443 "The current GauXC isolated-cell reference path is limited to GTH pseudopotentials. "// &
444 "Use non-periodic all-electron GAPW validation for molecular GAPW cases.")
445 END IF
446 END DO
447
448 END SUBROUTINE ensure_gauxc_periodic_reference_scope
449
450! **************************************************************************************************
451!> \brief adds a replicated GauXC energy gradient to the local CP2K force accumulator
452!> \param exc_grad ...
453!> \param force ...
454!> \param atomic_kind_set ...
455!> \param para_env ...
456! **************************************************************************************************
457 SUBROUTINE add_gauxc_gradient_to_force(exc_grad, force, atomic_kind_set, para_env)
458 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: exc_grad
459 TYPE(qs_force_type), DIMENSION(:), POINTER :: force
460 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
461 TYPE(mp_para_env_type), POINTER :: para_env
462
463 INTEGER :: ia, iatom, ikind, natom_kind
464 TYPE(atomic_kind_type), POINTER :: atomic_kind
465
466 cpassert(ASSOCIATED(force))
467 cpassert(ASSOCIATED(atomic_kind_set))
468
469 IF (para_env%mepos /= 0) RETURN
470
471 DO ikind = 1, SIZE(atomic_kind_set, 1)
472 atomic_kind => atomic_kind_set(ikind)
473 CALL get_atomic_kind(atomic_kind=atomic_kind, natom=natom_kind)
474 DO ia = 1, natom_kind
475 iatom = atomic_kind%atom_list(ia)
476 force(ikind)%rho_elec(:, ia) = force(ikind)%rho_elec(:, ia) + &
477 exc_grad(3*iatom - 2:3*iatom)
478 END DO
479 END DO
480
481 END SUBROUTINE add_gauxc_gradient_to_force
482
483! **************************************************************************************************
484!> \brief compute a GauXC XC energy for diagnostic finite differences
485!> \param particle_set_eval ...
486!> \param qs_kind_set ...
487!> \param density_scalar ...
488!> \param nspins ...
489!> \param model_name ...
490!> \param xc_fun_name ...
491!> \param grid_type ...
492!> \param radial_quadrature ...
493!> \param pruning_scheme ...
494!> \param lb_exec_space ...
495!> \param int_exec_space ...
496!> \param batch_size ...
497!> \param device_runtime_fill_fraction ...
498!> \param exc ...
499!> \param density_zeta ...
500! **************************************************************************************************
501 SUBROUTINE gauxc_xc_energy_for_particles( &
502 particle_set_eval, qs_kind_set, density_scalar, nspins, model_name, &
503 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
504 int_exec_space, batch_size, device_runtime_fill_fraction, exc, density_zeta)
505 TYPE(particle_type), DIMENSION(:), INTENT(IN) :: particle_set_eval
506 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
507 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: density_scalar
508 INTEGER, INTENT(IN) :: nspins
509 CHARACTER(len=*), INTENT(IN) :: model_name, xc_fun_name, grid_type, &
510 radial_quadrature, pruning_scheme, &
511 lb_exec_space, int_exec_space
512 INTEGER, INTENT(IN) :: batch_size
513 REAL(kind=dp), INTENT(IN) :: device_runtime_fill_fraction
514 REAL(kind=dp), INTENT(OUT) :: exc
515 REAL(kind=dp), DIMENSION(:, :), INTENT(IN), &
516 OPTIONAL :: density_zeta
517
518 TYPE(cp_gauxc_basisset_type) :: gauxc_basis_fd
519 TYPE(cp_gauxc_grid_type) :: gauxc_grid_fd
520 TYPE(cp_gauxc_integrator_type) :: gauxc_integrator_fd
521 TYPE(cp_gauxc_molecule_type) :: gauxc_mol_fd
522 TYPE(cp_gauxc_status_type) :: gauxc_status
523 TYPE(cp_gauxc_xc_type) :: gauxc_xc_result
524
525 gauxc_mol_fd = gauxc_create_molecule(particle_set_eval, gauxc_status)
526 CALL gauxc_check_status(gauxc_status)
527 gauxc_basis_fd = gauxc_create_basisset(qs_kind_set, particle_set_eval, gauxc_status)
528 CALL gauxc_check_status(gauxc_status)
529 gauxc_grid_fd = gauxc_create_grid( &
530 gauxc_mol_fd, &
531 gauxc_basis_fd, &
532 grid_type, &
533 radial_quadrature, &
534 pruning_scheme, &
535 lb_exec_space, &
536 batch_size, &
537 device_runtime_fill_fraction, &
538 gauxc_status, &
539 mpi_comm=mp_comm_self%get_handle())
540 CALL gauxc_check_status(gauxc_status)
541 gauxc_integrator_fd = gauxc_create_integrator( &
542 trim(xc_fun_name), &
543 gauxc_grid_fd, &
544 int_exec_space, &
545 nspins, &
546 gauxc_status)
547 CALL gauxc_check_status(gauxc_status)
548
549 IF (nspins == 1) THEN
550 gauxc_xc_result = gauxc_compute_xc( &
551 gauxc_integrator_fd, &
552 density_scalar, &
553 nspins=nspins, &
554 status=gauxc_status, &
555 model=trim(model_name))
556 ELSE
557 cpassert(nspins == 2)
558 cpassert(PRESENT(density_zeta))
559 gauxc_xc_result = gauxc_compute_xc( &
560 gauxc_integrator_fd, &
561 density_scalar, &
562 density_zeta, &
563 nspins, &
564 gauxc_status, &
565 model=trim(model_name))
566 END IF
567 CALL gauxc_check_status(gauxc_status)
568 exc = gauxc_xc_result%exc
569
570 IF (ALLOCATED(gauxc_xc_result%vxc_scalar)) DEALLOCATE (gauxc_xc_result%vxc_scalar)
571 IF (ALLOCATED(gauxc_xc_result%vxc_zeta)) DEALLOCATE (gauxc_xc_result%vxc_zeta)
572
573 CALL gauxc_destroy_integrator(gauxc_integrator_fd, gauxc_status)
574 CALL gauxc_check_status(gauxc_status)
575 CALL gauxc_destroy_grid(gauxc_grid_fd, gauxc_status)
576 CALL gauxc_check_status(gauxc_status)
577 CALL gauxc_destroy_basisset(gauxc_basis_fd, gauxc_status)
578 CALL gauxc_check_status(gauxc_status)
579 CALL gauxc_destroy_molecule(gauxc_mol_fd, gauxc_status)
580 CALL gauxc_check_status(gauxc_status)
581
582 END SUBROUTINE gauxc_xc_energy_for_particles
583
584! **************************************************************************************************
585!> \brief compute a finite-difference GauXC XC nuclear gradient at fixed density
586!> \param particle_set ...
587!> \param qs_kind_set ...
588!> \param density_scalar ...
589!> \param nspins ...
590!> \param model_name ...
591!> \param xc_fun_name ...
592!> \param grid_type ...
593!> \param radial_quadrature ...
594!> \param pruning_scheme ...
595!> \param lb_exec_space ...
596!> \param int_exec_space ...
597!> \param batch_size ...
598!> \param device_runtime_fill_fraction ...
599!> \param dx ...
600!> \param para_env ...
601!> \param exc_grad ...
602!> \param density_zeta ...
603! **************************************************************************************************
604 SUBROUTINE gauxc_xc_gradient_fd( &
605 particle_set, qs_kind_set, density_scalar, nspins, model_name, &
606 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
607 int_exec_space, batch_size, device_runtime_fill_fraction, dx, para_env, exc_grad, &
608 density_zeta)
609 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
610 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
611 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: density_scalar
612 INTEGER, INTENT(IN) :: nspins
613 CHARACTER(len=*), INTENT(IN) :: model_name, xc_fun_name, grid_type, &
614 radial_quadrature, pruning_scheme, &
615 lb_exec_space, int_exec_space
616 INTEGER, INTENT(IN) :: batch_size
617 REAL(kind=dp), INTENT(IN) :: device_runtime_fill_fraction, dx
618 TYPE(mp_para_env_type), POINTER :: para_env
619 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), &
620 INTENT(OUT) :: exc_grad
621 REAL(kind=dp), DIMENSION(:, :), INTENT(IN), &
622 OPTIONAL :: density_zeta
623
624 INTEGER :: iatom, idir
625 REAL(kind=dp) :: xc_minus, xc_plus
626 TYPE(particle_type), ALLOCATABLE, DIMENSION(:) :: particle_set_minus, particle_set_plus
627
628 cpassert(ASSOCIATED(particle_set))
629 cpassert(dx > 0.0_dp)
630
631 ALLOCATE (exc_grad(3*SIZE(particle_set)))
632 exc_grad = 0.0_dp
633
634 IF (para_env%mepos == 0) THEN
635 ALLOCATE (particle_set_minus(SIZE(particle_set)), particle_set_plus(SIZE(particle_set)))
636
637 DO iatom = 1, SIZE(particle_set)
638 DO idir = 1, 3
639 particle_set_minus = particle_set
640 particle_set_plus = particle_set
641 particle_set_minus(iatom)%r(idir) = particle_set_minus(iatom)%r(idir) - dx
642 particle_set_plus(iatom)%r(idir) = particle_set_plus(iatom)%r(idir) + dx
643 IF (PRESENT(density_zeta)) THEN
644 CALL gauxc_xc_energy_for_particles( &
645 particle_set_plus, qs_kind_set, density_scalar, nspins, model_name, &
646 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
647 int_exec_space, batch_size, device_runtime_fill_fraction, xc_plus, &
648 density_zeta=density_zeta)
649 CALL gauxc_xc_energy_for_particles( &
650 particle_set_minus, qs_kind_set, density_scalar, nspins, model_name, &
651 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
652 int_exec_space, batch_size, device_runtime_fill_fraction, xc_minus, &
653 density_zeta=density_zeta)
654 ELSE
655 CALL gauxc_xc_energy_for_particles( &
656 particle_set_plus, qs_kind_set, density_scalar, nspins, model_name, &
657 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
658 int_exec_space, batch_size, device_runtime_fill_fraction, xc_plus)
659 CALL gauxc_xc_energy_for_particles( &
660 particle_set_minus, qs_kind_set, density_scalar, nspins, model_name, &
661 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
662 int_exec_space, batch_size, device_runtime_fill_fraction, xc_minus)
663 END IF
664 exc_grad(3*iatom - 3 + idir) = (xc_plus - xc_minus)/(2.0_dp*dx)
665 END DO
666 END DO
667
668 DEALLOCATE (particle_set_minus, particle_set_plus)
669 END IF
670
671 CALL para_env%bcast(exc_grad, 0)
672
673 END SUBROUTINE gauxc_xc_gradient_fd
674
675! **************************************************************************************************
676!> \brief finite-difference check of the molecular GauXC XC virial diagnostic
677!> \param exc_grad ...
678!> \param particle_set ...
679!> \param qs_kind_set ...
680!> \param density_scalar ...
681!> \param nspins ...
682!> \param model_name ...
683!> \param xc_fun_name ...
684!> \param grid_type ...
685!> \param radial_quadrature ...
686!> \param pruning_scheme ...
687!> \param lb_exec_space ...
688!> \param int_exec_space ...
689!> \param batch_size ...
690!> \param device_runtime_fill_fraction ...
691!> \param dx ...
692!> \param para_env ...
693!> \param density_zeta ...
694! **************************************************************************************************
695 SUBROUTINE debug_gauxc_molecular_virial( &
696 exc_grad, particle_set, qs_kind_set, density_scalar, nspins, model_name, &
697 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
698 int_exec_space, batch_size, device_runtime_fill_fraction, dx, para_env, density_zeta)
699 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: exc_grad
700 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
701 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
702 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: density_scalar
703 INTEGER, INTENT(IN) :: nspins
704 CHARACTER(len=*), INTENT(IN) :: model_name, xc_fun_name, grid_type, &
705 radial_quadrature, pruning_scheme, &
706 lb_exec_space, int_exec_space
707 INTEGER, INTENT(IN) :: batch_size
708 REAL(kind=dp), INTENT(IN) :: device_runtime_fill_fraction, dx
709 TYPE(mp_para_env_type), POINTER :: para_env
710 REAL(kind=dp), DIMENSION(:, :), INTENT(IN), &
711 OPTIONAL :: density_zeta
712
713 INTEGER :: iatom, iw
714 REAL(kind=dp) :: analytic_trace, diff_trace, &
715 numerical_trace, xc_minus, xc_plus
716 REAL(kind=dp), DIMENSION(3) :: center, displacement, grad
717 TYPE(particle_type), ALLOCATABLE, DIMENSION(:) :: particle_set_minus, particle_set_plus
718
719 cpassert(ASSOCIATED(particle_set))
720 cpassert(SIZE(exc_grad) == 3*SIZE(particle_set))
721
722 IF (para_env%mepos /= 0) RETURN
723
724 center = 0.0_dp
725 DO iatom = 1, SIZE(particle_set)
726 center = center + particle_set(iatom)%r
727 END DO
728 center = center/real(SIZE(particle_set), dp)
729
730 ALLOCATE (particle_set_minus(SIZE(particle_set)), particle_set_plus(SIZE(particle_set)))
731 particle_set_minus = particle_set
732 particle_set_plus = particle_set
733
734 analytic_trace = 0.0_dp
735 DO iatom = 1, SIZE(particle_set)
736 grad = exc_grad(3*iatom - 2:3*iatom)
737 displacement = particle_set(iatom)%r - center
738 analytic_trace = analytic_trace + dot_product(grad, displacement)
739 particle_set_minus(iatom)%r = center + (1.0_dp - dx)*displacement
740 particle_set_plus(iatom)%r = center + (1.0_dp + dx)*displacement
741 END DO
742 analytic_trace = analytic_trace/3.0_dp
743
744 IF (PRESENT(density_zeta)) THEN
745 CALL gauxc_xc_energy_for_particles( &
746 particle_set_plus, qs_kind_set, density_scalar, nspins, model_name, &
747 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
748 int_exec_space, batch_size, device_runtime_fill_fraction, xc_plus, &
749 density_zeta=density_zeta)
750 CALL gauxc_xc_energy_for_particles( &
751 particle_set_minus, qs_kind_set, density_scalar, nspins, model_name, &
752 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
753 int_exec_space, batch_size, device_runtime_fill_fraction, xc_minus, &
754 density_zeta=density_zeta)
755 ELSE
756 CALL gauxc_xc_energy_for_particles( &
757 particle_set_plus, qs_kind_set, density_scalar, nspins, model_name, &
758 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
759 int_exec_space, batch_size, device_runtime_fill_fraction, xc_plus)
760 CALL gauxc_xc_energy_for_particles( &
761 particle_set_minus, qs_kind_set, density_scalar, nspins, model_name, &
762 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, lb_exec_space, &
763 int_exec_space, batch_size, device_runtime_fill_fraction, xc_minus)
764 END IF
765
766 numerical_trace = (xc_plus - xc_minus)/(2.0_dp*dx)/3.0_dp
767 diff_trace = analytic_trace - numerical_trace
768
770 IF (iw > 0) THEN
771 WRITE (unit=iw, fmt="(/,T2,A,1X,ES11.4)") &
772 "GAUXC| Molecular XC virial finite-difference dx", dx
773 WRITE (unit=iw, fmt="(T2,A,3(1X,ES19.11))") &
774 "GAUXC| Molecular XC virial FD 1/3 Trace", &
775 analytic_trace, numerical_trace, diff_trace
776 END IF
777
778 DEALLOCATE (particle_set_minus, particle_set_plus)
779
780 END SUBROUTINE debug_gauxc_molecular_virial
781
782! **************************************************************************************************
783!> \brief prints a force-based molecular XC virial diagnostic from GauXC gradients
784!> \param exc_grad ...
785!> \param particle_set ...
786!> \param para_env ...
787! **************************************************************************************************
788 SUBROUTINE print_gauxc_molecular_virial(exc_grad, particle_set, para_env)
789 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: exc_grad
790 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
791 TYPE(mp_para_env_type), POINTER :: para_env
792
793 CHARACTER(len=1), DIMENSION(3), PARAMETER :: label = ["x", "y", "z"]
794
795 INTEGER :: i, iatom, iw, j
796 REAL(kind=dp), DIMENSION(3) :: center, displacement, grad, grad_sum
797 REAL(kind=dp), DIMENSION(3, 3) :: molecular_virial
798
799 cpassert(ASSOCIATED(particle_set))
800 cpassert(SIZE(exc_grad) == 3*SIZE(particle_set))
801
802 IF (para_env%mepos /= 0) RETURN
803
804 center = 0.0_dp
805 DO iatom = 1, SIZE(particle_set)
806 center = center + particle_set(iatom)%r
807 END DO
808 center = center/real(SIZE(particle_set), dp)
809
810 grad_sum = 0.0_dp
811 molecular_virial = 0.0_dp
812 DO iatom = 1, SIZE(particle_set)
813 grad = exc_grad(3*iatom - 2:3*iatom)
814 displacement = particle_set(iatom)%r - center
815 grad_sum = grad_sum + grad
816 DO i = 1, 3
817 DO j = 1, 3
818 molecular_virial(i, j) = molecular_virial(i, j) + grad(i)*displacement(j)
819 END DO
820 END DO
821 END DO
822
824 IF (iw <= 0) RETURN
825
826 WRITE (unit=iw, fmt="(/,T2,A)") &
827 "GAUXC| Molecular XC gradient virial diagnostic [a.u.]"
828 WRITE (unit=iw, fmt="(T2,A,T20,A,T40,A,T60,A)") "GAUXC|", "x", "y", "z"
829 DO i = 1, 3
830 WRITE (unit=iw, fmt="(T2,A,1X,A1,3(1X,ES19.11))") &
831 "GAUXC|", label(i), molecular_virial(i, :)
832 END DO
833 WRITE (unit=iw, fmt="(T2,A,1X,ES19.11)") &
834 "GAUXC| Molecular XC gradient virial 1/3 Trace", &
835 (molecular_virial(1, 1) + molecular_virial(2, 2) + molecular_virial(3, 3))/3.0_dp
836 WRITE (unit=iw, fmt="(T2,A,3(1X,ES19.11))") &
837 "GAUXC| Molecular XC gradient sum", grad_sum
838 WRITE (unit=iw, fmt="(T2,A)") &
839 "GAUXC| Diagnostic only; this is not an analytical periodic stress tensor."
840
841 END SUBROUTINE print_gauxc_molecular_virial
842
843! **************************************************************************************************
844!> \brief Return information about the Skala functional
845!> \param functional section containing the SKALA subsection
846!> \param lsd if you are using lsd or lda
847!> \param reference the reference to the article where the functional is explained
848!> \param shortform the short definition of the functional
849!> \param needs the flags corresponding to the inputs needed by this
850!> functional are set to true (the flags not needed aren't touched)
851!> \param max_deriv the maximal derivative available
852! **************************************************************************************************
853 SUBROUTINE skala_info(functional, lsd, reference, shortform, needs, max_deriv)
854 TYPE(section_vals_type), POINTER :: functional
855 LOGICAL, INTENT(in) :: lsd
856 CHARACTER(LEN=*), INTENT(OUT), OPTIONAL :: reference, shortform
857 TYPE(xc_rho_cflags_type), INTENT(inout), OPTIONAL :: needs
858 INTEGER, INTENT(out), OPTIONAL :: max_deriv
859
860 CHARACTER(len=default_path_length) :: model_key, model_name
861 CHARACTER(len=default_string_length) :: xc_fun_name
862 LOGICAL :: native_grid
863
864 CALL section_vals_val_get(functional, "FUNCTIONAL", c_val=xc_fun_name)
865 CALL section_vals_val_get(functional, "MODEL", c_val=model_name)
866 CALL section_vals_val_get(functional, "NATIVE_GRID", l_val=native_grid)
867 model_key = adjustl(model_name)
868 CALL uppercase(model_key)
869
870 IF (PRESENT(reference)) THEN
871 IF (trim(model_key) == "NONE" .OR. trim(model_key) == "") THEN
872 reference = "Functional computed by GauXC (underlying: "//trim(xc_fun_name)//")"
873 ELSE
874 reference = "Functional computed by GauXC OneDFT model "//trim(model_name)
875 END IF
876 END IF
877 IF (PRESENT(shortform)) THEN
878 IF (trim(model_key) == "NONE" .OR. trim(model_key) == "") THEN
879 shortform = "GAUXC ("//trim(xc_fun_name)//")"
880 ELSE
881 shortform = "GAUXC OneDFT"
882 END IF
883 END IF
884 IF (PRESENT(needs)) THEN
885 IF (native_grid .AND. trim(model_key) /= "NONE" .AND. trim(model_key) /= "") THEN
886 IF (lsd) THEN
887 needs%rho_spin = .true.
888 needs%drho_spin = .true.
889 needs%tau_spin = .true.
890 ELSE
891 needs%rho = .true.
892 needs%drho = .true.
893 needs%tau = .true.
894 END IF
895 ELSE
896 needs%rho = .true.
897 IF (lsd) THEN
898 needs%rho_spin = .true.
899 END IF
900 END IF
901 END IF
902 IF (PRESENT(max_deriv)) max_deriv = 1
903
904 END SUBROUTINE skala_info
905
906! GauXC uses replicated dense density and VXC matrices. The DBCSR density matrix
907! is distributed over MPI ranks, so apply_gauxc allreduces the dense copy before
908! passing it to GauXC.
909
910! **************************************************************************************************
911!> \brief ...
912!> \param qs_env ...
913!> \param xc_section ...
914!> \param calculate_forces ...
915! **************************************************************************************************
916 SUBROUTINE apply_gauxc(qs_env, xc_section, calculate_forces)
917 TYPE(qs_environment_type), INTENT(in), POINTER :: qs_env
918 TYPE(section_vals_type), INTENT(in), POINTER :: xc_section
919 LOGICAL, INTENT(IN) :: calculate_forces
920
921 CHARACTER(len=*), PARAMETER :: gapw_xc_abort_message = &
922 "GauXC with METHOD GAPW_XC is not supported yet. "// &
923 "The GAPW_XC one-center XC correction needs a dedicated GauXC design.", &
924 nonlocal_vdw_abort_message = &
925 "GauXC does not support non-local VDW_POTENTIAL corrections. "// &
926 "Use an additive PAIR_POTENTIAL dispersion correction or disable GauXC."
927 REAL(kind=dp), PARAMETER :: gapw_fd_gradient_dx = 1.0e-4_dp
928
929 CHARACTER(len=32) :: gradient_mpi_runtime_env
930 CHARACTER(len=default_path_length) :: model_key, model_name, output_path
931 CHARACTER(len=default_string_length) :: grid_key, grid_type, int_exec_space, lb_exec_space, &
932 pruning_key, pruning_scheme, radial_quadrature, skala_runtime, skala_runtime_key, &
933 xc_fun_name
934 INTEGER :: batch_size, env_status, img, ispin, &
935 natom, nimages, nspins, &
936 onedft_atom_chunk_size
937 LOGICAL :: do_kpoints, grid_explicit, hdf5_output, is_periodic, molecular_virial, &
938 molecular_virial_debug, onedft_atom_chunk_size_explicit, periodic_reference, &
939 pruning_explicit, use_fd_gradient, use_gradient_mpi_runtime, use_gradient_self_runtime, &
940 use_onedft, use_self_runtime, use_skala_model, write_hdf5_output
941 REAL(kind=dp) :: device_runtime_fill_fraction, &
942 molecular_virial_debug_dx
943 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: density_scalar, density_zeta
944 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
945 TYPE(cell_type), POINTER :: cell
946 TYPE(cp_gauxc_basisset_type) :: gauxc_basis
947 TYPE(cp_gauxc_grid_type) :: gauxc_gradient_grid_result, &
948 gauxc_grid_result
949 TYPE(cp_gauxc_integrator_type) :: gauxc_gradient_integrator_result, gauxc_integrator_result
950 TYPE(cp_gauxc_molecule_type) :: gauxc_mol
951 TYPE(cp_gauxc_status_type) :: gauxc_status
952 TYPE(cp_gauxc_xc_gradient_type) :: exc_grad
953 TYPE(cp_gauxc_xc_type) :: gauxc_xc_result
954 TYPE(dbcsr_p_type) :: vxc_zeta_tmp
955 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_vxc
956 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: rho_ao
957 TYPE(dft_control_type), POINTER :: dft_control
958 TYPE(mp_para_env_type), POINTER :: para_env
959 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
960 TYPE(qs_energy_type), POINTER :: energy
961 TYPE(qs_force_type), DIMENSION(:), POINTER :: force
962 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
963 TYPE(qs_ks_env_type), POINTER :: ks_env
964 TYPE(qs_rho_type), POINTER :: rho, rho_use, rho_xc
965 TYPE(qs_scf_env_type), POINTER :: scf_env
966 TYPE(section_vals_type), POINTER :: gauxc_functional_section
967
968 NULLIFY ( &
969 atomic_kind_set, &
970 cell, &
971 dft_control, &
972 energy, &
973 force, &
974 ks_env, &
975 matrix_vxc, &
976 para_env, &
977 particle_set, &
978 qs_kind_set, &
979 rho, &
980 rho_use, &
981 rho_xc, &
982 rho_ao, &
983 scf_env)
984
985 CALL get_qs_env( &
986 qs_env, &
987 cell=cell, &
988 dft_control=dft_control, &
989 do_kpoints=do_kpoints, &
990 energy=energy, &
991 ks_env=ks_env, &
992 matrix_vxc=matrix_vxc, &
993 natom=natom, &
994 atomic_kind_set=atomic_kind_set, &
995 force=force, &
996 para_env=para_env, &
997 particle_set=particle_set, &
998 qs_kind_set=qs_kind_set, &
999 rho=rho, &
1000 rho_xc=rho_xc, &
1001 scf_env=scf_env)
1002
1003 IF (dft_control%qs_control%gapw_xc) THEN
1004 cpabort(gapw_xc_abort_message)
1005 END IF
1006 IF (dft_control%qs_control%gapw) THEN
1007 CALL ensure_gauxc_gapw_all_electron(qs_kind_set)
1008 END IF
1009 cpassert(ASSOCIATED(rho))
1010 rho_use => rho
1011 CALL qs_rho_get( &
1012 rho_use, &
1013 rho_ao_kp=rho_ao)
1014
1015 nimages = dft_control%nimages
1016 nspins = dft_control%nspins
1017 is_periodic = .false.
1018 IF (ASSOCIATED(cell)) is_periodic = any(cell%perd /= 0)
1019
1020 IF (ASSOCIATED(qs_env%dispersion_env)) THEN
1021 IF (qs_env%dispersion_env%type == xc_vdw_fun_nonloc) THEN
1022 cpabort(nonlocal_vdw_abort_message)
1023 END IF
1024 END IF
1025 NULLIFY (vxc_zeta_tmp%matrix)
1026
1027 gauxc_functional_section => get_gauxc_functional(xc_section)
1028 CALL section_vals_val_get( &
1029 gauxc_functional_section, &
1030 "FUNCTIONAL", &
1031 c_val=xc_fun_name)
1032 CALL section_vals_val_get( &
1033 gauxc_functional_section, &
1034 "MODEL", &
1035 c_val=model_name)
1036 CALL section_vals_val_get( &
1037 gauxc_functional_section, &
1038 "GRID", &
1039 c_val=grid_type, &
1040 explicit=grid_explicit)
1041 CALL section_vals_val_get( &
1042 gauxc_functional_section, &
1043 "RADIAL_QUADRATURE", &
1044 c_val=radial_quadrature)
1045 CALL section_vals_val_get( &
1046 gauxc_functional_section, &
1047 "PRUNING_SCHEME", &
1048 c_val=pruning_scheme, &
1049 explicit=pruning_explicit)
1050 CALL section_vals_val_get( &
1051 gauxc_functional_section, &
1052 "BATCH_SIZE", &
1053 i_val=batch_size)
1054 CALL section_vals_val_get( &
1055 gauxc_functional_section, &
1056 "DEVICE_RUNTIME_FILL_FRACTION", &
1057 r_val=device_runtime_fill_fraction)
1058 CALL section_vals_val_get( &
1059 gauxc_functional_section, &
1060 "ONEDFT_ATOM_CHUNK_SIZE", &
1061 i_val=onedft_atom_chunk_size, &
1062 explicit=onedft_atom_chunk_size_explicit)
1063 CALL section_vals_val_get( &
1064 gauxc_functional_section, &
1065 "PERIODIC_REFERENCE", &
1066 l_val=periodic_reference)
1067 CALL section_vals_val_get( &
1068 gauxc_functional_section, &
1069 "MOLECULAR_VIRIAL", &
1070 l_val=molecular_virial)
1071 CALL section_vals_val_get( &
1072 gauxc_functional_section, &
1073 "MOLECULAR_VIRIAL_DEBUG", &
1074 l_val=molecular_virial_debug)
1075 CALL section_vals_val_get( &
1076 gauxc_functional_section, &
1077 "MOLECULAR_VIRIAL_DEBUG_DX", &
1078 r_val=molecular_virial_debug_dx)
1079 CALL section_vals_val_get( &
1080 gauxc_functional_section, &
1081 "LB_EXECUTION_SPACE", &
1082 c_val=lb_exec_space)
1083 CALL section_vals_val_get( &
1084 gauxc_functional_section, &
1085 "INT_EXECUTION_SPACE", &
1086 c_val=int_exec_space)
1087 CALL section_vals_val_get( &
1088 gauxc_functional_section, &
1089 "SKALA_RUNTIME", &
1090 c_val=skala_runtime)
1091 CALL section_vals_val_get( &
1092 gauxc_functional_section, &
1093 "OUTPUT_PATH", &
1094 c_val=output_path)
1095
1096 model_key = adjustl(model_name)
1097 CALL uppercase(model_key)
1098 skala_runtime_key = adjustl(skala_runtime)
1099 CALL uppercase(skala_runtime_key)
1100 use_onedft = (trim(model_key) /= "" .AND. trim(model_key) /= "NONE")
1101 use_skala_model = (index(trim(model_key), "SKALA") > 0)
1102 IF (device_runtime_fill_fraction <= 0.0_dp .OR. device_runtime_fill_fraction > 1.0_dp) THEN
1103 CALL cp_abort(__location__, &
1104 "GAUXC%DEVICE_RUNTIME_FILL_FRACTION must be > 0 and <= 1.")
1105 END IF
1106 IF (onedft_atom_chunk_size < -1) THEN
1107 CALL cp_abort(__location__, &
1108 "GAUXC%ONEDFT_ATOM_CHUNK_SIZE must be -1, zero, or positive.")
1109 END IF
1110 IF (molecular_virial_debug) THEN
1111 IF (molecular_virial_debug_dx <= 0.0_dp) THEN
1112 CALL cp_abort(__location__, &
1113 "GauXC MOLECULAR_VIRIAL_DEBUG_DX must be positive.")
1114 END IF
1115 molecular_virial = .true.
1116 END IF
1117 CALL ensure_gauxc_periodic_reference_scope( &
1118 dft_control, cell, qs_kind_set, do_kpoints, periodic_reference)
1119 IF (is_periodic .AND. periodic_reference .AND. para_env%mepos == 0) THEN
1120 IF (ASSOCIATED(scf_env)) THEN
1121 IF (scf_env%iter_count == 1) THEN
1122 CALL cp_warn( &
1123 __location__, &
1124 "GAUXC%PERIODIC_REFERENCE uses GauXC molecular quadrature for isolated validation "// &
1125 "cells. Compact periodic materials require a dedicated periodic GauXC interface.")
1126 END IF
1127 END IF
1128 END IF
1129 IF (dft_control%qs_control%gapw_xc .AND. use_onedft) THEN
1130 CALL cp_abort(__location__, &
1131 "GauXC OneDFT/SKALA with METHOD GAPW_XC is not implemented. "// &
1132 "The GAPW_XC one-center XC correction must be evaluated by the same GauXC model.")
1133 END IF
1134 IF (use_onedft) THEN
1135 IF (has_nlcc(qs_kind_set)) THEN
1136 CALL cp_abort(__location__, &
1137 "GauXC OneDFT/SKALA with NLCC pseudopotentials is not implemented. "// &
1138 "The frozen core density would need a SKALA-consistent feature definition.")
1139 END IF
1140 END IF
1141 IF (use_onedft) THEN
1142 CALL set_gauxc_onedft_atom_chunk_env( &
1143 onedft_atom_chunk_size, onedft_atom_chunk_size_explicit)
1144 IF (.NOT. grid_explicit) grid_type = "SUPERFINE"
1145 IF (.NOT. pruning_explicit) pruning_scheme = "UNPRUNED"
1146
1147 grid_key = adjustl(grid_type)
1148 pruning_key = adjustl(pruning_scheme)
1149 CALL uppercase(grid_key)
1150 CALL uppercase(pruning_key)
1151 IF (use_skala_model .AND. (calculate_forces .OR. molecular_virial) .AND. &
1152 (trim(grid_key) /= "SUPERFINE" .OR. trim(pruning_key) /= "UNPRUNED")) THEN
1153 CALL cp_warn( &
1154 __location__, &
1155 "GauXC OneDFT/SKALA nuclear gradients are sensitive to the GauXC molecular grid. "// &
1156 "Use GRID SUPERFINE and PRUNING_SCHEME UNPRUNED for quantitative force checks.")
1157 END IF
1158 IF (trim(model_key) == "SKALA") THEN
1159 CALL get_environment_variable("GAUXC_SKALA_MODEL", model_name, status=env_status)
1160 IF (env_status /= 0 .OR. len_trim(model_name) == 0) THEN
1161 cpabort("MODEL SKALA requires the GAUXC_SKALA_MODEL environment variable")
1162 END IF
1163 END IF
1164 END IF
1165 SELECT CASE (trim(skala_runtime_key))
1166 CASE ("AUTO")
1167 use_self_runtime = use_skala_model .AND. para_env%num_pe > 1 .AND. nspins > 1
1168 CASE ("MPI")
1169 use_self_runtime = .false.
1170 CASE ("SELF")
1171 use_self_runtime = use_skala_model .AND. para_env%num_pe > 1
1172 CASE DEFAULT
1173 CALL cp_abort(__location__, "Unknown GAUXC%SKALA_RUNTIME value.")
1174 END SELECT
1175 IF (.NOT. use_skala_model) use_self_runtime = .false.
1176 use_gradient_self_runtime = calculate_forces .AND. use_onedft .AND. &
1177 para_env%num_pe > 1 .AND. .NOT. use_self_runtime
1178 IF (use_skala_model .AND. para_env%num_pe > 1 .AND. .NOT. use_self_runtime .AND. &
1179 para_env%mepos == 0 .AND. ASSOCIATED(scf_env)) THEN
1180 IF (scf_env%iter_count == 1) THEN
1181 CALL cp_warn( &
1182 __location__, &
1183 "GAUXC%SKALA_RUNTIME uses the MPI communicator for energy/VXC. "// &
1184 "SKALA Torch atom chunks can be distributed across MPI ranks; "// &
1185 "set GAUXC_ONEDFT_DISTRIBUTED_TORCH=0 to force rank-0 Torch inference.")
1186 END IF
1187 END IF
1188 use_gradient_mpi_runtime = .false.
1189 CALL get_environment_variable("GAUXC_GRADIENT_MPI_RUNTIME", gradient_mpi_runtime_env, &
1190 status=env_status)
1191 IF (env_status == 0) THEN
1192 CALL uppercase(gradient_mpi_runtime_env)
1193 SELECT CASE (trim(gradient_mpi_runtime_env))
1194 CASE ("", "0", "FALSE", "F", "OFF", "AUTO")
1195 use_gradient_mpi_runtime = .false.
1196 CASE ("1", "TRUE", "T", "ON", "MPI")
1197 use_gradient_mpi_runtime = .true.
1198 CASE DEFAULT
1199 CALL cp_abort(__location__, &
1200 "GAUXC_GRADIENT_MPI_RUNTIME must be 0/1, TRUE/FALSE, ON/OFF, or MPI.")
1201 END SELECT
1202 END IF
1203
1204 gauxc_mol = gauxc_create_molecule( &
1205 particle_set, &
1206 gauxc_status)
1207 CALL gauxc_check_status(gauxc_status)
1208 gauxc_basis = gauxc_create_basisset( &
1209 qs_kind_set, &
1210 particle_set, &
1211 gauxc_status)
1212 CALL gauxc_check_status(gauxc_status)
1213 hdf5_output = (trim(output_path) /= "")
1214 write_hdf5_output = hdf5_output .AND. para_env%mepos == 0
1215 IF (write_hdf5_output .AND. ASSOCIATED(scf_env)) THEN
1216 write_hdf5_output = scf_env%iter_count == 1
1217 END IF
1218 IF (write_hdf5_output) THEN
1219 CALL gauxc_write_molecule_hdf5( &
1220 gauxc_mol, &
1221 output_path, &
1222 "molecule.h5", &
1223 "molecule", &
1224 gauxc_status)
1225 CALL gauxc_check_status(gauxc_status)
1226 CALL gauxc_write_basisset_hdf5( &
1227 gauxc_basis, &
1228 output_path, &
1229 "basisset.h5", &
1230 "basisset", &
1231 gauxc_status)
1232 CALL gauxc_check_status(gauxc_status)
1233 END IF
1234 use_fd_gradient = dft_control%qs_control%gapw .AND. calculate_forces .AND. &
1235 (use_skala_model .OR. gauxc_basis%max_l > 3)
1236 IF (use_gradient_mpi_runtime) THEN
1237 use_gradient_self_runtime = .false.
1238 END IF
1239 IF (use_fd_gradient .AND. para_env%mepos == 0) THEN
1240 CALL cp_warn( &
1241 __location__, &
1242 "Using finite-difference GauXC XC gradients for METHOD GAPW with SKALA or "// &
1243 "all-electron basis functions beyond f shells. The upstream analytical GauXC "// &
1244 "gradient path is not yet reliable for this case.")
1245 END IF
1246 IF (use_self_runtime) THEN
1247 ! SKALA currently needs a replicated molecular runtime for reproducible
1248 ! open-shell densities across CP2K MPI ranks.
1249 gauxc_grid_result = gauxc_create_grid( &
1250 gauxc_mol, &
1251 gauxc_basis, &
1252 grid_type, &
1253 radial_quadrature, &
1254 pruning_scheme, &
1255 lb_exec_space, &
1256 batch_size, &
1257 device_runtime_fill_fraction, &
1258 gauxc_status, &
1259 mpi_comm=mp_comm_self%get_handle())
1260 ELSE
1261 ! Use the QS force-evaluation communicator rather than GauXC's global
1262 ! runtime. In mixed CDFT the diabatic states can be built on disjoint
1263 ! MPI subgroups; using the global communicator would make GauXC wait
1264 ! for ranks that are working on another state.
1265 gauxc_grid_result = gauxc_create_grid( &
1266 gauxc_mol, &
1267 gauxc_basis, &
1268 grid_type, &
1269 radial_quadrature, &
1270 pruning_scheme, &
1271 lb_exec_space, &
1272 batch_size, &
1273 device_runtime_fill_fraction, &
1274 gauxc_status, &
1275 mpi_comm=para_env%get_handle())
1276 END IF
1277 CALL gauxc_check_status(gauxc_status)
1278 gauxc_integrator_result = gauxc_create_integrator( &
1279 trim(xc_fun_name), &
1280 gauxc_grid_result, &
1281 int_exec_space, &
1282 nspins, &
1283 gauxc_status)
1284 CALL gauxc_check_status(gauxc_status)
1285
1286 IF (use_gradient_self_runtime) THEN
1287 ! Upstream GauXC does not yet support OneDFT/SKALA nuclear gradients
1288 ! on an MPI runtime. Keep the energy/VXC path on the normal MPI
1289 ! runtime and use an isolated runtime only for the replicated gradient.
1290 gauxc_gradient_grid_result = gauxc_create_grid( &
1291 gauxc_mol, &
1292 gauxc_basis, &
1293 grid_type, &
1294 radial_quadrature, &
1295 pruning_scheme, &
1296 lb_exec_space, &
1297 batch_size, &
1298 device_runtime_fill_fraction, &
1299 gauxc_status, &
1300 mpi_comm=mp_comm_self%get_handle())
1301 CALL gauxc_check_status(gauxc_status)
1302 gauxc_gradient_integrator_result = gauxc_create_integrator( &
1303 trim(xc_fun_name), &
1304 gauxc_gradient_grid_result, &
1305 int_exec_space, &
1306 nspins, &
1307 gauxc_status)
1308 CALL gauxc_check_status(gauxc_status)
1309 END IF
1310
1311 IF (qs_env%run_rtp) THEN
1312 cpabort("GAUXC XC energy currently does not support real-time propagation")
1313 END IF
1314
1315 energy%exc = 0
1316
1317 IF (ASSOCIATED(matrix_vxc)) CALL dbcsr_deallocate_matrix_set(matrix_vxc)
1318 CALL dbcsr_allocate_matrix_set(matrix_vxc, nspins)
1319
1320 DO img = 1, nimages
1321 IF (img > 1) THEN
1322 cpabort("UNIMPLEMENTED: Handling nimg>1 in k-point integration")
1323 END IF
1324 CALL dbcsr_to_dense(rho_ao(1, img), density_scalar)
1325 CALL para_env%sum(density_scalar)
1326 IF (nspins == 1) THEN
1327 gauxc_xc_result = gauxc_compute_xc( &
1328 gauxc_integrator_result, &
1329 density_scalar, &
1330 nspins=nspins, &
1331 status=gauxc_status, &
1332 model=trim(model_name))
1333 CALL gauxc_check_status(gauxc_status)
1334 IF (calculate_forces) THEN
1335 IF (use_fd_gradient) THEN
1336 CALL gauxc_xc_gradient_fd( &
1337 particle_set, qs_kind_set, density_scalar, nspins, model_name, &
1338 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, &
1339 lb_exec_space, int_exec_space, batch_size, &
1340 device_runtime_fill_fraction, gapw_fd_gradient_dx, para_env, &
1341 exc_grad%exc_grad)
1342 ELSE IF (use_gradient_self_runtime) THEN
1343 exc_grad = gauxc_compute_xc_gradient( &
1344 gauxc_gradient_integrator_result, &
1345 density_scalar, &
1346 nspins=nspins, &
1347 natom=natom, &
1348 status=gauxc_status, &
1349 model=trim(model_name))
1350 ELSE
1351 exc_grad = gauxc_compute_xc_gradient( &
1352 gauxc_integrator_result, &
1353 density_scalar, &
1354 nspins=nspins, &
1355 natom=natom, &
1356 status=gauxc_status, &
1357 model=trim(model_name))
1358 END IF
1359 CALL gauxc_check_status(gauxc_status)
1360 CALL add_gauxc_gradient_to_force( &
1361 exc_grad%exc_grad, &
1362 force, &
1363 atomic_kind_set, &
1364 para_env)
1365 IF (molecular_virial) THEN
1366 CALL print_gauxc_molecular_virial(exc_grad%exc_grad, particle_set, para_env)
1367 END IF
1368 IF (molecular_virial_debug) THEN
1369 CALL debug_gauxc_molecular_virial( &
1370 exc_grad%exc_grad, particle_set, qs_kind_set, density_scalar, nspins, &
1371 model_name, xc_fun_name, grid_type, radial_quadrature, pruning_scheme, &
1372 lb_exec_space, int_exec_space, batch_size, device_runtime_fill_fraction, &
1373 molecular_virial_debug_dx, para_env)
1374 END IF
1375 DEALLOCATE (exc_grad%exc_grad)
1376 END IF
1377 ELSE
1378 cpassert(nspins == 2)
1379 ! In here:
1380 ! scalar <- rho_ao(1, :) + rho_ao(2, :)
1381 ! zeta <- rho_ao(1, :) - rho_ao(2, :)
1382 CALL dbcsr_to_dense(rho_ao(2, img), density_zeta)
1383 CALL para_env%sum(density_zeta)
1384 ! Do NOT reorder the following lines!
1385 density_scalar(:, :) = density_scalar(:, :) + density_zeta(:, :)
1386 ! Factor two because the next line is evaluated after the above line.
1387 ! We need to subtract density_zeta once to undo the above line and
1388 ! a second time because that is what UKS requires.
1389 ! This style lowers memory footprint.
1390 density_zeta(:, :) = density_scalar(:, :) - 2.0_dp*density_zeta(:, :)
1391 gauxc_xc_result = gauxc_compute_xc( &
1392 gauxc_integrator_result, &
1393 density_scalar, &
1394 density_zeta, &
1395 nspins, &
1396 gauxc_status, &
1397 model=trim(model_name))
1398 CALL gauxc_check_status(gauxc_status)
1399 IF (calculate_forces) THEN
1400 IF (use_fd_gradient) THEN
1401 CALL gauxc_xc_gradient_fd( &
1402 particle_set, qs_kind_set, density_scalar, nspins, model_name, &
1403 xc_fun_name, grid_type, radial_quadrature, pruning_scheme, &
1404 lb_exec_space, int_exec_space, batch_size, &
1405 device_runtime_fill_fraction, gapw_fd_gradient_dx, para_env, &
1406 exc_grad%exc_grad, density_zeta=density_zeta)
1407 ELSE IF (use_gradient_self_runtime) THEN
1408 exc_grad = gauxc_compute_xc_gradient( &
1409 gauxc_gradient_integrator_result, &
1410 density_scalar, &
1411 density_zeta, &
1412 nspins, &
1413 natom, &
1414 gauxc_status, &
1415 model=trim(model_name))
1416 ELSE
1417 exc_grad = gauxc_compute_xc_gradient( &
1418 gauxc_integrator_result, &
1419 density_scalar, &
1420 density_zeta, &
1421 nspins, &
1422 natom, &
1423 gauxc_status, &
1424 model=trim(model_name))
1425 END IF
1426 CALL gauxc_check_status(gauxc_status)
1427 CALL add_gauxc_gradient_to_force( &
1428 exc_grad%exc_grad, &
1429 force, &
1430 atomic_kind_set, &
1431 para_env)
1432 IF (molecular_virial) THEN
1433 CALL print_gauxc_molecular_virial(exc_grad%exc_grad, particle_set, para_env)
1434 END IF
1435 IF (molecular_virial_debug) THEN
1436 CALL debug_gauxc_molecular_virial( &
1437 exc_grad%exc_grad, particle_set, qs_kind_set, density_scalar, nspins, &
1438 model_name, xc_fun_name, grid_type, radial_quadrature, pruning_scheme, &
1439 lb_exec_space, int_exec_space, batch_size, device_runtime_fill_fraction, &
1440 molecular_virial_debug_dx, para_env, density_zeta=density_zeta)
1441 END IF
1442 DEALLOCATE (exc_grad%exc_grad)
1443 END IF
1444 END IF
1445
1446 energy%exc = energy%exc + gauxc_xc_result%exc
1447
1448 IF (nspins == 1) THEN
1449 IF (img == 1) THEN
1450 matrix_vxc(1) = dense_to_dbcsr(gauxc_xc_result%vxc_scalar, rho_ao(1, img))
1451 ELSE
1452 cpabort("UNIMPLEMENTED: Handling multiple result matrices in k-point integration")
1453 END IF
1454 ELSE
1455 cpassert(nspins == 2)
1456 ! Transform derivatives from total/spin density back to alpha/beta channels.
1457 vxc_zeta_tmp = dense_to_dbcsr(gauxc_xc_result%vxc_zeta, rho_ao(1, img))
1458 IF (img == 1) THEN
1459 DO ispin = 1, 2
1460 matrix_vxc(ispin) = dense_to_dbcsr(gauxc_xc_result%vxc_scalar, rho_ao(ispin, 1))
1461 CALL dbcsr_add( &
1462 matrix_vxc(ispin)%matrix, &
1463 vxc_zeta_tmp%matrix, &
1464 1.0_dp, &
1465 ! 1.0 for ispin==1, -1.0 for ispin==2
1466 1.0_dp - real(ispin - 1, dp)*2.0_dp)
1467 END DO
1468 ELSE
1469 cpabort("UNIMPLEMENTED: Handling multiple result matrices in k-point integration")
1470 END IF
1471 CALL dbcsr_release(vxc_zeta_tmp%matrix)
1472 DEALLOCATE (vxc_zeta_tmp%matrix)
1473 END IF
1474 END DO
1475
1476 DEALLOCATE (density_scalar)
1477 IF (ALLOCATED(density_zeta)) DEALLOCATE (density_zeta)
1478 DEALLOCATE (gauxc_xc_result%vxc_scalar)
1479 IF (ALLOCATED(gauxc_xc_result%vxc_zeta)) DEALLOCATE (gauxc_xc_result%vxc_zeta)
1480
1481 CALL set_ks_env(ks_env, matrix_vxc=matrix_vxc)
1482 DO ispin = 1, nspins
1483 CALL dbcsr_finalize(matrix_vxc(ispin)%matrix)
1484 END DO
1485
1486 IF (use_gradient_self_runtime) THEN
1487 CALL gauxc_destroy_integrator(gauxc_gradient_integrator_result, gauxc_status)
1488 CALL gauxc_check_status(gauxc_status)
1489 CALL gauxc_destroy_grid(gauxc_gradient_grid_result, gauxc_status)
1490 CALL gauxc_check_status(gauxc_status)
1491 END IF
1492 CALL gauxc_destroy_integrator(gauxc_integrator_result, gauxc_status)
1493 CALL gauxc_check_status(gauxc_status)
1494 CALL gauxc_destroy_grid(gauxc_grid_result, gauxc_status)
1495 CALL gauxc_check_status(gauxc_status)
1496 CALL gauxc_destroy_basisset(gauxc_basis, gauxc_status)
1497 CALL gauxc_check_status(gauxc_status)
1498 CALL gauxc_destroy_molecule(gauxc_mol, gauxc_status)
1499 CALL gauxc_check_status(gauxc_status)
1500
1501 END SUBROUTINE apply_gauxc
1502
1503END MODULE xc_gauxc_functional
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
Handles all functions related to the CELL.
Definition cell_types.F:15
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
character function, public dbcsr_get_matrix_type(matrix)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_get_stored_coordinates(matrix, row, column, processor)
...
subroutine, public dbcsr_init_p(matrix)
...
subroutine, public dbcsr_work_create(matrix, nblks_guess, sizedata_guess, n, work_mutable)
...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset, transposed)
...
subroutine, public dbcsr_finalize(matrix)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_release(matrix)
...
subroutine, public dbcsr_put_block(matrix, row, col, block, summation)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
subroutine, public dbcsr_distribution_get(dist, row_dist, col_dist, nrows, ncols, has_threads, group, mynode, numnodes, nprows, npcols, myprow, mypcol, pgrid, subgroups_defined, prow_group, pcol_group)
...
DBCSR operations in CP2K.
various routines to log and control the output. The idea is that decisions about where to log should ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
Definition of the atomic potential types.
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public xc_vdw_fun_nonloc
objects that represent the structure of input sections and the data contained in an input section
type(section_vals_type) function, pointer, public section_vals_get_subs_vals2(section_vals, i_section, i_rep_section)
returns the values of the n-th non default subsection (null if no such section exists (not so many no...
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
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
integer, parameter, public default_path_length
Definition kinds.F:58
Interface to the message passing library MPI.
type(mp_comm_type), parameter, public mp_comm_self
Define the data structure for the particle information.
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, mimic, 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, sab_cneo, 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, xcint_weights, 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, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_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, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
logical function, public has_nlcc(qs_kind_set)
finds if a given qs run needs to use nlcc
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, cneo_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, monovalent, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
subroutine, public set_ks_env(ks_env, v_hartree_rspace, s_mstruct_changed, rho_changed, exc_accint, 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, xcint_weights, 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, sab_cneo, task_list, task_list_soft, subsys, dft_control, dbcsr_dist, distribution_2d, pw_env, para_env, blacs_env)
...
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...
module that contains the definitions of the scf types
Utilities for string manipulations.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
logical function, public xc_section_uses_gauxc(xc_section)
...
subroutine, public skala_info(functional, lsd, reference, shortform, needs, max_deriv)
Return information about the Skala functional.
subroutine, public apply_gauxc(qs_env, xc_section, calculate_forces)
...
contains the structure
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:60
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.
contains a flag for each component of xc_rho_set, so that you can use it to tell which components you...