Octopus
states_elec_io.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
22
25 use box_oct_m
26 use comm_oct_m
27 use debug_oct_m
29 use global_oct_m
30 use grid_oct_m
31 use io_oct_m
32 use ions_oct_m
33 use, intrinsic :: iso_fortran_env
35 use mesh_oct_m
38 use mpi_oct_m
42 use parser_oct_m
43 use phase_oct_m
45 use smear_oct_m
46 use sort_oct_m
47 use space_oct_m
53 use unit_oct_m
55 use utils_oct_m
56
57 implicit none
58
59 private
60
61 public :: &
66
67contains
68
69 ! ---------------------------------------------------------
71 !
72 subroutine states_elec_write_eigenvalues(nst, st, space, kpoints, error, st_start, compact, iunit, namespace)
73 integer, intent(in) :: nst
74 type(states_elec_t), intent(in) :: st
75 class(space_t), intent(in) :: space
76 type(kpoints_t), intent(in) :: kpoints
77 real(real64), optional, intent(in) :: error(:,:)
78 integer, optional, intent(in) :: st_start
79 logical, optional, intent(in) :: compact
80 integer, optional, intent(in) :: iunit
81 type(namespace_t), optional, intent(in) :: namespace
82
83 integer :: ik, ikk, ist, ns, is, idir, st_start_, iflat, iqn, homo_index, not_printed
84 logical :: print_eigenval
85 real(real64) :: kpoint(space%dim), max_error
86 character(len=120) :: tmp_str(max(space%dim, 3)), cspin
87
88 real(real64), allocatable :: flat_eigenval(:)
89 integer, allocatable :: flat_indices(:, :)
90 integer, parameter :: print_range = 8
91
93
94 st_start_ = 1
95 if (present(st_start)) st_start_ = st_start
96 assert(nst <= st%nst)
97
98
99 if (.not. mpi_grp_is_root(mpi_world)) then
101 return
102 end if
103
104 if (.not. optional_default(compact, .false.)) then
105
106 ns = 1
107 if (st%d%nspin == 2) ns = 2
108
109 message(1) = 'Eigenvalues [' // trim(units_abbrev(units_out%energy)) // ']'
110 call messages_info(1, iunit=iunit, namespace=namespace)
111
112 if (st%d%ispin == spinors) then
113 write(message(1), '(a4,1x,a5,1x,a12,1x,a12,2x,a4,4x,a4,4x,a4)') &
114 '#st',' Spin',' Eigenvalue', 'Occupation ', '<Sx>', '<Sy>', '<Sz>'
115 else
116 write(message(1), '(a4,1x,a5,1x,a12,4x,a12)') &
117 '#st',' Spin',' Eigenvalue', 'Occupation'
118 end if
119 if (present(error)) then
120 write(message(1),'(a,a10)') trim(message(1)), ' Error'
121 end if
122 call messages_info(1, iunit=iunit, namespace=namespace)
123
124 do ik = 1, st%nik, ns
125 if (space%is_periodic()) then
126 ikk = st%d%get_kpoint_index(ik)
127 kpoint(1:space%dim) = kpoints%get_point(ikk, absolute_coordinates = .false.)
128 write(message(1), '(a,i8,a)') '#k =', ikk, ', k = ('
129 do idir = 1, space%dim
130 write(tmp_str(1), '(f10.6)') kpoint(idir)
131 message(1) = trim(message(1))//trim(tmp_str(1))
132 if (idir < space%dim) message(1) = trim(message(1))//','
133 end do
134 message(1) = trim(message(1))//')'
135 call messages_info(1, iunit=iunit, namespace=namespace)
136 end if
137
138 do ist = st_start_, nst
139 do is = 0, ns-1
140 if (is == 0) cspin = 'up'
141 if (is == 1) cspin = 'dn'
142 if (st%d%ispin == unpolarized .or. st%d%ispin == spinors) cspin = '--'
143
144 write(tmp_str(1), '(i4,3x,a2)') ist, trim(cspin)
145 if (st%d%ispin == spinors) then
146 write(tmp_str(2), '(1x,f12.6,5x,f5.2,3x,3f8.4)') &
147 units_from_atomic(units_out%energy, st%eigenval(ist, ik)), st%occ(ist, ik), st%spin(1:3, ist, ik)
148 if (present(error)) write(tmp_str(3), '(a3,es8.1,a1)')' (', error(ist, ik), ')'
149 else
150 write(tmp_str(2), '(1x,f12.6,3x,f12.6)') &
151 units_from_atomic(units_out%energy, st%eigenval(ist, ik+is)), st%occ(ist, ik+is)
152 if (present(error)) write(tmp_str(3), '(a7,es8.1,a1)')' (', error(ist, ik+is), ')'
153 end if
154 if (present(error)) then
155 message(1) = trim(tmp_str(1))//trim(tmp_str(2))//trim(tmp_str(3))
156 else
157 message(1) = trim(tmp_str(1))//trim(tmp_str(2))
158 end if
159 call messages_info(1, iunit=iunit, namespace=namespace)
160 end do
161 end do
162 end do
163
164 else
166 call messages_info(1, iunit=iunit, namespace=namespace)
167
168 safe_allocate(flat_eigenval(1:st%nik*nst))
169 safe_allocate(flat_indices(1:2, 1:st%nik*nst))
170
171 iflat = 1
172 do iqn = 1, st%nik
173 do ist = 1, nst
174
175 flat_eigenval(iflat) = st%eigenval(ist, iqn)
176 flat_indices(1:2, iflat) = (/iqn, ist/)
177
178 iflat = iflat + 1
179 end do
180 end do
181
182 call sort(flat_eigenval, flat_indices(:, :))
183
184 homo_index = st%nik*nst
185 do iflat = 1, st%nik*nst
186 iqn = flat_indices(1, iflat)
187 ist = flat_indices(2, iflat)
188 if (abs(st%occ(ist, iqn)) < 0.1_real64) then
189 homo_index = iflat - 1
190 exit
191 end if
192 end do
193
194 tmp_str(1) = '# State'
195
196 if (space%is_periodic()) tmp_str(1) = trim(tmp_str(1))//' KPoint'
197
198 if (st%d%ispin == spin_polarized) tmp_str(1) = trim(tmp_str(1))//' Spin'
199
200 tmp_str(1) = trim(tmp_str(1))//' Eigenvalue ['// trim(units_abbrev(units_out%energy)) // ']'
201
202 tmp_str(1) = trim(tmp_str(1))//' Occupation'
203
204 if (st%d%ispin == spinors) then
205 tmp_str(1) = trim(tmp_str(1))//' <Sx> <Sy> <Sz>'
206 end if
207
208 if (present(error)) tmp_str(1) = trim(tmp_str(1))//' Error'
209
210 call messages_write(tmp_str(1))
211 call messages_info(iunit=iunit, namespace=namespace)
212
213 not_printed = 0
214 max_error = 0.0_real64
215 do iflat = 1, st%nik*nst
216 iqn = flat_indices(1, iflat)
217 ist = flat_indices(2, iflat)
218 ik = st%d%get_kpoint_index(iqn)
219 is = st%d%get_spin_index(iqn)
220
221 print_eigenval = iflat <= print_range
222 print_eigenval = print_eigenval .or. st%nik*nst - iflat < print_range
223 print_eigenval = print_eigenval .or. abs(iflat - homo_index) <= print_range
224
225 if (print_eigenval) then
226
227 if (not_printed > 0) then
228 call messages_write('')
229 call messages_new_line()
230 call messages_write(' [output of ')
231 call messages_write(not_printed)
232 call messages_write(' eigenvalues skipped')
233 if (present(error)) then
234 call messages_write(': maximum error =')
235 call messages_write(max_error, fmt = '(es7.1)', align_left = .true.)
236 end if
237 call messages_write(']')
238 call messages_new_line()
239 call messages_write('')
240 call messages_info(iunit=iunit, namespace=namespace)
241
242 not_printed = 0
243 max_error = 0.0_real64
244
245 end if
246
247 write(tmp_str(1), '(i7)') ist
248
249 if (space%is_periodic()) then
250 write(tmp_str(1), '(2a,i7)') trim(tmp_str(1)), ' ', ik
251 end if
252
253 if (st%d%ispin == spin_polarized) then
254 if (is == 1) cspin = ' up'
255 if (is == 2) cspin = ' dn'
256 write(tmp_str(1), '(2a,a5)') trim(tmp_str(1)), ' ', cspin
257 end if
258
259
260 if (len(units_abbrev(units_out%energy)) == 1) then
261 write(tmp_str(1), '(2a,f14.6)') trim(tmp_str(1)), ' ', units_from_atomic(units_out%energy, st%eigenval(ist, iqn))
262 else
263 write(tmp_str(1), '(2a,f15.6)') trim(tmp_str(1)), ' ', units_from_atomic(units_out%energy, st%eigenval(ist, iqn))
264 end if
265
266 write(tmp_str(1), '(2a,f11.6)') trim(tmp_str(1)), ' ', st%occ(ist, iqn)
267
268 if (st%d%ispin == spinors) then
269 write(tmp_str(1), '(2a,3f9.4)') trim(tmp_str(1)), ' ', st%spin(1:3, ist, iqn)
270 end if
271
272 if (present(error)) then
273 write(tmp_str(1), '(2a,es8.1,a)') trim(tmp_str(1)), ' (', error(ist, iqn), ')'
274 end if
275
276 call messages_write(tmp_str(1))
277 call messages_info(iunit=iunit, namespace=namespace)
278
279 else
280
281 not_printed = not_printed + 1
282
283 if (present(error)) then
284 max_error = max(max_error, error(ist, iqn))
285 end if
286
287 end if
288
289 end do
290
291 if (nst*st%nik > 1) call print_dos()
292
293 safe_deallocate_a(flat_indices)
294 safe_deallocate_a(flat_eigenval)
295
296 end if
297
298 call smear_write_info(st%smear, namespace, iunit)
299
301
302 contains
303
304 subroutine print_dos()
305
306 integer, parameter :: ndiv = 70, height = 10
307 integer :: histogram(1:ndiv), iev, ien, iline, ife
308 real(real64) :: dhistogram(1:ndiv)
309 character(len=ndiv) :: line
310 real(real64) :: emin, emax, de, maxhist
311
313
314 emin = flat_eigenval(1)
315 emax = flat_eigenval(st%nik*nst)
316 de = (emax - emin)/(ndiv - m_one)
317
318 if (de < m_epsilon) then
320 return
321 end if
322
323 ife = nint((st%smear%e_fermi - emin)/de) + 1
324
325 histogram = 0
326 dhistogram = m_zero
327 do iev = 1, st%nik*nst
328 ien = nint((flat_eigenval(iev) - emin)/de) + 1
329 assert(ien >= 1)
330 assert(ien <= ndiv)
331 dhistogram(ien) = dhistogram(ien) + st%kweights(flat_indices(1, iev))*kpoints%full%npoints
332 end do
333
334 !normalize
335 maxhist = maxval(dhistogram)
336 do ien = 1, ndiv
337 histogram(ien) = nint(dhistogram(ien)*height/maxhist)
338 end do
339
340 call messages_new_line()
341 call messages_write('Density of states:')
342 call messages_new_line()
343 call messages_info(iunit=iunit, namespace=namespace)
344
345 !print histogram
346 do iline = height, 1, -1
347 do ien = 1, ndiv
348 if (histogram(ien) >= iline) then
349 call messages_write('%')
350 else
351 call messages_write('-')
352 end if
353 end do
354 call messages_info(iunit=iunit, namespace=namespace)
355 end do
356
357 line(1:ndiv) = ' '
358 line(ife:ife) = '^'
359 call messages_write(line)
360 call messages_new_line()
361 call messages_info(iunit=iunit, namespace=namespace)
362
364 end subroutine print_dos
365
366 end subroutine states_elec_write_eigenvalues
367
368 ! ---------------------------------------------------------
370 !
371 subroutine states_elec_write_gaps(iunit, st, space)
372 integer, intent(in) :: iunit
373 type(states_elec_t), intent(in) :: st
374 class(space_t), intent(in) :: space
375
376 integer :: ik, ist
377
378 real(real64) :: homo, lumo, egdir, egindir
379 integer :: homok, lumok, egdirk
380
381 push_sub(states_elec_write_gaps)
382
383 if (.not. mpi_grp_is_root(mpi_world) .or. .not. space%is_periodic() &
384 .or. st%smear%method /= smear_semiconductor) then
386 return
387 end if
388
389 !We do not want to compute the bandgap in case there is no extra states
390 if (ceiling(st%qtot/st%smear%el_per_state) == st%nst) then
392 return
393 end if
394
395 homo = -m_huge
396 homok = -1
397 lumo = m_huge
398 lumok = -1
399 egdir = m_huge
400 egdirk = -1
401
402 !TODO: This definition depends on the type of smearing
403 do ik = 1, st%nik
404 if (abs(st%kweights(ik)) < m_epsilon) cycle
405 do ist = 1,st%nst-1
406 if (st%occ(ist, ik) > m_epsilon .and. st%eigenval(ist, ik) > homo) then
407 homo = st%eigenval(ist, ik)
408 homok = ik
409 end if
410
411 if (st%occ(ist+1, ik) <= m_epsilon .and. st%eigenval(ist+1, ik) < lumo) then
412 lumo = st%eigenval(ist+1, ik)
413 lumok = ik
414 end if
415
416 if (st%occ(ist, ik) > m_epsilon .and. st%occ(ist+1, ik) <= m_epsilon &
417 .and. (st%eigenval(ist+1, ik) - st%eigenval(ist, ik)) < egdir) then
418 egdir = (st%eigenval(ist+1, ik) - st%eigenval(ist, ik))
419 egdirk = ik
420 end if
421 end do
422 end do
423
424 egindir = lumo-homo
425
426 if (lumok == -1 .or. egdir <= m_epsilon) then
427 write(message(1),'(a)') 'The system seems to have no gap.'
428 call messages_info(1, iunit)
429 else
430 write(message(1),'(a,i5,a,f7.4,a)') 'Direct gap at ik=', egdirk, ' of ', &
431 units_from_atomic(units_out%energy, egdir), ' ' // trim(units_abbrev(units_out%energy))
432 write(message(2),'(a,i5,a,i5,a,f7.4,a)') 'Indirect gap between ik=', homok, ' and ik=', lumok, &
433 ' of ', units_from_atomic(units_out%energy, egindir), ' ' // trim(units_abbrev(units_out%energy))
434 call messages_info(2, iunit)
435 end if
436
437
439 end subroutine states_elec_write_gaps
440
441
442 ! ---------------------------------------------------------
446 !
447 subroutine states_elec_write_tpa(dir, namespace, space, mesh, st)
448 character(len=*), intent(in) :: dir
449 type(namespace_t), intent(in) :: namespace
450 class(space_t), intent(in) :: space
451 class(mesh_t), intent(in) :: mesh
452 type(states_elec_t), intent(in) :: st
453
454 type(block_t) :: blk
455 integer :: ncols, icoord, ist, ik, tpa_initialst, tpa_initialk
456 integer :: iunit
457
458 real(real64), allocatable :: ff(:)
459 complex(real64), allocatable :: cff(:)
460 real(real64), allocatable :: osc(:)
461 real(real64) :: transition_energy, osc_strength, dsf
462
463 logical :: use_qvector
464 real(real64), allocatable :: qvector(:), psi_initial(:, :), psi_ist(:, :)
465
466 push_sub(states_elec_write_tpa)
467
468 use_qvector = .false.
469
470 ! find the orbital with half-occupation
471 tpa_initialst = -1
472 do ist = 1, st%nst
473 do ik = 1, st%nik
474 if (abs(st%occ(ist,ik)-m_half) < m_min_occ) then
475 tpa_initialst = ist
476 tpa_initialk = ik
477 end if
478 end do
479 end do
480
481 ! make sure that half-occupancy was found
482 if (tpa_initialst == -1) then
483 if (mpi_grp_is_root(mpi_world)) then
484
485 call messages_write('No orbital with half-occupancy found. TPA output is not written.')
486 call messages_warning(namespace=namespace)
487
488 end if
489 pop_sub(states_elec_write_tpa)
490 return
491 end if
492
493 !%Variable MomentumTransfer
494 !%Type block
495 !%Section Output
496 !%Description
497 !% Momentum-transfer vector <math>\vec{q}</math> to be used when calculating matrix elements
498 !% <math>\left< f \left| e^{i \vec{q} \cdot \vec{r}} \right| i \right></math>.
499 !% This enables the calculation of the dynamical structure factor,
500 !% which is closely related to generalized oscillator strengths.
501 !% If the vector is not given, but TPA output is requested (<tt>Output = TPA</tt>),
502 !% only the oscillator strengths are written in the output file.
503 !% For example, to use <math>\vec{q}</math> = (0.1, 0.2, 0.3), set
504 !%
505 !% <tt>%MomentumTransfer
506 !% <br>&nbsp;&nbsp; 0.1 | 0.2 | 0.3
507 !% <br>%</tt>
508 !%End
509 if (parse_block(namespace, 'MomentumTransfer', blk) == 0) then
510
511 ! check if input makes sense
512 ncols = parse_block_cols(blk, 0)
513
514 if (ncols /= space%dim) then ! wrong size
515
516 if (mpi_grp_is_root(mpi_world)) then
517 call messages_write('Inconsistent size of momentum-transfer vector. It will not be used in the TPA calculation.')
518 call messages_warning(namespace=namespace)
519 end if
520
521 else ! correct size
522
523 use_qvector = .true.
524 safe_allocate(qvector(1:space%dim))
525
526 do icoord = 1, space%dim !for x,y,z
527 call parse_block_float(blk, 0, icoord-1, qvector(icoord))
528 qvector(icoord) = units_to_atomic(unit_one / units_inp%length, qvector(icoord))
529 end do
530
531 end if
532
533 end if
534
535 ! calculate the matrix elements
536
537 safe_allocate(ff(1:mesh%np))
538 if (use_qvector) then
539 safe_allocate(cff(1:mesh%np))
540 end if
541 safe_allocate(osc(1:space%dim))
542
543 ! root writes output to file
544
545 if (mpi_grp_is_root(mpi_world)) then
546
547 iunit = io_open(trim(dir)//'/'//trim('tpa_xas'), namespace, action='write')
548
549 ! header
550 if (use_qvector) then
551 write (message(1),'(a1,a30,3(es14.5,1x),a1)') '#', ' momentum-transfer vector : (', &
552 (units_from_atomic(unit_one / units_out%length, qvector(icoord)), icoord=1, space%dim),')'
553 select case (space%dim)
554 case (1)
555 write(message(2), '(a1,4(a15,1x))') '#', 'E' , '<x>', '<f>', 'S(q,omega)'
556 case (2)
557 write(message(2), '(a1,5(a15,1x))') '#', 'E' , '<x>', '<y>', '<f>', 'S(q,omega)'
558 case (3)
559 write(message(2), '(a1,6(a15,1x))') '#', 'E' , '<x>', '<y>', '<z>', '<f>', 'S(q,omega)'
560 end select
561 call messages_info(2, iunit, namespace=namespace)
562 else
563 select case (space%dim)
564 case (1)
565 write(message(1), '(a1,3(a15,1x))') '#', 'E' , '<x>', '<f>'
566 case (2)
567 write(message(1), '(a1,4(a15,1x))') '#', 'E' , '<x>', '<y>', '<f>'
568 case (3)
569 write(message(1), '(a1,5(a15,1x))') '#', 'E' , '<x>', '<y>', '<z>', '<f>'
570 end select
571 call messages_info(1, iunit, namespace=namespace)
572 end if
573
574 end if
575
576 safe_allocate(psi_initial(1:mesh%np, 1:st%d%dim))
577 safe_allocate(psi_ist(1:mesh%np, 1:st%d%dim))
578
579 call states_elec_get_state(st, mesh, tpa_initialst, tpa_initialk, psi_initial)
580
581 ! loop through every state
582 do ist = 1, st%nst
583
584 call states_elec_get_state(st, mesh, ist, tpa_initialk, psi_ist)
585
586 ! final states are the unoccupied ones
587 if (abs(st%occ(ist,tpa_initialk)) < m_min_occ) then
588
589 osc_strength = m_zero
590 transition_energy = st%eigenval(ist, tpa_initialk) - st%eigenval(tpa_initialst, tpa_initialk)
591
592 ! dipole matrix elements <f|x|i> etc. -> oscillator strengths
593 do icoord = 1, space%dim ! for x,y,z
594
595 ff(1:mesh%np) = psi_initial(1:mesh%np, 1)*mesh%x(1:mesh%np, icoord)* &
596 psi_ist(1:mesh%np, 1)
597 osc(icoord) = dmf_integrate(mesh, ff)
598 osc_strength = osc_strength + m_two/real(space%dim, real64) *transition_energy*abs(osc(icoord))**2
599
600 end do
601
602 ! matrix elements <f|exp(iq.r)|i> -> dynamic structure factor
603 if (use_qvector) then
604
605 cff(1:mesh%np) = psi_initial(1:mesh%np, 1)*psi_ist(1:mesh%np, 1)
606
607 do icoord = 1, space%dim ! for x,y,z
608 cff(1:mesh%np) = cff(1:mesh%np)*exp(m_zi*mesh%x(1:mesh%np, icoord)*qvector(icoord))
609 end do
610
611 dsf = abs(zmf_integrate(mesh, cff))**2
612 end if
613
614 ! write oscillator strengths (+ dynamic structure factor if qvector is given) into file
615 if (mpi_grp_is_root(mpi_world)) then
616
617 if (use_qvector) then
618 write(message(1), '(1x,6(es15.8,1x))') units_from_atomic(units_out%energy, transition_energy), osc(:), osc_strength, &
620 else
621 write(message(1), '(1x,6(es15.8,1x))') units_from_atomic(units_out%energy, transition_energy), osc(:), osc_strength
622 end if
623
624 call messages_info(1, iunit, namespace=namespace)
625
626 end if
627
628 end if
629
630 end do
631
632 ! finally close the file
633 if (mpi_grp_is_root(mpi_world)) then
634 call io_close(iunit)
635 end if
636
637 safe_deallocate_a(psi_initial)
638 safe_deallocate_a(psi_ist)
639
640 safe_deallocate_a(ff)
641 if (use_qvector) then
642 safe_deallocate_a(cff)
643 end if
644 safe_deallocate_a(osc)
645 if (use_qvector) then
646 safe_deallocate_a(qvector)
647 end if
648
649 pop_sub(states_elec_write_tpa)
650
651 end subroutine states_elec_write_tpa
652
653
654 ! ---------------------------------------------------------
656 !
657 subroutine states_elec_write_bandstructure(dir, namespace, nst, st, ions, mesh, kpoints, &
658 phase, vec_pot, vec_pot_var)
659 character(len=*), intent(in) :: dir
660 type(namespace_t), intent(in) :: namespace
661 integer, intent(in) :: nst
662 type(states_elec_t), intent(in) :: st
663 type(ions_t), target, intent(in) :: ions
664 class(mesh_t), intent(in) :: mesh
665 type(kpoints_t), intent(in) :: kpoints
666 type(phase_t), intent(in) :: phase
667 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
668 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
670
671 integer :: idir, ist, ik, ns, is,npath
672 integer, allocatable :: iunit(:)
673 real(real64) :: red_kpoint(mesh%box%dim)
674 character(len=80) :: filename
675
676 logical :: projection
677 integer :: ii, ll, mm, nn, work, norb, work2
678 integer :: ia, iorb, idim, maxnorb
679 real(real64), allocatable :: dpsi(:,:), ddot(:,:)
680 complex(real64), allocatable :: zpsi(:,:), zdot(:,:)
681 real(real64), allocatable :: weight(:,:,:,:,:)
682 type(orbitalset_t) :: os
683
684
686
687 ! shortcuts
688 ns = 1
689 if (st%d%nspin == 2) ns = 2
690
691 !%Variable BandStructureComputeProjections
692 !%Type logical
693 !%Default false
694 !%Section Output
695 !%Description
696 !% Determines if projections of wavefunctions on the atomic orbitals
697 !% are computed or not for obtaining the orbital resolved band-structure.
698 !%End
699 call parse_variable(namespace, 'BandStructureComputeProjections', .false., projection)
700
701
702 if (mpi_grp_is_root(mpi_world)) then
703 safe_allocate(iunit(0:ns-1))
704
705 do is = 0, ns-1
706 if (ns > 1) then
707 write(filename, '(a,i1.1,a)') 'bandstructure-sp', is+1
708 else
709 write(filename, '(a)') 'bandstructure'
710 end if
711 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
712
713 ! write header
714 write(iunit(is),'(a)',advance='no') '# coord. '
715 do idir = 1, mesh%box%dim
716 write(iunit(is),'(3a)',advance='no') 'k', index2axis(idir), ' '
717 end do
718 if (.not. projection) then
719 write(iunit(is),'(a,i6,3a)') '(red. coord.), bands:', nst, ' [', trim(units_abbrev(units_out%energy)), ']'
720 else
721 write(iunit(is),'(a,i6,3a)',advance='no') '(red. coord.), bands:', nst, ' [', trim(units_abbrev(units_out%energy)), '] '
722 do ia = 1, ions%natoms
723 work = orbitalset_utils_count(ions%atom(ia)%species)
724 do norb = 1, work
725 work2 = orbitalset_utils_count(ions%atom(ia)%species, norb)
726 write(iunit(is),'(a, i3.3,a,i1.1,a)',advance='no') 'w(at=',ia,',os=',norb,') '
727 end do
728 end do
729 write(iunit(is),'(a)') ''
730 end if
731 end do
732 end if
733
734 npath = kpoints%nkpt_in_path()*ns
735
736
737 !We need to compute the projections of each wavefunctions on the localized basis
738 if (projection) then
739
740 if (states_are_real(st)) then
741 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
742 else
743 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
744 end if
745
746 maxnorb = 0
747 do ia = 1, ions%natoms
748 maxnorb = max(maxnorb, orbitalset_utils_count(ions%atom(ia)%species))
749 end do
751 safe_allocate(weight(1:st%nik,1:st%nst, 1:maxnorb, 1:max_l, 1:ions%natoms))
752 weight(1:st%nik,1:st%nst, 1:maxnorb, 1:max_l, 1:ions%natoms) = m_zero
753
754 do ia = 1, ions%natoms
755
756 !We first count how many orbital set we have
757 work = orbitalset_utils_count(ions%atom(ia)%species)
758
759 !We loop over the orbital sets of the atom ia
760 do norb = 1, work
761 call orbitalset_init(os)
762
763 !We count the orbitals
764 work2 = 0
765 do iorb = 1, ions%atom(ia)%species%get_niwfs()
766 call ions%atom(ia)%species%get_iwf_ilm(iorb, 1, ii, ll, mm)
767 call ions%atom(ia)%species%get_iwf_n(iorb, 1, nn)
768 if (ii == norb) then
769 os%ll = ll
770 os%nn = nn
771 os%ii = ii
772 os%radius = atomic_orbital_get_radius(ions%atom(ia)%species, mesh, iorb, 1, option__aotruncation__ao_full, 0.01_real64)
773 work2 = work2 + 1
774 end if
775 end do
776 os%norbs = work2
777 os%ndim = 1
778 os%use_submesh = .false.
779 os%spec => ions%atom(ia)%species
780
781 do iorb = 1, os%norbs
782 ! We obtain the orbital
783 if (states_are_real(st)) then
784 call dget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
785 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
786 os, iorb, os%radius, os%ndim, use_mesh=.not.os%use_submesh, &
787 normalize = .true.)
788 else
789 call zget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
790 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
791 os, iorb, os%radius, os%ndim, use_mesh=.not.os%use_submesh, &
792 normalize = .true.)
793 end if
794 end do !iorb
795
796 if (phase%is_allocated()) then
797 ! In case of complex wavefunction, we allocate the array for the phase correction
798 safe_allocate(os%phase(1:os%sphere%np, st%d%kpt%start:st%d%kpt%end))
799 os%phase(:,:) = m_zero
800 if (.not. os%use_submesh) then
801 safe_allocate(os%eorb_mesh(1:mesh%np, 1:os%norbs, 1:os%ndim, st%d%kpt%start:st%d%kpt%end))
802 os%eorb_mesh(:,:,:,:) = m_zero
803 else
804 safe_allocate(os%eorb_submesh(1:os%sphere%np, 1:os%ndim, 1:os%norbs, st%d%kpt%start:st%d%kpt%end))
805 os%eorb_submesh(:,:,:,:) = m_zero
806 end if
807 call orbitalset_update_phase(os, mesh%box%dim, st%d%kpt, kpoints, (st%d%ispin==spin_polarized), &
808 vec_pot, vec_pot_var)
809 end if
810
811 if (states_are_real(st)) then
812 safe_allocate(ddot(1:st%d%dim,1:os%norbs))
813 else
814 safe_allocate(zdot(1:st%d%dim,1:os%norbs))
815 end if
816
817 do ist = st%st_start, st%st_end
818 do ik = st%d%kpt%start, st%d%kpt%end
819 if (ik < st%nik - npath + 1) cycle ! We only want points inside the k-point path
820 if (states_are_real(st)) then
821 call states_elec_get_state(st, mesh, ist, ik, dpsi)
822 call dorbitalset_get_coefficients(os, st%d%dim, dpsi, ik, .false., ddot(1:st%d%dim,1:os%norbs))
823 do iorb = 1, os%norbs
824 do idim = 1, st%d%dim
825 weight(ik,ist,iorb,norb,ia) = weight(ik,ist,iorb,norb,ia) + abs(ddot(idim,iorb))**2
826 end do
827 end do
828 else
829 call states_elec_get_state(st, mesh, ist, ik, zpsi)
830 if (phase%is_allocated()) then
831 ! Apply the phase that contains both the k-point and vector-potential terms.
832 call phase%apply_to_single(zpsi, mesh%np, st%d%dim, ik, .false.)
833 end if
834 call zorbitalset_get_coefficients(os, st%d%dim, zpsi, ik, phase%is_allocated(), &
835 zdot(1:st%d%dim,1:os%norbs))
836 do iorb = 1, os%norbs
837 do idim = 1, st%d%dim
838 weight(ik,ist,iorb,norb,ia) = weight(ik,ist,iorb,norb,ia) + abs(zdot(idim,iorb))**2
839 end do
840 end do
841 end if
842 end do
843 end do
844
845 safe_deallocate_a(ddot)
846 safe_deallocate_a(zdot)
847
848 call orbitalset_end(os)
849 end do !norb
850
851 if (st%parallel_in_states .or. st%d%kpt%parallel) then
852 call comm_allreduce(st%st_kpt_mpi_grp, weight(1:st%nik,1:st%nst, 1:maxnorb, 1:max_l,ia))
853 end if
854 end do !ia
855
856 safe_deallocate_a(dpsi)
857 safe_deallocate_a(zpsi)
858
859 end if !projection
860
861 if (mpi_grp_is_root(mpi_world)) then
862 ! output bands
863 do ik = st%nik-npath+1, st%nik, ns
864 do is = 0, ns - 1
865 red_kpoint(1:mesh%box%dim) = kpoints%get_point(st%d%get_kpoint_index(ik + is), absolute_coordinates=.false.)
866 write(iunit(is),'(1x)',advance='no')
867 if (st%nik > npath) then
868 write(iunit(is),'(f14.8)',advance='no') kpoints_get_path_coord(kpoints, &
869 st%d%get_kpoint_index(ik + is)-st%d%get_kpoint_index(st%nik -npath))
870 else
871 write(iunit(is),'(f14.8)',advance='no') kpoints_get_path_coord(kpoints, &
872 st%d%get_kpoint_index(ik + is))
873 end if
874 do idir = 1, mesh%box%dim
875 write(iunit(is),'(f14.8)',advance='no') red_kpoint(idir)
876 end do
877 do ist = 1, nst
878 write(iunit(is),'(1x,f14.8)',advance='no') units_from_atomic(units_out%energy, st%eigenval(ist, ik + is))
879 end do
880 if (projection) then
881 do ia = 1, ions%natoms
882 work = orbitalset_utils_count(ions%atom(ia)%species)
883 do norb = 1, work
884 work2 = orbitalset_utils_count(ions%atom(ia)%species, norb)
885 do iorb = 1, work2
886 do ist = 1, nst
887 write(iunit(is),'(es15.8)',advance='no') weight(ik+is,ist,iorb,norb,ia)
888 end do
889 end do
890 end do
891 end do
892 end if
893 end do
894
895 do is = 0, ns-1
896 write(iunit(is), '(a)')
897 end do
898 end do
899
900 do is = 0, ns-1
901 call io_close(iunit(is))
902 end do
903 end if
904
905 if (projection) then
906 safe_deallocate_a(weight)
907 end if
908
909 safe_deallocate_a(iunit)
910
911
913
915
916end module states_elec_io_oct_m
917
918
919!! Local Variables:
920!! mode: f90
921!! coding: utf-8
922!! End:
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:149
double exp(double __x) __attribute__((__nothrow__
integer, parameter, public max_l
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.
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.
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_huge
Definition: global.F90:205
real(real64), parameter, public m_zero
Definition: global.F90:187
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
real(real64), parameter, public m_min_occ
Minimal occupation that is considered to be non-zero.
Definition: global.F90:210
This module implements the underlying real-space grid.
Definition: grid.F90:117
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
real(real64) pure function, public kpoints_get_path_coord(this, ind)
Definition: kpoints.F90:1656
This module defines various routines, operating on mesh functions.
complex(real64) function, public zmf_integrate(mesh, ff, mask, reduce)
Integrate a function on the mesh.
real(real64) function, public dmf_integrate(mesh, ff, mask, reduce)
Integrate a function on the mesh.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
subroutine, public messages_new_line()
Definition: messages.F90:1146
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
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)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public smear_write_info(this, namespace, iunit)
Definition: smear.F90:836
integer, parameter, public smear_semiconductor
Definition: smear.F90:171
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
This module defines routines to write information about states.
subroutine, public states_elec_write_tpa(dir, namespace, space, mesh, st)
calculate the transition potential and write to a file
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_gaps(iunit, st, space)
calculate gaps and write 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
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
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
character pure function, public index2axis(idir)
Definition: utils.F90:202
subroutine print_dos()
Describes mesh distribution to nodes.
Definition: mesh.F90:186
A container for the phase.
Definition: phase.F90:177
The states_elec_t class contains all electronic wave functions.
int true(void)