MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_cak_opacity.t
Go to the documentation of this file.
1!> This module reads in CAK line opacities in the Gayley (1995) notation
2!> (alpha, Qbar, Q0, kappae) from corresponding tables. Tabulated values assume
3!> LTE conditions and are a function of mass density (D) and temperature (T),
4!> which are both given in base 10 logarithm.
5!> The construction of the tables is outlined in Poniatowski+ (2021), A&A, 667.
6
8
9 implicit none
10
11 !> min and max indices for R,T-range in opacity table
12 integer, parameter :: idmin = 2, idmax = 16
13 integer, parameter :: itmin = 2, itmax = 51
14
15 !> The opacity tables are read once and stored globally
16 double precision, public :: alpha_vals(idmin:idmax,itmin:itmax)
17 double precision, public :: qbar_vals(idmin:idmax,itmin:itmax)
18 double precision, public :: q0_vals(idmin:idmax,itmin:itmax)
19 double precision, public :: kappae_vals(idmin:idmax,itmin:itmax)
20 double precision, public :: logd_list(idmin:idmax), logt_list(itmin:itmax)
21
22 public :: init_cak_table
23 public :: set_cak_opacity
24
25contains
26
27 !> This routine is called when the FLD radiation module is initialised.
28 subroutine init_cak_table(tabledir)
29
30 character(len=*), intent(in) :: tabledir
31
32 ! Local variables
33 character(len=256) :: amrvac_dir, path_table_dir
34
35 call get_environment_variable("AMRVAC_DIR", amrvac_dir)
36 path_table_dir = trim(amrvac_dir)//"src/tables/CAK_tables/"//trim(tabledir)
37
38 call read_table(logd_list, logt_list, alpha_vals, trim(path_table_dir)//"/al_TD")
39 call read_table(logd_list, logt_list, qbar_vals, trim(path_table_dir)//"/Qb_TD")
40 call read_table(logd_list, logt_list, q0_vals, trim(path_table_dir)//"/Q0_TD")
41 call read_table(logd_list, logt_list, kappae_vals, trim(path_table_dir)//"/Ke_TD")
42
43 end subroutine init_cak_table
44
45 !> This subroutine calculates the opacity for a given temperature-density
46 !> structure. Opacities are read from a table with given metalicity.
47 subroutine set_cak_opacity(rho, temp, alpha_output, Qbar_output, Q0_output, kappae_output)
48 double precision, intent(in) :: rho, temp
49 double precision, intent(out) :: alpha_output, qbar_output, q0_output, kappae_output
50 ! Local variables
51 double precision :: d_input, t_input, kappa_cak
52
53 d_input = log10(rho)
54 t_input = log10(temp)
55
56 d_input = min(-7.d0-1.d-5, d_input)
57 d_input = max(-16.d0+1.d-5, d_input)
58 t_input = min(5d0-1.d-5, t_input)
59 t_input = max(4d0+1.d-5, t_input)
60
62 logd_list, logt_list, d_input, t_input, &
63 alpha_output, qbar_output, q0_output, kappae_output)
64 end subroutine set_cak_opacity
65
66 !> This routine reads in 1-D (T,rho) values from a line opacity table and gives
67 !> back as output the 1-D (T,rho) values and 2-D line opacity value
68 subroutine read_table(D, T, K, filename)
69
70 character(*), intent(in) :: filename
71 double precision, intent(out) :: K(iDmin:iDmax,iTmin:iTmax), D(iDmin:iDmax), T(iTmin:iTmax)
72
73 ! Local variables
74 character :: dum
75 integer :: row, col
76
77 open(unit=1, status='old', file=trim(filename))
78
79 ! Read temperature
80 read(1,*) dum, t(itmin:itmax)
81
82 ! Read rho and kappa
83 do row = idmin,idmax
84 read(1,*) d(row), k(row,itmin:itmax)
85 enddo
86
87 close(1)
88
89 end subroutine read_table
90
91 !> This subroutine looks in the table for the four couples (T,rho) surrounding
92 !> a given input T and rho
93 subroutine get_val_comb(K1_vals,K2_vals,K3_vals,K4_vals, &
94 logD_list, logT_list, D, T, &
95 K1, K2, K3, K4)
96
97 double precision, intent(in) :: K1_vals(iDmin:iDmax,iTmin:iTmax)
98 double precision, intent(in) :: K2_vals(iDmin:iDmax,iTmin:iTmax)
99 double precision, intent(in) :: K3_vals(iDmin:iDmax,iTmin:iTmax)
100 double precision, intent(in) :: K4_vals(iDmin:iDmax,iTmin:iTmax)
101 double precision, intent(in) :: logD_list(iDmin:iDmax), logT_list(iTmin:iTmax)
102 double precision, intent(in) :: D, T
103 double precision, intent(out) :: K1, K2, K3, K4
104
105 ! Local variables
106 integer :: low_D_index, up_D_index, low_T_index, up_T_index
107
108 if (d .gt. maxval(logd_list)) then
109 ! print*, 'Extrapolating in logR'
110 low_d_index = idmax-1
111 up_d_index = idmax
112 elseif (d .lt. minval(logd_list)) then
113 ! print*, 'Extrapolating in logR'
114 low_d_index = idmin
115 up_d_index = idmin+1
116 else
117 call get_low_up_index(d, logd_list, idmin, idmax, low_d_index, up_d_index)
118 endif
119
120 if (t .gt. maxval(logt_list)) then
121 ! print*, 'Extrapolating in logT'
122 low_t_index = itmax-1
123 up_t_index = itmax
124 elseif (t .lt. minval(logt_list)) then
125 ! print*, 'Extrapolating in logT'
126 low_t_index = itmin
127 up_t_index = itmin+1
128 else
129 call get_low_up_index(t, logt_list, itmin, itmax, low_t_index, up_t_index)
130 endif
131
132 call interpolate_krt(low_d_index, up_d_index, low_t_index, up_t_index, &
133 logd_list, logt_list, k1_vals, d, t, k1)
134
135 call interpolate_krt(low_d_index, up_d_index, low_t_index, up_t_index, &
136 logd_list, logt_list, k2_vals, d, t, k2)
137
138 call interpolate_krt(low_d_index, up_d_index, low_t_index, up_t_index, &
139 logd_list, logt_list, k3_vals, d, t, k3)
140
141 call interpolate_krt(low_d_index, up_d_index, low_t_index, up_t_index, &
142 logd_list, logt_list, k4_vals, d, t, k4)
143
144 end subroutine get_val_comb
145
146 !> This subroutine finds the indices in rho and T arrays of the two values
147 !> surrounding the input rho and T
148 subroutine get_low_up_index(var_in, var_list, imin, imax, low_i, up_i)
149
150 integer, intent(in) :: imin, imax
151 double precision, intent(in) :: var_in, var_list(imin:imax)
152 integer, intent(out) :: low_i, up_i
153
154 ! Local variables
155 double precision :: low_val, up_val, res(imin:imax)
156
157 res = var_list - var_in
158
159 ! Find all bounding values for given input
160 up_val = minval(res, mask = res .ge. 0) + var_in
161 low_val = maxval(res, mask = res .le. 0) + var_in
162
163 ! Find all bounding indices
164 up_i = minloc(abs(var_list - up_val),1) + imin -1
165 low_i = minloc(abs(var_list - low_val),1) + imin -1
166
167 if (up_i .eq. low_i) low_i = low_i - 1
168
169 end subroutine get_low_up_index
170
171 !> This subroutine does a bilinear interpolation in the R,T-plane
172 subroutine interpolate_krt(low_r, up_r, low_t, up_t, logD_list, logT_list, Kappa_vals, D, T, kappa_interp)
173
174 integer, intent(in) :: low_r, up_r, low_t, up_t
175 double precision, intent(in) :: Kappa_vals(iDmin:iDmax,iTmin:iTmax)
176 double precision, intent(in) :: logD_list(iDmin:iDmax), logT_list(iTmin:iTmax)
177 double precision, intent(in) :: D, T
178 double precision, intent(out) :: kappa_interp
179
180 ! Local variables
181 double precision :: r1, r2, t1, t2, k1, k2, k3, k4, ka, kb
182
183 ! First interpolate twice in the T coord to get ka and kb, then interpolate
184 ! in the R coord to get ki
185
186 ! r_1 rho r_2
187 ! | |
188 ! | |
189 ! ----k1--------ka-----k2----- t_1
190 ! | | |
191 ! | | |
192 ! T | | |
193 ! | | |
194 ! | ki |
195 ! | | |
196 ! ----k3--------kb-----k4----- t_2
197 ! | |
198 ! | |
199
200 r1 = logd_list(low_r)
201 r2 = logd_list(up_r)
202 t1 = logt_list(low_t)
203 t2 = logt_list(up_t)
204
205 k1 = kappa_vals(low_r, low_t)
206 k2 = kappa_vals(low_r, up_t)
207 k3 = kappa_vals(up_r, low_t)
208 k4 = kappa_vals(up_r, up_t)
209
210 call interpolate1d(r1,r2,d,k1,k2,ka)
211 call interpolate1d(r1,r2,d,k3,k4,kb)
212 call interpolate1d(t1,t2,t,ka,kb,kappa_interp)
213
214 end subroutine interpolate_krt
215
216 !> Interpolation in one dimension
217 subroutine interpolate1d(x1, x2, x, y1, y2, y)
218
219 double precision, intent(in) :: x, x1, x2, y1, y2
220 double precision, intent(out) :: y
221
222 y = y1 + (x-x1)*(y2-y1)/(x2-x1)
223
224 end subroutine interpolate1d
225
226end module mod_cak_opacity
This module reads in CAK line opacities in the Gayley (1995) notation (alpha, Qbar,...
double precision, dimension(itmin:itmax), public logt_list
double precision, dimension(idmin:idmax, itmin:itmax), public q0_vals
double precision, dimension(idmin:idmax, itmin:itmax), public qbar_vals
integer, parameter idmin
min and max indices for R,T-range in opacity table
double precision, dimension(idmin:idmax, itmin:itmax), public alpha_vals
The opacity tables are read once and stored globally.
subroutine interpolate_krt(low_r, up_r, low_t, up_t, logd_list, logt_list, kappa_vals, d, t, kappa_interp)
This subroutine does a bilinear interpolation in the R,T-plane.
integer, parameter itmin
integer, parameter idmax
double precision, dimension(idmin:idmax), public logd_list
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 interpolate1d(x1, x2, x, y1, y2, y)
Interpolation in one dimension.
subroutine, public set_cak_opacity(rho, temp, alpha_output, qbar_output, q0_output, kappae_output)
This subroutine calculates the opacity for a given temperature-density structure. Opacities are read ...
subroutine, public init_cak_table(tabledir)
This routine is called when the FLD radiation module is initialised.
subroutine read_table(d, t, k, filename)
This routine reads in 1-D (T,rho) values from a line opacity table and gives back as output the 1-D (...
double precision, dimension(idmin:idmax, itmin:itmax), public kappae_vals
integer, parameter itmax
subroutine get_val_comb(k1_vals, k2_vals, k3_vals, k4_vals, logd_list, logt_list, d, t, k1, k2, k3, k4)
This subroutine looks in the table for the four couples (T,rho) surrounding a given input T and rho.