37 use,
intrinsic :: iso_fortran_env
85 integer :: nkpnts_start, nkpnts_end
87 integer :: nk_start, nk_end
90 integer :: nstepsthetak, nstepsphik
91 real(real64) :: thetak_rng(1:2)
92 real(real64) :: phik_rng(1:2)
93 integer :: nstepsomegak
94 integer :: nstepsomegak_start, nstepsomegak_end
95 real(real64) :: ktransf(1:3,1:3)
97 real(real64),
allocatable :: klinear(:,:)
102 real(real64),
allocatable :: kcoords_cub(:,:,:)
103 real(real64),
allocatable :: kcoords_sph(:,:,:)
107 integer,
public :: surf_shape
109 integer :: nsrfcpnts_start, nsrfcpnts_end
110 real(real64),
allocatable :: srfcnrml(:,:)
111 real(real64),
allocatable :: rcoords(:,:)
112 integer,
allocatable :: srfcpnt(:)
113 integer,
allocatable :: rankmin(:)
115 complex(real64),
allocatable :: ylm_r(:,:,:)
116 complex(real64),
allocatable :: ylm_k(:,:,:)
117 real(real64),
allocatable :: j_l(:,:)
118 real(real64) :: radius
119 integer,
allocatable :: face_idx_range(:,:)
122 complex(real64),
allocatable :: bvk_phase(:,:)
123 complex(real64),
allocatable :: expkr(:,:,:,:)
124 complex(real64),
allocatable :: expkr_perp(:,:)
125 real(real64),
allocatable :: LLr(:,:)
126 integer,
allocatable :: NN(:,:)
127 integer,
allocatable :: Lkpuvz_inv(:,:,:)
137 complex(real64),
allocatable :: wf(:,:,:,:,:)
138 complex(real64),
allocatable :: gwf(:,:,:,:,:,:)
139 real(real64),
allocatable :: veca(:,:)
140 complex(real64),
allocatable :: conjgphase_prev(:,:)
141 complex(real64),
allocatable :: spctramp_cub(:,:,:,:)
142 complex(real64),
allocatable :: spctramp_sph(:,:,:,:,:)
144 integer,
public :: ll(3)
147 type(mesh_interpolation_t) :: interp
149 logical :: parallel_in_momentum
150 logical :: arpes_grid
151 logical :: surf_interp
152 logical :: use_symmetry
153 logical :: runtime_output
154 logical :: anisotrpy_correction
156 integer :: par_strategy
162 integer,
parameter :: &
167 integer,
parameter :: &
176 subroutine pes_flux_init(this, namespace, space, mesh, st, ext_partners, kpoints, abb, save_iter, max_iter)
185 integer,
intent(in) :: save_iter
186 integer,
intent(in) :: max_iter
189 real(real64) :: border(space%dim)
190 real(real64) :: offset(space%dim)
191 integer :: stst, stend, kptst, kptend, sdim, mdim, pdim
196 integer :: nstepsphir, nstepsthetar
198 integer :: default_shape
199 real(real64) :: fc_ptdens
201 integer :: par_strategy
202 logical :: use_symmetry
210 kptst = st%d%kpt%start
211 kptend = st%d%kpt%end
214 pdim = space%periodic_dim
216 this%surf_interp = .false.
220 if(
associated(lasers))
then
221 do il = 1, lasers%no_lasers
223 message(1) =
't-surff only works in velocity gauge.'
229 message(1) =
'Info: Calculating PES using t-surff technique.'
254 if (space%is_periodic())
then
256 else if (mdim <= 2)
then
257 default_shape = pes_cubic
259 select type (box => mesh%box)
261 default_shape = pes_cubic
265 call parse_variable(namespace,
'PES_Flux_Shape', default_shape, this%surf_shape)
270 message(1) =
'Spherical grid works only in 3d.'
283 if (this%surf_shape == pes_cubic)
then
284 call parse_variable(namespace,
'PES_Flux_AnisotropyCorrection', .false., this%anisotrpy_correction)
287 this%anisotrpy_correction = .false.
302 if (
parse_block(namespace,
'PES_Flux_Offset', blk) == 0)
then
324 if (this%surf_shape == pes_cubic .or. this%surf_shape ==
pes_plane)
then
344 if (
parse_block(namespace,
'PES_Flux_Lsize', blk) == 0)
then
349 border(1:mdim) = int(border(1:mdim) / mesh%spacing(1:mdim)) * mesh%spacing(1:mdim)
353 border(mdim) = mesh%box%bounding_box_l(mdim) *
m_half
354 if (space%is_periodic())
then
356 border(1:pdim)= mesh%box%bounding_box_l(1:pdim) *
m_two
357 call parse_variable(namespace,
'PES_Flux_Lsize', border(mdim), border(mdim))
359 border(mdim) =
floor(border(mdim) / mesh%spacing(mdim)) * mesh%spacing(mdim)
361 call parse_variable(namespace,
'PES_Flux_Lsize', border(mdim), border(mdim))
362 border(1:mdim - 1) = border(mdim)
363 border(1:mdim) =
floor(border(1:mdim) / mesh%spacing(1:mdim)) * mesh%spacing(1:mdim)
366 select type (box => mesh%box)
370 border(1:mdim) = box%half_length(1:mdim) *
m_half
372 call messages_write(
'PES_Flux_Lsize not specified. No default values available for this box shape.')
374 call messages_write(
'Specify the location of the parallelepiped with block PES_Flux_Lsize.')
377 call messages_write(
'PES_Flux_Lsize not specified. Using default values.')
391 this%surf_interp = .
true.
398 this%surf_interp = .
true.
411 select type (box => mesh%box)
413 this%radius = box%radius
415 this%radius = minval(box%half_length(1:mdim))
417 message(1) =
'PES_Flux_Radius not specified. No default values available for this box shape.'
418 message(2) =
'Specify the radius of the sphere with variable PES_Flux_Radius.'
421 message(1) =
'PES_Flux_Radius not specified. Using default values.'
433 call parse_variable(namespace,
'PES_Flux_StepsThetaR', 2 * this%lmax + 1, nstepsthetar)
444 call parse_variable(namespace,
'PES_Flux_StepsPhiR', 2 * this%lmax + 1, nstepsphir)
456 if (this%surf_shape == pes_cubic .or. this%surf_shape ==
pes_plane)
then
491 this%par_strategy = option__pes_flux_parallelization__pf_none
492 if (mesh%parallel_in_domains)
then
495 this%par_strategy = option__pes_flux_parallelization__pf_time &
496 + option__pes_flux_parallelization__pf_surface
498 this%par_strategy = option__pes_flux_parallelization__pf_surface
499 if (space%dim == 1) this%par_strategy = option__pes_flux_parallelization__pf_time
501 par_strategy = this%par_strategy
502 call parse_variable(namespace,
'PES_Flux_Parallelization', par_strategy, this%par_strategy)
510 if (
bitand(this%par_strategy, option__pes_flux_parallelization__pf_surface) /= 0 .and. &
511 bitand(this%par_strategy, option__pes_flux_parallelization__pf_momentum) /= 0)
then
512 call messages_input_error(namespace,
'PES_Flux_Parallelization',
"Cannot combine pf_surface and pf_momentum")
515 write(
message(1),
'(a)')
'Input: [PES_Flux_Parallelization = '
516 if (
bitand(this%par_strategy, option__pes_flux_parallelization__pf_none) /= 0)
then
519 if (
bitand(this%par_strategy, option__pes_flux_parallelization__pf_time) /= 0)
then
522 if (
bitand(this%par_strategy, option__pes_flux_parallelization__pf_momentum) /= 0)
then
525 if (
bitand(this%par_strategy, option__pes_flux_parallelization__pf_surface) /= 0)
then
535 if (
bitand(this%par_strategy, option__pes_flux_parallelization__pf_surface) /= 0)
then
537 call pes_flux_distribute(1, this%nsrfcpnts, this%nsrfcpnts_start, this%nsrfcpnts_end, mesh%mpi_grp%comm)
550 this%nsrfcpnts_start = 1
551 this%nsrfcpnts_end = this%nsrfcpnts
578 use_symmetry = .false.
579 if ((this%surf_shape == pes_cubic .or. this%surf_shape ==
pes_plane) &
581 if (this%surf_shape ==
pes_spherical .and. this%kgrid == pes_polar) use_symmetry = .
true.
582 call parse_variable(namespace,
'PES_Flux_UseSymmetries', use_symmetry, this%use_symmetry)
589 safe_allocate(this%ylm_r(0:this%lmax, -this%lmax:this%lmax, this%nsrfcpnts_start:this%nsrfcpnts_end))
592 if (this%nsrfcpnts_start > 0)
then
593 do isp = this%nsrfcpnts_start, this%nsrfcpnts_end
596 call ylmr_cmplx(this%rcoords(1:3, isp), ll, mm, this%ylm_r(ll, mm, isp))
597 this%ylm_r(ll, mm, isp) = conjg(this%ylm_r(ll, mm, isp))
617 call parse_variable(namespace,
'PES_Flux_RuntimeOutput', .false., this%runtime_output)
625 this%max_iter = max_iter
626 this%save_iter = save_iter
630 if (mesh%parallel_in_domains .and.
bitand(this%par_strategy, option__pes_flux_parallelization__pf_time) /= 0)
then
631 this%tdsteps = mesh%mpi_grp%size
638 safe_allocate(this%wf(stst:stend, 1:sdim, kptst:kptend, 0:this%nsrfcpnts, 1:this%tdsteps))
641 safe_allocate(this%gwf(stst:stend, 1:sdim, kptst:kptend, 0:this%nsrfcpnts, 1:this%tdsteps, 1:mdim))
644 safe_allocate(this%veca(1:mdim, 1:this%tdsteps))
648 safe_allocate(this%spctramp_sph(stst:stend, 1:sdim, kptst:kptend, 1:this%nk, 1:this%nstepsomegak))
649 this%spctramp_sph =
m_z0
651 safe_allocate(this%conjgphase_prev(1:this%nk, 1:this%nstepsomegak))
654 safe_allocate(this%spctramp_cub(stst:stend, 1:sdim, kptst:kptend, 1:this%nkpnts))
655 this%spctramp_cub =
m_z0
657 safe_allocate(this%conjgphase_prev(1:this%nkpnts, kptst:kptend))
661 this%conjgphase_prev =
m_z1
681 safe_deallocate_a(this%kcoords_sph)
682 safe_deallocate_a(this%ylm_k)
683 safe_deallocate_a(this%j_l)
684 safe_deallocate_a(this%ylm_r)
685 safe_deallocate_a(this%conjgphase_prev)
686 safe_deallocate_a(this%spctramp_sph)
688 safe_deallocate_a(this%kcoords_cub)
689 safe_deallocate_a(this%conjgphase_prev)
690 safe_deallocate_a(this%spctramp_cub)
692 if (.not. this%surf_interp)
then
693 safe_deallocate_a(this%srfcpnt)
695 safe_deallocate_a(this%rankmin)
697 safe_deallocate_a(this%face_idx_range)
698 safe_deallocate_a(this%LLr)
699 safe_deallocate_a(this%NN)
701 safe_deallocate_a(this%expkr)
702 safe_deallocate_a(this%expkr_perp)
703 safe_deallocate_a(this%bvk_phase)
705 safe_deallocate_a(this%Lkpuvz_inv)
709 safe_deallocate_a(this%klinear)
711 safe_deallocate_a(this%srfcnrml)
712 safe_deallocate_a(this%rcoords)
714 safe_deallocate_a(this%wf)
715 safe_deallocate_a(this%gwf)
716 safe_deallocate_a(this%veca)
725 class(
space_t),
intent(in) :: space
728 type(mpi_comm),
intent(in) :: comm
729 logical,
optional,
intent(in) :: post
731 integer :: mdim, pdim
732 integer :: kptst, kptend
734 integer :: ll, mm, idim, idir, ispin, nspin
735 integer :: ikk, ith, iph, iomk,ie, ik1, ik2, ik3, kgrid_block_dim
736 real(real64) :: kmax(space%dim), kmin(space%dim), kact, thetak, phik
739 real(real64) :: emin, emax, de , kvec(1:3), xx(3)
740 integer :: nkp_out, nkmin, nkmax
743 real(real64),
allocatable :: gpoints(:,:), gpoints_reduced(:,:)
744 real(real64) :: dk(1:3), kpoint(1:3), dthetak, dphik
745 logical :: use_enegy_grid, arpes_grid
749 kptst = st%d%kpt%start
750 kptend = st%d%kpt%end
752 pdim = space%periodic_dim
774 select case (this%surf_shape)
781 if (mdim == 1) kgrid = pes_polar
788 call parse_variable(namespace,
'PES_Flux_Momenutum_Grid', kgrid, this%kgrid)
801 if (space%is_periodic())
then
809 if (this%surf_shape == pes_cubic)
then
810 if (space%is_periodic())
then
831 use_enegy_grid = .false.
833 if (
parse_block(namespace,
'PES_Flux_EnergyGrid', blk) == 0)
then
844 use_enegy_grid = .
true.
866 call parse_variable(namespace,
'PES_Flux_ARPES_grid', .false., this%arpes_grid)
867 if (.not. use_enegy_grid .or. this%arpes_grid)
then
878 if (
parse_block(namespace,
'PES_Flux_Kmax', blk) == 0)
then
883 kgrid_block_dim = mdim
885 message(1) =
'Wrong block format for PES_Flux_Kmax and non-cartesian grid'
902 if (
parse_block(namespace,
'PES_Flux_Kmin', blk) == 0)
then
907 kgrid_block_dim = mdim
909 message(1) =
'Wrong block format for PES_Flux_Kmin and non-cartesian grid'
925 call parse_variable(namespace,
'PES_Flux_DeltaK', 0.02_real64, this%dk)
930 do idim = 1, kgrid_block_dim
931 if (kgrid_block_dim == 1)
then
950 if (this%kgrid == pes_polar)
then
973 this%nstepsthetak = 45
975 if (
parse_block(namespace,
'PES_Flux_ThetaK', blk) == 0)
then
981 if (this%thetak_rng(idim) <
m_zero .or. this%thetak_rng(idim) >
m_pi)
then
995 call parse_variable(namespace,
'PES_Flux_StepsThetaK', 45, this%nstepsthetak)
1020 this%nstepsphik = 90
1022 if (mdim == 1) this%phik_rng(1:2) = (/
m_pi,
m_zero/)
1023 if (
parse_block(namespace,
'PES_Flux_PhiK', blk) == 0)
then
1029 if (this%phik_rng(idim) <
m_zero .or. this%phik_rng(idim) >
m_two *
m_pi)
then
1044 call parse_variable(namespace,
'PES_Flux_StepsPhiK', 90, this%nstepsphik)
1046 if (this%nstepsphik == 0) this%nstepsphik = 1
1051 if (mdim == 3) dthetak = abs(this%thetak_rng(2) - this%thetak_rng(1)) / (this%nstepsthetak)
1052 dphik = abs(this%phik_rng(2) - this%phik_rng(1)) / (this%nstepsphik)
1057 this%nstepsthetak = 0
1059 this%nstepsomegak = this%nstepsphik
1062 this%nstepsthetak = 0
1063 this%nstepsomegak = this%nstepsphik
1066 if (this%nstepsthetak <= 1)
then
1068 this%nstepsthetak = 1
1072 do ith = 0, this%nstepsthetak
1073 thetak = ith * dthetak + this%thetak_rng(1)
1074 do iph = 0, this%nstepsphik - 1
1075 phik = iph * dphik + this%phik_rng(1)
1080 this%nstepsomegak = iomk
1087 write(
message(1),
'(a)')
"Polar momentum grid:"
1090 write(
message(1),
'(a,f12.6,a,f12.6,a, i6)') &
1091 " Theta = (", this%thetak_rng(1),
",",this%thetak_rng(2), &
1092 "); n = ", this%nstepsthetak
1095 write(
message(1),
'(a,f12.6,a,f12.6,a, i6)') &
1096 " Phi = (", this%phik_rng(1),
",",this%phik_rng(2), &
1097 "); n = ", this%nstepsphik
1102 if (use_enegy_grid)
then
1103 this%nk = nint(abs(emax - emin) / de)
1105 this%nk = nint(abs(kmax(1) - kmin(1)) / this%dk)
1107 this%nkpnts = this%nstepsomegak * this%nk
1109 this%ll(1) = this%nk
1110 this%ll(2) = this%nstepsphik
1111 this%ll(3) = this%nstepsthetak
1112 this%ll(mdim+1:3) = 1
1114 safe_allocate(this%klinear(1:this%nk, 1))
1120 dk(1:mdim) =
m_one/kpoints%nik_axis(1:mdim)
1122 this%arpes_grid = .false.
1123 if (space%is_periodic())
then
1134 arpes_grid = kpoints%have_zero_weight_path()
1135 call parse_variable(namespace,
'PES_Flux_ARPES_grid', arpes_grid, this%arpes_grid)
1144 if (kpoints%have_zero_weight_path())
then
1146 if (this%arpes_grid)
then
1147 nkmax =
floor(emax / de)
1148 nkmin =
floor(emin / de)
1151 nkmax =
floor(kmax(mdim) / this%dk)
1152 nkmin =
floor(kmin(mdim) / this%dk)
1156 this%ll(mdim) = abs(nkmax - nkmin) + 1
1158 this%nk = this%ll(mdim)
1163 if (.not. this%arpes_grid)
then
1164 this%ll(1:mdim) =
floor(abs(kmax(1:mdim) - kmin(1:mdim))/this%dk) + 1
1165 this%nk = maxval(this%ll(1:mdim))
1169 nkmax =
floor(emax / de)
1170 nkmin =
floor(emin / de)
1171 this%nk = abs(nkmax - nkmin) + 1
1173 this%ll(1:pdim) =
floor(abs(kmax(1:pdim) - kmin(1:pdim))/this%dk) + 1
1174 this%ll(mdim) = this%nk
1178 safe_allocate(this%klinear(1:maxval(this%ll(1:mdim)), 1:mdim))
1184 this%nkpnts = product(this%ll(1:mdim))
1193 this%parallel_in_momentum = .false.
1196 select case (this%kgrid)
1200 if (use_enegy_grid)
then
1202 this%klinear(ie, 1) =
sqrt(
m_two * (ie * de + emin))
1206 this%klinear(ikk, 1) = ikk * this%dk + kmin(1)
1220 if ((this%nk_end - this%nk_start + 1) < this%nk) this%parallel_in_momentum = .
true.
1221 call pes_flux_distribute(1, this%nstepsomegak, this%nstepsomegak_start, this%nstepsomegak_end, comm)
1233 safe_allocate(this%j_l(0:this%lmax, this%nk_start:this%nk_end))
1236 safe_allocate(this%kcoords_sph(1:3, this%nk_start:this%nk_end, 1:this%nstepsomegak))
1237 this%kcoords_sph =
m_zero
1239 safe_allocate(this%ylm_k(0:this%lmax, - this%lmax:this%lmax, this%nstepsomegak_start:this%nstepsomegak_end))
1244 do ith = 0, this%nstepsthetak
1245 thetak = ith * dthetak + this%thetak_rng(1)
1246 do iph = 0, this%nstepsphik - 1
1247 phik = iph * dphik + this%phik_rng(1)
1249 if (iomk >= this%nstepsomegak_start .and. iomk <= this%nstepsomegak_end)
then
1250 do ll = 0, this%lmax
1252 xx(1) =
cos(phik) *
sin(thetak)
1253 xx(2) =
sin(phik) *
sin(thetak)
1255 call ylmr_cmplx(xx, ll, mm, this%ylm_k(ll, mm, iomk))
1259 if (this%nk_start > 0)
then
1260 this%kcoords_sph(1, this%nk_start:this%nk_end, iomk) =
cos(phik) *
sin(thetak)
1261 this%kcoords_sph(2, this%nk_start:this%nk_end, iomk) =
sin(phik) *
sin(thetak)
1262 this%kcoords_sph(3, this%nk_start:this%nk_end, iomk) =
cos(thetak)
1268 if (this%nk_start > 0)
then
1270 do ikk = this%nk_start, this%nk_end
1271 kact = this%klinear(ikk,1)
1272 do ll = 0, this%lmax
1276 this%kcoords_sph(:, ikk, :) = kact * this%kcoords_sph(:, ikk, :)
1284 this%nkpnts_start = 1
1285 this%nkpnts_end = this%nkpnts
1288 safe_allocate(this%kcoords_cub(1:mdim, this%nkpnts_start:this%nkpnts_end, 1))
1290 this%kcoords_cub =
m_zero
1295 do ikk = -this%nk, this%nk
1298 kact = sign(1,ikk) * this%klinear(abs(ikk), 1)
1299 this%kcoords_cub(1, ikp, 1) = kact
1307 do ith = 0, this%nstepsthetak
1308 if (mdim == 3) thetak = ith * dthetak + this%thetak_rng(1)
1309 do iph = 0, this%nstepsphik - 1
1311 phik = iph * dphik + this%phik_rng(1)
1312 kact = this%klinear(ikk,1)
1313 this%kcoords_cub(1, ikp,1) = kact *
cos(phik) *
sin(thetak)
1314 this%kcoords_cub(2, ikp,1) = kact *
sin(phik) *
sin(thetak)
1315 if (mdim == 3) this%kcoords_cub(3, ikp,1) = kact *
cos(thetak)
1326 this%nkpnts_start = 1
1327 this%nkpnts_end = this%nkpnts
1329 if (kpoints%have_zero_weight_path())
then
1334 safe_allocate(this%kcoords_cub(1:mdim, this%nkpnts_start:this%nkpnts_end, kptst:kptend))
1335 this%kcoords_cub =
m_zero
1344 do ikpt = (kptst-1)+ispin, kptend, nspin
1346 do ikk = nkmin, nkmax
1347 kvec(1:mdim) = - kpoints%get_point(st%d%get_kpoint_index(ikpt))
1355 safe_allocate(this%kcoords_cub(1:mdim, this%nkpnts_start:this%nkpnts_end, 1))
1357 safe_allocate(this%Lkpuvz_inv(1:this%ll(1), 1:this%ll(2), 1:this%ll(3)))
1360 do ikk = 1, this%ll(idir)
1361 this%klinear(ikk, idir) = ikk * this%dk + kmin(idir)
1366 if (.not. this%arpes_grid)
then
1373 do ik3 = 1, this%ll(3)
1374 if (mdim>2) kvec(3) = this%klinear(ik3, 3)
1375 do ik2 = 1, this%ll(2)
1376 if (mdim>1) kvec(2) = this%klinear(ik2, 2)
1377 do ik1 = 1, this%ll(1)
1379 kvec(1) = this%klinear(ik1, 1)
1380 this%kcoords_cub(1:mdim, ikp, ikpt) = kvec(1:mdim)
1381 this%Lkpuvz_inv(ik1, ik2, ik3) = ikp
1393 do ikk = nkmin, nkmax
1399 do ik1 = 1, this%ll(1)
1400 kvec(1) = this%klinear(ik1,1)
1401 kvec(1:pdim) = matmul(kpoints%latt%klattice_primitive(1:pdim, 1:pdim), kvec(1:pdim))
1407 do ik2 = 1, this%ll(2)
1408 do ik1 = 1, this%ll(1)
1409 kvec(1:2) = (/this%klinear(ik1,1), this%klinear(ik2,2)/)
1410 kvec(1:pdim) = matmul(kpoints%latt%klattice_primitive(1:pdim, 1:pdim), kvec(1:pdim))
1413 this%Lkpuvz_inv(ik1, ik2, ikk - nkmin + 1) = ikp
1466 if (this%arpes_grid)
then
1475 safe_deallocate_a(gpoints)
1476 safe_deallocate_a(gpoints_reduced)
1497 this%ktransf(:,:) =
m_zero
1499 this%ktransf(idim,idim) =
m_one
1502 if (
parse_block(namespace,
'PES_Flux_GridTransformMatrix', blk) == 0)
then
1510 write(
message(1),
'(a)')
'Momentum grid transformation matrix :'
1511 do idir = 1, space%dim
1512 write(
message(1 + idir),
'(9f12.6)') ( this%ktransf(idim, idir), idim = 1, mdim)
1521 do ith = 0, this%nstepsthetak
1522 do iph = 0, this%nstepsphik - 1
1524 do ikk = this%nk_start, this%nk_end
1525 kvec(1:mdim) = this%kcoords_sph(1:mdim, ikk, iomk)
1526 this%kcoords_sph(1:mdim, ikk, iomk) = matmul(this%ktransf(1:mdim, 1:mdim), kvec(1:mdim))
1528 if (ith == 0 .or. ith == this%nstepsthetak)
exit
1534 do ikpt = kptst, kptend + 1
1535 if (ikpt == kptend + 1)
then
1536 kpoint(1:space%dim) =
m_zero
1538 kpoint(1:space%dim) = kpoints%get_point(st%d%get_kpoint_index(ikpt))
1541 do ikp = 1, this%nkpnts
1543 kvec(1:mdim) = this%kcoords_cub(1:mdim, ikp, ikpt) - kpoint(1:mdim)
1544 this%kcoords_cub(1:mdim, ikp, ikpt) = matmul(this%ktransf(1:mdim, 1:mdim), kvec(1:mdim)) &
1566 real(real64) :: kpar(1:pdim), val
1571 if (ikk /= 0) sign= ikk / abs(ikk)
1573 kpar(1:pdim) = kvec(1:pdim)
1574 val = abs(ikk) * de *
m_two - sum(kpar(1:pdim)**2)
1576 kvec(mdim) = sign *
sqrt(val)
1579 kvec(mdim) =
sqrt(val)
1580 nkp_out = nkp_out + 1
1583 this%kcoords_cub(1:mdim, ikp, ikpt) = kvec(1:mdim)
1593 subroutine pes_flux_calc(this, space, mesh, st, der, ext_partners, kpoints, iter, dt, stopping)
1595 class(
space_t),
intent(in) :: space
1596 type(
mesh_t),
intent(in) :: mesh
1601 integer,
intent(in) :: iter
1602 real(real64),
intent(in) :: dt
1603 logical,
intent(in) :: stopping
1605 integer :: stst, stend, kptst, kptend, sdim, mdim
1606 integer :: ist, ik, isdim, imdim
1608 integer :: il, ip_local
1609 complex(real64),
allocatable :: gpsi(:,:), psi(:)
1610 complex(real64),
allocatable :: interp_values(:)
1611 logical :: contains_isp
1612 real(real64) :: kpoint(1:3)
1613 complex(real64),
allocatable :: tmp_array(:)
1623 kptst = st%d%kpt%start
1624 kptend = st%d%kpt%end
1628 safe_allocate(psi(1:mesh%np_part))
1629 safe_allocate(gpsi(1:mesh%np_part, 1:mdim))
1631 if (this%surf_interp)
then
1632 safe_allocate(interp_values(1:this%nsrfcpnts))
1635 this%itstep = this%itstep + 1
1639 if(
associated(lasers))
then
1640 do il = 1, lasers%no_lasers
1641 call laser_field(lasers%lasers(il), this%veca(1:mdim, this%itstep), iter*dt)
1644 this%veca(:, this%itstep) = - this%veca(:, this%itstep)
1649 if(mesh%parallel_in_domains)
then
1650 safe_allocate(tmp_array(1:this%nsrfcpnts))
1654 do ik = kptst, kptend
1655 do ist = stst, stend
1664 kpoint(1:mdim) = kpoints%get_point(st%d%get_kpoint_index(ik))
1667 do ip = 1, mesh%np_part
1668 psi(ip) =
exp(-
m_zi * sum(mesh%x(ip, 1:mdim) * kpoint(1:mdim))) * psi(ip)
1675 if (this%surf_interp)
then
1677 this%rcoords(1:mdim, 1:this%nsrfcpnts), interp_values(1:this%nsrfcpnts))
1678 this%wf(ist, isdim, ik, 1:this%nsrfcpnts, this%itstep) = st%occ(ist, ik) * interp_values(1:this%nsrfcpnts)
1681 this%rcoords(1:mdim, 1:this%nsrfcpnts), interp_values(1:this%nsrfcpnts))
1682 this%gwf(ist, isdim, ik, 1:this%nsrfcpnts, this%itstep, imdim) = st%occ(ist, ik) * interp_values(1:this%nsrfcpnts)
1687 contains_isp = .
true.
1688 do isp = 1, this%nsrfcpnts
1689 if (mesh%parallel_in_domains)
then
1690 if (mesh%mpi_grp%rank == this%rankmin(isp))
then
1691 contains_isp = .
true.
1693 contains_isp = .false.
1697 if (contains_isp)
then
1698 ip_local = this%srfcpnt(isp)
1699 this%wf(ist, isdim, ik, isp, this%itstep) = st%occ(ist, ik) * psi(ip_local)
1700 this%gwf(ist, isdim, ik, isp, this%itstep, 1:mdim) = &
1701 st%occ(ist, ik) * gpsi(ip_local, 1:mdim)
1704 if (mesh%parallel_in_domains)
then
1706 tmp_array = this%wf(ist, isdim, ik, 1:this%nsrfcpnts, this%itstep)
1707 call mesh%allreduce(tmp_array)
1708 this%wf(ist, isdim, ik, 1:this%nsrfcpnts, this%itstep) = tmp_array
1711 tmp_array = this%gwf(ist, isdim, ik, 1:this%nsrfcpnts, this%itstep, imdim)
1712 call mesh%allreduce(tmp_array)
1713 this%gwf(ist, isdim, ik, 1:this%nsrfcpnts, this%itstep, imdim) = tmp_array
1721 if(mesh%parallel_in_domains)
then
1722 safe_deallocate_a(tmp_array)
1725 safe_deallocate_a(psi)
1726 safe_deallocate_a(gpsi)
1728 if (this%itstep == this%tdsteps .or. mod(iter, this%save_iter) == 0 .or. iter == this%max_iter .or. stopping)
then
1729 if (this%surf_shape == pes_cubic .or. this%surf_shape ==
pes_plane)
then
1752 class(
space_t),
intent(in) :: space
1753 type(
mesh_t),
intent(in) :: mesh
1757 integer :: kptst, kptend, mdim
1758 integer :: ik, j1, j2, jvec(1:2)
1759 integer :: isp, ikp, ikp_start, ikp_end
1762 integer :: idir, pdim, nfaces, ifc, n_dir
1763 complex(real64) :: tmp
1764 real(real64) :: kpoint(1:3), vec(1:3), lvec(1:3)
1768 if (kpoints%have_zero_weight_path())
then
1769 kptst = st%d%kpt%start
1770 kptend = st%d%kpt%end
1777 pdim = space%periodic_dim
1779 ikp_start = this%nkpnts_start
1780 ikp_end = this%nkpnts_end
1786 this%srfcnrml(:, 1:this%nsrfcpnts) = this%srfcnrml(:, 1:this%nsrfcpnts) * mesh%coord_system%surface_element(3)
1789 if (.not. this%use_symmetry .or. kpoints%have_zero_weight_path())
then
1791 safe_allocate(this%expkr(1:this%nsrfcpnts,ikp_start:ikp_end,kptst:kptend,1))
1793 do ik = kptst, kptend
1794 do ikp = ikp_start, ikp_end
1795 do isp = 1, this%nsrfcpnts
1796 this%expkr(isp,ikp,ik,1) =
exp(-
m_zi * sum(this%rcoords(1:mdim,isp) &
1797 * this%kcoords_cub(1:mdim, ikp, ik))) &
1809 if (this%surf_shape ==
pes_plane) nfaces = 1
1812 safe_allocate(this%expkr(1:this%nsrfcpnts,maxval(this%ll(1:mdim)),kptst:kptend,1:mdim))
1813 this%expkr(:,:,:,:) =
m_z1
1815 do ik = kptst, kptend
1817 do ikp = 1, this%ll(idir)
1818 do isp = 1, this%nsrfcpnts
1819 this%expkr(isp,ikp,ik,idir) =
exp(-
m_zi * this%rcoords(idir,isp) &
1820 * this%klinear(ikp, idir)) &
1828 safe_allocate(this%expkr_perp(1:maxval(this%ll(1:mdim)), 1:nfaces))
1829 this%expkr_perp(:,:) =
m_z1
1832 if (this%face_idx_range(ifc, 1) < 0) cycle
1833 isp = this%face_idx_range(ifc, 1)
1835 if (abs(this%srfcnrml(idir, isp)) >=
m_epsilon) n_dir = idir
1838 do ikp = 1, this%ll(n_dir)
1839 this%expkr_perp(ikp, ifc) =
exp(-
m_zi * this%rcoords(n_dir,isp) &
1840 * this%klinear(ikp, n_dir)) &
1850 if (space%is_periodic())
then
1852 safe_allocate(this%bvk_phase(ikp_start:ikp_end,st%d%kpt%start:st%d%kpt%end))
1854 this%bvk_phase(:,:) =
m_z0
1857 do ik = st%d%kpt%start, st%d%kpt%end
1858 kpoint(1:mdim) = kpoints%get_point(st%d%get_kpoint_index(ik))
1859 if (kpoints%have_zero_weight_path())
then
1864 do ikp = ikp_start, ikp_end
1865 vec(1:pdim) = this%kcoords_cub(1:pdim, ikp, ik_map) + kpoint(1:pdim)
1866 do j1 = 0, kpoints%nik_axis(1) - 1
1867 do j2 = 0, kpoints%nik_axis(2) - 1
1868 jvec(1:2)=(/j1, j2/)
1869 lvec(1:pdim) = matmul(kpoints%latt%rlattice(1:pdim, 1:2), jvec(1:2))
1870 tmp = sum(lvec(1:pdim) * vec(1:pdim))
1871 this%bvk_phase(ikp, ik) = this%bvk_phase(ikp,ik) &
1879 this%bvk_phase(:,:) = this%bvk_phase(:,:) *
m_one / product(kpoints%nik_axis(1:pdim))
1901 pure function get_ikp(this, ikpu, ikpv, ikpz, n_dir)
result(ikp)
1903 integer,
intent(in) :: ikpu
1904 integer,
intent(in) :: ikpv
1905 integer,
intent(in) :: ikpz
1906 integer,
intent(in) :: n_dir
1911 ikp = this%Lkpuvz_inv(ikpz, ikpu, ikpv)
1913 ikp = this%Lkpuvz_inv(ikpu, ikpz, ikpv)
1915 ikp = this%Lkpuvz_inv(ikpu, ikpv, ikpz)
1927 type(
mesh_t),
intent(in) :: mesh
1930 integer :: mdim, nfaces, ifc, ifp_start, ifp_end
1937 if (this%surf_shape ==
pes_plane) nfaces = 1
1941 ifp_start = this%face_idx_range(ifc, 1)
1942 ifp_end = this%face_idx_range(ifc, 2)
1945 if (this%nsrfcpnts_start <= ifp_end)
then
1947 if (this%nsrfcpnts_start >= ifp_start)
then
1948 if (this%nsrfcpnts_start <= ifp_end)
then
1949 this%face_idx_range(ifc, 1) = this%nsrfcpnts_start
1951 this%face_idx_range(ifc, 1:2) = -1
1955 if (this%nsrfcpnts_end <= ifp_end)
then
1956 if (this%nsrfcpnts_end >= ifp_start)
then
1957 this%face_idx_range(ifc, 2) = this%nsrfcpnts_end
1959 this%face_idx_range(ifc, 1:2) = -1
1965 this%face_idx_range(ifc, 1:2) = -1
1984 class(
space_t),
intent(in) :: space
1985 type(
mesh_t),
intent(in) :: mesh
1988 real(real64),
intent(in) :: dt
1990 integer :: stst, stend, kptst, kptend, sdim, mdim
1991 integer :: ist, ik, isdim, imdim, ik_map
1992 integer :: ikp, itstep
1993 integer :: idir, n_dir, nfaces
1994 complex(real64),
allocatable :: jk_cub(:,:,:,:), spctramp_cub(:,:,:,:)
1995 integer :: ikp_start, ikp_end, isp_start, isp_end
1996 real(real64) :: vec, kpoint(1:3)
1998 complex(real64),
allocatable :: wfpw(:), gwfpw(:)
1999 complex(real64),
allocatable :: phase(:,:),vphase(:,:)
2001 integer :: tdstep_on_node
2005 integer :: ikpu, ikpv, ikpz, dir_on_face(1:2)
2006 complex(real64) :: face_int_gwf, face_int_wf
2015 call messages_write(
"Debug: calculating pes_flux cub surface integral (accelerated direct expression)")
2020 if (
bitand(this%par_strategy, option__pes_flux_parallelization__pf_time) /= 0)
then
2021 if (mesh%parallel_in_domains) tdstep_on_node = mesh%mpi_grp%rank + 1
2026 kptst = st%d%kpt%start
2027 kptend = st%d%kpt%end
2031 ikp_start = this%nkpnts_start
2032 ikp_end = this%nkpnts_end
2036 safe_allocate(jk_cub(stst:stend, 1:sdim, kptst:kptend, ikp_start:ikp_end))
2037 safe_allocate(spctramp_cub(stst:stend, 1:sdim, kptst:kptend, ikp_start:ikp_end))
2041 safe_allocate(vphase(ikp_start:ikp_end, kptst:kptend))
2042 safe_allocate(phase(ikp_start:ikp_end, kptst:kptend))
2048 if (this%surf_shape ==
pes_plane) nfaces = 1
2050 safe_allocate( wfpw(ikp_start:ikp_end))
2051 safe_allocate(gwfpw(ikp_start:ikp_end))
2056 if (this%face_idx_range(ifc, 1)<0) cycle
2059 isp_start = this%face_idx_range(ifc, 1)
2060 isp_end = this%face_idx_range(ifc, 2)
2063 nfp = isp_end - isp_start + 1
2071 if (abs(this%srfcnrml(idir, isp_start)) >=
m_epsilon)
then
2074 dir_on_face(imdim)=idir
2081 vphase(ikp_start:ikp_end,:) = this%conjgphase_prev(ikp_start:ikp_end,:)
2083 jk_cub(:, :, :, :) =
m_z0
2085 do itstep = 1, this%itstep
2087 do ik = kptst, kptend
2089 if (kpoints%have_zero_weight_path())
then
2095 kpoint(1:mdim) = kpoints%get_point(st%d%get_kpoint_index(ik))
2096 do ikp = ikp_start, ikp_end
2097 vec = sum((this%kcoords_cub(1:mdim, ikp, ik_map) - this%veca(1:mdim, itstep) /
p_c)**2)
2098 vphase(ikp, ik) = vphase(ikp, ik) *
exp(
m_zi * vec * dt /
m_two)
2101 if (space%is_periodic())
then
2102 phase(ikp, ik) = vphase(ikp, ik) * this%bvk_phase(ikp,ik)
2104 phase(ikp, ik) = vphase(ikp, ik)
2109 if (itstep /= tdstep_on_node) cycle
2111 do ist = stst, stend
2116 if (.not. this%use_symmetry .or. kpoints%have_zero_weight_path())
then
2118 do ikp = ikp_start , ikp_end
2121 sum(this%gwf(ist, isdim, ik, isp_start:isp_end, itstep, n_dir) &
2122 * this%expkr(isp_start:isp_end, ikp, ik_map, 1) &
2123 * this%srfcnrml(n_dir, isp_start:isp_end))
2127 sum(this%wf(ist, isdim, ik, isp_start:isp_end, itstep) &
2128 * this%expkr(isp_start:isp_end, ikp,ik_map, 1) &
2129 * this%srfcnrml(n_dir, isp_start:isp_end))
2134 do ikpu = 1, this%ll(dir_on_face(1))
2135 do ikpv = 1, this%ll(dir_on_face(2))
2138 face_int_gwf = sum(this%gwf(ist, isdim, ik, isp_start:isp_end, itstep, n_dir) &
2139 * this%expkr(isp_start:isp_end, ikpu, ik_map, dir_on_face(1)) &
2140 * this%expkr(isp_start:isp_end, ikpv, ik_map, dir_on_face(2)) &
2141 * this%srfcnrml(n_dir, isp_start:isp_end))
2143 face_int_wf = sum(this%wf(ist, isdim, ik, isp_start:isp_end, itstep) &
2144 * this%expkr(isp_start:isp_end, ikpu, ik_map, dir_on_face(1)) &
2145 * this%expkr(isp_start:isp_end, ikpv, ik_map, dir_on_face(2)) &
2146 * this%srfcnrml(n_dir, isp_start:isp_end))
2148 do ikpz = 1, this%ll(n_dir)
2150 gwfpw(
get_ikp(this, ikpu, ikpv, ikpz, n_dir)) = face_int_gwf &
2151 * this%expkr_perp(ikpz, ifc)
2152 wfpw(
get_ikp(this, ikpu, ikpv, ikpz, n_dir)) = face_int_wf &
2153 * this%expkr_perp(ikpz, ifc)
2165 jk_cub(ist, isdim, ik, ikp_start:ikp_end) = jk_cub(ist, isdim, ik, ikp_start:ikp_end) + &
2166 phase(ikp_start:ikp_end, ik) * (wfpw(ikp_start:ikp_end) * &
2167 (
m_two * this%veca(n_dir, itstep) /
p_c - this%kcoords_cub(n_dir, ikp_start:ikp_end, ik_map)) + &
2168 m_zi * gwfpw(ikp_start:ikp_end))
2178 spctramp_cub(:,:,:,:) = spctramp_cub(:,:,:,:) + jk_cub(:,:,:,:) *
m_half
2184 this%conjgphase_prev(ikp_start:ikp_end, :) = vphase(ikp_start:ikp_end, :)
2187 if (mesh%parallel_in_domains .and.(
bitand(this%par_strategy, option__pes_flux_parallelization__pf_time) /= 0 &
2188 .or.
bitand(this%par_strategy, option__pes_flux_parallelization__pf_surface) /= 0))
then
2189 call mesh%allreduce(spctramp_cub)
2194 this%spctramp_cub = this%spctramp_cub + spctramp_cub * dt
2196 safe_deallocate_a(gwfpw)
2197 safe_deallocate_a(wfpw)
2199 safe_deallocate_a(jk_cub)
2200 safe_deallocate_a(spctramp_cub)
2201 safe_deallocate_a(phase)
2202 safe_deallocate_a(vphase)
2214 type(
mesh_t),
intent(in) :: mesh
2216 real(real64),
intent(in) :: dt
2218 integer :: stst, stend, kptst, kptend, sdim
2219 integer :: ist, ik, isdim, imdim
2220 integer :: itstep, lmax
2222 complex(real64),
allocatable :: s1_node(:,:,:,:,:,:), s1_act(:,:,:)
2223 complex(real64),
allocatable :: s2_node(:,:,:,:,:), s2_act(:,:)
2224 complex(real64),
allocatable :: integ11_t(:,:), integ12_t(:,:,:)
2225 complex(real64),
allocatable :: integ21_t(:), integ22_t(:,:)
2226 complex(real64),
allocatable :: spctramp_sph(:,:,:,:,:)
2227 integer :: ikk, ikk_start, ikk_end
2228 integer :: isp_start, isp_end
2229 integer :: iomk, iomk_start, iomk_end
2230 complex(real64),
allocatable :: phase_act(:,:)
2232 integer :: tdstep_on_node
2238 call messages_write(
"Debug: calculating pes_flux sph surface integral")
2242 if (mesh%parallel_in_domains) tdstep_on_node = mesh%mpi_grp%rank + 1
2246 kptst = st%d%kpt%start
2247 kptend = st%d%kpt%end
2251 ikk_start = this%nk_start
2252 ikk_end = this%nk_end
2253 isp_start = this%nsrfcpnts_start
2254 isp_end = this%nsrfcpnts_end
2255 iomk_start = this%nstepsomegak_start
2256 iomk_end = this%nstepsomegak_end
2259 safe_allocate(s1_node(stst:stend, 1:sdim, kptst:kptend, 0:lmax, -lmax:lmax, 1:3))
2260 safe_allocate(s2_node(stst:stend, 1:sdim, kptst:kptend, 0:lmax, -lmax:lmax))
2262 safe_allocate(s1_act(0:lmax, -lmax:lmax, 1:3))
2263 safe_allocate(s2_act(0:lmax, -lmax:lmax))
2265 do itstep = 1, this%itstep
2267 do ik = kptst, kptend
2268 do ist = stst, stend
2277 s1_act(ll, mm, imdim) = &
2278 sum(this%ylm_r(ll, mm, isp_start:isp_end) * this%srfcnrml(imdim, isp_start:isp_end) * &
2279 this%wf(ist, isdim, ik, isp_start:isp_end, itstep))
2281 s2_act(ll, mm) = s2_act(ll, mm) + &
2282 sum(this%ylm_r(ll, mm, isp_start:isp_end) * this%srfcnrml(imdim, isp_start:isp_end) * &
2283 this%gwf(ist, isdim, ik, isp_start:isp_end, itstep, imdim))
2288 if (mesh%parallel_in_domains .and.
bitand(this%par_strategy, option__pes_flux_parallelization__pf_surface) /= 0)
then
2289 call mesh%allreduce(s1_act)
2290 call mesh%allreduce(s2_act)
2293 if (itstep == tdstep_on_node)
then
2294 s1_node(ist, isdim, ik, :, :, :) = s1_act(:, :, :)
2295 s2_node(ist, isdim, ik, :, :) = s2_act(:, :)
2304 safe_allocate(integ11_t(0:this%nstepsomegak, 1:3))
2305 safe_allocate(integ21_t(0:this%nstepsomegak))
2306 safe_allocate(integ12_t(ikk_start:ikk_end, 1:this%nstepsomegak, 1:3))
2307 safe_allocate(integ22_t(ikk_start:ikk_end, 1:this%nstepsomegak))
2309 safe_allocate(spctramp_sph(stst:stend, 1:sdim, kptst:kptend, 0:this%nk, 1:this%nstepsomegak))
2313 safe_allocate(phase_act(ikk_start:ikk_end, 1:this%nstepsomegak))
2314 if (ikk_start > 0)
then
2315 phase_act(ikk_start:ikk_end, :) = this%conjgphase_prev(ikk_start:ikk_end, :)
2317 phase_act(ikk_start:ikk_end, :) =
m_z0
2320 do itstep = 1, this%itstep
2322 do ikk = ikk_start, ikk_end
2323 do iomk = 1, this%nstepsomegak
2324 vec = sum((this%kcoords_sph(1:3, ikk, iomk) - this%veca(1:3, itstep) /
p_c)**
m_two)
2325 phase_act(ikk, iomk) = phase_act(ikk, iomk) *
exp(
m_zi * vec * dt /
m_two)
2329 do ik = kptst, kptend
2330 do ist = stst, stend
2334 if (itstep == tdstep_on_node)
then
2335 s1_act(:,:,:) = s1_node(ist, isdim, ik, :, :, :)
2336 s2_act(:,:) = s2_node(ist, isdim, ik, :, :)
2338 if (mesh%parallel_in_domains)
then
2339 call mesh%mpi_grp%bcast(s1_act, (lmax + 1) * (2 * lmax + 1) * 3, mpi_double_complex, itstep - 1)
2340 call mesh%mpi_grp%bcast(s2_act, (lmax + 1) * (2 * lmax + 1), mpi_double_complex, itstep - 1)
2354 integ11_t(iomk_start:iomk_end, imdim) = integ11_t(iomk_start:iomk_end, imdim) + &
2355 s1_act(ll, mm, imdim) * this%ylm_k(ll, mm, iomk_start:iomk_end)
2357 integ21_t(iomk_start:iomk_end) = integ21_t(iomk_start:iomk_end) + &
2358 s2_act(ll, mm) * this%ylm_k(ll, mm, iomk_start:iomk_end)
2361 call mesh%allreduce(integ11_t)
2362 call mesh%allreduce(integ21_t)
2365 do ikk = ikk_start, ikk_end
2366 integ12_t(ikk, 1:this%nstepsomegak, 1:3) = integ12_t(ikk, 1:this%nstepsomegak, 1:3) + &
2367 integ11_t(1:this%nstepsomegak, 1:3) * this%j_l(ll, ikk) * ( -
m_zi)**ll
2368 integ22_t(ikk, 1:this%nstepsomegak) = integ22_t(ikk, 1:this%nstepsomegak) + &
2369 integ21_t(1:this%nstepsomegak) * this%j_l(ll, ikk) * ( -
m_zi)**ll
2374 spctramp_sph(ist, isdim, ik, ikk_start:ikk_end, :) = &
2375 spctramp_sph(ist, isdim, ik, ikk_start:ikk_end, :) + &
2376 phase_act(ikk_start:ikk_end,:) * (integ12_t(ikk_start:ikk_end,:, imdim) * &
2377 (
m_two * this%veca(imdim, itstep) /
p_c - this%kcoords_sph(imdim, ikk_start:ikk_end,:)))
2379 spctramp_sph(ist, isdim, ik, ikk_start:ikk_end, :) = &
2380 spctramp_sph(ist, isdim, ik, ikk_start:ikk_end,:) + &
2381 phase_act(ikk_start:ikk_end,:) * integ22_t(ikk_start:ikk_end,:) *
m_zi
2387 safe_deallocate_a(s1_node)
2388 safe_deallocate_a(s2_node)
2389 safe_deallocate_a(s1_act)
2390 safe_deallocate_a(s2_act)
2392 safe_deallocate_a(integ11_t)
2393 safe_deallocate_a(integ12_t)
2394 safe_deallocate_a(integ21_t)
2395 safe_deallocate_a(integ22_t)
2398 this%conjgphase_prev =
m_z0
2400 if (ikk_start > 0)
then
2401 this%conjgphase_prev(ikk_start:ikk_end, :) = phase_act(ikk_start:ikk_end, :)
2403 safe_deallocate_a(phase_act)
2405 if (mesh%parallel_in_domains .and.
bitand(this%par_strategy, option__pes_flux_parallelization__pf_time) /= 0)
then
2406 call mesh%allreduce(this%conjgphase_prev)
2407 call mesh%allreduce(spctramp_sph)
2410 this%spctramp_sph(:, :, :, 1:this%nk, :) = this%spctramp_sph(:, :, :, 1:this%nk, :) &
2411 + spctramp_sph(:, :, :, 1:this%nk, :) * dt
2412 safe_deallocate_a(spctramp_sph)
2419 type(
mesh_t),
intent(in) :: mesh
2422 real(real64),
intent(in) :: border(:)
2423 real(real64),
intent(in) :: offset(:)
2424 real(real64),
intent(in) :: fc_ptdens
2426 integer,
allocatable :: which_surface(:)
2427 real(real64) :: xx(mesh%box%dim), dd, area, dS(mesh%box%dim, 2), factor
2428 integer :: mdim, imdim, idir, isp, pm,nface, idim, ndir, iu,iv, iuv(1:2)
2429 integer(int64) :: ip_global
2431 integer :: rankmin, nsurfaces
2433 integer :: ip_local, NN(mesh%box%dim, 2), idx(mesh%box%dim, 2)
2434 integer :: isp_end, isp_start, ifc, n_dir, nfaces, mindim
2435 real(real64) :: RSmax(2), RSmin(2), RS(2), dRS(2), weight
2447 safe_allocate(this%face_idx_range(1:mdim * 2, 1:2))
2448 safe_allocate(this%LLr(mdim, 1:2))
2449 safe_allocate(this%NN(mdim, 1:2))
2451 this%face_idx_range(:, :) = 0
2456 if (this%surf_interp)
then
2472 do ndir = mdim, mindim, -1
2475 if (idir == ndir) cycle
2476 area = area * border(idir) * factor
2479 npface = int(fc_ptdens * area)
2483 if (idir == ndir) cycle
2484 nn(ndir, idim) = int(
sqrt(npface * (border(idir) * factor)**2 / area))
2485 ds(ndir, idim) = border(idir) * factor / nn(ndir, idim)
2486 this%LLr(ndir, idim) = nn(ndir, idim) * ds(ndir, idim)
2487 idx(ndir, idim) = idir
2490 this%nsrfcpnts = this%nsrfcpnts + 2 * product(nn(ndir, 1:mdim - 1))
2493 this%NN(1:mdim, :) = nn(1:mdim, :)
2496 assert(this%nsrfcpnts > 0)
2499 safe_allocate(this%srfcnrml(1:mdim, 0:this%nsrfcpnts))
2500 safe_allocate(this%rcoords(1:mdim, 0:this%nsrfcpnts))
2502 this%srfcnrml(:, :) =
m_zero
2503 this%rcoords(:, :) =
m_zero
2508 do ndir = mdim, mindim, -1
2511 this%face_idx_range(ifc,1) = isp
2512 do iu = 1, nn(ndir,1)
2513 do iv = 1, nn(ndir,2)
2514 this%rcoords(ndir, isp) = border(ndir)
2515 this%srfcnrml(ndir, isp) = product(ds(ndir,1:mdim-1))
2518 this%rcoords(idx(ndir, idim), isp) = (-nn(ndir, idim) *
m_half -
m_half + iuv(idim)) * ds(ndir, idim)
2523 this%face_idx_range(ifc, 2) = isp-1
2528 this%face_idx_range(ifc, 1) = isp
2529 do iu = 1, nn(ndir, 1)
2530 do iv = 1, nn(ndir, 2)
2531 this%rcoords(ndir, isp) = -border(ndir)
2532 this%srfcnrml(ndir, isp) = -product(ds(ndir, 1:mdim - 1))
2535 this%rcoords(idx(ndir, idim), isp) = (-nn(ndir, idim) *
m_half -
m_half + iuv(idim)) * ds(ndir, idim)
2540 this%face_idx_range(ifc, 2) = isp-1
2545 do isp = 1, this%nsrfcpnts
2546 this%rcoords(1:mdim, isp) = mesh%coord_system%to_cartesian(this%rcoords(1:mdim, isp))
2553 if (this%surf_shape ==
pes_plane) nfaces = 2
2557 safe_allocate(which_surface(1:mesh%np_global))
2562 do ip_local = 1, mesh%np
2567 xx(1:mdim) = mesh%x(ip_local, 1:mdim) - offset(1:mdim)
2570 select case (abb%abtype)
2574 in_ab = abs(abb%mf(ip_local)) >
m_epsilon
2580 if (all(abs(xx(1:mdim)) <= border(1:mdim)) .and. .not. in_ab)
then
2584 dd = border(imdim) - xx(imdim)
2586 dd = border(imdim) - abs(xx(imdim))
2588 if (dd < mesh%spacing(imdim) /
m_two)
then
2589 nsurfaces = nsurfaces + 1
2590 which_surface(ip_global) = int(sign(
m_one, xx(imdim))) * imdim
2595 if (nsurfaces == 1)
then
2596 this%nsrfcpnts = this%nsrfcpnts + 1
2598 which_surface(ip_global) = 0
2604 if (mesh%parallel_in_domains)
then
2605 assert(mesh%np_global < huge(0_int32))
2606 call mesh%allreduce(this%nsrfcpnts)
2607 call mesh%allreduce(which_surface)
2610 safe_allocate(this%srfcpnt(1:this%nsrfcpnts))
2611 safe_allocate(this%srfcnrml(1:mdim, 0:this%nsrfcpnts))
2612 safe_allocate(this%rcoords(1:mdim, 0:this%nsrfcpnts))
2613 safe_allocate(this%rankmin(1:this%nsrfcpnts))
2618 this%face_idx_range(:,:) = 0
2622 do idir = mdim, 1, -1
2625 this%face_idx_range(nface, 1) = isp + 1
2627 do ip_global = 1, mesh%np_global
2628 if (abs(which_surface(ip_global)) == idir .and. sign(1, which_surface(ip_global)) == pm)
then
2632 this%rcoords(1:mdim, isp) = xx(1:mdim)
2635 this%rankmin(isp) = rankmin
2637 this%srfcnrml(idir, isp) = sign(1, which_surface(ip_global))
2640 if (imdim == idir) cycle
2641 this%srfcnrml(idir, isp) = this%srfcnrml(idir, isp) * mesh%spacing(imdim)
2646 this%face_idx_range(nface,2) = isp
2653 do ifc = 1, nint((nfaces + 0.5) / 2)
2654 isp_start = this%face_idx_range(ifc, 1)
2655 isp_end = this%face_idx_range(ifc, 2)
2660 if (abs(this%srfcnrml(idir, isp_start)) >=
m_epsilon) n_dir = idir
2667 if (idir == n_dir) cycle
2668 drs(idim)= mesh%spacing(idir)
2676 do isp = isp_start, isp_end
2678 xx = mesh%coord_system%from_cartesian(this%rcoords(1:mdim, isp))
2681 if (idir == n_dir) cycle
2683 if (rs(idim) < rsmin(idim)) rsmin(idim) = rs(idim)
2684 if (rs(idim) > rsmax(idim)) rsmax(idim) = rs(idim)
2690 do idir = 1, mdim - 1
2691 this%LLr(n_dir, idir) = rsmax(idir) - rsmin(idir) + drs(idir)
2692 if (drs(idir) >
m_zero) this%NN(n_dir, idir) = int(this%LLr(n_dir, idir) / drs(idir))
2700 if (this%anisotrpy_correction)
then
2703 do isp = 1, this%nsrfcpnts
2704 weight = sum(this%rcoords(1:mdim, isp) * this%srfcnrml(1:mdim, isp))
2705 weight = weight/sum(this%rcoords(1:mdim, isp)**2) / sum(this%srfcnrml(1:mdim, isp)**2)
2706 this%srfcnrml(1:mdim, isp) = this%srfcnrml(1:mdim, isp) * abs(weight)
2712 safe_deallocate_a(which_surface)
2720 type(
mesh_t),
intent(in) :: mesh
2721 integer,
intent(inout) :: nstepsthetar, nstepsphir
2724 integer :: isp, ith, iph
2725 real(real64) :: thetar, phir
2726 real(real64) :: weight, dthetar
2732 if (nstepsphir == 0) nstepsphir = 1
2734 if (nstepsthetar <= 1)
then
2738 dthetar =
m_pi / nstepsthetar
2739 this%nsrfcpnts = nstepsphir * (nstepsthetar - 1) + 2
2741 safe_allocate(this%srfcnrml(1:mdim, 0:this%nsrfcpnts))
2742 safe_allocate(this%rcoords(1:mdim, 0:this%nsrfcpnts))
2749 do ith = 0, nstepsthetar
2750 thetar = ith * dthetar
2752 if (ith == 0 .or. ith == nstepsthetar)
then
2755 weight = abs(
cos(thetar - dthetar /
m_two) -
cos(thetar + dthetar /
m_two)) &
2759 do iph = 0, nstepsphir - 1
2762 this%srfcnrml(1, isp) =
cos(phir) *
sin(thetar)
2763 if (mdim >= 2) this%srfcnrml(2, isp) =
sin(phir) *
sin(thetar)
2764 if (mdim == 3) this%srfcnrml(3, isp) =
cos(thetar)
2765 this%rcoords(1:mdim, isp) = this%radius * this%srfcnrml(1:mdim, isp)
2767 this%srfcnrml(1:mdim, isp) = this%radius**
m_two * weight * this%srfcnrml(1:mdim, isp)
2768 if (ith == 0 .or. ith == nstepsthetar)
exit
2780 integer,
intent(in) :: istart_global
2781 integer,
intent(in) :: iend_global
2782 integer,
intent(inout) :: istart
2783 integer,
intent(inout) :: iend
2784 type(mpi_comm),
intent(in) :: comm
2786#if defined(HAVE_MPI)
2787 integer,
allocatable :: dimrank(:)
2788 integer :: mpisize, mpirank, irest, irank
2794#if defined(HAVE_MPI)
2795 call mpi_comm_size(comm, mpisize,
mpi_err)
2796 call mpi_comm_rank(comm, mpirank,
mpi_err)
2798 safe_allocate(dimrank(0:mpisize-1))
2800 inumber = iend_global - istart_global + 1
2801 irest = mod(inumber, mpisize)
2802 do irank = 0, mpisize - 1
2804 dimrank(irank) = int(inumber / mpisize) + 1
2807 dimrank(irank) = int(inumber / mpisize)
2811 iend = istart_global + sum(dimrank(:mpirank)) - 1
2812 istart = iend - dimrank(mpirank) + 1
2814 if (istart > iend)
then
2819 safe_deallocate_a(dimrank)
2822 istart = istart_global
2830#include "pes_flux_out_inc.F90"
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
integer, parameter, public imaginary_absorbing
integer, parameter, public mask_absorbing
Module implementing boundary conditions in Octopus.
This module handles the calculation mode.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public zderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
integer, parameter, public spin_polarized
type(lasers_t) function, pointer, public list_get_lasers(partners)
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
complex(real64), parameter, public m_z0
complex(real64), parameter, public m_zi
real(real64), parameter, public m_epsilon
complex(real64), parameter, public m_z1
real(real64), parameter, public m_half
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
real(real64), parameter, public m_one
real(real64), parameter, public m_three
This module defines classes and functions for interaction partners.
integer pure function, public kpoints_number(this)
integer, parameter, public e_field_vector_potential
integer pure elemental function, public laser_kind(laser)
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...
This module is intended to contain "only mathematical" functions and procedures.
subroutine, public ylmr_cmplx(xx, li, mi, ylm)
Computes spherical harmonics ylm at position (x, y, z)
subroutine, public mesh_interpolation_init(this, mesh)
This module defines the meshes, which are used in Octopus.
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
integer(int64) function, public mesh_local2global(mesh, ip)
This function returns the global mesh index for a given local index.
real(real64) function, dimension(1:mesh%box%dim), public mesh_x_global(mesh, ipg)
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_new_line()
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
integer, public mpi_err
used to store return values of mpi calls
This module handles the communicators for the various parallelization strategies.
logical function, public parse_is_defined(namespace, name)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public pes_flux_map_from_states(this, restart, st, ll, pesP, krng, Lp, istin)
subroutine pes_flux_getcube(this, mesh, abb, border, offset, fc_ptdens)
subroutine, public pes_flux_calc(this, space, mesh, st, der, ext_partners, kpoints, iter, dt, stopping)
subroutine, public pes_flux_out_energy(this, pesK, file, namespace, ll, Ekin, dim)
subroutine pes_flux_getsphere(this, mesh, nstepsthetar, nstepsphir)
subroutine, public pes_flux_output(this, box, st, namespace)
subroutine, public pes_flux_pmesh(this, namespace, dim, kpoints, ll, pmesh, idxZero, krng, nspin, Lp, Ekin)
subroutine, public pes_flux_reciprocal_mesh_gen(this, namespace, space, st, kpoints, comm, post)
subroutine, public pes_flux_load(restart, this, st, ierr)
subroutine, public pes_flux_dump(restart, this, mesh, st, ierr)
subroutine pes_flux_distribute_facepnts_cub(this, mesh)
subroutine pes_flux_integrate_sph(this, mesh, st, dt)
integer, parameter pes_cartesian
integer, parameter pes_plane
subroutine pes_flux_distribute(istart_global, iend_global, istart, iend, comm)
subroutine pes_flux_integrate_cub(this, space, mesh, st, kpoints, dt)
subroutine, public pes_flux_end(this)
pure integer function get_ikp(this, ikpu, ikpv, ikpz, n_dir)
subroutine pes_flux_integrate_cub_tabulate(this, space, mesh, st, kpoints)
subroutine, public pes_flux_init(this, namespace, space, mesh, st, ext_partners, kpoints, abb, save_iter, max_iter)
integer, parameter pes_spherical
subroutine, public pes_flux_out_vmap(this, pesK, file, namespace, ll, pmesh, dim)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
This module handles spin dimensions of the states and the k-point distribution.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
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
This module is intended to contain simple general-purpose utility functions and procedures.
character pure function, public index2axis(idir)
subroutine fill_non_periodic_dimension(this)
Class implementing a parallelepiped box. Currently this is restricted to a rectangular cuboid (all th...
Class implementing a spherical box.
class representing derivatives
Describes mesh distribution to nodes.
The states_elec_t class contains all electronic wave functions.