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
115 real(real64), allocatable, public :: vec_pot(:,:)
116
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
194
195 type(block_t) :: blk
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
202
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
207
209
210 mask%mesh => mesh
212 if (space%is_periodic()) then
213 call messages_experimental("PES_mask with periodic dimensions", namespace=namespace)
214 end if
216
217 write(message(1),'(a,i1,a)') 'Info: Calculating PES using mask technique.'
218 call messages_info(1, namespace=namespace)
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
230
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)
260
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 ! print *,mpi_world%rank, "mask%mesh%mpi_grp%comm", mask%mesh%mpi_grp%comm, mask%mesh%mpi_grp%size
454 ! print *,mpi_world%rank, "mask%cube%mpi_grp%comm", mask%cube%mpi_grp%comm, mask%cube%mpi_grp%size
455
456
457! mask%ll(1) = mask%cube%fs_n(3)
458! mask%ll(2) = mask%cube%fs_n(1)
459! mask%ll(3) = mask%cube%fs_n(2)
460! Note: even if tempting, setting
461! mask%ll(1:3) = mask%cube%fs_n(1:3)
462! results in the wrong index mapping! (1->3, 2->1, 3->2)
463! Well.. very much against intuition it turns out to be
464 mask%ll(1:3) = mask%cube%rs_n(1:3)
465
466 mask%fft = mask%cube%fft
467 if (mask%mesh%parallel_in_domains .and. mask%cube%parallel_in_domains) then
468 call mesh_cube_parallel_map_init(mask%mesh_cube_map, mask%mesh, mask%cube)
469 end if
470
471 case (pw_map_fft)
472 call cube_init(mask%cube, mask%ll, namespace, space, &
473 mesh%spacing, mesh%coord_system, &
474 fft_type = fft_complex, fft_library = fftlib_fftw, &
475 nn_out = ll)
476 call cube_init_cube_map(mask%cube, mesh)
477 mask%ll = ll
478 mask%fft = mask%cube%fft
479
480
481 case (pw_map_nfft)
482
483 !NFFT initialization
484 ! we just add 2 points for the enlarged region
485 if (.not. is_close(mask%enlarge_2p(1), m_one)) mask%ll(1:space%dim) = mask%ll(1:space%dim) + 2
486
487 call cube_init(mask%cube, mask%ll, namespace, space, &
488 mesh%spacing, mesh%coord_system, &
489 fft_type = fft_complex, fft_library = fftlib_nfft, &
490 nn_out = ll, tp_enlarge = mask%enlarge_2p)
491 call cube_init_cube_map(mask%cube, mesh)
492
493 mask%ll = ll
494 mask%fft = mask%cube%fft
495
496
497 case (pw_map_pnfft)
498
499 if (.not. is_close(mask%enlarge_2p(1), m_one)) mask%ll(1:space%dim) = mask%ll(1:space%dim) + 2
500
501 call cube_init(mask%cube, mask%ll, namespace, space, &
502 mask%mesh%spacing, mask%mesh%coord_system, &
503 fft_type = fft_complex, fft_library = fftlib_pnfft, &
504 nn_out = ll, tp_enlarge = mask%enlarge_2p, &
505 mpi_grp = mask%mesh%mpi_grp, need_partition=.true.)
506 call cube_init_cube_map(mask%cube, mesh)
507
508 mask%ll(1:3) = mask%cube%fs_n(1:3)
509
510 mask%fft = mask%cube%fft
511 if (mask%mesh%parallel_in_domains .and. mask%cube%parallel_in_domains) then
512 call mesh_cube_parallel_map_init(mask%mesh_cube_map, mask%mesh, mask%cube)
513 end if
514
515 case default
516 !Program should die before coming here
517 write(message(1),'(a)') "PESMaskPlaneWaveProjection unrecognized option."
518 call messages_fatal(1, namespace=namespace)
519
520 end select
521
522 !Indices
523
524 if (mask%pw_map_how == pw_map_pfft) then
525 mask%fs_istart = mask%cube%rs_istart
526 mask%fs_n = mask%cube%rs_n
527 mask%fs_n_global = mask%cube%rs_n_global
528 else
529 mask%fs_istart = mask%cube%fs_istart
530 mask%fs_n = mask%cube%fs_n
531 mask%fs_n_global = mask%cube%fs_n_global
532 end if
533
534 if (debug%info) then
535 print *,mpi_world%rank, "mask%ll ", mask%ll(:)
536 print *,mpi_world%rank, "mask%cube%fs_n_global(:) ", mask%cube%fs_n_global(:)
537 print *,mpi_world%rank, "mask%cube%fs_n(:) ", mask%cube%fs_n(:)
538 print *,mpi_world%rank, "mask%cube%fs_istart(:) ", mask%cube%fs_istart(:)
539 print *,mpi_world%rank, "mask%cube%rs_n_global(:) ", mask%cube%rs_n_global(:)
540 print *,mpi_world%rank, "mask%cube%rs_n(:) ", mask%cube%rs_n(:)
541 print *,mpi_world%rank, "mask%cube%rs_istart(:) ", mask%cube%rs_istart(:)
542 end if
543
544! print *, mpi_world%rank, " mask%ll", mask%ll(1:3), "states -",st%st_start,st%st_end
545 !Allocations
546 call zcube_function_alloc_rs(mask%cube, mask%cM, force_alloc = .true.)
547
548 safe_allocate(mask%Lk(1:maxval(mask%fs_n_global(:)),1:3))
549 mask%Lk(:,:) = m_zero
550
551 st1 = st%st_start
552 st2 = st%st_end
553 k1 = st%d%kpt%start
554 k2 = st%d%kpt%end
555 safe_allocate(mask%k(1:mask%ll(1),1:mask%ll(2),1:mask%ll(3),1:st%d%dim,st1:st2,k1:k2))
556 mask%k = m_z0
557
558
559 ! generate the map between mesh and cube
560 call pes_mask_generate_lk(mask) ! generate the physical momentum vector
561
562
563
564!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
565!! Mask Function options
566!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
567
568 safe_allocate(mask%mask_R(1:2))
569
570 !%Variable PESMaskShape
571 !%Type integer
572 !%Default m_sin2
573 !%Section Time-Dependent::PhotoElectronSpectrum
574 !%Description
575 !% The mask function shape.
576 !%Option m_sin2 1
577 !% sin2 mask.
578 !%Option m_step 2
579 !%step function.
580 !%Option m_erf 3
581 !%Error function. Not Implemented.
582 !%End
583 call parse_variable(namespace, 'PESMaskShape', defaultmask, mask%shape)
584 if (.not. varinfo_valid_option('PESMaskShape', mask%shape)) call messages_input_error(namespace, 'PESMaskShape')
585 call messages_print_var_option("PESMaskShape", mask%shape, namespace=namespace)
586
587 !%Variable PESMaskSize
588 !%Type block
589 !%Section Time-Dependent::PhotoElectronSpectrum
590 !%Description
591 !% Set the size of the mask function.
592 !% Here you can set the inner (R1) and outer (R2) radius by setting
593 !% the block as follows:
594 !%
595 !% <tt>%PESMaskSize
596 !% <br>&nbsp;&nbsp; R1 | R2 | "user-defined"
597 !% <br>%</tt>
598 !%
599 !% The optional 3rd column is a user-defined expression for the mask
600 !% function. For example, <i>r</i> creates a spherical mask (which is the
601 !% default for <tt>BoxShape = sphere</tt>). Note, values R2 larger than
602 !% the box size may lead in this case to unexpected reflection
603 !% behaviours.
604 !%End
605 cols_pesmask_block = 0
606 if (parse_block(namespace, 'PESMaskSize', blk) == 0) then
607 cols_pesmask_block = parse_block_cols(blk, 0)
608 end if
609
610 mask%user_def = .false.
611
612 select case (cols_pesmask_block)
613 case (0)
614 select type (box => box)
615 type is (box_sphere_t)
616 mask%mask_R(1) = box%radius/m_two
617 mask%mask_R(2) = box%radius
618 message(1) = "Input: PESMaskSize R(1) and R(2) not specified. Using default values for spherical mask."
619 type is (box_parallelepiped_t)
620 mask%mask_R(1) = box%half_length(1)/m_two
621 mask%mask_R(2) = box%half_length(1)
622 message(1) = "Input: PESMaskSize R(1) and R(2) not specified. Using default values for cubic mask."
623 end select
624 call messages_info(1, namespace=namespace)
625 case (1)
626 call parse_block_float(blk, 0, 0, mask%mask_R(1), units_inp%length)
627 select type (box => box)
628 type is (box_sphere_t)
629 mask%mask_R(2) = box%radius
630 message(1) = "Input: PESMaskSize R(2) not specified. Using default value for spherical mask."
631 type is (box_parallelepiped_t)
632 mask%mask_R(2) = box%half_length(1)
633 message(1) = "Input: PESMaskSize R(2) not specified. Using default value for cubic mask."
634 end select
635 call messages_info(1, namespace=namespace)
636 case (2)
637 call parse_block_float(blk, 0, 0, mask%mask_R(1), units_inp%length)
638 call parse_block_float(blk, 0, 1, mask%mask_R(2), units_inp%length)
639
640 select type (box => box)
641 type is (box_sphere_t)
642 if (mask%mask_R(2) > box%radius) mask%mask_R(2) = box%radius
643 message(1) = "Info: using spherical mask."
644 type is (box_parallelepiped_t)
645 if (mask%mask_R(2) > box%half_length(1)) mask%mask_R(2) = box%half_length(1)
646 message(1) = "Info: using cubic mask."
647 end select
648 call messages_info(1, namespace=namespace)
649
650 case (3)
651 mask%user_def = .true.
652 safe_allocate(mask%ufn(1:mask%np))
653 mask%ufn = m_zero
654 call parse_block_float(blk, 0, 0, mask%mask_R(1), units_inp%length)
655 call parse_block_float(blk, 0, 1, mask%mask_R(2), units_inp%length)
656 call parse_block_string(blk, 0, 2, user_def_expr)
657 do ip = 1, mask%np
658 xx = m_zero
659 xx(1:space%dim) = mesh%x(ip, 1:space%dim)
660 r = units_from_atomic(units_inp%length, norm2(xx(1:space%dim)))
661 do idim = 1, space%dim
662 xx(idim) = units_from_atomic(units_inp%length, xx(idim))
663 end do
664 call parse_expression(ufn_re, ufn_im, space%dim, xx, r, m_zero, user_def_expr)
665 mask%ufn(ip) = ufn_re
666 end do
667 message(1) = "Input: using user-defined mask function from expression:"
668 write(message(2),'(a,a)') ' R = ', trim(user_def_expr)
669 call messages_info(2, namespace=namespace)
670 end select
671
672 call parse_block_end(blk)
673
674 write(message(1),'(a,es10.3,3a)') &
675 " R1 = ", units_from_atomic(units_inp%length, mask%mask_R(1)),&
676 ' [', trim(units_abbrev(units_inp%length)), ']'
677 write(message(2),'(a,es10.3,3a)') &
678 " R2 = ", units_from_atomic(units_inp%length, mask%mask_R(2)),&
679 ' [', trim(units_abbrev(units_inp%length)), ']'
680 call messages_info(2, namespace=namespace)
681
682
683 call pes_mask_generate_mask(mask, namespace, mesh)
684
685 !%Variable PESMaskFilterCutOff
686 !%Type float
687 !%Default -1
688 !%Section Time-Dependent::PhotoElectronSpectrum
689 !%Description
690 !% In calculation with <tt>PESMaskMode = fullmask_mode</tt> and NFFT, spurious frequencies
691 !% may lead to numerical instability of the algorithm. This option gives the possibility
692 !% to filter out the unwanted components by setting an energy cut-off.
693 !% If <tt>PESMaskFilterCutOff = -1</tt> no filter is applied.
694 !%End
695 call parse_variable(namespace, 'PESMaskFilterCutOff', -m_one, pcutoff, unit = units_inp%energy)
696
697 mask%filter_k = .false.
698
699 if (pcutoff > m_zero) then
700 call messages_print_var_value("PESMaskFilterCutOff", pcutoff, unit = units_out%energy, namespace=namespace)
701 mask%filter_k = .true.
702
703 safe_allocate(mask%Mk(1:mask%ll(1),1:mask%ll(2),1:mask%ll(3)))
704
705 call pes_mask_generate_filter(mask,pcutoff)
706 end if
707
708!! Output
709
710 !%Variable PESMaskIncludePsiA
711 !%Type logical
712 !%Default false
713 !%Section Time-Dependent::PhotoElectronSpectrum
714 !%Description
715 !% Add the contribution of <math>\Psi_A</math> in the mask region to the photo-electron spectrum.
716 !% Literally adds the Fourier components of:
717 !% <math>\Theta(r-R1) \Psi_A(r)</math>
718 !% with <math>\Theta</math> being the Heaviside step function.
719 !% With this option PES will contain all the contributions starting from the inner
720 !% radius <math>R1</math>. Use this option to improve convergence with respect to the box size
721 !% and total simulation time.
722 !% Note: Carefully choose <math>R1</math> in order to avoid contributions from returning electrons.
723 !%End
724 call parse_variable(namespace, 'PESMaskIncludePsiA', .false., mask%add_psia)
725 if (mask%add_psia) then
726 message(1)= "Input: Include contribution from Psi_A."
727 call messages_info(1, namespace=namespace)
728 end if
729
730
731 !%Variable PESMaskSpectEnergyMax
732 !%Type float
733 !%Default maxval(mask%Lk)<math>^2/2</math>
734 !%Section Time-Dependent::PhotoElectronSpectrum
735 !%Description
736 !% The maximum energy for the PES spectrum.
737 !%End
738 maxe = m_epsilon
739 do idim = 1, space%dim
740 tmp = maxval(mask%Lk(1:mask%ll(idim),1:space%dim))**m_two/m_two
741 if (tmp > maxe) maxe = tmp
742 end do
743 call parse_variable(namespace, 'PESMaskSpectEnergyMax', maxe, mask%energyMax, unit = units_inp%energy)
744 call messages_print_var_value("PESMaskSpectEnergyMax", mask%energyMax, unit = units_out%energy, namespace=namespace)
745
746 !%Variable PESMaskSpectEnergyStep
747 !%Type float
748 !%Section Time-Dependent::PhotoElectronSpectrum
749 !%Description
750 !% The PES spectrum energy step.
751 !%End
752 deltae = minval(mask%Lk(2,1:space%dim) - mask%Lk(1,1:space%dim))**m_two/m_two
753 call parse_variable(namespace, 'PESMaskSpectEnergyStep', deltae, mask%energyStep, unit = units_inp%energy)
754 call messages_print_var_value("PESMaskSpectEnergyStep", mask%energyStep, unit = units_out%energy, namespace=namespace)
755
756
757!! Set external fields
758
759 safe_allocate(mask%vec_pot(0:max_iter,1:3))
760 mask%vec_pot=m_zero
761
762 lasers => list_get_lasers(ext_partners)
763 if(associated(lasers)) then
764 do il = 1, lasers%no_lasers
765 select case (laser_kind(lasers%lasers(il)))
767 do it = 1, max_iter
768 field=m_zero
769 call laser_field(lasers%lasers(il), field, it*dt)
770 ! We must sum with a -1 sign to account for the
771 ! electron charge.
772 mask%vec_pot(it,:)= mask%vec_pot(it,:) - field(:)
773 end do
774
775 case default
776 write(message(1),'(a)') 'PESMask should work only with TDExternalFields = vector_potential.'
777 write(message(2),'(a)') 'Unless PESMaskMode = passive_mode the results are likely to be wrong. '
778 call messages_warning(2, namespace=namespace)
779
780 end select
781 end do
782 end if
783
784
785
786 ! Compensate the sign for the forward <-> backward FT inversion
787 ! used with NFFT. See pes_mask_X_to_K comment for more info
788 if (mask%pw_map_how == pw_map_pnfft .or. &
789 mask%pw_map_how == pw_map_nfft) mask%Lk = - mask%Lk
790
791
792 pop_sub(pes_mask_init)
793 end subroutine pes_mask_init
794
795
796 ! ---------------------------------------------------------
797 subroutine pes_mask_end(mask)
798 type(pes_mask_t), intent(inout) :: mask
799
800 push_sub(pes_mask_end)
801
802 safe_deallocate_a(mask%k)
803
804 safe_deallocate_a(mask%vec_pot)
805 safe_deallocate_a(mask%mask_R)
806 safe_deallocate_a(mask%Lk)
807
808
809 if (mask%filter_k) then
810 safe_deallocate_a(mask%Mk)
811 end if
812
813 if (mask%mesh%parallel_in_domains .and. mask%cube%parallel_in_domains) then
814 call mesh_cube_parallel_map_end(mask%mesh_cube_map)
815 end if
816
817 call zcube_function_free_rs(mask%cube, mask%cM)
818 call cube_end(mask%cube)
819
820 if (mask%user_def) then
821 safe_deallocate_a(mask%ufn)
822 end if
823
824 pop_sub(pes_mask_end)
825 end subroutine pes_mask_end
826
827
828 ! --------------------------------------------------------
829 subroutine pes_mask_generate_lk(mask)
830 type(pes_mask_t), intent(inout) :: mask
831
832 integer :: ii, dim
833
834 push_sub(pes_mask_generate_lk)
835
836 dim = mask%mesh%box%dim
837
838 do ii = 1, maxval(mask%ll(:))
839 mask%Lk(ii,1:dim)= matmul(mask%cube%latt%klattice_primitive(1:dim,1:dim), mask%cube%Lfs(ii,1:dim))
840 end do
841
842 pop_sub(pes_mask_generate_lk)
843 end subroutine pes_mask_generate_lk
844
845 ! ======================================
847 ! ======================================
848 subroutine pes_mask_generate_filter(mask,cutOff)
849 type(pes_mask_t), intent(inout) :: mask
850 real(real64), intent(in) ::cutoff
851
852
853 integer :: kx,ky,kz, power
854 real(real64) :: kk(3), ee, emax
855
857
858 mask%Mk = m_zero
859
860 emax = maxval(mask%Lk(:,:))**2 / m_two
861
862 power = 8
863
864 do kx = 1, mask%ll(1)
865 kk(1) = mask%Lk(kx,1)
866 do ky = 1, mask%ll(2)
867 kk(2) = mask%Lk(ky,2)
868 do kz = 1, mask%ll(3)
869 kk(3) = mask%Lk(kz,3)
870
871 ee = sum(kk(1:mask%mesh%box%dim)**2) / m_two
872
873 if (ee > cutoff .and. ee < emax) then
874 mask%Mk(kx,ky,kz) = m_one * sin((emax-ee) * m_pi / (m_two * (cutoff)))**power
875 else
876 if (ee <= cutoff) mask%Mk(kx,ky,kz) = m_one
877 end if
878
879 end do
880 end do
881 end do
882
883
885 end subroutine pes_mask_generate_filter
886
887
888 ! --------------------------------------------------------
891 ! ---------------------------------------------------------
892 subroutine pes_mask_generate_mask(mask, namespace, mesh)
893 type(pes_mask_t), intent(inout) :: mask
894 type(namespace_t), intent(in) :: namespace
895 type(mesh_t), intent(in) :: mesh
896
897 push_sub(pes_mask_generate_mask)
898
899 call pes_mask_generate_mask_function(mask, namespace, mesh, mask%shape, mask%mask_R)
902
903 end subroutine pes_mask_generate_mask
904
905 ! --------------------------------------------------------
908 ! ---------------------------------------------------------
909 subroutine pes_mask_generate_mask_function(mask, namespace, mesh, shape, R, mask_sq)
910 type(pes_mask_t), intent(inout) :: mask
911 type(namespace_t), intent(in) :: namespace
912 type(mesh_t), intent(in) :: mesh
913 integer, intent(in) :: shape
914 real(real64), intent(in) :: R(2)
915 real(real64), optional, intent(out) :: mask_sq(:,:,:)
916
917 integer :: ip, dir
918 real(real64) :: width
919 real(real64) :: xx(1:mesh%box%dim), rr, dd, ddv(1:mesh%box%dim), tmp(1:mesh%box%dim)
920 complex(real64), allocatable :: mask_fn(:)
921
923
924
925 ! generate the mask function on the mesh
926 safe_allocate(mask_fn(1:mask%np))
927
928 mask_fn = m_zero
929 width = r(2) - r(1)
930 xx = m_zero
931
932 select case (shape)
933 case (m_sin2)
934 do ip = 1, mask%np
935 call mesh_r(mesh, ip, rr, coords=xx)
936
937 if (mask%user_def) then
938 dd = mask%ufn(ip) - r(1)
939 if (dd > m_zero) then
940 if (mask%ufn(ip) < r(2)) then
941 mask_fn(ip) = m_one * sin(dd * m_pi / (m_two * (width)))**2
942 else
943 mask_fn(ip) = m_one
944 end if
945 end if
946
947 else ! mask%user_def == .false.
948
949 select type (box => mesh%box)
950 type is (box_sphere_t)
951 dd = rr - r(1)
952 if (dd > m_zero) then
953 if (dd < width) then
954 mask_fn(ip) = m_one * sin(dd * m_pi / (m_two * (width)))**2
955 else
956 mask_fn(ip) = m_one
957 end if
958 end if
959
960 type is (box_parallelepiped_t)
961
962 ! We are filling from the center opposite to the spherical case
963 tmp = m_one
964 mask_fn(ip) = m_one
965 ddv = abs(xx) - r(1)
966 do dir=1, mesh%box%dim
967 if (ddv(dir) > m_zero) then
968 if (ddv(dir) < width) then
969 tmp(dir) = m_one - sin(ddv(dir) * m_pi / (m_two * width))**2
970 else
971 tmp(dir) = m_zero
972 end if
973 end if
974 mask_fn(ip) = mask_fn(ip) * tmp(dir)
975 end do
976 mask_fn(ip) = m_one - mask_fn(ip)
977
978 end select
979 end if
980 end do
981
982 case (m_step)
983 do ip = 1, mask%np
984 call mesh_r(mesh, ip, rr, coords=xx)
985 dd = rr - r(1)
986 if (dd > m_zero) then
987 if (dd < width) then
988 mask_fn(ip) = m_one
989 else
990 mask_fn(ip) = m_zero
991 end if
992 end if
993 end do
994
995 case (m_erf)
996
997 case default
998 message(1)="PhotoElectronSpectrum = pes_mask. Unrecognized mask type."
999 call messages_fatal(1, namespace=namespace)
1000 end select
1001
1002
1003
1004 mask_fn(:) = m_one - mask_fn(:)
1005
1006
1007 call pes_mask_mesh_to_cube(mask, mask_fn, mask%cM)
1008
1009 if (present(mask_sq)) mask_sq = real(mask%cM%zRS, real64)
1010
1011
1012
1013 safe_deallocate_a(mask_fn)
1014
1016
1017 end subroutine pes_mask_generate_mask_function
1018
1019
1020 ! --------------------------------------------------------
1021 subroutine pes_mask_apply_mask(mask, st)
1022 type(states_elec_t), intent(inout) :: st
1023 type(pes_mask_t), intent(in) :: mask
1024
1025 integer :: ik, ist, idim
1026 complex(real64), allocatable :: mmask(:), psi(:)
1027
1028 push_sub(pes_mask_apply_mask)
1029 safe_allocate(mmask(1:mask%mesh%np))
1030 safe_allocate(psi(1:mask%mesh%np))
1031
1032 call pes_mask_cube_to_mesh(mask, mask%cM, mmask)
1033
1034 do ik = st%d%kpt%start, st%d%kpt%end
1035 do ist = st%st_start, st%st_end
1036 do idim = 1, st%d%dim
1037 call states_elec_get_state(st, mask%mesh, idim, ist, ik, psi)
1038 psi(1:mask%mesh%np) = psi(1:mask%mesh%np) * mmask(1:mask%mesh%np)
1039 call states_elec_set_state(st, mask%mesh, idim, ist, ik, psi)
1040 end do
1041 end do
1042 end do
1043
1044 safe_deallocate_a(mmask)
1045
1046 pop_sub(pes_mask_apply_mask)
1047 end subroutine pes_mask_apply_mask
1048
1049
1050
1051
1052
1053
1054 ! ---------------------------------------------------------
1066 ! ---------------------------------------------------------
1067 subroutine pes_mask_volkov_time_evolution_wf(mask, space, mesh, kpoints, dt, iter, wf, ikpoint)
1068 type(pes_mask_t), intent(in) :: mask
1069 class(space_t), intent(in) :: space
1070 type(mesh_t), intent(in) :: mesh
1071 type(kpoints_t), intent(in) :: kpoints
1072 real(real64), intent(in) :: dt
1073 integer, intent(in) :: iter
1074 complex(real64), intent(inout) :: wf(:,:,:)
1075 integer, intent(in) :: ikpoint
1076
1077 integer :: ix, iy, iz
1078 real(real64) :: vec
1079 real(real64) :: kk(1:3), kpoint(1:3)
1080
1082
1083 kpoint = m_zero
1084 if (space%is_periodic()) then
1085 kpoint(1:mesh%box%dim) = kpoints%get_point(ikpoint)
1086 end if
1087
1088 do ix = 1, mask%ll(1)
1089 kk(1) = mask%Lk(ix + mask%fs_istart(1) - 1, 1)
1090 do iy = 1, mask%ll(2)
1091 kk(2) = mask%Lk(iy + mask%fs_istart(2) - 1, 2)
1092 do iz = 1, mask%ll(3)
1093 kk(3) = mask%Lk(iz + mask%fs_istart(3) - 1, 3)
1094 ! The k-points have the same sign as the vector potential consistently
1095 ! with what is done to generate the phase (hm%phase) in hamiltonian_elec_update()
1096 vec = sum(( kk(1:mesh%box%dim) &
1097 - kpoint(1:mesh%box%dim) &
1098 - mask%vec_pot(iter,1:mesh%box%dim)/p_c)**2) / m_two
1099 wf(ix, iy, iz) = wf(ix, iy, iz) * exp(-m_zi * dt * vec)
1100 end do
1101 end do
1102 end do
1103
1104
1105! do ix = 1, mask%ll(1)
1106! KK(1) = mask%Lk(ix + mask%fs_istart(1) - 1)
1107! do iy = 1, mask%ll(2)
1108! KK(2) = mask%Lk(iy + mask%fs_istart(2) - 1)
1109! do iz = 1, mask%ll(3)
1110! KK(3) = mask%Lk(iz + mask%fs_istart(3) - 1)
1111! vec = sum(( KK(1:mesh%box%dim) - mask%vec_pot(iter,1:mesh%box%dim)/P_C)**2) / M_TWO
1112! wf(ix, iy, iz) = wf(ix, iy, iz) * exp(-M_zI * dt * vec)
1113! end do
1114! end do
1115! end do
1116
1119
1120
1121
1122
1123 !---------------------------------------------------------
1125 !---------------------------------------------------------
1126 subroutine pes_mask_x_to_k(mask, wfin, wfout)
1127 type(pes_mask_t), intent(in) :: mask
1128 complex(real64), intent(inout) :: wfin(:,:,:)
1129 complex(real64), intent(out) :: wfout(:,:,:)
1130
1131 type(cube_function_t) :: cf_tmp
1132 real(real64) :: norm
1133
1134 call profiling_in("PESMASK_X_to_K")
1135
1136 push_sub(pes_mask_x_to_k)
1137
1138 wfout=m_z0
1139
1140 select case (mask%pw_map_how)
1141
1142 case (pw_map_fft)
1143 call zfft_forward(mask%cube%fft, wfin, wfout)
1144
1145 case (pw_map_nfft)
1146 ! By definition NFFT forward transform generates a function defined on an
1147 ! unstructured grid from its Fourier components (defined on a cubic equispaced grid).
1148 ! This is not what we want since for us it is the real-space grid the one that can be
1149 ! unstructured.
1150 ! In order to preserve the possibility to have an unstructured rs-grid we use the
1151 ! backward transform and flip the sign of all the momenta (Lk = -Lk).
1152
1153 call zfft_backward(mask%cube%fft, wfin, wfout, norm)
1154 wfout = wfout * norm
1155
1156! call zfft_forward(mask%cube%fft, wfin, wfout)
1157
1158 case (pw_map_pfft)
1159 call zcube_function_alloc_rs(mask%cube, cf_tmp)
1160 call cube_function_alloc_fs(mask%cube, cf_tmp)
1161 cf_tmp%zRs = wfin
1162 call zfft_forward(mask%cube%fft, cf_tmp%zRs, cf_tmp%fs)
1163 wfout = cf_tmp%fs
1164 call zcube_function_free_rs(mask%cube, cf_tmp)
1165 call cube_function_free_fs(mask%cube, cf_tmp)
1166
1167 case (pw_map_pnfft)
1168 call zfft_backward(mask%cube%fft, wfin, wfout, norm)
1169 wfout = wfout * norm
1170
1171! call zfft_forward(mask%cube%fft, wfin, wfout)
1172
1173! call zfft_backward(mask%cube%fft, M_zI*conjg(wfin), wfout)
1174! wfout = M_zI*conjg(wfout)
1175
1176 case default
1177
1178 end select
1179
1180
1181 pop_sub(pes_mask_x_to_k)
1182 call profiling_out("PESMASK_X_to_K")
1183
1184 end subroutine pes_mask_x_to_k
1185
1186 ! ------------------------------------------------
1187 subroutine pes_mask_k_to_x(mask, wfin, wfout)
1188 type(pes_mask_t), intent(in) :: mask
1189 complex(real64), intent(inout) :: wfin(:,:,:)
1190 complex(real64), intent(out) :: wfout(:,:,:)
1191
1192 type(cube_function_t) :: cf_tmp
1193 real(real64) :: norm
1194
1195 call profiling_in("PESMASK_K_toX")
1196
1198
1199 wfout=m_z0
1200
1201 select case (mask%pw_map_how)
1202
1203 case (pw_map_fft)
1204 call zfft_backward(mask%cube%fft, wfin,wfout)
1205
1206 case (pw_map_nfft)
1207 call zfft_forward(mask%cube%fft, wfin, wfout, norm)
1208 wfout = wfout / norm
1209! call zfft_backward(mask%cube%fft, wfin,wfout)
1210
1211 case (pw_map_pfft)
1212
1213 call zcube_function_alloc_rs(mask%cube, cf_tmp)
1214 call cube_function_alloc_fs(mask%cube, cf_tmp)
1215 cf_tmp%fs = wfin
1216 call zfft_backward(mask%cube%fft, cf_tmp%fs, cf_tmp%zRs)
1217 wfout = cf_tmp%zRs
1218 call zcube_function_free_rs(mask%cube, cf_tmp)
1219 call cube_function_free_fs(mask%cube, cf_tmp)
1220
1221 case (pw_map_pnfft)
1222 call zfft_forward(mask%cube%fft, wfin, wfout, norm)
1223 wfout = wfout / norm
1224
1225! call zfft_backward(mask%cube%fft, wfin,wfout)
1226
1227! call zfft_forward(mask%cube%fft, M_zI*conjg(wfin),wfout)
1228! wfout = M_zI*conjg(wfout)
1229
1230
1231 case default
1232
1233 end select
1234
1235
1236 pop_sub(pes_mask_k_to_x)
1237
1238 call profiling_out("PESMASK_K_toX")
1239
1240 end subroutine pes_mask_k_to_x
1241
1242 !---------------------------------------------------------
1243 subroutine pes_mask_mesh_to_cube(mask, mf, cf)
1244 type(pes_mask_t), intent(in) :: mask
1245 complex(real64), contiguous, intent(in) :: mf(:)
1246 type(cube_function_t), intent(inout) :: cf
1247
1248 push_sub(pes_mask_mesh_to_cube)
1249
1250 if (mask%cube%parallel_in_domains) then
1251 call zmesh_to_cube_parallel(mask%mesh, mf, mask%cube, cf, mask%mesh_cube_map)
1252 else
1253 call zmesh_to_cube(mask%mesh, mf, mask%cube, cf)
1254 end if
1255
1256 pop_sub(pes_mask_mesh_to_cube)
1257 end subroutine pes_mask_mesh_to_cube
1259
1260 !---------------------------------------------------------
1261 subroutine pes_mask_cube_to_mesh(mask, cf, mf)
1262 type(pes_mask_t), intent(in) :: mask
1263 complex(real64), contiguous, intent(out):: mf(:)
1264 type(cube_function_t), intent(in) :: cf
1265
1266 push_sub(pes_mask_cube_to_mesh)
1267
1268 if (mask%cube%parallel_in_domains) then
1269 call zcube_to_mesh_parallel(mask%cube, cf, mask%mesh, mf, mask%mesh_cube_map)
1270 else
1271 call zcube_to_mesh(mask%cube, cf, mask%mesh, mf)
1272 end if
1273
1274 pop_sub(pes_mask_cube_to_mesh)
1275 end subroutine pes_mask_cube_to_mesh
1276
1277
1278 !---------------------------------------------------------
1279 !
1280 ! Performs all the dirty work
1281 !
1282 !---------------------------------------------------------
1283 subroutine pes_mask_calc(mask, namespace, space, mesh, st, kpoints, dt, iter)
1284 type(pes_mask_t), intent(inout) :: mask
1285 type(namespace_t), intent(in) :: namespace
1286 class(space_t), intent(in) :: space
1287 type(mesh_t), intent(in) :: mesh
1288 type(states_elec_t), intent(inout) :: st
1289 type(kpoints_t), intent(in) :: kpoints
1290 real(real64), intent(in) :: dt
1291 integer, intent(in) :: iter
1292
1293 integer :: idim, ist, ik
1294 type(cube_function_t):: cf1, cf2
1295 complex(real64), allocatable :: mf(:), psi(:)
1296
1297 real(real64) :: time
1298
1299
1300 call profiling_in("PESMASK_calc")
1301
1302 push_sub(pes_mask_calc)
1303
1304 time = iter *dt
1305
1306 if (time > mask%start_time) then ! record photoelectrons only after mask%start_time
1307
1308 call zcube_function_alloc_rs(mask%cube, cf1, force_alloc = .true.)
1309 call cube_function_alloc_fs(mask%cube, cf1, force_alloc = .true.)
1310 call zcube_function_alloc_rs(mask%cube, cf2, force_alloc = .true.)
1311 call cube_function_alloc_fs(mask%cube, cf2, force_alloc = .true.)
1312
1313 select case (mask%mode)
1314 case (pes_mask_mode_mask)
1315 if (mask%back_action .eqv. .true.) then
1316 safe_allocate(mf(1:mask%mesh%np_part))
1317 end if
1318 end select
1319
1320 safe_allocate(psi(1:mask%mesh%np_part))
1321
1322 do ik = st%d%kpt%start, st%d%kpt%end
1323 do ist = st%st_start, st%st_end
1324 do idim = 1, st%d%dim
1325
1326 call states_elec_get_state(st, mask%mesh, idim, ist, ik, psi)
1327
1328 cf1%zRs(:,:,:) = m_z0
1329 cf2%zRS(:,:,:) = m_z0
1330 cf1%Fs(:,:,:) = m_z0
1331 cf2%Fs(:,:,:) = m_z0
1333 call pes_mask_mesh_to_cube(mask, psi, cf1)
1334
1335 select case (mask%mode)
1336 !-----------------------------------------
1337 ! Mask Method
1338 !----------------------------------------
1339 case (pes_mask_mode_mask)
1340
1341 cf1%zRs = (m_one - mask%cM%zRs) * cf1%zRs ! cf1 =(1-M)*\Psi_A(x,t2)
1342 call pes_mask_x_to_k(mask, cf1%zRs, cf2%Fs) ! cf2 = \tilde{\Psi}_A(k,t2)
1343
1344
1345 if (mask%filter_k) then ! apply a filter to the Fourier transform to remove unwanted energies
1346 assert(allocated(mask%Mk))
1347 cf2%Fs = cf2%Fs * mask%Mk
1348 end if
1349
1350
1351 cf1%Fs(:,:,:) = mask%k(:,:,:, idim, ist, ik) ! cf1 = \Psi_B(k,t1)
1352 mask%k(:,:,:, idim, ist, ik) = cf2%Fs(:,:,:) ! mask%k = \tilde{\Psi}_A(k,t2)
1353 call pes_mask_volkov_time_evolution_wf(mask, space, mesh, kpoints, dt,iter-1,cf1%Fs, & ! cf1 = \tilde{\Psi}_B(k,t2)
1354 st%d%get_kpoint_index(ik))
1355
1356 mask%k(:,:,:, idim, ist, ik) = mask%k(:,:,:, idim, ist, ik)&
1357 + cf1%Fs(:,:,:) ! mask%k = \tilde{\Psi}_A(k,t2) + \tilde{\Psi}_B(k,t2)
1358
1359 if (mask%back_action .eqv. .true.) then
1360
1361 ! Apply Back-action to wavefunction in A
1362 call pes_mask_k_to_x(mask, cf1%Fs, cf2%zRs) ! cf2 = \Psi_B(x,t2)
1363 call pes_mask_cube_to_mesh(mask, cf2, mf)
1364 psi(1:mask%mesh%np) = psi(1:mask%mesh%np) + mf(1:mask%mesh%np)
1365 call states_elec_set_state(st, mask%mesh, idim, ist, ik, psi)
1366
1367 ! Apply correction to wavefunction in B
1368 cf2%zRs= (mask%cM%zRs) * cf2%zRs ! cf2 = M*\Psi_B(x,t1)
1369 call pes_mask_x_to_k(mask, cf2%zRs, cf1%Fs)
1370
1371 mask%k(:,:,:, idim, ist, ik) = mask%k(:,:,:, idim, ist, ik) - cf1%Fs
1372
1373 end if
1374
1375
1376 !-----------------------------------------
1377 ! Passive Mask method
1378 !----------------------------------------
1380
1381
1382 cf1%zRs = (m_one-mask%cM%zRs) * cf1%zRs
1383 call pes_mask_x_to_k(mask, cf1%zRs, cf2%Fs)
1384
1385 mask%k(:,:,:, idim, ist, ik) = cf2%Fs(:,:,:)
1386
1387
1388 case default
1389 !Program should die before coming here
1390 write(message(1),'(a)') "PhotoElectroSpectrum = pes_mask. Unrecognized calculation mode."
1391 call messages_fatal(1, namespace=namespace)
1392
1393 end select
1394
1395
1396 end do
1397 end do
1398 end do
1399
1400 safe_deallocate_a(psi)
1401
1402 call zcube_function_free_rs(mask%cube, cf1)
1403 call cube_function_free_fs(mask%cube, cf1)
1404 call zcube_function_free_rs(mask%cube, cf2)
1405 call cube_function_free_fs(mask%cube, cf2)
1406
1407 select case (mask%mode)
1408 case (pes_mask_mode_mask)
1409 if (mask%back_action .eqv. .true.) then
1410 safe_deallocate_a(mf)
1411 end if
1412 end select
1413
1414 end if ! time > mask%start_time
1415
1416
1417 if (mask%mode == pes_mask_mode_mask) call pes_mask_apply_mask(mask, st) !apply the mask to all the KS orbitals
1418
1419
1420
1421
1422 pop_sub(pes_mask_calc)
1423
1424 call profiling_out("PESMASK_calc")
1425 end subroutine pes_mask_calc
1426
1427#include "pes_mask_out_inc.F90"
1428
1429end module pes_mask_oct_m
1430
1431!! Local Variables:
1432!! mode: f90
1433!! coding: utf-8
1434!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:180
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)
The next two subroutines convert a function between the normal mesh and the cube.
subroutine, public zcube_to_mesh_parallel(cube, cf, mesh, mf, map)
subroutine, public zcube_to_mesh(cube, cf, mesh, mf)
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:202
subroutine, public cube_end(cube)
Definition: cube.F90:378
subroutine, public cube_init_cube_map(cube, mesh)
Definition: cube.F90:823
type(debug_t), save, public debug
Definition: debug.F90:154
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
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:118
integer, parameter, public fft_complex
Definition: fft.F90:177
integer, parameter, public fftlib_nfft
Definition: fft.F90:182
integer, parameter, public fftlib_pnfft
Definition: fft.F90:182
integer, parameter, public fftlib_pfft
Definition: fft.F90:182
integer, parameter, public fftlib_fftw
Definition: fft.F90:182
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:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
complex(real64), parameter, public m_z0
Definition: global.F90:197
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:223
real(real64), parameter, public m_one
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module defines classes and functions for interaction partners.
Definition: io.F90:114
integer, parameter, public e_field_vector_potential
Definition: lasers.F90:176
integer pure elemental function, public laser_kind(laser)
Definition: lasers.F90:711
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:1111
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
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:118
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:336
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
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
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
this module contains the output system
Definition: output_low.F90:115
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public pes_mask_calc(mask, namespace, space, mesh, st, kpoints, dt, iter)
Definition: pes_mask.F90:1355
integer, parameter, public pes_mask_mode_backaction
Definition: pes_mask.F90:242
integer, parameter, public integrate_r
Definition: pes_mask.F90:262
subroutine, public pes_mask_dump(mask, namespace, restart, st, ierr)
Definition: pes_mask.F90:3603
subroutine, public pes_mask_output_full_mapm(pesK, file, namespace, space, Lk, ll, how, mesh, pmesh)
Definition: pes_mask.F90:2276
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:1139
integer, parameter pw_map_pnfft
use PNFFT
Definition: pes_mask.F90:247
subroutine, public pes_mask_cube_to_mesh(mask, cf, mf)
Definition: pes_mask.F90:1333
subroutine, public pes_mask_output_full_mapm_cut(pesK, file, namespace, ll, dim, pol, dir, integrate, pos, Lk, pmesh)
Definition: pes_mask.F90:2398
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:3369
integer, parameter pw_map_nfft
non-equispaced fft (NFFT)
Definition: pes_mask.F90:247
subroutine pes_mask_apply_mask(mask, st)
Definition: pes_mask.F90:1093
integer, parameter, public integrate_kx
Definition: pes_mask.F90:262
subroutine, public pes_mask_init(mask, namespace, space, mesh, box, st, ext_partners, abs_boundaries, max_iter, dt)
Definition: pes_mask.F90:277
subroutine, public pes_mask_pmesh(namespace, dim, kpoints, ll, LG, pmesh, idxZero, krng, Lp)
Definition: pes_mask.F90:1557
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:981
subroutine pes_mask_generate_mask(mask, namespace, mesh)
Generate the mask function on the cubic mesh containing the simulation box.
Definition: pes_mask.F90:964
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 pes_mask_generate_filter(mask, cutOff)
Generate the momentum-space filter.
Definition: pes_mask.F90:920
subroutine, public pes_mask_output_ar_spherical_cut_m(pesK, file, namespace, Lk, ll, dim, dir, Emin, Emax, Estep)
Definition: pes_mask.F90:2904
subroutine, public pes_mask_load(mask, namespace, restart, st, ierr)
Definition: pes_mask.F90:3698
subroutine, public pes_mask_mesh_to_cube(mask, mf, cf)
Definition: pes_mask.F90:1315
subroutine, public pes_mask_x_to_k(mask, wfin, wfout)
Project the wavefunction on plane waves.
Definition: pes_mask.F90:1198
integer, parameter, public pes_mask_mode_passive
Definition: pes_mask.F90:242
integer, parameter m_step
Definition: pes_mask.F90:253
subroutine pes_mask_generate_lk(mask)
Definition: pes_mask.F90:901
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
integer, parameter pw_map_pfft
use PFFT
Definition: pes_mask.F90:247
integer, parameter out
Definition: pes_mask.F90:258
subroutine, public pes_mask_map_from_states(restart, st, ll, pesK, krng, Lp, istin)
Definition: pes_mask.F90:1797
subroutine, public pes_mask_end(mask)
Definition: pes_mask.F90:869
integer, parameter, public integrate_ky
Definition: pes_mask.F90:262
integer, parameter m_erf
Definition: pes_mask.F90:253
subroutine, public pes_mask_k_to_x(mask, wfin, wfout)
Definition: pes_mask.F90:1259
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 profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
Class implementing a parallelepiped box. Currently this is restricted to a rectangular cuboid (all th...
Class implementing a spherical box.
Definition: box_sphere.F90:132
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
int true(void)