MPI-AMRVAC 3.1
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 subroutine get_sticking(w, x, ixI^L, ixO^L, alpha_T, ptherm)
440 integer, intent(in) :: ixi^l, ixo^l
441 double precision, intent(in) :: x(ixi^s, 1:ndim)
442 double precision, intent(in) :: w(ixi^s, 1:nw)
443 double precision, intent(out) :: alpha_t(ixi^s, 1:dust_n_species)
444 double precision, intent(in) :: ptherm(ixi^s)
445 double precision :: tgas(ixi^s)
446 integer :: n
447
448 ! get the dust species T in K
449 call get_tdust(w, x, ixi^l, ixo^l, alpha_t)
450
451 ! convert dimensionless gas T to K
452 tgas(ixo^s) = (ptherm(ixo^s)/w(ixo^s, gas_rho_))*unit_temperature
453
454 do n = 1, dust_n_species
455 alpha_t(ixo^s,n) = 0.35d0 * dexp(-dsqrt((tgas(ixo^s) + &
456 alpha_t(ixo^s,n))/5.0d2))+0.1d0
457 end do
458
459 end subroutine get_sticking
460
461 !> Returns dust temperature (in K), either as constant or based on equ. 5.41,
462 !> 5.42 and 5.44 from Tielens (2005)
463 !>
464 !> Note that this calculation assumes cgs!!!!
465 !>
466 !> It takes as input the stellar luminosity in solar units in 'stellar' case
467 !> or a fixed dust input temperature in Kelvin when 'constant' or does case 'ism'
468 subroutine get_tdust(w, x, ixI^L, ixO^L, Td)
470 use mod_geometry
471
472 integer, intent(in) :: ixi^l, ixo^l
473 double precision, intent(in) :: x(ixi^s, 1:ndim)
474 double precision, intent(in) :: w(ixi^s, 1:nw)
475 double precision, intent(out) :: td(ixi^s, 1:dust_n_species)
476 double precision :: g0(ixo^s)
477 integer :: n
478
479 select case( trim(dust_temperature_type) )
480 case( 'constant' )
481 td(ixo^s, :) = dust_temperature
482 case( 'ism' )
483 select case( trim(dust_species) )
484 case( 'graphite' )
485 do n = 1, dust_n_species
486 td(ixo^s, n) = 15.8d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0)
487 end do
488 case( 'silicate' )
489 do n = 1, dust_n_species
490 td(ixo^s, n) = 13.6d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0)
491 end do
492 case default
493 call mpistop( "=== Dust species undetermined===" )
494 end select
495 case( 'stellar' )
496 select case(coordinate)
497 case(spherical)
498 g0(ixo^s) = max(x(ixo^s, 1)*unit_length, smalldouble)
499 !!!case(cylindrical) convert R,Z to spherical radial coordinate r here
500 !!! but only ok for 2D (R,Z) or 2.5D (R,Z) case
501 !!! G0(ixO^S) = max(dsqrt(sum(x(ixO^S,:)**2,dim=ndim+1))*unit_length, smalldouble)
502 case default
503 call mpistop('stellar case not available in this coordinate system')
504 end select
505
506 g0(ixo^s) = 2.1d4*(dust_stellar_luminosity/1.0d8)*((3.0857d17/g0(ixo^s))**2)
507
508 select case( trim(dust_species) )
509 case( 'graphite' )
510 do n = 1, dust_n_species
511 td(ixo^s, n) = 61.0d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0) &
512 *(g0(ixo^s)**(one/5.8d0))
513 end do
514 case( 'silicate' )
515 do n = 1, dust_n_species
516 td(ixo^s, n) = 50.0d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0) &
517 *(g0(ixo^s)**(one/6.0d0))
518 end do
519 case default
520 call mpistop( "=== Dust species undetermined===" )
521 end select
522 case default
523 call mpistop( "=== Dust temperature undetermined===" )
524 end select
525
526 end subroutine get_tdust
527
528 !> w[iw]= w[iw]+qdt*S[wCT, x] where S is the source based on wCT within ixO
529 subroutine dust_add_source(qdt, ixI^L, ixO^L, wCT, w, x, qsourcesplit, active)
531
532 integer, intent(in) :: ixi^l, ixo^l
533 double precision, intent(in) :: qdt
534 double precision, intent(in) :: wct(ixi^s, 1:nw), x(ixi^s, 1:ndim)
535 double precision, intent(inout) :: w(ixi^s, 1:nw)
536 logical, intent(in) :: qsourcesplit
537 logical, intent(inout) :: active
538
539 double precision :: ptherm(ixi^s), vgas(ixi^s, 1:ndir)
540 double precision :: fdrag(ixi^s, 1:ndir, 1:dust_n_species)
541 integer :: n, idir
542
543 select case( trim(dust_method) )
544 case( 'none' )
545 !do nothing here
546 case default !all regular dust methods here
547 if (qsourcesplit .eqv. dust_source_split) then
548 active = .true.
549
550 call phys_get_pthermal(wct, x, ixi^l, ixo^l, ptherm)
551 do idir=1,ndir
552 vgas(ixo^s,idir)=wct(ixo^s,gas_mom(idir))/wct(ixo^s,gas_rho_)
553 end do
554
555 call get_3d_dragforce(ixi^l, ixo^l, wct, x, fdrag, ptherm, vgas)
556 fdrag(ixo^s, 1:ndir, 1:dust_n_species) = fdrag(ixo^s, 1:ndir, 1:dust_n_species) * qdt
557
558 do idir = 1, ndir
559
560 do n = 1, dust_n_species
561 if (dust_backreaction) then
562 w(ixo^s, gas_mom(idir)) = w(ixo^s, gas_mom(idir)) + &
563 fdrag(ixo^s, idir, n)
564 if (gas_e_ > 0) then
565 w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + vgas(ixo^s, idir) &
566 * fdrag(ixo^s, idir, n)
567 end if
568 end if
569
570
571 w(ixo^s, dust_mom(idir, n)) = w(ixo^s, dust_mom(idir, n)) - &
572 fdrag(ixo^s, idir, n)
573 end do
574 end do
575
576 endif
577 end select
578
579 end subroutine dust_add_source
580
581 !> inplace update of psa==>F_im(psa)
582 subroutine dust_evaluate_implicit(qtC,psa)
584 type(state), target :: psa(max_blocks)
585 double precision, intent(in) :: qtc
586
587 integer :: iigrid, igrid, level
588
589 !dust_method = 'none' not used
590
591 !$OMP PARALLEL DO PRIVATE(igrid)
592 do iigrid=1,igridstail; igrid=igrids(iigrid);
593 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
594 block=>psa(igrid)
595 call dust_terms(ixg^ll,ixm^ll,psa(igrid)%w,psa(igrid)%x)
596 end do
597 !$OMP END PARALLEL DO
598
599 end subroutine dust_evaluate_implicit
600
601
602
603
604 subroutine dust_terms(ixI^L,ixO^L,w,x)
606 integer, intent(in) :: ixi^l, ixo^l
607 double precision, intent(inout) :: w(ixi^s, 1:nw)
608 double precision, intent(in) :: x(ixi^s,1:ndim)
609
610 double precision :: tmp(ixi^s), vgas(ixi^s, 1:ndir)
611 double precision :: alpha(ixi^s, 1:ndir, 1:dust_n_species)
612 integer :: n, idir
613
614 do idir=1,ndir
615 vgas(ixo^s,idir)=w(ixo^s,gas_mom(idir))/w(ixo^s,gas_rho_)
616 end do
617 call get_alpha_dust(ixi^l, ixo^l, w, vgas,x, alpha)
618 w(ixo^s, gas_e_)=0d0
619 do idir = 1, ndir
620
621 w(ixo^s, gas_mom(idir))=0d0
622 do n = 1, dust_n_species
623 ! contribution for gas momentum
624 tmp(ixo^s) = alpha(ixo^s, idir,n) * ( &
625 w(ixo^s,dust_rho(n)) * w(ixo^s, gas_mom(idir)) - &
626 w(ixo^s,gas_rho_) * w(ixo^s, dust_mom(idir, n)))
627 w(ixo^s, dust_mom(idir, n)) = -tmp(ixo^s)
628 if (dust_backreaction) then
629 w(ixo^s, gas_mom(idir)) = w(ixo^s, gas_mom(idir)) + tmp(ixo^s)
630 if (gas_e_ > 0) then
631 if(dust_backreaction_fh) then
632 where(w(ixo^s,dust_rho(n)) > 0d0)
633 w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + alpha(ixo^s, idir,n) * &
634 (w(ixo^s,gas_rho_) * (w(ixo^s, dust_mom(idir,n))**2/w(ixo^s,dust_rho(n))) - &
635 w(ixo^s,dust_rho(n)) * (w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_)))
636 elsewhere
637 w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + alpha(ixo^s, idir,n) * ( - &
638 w(ixo^s,dust_rho(n)) * (w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_)))
639 endwhere
640 else
641 w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + vgas(ixo^s, idir) &
642 * tmp(ixo^s)
643 end if
644 end if
645 end if
646 end do
647 end do
648 end subroutine dust_terms
649
650 !> Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
651 subroutine dust_implicit_update(dtfactor,qdt,qtC,psb,psa)
653 !use mod_ghostcells_update
654
655 type(state), target :: psa(max_blocks)
656 type(state), target :: psb(max_blocks)
657 double precision, intent(in) :: qdt
658 double precision, intent(in) :: qtc
659 double precision, intent(in) :: dtfactor
660
661 integer :: iigrid, igrid
662
663 !call getbc(global_time,0.d0,psa,1,nw)
664 !$OMP PARALLEL DO PRIVATE(igrid)
665 do iigrid=1,igridstail; igrid=igrids(iigrid);
666 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
667 block=>psa(igrid)
668 call dust_advance_implicit_grid(ixg^ll, ixg^ll, psa(igrid)%w, psb(igrid)%w, psa(igrid)%x, dtfactor,qdt)
669 end do
670 !$OMP END PARALLEL DO
671
672 end subroutine dust_implicit_update
673
674 subroutine dust_advance_implicit_grid(ixI^L, ixO^L, w, wout, x, dtfactor,qdt)
676 integer, intent(in) :: ixi^l, ixo^l
677 double precision, intent(in) :: qdt
678 double precision, intent(in) :: dtfactor
679 double precision, intent(in) :: w(ixi^s,1:nw)
680 double precision, intent(in) :: x(ixi^s,1:ndim)
681 double precision, intent(out) :: wout(ixi^s,1:nw)
682
683 integer :: n, m, idir
684 double precision :: alpha(ixi^s, 1:ndir, 1:dust_n_species)
685 double precision :: tmp(ixi^s),tmp2(ixi^s)
686 double precision :: tmp3(ixi^s)
687 double precision :: vgas(ixi^s, 1:ndir)
688
689
690 do idir = 1, ndir
691 vgas(ixo^s,idir)=w(ixo^s,gas_mom(idir))/w(ixo^s,gas_rho_)
692 end do
693 call get_alpha_dust(ixi^l, ixo^l, w, vgas, x, alpha)
694 !TODO this is still neeed?
695 wout(ixo^s,1:nw) = w(ixo^s,1:nw)
696
697 do idir = 1, ndir
698 ! d1 from Eq 6
699 tmp2(ixo^s) = 0d0
700 do n = 1, dust_n_species
701 tmp2(ixo^s) = tmp2(ixo^s) + alpha(ixo^s, idir,n) * &
702 (w(ixo^s,gas_rho_) + w(ixo^s,dust_rho(n)))
703
704 enddo
705 ! store D in tmp
706 tmp(ixo^s) = 1d0 + tmp2(ixo^s) * qdt
707 if(dust_implicit_second_order) then
708 ! d2 from Eq 6
709 tmp2(ixo^s) = 0d0
710 do n = 1, dust_n_species
711 do m = n+1, dust_n_species
712 tmp2(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) *&
713 (w(ixo^s,gas_rho_) + w(ixo^s,dust_rho(n))+w(ixo^s,dust_rho(m)))
714 enddo
715 enddo
716 ! multiplied at the end by rho_gas
717 tmp(ixo^s) = tmp(ixo^s) + w(ixo^s,gas_rho_)*tmp2(ixo^s) * (qdt**2)
718 endif
719
720
721
722 do n = 1, dust_n_species
723 ! ni1 from eq 7
724 tmp2(ixo^s) = alpha(ixo^s, idir,n) * ( &
725 w(ixo^s,dust_rho(n)) * w(ixo^s, gas_mom(idir)) - &
726 w(ixo^s,gas_rho_) * w(ixo^s, dust_mom(idir, n))) * qdt
727
728 if(dust_implicit_second_order) then
729 ! ni2 from eq 7
730 tmp3(ixo^s) = 0d0
731 do m = n+1, dust_n_species
732 tmp3(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) * &
733 ( w(ixo^s,dust_rho(n)) * (w(ixo^s, dust_mom(idir, n)) + w(ixo^s, gas_mom(idir))) - &
734 (w(ixo^s,gas_rho_) + w(ixo^s,dust_rho(m))) * w(ixo^s, dust_mom(idir, n)) )
735 enddo
736 ! tmp3 multiplied at the end by rho_gas
737 tmp2(ixo^s) = tmp2(ixo^s) + tmp3(ixo^s) * w(ixo^s,gas_rho_)* (qdt**2)
738 endif
739 tmp2(ixo^s) = tmp2(ixo^s)/tmp(ixo^s)
740 wout(ixo^s, dust_mom(idir,n)) = w(ixo^s, dust_mom(idir,n)) + tmp2(ixo^s)
741 enddo
742
743 if (dust_backreaction) then
744 tmp2(ixo^s) = 0d0
745 !n1 from eq 8
746 do n = 1, dust_n_species
747 tmp2(ixo^s) = tmp2(ixo^s) + alpha(ixo^s, idir,n) * &
748 (w(ixo^s,gas_rho_) * w(ixo^s, dust_mom(idir,n)) - &
749 w(ixo^s,dust_rho(n)) * w(ixo^s, gas_mom(idir)))
750
751 enddo
752 tmp2(ixo^s) = qdt * tmp2(ixo^s)
753 if(dust_implicit_second_order) then
754 !n2 from eq 8
755 tmp3(ixo^s) = 0d0
756 do n = 1, dust_n_species
757 do m = n+1, dust_n_species
758 tmp3(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) * &
759 (w(ixo^s,gas_rho_) * (w(ixo^s, dust_mom(idir, n)) + w(ixo^s, dust_mom(idir, m))) - &
760 (w(ixo^s,dust_rho(n)) + w(ixo^s,dust_rho(m)))* w(ixo^s, gas_mom(idir)))
761 enddo
762 enddo
763 ! tmp3 multiplied at the end by rho_gas
764 tmp2(ixo^s) = tmp2(ixo^s) + (qdt**2)*tmp3(ixo^s)* w(ixo^s,gas_rho_)
765 endif
766 ! store in tmp2 contribution to momentum
767 ! so that it is used when dust_backreaction_fh = .false.
768 tmp2(ixo^s) = tmp2(ixo^s) / tmp(ixo^s)
769 wout(ixo^s, gas_mom(idir)) = w(ixo^s, gas_mom(idir)) + tmp2(ixo^s)
770
771 ! kinetic energy update
772 if (gas_e_ > 0) then
773 if(dust_backreaction_fh) then
774 ! add work done by coll terms + FrictionalHeating
775 tmp2(ixo^s) = 0d0
776 do n = 1, dust_n_species
777 ! 2*dust kinetic energy: dust rho can be 0
778 where(w(ixo^s,dust_rho(n)) > 0d0)
779 tmp3(ixo^s)= w(ixo^s, dust_mom(idir,n))**2/w(ixo^s,dust_rho(n))
780 elsewhere
781 tmp3(ixo^s) = 0d0
782 endwhere
783 tmp2(ixo^s) = tmp2(ixo^s) + alpha(ixo^s, idir,n) * &
784 (w(ixo^s,gas_rho_) * tmp3(ixo^s) - &
785 w(ixo^s,dust_rho(n)) * (w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_)))
786
787 enddo
788 tmp2(ixo^s) = qdt * tmp2(ixo^s)
789 if(dust_implicit_second_order) then
790 tmp3(ixo^s) = 0d0
791 do n = 1, dust_n_species
792 do m = n+1, dust_n_species
793 tmp3(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) * &
794 (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))) - &
795 (w(ixo^s,dust_rho(n)) + w(ixo^s,dust_rho(m)))* w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_))
796 enddo
797 enddo
798 ! tmp3 multiplied at the end by rho_gas
799 tmp2(ixo^s) = tmp2(ixo^s) + (qdt**2)*tmp3(ixo^s)* w(ixo^s,gas_rho_)
800 endif
801 wout(ixo^s, gas_e_) = wout(ixo^s, gas_e_) + 0.5d0 * tmp2(ixo^s) / tmp(ixo^s)
802 else
803 ! dust_backreaction_fh = .false.
804 ! add only work done by coll term by multiplyting the contribution in mom eq. by velocity
805 wout(ixo^s, gas_e_) = wout(ixo^s, gas_e_) + vgas(ixo^s, idir) * tmp2(ixo^s)
806 endif
807 end if
808 end if
809 end do !1..ndir
810
811
812 end subroutine dust_advance_implicit_grid
813
814 ! copied from get_3d_dragforce subroutine
815 subroutine get_alpha_dust(ixI^L, ixO^L, w, vgas,x, alpha)
818 integer, intent(in) :: ixi^l, ixo^l
819 double precision, intent(in) :: x(ixi^s, 1:ndim)
820 double precision, intent(in) :: w(ixi^s, 1:nw)
821 double precision,intent(in) :: vgas(ixi^s, 1:ndir)
822 double precision, intent(out) :: &
823 alpha(ixi^s, 1:ndir, 1:dust_n_species)
824
825 double precision :: ptherm(ixi^s)
826 double precision, dimension(ixI^S) :: vt2, deltav, fd, vdust
827 double precision :: alpha_t(ixi^s, 1:dust_n_species)
828 integer :: n, idir
829
830 call phys_get_pthermal(w, x, ixi^l, ixo^l, ptherm)
831
832 vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/w(ixo^s, gas_rho_)
833
834 select case( trim(dust_method) )
835 case ('Kwok') ! assume sticking coefficient equals 0.25
836
837 do idir = 1, ndir
838 do n = 1, dust_n_species
839 where(w(ixo^s, dust_rho(n)) > 0.0d0)
840
841 ! 0.75 from sticking coefficient
842 fd(ixo^s) = 0.75d0 / (dust_density(n) * dust_size(n))
843
844 ! 0.75 from spherical grainvolume
845 vdust(ixo^s) = w(ixo^s, dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
846 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
847 fd(ixo^s) = fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
848 elsewhere
849 fd(ixo^s) = 0.0d0
850 end where
851 alpha(ixo^s, idir, n) = fd(ixo^s)
852 end do
853 end do
854
855 case ('sticking') ! Calculate sticking coefficient based on the gas and dust temperatures
856
857 call get_sticking(w, x, ixi^l, ixo^l, alpha_t, ptherm)
858
859 do idir = 1, ndir
860 do n = 1, dust_n_species
861 where(w(ixo^s, dust_rho(n))>0.0d0)
862 ! sticking
863 fd(ixo^s) = (one-alpha_t(ixo^s,n)) / (dust_density(n)*dust_size(n))
864 ! 0.75 from spherical grainvolume
865 vdust(ixo^s) = w(ixo^s,dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
866 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
867 fd(ixo^s) = fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
868 else where
869 fd(ixo^s) = 0.0d0
870 end where
871 alpha(ixo^s, idir,n) = fd(ixo^s)
872 end do
873 end do
874
875 case('linear') !linear with Deltav, for testing (see Laibe & Price 2011)
876 do idir = 1, ndir
877 do n = 1, dust_n_species
878 where(w(ixo^s, dust_rho(n))>0.0d0)
879 fd(ixo^s) = dust_k_lineardrag/(w(ixo^s,gas_rho_)*w(ixo^s, dust_rho(n)))
880 else where
881 fd(ixo^s) = 0.0d0
882 end where
883 alpha(ixo^s, idir,n) = fd(ixo^s)
884 end do
885 end do
886
887 case('none')
888 alpha(ixo^s, :, :) = 0.0d0
889 case default
890 call mpistop( "=== This dust method has not been implemented===" )
891 end select
892
893 end subroutine get_alpha_dust
894
895
896 !> Get dt related to dust and gas stopping time (Laibe 2011)
897 subroutine dust_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
900
901 integer, intent(in) :: ixi^l, ixo^l
902 double precision, intent(in) :: dx^d, x(ixi^s, 1:ndim)
903 double precision, intent(in) :: w(ixi^s, 1:nw)
904 double precision, intent(inout) :: dtnew
905
906 double precision :: ptherm(ixi^s), vgas(ixi^s, 1:ndir)
907 double precision, dimension(1:dust_n_species):: dtdust
908 double precision, dimension(ixI^S) :: vt2, deltav, tstop, vdust
909 double precision, dimension(ixI^S, 1:dust_n_species) :: alpha_t
910 integer :: n, idir
911
912 if(dust_dtpar .le. 0) return
913
914 call phys_get_pthermal(w, x, ixi^l, ixo^l, ptherm)
915 do idir = 1, ndir
916 vgas(ixo^s,idir)=w(ixo^s,gas_mom(idir))/w(ixo^s,gas_rho_)
917 end do
918
919 select case( trim(dust_method) )
920
921 case( 'Kwok' ) ! assume sticking coefficient equals 0.25
922 dtdust(:) = bigdouble
923
924 vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/w(ixo^s, gas_rho_)
925
926 do idir = 1, ndir
927 do n = 1, dust_n_species
928 where(w(ixo^s, dust_rho(n))>0.0d0)
929 vdust(ixo^s) = w(ixo^s,dust_mom(idir, n))/w(ixo^s, dust_rho(n))
930 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
931 tstop(ixo^s) = 4.0d0*(dust_density(n)*dust_size(n))/ &
932 (3.0d0*(0.75d0)*dsqrt(vt2(ixo^s) + &
933 deltav(ixo^s)**2)*(w(ixo^s, dust_rho(n)) + &
934 w(ixo^s, gas_rho_)))
935 else where
936 tstop(ixo^s) = bigdouble
937 end where
938
939 dtdust(n) = min(minval(tstop(ixo^s)), dtdust(n))
940 end do
941 end do
942
943 dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
944
945 case( 'sticking' ) ! Calculate sticking coefficient based on the gas temperature
946 dtdust(:) = bigdouble
947
948 vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/w(ixo^s, gas_rho_)
949
950 ! Sticking coefficient
951 call get_sticking(w, x, ixi^l, ixo^l, alpha_t, ptherm)
952
953 do idir = 1, ndir
954 do n = 1, dust_n_species
955 where(w(ixo^s, dust_rho(n))>0.0d0)
956 vdust(ixo^s) = w(ixo^s,dust_mom(idir, n))/w(ixo^s, dust_rho(n))
957 deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
958 tstop(ixo^s) = 4.0d0*(dust_density(n)*dust_size(n))/ &
959 (3.0d0*(one-alpha_t(ixo^s,n))*dsqrt(vt2(ixo^s) + &
960 deltav(ixo^s)**2)*(w(ixo^s, dust_rho(n)) + &
961 w(ixo^s, gas_rho_)))
962 else where
963 tstop(ixo^s) = bigdouble
964 end where
965
966 dtdust(n) = min(minval(tstop(ixo^s)), dtdust(n))
967 end do
968 end do
969
970 dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
971
972 case('linear') !linear with Deltav, for testing (see Laibe & Price 2011)
973 dtdust(:) = bigdouble
974
975 do n = 1, dust_n_species
976 where(w(ixo^s, dust_rho(n))>0.0d0)
977 tstop(ixo^s) = (w(ixo^s, dust_rho(n))*w(ixo^s, gas_rho_))/ &
978 (dust_k_lineardrag*(w(ixo^s, dust_rho(n)) + w(ixo^s, gas_rho_)))
979 else where
980 tstop(ixo^s) = bigdouble
981 end where
982
983 dtdust(n) = min(minval(tstop(ixo^s)), dtdust(n))
984 end do
985
986 dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
987 case('usr')
988 dtdust(:) = bigdouble
989 call usr_dust_get_dt(w, ixi^l, ixo^l, dtdust, dx^d, x, dust_n_species)
990 dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
991 case('none')
992 ! no dust timestep
993 case default
994 call mpistop( "=== This dust method has not been implemented===" )
995 end select
996
997 if (dtnew < dtmin) then
998 write(unitterm,*)"-------------------------------------"
999 write(unitterm,*)"Warning: found DUST related time step too small! dtnew=", dtnew
1000 write(unitterm,*)"on grid with index:", block%igrid," grid level=", block%level
1001 write(unitterm,*)"grid corners are=",{^d&rnode(rpxmin^d_, block%igrid), rnode(rpxmax^d_, block%igrid)}
1002 write(unitterm,*)" dtdust =", dtdust(:)
1003 write(unitterm,*)"on processor:", mype
1004 write(unitterm,*)"-------------------------------------"
1005 endif
1006
1007 end subroutine dust_get_dt
1008
1009 ! Note that cmax and cmin are assumed to be initialized
1010 subroutine dust_get_cmax(w, x, ixI^L, ixO^L, idim, cmax, cmin)
1012 use mod_variables
1013
1014 integer, intent(in) :: ixi^l, ixo^l, idim
1015 double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
1016 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
1017 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
1018 double precision :: vdust(ixo^s)
1019 integer :: n
1020
1021 do n = 1, dust_n_species
1022 vdust(ixo^s) = get_vdust(w, ixi^l, ixo^l, idim, n)
1023
1024 if (present(cmin)) then
1025 cmin(ixo^s,1) = min(cmin(ixo^s,1), vdust(ixo^s))
1026 cmax(ixo^s,1) = max(cmax(ixo^s,1), vdust(ixo^s))
1027 else
1028 cmax(ixo^s,1) = max(cmax(ixo^s,1), abs(vdust(ixo^s)))
1029 end if
1030 end do
1031 end subroutine dust_get_cmax
1032
1033 ! Note that cmax and cmin are assumed to be initialized
1034 subroutine dust_get_cmax_prim(w, x, ixI^L, ixO^L, idim, cmax, cmin)
1036
1037 integer, intent(in) :: ixi^l, ixo^l, idim
1038 double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
1039 double precision, intent(inout) :: cmax(ixi^s)
1040 double precision, intent(inout), optional :: cmin(ixi^s)
1041 double precision :: vdust(ixo^s)
1042 integer :: n
1043
1044 do n = 1, dust_n_species
1045 vdust(ixo^s) = get_vdust_prim(w, ixi^l, ixo^l, idim, n)
1046
1047 if (present(cmin)) then
1048 cmin(ixo^s) = min(cmin(ixo^s), vdust(ixo^s))
1049 cmax(ixo^s) = max(cmax(ixo^s), vdust(ixo^s))
1050 else
1051 cmax(ixo^s) = max(cmax(ixo^s), abs(vdust(ixo^s)))
1052 end if
1053 end do
1054 end subroutine dust_get_cmax_prim
1055
1056end 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:530
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:583
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_dt(w, ixil, ixol, dtnew, dxd, x)
Get dt related to dust and gas stopping time (Laibe 2011)
Definition mod_dust.t:898
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:1011
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:1035
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:652
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
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...