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_world%is_root()) 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_world%is_root()) 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_new_line()
182 call messages_write('or define the PhotonModes variable!')
183 call messages_fatal(namespace=namespace)
184 end if
185 end if
186
187 !%Variable PhotonModes
188 !%Type block
189 !%Section Hamiltonian::XC
190 !%Description
191 !% Each line of the block should specify one photon mode. The syntax is the following:
192 !%
193 !% %PhotonModes
194 !% omega1 | lambda1| PolX1 | PolY1 | PolZ1
195 !% ...
196 !% %
197 !%
198 !% The first column is the mode frequency, in units of energy.
199 !% The second column is the coupling strength, in units of energy.
200 !% The remaining columns specify the polarization direction of the mode.
201 !% If the polarization vector should be normalized to one. If that is not the case
202 !% the code will normalize it.
203 !%End
204
205 if (.not. file_exists) then
206 if (parse_block(namespace, 'PhotonModes', blk) == 0) then
207
208 this%nmodes = parse_block_n(blk)
209 safe_allocate(this%omega(1:this%nmodes))
210 safe_allocate(this%lambda(1:this%nmodes))
211 safe_allocate(this%pol(1:this%nmodes, 1:this%dim))
212
213 this%pol = m_zero
214
215 do ii = 1, this%nmodes
216 ncols = parse_block_cols(blk, ii-1)
217
218 ! Read line
219 call parse_block_float(blk, ii-1, 0, this%omega(ii), units_inp%energy) ! frequency
220 call parse_block_float(blk, ii-1, 1, this%lambda(ii), units_inp%energy) ! coupling strength
221 do idir = 1, this%dim
222 call parse_block_float(blk, ii-1, idir + 1, this%pol(ii, idir)) ! polarization vector components
223 end do
224
225 if (ncols > this%dim + 2) then
226 this%has_q0_p0 = .true.
227 if (.not. associated(this%pt_coord_q0)) then
228 safe_allocate(this%pt_coord_q0(1:this%nmodes))
229 end if
230 if (.not. associated(this%pt_momen_p0)) then
231 safe_allocate(this%pt_momen_p0(1:this%nmodes))
232 end if
233 call parse_block_float(blk, ii-1, this%dim+2, this%pt_coord_q0(ii)) !row, column
234 call parse_block_float(blk, ii-1, this%dim+3, this%pt_momen_p0(ii)) !row, column
235 end if
236 ! Sanity check
237 if ((ncols /= this%dim + 2) .and. (ncols /= this%dim + 4)) then
238 call messages_input_error(namespace, 'PhotonModes', 'Incorrect number of columns')
239 end if
240
241 ! Normalize polarization vector
242 this%pol(ii,:) = this%pol(ii,:)/norm2(this%pol(ii,:))
243
244 end do
245 call parse_block_end(blk)
246 else
247 call messages_write('You need to specify the photon modes!')
248 call messages_fatal(namespace=namespace)
249 end if
250 end if
251
252 !%Variable TDPhotonicTimeScale
253 !%Type float
254 !%Default 1.0
255 !%Section Time-Dependent::Propagation
256 !%Description
257 !% This variable defines the factor between the timescale of photonic
258 !% and electronic movement.
259 !% for more details see the documentation of TDIonicTimeScale
260 !% If you also use TDIonicTimeScale, we advise to set
261 !% TDPhotonicTimeScale = TDIonicTimeScale, in the case the
262 !% photon frequency is in a vibrational energy range.
263 !% Important: The electronic time step will be the value of
264 !% <tt>TDTimeStep</tt> divided by this variable, so if you have determined an
265 !% optimal electronic time step (that we can call <i>dte</i>), it is
266 !% recommended that you define your time step as:
267 !%
268 !% <tt>TDTimeStep</tt> = <i>dte</i> * <tt>TDPhotonicTimeScale</tt>
269 !%
270 !% so you will always use the optimal electronic time step
271 !% (<a href=http://arxiv.org/abs/0710.3321>more details</a>).
272 !%End
273
274 call parse_variable(namespace, 'TDPhotonicTimeScale', 1.0_real64, this%mu)
275
276 this%ex = m_zero
277 safe_allocate(this%number(1:this%nmodes))
278 this%number = m_zero
279
280 ! photon-free approach, which requires the dressed photon modes
281 this%use_photon_exchange = .false.
282
283 this%use_photon_exchange = optional_default(photon_free, .false.)
284 if (this%use_photon_exchange) then
285 safe_allocate(this%dressed_omega(1:this%nmodes))
286 safe_allocate(this%dressed_lambda(1:this%nmodes))
287 safe_allocate(this%dressed_pol(1:this%dim, 1:this%nmodes))
288 safe_allocate(this%rotated_omega(1:this%nmodes))
289 end if
290
291 pop_sub(photon_mode_init)
292 end subroutine photon_mode_init
293
294 ! ---------------------------------------------------------
295 subroutine photon_mode_end(this)
296 type(photon_mode_t), intent(inout) :: this
297
298 push_sub(photon_mode_end)
299
300 safe_deallocate_a(this%correlator)
301
302 safe_deallocate_a(this%omega)
303 safe_deallocate_a(this%lambda)
304 safe_deallocate_a(this%number)
305
306 safe_deallocate_a(this%pol)
307 safe_deallocate_a(this%pol_dipole)
308
309 if (this%has_q0_p0) then
310 safe_deallocate_p(this%pt_coord_q0)
311 safe_deallocate_p(this%pt_momen_p0)
312 end if
313
314 if (this%use_photon_exchange) then
315 safe_deallocate_a(this%dressed_omega)
316 safe_deallocate_a(this%dressed_lambda)
317 safe_deallocate_a(this%dressed_pol)
318 safe_deallocate_a(this%rotated_omega)
319 end if
320
321 pop_sub(photon_mode_end)
322 end subroutine photon_mode_end
323
324 !-----------------------------------------------------------------
325 subroutine photon_mode_write_info(this, iunit, namespace)
326 type(photon_mode_t), intent(in) :: this
327 integer, optional, intent(in) :: iunit
328 type(namespace_t), optional, intent(in) :: namespace
329
330 integer :: im, idir
331
332 push_sub(photon_mode_write_info)
333
334 call messages_print_with_emphasis(msg="Photon Modes", iunit=iunit, namespace=namespace)
335 write(iunit, '(6x,a,10x,a,3x)', advance='no') 'Omega', 'Lambda'
336 write(iunit, '(1x,a,i1,a)') ('Pol.(', idir, ')', idir = 1, this%dim)
337 do im = 1, this%nmodes
338 write(iunit, '(1x,f14.12)', advance='no') this%omega(im)
339 write(iunit, '(1x,f14.12)', advance='no') this%lambda(im)
340 write(iunit, '(2X,f5.3,1X)') (this%pol(im, idir), idir = 1, this%dim)
341 end do
342 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
343
345 end subroutine photon_mode_write_info
346
347 !-----------------------------------------------------------------
348 subroutine photon_mode_set_n_electrons(this, qtot)
349 type(photon_mode_t), intent(inout) :: this
350 real(real64), intent(in) :: qtot
351
353
354 this%n_electrons = qtot
355
357 end subroutine photon_mode_set_n_electrons
358
359 !-----------------------------------------------------------------
361 subroutine photon_mode_compute_dipoles(this, mesh)
362 type(photon_mode_t), intent(inout) :: this
363 class(mesh_t), intent(in) :: mesh
364
365 integer :: i, ip
366
368
369 safe_allocate(this%pol_dipole(1:mesh%np, 1:this%nmodes))
370
371 do i = 1, this%nmodes
372 ! Calculate polarization dipole
373 do ip = 1, mesh%np
374 this%pol_dipole(ip, i) = dot_product(this%pol(i, 1:this%dim), mesh%x(ip, 1:this%dim))
375 end do
376 end do
377
379 end subroutine photon_mode_compute_dipoles
380
381 !-----------------------------------------------------------------
382 subroutine photon_mode_add_poisson_terms(this, mesh, rho, pot)
383 type(photon_mode_t), intent(in) :: this
384 type(mesh_t), intent(in) :: mesh
385 real(real64), intent(in) :: rho(:)
386 real(real64), intent(inout) :: pot(:)
387
388 integer :: ip
389 real(real64) :: lx, ld, dipole(mesh%box%dim)
392
393 ! Currently this only works with one photon mode
394 assert(this%nmodes == 1)
395 assert(this%n_electrons > m_epsilon)
396
397 call dmf_dipole(mesh, rho, dipole)
398 ld = dot_product(this%pol(1, 1:this%dim), dipole(1:this%dim))*this%lambda(1)
399
400 do ip = 1, mesh%np
401 lx = this%pol_dipole(ip, 1)*this%lambda(1)
402 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
403 end do
404
406 end subroutine photon_mode_add_poisson_terms
407
408 !-----------------------------------------------------------------
409 subroutine photon_mode_dressed(this)
410 type(photon_mode_t), intent(inout) :: this
411
412 integer :: ia, iap
413 real(real64) :: eps_epsp, wmat_off
414 real(real64), allocatable :: wmat(:,:)
415 real(real64), allocatable :: dressed_omega_sq(:)
416 real(real64), allocatable :: lambda_pol(:,:)
417 real(real64), allocatable :: dressed_lambda_pol(:,:)
418
419 push_sub(photon_mode_dressed)
421 safe_allocate(wmat(1:this%nmodes,1:this%nmodes))
422 safe_allocate(dressed_omega_sq(1:this%nmodes))
423 safe_allocate(lambda_pol(1:this%nmodes,1:this%dim))
424 safe_allocate(dressed_lambda_pol(1:this%dim, 1:this%nmodes))
425
426 ! compute \lambda**2 * polarization array
427 do ia = 1, this%nmodes
428 lambda_pol(ia, 1:this%dim) = this%lambda(ia) * this%pol(ia, 1:this%dim)
429 end do
430
431 ! compute the photon mode matrix, the W matrix
432 do ia = 1, this%nmodes
433 do iap = ia, this%nmodes
434 eps_epsp = dot_product(this%pol(ia, 1:this%dim), this%pol(iap, 1:this%dim))
435 ! TODO: n_electrons should be the total number of electrons including the core levels
436 wmat_off = this%n_electrons * this%lambda(ia) * this%lambda(iap) * eps_epsp
437 wmat(ia,iap) = wmat_off
438 wmat(iap,ia) = wmat_off
439 end do
440 wmat(ia,ia) = this%omega(ia)**2 + wmat(ia,ia)
441 end do
442
443 ! diagonalize the W matrix
444 call lalg_eigensolve(this%nmodes, wmat, dressed_omega_sq)
445
446 do ia = 1, this%nmodes
447 this%dressed_omega(ia) = sqrt(dressed_omega_sq(ia))
448 end do
449
450 ! rotate the ( coupling times polarization )
451 call dgemm('T','N', this%dim, this%nmodes, this%nmodes,&
452 m_one, lambda_pol, this%nmodes, &
453 wmat,this%nmodes, &
454 m_zero,dressed_lambda_pol,this%dim)
455
456 ! compute the rotated coupling and polarization direction
457 do ia = 1, this%nmodes
458 this%dressed_lambda(ia) = norm2(dressed_lambda_pol(1:this%dim, ia))
459 !
460 if (this%dressed_lambda(ia) .lt. m_epsilon) THEN
461 ! to avoid dividing by zero
462 this%dressed_pol(1:this%dim, ia) = m_zero
463 else
464 this%dressed_pol(1:this%dim, ia) = dressed_lambda_pol(1:this%dim, ia)/this%dressed_lambda(ia)
465 end if
466 end do
467
468 ! rotate the bare frequency
469 this%rotated_omega = m_zero
470 ! here we define the rotated bare frequency from the dressed frequency and coupling
471 do ia = 1, this%nmodes
472 this%rotated_omega(ia) = sqrt(this%dressed_omega(ia)**2 - this%n_electrons * this%dressed_lambda(ia)**2)
473 end do
474
475 safe_deallocate_a(wmat)
476 safe_deallocate_a(dressed_omega_sq)
477 safe_deallocate_a(lambda_pol)
478 safe_deallocate_a(dressed_lambda_pol)
479
480 pop_sub(photon_mode_dressed)
481
482 end subroutine photon_mode_dressed
483
484end module photon_mode_oct_m
485
486!! Local Variables:
487!! mode: f90
488!! coding: utf-8
489!! End:
real(real64), parameter, public m_zero
Definition: global.F90:190
real(real64), parameter, public m_epsilon
Definition: global.F90:206
real(real64), parameter, public m_one
Definition: global.F90:191
Definition: io.F90:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:466
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:401
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:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:904
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_new_line()
Definition: messages.F90:1118
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:416
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:697
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:600
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:505
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:621
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:134
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:188
int true(void)