63 subroutine unocc_run(system, from_scratch)
64 class(*),
intent(inout) :: system
65 logical,
intent(in) :: from_scratch
71 message(1) =
"CalculationMode = unocc not implemented for multi-system calculations"
82 type(electrons_t),
intent(inout) :: sys
83 logical,
intent(in) :: fromscratch
85 type(eigensolver_t) :: eigens
86 integer :: iunit, ierr, iter, ierr_rho, ik
87 integer(int64) :: what_it
88 logical :: read_gs, converged, forced_finish, showoccstates, is_orbital_dependent, occ_missing
89 integer :: max_iter, nst_calculated, showstart
90 integer :: n_filled, n_partially_filled, n_half_filled
91 integer,
allocatable :: lowest_missing(:, :), occ_states(:)
92 character(len=10) :: dirname
93 type(restart_t) :: restart_load_unocc, restart_load_gs, restart_dump
94 logical :: write_density, bandstructure_mode, read_td_states, output_iter
98 if (sys%hm%pcm%run_pcm)
then
105 if (max_iter < 0) max_iter = huge(max_iter)
116 call parse_variable(sys%namespace,
'UnoccShowOccStates', .false., showoccstates)
118 bandstructure_mode = .false.
119 if (sys%space%is_periodic() .and. sys%kpoints%get_kpoint_method() ==
kpoints_path)
then
120 bandstructure_mode = .
true.
123 call init_(sys%gr, sys%st)
126 read_td_states = .false.
127 if (bandstructure_mode)
then
136 call parse_variable(sys%namespace,
'UnoccUseTD', .false., read_td_states)
140 safe_allocate(lowest_missing(1:sys%st%d%dim, 1:sys%st%nik))
143 lowest_missing(:,:) = 0
146 if (.not. fromscratch)
then
148 mesh = sys%gr, exact = .
true.)
151 call states_elec_load(restart_load_unocc, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
152 ierr, lowest_missing = lowest_missing)
162 if (read_td_states)
then
170 if (ierr_rho == 0)
then
172 call states_elec_load(restart_load_gs, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
173 ierr, lowest_missing = lowest_missing)
176 call lda_u_load(restart_load_gs, sys%hm%lda_u, sys%st, sys%hm%energy%dft_u, ierr)
178 message(1) =
"Unable to read DFT+U information. DFT+U data will be calculated from states."
187 write_density = .
true.
190 safe_allocate(occ_states(1:sys%st%nik))
191 do ik = 1, sys%st%nik
192 call occupied_states(sys%st, sys%namespace, ik, n_filled, n_partially_filled, n_half_filled)
193 occ_states(ik) = n_filled + n_partially_filled + n_half_filled
196 is_orbital_dependent = (sys%ks%theory_level ==
hartree .or. sys%ks%theory_level ==
hartree_fock .or. &
200 if (is_orbital_dependent)
then
201 message(1) =
"Be sure your gs run is well converged since you have an orbital-dependent functional."
202 message(2) =
"Otherwise, the occupied states may change in CalculationMode = unocc, and your"
203 message(3) =
"unoccupied states will not be consistent with the gs run."
207 if (ierr_rho /= 0 .or. is_orbital_dependent)
then
208 occ_missing = .false.
209 do ik = 1, sys%st%nik
210 if (any(lowest_missing(1:sys%st%d%dim, ik) <= occ_states(ik)))
then
215 if (occ_missing)
then
216 if (is_orbital_dependent)
then
217 message(1) =
"For an orbital-dependent functional, all occupied orbitals must be provided."
218 else if (ierr_rho /= 0)
then
219 message(1) =
"Since density could not be read, all occupied orbitals must be provided."
222 message(2) =
"Not all the occupied orbitals could be read."
223 message(3) =
"Please run a ground-state calculation first!"
227 message(1) =
"Unable to read density: Building density from wavefunctions."
235 if (fromscratch .or. ierr /= 0)
then
236 if (fromscratch)
then
238 nst_calculated = min(maxval(occ_states), minval(lowest_missing) - 1)
241 nst_calculated = minval(lowest_missing) - 1
243 showstart = max(nst_calculated + 1, 1)
244 call lcao_run(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, &
245 sys%hm, st_start = showstart)
248 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, sys%st, sys%ions, sys%ext_partners, &
249 calc_eigenval = .false.)
250 showstart = minval(occ_states(:)) + 1
257 if (showstart > sys%st%nst) showstart = 1
259 safe_deallocate_a(lowest_missing)
261 if (showoccstates) showstart = 1
265 if (bandstructure_mode)
then
266 message(1) =
"Info: The code will run in band structure mode."
267 message(2) =
" No restart information will be printed."
271 if (.not. bandstructure_mode)
then
276 if (write_density)
then
281 message(1) =
"Info: Starting calculation of unoccupied states."
285 eigens%converged(:) = 0
292 if (sys%st%pack_states .and. sys%hm%apply_packed())
call sys%st%pack()
294 do iter = 1, max_iter
295 output_iter = .false.
296 call eigens%run(sys%namespace, sys%gr, sys%st, sys%hm, 1, converged, sys%st%nst_conv)
311 write(iunit,
'(a)')
'All states converged.'
313 write(iunit,
'(a)')
'Some of the states are not fully converged!'
315 write(iunit,
'(a, e17.6)')
'Criterion = ', eigens%tolerance
321 forced_finish =
clean_stop(sys%mc%master_comm)
323 if (.not. bandstructure_mode)
then
325 if (converged .or. (modulo(iter, sys%outp%restart_write_interval) == 0) &
326 .or. iter == max_iter .or. forced_finish)
then
327 call states_elec_dump(restart_dump, sys%space, sys%st, sys%gr, sys%kpoints, ierr, iter=iter)
329 message(1) =
"Unable to write states wavefunctions."
335 do what_it = lbound(sys%outp%output_interval, 1), ubound(sys%outp%output_interval, 1)
336 if (sys%outp%what_now(what_it, iter))
then
342 if (output_iter .and. sys%outp%duringscf)
then
343 write(dirname,
'(a,i4.4)')
"unocc.",iter
344 call output_all(sys%outp, sys%namespace, sys%space, dirname, sys%gr, sys%ions, iter, sys%st, sys%hm, sys%ks)
345 call output_modelmb(sys%outp, sys%namespace, sys%space, dirname, sys%gr, sys%ions, iter, sys%st)
348 if (converged .or. forced_finish)
exit
352 if (.not. bandstructure_mode)
call restart_end(restart_dump)
354 if (sys%st%pack_states .and. sys%hm%apply_packed())
then
358 if (any(eigens%converged(:) < occ_states(:)))
then
359 write(
message(1),
'(a)')
'Some of the occupied states are not fully converged!'
363 safe_deallocate_a(occ_states)
365 if (.not. converged)
then
366 write(
message(1),
'(a)')
'Some of the unoccupied states are not fully converged!'
370 if (sys%space%is_periodic().and. sys%st%nik > sys%st%d%nspin)
then
373 sys%ions, sys%gr, sys%kpoints, &
374 sys%hm%phase, vec_pot = sys%hm%hm_base%uniform_vector_potential, &
375 vec_pot_var = sys%hm%hm_base%vector_potential)
380 call output_all(sys%outp, sys%namespace, sys%space,
static_dir, sys%gr, sys%ions, -1, sys%st, sys%hm, sys%ks)
389 subroutine init_(mesh, st)
390 class(
mesh_t),
intent(in) :: mesh
403 message(1) =
"With the RMMDIIS eigensolver for unocc, you will need to stop the calculation"
404 message(2) =
"by hand, since the highest states will probably never converge."
426 character(len=50) :: str
430 write(str,
'(a,i5)')
'Unoccupied states iteration #', iter
433 write(
message(1),
'(a,i6,a,i6)')
'Converged states: ', minval(eigens%converged(1:st%nik))
437 eigens%diff, st_start = showstart, compact = .
true., namespace=sys%namespace)
subroutine init_(fromscratch)
This module implements a calculator for the density and defines related functions.
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
subroutine, public eigensolver_init(eigens, namespace, gr, st, mc, space)
subroutine, public eigensolver_end(eigens)
integer, parameter, public rs_rmmdiis
character(len= *), parameter, public static_dir
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public hartree
integer, parameter, public hartree_fock
integer, parameter, public kohn_sham_dft
subroutine, public io_close(iunit, grp)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
integer, parameter, public kpoints_path
subroutine, public lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r)
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
integer, parameter, public dft_u_none
This module defines the meshes, which are used in Octopus.
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
subroutine, public messages_not_implemented(feature, namespace)
character(len=512), private msg
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
type(mpi_grp_t), public mpi_world
This module implements the basic mulsisystem class, a container system for other systems.
subroutine, public output_modelmb(outp, namespace, space, dir, gr, ions, iter, st)
this module contains the output system
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
logical function, public clean_stop(comm)
returns true if a file named stop exists
integer, parameter, public restart_gs
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
integer, parameter, public restart_type_dump
logical pure function, public restart_has_map(restart)
Returns true if the restart was from a different order of mesh points.
integer, parameter, public restart_td
integer, parameter, public restart_type_load
integer, parameter, public restart_unocc
logical pure function, public restart_are_basedirs_equal(type1, type2)
Returns true if...
subroutine, public restart_end(restart)
subroutine, public scf_state_info(namespace, st)
subroutine, public scf_print_mem_use(namespace)
This module defines routines to write information about states.
subroutine, public states_elec_write_eigenvalues(nst, st, space, kpoints, error, st_start, compact, iunit, namespace)
write the eigenvalues for some states to a file.
subroutine, public states_elec_write_bandstructure(dir, namespace, nst, st, ions, mesh, kpoints, phase, vec_pot, vec_pot_var)
calculate and write the bandstructure
subroutine, public states_elec_fermi(st, namespace, mesh, compute_spin)
calculate the Fermi level for the states in this object
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public occupied_states(st, namespace, ik, n_filled, n_partially_filled, n_half_filled, filled, partially_filled, half_filled)
return information about occupied orbitals in many-body state
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_dump(restart, space, st, mesh, kpoints, ierr, iter, lr, st_start_writing, verbose)
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, ierr, iter, lr, lowest_missing, label, verbose, skip)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
subroutine, public states_elec_load_rho(restart, space, st, mesh, ierr)
subroutine, public states_elec_dump_rho(restart, space, st, mesh, ierr, iter)
subroutine, public unocc_run(system, from_scratch)
subroutine unocc_run_legacy(sys, fromscratch)
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
logical pure function, public xc_is_orbital_dependent(xcs)
Is the xc family orbital dependent.
Class describing the electron system.
Describes mesh distribution to nodes.
Container class for lists of system_oct_m::system_t.
The states_elec_t class contains all electronic wave functions.
subroutine write_iter_(namespace, st)