Octopus
photon_mode.F90
Go to the documentation of this file.
1!! Copyright (C) 2017 Johannes Flick
2!! Copyright (C) 2019 F. Buchholz, M. Oliveira
3!! Copyright (C) 2020 S. Ohlmann
4!!
5!! This program is free software; you can redistribute it and/or modify
6!! it under the terms of the GNU General Public License as published by
7!! the Free Software Foundation; either version 2, or (at your option)
8!! any later version.
9!!
10!! This program is distributed in the hope that it will be useful,
11!! but WITHOUT ANY WARRANTY; without even the implied warranty of
12!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13!! GNU General Public License for more details.
14!!
15!! You should have received a copy of the GNU General Public License
16!! along with this program; if not, write to the Free Software
17!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
18!! 02110-1301, USA.
19!!
20
21#include "global.h"
22
24 use comm_oct_m
25 use debug_oct_m
26 use global_oct_m
27 use io_oct_m
29 use, intrinsic :: iso_fortran_env
30 use mesh_oct_m
33 use mpi_oct_m
35 use parser_oct_m
37 use unit_oct_m
39
40 implicit none
41
42 private
43 public :: &
52
53 type photon_mode_t
54 ! All components are public by default
55 integer :: nmodes
56 integer :: dim
57 real(real64), allocatable :: omega(:)
58 real(real64), allocatable :: lambda(:)
59 real(real64), allocatable :: pol(:,:)
60 real(real64), allocatable :: pol_dipole(:,:)
61 real(real64) :: ex
62 real(real64), allocatable :: number(:)
63 real(real64), allocatable :: correlator(:,:)
64 real(real64) :: n_electrons
65 real(real64), pointer :: pt_coord_q0(:) => null()
66 real(real64), pointer :: pt_momen_p0(:) => null()
67 real(real64) :: mu
68 ! dressed photons for photon-free approach
69 logical :: use_photon_exchange
70 real(real64), allocatable :: dressed_omega(:)
71 real(real64), allocatable :: dressed_lambda(:)
72 real(real64), allocatable :: dressed_pol(:,:)
73 ! !! dimension (1:this%dim, 1:this%nmodes)
74 real(real64), allocatable :: rotated_omega(:)
75
76 logical :: has_q0_p0 = .false.
77 end type photon_mode_t
78
79contains
80
81 ! ---------------------------------------------------------
82 subroutine photon_mode_init(this, namespace, dim, photon_free)
83 type(photon_mode_t), intent(out) :: this
84 type(namespace_t), intent(in) :: namespace
85 integer, intent(in) :: dim
86 logical, optional, intent(in) :: photon_free
87
88 type(block_t) :: blk
89 integer :: ii, idir, iunit, ncols
90 logical :: file_exists
91 character(256) :: filename
92
93 push_sub(photon_mode_init)
94
95 this%nmodes = 0
96
97 this%dim = dim
98 this%has_q0_p0 = .false.
99
100 !%Variable PhotonmodesFilename
101 !%Type string
102 !%Default "photonmodes"
103 !%Section Linear Response::Casida
104 !%Description
105 !% Filename for photon modes in text format
106 !% - first line contains 2 integers: number of photon modes and number of
107 !% columns
108 !% - each further line contains the given number of floats for one photon
109 !% mode
110 !%End
111 call parse_variable(namespace, 'PhotonmodesFilename', 'photonmodes', filename)
112 inquire(file=trim(filename), exist=file_exists)
113 if (file_exists) then
114 if (mpi_grp_is_root(mpi_world)) then
115 message(1) = 'Opening '//trim(filename)
116 call messages_info(1, namespace=namespace)
117 ! open file on root
118 iunit = io_open(trim(filename), namespace, action='read', form='formatted')
119
120 ! get dimensions from first line
121 read(iunit, *) this%nmodes, ncols
122
123 write(message(1), '(3a,i7,a,i3,a)') 'Reading file ', trim(filename), ' with ', &
124 this%nmodes, ' photon modes and ', ncols, ' columns.'
125 call messages_info(1, namespace=namespace)
126
127 safe_allocate(this%omega(1:this%nmodes))
128 safe_allocate(this%lambda(1:this%nmodes))
129 safe_allocate(this%pol(1:this%nmodes,1:3))
130
131 ! now read in all modes
132 do ii = 1, this%nmodes
133 if (ncols == 5) then
134 read(iunit, *) this%omega(ii), this%lambda(ii), &
135 this%pol(ii,1), this%pol(ii,2), this%pol(ii,3)
136 else if (ncols == 7) then
137 this%has_q0_p0 = .true.
138 if (.not. associated(this%pt_coord_q0)) then
139 safe_allocate(this%pt_coord_q0(1:this%nmodes))
140 end if
141 if (.not. associated(this%pt_momen_p0)) then
142 safe_allocate(this%pt_momen_p0(1:this%nmodes))
143 end if
144 read(iunit, *) this%omega(ii), this%lambda(ii), &
145 this%pol(ii,1), this%pol(ii,2), this%pol(ii,3), &
146 this%pt_coord_q0(ii), this%pt_momen_p0(ii)
147 else
148 ! error if not 5 columns
149 message(1) = 'Error: unexpected number of columns in file:'
150 message(2) = filename
151 call messages_fatal(2, namespace=namespace)
152 end if
154 ! Normalize polarization vector
155 this%pol(ii,:) = this%pol(ii,:)/norm2(this%pol(ii,:))
156 end do
157 call io_close(iunit)
158 end if
159 ! broadcast first array dimensions, then allocate and broadcast arrays
160 call mpi_world%bcast(this%nmodes, 1, mpi_integer, 0)
161 call mpi_world%bcast(ncols, 1, mpi_integer, 0)
162 if (.not. mpi_grp_is_root(mpi_world)) then
163 safe_allocate(this%omega(1:this%nmodes))
164 safe_allocate(this%lambda(1:this%nmodes))
165 safe_allocate(this%pol(1:this%nmodes,1:3))
166 if (ncols == 7) then
167 safe_allocate(this%pt_coord_q0(1:this%nmodes))
168 safe_allocate(this%pt_momen_p0(1:this%nmodes))
169 end if
170 end if
171 call mpi_world%bcast(this%omega(1), this%nmodes, mpi_double_precision, 0)
172 call mpi_world%bcast(this%lambda(1), this%nmodes, mpi_double_precision, 0)
173 call mpi_world%bcast(this%pol(1,1), this%nmodes*3, mpi_double_precision, 0)
174 if (ncols == 7) then
175 call mpi_world%bcast(this%pt_coord_q0(1), this%nmodes, mpi_double_precision, 0)
176 call mpi_world%bcast(this%pt_momen_p0(1), this%nmodes, mpi_double_precision, 0)
177 end if
178 else
179 if (.not. parse_is_defined(namespace, 'PhotonModes')) then
180 call messages_write('You need to specify the correct external photon modes file')
181 call messages_write('or define the PhotonModes variable!')
182 call messages_fatal(namespace=namespace)
183 end if
184 end if
185
186 !%Variable PhotonModes
187 !%Type block
188 !%Section Hamiltonian::XC
189 !%Description
190 !% Each line of the block should specify one photon mode. The syntax is the following:
191 !%
192 !% %PhotonModes
193 !% omega1 | lambda1| PolX1 | PolY1 | PolZ1
194 !% ...
195 !% %
196 !%
197 !% The first column is the mode frequency, in units of energy.
198 !% The second column is the coupling strength, in units of energy.
199 !% The remaining columns specify the polarization direction of the mode.
200 !% If the polarization vector should be normalized to one. If that is not the case
201 !% the code will normalize it.
202 !%End
203
204 if (.not. file_exists) then
205 if (parse_block(namespace, 'PhotonModes', blk) == 0) then
206
207 this%nmodes = parse_block_n(blk)
208 safe_allocate(this%omega(1:this%nmodes))
209 safe_allocate(this%lambda(1:this%nmodes))
210 safe_allocate(this%pol(1:this%nmodes, 1:this%dim))
211
212 this%pol = m_zero
213
214 do ii = 1, this%nmodes
215 ncols = parse_block_cols(blk, ii-1)
216
217 ! Read line
218 call parse_block_float(blk, ii-1, 0, this%omega(ii), units_inp%energy) ! frequency
219 call parse_block_float(blk, ii-1, 1, this%lambda(ii), units_inp%energy) ! coupling strength
220 do idir = 1, this%dim
221 call parse_block_float(blk, ii-1, idir + 1, this%pol(ii, idir)) ! polarization vector components
222 end do
223
224 if (ncols > this%dim + 2) then
225 this%has_q0_p0 = .true.
226 if (.not. associated(this%pt_coord_q0)) then
227 safe_allocate(this%pt_coord_q0(1:this%nmodes))
228 end if
229 if (.not. associated(this%pt_momen_p0)) then
230 safe_allocate(this%pt_momen_p0(1:this%nmodes))
231 end if
232 call parse_block_float(blk, ii-1, this%dim+2, this%pt_coord_q0(ii)) !row, column
233 call parse_block_float(blk, ii-1, this%dim+3, this%pt_momen_p0(ii)) !row, column
234 end if
235 ! Sanity check
236 if ((ncols /= this%dim + 2) .and. (ncols /= this%dim + 4)) then
237 call messages_input_error(namespace, 'PhotonModes', 'Incorrect number of columns')
238 end if
239
240 ! Normalize polarization vector
241 this%pol(ii,:) = this%pol(ii,:)/norm2(this%pol(ii,:))
242
243 end do
244 call parse_block_end(blk)
245 else
246 call messages_write('You need to specify the photon modes!')
247 call messages_fatal(namespace=namespace)
248 end if
249 end if
250
251 !%Variable TDPhotonicTimeScale
252 !%Type float
253 !%Default 1.0
254 !%Section Time-Dependent::Propagation
255 !%Description
256 !% This variable defines the factor between the timescale of photonic
257 !% and electronic movement.
258 !% for more details see the documentation of TDIonicTimeScale
259 !% If you also use TDIonicTimeScale, we advise to set
260 !% TDPhotonicTimeScale = TDIonicTimeScale, in the case the
261 !% photon frequency is in a vibrational energy range.
262 !% Important: The electronic time step will be the value of
263 !% <tt>TDTimeStep</tt> divided by this variable, so if you have determined an
264 !% optimal electronic time step (that we can call <i>dte</i>), it is
265 !% recommended that you define your time step as:
266 !%
267 !% <tt>TDTimeStep</tt> = <i>dte</i> * <tt>TDPhotonicTimeScale</tt>
268 !%
269 !% so you will always use the optimal electronic time step
270 !% (<a href=http://arxiv.org/abs/0710.3321>more details</a>).
271 !%End
272
273 call parse_variable(namespace, 'TDPhotonicTimeScale', 1.0_real64, this%mu)
274
275 this%ex = m_zero
276 safe_allocate(this%number(1:this%nmodes))
277 this%number = m_zero
278
279 ! photon-free approach, which requires the dressed photon modes
280 this%use_photon_exchange = .false.
281
282 this%use_photon_exchange = optional_default(photon_free, .false.)
283 if (this%use_photon_exchange) then
284 safe_allocate(this%dressed_omega(1:this%nmodes))
285 safe_allocate(this%dressed_lambda(1:this%nmodes))
286 safe_allocate(this%dressed_pol(1:this%dim, 1:this%nmodes))
287 safe_allocate(this%rotated_omega(1:this%nmodes))
288 end if
289
290 pop_sub(photon_mode_init)
291 end subroutine photon_mode_init
292
293 ! ---------------------------------------------------------
294 subroutine photon_mode_end(this)
295 type(photon_mode_t), intent(inout) :: this
296
297 push_sub(photon_mode_end)
298
299 safe_deallocate_a(this%correlator)
300
301 safe_deallocate_a(this%omega)
302 safe_deallocate_a(this%lambda)
303 safe_deallocate_a(this%number)
304
305 safe_deallocate_a(this%pol)
306 safe_deallocate_a(this%pol_dipole)
307
308 if (this%has_q0_p0) then
309 safe_deallocate_p(this%pt_coord_q0)
310 safe_deallocate_p(this%pt_momen_p0)
311 end if
312
313 if (this%use_photon_exchange) then
314 safe_deallocate_a(this%dressed_omega)
315 safe_deallocate_a(this%dressed_lambda)
316 safe_deallocate_a(this%dressed_pol)
317 safe_deallocate_a(this%rotated_omega)
318 end if
319
320 pop_sub(photon_mode_end)
321 end subroutine photon_mode_end
322
323 !-----------------------------------------------------------------
324 subroutine photon_mode_write_info(this, iunit, namespace)
325 type(photon_mode_t), intent(in) :: this
326 integer, optional, intent(in) :: iunit
327 type(namespace_t), optional, intent(in) :: namespace
328
329 integer :: im, idir
330
331 push_sub(photon_mode_write_info)
332
333 call messages_print_with_emphasis(msg="Photon Modes", iunit=iunit, namespace=namespace)
334 write(iunit, '(6x,a,10x,a,3x)', advance='no') 'Omega', 'Lambda'
335 write(iunit, '(1x,a,i1,a)') ('Pol.(', idir, ')', idir = 1, this%dim)
336 do im = 1, this%nmodes
337 write(iunit, '(1x,f14.12)', advance='no') this%omega(im)
338 write(iunit, '(1x,f14.12)', advance='no') this%lambda(im)
339 write(iunit, '(2X,f5.3,1X)') (this%pol(im, idir), idir = 1, this%dim)
340 end do
341 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
342
344 end subroutine photon_mode_write_info
345
346 !-----------------------------------------------------------------
347 subroutine photon_mode_set_n_electrons(this, qtot)
348 type(photon_mode_t), intent(inout) :: this
349 real(real64), intent(in) :: qtot
350
352
353 this%n_electrons = qtot
354
356 end subroutine photon_mode_set_n_electrons
357
358 !-----------------------------------------------------------------
360 subroutine photon_mode_compute_dipoles(this, mesh)
361 type(photon_mode_t), intent(inout) :: this
362 class(mesh_t), intent(in) :: mesh
363
364 integer :: i, ip
365
367
368 safe_allocate(this%pol_dipole(1:mesh%np, 1:this%nmodes))
369
370 do i = 1, this%nmodes
371 ! Calculate polarization dipole
372 do ip = 1, mesh%np
373 this%pol_dipole(ip, i) = dot_product(this%pol(i, 1:this%dim), mesh%x(ip, 1:this%dim))
374 end do
375 end do
376
378 end subroutine photon_mode_compute_dipoles
379
380 !-----------------------------------------------------------------
381 subroutine photon_mode_add_poisson_terms(this, mesh, rho, pot)
382 type(photon_mode_t), intent(in) :: this
383 type(mesh_t), intent(in) :: mesh
384 real(real64), intent(in) :: rho(:)
385 real(real64), intent(inout) :: pot(:)
386
387 integer :: ip
388 real(real64) :: lx, ld, dipole(mesh%box%dim)
389
391
392 ! Currently this only works with one photon mode
393 assert(this%nmodes == 1)
394 assert(this%n_electrons > m_epsilon)
395
396 call dmf_dipole(mesh, rho, dipole)
397 ld = dot_product(this%pol(1, 1:this%dim), dipole(1:this%dim))*this%lambda(1)
398
399 do ip = 1, mesh%np
400 lx = this%pol_dipole(ip, 1)*this%lambda(1)
401 pot(ip) = pot(ip) - this%omega(1)/sqrt(this%n_electrons)*(mesh%x(ip, this%dim + 1)*ld + lx*dipole(this%dim + 1)) + lx*ld
402 end do
403
405 end subroutine photon_mode_add_poisson_terms
406
407 !-----------------------------------------------------------------
408 subroutine photon_mode_dressed(this)
409 type(photon_mode_t), intent(inout) :: this
410
411 integer :: ia, iap
412 real(real64) :: eps_epsp, wmat_off
413 real(real64), allocatable :: wmat(:,:)
414 real(real64), allocatable :: dressed_omega_sq(:)
415 real(real64), allocatable :: lambda_pol(:,:)
416 real(real64), allocatable :: dressed_lambda_pol(:,:)
418 push_sub(photon_mode_dressed)
419
420 safe_allocate(wmat(1:this%nmodes,1:this%nmodes))
421 safe_allocate(dressed_omega_sq(1:this%nmodes))
422 safe_allocate(lambda_pol(1:this%nmodes,1:this%dim))
423 safe_allocate(dressed_lambda_pol(1:this%dim, 1:this%nmodes))
424
425 ! compute \lambda**2 * polarization array
426 do ia = 1, this%nmodes
427 lambda_pol(ia, 1:this%dim) = this%lambda(ia) * this%pol(ia, 1:this%dim)
428 end do
429
430 ! compute the photon mode matrix, the W matrix
431 do ia = 1, this%nmodes
432 do iap = ia, this%nmodes
433 eps_epsp = dot_product(this%pol(ia, 1:this%dim), this%pol(iap, 1:this%dim))
434 ! TODO: n_electrons should be the total number of electrons including the core levels
435 wmat_off = this%n_electrons * this%lambda(ia) * this%lambda(iap) * eps_epsp
436 wmat(ia,iap) = wmat_off
437 wmat(iap,ia) = wmat_off
438 end do
439 wmat(ia,ia) = this%omega(ia)**2 + wmat(ia,ia)
440 end do
441
442 ! diagonalize the W matrix
443 call lalg_eigensolve(this%nmodes, wmat, dressed_omega_sq)
444
445 do ia = 1, this%nmodes
446 this%dressed_omega(ia) = sqrt(dressed_omega_sq(ia))
447 end do
448
449 ! rotate the ( coupling times polarization )
450 call dgemm('T','N', this%dim, this%nmodes, this%nmodes,&
451 m_one, lambda_pol, this%nmodes, &
452 wmat,this%nmodes, &
453 m_zero,dressed_lambda_pol,this%dim)
454
455 ! compute the rotated coupling and polarization direction
456 do ia = 1, this%nmodes
457 this%dressed_lambda(ia) = norm2(dressed_lambda_pol(1:this%dim, ia))
458 !
459 if (this%dressed_lambda(ia) .lt. m_epsilon) THEN
460 ! to avoid dividing by zero
461 this%dressed_pol(1:this%dim, ia) = m_zero
462 else
463 this%dressed_pol(1:this%dim, ia) = dressed_lambda_pol(1:this%dim, ia)/this%dressed_lambda(ia)
464 end if
465 end do
466
467 ! rotate the bare frequency
468 this%rotated_omega = m_zero
469 ! here we define the rotated bare frequency from the dressed frequency and coupling
470 do ia = 1, this%nmodes
471 this%rotated_omega(ia) = sqrt(this%dressed_omega(ia)**2 - this%n_electrons * this%dressed_lambda(ia)**2)
472 end do
473
474 safe_deallocate_a(wmat)
475 safe_deallocate_a(dressed_omega_sq)
476 safe_deallocate_a(lambda_pol)
477 safe_deallocate_a(dressed_lambda_pol)
478
479 pop_sub(photon_mode_dressed)
480
481 end subroutine photon_mode_dressed
482
483end module photon_mode_oct_m
484
485!! Local Variables:
486!! mode: f90
487!! coding: utf-8
488!! End:
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_one
Definition: global.F90:189
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
This module defines various routines, operating on mesh functions.
subroutine, public dmf_dipole(mesh, ff, dipole, mask)
This routine calculates the dipole of a function ff, for arbitrary dimensions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:920
character(len=512), private msg
Definition: messages.F90:165
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:414
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:434
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:270
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public photon_mode_compute_dipoles(this, mesh)
Computes the polarization dipole.
subroutine, public photon_mode_add_poisson_terms(this, mesh, rho, pot)
subroutine, public photon_mode_write_info(this, iunit, namespace)
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)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
Describes mesh distribution to nodes.
Definition: mesh.F90:186
int true(void)