33#include "./base/base_uses.f90"
39 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_mo_occupation'
44 MODULE PROCEDURE set_mo_occupation_1, set_mo_occupation_2
59 SUBROUTINE set_mo_occupation_3(mo_array, smear)
61 TYPE(mo_set_type),
DIMENSION(2),
INTENT(INOUT) :: mo_array
62 TYPE(smear_type),
POINTER :: smear
64 CHARACTER(LEN=*),
PARAMETER :: routineN =
'set_mo_occupation_3'
66 INTEGER :: all_nmo, handle, homo_a, homo_b, i, &
67 lfomo_a, lfomo_b, nmo_a, nmo_b, &
69 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: all_index
71 REAL(KIND=
dp) :: all_nelec, kts, mu, nelec_a, nelec_b, &
73 REAL(KIND=
dp),
ALLOCATABLE,
DIMENSION(:) :: all_eigval, all_occ
74 REAL(KIND=
dp),
DIMENSION(:),
POINTER :: eigval_a, eigval_b, occ_a, occ_b
76 CALL timeset(routinen, handle)
78 NULLIFY (eigval_a, eigval_b, occ_a, occ_b)
79 CALL get_mo_set(mo_set=mo_array(1), nmo=nmo_a, eigenvalues=eigval_a, &
80 occupation_numbers=occ_a)
81 CALL get_mo_set(mo_set=mo_array(2), nmo=nmo_b, eigenvalues=eigval_b, &
82 occupation_numbers=occ_b)
83 all_nmo = nmo_a + nmo_b
84 ALLOCATE (all_eigval(all_nmo))
85 ALLOCATE (all_occ(all_nmo))
86 ALLOCATE (all_index(all_nmo))
88 all_eigval(1:nmo_a) = eigval_a(1:nmo_a)
89 all_eigval(nmo_a + 1:all_nmo) = eigval_b(1:nmo_b)
91 CALL sort(all_eigval, all_nmo, all_index)
101 all_nelec = nelec_a + nelec_b
104 IF (all_index(i) <= nmo_a)
THEN
105 all_occ(i) = occ_a(all_index(i))
107 all_occ(i) = occ_b(all_index(i) - nmo_a)
111 CALL fermifixed(all_occ, mu, kts, all_eigval, all_nelec, &
112 smear%electronic_temperature, 1._dp, xas_estate, occ_estate)
114 is_large = abs(maxval(all_occ) - 1.0_dp) > smear%eps_fermi_dirac
116 cpwarn_if(is_large,
"Fermi-Dirac smearing includes the first MO")
118 is_large = abs(minval(all_occ)) > smear%eps_fermi_dirac
120 CALL cp_warn(__location__, &
121 "Fermi-Dirac smearing includes the last MO => "// &
122 "Add more MOs for proper smearing.")
125 is_large = (abs(all_nelec -
accurate_sum(all_occ(:))) > smear%eps_fermi_dirac*all_nelec)
126 cpwarn_if(is_large,
"Total number of electrons is not accurate")
129 IF (all_index(i) <= nmo_a)
THEN
130 occ_a(all_index(i)) = all_occ(i)
131 eigval_a(all_index(i)) = all_eigval(i)
133 occ_b(all_index(i) - nmo_a) = all_occ(i)
134 eigval_b(all_index(i) - nmo_a) = all_eigval(i)
142 IF (occ_a(i) < 1.0_dp)
THEN
148 IF (occ_b(i) < 1.0_dp)
THEN
154 DO i = nmo_a, lfomo_a, -1
155 IF (occ_a(i) > smear%eps_fermi_dirac)
THEN
161 DO i = nmo_b, lfomo_b, -1
162 IF (occ_b(i) > smear%eps_fermi_dirac)
THEN
168 CALL set_mo_set(mo_set=mo_array(1), kts=kts/2.0_dp, mu=mu, n_el_f=nelec_a, &
169 lfomo=lfomo_a, homo=homo_a, uniform_occupation=.false.)
170 CALL set_mo_set(mo_set=mo_array(2), kts=kts/2.0_dp, mu=mu, n_el_f=nelec_b, &
171 lfomo=lfomo_b, homo=homo_b, uniform_occupation=.false.)
173 CALL timestop(handle)
175 END SUBROUTINE set_mo_occupation_3
194 SUBROUTINE set_mo_occupation_2(mo_array, smear, eval_deriv, tot_zeff_corr)
196 TYPE(
mo_set_type),
DIMENSION(:),
INTENT(INOUT) :: mo_array
198 REAL(KIND=
dp),
DIMENSION(:),
OPTIONAL,
POINTER :: eval_deriv
199 REAL(KIND=
dp),
OPTIONAL :: tot_zeff_corr
201 CHARACTER(LEN=*),
PARAMETER :: routineN =
'set_mo_occupation_2'
203 INTEGER :: handle, i, lumo_a, lumo_b, &
204 multiplicity_new, multiplicity_old, &
206 REAL(KIND=
dp) :: nelec_f, threshold
207 REAL(KIND=
dp),
DIMENSION(:),
POINTER :: eigval_a, eigval_b
209 CALL timeset(routinen, handle)
212 IF (
SIZE(mo_array) == 1)
THEN
213 IF (
PRESENT(eval_deriv))
THEN
215 CALL set_mo_occupation_1(mo_array(1), smear=smear, eval_deriv=eval_deriv)
217 IF (
PRESENT(tot_zeff_corr))
THEN
218 CALL set_mo_occupation_1(mo_array(1), smear=smear, tot_zeff_corr=tot_zeff_corr)
220 CALL set_mo_occupation_1(mo_array(1), smear=smear)
223 CALL timestop(handle)
227 IF (smear%do_smear)
THEN
228 IF (smear%fixed_mag_mom < 0.0_dp)
THEN
229 IF (
PRESENT(tot_zeff_corr))
THEN
230 CALL cp_warn(__location__, &
231 "CORE_CORRECTION /= 0.0 might cause the cell to charge up "// &
232 "that will lead to application of different background "// &
233 "correction compared to the reference system. "// &
234 "Use FIXED_MAGNETIC_MOMENT >= 0.0 if using SMEAR keyword "// &
235 "to correct the electron density")
237 IF (smear%fixed_mag_mom /= -1.0_dp)
THEN
238 cpassert(.NOT. (
PRESENT(eval_deriv)))
239 CALL set_mo_occupation_3(mo_array, smear=smear)
240 CALL timestop(handle)
244 nelec_f = mo_array(1)%n_el_f + mo_array(2)%n_el_f
245 IF (abs((mo_array(1)%n_el_f - mo_array(2)%n_el_f) - smear%fixed_mag_mom) > smear%eps_fermi_dirac*nelec_f)
THEN
246 mo_array(1)%n_el_f = nelec_f/2.0_dp + smear%fixed_mag_mom/2.0_dp
247 mo_array(2)%n_el_f = nelec_f/2.0_dp - smear%fixed_mag_mom/2.0_dp
249 cpassert(.NOT. (
PRESENT(eval_deriv)))
250 IF (
PRESENT(tot_zeff_corr))
THEN
251 CALL set_mo_occupation_1(mo_array(1), smear=smear, tot_zeff_corr=tot_zeff_corr)
252 CALL set_mo_occupation_1(mo_array(2), smear=smear, tot_zeff_corr=tot_zeff_corr)
254 CALL set_mo_occupation_1(mo_array(1), smear=smear)
255 CALL set_mo_occupation_1(mo_array(2), smear=smear)
260 IF (.NOT. ((mo_array(1)%flexible_electron_count > 0.0_dp) .AND. &
261 (mo_array(2)%flexible_electron_count > 0.0_dp)))
THEN
262 IF (
PRESENT(eval_deriv))
THEN
263 CALL set_mo_occupation_1(mo_array(1), smear=smear, eval_deriv=eval_deriv)
264 CALL set_mo_occupation_1(mo_array(2), smear=smear, eval_deriv=eval_deriv)
266 IF (
PRESENT(tot_zeff_corr))
THEN
267 CALL set_mo_occupation_1(mo_array(1), smear=smear, tot_zeff_corr=tot_zeff_corr)
268 CALL set_mo_occupation_1(mo_array(2), smear=smear, tot_zeff_corr=tot_zeff_corr)
270 CALL set_mo_occupation_1(mo_array(1), smear=smear)
271 CALL set_mo_occupation_1(mo_array(2), smear=smear)
274 CALL timestop(handle)
278 nelec = mo_array(1)%nelectron + mo_array(2)%nelectron
280 multiplicity_old = mo_array(1)%nelectron - mo_array(2)%nelectron + 1
282 IF (mo_array(1)%nelectron >= mo_array(1)%nmo) &
283 CALL cp_warn(__location__, &
284 "All alpha MOs are occupied. Add more alpha MOs to "// &
285 "allow for a higher multiplicity")
286 IF ((mo_array(2)%nelectron >= mo_array(2)%nmo) .AND. (mo_array(2)%nelectron /= mo_array(1)%nelectron)) &
287 CALL cp_warn(__location__,
"All beta MOs are occupied. Add more beta MOs to "// &
288 "allow for a lower multiplicity")
290 eigval_a => mo_array(1)%eigenvalues
291 eigval_b => mo_array(2)%eigenvalues
300 threshold = max(mo_array(1)%flexible_electron_count, mo_array(2)%flexible_electron_count)
301 IF ((eigval_a(lumo_a) - threshold) < eigval_b(lumo_b))
THEN
306 IF (lumo_a > mo_array(1)%nmo)
THEN
308 CALL cp_warn(__location__, &
309 "All alpha MOs are occupied. Add more alpha MOs to "// &
310 "allow for a higher multiplicity")
316 IF (lumo_b > mo_array(2)%nmo)
THEN
317 IF (lumo_b < lumo_a) &
318 CALL cp_warn(__location__, &
319 "All beta MOs are occupied. Add more beta MOs to "// &
320 "allow for a lower multiplicity")
328 mo_array(1)%homo = lumo_a - 1
329 mo_array(2)%homo = lumo_b - 1
331 IF (mo_array(2)%homo > mo_array(1)%homo)
THEN
332 CALL cp_warn(__location__, &
337 ") MOs are occupied. Resorting to low spin state")
338 mo_array(1)%homo = nelec/2 +
modulo(nelec, 2)
339 mo_array(2)%homo = nelec/2
342 mo_array(1)%nelectron = mo_array(1)%homo
343 mo_array(2)%nelectron = mo_array(2)%homo
344 multiplicity_new = mo_array(1)%nelectron - mo_array(2)%nelectron + 1
346 IF (multiplicity_new /= multiplicity_old) &
347 CALL cp_warn(__location__, &
348 "Multiplicity changed from "// &
349 trim(adjustl(
cp_to_string(multiplicity_old)))//
" to "// &
352 IF (
PRESENT(eval_deriv))
THEN
353 CALL set_mo_occupation_1(mo_array(1), smear=smear, eval_deriv=eval_deriv)
354 CALL set_mo_occupation_1(mo_array(2), smear=smear, eval_deriv=eval_deriv)
356 IF (
PRESENT(tot_zeff_corr))
THEN
357 CALL set_mo_occupation_1(mo_array(1), smear=smear, tot_zeff_corr=tot_zeff_corr)
358 CALL set_mo_occupation_1(mo_array(2), smear=smear, tot_zeff_corr=tot_zeff_corr)
360 CALL set_mo_occupation_1(mo_array(1), smear=smear)
361 CALL set_mo_occupation_1(mo_array(2), smear=smear)
365 CALL timestop(handle)
367 END SUBROUTINE set_mo_occupation_2
387 SUBROUTINE set_mo_occupation_1(mo_set, smear, eval_deriv, xas_env, tot_zeff_corr)
391 REAL(KIND=
dp),
DIMENSION(:),
OPTIONAL,
POINTER :: eval_deriv
393 REAL(kind=
dp),
OPTIONAL :: tot_zeff_corr
395 CHARACTER(LEN=*),
PARAMETER :: routinen =
'set_mo_occupation_1'
397 INTEGER :: handle, i_first, imo, ir, irmo, nmo, &
399 LOGICAL :: equal_size, is_large
400 REAL(kind=
dp) :: delectron, e1, e2, edelta, edist, &
401 el_count, lengthscale, nelec, &
402 occ_estate, total_zeff_corr, &
404 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: dfde
406 CALL timeset(routinen, handle)
408 cpassert(
ASSOCIATED(mo_set%eigenvalues))
409 cpassert(
ASSOCIATED(mo_set%occupation_numbers))
410 mo_set%occupation_numbers(:) = 0.0_dp
413 IF (mo_set%nelectron == 0)
THEN
414 CALL timestop(handle)
420 IF (
PRESENT(xas_env))
THEN
421 CALL get_xas_env(xas_env=xas_env, xas_nelectron=xas_nelectron, occ_estate=occ_estate, xas_estate=xas_estate)
422 nomo = ceiling(xas_nelectron + 1.0 - occ_estate - epsilon(0.0_dp))
424 mo_set%occupation_numbers(1:nomo) = mo_set%maxocc
425 IF (xas_estate > 0) mo_set%occupation_numbers(xas_estate) = occ_estate
426 el_count = sum(mo_set%occupation_numbers(1:nomo))
427 IF (el_count > xas_nelectron) &
428 mo_set%occupation_numbers(nomo) = mo_set%occupation_numbers(nomo) - (el_count - xas_nelectron)
429 el_count = sum(mo_set%occupation_numbers(1:nomo))
430 is_large = abs(el_count - xas_nelectron) > xas_nelectron*epsilon(el_count)
431 cpassert(.NOT. is_large)
433 IF (
modulo(mo_set%nelectron, int(mo_set%maxocc)) == 0)
THEN
434 nomo = nint(mo_set%nelectron/mo_set%maxocc)
436 mo_set%occupation_numbers(1:nomo) = mo_set%maxocc
438 nomo = int(mo_set%nelectron/mo_set%maxocc) + 1
440 mo_set%occupation_numbers(1:nomo - 1) = mo_set%maxocc
441 mo_set%occupation_numbers(nomo) = mo_set%nelectron - (nomo - 1)*mo_set%maxocc
448 IF (
PRESENT(tot_zeff_corr))
THEN
450 total_zeff_corr = tot_zeff_corr
451 IF (int(mo_set%maxocc) == 1) total_zeff_corr = total_zeff_corr/2.0_dp
453 IF (total_zeff_corr < 0.0_dp)
THEN
455 delectron = abs(total_zeff_corr) - real(mo_set%maxocc, kind=
dp)
456 IF (delectron > 0.0_dp)
THEN
457 mo_set%occupation_numbers(nomo) = 0.0_dp
458 irmo = ceiling(delectron/real(mo_set%maxocc, kind=
dp))
460 delectron = delectron - real(mo_set%maxocc, kind=
dp)
461 IF (delectron < 0.0_dp)
THEN
462 mo_set%occupation_numbers(nomo - ir) = -delectron
464 mo_set%occupation_numbers(nomo - ir) = 0.0_dp
468 IF (mo_set%occupation_numbers(nomo) == 0.0_dp) nomo = nomo - 1
469 ELSEIF (delectron < 0.0_dp)
THEN
470 mo_set%occupation_numbers(nomo) = -delectron
472 mo_set%occupation_numbers(nomo) = 0.0_dp
475 ELSEIF (total_zeff_corr > 0.0_dp)
THEN
477 delectron = total_zeff_corr - real(mo_set%maxocc, kind=
dp)
478 IF (delectron > 0.0_dp)
THEN
479 mo_set%occupation_numbers(nomo + 1) = real(mo_set%maxocc, kind=
dp)
481 irmo = ceiling(delectron/real(mo_set%maxocc, kind=
dp))
483 delectron = delectron - real(mo_set%maxocc, kind=
dp)
484 IF (delectron < 0.0_dp)
THEN
485 mo_set%occupation_numbers(nomo + ir) = delectron + real(mo_set%maxocc, kind=
dp)
487 mo_set%occupation_numbers(nomo + ir) = real(mo_set%maxocc, kind=
dp)
492 mo_set%occupation_numbers(nomo + 1) = total_zeff_corr
498 nmo =
SIZE(mo_set%eigenvalues)
500 cpassert(nmo >= nomo)
501 cpassert((
SIZE(mo_set%occupation_numbers) == nmo))
504 mo_set%lfomo = nomo + 1
505 mo_set%mu = mo_set%eigenvalues(nomo)
508 IF (
PRESENT(eval_deriv))
THEN
509 equal_size = (
SIZE(mo_set%occupation_numbers, 1) ==
SIZE(eval_deriv, 1))
514 IF (.NOT.
PRESENT(smear))
THEN
516 mo_set%uniform_occupation = .true.
517 IF (
PRESENT(eval_deriv))
THEN
520 CALL timestop(handle)
526 IF ((abs(mo_set%eigenvalues(1)) < 1.0e-12_dp) .AND. &
527 (abs(mo_set%eigenvalues(nmo)) < 1.0e-12_dp))
THEN
528 CALL timestop(handle)
534 IF (smear%do_smear)
THEN
535 IF (
PRESENT(xas_env))
THEN
536 i_first = xas_estate + 1
537 nelec = xas_nelectron
540 IF (smear%fixed_mag_mom == -1.0_dp)
THEN
541 nelec = real(mo_set%nelectron,
dp)
543 nelec = mo_set%n_el_f
546 SELECT CASE (smear%method)
548 IF (.NOT.
PRESENT(eval_deriv))
THEN
549 CALL fermifixed(mo_set%occupation_numbers(1:mo_set%nmo), mo_set%mu, mo_set%kTS, &
550 mo_set%eigenvalues(1:mo_set%nmo), nelec, &
551 smear%electronic_temperature, mo_set%maxocc, xas_estate, occ_estate)
555 ALLOCATE (dfde(nmo, nmo))
557 lengthscale = 10*smear%electronic_temperature
559 CALL fermifixedderiv(dfde, mo_set%occupation_numbers(1:mo_set%nmo), mo_set%mu, &
560 mo_set%kTS, mo_set%eigenvalues(1:mo_set%nmo), nelec, &
561 smear%electronic_temperature, mo_set%maxocc, lengthscale, xas_estate, occ_estate)
564 eval_deriv = eval_deriv - mo_set%eigenvalues + mo_set%mu
566 eval_deriv = matmul(transpose(dfde), eval_deriv)
572 DO imo = i_first, nmo
573 IF (mo_set%occupation_numbers(imo) < mo_set%maxocc)
THEN
578 is_large = abs(maxval(mo_set%occupation_numbers) - mo_set%maxocc) > smear%eps_fermi_dirac
580 cpwarn_if(is_large,
"Fermi-Dirac smearing includes the first MO")
583 DO imo = nmo, mo_set%lfomo, -1
584 IF (mo_set%occupation_numbers(imo) > smear%eps_fermi_dirac)
THEN
589 is_large = abs(minval(mo_set%occupation_numbers)) > smear%eps_fermi_dirac
591 CALL cp_warn(__location__, &
592 "Fermi-Dirac smearing includes the last MO => "// &
593 "Add more MOs for proper smearing.")
596 is_large = (abs(nelec -
accurate_sum(mo_set%occupation_numbers(:))) > smear%eps_fermi_dirac*nelec)
597 cpwarn_if(is_large,
"Total number of electrons is not accurate")
601 cpassert(.NOT.
PRESENT(eval_deriv))
604 e1 = mo_set%eigenvalues(mo_set%homo) - 0.5_dp*smear%window_size
605 IF (e1 <= mo_set%eigenvalues(1))
THEN
606 cpwarn(
"Energy window for smearing includes the first MO")
609 e2 = mo_set%eigenvalues(mo_set%homo) + 0.5_dp*smear%window_size
610 IF (e2 >= mo_set%eigenvalues(nmo)) &
611 CALL cp_warn(__location__, &
612 "Energy window for smearing includes the last MO => "// &
613 "Add more MOs for proper smearing.")
616 DO imo = i_first, nomo
617 IF (mo_set%eigenvalues(imo) > e1)
THEN
624 DO imo = nmo, nomo, -1
625 IF (mo_set%eigenvalues(imo) < e2)
THEN
635 DO imo = mo_set%lfomo, mo_set%homo
636 nelec = nelec + mo_set%occupation_numbers(imo)
637 edist = edist + abs(e2 - mo_set%eigenvalues(imo))
641 DO imo = mo_set%lfomo, mo_set%homo
642 edelta = abs(e2 - mo_set%eigenvalues(imo))
643 mo_set%occupation_numbers(imo) = min(mo_set%maxocc, nelec*edelta/edist)
644 nelec = nelec - mo_set%occupation_numbers(imo)
645 edist = edist - edelta
649 equal_size =
SIZE(mo_set%occupation_numbers, 1) ==
SIZE(smear%list, 1)
651 mo_set%occupation_numbers = smear%list
653 IF (
PRESENT(eval_deriv))
THEN
662 IF (mo_set%lfomo == mo_set%homo)
THEN
664 mo_set%lfomo = nomo + 1
666 mo_set%uniform_occupation = .false.
674 CALL timestop(handle)
676 END SUBROUTINE set_mo_occupation_1
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
various routines to log and control the output. The idea is that decisions about where to log should ...
deal with the Fermi distribution, compute it, fix mu, get derivs
subroutine, public fermifixedderiv(dfde, f, mu, kts, e, n, t, maxocc, l, estate, festate)
returns f and dfde for a given set of energies and number of electrons it is a numerical derivative,...
subroutine, public fermifixed(f, mu, kts, e, n, t, maxocc, estate, festate)
returns occupations according to Fermi-Dirac statistics for a given set of energies and number of ele...
sums arrays of real/complex numbers with much reduced round-off as compared to a naive implementation...
Defines the basic variable types.
integer, parameter, public dp
Set occupation of molecular orbitals.
Definition and initialisation of the mo data type.
subroutine, public set_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, uniform_occupation, kts, mu, flexible_electron_count)
Set the components of a MO set data structure.
logical function, public has_uniform_occupation(mo_set, first_mo, last_mo, occupation, tolerance)
Check if the set of MOs in mo_set specifed by the MO index range [first_mo,last_mo] an integer occupa...
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
parameters that control an scf iteration
All kind of helpful little routines.
define create destroy get and put information in xas_env to calculate the x-ray absorption spectra
subroutine, public get_xas_env(xas_env, exc_state, nao, nvirtual, nvirtual2, centers_wfn, atom_of_state, exc_atoms, nexc_states, type_of_state, mykind_of_atom, mykind_of_kind, state_of_atom, spectrum, groundstate_coeff, ostrength_sm, dip_fm_set, excvec_coeff, excvec_overlap, unoccupied_orbs, unoccupied_evals, unoccupied_max_iter, unoccupied_eps, all_vectors, all_evals, my_gto_basis, qs_loc_env, stogto_overlap, occ_estate, xas_nelectron, xas_estate, nexc_atoms, nexc_search, spin_channel, scf_env, scf_control)
...
contains the parameters needed by a scf run