16 USE iso_c_binding,
ONLY: c_double,&
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"
88 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'sirius_interface'
107 CALL sirius_initialize(.false.)
126 CALL sirius_finalize(.false., .false., .false.)
138 TYPE(pwdft_environment_type),
POINTER :: pwdft_env
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, &
174 TYPE(sirius_context_handler) :: sctx
175 TYPE(sirius_ground_state_handler) :: gs_handler
176 TYPE(sirius_kpoint_set_handler) :: ks_handler
178 cpassert(
ASSOCIATED(pwdft_env))
183 sirius_mpi_comm = para_env%get_handle()
184 CALL sirius_create_context(sirius_mpi_comm, sctx)
188 CALL pwdft_env_get(pwdft_env=pwdft_env, pwdft_input=pwdft_section, xc_input=xc_section)
191 l_val=pwdft_env%ignore_convergence_failure)
195 IF (
ASSOCIATED(xc_section))
THEN
200 IF (.NOT.
ASSOCIATED(xc_fun))
EXIT
203 CALL sirius_add_xc_functional(sctx,
"XC_"//trim(xc_fun%section%name))
209 IF (
ASSOCIATED(pwdft_sub_section))
THEN
210 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section,
"control")
217 IF (
ASSOCIATED(pwdft_sub_section))
THEN
218 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section,
"parameters")
236 SELECT CASE (vdw_func)
238 CALL sirius_add_xc_functional(sctx,
"XC_FUNC_VDWDF")
240 CALL sirius_add_xc_functional(sctx,
"XC_FUNC_VDWDF2")
242 CALL sirius_add_xc_functional(sctx,
"XC_FUNC_VDWDFCX")
251 IF (
ASSOCIATED(pwdft_sub_section))
THEN
252 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section,
"mixer")
257 IF (
ASSOCIATED(pwdft_sub_section))
THEN
258 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section,
"settings")
263 IF (
ASSOCIATED(pwdft_sub_section))
THEN
264 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section,
"iterative_solver")
267#if defined(__SIRIUS_DFTD3)
270 IF (
ASSOCIATED(pwdft_sub_section))
THEN
271 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section,
"dftd3")
274#if defined(__SIRIUS_DFTD4)
277 IF (
ASSOCIATED(pwdft_sub_section))
THEN
278 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section,
"dftd4")
282#if defined(__SIRIUS_NLCG)
285 IF (
ASSOCIATED(pwdft_sub_section))
THEN
286 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section,
"nlcg")
290#if defined(__SIRIUS_VCSQNM)
292 IF (
ASSOCIATED(pwdft_sub_section))
THEN
293 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section,
"vcsqnm")
298 CALL sirius_import_parameters(sctx,
'{}')
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))
308 IF (use_ref_cell)
THEN
309 cpwarn(
"SIRIUS| The specified CELL_REF will be ignored for PW_DFT calculations")
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)
320 name=label, element_symbol=element_symbol, mass=mass)
322 NULLIFY (upf_pot, gth_potential)
323 CALL get_qs_kind(qs_kind_set(ikind), upf_potential=upf_pot, gth_potential=gth_potential)
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))
330 ELSEIF (
ASSOCIATED(gth_potential))
THEN
337 ALLOCATE (rp(nmesh), fun(nmesh))
338 IF (atom_grid%rad(1) < atom_grid%rad(nmesh))
THEN
344 rp(1:nmesh) = atom_grid%rad(1:nmesh)
347 rp(i) = atom_grid%rad(nmesh - i + 1)
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)
357 ALLOCATE (gth_atompot)
360 fun(1:nmesh) = rp(1:nmesh)
361 CALL sirius_set_atom_type_radial_grid(sctx, label, nmesh, fun(1))
364 ALLOCATE (ef(nmesh), beta(nmesh))
367 IF (gth_atompot%nl(l) == 0) cycle
368 rl = gth_atompot%rcnl(l)
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))
374 pf = sqrt(2._dp)/(pf*sqrt(
gamma1(j)))
375 beta(:) = pf*rp**(l + 2*i - 2)*ef
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)
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)
387 DEALLOCATE (ef, beta)
391 IF (gth_atompot%soc)
THEN
392 nbs = 2*nbeta - gth_atompot%nl(0)
393 ALLOCATE (dion(nbs, nbs))
395 ALLOCATE (dion(nbeta, nbeta))
398 IF (gth_atompot%soc)
THEN
399 ns = gth_atompot%nl(0)
401 dion(1:ns, 1:ns) = gth_atompot%hnl(1:ns, 1:ns, 0)
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)
414 dion(ii + 1, jj + 1) = vnlm
418 CALL sirius_set_atom_type_dion(sctx, label, nbs, dion(1, 1))
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)
426 CALL sirius_set_atom_type_dion(sctx, label, nbeta, dion(1, 1))
432 IF (gth_atompot%nlcc)
THEN
433 ALLOCATE (corden(nmesh), fe(nmesh), rc(nmesh))
435 n = gth_atompot%nexp_nlcc
437 al = gth_atompot%alpha_nlcc(i)
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
445 fun(1:nmesh) = corden(1:nmesh)*rp(1:nmesh)
446 CALL sirius_add_atom_type_radial_function(sctx, label,
"ps_rho_core", &
448 DEALLOCATE (corden, fe, rc)
452 ALLOCATE (locpot(nmesh))
455 fun(1:nmesh) = locpot(1:nmesh)
456 CALL sirius_add_atom_type_radial_function(sctx, label,
"vloc", &
460 NULLIFY (density, wavefunction, wfninfo)
462 density=density, wavefunction=wavefunction, &
463 wfninfo=wfninfo, agrid=atom_grid)
466 DO iwf = 1,
SIZE(wavefunction, 2)
467 focc = wfninfo(1, iwf)
468 l = nint(wfninfo(2, iwf))
472 fun(1:nmesh) = wavefunction(1:nmesh, iwf)*rp(i)
475 fun(i) = wavefunction(nmesh - i + 1, iwf)*rp(i)
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)
485 fun(1:nmesh) =
fourpi*density(1:nmesh)*atom_grid%rad(1:nmesh)**2
488 fun(i) =
fourpi*density(nmesh - i + 1)*atom_grid%rad(nmesh - i + 1)**2
491 CALL sirius_add_atom_type_radial_function(sctx, label,
"ps_rho_total", &
494 IF (
ASSOCIATED(density))
DEALLOCATE (density)
495 IF (
ASSOCIATED(wavefunction))
DEALLOCATE (wavefunction)
496 IF (
ASSOCIATED(wfninfo))
DEALLOCATE (wfninfo)
500 DEALLOCATE (gth_atompot)
503 CALL cp_abort(__location__, &
504 "CP2K/SIRIUS: atomic kind needs UPF or GTH potential definition")
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)
519 IF (dft_plus_u_atom)
THEN
521 cpabort(
"CP2K/SIRIUS (hubbard): principal quantum number not specified")
525 cpabort(
"CP2K/SIRIUS (hubbard): l can not be negative.")
528 IF (occ_u < 0.0)
THEN
529 cpabort(
"CP2K/SIRIUS (hubbard): the occupation number can not be negative.")
533 IF (abs(u_minus_j) < 1e-8)
THEN
535 CALL sirius_set_atom_type_hubbard(sctx, label, lu, nu, &
536 occ_u, u_u, j_t, alpha_u, beta_u, j0_u)
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)
547 natom =
SIZE(particle_set)
549 vr(1:3) = particle_set(iatom)%r(1:3)
551 atomic_kind => particle_set(iatom)%atomic_kind
552 ikind = atomic_kind%kind_number
554 CALL get_qs_kind(qs_kind_set(ikind), zeff=zeff, magnetization=magnetization)
558 IF (num_mag_dims == 3)
THEN
561 v1(1) = magnetization*sin(angle1)*cos(angle2)
562 v1(2) = magnetization*sin(angle1)*sin(angle2)
563 v1(3) = magnetization*cos(angle1)
566 v1(3) = magnetization
569 CALL sirius_add_atom(sctx, label, v2(1), v1(1))
572 CALL sirius_set_mpi_grid_dims(sctx, 2, mpi_grid_dims)
575 CALL sirius_initialize_context(sctx)
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)
582 CALL sirius_create_kset_from_grid(sctx, k_grid(1), k_shift(1), use_symmetry=.false., kset_handler=ks_handler)
585 CALL sirius_create_ground_state(ks_handler, gs_handler)
587 CALL pwdft_env_set(pwdft_env, sctx=sctx, gs_handler=gs_handler, ks_handler=ks_handler)
600 TYPE(pwdft_environment_type),
POINTER :: pwdft_env
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
611 cpassert(
ASSOCIATED(pwdft_env))
612 CALL pwdft_env_get(pwdft_env, sctx=sctx, gs_handler=gs_handler)
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))
626 natom =
SIZE(particle_set)
628 vr(1:3) = particle_set(iatom)%r(1:3)
631 CALL sirius_set_atom_position(sctx, iatom, v2(1))
635 CALL sirius_update_ground_state(gs_handler)
637 CALL pwdft_env_set(pwdft_env, sctx=sctx, gs_handler=gs_handler)
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
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
661 LOGICAL,
DIMENSION(:),
POINTER :: lvals
662 LOGICAL,
TARGET :: found, lval
663 REAL(kind=
dp),
DIMENSION(:),
POINTER :: rvals
664 REAL(kind=
dp),
TARGET :: rval
668 CALL sirius_option_get_section_length(section_name, number_of_options)
670 DO elem = 1, number_of_options
672 CALL sirius_option_get_info(section_name, &
677 num_possible_values, &
685 IF (((section_name ==
'dftd3') .OR. (section_name ==
'dftd4')) .AND. (option_name ==
'parameters'))
THEN
689 IF ((option_name /=
'memory_usage') .AND. (option_name /=
'xc_functionals') .AND. (option_name /=
'vk'))
THEN
693 CASE (sirius_integer_type)
695 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(ival))
696 CASE (sirius_number_type)
698 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(rval))
699 CASE (sirius_logical_type)
701 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(lval))
702 CASE (sirius_string_type)
705 str = trim(adjustl(str))
708 IF (ic >= 65 .AND. ic < 90) str(j:j) = char(ic + 32)
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)
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)
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)
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)
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.)
738 END SUBROUTINE cp_sirius_fill_in_section
751 TYPE(pwdft_environment_type),
INTENT(INOUT), &
753 LOGICAL,
INTENT(IN) :: calculate_forces, calculate_stress_tensor
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
767 cpassert(
ASSOCIATED(pwdft_env))
773 CALL pwdft_env_get(pwdft_env=pwdft_env, gs_handler=gs_handler)
774 CALL sirius_find_ground_state(gs_handler, converged=gs_converged)
776 IF (gs_converged)
THEN
777 IF (iw > 0)
WRITE (iw,
'(A)')
"CP2K/SIRIUS: ground state is converged"
779 IF (pwdft_env%ignore_convergence_failure)
THEN
780 IF (iw > 0)
WRITE (iw,
'(A)')
"CP2K/SIRIUS Warning: ground state is not converged"
782 cpabort(
"CP2K/SIRIUS (ground state): SIRIUS did not converge.")
788 etotal = 0.0_c_double
790 CALL sirius_get_energy(gs_handler,
'band-gap', etotal)
791 energy%band_gap = etotal
793 etotal = 0.0_c_double
794 CALL sirius_get_energy(gs_handler,
'total', etotal)
795 energy%etotal = etotal
799 etotal = 0.0_c_double
800 CALL sirius_get_energy(gs_handler,
'demet', etotal)
801 energy%entropy = -etotal
803 IF (calculate_forces)
THEN
808 ALLOCATE (cforces(n2, n1))
809 cforces = 0.0_c_double
810 CALL sirius_get_forces(gs_handler,
'total', cforces)
816 forces = -transpose(cforces(:, :))
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)
827 CALL pwdft_env_get(pwdft_env=pwdft_env, pwdft_input=pwdft_input)
831 CALL cp_sirius_print_results(pwdft_env, print_section)
844 SUBROUTINE cp_sirius_print_results(pwdft_env, print_section)
845 TYPE(pwdft_environment_type),
INTENT(INOUT), &
847 TYPE(section_vals_type),
POINTER :: print_section
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
867 ionode = logger%para_env%is_source()
880 CALL sirius_get_num_kpoints(ks_handler, cint)
882 CALL sirius_get_parameters(sctx, num_bands=cint)
884 CALL sirius_get_parameters(sctx, num_spins=cint)
887 ALLOCATE (energies(nbands, nspins, nkpts))
889 ALLOCATE (occupations(nbands, nspins, nkpts))
891 ALLOCATE (wkpt(nkpts))
892 ALLOCATE (slist(nbands))
894 CALL sirius_get_kpoint_properties(ks_handler, ik, creal)
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)
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))
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)
924 hist = hist/real(nbands, kind=
dp)
926 ehist(i, 1:nspins) = emin + (i - 1)*de
929 iterstep = logger%iter_info%iteration(logger%iter_info%n_rlevel)
931 IF (append .AND. iterstep > 1)
THEN
938 extension=
".dos", file_position=my_pos, file_action=my_act, &
939 file_form=
"FORMATTED")
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"
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"
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)
957 WRITE (unit=iw, fmt=
"(F15.8,2F15.4)") eval, hist(i, 1), occval(i, 1)
963 DEALLOCATE (energies, occupations, wkpt, slist)
964 DEALLOCATE (hist, occval, ehist)
966 END SUBROUTINE cp_sirius_print_results
977#include "./base/base_uses.f90"
1021 mark_used(pwdft_env)
1022 cpabort(
"Sirius library is missing")
1033 LOGICAL :: calculate_forces, calculate_stress
1035 mark_used(pwdft_env)
1036 mark_used(calculate_forces)
1037 mark_used(calculate_stress)
1038 cpabort(
"Sirius library is missing")
1048 mark_used(pwdft_env)
1049 cpabort(
"Sirius library is missing")
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.
Routines that process Quantum Espresso UPF files.
Some basic routines for atomic calculations.
pure subroutine, public atom_local_potential(locpot, gthpot, rr)
...
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.
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
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.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
Machine interface based on Fortran 2003 and POSIX.
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
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:
real(kind=dp), parameter, public massunit
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.
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
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
The PWDFT environment type.
Provides all information about a quickstep kind.