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