Octopus
xc_photons.F90
Go to the documentation of this file.
1!! Copyright (C) 2022 I.-T Lu
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
27
29 use debug_oct_m
31 use epot_oct_m
32 use global_oct_m
33 use grid_oct_m
34 use iso_c_binding
38 use mesh_oct_m
45 use space_oct_m
47 use parser_oct_m
50
51 implicit none
52
53 private
54 public :: xc_photons_t
55
64 !
65 type xc_photons_t
66 private
67 integer :: method = 0
68 real(real64), allocatable, public :: vpx(:)
69 real(real64), public :: ex
70 type(photon_mode_t) :: pt
71 real(real64) :: pxlda_kappa
72 real(real64) :: eta_c
74 integer :: energy_method = 0
75 logical :: lcorrelations = .false.
76 logical :: llamb_re_mass = .false.
77 logical :: llamb_freespace =.false.
78 real(real64) :: lamb_omega
79
80 logical, public :: lpfmf = .false.
81 real(real64), allocatable, public :: mf_vector_potential(:)
82 real(real64), allocatable :: jp_proj_eo(:,:)
85
86 contains
87 procedure :: init => xc_photons_init
88 procedure :: end => xc_photons_end
89 procedure :: wants_to_renormalize_mass => xc_photons_wants_to_renormalize_mass
90 procedure :: get_renormalized_mass => xc_photons_get_renormalized_emass
91 procedure :: mf_dump => xc_photons_mf_dump
92 procedure :: mf_load => xc_photons_mf_load
93 procedure :: v_ks => xc_photons_v_ks
94 procedure :: add_mean_field => xc_photons_add_mean_field
95 end type xc_photons_t
96
97 ! the PhotonXCXCMethod
98 integer, private, parameter :: &
99 XC_PHOTONS_NONE = 0, &
100 xc_photons_lda = 1, &
102
103contains
104
105 ! ---------------------------------------------------------
107 !
108 subroutine xc_photons_init(xc_photons, namespace, xc_photon, space, gr, st)
109 class(xc_photons_t), intent(out) :: xc_photons
110 type(namespace_t), intent(in) :: namespace
111 integer, intent(in) :: xc_photon
112 class(space_t), intent(in) :: space
113 type(grid_t), intent(in) :: gr
114 type(states_elec_t), intent(in) :: st
115
116 push_sub(xc_photons_init)
117
118 xc_photons%lpfmf = .false.
119
120 call messages_experimental("XCPhotonFunctional /= none")
121
122 call photon_mode_init(xc_photons%pt, namespace, space%dim, .true.)
123 call photon_mode_set_n_electrons(xc_photons%pt, st%qtot)
124
125 select case(xc_photon)
126 case(option__xcphotonfunctional__photon_x_lda)
127 xc_photons%method = xc_photons_lda
128 xc_photons%lcorrelations = .false.
129 case(option__xcphotonfunctional__photon_xc_lda)
130 xc_photons%method = xc_photons_lda
131 xc_photons%lcorrelations = .true.
132 case(option__xcphotonfunctional__photon_x_wfn)
133 xc_photons%method = xc_photons_wfs
134 xc_photons%lcorrelations = .false.
135 case(option__xcphotonfunctional__photon_xc_wfn)
136 xc_photons%method = xc_photons_wfs
137 xc_photons%lcorrelations = .true.
138 case default
139 xc_photons%method = xc_photons_none
140 return
141 end select
142
143
144 if (xc_photons%method == xc_photons_lda) then
145
146 !%Variable PhotonXCLDAKappa
147 !%Type float
148 !%Default 1.0
149 !%Section Hamiltonian::XC
150 !%Description
151 !% the scaling factor for px-LDA potential
152 !%End
153 call parse_variable(namespace, 'PhotonXCLDAKappa', m_one, xc_photons%pxlda_kappa)
154
155 end if
156
157 !%Variable PhotonXCEnergyMethod
158 !%Type integer
159 !%Default 1
160 !%Section Hamiltonian::XC
161 !%Description
162 !% There are different ways to calculate the energy,
163 !%Option virial 1
164 !% (modified) virial approach</br>
165 !% <math>
166 !% (E_{\rm{px}}^{\rm{virial}} = \frac{1}{2}\int d\mathbf{r}\ \mathbf{r}\cdot[
167 !% -\rho(\mathbf{r})\nabla v_{\rm{px}}(\mathbf{r})])
168 !% </math></br>
169 !%Option expectation_value 2
170 !% expectation value w.tr.t. the wave functions (valid only for 1 electron)</br>
171 !% <math>
172 !% E_{\rm{px}}[\rho] = -\sum_{\alpha=1}^{M_{p}}\frac{\tilde{\lambda}_{\alpha}^{2}}{2\tilde{\omega}_{\alpha}^{2}}
173 !% \langle (\tilde{\mathbf{{\varepsilon}}}_{\alpha}\cdot\hat{\mathbf{J}}_{\rm{p}})\Phi[\rho]
174 !% | (\tilde{\mathbf{{\varepsilon}}}_{\alpha}\cdot\hat{\mathbf{J}}_{\rm{p}})\Phi[\rho] \rangle
175 !% </math></br>
176 !% This option only works for the wave function based electron-photon functionals
177 !%Option LDA 3
178 !% energy from electron density</br>
179 !% <math>
180 !% E_{\rm pxLDA}[\rho] = \frac{-2\pi^{2}}{(d+2)({2V_{d}})^{\frac{2}{d}}}
181 !% \sum_{\alpha=1}^{M_{p}}\frac{\tilde{\lambda}_{\alpha}^{2}}{\tilde{\omega}_{\alpha}^{2}}
182 !% \int d\mathbf{r}\ \rho^{\frac{2+d}{d}}(\mathbf{r})
183 !% </math></br>
184 !% This option only works with LDA electron-photon functionals.
185 !%End
187 call parse_variable(namespace, 'PhotonXCEnergyMethod', 1, xc_photons%energy_method)
189 if( xc_photons%method == xc_photons_wfs .and. xc_photons%energy_method == option__photonxcenergymethod__lda ) then
190 message(1) = "Calculating the electron-photon energy from the LDA expression"
191 message(2) = "is not implemented for wave function based electron-photon functionals"
192 call messages_fatal(2, namespace=namespace)
193 end if
194
195
196 if (xc_photons%lcorrelations) then
197
198 !%Variable PhotonXCEtaC
199 !%Type float
200 !%Default 1.0
201 !%Section Hamiltonian::XC
202 !%Description
203 !% The scaling factor for the px potential to reduce the weak coupling perturbation regime
204 !%End
205
206 if (parse_is_defined(namespace, 'PhotonXCEtaC')) then
207 call parse_variable(namespace, 'PhotonXCEtaC', m_one, xc_photons%eta_c)
208 else
209 message(1) = "Defining PhotonXCEtaC is required for photon functionals containing correlation."
210 call messages_fatal(1, namespace=namespace)
211 end if
212
213 else
214
215 xc_photons%eta_c = m_one
216
217 end if
218
219 ! This variable will keep vpx across iterations
220 safe_allocate(xc_photons%vpx(1:gr%np_part))
221
222 xc_photons%vpx = m_zero
223
224 !%Variable PhotonXCLambShift
225 !%Type logical
226 !%Default .false.
227 !%Section Hamiltonian::XC
228 !%Description
229 !% to deal with the photon free exchange potential for continuum mode in free space
230 !%End
231
232 call parse_variable(namespace, 'PhotonXCLambShift', .false., xc_photons%llamb_freespace)
233 call messages_experimental("PhotonXCLambShift", namespace=namespace)
234
235 if (xc_photons%llamb_freespace) then
236
237 !%Variable PhotonXCLambShiftOmegaCutoff
238 !%Type float
239 !%Default 0.0
240 !%Section Hamiltonian::XC
241 !%Description
242 !% the cutoff frequency (Ha) for Lamb shift
243 !%End
244
245 call parse_variable(namespace, 'PhotonXCLambShiftOmegaCutoff', m_zero, xc_photons%lamb_omega)
246
247 !%Variable PhotonXCLambShiftRenormalizeMass
248 !%Type logical
249 !%Default .false.
250 !%Section Hamiltonian::XC
251 !%Description
252 !% to deal with the photon free exchange potential for continuum mode in free space
253 !%End
254
255 call parse_variable(namespace, 'PhotonXCLambShiftRenormalizeMass', .false., xc_photons%llamb_re_mass)
256
257 end if
258
259 ! compute the dressed photon modes
260 call photon_mode_dressed(xc_photons%pt)
261
262
263 pop_sub(xc_photons_init)
264
265 end subroutine xc_photons_init
266
267 ! ---------------------------------------------------------
268 subroutine xc_photons_end(this)
269 class(xc_photons_t), intent(inout) :: this
270
271 push_sub(xc_photons_end)
272
273 call photon_mode_end(this%pt)
274 safe_deallocate_a(this%vpx)
275
276 if (allocated(this%mf_vector_potential)) then
277 safe_deallocate_a(this%mf_vector_potential)
278 end if
279 if (allocated(this%jp_proj_eo)) then
280 safe_deallocate_a(this%jp_proj_eo)
281 end if
282
283 pop_sub(xc_photons_end)
284 end subroutine xc_photons_end
285
291 !
292 subroutine xc_photons_v_ks(xc_photons, namespace, total_density, gr, space, psolver, st)
293 class(xc_photons_t), intent(inout) :: xc_photons
294 type(namespace_t), intent(in) :: namespace
295 real(real64), pointer, contiguous, intent(in) :: total_density(:)
296 class(grid_t), intent(in) :: gr
297 type(space_t), intent(in) :: space
298 type(poisson_t), intent(in) :: psolver
299 type(states_elec_t), intent(inout) :: st
300
301 integer :: ia
302
303 push_sub(xc_photons_v_ks)
304
305 xc_photons%lpfmf = xc_photons%method > 0
306
307 xc_photons%vpx = m_zero
308 xc_photons%ex = m_zero
309
310 if ( .not. allocated(xc_photons%mf_vector_potential) ) then
311 safe_allocate(xc_photons%mf_vector_potential(1:space%dim))
312 xc_photons%mf_vector_potential = m_zero
313 end if
314 if ( .not. allocated(xc_photons%jp_proj_eo)) then
315 safe_allocate(xc_photons%jp_proj_eo(1:xc_photons%pt%nmodes,1:2))
316 xc_photons%jp_proj_eo = m_zero
317 end if
318
319
320 select case(xc_photons%method)
321 case(xc_photons_lda) ! LDA approximation for px potential
322 call photon_free_vpx_lda(namespace, xc_photons, total_density, gr, space, psolver)
323 case(xc_photons_wfs) ! wave function approxmation for px potential
324 call photon_free_vpx_wfc(namespace, xc_photons, total_density, gr, space, st)
325 case(xc_photons_none) ! no photon-exchange potential
326 call messages_write('Photon-free px potential is not computed', new_line = .true.)
327 call messages_info()
328 case default
329 assert(.false.)
330 end select
331
332
333 if (.not. xc_photons%llamb_freespace) then
334 ! add the constant energy shift
335 do ia = 1, xc_photons%pt%nmodes
336 xc_photons%ex = xc_photons%ex + 0.5_real64 * (xc_photons%pt%dressed_omega(ia)-xc_photons%pt%omega(ia))
337 end do
338 end if
339
340 pop_sub(xc_photons_v_ks)
341 end subroutine xc_photons_v_ks
342 ! ---------------------------------------------------------
343
344 ! ---------------------------------------------------------
345 ! ---------------------------------------------------------
346 !
378 ! The Lamb shift code is experimental and untested
379 subroutine photon_free_vpx_lda(namespace, xc_photons, total_density, gr, space, psolver)
380 type(namespace_t), intent(in) :: namespace
381 type(xc_photons_t), intent(inout) :: xc_photons
382 real(real64), pointer, contiguous, intent(in) :: total_density(:)
383 type(grid_t), target, intent(in) :: gr
384 type(space_t), intent(in) :: space
385 type(poisson_t), intent(in) :: psolver
386
387 integer :: ia, ip, iter
388 real(real64) :: unit_volume, r, res, presum, prefact
389 real(real64) :: xx(space%dim), prefactor_lamb
390 real(real64), allocatable :: prefactor(:)
391 real(real64), allocatable :: rho_aux(:)
392 real(real64), allocatable :: grad_rho_aux(:,:)
393 real(real64), allocatable :: px_source(:)
394 real(real64), allocatable :: tmp1(:)
395 real(real64), allocatable :: tmp2(:,:)
396 real(real64), allocatable :: tmp3(:)
397 real(real64), allocatable :: grad_vpx(:,:)
398 real(real64), allocatable :: epsgrad_epsgrad_rho_aux(:)
399 real(real64), allocatable :: epx_force_module(:)
400
401 real(real64), parameter :: threshold = 1e-7_real64
402
403 push_sub(photon_free_vpx_lda)
404
405 if (xc_photons%energy_method == 2 .and. xc_photons%pt%n_electrons >1) then
406 call messages_not_implemented("expectation value for energy for pxLDA more than 1 electron", namespace=namespace)
407 end if
408
409 xc_photons%vpx = m_zero
410 xc_photons%ex = m_zero
411
412 safe_allocate(prefactor(1:xc_photons%pt%nmodes))
413 prefactor = m_zero
414 ! here we will use only one spin channel
415 safe_allocate(rho_aux(1:gr%np_part))
416 safe_allocate(grad_rho_aux(1:gr%np, 1:xc_photons%pt%dim))
417 safe_allocate(px_source(1:gr%np_part))
418 safe_allocate(tmp1(1:gr%np_part))
419 safe_allocate(tmp2(1:gr%np, 1:xc_photons%pt%dim))
420 safe_allocate(tmp3(1:gr%np_part))
421 safe_allocate(grad_vpx(1:gr%np, 1:xc_photons%pt%dim))
422 grad_vpx = m_zero
423 safe_allocate(epx_force_module(1:gr%np_part))
424 epx_force_module = m_zero
425
426 select case(xc_photons%pt%dim)
427 case(1)
428 unit_volume = m_two
429 case(2)
430 unit_volume = m_pi
431 case(3)
432 unit_volume = m_four*m_pi/m_three
433 case default
434 call messages_not_implemented("LDA px more than 3 dimension", namespace=namespace)
435 end select
436
437 !$omp parallel do
438 do ip=1, gr%np
439 rho_aux(ip) = ( abs(total_density(ip))/(m_two*unit_volume) )**(m_two/(xc_photons%pt%dim*m_one))
440 end do
441 !$omp end parallel do
442
443
444 ! compute the electron-photon exchange potential
445
446 if (xc_photons%llamb_freespace) then
447
448 ! Note: The Lamb shift part is currently experimental and untested!
449 ! compute the electron-photon exchange potential
450
451 prefactor_lamb = -(8.0_real64*m_pi*m_third) * xc_photons%lamb_omega / (p_c**3)
452
453 !$OMP parallel do
454 do ip=1,gr%np
455 xc_photons%vpx(ip) = prefactor_lamb*rho_aux(ip)
456 end do
457 !$OMP end parallel do
458
459 else
460
461 do ia = 1, xc_photons%pt%nmodes
462 prefactor(ia) = -m_two*(m_pi * xc_photons%pt%dressed_lambda(ia) / xc_photons%pt%dressed_omega(ia))**2
463 end do
464
465 select case(xc_photons%pt%dim)
466 case(1)
467 ! solve the pxLDA potential using the analytical form in 1D
468 px_source = m_zero
469 do ia = 1, xc_photons%pt%nmodes
470 !$OMP parallel do
471 do ip=1,gr%np_part
472 px_source(ip) = px_source(ip) + prefactor(ia)
473 end do
474 !$OMP end parallel do
475 end do
476
477 !$OMP parallel do
478 do ip=1,gr%np
479 xc_photons%vpx(ip) = px_source(ip)*rho_aux(ip)
480 end do
481 !$OMP end parallel do
482 case(2)
483 call get_px_source(px_source)
484 ! for 2D we solve the Poisson equation using the conjugate gradient method
485
486 ! Note that we need to solve -\Delta f = v, as CG requires SPD operator
487 call lalg_scal(gr%np, -m_one, px_source)
488
489 mesh_aux => gr%der%mesh
490 iter = 1000
491 call dconjugate_gradients(gr%np, xc_photons%vpx(:), px_source, dlaplacian_op, dmf_dotp_aux, iter, res, threshold, &
492 userdata=[c_loc(gr%der)])
493 write(message(1),'(a,i6,a)') "Info: CG converged with ", iter, " iterations."
494 write(message(2),'(a,e14.6)') "Info: The residue is ", res
495 call messages_info(2, namespace=namespace)
496
497 case(3)
498 ! for 3D we use thepoisson solver including the prefactor (-4*pi)
499 ! therefore, we need to remove the factor in advance
500 call get_px_source(px_source)
501
502 call lalg_scal(gr%np, m_one/(-m_four*m_pi), px_source)
503
504 ! solve the Poisson equation
505 call dpoisson_solve(psolver, namespace, xc_photons%vpx(:), px_source)
506 case default
507 assert(.false.)
508 end select
509
510 end if
511
512 ! scaling the potential
513 call lalg_scal(gr%np, (xc_photons%eta_c * xc_photons%pxlda_kappa), xc_photons%vpx)
514
515 ! compute electron-photon energy
516
517 select case (xc_photons%energy_method)
518 case(1) ! compute the epx energy from the virial relation
519
520 do ia = 1, xc_photons%pt%nmodes
521
522 ! compute the electron-photon force
523 !$omp parallel do
524 do ip = 1, gr%np
525 epx_force_module(ip) = -prefactor(ia)*m_two*abs(total_density(ip))*rho_aux(ip)/(xc_photons%pt%dim*m_one+m_two)
526 end do
527 !$omp end parallel do
528
529 call dderivatives_grad(gr%der, epx_force_module(:), tmp2)
530 call lalg_gemv(gr%np, xc_photons%pt%dim, m_one, tmp2, xc_photons%pt%dressed_pol(1:xc_photons%pt%dim, ia), m_zero, tmp1)
531
532 !$omp parallel do private(r, xx)
533 do ip = 1, gr%np
534 call mesh_r(gr, ip, r, coords=xx)
535 tmp3(ip) = tmp1(ip)*dot_product(xx(1:xc_photons%pt%dim), xc_photons%pt%dressed_pol(1:xc_photons%pt%dim, ia))
536 end do
537 !$omp end parallel do
538
539 xc_photons%ex = xc_photons%ex + m_half*dmf_integrate(gr, tmp3)
540 end do
541
542 xc_photons%ex = xc_photons%eta_c * xc_photons%pxlda_kappa * xc_photons%ex
543
544 case(2) ! compute the energy as expetation value wrt to the wave functions
545
546 rho_aux(1:gr%np) = sqrt(abs( total_density(1:gr%np)))
547
548 safe_allocate(epsgrad_epsgrad_rho_aux(1:gr%np))
549
550 call dderivatives_grad(gr%der, rho_aux(1:gr%np_part), grad_rho_aux)
551
552 do ia = 1, xc_photons%pt%nmodes
553 prefact = (xc_photons%pt%dressed_lambda(ia)**2) / (m_two*xc_photons%pt%dressed_omega(ia)**2)
554
555 call lalg_gemv(gr%np, xc_photons%pt%dim, m_one, grad_rho_aux, xc_photons%pt%dressed_pol(:, ia), m_zero, tmp1)
556 call dderivatives_grad(gr%der, tmp1, tmp2)
557 call lalg_gemv(gr%np, xc_photons%pt%dim, m_one, tmp2, xc_photons%pt%dressed_pol(1:xc_photons%pt%dim, ia), m_zero, tmp1)
558
559 !$OMP parallel do
560 do ip=1, gr%np
561 epsgrad_epsgrad_rho_aux(ip) = epsgrad_epsgrad_rho_aux(ip) + prefact*tmp1(ip)
562 end do
563 !$OMP end parallel do
564
565 end do
566
567 xc_photons%ex = xc_photons%eta_c * dmf_dotp(gr, rho_aux(1:gr%np), epsgrad_epsgrad_rho_aux(1:gr%np))
568
569 safe_deallocate_a(epsgrad_epsgrad_rho_aux)
570
571 case(3) ! integrate the aux electron density over the volume
572
573 !$OMP parallel do
574 do ip=1,gr%np
575 tmp1(ip) = abs( total_density(ip))**((m_one*xc_photons%pt%dim+m_two)/(xc_photons%pt%dim*m_one))
576 end do
577 !$OMP end parallel do
578
579 ! sum over the prefactors
580 presum = m_zero
581 do ia = 1, xc_photons%pt%nmodes
582 presum = presum + prefactor(ia)
583 end do
584 presum = presum * (m_one/(m_two*unit_volume))**(m_two/(xc_photons%pt%dim*m_one)) / (xc_photons%pt%dim*m_one+m_two)
585 presum = presum * xc_photons%eta_c * xc_photons%pxlda_kappa
586
587 xc_photons%ex = xc_photons%eta_c * xc_photons%pxlda_kappa * presum * dmf_integrate(gr, tmp1)
588
589 end select
590
591 safe_deallocate_a(prefactor)
592 safe_deallocate_a(rho_aux)
593 safe_deallocate_a(grad_rho_aux)
594 safe_deallocate_a(px_source)
595 safe_deallocate_a(tmp1)
596 safe_deallocate_a(tmp2)
597 safe_deallocate_a(tmp3)
598 safe_deallocate_a(grad_vpx)
599 safe_deallocate_a(epx_force_module)
600
601 pop_sub(photon_free_vpx_lda)
602
603 contains
604
605 subroutine get_px_source(px_source)
606 real(real64), contiguous, intent(out) :: px_source(:)
607
608 call dderivatives_grad(gr%der, rho_aux(1:gr%np_part), grad_rho_aux)
609
610 px_source = m_zero
611 do ia = 1, xc_photons%pt%nmodes
612 call lalg_gemv(gr%np, xc_photons%pt%dim, m_one, grad_rho_aux, xc_photons%pt%dressed_pol(:, ia), m_zero, tmp1)
613 call dderivatives_grad(gr%der, tmp1, tmp2)
614 call lalg_gemv(gr%np, xc_photons%pt%dim, m_one, tmp2, xc_photons%pt%dressed_pol(1:xc_photons%pt%dim, ia), m_zero, tmp1)
615 call lalg_axpy(gr%np, prefactor(ia), tmp1, px_source)
616 end do
617
618 end subroutine get_px_source
619
620 end subroutine photon_free_vpx_lda
621 ! ---------------------------------------------------------
622
623 ! ---------------------------------------------------------
624 !
650 !
651 ! The Lamb shift code is experimental and untested
652
653 subroutine photon_free_vpx_wfc(namespace, xc_photons, total_density, gr, space, st)
654 type(namespace_t), intent(in) :: namespace
655 type(xc_photons_t), intent(inout) :: xc_photons
656 real(real64), pointer, contiguous, intent(in) :: total_density(:)
657 class(grid_t), intent(in) :: gr
658 type(space_t), intent(in) :: space
659 type(states_elec_t), intent(in) :: st
660
661 integer :: ia, ip
662 real(real64) :: prefactor_lamb
663 real(real64) :: xx(space%dim), r
664 real(real64), allocatable :: prefactor(:)
665 real(real64), allocatable :: rho_aux(:)
666 real(real64), allocatable :: grad_rho_aux(:,:)
667 real(real64), allocatable :: grad_vpx(:,:)
668 real(real64), allocatable :: epsgrad_epsgrad_rho_aux(:)
669 real(real64), allocatable :: tmp1(:)
670 real(real64), allocatable :: tmp2(:,:)
671 real(real64) :: shift
672
673 push_sub(photon_free_vpx_wfc)
674
675 if (st%d%nspin >1) then
676 call messages_not_implemented("PhotonXCXCMethod = wavefunction for polarized and spinor cases", namespace=namespace)
677 end if
678
679 if (xc_photons%pt%n_electrons >1) then
680 call messages_not_implemented("PhotonXCXCMethod = wavefunction for more than 1 electron", namespace=namespace)
681 end if
682
683 xc_photons%vpx = m_zero
684 xc_photons%ex = m_zero
685
686 safe_allocate(prefactor(1:xc_photons%pt%nmodes))
687 prefactor = m_zero
688 safe_allocate(rho_aux(1:gr%np_part))
689 rho_aux = m_zero
690 safe_allocate(grad_rho_aux(1:gr%np, 1:xc_photons%pt%dim))
691 grad_rho_aux = m_zero
692 safe_allocate(epsgrad_epsgrad_rho_aux(1:gr%np))
693 safe_allocate(grad_vpx(1:gr%np, 1:xc_photons%pt%dim))
694 safe_allocate(tmp1(1:gr%np_part))
695 safe_allocate(tmp2(1:gr%np_part, 1:xc_photons%pt%dim))
696
697
698 rho_aux(1:gr%np) = sqrt(abs( total_density(1:gr%np)))
699 !$OMP parallel do
700 do ip = 1, gr%np
701 rho_aux(ip) = safe_tol(rho_aux(ip),1e-18_real64)
702 end do
703 !$OMP end parallel do
704
705 if (xc_photons%llamb_freespace) then
706
707 ! experimental Lamb shift mode
708
709 prefactor_lamb = ( m_two/(m_three*m_pi) ) * xc_photons%lamb_omega / p_c**3
710 call dderivatives_lapl(gr%der, rho_aux(1:gr%np_part), epsgrad_epsgrad_rho_aux)
711 call lalg_scal(gr%np, prefactor_lamb, epsgrad_epsgrad_rho_aux)
712
713 else
714
715 do ia = 1, xc_photons%pt%nmodes
716 prefactor(ia) = (xc_photons%pt%dressed_lambda(ia)**2) / (m_two*xc_photons%pt%dressed_omega(ia)**2)
717 end do
718
719 call dderivatives_grad(gr%der, rho_aux, grad_rho_aux)
720
721 epsgrad_epsgrad_rho_aux = m_zero
722 do ia = 1, xc_photons%pt%nmodes
723 tmp1 = m_zero
724 call lalg_gemv(gr%np, xc_photons%pt%dim, m_one, grad_rho_aux, xc_photons%pt%dressed_pol(:, ia), m_zero, tmp1)
725 call dderivatives_grad(gr%der, tmp1, tmp2)
726 call lalg_gemv(gr%np, xc_photons%pt%dim, prefactor(ia), tmp2, xc_photons%pt%dressed_pol(:, ia), &
727 m_one, epsgrad_epsgrad_rho_aux)
728 end do
729
730 end if
731
732 !$OMP parallel do
733 do ip = 1, gr%np
734 xc_photons%vpx(ip) = epsgrad_epsgrad_rho_aux(ip)/rho_aux(ip)
735 end do
736 !$OMP end parallel do
737
738 if(st%eigenval(1,1) < m_huge .and. .not. space%is_periodic()) then
739 shift = m_two * st%eigenval(1,1) * prefactor(1)
740 !$OMP parallel do
741 do ip=1, gr%np
742 xc_photons%vpx(ip) = xc_photons%vpx(ip) + shift
743 end do
744 !$OMP end parallel do
745 end if
746
747 call lalg_scal(gr%np, xc_photons%eta_c, xc_photons%vpx(:))
749 select case (xc_photons%energy_method)
750 case(1) ! virial
751
752 call dderivatives_grad(gr%der, xc_photons%vpx(:), grad_vpx)
753 !$OMP parallel do private(r, xx)
754 do ip = 1, gr%np
755 call mesh_r(gr, ip, r, coords=xx)
756 tmp1(ip) = - total_density(ip)*dot_product(xx(1:xc_photons%pt%dim), grad_vpx(ip,1:xc_photons%pt%dim))
757 end do
758 !$OMP end parallel do
759
760 xc_photons%ex = m_half*dmf_integrate(gr, tmp1)
761
762 case(2) ! expectation_value
763 xc_photons%ex = xc_photons%eta_c * dmf_dotp(gr, rho_aux(1:gr%np), epsgrad_epsgrad_rho_aux)
764
765 case default
766 assert(.false.)
767 end select
768
769 safe_deallocate_a(prefactor)
770 safe_deallocate_a(rho_aux)
771 safe_deallocate_a(grad_rho_aux)
772 safe_deallocate_a(grad_vpx)
773 safe_deallocate_a(epsgrad_epsgrad_rho_aux)
774 safe_deallocate_a(tmp1)
775 safe_deallocate_a(tmp2)
776
777 pop_sub(photon_free_vpx_wfc)
778
779 end subroutine photon_free_vpx_wfc
780 ! ---------------------------------------------------------
781
782
783 ! ---------------------------------------------------------
793 !
794 subroutine xc_photons_add_mean_field(xc_photons, gr, space, kpoints, st, time, dt)
795 class(xc_photons_t), intent(inout) :: xc_photons
796 class(grid_t), intent(in) :: gr
797 class(space_t), intent(in) :: space
798 type(kpoints_t), intent(in) :: kpoints
799 type(states_elec_t), intent(inout) :: st
800 real(real64), intent(in) :: time
801 real(real64), intent(in) :: dt
802
803
804 integer :: ia, idir, ispin
805 real(real64) :: pol_dot_jp
806 real(real64), allocatable :: current(:,:,:)
807 real(real64), allocatable :: jp(:)
808
810
811 ! compute the dressed photon-free vector potential
812 xc_photons%mf_vector_potential = m_zero
813
814 safe_allocate(current(1:gr%np_part, 1:space%dim, 1:st%d%nspin))
815 current = m_zero
816 ! here we use the paramagnetic current; note that the physical current here
817 ! only contains the paramagnetic current.
818 call states_elec_calc_quantities(gr, st, kpoints, .false., paramagnetic_current = current)
819
820 ! compute the paramagnetic current
821 safe_allocate(jp(1:space%dim))
822 do idir = 1, space%dim
823 jp(idir) = m_zero
824 do ispin = 1, st%d%spin_channels
825 jp(idir) = jp(idir) + dmf_integrate(gr%der%mesh, current(:,idir,ispin))
826 end do
827 end do
828
829 ! update the 'projected' current
830 do ia = 1, xc_photons%pt%nmodes
831 pol_dot_jp = dot_product(xc_photons%pt%dressed_pol(1:space%dim, ia),jp(1:space%dim))
832 xc_photons%jp_proj_eo(ia,1) = xc_photons%jp_proj_eo(ia,1) + &
833 cos(xc_photons%pt%dressed_omega(ia)*( time-dt))*pol_dot_jp*dt
834 xc_photons%jp_proj_eo(ia,2) = xc_photons%jp_proj_eo(ia,2) + &
835 sin(xc_photons%pt%dressed_omega(ia)*( time-dt))*pol_dot_jp*dt
836 end do
837
838 do ia = 1, xc_photons%pt%nmodes
839 xc_photons%mf_vector_potential(1:xc_photons%pt%dim) = xc_photons%mf_vector_potential(1:xc_photons%pt%dim) &
840 + (-p_c*(xc_photons%pt%dressed_lambda(ia)**2) / xc_photons%pt%dressed_omega(ia)) &
841 * (xc_photons%jp_proj_eo(ia,1)*sin(xc_photons%pt%dressed_omega(ia)* time) &
842 - xc_photons%jp_proj_eo(ia,2)*cos(xc_photons%pt%dressed_omega(ia)* time)) &
843 * xc_photons%pt%dressed_pol(1:xc_photons%pt%dim, ia)
844 end do
845
846 safe_deallocate_a(current)
847 safe_deallocate_a(jp)
848
850
851 end subroutine xc_photons_add_mean_field
852 ! ---------------------------------------------------------
853
854
858 !
859 logical pure function xc_photons_wants_to_renormalize_mass(xc_photons) result (renorm)
860 class(xc_photons_t), intent(in) :: xc_photons
861
862 renorm = (xc_photons%method>0) .and. xc_photons%llamb_freespace .and. xc_photons%llamb_re_mass
863
865
867 real(real64) pure function xc_photons_get_renormalized_emass(xc_photons) result(mass)
868 class(xc_photons_t), intent(in) :: xc_photons
869
870 mass = m_one - (m_four*xc_photons%lamb_omega) / (3.0*m_pi * p_c**3)
871
873
875 !
876 subroutine xc_photons_mf_dump(xc_photons, restart, ierr)
877 class(xc_photons_t), intent(in) :: xc_photons
878 type(restart_t), intent(in) :: restart
879 integer, intent(out) :: ierr
880
881 character(len=80), allocatable :: lines(:)
882 integer :: iunit, err, jj, nmodes
883 integer :: pt_dim
884
885 push_sub(photon_free_mf_dump)
886 nmodes = xc_photons%pt%nmodes
887 pt_dim = xc_photons%pt%dim
888
889 safe_allocate(lines(1:nmodes+pt_dim))
890
891 ierr = 0
892
893 iunit = restart%open('photon_free_mf')
894
895 do jj = 1, pt_dim
896 write(lines(jj), '(2x, es19.12)') xc_photons%mf_vector_potential(jj)
897 end do
898
899 do jj = 1, nmodes
900 write(lines(jj+pt_dim), '(a10,1x,I8,a1,2x,2(es19.12,2x))') 'Mode ', jj, ":", xc_photons%jp_proj_eo(jj,1:2)
901 end do
902
903 call restart%write(iunit, lines, nmodes+pt_dim, err)
904 if (err /= 0) ierr = ierr + 1
905 call restart%close(iunit)
906
907 safe_deallocate_a(lines)
908
909 pop_sub(photon_free_mf_dump)
910 end subroutine xc_photons_mf_dump
911
912! ---------------------------------------------------------
913
915 !
916 subroutine xc_photons_mf_load(xc_photons, restart, space, ierr)
917 class(xc_photons_t), intent(inout) :: xc_photons
918 type(restart_t), intent(in) :: restart
919 class(space_t), intent(in) :: space
920 integer, intent(out) :: ierr
921
922 character(len=80), allocatable :: lines(:)
923 character(len=7) :: sdummy
924 integer :: idummy
925 integer :: iunit, err, jj, nmodes
926 integer :: pt_dim
927
928 push_sub(photon_free_mf_load)
929
930 ierr = 0
931 nmodes = xc_photons%pt%nmodes
932 pt_dim = xc_photons%pt%dim
933
934 if (restart%skip()) then
935 ierr = -1
936 pop_sub(photon_free_mf_load)
937 return
938 end if
939
940 message(1) = "Debug: Reading Photon-Free Photons restart."
941 call messages_info(1, namespace=restart%namespace, debug_only=.true.)
942
943 if ( .not. allocated(xc_photons%jp_proj_eo)) then
944 safe_allocate(xc_photons%jp_proj_eo(1:xc_photons%pt%nmodes, 1:2))
945 xc_photons%jp_proj_eo = m_zero
946 end if
947 if ( .not. allocated(xc_photons%mf_vector_potential)) then
948 safe_allocate(xc_photons%mf_vector_potential(1:space%dim))
949 xc_photons%mf_vector_potential = m_zero
950 end if
951
952 safe_allocate(lines(1:nmodes+pt_dim))
953 iunit = restart%open('photon_free_mf')
954 call restart%read(iunit, lines, nmodes+pt_dim, err)
955 if (err /= 0) then
956 ierr = ierr + 1
957 else
958 do jj = 1, pt_dim
959 read(lines(jj),'(2x, es19.12)') xc_photons%mf_vector_potential(jj)
960 end do
961
962 do jj = 1, nmodes
963 read(lines(jj+pt_dim), '(a10,1x,I8,a1,2x,2(es19.12,2x))') sdummy, idummy, sdummy, xc_photons%jp_proj_eo(jj,1:2)
964 end do
965 end if
966 call restart%close(iunit)
967
968 message(1) = "Debug: Reading Photons restart done."
969 call messages_info(1, namespace=restart%namespace, debug_only=.true.)
970
971 safe_deallocate_a(lines)
972
973 pop_sub(photon_free_mf_load)
974 end subroutine xc_photons_mf_load
975
976end module xc_photons_oct_m
977
978!! Local Variables:
979!! mode: f90
980!! coding: utf-8
981!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:173
scales a vector by a constant
Definition: lalg_basic.F90:159
double sin(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine, public dderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_huge
Definition: global.F90:209
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_third
Definition: global.F90:198
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:233
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_three
Definition: global.F90:194
This module implements the underlying real-space grid.
Definition: grid.F90:119
Computes and , suitable as an operator callback for iterative solvers (CG, QMR, etc....
subroutine, public dlaplacian_op(x, hx, userdata)
Computes the negative Laplacian operator action: .
This module defines various routines, operating on mesh functions.
class(mesh_t), pointer, public mesh_aux
Globally-scoped pointer to the mesh instance.
real(real64) function, public dmf_dotp_aux(f1, f2)
dot product between two vectors (mesh functions)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:342
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1091
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:455
subroutine, public photon_mode_end(this)
subroutine, public photon_mode_dressed(this)
subroutine, public photon_mode_set_n_electrons(this, qtot)
subroutine, public photon_mode_init(this, namespace, dim, photon_free)
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:875
This module is intended to contain "only mathematical" functions and procedures.
Definition: solvers.F90:117
subroutine, public states_elec_calc_quantities(gr, st, kpoints, nlcc, kinetic_energy_density, paramagnetic_current, density_gradient, density_laplacian, gi_kinetic_energy_density, st_end)
calculated selected quantities
This module implements the "photon-free" electron-photon exchange-correlation functional.
Definition: xc_photons.F90:123
subroutine photon_free_vpx_wfc(namespace, xc_photons, total_density, gr, space, st)
compute the electron-photon exchange potential based on wave functions
Definition: xc_photons.F90:749
logical pure function xc_photons_wants_to_renormalize_mass(xc_photons)
indicate whether the photon-exchange requires a renormalized electron mass
Definition: xc_photons.F90:955
subroutine xc_photons_init(xc_photons, namespace, xc_photon, space, gr, st)
initialize the photon-exchange functional
Definition: xc_photons.F90:204
subroutine xc_photons_add_mean_field(xc_photons, gr, space, kpoints, st, time, dt)
accumulate the results of time integral the paramagnetic current.
Definition: xc_photons.F90:890
subroutine xc_photons_mf_dump(xc_photons, restart, ierr)
write restart information
Definition: xc_photons.F90:972
subroutine xc_photons_v_ks(xc_photons, namespace, total_density, gr, space, psolver, st)
evaluate the KS potential and energy for the given functional
Definition: xc_photons.F90:388
real(real64) pure function xc_photons_get_renormalized_emass(xc_photons)
return the renormalized electron mass for the electron-photon exhange
Definition: xc_photons.F90:963
integer, parameter, private xc_photons_lda
Definition: xc_photons.F90:193
subroutine photon_free_vpx_lda(namespace, xc_photons, total_density, gr, space, psolver)
compute the electron-photon exchange potential within the LDA
Definition: xc_photons.F90:475
subroutine xc_photons_end(this)
Definition: xc_photons.F90:364
integer, parameter, private xc_photons_wfs
Definition: xc_photons.F90:193
subroutine xc_photons_mf_load(xc_photons, restart, space, ierr)
load restart information
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
The states_elec_t class contains all electronic wave functions.
This class described the 'photon-exchange' electron-photon xc functionals, based on QEDFT.
Definition: xc_photons.F90:160
int true(void)
subroutine get_px_source(px_source)
Definition: xc_photons.F90:701