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