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