Octopus
lasers.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 N. Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22module lasers_oct_m
23 use clock_oct_m
24 use debug_oct_m
26 use global_oct_m
30 use, intrinsic :: iso_fortran_env
34 use math_oct_m
35 use mpi_oct_m
36 use mesh_oct_m
39 use parser_oct_m
42 use space_oct_m
45 use unit_oct_m
48
49 implicit none
50
51 private
52 public :: &
53 lasers_t, &
58 laser_t, &
69 laser_kind, &
80
81
82 ! TODO: (Micael, Alex) Issue 836. Remove the following paramaters and use the
83 ! corresponding quantities defined in the multisystem (E_FIELD, B_FIELD, etc)
84 integer, public, parameter :: &
85 E_FIELD_NONE = 0, &
86 e_field_electric = 1, &
87 e_field_magnetic = 2, &
90
91 type laser_t
92 private
93 integer :: field = e_field_none
94 complex(real64), allocatable :: pol(:)
95 real(real64), allocatable :: prop(:)
96 type(tdf_t) :: f
97 type(tdf_t) :: phi
98 real(real64) :: omega = m_zero
99
100 real(real64), allocatable :: v(:)
101 real(real64), allocatable :: a(:, :)
102 end type laser_t
103
104 type, extends(interaction_partner_t) :: lasers_t
105 private
106
107 integer, public :: no_lasers
108 type(laser_t), allocatable, public :: lasers(:)
109 character(len=200) :: scalar_pot_expression
110 character(len=200) :: envelope_expression
111 character(len=200) :: phase_expression
112
113 real(real64), allocatable :: e(:)
114 real(real64), allocatable :: b(:)
115 real(real64), allocatable :: integrated_nondipole_afield(:)
116 real(real64) :: nd_integration_time
117 real(real64) :: nd_integration_step
118 contains
119 procedure :: init_interaction_as_partner => lasers_init_interaction_as_partner
120 procedure :: update_quantity => lasers_update_quantity
121 procedure :: copy_quantities_to_interaction => lasers_copy_quantities_to_interaction
122 final :: lasers_finalize
123 end type
124
125 interface lasers_t
126 module procedure lasers_constructor
127 end interface lasers_t
128
129
130contains
131
132 function lasers_constructor(namespace) result(this)
133 class(lasers_t), pointer :: this
134 type(namespace_t), intent(in) :: namespace
135
136 push_sub(lasers_constructor)
137
138 safe_allocate(this)
139
140 this%namespace = namespace_t("Lasers", parent=namespace)
141
142 safe_allocate(this%e(1:3))
143 safe_allocate(this%b(1:3))
144 this%e = m_zero
145 this%b = m_zero
146
147 pop_sub(lasers_constructor)
148 end function lasers_constructor
149
150 subroutine lasers_parse_external_fields(this)
151 class(lasers_t), intent(inout) :: this
152
153 type(block_t) :: blk
154 integer :: il, jj, ierr, k
155 real(real64) :: omega0
156 complex(real64) :: cprop(3)
157
159
160 call messages_obsolete_variable(this%namespace, "TDLasers", "TDExternalFields")
161 !%Variable TDExternalFields
162 !%Type block
163 !%Section Time-Dependent
164 !%Description
165 !% The block <tt>TDExternalFields</tt> describes the type and shape of time-dependent
166 !% external perturbations that are applied to the system, in the form
167 !% <math>f(x,y,z) \cos(\omega t + \phi (t)) g(t)</math>, where <math>f(x,y,z)</math> is defined by
168 !% by a field type and polarization or a scalar potential, as below; <math>\omega</math>
169 !% is defined by <tt>omega</tt>; <math>g(t)</math> is defined by
170 !% <tt>envelope_function_name</tt>; and <math>\phi(t)</math> is the (time-dependent) phase from <tt>phase</tt>.
171 !%
172 !% These perturbations are only applied for time-dependent runs. If
173 !% you want the value of the perturbation at time zero to be
174 !% applied for time-independent runs, use <tt>TimeZero = yes</tt>.
175 !%
176 !% Each line of the block describes an external field; this way you can actually have more
177 !% than one laser (<i>e.g.</i> a "pump" and a "probe").
178 !%
179 !% There are two ways to specify <math>f(x,y,z)</math> but both use the same <tt>omega | envelope_function_name [| phase]</tt>
180 !% for the time-dependence.
181 !% The float <tt>omega</tt> will be the carrier frequency of the
182 !% pulse (in energy units). The envelope of the field is a time-dependent function whose definition
183 !% must be given in a <tt>TDFunctions</tt> block. <tt>envelope_function_name</tt> is a string (and therefore
184 !% it must be surrounded by quotation marks) that must match one of the function names
185 !% given in the first column of the <tt>TDFunctions</tt> block.
186 !% <tt>phase</tt> is optional and is taken to be zero if not provided, and is also a string specifying
187 !% a time-dependent function.
188 !%
189 !% (A) type = <tt>electric field, magnetic field, vector_potential</tt>
190 !%
191 !% For these cases, the syntax is:
192 !%
193 !% <tt>%TDExternalFields
194 !% <br>&nbsp;&nbsp; type | nx | ny | nz | omega | envelope_function_name | phase (| px | py | pz )
195 !% <br>%</tt>
196 !%
197 !% The <tt>vector_potential</tt> option (constant in space) permits us to describe
198 !% an electric perturbation in the velocity gauge.
199 !% The three (possibly complex) numbers (<tt>nx</tt>, <tt>ny</tt>, <tt>nz</tt>) mark the polarization
200 !% direction of the field.
201 !% By default, (<tt>nx</tt>, <tt>ny</tt>, <tt>nz</tt>) are defined in Cartesian space.
202 !% However, it is possible for solids to define them using the Miller indices.
203 !% This can be achieved by defining the block <tt>MillerIndicesBasis</tt>. It is also possible to activate
204 !% the nondipole SFA correction by inserting the laser propagation direction px | py | pz as prescribed in
205 !% Phys. Rev. A 101, 043408 (2020).
206 !%
207 !% (B) type = <tt>scalar_potential</tt>
208 !%
209 !% <tt>%TDExternalFields
210 !% <br>&nbsp;&nbsp; scalar_potential | "spatial_expression" | omega | envelope_function_name | phase
211 !% <br>%</tt>
212 !%
213 !% The scalar potential is any expression of the spatial coordinates given by the string
214 !% "spatial_expression", allowing a field beyond the dipole approximation.
215 !%
216 !% For DFTB runs, only fields of type type = <tt>electric field</tt> are allowed for the moment, and the
217 !% <tt>type</tt> keyword is omitted.
218 !%
219 !% A NOTE ON UNITS:
220 !%
221 !% It is very common to describe the strength of a laser field by its intensity, rather
222 !% than using the electric-field amplitude. In atomic units (or, more precisely, in any
223 !% Gaussian system of units), the relationship between instantaneous electric field
224 !% and intensity is:
225 !% <math> I(t) = \frac{c}{8\pi} E^2(t) </math>.
226 !%
227 !% It is common to read intensities in W/cm<math>^2</math>. The dimensions of intensities are
228 !% [W]/(L<math>^2</math>T), where [W] are the dimensions of energy. The relevant conversion factors
229 !% are:
230 !%
231 !% Hartree / (<math>a_0^2</math> atomic_time) = <math>6.4364086 \times 10^{15} \mathrm{W/cm}^2</math>
232 !%
233 !% eV / ( &Aring;<math>^2 (\hbar</math>/eV) ) = <math>2.4341348 \times 10^{12} \mathrm{W/cm}^2</math>
234 !%
235 !% If, in atomic units, we set the electric-field amplitude to <math>E_0</math>,
236 !% then the intensity is:
237 !%
238 !% <math> I_0 = 3.51 \times 10^{16} \mathrm{W/cm}^2 (E_0^2) </math>
239 !%
240 !% If, working with <tt>Units = ev_angstrom</tt>, we set <math>E_0</math>, then the intensity is:
241 !%
242 !% <math> I_0 = 1.327 \times 10^{13} (E_0^2) \mathrm{W/cm}^2 </math>
243 !%
244 !%Option electric_field 1
245 !% The external field is an electric field, the usual case when we want to describe a
246 !% laser in the length gauge.
247 !%Option magnetic_field 2
248 !% The external field is a (homogeneous) time-dependent magnetic field.
249 !%Option vector_potential 3
250 !% The external field is a time-dependent homogeneous vector potential, which may describe
251 !% a laser field in the velocity gauge.
252 !%Option scalar_potential 4
253 !% The external field is an arbitrary scalar potential, which may describe an
254 !% inhomogeneous electrical field.
255 !%End
256
257 this%no_lasers = 0
258 if (parse_block(this%namespace, 'TDExternalFields', blk) == 0) then
259 this%no_lasers = parse_block_n(blk)
260 safe_allocate(this%lasers(1:this%no_lasers))
261
262 do il = 1, this%no_lasers
263 safe_allocate(this%lasers(il)%pol(1:3))
264 this%lasers(il)%pol = m_z0
265
266 call parse_block_integer(blk, il-1, 0, this%lasers(il)%field)
267
268 select case (this%lasers(il)%field)
270 call parse_block_string(blk, il-1, 1, this%scalar_pot_expression)
271 jj = 1
272 this%lasers(il)%pol = m_z1
273 case default
274 call parse_block_cmplx(blk, il-1, 1, this%lasers(il)%pol(1))
275 call parse_block_cmplx(blk, il-1, 2, this%lasers(il)%pol(2))
276 call parse_block_cmplx(blk, il-1, 3, this%lasers(il)%pol(3))
277 jj = 3
278 end select
279
280 call parse_block_float(blk, il-1, jj+1, omega0)
281 omega0 = units_to_atomic(units_inp%energy, omega0)
282
283 this%lasers(il)%omega = omega0
284
285 call parse_block_string(blk, il-1, jj+2, this%envelope_expression)
286 call tdf_read(this%lasers(il)%f, this%namespace, trim(this%envelope_expression), ierr)
287
288 ! Check if there is a phase.
289 if (parse_block_cols(blk, il-1) > jj+3) then
290 call parse_block_string(blk, il-1, jj+3, this%phase_expression)
291 call tdf_read(this%lasers(il)%phi, this%namespace, trim(this%phase_expression), ierr)
292 if (ierr /= 0) then
293 write(message(1),'(3A)') 'Error in the "', trim(this%envelope_expression), &
294 '" field defined in the TDExternalFields block:'
295 write(message(2),'(3A)') 'Time-dependent phase function "', trim(this%phase_expression), &
296 '" not found.'
297 call messages_warning(2, namespace=this%namespace)
298 end if
299 else
300 call tdf_init(this%lasers(il)%phi)
301 end if
302
303 ! Check if there is nondipoole fields to activate the nondipole SFA formalism of Phys. Rev. A 101, 043408 (2020)
304 if (parse_block_cols(blk, il-1) > jj+4) then
305 safe_allocate(this%lasers(il)%prop(1:size(cprop)))
306 do k = 1, size(cprop)
307 call parse_block_cmplx(blk, il-1, jj+size(cprop)+k, cprop(k))
308
309 end do
310 if (any(abs(aimag(cprop)) > m_epsilon)) then
311 write(message(1),'(3A)') 'Error in the "', trim(this%envelope_expression), &
312 '" field defined in the TDExternalFields block:'
313 write(message(2),'(A)') 'Propagation direction cannot be complex.'
314 call messages_fatal(2, namespace=this%namespace)
315 end if
316 this%lasers(il)%prop(:) = real(cprop(:), real64)
317 if (.not.allocated(this%integrated_nondipole_afield)) then
318 safe_allocate(this%integrated_nondipole_afield(1:size(cprop)))
319 this%integrated_nondipole_afield(1:3)=m_zero
320 this%nd_integration_time=m_zero
321 !%Variable NDSFATimeIntegrationStep
322 !%Type float
323 !%Default 0.01
324 !%Section Hamiltonian
325 !%Description
326 !% Timestep for the integration of the NDSFA first order
327 !% response numerically as prescribed in Phys. Rev. A 101, 043408 (2020).
328 !% Typically needs to be an order of magniude less than
329 !% the TDtimestep. However, for writing "laser" the TDtimestep
330 !% is used instead.
331 !%End
332 call parse_variable(this%namespace, 'NDSFATimeIntegrationStep', 0.01_real64, this%nd_integration_step)
333
334 if (is_close(this%nd_integration_step, 0.01_real64)) then
335 message(1) = "The default timestep of 0.01 is utilized for the nondipole SFA integration."
336 message(2) = "Be aware that this should be at least an order of magniude less than the TDtimestep"
337 call messages_info(2, namespace=this%namespace)
338 end if
339 end if
340 end if
341
342 end do
343
344 call parse_block_end(blk)
345 end if
346
348 end subroutine lasers_parse_external_fields
349
350 subroutine lasers_generate_potentials(this, mesh, space, latt)
351 class(lasers_t), intent(inout) :: this
352 class(mesh_t), intent(in) :: mesh
353 class(space_t), intent(in) :: space
354 type(lattice_vectors_t), intent(in) :: latt
355
356 type(block_t) :: blk2
357 integer :: il, ip, idir, idir2
358 real(real64) :: rr, pot_re, pot_im, xx(3)
359 real(real64) :: miller(space%dim,space%dim), miller_red(space%dim,space%dim)
360
362
363 do il = 1, this%no_lasers
364 ! For periodic systems, we might want to define the polarization direction in
365 ! terms of Miller indices
366 ! In this case, the use must provide the coordinates of the X, Y, and Z high symmetry points
367 ! such that the code applying a laser along say the [100] direction converts it to
368 ! the proper coordinates in Cartesian space.
369 ! This is usefull for arbitrarily rotated crystals.
370
371 !%Variable MillerIndicesBasis
372 !%Type block
373 !%Section Time-Dependent
374 !%Description
375 !% When this block is given, the polarisation of the TDExternalFields is
376 !% understood to be defined in terms of Miller indices.
377 !% This block define the corresponding basis, by defining the reduced coordinates
378 !% of the X, Y, and Z high symmetry points, such that the code can do the corresponding
379 !% transformation.
380 !%
381 !% For example, in an FCC crystal with the conventional primitive cell,
382 !% the following input allows to define the polarization in terms
383 !% message(1) = "The lasers break (at least) one of the symmetries used to reduce the k-points ."
384 !%
385 !% <tt>%MillerIndicesBasis
386 !% <br> 0.0 | 0.5 | 0.5
387 !% <br> 0.5 | 0.0 | 0.5
388 !% <br> 0.5 | 0.5 | 0.0
389 !% <br>%</tt>
390 !%
391 !% Indeed, in this case, the reciprocal lattice vectors are (-1, 1, 1), (1, -1, 1),
392 !% and (1, 1, -1) in units of 2*pi/a.
393 !% This directly gives that the [100] direction correspond to the x direction, [111]
394 !% gives the vector (1,1,1), etc.
395 !%
396 !%End
397 if (parse_block(this%namespace, 'MillerIndicesBasis', blk2) == 0) then
398 if(.not. space%is_periodic()) then
399 write(message(1),'(a)') 'MillerIndicesBasis can only be used for periodic systems.'
400 call messages_fatal(1, namespace=this%namespace)
401 end if
402
403 do idir = 1, space%dim
404 do idir2 = 1, space%dim
405 call parse_block_float(blk2, idir-1, idir2-1, miller_red(idir2, idir))
406 end do
407 call kpoints_to_absolute(latt, miller_red(:, idir), miller(:, idir))
408 end do
409
410
411 this%lasers(il)%pol = matmul(miller, this%lasers(il)%pol)
412 call parse_block_end(blk2)
413 end if
414
415 this%lasers(il)%pol(:) = this%lasers(il)%pol(:)/sqrt(sum(abs(this%lasers(il)%pol(:))**2))
416
417 select case (this%lasers(il)%field)
419 safe_allocate(this%lasers(il)%v(1:mesh%np_part))
420 this%lasers(il)%v = m_zero
421 do ip = 1, mesh%np
422 call mesh_r(mesh, ip, rr, coords = xx(1:space%dim))
423 xx(space%dim+1:3) = m_zero
424 call parse_expression(pot_re, pot_im, 3, xx, rr, m_zero, trim(this%scalar_pot_expression))
425 this%lasers(il)%v(ip) = pot_re
426 end do
427
428 case (e_field_magnetic)
429 ! \warning: note that for the moment we are ignoring the possibility of a complex
430 ! polarizability vector for the td magnetic-field case.
431 safe_allocate(this%lasers(il)%a(1:mesh%np_part, 1:3))
432 this%lasers(il)%a = m_zero
433 do ip = 1, mesh%np
434 xx(1:space%dim) = mesh%x(ip, :)
435 xx(space%dim+1:3) = m_zero
436 ! Compute the vector potential from a uniform B field.
437 ! The sign is determined by the relation $\vec{B} = \nabla \times \vec{A}$.
438 ! This leads to $\vec{A} = -\frac{1}{2}\vec{r}\times\vec{B}$.
439 select case (space%dim)
440 case (2)
441 this%lasers(il)%a(ip, :) = (/xx(2), -xx(1)/) * sign(m_one, real(this%lasers(il)%pol(3)))
442 case (3)
443 this%lasers(il)%a(ip, :) = (/ xx(2)*real(this%lasers(il)%pol(3)) - xx(3)*real(this%lasers(il)%pol(2)), &
444 xx(3)*real(this%lasers(il)%pol(1)) - xx(1)*real(this%lasers(il)%pol(3)), &
445 xx(1)*real(this%lasers(il)%pol(2)) - xx(2)*real(this%lasers(il)%pol(1)) /)
446 case default
447 message(1) = "Magnetic fields only allowed in 2 or 3D."
448 call messages_fatal(1, namespace=this%namespace)
449 end select
450 end do
451 this%lasers(il)%a = -m_half * this%lasers(il)%a
452
453 end select
454
455 end do
456
458 end subroutine lasers_generate_potentials
459
460
461 ! ---------------------------------------------------------
462 subroutine lasers_check_symmetries(this, kpoints)
463 type(lasers_t), intent(in) :: this
464 type(kpoints_t), intent(in) :: kpoints
465
466 integer :: iop, il
467
469
470 if (kpoints%use_symmetries) then
471 do iop = 1, symmetries_number(kpoints%symm)
472 if (iop == symmetries_identity_index(kpoints%symm)) cycle
473 do il = 1, this%no_lasers
474 if (.not. symm_op_invariant_cart(kpoints%symm%ops(iop), this%lasers(il)%pol(:), symprec)) then
475 message(1) = "The lasers break (at least) one of the symmetries used to reduce the k-points ."
476 message(2) = "Set SymmetryBreakDir accordingly to your laser fields."
477 call messages_fatal(2, namespace=this%namespace)
478 end if
479 end do
480 end do
481 end if
482
484 end subroutine lasers_check_symmetries
485
486 ! ---------------------------------------------------------
487 subroutine lasers_finalize(this)
488 type(lasers_t), intent(inout) :: this
489
490 push_sub(lasers_finalize)
491
492 call lasers_deallocate(this)
493
494 pop_sub(lasers_finalize)
495 end subroutine lasers_finalize
496
497 ! ---------------------------------------------------------
498 subroutine lasers_deallocate(this)
499 class(lasers_t), intent(inout) :: this
500
501 integer :: il
502
503 push_sub(lasers_deallocate)
504
505 safe_deallocate_a(this%e)
506 safe_deallocate_a(this%b)
507 safe_deallocate_a(this%integrated_nondipole_afield)
508
509 do il = 1, this%no_lasers
510 safe_deallocate_a(this%lasers(il)%pol)
511 call tdf_end(this%lasers(il)%f)
512 call tdf_end(this%lasers(il)%phi)
513 select case (this%lasers(il)%field)
515 safe_deallocate_a(this%lasers(il)%v)
516 case (e_field_magnetic)
517 safe_deallocate_a(this%lasers(il)%a)
518 end select
519 safe_deallocate_a(this%lasers(il)%prop)
520 end do
521 safe_deallocate_a(this%lasers)
522
523 pop_sub(lasers_deallocate)
524 end subroutine lasers_deallocate
525
526
527 ! ---------------------------------------------------------
528 real(real64) function laser_carrier_frequency(laser) result(w0)
529 type(laser_t), intent(in) :: laser
530
532
533 w0 = laser%omega
534
536 end function laser_carrier_frequency
537
538 ! ---------------------------------------------------------
539 subroutine lasers_init_interaction_as_partner(partner, interaction)
540 class(lasers_t), intent(in) :: partner
541 class(interaction_surrogate_t), intent(inout) :: interaction
542
544
545 select type (interaction)
546 type is (lorentz_force_t)
547 ! Nothing to be initialized for the Lorentz force.
548 class default
549 message(1) = "Unsupported interaction."
550 call messages_fatal(1, namespace=partner%namespace)
551 end select
552
556 ! ---------------------------------------------------------
557 subroutine lasers_update_quantity(this, label)
558 class(lasers_t), intent(inout) :: this
559 character(len=*), intent(in) :: label
560
561 type(quantity_t), pointer :: quantity
562 integer :: il
563
564 push_sub(lasers_update_quantity)
565
566 if (any(laser_kind(this%lasers) == e_field_vector_potential) .or. &
567 any(laser_kind(this%lasers) == e_field_scalar_potential)) then
568 call messages_not_implemented("Laser vector potentials and scalar potentials in multi-system framework", &
569 namespace=this%namespace)
570 end if
571
572 quantity => this%quantities%get(label)
573 select case (label)
574 case ("E field")
575 this%e = m_zero
576 do il = 1, this%no_lasers
577 if (laser_kind(this%lasers(il)) == e_field_electric) then
578 call laser_field(this%lasers(il), this%e, quantity%iteration%value())
579 end if
580 end do
581
582 case ("B field")
583 this%b = m_zero
584 do il = 1, this%no_lasers
585 if (laser_kind(this%lasers(il)) == e_field_magnetic) then
586 call laser_field(this%lasers(il), this%b, quantity%iteration%value())
587 end if
588 end do
589
590 case default
591 message(1) = "Incompatible quantity."
592 call messages_fatal(1, namespace=this%namespace)
593 end select
594
596 end subroutine lasers_update_quantity
597
598 ! ---------------------------------------------------------
599 subroutine lasers_copy_quantities_to_interaction(partner, interaction)
600 class(lasers_t), intent(inout) :: partner
601 class(interaction_surrogate_t), intent(inout) :: interaction
602
603 integer :: ip
604
606
607 select type (interaction)
608 type is (lorentz_force_t)
609 do ip = 1, interaction%system_np
610 interaction%partner_e_field(:, ip) = partner%e
611 interaction%partner_b_field(:, ip) = partner%b
612 end do
613 class default
614 message(1) = "Unsupported interaction."
615 call messages_fatal(1, namespace=partner%namespace)
616 end select
617
620
621 ! ---------------------------------------------------------
622 integer pure elemental function laser_kind(laser)
623 type(laser_t), intent(in) :: laser
624
625 ! no push_sub allowed in pure function
626 laser_kind = laser%field
627
628 end function laser_kind
629 ! ---------------------------------------------------------
630
631
632 ! ---------------------------------------------------------
633 function laser_polarization(laser) result(pol)
634 type(laser_t), intent(in) :: laser
635 complex(real64) :: pol(3)
636
637 push_sub(laser_polarization)
638
639 pol = laser%pol
640
641 pop_sub(laser_polarization)
642 end function laser_polarization
643 ! ---------------------------------------------------------
644
646 function lasers_with_nondipole_field(lasers) result(isnondipole)
647 type(lasers_t), intent(in) :: lasers
648 logical :: isnondipole
649
651 isnondipole = .false.
652 if(allocated(lasers%integrated_nondipole_afield)) isnondipole = .true.
654 end function lasers_with_nondipole_field
655 ! ---------------------------------------------------------
656
657
659 subroutine lasers_set_nondipole_parameters(this, ndfield, nd_integration_time)
660 type(lasers_t), intent(inout) :: this
661 real(real64), intent(in) :: ndfield(:)
662 real(real64), intent(in) :: nd_integration_time
663 integer :: dim
664 dim = size(ndfield)
665
667 this%integrated_nondipole_afield(1:dim) = ndfield(1:dim)
668 this%nd_integration_time = nd_integration_time
671
672
673 ! ---------------------------------------------------------
674 subroutine laser_get_f(laser, ff)
675 type(laser_t), intent(in) :: laser
676 type(tdf_t), intent(inout) :: ff
677
678 push_sub(laser_get_f)
679 call tdf_copy(ff, laser%f)
680
681 pop_sub(laser_get_f)
682 end subroutine laser_get_f
683 ! ---------------------------------------------------------
684
685
686 ! ---------------------------------------------------------
687 subroutine laser_set_f(laser, ff)
688 type(laser_t), intent(inout) :: laser
689 type(tdf_t), intent(inout) :: ff
690
691 push_sub(laser_set_f)
693 call tdf_end(laser%f)
694 call tdf_copy(laser%f, ff)
695
696 pop_sub(laser_set_f)
697 end subroutine laser_set_f
698 ! ---------------------------------------------------------
699
700
701 ! ---------------------------------------------------------
702 subroutine laser_get_phi(laser, phi)
703 type(laser_t), intent(in) :: laser
704 type(tdf_t), intent(inout) :: phi
705
706 push_sub(laser_get_phi)
707 call tdf_copy(phi, laser%phi)
708
709 pop_sub(laser_get_phi)
710 end subroutine laser_get_phi
711 ! ---------------------------------------------------------
712
713
714 ! ---------------------------------------------------------
715 subroutine laser_set_phi(laser, phi)
716 type(laser_t), intent(inout) :: laser
717 type(tdf_t), intent(inout) :: phi
718
719 push_sub(laser_set_phi)
720
721 call tdf_end(laser%phi)
722 call tdf_copy(laser%phi, phi)
723
724 pop_sub(laser_set_phi)
725 end subroutine laser_set_phi
726 ! ---------------------------------------------------------
727
728
729 ! ---------------------------------------------------------
730 subroutine laser_set_empty_phi(laser)
731 type(laser_t), intent(inout) :: laser
732
733 push_sub(laser_set_empty_phi)
734
735 call tdf_init(laser%phi)
736
737 pop_sub(laser_set_empty_phi)
738 end subroutine laser_set_empty_phi
739 ! ---------------------------------------------------------
740
741 ! ---------------------------------------------------------
742 subroutine laser_set_f_value(laser, ii, xx)
743 type(laser_t), intent(inout) :: laser
744 integer, intent(in) :: ii
745 real(real64), intent(in) :: xx
746
747 push_sub(laser_set_f_value)
748 call tdf_set_numerical(laser%f, ii, xx)
749
750 pop_sub(laser_set_f_value)
751 end subroutine laser_set_f_value
752 ! ---------------------------------------------------------
753
754
755 ! ---------------------------------------------------------
756 subroutine laser_set_frequency(laser, omega)
757 type(laser_t), intent(inout) :: laser
758 real(real64), intent(in) :: omega
759
760 push_sub(laser_set_frequency)
761 laser%omega = omega
762
763 pop_sub(laser_set_frequency)
764 end subroutine laser_set_frequency
765 ! ---------------------------------------------------------
766
768 ! ---------------------------------------------------------
769 subroutine laser_set_polarization(laser, pol)
770 type(laser_t), intent(inout) :: laser
771 complex(real64), intent(in) :: pol(:)
772
773 push_sub(laser_set_polarization)
774
775 laser%pol = pol
776
778 end subroutine laser_set_polarization
779 ! ---------------------------------------------------------
781
782 ! ---------------------------------------------------------
789 ! ---------------------------------------------------------
790 subroutine laser_to_numerical_all(laser, dt, max_iter, omegamax)
791 type(laser_t), intent(inout) :: laser
792 real(real64), intent(in) :: dt
793 integer, intent(in) :: max_iter
794 real(real64), intent(in) :: omegamax
796 integer :: iter
797 real(real64) :: tt, fj, phi
798
799 push_sub(lasers_to_numerical_all)
800
801 call tdf_to_numerical(laser%f, max_iter, dt, omegamax)
802 do iter = 1, max_iter + 1
803 tt = (iter-1)*dt
804 fj = tdf(laser%f, iter)
805 phi = tdf(laser%phi, tt)
806 call tdf_set_numerical(laser%f, iter, fj*cos(laser%omega*tt+phi))
807 end do
808 call tdf_end(laser%phi)
809 call tdf_init_cw(laser%phi, m_zero, m_zero)
810 laser%omega = m_zero
811
812 pop_sub(lasers_to_numerical_all)
813 end subroutine laser_to_numerical_all
814 ! ---------------------------------------------------------
815
816
817 ! ---------------------------------------------------------
820 ! ---------------------------------------------------------
821 subroutine laser_to_numerical(laser, dt, max_iter, omegamax)
822 type(laser_t), intent(inout) :: laser
823 real(real64), intent(in) :: dt
824 integer, intent(in) :: max_iter
825 real(real64), intent(in) :: omegamax
826
827 push_sub(lasers_to_numerical)
828
829 call tdf_to_numerical(laser%f, max_iter, dt, omegamax)
830 call tdf_to_numerical(laser%phi, max_iter, dt, omegamax)
831
832 pop_sub(lasers_to_numerical)
833 end subroutine laser_to_numerical
834 ! ---------------------------------------------------------
836 ! ---------------------------------------------------------
837 subroutine laser_write_info(lasers, namespace, dt, max_iter, iunit)
838 type(laser_t), intent(in) :: lasers(:)
839 type(namespace_t), intent(in) :: namespace
840 real(real64), optional, intent(in) :: dt
841 integer, optional, intent(in) :: max_iter
842 integer, optional, intent(in) :: iunit
843
844 real(real64) :: tt, fluence, max_intensity, intensity, dt_, field(3), up, maxfield,tmp
845 integer :: il, iter, no_l, max_iter_
846
847 push_sub(laser_write_info)
848
849 no_l = size(lasers)
850
851 do il = 1, no_l
852
853 if (present(dt)) then
854 dt_ = dt
855 else
856 dt_ = tdf_dt(lasers(il)%f)
857 end if
858 if (present(max_iter)) then
859 max_iter_ = max_iter
860 else
861 max_iter_ = tdf_niter(lasers(il)%f)
862 end if
863
864 write(message(1),'(i2,a)') il, ':'
865 select case (lasers(il)%field)
866 case (e_field_electric)
867 message(2) = ' Electric Field.'
868 case (e_field_magnetic)
869 message(2) = ' Magnetic Field.'
871 message(2) = ' Vector Potential.'
873 message(2) = ' Scalar Potential.'
874 end select
875 call messages_info(2, iunit=iunit, namespace=namespace)
876
877 if (lasers(il)%field /= e_field_scalar_potential) then
878 write(message(1),'(3x,a,3(a1,f7.4,a1,f7.4,a1))') 'Polarization: ', &
879 '(', real(lasers(il)%pol(1), real64), ',', aimag(lasers(il)%pol(1)), '), ', &
880 '(', real(lasers(il)%pol(2), real64), ',', aimag(lasers(il)%pol(2)), '), ', &
881 '(', real(lasers(il)%pol(3), real64), ',', aimag(lasers(il)%pol(3)), ')'
882 call messages_info(1, iunit=iunit, namespace=namespace)
883 end if
884
885 write(message(1),'(3x,a,f14.8,3a)') 'Carrier frequency = ', &
886 units_from_atomic(units_out%energy, lasers(il)%omega), &
887 ' [', trim(units_abbrev(units_out%energy)), ']'
888 message(2) = ' Envelope: '
889 call messages_info(2, iunit=iunit, namespace=namespace)
890 call tdf_write(lasers(il)%f, iunit)
891
892 if (.not. tdf_is_empty(lasers(il)%phi)) then
893 message(1) = ' Phase: '
894 call messages_info(1, iunit=iunit, namespace=namespace)
895 call tdf_write(lasers(il)%phi, iunit)
896 end if
897
898 ! 1 atomic unit of intensity = 3.5094448e+16 W / cm^2
899 ! In a Gaussian system of units,
900 ! I(t) = (1/(8\pi)) * c * E(t)^2
901 ! (1/(8\pi)) * c = 5.4525289841210 a.u.
902 if (lasers(il)%field == e_field_electric .or. lasers(il)%field == e_field_vector_potential) then
903 fluence = m_zero
904 max_intensity = m_zero
905 maxfield= m_zero
906 do iter = 1, max_iter_
907 tt = iter * dt_
908 call laser_electric_field(lasers(il), field, tt, dt_)
909 intensity = 5.4525289841210_real64*sum(field**2)
910 fluence = fluence + intensity
911 if (intensity > max_intensity) max_intensity = intensity
912
913 tmp = sum(field(:)**2)
914 if (tmp > maxfield) maxfield = tmp
915 end do
916 fluence = fluence * dt_
917
918 write(message(1),'(a,es17.6,3a)') ' Peak intensity = ', max_intensity, ' [a.u]'
919 write(message(2),'(a,es17.6,3a)') ' = ', &
920 max_intensity * 6.4364086e+15_real64, ' [W/cm^2]'
921 write(message(3),'(a,es17.6,a)') ' Int. intensity = ', fluence, ' [a.u]'
922 write(message(4),'(a,es17.6,a)') ' Fluence = ', &
923 fluence / 5.4525289841210_real64 , ' [a.u]'
924 call messages_info(4, iunit=iunit, namespace=namespace)
925
926 if (abs(lasers(il)%omega) > m_epsilon) then
927 ! Ponderomotive Energy is the cycle-averaged kinetic energy of
928 ! a free electron quivering in the field
929 ! Up = E^2/(4*\omega^2)
931 ! subroutine laser_to_numerical_all sets lasers%omega to zero
932 !
933 up = maxfield/(4*lasers(il)%omega**2)
934
935 write(message(1),'(a,es17.6,3a)') ' Ponderomotive energy = ', &
936 units_from_atomic(units_out%energy, up) ,&
937 ' [', trim(units_abbrev(units_out%energy)), ']'
938 call messages_info(1, iunit=iunit, namespace=namespace)
939 end if
940 end if
941
942 end do
943
944 pop_sub(laser_write_info)
945 end subroutine laser_write_info
946 ! ---------------------------------------------------------
947
948
949 ! ---------------------------------------------------------
950 subroutine laser_potential(laser, mesh, pot, time)
951 type(laser_t), intent(in) :: laser
952 class(mesh_t), intent(in) :: mesh
953 real(real64), intent(inout) :: pot(:)
954 real(real64), optional, intent(in) :: time
955
956 complex(real64) :: amp
957 integer :: ip
958 real(real64) :: field(3)
959
960 push_sub(laser_potential)
961
962 if (present(time)) then
963 amp = tdf(laser%f, time) * exp(m_zi * (laser%omega * time + tdf(laser%phi, time)))
964 else
965 amp = m_z1
966 end if
967
968 select case (laser%field)
970 pot(1:mesh%np) = pot(1:mesh%np) + real(amp, real64) * laser%v(1:mesh%np)
971 case default
972 field(:) = real(amp * laser%pol(:), real64)
973 do ip = 1, mesh%np
974 ! The -1 sign is missing here. Check epot.F90 for the explanation.
975 pot(ip) = pot(ip) + sum(field(1:mesh%box%dim) * mesh%x(ip, 1:mesh%box%dim))
976 end do
977 end select
978
979 pop_sub(laser_potential)
980 end subroutine laser_potential
981 ! ---------------------------------------------------------
982
983
984 ! ---------------------------------------------------------
985 subroutine laser_vector_potential(laser, mesh, aa, time)
986 type(laser_t), intent(in) :: laser
987 type(mesh_t), intent(in) :: mesh
988 real(real64), intent(inout) :: aa(:, :)
989 real(real64), optional, intent(in) :: time
990
991 real(real64) :: amp
992 integer :: ip, idir
993
994 push_sub(laser_vector_potential)
995
996 if (present(time)) then
997 amp = tdf(laser%f, time)*cos((laser%omega*time + tdf(laser%phi, time)))
998 do idir = 1, mesh%box%dim
999 do ip = 1, mesh%np
1000 aa(ip, idir) = aa(ip, idir) + amp*laser%a(ip, idir)
1001 end do
1002 end do
1003 else
1004 do idir = 1, mesh%box%dim
1005 do ip = 1, mesh%np
1006 aa(ip, idir) = aa(ip, idir) + laser%a(ip, idir)
1007 end do
1008 end do
1009 end if
1010
1011 pop_sub(laser_vector_potential)
1012 end subroutine laser_vector_potential
1013 ! ---------------------------------------------------------
1014
1015
1016 ! ---------------------------------------------------------
1022 subroutine laser_field(laser, field, time)
1023 type(laser_t), intent(in) :: laser
1024 real(real64), intent(inout) :: field(:)
1025 real(real64), optional, intent(in) :: time
1026
1027 integer :: dim
1028 complex(real64) :: amp
1029
1030 !no PUSH SUB, called too often
1031
1032 dim = size(field)
1033
1034 if (present(time)) then
1035 amp = tdf(laser%f, time) * exp(m_zi * (laser%omega * time + tdf(laser%phi, time)))
1036 else
1037 amp = m_z1
1038 end if
1039 if (laser%field == e_field_scalar_potential) then
1040 ! In this case we will just return the value of the time function. The "field", in fact,
1041 ! should be a function of the position in space (thus, a true "field"), given by the
1042 ! gradient of the scalar potential.
1043 field(1) = field(1) + real(amp, real64)
1044 else
1045 field(1:dim) = field(1:dim) + real(amp*laser%pol(1:dim), real64)
1046 end if
1047 end subroutine laser_field
1048
1049! ---------------------------------------------------------
1057 subroutine lasers_nondipole_laser_field_step(this, field, time)
1058 type(lasers_t), intent(in) :: this
1059 real(real64), intent(out) :: field(:)
1060 real(real64), intent(in) :: time
1061
1062 real(real64) :: a0(3)
1063 real(real64) :: e0(3)
1064 integer :: ilaser
1065 integer :: jlaser
1066 integer :: dim
1067 integer :: iter
1068 dim = size(field)
1069
1070 field(1:dim) = this%integrated_nondipole_afield(1:dim)
1071 if (time - this%nd_integration_time < 0) then
1072 return
1073 endif
1074 do iter = 1, nint((time-this%nd_integration_time)/this%nd_integration_step)
1075 do ilaser = 1, this%no_lasers
1076 if(laser_kind(this%lasers(ilaser)) /= e_field_vector_potential) cycle
1077 do jlaser = 1, this%no_lasers
1078 if(laser_kind(this%lasers(jlaser)) /= e_field_vector_potential) cycle
1079 if(allocated(this%lasers(jlaser)%prop)) then
1080 a0 = m_zero
1081 call laser_field(this%lasers(ilaser), a0, &
1082 iter*this%nd_integration_step+this%nd_integration_time)
1083 call laser_electric_field(this%lasers(jlaser), e0, &
1084 iter*this%nd_integration_step+this%nd_integration_time, this%nd_integration_step)
1085 ! We have front factor q/(mc) with q=-abs(e) =-1 a.u. and m = 1 a.u.
1086 field(1:dim) = field(1:dim) - m_one/(p_c) * this%nd_integration_step &
1087 * dot_product(a0,e0) * this%lasers(jlaser)%prop(1:dim)
1088 end if
1089 end do
1090 end do
1091 end do
1093
1094
1095 ! ---------------------------------------------------------
1098 subroutine laser_electric_field(laser, field, time, dt)
1099 type(laser_t), intent(in) :: laser
1100 real(real64), intent(out) :: field(:)
1101 real(real64), intent(in) :: time
1102 real(real64), intent(in) :: dt
1103
1104 integer :: dim
1105 real(real64), allocatable :: field1(:), field2(:)
1106
1107 !no PUSH SUB, called too often
1108
1109 dim = size(field)
1110
1111 select case (laser%field)
1112 case (e_field_electric)
1113 field = m_zero
1114 call laser_field(laser, field(1:dim), time)
1116 safe_allocate(field1(1:dim))
1117 safe_allocate(field2(1:dim))
1118 field1 = m_zero
1119 field2 = m_zero
1120 call laser_field(laser, field1(1:dim), time - dt)
1121 call laser_field(laser, field2(1:dim), time + dt)
1122 field = - (field2 - field1) / (m_two * p_c * dt)
1123 safe_deallocate_a(field1)
1124 safe_deallocate_a(field2)
1125 case default
1126 field = m_zero
1127 end select
1128
1129 end subroutine laser_electric_field
1130 ! ---------------------------------------------------------
1131
1132 ! ---------------------------------------------------------
1133 ! Loading of the lasers for the multisystem framework
1134 ! Ultimately this should be done in laser_init
1135 subroutine load_lasers(partners, namespace)
1136 class(partner_list_t), intent(inout) :: partners
1137 type(namespace_t), intent(in) :: namespace
1138
1139 class(lasers_t), pointer :: lasers
1140 integer :: il
1141
1142 push_sub(load_lasers)
1143
1144 lasers => lasers_t(namespace)
1145
1146 call lasers_parse_external_fields(lasers)
1147
1148 ! TODO: see how to do this in the multisystem framework
1149 if (parse_is_defined(namespace, 'MillerIndicesBasis')) then
1150 call messages_not_implemented("MillerIndicesBasis with load_lasers routine")
1151 end if
1152
1153 ! This is done normally in lasers_generate_potential, but we cannot call this routine here,
1154 ! so we do it here
1155 do il = 1, lasers%no_lasers
1156 lasers%lasers(il)%pol(:) = lasers%lasers(il)%pol(:)/sqrt(sum(abs(lasers%lasers(il)%pol(:))**2))
1157 end do
1158
1159 call lasers%quantities%add(quantity_t("E field", always_available = .true., updated_on_demand = .true., iteration = clock_t()))
1160 call lasers%quantities%add(quantity_t("B field", always_available = .true., updated_on_demand = .true., iteration = clock_t()))
1161
1162 lasers%supported_interactions_as_partner = [lorentz_force]
1163
1164 if (lasers%no_lasers > 0) then
1165 call partners%add(lasers)
1166 else
1167 safe_deallocate_p(lasers)
1168 end if
1169
1170 pop_sub(load_lasers)
1171 end subroutine load_lasers
1172
1173end module lasers_oct_m
1174
1175!! Local Variables:
1176!! mode: f90
1177!! coding: utf-8
1178!! End:
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_zero
Definition: global.F90:188
complex(real64), parameter, public m_z0
Definition: global.F90:198
real(real64), parameter, public m_epsilon
Definition: global.F90:204
complex(real64), parameter, public m_z1
Definition: global.F90:199
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
This module defines classes and functions for interaction partners.
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1031
subroutine, public load_lasers(partners, namespace)
Definition: lasers.F90:1229
complex(real64) function, dimension(3), public laser_polarization(laser)
Definition: lasers.F90:727
subroutine, public laser_set_phi(laser, phi)
Definition: lasers.F90:809
subroutine, public lasers_check_symmetries(this, kpoints)
Definition: lasers.F90:556
subroutine lasers_copy_quantities_to_interaction(partner, interaction)
Definition: lasers.F90:693
subroutine, public laser_to_numerical_all(laser, dt, max_iter, omegamax)
The td functions that describe the laser field are transformed to a "numerical" representation (i....
Definition: lasers.F90:884
subroutine, public laser_vector_potential(laser, mesh, aa, time)
Definition: lasers.F90:1079
subroutine, public lasers_nondipole_laser_field_step(this, field, time)
Retrieves the NDSFA vector_potential correction. The nondipole field is obtained for consecutive time...
Definition: lasers.F90:1151
subroutine, public lasers_parse_external_fields(this)
Definition: lasers.F90:244
subroutine, public lasers_set_nondipole_parameters(this, ndfield, nd_integration_time)
Set parameters for nondipole SFA calculation.
Definition: lasers.F90:753
subroutine lasers_update_quantity(this, label)
Definition: lasers.F90:651
subroutine, public laser_get_f(laser, ff)
Definition: lasers.F90:768
subroutine, public laser_set_f(laser, ff)
Definition: lasers.F90:781
logical function, public lasers_with_nondipole_field(lasers)
Check if a nondipole SFA correction should be computed for the given laser.
Definition: lasers.F90:740
subroutine, public laser_write_info(lasers, namespace, dt, max_iter, iunit)
Definition: lasers.F90:931
subroutine, public laser_set_empty_phi(laser)
Definition: lasers.F90:824
real(real64) function, public laser_carrier_frequency(laser)
Definition: lasers.F90:622
integer, parameter, public e_field_electric
Definition: lasers.F90:177
integer, parameter, public e_field_vector_potential
Definition: lasers.F90:177
subroutine, public laser_to_numerical(laser, dt, max_iter, omegamax)
The td functions that describe the laser field are transformed to a "numerical" representation (i....
Definition: lasers.F90:915
subroutine, public laser_electric_field(laser, field, time, dt)
Returns a vector with the electric field, no matter whether the laser is described directly as an ele...
Definition: lasers.F90:1192
subroutine, public lasers_generate_potentials(this, mesh, space, latt)
Definition: lasers.F90:444
subroutine lasers_init_interaction_as_partner(partner, interaction)
Definition: lasers.F90:633
subroutine lasers_finalize(this)
Definition: lasers.F90:581
subroutine, public laser_potential(laser, mesh, pot, time)
Definition: lasers.F90:1044
class(lasers_t) function, pointer lasers_constructor(namespace)
Definition: lasers.F90:226
integer, parameter, public e_field_scalar_potential
Definition: lasers.F90:177
integer pure elemental function, public laser_kind(laser)
Definition: lasers.F90:716
subroutine, public laser_set_frequency(laser, omega)
Definition: lasers.F90:850
subroutine, public laser_field(laser, field, time)
Retrieves the value of either the electric or the magnetic field. If the laser is given by a scalar p...
Definition: lasers.F90:1116
subroutine, public laser_set_polarization(laser, pol)
Definition: lasers.F90:863
subroutine, public laser_set_f_value(laser, ii, xx)
Definition: lasers.F90:836
subroutine, public laser_get_phi(laser, phi)
Definition: lasers.F90:796
subroutine lasers_deallocate(this)
Definition: lasers.F90:592
integer, parameter, public e_field_magnetic
Definition: lasers.F90:177
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:336
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1045
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:138
integer pure function, public symmetries_identity_index(this)
Definition: symmetries.F90:589
real(real64), public symprec
Definition: symmetries.F90:169
integer pure function, public symmetries_number(this)
Definition: symmetries.F90:543
subroutine, public tdf_end(f)
Definition: tdfunction.F90:974
subroutine, public tdf_init(f)
Definition: tdfunction.F90:388
subroutine, public tdf_read(f, namespace, function_name, ierr)
This function initializes "f" from the TDFunctions block.
Definition: tdfunction.F90:218
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
__clock_t clock_t
Definition: recipes.c:399
abstract class for general interaction partners
Describes mesh distribution to nodes.
Definition: mesh.F90:186
int true(void)