30 use,
intrinsic :: iso_fortran_env
84 integer,
public,
parameter :: &
93 integer :: field = e_field_none
94 complex(real64),
allocatable :: pol(:)
95 real(real64),
allocatable :: prop(:)
98 real(real64) :: omega =
m_zero
100 real(real64),
allocatable :: v(:)
101 real(real64),
allocatable :: a(:, :)
102 character(len=200) :: scalar_pot_expression
103 character(len=200) :: phase_expression
104 character(len=200) :: envelope_expression
110 integer,
public :: no_lasers
111 type(laser_t),
allocatable,
public :: lasers(:)
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
133 class(lasers_t),
pointer :: this
134 type(namespace_t),
intent(in) :: namespace
140 this%namespace =
namespace_t(
"Lasers", parent=namespace)
142 safe_allocate(this%e(1:3))
143 safe_allocate(this%b(1:3))
151 class(lasers_t),
intent(inout) :: this
154 integer :: il, jj, ierr, k
155 real(real64) :: omega0
156 complex(real64) :: cprop(3)
258 if (
parse_block(this%namespace,
'TDExternalFields', blk) == 0)
then
260 safe_allocate(this%lasers(1:this%no_lasers))
262 do il = 1, this%no_lasers
263 safe_allocate(this%lasers(il)%pol(1:3))
264 this%lasers(il)%pol =
m_z0
268 select case (this%lasers(il)%field)
272 this%lasers(il)%pol =
m_z1
283 this%lasers(il)%omega = omega0
286 call tdf_read(this%lasers(il)%f, this%namespace, trim(this%lasers(il)%envelope_expression), ierr)
291 call tdf_read(this%lasers(il)%phi, this%namespace, trim(this%lasers(il)%phase_expression), ierr)
293 write(
message(1),
'(3A)')
'Error in the "', trim(this%lasers(il)%envelope_expression), &
294 '" field defined in the TDExternalFields block:'
295 write(
message(2),
'(3A)')
'Time-dependent phase function "', trim(this%lasers(il)%phase_expression), &
305 safe_allocate(this%lasers(il)%prop(1:
size(cprop)))
306 do k = 1,
size(cprop)
310 if (any(abs(aimag(cprop)) >
m_epsilon))
then
311 write(
message(1),
'(3A)')
'Error in the "', trim(this%lasers(il)%envelope_expression), &
312 '" field defined in the TDExternalFields block:'
313 write(
message(2),
'(A)')
'Propagation direction cannot be complex.'
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
332 call parse_variable(this%namespace,
'NDSFATimeIntegrationStep', 0.01_real64, this%nd_integration_step)
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"
351 class(
lasers_t),
intent(inout) :: this
352 class(
mesh_t),
intent(in) :: mesh
353 class(
space_t),
intent(in) :: space
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)
363 do il = 1, this%no_lasers
396 if (
parse_block(this%namespace,
'MillerIndicesBasis', blk2) == 0)
then
397 if(.not. space%is_periodic())
then
398 write(
message(1),
'(a)')
'MillerIndicesBasis can only be used for periodic systems.'
402 do idir = 1, space%dim
403 do idir2 = 1, space%dim
410 this%lasers(il)%pol = matmul(miller, this%lasers(il)%pol)
414 this%lasers(il)%pol(:) = this%lasers(il)%pol(:)/
sqrt(sum(abs(this%lasers(il)%pol(:))**2))
416 select case (this%lasers(il)%field)
418 safe_allocate(this%lasers(il)%v(1:mesh%np_part))
419 this%lasers(il)%v =
m_zero
421 call mesh_r(mesh, ip, rr, coords = xx(1:space%dim))
422 xx(space%dim+1:3) =
m_zero
424 this%lasers(il)%v(ip) = pot_re
430 safe_allocate(this%lasers(il)%a(1:mesh%np_part, 1:3))
431 this%lasers(il)%a =
m_zero
433 xx(1:space%dim) = mesh%x(ip, :)
434 xx(space%dim+1:3) =
m_zero
438 select case (space%dim)
440 this%lasers(il)%a(ip, 1:2) = (/xx(2), -xx(1)/) * sign(
m_one, real(this%lasers(il)%pol(3)))
442 this%lasers(il)%a(ip, :) = (/ xx(2)*real(this%lasers(il)%pol(3)) - xx(3)*real(this%lasers(il)%pol(2)), &
443 xx(3)*real(this%lasers(il)%pol(1)) - xx(1)*real(this%lasers(il)%pol(3)), &
444 xx(1)*real(this%lasers(il)%pol(2)) - xx(2)*real(this%lasers(il)%pol(1)) /)
446 message(1) =
"Magnetic fields only allowed in 2 or 3D."
450 this%lasers(il)%a = -
m_half * this%lasers(il)%a
469 if (kpoints%use_symmetries)
then
472 do il = 1, this%no_lasers
474 message(1) =
"The lasers break (at least) one of the symmetries used to reduce the k-points ."
475 message(2) =
"Set SymmetryBreakDir accordingly to your laser fields."
487 type(
lasers_t),
intent(inout) :: this
498 class(
lasers_t),
intent(inout) :: this
504 safe_deallocate_a(this%e)
505 safe_deallocate_a(this%b)
506 safe_deallocate_a(this%integrated_nondipole_afield)
508 do il = 1, this%no_lasers
509 safe_deallocate_a(this%lasers(il)%pol)
510 call tdf_end(this%lasers(il)%f)
511 call tdf_end(this%lasers(il)%phi)
512 select case (this%lasers(il)%field)
514 safe_deallocate_a(this%lasers(il)%v)
516 safe_deallocate_a(this%lasers(il)%a)
518 safe_deallocate_a(this%lasers(il)%prop)
520 safe_deallocate_a(this%lasers)
528 type(
laser_t),
intent(in) :: laser
539 class(
lasers_t),
intent(in) :: partner
540 class(interaction_surrogate_t),
intent(inout) :: interaction
544 select type (interaction)
545 type is (lorentz_force_t)
548 message(1) =
"Unsupported interaction."
549 call messages_fatal(1, namespace=partner%namespace)
557 class(
lasers_t),
intent(inout) :: this
558 character(len=*),
intent(in) :: label
560 type(quantity_t),
pointer :: quantity
565 if (
allocated(this%lasers))
then
568 call messages_not_implemented(
"Laser vector potentials and scalar potentials in multi-system framework", &
569 namespace=this%namespace)
573 quantity => this%quantities%get(label)
577 do il = 1, this%no_lasers
579 call laser_field(this%lasers(il), this%e, quantity%iteration%value())
585 do il = 1, this%no_lasers
587 call laser_field(this%lasers(il), this%b, quantity%iteration%value())
592 message(1) =
"Incompatible quantity."
593 call messages_fatal(1, namespace=this%namespace)
601 class(
lasers_t),
intent(inout) :: partner
602 class(interaction_surrogate_t),
intent(inout) :: interaction
608 select type (interaction)
609 type is (lorentz_force_t)
610 do ip = 1, interaction%system_np
611 interaction%partner_e_field(:, ip) = partner%e
612 interaction%partner_b_field(:, ip) = partner%b
615 message(1) =
"Unsupported interaction."
616 call messages_fatal(1, namespace=partner%namespace)
623 integer pure elemental function
laser_kind(laser)
624 type(
laser_t),
intent(in) :: laser
635 type(
laser_t),
intent(in) :: laser
636 complex(real64) :: pol(3)
648 type(
lasers_t),
intent(in) :: lasers
649 logical :: isnondipole
652 isnondipole = .false.
653 if(
allocated(lasers%integrated_nondipole_afield)) isnondipole = .
true.
661 type(
lasers_t),
intent(inout) :: this
662 real(real64),
intent(in) :: ndfield(:)
663 real(real64),
intent(in) :: nd_integration_time
669 this%integrated_nondipole_afield(1:dim) = ndfield(1:dim)
670 this%nd_integration_time = nd_integration_time
678 type(
laser_t),
intent(in) :: laser
679 type(tdf_t),
intent(inout) :: ff
682 call tdf_copy(ff, laser%f)
691 type(
laser_t),
intent(inout) :: laser
692 type(tdf_t),
intent(inout) :: ff
696 call tdf_end(laser%f)
697 call tdf_copy(laser%f, ff)
706 type(
laser_t),
intent(in) :: laser
707 type(tdf_t),
intent(inout) :: phi
710 call tdf_copy(phi, laser%phi)
719 type(
laser_t),
intent(inout) :: laser
720 type(tdf_t),
intent(inout) :: phi
724 call tdf_end(laser%phi)
725 call tdf_copy(laser%phi, phi)
734 type(
laser_t),
intent(inout) :: laser
738 call tdf_init(laser%phi)
746 type(
laser_t),
intent(inout) :: laser
747 integer,
intent(in) :: ii
748 real(real64),
intent(in) :: xx
751 call tdf_set_numerical(laser%f, ii, xx)
760 type(
laser_t),
intent(inout) :: laser
761 real(real64),
intent(in) :: omega
773 type(
laser_t),
intent(inout) :: laser
774 complex(real64),
intent(in) :: pol(:)
794 type(
laser_t),
intent(inout) :: laser
795 real(real64),
intent(in) :: dt
796 integer,
intent(in) :: max_iter
797 real(real64),
intent(in) :: omegamax
800 real(real64) :: tt, fj, phi
802 push_sub(lasers_to_numerical_all)
804 call tdf_to_numerical(laser%f, max_iter, dt, omegamax)
805 do iter = 1, max_iter + 1
807 fj = tdf(laser%f, iter)
808 phi = tdf(laser%phi, tt)
809 call tdf_set_numerical(laser%f, iter, fj*
cos(laser%omega*tt+phi))
811 call tdf_end(laser%phi)
812 call tdf_init_cw(laser%phi, m_zero, m_zero)
815 pop_sub(lasers_to_numerical_all)
825 type(
laser_t),
intent(inout) :: laser
826 real(real64),
intent(in) :: dt
827 integer,
intent(in) :: max_iter
828 real(real64),
intent(in) :: omegamax
830 push_sub(lasers_to_numerical)
832 call tdf_to_numerical(laser%f, max_iter, dt, omegamax)
833 call tdf_to_numerical(laser%phi, max_iter, dt, omegamax)
835 pop_sub(lasers_to_numerical)
841 type(
laser_t),
intent(in) :: lasers(:)
842 type(namespace_t),
intent(in) :: namespace
843 real(real64),
optional,
intent(in) :: dt
844 integer,
optional,
intent(in) :: max_iter
845 integer,
optional,
intent(in) :: iunit
847 real(real64) :: tt, fluence, max_intensity, intensity, dt_, field(3), up, maxfield,tmp
848 integer :: il, iter, no_l, max_iter_
856 if (
present(dt))
then
859 dt_ = tdf_dt(lasers(il)%f)
861 if (
present(max_iter))
then
864 max_iter_ = tdf_niter(lasers(il)%f)
867 write(message(1),
'(i2,a)') il,
':'
868 select case (lasers(il)%field)
870 message(2) =
' Electric Field.'
872 message(2) =
' Magnetic Field.'
874 message(2) =
' Vector Potential.'
876 message(2) =
' Scalar Potential.'
878 call messages_info(2, iunit=iunit, namespace=namespace)
881 write(message(1),
'(3x,a,3(a1,f7.4,a1,f7.4,a1))')
'Polarization: ', &
882 '(', real(lasers(il)%pol(1), real64),
',', aimag(lasers(il)%pol(1)),
'), ', &
883 '(', real(lasers(il)%pol(2), real64),
',', aimag(lasers(il)%pol(2)),
'), ', &
884 '(', real(lasers(il)%pol(3), real64),
',', aimag(lasers(il)%pol(3)),
')'
885 call messages_info(1, iunit=iunit, namespace=namespace)
888 write(message(1),
'(3x,a,f14.8,3a)')
'Carrier frequency = ', &
889 units_from_atomic(units_out%energy, lasers(il)%omega), &
890 ' [', trim(units_abbrev(units_out%energy)),
']'
891 message(2) =
' Envelope: '
892 call messages_info(2, iunit=iunit, namespace=namespace)
893 call tdf_write(lasers(il)%f, iunit)
895 if (.not. tdf_is_empty(lasers(il)%phi))
then
896 message(1) =
' Phase: '
897 call messages_info(1, iunit=iunit, namespace=namespace)
898 call tdf_write(lasers(il)%phi, iunit)
907 max_intensity = m_zero
909 do iter = 1, max_iter_
912 intensity = 5.4525289841210_real64*sum(field**2)
913 fluence = fluence + intensity
914 if (intensity > max_intensity) max_intensity = intensity
916 tmp = sum(field(:)**2)
917 if (tmp > maxfield) maxfield = tmp
919 fluence = fluence * dt_
921 write(message(1),
'(a,es17.6,3a)')
' Peak intensity = ', max_intensity,
' [a.u]'
922 write(message(2),
'(a,es17.6,3a)')
' = ', &
923 max_intensity * 6.4364086e+15_real64,
' [W/cm^2]'
924 write(message(3),
'(a,es17.6,a)')
' Int. intensity = ', fluence,
' [a.u]'
925 write(message(4),
'(a,es17.6,a)')
' Fluence = ', &
926 fluence / 5.4525289841210_real64 ,
' [a.u]'
927 call messages_info(4, iunit=iunit, namespace=namespace)
929 if (abs(lasers(il)%omega) > m_epsilon)
then
936 up = maxfield/(4*lasers(il)%omega**2)
938 write(message(1),
'(a,es17.6,3a)')
' Ponderomotive energy = ', &
939 units_from_atomic(units_out%energy, up) ,&
940 ' [', trim(units_abbrev(units_out%energy)),
']'
941 call messages_info(1, iunit=iunit, namespace=namespace)
954 type(
laser_t),
intent(in) :: laser
955 class(mesh_t),
intent(in) :: mesh
956 real(real64),
intent(inout) :: pot(:)
957 real(real64),
optional,
intent(in) :: time
959 complex(real64) :: amp
961 real(real64) :: field(3)
965 if (
present(time))
then
966 amp = tdf(laser%f, time) *
exp(m_zi * (laser%omega * time + tdf(laser%phi, time)))
971 select case (laser%field)
973 pot(1:mesh%np) = pot(1:mesh%np) + real(amp, real64) * laser%v(1:mesh%np)
975 field(:) = real(amp * laser%pol(:), real64)
978 pot(ip) = pot(ip) + sum(field(1:mesh%box%dim) * mesh%x(ip, 1:mesh%box%dim))
989 type(
laser_t),
intent(in) :: laser
990 type(mesh_t),
intent(in) :: mesh
991 real(real64),
intent(inout) :: aa(:, :)
992 real(real64),
optional,
intent(in) :: time
999 if (
present(time))
then
1000 amp = tdf(laser%f, time)*
cos((laser%omega*time + tdf(laser%phi, time)))
1001 do idir = 1, mesh%box%dim
1003 aa(ip, idir) = aa(ip, idir) + amp*laser%a(ip, idir)
1007 do idir = 1, mesh%box%dim
1009 aa(ip, idir) = aa(ip, idir) + laser%a(ip, idir)
1026 type(
laser_t),
intent(in) :: laser
1027 real(real64),
intent(inout) :: field(:)
1028 real(real64),
optional,
intent(in) :: time
1031 complex(real64) :: amp
1037 if (
present(time))
then
1038 amp = tdf(laser%f, time) *
exp(m_zi * (laser%omega * time + tdf(laser%phi, time)))
1046 field(1) = field(1) + real(amp, real64)
1048 field(1:dim) = field(1:dim) + real(amp*laser%pol(1:dim), real64)
1062 real(real64),
intent(out) :: field(:)
1063 real(real64),
intent(in) :: time
1065 real(real64) :: a0(3)
1066 real(real64) :: e0(3)
1073 field(1:dim) = this%integrated_nondipole_afield(1:dim)
1074 if (time - this%nd_integration_time < 0)
then
1077 do iter = 1, nint((time-this%nd_integration_time)/this%nd_integration_step)
1078 do ilaser = 1, this%no_lasers
1080 do jlaser = 1, this%no_lasers
1082 if(
allocated(this%lasers(jlaser)%prop))
then
1085 iter*this%nd_integration_step+this%nd_integration_time)
1087 iter*this%nd_integration_step+this%nd_integration_time, this%nd_integration_step)
1089 field(1:dim) = field(1:dim) - m_one/(p_c) * this%nd_integration_step &
1090 * dot_product(a0,e0) * this%lasers(jlaser)%prop(1:dim)
1102 type(
laser_t),
intent(in) :: laser
1103 real(real64),
intent(out) :: field(:)
1104 real(real64),
intent(in) :: time
1105 real(real64),
intent(in) :: dt
1108 real(real64),
allocatable :: field1(:), field2(:)
1114 select case (laser%field)
1119 safe_allocate(field1(1:dim))
1120 safe_allocate(field2(1:dim))
1125 field = - (field2 - field1) / (m_two * p_c * dt)
1126 safe_deallocate_a(field1)
1127 safe_deallocate_a(field2)
1139 class(partner_list_t),
intent(inout) :: partners
1140 type(namespace_t),
intent(in) :: namespace
1152 if (parse_is_defined(namespace,
'MillerIndicesBasis'))
then
1153 call messages_not_implemented(
"MillerIndicesBasis with load_lasers routine")
1158 do il = 1, lasers%no_lasers
1159 lasers%lasers(il)%pol(:) = lasers%lasers(il)%pol(:)/
sqrt(sum(abs(lasers%lasers(il)%pol(:))**2))
1162 call lasers%quantities%add(quantity_t(
"E field", always_available = .
true., updated_on_demand = .
true., iteration =
clock_t()))
1163 call lasers%quantities%add(quantity_t(
"B field", always_available = .
true., updated_on_demand = .
true., iteration =
clock_t()))
1165 lasers%supported_interactions_as_partner = [lorentz_force]
1167 if (lasers%no_lasers > 0)
then
1168 call partners%add(lasers)
1170 safe_deallocate_p(lasers)
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_zero
complex(real64), parameter, public m_z0
real(real64), parameter, public m_epsilon
complex(real64), parameter, public m_z1
real(real64), parameter, public m_half
real(real64), parameter, public m_one
This module defines classes and functions for interaction partners.
subroutine, public kpoints_to_absolute(latt, kin, kout)
subroutine, public load_lasers(partners, namespace)
complex(real64) function, dimension(3), public laser_polarization(laser)
subroutine, public laser_set_phi(laser, phi)
subroutine, public lasers_check_symmetries(this, kpoints)
subroutine lasers_copy_quantities_to_interaction(partner, interaction)
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....
subroutine, public laser_vector_potential(laser, mesh, aa, time)
subroutine, public lasers_nondipole_laser_field_step(this, field, time)
Retrieves the NDSFA vector_potential correction. The nondipole field is obtained for consecutive time...
subroutine, public lasers_parse_external_fields(this)
subroutine, public lasers_set_nondipole_parameters(this, ndfield, nd_integration_time)
Set parameters for nondipole SFA calculation.
subroutine lasers_update_quantity(this, label)
subroutine, public laser_get_f(laser, ff)
subroutine, public laser_set_f(laser, ff)
logical function, public lasers_with_nondipole_field(lasers)
Check if a nondipole SFA correction should be computed for the given laser.
subroutine, public laser_write_info(lasers, namespace, dt, max_iter, iunit)
subroutine, public laser_set_empty_phi(laser)
real(real64) function, public laser_carrier_frequency(laser)
integer, parameter, public e_field_electric
integer, parameter, public e_field_vector_potential
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....
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...
subroutine, public lasers_generate_potentials(this, mesh, space, latt)
subroutine lasers_init_interaction_as_partner(partner, interaction)
subroutine lasers_finalize(this)
subroutine, public laser_potential(laser, mesh, pot, time)
class(lasers_t) function, pointer lasers_constructor(namespace)
integer, parameter, public e_field_scalar_potential
integer pure elemental function, public laser_kind(laser)
subroutine, public laser_set_frequency(laser, omega)
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...
subroutine, public laser_set_polarization(laser, pol)
subroutine, public laser_set_f_value(laser, ii, xx)
subroutine, public laser_get_phi(laser, phi)
subroutine lasers_deallocate(this)
integer, parameter, public e_field_magnetic
This module is intended to contain "only mathematical" functions and procedures.
This module defines the meshes, which are used in Octopus.
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
integer pure function, public symmetries_identity_index(this)
real(real64), public symprec
integer pure function, public symmetries_number(this)
subroutine, public tdf_end(f)
subroutine, public tdf_init(f)
subroutine, public tdf_read(f, namespace, function_name, ierr)
This function initializes "f" from the TDFunctions block.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
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.