Octopus
unocc.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module unocc_oct_m
22 use debug_oct_m
26 use global_oct_m
27 use output_oct_m
30 use io_oct_m
32 use lcao_oct_m
33 use lda_u_oct_m
35 use mesh_oct_m
37 use mpi_oct_m
40 use parser_oct_m
43 use scf_oct_m
44 use space_oct_m
50 use v_ks_oct_m
51 use xc_oct_m
52
53 implicit none
54
55 private
56 public :: &
58
59
60contains
61
62 ! ---------------------------------------------------------
63 subroutine unocc_run(system, from_scratch)
64 class(*), intent(inout) :: system
65 logical, intent(in) :: from_scratch
66
67 push_sub(unocc_run)
68
69 select type (system)
70 class is (multisystem_basic_t)
71 message(1) = "CalculationMode = unocc not implemented for multi-system calculations"
72 call messages_fatal(1, namespace=system%namespace)
73 type is (electrons_t)
74 call unocc_run_legacy(system, from_scratch)
75 end select
76
77 pop_sub(unocc_run)
78 end subroutine unocc_run
79
80 ! ---------------------------------------------------------
81 subroutine unocc_run_legacy(sys, fromscratch)
82 type(electrons_t), intent(inout) :: sys
83 logical, intent(in) :: fromscratch
84
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
95
96 push_sub(unocc_run_legacy)
97
98 if (sys%hm%pcm%run_pcm) then
99 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
100 end if
101
102 ! This variable is documented in scf.F90
103 call parse_variable(sys%namespace, 'MaximumIter', 50, max_iter)
104 call messages_obsolete_variable(sys%namespace, 'UnoccMaximumIter', 'MaximumIter')
105 if (max_iter < 0) max_iter = huge(max_iter)
106
107 !%Variable UnoccShowOccStates
108 !%Type logical
109 !%Default false
110 !%Section Calculation Modes::Unoccupied States
111 !%Description
112 !% If true, the convergence for the occupied states will be shown too in the output.
113 !% This is useful for testing, or if the occupied states fail to converge.
114 !% It will be enabled automatically if only occupied states are being calculated.
115 !%End
116 call parse_variable(sys%namespace, 'UnoccShowOccStates', .false., showoccstates)
117
118 bandstructure_mode = .false.
119 if (sys%space%is_periodic() .and. sys%kpoints%get_kpoint_method() == kpoints_path) then
120 bandstructure_mode = .true.
121 end if
122
123 call init_(sys%gr, sys%st)
124 converged = .false.
125
126 read_td_states = .false.
127 if (bandstructure_mode) then
128 !%Variable UnoccUseTD
129 !%Type logical
130 !%Default no
131 !%Section Calculation Modes::Unoccupied States
132 !%Description
133 !% If true, Octopus will use the density and states from the restart/td folder to compute
134 !% the bandstructure, instead of the restart/gs ones.
135 !%End
136 call parse_variable(sys%namespace, 'UnoccUseTD', .false., read_td_states)
137 end if
138
139
140 safe_allocate(lowest_missing(1:sys%st%d%dim, 1:sys%st%nik))
141 ! if there is no restart info to read, this will not get set otherwise
142 ! setting to zero means everything is missing.
143 lowest_missing(:,:) = 0
144
145 read_gs = .true.
146 if (.not. fromscratch) then
147 call restart_init(restart_load_unocc, sys%namespace, restart_unocc, restart_type_load, sys%mc, ierr, &
148 mesh = sys%gr, exact = .true.)
149
150 if (ierr == 0) then
151 call states_elec_load(restart_load_unocc, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
152 ierr, lowest_missing = lowest_missing)
153 call restart_end(restart_load_unocc)
154 end if
155
156 ! If successfully read states from unocc, do not read from gs.
157 ! If RESTART_GS and RESTART_UNOCC have the same directory (the default), and we tried RESTART_UNOCC
158 ! already and failed, it is a waste of time to try to read again.
159 if (ierr == 0 .or. restart_are_basedirs_equal(restart_gs, restart_unocc)) read_gs = .false.
160 end if
161
162 if (read_td_states) then
163 call restart_init(restart_load_gs, sys%namespace, restart_td, restart_type_load, sys%mc, ierr_rho, mesh=sys%gr, &
164 exact=.true.)
165 else
166 call restart_init(restart_load_gs, sys%namespace, restart_gs, restart_type_load, sys%mc, ierr_rho, mesh=sys%gr, &
167 exact=.true.)
168 end if
169
170 if (ierr_rho == 0) then
171 if (read_gs) 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)
174 end if
175 if (sys%hm%lda_u_level /= dft_u_none) then
176 call lda_u_load(restart_load_gs, sys%hm%lda_u, sys%st, sys%hm%energy%dft_u, ierr)
177 if (ierr /= 0) then
178 message(1) = "Unable to read DFT+U information. DFT+U data will be calculated from states."
179 call messages_warning(1, namespace=sys%namespace)
180 end if
181 end if
182
183 call states_elec_load_rho(restart_load_gs, sys%space, sys%st, sys%gr, ierr_rho)
184 write_density = restart_has_map(restart_load_gs)
185 call restart_end(restart_load_gs)
186 else
187 write_density = .true.
188 end if
189
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
194 end do
195
196 is_orbital_dependent = (sys%ks%theory_level == hartree .or. sys%ks%theory_level == hartree_fock .or. &
197 (sys%ks%theory_level == kohn_sham_dft .and. xc_is_orbital_dependent(sys%ks%xc)) &
198 .or. (sys%ks%theory_level == generalized_kohn_sham_dft .and. xc_is_orbital_dependent(sys%ks%xc)))
199
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."
204 call messages_warning(3, namespace=sys%namespace)
205 end if
206
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
211 occ_missing = .true.
212 end if
213 end do
214
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."
220 end if
221
222 message(2) = "Not all the occupied orbitals could be read."
223 message(3) = "Please run a ground-state calculation first!"
224 call messages_fatal(3, only_root_writes = .true., namespace=sys%namespace)
225 end if
226
227 message(1) = "Unable to read density: Building density from wavefunctions."
228 call messages_info(1, namespace=sys%namespace)
229
230 call density_calc(sys%st, sys%gr, sys%st%rho)
231 end if
232
233 call scf_state_info(sys%namespace, sys%st)
234
235 if (fromscratch .or. ierr /= 0) then
236 if (fromscratch) then
237 ! do not use previously calculated occupied states
238 nst_calculated = min(maxval(occ_states), minval(lowest_missing) - 1)
239 else
240 ! or, use as many states as have been calculated
241 nst_calculated = minval(lowest_missing) - 1
242 end if
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)
246 else
247 ! we successfully read all the states and are planning to use them, no need for LCAO
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
251 end if
252
253
254
255 ! it is strange and useless to see no eigenvalues written if you are only calculating
256 ! occupied states, on a different k-point.
257 if (showstart > sys%st%nst) showstart = 1
258
259 safe_deallocate_a(lowest_missing)
260
261 if (showoccstates) showstart = 1
262
263 ! In the case of someone using KPointsPath, the code assume that this is only for plotting a
264 ! bandstructure. This mode ensure that no restart information will be written for the new grid
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."
268 call messages_info(2, namespace=sys%namespace)
269 end if
270
271 if (.not. bandstructure_mode) then
272 ! Restart dump should be initialized after restart_load, as the mesh might have changed
273 call restart_init(restart_dump, sys%namespace, restart_unocc, restart_type_dump, sys%mc, ierr, mesh=sys%gr)
274
275 ! make sure the density is defined on the same mesh as the wavefunctions that will be written
276 if (write_density) then
277 call states_elec_dump_rho(restart_dump, sys%space, sys%st, sys%gr, ierr_rho)
278 end if
279 end if
280
281 message(1) = "Info: Starting calculation of unoccupied states."
282 call messages_info(1, namespace=sys%namespace)
283
284 ! reset this variable, so that the eigensolver passes through all states
285 eigens%converged(:) = 0
286
287 ! If not all gs wavefunctions were read when starting, in particular for nscf with different k-points,
288 ! the occupations must be recalculated each time, though they do not affect the result of course.
289 ! FIXME: This is wrong for metals where we must use the Fermi level from the original calculation!
290 call states_elec_fermi(sys%st, sys%namespace, sys%gr)
291
292 if (sys%st%pack_states .and. sys%hm%apply_packed()) call sys%st%pack()
293
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)
297
298 ! If not all gs wavefunctions were read when starting, in particular for nscf with different k-points,
299 ! the occupations must be recalculated each time, though they do not affect the result of course.
300 ! FIXME: This is wrong for metals where we must use the Fermi level from the original calculation!
301 call states_elec_fermi(sys%st, sys%namespace, sys%gr)
302
303 call write_iter_(sys%namespace, sys%st)
304
305 ! write output file
306 if (mpi_grp_is_root(mpi_world)) then
307 call io_mkdir(static_dir, sys%namespace)
308 iunit = io_open(static_dir//'/eigenvalues', sys%namespace, action='write')
309
310 if (converged) then
311 write(iunit,'(a)') 'All states converged.'
312 else
313 write(iunit,'(a)') 'Some of the states are not fully converged!'
314 end if
315 write(iunit,'(a, e17.6)') 'Criterion = ', eigens%tolerance
316 write(iunit,'(1x)')
317 call states_elec_write_eigenvalues(sys%st%nst, sys%st, sys%space, sys%kpoints, eigens%diff, iunit=iunit)
318 call io_close(iunit)
319 end if
320
321 forced_finish = clean_stop(sys%mc%master_comm)
322
323 if (.not. bandstructure_mode) then
324 ! write restart information.
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)
328 if (ierr /= 0) then
329 message(1) = "Unable to write states wavefunctions."
330 call messages_warning(1, namespace=sys%namespace)
331 end if
332 end if
333 end if
334
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
337 output_iter = .true.
338 exit
339 end if
340 end do
341
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)
346 end if
347
348 if (converged .or. forced_finish) exit
349
350 end do
351
352 if (.not. bandstructure_mode) call restart_end(restart_dump)
353
354 if (sys%st%pack_states .and. sys%hm%apply_packed()) then
355 call sys%st%unpack()
356 end if
357
358 if (any(eigens%converged(:) < occ_states(:))) then
359 write(message(1),'(a)') 'Some of the occupied states are not fully converged!'
360 call messages_warning(1, namespace=sys%namespace)
361 end if
362
363 safe_deallocate_a(occ_states)
364
365 if (.not. converged) then
366 write(message(1),'(a)') 'Some of the unoccupied states are not fully converged!'
367 call messages_warning(1, namespace=sys%namespace)
368 end if
369
370 if (sys%space%is_periodic().and. sys%st%nik > sys%st%d%nspin) then
371 if (bitand(sys%kpoints%method, kpoints_path) /= 0) then
372 call states_elec_write_bandstructure(static_dir, sys%namespace, sys%st%nst, sys%st, &
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)
376 end if
377 end if
378
379
380 call output_all(sys%outp, sys%namespace, sys%space, static_dir, sys%gr, sys%ions, -1, sys%st, sys%hm, sys%ks)
381 call output_modelmb(sys%outp, sys%namespace, sys%space, static_dir, sys%gr, sys%ions, -1, sys%st)
382
383 call end_()
384 pop_sub(unocc_run_legacy)
385
386 contains
387
388 ! ---------------------------------------------------------
389 subroutine init_(mesh, st)
390 class(mesh_t), intent(in) :: mesh
391 type(states_elec_t), intent(inout) :: st
392
393 push_sub(unocc_run_legacy.init_)
394
395 call messages_obsolete_variable(sys%namespace, "NumberUnoccStates", "ExtraStates")
396
397 call states_elec_allocate_wfns(st, mesh, packed=.true.)
398
399 ! now the eigensolver stuff
400 call eigensolver_init(eigens, sys%namespace, sys%gr, st, sys%mc, sys%space)
401
402 if (eigens%es_type == rs_rmmdiis) then
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."
405 call messages_warning(2, namespace=sys%namespace)
406 end if
407
408 pop_sub(unocc_run_legacy.init_)
409 end subroutine init_
410
411
412 ! ---------------------------------------------------------
413 subroutine end_()
414 push_sub(unocc_run_legacy.end_)
415
416 call eigensolver_end(eigens)
417
418 pop_sub(unocc_run_legacy.end_)
419 end subroutine end_
420
421 ! ---------------------------------------------------------
422 subroutine write_iter_(namespace, st)
423 type(namespace_t), intent(in) :: namespace
424 type(states_elec_t), intent(in) :: st
425
426 character(len=50) :: str
427
429
430 write(str, '(a,i5)') 'Unoccupied states iteration #', iter
431 call messages_print_with_emphasis(msg=trim(str), namespace=sys%namespace)
432
433 write(message(1),'(a,i6,a,i6)') 'Converged states: ', minval(eigens%converged(1:st%nik))
434 call messages_info(1, namespace=namespace)
435
436 call states_elec_write_eigenvalues(sys%st%nst, sys%st, sys%space, sys%kpoints, &
437 eigens%diff, st_start = showstart, compact = .true., namespace=sys%namespace)
438
439 call scf_print_mem_use(namespace)
440
441 call messages_print_with_emphasis(namespace=sys%namespace)
442
444 end subroutine write_iter_
445
446
447 end subroutine unocc_run_legacy
448
449
450end module unocc_oct_m
451
452!! Local Variables:
453!! mode: f90
454!! coding: utf-8
455!! End:
subroutine init_(fromscratch)
Definition: geom_opt.F90:354
subroutine end_()
Definition: geom_opt.F90:804
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:608
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
Definition: global.F90:251
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public hartree
integer, parameter, public hartree_fock
integer, parameter, public kohn_sham_dft
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:354
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
integer, parameter, public kpoints_path
Definition: kpoints.F90:209
subroutine, public lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r)
Definition: lcao.F90:799
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
Definition: lda_u_io.F90:731
integer, parameter, public dft_u_none
Definition: lda_u.F90:200
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:930
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
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
Definition: output.F90:115
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
Definition: output.F90:491
logical function, public clean_stop(comm)
returns true if a file named stop exists
Definition: restart.F90:274
integer, parameter, public restart_gs
Definition: restart.F90:229
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:514
integer, parameter, public restart_type_dump
Definition: restart.F90:225
logical pure function, public restart_has_map(restart)
Returns true if the restart was from a different order of mesh points.
Definition: restart.F90:988
integer, parameter, public restart_td
Definition: restart.F90:229
integer, parameter, public restart_type_load
Definition: restart.F90:225
integer, parameter, public restart_unocc
Definition: restart.F90:229
logical pure function, public restart_are_basedirs_equal(type1, type2)
Returns true if...
Definition: restart.F90:998
subroutine, public restart_end(restart)
Definition: restart.F90:720
subroutine, public scf_state_info(namespace, st)
Definition: scf.F90:1428
subroutine, public scf_print_mem_use(namespace)
Definition: scf.F90:1446
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)
Definition: unocc.F90:157
subroutine unocc_run_legacy(sys, fromscratch)
Definition: unocc.F90:175
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:732
Definition: xc.F90:114
logical pure function, public xc_is_orbital_dependent(xcs)
Is the xc family orbital dependent.
Definition: xc.F90:541
Class describing the electron system.
Definition: electrons.F90:214
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Container class for lists of system_oct_m::system_t.
The states_elec_t class contains all electronic wave functions.
int true(void)
subroutine write_iter_(namespace, st)
Definition: unocc.F90:516