Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2017-2019 Josiah Schwab & The MESA Team
4 : !
5 : ! This program is free software: you can redistribute it and/or modify
6 : ! it under the terms of the GNU Lesser General Public License
7 : ! as published by the Free Software Foundation,
8 : ! either version 3 of the License, or (at your option) any later version.
9 : !
10 : ! This program is distributed in the hope that it will be useful,
11 : ! but WITHOUT ANY WARRANTY; without even the implied warranty of
12 : ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
13 : ! See the GNU Lesser General Public License for more details.
14 : !
15 : ! You should have received a copy of the GNU Lesser General Public License
16 : ! along with this program. If not, see <https://www.gnu.org/licenses/>.
17 : !
18 : ! ***********************************************************************
19 :
20 :
21 : module weaklib_tables
22 :
23 : use const_def, only: dp, ln10
24 : use math_lib
25 : use utils_lib, only: mesa_error, is_bad
26 : use rates_def
27 :
28 : implicit none
29 :
30 : type, extends(weak_rate_table) :: weaklib_rate_table
31 : integer :: i_ldecay = 1
32 : integer :: i_lcapture = 2
33 : integer :: i_lneutrino = 3
34 :
35 : contains
36 :
37 : procedure :: setup => setup_weaklib_table
38 : procedure :: interpolate => interpolate_weaklib_table
39 :
40 : final :: deallocate_weaklib_rate_table, deallocate_weaklib_rate_table_array
41 : end type weaklib_rate_table
42 :
43 : interface weaklib_rate_table
44 : module procedure new_weaklib_rate_table
45 : end interface weaklib_rate_table
46 :
47 : contains
48 :
49 :
50 2313 : function new_weaklib_rate_table(T9s, lYeRhos)
51 : real(dp), intent(in), dimension(:) :: T9s, lYeRhos
52 : type(weaklib_rate_table) :: new_weaklib_rate_table
53 :
54 2313 : new_weaklib_rate_table% num_T9 = size(T9s)
55 2313 : allocate(new_weaklib_rate_table% T9s(new_weaklib_rate_table% num_T9))
56 32499 : new_weaklib_rate_table% T9s = T9s
57 :
58 2313 : new_weaklib_rate_table% num_lYeRho = size(lYeRhos)
59 2313 : allocate(new_weaklib_rate_table% lYeRhos(new_weaklib_rate_table% num_lYeRho))
60 30252 : new_weaklib_rate_table% lYeRhos = lYeRhos
61 :
62 : allocate(new_weaklib_rate_table% data(4, &
63 : new_weaklib_rate_table% num_T9, &
64 9252 : new_weaklib_rate_table% num_lYeRho, 3))
65 2313 : end function new_weaklib_rate_table
66 :
67 :
68 2313 : subroutine setup_weaklib_table(table, ierr)
69 : class(weaklib_rate_table), intent(inout) :: table
70 : integer, intent(out) :: ierr
71 :
72 : integer :: ii
73 9252 : do ii = 1, size(table% data, dim=4)
74 9252 : if (weak_bicubic) then
75 0 : call create_bicubic_interpolant(table,ii,ierr)
76 : else
77 6939 : call create_pm_T9_interpolant(table,ii,ierr)
78 : end if
79 : end do
80 :
81 : contains
82 :
83 6939 : subroutine create_pm_T9_interpolant(table,ii,ierr)
84 : ! piecewise monotonic interpolation in T9 for each lYeRho in table
85 : use interp_1d_def
86 : use interp_1d_lib
87 : class(weaklib_rate_table), intent(inout) :: table
88 : integer, intent(in) :: ii
89 : integer, intent(out) :: ierr
90 : integer :: j, m, n
91 :
92 : integer :: nx ! length of x vector (>= 2)
93 : real(dp), pointer :: x(:) ! (nx) ! junction points, strictly monotonic
94 6939 : real(dp), pointer :: f1(:), f(:,:) ! (4,nx) ! data & interpolation coefficients
95 : integer, parameter :: nwork = pm_work_size
96 6939 : real(dp), pointer :: work(:) ! =(nx,nwork)
97 :
98 6939 : ierr = 0
99 :
100 6939 : nx = table % num_T9
101 6939 : allocate(x(nx), f1(4*nx), work(nx*nwork), stat=ierr)
102 6939 : if (ierr /= 0) return
103 :
104 6939 : f(1:4,1:nx) => f1(1:4*nx)
105 :
106 90558 : x = table % T9s
107 :
108 83817 : do j=1,table % num_lYeRho
109 1024686 : do m=1,nx
110 1024686 : f(1,m) = table % data(1,m,j,ii)
111 : end do
112 76878 : call interp_pm(x, nx, f1, nwork, work, 'create_pm_T9_interpolant', ierr)
113 76878 : if (ierr /= 0) return
114 1031625 : do n=1,nx
115 4815918 : do m=1,4
116 4739040 : table % data(m,n,j,ii) = real(f(m,n),kind=dp)
117 : end do
118 : end do
119 : end do
120 :
121 6939 : deallocate(x, f1, work)
122 :
123 13878 : end subroutine create_pm_T9_interpolant
124 :
125 :
126 0 : subroutine create_bicubic_interpolant(table,ii,ierr)
127 6939 : use interp_2d_lib_db
128 : class(weaklib_rate_table), intent(inout) :: table
129 : integer, intent(in) :: ii
130 : integer, intent(out) :: ierr
131 : integer :: ibcxmin ! bc flag for x=xmin
132 : real(dp), allocatable :: bcxmin(:) ! bc data vs. y at x=xmin
133 : integer :: ibcxmax ! bc flag for x=xmax
134 : real(dp), allocatable :: bcxmax(:) ! bc data vs. y at x=xmax
135 : integer :: ibcymin ! bc flag for y=ymin
136 : real(dp), allocatable :: bcymin(:) ! bc data vs. x at y=ymin
137 : integer :: ibcymax ! bc flag for y=ymax
138 : real(dp), allocatable :: bcymax(:) ! bc data vs. x at y=ymax
139 : integer :: il1, il2, j, k, m
140 :
141 : integer :: num_T9, num_lYeRho
142 0 : real(dp), pointer :: f1(:), f3(:,:,:)
143 :
144 0 : num_T9 = table % num_T9
145 0 : num_lYeRho = table % num_lYeRho
146 :
147 0 : allocate(bcxmin(num_lYeRho), bcxmax(num_lYeRho))
148 0 : allocate(bcymin(num_T9), bcymax(num_T9))
149 :
150 : ! just use "not a knot" bc's at edges of tables
151 0 : ibcxmin = 0; bcxmin(:) = 0
152 0 : ibcxmax = 0; bcxmax(:) = 0
153 0 : ibcymin = 0; bcymin(:) = 0
154 0 : ibcymax = 0; bcymax(:) = 0
155 :
156 0 : allocate(f1(4*num_T9*num_lYeRho))
157 :
158 0 : f3(1:4,1:num_T9,1:num_lYeRho) => f1(1:4*num_T9*num_lYeRho)
159 0 : do k = 1,num_T9
160 0 : do j = 1,4
161 0 : do m = 1,num_lYeRho
162 0 : f3(j,k,m) = table % data(j,k,m,ii)
163 : end do
164 : end do
165 : end do
166 : call interp_mkbicub_db( &
167 : table % T9s, num_T9, &
168 : table % lYeRhos, num_lYeRho, &
169 : f1, num_T9, &
170 : ibcxmin, bcxmin, ibcxmax, bcxmax, &
171 : ibcymin, bcymin, ibcymax, bcymax, &
172 0 : il1, il2, ierr)
173 0 : do k = 1,num_T9
174 0 : do j = 1,4
175 0 : do m = 1,num_lYeRho
176 0 : table % data(j,k,m,ii) = f3(j,k,m)
177 : end do
178 : end do
179 : end do
180 0 : deallocate(f1, bcxmin, bcxmax, bcymin, bcymax)
181 0 : end subroutine create_bicubic_interpolant
182 :
183 : end subroutine setup_weaklib_table
184 :
185 10060 : subroutine interpolate_weaklib_table(table, T9, lYeRho, &
186 : lambda, dlambda_dlnT, dlambda_dlnRho, &
187 : Qneu, dQneu_dlnT, dQneu_dlnRho, ierr)
188 : class(weaklib_rate_table), intent(inout) :: table
189 : real(dp), intent(in) :: T9, lYeRho
190 : real(dp), intent(out) :: lambda, dlambda_dlnT, dlambda_dlnRho
191 : real(dp), intent(out) :: Qneu, dQneu_dlnT, dQneu_dlnRho
192 : integer, intent(out) :: ierr
193 :
194 : integer :: ix, jy ! target cell in the spline data
195 10060 : real(dp) :: x0, xget, x1 ! x0 <= xget <= x1; x0 = xs(ix), x1 = xs(ix+1)
196 10060 : real(dp) :: y0, yget, y1 ! y0 <= yget <= y1; y0 = ys(jy), y1 = ys(jy+1)
197 10060 : real(dp) :: xp, xpi, xp2, xpi2, cx, cxi, hx2, cxd, cxdi, hx, hxi
198 10060 : real(dp) :: yp, ypi, yp2, ypi2, cy, cyi, hy2, cyd, cydi, hy, hyi
199 :
200 10060 : real(dp) :: delta_T9, dT9, dlYeRho, delta_lYeRho, y_alfa, y_beta, x_alfa, x_beta
201 :
202 10060 : real(dp) :: ldecay, d_ldecay_dT9, d_ldecay_dlYeRho, &
203 10060 : lcapture, d_lcapture_dT9, d_lcapture_dlYeRho, &
204 10060 : lneutrino, d_lneutrino_dT9, d_lneutrino_dlYeRho
205 :
206 10060 : real(dp) :: decay, capture
207 :
208 : logical, parameter :: dbg = .false.
209 :
210 10060 : xget = T9
211 10060 : yget = lYeRho
212 :
213 10060 : if (weak_bicubic) then
214 0 : call setup_for_bicubic_interpolations
215 : else
216 10060 : call setup_for_linear_interp
217 : end if
218 :
219 10060 : if (weak_bicubic) then
220 :
221 : call do_bicubic_interpolations( &
222 : table % data(1:4,1:table%num_T9,1:table%num_lYeRho,table%i_ldecay), &
223 0 : ldecay, d_ldecay_dT9, d_ldecay_dlYeRho, ierr)
224 :
225 : call do_bicubic_interpolations( &
226 : table % data(1:4,1:table%num_T9,1:table%num_lYeRho,table%i_lcapture), &
227 0 : lcapture, d_lcapture_dT9, d_lcapture_dlYeRho, ierr)
228 :
229 : call do_bicubic_interpolations( &
230 : table % data(1:4,1:table%num_T9,1:table%num_lYeRho,table%i_lneutrino), &
231 0 : lneutrino, d_lneutrino_dT9, d_lneutrino_dlYeRho, ierr)
232 :
233 : else
234 :
235 : call do_linear_interp( &
236 : table % data(1:4,1:table%num_T9,1:table%num_lYeRho,table%i_ldecay), &
237 10060 : ldecay, d_ldecay_dT9, d_ldecay_dlYeRho, ierr)
238 :
239 : call do_linear_interp( &
240 : table % data(1:4,1:table%num_T9,1:table%num_lYeRho,table%i_lcapture), &
241 10060 : lcapture, d_lcapture_dT9, d_lcapture_dlYeRho, ierr)
242 :
243 : call do_linear_interp( &
244 : table % data(1:4,1:table%num_T9,1:table%num_lYeRho,table%i_lneutrino), &
245 10060 : lneutrino, d_lneutrino_dT9, d_lneutrino_dlYeRho, ierr)
246 :
247 : end if
248 :
249 10060 : decay = exp10(ldecay)
250 10060 : capture = exp10(lcapture)
251 10060 : lambda = decay + capture
252 :
253 : ! lrates are log10
254 : ! T9 is linear
255 : ! lYeRho is log10
256 :
257 : ! drate_dlnT = T9*drate_dT9 = ln10 * rate * d_lrate_dT9
258 10060 : dlambda_dlnT = ln10*T9*(decay*d_ldecay_dT9 + capture*d_lcapture_dT9)
259 : ! drate_dlnRho = drate_dlnYeRho (assuming Ye held fixed)
260 : ! = rate * dlnrate_dlnYeRho = rate * dlrate_dlYeRho
261 10060 : dlambda_dlnRho = decay*d_ldecay_dlYeRho + capture*d_lcapture_dlYeRho
262 10060 : Qneu = exp10(lneutrino)/lambda
263 10060 : dQneu_dlnT = ln10*T9*Qneu*d_lneutrino_dT9 - dlambda_dlnT*Qneu/lambda
264 10060 : dQneu_dlnRho = Qneu*d_lneutrino_dlYeRho - dlambda_dlnRho*Qneu/lambda
265 :
266 : contains
267 :
268 40240 : subroutine find_location ! set ix, jy; x is T9; y is lYeRho
269 : integer :: i, j
270 : include 'formats'
271 : ! x0 <= T9 <= x1
272 10060 : ix = table % num_T9-1 ! since weak_num_T9 is small, just do a linear search
273 67138 : do i = 2, table % num_T9-1
274 67136 : if (T9 > table% T9s(i)) cycle
275 10058 : ix = i-1
276 67138 : exit
277 : end do
278 :
279 : ! y0 <= lYeRho <= y1
280 10060 : jy = table % num_lYeRho-1 ! since weak_num_lYeRho is small, just do a linear search
281 61526 : do j = 2, table % num_lYeRho-1
282 61524 : if (lYeRho > table % lYeRhos(j)) cycle
283 10058 : jy = j-1
284 61526 : exit
285 : end do
286 :
287 10060 : x0 = table % T9s(ix)
288 10060 : x1 = table % T9s(ix+1)
289 10060 : y0 = table % lYeRhos(jy)
290 10060 : y1 = table % lYeRhos(jy+1)
291 :
292 10060 : end subroutine find_location
293 :
294 0 : subroutine setup_for_bicubic_interpolations
295 :
296 : include 'formats'
297 :
298 0 : call find_location
299 :
300 : ! set factors for interpolation
301 :
302 0 : hx=x1-x0
303 0 : hxi=1.0d0/hx
304 0 : hx2=hx*hx
305 :
306 0 : xp=(xget-x0)*hxi
307 0 : xpi=1.0d0-xp
308 0 : xp2=xp*xp
309 0 : xpi2=xpi*xpi
310 :
311 0 : cx=xp*(xp2-1.0d0)
312 0 : cxi=xpi*(xpi2-1.0d0)
313 0 : cxd=3.0d0*xp2-1.0d0
314 0 : cxdi=-3.0d0*xpi2+1.0d0
315 :
316 0 : hy=y1-y0
317 0 : hyi=1.0d0/hy
318 0 : hy2=hy*hy
319 :
320 0 : yp=(yget-y0)*hyi
321 0 : ypi=1.0d0-yp
322 0 : yp2=yp*yp
323 0 : ypi2=ypi*ypi
324 :
325 0 : cy=yp*(yp2-1.0d0)
326 0 : cyi=ypi*(ypi2-1.0d0)
327 0 : cyd=3.0d0*yp2-1.0d0
328 0 : cydi=-3.0d0*ypi2+1.0d0
329 :
330 : if (dbg) then
331 : write(*,2) 'T9', ix, x0, T9, x1
332 : write(*,2) 'lYeRho', jy, y0, lYeRho, y1
333 : write(*,1) 'xpi', xpi
334 : write(*,1) 'ypi', ypi
335 : write(*,'(A)')
336 : end if
337 :
338 0 : end subroutine setup_for_bicubic_interpolations
339 :
340 0 : subroutine do_bicubic_interpolations(fin, fval, df_dx, df_dy, ierr)
341 : ! derived from routines in the PSPLINE package written by Doug McCune
342 : real(dp), dimension(:,:,:) :: fin ! the spline data array, dimensions (4, nx, ny)
343 : real(dp), intent(out) :: fval, df_dx, df_dy
344 : integer, intent(out) :: ierr
345 :
346 : real(dp), parameter :: one_sixth = 1d0/6d0
347 : real(dp), parameter :: z36th = 1d0/36d0
348 :
349 0 : ierr = 0
350 :
351 : ! bicubic spline interpolation
352 : fval = &
353 : xpi*( &
354 : ypi*fin(1,ix,jy) +yp*fin(1,ix,jy+1)) &
355 : +xp*(ypi*fin(1,ix+1,jy)+yp*fin(1,ix+1,jy+1)) &
356 : +one_sixth*hx2*( &
357 : cxi*(ypi*fin(2,ix,jy) +yp*fin(2,ix,jy+1))+ &
358 : cx*(ypi*fin(2,ix+1,jy)+yp*fin(2,ix+1,jy+1))) &
359 : +one_sixth*hy2*( &
360 : xpi*(cyi*fin(3,ix,jy) +cy*fin(3,ix,jy+1))+ &
361 : xp*(cyi*fin(3,ix+1,jy)+cy*fin(3,ix+1,jy+1))) &
362 : +z36th*hx2*hy2*( &
363 : cxi*(cyi*fin(4,ix,jy) +cy*fin(4,ix,jy+1))+ &
364 0 : cx*(cyi*fin(4,ix+1,jy)+cy*fin(4,ix+1,jy+1)))
365 :
366 : include 'formats'
367 : if (dbg) then
368 : write(*,1) 'fin(1,ix,jy)', fin(1,ix,jy)
369 : write(*,1) 'fin(1,ix,jy+1)', fin(1,ix,jy+1)
370 : write(*,1) 'fin(1,ix+1,jy)', fin(1,ix+1,jy)
371 : write(*,1) 'fin(1,ix+1,jy+1)', fin(1,ix+1,jy+1)
372 : write(*,1) 'fval', fval
373 :
374 : write(*,'(A)')
375 : call mesa_error(__FILE__,__LINE__,'debug: do_bicubic_interpolations')
376 : end if
377 :
378 : ! derivatives of bicubic splines
379 : df_dx = &
380 : hxi*( &
381 : -(ypi*fin(1,ix,jy) +yp*fin(1,ix,jy+1)) &
382 : +(ypi*fin(1,ix+1,jy)+yp*fin(1,ix+1,jy+1))) &
383 : +one_sixth*hx*( &
384 : cxdi*(ypi*fin(2,ix,jy) +yp*fin(2,ix,jy+1))+ &
385 : cxd*(ypi*fin(2,ix+1,jy)+yp*fin(2,ix+1,jy+1))) &
386 : +one_sixth*hxi*hy2*( &
387 : -(cyi*fin(3,ix,jy) +cy*fin(3,ix,jy+1)) &
388 : +(cyi*fin(3,ix+1,jy)+cy*fin(3,ix+1,jy+1))) &
389 : +z36th*hx*hy2*( &
390 : cxdi*(cyi*fin(4,ix,jy) +cy*fin(4,ix,jy+1))+ &
391 0 : cxd*(cyi*fin(4,ix+1,jy)+cy*fin(4,ix+1,jy+1)))
392 :
393 : df_dy = &
394 : hyi*( &
395 : xpi*(-fin(1,ix,jy) +fin(1,ix,jy+1))+ &
396 : xp*(-fin(1,ix+1,jy)+fin(1,ix+1,jy+1))) &
397 : +one_sixth*hx2*hyi*( &
398 : cxi*(-fin(2,ix,jy) +fin(2,ix,jy+1))+ &
399 : cx*(-fin(2,ix+1,jy)+fin(2,ix+1,jy+1))) &
400 : +one_sixth*hy*( &
401 : xpi*(cydi*fin(3,ix,jy) +cyd*fin(3,ix,jy+1))+ &
402 : xp*(cydi*fin(3,ix+1,jy)+cyd*fin(3,ix+1,jy+1))) &
403 : +z36th*hx2*hy*( &
404 : cxi*(cydi*fin(4,ix,jy) +cyd*fin(4,ix,jy+1))+ &
405 0 : cx*(cydi*fin(4,ix+1,jy)+cyd*fin(4,ix+1,jy+1)))
406 :
407 0 : end subroutine do_bicubic_interpolations
408 :
409 10060 : subroutine setup_for_linear_interp
410 : include 'formats'
411 :
412 10060 : call find_location
413 :
414 10060 : dT9 = T9 - x0
415 10060 : delta_T9 = x1 - x0
416 10060 : x_beta = dT9 / delta_T9 ! fraction of x1 result
417 10060 : x_alfa = 1.0d0 - x_beta ! fraction of x0 result
418 10060 : if (x_alfa < 0 .or. x_alfa > 1) then
419 0 : write(*,1) 'weaklib: x_alfa', x_alfa
420 0 : write(*,1) 'T9', T9
421 0 : write(*,1) 'x0', x0
422 0 : write(*,1) 'x1', x1
423 0 : call mesa_error(__FILE__,__LINE__)
424 : end if
425 :
426 10060 : dlYeRho = lYeRho - y0
427 10060 : delta_lYeRho = y1 - y0
428 10060 : y_beta = dlYeRho / delta_lYeRho ! fraction of y1 result
429 10060 : y_alfa = 1 - y_beta ! fraction of y0 result
430 10060 : if (is_bad(y_alfa) .or. y_alfa < 0 .or. y_alfa > 1) then
431 0 : write(*,1) 'weaklib: y_alfa', y_alfa
432 0 : write(*,1) 'T9', T9
433 0 : write(*,1) 'x0', x0
434 0 : write(*,1) 'dT9', dT9
435 0 : write(*,1) 'delta_T9', delta_T9
436 0 : write(*,1) 'lYeRho', lYeRho
437 0 : write(*,1) 'y0', y0
438 0 : write(*,1) 'dlYeRho', dlYeRho
439 0 : write(*,1) 'y1', y1
440 0 : write(*,1) 'delta_lYeRho', delta_lYeRho
441 0 : write(*,1) 'y_beta', y_beta
442 : !call mesa_error(__FILE__,__LINE__,'weak setup_for_linear_interp')
443 : end if
444 :
445 : if (dbg) then
446 : write(*,2) 'T9', ix, x0, T9, x1
447 : write(*,2) 'lYeRho', jy, y0, lYeRho, y1
448 : write(*,1) 'x_alfa, x_beta', x_alfa, x_beta
449 : write(*,1) 'y_alfa, y_beta', y_alfa, y_beta
450 : write(*,'(A)')
451 : end if
452 :
453 10060 : end subroutine setup_for_linear_interp
454 :
455 30180 : subroutine do_linear_interp(f, fval, df_dx, df_dy, ierr)
456 : use interp_1d_lib
457 : use utils_lib, only: is_bad
458 : real(dp), dimension(:,:,:) :: f ! (4, nx, ny)
459 : real(dp), intent(out) :: fval, df_dx, df_dy
460 : integer, intent(out) :: ierr
461 :
462 30180 : real(dp) :: fx0, fx1, fy0, fy1
463 :
464 : include 'formats'
465 :
466 30180 : ierr = 0
467 :
468 30180 : fx0 = y_alfa*f(1,ix,jy) + y_beta*f(1,ix,jy+1)
469 30180 : fx1 = y_alfa*f(1,ix+1,jy) + y_beta*f(1,ix+1,jy+1)
470 :
471 30180 : fy0 = x_alfa*f(1,ix,jy) + x_beta*f(1,ix+1,jy)
472 30180 : fy1 = x_alfa*f(1,ix,jy+1) + x_beta*f(1,ix+1,jy+1)
473 :
474 30180 : fval = x_alfa*fx0 + x_beta*fx1
475 30180 : df_dx = (fx1 - fx0)/(x1 - x0)
476 30180 : df_dy = (fy1 - fy0)/(y1 - y0)
477 :
478 30180 : if (is_bad(fval)) then
479 0 : ierr = -1
480 : return
481 :
482 : write(*,1) 'x_alfa', x_alfa
483 : write(*,1) 'x_beta', x_beta
484 : write(*,1) 'fx0', fx0
485 : write(*,1) 'fx1', fx1
486 : write(*,1) 'y_alfa', y_alfa
487 : write(*,1) 'y_beta', y_beta
488 : write(*,1) 'f(1,ix,jy)', f(1,ix,jy)
489 : write(*,1) 'f(1,ix,jy+1)', f(1,ix,jy+1)
490 : !call mesa_error(__FILE__,__LINE__,'weak do_linear_interp')
491 : end if
492 :
493 30180 : end subroutine do_linear_interp
494 :
495 : end subroutine interpolate_weaklib_table
496 :
497 :
498 6031 : subroutine deallocate_weaklib_rate_table(self)
499 : type(weaklib_rate_table), intent(inout) :: self
500 6031 : if(allocated(self % T9s)) deallocate(self % T9s)
501 6031 : if(allocated(self % lYeRhos)) deallocate(self % lYeRhos)
502 6031 : if(allocated(self % data)) deallocate(self % data)
503 6031 : end subroutine deallocate_weaklib_rate_table
504 :
505 0 : subroutine deallocate_weaklib_rate_table_array(self)
506 : type(weaklib_rate_table), intent(inout) :: self(:)
507 : integer :: i
508 0 : do i = 1, size(self)
509 0 : call deallocate_weaklib_rate_table(self(i))
510 : end do
511 0 : end subroutine deallocate_weaklib_rate_table_array
512 :
513 :
514 20406 : end module weaklib_tables
|