Octopus
unfold.F90
Go to the documentation of this file.
1!! Copyright (C) 2018-2019 M. S. Mrudul, N. Tancogne-Dejean
2!!
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20
21#include "global.h"
22
23
24program oct_unfold
25 use batch_oct_m
28 use comm_oct_m
30 use cube_oct_m
32 use debug_oct_m
34 use fft_oct_m
38 use global_oct_m
39 use grid_oct_m
40 use io_oct_m
43 use, intrinsic :: iso_fortran_env
47 use loct_oct_m
48 use math_oct_m
49 use mesh_oct_m
52 use mpi_oct_m
56 use parser_oct_m
60 use space_oct_m
64 use unit_oct_m
66 use utils_oct_m
68 use v_ks_oct_m
69 use xc_oct_m
70
71 implicit none
72
73 type(electrons_t), pointer :: sys
74 integer :: ik, nkpoints
75 type(restart_t) :: restart
76 type(cube_t) :: zcube
77 type(cube_function_t) :: cf
78
79 type(lattice_vectors_t) :: pc
80 integer :: ierr, run_mode, file_gvec
81 type(block_t) :: blk
82 integer :: nhighsympoints, nsegments
83 integer :: icol, idir, ncols
84
85 integer, allocatable :: resolution(:)
86 real(real64), allocatable :: highsympoints(:,:), coord_along_path(:)
87 type(kpoints_grid_t) :: path_kpoints_grid
88
89
90 ! the usual initializations
91 call global_init()
92 call parser_init()
93
94 call messages_init()
95
96 call io_init()
98
99 call print_header()
100 call messages_print_with_emphasis(msg="Unfolding Band-structure", namespace=global_namespace)
102
103 call messages_experimental("oct-unfold utility")
107
108 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
110 call sys%init_parallelization(mpi_world)
111
112 if (sys%space%periodic_dim == 0) then
113 message(1) = "oct-unfold can only be used for periodic systems."
114 call messages_fatal(1)
115 end if
116
117 if (sys%st%parallel_in_states) then
118 call messages_not_implemented("oct-unfold with states parallelization")
119 end if
120
121 if (sys%st%d%ispin == spinors) then
122 call messages_not_implemented("oct-unfold for spinors")
123 end if
124
125
126 !%Variable UnfoldMode
127 !%Type flag
128 !%Default none
129 !%Section Utilities::oct-unfold
130 !%Description
131 !% Specifies which stage of the unfolding tool to use
132 !%Option unfold_setup bit(1)
133 !% Writes the list of k-points corresponding to the path specified by <tt>UnfoldKPointPath</tt>.
134 !% This list of k-point (unfold_kpt.dat) must be used for an unocc calculation of the supercell,
135 !% adding the line "include 'unfold_kpt.dat'" to the inp file and removing the KPointGrid information.
136 !%Option unfold_run bit(2)
137 !% Perform the actual unfolding, based on the states obtained from the previous unocc run.
138 !%End
139 call parse_variable(global_namespace, 'UnfoldMode', 0, run_mode)
140 if (.not. varinfo_valid_option('UnfoldMode', run_mode)) then
141 call messages_input_error(global_namespace, "UnfoldMode must be set to a value different from 0.")
142 end if
143
144 !%Variable UnfoldLatticeParameters
145 !%Type block
146 !%Section Utilities::oct-unfold
147 !%Description
148 !% The lattice parameters of the primitive cell, on which unfolding is performed.
149 !% See the LatticeParameters variable for a more detailed description.
150 !%End
151
152 !%Variable UnfoldLatticeVectors
153 !%Type block
154 !%Default simple cubic
155 !%Section Utilities::oct-unfold
156 !%Description
157 !% Lattice vectors of the primitive cell on which the unfolding is performed.
158 !% See the LatticeVectors variable for a more detailed description.
159 !%End
160 pc = lattice_vectors_t(global_namespace, sys%space, variable_prefix='Unfold')
161
162 !%Variable UnfoldKPointsPath
163 !%Type block
164 !%Section Utilities::oct-unfold
165 !%Description
166 !% Specifies the k-point path for which the unfolding need to be done.
167 !% The syntax is identical to <tt>KPointsPath</tt>.
168 !%End
169 if (parse_block(global_namespace, 'UnfoldKPointsPath', blk) /= 0) then
170 write(message(1),'(a)') 'Error while reading UnfoldPointsPath.'
171 call messages_fatal(1)
172 end if
173
174 ! There is one high symmetry k-point per line
175 nsegments = parse_block_cols(blk, 0)
176 nhighsympoints = parse_block_n(blk) - 1
177 if (nhighsympoints /= nsegments+1) then
178 write(message(1),'(a,i3,a,i3)') 'The first row of UnfoldPointsPath is not compatible with the number of specified k-points.'
179 call messages_fatal(1)
180 end if
181
182 safe_allocate(resolution(1:nsegments))
183 do icol = 1, nsegments
184 call parse_block_integer(blk, 0, icol-1, resolution(icol))
185 end do
186 !Total number of points in the segment
187 nkpoints = sum(resolution) + 1
188
189 safe_allocate(highsympoints(1:sys%space%dim, 1:nhighsympoints))
190 do ik = 1, nhighsympoints
191 !Sanity check
192 ncols = parse_block_cols(blk, ik)
193 if (ncols /= sys%space%dim) then
194 write(message(1),'(a,i8,a,i3)') 'UnfoldPointsPath row ', ik, ' has ', ncols, ' columns but must have ', sys%space%dim
195 call messages_fatal(1)
196 end if
197
198 do idir = 1, sys%space%dim
199 call parse_block_float(blk, ik, idir-1, highsympoints(idir, ik))
200 end do
201 end do
202
203 call parse_block_end(blk)
204
205 call kpoints_grid_init(sys%space%dim, path_kpoints_grid, nkpoints, 1)
206 ! For the output of band-structures
207 safe_allocate(coord_along_path(1:nkpoints))
208
209 call kpoints_path_generate(sys%space%dim, pc, nkpoints, nsegments, resolution, &
210 highsympoints, path_kpoints_grid%point, coord_along_path)
211
212 safe_deallocate_a(resolution)
213 safe_deallocate_a(highsympoints)
214
215 !We convert the k-point to the reduced coordinate of the supercell
216 do ik = 1, path_kpoints_grid%npoints
217 call kpoints_to_reduced(sys%kpoints%latt, path_kpoints_grid%point(:, ik), path_kpoints_grid%red_point(:, ik))
218 end do
219
220 call kpoints_fold_to_1bz(path_kpoints_grid, pc)
221
222 if (run_mode == option__unfoldmode__unfold_setup) then
223
224 call unfold_setup()
225
226 else if (run_mode == option__unfoldmode__unfold_run) then
227
228 !Sanity check
229 file_gvec = io_open('unfold_gvec.dat', global_namespace, action='read')
230 read(file_gvec, *)
231 read(file_gvec, *) ik
232 if (ik /= path_kpoints_grid%npoints) then
233 message(1) = 'There is an inconsistency between unfold_gvec.dat and the input file'
234 call messages_fatal(1)
235 end if
236 call io_close(file_gvec)
237
238 call states_elec_allocate_wfns(sys%st, sys%gr)
239
240 call restart_init(restart, global_namespace, restart_unocc, restart_type_load, sys%mc, ierr, mesh=sys%gr, exact=.true.)
241 if (ierr == 0) then
242 call states_elec_load(restart, global_namespace, sys%space, sys%st, sys%gr, sys%kpoints, fixed_occ=.true., &
243 ierr=ierr, label = ": unfold")
244 end if
245 if (ierr /= 0) then
246 message(1) = 'Unable to read unocc wavefunctions.'
247 call messages_fatal(1)
248 end if
249 call restart_end(restart)
250
251 call cube_init(zcube, sys%gr%idx%ll, global_namespace, sys%space, sys%gr%spacing, &
252 sys%gr%coord_system, fft_type = fft_complex, dont_optimize = .true.)
253 call cube_init_cube_map(zcube, sys%gr)
254
255 call zcube_function_alloc_rs(zcube, cf)
256 call cube_function_alloc_fs(zcube, cf)
257
258 call wfs_extract_spec_fn(sys%space, sys%st, sys%gr, zcube, cf)
259
260 call cube_function_free_fs(zcube, cf)
261 call zcube_function_free_rs(zcube, cf)
262 call cube_end(zcube)
263
264 else
265 message(1) = "Unsupported or incorrect value of UnfoldMode."
266 call messages_fatal(1)
267 end if
268
269 safe_deallocate_a(coord_along_path)
270
271 call kpoints_grid_end(path_kpoints_grid)
272
273 safe_deallocate_p(sys)
274 call fft_all_end()
276 call io_end()
277 call print_date("Calculation ended on ")
278 call messages_end()
279 call parser_end()
280 call global_end()
281
282
283
284contains
285 !-----------------------------------------------------------------
286 subroutine unfold_setup()
287 integer :: file_gvec, file_kpts, idir
288 integer :: gvec(sys%space%dim)
289
290 push_sub(unfold_setup)
291
292 if (mpi_grp_is_root(mpi_world)) then
293 file_gvec = io_open('unfold_gvec.dat', global_namespace, action='write')
294 file_kpts = io_open('unfold_kpt.dat', global_namespace, action='write')
295 write(file_kpts,'(a)') '%KpointsReduced'
296 write(file_gvec,'(a)') '#Created by oct-unfold'
297 write(file_gvec,'(i5)') path_kpoints_grid%npoints
298
299 !We convert the k-point to the reduce coordinate of the supercell
300 do ik = 1, path_kpoints_grid%npoints
301 gvec(:) = nint(path_kpoints_grid%red_point(:, ik) + m_half * 1e-7_real64)
302 write(file_kpts,'(a3)', advance='no') ' 1.'
303 write(file_kpts,'(*(a3,f12.8))') &
304 (' | ', path_kpoints_grid%red_point(idir, ik) - gvec(idir), idir = 1, sys%space%dim)
305 write(file_gvec,'(*(i5))') (gvec(idir), idir = 1, sys%space%dim)
306 end do
307 write(file_kpts,'(a)') '%'
308 call io_close(file_gvec)
309 call io_close(file_kpts)
310 end if
311
312 pop_sub(unfold_setup)
313 end subroutine unfold_setup
314
315 !--------------------------------------------------------------------
316 subroutine wfs_extract_spec_fn(space, st, gr, zcube, cf)
317 class(space_t), intent(in) :: space
318 type(states_elec_t), intent(in) :: st
319 type(grid_t), intent(in) :: gr
320 type(cube_t), intent(inout) :: zcube
321 type(cube_function_t), intent(inout) :: cf
322
323 real(real64), allocatable :: pkm(:,:), ake(:,:), eigs(:)
324 complex(real64), allocatable :: zpsi(:), field_g(:)
325 integer :: file_ake, iq, ist, idim, nenergy
326 integer :: ig, ix, iy, iz, ik, ie, gmin, gmax
327 real(real64) :: eigmin, eigmax, de, norm
328 real(real64), parameter :: tol = 1e-7_real64
329 integer, parameter :: nextend = 10
330 real(real64) :: vec_pc(space%dim), vec_sc(space%dim)
331 type(fourier_shell_t) :: shell
332 character(len=MAX_PATH_LEN) :: filename
333 real(real64), allocatable :: gvec_abs(:,:)
334 logical, allocatable :: g_select(:)
335
336 push_sub(wfs_extract_spec_fn)
337
338 safe_allocate(zpsi(1:gr%np))
339
340 !%Variable UnfoldEnergyStep
341 !%Type float
342 !%Default 0
343 !%Section Utilities::oct-unfold
344 !%Description
345 !% Specifies the energy resolution for the unfolded band structure.
346 !% If you specify 0, the resolution will be set to be 1/1000 points between <tt>UnfoldMinEnergy</tt>
347 !% and <tt>UnfoldMaxEnergy</tt>
348 !%End
349 call parse_variable(global_namespace, 'UnfoldEnergyStep', m_zero, de)
350 if (de < m_zero) then
351 message(1) = "UnfoldEnergyStep must be positive"
352 call messages_fatal(1)
353 end if
354
355 !%Variable UnfoldMinEnergy
356 !%Type float
357 !%Section Utilities::oct-unfold
358 !%Description
359 !% Specifies the start of the energy range for the unfolded band structure.
360 !% The default value correspond to the samllest eigenvalue.
361 !%End
362 call parse_variable(global_namespace, 'UnfoldMinEnergy', minval(st%eigenval(:, :)), eigmin)
363
364 !%Variable UnfoldMaxEnergy
365 !%Type float
366 !%Section Utilities::oct-unfold
367 !%Description
368 !% Specifies the end of the energy range for the unfolded band structure.
369 !% The default value correspond to the largest eigenvalue.
370 !%End
371 call parse_variable(global_namespace, 'UnfoldMaxEnergy', maxval(st%eigenval(:, :)), eigmax)
372
373 if (abs(de) <= m_epsilon) then
374 de = (eigmax - eigmin) / 1000_real64
375 end if
376
377 !We increase a bit the energy range
378 nenergy = nint((eigmax - eigmin + 2 * nextend * de) / de)
379 safe_allocate(eigs(1:nenergy))
380 do ie = 1, nenergy
381 eigs(ie) = eigmin - nextend * de + (ie - 1) * de
382 end do
383
384 safe_allocate(gvec_abs(1:sys%space%periodic_dim, 1:sys%kpoints%reduced%npoints))
385 gvec_abs = 0
386 file_gvec = io_open('./unfold_gvec.dat', global_namespace, action='read')
387 read(file_gvec,*)
388 read(file_gvec,*)
389 do ik = 1, sys%kpoints%reduced%npoints
390 read(file_gvec,*) vec_sc(1:space%dim)
391 call kpoints_to_absolute(sys%kpoints%latt, vec_sc, gvec_abs(:, ik))
392 end do
393 call io_close(file_gvec)
394
395 if (mpi_grp_is_root(mpi_world)) call loct_progress_bar(-1, (st%d%kpt%end - st%d%kpt%start + 1) * st%nst)
396
397 safe_allocate(ake(1:nenergy, 1:st%nik))
398 ake(:, :) = m_zero
399
400 safe_allocate(pkm(st%d%kpt%start:st%d%kpt%end, 1:st%nst))
401 pkm(:, :) = m_zero
402 do ik = st%d%kpt%start, st%d%kpt%end
403 iq = st%d%get_kpoint_index(ik)
404
405 call fourier_shell_init(shell, global_namespace, space, zcube, gr, kk = sys%kpoints%reduced%red_point(:, iq))
406
407 gmin = minval(shell%red_gvec(:,:))
408 gmax = maxval(shell%red_gvec(:,:))
410 safe_allocate(g_select(1:shell%ngvectors))
411 g_select(:) = .false.
412
413 select case (sys%space%periodic_dim)
414 case (3)
415 do ig = 1, shell%ngvectors
416 call kpoints_to_absolute(sys%kpoints%latt, real(shell%red_gvec(:,ig), real64), vec_sc(:))
417 do ix = gmin, gmax
418 do iy = gmin, gmax
419 do iz = gmin, gmax
420 vec_pc(1:3) = ix * pc%klattice(1:3,1) + iy * pc%klattice(1:3,2) + iz * pc%klattice(1:3,3)
421 if (abs(vec_sc(1) - vec_pc(1) - gvec_abs(1, iq)) < tol &
422 .and. abs(vec_sc(2) - vec_pc(2)-gvec_abs(2, iq)) < tol &
423 .and. abs(vec_sc(3) - vec_pc(3)-gvec_abs(3, iq)) < tol) then
424 g_select(ig) = .true.
425 end if
426 end do !iz
427 end do !iy
428 end do !ix
429 end do !ig
430
431 case (2)
432
433 do ig = 1, shell%ngvectors
434 call kpoints_to_absolute(sys%kpoints%latt, real(shell%red_gvec(:,ig), real64), vec_sc(:))
435 do ix = gmin, gmax
436 do iy = gmin, gmax
437 vec_pc(1:2) = ix * pc%klattice(1:2,1) + iy * pc%klattice(1:2,2)
438 if (abs(vec_sc(1) - vec_pc(1) - gvec_abs(1, iq)) < tol &
439 .and. abs(vec_sc(2) - vec_pc(2) - gvec_abs(2, iq)) < tol) then
440 g_select(ig) = .true.
441 end if
442 end do !ix
443 end do !ix
444 end do !ig
445
446 case default
447 call messages_not_implemented("Unfolding for periodic dimensions other than 2 or 3")
448 end select
449
450 if (mpi_grp_is_root(gr%mpi_grp)) then
451 write(filename,"(a,i3.3,a4)") trim(adjustl(static_dir))//"ake_",ik,".dat"
452 print *, '-'// filename
453 file_ake = io_open(trim(filename), global_namespace, action='write')
454 write(file_ake, '(a)') '#Energy Ak(E)'
455 write(file_ake, '(a, i5)') '#Number of points in energy window ', nenergy
456 end if
457
458 do ist = 1, st%nst
459 !loop over states
460 do idim = 1, st%d%dim
461 ! Getting wavefunctions
462 ! for the moment we treat all functions as complex
463 call states_elec_get_state(st, gr, idim, ist, ik, zpsi)
464
465 call zmesh_to_cube(gr, zpsi, zcube, cf)
466
467 !Fourier transform from real-space to fourier space
468 call zcube_function_rs2fs(zcube, cf)
469
470 ! Normalisation
471 safe_allocate(field_g(1:shell%ngvectors))
472 norm = m_zero
473 do ig = 1, shell%ngvectors
474 field_g(ig) = cf%fs(shell%coords(1, ig), shell%coords(2, ig), shell%coords(3, ig))
475 norm = norm + abs(field_g(ig))**2
476 end do
477 field_g(:) = field_g(:) / sqrt(norm)
478
479 !Finding sub-g and calculating the Projection Pkm
480 do ig = 1, shell%ngvectors
481 if (.not. g_select(ig)) cycle
482 pkm(ik,ist) = pkm(ik,ist) + abs(field_g(ig))**2
483 end do
484
485 safe_deallocate_a(field_g)
486 end do !idim
487
488 !ist loop end
489 if (mpi_grp_is_root(mpi_world)) then
490 call loct_progress_bar((ik - st%d%kpt%start) * st%nst + ist, &
491 (st%d%kpt%end - st%d%kpt%start + 1) * st%nst)
492 end if
493 end do !ist
494
495 ! Calculating the spectral function
496 !TODO: We could implement here a different broadening
497 do ist = 1, st%nst
498 do ie = 1, nenergy
499 ake(ie, ik) = ake(ie, ik) + pkm(ik, ist) * (m_three * de / m_pi) / &
500 ((eigs(ie) - st%eigenval(ist, ik))**2 + (m_three * de)**2)
501 end do
502 end do
503
504 ! writing Spectral-Function
505 if (mpi_grp_is_root(gr%mpi_grp)) then
506 do ie = 1, nenergy
507 write(file_ake, '(1es19.12,1x,1es19.12)') eigs(ie), ake(ie, ik)
508 end do
509
510 call io_close(file_ake)
511 end if
512
513 call fourier_shell_end(shell)
514 safe_deallocate_a(g_select)
515
516 end do !ik
517
518 if (st%d%kpt%parallel) then
519 call comm_allreduce(st%st_kpt_mpi_grp, ake)
520 end if
521
522 if (mpi_grp_is_root(mpi_world)) then
523 file_ake = io_open(trim(adjustl(static_dir))//"ake.dat", global_namespace, action='write')
524 write(file_ake, '(a)') '#Energy Ak(E)'
525 write(file_ake, '(a, i5)') '#Number of points in energy window ', nenergy
526 do ik = 1, nkpoints
527 do ie = 1, nenergy
528 write(file_ake,fmt ='(1es19.12,1x,1es19.12,1x,1es19.12)') coord_along_path(ik), &
529 eigs(ie), ake(ie, ik)
530 end do
531 end do
532
533 call io_close(file_ake)
534 end if
535
536 safe_deallocate_a(eigs)
537 safe_deallocate_a(ake)
538
539 safe_deallocate_a(gvec_abs)
540 safe_deallocate_a(pkm)
541 safe_deallocate_a(zpsi)
542 pop_sub(wfs_extract_spec_fn)
543
544 end subroutine wfs_extract_spec_fn
545 !------------------------------------------------------------
546
547end program oct_unfold
548!! Local Variables:
549!! mode: f90
550!! coding: utf-8
551!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
subroutine, public zmesh_to_cube(mesh, mf, cube, cf)
Convert a function from the mesh to the cube.
subroutine, public zcube_function_free_rs(cube, cf)
Deallocates the real space grid.
subroutine, public zcube_function_alloc_rs(cube, cf, in_device, force_alloc)
Allocates locally the real space grid, if PFFT library is not used. Otherwise, it assigns the PFFT re...
subroutine, public cube_init(cube, nn, namespace, space, spacing, coord_system, fft_type, fft_library, dont_optimize, nn_out, mpi_grp, need_partition, tp_enlarge, blocksize)
Definition: cube.F90:202
subroutine, public cube_end(cube)
Definition: cube.F90:378
subroutine, public cube_init_cube_map(cube, mesh)
Definition: cube.F90:823
integer, parameter, public spinors
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:118
subroutine, public fft_all_init(namespace)
initialize the table
Definition: fft.F90:278
subroutine, public fft_all_end()
delete all plans
Definition: fft.F90:391
integer, parameter, public fft_complex
Definition: fft.F90:178
subroutine, public fourier_shell_init(this, namespace, space, cube, mesh, kk)
subroutine, public fourier_shell_end(this)
subroutine, public cube_function_free_fs(cube, cf)
Deallocates the Fourier space grid.
subroutine, public zcube_function_rs2fs(cube, cf)
The following routines convert the function between real space and Fourier space Note that the dimens...
subroutine, public cube_function_alloc_fs(cube, cf, force_alloc)
Allocates locally the Fourier space grid, if PFFT library is not used. Otherwise, it assigns the PFFT...
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:382
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
real(real64), parameter, public m_epsilon
Definition: global.F90:204
character(len= *), parameter, public static_dir
Definition: global.F90:252
subroutine, public global_init(communicator)
Initialise Octopus.
Definition: global.F90:325
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_three
Definition: global.F90:191
This module implements the underlying real-space grid.
Definition: grid.F90:117
Definition: io.F90:114
subroutine, public io_init(defaults)
If the argument defaults is present and set to true, then the routine will not try to read anything f...
Definition: io.F90:161
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_end()
Definition: io.F90:258
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
subroutine, public kpoints_path_generate(dim, latt, nkpoints, nsegments, resolution, highsympoints, kpoints, coord)
Generate the k-point along a path.
Definition: kpoints.F90:1243
subroutine, public kpoints_fold_to_1bz(grid, latt)
Definition: kpoints.F90:1420
subroutine, public kpoints_grid_end(this)
Definition: kpoints.F90:240
subroutine, public kpoints_to_reduced(latt, kin, kout)
Definition: kpoints.F90:1044
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1031
subroutine, public kpoints_grid_init(dim, this, npoints, nshifts)
Definition: kpoints.F90:219
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_end()
Definition: messages.F90:278
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:921
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_init(output_dir)
Definition: messages.F90:225
subroutine, public print_date(str)
Definition: messages.F90:1006
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1086
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:115
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 handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
type(namespace_t), public global_namespace
Definition: namespace.F90:132
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:450
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:481
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public profiling_end(namespace)
Definition: profiling.F90:413
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:255
subroutine, public restart_module_init(namespace)
Definition: restart.F90:308
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_load
Definition: restart.F90:245
integer, parameter, public restart_unocc
Definition: restart.F90:200
subroutine, public restart_end(restart)
Definition: restart.F90:722
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, 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...
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
subroutine, public unit_system_init(namespace)
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
Definition: utils.F90:301
Definition: xc.F90:114
Class describing the electron system.
Definition: electrons.F90:218
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
The states_elec_t class contains all electronic wave functions.
int true(void)
subroutine unfold_setup()
Definition: unfold.F90:380
subroutine wfs_extract_spec_fn(space, st, gr, zcube, cf)
Definition: unfold.F90:410
program oct_unfold
Definition: unfold.F90:117