File: writeafpdt.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (51 lines) | stat: -rw-r--r-- 1,399 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51

! Copyright (C) 2014 K. Krieger, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine writeafpdt
use modmain
use modtddft
implicit none
! local variables
integer its,i
real(8) ed,t1
! conversion factor of energy density to J/cm^2
real(8), parameter :: ced=ha_si/(100.d0*br_si)**2
! allocatable arrays
real(8), allocatable :: f(:),g(:),pd(:)
! external functions
real(8) splint
external splint
! allocate local arrays
allocate(f(ntimes),g(ntimes),pd(ntimes))
! compute the power density at each time step
pd(:)=0.d0
do i=1,3
  f(:)=afieldt(i,:)
  call fderiv(1,ntimes,times,f,g)
  pd(:)=pd(:)+g(:)**2
end do
t1=1.d0/(8.d0*pi*solsc)
pd(:)=t1*pd(:)
! write the power density to file
open(50,file='AFPDT.OUT',form='FORMATTED')
do its=1,ntimes
  write(50,'(2G18.10)') times(its),pd(its)
end do
close(50)
! integrate power density to find the total energy density
ed=splint(ntimes,times,pd)
open(50,file='AFTED.OUT',form='FORMATTED')
write(50,*)
write(50,'("Total energy density : ",G18.10)') ed
write(50,'(" in J/cm^2           : ",G18.10)') ed*ced
close(50)
write(*,*)
write(*,'("Info(writeafpdt):")')
write(*,'(" Power density of A-field written to AFPDT.OUT")')
write(*,'(" Total energy density of A-field written to AFTED.OUT")')
deallocate(f,g,pd)
return
end subroutine