(git:41cb813)
Loading...
Searching...
No Matches
sirius_interface.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
8!***************************************************************************************************
9!> \brief Interface to the SIRIUS Library
10!> \par History
11!> 07.2018 initial create
12!> \author JHU
13!***************************************************************************************************
14#if defined(__SIRIUS)
16 USE iso_c_binding, ONLY: c_double,&
17 c_int,&
18 c_loc
22 USE atom_upf, ONLY: atom_upfpot_type
26 USE cell_types, ONLY: cell_type,&
31 USE cp_output_handling, ONLY: cp_p_file,&
45 USE kinds, ONLY: default_string_length,&
46 dp
47 USE machine, ONLY: m_flush
48 USE mathconstants, ONLY: fourpi,&
49 gamma1
52 USE physcon, ONLY: massunit
61 USE qs_kind_types, ONLY: get_qs_kind,&
65 USE sirius, ONLY: &
66 sirius_integer_array_type, sirius_integer_type, sirius_logical_array_type, &
67 sirius_logical_type, sirius_number_array_type, sirius_number_type, &
68 sirius_string_array_type, sirius_string_type, sirius_add_atom, sirius_add_atom_type, &
69 sirius_add_atom_type_radial_function, sirius_add_xc_functional, sirius_context_handler, &
70 sirius_create_context, sirius_create_ground_state, sirius_create_kset_from_grid, &
71 sirius_finalize, sirius_find_ground_state, sirius_get_band_energies, &
72 sirius_get_band_occupancies, sirius_get_energy, sirius_get_forces, &
73 sirius_get_kpoint_properties, sirius_get_num_kpoints, sirius_get_parameters, &
74 sirius_get_stress_tensor, sirius_ground_state_handler, sirius_import_parameters, &
75 sirius_initialize, sirius_initialize_context, sirius_is_initialized, &
76 sirius_kpoint_set_handler, sirius_option_get_info, sirius_option_get_section_length, &
77 sirius_option_set, sirius_set_atom_position, sirius_set_atom_type_dion, &
78 sirius_set_atom_type_hubbard, sirius_set_atom_type_radial_grid, &
79 sirius_set_lattice_vectors, sirius_set_mpi_grid_dims, sirius_update_ground_state
80#include "./base/base_uses.f90"
81
82 IMPLICIT NONE
83
84 PRIVATE
85
86 ! Global parameters
87
88 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'sirius_interface'
89
90 ! Public subroutines
91
92 PUBLIC :: cp_sirius_create_env, &
98
99CONTAINS
100
101!***************************************************************************************************
102!> \brief Initialize the Sirius library
103!> \par Creation (07.2018, JHU)
104!> \author JHU
105! **************************************************************************************************
106 SUBROUTINE cp_sirius_init()
107 CALL sirius_initialize(.false.)
108 END SUBROUTINE cp_sirius_init
109
110!***************************************************************************************************
111!> \brief Check the initialisation status of the Sirius library
112!> \return Return the initialisation status of the Sirius library as boolean
113!> \par Creation (03.12.2025, MK)
114!> \author Matthias Krack
115! **************************************************************************************************
116 LOGICAL FUNCTION cp_sirius_is_initialized()
117 CALL sirius_is_initialized(cp_sirius_is_initialized)
118 END FUNCTION cp_sirius_is_initialized
119
120!***************************************************************************************************
121!> \brief Finalize the Sirius library
122!> \par Creation (07.2018, JHU)
123!> \author JHU
124! **************************************************************************************************
125 SUBROUTINE cp_sirius_finalize()
126 CALL sirius_finalize(.false., .false., .false.)
127 END SUBROUTINE cp_sirius_finalize
128
129!***************************************************************************************************
130!> \brief ...
131!> \param pwdft_env ...
132!> \param
133!> \par History
134!> 07.2018 Create the Sirius environment
135!> \author JHU
136! **************************************************************************************************
137 SUBROUTINE cp_sirius_create_env(pwdft_env)
138 TYPE(pwdft_environment_type), POINTER :: pwdft_env
139#if defined(__SIRIUS)
140
141 CHARACTER(len=2) :: element_symbol
142 CHARACTER(len=default_string_length) :: label
143 INTEGER :: i, ii, jj, iatom, ibeta, ifun, ikind, iwf, j, l, &
144 n, ns, natom, nbeta, nbs, nkind, nmesh, &
145 num_mag_dims, sirius_mpi_comm, vdw_func, nu, lu, output_unit
146 INTEGER, DIMENSION(:), POINTER :: mpi_grid_dims
147 INTEGER(KIND=C_INT), DIMENSION(3) :: k_grid, k_shift
148 INTEGER, DIMENSION(:), POINTER :: kk
149 LOGICAL :: up, use_ref_cell
150 LOGICAL(4) :: use_so, use_symmetry, dft_plus_u_atom
151 REAL(KIND=c_double), ALLOCATABLE, DIMENSION(:) :: fun
152 REAL(KIND=c_double), ALLOCATABLE, DIMENSION(:, :) :: dion
153 REAL(KIND=c_double), DIMENSION(3) :: a1, a2, a3, v1, v2
154 REAL(KIND=dp) :: al, angle1, angle2, cval, focc, &
155 magnetization, mass, pf, rl, zeff, alpha_u, beta_u, &
156 j0_u, j_u, u_u, occ_u, u_minus_j, vnlp, vnlm
157 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: beta, corden, ef, fe, locpot, rc, rp
158 REAL(KIND=dp), DIMENSION(3) :: vr, vs, j_t
159 REAL(KIND=dp), DIMENSION(:), POINTER :: density
160 REAL(KIND=dp), DIMENSION(:, :), POINTER :: wavefunction, wfninfo
161 TYPE(atom_gthpot_type), POINTER :: gth_atompot
162 TYPE(atom_upfpot_type), POINTER :: upf_pot
163 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
164 TYPE(atomic_kind_type), POINTER :: atomic_kind
165 TYPE(cell_type), POINTER :: my_cell
166 TYPE(mp_para_env_type), POINTER :: para_env
167 TYPE(grid_atom_type), POINTER :: atom_grid
168 TYPE(gth_potential_type), POINTER :: gth_potential
169 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
170 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
171 TYPE(qs_subsys_type), POINTER :: qs_subsys
172 TYPE(section_vals_type), POINTER :: pwdft_section, pwdft_sub_section, &
173 xc_fun, xc_section
174 TYPE(sirius_context_handler) :: sctx
175 TYPE(sirius_ground_state_handler) :: gs_handler
176 TYPE(sirius_kpoint_set_handler) :: ks_handler
177
178 cpassert(ASSOCIATED(pwdft_env))
179
180 output_unit = cp_logger_get_default_io_unit()
181 ! create context of simulation
182 CALL pwdft_env_get(pwdft_env, para_env=para_env)
183 sirius_mpi_comm = para_env%get_handle()
184 CALL sirius_create_context(sirius_mpi_comm, sctx)
185
186! the "fun" starts.
187
188 CALL pwdft_env_get(pwdft_env=pwdft_env, pwdft_input=pwdft_section, xc_input=xc_section)
189
190 CALL section_vals_val_get(pwdft_section, "ignore_convergence_failure", &
191 l_val=pwdft_env%ignore_convergence_failure)
192 ! cp2k should *have* a function that return all xc_functionals. Doing
193 ! manually is prone to errors
194
195 IF (ASSOCIATED(xc_section)) THEN
196 ifun = 0
197 DO
198 ifun = ifun + 1
199 xc_fun => section_vals_get_subs_vals2(xc_section, i_section=ifun)
200 IF (.NOT. ASSOCIATED(xc_fun)) EXIT
201 ! Here, we do not have to check whether the functional name starts with XC_
202 ! because we only allow the shorter form w/o XC_
203 CALL sirius_add_xc_functional(sctx, "XC_"//trim(xc_fun%section%name))
204 END DO
205 END IF
206
207 ! import control section
208 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "control")
209 IF (ASSOCIATED(pwdft_sub_section)) THEN
210 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "control")
211 CALL section_vals_val_get(pwdft_sub_section, "mpi_grid_dims", i_vals=mpi_grid_dims)
212 END IF
213
214! import parameters section
215 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "parameters")
216
217 IF (ASSOCIATED(pwdft_sub_section)) THEN
218 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "parameters")
219 CALL section_vals_val_get(pwdft_sub_section, "ngridk", i_vals=kk)
220 k_grid(1) = kk(1)
221 k_grid(2) = kk(2)
222 k_grid(3) = kk(3)
223
224 CALL section_vals_val_get(pwdft_sub_section, "shiftk", i_vals=kk)
225 k_shift(1) = kk(1)
226 k_shift(2) = kk(2)
227 k_shift(3) = kk(3)
228 CALL section_vals_val_get(pwdft_sub_section, "num_mag_dims", i_val=num_mag_dims)
229 CALL section_vals_val_get(pwdft_sub_section, "use_symmetry", l_val=use_symmetry)
230 CALL section_vals_val_get(pwdft_sub_section, "so_correction", l_val=use_so)
231
232! now check if van der walls corrections are needed
233 vdw_func = -1
234#ifdef __LIBVDWXC
235 CALL section_vals_val_get(pwdft_sub_section, "vdw_functional", i_val=vdw_func)
236 SELECT CASE (vdw_func)
237 CASE (sirius_func_vdwdf)
238 CALL sirius_add_xc_functional(sctx, "XC_FUNC_VDWDF")
239 CASE (sirius_func_vdwdf2)
240 CALL sirius_add_xc_functional(sctx, "XC_FUNC_VDWDF2")
242 CALL sirius_add_xc_functional(sctx, "XC_FUNC_VDWDFCX")
243 CASE default
244 END SELECT
245#endif
246
247 END IF
248
249! import mixer section
250 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "mixer")
251 IF (ASSOCIATED(pwdft_sub_section)) THEN
252 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "mixer")
253 END IF
254
255! import settings section
256 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "settings")
257 IF (ASSOCIATED(pwdft_sub_section)) THEN
258 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "settings")
259 END IF
260
261 ! import solver section
262 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "iterative_solver")
263 IF (ASSOCIATED(pwdft_sub_section)) THEN
264 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "iterative_solver")
265 END IF
266
267#if defined(__SIRIUS_DFTD3)
268 ! import dftd3 section
269 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "dftd3")
270 IF (ASSOCIATED(pwdft_sub_section)) THEN
271 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "dftd3")
272 END IF
273#endif
274#if defined(__SIRIUS_DFTD4)
275 ! import dftd4 section
276 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "dftd4")
277 IF (ASSOCIATED(pwdft_sub_section)) THEN
278 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "dftd4")
279 END IF
280#endif
281
282#if defined(__SIRIUS_NLCG)
283 ! import nlcg section
284 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "nlcg")
285 IF (ASSOCIATED(pwdft_sub_section)) THEN
286 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "nlcg")
287 END IF
288#endif
289
290#if defined(__SIRIUS_VCSQNM)
291 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "vcsqnm")
292 IF (ASSOCIATED(pwdft_sub_section)) THEN
293 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "vcsqnm")
294 END IF
295#endif
296
297 !CALL sirius_dump_runtime_setup(sctx, "runtime.json")
298 CALL sirius_import_parameters(sctx, '{}')
299
300! lattice vectors of the unit cell should be in [a.u.] (length is in [a.u.])
301 CALL pwdft_env_get(pwdft_env=pwdft_env, qs_subsys=qs_subsys)
302 CALL qs_subsys_get(qs_subsys, cell=my_cell, use_ref_cell=use_ref_cell)
303 a1(:) = my_cell%hmat(:, 1)
304 a2(:) = my_cell%hmat(:, 2)
305 a3(:) = my_cell%hmat(:, 3)
306 CALL sirius_set_lattice_vectors(sctx, a1(1), a2(1), a3(1))
307
308 IF (use_ref_cell) THEN
309 cpwarn("SIRIUS| The specified CELL_REF will be ignored for PW_DFT calculations")
310 END IF
311
312! set up the atomic type definitions
313 CALL qs_subsys_get(qs_subsys, &
314 atomic_kind_set=atomic_kind_set, &
315 qs_kind_set=qs_kind_set, &
316 particle_set=particle_set)
317 nkind = SIZE(atomic_kind_set)
318 DO ikind = 1, nkind
319 CALL get_atomic_kind(atomic_kind_set(ikind), &
320 name=label, element_symbol=element_symbol, mass=mass)
321 CALL get_qs_kind(qs_kind_set(ikind), zeff=zeff)
322 NULLIFY (upf_pot, gth_potential)
323 CALL get_qs_kind(qs_kind_set(ikind), upf_potential=upf_pot, gth_potential=gth_potential)
324
325 IF (ASSOCIATED(upf_pot)) THEN
326 CALL sirius_add_atom_type(sctx, label, fname=upf_pot%filename, &
327 symbol=element_symbol, &
328 mass=real(mass/massunit, kind=c_double))
329
330 ELSEIF (ASSOCIATED(gth_potential)) THEN
331!
332 NULLIFY (atom_grid)
333 CALL allocate_grid_atom(atom_grid)
334 nmesh = 929
335 atom_grid%nr = nmesh
336 CALL create_grid_atom(atom_grid, nmesh, 1, 1, 0, do_gapw_log)
337 ALLOCATE (rp(nmesh), fun(nmesh))
338 IF (atom_grid%rad(1) < atom_grid%rad(nmesh)) THEN
339 up = .true.
340 ELSE
341 up = .false.
342 END IF
343 IF (up) THEN
344 rp(1:nmesh) = atom_grid%rad(1:nmesh)
345 ELSE
346 DO i = 1, nmesh
347 rp(i) = atom_grid%rad(nmesh - i + 1)
348 END DO
349 END IF
350! add new atom type
351 CALL sirius_add_atom_type(sctx, label, &
352 zn=nint(zeff + 0.001d0), &
353 symbol=element_symbol, &
354 mass=real(mass/massunit, kind=c_double), &
355 spin_orbit=gth_potential%soc)
356!
357 ALLOCATE (gth_atompot)
358 CALL gth_potential_conversion(gth_potential, gth_atompot)
359! set radial grid
360 fun(1:nmesh) = rp(1:nmesh)
361 CALL sirius_set_atom_type_radial_grid(sctx, label, nmesh, fun(1))
362! set beta-projectors
363! GTH SOC uses the same projectors, SIRIUS can use the same or different projectors for l+1/2, l-1/2 (l > 0 l+1/2 l < 0 l-/2 )
364 ALLOCATE (ef(nmesh), beta(nmesh))
365 ibeta = 0
366 DO l = 0, 3
367 IF (gth_atompot%nl(l) == 0) cycle
368 rl = gth_atompot%rcnl(l)
369! we need to multiply by r so that data transferred to sirius are r \beta(r) not beta(r)
370 ef(1:nmesh) = exp(-0.5_dp*rp(1:nmesh)*rp(1:nmesh)/(rl*rl))
371 DO i = 1, gth_atompot%nl(l)
372 pf = rl**(l + 0.5_dp*(4._dp*i - 1._dp))
373 j = l + 2*i - 1
374 pf = sqrt(2._dp)/(pf*sqrt(gamma1(j)))
375 beta(:) = pf*rp**(l + 2*i - 2)*ef
376 ibeta = ibeta + 1
377 fun(1:nmesh) = beta(1:nmesh)*rp(1:nmesh)
378 CALL sirius_add_atom_type_radial_function(sctx, label, &
379 "beta", fun(1), nmesh, l=l)
380 ! we double the number of beta projectors for SO and l>0
381 IF (gth_atompot%soc .AND. l /= 0) THEN
382 CALL sirius_add_atom_type_radial_function(sctx, label, &
383 "beta", fun(1), nmesh, l=-l)
384 END IF
385 END DO
386 END DO
387 DEALLOCATE (ef, beta)
388 nbeta = ibeta
389
390! nonlocal PP matrix elements
391 IF (gth_atompot%soc) THEN
392 nbs = 2*nbeta - gth_atompot%nl(0)
393 ALLOCATE (dion(nbs, nbs))
394 ELSE
395 ALLOCATE (dion(nbeta, nbeta))
396 END IF
397 dion = 0.0_dp
398 IF (gth_atompot%soc) THEN
399 ns = gth_atompot%nl(0)
400 IF (ns /= 0) THEN
401 dion(1:ns, 1:ns) = gth_atompot%hnl(1:ns, 1:ns, 0)
402 END IF
403 DO l = 1, 3
404 IF (gth_atompot%nl(l) == 0) cycle
405 DO i = 1, gth_atompot%nl(l)
406 ii = ns + 2*sum(gth_atompot%nl(1:l - 1))
407 ii = ii + 2*(i - 1) + 1
408 DO j = 1, gth_atompot%nl(l)
409 jj = ns + 2*sum(gth_atompot%nl(1:l - 1))
410 jj = jj + 2*(j - 1) + 1
411 vnlp = gth_atompot%hnl(i, j, l) + 0.5_dp*l*gth_atompot%knl(i, j, l)
412 vnlm = gth_atompot%hnl(i, j, l) - 0.5_dp*(l + 1)*gth_atompot%knl(i, j, l)
413 dion(ii, jj) = vnlp
414 dion(ii + 1, jj + 1) = vnlm
415 END DO
416 END DO
417 END DO
418 CALL sirius_set_atom_type_dion(sctx, label, nbs, dion(1, 1))
419 ELSE
420 DO l = 0, 3
421 IF (gth_atompot%nl(l) == 0) cycle
422 ibeta = sum(gth_atompot%nl(0:l - 1)) + 1
423 i = ibeta + gth_atompot%nl(l) - 1
424 dion(ibeta:i, ibeta:i) = gth_atompot%hnl(1:gth_atompot%nl(l), 1:gth_atompot%nl(l), l)
425 END DO
426 CALL sirius_set_atom_type_dion(sctx, label, nbeta, dion(1, 1))
427 END IF
428
429 DEALLOCATE (dion)
430
431! set non-linear core correction
432 IF (gth_atompot%nlcc) THEN
433 ALLOCATE (corden(nmesh), fe(nmesh), rc(nmesh))
434 corden(:) = 0.0_dp
435 n = gth_atompot%nexp_nlcc
436 DO i = 1, n
437 al = gth_atompot%alpha_nlcc(i)
438 rc(:) = rp(:)/al
439 fe(:) = exp(-0.5_dp*rc(:)*rc(:))
440 DO j = 1, gth_atompot%nct_nlcc(i)
441 cval = gth_atompot%cval_nlcc(j, i)
442 corden(:) = corden(:) + fe(:)*rc(:)**(2*j - 2)*cval
443 END DO
444 END DO
445 fun(1:nmesh) = corden(1:nmesh)*rp(1:nmesh)
446 CALL sirius_add_atom_type_radial_function(sctx, label, "ps_rho_core", &
447 fun(1), nmesh)
448 DEALLOCATE (corden, fe, rc)
449 END IF
450
451! local potential
452 ALLOCATE (locpot(nmesh))
453 locpot(:) = 0.0_dp
454 CALL atom_local_potential(locpot, gth_atompot, rp)
455 fun(1:nmesh) = locpot(1:nmesh)
456 CALL sirius_add_atom_type_radial_function(sctx, label, "vloc", &
457 fun(1), nmesh)
458 DEALLOCATE (locpot)
459!
460 NULLIFY (density, wavefunction, wfninfo)
461 CALL calculate_atomic_orbitals(atomic_kind_set(ikind), qs_kind_set(ikind), &
462 density=density, wavefunction=wavefunction, &
463 wfninfo=wfninfo, agrid=atom_grid)
464
465! set the atomic radial functions
466 DO iwf = 1, SIZE(wavefunction, 2)
467 focc = wfninfo(1, iwf)
468 l = nint(wfninfo(2, iwf))
469! we can not easily get the principal quantum number
470 nu = -1
471 IF (up) THEN
472 fun(1:nmesh) = wavefunction(1:nmesh, iwf)*rp(i)
473 ELSE
474 DO i = 1, nmesh
475 fun(i) = wavefunction(nmesh - i + 1, iwf)*rp(i)
476 END DO
477 END IF
478 CALL sirius_add_atom_type_radial_function(sctx, &
479 label, "ps_atomic_wf", &
480 fun(1), nmesh, l=l, occ=real(focc, kind=c_double), n=nu)
481 END DO
482
483! set total charge density of a free atom (to compute initial rho(r))
484 IF (up) THEN
485 fun(1:nmesh) = fourpi*density(1:nmesh)*atom_grid%rad(1:nmesh)**2
486 ELSE
487 DO i = 1, nmesh
488 fun(i) = fourpi*density(nmesh - i + 1)*atom_grid%rad(nmesh - i + 1)**2
489 END DO
490 END IF
491 CALL sirius_add_atom_type_radial_function(sctx, label, "ps_rho_total", &
492 fun(1), nmesh)
493
494 IF (ASSOCIATED(density)) DEALLOCATE (density)
495 IF (ASSOCIATED(wavefunction)) DEALLOCATE (wavefunction)
496 IF (ASSOCIATED(wfninfo)) DEALLOCATE (wfninfo)
497
498 CALL deallocate_grid_atom(atom_grid)
499 DEALLOCATE (rp, fun)
500 DEALLOCATE (gth_atompot)
501!
502 ELSE
503 CALL cp_abort(__location__, &
504 "CP2K/SIRIUS: atomic kind needs UPF or GTH potential definition")
505 END IF
506
507 CALL get_qs_kind(qs_kind_set(ikind), &
508 dft_plus_u_atom=dft_plus_u_atom, &
509 l_of_dft_plus_u=lu, &
510 n_of_dft_plus_u=nu, &
511 u_minus_j_target=u_minus_j, &
512 u_of_dft_plus_u=u_u, &
513 j_of_dft_plus_u=j_u, &
514 alpha_of_dft_plus_u=alpha_u, &
515 beta_of_dft_plus_u=beta_u, &
516 j0_of_dft_plus_u=j0_u, &
517 occupation_of_dft_plus_u=occ_u)
518
519 IF (dft_plus_u_atom) THEN
520 IF (nu < 1) THEN
521 cpabort("CP2K/SIRIUS (hubbard): principal quantum number not specified")
522 END IF
523
524 IF (lu < 0) THEN
525 cpabort("CP2K/SIRIUS (hubbard): l can not be negative.")
526 END IF
527
528 IF (occ_u < 0.0) THEN
529 cpabort("CP2K/SIRIUS (hubbard): the occupation number can not be negative.")
530 END IF
531
532 j_t(:) = 0.0
533 IF (abs(u_minus_j) < 1e-8) THEN
534 j_t(1) = j_u
535 CALL sirius_set_atom_type_hubbard(sctx, label, lu, nu, &
536 occ_u, u_u, j_t, alpha_u, beta_u, j0_u)
537 ELSE
538 CALL sirius_set_atom_type_hubbard(sctx, label, lu, nu, &
539 occ_u, u_minus_j, j_t, alpha_u, beta_u, j0_u)
540 END IF
541 END IF
542
543 END DO
544
545! add atoms to the unit cell
546! WARNING: sirius accepts only fractional coordinates;
547 natom = SIZE(particle_set)
548 DO iatom = 1, natom
549 vr(1:3) = particle_set(iatom)%r(1:3)
550 CALL real_to_scaled(vs, vr, my_cell)
551 atomic_kind => particle_set(iatom)%atomic_kind
552 ikind = atomic_kind%kind_number
553 CALL get_atomic_kind(atomic_kind, name=label)
554 CALL get_qs_kind(qs_kind_set(ikind), zeff=zeff, magnetization=magnetization)
555! angle of magnetization might come from input Atom x y z mx my mz
556! or as an angle?
557! Answer : SIRIUS only accept the magnetization as mx, my, mz
558 IF (num_mag_dims == 3) THEN
559 angle1 = 0.0_dp
560 angle2 = 0.0_dp
561 v1(1) = magnetization*sin(angle1)*cos(angle2)
562 v1(2) = magnetization*sin(angle1)*sin(angle2)
563 v1(3) = magnetization*cos(angle1)
564 ELSE
565 v1 = 0._dp
566 v1(3) = magnetization
567 END IF
568 v2(1:3) = vs(1:3)
569 CALL sirius_add_atom(sctx, label, v2(1), v1(1))
570 END DO
571
572 CALL sirius_set_mpi_grid_dims(sctx, 2, mpi_grid_dims)
573
574! initialize global variables/indices/arrays/etc. of the simulation
575 CALL sirius_initialize_context(sctx)
576
577 ! strictly speaking the parameter use_symmetry is initialized at the
578 ! beginning but it does no harm to do it that way
579 IF (use_symmetry) THEN
580 CALL sirius_create_kset_from_grid(sctx, k_grid(1), k_shift(1), use_symmetry=.true., kset_handler=ks_handler)
581 ELSE
582 CALL sirius_create_kset_from_grid(sctx, k_grid(1), k_shift(1), use_symmetry=.false., kset_handler=ks_handler)
583 END IF
584! create ground-state class
585 CALL sirius_create_ground_state(ks_handler, gs_handler)
586
587 CALL pwdft_env_set(pwdft_env, sctx=sctx, gs_handler=gs_handler, ks_handler=ks_handler)
588#endif
589 END SUBROUTINE cp_sirius_create_env
590
591!***************************************************************************************************
592!> \brief ...
593!> \param pwdft_env ...
594!> \param
595!> \par History
596!> 07.2018 Update the Sirius environment
597!> \author JHU
598! **************************************************************************************************
599 SUBROUTINE cp_sirius_update_context(pwdft_env)
600 TYPE(pwdft_environment_type), POINTER :: pwdft_env
601
602 INTEGER :: iatom, natom
603 REAL(KIND=c_double), DIMENSION(3) :: a1, a2, a3, v2
604 REAL(KIND=dp), DIMENSION(3) :: vr, vs
605 TYPE(cell_type), POINTER :: my_cell
606 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
607 TYPE(qs_subsys_type), POINTER :: qs_subsys
608 TYPE(sirius_context_handler) :: sctx
609 TYPE(sirius_ground_state_handler) :: gs_handler
610
611 cpassert(ASSOCIATED(pwdft_env))
612 CALL pwdft_env_get(pwdft_env, sctx=sctx, gs_handler=gs_handler)
613
614! get current positions and lattice vectors
615 CALL pwdft_env_get(pwdft_env=pwdft_env, qs_subsys=qs_subsys)
616
617! lattice vectors of the unit cell should be in [a.u.] (length is in [a.u.])
618 CALL qs_subsys_get(qs_subsys, cell=my_cell)
619 a1(:) = my_cell%hmat(:, 1)
620 a2(:) = my_cell%hmat(:, 2)
621 a3(:) = my_cell%hmat(:, 3)
622 CALL sirius_set_lattice_vectors(sctx, a1(1), a2(1), a3(1))
623
624! new atomic positions
625 CALL qs_subsys_get(qs_subsys, particle_set=particle_set)
626 natom = SIZE(particle_set)
627 DO iatom = 1, natom
628 vr(1:3) = particle_set(iatom)%r(1:3)
629 CALL real_to_scaled(vs, vr, my_cell)
630 v2(1:3) = vs(1:3)
631 CALL sirius_set_atom_position(sctx, iatom, v2(1))
632 END DO
633
634! update ground-state class
635 CALL sirius_update_ground_state(gs_handler)
636
637 CALL pwdft_env_set(pwdft_env, sctx=sctx, gs_handler=gs_handler)
638
639 END SUBROUTINE cp_sirius_update_context
640
641! **************************************************************************************************
642!> \brief ...
643!> \param sctx ...
644!> \param section ...
645!> \param section_name ...
646! **************************************************************************************************
647 SUBROUTINE cp_sirius_fill_in_section(sctx, section, section_name)
648 TYPE(sirius_context_handler), INTENT(INOUT) :: sctx
649 TYPE(section_vals_type), POINTER :: section
650 CHARACTER(*), INTENT(in) :: section_name
651
652 CHARACTER(len=256), TARGET :: option_name
653 CHARACTER(len=4096) :: description, usage
654 CHARACTER(len=80), DIMENSION(:), POINTER :: tmp
655 CHARACTER(len=80), TARGET :: str
656 INTEGER :: ctype, elem, ic, j
657 INTEGER, DIMENSION(:), POINTER :: ivals
658 INTEGER, TARGET :: enum_length, ival, length, &
659 num_possible_values, number_of_options
660 LOGICAL :: explicit
661 LOGICAL, DIMENSION(:), POINTER :: lvals
662 LOGICAL, TARGET :: found, lval
663 REAL(kind=dp), DIMENSION(:), POINTER :: rvals
664 REAL(kind=dp), TARGET :: rval
665
666 NULLIFY (rvals)
667 NULLIFY (ivals)
668 CALL sirius_option_get_section_length(section_name, number_of_options)
669
670 DO elem = 1, number_of_options
671 option_name = ''
672 CALL sirius_option_get_info(section_name, &
673 elem, &
674 option_name, &
675 256, &
676 ctype, &
677 num_possible_values, &
678 enum_length, &
679 description, &
680 4096, &
681 usage, &
682 4096)
683
684 ! ignore the keyword parametes for sections dftd3 and dftd4.
685 IF (((section_name == 'dftd3') .OR. (section_name == 'dftd4')) .AND. (option_name == 'parameters')) THEN
686 cycle
687 END IF
688
689 IF ((option_name /= 'memory_usage') .AND. (option_name /= 'xc_functionals') .AND. (option_name /= 'vk')) THEN
690 CALL section_vals_val_get(section, option_name, explicit=found)
691 IF (found) THEN
692 SELECT CASE (ctype)
693 CASE (sirius_integer_type)
694 CALL section_vals_val_get(section, option_name, i_val=ival)
695 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(ival))
696 CASE (sirius_number_type)
697 CALL section_vals_val_get(section, option_name, r_val=rval)
698 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(rval))
699 CASE (sirius_logical_type)
700 CALL section_vals_val_get(section, option_name, l_val=lval)
701 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(lval))
702 CASE (sirius_string_type) ! string nightmare
703 str = ''
704 CALL section_vals_val_get(section, option_name, explicit=explicit, c_val=str)
705 str = trim(adjustl(str))
706 DO j = 1, len(str)
707 ic = ichar(str(j:j))
708 IF (ic >= 65 .AND. ic < 90) str(j:j) = char(ic + 32)
709 END DO
710
711 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(str), max_length=len_trim(str))
712 CASE (sirius_integer_array_type)
713 CALL section_vals_val_get(section, option_name, i_vals=ivals)
714 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(ivals(1)), &
715 max_length=num_possible_values)
716 CASE (sirius_number_array_type)
717 CALL section_vals_val_get(section, option_name, r_vals=rvals)
718 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(rvals(1)), &
719 max_length=num_possible_values)
720 CASE (sirius_logical_array_type)
721 CALL section_vals_val_get(section, option_name, l_vals=lvals)
722 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(lvals(1)), &
723 max_length=num_possible_values)
724 CASE (sirius_string_array_type)
725 CALL section_vals_val_get(section, option_name, explicit=explicit, n_rep_val=length)
726 DO j = 1, length
727 str = ''
728 CALL section_vals_val_get(section, option_name, i_rep_val=j, explicit=explicit, c_vals=tmp)
729 str = trim(adjustl(tmp(j)))
730 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(str), &
731 max_length=len_trim(str), append=.true.)
732 END DO
733 CASE DEFAULT
734 END SELECT
735 END IF
736 END IF
737 END DO
738 END SUBROUTINE cp_sirius_fill_in_section
739
740!***************************************************************************************************
741!> \brief ...
742!> \param pwdft_env ...
743!> \param calculate_forces ...
744!> \param calculate_stress_tensor ...
745!> \param
746!> \par History
747!> 07.2018 start the Sirius library
748!> \author JHU
749! **************************************************************************************************
750 SUBROUTINE cp_sirius_energy_force(pwdft_env, calculate_forces, calculate_stress_tensor)
751 TYPE(pwdft_environment_type), INTENT(INOUT), &
752 POINTER :: pwdft_env
753 LOGICAL, INTENT(IN) :: calculate_forces, calculate_stress_tensor
754
755 INTEGER :: iw, n1, n2
756 LOGICAL :: do_print, gs_converged
757 REAL(KIND=c_double) :: etotal
758 REAL(KIND=c_double), ALLOCATABLE, DIMENSION(:, :) :: cforces
759 REAL(KIND=c_double), DIMENSION(3, 3) :: cstress
760 REAL(KIND=dp), DIMENSION(3, 3) :: stress
761 REAL(KIND=dp), DIMENSION(:, :), POINTER :: forces
762 TYPE(cp_logger_type), POINTER :: logger
763 TYPE(pwdft_energy_type), POINTER :: energy
764 TYPE(section_vals_type), POINTER :: print_section, pwdft_input
765 TYPE(sirius_ground_state_handler) :: gs_handler
766
767 cpassert(ASSOCIATED(pwdft_env))
768
769 NULLIFY (logger)
770 logger => cp_get_default_logger()
772
773 CALL pwdft_env_get(pwdft_env=pwdft_env, gs_handler=gs_handler)
774 CALL sirius_find_ground_state(gs_handler, converged=gs_converged)
775
776 IF (gs_converged) THEN
777 IF (iw > 0) WRITE (iw, '(A)') "CP2K/SIRIUS: ground state is converged"
778 ELSE
779 IF (pwdft_env%ignore_convergence_failure) THEN
780 IF (iw > 0) WRITE (iw, '(A)') "CP2K/SIRIUS Warning: ground state is not converged"
781 ELSE
782 cpabort("CP2K/SIRIUS (ground state): SIRIUS did not converge.")
783 END IF
784 END IF
785 IF (iw > 0) CALL m_flush(iw)
786
787 CALL pwdft_env_get(pwdft_env=pwdft_env, energy=energy)
788 etotal = 0.0_c_double
789
790 CALL sirius_get_energy(gs_handler, 'band-gap', etotal)
791 energy%band_gap = etotal
792
793 etotal = 0.0_c_double
794 CALL sirius_get_energy(gs_handler, 'total', etotal)
795 energy%etotal = etotal
796
797 ! extract entropy (TS returned by sirius is always negative, sign
798 ! convention in QE)
799 etotal = 0.0_c_double
800 CALL sirius_get_energy(gs_handler, 'demet', etotal)
801 energy%entropy = -etotal
802
803 IF (calculate_forces) THEN
804 CALL pwdft_env_get(pwdft_env=pwdft_env, forces=forces)
805 n1 = SIZE(forces, 1)
806 n2 = SIZE(forces, 2)
807
808 ALLOCATE (cforces(n2, n1))
809 cforces = 0.0_c_double
810 CALL sirius_get_forces(gs_handler, 'total', cforces)
811 ! Sirius computes the forces but cp2k use the gradient everywhere
812 ! so a minus sign is needed.
813 ! note also that sirius and cp2k store the forces transpose to each other
814 ! sirius : forces(coordinates, atoms)
815 ! cp2k : forces(atoms, coordinates)
816 forces = -transpose(cforces(:, :))
817 DEALLOCATE (cforces)
818 END IF
819
820 IF (calculate_stress_tensor) THEN
821 cstress = 0.0_c_double
822 CALL sirius_get_stress_tensor(gs_handler, 'total', cstress)
823 stress(1:3, 1:3) = cstress(1:3, 1:3)
824 CALL pwdft_env_set(pwdft_env=pwdft_env, stress=stress)
825 END IF
826
827 CALL pwdft_env_get(pwdft_env=pwdft_env, pwdft_input=pwdft_input)
828 print_section => section_vals_get_subs_vals(pwdft_input, "PRINT")
829 CALL section_vals_get(print_section, explicit=do_print)
830 IF (do_print) THEN
831 CALL cp_sirius_print_results(pwdft_env, print_section)
832 END IF
833 END SUBROUTINE cp_sirius_energy_force
834
835!***************************************************************************************************
836!> \brief ...
837!> \param pwdft_env ...
838!> \param print_section ...
839!> \param
840!> \par History
841!> 12.2019 init
842!> \author JHU
843! **************************************************************************************************
844 SUBROUTINE cp_sirius_print_results(pwdft_env, print_section)
845 TYPE(pwdft_environment_type), INTENT(INOUT), &
846 POINTER :: pwdft_env
847 TYPE(section_vals_type), POINTER :: print_section
848
849 CHARACTER(LEN=default_string_length) :: my_act, my_pos
850 INTEGER :: i, ik, iounit, ispn, iterstep, iv, iw, &
851 nbands, nhist, nkpts, nspins
852 INTEGER(KIND=C_INT) :: cint
853 LOGICAL :: append, dos, ionode
854 REAL(KIND=c_double) :: creal
855 REAL(KIND=c_double), ALLOCATABLE, DIMENSION(:) :: slist
856 REAL(KIND=dp) :: de, e_fermi(2), emax, emin, eval
857 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: wkpt
858 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: ehist, hist, occval
859 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: energies, occupations
860 TYPE(cp_logger_type), POINTER :: logger
861 TYPE(sirius_context_handler) :: sctx
862 TYPE(sirius_ground_state_handler) :: gs_handler
863 TYPE(sirius_kpoint_set_handler) :: ks_handler
864
865 NULLIFY (logger)
866 logger => cp_get_default_logger()
867 ionode = logger%para_env%is_source()
868 iounit = cp_logger_get_default_io_unit(logger)
869
870 ! Density of States
871 dos = btest(cp_print_key_should_output(logger%iter_info, print_section, "DOS"), cp_p_file)
872 IF (dos) THEN
873 CALL pwdft_env_get(pwdft_env, ks_handler=ks_handler)
874 CALL pwdft_env_get(pwdft_env, gs_handler=gs_handler)
875 CALL pwdft_env_get(pwdft_env, sctx=sctx)
876
877 CALL section_vals_val_get(print_section, "DOS%DELTA_E", r_val=de)
878 CALL section_vals_val_get(print_section, "DOS%APPEND", l_val=append)
879
880 CALL sirius_get_num_kpoints(ks_handler, cint)
881 nkpts = cint
882 CALL sirius_get_parameters(sctx, num_bands=cint)
883 nbands = cint
884 CALL sirius_get_parameters(sctx, num_spins=cint)
885 nspins = cint
886 e_fermi(:) = 0.0_dp
887 ALLOCATE (energies(nbands, nspins, nkpts))
888 energies = 0.0_dp
889 ALLOCATE (occupations(nbands, nspins, nkpts))
890 occupations = 0.0_dp
891 ALLOCATE (wkpt(nkpts))
892 ALLOCATE (slist(nbands))
893 DO ik = 1, nkpts
894 CALL sirius_get_kpoint_properties(ks_handler, ik, creal)
895 wkpt(ik) = creal
896 END DO
897 DO ik = 1, nkpts
898 DO ispn = 1, nspins
899 CALL sirius_get_band_energies(ks_handler, ik, ispn, slist)
900 energies(1:nbands, ispn, ik) = slist(1:nbands)
901 CALL sirius_get_band_occupancies(ks_handler, ik, ispn, slist)
902 occupations(1:nbands, ispn, ik) = slist(1:nbands)
903 END DO
904 END DO
905 emin = minval(energies)
906 emax = maxval(energies)
907 nhist = nint((emax - emin)/de) + 1
908 ALLOCATE (hist(nhist, nspins), occval(nhist, nspins), ehist(nhist, nspins))
909 hist = 0.0_dp
910 occval = 0.0_dp
911 ehist = 0.0_dp
912
913 DO ik = 1, nkpts
914 DO ispn = 1, nspins
915 DO i = 1, nbands
916 eval = energies(i, ispn, ik) - emin
917 iv = nint(eval/de) + 1
918 cpassert((iv > 0) .AND. (iv <= nhist))
919 hist(iv, ispn) = hist(iv, ispn) + wkpt(ik)
920 occval(iv, ispn) = occval(iv, ispn) + wkpt(ik)*occupations(i, ispn, ik)
921 END DO
922 END DO
923 END DO
924 hist = hist/real(nbands, kind=dp)
925 DO i = 1, nhist
926 ehist(i, 1:nspins) = emin + (i - 1)*de
927 END DO
928
929 iterstep = logger%iter_info%iteration(logger%iter_info%n_rlevel)
930 my_act = "WRITE"
931 IF (append .AND. iterstep > 1) THEN
932 my_pos = "APPEND"
933 ELSE
934 my_pos = "REWIND"
935 END IF
936
937 iw = cp_print_key_unit_nr(logger, print_section, "DOS", &
938 extension=".dos", file_position=my_pos, file_action=my_act, &
939 file_form="FORMATTED")
940 IF (iw > 0) THEN
941 IF (nspins == 2) THEN
942 WRITE (unit=iw, fmt="(T2,A,I0,A,2F12.6)") &
943 "# DOS at iteration step i = ", iterstep, ", E_Fermi[a.u.] = ", e_fermi(1:2)
944 WRITE (unit=iw, fmt="(T2,A, A)") " Energy[a.u.] Alpha_Density Occupation", &
945 " Beta_Density Occupation"
946 ELSE
947 WRITE (unit=iw, fmt="(T2,A,I0,A,F12.6)") &
948 "# DOS at iteration step i = ", iterstep, ", E_Fermi[a.u.] = ", e_fermi(1)
949 WRITE (unit=iw, fmt="(T2,A)") " Energy[a.u.] Density Occupation"
950 END IF
951 DO i = 1, nhist
952 eval = emin + (i - 1)*de
953 IF (nspins == 2) THEN
954 WRITE (unit=iw, fmt="(F15.8,4F15.4)") eval, hist(i, 1), occval(i, 1), &
955 hist(i, 2), occval(i, 2)
956 ELSE
957 WRITE (unit=iw, fmt="(F15.8,2F15.4)") eval, hist(i, 1), occval(i, 1)
958 END IF
959 END DO
960 END IF
961 CALL cp_print_key_finished_output(iw, logger, print_section, "DOS")
962
963 DEALLOCATE (energies, occupations, wkpt, slist)
964 DEALLOCATE (hist, occval, ehist)
965 END IF
966 END SUBROUTINE cp_sirius_print_results
967
968END MODULE sirius_interface
969
970#else
971
972!***************************************************************************************************
973!> \brief Empty implementation in case SIRIUS is not compiled in.
974!***************************************************************************************************
977#include "./base/base_uses.f90"
978
979 IMPLICIT NONE
980
981 PRIVATE
982
983 ! Public subroutines
984
985 PUBLIC :: cp_sirius_create_env, &
991
992CONTAINS
993
994! **************************************************************************************************
995!> \brief Empty implementation in case SIRIUS is not compiled in.
996! **************************************************************************************************
997 SUBROUTINE cp_sirius_init()
998 END SUBROUTINE cp_sirius_init
999
1000! **************************************************************************************************
1001!> \brief Return always .FALSE. because the Sirius library is not compiled in.
1002!> \return Return the initialisation status of the Sirius library as boolean
1003! **************************************************************************************************
1004 LOGICAL FUNCTION cp_sirius_is_initialized()
1005 cp_sirius_is_initialized = .false.
1006 END FUNCTION cp_sirius_is_initialized
1007
1008! **************************************************************************************************
1009!> \brief Empty implementation in case SIRIUS is not compiled in.
1010! **************************************************************************************************
1012 END SUBROUTINE cp_sirius_finalize
1013
1014! **************************************************************************************************
1015!> \brief Empty implementation in case SIRIUS is not compiled in.
1016!> \param pwdft_env ...
1017! **************************************************************************************************
1018 SUBROUTINE cp_sirius_create_env(pwdft_env)
1019 TYPE(pwdft_environment_type), POINTER :: pwdft_env
1020
1021 mark_used(pwdft_env)
1022 cpabort("Sirius library is missing")
1023 END SUBROUTINE cp_sirius_create_env
1024
1025! **************************************************************************************************
1026!> \brief Empty implementation in case SIRIUS is not compiled in.
1027!> \param pwdft_env ...
1028!> \param calculate_forces ...
1029!> \param calculate_stress ...
1030! **************************************************************************************************
1031 SUBROUTINE cp_sirius_energy_force(pwdft_env, calculate_forces, calculate_stress)
1032 TYPE(pwdft_environment_type), POINTER :: pwdft_env
1033 LOGICAL :: calculate_forces, calculate_stress
1034
1035 mark_used(pwdft_env)
1036 mark_used(calculate_forces)
1037 mark_used(calculate_stress)
1038 cpabort("Sirius library is missing")
1039 END SUBROUTINE cp_sirius_energy_force
1040
1041! **************************************************************************************************
1042!> \brief Empty implementation in case SIRIUS is not compiled in.
1043!> \param pwdft_env ...
1044! **************************************************************************************************
1045 SUBROUTINE cp_sirius_update_context(pwdft_env)
1046 TYPE(pwdft_environment_type), POINTER :: pwdft_env
1047
1048 mark_used(pwdft_env)
1049 cpabort("Sirius library is missing")
1050 END SUBROUTINE cp_sirius_update_context
1051
1052END MODULE sirius_interface
1053
1054#endif
calculate the orbitals for a given atomic kind type
subroutine, public calculate_atomic_orbitals(atomic_kind, qs_kind, agrid, iunit, pmat, fmat, density, wavefunction, wfninfo, confine, xc_section, nocc)
...
subroutine, public gth_potential_conversion(gth_potential, gth_atompot)
...
Define the atom type and its sub types.
Definition atom_types.F:15
Routines that process Quantum Espresso UPF files.
Definition atom_upf.F:14
Some basic routines for atomic calculations.
Definition atom_utils.F:15
pure subroutine, public atom_local_potential(locpot, gthpot, rr)
...
Definition atom_utils.F:799
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
Definition cell_types.F:491
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...
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...
Definition of the atomic potential types.
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_gapw_log
integer, parameter, public sirius_func_vdwdfcx
integer, parameter, public sirius_func_vdwdf2
integer, parameter, public sirius_func_vdwdf
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_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
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
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:136
Definition of mathematical constants and functions.
real(kind=dp), dimension(0:maxfac), parameter, public gamma1
real(kind=dp), parameter, public fourpi
Interface to the message passing library MPI.
Define the data structure for the particle information.
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public massunit
Definition physcon.F:141
The type definitions for the PWDFT environment.
subroutine, public pwdft_env_get(pwdft_env, pwdft_input, force_env_input, xc_input, cp_subsys, qs_subsys, para_env, energy, forces, stress, sctx, gs_handler, ks_handler)
Returns various attributes of the pwdft environment.
subroutine, public pwdft_env_set(pwdft_env, pwdft_input, force_env_input, xc_input, qs_subsys, cp_subsys, para_env, energy, forces, stress, sctx, gs_handler, ks_handler)
Sets various attributes of the pwdft environment.
subroutine, public deallocate_grid_atom(grid_atom)
Deallocate a Gaussian-type orbital (GTO) basis set data set.
subroutine, public allocate_grid_atom(grid_atom)
Initialize components of the grid_atom_type structure.
subroutine, public create_grid_atom(grid_atom, nr, na, llmax, ll, quadrature)
...
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, 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.
types that represent a quickstep subsys
subroutine, public qs_subsys_get(subsys, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell, cell_ref, use_ref_cell, energy, force, qs_kind_set, cp_subsys, nelectron_total, nelectron_spin)
...
Interface to the SIRIUS Library.
subroutine, public cp_sirius_update_context(pwdft_env)
Empty implementation in case SIRIUS is not compiled in.
subroutine, public cp_sirius_init()
Empty implementation in case SIRIUS is not compiled in.
subroutine, public cp_sirius_energy_force(pwdft_env, calculate_forces, calculate_stress)
Empty implementation in case SIRIUS is not compiled in.
logical function, public cp_sirius_is_initialized()
Return always .FALSE. because the Sirius library is not compiled in.
subroutine, public cp_sirius_finalize()
Empty implementation in case SIRIUS is not compiled in.
subroutine, public cp_sirius_create_env(pwdft_env)
Empty implementation in case SIRIUS is not compiled in.
Provides all information about a pseudopotential.
Definition atom_types.F:98
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:60
type of a logger, at the moment it contains just a print level starting at which level it should be l...
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.