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 ierr = 0
151 if (.not. fromscratch) then
152 call restart_load_unocc%init(sys%namespace, restart_unocc, restart_type_load, sys%mc, ierr, &
153 mesh = sys%gr, exact = .true.)
154
155 if (ierr == 0) then
156 call states_elec_load(restart_load_unocc, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
157 ierr, lowest_missing = lowest_missing)
158 call restart_load_unocc%end()
159 end if
160
161 ! If successfully read states from unocc, do not read from gs.
162 ! If RESTART_GS and RESTART_UNOCC have the same directory (the default), and we tried RESTART_UNOCC
163 ! already and failed, it is a waste of time to try to read again.
164 if (ierr == 0 .or. restart_load_gs%basedir == restart_load_unocc%basedir) read_gs = .false.
165 end if
166
167 if (read_td_states) then
168 call restart_load_gs%init(sys%namespace, restart_td, restart_type_load, sys%mc, ierr_rho, mesh=sys%gr, &
169 exact=.true.)
170 else
171 call restart_load_gs%init(sys%namespace, restart_gs, restart_type_load, sys%mc, ierr_rho, mesh=sys%gr, &
172 exact=.true.)
173 end if
174
175 if (ierr_rho == 0) then
176 if (read_gs) then
177 call states_elec_load(restart_load_gs, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
178 ierr, lowest_missing = lowest_missing)
179 end if
180 if (sys%hm%lda_u_level /= dft_u_none) then
181 call lda_u_load(restart_load_gs, sys%hm%lda_u, sys%st, sys%hm%energy%dft_u, ierr)
182 if (ierr /= 0) then
183 message(1) = "Unable to read DFT+U information. DFT+U data will be calculated from states."
184 call messages_warning(1, namespace=sys%namespace)
185 end if
186 end if
187
188 call states_elec_load_rho(restart_load_gs, sys%st, sys%gr, ierr_rho)
189 write_density = restart_load_gs%has_map()
190 call restart_load_gs%end()
191 else
192 write_density = .true.
193 end if
194
195 safe_allocate(occ_states(1:sys%st%nik))
196 do ik = 1, sys%st%nik
197 call occupied_states(sys%st, sys%namespace, ik, n_filled, n_partially_filled, n_half_filled)
198 occ_states(ik) = n_filled + n_partially_filled + n_half_filled
199 end do
200
201 is_orbital_dependent = (sys%ks%theory_level == hartree .or. sys%ks%theory_level == hartree_fock .or. &
202 (sys%ks%theory_level == kohn_sham_dft .and. xc_is_orbital_dependent(sys%ks%xc)) &
203 .or. (sys%ks%theory_level == generalized_kohn_sham_dft .and. xc_is_orbital_dependent(sys%ks%xc)))
204
205 if (is_orbital_dependent) then
206 message(1) = "Be sure your gs run is well converged since you have an orbital-dependent functional."
207 message(2) = "Otherwise, the occupied states may change in CalculationMode = unocc, and your"
208 message(3) = "unoccupied states will not be consistent with the gs run."
209 call messages_warning(3, namespace=sys%namespace)
210 end if
211
212 if (ierr_rho /= 0 .or. is_orbital_dependent) then
213 occ_missing = .false.
214 do ik = 1, sys%st%nik
215 if (any(lowest_missing(1:sys%st%d%dim, ik) <= occ_states(ik))) then
216 occ_missing = .true.
217 end if
218 end do
219
220 if (occ_missing) then
221 if (is_orbital_dependent) then
222 message(1) = "For an orbital-dependent functional, all occupied orbitals must be provided."
223 else if (ierr_rho /= 0) then
224 message(1) = "Since density could not be read, all occupied orbitals must be provided."
225 end if
226
227 message(2) = "Not all the occupied orbitals could be read."
228 message(3) = "Please run a ground-state calculation first!"
229 call messages_fatal(3, only_root_writes = .true., namespace=sys%namespace)
230 end if
231
232 message(1) = "Unable to read density: Building density from wavefunctions."
233 call messages_info(1, namespace=sys%namespace)
234
235 call density_calc(sys%st, sys%gr, sys%st%rho)
236 end if
237
238 call scf_state_info(sys%namespace, sys%st)
239
240 if (fromscratch .or. ierr /= 0) then
241 if (fromscratch) then
242 ! do not use previously calculated occupied states
243 nst_calculated = min(maxval(occ_states), minval(lowest_missing) - 1)
244 else
245 ! or, use as many states as have been calculated
246 nst_calculated = minval(lowest_missing) - 1
247 end if
248 showstart = max(nst_calculated + 1, 1)
249 call lcao_run(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, &
250 sys%hm, st_start = showstart, known_lower_bound=known_lower_bound)
251 call eigens%set_lower_bound_is_known(known_lower_bound)
252 else
253 ! we successfully read all the states and are planning to use them, no need for LCAO
254 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, sys%st, sys%ions, sys%ext_partners, &
255 calc_eigenval = .false.)
256 showstart = minval(occ_states(:)) + 1
257
258 call eigens%set_lower_bound_is_known(.true.) ! Lower bound is known here
259 end if
260
261
262
263 ! it is strange and useless to see no eigenvalues written if you are only calculating
264 ! occupied states, on a different k-point.
265 if (showstart > sys%st%nst) showstart = 1
266
267 safe_deallocate_a(lowest_missing)
268
269 if (showoccstates) showstart = 1
270
271 ! In the case of someone using KPointsPath, the code assume that this is only for plotting a
272 ! bandstructure. This mode ensure that no restart information will be written for the new grid
273 if (bandstructure_mode) then
274 message(1) = "Info: The code will run in band structure mode."
275 message(2) = " No restart information will be printed."
276 call messages_info(2, namespace=sys%namespace)
277 end if
278
279 if (.not. bandstructure_mode) then
280 ! Restart dump should be initialized after restart_load, as the mesh might have changed
281 call restart_dump%init(sys%namespace, restart_unocc, restart_type_dump, sys%mc, ierr, mesh=sys%gr)
282
283 ! make sure the density is defined on the same mesh as the wavefunctions that will be written
284 if (write_density) then
285 call states_elec_dump_rho(restart_dump, sys%st, sys%gr, ierr_rho)
286 end if
287 end if
288
289 message(1) = "Info: Starting calculation of unoccupied states."
290 call messages_info(1, namespace=sys%namespace)
291
292 ! reset this variable, so that the eigensolver passes through all states
293 eigens%converged(:) = 0
294
295 ! If not all gs wavefunctions were read when starting, in particular for nscf with different k-points,
296 ! the occupations must be recalculated each time, though they do not affect the result of course.
297 ! FIXME: This is wrong for metals where we must use the Fermi level from the original calculation!
298 call states_elec_fermi(sys%st, sys%namespace, sys%gr)
299
300 if(output_needs_current(sys%outp, states_are_real(sys%st))) then
301 call current_init(current_calculator, sys%namespace)
302 end if
303
304 if (sys%st%pack_states .and. sys%hm%apply_packed()) call sys%st%pack()
305
306 do iter = 1, max_iter
307 output_iter = .false.
308 call eigens%run(sys%namespace, sys%gr, sys%st, sys%hm, sys%space, sys%ext_partners, 1, converged, sys%st%nst_conv)
309
310 ! If not all gs wavefunctions were read when starting, in particular for nscf with different k-points,
311 ! the occupations must be recalculated each time, though they do not affect the result of course.
312 ! FIXME: This is wrong for metals where we must use the Fermi level from the original calculation!
313 call states_elec_fermi(sys%st, sys%namespace, sys%gr)
314
315 call write_iter_(sys%namespace, sys%st)
316
317 ! write output file
318 if (sys%st%system_grp%is_root()) then
319 call io_mkdir(static_dir, sys%namespace)
320 iunit = io_open(static_dir//'/eigenvalues', sys%namespace, action='write')
321
322 if (converged) then
323 write(iunit,'(a)') 'All states converged.'
324 else
325 write(iunit,'(a)') 'Some of the states are not fully converged!'
326 end if
327 write(iunit,'(a, e17.6)') 'Criterion = ', eigens%tolerance
328 write(iunit,'(1x)')
329 call states_elec_write_eigenvalues(sys%st%nst, sys%st, sys%space, sys%kpoints, eigens%diff, iunit=iunit)
330 call io_close(iunit)
331 end if
332
333 forced_finish = clean_stop(sys%mc%master_comm)
334
335 if (.not. bandstructure_mode) then
336 ! write restart information.
337 if (converged .or. restart_walltime_period_alarm(sys%mc%master_comm) &
338 .or. iter == max_iter .or. forced_finish) then
339 call states_elec_dump(restart_dump, sys%space, sys%st, sys%gr, sys%kpoints, ierr, iter=iter)
340 if (ierr /= 0) then
341 message(1) = "Unable to write states wavefunctions."
342 call messages_warning(1, namespace=sys%namespace)
343 end if
344 end if
345 end if
346
347 do what_it = lbound(sys%outp%output_interval, 1), ubound(sys%outp%output_interval, 1)
348 if (sys%outp%what_now(what_it, iter)) then
349 output_iter = .true.
350 exit
351 end if
352 end do
353
354 if (output_iter .and. sys%outp%duringscf) then
355 write(dirname,'(a,i4.4)') "unocc.",iter
356 if(output_needs_current(sys%outp, states_are_real(sys%st))) then
357 call states_elec_allocate_current(sys%st, sys%space, sys%gr)
358 call current_calculate(current_calculator, sys%namespace, sys%gr, sys%hm, sys%space, sys%st)
359 end if
360 call output_all(sys%outp, sys%namespace, sys%space, dirname, sys%gr, sys%ions, iter, sys%st, sys%hm, sys%ks)
361 call output_modelmb(sys%outp, sys%namespace, sys%space, dirname, sys%gr, sys%ions, iter, sys%st)
362 end if
363
364 if (converged .or. forced_finish) exit
365
366 end do
367
368 if (.not. bandstructure_mode) call restart_dump%end()
369
370 if (sys%st%pack_states .and. sys%hm%apply_packed()) then
371 call sys%st%unpack()
372 end if
373
374 if (any(eigens%converged(:) < occ_states(:))) then
375 write(message(1),'(a)') 'Some of the occupied states are not fully converged!'
376 call messages_warning(1, namespace=sys%namespace)
377 end if
378
379 safe_deallocate_a(occ_states)
380
381 if (.not. converged) then
382 write(message(1),'(a)') 'Some of the unoccupied states are not fully converged!'
383 call messages_warning(1, namespace=sys%namespace)
384 end if
385
386 if (sys%space%is_periodic().and. sys%st%nik > sys%st%d%nspin) then
387 if (bitand(sys%kpoints%method, kpoints_path) /= 0) then
388 call states_elec_write_bandstructure(static_dir, sys%namespace, sys%st%nst, sys%st, &
389 sys%ions, sys%gr, sys%kpoints, &
390 sys%hm%phase, vec_pot = sys%hm%hm_base%uniform_vector_potential, &
391 vec_pot_var = sys%hm%hm_base%vector_potential)
392 end if
393 end if
394
395 if(output_needs_current(sys%outp, states_are_real(sys%st))) then
396 call states_elec_allocate_current(sys%st, sys%space, sys%gr)
397 call current_calculate(current_calculator, sys%namespace, sys%gr, sys%hm, sys%space, sys%st)
398 end if
399 call output_all(sys%outp, sys%namespace, sys%space, static_dir, sys%gr, sys%ions, -1, sys%st, sys%hm, sys%ks)
400 call output_modelmb(sys%outp, sys%namespace, sys%space, static_dir, sys%gr, sys%ions, -1, sys%st)
401
402 call end_()
403 pop_sub(unocc_run_legacy)
404
405 contains
406
407 ! ---------------------------------------------------------
408 subroutine init_(mesh, st)
409 class(mesh_t), intent(in) :: mesh
410 type(states_elec_t), intent(inout) :: st
411
412 push_sub(unocc_run_legacy.init_)
413
414 call messages_obsolete_variable(sys%namespace, "NumberUnoccStates", "ExtraStates")
415
416 call states_elec_allocate_wfns(st, mesh, packed=.true.)
417
418 ! now the eigensolver stuff
419 call eigensolver_init(eigens, sys%namespace, sys%gr, st, sys%hm, sys%mc, sys%space, deactivate_oracle=.true.)
420
421 if (eigens%es_type == rs_rmmdiis) then
422 message(1) = "With the RMMDIIS eigensolver for unocc, you will need to stop the calculation"
423 message(2) = "by hand, since the highest states will probably never converge."
424 call messages_warning(2, namespace=sys%namespace)
425 end if
426
427 pop_sub(unocc_run_legacy.init_)
428 end subroutine init_
429
430
431 ! ---------------------------------------------------------
432 subroutine end_()
433 push_sub(unocc_run_legacy.end_)
434
435 call eigensolver_end(eigens)
436
437 pop_sub(unocc_run_legacy.end_)
438 end subroutine end_
439
440 ! ---------------------------------------------------------
441 subroutine write_iter_(namespace, st)
442 type(namespace_t), intent(in) :: namespace
443 type(states_elec_t), intent(in) :: st
444
445 character(len=50) :: str
446
448
449 write(str, '(a,i5)') 'Unoccupied states iteration #', iter
450 call messages_print_with_emphasis(msg=trim(str), namespace=namespace)
451
452 write(message(1),'(a,i6)') 'Converged eigenvectors: ', minval(eigens%converged(1:st%nik))
453 call messages_info(1, namespace=namespace)
454
455 call states_elec_write_eigenvalues(sys%st%nst, sys%st, sys%space, sys%kpoints, &
456 eigens%diff, st_start = showstart, compact = .true., namespace=namespace)
457
458 call scf_print_mem_use(namespace)
459
460 call messages_print_with_emphasis(namespace=namespace)
461
463 end subroutine write_iter_
464
465
466 end subroutine unocc_run_legacy
467
468
469end module unocc_oct_m
470
471!! Local Variables:
472!! mode: f90
473!! coding: utf-8
474!! 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:372
subroutine, public current_init(this, namespace)
Definition: current.F90:180
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:619
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:250
integer, parameter, public generalized_kohn_sham_dft
Definition: global.F90:250
integer, parameter, public kohn_sham_dft
Definition: global.F90:250
character(len= *), parameter, public static_dir
Definition: global.F90:279
integer, parameter, public hartree
Definition: global.F90:250
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:223
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:769
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
Definition: lda_u_io.F90:723
integer, parameter, public dft_u_none
Definition: lda_u.F90:205
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:1068
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:1000
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
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:984
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
Definition: output.F90:498
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:1540
subroutine, public scf_print_mem_use(namespace)
Definition: scf.F90:1558
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:754
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:537