MPI-AMRVAC 3.2
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 use mod_comm_lib, only: mpistop
11
12 implicit none
13 private
14
15 !> min and max indices for R,T-range in opacity table
16 integer, parameter :: iDmin = 2, idmax = 51
17 integer, parameter :: iTmin = 2, itmax = 21
18
19 !> If user wants to read tables from different location than AMRVAC source
20 !> NOTE: the tables should have same size as AMRVAC tables
21 logical :: set_custom_tabledir = .false.
22
23 !> The opacity tables are read once and stored globally
24 double precision, public :: alpha_vals(idmin:idmax,itmin:itmax)
25 double precision, public :: qbar_vals(idmin:idmax,itmin:itmax)
26 double precision, public :: q0_vals(idmin:idmax,itmin:itmax)
27 double precision, public :: kappae_vals(idmin:idmax,itmin:itmax)
28 double precision, public :: logd_list(idmin:idmax), logt_list(itmin:itmax)
29
30 public :: init_cak_table
31 public :: set_cak_opacity
32
33contains
34
35 !> This routine is called when the FLD radiation module is initialised.
36 subroutine init_cak_table(tabledir,set_custom_tabledir)
37
38 character(len=*), intent(in) :: tabledir
39 logical, optional, intent(in) :: set_custom_tabledir
40
41 ! Local variables
42 character(len=256) :: amrvac_dir, path_table_dir
43
44 if (present(set_custom_tabledir) .and. set_custom_tabledir) then
45 path_table_dir = trim(tabledir)
46 else
47 call get_environment_variable("AMRVAC_DIR", amrvac_dir)
48 path_table_dir = trim(amrvac_dir)//"/src/tables/CAK_tables/"//trim(tabledir)
49 endif
50
51 call read_table(logd_list, logt_list, alpha_vals, trim(path_table_dir)//"/al_TD")
52 call read_table(logd_list, logt_list, qbar_vals, trim(path_table_dir)//"/Qb_TD")
53 call read_table(logd_list, logt_list, q0_vals, trim(path_table_dir)//"/Q0_TD")
54 call read_table(logd_list, logt_list, kappae_vals, trim(path_table_dir)//"/Ke_TD")
55
56 if(mype==0)then
57 write(*,*)'Initialized CAK tables, read in from',path_table_dir
58 endif
59 end subroutine init_cak_table
60
61 !> This subroutine calculates the opacity for a given temperature-density
62 !> structure. Opacities are read from a table with given metalicity.
63 subroutine set_cak_opacity(rho, temp, alpha_output, Qbar_output, Q0_output, kappae_output)
64
65 double precision, intent(in) :: rho, temp
66 double precision, intent(out) :: alpha_output, qbar_output, q0_output, kappae_output
67
68 ! Local variables
69 double precision :: d_input, t_input, kappa_cak
70
71 if (temp <= 0.0d0) call mpistop('Input temperature < 0 in set_cak_opacity.')
72 if (rho <= 0.0d0) call mpistop('Input density < 0 in set_cak_opacity.')
73
74 d_input = log10(rho)
75 t_input = log10(temp)
76
77 call get_val_comb(alpha_vals,qbar_vals,q0_vals,kappae_vals, &
78 logd_list, logt_list, d_input, t_input, &
79 alpha_output, qbar_output, q0_output, kappae_output)
80
81 end subroutine set_cak_opacity
82
83 !> This routine reads in 1-D (T,rho) values from a line opacity table and gives
84 !> back as output the 1-D (T,rho) values and 2-D line opacity value
85 subroutine read_table(D, T, K, filename)
86
87 character(*), intent(in) :: filename
88 double precision, intent(out) :: k(idmin:idmax,itmin:itmax), d(idmin:idmax), t(itmin:itmax)
89
90 ! Local variables
91 character :: dum
92 integer :: row, col, funit=99
93 logical :: alive
94
95 inquire(file=trim(filename), exist=alive)
96
97 if (alive) then
98 open(funit, file=trim(filename), status='unknown')
99 else
100 call mpistop('Table file you want to use cannot be found: '//filename)
101 endif
102
103 ! Read temperature
104 read(funit,*) dum, t(itmin:itmax)
105
106 ! Read rho and kappa
107 do row = idmin,idmax
108 read(funit,*) d(row), k(row,itmin:itmax)
109 enddo
110
111 close(funit)
112
113 end subroutine read_table
114
115 !> This subroutine looks in the table for the four couples (T,rho) surrounding
116 !> a given input T and rho
117 subroutine get_val_comb(K1_vals,K2_vals,K3_vals,K4_vals, &
118 logD_list, logT_list, D, T, K1, K2, K3, K4)
119
120 double precision, intent(in) :: k1_vals(idmin:idmax,itmin:itmax)
121 double precision, intent(in) :: k2_vals(idmin:idmax,itmin:itmax)
122 double precision, intent(in) :: k3_vals(idmin:idmax,itmin:itmax)
123 double precision, intent(in) :: k4_vals(idmin:idmax,itmin:itmax)
124 double precision, intent(in) :: logd_list(idmin:idmax), logt_list(itmin:itmax)
125 double precision, intent(in) :: d, t
126 double precision, intent(out) :: k1, k2, k3, k4
127
128 ! Local variables
129 integer :: low_d_index, up_d_index, low_t_index, up_t_index
130
131 if (d .ge. maxval(logd_list)) then
132 ! print*, 'Extrapolating in logR'
133 low_d_index = idmax-1
134 up_d_index = idmax
135 elseif (d .le. minval(logd_list)) then
136 ! print*, 'Extrapolating in logR'
137 low_d_index = idmin
138 up_d_index = idmin+1
139 else
140 call get_low_up_index(d, logd_list, idmin, idmax, low_d_index, up_d_index)
141 endif
142
143 if (t .ge. maxval(logt_list)) then
144 ! print*, 'Extrapolating in logT'
145 low_t_index = itmax-1
146 up_t_index = itmax
147 elseif (t .le. minval(logt_list)) then
148 ! print*, 'Extrapolating in logT'
149 low_t_index = itmin
150 up_t_index = itmin+1
151 else
152 call get_low_up_index(t, logt_list, itmin, itmax, low_t_index, up_t_index)
153 endif
154
155 call interpolate_krt(low_d_index, up_d_index, low_t_index, up_t_index, &
156 logd_list, logt_list, k1_vals, d, t, k1)
157
158 call interpolate_krt(low_d_index, up_d_index, low_t_index, up_t_index, &
159 logd_list, logt_list, k2_vals, d, t, k2)
160
161 call interpolate_krt(low_d_index, up_d_index, low_t_index, up_t_index, &
162 logd_list, logt_list, k3_vals, d, t, k3)
163
164 call interpolate_krt(low_d_index, up_d_index, low_t_index, up_t_index, &
165 logd_list, logt_list, k4_vals, d, t, k4)
166
167 end subroutine get_val_comb
168
169 !> This subroutine finds the indices in rho and T arrays of the two values
170 !> surrounding the input rho and T
171 subroutine get_low_up_index(var_in, var_list, imin, imax, low_i, up_i)
172
173 integer, intent(in) :: imin, imax
174 double precision, intent(in) :: var_in, var_list(imin:imax)
175 integer, intent(out) :: low_i, up_i
176
177 ! Local variables
178 double precision :: low_val, up_val, res(imin:imax)
179
180 res = var_list - var_in
181
182 ! Find all bounding values for given input
183 up_val = minval(res, mask = res >= 0) + var_in
184 low_val = maxval(res, mask = res <= 0) + var_in
185
186 ! Find all bounding indices
187 up_i = minloc(abs(var_list - up_val),1) + imin -1
188 low_i = minloc(abs(var_list - low_val),1) + imin -1
189
190 if (up_i == low_i) low_i = low_i - 1
191
192 end subroutine get_low_up_index
193
194 !> This subroutine does a bilinear interpolation in the R,T-plane
195 subroutine interpolate_krt(low_r, up_r, low_t, up_t, logD_list, logT_list, Kappa_vals, D, T, kappa_interp)
196
197 integer, intent(in) :: low_r, up_r, low_t, up_t
198 double precision, intent(in) :: kappa_vals(idmin:idmax,itmin:itmax)
199 double precision, intent(in) :: logd_list(idmin:idmax), logt_list(itmin:itmax)
200 double precision, intent(in) :: d, t
201 double precision, intent(out) :: kappa_interp
202
203 ! Local variables
204 double precision :: r1, r2, t1, t2, k1, k2, k3, k4, ka, kb
205
206 ! First interpolate twice in the T coord to get ka and kb, then interpolate
207 ! in the R coord to get ki
208
209 ! r_1 rho r_2
210 ! | |
211 ! | |
212 ! ----k1--------ka-----k2----- t_1
213 ! | | |
214 ! | | |
215 ! T | | |
216 ! | | |
217 ! | ki |
218 ! | | |
219 ! ----k3--------kb-----k4----- t_2
220 ! | |
221 ! | |
222
223 r1 = logd_list(low_r)
224 r2 = logd_list(up_r)
225 t1 = logt_list(low_t)
226 t2 = logt_list(up_t)
227
228 k1 = kappa_vals(low_r, low_t)
229 k2 = kappa_vals(low_r, up_t)
230 k3 = kappa_vals(up_r, low_t)
231 k4 = kappa_vals(up_r, up_t)
232
233 call interpolate1d(r1,r2,d,k1,k2,ka)
234 call interpolate1d(r1,r2,d,k3,k4,kb)
235 call interpolate1d(t1,t2,t,ka,kb,kappa_interp)
236
237 end subroutine interpolate_krt
238
239 !> Interpolation in one dimension
240 subroutine interpolate1d(x1, x2, x, y1, y2, y)
241
242 double precision, intent(in) :: x, x1, x2, y1, y2
243 double precision, intent(out) :: y
244
245 y = y1 + (x-x1)*(y2-y1)/(x2-x1)
246
247 end subroutine interpolate1d
248
249end 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
subroutine, public init_cak_table(tabledir, set_custom_tabledir)
This routine is called when the FLD radiation module is initialised.
double precision, dimension(idmin:idmax, itmin:itmax), public alpha_vals
The opacity tables are read once and stored globally.
double precision, dimension(idmin:idmax), public logd_list
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 ...
double precision, dimension(idmin:idmax, itmin:itmax), public kappae_vals
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d