Shamrock 2025.10.0
Astrophysical Code
Loading...
Searching...
No Matches
generate_tabflux_tabintflux.py
1import sys
2
3import numpy as np
4
5# from progressbar import Bar,Percentage,ProgressBar
6# from numba import jit, njit, prange
8from .utils_polynomials import *
9
10
11def compute_coagtabflux_k0(
12 kernel, K0, Q, vecnodes, vecweights, nbins, massgrid, mat_coeffs_leg, tensor_tabflux_coag
13):
14 """
15 Function to precompute array depending only on massgrid to evaluate the coagulation flux
16 DG scheme with piecewise constant approximation
17
18 Parameters
19 ----------
20 kernel : scalar, type -> integer
21 select the collisional kernel function
22 K0 : scalar, type -> float
23 constant value of the kernel function (used to adapt to code unit)
24 Q : scalar, type -> integer
25 number of points for Gauss-Legendre quadrature
26 vecnodes : 1D array (dim = Q), type -> float
27 nodes of the Legendre polynomials
28 vecweights : 1D array (dim = Q), type -> float
29 weights coefficients for the Gauss-Legendre polynomials
30 nbins : scalar, type -> integer
31 number of dust bins
32 massgrid : 1D array (dim = nbins+1), type -> float
33 grid of masses given borders value of mass bins
34 mat_coeffs_leg : 2D array (dmin = (kpol+1,kpol+1)), type -> float
35 array containing on each line Legendre polynomial coefficients from degree 0 to kpol inclusive
36 on each line coefficients are ordered from low to high orders
37 tensor_tabflux_coag : 3D array (dim = (nbins,nbins,nbins)), type -> float
38 array to evaluate coagulation flux
39
40 Returns
41 -------
42 filled array tensor_tabflux_coag
43
44 """
45
46 # display progress bar
47 # bar = ProgressBar(widgets=[Percentage(), Bar()], maxval=nbins).start()
48
49 for j in range(nbins):
50 for lp in range(j + 1):
51 for l in range(nbins):
52 res = coagfluxfunction(
53 kernel,
54 K0,
55 Q,
56 vecnodes,
57 vecweights,
58 nbins,
59 massgrid,
60 mat_coeffs_leg,
61 j,
62 lp,
63 l,
64 0,
65 0,
66 )
67
68 if res != 0.0:
69 tensor_tabflux_coag[j, lp, l] = res
70
71 # bar.update(j+1)
72
73 # bar.finish()
74
75
76
77def compute_coagtabflux_k0_numba(
78 kernel, K0, Q, vecnodes, vecweights, nbins, massgrid, mat_coeffs_leg, tensor_tabflux_coag
79): # ,progress):
80 """
81 Function to precompute array depending only on massgrid to evaluate the coagulation flux
82 DG scheme with piecewise constant approximation
83 Numba formalism
84
85 Parameters
86 ----------
87 kernel : scalar, type -> integer
88 select the collisional kernel function
89 K0 : scalar, type -> float
90 constant value of the kernel function (used to adapt to code unit)
91 Q : scalar, type -> integer
92 number of points for Gauss-Legendre quadrature
93 vecnodes : 1D array (dim = Q), type -> float
94 nodes of the Legendre polynomials
95 vecweights : 1D array (dim = Q), type -> float
96 weights coefficients for the Gauss-Legendre polynomials
97 nbins : scalar, type -> integer
98 number of dust bins
99 massgrid : 1D array (dim = nbins+1), type -> float
100 grid of masses given borders value of mass bins
101 mat_coeffs_leg : 2D array (dmin = (kpol+1,kpol+1)), type -> float
102 array containing on each line Legendre polynomial coefficients from degree 0 to kpol inclusive
103 on each line coefficients are ordered from low to high orders
104 tensor_tabflux_coag : 3D array (dim = (nbins,nbins,nbins)), type -> float
105 array to evaluate coagulation flux
106 progress : used for numba progress bar
107
108 Returns
109 -------
110 filled array tensor_tabflux_coag
111
112 """
113
114 for j in range(nbins):
115 # progress.update(1) # update on each outer loop iteration
116 for lp in range(j + 1):
117 for l in range(nbins):
118 res = coagfluxfunction_numba(
119 kernel,
120 K0,
121 Q,
122 vecnodes,
123 vecweights,
124 nbins,
125 massgrid,
126 mat_coeffs_leg,
127 j,
128 lp,
129 l,
130 0,
131 0,
132 )
133
134 if res != 0.0:
135 tensor_tabflux_coag[j, lp, l] = res
136
137
138def compute_coagtabflux(
139 kernel, K0, Q, vecnodes, vecweights, nbins, kpol, massgrid, mat_coeffs_leg, tensor_tabflux_coag
140):
141 """
142 Function to precompute array depending only on massgrid to evaluate the coagulation flux
143 DG scheme with piecewise polynomial approximation
144
145 Parameters
146 ----------
147 kernel : scalar, type -> integer
148 select the collisional kernel function
149 K0 : scalar, type -> float
150 constant value of the kernel function (used to adapt to code unit)
151 Q : scalar, type -> integer
152 number of points for Gauss-Legendre quadrature
153 vecnodes : 1D array (dim = Q), type -> float
154 nodes of the Legendre polynomials
155 vecweights : 1D array (dim = Q), type -> float
156 weights coefficients for the Gauss-Legendre polynomials
157 nbins : scalar, type -> integer
158 number of dust bins
159 kpol : scalar, type -> integer
160 degree of polynomials for approximation
161 massgrid : 1D array (dim = nbins+1), type -> float
162 grid of masses given borders value of mass bins
163 mat_coeffs_leg : 2D array (dmin = (kpol+1,kpol+1)), type -> float
164 array containing on each line Legendre polynomial coefficients from degree 0 to kpol inclusive
165 on each line coefficients are ordered from low to high orders
166 tensor_tabflux_coag : 5D array (dim = (nbins,nbins,nbins,kpol+1,kpol+1)), type -> float
167 array to evaluate coagulation flux
168
169
170 Returns
171 -------
172 filled array tensor_tabflux_coag
173
174 """
175
176 # display progress bar
177 # bar = ProgressBar(widgets=[Percentage(), Bar()], maxval=nbins).start()
178
179 for j in range(nbins):
180 for lp in range(j + 1):
181 for l in range(nbins):
182 for ip in range(kpol + 1):
183 for i in range(kpol + 1):
184 res = coagfluxfunction(
185 kernel,
186 K0,
187 Q,
188 vecnodes,
189 vecweights,
190 nbins,
191 massgrid,
192 mat_coeffs_leg,
193 j,
194 lp,
195 l,
196 ip,
197 i,
198 )
199
200 if res != 0.0:
201 tensor_tabflux_coag[j, lp, l, ip, i] = res
202
203 # bar.update(j+1)
204
205 # bar.finish()
206
207
208# @njit
209def compute_coagtabflux_numba(
210 kernel, K0, Q, vecnodes, vecweights, nbins, kpol, massgrid, mat_coeffs_leg, tensor_tabflux_coag
211): # ,progress):
212 """
213 Function to precompute array depending only on massgrid to evaluate the coagulation flux
214 DG scheme with piecewise polynomial approximation
215 Numba formalism
216
217 Parameters
218 ----------
219 kernel : scalar, type -> integer
220 select the collisional kernel function
221 K0 : scalar, type -> float
222 constant value of the kernel function (used to adapt to code unit)
223 Q : scalar, type -> integer
224 number of points for Gauss-Legendre quadrature
225 vecnodes : 1D array (dim = Q), type -> float
226 nodes of the Legendre polynomials
227 vecweights : 1D array (dim = Q), type -> float
228 weights coefficients for the Gauss-Legendre polynomials
229 nbins : scalar, type -> integer
230 number of dust bins
231 kpol : scalar, type -> integer
232 degree of polynomials for approximation
233 massgrid : 1D array (dim = nbins+1), type -> float
234 grid of masses given borders value of mass bins
235 mat_coeffs_leg : 2D array (dmin = (kpol+1,kpol+1)), type -> float
236 array containing on each line Legendre polynomial coefficients from degree 0 to kpol inclusive
237 on each line coefficients are ordered from low to high orders
238 tensor_tabflux_coag : 5D array (dim = (nbins,nbins,nbins,kpol+1,kpol+1)), type -> float
239 array to evaluate coagulation flux
240
241
242 Returns
243 -------
244 filled array tensor_tabflux_coag
245
246 """
247
248 for j in range(nbins):
249 # progress.update(1) # update on each outer loop iteration
250 for lp in range(j + 1):
251 for l in range(nbins):
252 for ip in range(kpol + 1):
253 for i in range(kpol + 1):
254 res = coagfluxfunction_numba(
255 kernel,
256 K0,
257 Q,
258 vecnodes,
259 vecweights,
260 nbins,
261 massgrid,
262 mat_coeffs_leg,
263 j,
264 lp,
265 l,
266 ip,
267 i,
268 )
269
270 if res != 0.0:
271 tensor_tabflux_coag[j, lp, l, ip, i] = res
272
273
274def compute_coagtabintflux(
275 kernel,
276 K0,
277 Q,
278 vecnodes,
279 vecweights,
280 nbins,
281 kpol,
282 massgrid,
283 mat_coeffs_leg,
284 tensor_tabintflux_coag,
285):
286 """
287 Function to precompute array depending only on massgrid to evaluate the term including the integral of coagulation flux
288 DG scheme with piecewise polynomial approximation
289
290 Parameters
291 ----------
292 kernel : scalar, type -> integer
293 select the collisional kernel function
294 K0 : scalar, type -> float
295 constant value of the kernel function (used to adapt to code unit)
296 Q : scalar, type -> integer
297 number of points for Gauss-Legendre quadrature
298 vecnodes : 1D array (dim = Q), type -> float
299 nodes of the Legendre polynomials
300 vecweights : 1D array (dim = Q), type -> float
301 weights coefficients for the Gauss-Legendre polynomials
302 nbins : scalar, type -> integer
303 number of dust bins
304 kpol : scalar, type -> integer
305 degree of polynomials for approximation
306 massgrid : 1D array (dim = nbins+1), type -> float
307 grid of masses given borders value of mass bins
308 mat_coeffs_leg : 2D array (dmin = (kpol+1,kpol+1)), type -> float
309 array containing on each line Legendre polynomial coefficients from degree 0 to kpol inclusive
310 on each line coefficients are ordered from low to high orders
311 tensor_tabintflux_coag : 6D array (dim = (nbins,kpol+1,nbins,nbins,kpol+1,kpol+1)), type -> float
312 array to evaluate the term including the integral of the coagulation flux
313
314
315 Returns
316 -------
317 filled array tensor_tabintflux_coag
318
319 """
320
321 # bar = ProgressBar(widgets=[Percentage(), Bar()], maxval=nbins*(kpol+1)).start()
322
323 for j in range(nbins):
324 for k in range(1, kpol + 1):
325 for lp in range(j + 1):
326 for l in range(nbins):
327 for ip in range(kpol + 1):
328 for i in range(kpol + 1):
329 res = coagintfluxfunction(
330 kernel,
331 K0,
332 Q,
333 vecnodes,
334 vecweights,
335 nbins,
336 massgrid,
337 mat_coeffs_leg,
338 j,
339 k,
340 lp,
341 l,
342 ip,
343 i,
344 )
345
346 if res != 0.0:
347 tensor_tabintflux_coag[j, k, lp, l, ip, i] = res
348
349 # bar.update(j+1)
350
351 # bar.finish()
352
353
354# @njit
355def compute_coagtabintflux_numba(
356 kernel,
357 K0,
358 Q,
359 vecnodes,
360 vecweights,
361 nbins,
362 kpol,
363 massgrid,
364 mat_coeffs_leg,
365 tensor_tabintflux_coag,
366): # ,progress):
367 """
368 Function to precompute array depending only on massgrid to evaluate the term including the integral of coagulation flux
369 DG scheme with piecewise polynomial approximation
370 Numba formalism
371
372 Parameters
373 ----------
374 kernel : scalar, type -> integer
375 select the collisional kernel function
376 K0 : scalar, type -> float
377 constant value of the kernel function (used to adapt to code unit)
378 Q : scalar, type -> integer
379 number of points for Gauss-Legendre quadrature
380 vecnodes : 1D array (dim = Q), type -> float
381 nodes of the Legendre polynomials
382 vecweights : 1D array (dim = Q), type -> float
383 weights coefficients for the Gauss-Legendre polynomials
384 nbins : scalar, type -> integer
385 number of dust bins
386 kpol : scalar, type -> integer
387 degree of polynomials for approximation
388 massgrid : 1D array (dim = nbins+1), type -> float
389 grid of masses given borders value of mass bins
390 mat_coeffs_leg : 2D array (dmin = (kpol+1,kpol+1)), type -> float
391 array containing on each line Legendre polynomial coefficients from degree 0 to kpol inclusive
392 on each line coefficients are ordered from low to high orders
393 tensor_tabintflux_coag : 6D array (dim = (nbins,kpol+1,nbins,nbins,kpol+1,kpol+1)), type -> float
394 array to evaluate the term including the integral of the coagulation flux
395 progress : used for numba progress bar
396
397
398 Returns
399 -------
400 filled array tensor_tabintflux_coag
401
402
403 """
404
405 for j in range(nbins):
406 # progress.update(1) # update on each outer loop iteration
407 for k in range(1, kpol + 1):
408 for lp in range(j + 1):
409 for l in range(nbins):
410 for ip in range(kpol + 1):
411 for i in range(kpol + 1):
412 res = coagintfluxfunction_numba(
413 kernel,
414 K0,
415 Q,
416 vecnodes,
417 vecweights,
418 nbins,
419 massgrid,
420 mat_coeffs_leg,
421 j,
422 k,
423 lp,
424 l,
425 ip,
426 i,
427 )
428
429 if res != 0.0:
430 tensor_tabintflux_coag[j, k, lp, l, ip, i] = res