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