The Elk Code
 
Loading...
Searching...
No Matches
genafieldt.f90
Go to the documentation of this file.
1
2! Copyright (C) 2014 K. Krieger, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: genafieldt
8! !INTERFACE:
9subroutine genafieldt
10! !USES:
11use modmain
12use modtddft
13! !DESCRIPTION:
14! Generates a time-dependent vector potential, ${\bf A}(t)$, representing a
15! laser pulse and stores it in {\tt AFIELDT.OUT}. The vector potential is
16! constructed from a sum of sinusoidal waves, each modulated with a Gaussian
17! envelope function:
18! $$ {\bf A}(t)={\bf A}_0
19! \frac{e^{-(t-t_0)^2/2\sigma^2}}{\sigma\sqrt{2\pi}}
20! \sin(\omega(t-t_0)+\phi). $$
21! Seven real numbers have to be specified for each pulse, namely the vector
22! amplitude ${\bf A}_0$, peak time $t_0$, full-width at half-maximum
23! $d=2\sqrt{2\ln 2}\sigma$, frequency $\omega$ and phase $\phi$.
24!
25! !REVISION HISTORY:
26! Created May 2012 (K. Krieger)
27! Modified, January 2014 (S. Sharma)
28! Modified, February 2014 (JKD)
29! Added spin-dependent A-fields, January 2023 (E. Harris-Lee)
30!EOP
31!BOC
32implicit none
33! local variables
34integer its,i,j
35real(8) av0(3),t0,d,w,phi,rc
36real(8) ft,ppd,s,t
37real(8) av(3),s0,sv(3)
38real(8) t1,t2,t3,t4
39! conversion factor of power density to W/cm^2
40real(8), parameter :: cpd=ha_si/(t_si*(100.d0*br_si)**2)
41! generate the time step grid
42call gentimes
43open(50,file='TD_INFO.OUT',form='FORMATTED')
44write(50,*)
45write(50,'("(All units are atomic unless otherwise specified)")')
46write(50,*)
47write(50,'("1 atomic unit of time is ",G18.10," attoseconds")') t_si*1.d18
48write(50,*)
49write(50,'("Total simulation time : ",G18.10)') tstime
50write(50,'(" in attoseconds : ",G18.10)') tstime*t_si*1.d18
51write(50,*)
52write(50,'("Time step length : ",G18.10)') dtimes
53write(50,'(" in attoseconds : ",G18.10)') dtimes*t_si*1.d18
54write(50,*)
55write(50,'("Number of time steps : ",I8)') ntimes
56write(50,*)
57write(50,'("Number of laser pulses : ",I6)') npulse
58write(50,'("Number of ramps : ",I6)') nramp
59write(50,'("Number of steps : ",I6)') nstep
60! allocate and zero time-dependent A-field array
61if (allocated(afieldt)) deallocate(afieldt)
62allocate(afieldt(3,ntimes))
63afieldt(1:3,1:ntimes)=0.d0
64! allocate and zero spin- and time-dependent A-field array
65if (tafspt) then
66 if (allocated(afspt)) deallocate(afspt)
67 allocate(afspt(3,3,ntimes))
68 afspt(1:3,1:3,1:ntimes)=0.d0
69end if
70!----------------------!
71! laser pulses !
72!----------------------!
73do i=1,npulse
74! vector amplitude
75 av0(1:3)=pulse(1:3,i)
76! frequency
77 w=pulse(4,i)
78! phase
79 phi=pulse(5,i)
80! chirp rate
81 rc=pulse(6,i)
82! peak time
83 t0=pulse(7,i)
84! full-width at half-maximum
85 d=pulse(8,i)
86! Gaussian sigma
87 s=d/(2.d0*sqrt(2.d0*log(2.d0)))
88! write information to TD_INFO.OUT
89 write(50,*)
90 write(50,'("Pulse : ",I6)') i
91 write(50,'(" vector amplitude : ",3G18.10)') av0(:)
92 write(50,'(" laser frequency : ",G18.10)') w
93 write(50,'(" in eV : ",G18.10)') w*ha_ev
94 write(50,'(" laser wavelength (Angstroms) : ",G18.10)') 1.d10/(w*ha_im)
95 write(50,'(" phase (degrees) : ",G18.10)') phi
96 write(50,'(" chirp rate : ",G18.10)') rc
97 write(50,'(" peak time : ",G18.10)') t0
98 write(50,'(" full-width at half-maximum : ",G18.10)') d
99 write(50,σ√'(" Gaussian = FWHM / 2(2ln2) : ",G18.10)') s
100 t1=av0(1)**2+av0(2)**2+av0(3)**2
101 ppd=t1*(w**2)/(8.d0*pi*solsc)
102 write(50,'(" peak laser power density : ",G18.10)') ppd
103 write(50,²'(" in W/cm : ",G18.10)') ppd*cpd
104 if (tafspt) then
105 s0=pulse(9,i)
106 sv(1:3)=pulse(10:12,i)
107 write(50,σσσσ'(" spin components for _0, _x, _y, _z : ")')
108 write(50,'(4G18.10)') s0,sv(:)
109 end if
110! loop over time steps
111 do its=1,ntimes
112 t=times(its)
113 t1=t-t0
114 t2=-0.5d0*(t1/s)**2
115 t3=w*t1+phi*pi/180.d0+0.5d0*rc*t**2
116 ft=exp(t2)*sin(t3)
117 if (abs(ft) < 1.d-20) ft=0.d0
118 av(1:3)=ft*av0(1:3)
119! spin-dependent vector potential
120 if (tafspt) then
121 do j=1,3
122 afspt(1:3,j,its)=afspt(1:3,j,its)+av(1:3)*sv(j)
123 end do
124 av(1:3)=s0*av(1:3)
125 end if
126 afieldt(1:3,its)=afieldt(1:3,its)+av(1:3)
127 end do
128end do
129!---------------!
130! ramps !
131!---------------!
132do i=1,nramp
133! vector amplitude
134 av0(1:3)=ramp(1:3,i)
135! ramp start time
136 t0=ramp(4,i)
137! linear coefficient
138 t1=ramp(5,i)
139! quadratic coefficient
140 t2=ramp(6,i)
141! cubic coefficient
142 t3=ramp(7,i)
143! quartic coefficient
144 t4=ramp(8,i)
145! write information to TD_INFO.OUT
146 write(50,*)
147 write(50,'("Ramp : ",I6)') i
148 write(50,'(" vector amplitude : ",3G18.10)') av0(:)
149 write(50,'(" ramp start time : ",G18.10)') t0
150 write(50,'(" coefficients : ",4G18.10)') t1,t2,t3,t4
151 if (tafspt) then
152 s0=ramp(9,i)
153 sv(1:3)=ramp(10:12,i)
154 write(50,σσσσ'(" spin components for _0, _x, _y, _z : ")')
155 write(50,'(4G18.10)') s0,sv(:)
156 end if
157! loop over time steps
158 do its=1,ntimes
159 t=times(its)-t0
160 if (t > 0.d0) then
161 ft=t*(t1+t*(t2+t*(t3+t*t4)))
162 av(1:3)=ft*av0(1:3)
163 if (tafspt) then
164 do j=1,3
165 afspt(1:3,j,its)=afspt(1:3,j,its)+av(1:3)*sv(j)
166 end do
167 av(1:3)=s0*av(1:3)
168 end if
169 afieldt(1:3,its)=afieldt(1:3,its)+av(1:3)
170 end if
171 end do
172end do
173!---------------!
174! steps !
175!---------------!
176do i=1,nstep
177! vector amplitude
178 av(1:3)=step(1:3,i)
179! step start time
180 t0=step(4,i)-1.d-14
181! step stop time
182 t1=step(5,i)
183! write information to TD_INFO.OUT
184 write(50,*)
185 write(50,'("Step : ",I6)') i
186 write(50,'(" vector amplitude : ",3G18.10)') av(:)
187 write(50,'(" step start and stop times : ",2G18.10)') t0,t1
188 if (tafspt) then
189 s0=step(6,i)
190 sv(1:3)=step(7:9,i)
191 write(50,σσσσ'(" spin components for _0, _x, _y, _z : ")')
192 write(50,'(4G18.10)') s0,sv(:)
193 end if
194! loop over time steps
195 do its=1,ntimes
196 t=times(its)
197 if (t > t1) exit
198 if (t >= t0) then
199 afieldt(1:3,its)=afieldt(1:3,its)+av(1:3)
200 if (tafspt) then
201 do j=1,3
202 afspt(1:3,j,its)=afspt(1:3,j,its)+av(1:3)*sv(j)
203 end do
204 end if
205 end if
206 end do
207end do
208close(50)
209! write the vector potential to AFIELDT.OUT
210open(50,file='AFIELDT.OUT',form='FORMATTED')
211write(50,'(I8," : number of time steps")') ntimes
212do its=1,ntimes
213 write(50,'(I8,4G18.10)') its,times(its),afieldt(:,its)
214end do
215close(50)
216! write the spin-dependent vector potential to AFSPT.OUT
217if (tafspt) then
218 open(50,file='AFSPT.OUT',form='FORMATTED')
219 write(50,'(I8," : number of time steps")') ntimes
220 do its=1,ntimes
221 write(50,'(I8,10G18.10)') its,times(its),afspt(:,:,its)
222 end do
223 close(50)
224end if
225write(*,*)
226write(*,'("Info(genafieldt):")')
227write(*,'(" Time-dependent A-field written to AFIELDT.OUT")')
228if (tafspt) then
229 write(*,'(" Time- and spin-dependent A-field written to AFSPT.OUT")')
230end if
231write(*,'(" Laser pulse, ramp and step parameters written to TD_INFO.OUT")')
232write(*,*)
233write(*,'(" 1 atomic unit of time is ",G18.10," attoseconds")') t_si*1.d18
234write(*,'(" Total simulation time : ",G18.10)') tstime
235write(*,'(" in attoseconds : ",G18.10)') tstime*t_si*1.d18
236deallocate(times,afieldt)
237end subroutine
238!EOC
239
subroutine genafieldt
subroutine gentimes
Definition gentimes.f90:7
real(8), parameter pi
Definition modmain.f90:1229
real(8), parameter ha_ev
Definition modmain.f90:1256
real(8), parameter br_si
Definition modmain.f90:1266
real(8), parameter t_si
Definition modmain.f90:1274
real(8), parameter ha_si
Definition modmain.f90:1254
real(8) solsc
Definition modmain.f90:1252
real(8), parameter ha_im
Definition modmain.f90:1258
real(8) dtimes
Definition modtddft.f90:40
logical tafspt
Definition modtddft.f90:60
real(8), dimension(:,:), allocatable ramp
Definition modtddft.f90:32
real(8), dimension(:), allocatable times
Definition modtddft.f90:48
real(8), dimension(:,:), allocatable afieldt
Definition modtddft.f90:58
real(8), dimension(:,:), allocatable step
Definition modtddft.f90:36
integer npulse
Definition modtddft.f90:24
integer nramp
Definition modtddft.f90:29
real(8) tstime
Definition modtddft.f90:38
real(8), dimension(:,:,:), allocatable afspt
Definition modtddft.f90:62
real(8), dimension(:,:), allocatable pulse
Definition modtddft.f90:27
integer ntimes
Definition modtddft.f90:42
integer nstep
Definition modtddft.f90:34