File: hydrogen.f90

package info (click to toggle)
neat 2.3.2-2
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, forky, sid, trixie
  • size: 4,108 kB
  • sloc: f90: 7,385; python: 211; makefile: 78; sh: 64
file content (272 lines) | stat: -rw-r--r-- 7,929 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
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
!hydrogen.f90, calculations of balmer and paschen densities and temperatures
!(C) Roger Wesson
module mod_hydrogen
use mod_types
use mod_globals

implicit none
real(kind=dp), dimension(:,:,:,:), allocatable :: hidata
real(kind=dp), dimension(:), allocatable :: temperatures
integer :: ntemps, ndens, nlevs

contains

subroutine read_hydrogen
!get Case B emissivities from file, values are from Storey and Hummer 1995.
implicit none
character(len=1) :: junk
character(len=20), dimension(8) :: invar
integer :: i,j,k,l !counters

!debugging
#ifdef CO
        print *,"subroutine: read_hydrogen"
#endif

open (unit=357, file=trim(PREFIX)//"/share/neat/RHi.dat")
read (357,*) junk !first line is a comment
read (357,"(I3,I3)") ntemps, ndens !second line has number of temperatures and densities
nlevs=25 ! maximum number of levels shown in intrat data file

!allocate the emissivities array, dimensions are temperature, density, level 1, level 2

allocate(hidata(ntemps, ndens, nlevs, nlevs))
allocate(temperatures(ntemps))

hidata(:,:,:,:)=0.d0

!now, loop through the temperatures and densities and read in the emissivities

do i=1,ntemps
  do j=1,ndens
    read (357,*) invar(1:6) !line at top of each block has density, charge, temperature, case, maximum level calculated and maximum level displayed
    if (j.eq.1) then
      read(invar(3),"(E9.2)") temperatures(i) ! get the temperatures for later use in interpolating
    endif
    do k=nlevs,1,-1
      do l=1,k-1
        read (357,"(E10.3)", advance="no") hidata(i,j,k,l)
      enddo
    enddo
  enddo
enddo
close (357)

end subroutine read_hydrogen

subroutine balmer_densities(linelist,H_Balmer,medtemp,density)

implicit none
integer, dimension(3:40), intent(in) :: H_Balmer
type(line), dimension(:), intent(in) :: linelist
real(kind=dp), dimension(10:25) :: ratios,densities
real(kind=dp), dimension(:,:,:), allocatable :: searcharray
real(kind=dp) :: density, interp_t, weight
real(kind=dp), intent(in) :: medtemp
integer :: i,j !counters
real(kind=dp) :: r1, r2
integer :: H ! temporary, which line to look at

!debugging
#ifdef CO
        print *,"subroutine: balmer_densities"
#endif

density=0.d0
weight=0.d0

!H_Balmer indexing is such that entry 3 contains H3, etc
!allocate arrays to store values, high order lines (n>10) only

where (H_balmer(10:25) .gt. 0)
  ratios = linelist(H_balmer(10:25))%int_dered/linelist(H_Balmer(4))%int_dered
elsewhere
  ratios = 0.d0
endwhere

densities=0.d0

!allocate the search array.  emissivities array has dimensions of temperature, density, level 1, level 2
!search array just has dimensions of density, level 1, level 2

allocate(searcharray(size(hidata,2),size(hidata,3),size(hidata,4)))

!now find the temperature

do i=1,ntemps
  if (medtemp .lt. temperatures(i)) then
    exit
  endif
enddo

if (i>ntemps) i=ntemps ! because in do i=1,10, final value of i is 11

!interpolate linearly between temperatures to get line ratios with density at specified temperature
!if temperature is outside data limits, it is just set to the limit
!lower limit is 500K, upper is 30000K

if (medtemp .lt. temperatures(1) .or. medtemp .gt. temperatures(size(temperatures))) then
  searcharray(:,:,:)=hidata(i,:,:,:)
else
  interp_t=(medtemp-temperatures(i-1))/(temperatures(i)-temperatures(i-1))
  searcharray(:,:,:)=hidata(i-1,:,:,:)+interp_t*(hidata(i,:,:,:)-hidata(i,:,:,:))
endif

!now go through H_Balmer from 10 to 25, calculate the density for each line ratio, store it
!intrat data only goes up to n=25

do H=10,25
  if (H_Balmer(H) .gt. 0) then !line is present, calculate a density
    if (ratios(H) .lt. (searcharray(1,H,2)/searcharray(1,4,2))) then
      densities(H) = 2.
    elseif (ratios(H) .gt. (searcharray(ndens,H,2)/searcharray(ndens,4,2))) then
      densities(H) = 8.
    else !search for the value
      do j=1,ndens-1
        r1=searcharray(j,H,2)/searcharray(j,4,2)
        r2=searcharray(j+1,H,2)/searcharray(j+1,4,2)
        if (ratios(H) .gt. r1 .and. ratios(H) .lt. r2) then
          densities(H)=j+1+(ratios(H)-r1)/(r2-r1)
        endif
      enddo
    endif
  endif
enddo

!now we have an array with all the densities implied by the ratios.  Derive a flux weighted value

do i=10,25
  if (densities(i).gt.2.d0) then
    density=density+linelist(H_Balmer(i))%weight*densities(i)
    weight=weight+linelist(H_Balmer(i))%weight
  endif
enddo

if (weight .gt. 0.d0) then
  density=10**(density/weight)
else
  density=0.d0
endif

!check if lower limit, Storey and Hummer go down to 100cm-3

if (density .lt. 100. .and. density .ne. 0.d0) then
  density = 100.d0
endif

!deallocate search array

deallocate(searcharray)

end subroutine balmer_densities

subroutine paschen_densities(linelist,H_Paschen,medtemp,density)
! computed using ratios of Paschen lines to Hbeta
! will change to P9, checking if it is observed
implicit none
integer, dimension(4:39),intent(in) :: H_Paschen
type(line), dimension(:), intent(in) :: linelist
real(kind=dp), dimension(10:25) :: ratios,densities
real(kind=dp), dimension(:,:,:), allocatable :: searcharray
real(kind=dp) :: density, interp_t, weight
real(kind=dp), intent(in) :: medtemp
integer :: i,j !counters
real(kind=dp) :: r1, r2
integer :: H ! temporary, which line to look at

!debugging
#ifdef CO
        print *,"subroutine: paschendensities"
#endif

density=0.d0
weight=0.d0

!H_Paschen indexing is such that entry 4 contains P4, etc
!allocate arrays to store values, high order lines (n>10) only
!todo: this will be wrong if spectrum is not normalised to Hb=100
!should check for this

where (H_paschen(10:25) .gt. 0)
  ratios = linelist(H_paschen(10:25))%int_dered/100.d0
elsewhere
  ratios = 0.d0
endwhere

densities=0.d0

!allocate the search array.  emissivities array has dimensions of temperature, density, level 1, level 2
!search array just has dimensions of density, level 1, level 2

allocate(searcharray(size(hidata,2),size(hidata,3),size(hidata,4)))

!now find the temperature

do i=1,ntemps
  if (medtemp .lt. temperatures(i)) then
    exit
  endif
enddo

if (i>ntemps) i=ntemps ! because in do i=1,10, final value of i is 11

!interpolate linearly between temperatures to get line ratios with density at specified temperature
!if temperature is outside data limits, it is just set to the limit
!lower limit is 500K, upper is 30000K

if (medtemp .lt. temperatures(1) .or. medtemp .gt. temperatures(size(temperatures))) then
  searcharray(:,:,:)=hidata(i,:,:,:)
else
  interp_t=(medtemp-temperatures(i-1))/(temperatures(i)-temperatures(i-1))
  searcharray(:,:,:)=hidata(i-1,:,:,:)+interp_t*(hidata(i,:,:,:)-hidata(i,:,:,:))
endif

!now go through H_Paschen from 10 to 25, calculate the density for each line ratio, store it
!intrat data only goes up to n=25

do H=10,25
  if (H_Paschen(H) .gt. 0) then !line is present, calculate a density
    if (ratios(H) .lt. (searcharray(1,H,3)/searcharray(1,4,2))) then
      densities(H) = 2.
    elseif (ratios(H) .gt. (searcharray(ndens,H,3)/searcharray(ndens,4,2))) then
      densities(H) = 8.
    else !search for the value
      do j=1,ndens-1
        r1=searcharray(j,H,3)/searcharray(j,4,2)
        r2=searcharray(j+1,H,3)/searcharray(j+1,4,2)
        if (ratios(H) .gt. r1 .and. ratios(H) .lt. r2) then
          densities(H)=j+1+(ratios(H)-r1)/(r2-r1)
        endif
      enddo
    endif
  endif
enddo

!now we have an array with all the densities implied by the ratios.  Derive a flux weighted value

do i=10,25
  if (densities(i).gt.2.d0) then
    density=density+linelist(H_Paschen(i))%weight*densities(i)
    weight=weight+linelist(H_Paschen(i))%weight
  endif
enddo

if (weight .gt. 0.d0) then
  density=10**(density/weight)
else
  density=0.d0
endif

!check if lower limit, Storey and Hummer go down to 100cm-3

if (density .lt. 100. .and. density .ne. 0.d0) then
  density = 100.d0
endif

!deallocate search array

deallocate(searcharray)

end subroutine paschen_densities

end module mod_hydrogen