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
23 use debug_oct_m
27 use global_oct_m
28 use output_oct_m
31 use io_oct_m
34 use lcao_oct_m
35 use lda_u_oct_m
37 use mesh_oct_m
39 use mpi_oct_m
42 use parser_oct_m
45 use scf_oct_m
46 use space_oct_m
52 use v_ks_oct_m
54 use xc_oct_m
55
56 implicit none
57
58 private
59 public :: &
61
62
63contains
64
65 ! ---------------------------------------------------------
66 subroutine unocc_run(system, from_scratch)
67 class(*), intent(inout) :: system
68 logical, intent(in) :: from_scratch
69
70 push_sub(unocc_run)
71
72 select type (system)
73 class is (multisystem_basic_t)
74 message(1) = "CalculationMode = unocc not implemented for multi-system calculations"
75 call messages_fatal(1, namespace=system%namespace)
76 type is (electrons_t)
77 call unocc_run_legacy(system, from_scratch)
78 end select
79
80 pop_sub(unocc_run)
81 end subroutine unocc_run
82
83 ! ---------------------------------------------------------
84 subroutine unocc_run_legacy(sys, fromscratch)
85 type(electrons_t), intent(inout) :: sys
86 logical, intent(in) :: fromscratch
87
88 type(eigensolver_t) :: eigens
89 integer :: iunit, ierr, iter, ierr_rho, ik
90 integer(int64) :: what_it
91 logical :: read_gs, converged, forced_finish, showoccstates, is_orbital_dependent, occ_missing
92 integer :: max_iter, nst_calculated, showstart
93 integer :: n_filled, n_partially_filled, n_half_filled
94 integer, allocatable :: lowest_missing(:, :), occ_states(:)
95 character(len=10) :: dirname
96 type(restart_t) :: restart_load_unocc, restart_load_gs, restart_dump
97 logical :: write_density, bandstructure_mode, read_td_states, output_iter, known_lower_bound
98 type(current_t) :: current_calculator
99
100 push_sub(unocc_run_legacy)
101
102 if (sys%hm%pcm%run_pcm) then
103 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
104 end if
105
106 ! This variable is documented in scf.F90
107 call parse_variable(sys%namespace, 'MaximumIter', 50, max_iter)
108 call messages_obsolete_variable(sys%namespace, 'UnoccMaximumIter', 'MaximumIter')
109 if (max_iter < 0) max_iter = huge(max_iter)
110
111 !%Variable UnoccShowOccStates
112 !%Type logical
113 !%Default false
114 !%Section Calculation Modes::Unoccupied States
115 !%Description
116 !% If true, the convergence for the occupied states will be shown too in the output.
117 !% This is useful for testing, or if the occupied states fail to converge.
118 !% It will be enabled automatically if only occupied states are being calculated.
119 !%End
120 call parse_variable(sys%namespace, 'UnoccShowOccStates', .false., showoccstates)
121
122 bandstructure_mode = .false.
123 if (sys%space%is_periodic() .and. sys%kpoints%get_kpoint_method() == kpoints_path) then
124 bandstructure_mode = .true.
125 end if
126
127 call init_(sys%gr, sys%st)
128 converged = .false.
129
130 read_td_states = .false.
131 if (bandstructure_mode) then
132 !%Variable UnoccUseTD
133 !%Type logical
134 !%Default no
135 !%Section Calculation Modes::Unoccupied States
136 !%Description
137 !% If true, Octopus will use the density and states from the restart/td folder to compute
138 !% the bandstructure, instead of the restart/gs ones.
139 !%End
140 call parse_variable(sys%namespace, 'UnoccUseTD', .false., read_td_states)
141 end if
142
143
144 safe_allocate(lowest_missing(1:sys%st%d%dim, 1:sys%st%nik))
145 ! if there is no restart info to read, this will not get set otherwise
146 ! setting to zero means everything is missing.
147 lowest_missing(:,:) = 0
148
149 read_gs = .true.
150 if (.not. fromscratch) then
151 call restart_load_unocc%init(sys%namespace, restart_unocc, restart_type_load, sys%mc, ierr, &
152 mesh = sys%gr, exact = .true.)
153
154 if (ierr == 0) then
155 call states_elec_load(restart_load_unocc, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
156 ierr, lowest_missing = lowest_missing)
157 call restart_load_unocc%end()
158 end if
159
160 ! If successfully read states from unocc, do not read from gs.
161 ! If RESTART_GS and RESTART_UNOCC have the same directory (the default), and we tried RESTART_UNOCC
162 ! already and failed, it is a waste of time to try to read again.
163 if (ierr == 0 .or. restart_load_gs%basedir == restart_load_unocc%basedir) read_gs = .false.
164 end if
165
166 if (read_td_states) then
167 call restart_load_gs%init(sys%namespace, restart_td, restart_type_load, sys%mc, ierr_rho, mesh=sys%gr, &
168 exact=.true.)
169 else
170 call restart_load_gs%init(sys%namespace, restart_gs, restart_type_load, sys%mc, ierr_rho, mesh=sys%gr, &
171 exact=.true.)
172 end if
173
174 if (ierr_rho == 0) then
175 if (read_gs) then
176 call states_elec_load(restart_load_gs, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
177 ierr, lowest_missing = lowest_missing)
178 end if
179 if (sys%hm%lda_u_level /= dft_u_none) then
180 call lda_u_load(restart_load_gs, sys%hm%lda_u, sys%st, sys%hm%energy%dft_u, ierr)
181 if (ierr /= 0) then
182 message(1) = "Unable to read DFT+U information. DFT+U data will be calculated from states."
183 call messages_warning(1, namespace=sys%namespace)
184 end if
185 end if
186
187 call states_elec_load_rho(restart_load_gs, sys%st, sys%gr, ierr_rho)
188 write_density = restart_load_gs%has_map()
189 call restart_load_gs%end()
190 else
191 write_density = .true.
192 end if
193
194 safe_allocate(occ_states(1:sys%st%nik))
195 do ik = 1, sys%st%nik
196 call occupied_states(sys%st, sys%namespace, ik, n_filled, n_partially_filled, n_half_filled)
197 occ_states(ik) = n_filled + n_partially_filled + n_half_filled
198 end do
199
200 is_orbital_dependent = (sys%ks%theory_level == hartree .or. sys%ks%theory_level == hartree_fock .or. &
201 (sys%ks%theory_level == kohn_sham_dft .and. xc_is_orbital_dependent(sys%ks%xc)) &
202 .or. (sys%ks%theory_level == generalized_kohn_sham_dft .and. xc_is_orbital_dependent(sys%ks%xc)))
203
204 if (is_orbital_dependent) then
205 message(1) = "Be sure your gs run is well converged since you have an orbital-dependent functional."
206 message(2) = "Otherwise, the occupied states may change in CalculationMode = unocc, and your"
207 message(3) = "unoccupied states will not be consistent with the gs run."
208 call messages_warning(3, namespace=sys%namespace)
209 end if
210
211 if (ierr_rho /= 0 .or. is_orbital_dependent) then
212 occ_missing = .false.
213 do ik = 1, sys%st%nik
214 if (any(lowest_missing(1:sys%st%d%dim, ik) <= occ_states(ik))) then
215 occ_missing = .true.
216 end if
217 end do
218
219 if (occ_missing) then
220 if (is_orbital_dependent) then
221 message(1) = "For an orbital-dependent functional, all occupied orbitals must be provided."
222 else if (ierr_rho /= 0) then
223 message(1) = "Since density could not be read, all occupied orbitals must be provided."
224 end if
225
226 message(2) = "Not all the occupied orbitals could be read."
227 message(3) = "Please run a ground-state calculation first!"
228 call messages_fatal(3, only_root_writes = .true., namespace=sys%namespace)
229 end if
230
231 message(1) = "Unable to read density: Building density from wavefunctions."
232 call messages_info(1, namespace=sys%namespace)
233
234 call density_calc(sys%st, sys%gr, sys%st%rho)
235 end if
236
237 call scf_state_info(sys%namespace, sys%st)
238
239 if (fromscratch .or. ierr /= 0) then
240 if (fromscratch) then
241 ! do not use previously calculated occupied states
242 nst_calculated = min(maxval(occ_states), minval(lowest_missing) - 1)
243 else
244 ! or, use as many states as have been calculated
245 nst_calculated = minval(lowest_missing) - 1
246 end if
247 showstart = max(nst_calculated + 1, 1)
248 call lcao_run(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, &
249 sys%hm, st_start = showstart, known_lower_bound=known_lower_bound)
250 call eigens%set_lower_bound_is_known(known_lower_bound)
251 else
252 ! we successfully read all the states and are planning to use them, no need for LCAO
253 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, sys%st, sys%ions, sys%ext_partners, &
254 calc_eigenval = .false.)
255 showstart = minval(occ_states(:)) + 1
256
257 call eigens%set_lower_bound_is_known(.true.) ! Lower bound is known here
258 end if
259
260
261
262 ! it is strange and useless to see no eigenvalues written if you are only calculating
263 ! occupied states, on a different k-point.
264 if (showstart > sys%st%nst) showstart = 1
265
266 safe_deallocate_a(lowest_missing)
267
268 if (showoccstates) showstart = 1
269
270 ! In the case of someone using KPointsPath, the code assume that this is only for plotting a
271 ! bandstructure. This mode ensure that no restart information will be written for the new grid
272 if (bandstructure_mode) then
273 message(1) = "Info: The code will run in band structure mode."
274 message(2) = " No restart information will be printed."
275 call messages_info(2, namespace=sys%namespace)
276 end if
277
278 if (.not. bandstructure_mode) then
279 ! Restart dump should be initialized after restart_load, as the mesh might have changed
280 call restart_dump%init(sys%namespace, restart_unocc, restart_type_dump, sys%mc, ierr, mesh=sys%gr)
281
282 ! make sure the density is defined on the same mesh as the wavefunctions that will be written
283 if (write_density) then
284 call states_elec_dump_rho(restart_dump, sys%st, sys%gr, ierr_rho)
285 end if
286 end if
287
288 message(1) = "Info: Starting calculation of unoccupied states."
289 call messages_info(1, namespace=sys%namespace)
290
291 ! reset this variable, so that the eigensolver passes through all states
292 eigens%converged(:) = 0
293
294 ! If not all gs wavefunctions were read when starting, in particular for nscf with different k-points,
295 ! the occupations must be recalculated each time, though they do not affect the result of course.
296 ! FIXME: This is wrong for metals where we must use the Fermi level from the original calculation!
297 call states_elec_fermi(sys%st, sys%namespace, sys%gr)
298
299 if(output_needs_current(sys%outp, states_are_real(sys%st))) then
300 call current_init(current_calculator, sys%namespace)
301 end if
302
303 if (sys%st%pack_states .and. sys%hm%apply_packed()) call sys%st%pack()
304
305 do iter = 1, max_iter
306 output_iter = .false.
307 call eigens%run(sys%namespace, sys%gr, sys%st, sys%hm, sys%space, sys%ext_partners, 1, converged, sys%st%nst_conv)
308
309 ! If not all gs wavefunctions were read when starting, in particular for nscf with different k-points,
310 ! the occupations must be recalculated each time, though they do not affect the result of course.
311 ! FIXME: This is wrong for metals where we must use the Fermi level from the original calculation!
312 call states_elec_fermi(sys%st, sys%namespace, sys%gr)
313
314 call write_iter_(sys%namespace, sys%st)
315
316 ! write output file
317 if (mpi_world%is_root()) then
318 call io_mkdir(static_dir, sys%namespace)
319 iunit = io_open(static_dir//'/eigenvalues', sys%namespace, action='write')
320
321 if (converged) then
322 write(iunit,'(a)') 'All states converged.'
323 else
324 write(iunit,'(a)') 'Some of the states are not fully converged!'
325 end if
326 write(iunit,'(a, e17.6)') 'Criterion = ', eigens%tolerance
327 write(iunit,'(1x)')
328 call states_elec_write_eigenvalues(sys%st%nst, sys%st, sys%space, sys%kpoints, eigens%diff, iunit=iunit)
329 call io_close(iunit)
330 end if
331
332 forced_finish = clean_stop(sys%mc%master_comm)
333
334 if (.not. bandstructure_mode) then
335 ! write restart information.
336 if (converged .or. restart_walltime_period_alarm(sys%mc%master_comm) &
337 .or. iter == max_iter .or. forced_finish) then
338 call states_elec_dump(restart_dump, sys%space, sys%st, sys%gr, sys%kpoints, ierr, iter=iter)
339 if (ierr /= 0) then
340 message(1) = "Unable to write states wavefunctions."
341 call messages_warning(1, namespace=sys%namespace)
342 end if
343 end if
344 end if
345
346 do what_it = lbound(sys%outp%output_interval, 1), ubound(sys%outp%output_interval, 1)
347 if (sys%outp%what_now(what_it, iter)) then
348 output_iter = .true.
349 exit
350 end if
351 end do
352
353 if (output_iter .and. sys%outp%duringscf) then
354 write(dirname,'(a,i4.4)') "unocc.",iter
355 if(output_needs_current(sys%outp, states_are_real(sys%st))) then
356 call states_elec_allocate_current(sys%st, sys%space, sys%gr)
357 call current_calculate(current_calculator, sys%namespace, sys%gr, sys%hm, sys%space, sys%st)
358 end if
359 call output_all(sys%outp, sys%namespace, sys%space, dirname, sys%gr, sys%ions, iter, sys%st, sys%hm, sys%ks)
360 call output_modelmb(sys%outp, sys%namespace, sys%space, dirname, sys%gr, sys%ions, iter, sys%st)
361 end if
362
363 if (converged .or. forced_finish) exit
364
365 end do
366
367 if (.not. bandstructure_mode) call restart_dump%end()
368
369 if (sys%st%pack_states .and. sys%hm%apply_packed()) then
370 call sys%st%unpack()
371 end if
372
373 if (any(eigens%converged(:) < occ_states(:))) then
374 write(message(1),'(a)') 'Some of the occupied states are not fully converged!'
375 call messages_warning(1, namespace=sys%namespace)
376 end if
377
378 safe_deallocate_a(occ_states)
379
380 if (.not. converged) then
381 write(message(1),'(a)') 'Some of the unoccupied states are not fully converged!'
382 call messages_warning(1, namespace=sys%namespace)
383 end if
384
385 if (sys%space%is_periodic().and. sys%st%nik > sys%st%d%nspin) then
386 if (bitand(sys%kpoints%method, kpoints_path) /= 0) then
387 call states_elec_write_bandstructure(static_dir, sys%namespace, sys%st%nst, sys%st, &
388 sys%ions, sys%gr, sys%kpoints, &
389 sys%hm%phase, vec_pot = sys%hm%hm_base%uniform_vector_potential, &
390 vec_pot_var = sys%hm%hm_base%vector_potential)
391 end if
392 end if
393
394 if(output_needs_current(sys%outp, states_are_real(sys%st))) then
395 call states_elec_allocate_current(sys%st, sys%space, sys%gr)
396 call current_calculate(current_calculator, sys%namespace, sys%gr, sys%hm, sys%space, sys%st)
397 end if
398 call output_all(sys%outp, sys%namespace, sys%space, static_dir, sys%gr, sys%ions, -1, sys%st, sys%hm, sys%ks)
399 call output_modelmb(sys%outp, sys%namespace, sys%space, static_dir, sys%gr, sys%ions, -1, sys%st)
400
401 call end_()
402 pop_sub(unocc_run_legacy)
403
404 contains
405
406 ! ---------------------------------------------------------
407 subroutine init_(mesh, st)
408 class(mesh_t), intent(in) :: mesh
409 type(states_elec_t), intent(inout) :: st
410
411 push_sub(unocc_run_legacy.init_)
412
413 call messages_obsolete_variable(sys%namespace, "NumberUnoccStates", "ExtraStates")
414
415 call states_elec_allocate_wfns(st, mesh, packed=.true.)
416
417 ! now the eigensolver stuff
418 call eigensolver_init(eigens, sys%namespace, sys%gr, st, sys%hm, sys%mc, sys%space, deactivate_oracle=.true.)
419
420 if (eigens%es_type == rs_rmmdiis) then
421 message(1) = "With the RMMDIIS eigensolver for unocc, you will need to stop the calculation"
422 message(2) = "by hand, since the highest states will probably never converge."
423 call messages_warning(2, namespace=sys%namespace)
424 end if
425
426 pop_sub(unocc_run_legacy.init_)
427 end subroutine init_
428
429
430 ! ---------------------------------------------------------
431 subroutine end_()
432 push_sub(unocc_run_legacy.end_)
433
434 call eigensolver_end(eigens)
435
436 pop_sub(unocc_run_legacy.end_)
437 end subroutine end_
438
439 ! ---------------------------------------------------------
440 subroutine write_iter_(namespace, st)
441 type(namespace_t), intent(in) :: namespace
442 type(states_elec_t), intent(in) :: st
443
444 character(len=50) :: str
445
447
448 write(str, '(a,i5)') 'Unoccupied states iteration #', iter
449 call messages_print_with_emphasis(msg=trim(str), namespace=sys%namespace)
450
451 write(message(1),'(a,i6,a,i6)') 'Converged eigenvectors: ', minval(eigens%converged(1:st%nik))
452 call messages_info(1, namespace=namespace)
453
454 call states_elec_write_eigenvalues(sys%st%nst, sys%st, sys%space, sys%kpoints, &
455 eigens%diff, st_start = showstart, compact = .true., namespace=sys%namespace)
456
457 call scf_print_mem_use(namespace)
458
459 call messages_print_with_emphasis(namespace=sys%namespace)
460
462 end subroutine write_iter_
463
464
465 end subroutine unocc_run_legacy
466
467
468end module unocc_oct_m
469
470!! Local Variables:
471!! mode: f90
472!! coding: utf-8
473!! End:
subroutine init_(fromscratch)
Definition: geom_opt.F90:364
subroutine end_()
Definition: geom_opt.F90:825
subroutine, public current_calculate(this, namespace, gr, hm, space, st)
Compute total electronic current density.
Definition: current.F90:371
subroutine, public current_init(this, namespace)
Definition: current.F90:181
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:612
subroutine, public eigensolver_init(eigens, namespace, gr, st, hm, mc, space, deactivate_oracle)
subroutine, public eigensolver_end(eigens)
integer, parameter, public rs_rmmdiis
integer, parameter, public hartree_fock
Definition: global.F90:237
integer, parameter, public generalized_kohn_sham_dft
Definition: global.F90:237
integer, parameter, public kohn_sham_dft
Definition: global.F90:237
character(len= *), parameter, public static_dir
Definition: global.F90:266
integer, parameter, public hartree
Definition: global.F90:237
Definition: io.F90:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:467
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:361
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
integer, parameter, public kpoints_path
Definition: kpoints.F90:211
A module to handle KS potential, without the external potential.
subroutine, public lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r, known_lower_bound)
Definition: lcao.F90:770
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
Definition: lda_u_io.F90:729
integer, parameter, public dft_u_none
Definition: lda_u.F90:203
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1091
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1023
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
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:117
logical function, public output_needs_current(outp, states_are_real)
Definition: output.F90:983
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
Definition: output.F90:497
logical function, public clean_stop(comm)
returns true if a file named stop exists
Definition: restart.F90:335
integer, parameter, public restart_gs
Definition: restart.F90:156
integer, parameter, public restart_type_dump
Definition: restart.F90:183
integer, parameter, public restart_td
Definition: restart.F90:156
integer, parameter, public restart_type_load
Definition: restart.F90:183
integer, parameter, public restart_unocc
Definition: restart.F90:156
subroutine, public scf_state_info(namespace, st)
Definition: scf.F90:1507
subroutine, public scf_print_mem_use(namespace)
Definition: scf.F90:1525
pure logical function, public states_are_real(st)
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
subroutine, public states_elec_allocate_current(st, space, mesh)
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, 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, st, mesh, ierr)
subroutine, public states_elec_dump_rho(restart, st, mesh, ierr, iter)
subroutine, public unocc_run(system, from_scratch)
Definition: unocc.F90:162
subroutine unocc_run_legacy(sys, fromscratch)
Definition: unocc.F90:180
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:752
This module provices a simple timer class which can be used to trigger the writing of a restart file ...
Definition: walltimer.F90:123
logical function, public restart_walltime_period_alarm(comm)
Definition: walltimer.F90:375
Definition: xc.F90:116
logical pure function, public xc_is_orbital_dependent(xcs)
Is the xc family orbital dependent.
Definition: xc.F90:542
Class describing the electron system.
Definition: electrons.F90:220
Describes mesh distribution to nodes.
Definition: mesh.F90:187
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:536