Octopus
pes_mask.F90
Go to the documentation of this file.
1!! Copyright (C) 2006-2011 U. De Giovannini, M. Marques
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module pes_mask_oct_m
23 use box_oct_m
26 use comm_oct_m
28 use cube_oct_m
29 use debug_oct_m
33 use fft_oct_m
35 use global_oct_m
36 use grid_oct_m
40 use io_oct_m
41 use ions_oct_m
43 use lasers_oct_m
44 use loct_oct_m
45 use math_oct_m
47 use mesh_oct_m
49 use mpi_oct_m
51#if defined(HAVE_NETCDF)
52 use netcdf
53#endif
55 use parser_oct_m
57 use qshep_oct_m
59 use space_oct_m
60 use sort_oct_m
62 use string_oct_m
63 use unit_oct_m
66 use vtk_oct_m
68
69 implicit none
70
71 private
72
73 public :: &
74 pes_mask_t, &
95
96
97 type pes_mask_t
98 private
99 complex(real64), allocatable, public :: k(:,:,:,:,:,:)
100 ! !< mask%k(ll(1),ll(2),ll(3),st%d%dim, st%nst, st%nik)
101
102 ! mesh- and cube-related stuff
103 integer :: np
105 integer :: ll(3)
106 integer :: fs_n_global(1:3)
107 integer :: fs_n(1:3)
108 integer :: fs_istart(1:3)
109
110 real(real64) :: spacing(3)
111
112 type(mesh_t), pointer, public :: mesh => null()
113 type(cube_t) :: cube
114
115 real(real64), allocatable, public :: vec_pot(:,:)
117 real(real64), allocatable, public :: Mk(:,:,:)
118 type(cube_function_t) :: cM
119 real(real64), allocatable, public :: mask_R(:)
120 integer :: shape
121 real(real64), allocatable :: ufn(:)
122 logical :: user_def
123
124 real(real64), allocatable, public :: Lk(:,:)
125
126 real(real64) :: enlarge(3)
127 real(real64) :: enlarge_2p(3)
128
129 real(real64) :: start_time
130 real(real64) :: energyMax
131 real(real64) :: energyStep
132
133 integer :: sw_evolve
134 logical :: back_action
135 logical :: add_psia
136 logical :: filter_k
137
138 integer :: mode
139 integer :: pw_map_how
140
141 type(fft_t) :: fft
142
143 type(mesh_cube_parallel_map_t) :: mesh_cube_map
144
145
146 end type pes_mask_t
147
148
149 integer, public, parameter :: &
150 PES_MASK_MODE_MASK = 1, &
153
154 integer, parameter :: &
155 PW_MAP_FFT = 3, & !< FFT - normally from fftw3
156 pw_map_nfft = 5, &
157 pw_map_pfft = 6, &
158 pw_map_pnfft = 7
159
160 integer, parameter :: &
161 M_SIN2 = 1, &
162 m_step = 2, &
163 m_erf = 3
164
165 integer, parameter :: &
166 IN = 1, &
167 out = 2
168
169 integer, public, parameter :: &
170 INTEGRATE_NONE = -1, &
171 integrate_phi = 1, &
172 integrate_theta = 2, &
173 integrate_r = 3, &
174 integrate_kx = 4, &
175 integrate_ky = 5, &
176 integrate_kz = 6
177
178
179contains
180
181
182 ! ---------------------------------------------------------
183 subroutine pes_mask_init(mask, namespace, space, mesh, box, st, ext_partners, abs_boundaries, max_iter,dt)
184 type(pes_mask_t), intent(out) :: mask
185 type(namespace_t), intent(in) :: namespace
186 class(space_t), intent(in) :: space
187 type(mesh_t), target, intent(in) :: mesh
188 class(box_t), intent(in) :: box
189 type(states_elec_t), intent(in) :: st
190 type(partner_list_t), intent(in) :: ext_partners
191 type(absorbing_boundaries_t), intent(in) :: abs_boundaries
192 integer, intent(in) :: max_iter
193 real(real64), intent(in) :: dt
195 type(block_t) :: blk
196
197 integer :: il, it, ll(3)
198 real(real64) :: field(3)
199 real(real64) :: deltae, maxe, pcutoff, tmp
200 integer :: defaultmask,k1,k2,st1,st2
201 integer :: cols_pesmask_block, idim, ip
203 real(real64) :: xx(space%dim), r
204 real(real64) :: ufn_re, ufn_im
205 character(len=1024) :: user_def_expr
206 type(lasers_t), pointer :: lasers
209
210 mask%mesh => mesh
211
212 if (space%is_periodic()) then
213 call messages_experimental("PES_mask with periodic dimensions", namespace=namespace)
214 end if
217 write(message(1),'(a,i1,a)') 'Info: Calculating PES using mask technique.'
218 call messages_info(1, namespace=namespace)
220
221 select type (box => box)
222 type is (box_sphere_t)
223 class default
224 if (.not. space%is_periodic()) then
225 message(1) = 'PhotoElectronSpectrum = pes_mask usually requires BoxShape = sphere.'
226 message(2) = 'Unless you know what you are doing modify this parameter and rerun.'
227 call messages_warning(2, namespace=namespace)
228 end if
229 end select
231 if (abs_boundaries%abtype /= not_absorbing) then
232 message(1) = 'PhotoElectronSpectrum = pes_mask already contains absorbing boundaries.'
233 message(2) = 'Set AbsorbingBoundaries = no and rerun.'
234 call messages_fatal(2, namespace=namespace)
235 end if
237
238!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
239!! Calculation mode
240!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
241
242 !%Variable PESMaskMode
243 !%Type integer
244 !%Default mask_mode
245 !%Section Time-Dependent::PhotoElectronSpectrum
246 !%Description
247 !% PES calculation mode.
248 !%Option mask_mode 1
249 !% Mask method.
250 !%Option fullmask_mode 2
251 !% Full mask method. This includes a back action of the momentum-space states on the
252 !% interaction region. This enables electrons to come back from the continuum.
253 !%Option passive_mode 3
254 !% Passive analysis of the wf. Simply analyze the plane-wave components of the
255 !% wavefunctions on the region <i>r</i> > <i>R1</i>. This mode employs a step masking function by default.
256 !%End
257 call parse_variable(namespace, 'PESMaskMode', pes_mask_mode_mask, mask%mode)
258 if (.not. varinfo_valid_option('PESMaskMode', mask%mode)) call messages_input_error(namespace, 'PESMaskMode')
259 call messages_print_var_option("PESMaskMode", mask%mode, namespace=namespace)
261 select case (mask%mode)
263 defaultmask = m_step
264 mask%back_action = .false.
265
267 defaultmask = m_sin2
268 mask%back_action = .true.
269 mask%mode = pes_mask_mode_mask
270 case default
271 defaultmask = m_sin2
272 mask%back_action = .false.
273 end select
274
275
276 !%Variable PESMaskStartTime
277 !%Type float
278 !%Default -1.0
279 !%Section Time-Dependent::PhotoElectronSpectrum
280 !%Description
281 !% The time photoelectrons start to be recorded. In pump-probe simulations, this allows
282 !% getting rid of an unwanted ionization signal coming from the pump.
283 !% NOTE: This will enforce the mask boundary conditions for all times.
284 !%End
285 call parse_variable(namespace, 'PESMaskStartTime', -m_one, mask%start_time, unit = units_inp%time)
286
287 !%Variable PESMaskPlaneWaveProjection
288 !%Type integer
289 !%Default fft_map
290 !%Section Time-Dependent::PhotoElectronSpectrum
291 !%Description
292 !% With the mask method, wavefunctions in the continuum are treated as plane waves.
293 !% This variable sets how to calculate the plane-wave projection in the buffer
294 !% region. We perform discrete Fourier transforms (DFT) in order to approximate
295 !% a continuous Fourier transform. The major drawback of this approach is the built-in
296 !% periodic boundary condition of DFT. Choosing an appropriate plane-wave projection
297 !% for a given simulation in addition to <tt>PESMaskEnlargeFactor</tt> and
298 !% <tt>PESMask2PEnlargeFactor</tt>will help to converge the results.
299 !%
300 !% NOTE: depending on the value of <tt>PESMaskMode</tt> <tt>PESMaskPlaneWaveProjection</tt>,
301 !% may affect not only performance but also the time evolution of the density.
302 !%Option fft_out 2
303 !% FFT filtered in order to keep only outgoing waves. 1D only.
304 !%Option fft_map 3
305 !% FFT transform.
306 !%Option nfft_map 5
307 !% Non-equispaced FFT map.
308 !%Option pfft_map 6
309 !% Use PFFT library.
310 !%Option pnfft_map 7
311 !% Use PNFFT library.
312 !%End
313 call parse_variable(namespace, 'PESMaskPlaneWaveProjection', pw_map_fft, mask%pw_map_how)
314
315 if (.not. varinfo_valid_option('PESMaskPlaneWaveProjection', mask%pw_map_how)) then
316 call messages_input_error(namespace, 'PESMaskPlaneWaveProjection')
317 end if
318
319 call messages_print_var_option("PESMaskPlaneWaveProjection", mask%pw_map_how, namespace=namespace)
320
321 if (mask%pw_map_how == pw_map_pfft .and. (.not. mask%mesh%parallel_in_domains)) then
322 message(1)= "Trying to use PESMaskPlaneWaveProjection = pfft_map with no domain parallelization."
323 message(2)= "Projection method changed to more efficient fft_map."
324 call messages_warning(2, namespace=namespace)
325 mask%pw_map_how = pw_map_fft
326 end if
327
328 if (mask%pw_map_how == pw_map_pnfft .and. (.not. mask%mesh%parallel_in_domains)) then
329 message(1)= "Trying to use PESMaskPlaneWaveProjection = pnfft_map with no domain parallelization."
330 message(2)= "Projection method changed to more efficient nfft_map."
331 call messages_warning(2, namespace=namespace)
332 mask%pw_map_how = pw_map_nfft
333 end if
334
335#if !defined(HAVE_NFFT)
336 if (mask%pw_map_how == pw_map_nfft) then
337 message(1) = "PESMaskPlaneWaveProjection = nfft_map requires NFFT but that library was not linked."
338 call messages_fatal(1, namespace=namespace)
339 end if
340#endif
341
342#if !defined(HAVE_PFFT)
343 if (mask%pw_map_how == pw_map_pfft) then
344 message(1) = "PESMaskPlaneWaveProjection = pfft_map requires PFFT but that library was not linked."
345 call messages_fatal(1, namespace=namespace)
346 end if
347#endif
348
349#if !defined(HAVE_PNFFT)
350 if (mask%pw_map_how == pw_map_pnfft) then
351 message(1) = "PESMaskPlaneWaveProjection = pnfft_map requires PNFFT but that library was not linked."
352 call messages_fatal(1, namespace=namespace)
353 end if
354#endif
355
356
357 !%Variable PESMaskEnlargeFactor
358 !%Type float
359 !%Default 1
360 !%Section Time-Dependent::PhotoElectronSpectrum
361 !%Description
362 !% Mask box enlargement level. Enlarges the mask bounding box by a <tt>PESMaskEnlargeFactor</tt>.
363 !% This helps to avoid wavefunction wrapping at the boundaries.
364 !%End
365
366 mask%enlarge = m_one
367 call parse_variable(namespace, 'PESMaskEnlargeFactor', m_one, mask%enlarge(1))
368
369 if (.not. is_close(mask%enlarge(1), m_one)) then
370
371 mask%enlarge(space%periodic_dim + 1:space%dim) = mask%enlarge(1)
372 mask%enlarge(1:space%periodic_dim) = m_one
373
374 if (space%is_periodic()) then
375 call messages_print_var_value("PESMaskEnlargeFactor", mask%enlarge(1:space%dim), namespace=namespace)
376 else
377 call messages_print_var_value("PESMaskEnlargeFactor", mask%enlarge(1), namespace=namespace)
378 end if
379
380 end if
381 if (mask%enlarge(1) < m_one) then
382 message(1) = "PESMaskEnlargeFactor must be bigger than one."
383 call messages_fatal(1, namespace=namespace)
384 end if
385
386 call messages_obsolete_variable(namespace, 'PESMaskEnlargeLev', 'PESMaskEnlargeFactor')
387
388 !%Variable PESMask2PEnlargeFactor
389 !%Type float
390 !%Default 1.0
391 !%Section Time-Dependent::PhotoElectronSpectrum
392 !%Description
393 !% Mask two points enlargement factor. Enlarges the mask box by adding two
394 !% points at the edges of the box in each direction (x,y,z) at a distance
395 !% L=Lb*<tt>PESMask2PEnlargeFactor</tt> where <i>Lb</i> is the box size.
396 !% This allows to run simulations with an additional void space at a price of
397 !% adding few points. The Fourier space associated with the new box is restricted
398 !% by the same factor.
399 !%
400 !% Note: needs <tt> PESMaskPlaneWaveProjection = nfft_map or pnfft_map </tt>.
401 !%End
402
403 mask%enlarge_2p = m_one
404 call parse_variable(namespace, 'PESMask2PEnlargeFactor', m_one, mask%enlarge_2p(1))
405
406
407 if (.not. is_close(mask%enlarge_2p(1), m_one)) then
408
409 mask%enlarge_2p(space%periodic_dim + 1:space%dim) = mask%enlarge_2p(1)
410 mask%enlarge_2p(1:space%periodic_dim) = m_one
411
412 if (space%is_periodic()) then
413 call messages_print_var_value("PESMask2PEnlargeFactor", mask%enlarge_2p(1:space%dim), namespace=namespace)
414 else
415 call messages_print_var_value("PESMask2PEnlargeFactor", mask%enlarge_2p(1), namespace=namespace)
416 end if
417
418 if (mask%pw_map_how /= pw_map_nfft .and. mask%pw_map_how /= pw_map_pnfft) then
419 message(1) = "PESMask2PEnlargeFactor requires PESMaskPlaneWaveProjection = nfft_map"
420 message(2) = "or pnfft_map in order to run properly."
421 call messages_warning(2, namespace=namespace)
422 end if
423 end if
424
425 if (mask%enlarge_2p(1) < m_one) then
426 message(1) = "PESMask2PEnlargeFactor must be bigger than one."
427 call messages_fatal(1, namespace=namespace)
428 end if
429
430 call messages_obsolete_variable(namespace, 'PESMaskNFFTEnlargeLev', 'PESMask2PEnlargeFactor')
431
432
433 mask%ll = 1
434 mask%spacing = -m_one
435
436 mask%spacing(1:space%dim) = mesh%spacing(1:space%dim)
437 mask%ll(1:space%dim) = mesh%idx%ll(1:space%dim)
438
439 !Enlarge the cube region
440 mask%ll(1:space%dim) = int(mask%ll(1:space%dim) * mask%enlarge(1:space%dim))
441
442 mask%np = mesh%np_part ! the mask is local
443
444 select case (mask%pw_map_how)
445
446 case (pw_map_pfft)
447 assert(mask%mesh%parallel_in_domains)
448 call cube_init(mask%cube, mask%ll, namespace, space, &
449 mesh%spacing, mesh%coord_system, &
450 fft_type = fft_complex, fft_library = fftlib_pfft, nn_out = ll, &
451 mpi_grp = mask%mesh%mpi_grp, need_partition=.true.)
452 call cube_init_cube_map(mask%cube, mesh)
453
454
455! mask%ll(1) = mask%cube%fs_n(3)
456! mask%ll(2) = mask%cube%fs_n(1)
457! mask%ll(3) = mask%cube%fs_n(2)
458! Note: even if tempting, setting
459! mask%ll(1:3) = mask%cube%fs_n(1:3)
460! results in the wrong index mapping! (1->3, 2->1, 3->2)
461! Well.. very much against intuition it turns out to be
462 mask%ll(1:3) = mask%cube%rs_n(1:3)
463
464 mask%fft = mask%cube%fft
465 if (mask%mesh%parallel_in_domains .and. mask%cube%parallel_in_domains) then
466 call mesh_cube_parallel_map_init(mask%mesh_cube_map, mask%mesh, mask%cube)
467 end if
468
469 case (pw_map_fft)
470 call cube_init(mask%cube, mask%ll, namespace, space, &
471 mesh%spacing, mesh%coord_system, &
472 fft_type = fft_complex, fft_library = fftlib_fftw, &
473 nn_out = ll)
474 call cube_init_cube_map(mask%cube, mesh)
475 mask%ll = ll
476 mask%fft = mask%cube%fft
477
478
479 case (pw_map_nfft)
480
481 !NFFT initialization
482 ! we just add 2 points for the enlarged region
483 if (.not. is_close(mask%enlarge_2p(1), m_one)) mask%ll(1:space%dim) = mask%ll(1:space%dim) + 2
484
485 call cube_init(mask%cube, mask%ll, namespace, space, &
486 mesh%spacing, mesh%coord_system, &
487 fft_type = fft_complex, fft_library = fftlib_nfft, &
488 nn_out = ll, tp_enlarge = mask%enlarge_2p)
489 call cube_init_cube_map(mask%cube, mesh)
490
491 mask%ll = ll
492 mask%fft = mask%cube%fft
493
494
495 case (pw_map_pnfft)
496
497 if (.not. is_close(mask%enlarge_2p(1), m_one)) mask%ll(1:space%dim) = mask%ll(1:space%dim) + 2
498
499 call cube_init(mask%cube, mask%ll, namespace, space, &
500 mask%mesh%spacing, mask%mesh%coord_system, &
501 fft_type = fft_complex, fft_library = fftlib_pnfft, &
502 nn_out = ll, tp_enlarge = mask%enlarge_2p, &
503 mpi_grp = mask%mesh%mpi_grp, need_partition=.true.)
504 call cube_init_cube_map(mask%cube, mesh)
505
506 mask%ll(1:3) = mask%cube%fs_n(1:3)
507
508 mask%fft = mask%cube%fft
509 if (mask%mesh%parallel_in_domains .and. mask%cube%parallel_in_domains) then
510 call mesh_cube_parallel_map_init(mask%mesh_cube_map, mask%mesh, mask%cube)
511 end if
512
513 case default
514 !Program should die before coming here
515 write(message(1),'(a)') "PESMaskPlaneWaveProjection unrecognized option."
516 call messages_fatal(1, namespace=namespace)
517
518 end select
519
520 !Indices
521
522 if (mask%pw_map_how == pw_map_pfft) then
523 mask%fs_istart = mask%cube%rs_istart
524 mask%fs_n = mask%cube%rs_n
525 mask%fs_n_global = mask%cube%rs_n_global
526 else
527 mask%fs_istart = mask%cube%fs_istart
528 mask%fs_n = mask%cube%fs_n
529 mask%fs_n_global = mask%cube%fs_n_global
530 end if
531
532 !Allocations
533 call zcube_function_alloc_rs(mask%cube, mask%cM, force_alloc = .true.)
534
535 safe_allocate(mask%Lk(1:maxval(mask%fs_n_global(:)),1:3))
536 mask%Lk(:,:) = m_zero
537
538 st1 = st%st_start
539 st2 = st%st_end
540 k1 = st%d%kpt%start
541 k2 = st%d%kpt%end
542 safe_allocate(mask%k(1:mask%ll(1),1:mask%ll(2),1:mask%ll(3),1:st%d%dim,st1:st2,k1:k2))
543 mask%k = m_z0
544
545
546 ! generate the map between mesh and cube
547 call pes_mask_generate_lk(mask) ! generate the physical momentum vector
548
549
550
551!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
552!! Mask Function options
553!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
554
555 safe_allocate(mask%mask_R(1:2))
556
557 !%Variable PESMaskShape
558 !%Type integer
559 !%Default m_sin2
560 !%Section Time-Dependent::PhotoElectronSpectrum
561 !%Description
562 !% The mask function shape.
563 !%Option m_sin2 1
564 !% sin2 mask.
565 !%Option m_step 2
566 !%step function.
567 !%Option m_erf 3
568 !%Error function. Not Implemented.
569 !%End
570 call parse_variable(namespace, 'PESMaskShape', defaultmask, mask%shape)
571 if (.not. varinfo_valid_option('PESMaskShape', mask%shape)) call messages_input_error(namespace, 'PESMaskShape')
572 call messages_print_var_option("PESMaskShape", mask%shape, namespace=namespace)
573
574 !%Variable PESMaskSize
575 !%Type block
576 !%Section Time-Dependent::PhotoElectronSpectrum
577 !%Description
578 !% Set the size of the mask function.
579 !% Here you can set the inner (R1) and outer (R2) radius by setting
580 !% the block as follows:
581 !%
582 !% <tt>%PESMaskSize
583 !% <br>&nbsp;&nbsp; R1 | R2 | "user-defined"
584 !% <br>%</tt>
585 !%
586 !% The optional 3rd column is a user-defined expression for the mask
587 !% function. For example, <i>r</i> creates a spherical mask (which is the
588 !% default for <tt>BoxShape = sphere</tt>). Note, values R2 larger than
589 !% the box size may lead in this case to unexpected reflection
590 !% behaviours.
591 !%End
592 cols_pesmask_block = 0
593 if (parse_block(namespace, 'PESMaskSize', blk) == 0) then
594 cols_pesmask_block = parse_block_cols(blk, 0)
595 end if
596
597 mask%user_def = .false.
598
599 select case (cols_pesmask_block)
600 case (0)
601 select type (box => box)
602 type is (box_sphere_t)
603 mask%mask_R(1) = box%radius/m_two
604 mask%mask_R(2) = box%radius
605 message(1) = "Input: PESMaskSize R(1) and R(2) not specified. Using default values for spherical mask."
606 type is (box_parallelepiped_t)
607 mask%mask_R(1) = box%half_length(1)/m_two
608 mask%mask_R(2) = box%half_length(1)
609 message(1) = "Input: PESMaskSize R(1) and R(2) not specified. Using default values for cubic mask."
610 end select
611 call messages_info(1, namespace=namespace)
612 case (1)
613 call parse_block_float(blk, 0, 0, mask%mask_R(1), units_inp%length)
614 select type (box => box)
615 type is (box_sphere_t)
616 mask%mask_R(2) = box%radius
617 message(1) = "Input: PESMaskSize R(2) not specified. Using default value for spherical mask."
618 type is (box_parallelepiped_t)
619 mask%mask_R(2) = box%half_length(1)
620 message(1) = "Input: PESMaskSize R(2) not specified. Using default value for cubic mask."
621 end select
622 call messages_info(1, namespace=namespace)
623 case (2)
624 call parse_block_float(blk, 0, 0, mask%mask_R(1), units_inp%length)
625 call parse_block_float(blk, 0, 1, mask%mask_R(2), units_inp%length)
626
627 select type (box => box)
628 type is (box_sphere_t)
629 if (mask%mask_R(2) > box%radius) mask%mask_R(2) = box%radius
630 message(1) = "Info: using spherical mask."
631 type is (box_parallelepiped_t)
632 if (mask%mask_R(2) > box%half_length(1)) mask%mask_R(2) = box%half_length(1)
633 message(1) = "Info: using cubic mask."
634 end select
635 call messages_info(1, namespace=namespace)
636
637 case (3)
638 mask%user_def = .true.
639 safe_allocate(mask%ufn(1:mask%np))
640 mask%ufn = m_zero
641 call parse_block_float(blk, 0, 0, mask%mask_R(1), units_inp%length)
642 call parse_block_float(blk, 0, 1, mask%mask_R(2), units_inp%length)
643 call parse_block_string(blk, 0, 2, user_def_expr)
644 do ip = 1, mask%np
645 xx(1:space%dim) = mesh%x(1:space%dim, ip)
646 r = units_from_atomic(units_inp%length, norm2(xx(1:space%dim)))
647 do idim = 1, space%dim
648 xx(idim) = units_from_atomic(units_inp%length, xx(idim))
649 end do
650 call parse_expression(ufn_re, ufn_im, space%dim, xx, r, m_zero, user_def_expr)
651 mask%ufn(ip) = ufn_re
652 end do
653 message(1) = "Input: using user-defined mask function from expression:"
654 write(message(2),'(a,a)') ' R = ', trim(user_def_expr)
655 call messages_info(2, namespace=namespace)
656 end select
657
658 call parse_block_end(blk)
659
660 write(message(1),'(a,es10.3,3a)') &
661 " R1 = ", units_from_atomic(units_inp%length, mask%mask_R(1)),&
662 ' [', trim(units_abbrev(units_inp%length)), ']'
663 write(message(2),'(a,es10.3,3a)') &
664 " R2 = ", units_from_atomic(units_inp%length, mask%mask_R(2)),&
665 ' [', trim(units_abbrev(units_inp%length)), ']'
666 call messages_info(2, namespace=namespace)
667
668
669 call pes_mask_generate_mask(mask, namespace, mesh)
670
671 !%Variable PESMaskFilterCutOff
672 !%Type float
673 !%Default -1
674 !%Section Time-Dependent::PhotoElectronSpectrum
675 !%Description
676 !% In calculation with <tt>PESMaskMode = fullmask_mode</tt> and NFFT, spurious frequencies
677 !% may lead to numerical instability of the algorithm. This option gives the possibility
678 !% to filter out the unwanted components by setting an energy cut-off.
679 !% If <tt>PESMaskFilterCutOff = -1</tt> no filter is applied.
680 !%End
681 call parse_variable(namespace, 'PESMaskFilterCutOff', -m_one, pcutoff, unit = units_inp%energy)
682
683 mask%filter_k = .false.
684
685 if (pcutoff > m_zero) then
686 call messages_print_var_value("PESMaskFilterCutOff", pcutoff, unit = units_out%energy, namespace=namespace)
687 mask%filter_k = .true.
688
689 safe_allocate(mask%Mk(1:mask%ll(1),1:mask%ll(2),1:mask%ll(3)))
690
691 call pes_mask_generate_filter(mask,pcutoff)
692 end if
693
694!! Output
695
696 !%Variable PESMaskIncludePsiA
697 !%Type logical
698 !%Default false
699 !%Section Time-Dependent::PhotoElectronSpectrum
700 !%Description
701 !% Add the contribution of <math>\Psi_A</math> in the mask region to the photo-electron spectrum.
702 !% Literally adds the Fourier components of:
703 !% <math>\Theta(r-R1) \Psi_A(r)</math>
704 !% with <math>\Theta</math> being the Heaviside step function.
705 !% With this option PES will contain all the contributions starting from the inner
706 !% radius <math>R1</math>. Use this option to improve convergence with respect to the box size
707 !% and total simulation time.
708 !% Note: Carefully choose <math>R1</math> in order to avoid contributions from returning electrons.
709 !%End
710 call parse_variable(namespace, 'PESMaskIncludePsiA', .false., mask%add_psia)
711 if (mask%add_psia) then
712 message(1)= "Input: Include contribution from Psi_A."
713 call messages_info(1, namespace=namespace)
714 end if
715
716
717 !%Variable PESMaskSpectEnergyMax
718 !%Type float
719 !%Default maxval(mask%Lk)<math>^2/2</math>
720 !%Section Time-Dependent::PhotoElectronSpectrum
721 !%Description
722 !% The maximum energy for the PES spectrum.
723 !%End
724 maxe = m_epsilon
725 do idim = 1, space%dim
726 tmp = maxval(mask%Lk(1:mask%ll(idim),1:space%dim))**m_two/m_two
727 if (tmp > maxe) maxe = tmp
728 end do
729 call parse_variable(namespace, 'PESMaskSpectEnergyMax', maxe, mask%energyMax, unit = units_inp%energy)
730 call messages_print_var_value("PESMaskSpectEnergyMax", mask%energyMax, unit = units_out%energy, namespace=namespace)
731
732 !%Variable PESMaskSpectEnergyStep
733 !%Type float
734 !%Section Time-Dependent::PhotoElectronSpectrum
735 !%Description
736 !% The PES spectrum energy step.
737 !%End
738 deltae = minval(mask%Lk(2,1:space%dim) - mask%Lk(1,1:space%dim))**m_two/m_two
739 call parse_variable(namespace, 'PESMaskSpectEnergyStep', deltae, mask%energyStep, unit = units_inp%energy)
740 call messages_print_var_value("PESMaskSpectEnergyStep", mask%energyStep, unit = units_out%energy, namespace=namespace)
741
742
743!! Set external fields
744
745 safe_allocate(mask%vec_pot(0:max_iter,1:3))
746 mask%vec_pot=m_zero
747
748 lasers => list_get_lasers(ext_partners)
749 if(associated(lasers)) then
750 do il = 1, lasers%no_lasers
751 select case (laser_kind(lasers%lasers(il)))
753 do it = 1, max_iter
754 field=m_zero
755 call laser_field(lasers%lasers(il), field, it*dt)
756 ! We must sum with a -1 sign to account for the
757 ! electron charge.
758 mask%vec_pot(it,:)= mask%vec_pot(it,:) - field(:)
759 end do
760
761 case default
762 write(message(1),'(a)') 'PESMask should work only with TDExternalFields = vector_potential.'
763 write(message(2),'(a)') 'Unless PESMaskMode = passive_mode the results are likely to be wrong. '
764 call messages_warning(2, namespace=namespace)
765
766 end select
767 end do
768 end if
769
770
771
772 ! Compensate the sign for the forward <-> backward FT inversion
773 ! used with NFFT. See pes_mask_X_to_K comment for more info
774 if (mask%pw_map_how == pw_map_pnfft .or. &
775 mask%pw_map_how == pw_map_nfft) mask%Lk = - mask%Lk
776
777
778 pop_sub(pes_mask_init)
779 end subroutine pes_mask_init
780
781
782 ! ---------------------------------------------------------
783 subroutine pes_mask_end(mask)
784 type(pes_mask_t), intent(inout) :: mask
785
786 push_sub(pes_mask_end)
787
788 safe_deallocate_a(mask%k)
789
790 safe_deallocate_a(mask%vec_pot)
791 safe_deallocate_a(mask%mask_R)
792 safe_deallocate_a(mask%Lk)
793
794
795 if (mask%filter_k) then
796 safe_deallocate_a(mask%Mk)
797 end if
798
799 if (mask%mesh%parallel_in_domains .and. mask%cube%parallel_in_domains) then
800 call mesh_cube_parallel_map_end(mask%mesh_cube_map)
801 end if
802
803 call zcube_function_free_rs(mask%cube, mask%cM)
804 call cube_end(mask%cube)
805
806 if (mask%user_def) then
807 safe_deallocate_a(mask%ufn)
808 end if
809
810 pop_sub(pes_mask_end)
811 end subroutine pes_mask_end
812
813
814 ! --------------------------------------------------------
815 subroutine pes_mask_generate_lk(mask)
816 type(pes_mask_t), intent(inout) :: mask
817
818 integer :: ii, dim
819
820 push_sub(pes_mask_generate_lk)
821
822 dim = mask%mesh%box%dim
823
824 do ii = 1, maxval(mask%ll(:))
825 mask%Lk(ii,1:dim)= matmul(mask%cube%latt%klattice_primitive(1:dim,1:dim), mask%cube%Lfs(ii,1:dim))
826 end do
827
828 pop_sub(pes_mask_generate_lk)
829 end subroutine pes_mask_generate_lk
830
831 ! ======================================
833 ! ======================================
834 subroutine pes_mask_generate_filter(mask,cutOff)
835 type(pes_mask_t), intent(inout) :: mask
836 real(real64), intent(in) ::cutoff
837
838
839 integer :: kx,ky,kz, power
840 real(real64) :: kk(3), ee, emax
841
843
844 mask%Mk = m_zero
845
846 emax = maxval(mask%Lk(:,:))**2 / m_two
847
848 power = 8
849
850 do kx = 1, mask%ll(1)
851 kk(1) = mask%Lk(kx,1)
852 do ky = 1, mask%ll(2)
853 kk(2) = mask%Lk(ky,2)
854 do kz = 1, mask%ll(3)
855 kk(3) = mask%Lk(kz,3)
856
857 ee = sum(kk(1:mask%mesh%box%dim)**2) / m_two
858
859 if (ee > cutoff .and. ee < emax) then
860 mask%Mk(kx,ky,kz) = m_one * sin((emax-ee) * m_pi / (m_two * (cutoff)))**power
861 else
862 if (ee <= cutoff) mask%Mk(kx,ky,kz) = m_one
863 end if
864
865 end do
866 end do
867 end do
868
869
871 end subroutine pes_mask_generate_filter
872
873
874 ! --------------------------------------------------------
877 ! ---------------------------------------------------------
878 subroutine pes_mask_generate_mask(mask, namespace, mesh)
879 type(pes_mask_t), intent(inout) :: mask
880 type(namespace_t), intent(in) :: namespace
881 type(mesh_t), intent(in) :: mesh
882
883 push_sub(pes_mask_generate_mask)
884
885 call pes_mask_generate_mask_function(mask, namespace, mesh, mask%shape, mask%mask_R)
886
888
889 end subroutine pes_mask_generate_mask
890
891 ! --------------------------------------------------------
894 ! ---------------------------------------------------------
895 subroutine pes_mask_generate_mask_function(mask, namespace, mesh, shape, R, mask_sq)
896 type(pes_mask_t), intent(inout) :: mask
897 type(namespace_t), intent(in) :: namespace
898 type(mesh_t), intent(in) :: mesh
899 integer, intent(in) :: shape
900 real(real64), intent(in) :: r(2)
901 real(real64), optional, intent(out) :: mask_sq(:,:,:)
902
903 integer :: ip, dir
904 real(real64) :: width
905 real(real64) :: xx(1:mesh%box%dim), rr, dd, ddv(1:mesh%box%dim), tmp(1:mesh%box%dim)
906 complex(real64), allocatable :: mask_fn(:)
907
909
911 ! generate the mask function on the mesh
912 safe_allocate(mask_fn(1:mask%np))
913
914 mask_fn = m_zero
915 width = r(2) - r(1)
916 xx = m_zero
917
918 select case (shape)
919 case (m_sin2)
920 do ip = 1, mask%np
921 call mesh_r(mesh, ip, rr, coords=xx)
922
923 if (mask%user_def) then
924 dd = mask%ufn(ip) - r(1)
925 if (dd > m_zero) then
926 if (mask%ufn(ip) < r(2)) then
927 mask_fn(ip) = m_one * sin(dd * m_pi / (m_two * (width)))**2
928 else
929 mask_fn(ip) = m_one
930 end if
931 end if
932
933 else ! mask%user_def == .false.
934
935 select type (box => mesh%box)
936 type is (box_sphere_t)
937 dd = rr - r(1)
938 if (dd > m_zero) then
939 if (dd < width) then
940 mask_fn(ip) = m_one * sin(dd * m_pi / (m_two * (width)))**2
941 else
942 mask_fn(ip) = m_one
943 end if
944 end if
945
946 type is (box_parallelepiped_t)
947
948 ! We are filling from the center opposite to the spherical case
949 tmp = m_one
950 mask_fn(ip) = m_one
951 ddv = abs(xx) - r(1)
952 do dir=1, mesh%box%dim
953 if (ddv(dir) > m_zero) then
954 if (ddv(dir) < width) then
955 tmp(dir) = m_one - sin(ddv(dir) * m_pi / (m_two * width))**2
956 else
957 tmp(dir) = m_zero
958 end if
959 end if
960 mask_fn(ip) = mask_fn(ip) * tmp(dir)
961 end do
962 mask_fn(ip) = m_one - mask_fn(ip)
963
964 end select
965 end if
966 end do
967
968 case (m_step)
969 do ip = 1, mask%np
970 call mesh_r(mesh, ip, rr, coords=xx)
971 dd = rr - r(1)
972 if (dd > m_zero) then
973 if (dd < width) then
974 mask_fn(ip) = m_one
975 else
976 mask_fn(ip) = m_zero
977 end if
978 end if
979 end do
980
981 case (m_erf)
982
983 case default
984 message(1)="PhotoElectronSpectrum = pes_mask. Unrecognized mask type."
985 call messages_fatal(1, namespace=namespace)
986 end select
987
988
989
990 mask_fn(:) = m_one - mask_fn(:)
991
992
993 call pes_mask_mesh_to_cube(mask, mask_fn, mask%cM)
994
995 if (present(mask_sq)) mask_sq = real(mask%cM%zRS, real64)
996
997
998
999 safe_deallocate_a(mask_fn)
1000
1002
1003 end subroutine pes_mask_generate_mask_function
1004
1005
1006 ! --------------------------------------------------------
1007 subroutine pes_mask_apply_mask(mask, st)
1008 type(states_elec_t), intent(inout) :: st
1009 type(pes_mask_t), intent(in) :: mask
1010
1011 integer :: ik, ist, idim
1012 complex(real64), allocatable :: mmask(:), psi(:)
1013
1014 push_sub(pes_mask_apply_mask)
1015 safe_allocate(mmask(1:mask%mesh%np))
1016 safe_allocate(psi(1:mask%mesh%np))
1017
1018 call pes_mask_cube_to_mesh(mask, mask%cM, mmask)
1019
1020 do ik = st%d%kpt%start, st%d%kpt%end
1021 do ist = st%st_start, st%st_end
1022 do idim = 1, st%d%dim
1023 call states_elec_get_state(st, mask%mesh, idim, ist, ik, psi)
1024 psi(1:mask%mesh%np) = psi(1:mask%mesh%np) * mmask(1:mask%mesh%np)
1025 call states_elec_set_state(st, mask%mesh, idim, ist, ik, psi)
1026 end do
1027 end do
1028 end do
1029
1030 safe_deallocate_a(mmask)
1031
1032 pop_sub(pes_mask_apply_mask)
1033 end subroutine pes_mask_apply_mask
1034
1035
1036
1037
1038
1039
1040 ! ---------------------------------------------------------
1052 ! ---------------------------------------------------------
1053 subroutine pes_mask_volkov_time_evolution_wf(mask, space, mesh, kpoints, dt, iter, wf, ikpoint)
1054 type(pes_mask_t), intent(in) :: mask
1055 class(space_t), intent(in) :: space
1056 type(mesh_t), intent(in) :: mesh
1057 type(kpoints_t), intent(in) :: kpoints
1058 real(real64), intent(in) :: dt
1059 integer, intent(in) :: iter
1060 complex(real64), intent(inout) :: wf(:,:,:)
1061 integer, intent(in) :: ikpoint
1062
1063 integer :: ix, iy, iz
1064 real(real64) :: vec
1065 real(real64) :: kk(1:3), kpoint(1:3)
1066
1068
1069 kpoint = m_zero
1070 if (space%is_periodic()) then
1071 kpoint(1:mesh%box%dim) = kpoints%get_point(ikpoint)
1072 end if
1073
1074 do ix = 1, mask%ll(1)
1075 kk(1) = mask%Lk(ix + mask%fs_istart(1) - 1, 1)
1076 do iy = 1, mask%ll(2)
1077 kk(2) = mask%Lk(iy + mask%fs_istart(2) - 1, 2)
1078 do iz = 1, mask%ll(3)
1079 kk(3) = mask%Lk(iz + mask%fs_istart(3) - 1, 3)
1080 ! The k-points have the same sign as the vector potential consistently
1081 ! with what is done to generate the phase (hm%phase) in hamiltonian_elec_update()
1082 vec = sum(( kk(1:mesh%box%dim) &
1083 - kpoint(1:mesh%box%dim) &
1084 - mask%vec_pot(iter,1:mesh%box%dim)/p_c)**2) / m_two
1085 wf(ix, iy, iz) = wf(ix, iy, iz) * exp(-m_zi * dt * vec)
1086 end do
1087 end do
1088 end do
1089
1090
1091! do ix = 1, mask%ll(1)
1092! KK(1) = mask%Lk(ix + mask%fs_istart(1) - 1)
1093! do iy = 1, mask%ll(2)
1094! KK(2) = mask%Lk(iy + mask%fs_istart(2) - 1)
1095! do iz = 1, mask%ll(3)
1096! KK(3) = mask%Lk(iz + mask%fs_istart(3) - 1)
1097! vec = sum(( KK(1:mesh%box%dim) - mask%vec_pot(iter,1:mesh%box%dim)/P_C)**2) / M_TWO
1098! wf(ix, iy, iz) = wf(ix, iy, iz) * exp(-M_zI * dt * vec)
1099! end do
1100! end do
1101! end do
1105
1106
1107
1108
1109 !---------------------------------------------------------
1111 !---------------------------------------------------------
1112 subroutine pes_mask_x_to_k(mask, wfin, wfout)
1113 type(pes_mask_t), intent(in) :: mask
1114 complex(real64), intent(inout) :: wfin(:,:,:)
1115 complex(real64), intent(out) :: wfout(:,:,:)
1116
1117 type(cube_function_t) :: cf_tmp
1118 real(real64) :: norm
1119
1120 call profiling_in("PESMASK_X_to_K")
1121
1122 push_sub(pes_mask_x_to_k)
1123
1124 wfout=m_z0
1125
1126 select case (mask%pw_map_how)
1127
1128 case (pw_map_fft)
1129 call zfft_forward(mask%cube%fft, wfin, wfout)
1130
1131 case (pw_map_nfft)
1132 ! By definition NFFT forward transform generates a function defined on an
1133 ! unstructured grid from its Fourier components (defined on a cubic equispaced grid).
1134 ! This is not what we want since for us it is the real-space grid the one that can be
1135 ! unstructured.
1136 ! In order to preserve the possibility to have an unstructured rs-grid we use the
1137 ! backward transform and flip the sign of all the momenta (Lk = -Lk).
1138
1139 call zfft_backward(mask%cube%fft, wfin, wfout, norm)
1140 wfout = wfout * norm
1141
1142! call zfft_forward(mask%cube%fft, wfin, wfout)
1143
1144 case (pw_map_pfft)
1145 call zcube_function_alloc_rs(mask%cube, cf_tmp)
1146 call cube_function_alloc_fs(mask%cube, cf_tmp)
1147 cf_tmp%zRs = wfin
1148 call zfft_forward(mask%cube%fft, cf_tmp%zRs, cf_tmp%fs)
1149 wfout = cf_tmp%fs
1150 call zcube_function_free_rs(mask%cube, cf_tmp)
1151 call cube_function_free_fs(mask%cube, cf_tmp)
1152
1153 case (pw_map_pnfft)
1154 call zfft_backward(mask%cube%fft, wfin, wfout, norm)
1155 wfout = wfout * norm
1156
1157! call zfft_forward(mask%cube%fft, wfin, wfout)
1158
1159! call zfft_backward(mask%cube%fft, M_zI*conjg(wfin), wfout)
1160! wfout = M_zI*conjg(wfout)
1161
1162 case default
1163
1164 end select
1165
1166
1167 pop_sub(pes_mask_x_to_k)
1168 call profiling_out("PESMASK_X_to_K")
1169
1170 end subroutine pes_mask_x_to_k
1171
1172 ! ------------------------------------------------
1173 subroutine pes_mask_k_to_x(mask, wfin, wfout)
1174 type(pes_mask_t), intent(in) :: mask
1175 complex(real64), intent(inout) :: wfin(:,:,:)
1176 complex(real64), intent(out) :: wfout(:,:,:)
1177
1178 type(cube_function_t) :: cf_tmp
1179 real(real64) :: norm
1180
1181 call profiling_in("PESMASK_K_toX")
1182
1183 push_sub(pes_mask_k_to_x)
1184
1185 wfout=m_z0
1186
1187 select case (mask%pw_map_how)
1188
1189 case (pw_map_fft)
1190 call zfft_backward(mask%cube%fft, wfin,wfout)
1191
1192 case (pw_map_nfft)
1193 call zfft_forward(mask%cube%fft, wfin, wfout, norm)
1194 wfout = wfout / norm
1195! call zfft_backward(mask%cube%fft, wfin,wfout)
1196
1197 case (pw_map_pfft)
1198
1199 call zcube_function_alloc_rs(mask%cube, cf_tmp)
1200 call cube_function_alloc_fs(mask%cube, cf_tmp)
1201 cf_tmp%fs = wfin
1202 call zfft_backward(mask%cube%fft, cf_tmp%fs, cf_tmp%zRs)
1203 wfout = cf_tmp%zRs
1204 call zcube_function_free_rs(mask%cube, cf_tmp)
1205 call cube_function_free_fs(mask%cube, cf_tmp)
1206
1208 call zfft_forward(mask%cube%fft, wfin, wfout, norm)
1209 wfout = wfout / norm
1210
1211! call zfft_backward(mask%cube%fft, wfin,wfout)
1212
1213! call zfft_forward(mask%cube%fft, M_zI*conjg(wfin),wfout)
1214! wfout = M_zI*conjg(wfout)
1215
1216
1217 case default
1218
1219 end select
1220
1221
1222 pop_sub(pes_mask_k_to_x)
1223
1224 call profiling_out("PESMASK_K_toX")
1225
1226 end subroutine pes_mask_k_to_x
1227
1228 !---------------------------------------------------------
1229 subroutine pes_mask_mesh_to_cube(mask, mf, cf)
1230 type(pes_mask_t), intent(in) :: mask
1231 complex(real64), contiguous, intent(in) :: mf(:)
1232 type(cube_function_t), intent(inout) :: cf
1233
1234 push_sub(pes_mask_mesh_to_cube)
1235
1236 if (mask%cube%parallel_in_domains) then
1237 call zmesh_to_cube_parallel(mask%mesh, mf, mask%cube, cf, mask%mesh_cube_map)
1238 else
1239 call zmesh_to_cube(mask%mesh, mf, mask%cube, cf)
1240 end if
1241
1242 pop_sub(pes_mask_mesh_to_cube)
1243 end subroutine pes_mask_mesh_to_cube
1244
1245
1246 !---------------------------------------------------------
1247 subroutine pes_mask_cube_to_mesh(mask, cf, mf)
1248 type(pes_mask_t), intent(in) :: mask
1249 complex(real64), contiguous, intent(out):: mf(:)
1250 type(cube_function_t), intent(in) :: cf
1251
1252 push_sub(pes_mask_cube_to_mesh)
1253
1254 if (mask%cube%parallel_in_domains) then
1255 call zcube_to_mesh_parallel(mask%cube, cf, mask%mesh, mf, mask%mesh_cube_map)
1256 else
1257 call zcube_to_mesh(mask%cube, cf, mask%mesh, mf)
1258 end if
1259
1260 pop_sub(pes_mask_cube_to_mesh)
1261 end subroutine pes_mask_cube_to_mesh
1262
1263
1264 !---------------------------------------------------------
1265 !
1266 ! Performs all the dirty work
1267 !
1268 !---------------------------------------------------------
1269 subroutine pes_mask_calc(mask, namespace, space, mesh, st, kpoints, dt, iter)
1270 type(pes_mask_t), intent(inout) :: mask
1271 type(namespace_t), intent(in) :: namespace
1272 class(space_t), intent(in) :: space
1273 type(mesh_t), intent(in) :: mesh
1274 type(states_elec_t), intent(inout) :: st
1275 type(kpoints_t), intent(in) :: kpoints
1276 real(real64), intent(in) :: dt
1277 integer, intent(in) :: iter
1278
1279 integer :: idim, ist, ik
1280 type(cube_function_t):: cf1, cf2
1281 complex(real64), allocatable :: mf(:), psi(:)
1282
1283 real(real64) :: time
1284
1285
1286 call profiling_in("PESMASK_calc")
1287
1288 push_sub(pes_mask_calc)
1289
1290 time = iter *dt
1291
1292 if (time > mask%start_time) then ! record photoelectrons only after mask%start_time
1293
1294 call zcube_function_alloc_rs(mask%cube, cf1, force_alloc = .true.)
1295 call cube_function_alloc_fs(mask%cube, cf1, force_alloc = .true.)
1296 call zcube_function_alloc_rs(mask%cube, cf2, force_alloc = .true.)
1297 call cube_function_alloc_fs(mask%cube, cf2, force_alloc = .true.)
1298
1299 select case (mask%mode)
1300 case (pes_mask_mode_mask)
1301 if (mask%back_action .eqv. .true.) then
1302 safe_allocate(mf(1:mask%mesh%np_part))
1303 end if
1304 end select
1305
1306 safe_allocate(psi(1:mask%mesh%np_part))
1307
1308 do ik = st%d%kpt%start, st%d%kpt%end
1309 do ist = st%st_start, st%st_end
1310 do idim = 1, st%d%dim
1311
1312 call states_elec_get_state(st, mask%mesh, idim, ist, ik, psi)
1313
1314 cf1%zRs(:,:,:) = m_z0
1315 cf2%zRS(:,:,:) = m_z0
1316 cf1%Fs(:,:,:) = m_z0
1317 cf2%Fs(:,:,:) = m_z0
1318
1319 call pes_mask_mesh_to_cube(mask, psi, cf1)
1320
1321 select case (mask%mode)
1322 !-----------------------------------------
1323 ! Mask Method
1324 !----------------------------------------
1325 case (pes_mask_mode_mask)
1326
1327 cf1%zRs = (m_one - mask%cM%zRs) * cf1%zRs ! cf1 =(1-M)*\Psi_A(x,t2)
1328 call pes_mask_x_to_k(mask, cf1%zRs, cf2%Fs) ! cf2 = \tilde{\Psi}_A(k,t2)
1329
1330
1331 if (mask%filter_k) then ! apply a filter to the Fourier transform to remove unwanted energies
1332 assert(allocated(mask%Mk))
1333 cf2%Fs = cf2%Fs * mask%Mk
1334 end if
1335
1336
1337 cf1%Fs(:,:,:) = mask%k(:,:,:, idim, ist, ik) ! cf1 = \Psi_B(k,t1)
1338 mask%k(:,:,:, idim, ist, ik) = cf2%Fs(:,:,:) ! mask%k = \tilde{\Psi}_A(k,t2)
1339 call pes_mask_volkov_time_evolution_wf(mask, space, mesh, kpoints, dt,iter-1,cf1%Fs, & ! cf1 = \tilde{\Psi}_B(k,t2)
1340 st%d%get_kpoint_index(ik))
1341
1342 mask%k(:,:,:, idim, ist, ik) = mask%k(:,:,:, idim, ist, ik)&
1343 + cf1%Fs(:,:,:) ! mask%k = \tilde{\Psi}_A(k,t2) + \tilde{\Psi}_B(k,t2)
1344
1345 if (mask%back_action .eqv. .true.) then
1346
1347 ! Apply Back-action to wavefunction in A
1348 call pes_mask_k_to_x(mask, cf1%Fs, cf2%zRs) ! cf2 = \Psi_B(x,t2)
1349 call pes_mask_cube_to_mesh(mask, cf2, mf)
1350 psi(1:mask%mesh%np) = psi(1:mask%mesh%np) + mf(1:mask%mesh%np)
1351 call states_elec_set_state(st, mask%mesh, idim, ist, ik, psi)
1352
1353 ! Apply correction to wavefunction in B
1354 cf2%zRs= (mask%cM%zRs) * cf2%zRs ! cf2 = M*\Psi_B(x,t1)
1355 call pes_mask_x_to_k(mask, cf2%zRs, cf1%Fs)
1356
1357 mask%k(:,:,:, idim, ist, ik) = mask%k(:,:,:, idim, ist, ik) - cf1%Fs
1358
1359 end if
1360
1361
1362 !-----------------------------------------
1363 ! Passive Mask method
1364 !----------------------------------------
1366
1367
1368 cf1%zRs = (m_one-mask%cM%zRs) * cf1%zRs
1369 call pes_mask_x_to_k(mask, cf1%zRs, cf2%Fs)
1370
1371 mask%k(:,:,:, idim, ist, ik) = cf2%Fs(:,:,:)
1372
1373
1374 case default
1375 !Program should die before coming here
1376 write(message(1),'(a)') "PhotoElectroSpectrum = pes_mask. Unrecognized calculation mode."
1377 call messages_fatal(1, namespace=namespace)
1378
1379 end select
1380
1381
1382 end do
1383 end do
1384 end do
1385
1386 safe_deallocate_a(psi)
1387
1388 call zcube_function_free_rs(mask%cube, cf1)
1389 call cube_function_free_fs(mask%cube, cf1)
1390 call zcube_function_free_rs(mask%cube, cf2)
1391 call cube_function_free_fs(mask%cube, cf2)
1392
1393 select case (mask%mode)
1394 case (pes_mask_mode_mask)
1395 if (mask%back_action .eqv. .true.) then
1396 safe_deallocate_a(mf)
1397 end if
1398 end select
1399
1400 end if ! time > mask%start_time
1401
1402
1403 if (mask%mode == pes_mask_mode_mask) call pes_mask_apply_mask(mask, st) !apply the mask to all the KS orbitals
1404
1405
1406
1407
1408 pop_sub(pes_mask_calc)
1409
1410 call profiling_out("PESMASK_calc")
1411 end subroutine pes_mask_calc
1412
1413#include "pes_mask_out_inc.F90"
1414
1415end module pes_mask_oct_m
1416
1417!! Local Variables:
1418!! mode: f90
1419!! coding: utf-8
1420!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:182
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
integer, parameter, public not_absorbing
subroutine, public zmesh_to_cube(mesh, mf, cube, cf)
Convert a function from the mesh to the cube.
subroutine, public zcube_to_mesh_parallel(cube, cf, mesh, mf, map)
subroutine, public zcube_to_mesh(cube, cf, mesh, mf)
Convert a function from the cube to the mesh.
subroutine, public zmesh_to_cube_parallel(mesh, mf, cube, cf, map)
The next two subroutines convert a function between the normal mesh and the cube in parallel.
subroutine, public zcube_function_free_rs(cube, cf)
Deallocates the real space grid.
subroutine, public zcube_function_alloc_rs(cube, cf, in_device, force_alloc)
Allocates locally the real space grid, if PFFT library is not used. Otherwise, it assigns the PFFT re...
subroutine, public cube_init(cube, nn, namespace, space, spacing, coord_system, fft_type, fft_library, dont_optimize, nn_out, mpi_grp, need_partition, tp_enlarge, blocksize)
Definition: cube.F90:204
subroutine, public cube_end(cube)
Definition: cube.F90:387
subroutine, public cube_init_cube_map(cube, mesh)
Definition: cube.F90:824
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
type(lasers_t) function, pointer, public list_get_lasers(partners)
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:120
integer, parameter, public fft_complex
Definition: fft.F90:174
integer, parameter, public fftlib_nfft
Definition: fft.F90:179
integer, parameter, public fftlib_pnfft
Definition: fft.F90:179
integer, parameter, public fftlib_pfft
Definition: fft.F90:179
integer, parameter, public fftlib_fftw
Definition: fft.F90:179
subroutine, public cube_function_free_fs(cube, cf)
Deallocates the Fourier space grid.
subroutine, public cube_function_alloc_fs(cube, cf, force_alloc)
Allocates locally the Fourier space grid, if PFFT library is not used. Otherwise, it assigns the PFFT...
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
complex(real64), parameter, public m_z0
Definition: global.F90:201
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:233
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
This module defines classes and functions for interaction partners.
Definition: io.F90:116
integer, parameter, public e_field_vector_potential
Definition: lasers.F90:179
integer pure elemental function, public laser_kind(laser)
Definition: lasers.F90:717
subroutine, public laser_field(laser, field, time)
Retrieves the value of either the electric or the magnetic field. If the laser is given by a scalar p...
Definition: lasers.F90:1117
System information (time, memory, sysname)
Definition: loct.F90:117
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
subroutine, public mesh_cube_parallel_map_end(this)
subroutine, public mesh_cube_parallel_map_init(this, mesh, cube)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:342
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1023
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
this module contains the low-level part of the output system
Definition: output_low.F90:117
subroutine, public parse_block_string(blk, l, c, res, convert_to_c)
Definition: parser.F90:810
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
subroutine, public pes_mask_calc(mask, namespace, space, mesh, st, kpoints, dt, iter)
Definition: pes_mask.F90:1365
integer, parameter, public pes_mask_mode_backaction
Definition: pes_mask.F90:244
integer, parameter, public integrate_r
Definition: pes_mask.F90:264
subroutine, public pes_mask_dump(mask, namespace, restart, st, ierr)
Definition: pes_mask.F90:3600
subroutine, public pes_mask_output_full_mapm(pesK, file, namespace, space, Lk, ll, how, mesh, pmesh)
Definition: pes_mask.F90:2286
subroutine pes_mask_volkov_time_evolution_wf(mask, space, mesh, kpoints, dt, iter, wf, ikpoint)
Propagate in time a wavefunction in momentum space with the Volkov Hamiltonian.
Definition: pes_mask.F90:1149
integer, parameter pw_map_pnfft
use PNFFT
Definition: pes_mask.F90:249
subroutine, public pes_mask_cube_to_mesh(mask, cf, mf)
Definition: pes_mask.F90:1343
subroutine, public pes_mask_output_full_mapm_cut(pesK, file, namespace, ll, dim, pol, dir, integrate, pos, Lk, pmesh)
Definition: pes_mask.F90:2395
subroutine, public pes_mask_output(mask, gr, st, outp, namespace, space, file, ions, iter)
This routine is the main routine dedicated to the output of PES data.
Definition: pes_mask.F90:3366
integer, parameter pw_map_nfft
non-equispaced fft (NFFT)
Definition: pes_mask.F90:249
subroutine pes_mask_apply_mask(mask, st)
Definition: pes_mask.F90:1103
integer, parameter, public integrate_kx
Definition: pes_mask.F90:264
subroutine, public pes_mask_init(mask, namespace, space, mesh, box, st, ext_partners, abs_boundaries, max_iter, dt)
Definition: pes_mask.F90:279
subroutine, public pes_mask_pmesh(namespace, dim, kpoints, ll, LG, pmesh, idxZero, krng, Lp)
Definition: pes_mask.F90:1567
subroutine, public pes_mask_generate_mask_function(mask, namespace, mesh, shape, R, mask_sq)
Generate the mask function on the cubic mesh containing the simulation box.
Definition: pes_mask.F90:991
subroutine pes_mask_generate_mask(mask, namespace, mesh)
Generate the mask function on the cubic mesh containing the simulation box.
Definition: pes_mask.F90:974
subroutine, public pes_mask_output_power_totalm(pesK, file, namespace, Lk, ll, dim, Emax, Estep, interpolate)
Definition: pes_mask.F90:3172
integer, parameter, public integrate_phi
Definition: pes_mask.F90:264
subroutine pes_mask_generate_filter(mask, cutOff)
Generate the momentum-space filter.
Definition: pes_mask.F90:930
subroutine, public pes_mask_output_ar_spherical_cut_m(pesK, file, namespace, Lk, ll, dim, dir, Emin, Emax, Estep)
Definition: pes_mask.F90:2901
subroutine, public pes_mask_load(mask, namespace, restart, st, ierr)
Definition: pes_mask.F90:3691
subroutine, public pes_mask_mesh_to_cube(mask, mf, cf)
Definition: pes_mask.F90:1325
subroutine, public pes_mask_x_to_k(mask, wfin, wfout)
Project the wavefunction on plane waves.
Definition: pes_mask.F90:1208
integer, parameter, public pes_mask_mode_passive
Definition: pes_mask.F90:244
integer, parameter m_step
Definition: pes_mask.F90:255
subroutine pes_mask_generate_lk(mask)
Definition: pes_mask.F90:911
integer, parameter, public integrate_theta
Definition: pes_mask.F90:264
subroutine, public pes_mask_output_ar_plane_m(pesK, file, namespace, Lk, ll, dim, dir, Emax, Estep)
Definition: pes_mask.F90:2771
integer, parameter pw_map_pfft
use PFFT
Definition: pes_mask.F90:249
integer, parameter out
Definition: pes_mask.F90:260
subroutine, public pes_mask_map_from_states(restart, st, ll, pesK, krng, Lp, istin)
Definition: pes_mask.F90:1807
subroutine, public pes_mask_end(mask)
Definition: pes_mask.F90:879
integer, parameter, public integrate_ky
Definition: pes_mask.F90:264
integer, parameter m_erf
Definition: pes_mask.F90:255
subroutine, public pes_mask_k_to_x(mask, wfin, wfout)
Definition: pes_mask.F90:1269
integer, parameter, public integrate_kz
Definition: pes_mask.F90:264
subroutine, public pes_mask_read_info(dir, namespace, dim, Emax, Estep, ll, Lk, RR)
Read pes info.
Definition: pes_mask.F90:3496
subroutine, public pes_mask_output_ar_polar_m(pesK, file, namespace, Lk, ll, dim, dir, Emax, Estep)
Definition: pes_mask.F90:2652
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:225
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
Class implementing a parallelepiped box. Currently this is restricted to a rectangular cuboid (all th...
Class implementing a spherical box.
Definition: box_sphere.F90:134
Describes mesh distribution to nodes.
Definition: mesh.F90:187
The states_elec_t class contains all electronic wave functions.
int true(void)