Octopus
photoelectron_spectrum.F90
Go to the documentation of this file.
1!! Copyright (C) 2011 U. De Giovannini
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
27 use debug_oct_m
29 use global_oct_m
32 use io_oct_m
33 use ions_oct_m
37 use parser_oct_m
43 use sort_oct_m
44 use space_oct_m
45 use string_oct_m
48 use unit_oct_m
50 use utils_oct_m
51 use mpi_oct_m
52
53 implicit none
54
55 type pesoutput_t
56 logical :: what(MAX_OUTPUT_TYPES)
57 integer(int64) :: how(0:MAX_OUTPUT_TYPES)
58 integer :: output_interval(0:MAX_OUTPUT_TYPES)
59 real(real64) :: pol(3)
60 real(real64) :: pvec(3)
61 end type pesoutput_t
62
63 integer :: ierr, integrate
64 integer :: dim, dir, idim, pdim
65 logical :: what(MAX_OUTPUT_TYPES)
66 integer :: llp(3), llpp(3)
67 real(real64) :: Emax, Emin, Estep, uEstep,uEspan(2), pol(3)
68 real(real64) :: uThstep, uThspan(2), uPhstep, uPhspan(2), pvec(3)
69 real(real64) :: center(3)
70 real(real64), allocatable :: Lg(:,:), RR(:)
71 real(real64), allocatable :: pmesh(:,:,:,:)
72 integer, allocatable :: Lp(:,:,:,:,:)
73 logical :: need_pmesh, resolve_states
74 integer :: ii, i1,i2,i3, idxZero(1:3), st_range(2)
75 type(block_t) :: blk
76
77 type(electron_space_t) :: space
78 type(ions_t), pointer :: ions
79 type(states_elec_t) :: st
80 type(kpoints_t) :: kpoints
81 type(restart_t) :: restart
82
83 character(len=512) :: filename
84
85 logical :: have_zweight_path, use_zweight_path
86 integer :: krng(2), nkpt, kpth_dir
87
88 type(pesoutput_t) :: pesout
89
90 integer :: ist, ispin
91 real(real64), pointer :: pesP_out(:,:,:)
92 real(real64), allocatable, target :: pesP(:,:,:,:)
93 real(real64), allocatable :: Ekin(:,:,:)
94
95 type(pes_flux_t) :: pflux
96 integer :: pes_method, option
97
98 type(calc_mode_par_t) :: calc_mode
99 type(multicomm_t) :: mc
100 integer(int64) :: index_range(4)
101 integer :: nspin
102
103 class(coordinate_system_t), pointer :: coord_system
104
105 call getopt_init(ierr)
106 if (ierr /= 0) then
107 message(1) = "Your Fortran compiler doesn't support command-line arguments;"
108 message(2) = "the oct-photoelectron-spectrum command is not available."
109 call messages_fatal(2)
110 end if
111
113
115
116 call messages_init()
117 call io_init()
119
120 !* In order to initialize k-points
122
124 ions => ions_t(global_namespace)
125
126 ! we need k-points for periodic systems
127 call kpoints_init(kpoints, global_namespace, ions%symm, space%dim, space%periodic_dim, ions%latt)
128 call states_elec_init(st, global_namespace, space, ions%val_charge(), kpoints)
129 !*
130
131 !Initialize variables
132 llp(:) = 1
133 llpp(:) = 1
134
135 need_pmesh = .false.
136 dim = space%dim ! The dimensionality dim = [1,2,3]
137 pdim = space%periodic_dim
138
139 call messages_print_with_emphasis(msg="Postprocessing", namespace=global_namespace)
140
141 !Figure out which method has been used to calculate the photoelectron data
142 call parse_variable(global_namespace, 'PhotoElectronSpectrum', option__photoelectronspectrum__none, pes_method)
143
144 select case (pes_method)
145 case (option__photoelectronspectrum__pes_mask)
146 call messages_write('Will process mask-method data.')
147 call messages_new_line()
149 ! Note that Lg(:,:) is allocated inside pes_mask_read_info
150 call pes_mask_read_info("td.general/", global_namespace, dim, emax, estep, llp(:), lg, rr)
151 ! Keep a copy the original dimensions vector
152 llpp(1:dim) = llp(1:dim)
154 call messages_write('Read PES_MASK info file.')
155 call messages_info()
156
157 need_pmesh = space%is_periodic()
158
159
160 case (option__photoelectronspectrum__pes_flux)
161 call messages_write('Will process flux-method data.')
162 call messages_new_line()
163 call messages_info()
164
165 option = option__pes_flux_shape__sph
166 if (dim <= 2) option = option__pes_flux_shape__cub
167 if (space%is_periodic()) option = option__pes_flux_shape__pln
168
169 call parse_variable(global_namespace, 'PES_Flux_Shape', option, pflux%surf_shape)
170 call pes_flux_reciprocal_mesh_gen(pflux, global_namespace, space, st, kpoints, mpi_comm_null, post = .true.)
171
172 llpp(1:dim) = pflux%ll(1:dim)
173 need_pmesh = .true.
174
175
176 case (option__photoelectronspectrum__pes_spm)
177 call messages_not_implemented('Postprocessing SPM data')
178 call messages_fatal()
179
180 case default
181 call messages_write('Could not find any photoelectron data')
182 call messages_fatal()
183
184 end select
185
186
187 !set default values
188
189
190 integrate = integrate_none
191 uestep = -1
192 uespan = (/-1,-1/)
193 uthstep = -1
194 uthspan = (/-1,-1/)
195 uphstep = -1
196 uphspan = (/-1,-1/)
197 center = (/0,0,0/)
198 pvec = (/1,0,0/)
199
200 what(option__photoelectronspectrumoutput__energy_tot) = .true.
201 pesout%pvec = pvec
202
203 have_zweight_path = kpoints%have_zero_weight_path()
204 use_zweight_path = have_zweight_path
205
206 call get_laser_polarization(pol)
207 ! if there is no laser set pol along the z-axis
208 if (sum(pol(1:3)**2) <= m_epsilon) pol = (/0,0,1/)
209
210 ! more defaults
211 if (space%is_periodic()) then
212 what(:) = .false.
213
214 if (dim == 2) then
215 ! write the velocity map on plane pz=0 as it contains all the informations
216 pol = (/0,1,0/)
217 pvec = (/0,0,1/)
218
219 what(option__photoelectronspectrumoutput__velocity_map_cut) = .true.
220 pesout%pol = (/0,1,0/)
221 pesout%pvec = (/0,0,1/)
222 end if
223
224 if (dim == 3) then
225 ! write the full ARPES in vtk format (this could be a big file)
226 pol = (/0,0,1/)
227
228 what(option__photoelectronspectrumoutput__arpes) = .true.
229 pesout%pol = (/0,0,1/)
230 end if
231
232 if (have_zweight_path) then
233 ! In this case the output is well defined only on a path in reciprocal space
234 ! so we are going to have only a 2D slice regardless of dim=2 or 3
235 pol = (/0,0,1/)
236 pvec = (/0,1,0/)
237
238 what(:) = .false.
239 what(option__photoelectronspectrumoutput__arpes_cut) = .true.
240 pesout%pol = (/0,0,1/)
241 pesout%pvec = (/0,1,0/)
242 end if
243
244 end if
245
246
247 ! Intialize mc
248 call calc_mode%unset_parallelization(p_strategy_domains)
249 call multicomm_init(mc, global_namespace, mpi_world, calc_mode, &
250 mpi_world%size, index_range, (/ 5000, 1, 1, 1 /))
251 index_range(:) = 100000
252
253
256 if (ierr /= 0) then
257 message(1) = "Unable to read time-dependent restart information."
258 call messages_fatal(1)
259 end if
260
261 !%Variable PhotoelectronSpectrumOutput
262 !%Type block
263 !%Default none
264 !%Section Utilities::oct-photoelectron_spectrum
265 !%Description
266 !% Specifies what to output extracting the photoelectron cross-section informations.
267 !% When we use polar coordinates the zenith axis is set by vec (default is the first
268 !% laser field polarization vector), theta is the inclination angle measured from
269 !% vec (from 0 to \pi), and phi is the azimuthal angle on a plane perpendicular to
270 !% vec (from 0 to 2\pi).
271 !% Each option must be in a separate row. Optionally individual output formats can be defined
272 !% for each row or they can be read separately from <tt>OutputFormat</tt> variable
273 !% in the input file.
274 !%
275 !% Example (minimal):
276 !% <br><br><tt>%PhotoelectronSpectrumOutput
277 !% <br>&nbsp;&nbsp;energy_tot
278 !% <br>&nbsp;&nbsp;velocity_map
279 !% <br>%<br></tt>
280 !%
281 !% Example (with OutputFormat):
282 !% <br><br><tt>%PhotoelectronSpectrumOutput
283 !% <br>&nbsp;&nbsp;arpes | vtk
284 !% <br>&nbsp;&nbsp;velocity_map | ncdf
285 !% <br>%<br></tt>
286 !%
287 !%Option energy_tot 1
288 !% Output the energy-resolved photoelectron spectrum: E.
289 !%Option energy_angle 2
290 !% Output the energy and angle resolved spectrum: (theta, E)
291 !% The result is integrated over phi.
292 !%Option velocity_map_cut 3
293 !% Velocity map on a plane orthogonal to pvec: (px, py). The allowed cutting planes
294 !% (pvec) can only be parallel to the x,y,z=0 planes.
295 !% Space is oriented so that the z-axis is along vec. Supports the -I option.
296 !%Option energy_xy 4
297 !% Angle and energy-resolved spectrum on the inclination plane: (Ex, Ey).
298 !% The result is integrated over ph;
299 !%Option energy_th_ph 5
300 !% Ionization probability integrated on spherical cuts: (theta, phi).
301 !%Option velocity_map 6
302 !% Full momentum-resolved ionization probability: (px, py, pz).
303 !% The output format can be controlled with <tt>OutputHow</tt> and can be vtk, ncdf or ascii.
304 !%Option arpes 7
305 !% Full ARPES for semi-periodic systems (vtk).
306 !%Option arpes_cut 8
307 !% ARPES cut on a plane following a zero-weight path in reciprocal space.
308 !%End
309 pesout%what = what
310 call io_function_read_what_how_when(global_namespace, space, pesout%what, pesout%how, pesout%output_interval, &
311 what_tag_in = 'PhotoelectronSpectrumOutput', ignore_error = .true.)
312
313 ! TODO: I think it would be better to move these options in the
314 ! input file to have more flexibility to combine and to keep
315 ! track of them. UDG
316 ! Read options from command line
317 call getopt_photoelectron_spectrum(uestep, uespan,&
318 uthstep, uthspan, uphstep, &
319 uphspan, pol, center, pvec, integrate)
320
321 call getopt_end()
322
323
324 !! set user values
325 if (uestep > 0 .and. uestep > estep) estep = uestep
326 if (uespan(1) > 0) emin = uespan(1)
327 if (uespan(2) > 0) emax = uespan(2)
328
329
330
331 !%Variable PhotoelectronSpectrumResolveStates
332 !%Type block
333 !%Default
334 !%Section Utilities::oct-photoelectron_spectrum
335 !%Description
336 !% If <tt>yes</tt> calculate the photoelectron spectrum resolved in each K.S. state.
337 !% Optionally a range of states can be given as two slot block where the
338 !% first slot is the lower state index and the second is the highest one.
339 !% For example to calculate the spectra from state i to state j:
340 !%
341 !% <tt>%PhotoelectronSpectrumResolveStates
342 !% <br> i | j
343 !% <br>%</tt>
344 !%End
345 st_range(1:2) = (/1, st%nst/)
346 resolve_states = .false.
347 if (parse_block(global_namespace, 'PhotoelectronSpectrumResolveStates', blk) == 0) then
348 if (parse_block_cols(blk,0) < 2) call messages_input_error(global_namespace, 'PhotoelectronSpectrumResolveStates')
349 do idim = 1, 2
350 call parse_block_integer(blk, 0, idim - 1, st_range(idim))
351 end do
352 call parse_block_end(blk)
353 if (abs(st_range(2)-st_range(1)) > 0)resolve_states = .true.
354 else
355 call parse_variable(global_namespace, 'PhotoelectronSpectrumResolveStates', .false., resolve_states)
356 end if
357
358
359 ! Determine the range of k-points used in the PES calculation
360 krng(1) = 1
361 krng(2) = kpoints_number(kpoints)
362
363 if (have_zweight_path) then
364 if (pesout%what(option__photoelectronspectrumoutput__arpes_cut)) then
365 !Use the path only when asked for ARPES on a cutting curve in reciprocal space(the path)
366 use_zweight_path = .true.
367
368 krng(1) = kpoints_number(kpoints) - kpoints%nik_skip + 1
369
370 call messages_print_with_emphasis(msg="Kpoint selection", namespace=global_namespace)
371 write(message(1), '(a)') 'Will use a zero-weight path in reciprocal space with the following points'
372 call messages_info(1)
373
374 kpth_dir = 1
375 pvec = (/0,1,0/)
376
377
378 else
379 use_zweight_path = .false.
380 krng(2) = kpoints_number(kpoints) - kpoints%nik_skip
381
382 end if
383 end if
384
385 call write_kpoints_info(kpoints, krng(1), krng(2))
387
388 nkpt = krng(2) - krng(1) + 1
389
390
391
392 if (use_zweight_path) then
393 llp(1:dim) = llpp(1:dim)
394 llp(kpth_dir) = llpp(kpth_dir) * nkpt
395 else
396 llp(1:dim) = llpp(1:dim)
397 end if
398
399 write(message(1),'(a,i4,i4,i4)') 'Debug : llp = ', llp(1:3)
400 write(message(2),'(a,i4,i4,i4)') 'Debug : llpp = ', llpp(1:3)
401 call messages_info(2, debug_only=.true.)
402
403 safe_allocate(pmesh(1:llp(1), 1:llp(2), 1:llp(3), 1:3 + 1))
404 safe_allocate( pesp(1:llp(1), 1:llp(2), 1:llp(3), 1:st%d%nspin))
405
406 select case (pes_method)
407 case (option__photoelectronspectrum__pes_mask)
408 safe_allocate(lp(1:llpp(1), 1:llpp(2), 1:llpp(3), krng(1):krng(2), 1:3))
409 call pes_mask_pmesh(global_namespace, dim, kpoints, llpp, lg, pmesh, idxzero, krng, lp)
410
411 case (option__photoelectronspectrum__pes_flux)
412 ! Lp is allocated inside pes_flux_pmesh to comply with the
413 ! declinations of the different surfaces
414 safe_allocate(ekin(1:llp(1), 1:llp(2), 1:llp(3)))
415 ekin = m_zero
416 nspin = 1
417 if (st%d%ispin == spin_polarized) nspin = 2
418 call pes_flux_pmesh(pflux, global_namespace, dim, kpoints, llpp, pmesh, idxzero, krng, nspin, lp, ekin)
419 end select
420
421
422 if (.not. need_pmesh) then
423 ! There is no need to use pmesh we just need to sort Lg in place
424 ! in order to have a coordinate ordering coherent with pesP
425 ! NOTE: this works only for the mask_method since Lg is well-defined
426 ! only in this case
427 do idim = 1, dim
428 call sort(lg(1:llp(idim), idim))
429 end do
430 end if
431
432
433 call messages_write('Read PES restart files.')
434 call messages_info()
435
436
437
439
440 write(message(1),'(a,f10.2,a2,f10.2,a2,f10.2,a1)') &
441 "Zenith axis: (",pol(1),", ",pol(2),", ",pol(3),")"
442 call messages_info(1)
443
444
445 ! Convert the grid units
446 if (need_pmesh) then
447 do ii = 1,3
448 do i3 = 1, llp(3)
449 do i2 = 1, llp(2)
450 do i1 = 1, llp(1)
451 pmesh(i1, i2, i3, ii) = units_from_atomic(sqrt(units_out%energy), pmesh(i1, i2, i3, ii))
452 end do
453 end do
454 end do
455 end do
456 end if
457
458 if (resolve_states) then
459 do ist = st_range(1), st_range(2)
460
461 select case (pes_method)
462 case (option__photoelectronspectrum__pes_mask)
463 call pes_mask_map_from_states(restart, st, llpp, pesp, krng, lp, ist)
464 case (option__photoelectronspectrum__pes_flux)
465 call pes_flux_map_from_states(pflux, restart, st, llpp, pesp, krng, lp, ist)
466 end select
467
468 call output_spin_pes()
469 end do
470
471 else
472 ! Read the data
473 ist = 0
474
475 select case (pes_method)
476 case (option__photoelectronspectrum__pes_mask)
477 call pes_mask_map_from_states(restart, st, llpp, pesp, krng, lp)
478 case (option__photoelectronspectrum__pes_flux)
479 call pes_flux_map_from_states(pflux, restart, st, llpp, pesp, krng, lp)
480 end select
481
482 call output_spin_pes()
483
484 end if
485
486
487
488 write(message(1), '(a)') 'Done'
489 call messages_info(1)
490
492
493 call restart_end(restart)
494
495 call states_elec_end(st)
496
497 safe_deallocate_p(ions)
498 call kpoints_end(kpoints)
499
501 call io_end()
502 call messages_end()
503
504 call parser_end()
505 call global_end()
506
507 safe_deallocate_a(pesp)
508 safe_deallocate_a(pmesh)
509 safe_deallocate_a(lp)
510 safe_deallocate_a(ekin)
511 if (.not. need_pmesh .or. pes_method == option__photoelectronspectrum__pes_mask) then
512 safe_deallocate_a(lg)
513 end if
514contains
515
516 function outfile(name, ist, ispin, extension) result(fname)
517 character(len=*), intent(in) :: name
518 integer, intent(in) :: ist
519 integer, intent(in) :: ispin
520 character(len=*), optional, intent(in) :: extension
521 character(len=512) :: fname
522
523 if (ist == 0) then
524 fname = trim(name)
525 else
526 write(fname, '(a,a,i4.4)') trim(name), '-st', ist
527 end if
528
529 if (ispin >0) then
530 write(fname, '(a,a,i1)') trim(fname),'-sp', ispin
531 end if
532
533 if (present(extension) .and. len(extension)>0) then
534 fname = trim(fname)//'.'//trim(extension)
535 end if
536
537 end function outfile
538
539 subroutine output_spin_pes()
540
541 ! Write total quantities (summed over spin)
542 ispin = 0
543 if (st%d%ispin == unpolarized) then
544 pesp_out => pesp(:,:,:,1)
545 call output_pes()
546 else
547 ! Write total quantities (summed over spin)
548 safe_allocate(pesp_out(1:llp(1), 1:llp(2), 1:llp(3)))
549 pesp_out(:,:,:) = pesp(:,:,:,1) + pesp(:,:,:,2)
550
551 call output_pes()
552
553 safe_deallocate_p(pesp_out)
554
555 ! spin-resolved
556 do ispin = 1, st%d%nspin
557 pesp_out => pesp(:,:,:,ispin)
558 call output_pes()
559 end do
560 end if
561 end subroutine output_spin_pes
562
563 subroutine output_pes()
564
565 ! choose what to calculate
566 ! these functions are defined in pes_mask_out_inc.F90
567
568 if (st%d%ispin /= unpolarized .or. ist>0) then
570 end if
571
572 if (ist > 0) then
573 write(message(1), '(a,i4)') 'State = ', ist
574 call messages_info(1)
575 end if
576
577 if (st%d%ispin /= unpolarized) then
578 if (ispin > 0) then
579 write(message(1), '(a,i1)') 'Spin component= ', ispin
580 call messages_info(1)
581 else
582 write(message(1), '(a)') 'Spinless'
583 call messages_info(1)
584 end if
585 end if
586
587 if (ions%latt%nonorthogonal) then
588 coord_system => affine_coordinates_t(global_namespace, space%dim, ions%latt%rlattice_primitive)
589 else
590 coord_system => cartesian_t(global_namespace, space%dim)
591 end if
592
593
594 if (pesout%what(option__photoelectronspectrumoutput__energy_tot)) then
595 call messages_print_with_emphasis(msg="Energy-resolved PES", namespace=global_namespace)
596
597 select case (pes_method)
598 case (option__photoelectronspectrum__pes_mask)
599 call pes_mask_output_power_totalm(pesp_out,outfile('./PES_energy',ist, ispin, 'sum'), &
600 global_namespace, lg, llp, dim, emax, estep, interpolate = .true.)
601 case (option__photoelectronspectrum__pes_flux)
602 call pes_flux_out_energy(pflux, pesp_out, outfile('./PES_energy',ist, ispin, 'sum'), global_namespace, llp, ekin, dim)
603 end select
604
605 end if
606
607 if (pesout%what(option__photoelectronspectrumoutput__energy_angle)) then
608 call messages_print_with_emphasis(msg="Angle- and energy-resolved PES", namespace=global_namespace)
610 select case (pes_method)
611 case (option__photoelectronspectrum__pes_mask)
612 call pes_mask_output_ar_polar_m(pesp_out,outfile('./PES_angle_energy',ist, ispin, 'map'), &
613 global_namespace, lg, llp, dim, pol, emax, estep)
614 case (option__photoelectronspectrum__pes_flux)
615 call messages_not_implemented("Angle- and energy-resolved PES for the flux method")
616 end select
617
618
619 end if
620
621 if (pesout%what(option__photoelectronspectrumoutput__velocity_map_cut)) then
622 call messages_print_with_emphasis(msg="Velocity map on a plane", namespace=global_namespace)
623 dir = -1
624 if (sum((pvec-(/1 ,0 ,0/))**2) <= m_epsilon) dir = 1
625 if (sum((pvec-(/0 ,1 ,0/))**2) <= m_epsilon) dir = 2
626 if (sum((pvec-(/0 ,0 ,1/))**2) <= m_epsilon) dir = 3
627
628 if (use_zweight_path) then
629 filename = outfile('PES_velocity_map', ist, ispin, 'path')
630 else
631 filename = outfile('PES_velocity_map', ist, ispin, 'p'//index2axis(dir)//'=0')
632 end if
633
634 if (dir == -1) then
635 write(message(1), '(a)') 'Unrecognized plane. Use -u to change.'
636 call messages_fatal(1)
637 else
638 write(message(1), '(a)') 'Save velocity map on plane: '//index2axis(dir)//" = 0"
639 call messages_info(1)
640 end if
641
642 if (integrate /= integrate_none) then
643 write(message(1), '(a)') 'Integrate on: '//index2var(integrate)
644 call messages_info(1)
645 filename = trim(filename)//'.i_'//trim(index2var(integrate))
646 end if
647
648 if (need_pmesh) then
649 call pes_out_velocity_map_cut(global_namespace, pesp_out, filename, llp, dim, pol, dir, integrate, &
650 pos = idxzero, pmesh = pmesh)
651 else
652 call pes_out_velocity_map_cut(global_namespace, pesp_out, filename, llp, dim, pol, dir, integrate, &
653 pos = idxzero, lk = lg)
654 end if
655 end if
657 if (pesout%what(option__photoelectronspectrumoutput__energy_xy)) then
658 call messages_print_with_emphasis(msg="Angle and energy-resolved on a plane", namespace=global_namespace)
659 if (uestep > 0 .and. uestep > estep) then
660 estep = uestep
661 else
662 estep = emax/size(lg,1)
663 end if
664
665 select case (pes_method)
666 case (option__photoelectronspectrum__pes_mask)
667 call pes_mask_output_ar_plane_m(pesp_out,outfile('./PES_energy', ist, ispin, 'map'), &
668 global_namespace, lg, llp, dim, pol, emax, estep)
669 case (option__photoelectronspectrum__pes_flux)
670 call messages_not_implemented("Angle and energy-resolved on a plane for the flux method")
671 end select
672
673 end if
674
675 if (pesout%what(option__photoelectronspectrumoutput__energy_th_ph)) then
676 call messages_print_with_emphasis(msg="PES on spherical cuts", namespace=global_namespace)
677
678 write(message(1), '(a,es19.12,a2,es19.12,2x,a19)') &
679 'Save PES on a spherical cut at E= ',emin,", ",emax, &
680 str_center('['//trim(units_abbrev(units_out%energy)) // ']', 19)
681 call messages_info(1)
682
683 if (uestep > 0 .and. uestep > estep) then
684 estep = uestep
685 else
686 estep = emax/size(lg,1)
687 end if
688
689 select case (pes_method)
690 case (option__photoelectronspectrum__pes_mask)
691 call pes_mask_output_ar_spherical_cut_m(pesp_out,outfile('./PES_sphere', ist, ispin, 'map'), &
692 global_namespace, lg, llp, dim, pol, emin, emax, estep)
693
694 case (option__photoelectronspectrum__pes_flux)
695 call messages_not_implemented("PES on spherical cuts for the flux method")
696 end select
697
698 end if
699
700 if (pesout%what(option__photoelectronspectrumoutput__velocity_map)) then
701
702 call messages_print_with_emphasis(msg="Full velocity map", namespace=global_namespace)
703
704 if (.not. (bitand(pesout%how(option__photoelectronspectrumoutput__velocity_map), option__outputformat__netcdf) /= 0) .and. &
705 .not. (bitand(pesout%how(option__photoelectronspectrumoutput__velocity_map), option__outputformat__vtk) /= 0) .and. &
706 .not. (bitand(pesout%how(option__photoelectronspectrumoutput__velocity_map), option__outputformat__ascii) /= 0)) then
707 message(1) = 'User must specify the format with "OutputFormat".'
708 message(2) = 'Available options are: necdf, vtk, ascii.'
709 call messages_fatal(2)
710
711 end if
712
713
714 filename = outfile('./PES_velocity_map', ist, ispin)
715 if (need_pmesh) then
716 !force vtk output
717! how = io_function_fill_how("VTK")
718 if (bitand(pesout%how(option__photoelectronspectrumoutput__velocity_map), option__outputformat__ascii) /= 0) then
719 call pes_flux_out_vmap(pflux, pesp_out, filename, global_namespace, llp, pmesh, space%dim)
720 else
721 call pes_out_velocity_map(pesp_out, filename, global_namespace, space, lg, llp, &
722 pesout%how(option__photoelectronspectrumoutput__velocity_map), &
723 (/m_one, m_one, m_one/), coord_system, pmesh)
724 end if
725 else
726 call pes_out_velocity_map(pesp_out, filename, global_namespace, space, lg, llp, &
727 pesout%how(option__photoelectronspectrumoutput__velocity_map), &
728 (/m_one, m_one, m_one/), coord_system)
729 end if
730
731 end if
732
733 if (pesout%what(option__photoelectronspectrumoutput__arpes)) then
734 call messages_print_with_emphasis(msg="ARPES", namespace=global_namespace)
735
736 do i3 = 1, llp(3)
737 do i2 = 1, llp(2)
738 do i1 = 1, llp(1)
739 pmesh(i1, i2, i3, dim) = units_from_atomic(units_out%energy, &
740 sign(m_one,pmesh( i1, i2, i3, dim)) * sum(pmesh(i1, i2, i3, 1:dim)**2) / m_two)
741 end do
742 end do
743 end do
744
745 pesout%how(option__photoelectronspectrumoutput__arpes) = io_function_fill_how("VTK")
746
747 call pes_out_velocity_map(pesp_out, outfile('./PES_ARPES', ist, ispin), &
748 global_namespace, space, lg, llp, pesout%how(option__photoelectronspectrumoutput__arpes), &
749 (/m_one, m_one, m_one/), coord_system, pmesh)
750 end if
751
752
753 if (pesout%what(option__photoelectronspectrumoutput__arpes_cut)) then
754 call messages_print_with_emphasis(msg="ARPES cut on reciprocal space path", namespace=global_namespace)
755
756 filename = outfile('./PES_ARPES', ist, ispin, "path")
757 call pes_out_arpes_cut(global_namespace, pesp_out, filename, dim, llp, pmesh, ekin)
758
759 end if
760
761 safe_deallocate_p(coord_system)
762
763 end subroutine output_pes
764
765
766 subroutine write_kpoints_info(kpoints, ikstart, ikend)
767 type(kpoints_t), intent(in) :: kpoints
768 integer, intent(in) :: ikstart
769 integer, intent(in) :: ikend
770
771 integer :: ik, idir
772 character(len=100) :: str_tmp
773
774 push_sub(write_kpoints_info)
775
776 do ik = ikstart, ikend
777 write(message(1),'(i8,1x)') ik
778 write(str_tmp,'(f12.6)') kpoints%get_weight(ik)
779 message(1) = trim(message(1)) // trim(str_tmp)//' |'
780 do idir = 1, kpoints%full%dim
781 write(str_tmp,'(f12.6)') kpoints%reduced%red_point(idir, ik)
782 message(1) = trim(message(1)) // trim(str_tmp)
783 end do
784 message(1) = trim(message(1)) //' |'
785 call messages_info(1)
786 end do
787
788 call messages_info(1)
789
790 pop_sub(write_kpoints_info)
791 end subroutine write_kpoints_info
792
793 subroutine get_laser_polarization(lPol)
794 real(real64), intent(out) :: lPol(:)
795
796 type(block_t) :: blk
797 integer :: no_l
798 complex(real64) :: cPol(1:3)
799
800 push_sub(get_laser_polarization)
801
802 cpol = m_zero
803
804 no_l = 0
805 if (parse_block(global_namespace, 'TDExternalFields', blk) == 0) then
806 no_l = parse_block_n(blk)
807
808 call parse_block_cmplx(blk, 0, 1, cpol(1))
809 call parse_block_cmplx(blk, 0, 2, cpol(2))
810 call parse_block_cmplx(blk, 0, 3, cpol(3))
811
812
813 call parse_block_end(blk)
814 end if
815
816 lpol(:) = abs(cpol)
817
818 if (no_l > 1) then
819 message(1) = "There is more than one external field. Polarization will be selected"
820 message(2) = "from the first field. Use -V to change axis."
821 call messages_info(2)
822 end if
823
825 end subroutine get_laser_polarization
826
827 character(5) pure function index2var(ivar) result(ch)
828 integer, intent(in) :: ivar
829
830 select case (ivar)
831 case (integrate_phi)
832 ch = 'phi'
833 case (integrate_theta)
834 ch = 'theta'
835 case (integrate_r)
836 ch = 'r'
837 case (integrate_kx)
838 ch = 'kx'
839 case (integrate_ky)
840 ch = 'ky'
841 case (integrate_kz)
842 ch = 'kz'
843 case default
844 write(ch,'(i1)') ivar
845 end select
846 end function index2var
847
848
849
850
851
852
853end program photoelectron_spectrum
854
855!! Local Variables:
856!! mode: f90
857!! coding: utf-8
858!! End:
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:149
This module handles the calculation mode.
integer, parameter, public p_strategy_domains
parallelization in domains
subroutine, public getopt_init(ierr)
Initializes the getopt machinery. Must be called before attempting to parse the options....
subroutine, public getopt_end
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
Definition: global.F90:190
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:382
real(real64), parameter, public m_zero
Definition: global.F90:188
type(mpi_comm), parameter, public serial_dummy_comm
Alias MPI_COMM_UNDEFINED for the specific use case of initialising Octopus utilities with no MPI supp...
Definition: global.F90:265
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
Definition: global.F90:354
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_one
Definition: global.F90:189
subroutine, public io_function_read_what_how_when(namespace, space, what, how, output_interval, what_tag_in, how_tag_in, output_interval_tag_in, ignore_error)
integer(int64) function, public io_function_fill_how(where)
Use this function to quickly plot functions for debugging purposes: call dio_function_output(io_funct...
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_end()
Definition: io.F90:258
subroutine, public kpoints_end(this)
Definition: kpoints.F90:1012
subroutine, public kpoints_init(this, namespace, symm, dim, periodic_dim, latt)
Definition: kpoints.F90:322
integer pure function, public kpoints_number(this)
Definition: kpoints.F90:1099
subroutine, public messages_end()
Definition: messages.F90:277
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:903
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1096
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_init(output_dir)
Definition: messages.F90:224
subroutine, public messages_new_line()
Definition: messages.F90:1117
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_input_error(namespace, var, details, row, column)
Definition: messages.F90:696
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:599
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:270
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:264
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 pes_flux_map_from_states(this, restart, st, ll, pesP, krng, Lp, istin)
Definition: pes_flux.F90:2992
subroutine, public pes_flux_out_energy(this, pesK, file, namespace, ll, Ekin, dim)
Definition: pes_flux.F90:3649
subroutine, public pes_flux_pmesh(this, namespace, dim, kpoints, ll, pmesh, idxZero, krng, nspin, Lp, Ekin)
Definition: pes_flux.F90:2951
subroutine, public pes_flux_reciprocal_mesh_gen(this, namespace, space, st, kpoints, comm, post)
Definition: pes_flux.F90:819
subroutine, public pes_flux_out_vmap(this, pesK, file, namespace, ll, pmesh, dim)
Definition: pes_flux.F90:3739
integer, parameter, public integrate_r
Definition: pes_mask.F90:262
integer, parameter, public integrate_kx
Definition: pes_mask.F90:262
subroutine, public pes_mask_pmesh(namespace, dim, kpoints, ll, LG, pmesh, idxZero, krng, Lp)
Definition: pes_mask.F90:1579
subroutine, public pes_mask_output_power_totalm(pesK, file, namespace, Lk, ll, dim, Emax, Estep, interpolate)
Definition: pes_mask.F90:3184
integer, parameter, public integrate_phi
Definition: pes_mask.F90:262
subroutine, public pes_mask_output_ar_spherical_cut_m(pesK, file, namespace, Lk, ll, dim, dir, Emin, Emax, Estep)
Definition: pes_mask.F90:2913
integer, parameter, public integrate_none
Definition: pes_mask.F90:262
integer, parameter, public integrate_theta
Definition: pes_mask.F90:262
subroutine, public pes_mask_output_ar_plane_m(pesK, file, namespace, Lk, ll, dim, dir, Emax, Estep)
Definition: pes_mask.F90:2783
subroutine, public pes_mask_map_from_states(restart, st, ll, pesK, krng, Lp, istin)
Definition: pes_mask.F90:1819
integer, parameter, public integrate_ky
Definition: pes_mask.F90:262
integer, parameter, public integrate_kz
Definition: pes_mask.F90:262
subroutine, public pes_mask_read_info(dir, namespace, dim, Emax, Estep, ll, Lk, RR)
Read pes info.
Definition: pes_mask.F90:3508
subroutine, public pes_mask_output_ar_polar_m(pesK, file, namespace, Lk, ll, dim, dir, Emax, Estep)
Definition: pes_mask.F90:2664
subroutine, public pes_out_velocity_map_cut(namespace, pesK, file, ll, dim, pol, dir, integrate, pos, Lk, pmesh)
Definition: pes_out.F90:292
subroutine, public pes_out_velocity_map(pesK, file, namespace, space, Lk, ll, how, spacing, coord_system, pmesh)
Definition: pes_out.F90:162
subroutine, public pes_out_arpes_cut(namespace, arpes, file, dim, ll, pmesh, Ekin)
Definition: pes_out.F90:216
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:309
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:517
integer, parameter, public restart_td
Definition: restart.F90:200
integer, parameter, public restart_type_load
Definition: restart.F90:246
subroutine, public restart_end(restart)
Definition: restart.F90:723
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_init(st, namespace, space, valence_charge, kpoints)
Initialize a new states_elec_t object.
character(len=80) function, public str_center(s_in, l_in)
puts space around string, so that it is centered
Definition: string.F90:197
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
subroutine, public unit_system_init(namespace)
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 write_kpoints_info(kpoints, ikstart, ikend)
subroutine get_laser_polarization(lPol)
subroutine output_spin_pes()
character(5) pure function index2var(ivar)
program photoelectron_spectrum
subroutine output_pes()
character(len=512) function outfile(name, ist, ispin, extension)
Extension of space that contains the knowledge of the spin dimension.
int true(void)