The Elk Code
writeafpdt.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 subroutine writeafpdt
7 use modmain
8 use modtddft
9 implicit none
10 ! local variables
11 integer its,i
12 real(8) ed,t1
13 ! conversion factor of energy density to J/cm²
14 real(8), parameter :: ced=ha_si/(100.d0*br_si)**2
15 ! allocatable arrays
16 real(8), allocatable :: f(:),g(:),pd(:)
17 ! external functions
18 real(8), external :: splint
19 ! read time-dependent A-field from file
20 call readafieldt
21 ! allocate local arrays
22 allocate(f(ntimes),g(ntimes),pd(ntimes))
23 ! compute the power density at each time step
24 pd(:)=0.d0
25 do i=1,3
26  f(:)=afieldt(i,:)
27  call fderiv(1,ntimes,times,f,g)
28  pd(:)=pd(:)+g(:)**2
29 end do
30 t1=1.d0/(8.d0*pi*solsc)
31 pd(:)=t1*pd(:)
32 ! write the power density to file
33 open(50,file='AFPDT.OUT',form='FORMATTED')
34 do its=1,ntimes
35  write(50,'(2G18.10)') times(its),pd(its)
36 end do
37 close(50)
38 ! integrate power density to find the total energy density
39 ed=splint(ntimes,times,pd)
40 open(50,file='AFTED.OUT',form='FORMATTED')
41 write(50,*)
42 write(50,'("Total energy density : ",G18.10)') ed
43 write(50,'(" in J/cm² : ",G18.10)') ed*ced
44 close(50)
45 write(*,*)
46 write(*,'("Info(writeafpdt):")')
47 write(*,'(" Power density of A-field written to AFPDT.OUT")')
48 write(*,'(" Total energy density of A-field written to AFTED.OUT")')
49 deallocate(f,g,pd)
50 end subroutine
51 
real(8), parameter ha_si
Definition: modmain.f90:1255
integer ntimes
Definition: modtddft.f90:42
subroutine fderiv(m, n, x, f, g)
Definition: fderiv.f90:10
real(8), parameter pi
Definition: modmain.f90:1232
subroutine readafieldt
Definition: readafieldt.f90:7
real(8) solsc
Definition: modmain.f90:1253
real(8), dimension(:,:), allocatable afieldt
Definition: modtddft.f90:58
real(8), dimension(:), allocatable times
Definition: modtddft.f90:48
subroutine writeafpdt
Definition: writeafpdt.f90:7
real(8), parameter br_si
Definition: modmain.f90:1267