MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_dust.t
Go to the documentation of this file.
1!> Module for including dust species, which interact with the gas through a drag
2!> force
3module mod_dust
4 use mod_global_parameters, only: std_len
6 use mod_comm_lib, only: mpistop
7
8 implicit none
9 private
10
11 !> Reduction of stopping time timestep limit
12 double precision :: dust_dtpar = 0.5d0
13
14 !> Factor used in squared thermal velocity
15 double precision :: gas_vtherm_factor = 3.0d0
16
17 !> Dust temperature in K (if dust_temperature_type is constant)
18 double precision :: dust_temperature = -1.0d0
19
20 !> Dust drag coefficient for linear drag (for testing dust_method=linear)
21 double precision :: dust_K_lineardrag = -1.0d0
22
23 !> If dust_temperature_type is stellar, it will be calculated according to Tielens (2005),
24 !> eqn. 5.44 using an input stellar luminosity in solar luminosities
25 double precision :: dust_stellar_luminosity = -1.0d0
26
27 !> Minimum dust density as used when dust_small_to_zero=T
28 double precision, public, protected :: dust_min_rho = -1.0d0
29
30 !> Size of each dust species, dimensionless expression
31 double precision, allocatable, public :: dust_size(:)
32
33 !> Internal density of each dust species, dimensionless expression
34 double precision, allocatable, public :: dust_density(:)
35
36 !> The number of dust species
37 integer, public, protected :: dust_n_species = 0
38
39 integer, protected :: gas_rho_ = -1
40 integer, allocatable, protected :: gas_mom(:)
41 integer, protected :: gas_e_ = -1
42
43 !> Indices of the dust densities
44 integer, allocatable, public, protected :: dust_rho(:)
45
46 !> Indices of the dust momentum densities
47 integer, allocatable, public, protected :: dust_mom(:, :)
48
49 !> Set small dust densities to zero to avoid numerical problems
50 logical, public, protected :: dust_small_to_zero = .false.
51
52 !> Adding dust in sourcesplit manner or not
53 logical :: dust_source_split = .false.
54
55 !> This can be turned off for testing purposes, if F then gas uncouples from dust
56 logical :: dust_backreaction = .true.
57
58 !> whether second order terms (relevant only when dust_n_species >=2) are included
59 !> there are the terms n2, ni2, d2 in Eqs 6,7,8 in amrvac 3.0 paper
60 logical :: dust_implicit_second_order = .true.
61
62 !> whether fh is added for gas energy: is only added in the impliict implementation, the explicit one was left as before
63 logical :: dust_backreaction_fh = .false.
64
65 !> What type of dust drag force to use. Can be 'Kwok', 'sticking', 'linear', 'usr' or 'none'.
66 character(len=std_len), public, protected :: dust_method = 'Kwok'
67
68 !> Can be 'graphite' or 'silicate', affects the dust temperature
69 character(len=std_len) :: dust_species = 'graphite'
70
71 !> Determines the dust temperature, can be 'constant', 'ism', or 'stellar'
72 character(len=std_len) :: dust_temperature_type = 'constant'
73
74
75 ! Public methods
76 public :: dust_init
77 public :: dust_get_dt
78 public :: dust_get_flux
79 public :: dust_get_cmax
80 public :: dust_get_flux_prim
81 public :: dust_get_cmax_prim
82 public :: dust_add_source
83 public :: dust_to_conserved
84 public :: dust_to_primitive
85 public :: dust_check_params
86 public :: dust_check_w
87 public :: set_dusttozero
88 public :: dust_implicit_update
90
91
92contains
93
94 subroutine dust_init(g_rho, g_mom, g_energy)
96
97 integer, intent(in) :: g_rho
98 integer, intent(in) :: g_mom(ndir)
99 integer, intent(in) :: g_energy ! Negative value if not present
100 integer :: n, idir
101 character(len=2) :: dim
102
103 call dust_read_params(par_files)
104
105 if(dust_source_split) any_source_split=.true.
106 allocate(gas_mom(ndir))
107 gas_rho_ = g_rho
108 gas_mom = g_mom
109 gas_e_ = g_energy
110
111 allocate(dust_size(dust_n_species))
113 dust_size(:) = -1.0d0
114 dust_density(:) = -1.0d0
115
116 allocate(dust_rho(dust_n_species))
117 allocate(dust_mom(ndir, dust_n_species))
118
119 ! Set index of dust densities
120 do n = 1, dust_n_species
121 dust_rho(n) = var_set_fluxvar("rhod", "rhod", n)
122 end do
123
124 ! Dust momentum
125 do idir = 1, ndir
126 write(dim, "(I0,A)") idir, "d"
127 do n = 1, dust_n_species
128 dust_mom(idir, n) = var_set_fluxvar("m"//dim, "v"//dim, n)
129 end do
130 end do
131
132 end subroutine dust_init
133
134 !> Read dust_list module parameters from a file
135 subroutine dust_read_params(files)
137 character(len=*), intent(in) :: files(:)
138 integer :: n
139
140 namelist /dust_list/ dust_n_species, dust_min_rho, dust_method, &
141 dust_k_lineardrag, dust_small_to_zero, dust_source_split, dust_temperature, &
142 dust_temperature_type, dust_backreaction, dust_dtpar, gas_vtherm_factor, dust_stellar_luminosity,&
143 dust_implicit_second_order, dust_backreaction_fh
144
145 do n = 1, size(files)
146 open(unitpar, file=trim(files(n)), status="old")
147 read(unitpar, dust_list, end=111)
148111 close(unitpar)
149 end do
150
151 end subroutine dust_read_params
152
153 subroutine dust_check_params()
156
157 if (dust_method == 'sticking') then
158 if (si_unit) call mpistop("Dust error: sticking assumes cgs units")
159 if (dust_temperature_type == "constant") then
160 if (dust_temperature < 0.0d0) then
161 call mpistop("Dust error: dust_temperature (in K) < 0 or not set")
162 end if
163 else if (dust_temperature_type == "stellar") then
164 if (dust_stellar_luminosity < 0.0d0) then
165 call mpistop("Dust error: dust_stellar_luminosity (in solar) < 0 or not set")
166 end if
167 end if
168 end if
169
170 if (dust_method == 'linear') then
171 if(dust_k_lineardrag<0.0d0) then
172 call mpistop("With dust_method=='linear', you must set a positive dust_K_lineardrag")
173 end if
174 end if
175
176 if (any(dust_size < 0.0d0)) &
177 call mpistop("Dust error: any(dust_size < 0) or not set")
178 if (any(dust_density < 0.0d0)) &
179 call mpistop("Dust error: any(dust_density < 0) or not set")
180
181 if (dust_method == 'usr') then
182 if (.not. associated(usr_get_3d_dragforce) .or. .not. associated(usr_dust_get_dt)) &
183 call mpistop("Dust error:usr_get_3d_dragforce and usr_dust_get_dt not defined")
184 end if
185
186 if(.not. use_imex_scheme .and. ((dust_dtpar .ge. 1d0).or.(dust_dtpar.le.0))) then
187 if(mype .eq. 0) print*, "EXPLICIT source for dust requires 0<dt_dustpar < 1, set to 0.8"
188 dust_dtpar = 0.8
189 endif
190
191 end subroutine dust_check_params
192
193 subroutine dust_check_w(ixI^L,ixO^L,w,flag)
195
196 integer, intent(in) :: ixi^l,ixo^l
197 double precision, intent(in):: w(ixi^s,1:nw)
198 logical, intent(inout) :: flag(ixi^s,1:nw)
199 integer :: n
200
201 do n = 1, dust_n_species
202 flag(ixo^s,dust_rho(n))=(w(ixo^s,dust_rho(n))<0.0d0)
203 enddo
204
205 end subroutine dust_check_w
206
207 subroutine dust_to_conserved(ixI^L, ixO^L, w, x)
209
210 integer, intent(in) :: ixi^l, ixo^l
211 double precision, intent(inout) :: w(ixi^s, 1:nw)
212 double precision, intent(in) :: x(ixi^s, 1:ndim)
213 integer :: n, idir
214
215 if(fix_small_values .and. dust_small_to_zero) call set_dusttozero(ixi^l, ixo^l, w, x)
216
217 do n = 1, dust_n_species
218 ! Convert velocity to momentum
219 do idir = 1, ndir
220 w(ixo^s, dust_mom(idir, n)) = w(ixo^s, dust_rho(n)) * &
221 w(ixo^s, dust_mom(idir, n))
222 end do
223 end do
224
225 end subroutine dust_to_conserved
226
227 subroutine dust_to_primitive(ixI^L, ixO^L, w, x)
229
230 integer, intent(in) :: ixi^l, ixo^l
231 double precision, intent(inout) :: w(ixi^s, 1:nw)
232 double precision, intent(in) :: x(ixi^s, 1:ndim)
233 integer :: n, idir
234
235 do n = 1, dust_n_species
236 ! Convert momentum to velocity
237 do idir = 1, ndir
238 where (w(ixo^s, dust_rho(n)) > 0.0d0)
239 w(ixo^s, dust_mom(idir, n)) = w(ixo^s, dust_mom(idir, n)) / &
240 w(ixo^s, dust_rho(n))
241 elsewhere
242 w(ixo^s, dust_mom(idir, n)) = 0.0d0
243 end where
244 end do
245 end do
246
247 if(fix_small_values .and. dust_small_to_zero) call set_dusttozero(ixi^l, ixo^l, w, x)
248
249 end subroutine dust_to_primitive
250
251 subroutine dust_get_flux(w, x, ixI^L, ixO^L, idim, f)
253
254 integer, intent(in) :: ixi^l, ixo^l, idim
255 double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
256 double precision, intent(inout) :: f(ixi^s, nwflux)
257 integer :: n, idir
258
259 do n = 1, dust_n_species
260 where (w(ixo^s, dust_rho(n)) > 0.0d0)
261 f(ixo^s, dust_rho(n)) = w(ixo^s, dust_mom(idim, n))
262 elsewhere
263 f(ixo^s, dust_rho(n)) = 0.0d0
264 end where
265
266 do idir = 1, ndir
267 f(ixo^s, dust_mom(idir, n)) = w(ixo^s, dust_mom(idir, n)) * &
268 get_vdust(w, ixi^l, ixo^l, idim, n)
269 end do
270 end do
271
272 end subroutine dust_get_flux
273
274 subroutine dust_get_flux_prim(w, x, ixI^L, ixO^L, idim, f)
276
277 integer, intent(in) :: ixi^l, ixo^l, idim
278 double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
279 double precision, intent(inout) :: f(ixi^s, nwflux)
280 integer :: n, idir
281
282 do n = 1, dust_n_species
283 where (w(ixo^s, dust_rho(n)) > 0.0d0)
284 f(ixo^s, dust_rho(n)) = w(ixo^s, dust_mom(idim, n))*w(ixo^s, dust_rho(n))
285 elsewhere
286 f(ixo^s, dust_rho(n)) = 0.0d0
287 end where
288
289 do idir = 1, ndir
290 f(ixo^s, dust_mom(idir, n)) = w(ixo^s, dust_mom(idir, n)) * &
291 w(ixo^s, dust_rho(n)) * get_vdust_prim(w, ixi^l, ixo^l, idim, n)
292 end do
293 end do
294
295 end subroutine dust_get_flux_prim
296
297 function get_vdust(w, ixI^L, ixO^L, idim, n) result(vdust)
298 use mod_global_parameters, only: nw
299 integer, intent(in) :: ixi^l, ixo^l, idim, n
300 double precision, intent(in) :: w(ixi^s, nw)
301 double precision :: vdust(ixo^s)
302
303 where (w(ixo^s, dust_rho(n)) > 0.0d0)
304 vdust(ixo^s) = w(ixo^s, dust_mom(idim, n)) / w(ixo^s, dust_rho(n))
305 elsewhere
306 vdust(ixo^s) = 0.0d0
307 end where
308
309 end function get_vdust
310
311 function get_vdust_prim(w, ixI^L, ixO^L, idim, n) result(vdust)
312 use mod_global_parameters, only: nw
313 integer, intent(in) :: ixi^l, ixo^l, idim, n
314 double precision, intent(in) :: w(ixi^s, nw)
315 double precision :: vdust(ixo^s)
316
317 where (w(ixo^s, dust_rho(n)) > 0.0d0)
318 vdust(ixo^s) = w(ixo^s, dust_mom(idim, n))
319 elsewhere
320 vdust(ixo^s) = 0.0d0
321 end where
322
323 end function get_vdust_prim
324
325 ! Force dust density to zero if dust_rho <= dust_min_rho
326 subroutine set_dusttozero(ixI^L, ixO^L, w, x)
328
329 integer, intent(in) :: ixi^l, ixo^l
330 double precision, intent(in) :: x(ixi^s, 1:ndim)
331 double precision, intent(inout) :: w(ixi^s, 1:nw)
332
333 integer :: n, idir
334 logical :: flag(ixi^s)
335
336 do n = 1, dust_n_species
337 flag(ixo^s)=(w(ixo^s, dust_rho(n)) <= dust_min_rho)
338 where (flag(ixo^s))
339 w(ixo^s, dust_rho(n)) = 0.0d0
340 end where
341 do idir = 1, ndir
342 where (flag(ixo^s))
343 w(ixo^s, dust_mom(idir, n)) = 0.0d0
344 end where
345 end do
346 end do
347
348 end subroutine set_dusttozero
349
350 ! Calculate drag force based on Epstein's law
351 ! From Kwok 1975, page 584 (between eqn 8 and 9)
352 subroutine get_3d_dragforce(ixI^L, ixO^L, w, x, fdrag, ptherm, vgas)
355 integer, intent(in) :: ixi^l, ixo^l
356 double precision, intent(in) :: x(ixi^s, 1:ndim)
357 double precision, intent(in) :: w(ixi^s, 1:nw)
358 double precision, intent(out) :: &
359 fdrag(ixi^s, 1:ndir, 1:dust_n_species)
360 double precision, intent(in) :: ptherm(ixi^s), vgas(ixi^s, 1:ndir)
361
362 double precision, dimension(ixI^S) :: vt2, deltav, fd, vdust
363 double precision :: alpha_t(ixi^s, 1:dust_n_species)
364 integer :: n, idir
365
366 vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/w(ixo^s, gas_rho_)
367
368 select case( trim(dust_method) )
369 case ('Kwok') ! assume sticking coefficient equals 0.25
370
371 do idir = 1, ndir
372 do n = 1, dust_n_species
373 where(w(ixo^s, dust_rho(n)) > 0.0d0)
374 vdust(ixo^s) = w(ixo^s, dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
375 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
376
377 ! 0.75 from sticking coefficient
378 fd(ixo^s) = 0.75d0*w(ixo^s, dust_rho(n))*w(ixo^s, gas_rho_)*deltav(ixo^s) &
379 / (dust_density(n) * dust_size(n))
380
381 ! 0.75 from spherical grainvolume
382 fd(ixo^s) = -fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
383 elsewhere
384 fd(ixo^s) = 0.0d0
385 end where
386 fdrag(ixo^s, idir, n) = fd(ixo^s)
387 end do
388 end do
389
390 case ('sticking') ! Calculate sticking coefficient based on the gas and dust temperatures
391
392 call get_sticking(w, x, ixi^l, ixo^l, alpha_t, ptherm)
393
394 do idir = 1, ndir
395 do n = 1, dust_n_species
396 where(w(ixo^s, dust_rho(n))>0.0d0)
397 vdust(ixo^s) = w(ixo^s,dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
398 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
399 fd(ixo^s) = (one-alpha_t(ixo^s,n)) * w(ixo^s, dust_rho(n))*w(ixo^s, gas_rho_) * &
400 deltav(ixo^s) / (dust_density(n)*dust_size(n))
401 fd(ixo^s) = -fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
402 else where
403 fd(ixo^s) = 0.0d0
404 end where
405 fdrag(ixo^s, idir,n) = fd(ixo^s)
406 end do
407 end do
408
409 case('linear') !linear with Deltav, for testing (see Laibe & Price 2011)
410 do idir = 1, ndir
411 do n = 1, dust_n_species
412 where(w(ixo^s, dust_rho(n))>0.0d0)
413 vdust(ixo^s) = w(ixo^s,dust_mom(idir, n))/w(ixo^s, dust_rho(n))
414 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
415
416 fd(ixo^s) = -dust_k_lineardrag*deltav(ixo^s)
417 else where
418 fd(ixo^s) = 0.0d0
419 end where
420 fdrag(ixo^s, idir,n) = fd(ixo^s)
421 end do
422 end do
423
424 case('usr')
425 call usr_get_3d_dragforce(ixi^l, ixo^l, w, x, fdrag, ptherm, vgas, dust_n_species)
426 case('none')
427 fdrag(ixo^s, :, :) = 0.0d0
428 case default
429 call mpistop( "=== This dust method has not been implemented===" )
430 end select
431
432 end subroutine get_3d_dragforce
433
434 !> Get sticking coefficient alpha_T (always between 0 and 1)
435 !>
436 !> Uses Temperatures in K
437 !> Equation from Decin et al. 2006
438 !> NOTE: w can also be wprim on entry, and hence only identical prim-cons variables must be used
439 subroutine get_sticking(w, x, ixI^L, ixO^L, alpha_T, ptherm)
441 integer, intent(in) :: ixi^l, ixo^l
442 double precision, intent(in) :: x(ixi^s, 1:ndim)
443 double precision, intent(in) :: w(ixi^s, 1:nw)
444 double precision, intent(out) :: alpha_t(ixi^s, 1:dust_n_species)
445 double precision, intent(in) :: ptherm(ixi^s)
446 double precision :: tgas(ixi^s)
447 integer :: n
448
449 ! get the dust species T in K
450 call get_tdust(w, x, ixi^l, ixo^l, alpha_t)
451
452 ! convert dimensionless gas T to K
453 tgas(ixo^s) = (ptherm(ixo^s)/w(ixo^s, gas_rho_))*unit_temperature
454
455 do n = 1, dust_n_species
456 alpha_t(ixo^s,n) = 0.35d0 * dexp(-dsqrt((tgas(ixo^s) + &
457 alpha_t(ixo^s,n))/5.0d2))+0.1d0
458 end do
459
460 end subroutine get_sticking
461
462 !> Returns dust temperature (in K), either as constant or based on equ. 5.41,
463 !> 5.42 and 5.44 from Tielens (2005)
464 !>
465 !> Note that this calculation assumes cgs!!!!
466 !>
467 !> It takes as input the stellar luminosity in solar units in 'stellar' case
468 !> or a fixed dust input temperature in Kelvin when 'constant' or does case 'ism'
469 !> NOTE: w can also be wprim on entry, and hence only identical prim-cons variables must be used
470 subroutine get_tdust(w, x, ixI^L, ixO^L, Td)
472 use mod_geometry
473
474 integer, intent(in) :: ixi^l, ixo^l
475 double precision, intent(in) :: x(ixi^s, 1:ndim)
476 double precision, intent(in) :: w(ixi^s, 1:nw)
477 double precision, intent(out) :: td(ixi^s, 1:dust_n_species)
478 double precision :: g0(ixo^s)
479 integer :: n
480
481 select case( trim(dust_temperature_type) )
482 case( 'constant' )
483 td(ixo^s, :) = dust_temperature
484 case( 'ism' )
485 select case( trim(dust_species) )
486 case( 'graphite' )
487 do n = 1, dust_n_species
488 td(ixo^s, n) = 15.8d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0)
489 end do
490 case( 'silicate' )
491 do n = 1, dust_n_species
492 td(ixo^s, n) = 13.6d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0)
493 end do
494 case default
495 call mpistop( "=== Dust species undetermined===" )
496 end select
497 case( 'stellar' )
498 select case(coordinate)
499 case(spherical)
500 g0(ixo^s) = max(x(ixo^s, 1)*unit_length, smalldouble)
501 !!!case(cylindrical) convert R,Z to spherical radial coordinate r here
502 !!! but only ok for 2D (R,Z) or 2.5D (R,Z) case
503 !!! G0(ixO^S) = max(dsqrt(sum(x(ixO^S,:)**2,dim=ndim+1))*unit_length, smalldouble)
504 case default
505 call mpistop('stellar case not available in this coordinate system')
506 end select
507
508 g0(ixo^s) = 2.1d4*(dust_stellar_luminosity/1.0d8)*((3.0857d17/g0(ixo^s))**2)
509
510 select case( trim(dust_species) )
511 case( 'graphite' )
512 do n = 1, dust_n_species
513 td(ixo^s, n) = 61.0d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0) &
514 *(g0(ixo^s)**(one/5.8d0))
515 end do
516 case( 'silicate' )
517 do n = 1, dust_n_species
518 td(ixo^s, n) = 50.0d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0) &
519 *(g0(ixo^s)**(one/6.0d0))
520 end do
521 case default
522 call mpistop( "=== Dust species undetermined===" )
523 end select
524 case default
525 call mpistop( "=== Dust temperature undetermined===" )
526 end select
527
528 end subroutine get_tdust
529
530 !> w[iw]= w[iw]+qdt*S[wCT, x] where S is the source based on wCT within ixO
531 subroutine dust_add_source(qdt, ixI^L, ixO^L, wCT, w, x, qsourcesplit, active)
533
534 integer, intent(in) :: ixi^l, ixo^l
535 double precision, intent(in) :: qdt
536 double precision, intent(in) :: wct(ixi^s, 1:nw), x(ixi^s, 1:ndim)
537 double precision, intent(inout) :: w(ixi^s, 1:nw)
538 logical, intent(in) :: qsourcesplit
539 logical, intent(inout) :: active
540
541 double precision :: ptherm(ixi^s), vgas(ixi^s, 1:ndir)
542 double precision :: fdrag(ixi^s, 1:ndir, 1:dust_n_species)
543 integer :: n, idir
544
545 select case( trim(dust_method) )
546 case( 'none' )
547 !do nothing here
548 case default !all regular dust methods here
549 if (qsourcesplit .eqv. dust_source_split) then
550 active = .true.
551
552 call phys_get_pthermal(wct, x, ixi^l, ixo^l, ptherm)
553 do idir=1,ndir
554 vgas(ixo^s,idir)=wct(ixo^s,gas_mom(idir))/wct(ixo^s,gas_rho_)
555 end do
556
557 call get_3d_dragforce(ixi^l, ixo^l, wct, x, fdrag, ptherm, vgas)
558 fdrag(ixo^s, 1:ndir, 1:dust_n_species) = fdrag(ixo^s, 1:ndir, 1:dust_n_species) * qdt
559
560 do idir = 1, ndir
561
562 do n = 1, dust_n_species
563 if (dust_backreaction) then
564 w(ixo^s, gas_mom(idir)) = w(ixo^s, gas_mom(idir)) + &
565 fdrag(ixo^s, idir, n)
566 if (gas_e_ > 0) then
567 w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + vgas(ixo^s, idir) &
568 * fdrag(ixo^s, idir, n)
569 end if
570 end if
571
572
573 w(ixo^s, dust_mom(idir, n)) = w(ixo^s, dust_mom(idir, n)) - &
574 fdrag(ixo^s, idir, n)
575 end do
576 end do
577
578 endif
579 end select
580
581 end subroutine dust_add_source
582
583 !> inplace update of psa==>F_im(psa)
584 subroutine dust_evaluate_implicit(qtC,psa)
586 type(state), target :: psa(max_blocks)
587 double precision, intent(in) :: qtc
588
589 integer :: iigrid, igrid, level
590
591 !dust_method = 'none' not used
592
593 !$OMP PARALLEL DO PRIVATE(igrid)
594 do iigrid=1,igridstail; igrid=igrids(iigrid);
595 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
596 block=>psa(igrid)
597 call dust_terms(ixg^ll,ixm^ll,psa(igrid)%w,psa(igrid)%x)
598 end do
599 !$OMP END PARALLEL DO
600
601 end subroutine dust_evaluate_implicit
602
603
604
605
606 subroutine dust_terms(ixI^L,ixO^L,w,x)
608 integer, intent(in) :: ixi^l, ixo^l
609 double precision, intent(inout) :: w(ixi^s, 1:nw)
610 double precision, intent(in) :: x(ixi^s,1:ndim)
611
612 double precision :: tmp(ixi^s), vgas(ixi^s, 1:ndir)
613 double precision :: alpha(ixi^s, 1:ndir, 1:dust_n_species)
614 integer :: n, idir
615
616 do idir=1,ndir
617 vgas(ixo^s,idir)=w(ixo^s,gas_mom(idir))/w(ixo^s,gas_rho_)
618 end do
619 call get_alpha_dust(ixi^l, ixo^l, w, vgas,x, alpha)
620 w(ixo^s, gas_e_)=0d0
621 do idir = 1, ndir
622
623 w(ixo^s, gas_mom(idir))=0d0
624 do n = 1, dust_n_species
625 ! contribution for gas momentum
626 tmp(ixo^s) = alpha(ixo^s, idir,n) * ( &
627 w(ixo^s,dust_rho(n)) * w(ixo^s, gas_mom(idir)) - &
628 w(ixo^s,gas_rho_) * w(ixo^s, dust_mom(idir, n)))
629 w(ixo^s, dust_mom(idir, n)) = -tmp(ixo^s)
630 if (dust_backreaction) then
631 w(ixo^s, gas_mom(idir)) = w(ixo^s, gas_mom(idir)) + tmp(ixo^s)
632 if (gas_e_ > 0) then
633 if(dust_backreaction_fh) then
634 where(w(ixo^s,dust_rho(n)) > 0d0)
635 w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + alpha(ixo^s, idir,n) * &
636 (w(ixo^s,gas_rho_) * (w(ixo^s, dust_mom(idir,n))**2/w(ixo^s,dust_rho(n))) - &
637 w(ixo^s,dust_rho(n)) * (w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_)))
638 elsewhere
639 w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + alpha(ixo^s, idir,n) * ( - &
640 w(ixo^s,dust_rho(n)) * (w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_)))
641 endwhere
642 else
643 w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + vgas(ixo^s, idir) &
644 * tmp(ixo^s)
645 end if
646 end if
647 end if
648 end do
649 end do
650 end subroutine dust_terms
651
652 !> Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
653 subroutine dust_implicit_update(dtfactor,qdt,qtC,psb,psa)
655 !use mod_ghostcells_update
656
657 type(state), target :: psa(max_blocks)
658 type(state), target :: psb(max_blocks)
659 double precision, intent(in) :: qdt
660 double precision, intent(in) :: qtc
661 double precision, intent(in) :: dtfactor
662
663 integer :: iigrid, igrid
664
665 !call getbc(global_time,0.d0,psa,1,nw)
666 !$OMP PARALLEL DO PRIVATE(igrid)
667 do iigrid=1,igridstail; igrid=igrids(iigrid);
668 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
669 block=>psa(igrid)
670 call dust_advance_implicit_grid(ixg^ll, ixg^ll, psa(igrid)%w, psb(igrid)%w, psa(igrid)%x, dtfactor,qdt)
671 end do
672 !$OMP END PARALLEL DO
673
674 end subroutine dust_implicit_update
675
676 subroutine dust_advance_implicit_grid(ixI^L, ixO^L, w, wout, x, dtfactor,qdt)
678 integer, intent(in) :: ixi^l, ixo^l
679 double precision, intent(in) :: qdt
680 double precision, intent(in) :: dtfactor
681 double precision, intent(in) :: w(ixi^s,1:nw)
682 double precision, intent(in) :: x(ixi^s,1:ndim)
683 double precision, intent(out) :: wout(ixi^s,1:nw)
684
685 integer :: n, m, idir
686 double precision :: alpha(ixi^s, 1:ndir, 1:dust_n_species)
687 double precision :: tmp(ixi^s),tmp2(ixi^s)
688 double precision :: tmp3(ixi^s)
689 double precision :: vgas(ixi^s, 1:ndir)
690
691
692 do idir = 1, ndir
693 vgas(ixo^s,idir)=w(ixo^s,gas_mom(idir))/w(ixo^s,gas_rho_)
694 end do
695 call get_alpha_dust(ixi^l, ixo^l, w, vgas, x, alpha)
696 !TODO this is still neeed?
697 wout(ixo^s,1:nw) = w(ixo^s,1:nw)
698
699 do idir = 1, ndir
700 ! d1 from Eq 6
701 tmp2(ixo^s) = 0d0
702 do n = 1, dust_n_species
703 tmp2(ixo^s) = tmp2(ixo^s) + alpha(ixo^s, idir,n) * &
704 (w(ixo^s,gas_rho_) + w(ixo^s,dust_rho(n)))
705
706 enddo
707 ! store D in tmp
708 tmp(ixo^s) = 1d0 + tmp2(ixo^s) * qdt
709 if(dust_implicit_second_order) then
710 ! d2 from Eq 6
711 tmp2(ixo^s) = 0d0
712 do n = 1, dust_n_species
713 do m = n+1, dust_n_species
714 tmp2(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) *&
715 (w(ixo^s,gas_rho_) + w(ixo^s,dust_rho(n))+w(ixo^s,dust_rho(m)))
716 enddo
717 enddo
718 ! multiplied at the end by rho_gas
719 tmp(ixo^s) = tmp(ixo^s) + w(ixo^s,gas_rho_)*tmp2(ixo^s) * (qdt**2)
720 endif
721
722
723
724 do n = 1, dust_n_species
725 ! ni1 from eq 7
726 tmp2(ixo^s) = alpha(ixo^s, idir,n) * ( &
727 w(ixo^s,dust_rho(n)) * w(ixo^s, gas_mom(idir)) - &
728 w(ixo^s,gas_rho_) * w(ixo^s, dust_mom(idir, n))) * qdt
729
730 if(dust_implicit_second_order) then
731 ! ni2 from eq 7
732 tmp3(ixo^s) = 0d0
733 do m = n+1, dust_n_species
734 tmp3(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) * &
735 ( w(ixo^s,dust_rho(n)) * (w(ixo^s, dust_mom(idir, n)) + w(ixo^s, gas_mom(idir))) - &
736 (w(ixo^s,gas_rho_) + w(ixo^s,dust_rho(m))) * w(ixo^s, dust_mom(idir, n)) )
737 enddo
738 ! tmp3 multiplied at the end by rho_gas
739 tmp2(ixo^s) = tmp2(ixo^s) + tmp3(ixo^s) * w(ixo^s,gas_rho_)* (qdt**2)
740 endif
741 tmp2(ixo^s) = tmp2(ixo^s)/tmp(ixo^s)
742 wout(ixo^s, dust_mom(idir,n)) = w(ixo^s, dust_mom(idir,n)) + tmp2(ixo^s)
743 enddo
744
745 if (dust_backreaction) then
746 tmp2(ixo^s) = 0d0
747 !n1 from eq 8
748 do n = 1, dust_n_species
749 tmp2(ixo^s) = tmp2(ixo^s) + alpha(ixo^s, idir,n) * &
750 (w(ixo^s,gas_rho_) * w(ixo^s, dust_mom(idir,n)) - &
751 w(ixo^s,dust_rho(n)) * w(ixo^s, gas_mom(idir)))
752
753 enddo
754 tmp2(ixo^s) = qdt * tmp2(ixo^s)
755 if(dust_implicit_second_order) then
756 !n2 from eq 8
757 tmp3(ixo^s) = 0d0
758 do n = 1, dust_n_species
759 do m = n+1, dust_n_species
760 tmp3(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) * &
761 (w(ixo^s,gas_rho_) * (w(ixo^s, dust_mom(idir, n)) + w(ixo^s, dust_mom(idir, m))) - &
762 (w(ixo^s,dust_rho(n)) + w(ixo^s,dust_rho(m)))* w(ixo^s, gas_mom(idir)))
763 enddo
764 enddo
765 ! tmp3 multiplied at the end by rho_gas
766 tmp2(ixo^s) = tmp2(ixo^s) + (qdt**2)*tmp3(ixo^s)* w(ixo^s,gas_rho_)
767 endif
768 ! store in tmp2 contribution to momentum
769 ! so that it is used when dust_backreaction_fh = .false.
770 tmp2(ixo^s) = tmp2(ixo^s) / tmp(ixo^s)
771 wout(ixo^s, gas_mom(idir)) = w(ixo^s, gas_mom(idir)) + tmp2(ixo^s)
772
773 ! kinetic energy update
774 if (gas_e_ > 0) then
775 if(dust_backreaction_fh) then
776 ! add work done by coll terms + FrictionalHeating
777 tmp2(ixo^s) = 0d0
778 do n = 1, dust_n_species
779 ! 2*dust kinetic energy: dust rho can be 0
780 where(w(ixo^s,dust_rho(n)) > 0d0)
781 tmp3(ixo^s)= w(ixo^s, dust_mom(idir,n))**2/w(ixo^s,dust_rho(n))
782 elsewhere
783 tmp3(ixo^s) = 0d0
784 endwhere
785 tmp2(ixo^s) = tmp2(ixo^s) + alpha(ixo^s, idir,n) * &
786 (w(ixo^s,gas_rho_) * tmp3(ixo^s) - &
787 w(ixo^s,dust_rho(n)) * (w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_)))
788
789 enddo
790 tmp2(ixo^s) = qdt * tmp2(ixo^s)
791 if(dust_implicit_second_order) then
792 tmp3(ixo^s) = 0d0
793 do n = 1, dust_n_species
794 do m = n+1, dust_n_species
795 tmp3(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) * &
796 (w(ixo^s,gas_rho_) * (w(ixo^s, dust_mom(idir, n))**2/w(ixo^s,dust_rho(n)) + w(ixo^s, dust_mom(idir,m))**2/w(ixo^s,dust_rho(m))) - &
797 (w(ixo^s,dust_rho(n)) + w(ixo^s,dust_rho(m)))* w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_))
798 enddo
799 enddo
800 ! tmp3 multiplied at the end by rho_gas
801 tmp2(ixo^s) = tmp2(ixo^s) + (qdt**2)*tmp3(ixo^s)* w(ixo^s,gas_rho_)
802 endif
803 wout(ixo^s, gas_e_) = wout(ixo^s, gas_e_) + 0.5d0 * tmp2(ixo^s) / tmp(ixo^s)
804 else
805 ! dust_backreaction_fh = .false.
806 ! add only work done by coll term by multiplyting the contribution in mom eq. by velocity
807 wout(ixo^s, gas_e_) = wout(ixo^s, gas_e_) + vgas(ixo^s, idir) * tmp2(ixo^s)
808 endif
809 end if
810 end if
811 end do !1..ndir
812
813
814 end subroutine dust_advance_implicit_grid
815
816 ! copied from get_3d_dragforce subroutine
817 subroutine get_alpha_dust(ixI^L, ixO^L, w, vgas,x, alpha)
820 integer, intent(in) :: ixi^l, ixo^l
821 double precision, intent(in) :: x(ixi^s, 1:ndim)
822 double precision, intent(in) :: w(ixi^s, 1:nw)
823 double precision,intent(in) :: vgas(ixi^s, 1:ndir)
824 double precision, intent(out) :: &
825 alpha(ixi^s, 1:ndir, 1:dust_n_species)
826
827 double precision :: ptherm(ixi^s)
828 double precision, dimension(ixI^S) :: vt2, deltav, fd, vdust
829 double precision :: alpha_t(ixi^s, 1:dust_n_species)
830 integer :: n, idir
831
832 call phys_get_pthermal(w, x, ixi^l, ixo^l, ptherm)
833
834 vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/w(ixo^s, gas_rho_)
835
836 select case( trim(dust_method) )
837 case ('Kwok') ! assume sticking coefficient equals 0.25
838
839 do idir = 1, ndir
840 do n = 1, dust_n_species
841 where(w(ixo^s, dust_rho(n)) > 0.0d0)
842
843 ! 0.75 from sticking coefficient
844 fd(ixo^s) = 0.75d0 / (dust_density(n) * dust_size(n))
845
846 ! 0.75 from spherical grainvolume
847 vdust(ixo^s) = w(ixo^s, dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
848 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
849 fd(ixo^s) = fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
850 elsewhere
851 fd(ixo^s) = 0.0d0
852 end where
853 alpha(ixo^s, idir, n) = fd(ixo^s)
854 end do
855 end do
856
857 case ('sticking') ! Calculate sticking coefficient based on the gas and dust temperatures
858
859 call get_sticking(w, x, ixi^l, ixo^l, alpha_t, ptherm)
860
861 do idir = 1, ndir
862 do n = 1, dust_n_species
863 where(w(ixo^s, dust_rho(n))>0.0d0)
864 ! sticking
865 fd(ixo^s) = (one-alpha_t(ixo^s,n)) / (dust_density(n)*dust_size(n))
866 ! 0.75 from spherical grainvolume
867 vdust(ixo^s) = w(ixo^s,dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
868 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
869 fd(ixo^s) = fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
870 else where
871 fd(ixo^s) = 0.0d0
872 end where
873 alpha(ixo^s, idir,n) = fd(ixo^s)
874 end do
875 end do
876
877 case('linear') !linear with Deltav, for testing (see Laibe & Price 2011)
878 do idir = 1, ndir
879 do n = 1, dust_n_species
880 where(w(ixo^s, dust_rho(n))>0.0d0)
881 fd(ixo^s) = dust_k_lineardrag/(w(ixo^s,gas_rho_)*w(ixo^s, dust_rho(n)))
882 else where
883 fd(ixo^s) = 0.0d0
884 end where
885 alpha(ixo^s, idir,n) = fd(ixo^s)
886 end do
887 end do
888
889 case('none')
890 alpha(ixo^s, :, :) = 0.0d0
891 case default
892 call mpistop( "=== This dust method has not been implemented===" )
893 end select
894
895 end subroutine get_alpha_dust
896
897
898 !> Get dt related to dust and gas stopping time (Laibe 2011)
899 subroutine dust_get_dt(wprim, ixI^L, ixO^L, dtnew, dx^D, x)
902
903 integer, intent(in) :: ixi^l, ixo^l
904 double precision, intent(in) :: dx^d, x(ixi^s, 1:ndim)
905 double precision, intent(in) :: wprim(ixi^s, 1:nw)
906 double precision, intent(inout) :: dtnew
907
908 double precision :: ptherm(ixi^s), vgas(ixi^s, 1:ndir)
909 double precision, dimension(1:dust_n_species):: dtdust
910 double precision, dimension(ixI^S) :: vt2, deltav, tstop, vdust
911 double precision, dimension(ixI^S, 1:dust_n_species) :: alpha_t
912 integer :: n, idir
913
914 if(dust_dtpar .le. 0) return
915
916 ! hydro has no splitting of pressure and primitives on entry
917 if(gas_e_>0) then
918 ptherm(ixo^s)=wprim(ixo^s,gas_e_)
919 else
920 call mpistop("adjust dust module for no energy for gas")
921 endif
922 do idir = 1, ndir
923 vgas(ixo^s,idir)=wprim(ixo^s,gas_mom(idir))
924 end do
925
926 select case( trim(dust_method) )
927
928 case( 'Kwok' ) ! assume sticking coefficient equals 0.25
929 dtdust(:) = bigdouble
930
931 vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/wprim(ixo^s, gas_rho_)
932
933 do idir = 1, ndir
934 do n = 1, dust_n_species
935 where(wprim(ixo^s, dust_rho(n))>0.0d0)
936 vdust(ixo^s) = wprim(ixo^s,dust_mom(idir, n))
937 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
938 tstop(ixo^s) = 4.0d0*(dust_density(n)*dust_size(n))/ &
939 (3.0d0*(0.75d0)*dsqrt(vt2(ixo^s) + &
940 deltav(ixo^s)**2)*(wprim(ixo^s, dust_rho(n)) + &
941 wprim(ixo^s, gas_rho_)))
942 else where
943 tstop(ixo^s) = bigdouble
944 end where
945
946 dtdust(n) = min(minval(tstop(ixo^s)), dtdust(n))
947 end do
948 end do
949
950 dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
951
952 case( 'sticking' ) ! Calculate sticking coefficient based on the gas temperature
953 dtdust(:) = bigdouble
954
955 vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/wprim(ixo^s, gas_rho_)
956
957 ! Sticking coefficient
958 call get_sticking(wprim, x, ixi^l, ixo^l, alpha_t, ptherm)
959
960 do idir = 1, ndir
961 do n = 1, dust_n_species
962 where(wprim(ixo^s, dust_rho(n))>0.0d0)
963 vdust(ixo^s) = wprim(ixo^s,dust_mom(idir, n))
964 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
965 tstop(ixo^s) = 4.0d0*(dust_density(n)*dust_size(n))/ &
966 (3.0d0*(one-alpha_t(ixo^s,n))*dsqrt(vt2(ixo^s) + &
967 deltav(ixo^s)**2)*(wprim(ixo^s, dust_rho(n)) + &
968 wprim(ixo^s, gas_rho_)))
969 else where
970 tstop(ixo^s) = bigdouble
971 end where
972
973 dtdust(n) = min(minval(tstop(ixo^s)), dtdust(n))
974 end do
975 end do
976
977 dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
978
979 case('linear') !linear with Deltav, for testing (see Laibe & Price 2011)
980 dtdust(:) = bigdouble
981
982 do n = 1, dust_n_species
983 where(wprim(ixo^s, dust_rho(n))>0.0d0)
984 tstop(ixo^s) = (wprim(ixo^s, dust_rho(n))*wprim(ixo^s, gas_rho_))/ &
985 (dust_k_lineardrag*(wprim(ixo^s, dust_rho(n)) + wprim(ixo^s, gas_rho_)))
986 else where
987 tstop(ixo^s) = bigdouble
988 end where
989
990 dtdust(n) = min(minval(tstop(ixo^s)), dtdust(n))
991 end do
992
993 dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
994 case('usr')
995 dtdust(:) = bigdouble
996 call usr_dust_get_dt(wprim, ixi^l, ixo^l, dtdust, dx^d, x, dust_n_species)
997 dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
998 case('none')
999 ! no dust timestep
1000 case default
1001 call mpistop( "=== This dust method has not been implemented===" )
1002 end select
1003
1004 if (dtnew < dtmin) then
1005 write(unitterm,*)"-------------------------------------"
1006 write(unitterm,*)"Warning: found DUST related time step too small! dtnew=", dtnew
1007 write(unitterm,*)"on grid with index:", block%igrid," grid level=", block%level
1008 write(unitterm,*)"grid corners are=",{^d&rnode(rpxmin^d_, block%igrid), rnode(rpxmax^d_, block%igrid)}
1009 write(unitterm,*)" dtdust =", dtdust(:)
1010 write(unitterm,*)"on processor:", mype
1011 write(unitterm,*)"-------------------------------------"
1012 endif
1013
1014 end subroutine dust_get_dt
1015
1016 ! Note that cmax and cmin are assumed to be initialized
1017 subroutine dust_get_cmax(w, x, ixI^L, ixO^L, idim, cmax, cmin)
1019 use mod_variables
1020
1021 integer, intent(in) :: ixi^l, ixo^l, idim
1022 double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
1023 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
1024 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
1025 double precision :: vdust(ixo^s)
1026 integer :: n
1027
1028 do n = 1, dust_n_species
1029 vdust(ixo^s) = get_vdust(w, ixi^l, ixo^l, idim, n)
1030
1031 if (present(cmin)) then
1032 cmin(ixo^s,1) = min(cmin(ixo^s,1), vdust(ixo^s))
1033 cmax(ixo^s,1) = max(cmax(ixo^s,1), vdust(ixo^s))
1034 else
1035 cmax(ixo^s,1) = max(cmax(ixo^s,1), abs(vdust(ixo^s)))
1036 end if
1037 end do
1038 end subroutine dust_get_cmax
1039
1040 ! Note that cmax and cmin are assumed to be initialized
1041 subroutine dust_get_cmax_prim(w, x, ixI^L, ixO^L, idim, cmax, cmin)
1043
1044 integer, intent(in) :: ixi^l, ixo^l, idim
1045 double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
1046 double precision, intent(inout) :: cmax(ixi^s)
1047 double precision, intent(inout), optional :: cmin(ixi^s)
1048 double precision :: vdust(ixo^s)
1049 integer :: n
1050
1051 do n = 1, dust_n_species
1052 vdust(ixo^s) = get_vdust_prim(w, ixi^l, ixo^l, idim, n)
1053
1054 if (present(cmin)) then
1055 cmin(ixo^s) = min(cmin(ixo^s), vdust(ixo^s))
1056 cmax(ixo^s) = max(cmax(ixo^s), vdust(ixo^s))
1057 else
1058 cmax(ixo^s) = max(cmax(ixo^s), abs(vdust(ixo^s)))
1059 end if
1060 end do
1061 end subroutine dust_get_cmax_prim
1062
1063end module mod_dust
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for including dust species, which interact with the gas through a drag force.
Definition mod_dust.t:3
double precision, public, protected dust_min_rho
Minimum dust density as used when dust_small_to_zero=T.
Definition mod_dust.t:28
subroutine, public dust_add_source(qdt, ixil, ixol, wct, w, x, qsourcesplit, active)
w[iw]= w[iw]+qdt*S[wCT, x] where S is the source based on wCT within ixO
Definition mod_dust.t:532
double precision, dimension(:), allocatable, public dust_size
Size of each dust species, dimensionless expression.
Definition mod_dust.t:31
subroutine, public dust_evaluate_implicit(qtc, psa)
inplace update of psa==>F_im(psa)
Definition mod_dust.t:585
character(len=std_len), public, protected dust_method
What type of dust drag force to use. Can be 'Kwok', 'sticking', 'linear', 'usr' or 'none'.
Definition mod_dust.t:66
subroutine, public dust_to_primitive(ixil, ixol, w, x)
Definition mod_dust.t:228
subroutine, public dust_get_flux(w, x, ixil, ixol, idim, f)
Definition mod_dust.t:252
integer, dimension(:, :), allocatable, public, protected dust_mom
Indices of the dust momentum densities.
Definition mod_dust.t:47
subroutine, public set_dusttozero(ixil, ixol, w, x)
Definition mod_dust.t:327
logical, public, protected dust_small_to_zero
Set small dust densities to zero to avoid numerical problems.
Definition mod_dust.t:50
subroutine, public dust_to_conserved(ixil, ixol, w, x)
Definition mod_dust.t:208
integer, public, protected dust_n_species
The number of dust species.
Definition mod_dust.t:37
subroutine, public dust_get_flux_prim(w, x, ixil, ixol, idim, f)
Definition mod_dust.t:275
subroutine, public dust_check_w(ixil, ixol, w, flag)
Definition mod_dust.t:194
integer, dimension(:), allocatable, public, protected dust_rho
Indices of the dust densities.
Definition mod_dust.t:44
subroutine, public dust_get_cmax(w, x, ixil, ixol, idim, cmax, cmin)
Definition mod_dust.t:1018
subroutine, public dust_check_params()
Definition mod_dust.t:154
subroutine, public dust_get_cmax_prim(w, x, ixil, ixol, idim, cmax, cmin)
Definition mod_dust.t:1042
subroutine, public dust_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
Get dt related to dust and gas stopping time (Laibe 2011)
Definition mod_dust.t:900
subroutine, public dust_init(g_rho, g_mom, g_energy)
Definition mod_dust.t:95
subroutine, public dust_implicit_update(dtfactor, qdt, qtc, psb, psa)
Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
Definition mod_dust.t:654
double precision, dimension(:), allocatable, public dust_density
Internal density of each dust species, dimensionless expression.
Definition mod_dust.t:34
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
integer coordinate
Definition mod_geometry.t:7
integer, parameter spherical
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
integer, parameter unitpar
file handle for IO
logical any_source_split
if any normal source term is added in split fasion
logical use_imex_scheme
whether IMEX in use or not
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
integer ixm
the mesh range of a physical block without ghost cells
integer, parameter unitterm
Unit for standard output.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision unit_temperature
Physical scaling factor for temperature.
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
logical fix_small_values
fix small values with average or replace methods
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
double precision dtmin
Stop the simulation when the time step becomes smaller than this value.
integer max_blocks
The maximum number of grid blocks in a processor.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
Module with all the methods that users can customize in AMRVAC.
procedure(phys_dust_get_dt), pointer usr_dust_get_dt
procedure(phys_dust_get_3d_dragforce), pointer usr_get_3d_dragforce
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...