MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_opal_opacity.t
Go to the documentation of this file.
1!> This module reads in Rosseland-mean opacities from OPAL tables. Table opacity
2!> values are given in base 10 logarithm and are a function of mass density (R)
3!> and temperature (T), which are also both given in base 10 logarithm.
4
6 implicit none
7
8 !> min and max indices for R,T-range in opacity table
9 integer, parameter :: irmin = 2, irmax = 20
10 integer, parameter :: itmin = 7, itmax = 76
11
12 !> The opacity tables are read once and stored globally in Kappa_vals
13 double precision, public :: kappa_vals(itmin:itmax,irmin:irmax)
14 double precision, public :: logr_list(irmin:irmax), logt_list(itmin:itmax)
15
16 public :: init_opal_table
17 public :: set_opal_opacity
18
19contains
20
21 !> This subroutine is called when the FLD radiation module is initialised.
22 !> The OPAL tables for different helium abundances are read and interpolated.
23 subroutine init_opal_table(tablename)
24
25 character(len=*), intent(in) :: tablename
26
27 ! Local variables
28 character(len=256) :: amrvac_dir, path_table_dir
29
30 call get_environment_variable("AMRVAC_DIR", amrvac_dir)
31
32 ! Absolute path to OPAL table location
33 path_table_dir = trim(amrvac_dir)//"src/tables/OPAL_tables/"//trim(tablename)
34 ! Read in the OPAL table
35 call read_table(logr_list, logt_list, kappa_vals, path_table_dir)
36
37 end subroutine init_opal_table
38
39 !> This subroutine calculates the opacity for a given temperature-density
40 !> structure. Opacities are read from a table with given metalicity.
41 subroutine set_opal_opacity(rho, temp, kappa)
42 double precision, intent(in) :: rho, temp
43 double precision, intent(out) :: kappa
44 ! Local variables
45 double precision :: logr_in, logt_in, logkappa_out
46
47 ! Convert input density and temperature to OPAL table convention
48 logr_in = log10(rho/(temp*1d-6)**3)
49 logt_in = log10(temp)
50
51 call get_kappa(kappa_vals, logr_list, logt_list, logr_in, logt_in, logkappa_out)
52
53 ! If the outcome is 9.999 (no table entry), look right in the table
54 do while (logkappa_out .gt. 9.0d0)
55 logr_in = logr_in + 0.5d0
56 call get_kappa(kappa_vals, logr_list, logt_list, logr_in, logt_in, logkappa_out)
57 enddo
58
59 ! If the outcome is NaN, look left in the table
60 do while (logkappa_out .eq. 0.0d0)
61 logr_in = logr_in - 0.5d0
62 call get_kappa(kappa_vals, logr_list, logt_list, logr_in, logt_in, logkappa_out)
63 enddo
64
65 ! Convert OPAL opacity for output
66 kappa = 10d0**logkappa_out
67 end subroutine set_opal_opacity
68
69 !> This routine reads in 1-D (rho,T) values from an opacity table and gives
70 !> back as output the 1-D (rho,T) values and 2-D opacity
71 subroutine read_table(R, T, K, filename)
72
73 character(*), intent(in) :: filename
74 double precision, intent(out) :: K(iTmin:iTmax,iRmin:iRmax), R(iRmin:iRmax), T(iTmin:iTmax)
75
76 ! Local variables
77 character :: dum
78 integer :: row, col
79
80 open(unit=1, status='old', file=trim(filename))
81
82 ! Skip first 4 lines
83 do row = 1,4
84 read(1,*)
85 enddo
86
87 ! Read rho
88 read(1,*) dum, r(irmin:irmax)
89 read(1,*)
90
91 ! Read T and kappa
92 do row = itmin,itmax
93 read(1,'(f4.2,19f7.3)') t(row), k(row,irmin:irmax)
94 enddo
95
96 close(1)
97
98 end subroutine read_table
99
100 !> This subroutine looks in the table for the four couples (T,rho) surrounding
101 !> a given input T and rho
102 subroutine get_kappa(Kappa_vals, logR_list, logT_list, R, T, K)
103
104 double precision, intent(in) :: Kappa_vals(iTmin:iTmax,iRmin:iRmax)
105 double precision, intent(in) :: logR_list(iRmin:iRmax), logT_list(iTmin:iTmax)
106 double precision, intent(in) :: R, T
107 double precision, intent(out) :: K
108
109 ! Local variables
110 integer :: low_R_index, up_R_index
111 integer :: low_T_index, up_T_index
112
113 if (r .gt. maxval(logr_list)) then
114 ! print*, 'Extrapolating in logR'
115 low_r_index = irmax-1
116 up_r_index = irmax
117 elseif (r .lt. minval(logr_list)) then
118 ! print*, 'Extrapolating in logR'
119 low_r_index = irmin
120 up_r_index = irmin+1
121 else
122 call get_low_up_index(r, logr_list, irmin, irmax, low_r_index, up_r_index)
123 endif
124
125 if (t .gt. maxval(logt_list)) then
126 ! print*, 'Extrapolating in logT'
127 low_t_index = itmax-1
128 up_t_index = itmax
129 elseif ( t .lt. minval(logt_list)) then
130 ! print*, 'Extrapolating in logT'
131 low_t_index = itmin
132 up_t_index = itmin+1
133 else
134 call get_low_up_index(t, logt_list, itmin, itmax, low_t_index, up_t_index)
135 endif
136
137 call interpolate_krt(low_r_index, up_r_index, low_t_index, up_t_index, logr_list, logt_list, kappa_vals, r, t, k)
138
139 end subroutine get_kappa
140
141 !> This subroutine finds the indices in rho and T arrays of the two values
142 !> surrounding the input rho and T
143 subroutine get_low_up_index(var_in, var_list, imin, imax, low_i, up_i)
144
145 integer, intent(in) :: imin, imax
146 double precision, intent(in) :: var_in, var_list(imin:imax)
147 integer, intent(out) :: low_i, up_i
148
149 ! Local variables
150 double precision :: low_val, up_val, res(imin:imax)
151
152 res = var_list - var_in
153
154 ! Find all bounding values for given input
155 up_val = minval(res, mask = res .ge. 0) + var_in
156 low_val = maxval(res, mask = res .le. 0) + var_in
157
158 ! Find all bounding indices
159 up_i = minloc(abs(var_list - up_val),1) + imin -1
160 low_i = minloc(abs(var_list - low_val),1) + imin -1
161
162 if (up_i .eq. low_i) low_i = low_i - 1
163
164 end subroutine get_low_up_index
165
166 !> This subroutine does a bilinear interpolation in the R,T-plane
167 subroutine interpolate_krt(low_r, up_r, low_t, up_t, logR_list, logT_list, Kappa_vals, R, T, kappa_interp)
168
169 integer, intent(in) :: low_r, up_r, low_t, up_t
170 double precision, intent(in) :: Kappa_vals(iTmin:iTmax,iRmin:iRmax)
171 double precision, intent(in) :: logR_list(iRmin:iRmax), logT_list(iTmin:iTmax)
172 double precision, intent(in) :: R, T
173 double precision, intent(out) :: kappa_interp
174
175 ! Local variables
176 double precision :: r1, r2, t1, t2, k1, k2, k3, k4, ka, kb
177
178 ! First interpolate twice in the T coord to get ka and kb, then interpolate
179 ! in the R coord to get ki
180
181 ! r_1 rho r_2
182 ! | |
183 ! | |
184 ! ----k1--------ka-----k2----- t_1
185 ! | | |
186 ! | | |
187 ! T | | |
188 ! | | |
189 ! | ki |
190 ! | | |
191 ! ----k3--------kb-----k4----- t_2
192 ! | |
193 ! | |
194
195 r1 = logr_list(low_r)
196 r2 = logr_list(up_r)
197 t1 = logt_list(low_t)
198 t2 = logt_list(up_t)
199
200 k1 = kappa_vals(low_t, low_r)
201 k2 = kappa_vals(low_t, up_r)
202 k3 = kappa_vals(up_t, low_r)
203 k4 = kappa_vals(up_t, up_r)
204
205 call interpolate1d(r1,r2,r,k1,k2,ka)
206 call interpolate1d(r1,r2,r,k3,k4,kb)
207 call interpolate1d(t1,t2,t,ka,kb,kappa_interp)
208
209 end subroutine interpolate_krt
210
211 !> Interpolation in one dimension
212 subroutine interpolate1d(x1, x2, x, y1, y2, y)
213
214 double precision, intent(in) :: x, x1, x2, y1, y2
215 double precision, intent(out) :: y
216
217 y = y1 + (x-x1)*(y2-y1)/(x2-x1)
218
219 end subroutine interpolate1d
220
221end module mod_opal_opacity
This module reads in Rosseland-mean opacities from OPAL tables. Table opacity values are given in bas...
integer, parameter irmin
min and max indices for R,T-range in opacity table
subroutine read_table(r, t, k, filename)
This routine reads in 1-D (rho,T) values from an opacity table and gives back as output the 1-D (rho,...
double precision, dimension(itmin:itmax), public logt_list
subroutine, public init_opal_table(tablename)
This subroutine is called when the FLD radiation module is initialised. The OPAL tables for different...
double precision, dimension(itmin:itmax, irmin:irmax), public kappa_vals
The opacity tables are read once and stored globally in Kappa_vals.
subroutine, public set_opal_opacity(rho, temp, kappa)
This subroutine calculates the opacity for a given temperature-density structure. Opacities are read ...
subroutine get_low_up_index(var_in, var_list, imin, imax, low_i, up_i)
This subroutine finds the indices in rho and T arrays of the two values surrounding the input rho and...
subroutine get_kappa(kappa_vals, logr_list, logt_list, r, t, k)
This subroutine looks in the table for the four couples (T,rho) surrounding a given input T and rho.
subroutine interpolate1d(x1, x2, x, y1, y2, y)
Interpolation in one dimension.
integer, parameter irmax
integer, parameter itmin
double precision, dimension(irmin:irmax), public logr_list
integer, parameter itmax
subroutine interpolate_krt(low_r, up_r, low_t, up_t, logr_list, logt_list, kappa_vals, r, t, kappa_interp)
This subroutine does a bilinear interpolation in the R,T-plane.