Octopus
dos.F90
Go to the documentation of this file.
1!! Copyright (C) 2017 N. Tancogne-Dejean
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 dos_oct_m
22 use accel_oct_m
24 use box_oct_m
25 use comm_oct_m
26 use debug_oct_m
28 use global_oct_m
30 use io_oct_m
31 use ions_oct_m
32 use, intrinsic :: iso_fortran_env
34 use math_oct_m
35 use mesh_oct_m
37 use mpi_oct_m
41 use parser_oct_m
47 use types_oct_m
48 use unit_oct_m
50
51 implicit none
52
53 private
54 public :: &
55 dos_t, &
56 dos_init, &
59
60 type dos_t
61 private
62 real(real64) :: emin
63 real(real64) :: emax
64 integer :: epoints
65 real(real64) :: gamma
66 real(real64) :: de
67 logical :: computepdos
68 end type dos_t
69
70contains
71
72 subroutine dos_init(this, namespace, st, kpoints)
73 type(dos_t), intent(out) :: this
74 type(namespace_t), intent(in) :: namespace
75 type(states_elec_t), intent(in) :: st
76 type(kpoints_t), intent(in) :: kpoints
77
78 real(real64) :: evalmin, evalmax, eextend
79 integer :: npath
80
81 push_sub(dos_init)
82
83 !The range of the dos is only calculated for physical points,
84 !without the one from a k-point path
85 npath = kpoints%nkpt_in_path()
86 if (st%nik > npath) then
87 evalmin = minval(st%eigenval(1:st%nst, 1:(st%nik-npath)))
88 evalmax = maxval(st%eigenval(1:st%nst, 1:(st%nik-npath)))
89 else !In case we only have a path, e.g., a bandstructure calculation
90 evalmin = minval(st%eigenval(1:st%nst, 1:st%nik))
91 evalmax = maxval(st%eigenval(1:st%nst, 1:st%nik))
92 end if
93 ! we extend the energy mesh by this amount
94 eextend = (evalmax - evalmin) / m_four
95
96 !%Variable DOSEnergyMin
97 !%Type float
98 !%Section Output
99 !%Description
100 !% Lower bound for the energy mesh of the DOS.
101 !% The default is the lowest eigenvalue, minus a quarter of the total range of eigenvalues.
102 !% This is ignored for the joint density of states, and the minimal energy is always set to zero.
103 !%End
104 call parse_variable(namespace, 'DOSEnergyMin', evalmin - eextend, this%emin, units_inp%energy)
105
106 !%Variable DOSEnergyMax
107 !%Type float
108 !%Section Output
109 !%Description
110 !% Upper bound for the energy mesh of the DOS.
111 !% The default is the highest eigenvalue, plus a quarter of the total range of eigenvalues.
112 !%End
113 call parse_variable(namespace, 'DOSEnergyMax', evalmax + eextend, this%emax, units_inp%energy)
115 !%Variable DOSEnergyPoints
116 !%Type integer
117 !%Default 500
118 !%Section Output
119 !%Description
120 !% Determines how many energy points <tt>Octopus</tt> should use for
121 !% the DOS energy grid.
122 !%End
123 call parse_variable(namespace, 'DOSEnergyPoints', 500, this%epoints)
124
125 !%Variable DOSGamma
126 !%Type float
127 !%Default 0.008 Ha
128 !%Section Output
129 !%Description
130 !% Determines the width of the Lorentzian which is used for the DOS sum.
131 !%End
132 call parse_variable(namespace, 'DOSGamma', units_from_atomic(units_inp%energy, 0.008_real64), this%gamma)
133 this%gamma = units_to_atomic(units_inp%energy, this%gamma)
134
135 !%Variable DOSComputePDOS
136 !%Type logical
137 !%Default false
138 !%Section Output
139 !%Description
140 !% Determines if projected dos are computed or not.
141 !% At the moment, the PDOS is computed from the bare pseudo-atomic orbitals, directly taken from
142 !% the pseudopotentials. The orbitals are not orthonormalized, in order to preserve their
143 !% atomic orbitals character. As a consequence, the sum of the different PDOS does not integrate
144 !% to the total DOS.
145 !%
146 !% The radii of the orbitals are controled by the threshold defined by <tt>AOThreshold</tt>,
147 !% and the fact that they are normalized or not by <tt>AONormalize</tt>.
148 !%End
149 call parse_variable(namespace, 'DOSComputePDOS', .false., this%computepdos)
150
151 ! spacing for energy mesh
152 this%de = (this%emax - this%emin) / (this%epoints - 1)
154 pop_sub(dos_init)
155 end subroutine dos_init
157 ! ---------------------------------------------------------
158 subroutine dos_write_dos(this, dir, st, box, ions, mesh, hm, namespace)
159 type(dos_t), intent(in) :: this
160 character(len=*), intent(in) :: dir
161 type(states_elec_t), intent(in) :: st
162 class(box_t), intent(in) :: box
163 type(ions_t), target, intent(in) :: ions
164 class(mesh_t), intent(in) :: mesh
165 type(hamiltonian_elec_t), intent(in) :: hm
166 type(namespace_t), intent(in) :: namespace
167
168 integer :: ie, ik, ist, is, ns, maxdos
169 integer, allocatable :: iunit(:)
170 real(real64) :: energy
171 real(real64), allocatable :: tdos(:)
172 real(real64), allocatable :: dos(:,:,:)
173 character(len=64) :: filename,format_str
174 logical :: normalize
175
176 integer :: ii, ll, mm, nn, work, norb, work2
177 integer :: ia, iorb, idim
178 real(real64) :: threshold
179 real(real64), allocatable :: dpsi(:,:), ddot(:,:)
180 complex(real64), allocatable :: zpsi(:,:), zdot(:,:)
181 real(real64), allocatable :: weight(:,:,:)
182 type(orbitalset_t) :: os
183
184 push_sub(dos_write_dos)
185
186 ! shortcuts
187 ns = 1
188 if (st%d%nspin == 2) ns = 2
189
190 if (mpi_grp_is_root(mpi_world)) then
191 ! space for state-dependent DOS
192 safe_allocate(dos(1:this%epoints, 1:st%nst, 0:ns-1))
193 safe_allocate(iunit(0:ns-1))
194
195 ! compute band/spin-resolved density of states
196 do ist = 1, st%nst
197
198 do is = 0, ns-1
199 if (ns > 1) then
200 write(filename, '(a,i4.4,a,i1.1,a)') 'dos-', ist, '-', is+1,'.dat'
201 else
202 write(filename, '(a,i4.4,a)') 'dos-', ist, '.dat'
203 end if
204 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
205 ! write header
206 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], band-resolved DOS'
207 end do
208
209 do ie = 1, this%epoints
210 energy = this%emin + (ie - 1) * this%de
211 dos(ie, ist, :) = m_zero
212 ! sum up Lorentzians
213 do ik = 1, st%nik, ns
214 do is = 0, ns-1
215 dos(ie, ist, is) = dos(ie, ist, is) + st%kweights(ik+is) * m_one/m_pi * &
216 this%gamma / ((energy - st%eigenval(ist, ik+is))**2 + this%gamma**2)
217 end do
218 end do
219 do is = 0, ns-1
220 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
221 units_from_atomic(unit_one / units_out%energy, dos(ie, ist, is))
222 call messages_info(1, iunit(is))
223 end do
224 end do
225
226 do is = 0, ns-1
227 call io_close(iunit(is))
228 end do
229 end do
230
231 safe_allocate(tdos(1))
232
233 ! for spin-polarized calculations also output spin-resolved tDOS
234 if (st%d%nspin > 1) then
235 do is = 0, ns-1
236 write(filename, '(a,i1.1,a)') 'total-dos-', is+1,'.dat'
237 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
238 ! write header
239 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS (spin-resolved)'
240
241 do ie = 1, this%epoints
242 energy = this%emin + (ie - 1) * this%de
243 tdos(1) = m_zero
244 do ist = 1, st%nst
245 tdos(1) = tdos(1) + dos(ie, ist, is)
246 end do
247 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
248 units_from_atomic(unit_one / units_out%energy, tdos(1))
249 call messages_info(1, iunit(is))
250 end do
252 call io_close(iunit(is))
253 end do
254 end if
255
256
257 iunit(0) = io_open(trim(dir)//'/'//'total-dos.dat', namespace, action='write')
258 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS'
259
260 ! compute total density of states
261 do ie = 1, this%epoints
262 energy = this%emin + (ie - 1) * this%de
263 tdos(1) = m_zero
264 do ist = 1, st%nst
265 do is = 0, ns-1
266 tdos(1) = tdos(1) + dos(ie, ist, is)
267 end do
268 end do
269 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
270 units_from_atomic(unit_one / units_out%energy, tdos(1))
271 call messages_info(1, iunit(0))
272 end do
273
274 call io_close(iunit(0))
275
276 safe_deallocate_a(tdos)
277
278
279 ! write Fermi file
280 iunit(0) = io_open(trim(dir)//'/'//'total-dos-efermi.dat', namespace, action='write')
281 write(message(1), '(3a)') '# Fermi energy [', trim(units_abbrev(units_out%energy)), &
282 '] in a format compatible with total-dos.dat'
283
284 ! this is the maximum that tdos can reach
285 maxdos = st%smear%el_per_state * st%nst
286
287 write(message(2), '(2f12.6)') units_from_atomic(units_out%energy, st%smear%e_fermi), m_zero
288 write(message(3), '(f12.6,i6)') units_from_atomic(units_out%energy, st%smear%e_fermi), maxdos
289
290 call messages_info(3, iunit(0))
291 call io_close(iunit(0))
292
293 end if
294
295 if (this%computepdos) then
296
297 !These variables are defined in basis_set/orbitalbasis.F90
298 call parse_variable(namespace, 'AOThreshold', 0.01_real64, threshold)
299 call parse_variable(namespace, 'AONormalize', .true., normalize)
300
301 if (states_are_real(st)) then
302 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
303 else
304 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
305 end if
306
307
308 do ia = 1, ions%natoms
309 !We first count how many orbital set we have
310 work = orbitalset_utils_count(ions%atom(ia)%species)
311
312 !We loop over the orbital sets of the atom ia
313 do norb = 1, work
314 call orbitalset_init(os)
315 os%spec => ions%atom(ia)%species
316
317 !We count the orbitals
318 work2 = 0
319 do iorb = 1, os%spec%get_niwfs()
320 call os%spec%get_iwf_ilm(iorb, 1, ii, ll, mm)
321 call os%spec%get_iwf_n(iorb, 1, nn)
322 if (ii == norb) then
323 os%ll = ll
324 os%nn = nn
325 os%ii = ii
326 os%radius = atomic_orbital_get_radius(os%spec, mesh, iorb, 1, &
327 option__aotruncation__ao_full, threshold)
328 work2 = work2 + 1
329 end if
330 end do
331 os%norbs = work2
332 os%ndim = 1
333 os%use_submesh = .false.
334
335 do work = 1, os%norbs
336 ! We obtain the orbital
337 if (states_are_real(st)) then
338 call dget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
339 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
340 os, work, os%radius, os%ndim, use_mesh=.not.os%use_submesh, &
341 normalize = normalize)
342 else
343 call zget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
344 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
345 os, work, os%radius, os%ndim, &
346 use_mesh = .not. hm%phase%is_allocated() .and. .not. os%use_submesh, &
347 normalize = normalize)
348 end if
349 end do
350
351 if (hm%phase%is_allocated()) then
352 ! In case of complex wavefunction, we allocate the array for the phase correction
353 safe_allocate(os%phase(1:os%sphere%np, st%d%kpt%start:st%d%kpt%end))
354 os%phase(:,:) = m_zero
355
356 if (.not. os%use_submesh) then
357 safe_allocate(os%eorb_mesh(1:mesh%np, 1:os%norbs, 1:os%ndim, st%d%kpt%start:st%d%kpt%end))
358 os%eorb_mesh(:,:,:,:) = m_zero
359 else
360 safe_allocate(os%eorb_submesh(1:os%sphere%np, 1:os%ndim, 1:os%norbs, st%d%kpt%start:st%d%kpt%end))
361 os%eorb_submesh(:,:,:,:) = m_zero
362 end if
363
364 if (accel_is_enabled() .and. st%d%dim == 1) then
365 os%ldorbs = max(pad_pow2(os%sphere%np), 1)
366 if(.not. os%use_submesh) os%ldorbs = max(pad_pow2(os%sphere%mesh%np), 1)
367
368 safe_allocate(os%buff_eorb(st%d%kpt%start:st%d%kpt%end))
369 do ik= st%d%kpt%start, st%d%kpt%end
370 call accel_create_buffer(os%buff_eorb(ik), accel_mem_read_only, type_cmplx, os%ldorbs*os%norbs)
371 end do
372 end if
373
374 call orbitalset_update_phase(os, box%dim, st%d%kpt, hm%kpoints, (st%d%ispin==spin_polarized), &
375 vec_pot = hm%hm_base%uniform_vector_potential, &
376 vec_pot_var = hm%hm_base%vector_potential)
377 end if
378
379 if (mpi_grp_is_root(mpi_world)) then
380 if (os%nn /= 0) then
381 write(filename,'(a, i3.3, a1, a, i1.1, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
382 os%nn, l_notation(os%ll), '.dat'
383 else
384 write(filename,'(a, i3.3, a1, a, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
385 l_notation(os%ll), '.dat'
386 end if
387
388 iunit(0) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
389 ! write header
390 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), &
391 '], projected DOS (total and orbital resolved)'
392 end if
393
394 if (states_are_real(st)) then
395 safe_allocate(ddot(1:st%d%dim,1:os%norbs))
396 else
397 safe_allocate(zdot(1:st%d%dim,1:os%norbs))
398 end if
399
400 safe_allocate(weight(1:os%norbs,1:st%nik,1:st%nst))
401 weight(1:os%norbs,1:st%nik,1:st%nst) = m_zero
402
403 do ist = st%st_start, st%st_end
404 do ik = st%d%kpt%start, st%d%kpt%end
405 if (abs(st%kweights(ik)) <= m_epsilon) cycle
406 if (states_are_real(st)) then
407 call states_elec_get_state(st, mesh, ist, ik, dpsi)
408 call dorbitalset_get_coefficients(os, st%d%dim, dpsi, ik, .false., ddot(:,1:os%norbs))
409 do iorb = 1, os%norbs
410 do idim = 1, st%d%dim
411 weight(iorb,ik,ist) = weight(iorb,ik,ist) + st%kweights(ik)*abs(ddot(idim,iorb))**2
412 end do
413 end do
414 else
415 call states_elec_get_state(st, mesh, ist, ik, zpsi)
416 if (hm%phase%is_allocated()) then
417 ! Apply the phase that contains both the k-point and vector-potential terms.
418 call hm%phase%apply_to_single(zpsi, mesh%np, st%d%dim, ik, .false.)
419 end if
420 call zorbitalset_get_coefficients(os, st%d%dim, zpsi, ik, hm%phase%is_allocated(), &
421 zdot(:,1:os%norbs))
422
423 do iorb = 1, os%norbs
424 do idim = 1, st%d%dim
425 weight(iorb,ik,ist) = weight(iorb,ik,ist) + st%kweights(ik)*abs(zdot(idim,iorb))**2
426 end do
427 end do
428 end if
429 end do
430 end do
431
432 if (st%parallel_in_states .or. st%d%kpt%parallel) then
433 call comm_allreduce(st%st_kpt_mpi_grp, weight)
434 end if
435
436 safe_deallocate_a(ddot)
437 safe_deallocate_a(zdot)
438
439 if (mpi_grp_is_root(mpi_world)) then
440 write(format_str,'(a,i2,a)') '(', os%norbs+2, 'f12.6)'
441 safe_allocate(tdos(1:os%norbs))
442 do ie = 1, this%epoints
443 energy = this%emin + (ie - 1) * this%de
444 do iorb = 1, os%norbs
445 tdos(iorb) = m_zero
446 do ist = 1, st%nst
447 do ik = 1, st%nik
448 tdos(iorb) = tdos(iorb) + weight(iorb,ik,ist) * m_one/m_pi * &
449 this%gamma / ((energy - st%eigenval(ist, ik))**2 + this%gamma**2)
450 end do
451 end do
452 end do
453
454 write(iunit(0), trim(format_str)) units_from_atomic(units_out%energy, energy), &
455 units_from_atomic(unit_one / units_out%energy, sum(tdos)), &
456 (units_from_atomic(unit_one / units_out%energy, tdos(iorb)), iorb=1,os%norbs)
457 end do
458 safe_deallocate_a(tdos)
459 call io_close(iunit(0))
460 end if
461
462 call orbitalset_end(os)
463 safe_deallocate_a(weight)
464
465 end do
466
467 end do
468
469 safe_deallocate_a(dpsi)
470 safe_deallocate_a(zpsi)
471 end if
472
473 safe_deallocate_a(iunit)
474 safe_deallocate_a(dos)
475
476 pop_sub(dos_write_dos)
477 end subroutine dos_write_dos
478
479 ! ---------------------------------------------------------
480 subroutine dos_write_jdos(this, dir, st, box, ions, mesh, hm, namespace)
481 type(dos_t), intent(in) :: this
482 character(len=*), intent(in) :: dir
483 type(states_elec_t), intent(in) :: st
484 class(box_t), intent(in) :: box
485 type(ions_t), target, intent(in) :: ions
486 class(mesh_t), intent(in) :: mesh
487 type(hamiltonian_elec_t), intent(in) :: hm
488 type(namespace_t), intent(in) :: namespace
489
490 integer :: ie, ik, val, cond, is, ns
491 integer, allocatable :: iunit(:)
492 real(real64) :: energy
493 real(real64) :: tjdos(1)
494 real(real64), allocatable :: jdos(:,:)
495 character(len=64) :: filename
496
497 push_sub(dos_write_jdos)
498
499 ! shortcuts
500 ns = 1
501 if (st%d%nspin == 2) ns = 2
502
503 if (mpi_grp_is_root(mpi_world)) then
504 ! space for state-dependent DOS
505 safe_allocate(jdos(1:this%epoints, 0:ns-1))
506 safe_allocate(iunit(0:ns-1))
507 jdos = m_zero
508
509 ! compute band/spin-resolved density of states
510 do val = 1, st%nst
511 do cond = val, st%nst
512 do ik = 1, st%nik, ns
513 do is = 0, ns-1
514 if(st%occ(val, ik+is) < m_epsilon) cycle
515 if(st%occ(cond, ik+is) > m_epsilon) cycle
516 do ie = 1, this%epoints
517 energy = (ie - 1) * this%de
518 ! sum up Lorentzians
519 jdos(ie, is) = jdos(ie, is) + st%kweights(ik+is) * m_one/m_pi * &
520 this%gamma / ((energy - (st%eigenval(cond, ik+is)-st%eigenval(val, ik+is)))**2 + this%gamma**2)
521 end do
522 end do
523 end do
524 end do
525 end do
526
527 ! for spin-polarized calculations also output spin-resolved tDOS
528 if (st%d%nspin > 1) then
529 do is = 0, ns-1
530 write(filename, '(a,i1.1,a)') 'total-jdos-', is+1,'.dat'
531 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
532 ! write header
533 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS (spin-resolved)'
534
535 do ie = 1, this%epoints
536 energy = (ie - 1) * this%de
537 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
538 units_from_atomic(unit_one / units_out%energy, jdos(ie, is))
539 call messages_info(1, iunit(is))
540 end do
541
542 call io_close(iunit(is))
543 end do
544 end if
545
546
547 iunit(0) = io_open(trim(dir)//'/'//'total-jdos.dat', namespace, action='write')
548 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS'
549
550 ! compute total joint density of states
551 do ie = 1, this%epoints
552 energy = (ie - 1) * this%de
553 tjdos(1) = m_zero
554 do is = 0, ns-1
555 tjdos(1) = tjdos(1) + jdos(ie, is)
556 end do
557 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
558 units_from_atomic(unit_one / units_out%energy, tjdos(1))
559 call messages_info(1, iunit(0))
560 end do
561
562 call io_close(iunit(0))
563 end if
564
565 safe_deallocate_a(iunit)
566 safe_deallocate_a(jdos)
567
568 pop_sub(dos_write_jdos)
569 end subroutine dos_write_jdos
570
571
572
574end module dos_oct_m
575
576!! Local Variables:
577!! mode: f90
578!! coding: utf-8
579!! End:
pure logical function, public accel_is_enabled()
Definition: accel.F90:400
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
subroutine, public zget_atomic_orbital(namespace, space, latt, pos, species, mesh, sm, ii, ll, jj, os, orbind, radius, d_dim, use_mesh, normalize, index_shift)
This routine returns the atomic orbital basis – provided by the pseudopotential structure in geo.
character(len=1), dimension(0:3), parameter, public l_notation
real(real64) function, public atomic_orbital_get_radius(species, mesh, iorb, ispin, truncation, threshold)
subroutine, public dget_atomic_orbital(namespace, space, latt, pos, species, mesh, sm, ii, ll, jj, os, orbind, radius, d_dim, use_mesh, normalize, index_shift)
This routine returns the atomic orbital basis – provided by the pseudopotential structure in geo.
subroutine, public dos_write_dos(this, dir, st, box, ions, mesh, hm, namespace)
Definition: dos.F90:252
subroutine, public dos_init(this, namespace, st, kpoints)
Definition: dos.F90:166
subroutine, public dos_write_jdos(this, dir, st, box, ions, mesh, hm, namespace)
Definition: dos.F90:574
integer, parameter, public spin_polarized
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_four
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_one
Definition: global.F90:188
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
integer pure function, public pad_pow2(size)
create array size, which is padded to powers of 2
Definition: math.F90:886
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:624
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
subroutine, public orbitalset_init(this)
Definition: orbitalset.F90:203
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:229
subroutine, public orbitalset_update_phase(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
Build the phase correction to the global phase in case the orbital crosses the border of the simulato...
Definition: orbitalset.F90:279
subroutine, public dorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
Definition: orbitalset.F90:521
subroutine, public zorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
integer function, public orbitalset_utils_count(species, iselect)
pure logical function, public states_are_real(st)
type(type_t), public type_cmplx
Definition: types.F90:134
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
class to tell whether a point is inside or outside
Definition: box.F90:141
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
int true(void)