MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
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
3 module mod_dust
4  use mod_global_parameters, only: std_len
5  use mod_physics
6  use mod_comm_lib, only: mpistop
7 
8  implicit none
9  private
10 
11  !> The number of dust species
12  integer, public, protected :: dust_n_species = 0
13 
14  integer, protected :: gas_rho_ = -1
15  integer, allocatable, protected :: gas_mom(:)
16  integer, protected :: gas_e_ = -1
17 
18  !> Indices of the dust densities
19  integer, allocatable, public, protected :: dust_rho(:)
20 
21  !> Indices of the dust momentum densities
22  integer, allocatable, public, protected :: dust_mom(:, :)
23 
24  !> Size of each dust species, dimensionless expression
25  double precision, allocatable, public :: dust_size(:)
26 
27  !> Internal density of each dust species, dimensionless expression
28  double precision, allocatable, public :: dust_density(:)
29 
30  !> Reduction of stopping time timestep limit
31  double precision :: dust_dtpar = 0.5d0
32 
33  !> Factor used in squared thermal velocity
34  double precision :: gas_vtherm_factor = 3.0d0
35 
36  !> Dust temperature in K (if dust_temperature_type is constant)
37  double precision :: dust_temperature = -1.0d0
38 
39  !> Dust drag coefficient for linear drag (for testing dust_method=linear)
40  double precision :: dust_k_lineardrag = -1.0d0
41 
42  !> If dust_temperature_type is stellar, it will be calculated according to Tielens (2005),
43  !> eqn. 5.44 using an input stellar luminosity in solar luminosities
44  double precision :: dust_stellar_luminosity = -1.0d0
45 
46  !> Set small dust densities to zero to avoid numerical problems
47  logical, public, protected :: dust_small_to_zero = .false.
48 
49  !> Minimum dust density as used when dust_small_to_zero=T
50  double precision, public, protected :: dust_min_rho = -1.0d0
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  !> What type of dust drag force to use. Can be 'Kwok', 'sticking', 'linear', 'usr' or 'none'.
59  character(len=std_len), public, protected :: dust_method = 'Kwok'
60 
61  !> Can be 'graphite' or 'silicate', affects the dust temperature
62  character(len=std_len) :: dust_species = 'graphite'
63 
64  !> Determines the dust temperature, can be 'constant', 'ism', or 'stellar'
65  character(len=std_len) :: dust_temperature_type = 'constant'
66 
67  !> whether second order terms (relevant only when dust_n_species >=2) are included
68  !> there are the terms n2, ni2, d2 in Eqs 6,7,8 in amrvac 3.0 paper
69  logical :: dust_implicit_second_order = .true.
70 
71  !> whether fh is added for gas energy: is only added in the impliict implementation, the explicit one was left as before
72  logical :: dust_backreaction_fh = .false.
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
89  public :: dust_evaluate_implicit
90 
91 
92 contains
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 
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))
112  allocate(dust_density(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)
136  use mod_global_parameters, only: unitpar
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)
148 111 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  logical :: flag(ixi^s)
333  integer :: n, idir
334 
335  do n = 1, dust_n_species
336  flag(ixo^s)=(w(ixo^s, dust_rho(n)) <= dust_min_rho)
337  where (flag(ixo^s))
338  w(ixo^s, dust_rho(n)) = 0.0d0
339  end where
340  do idir = 1, ndir
341  where (flag(ixo^s))
342  w(ixo^s, dust_mom(idir, n)) = 0.0d0
343  end where
344  end do
345  end do
346 
347  end subroutine set_dusttozero
348 
349  ! Calculate drag force based on Epstein's law
350  ! From Kwok 1975, page 584 (between eqn 8 and 9)
351  subroutine get_3d_dragforce(ixI^L, ixO^L, w, x, fdrag, ptherm, vgas)
354  integer, intent(in) :: ixI^L, ixO^L
355  double precision, intent(in) :: x(ixI^S, 1:ndim)
356  double precision, intent(in) :: w(ixI^S, 1:nw)
357  double precision, intent(out) :: &
358  fdrag(ixI^S, 1:ndir, 1:dust_n_species)
359  double precision, intent(in) :: ptherm(ixI^S), vgas(ixI^S, 1:ndir)
360 
361  double precision, dimension(ixI^S) :: vt2, deltav, fd, vdust
362  double precision :: alpha_T(ixI^S, 1:dust_n_species)
363  integer :: n, idir
364 
365  vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/w(ixo^s, gas_rho_)
366 
367  select case( trim(dust_method) )
368  case ('Kwok') ! assume sticking coefficient equals 0.25
369 
370  do idir = 1, ndir
371  do n = 1, dust_n_species
372  where(w(ixo^s, dust_rho(n)) > 0.0d0)
373  vdust(ixo^s) = w(ixo^s, dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
374  deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
375 
376  ! 0.75 from sticking coefficient
377  fd(ixo^s) = 0.75d0*w(ixo^s, dust_rho(n))*w(ixo^s, gas_rho_)*deltav(ixo^s) &
378  / (dust_density(n) * dust_size(n))
379 
380  ! 0.75 from spherical grainvolume
381  fd(ixo^s) = -fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
382  elsewhere
383  fd(ixo^s) = 0.0d0
384  end where
385  fdrag(ixo^s, idir, n) = fd(ixo^s)
386  end do
387  end do
388 
389  case ('sticking') ! Calculate sticking coefficient based on the gas and dust temperatures
390 
391  call get_sticking(w, x, ixi^l, ixo^l, alpha_t, ptherm)
392 
393  do idir = 1, ndir
394  do n = 1, dust_n_species
395  where(w(ixo^s, dust_rho(n))>0.0d0)
396  vdust(ixo^s) = w(ixo^s,dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
397  deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
398  fd(ixo^s) = (one-alpha_t(ixo^s,n)) * w(ixo^s, dust_rho(n))*w(ixo^s, gas_rho_) * &
399  deltav(ixo^s) / (dust_density(n)*dust_size(n))
400  fd(ixo^s) = -fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
401  else where
402  fd(ixo^s) = 0.0d0
403  end where
404  fdrag(ixo^s, idir,n) = fd(ixo^s)
405  end do
406  end do
407 
408  case('linear') !linear with Deltav, for testing (see Laibe & Price 2011)
409  do idir = 1, ndir
410  do n = 1, dust_n_species
411  where(w(ixo^s, dust_rho(n))>0.0d0)
412  vdust(ixo^s) = w(ixo^s,dust_mom(idir, n))/w(ixo^s, dust_rho(n))
413  deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
414 
415  fd(ixo^s) = -dust_k_lineardrag*deltav(ixo^s)
416  else where
417  fd(ixo^s) = 0.0d0
418  end where
419  fdrag(ixo^s, idir,n) = fd(ixo^s)
420  end do
421  end do
422 
423  case('usr')
424  call usr_get_3d_dragforce(ixi^l, ixo^l, w, x, fdrag, ptherm, vgas, dust_n_species)
425  case('none')
426  fdrag(ixo^s, :, :) = 0.0d0
427  case default
428  call mpistop( "=== This dust method has not been implemented===" )
429  end select
430 
431  end subroutine get_3d_dragforce
432 
433  !> Get sticking coefficient alpha_T (always between 0 and 1)
434  !>
435  !> Uses Temperatures in K
436  !> Equation from Decin et al. 2006
437  subroutine get_sticking(w, x, ixI^L, ixO^L, alpha_T, ptherm)
439  integer, intent(in) :: ixI^L, ixO^L
440  double precision, intent(in) :: x(ixI^S, 1:ndim)
441  double precision, intent(in) :: w(ixI^S, 1:nw)
442  double precision, intent(out) :: alpha_T(ixI^S, 1:dust_n_species)
443  double precision, intent(in) :: ptherm(ixI^S)
444  double precision :: Tgas(ixI^S)
445  integer :: n
446 
447  ! get the dust species T in K
448  call get_tdust(w, x, ixi^l, ixo^l, alpha_t)
449 
450  ! convert dimensionless gas T to K
451  tgas(ixo^s) = (ptherm(ixo^s)/w(ixo^s, gas_rho_))*unit_temperature
452 
453  do n = 1, dust_n_species
454  alpha_t(ixo^s,n) = 0.35d0 * dexp(-dsqrt((tgas(ixo^s) + &
455  alpha_t(ixo^s,n))/5.0d2))+0.1d0
456  end do
457 
458  end subroutine get_sticking
459 
460  !> Returns dust temperature (in K), either as constant or based on equ. 5.41,
461  !> 5.42 and 5.44 from Tielens (2005)
462  !>
463  !> Note that this calculation assumes cgs!!!!
464  !>
465  !> It takes as input the stellar luminosity in solar units in 'stellar' case
466  !> or a fixed dust input temperature in Kelvin when 'constant' or does case 'ism'
467  subroutine get_tdust(w, x, ixI^L, ixO^L, Td)
469  use mod_geometry
470 
471  integer, intent(in) :: ixI^L, ixO^L
472  double precision, intent(in) :: x(ixI^S, 1:ndim)
473  double precision, intent(in) :: w(ixI^S, 1:nw)
474  double precision, intent(out) :: Td(ixI^S, 1:dust_n_species)
475  double precision :: G0(ixO^S)
476  integer :: n
477 
478  select case( trim(dust_temperature_type) )
479  case( 'constant' )
480  td(ixo^s, :) = dust_temperature
481  case( 'ism' )
482  select case( trim(dust_species) )
483  case( 'graphite' )
484  do n = 1, dust_n_species
485  td(ixo^s, n) = 15.8d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0)
486  end do
487  case( 'silicate' )
488  do n = 1, dust_n_species
489  td(ixo^s, n) = 13.6d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0)
490  end do
491  case default
492  call mpistop( "=== Dust species undetermined===" )
493  end select
494  case( 'stellar' )
495  select case(coordinate)
496  case(spherical)
497  g0(ixo^s) = max(x(ixo^s, 1)*unit_length, smalldouble)
498  !!!case(cylindrical) convert R,Z to spherical radial coordinate r here
499  !!! but only ok for 2D (R,Z) or 2.5D (R,Z) case
500  !!! G0(ixO^S) = max(dsqrt(sum(x(ixO^S,:)**2,dim=ndim+1))*unit_length, smalldouble)
501  case default
502  call mpistop('stellar case not available in this coordinate system')
503  end select
504 
505  g0(ixo^s) = 2.1d4*(dust_stellar_luminosity/1.0d8)*((3.0857d17/g0(ixo^s))**2)
506 
507  select case( trim(dust_species) )
508  case( 'graphite' )
509  do n = 1, dust_n_species
510  td(ixo^s, n) = 61.0d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0) &
511  *(g0(ixo^s)**(one/5.8d0))
512  end do
513  case( 'silicate' )
514  do n = 1, dust_n_species
515  td(ixo^s, n) = 50.0d0*((0.0001d0/(dust_size(n)*unit_length))**0.06d0) &
516  *(g0(ixo^s)**(one/6.0d0))
517  end do
518  case default
519  call mpistop( "=== Dust species undetermined===" )
520  end select
521  case default
522  call mpistop( "=== Dust temperature undetermined===" )
523  end select
524 
525  end subroutine get_tdust
526 
527  !> w[iw]= w[iw]+qdt*S[wCT, x] where S is the source based on wCT within ixO
528  subroutine dust_add_source(qdt, ixI^L, ixO^L, wCT, w, x, qsourcesplit, active)
530 
531  integer, intent(in) :: ixi^l, ixo^l
532  double precision, intent(in) :: qdt
533  double precision, intent(in) :: wct(ixi^s, 1:nw), x(ixi^s, 1:ndim)
534  double precision, intent(inout) :: w(ixi^s, 1:nw)
535  logical, intent(in) :: qsourcesplit
536  logical, intent(inout) :: active
537 
538  double precision :: ptherm(ixi^s), vgas(ixi^s, 1:ndir)
539  double precision :: fdrag(ixi^s, 1:ndir, 1:dust_n_species)
540  integer :: n, idir
541 
542  select case( trim(dust_method) )
543  case( 'none' )
544  !do nothing here
545  case default !all regular dust methods here
546  if (qsourcesplit .eqv. dust_source_split) then
547  active = .true.
548 
549  call phys_get_pthermal(wct, x, ixi^l, ixo^l, ptherm)
550  do idir=1,ndir
551  vgas(ixo^s,idir)=wct(ixo^s,gas_mom(idir))/wct(ixo^s,gas_rho_)
552  end do
553 
554  call get_3d_dragforce(ixi^l, ixo^l, wct, x, fdrag, ptherm, vgas)
555  fdrag(ixo^s, 1:ndir, 1:dust_n_species) = fdrag(ixo^s, 1:ndir, 1:dust_n_species) * qdt
556 
557  do idir = 1, ndir
558 
559  do n = 1, dust_n_species
560  if (dust_backreaction) then
561  w(ixo^s, gas_mom(idir)) = w(ixo^s, gas_mom(idir)) + &
562  fdrag(ixo^s, idir, n)
563  if (gas_e_ > 0) then
564  w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + vgas(ixo^s, idir) &
565  * fdrag(ixo^s, idir, n)
566  end if
567  end if
568 
569 
570  w(ixo^s, dust_mom(idir, n)) = w(ixo^s, dust_mom(idir, n)) - &
571  fdrag(ixo^s, idir, n)
572  end do
573  end do
574 
575  endif
576  end select
577 
578  end subroutine dust_add_source
579 
580  !> inplace update of psa==>F_im(psa)
581  subroutine dust_evaluate_implicit(qtC,psa)
583  type(state), target :: psa(max_blocks)
584  double precision, intent(in) :: qtc
585 
586  integer :: iigrid, igrid, level
587 
588  !dust_method = 'none' not used
589 
590  !$OMP PARALLEL DO PRIVATE(igrid)
591  do iigrid=1,igridstail; igrid=igrids(iigrid);
592  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
593  block=>psa(igrid)
594  call dust_terms(ixg^ll,ixm^ll,psa(igrid)%w,psa(igrid)%x)
595  end do
596  !$OMP END PARALLEL DO
597 
598  end subroutine dust_evaluate_implicit
599 
600 
601 
602 
603  subroutine dust_terms(ixI^L,ixO^L,w,x)
605  integer, intent(in) :: ixI^L, ixO^L
606  double precision, intent(inout) :: w(ixI^S, 1:nw)
607  double precision, intent(in) :: x(ixI^S,1:ndim)
608 
609  double precision :: tmp(ixI^S), vgas(ixI^S, 1:ndir)
610  double precision :: alpha(ixI^S, 1:ndir, 1:dust_n_species)
611  integer :: n, idir
612 
613  do idir=1,ndir
614  vgas(ixo^s,idir)=w(ixo^s,gas_mom(idir))/w(ixo^s,gas_rho_)
615  end do
616  call get_alpha_dust(ixi^l, ixo^l, w, vgas,x, alpha)
617  w(ixo^s, gas_e_)=0d0
618  do idir = 1, ndir
619 
620  w(ixo^s, gas_mom(idir))=0d0
621  do n = 1, dust_n_species
622  ! contribution for gas momentum
623  tmp(ixo^s) = alpha(ixo^s, idir,n) * ( &
624  w(ixo^s,dust_rho(n)) * w(ixo^s, gas_mom(idir)) - &
625  w(ixo^s,gas_rho_) * w(ixo^s, dust_mom(idir, n)))
626  w(ixo^s, dust_mom(idir, n)) = -tmp(ixo^s)
627  if (dust_backreaction) then
628  w(ixo^s, gas_mom(idir)) = w(ixo^s, gas_mom(idir)) + tmp(ixo^s)
629  if (gas_e_ > 0) then
630  if(dust_backreaction_fh) then
631  where(w(ixo^s,dust_rho(n)) > 0d0)
632  w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + alpha(ixo^s, idir,n) * &
633  (w(ixo^s,gas_rho_) * (w(ixo^s, dust_mom(idir,n))**2/w(ixo^s,dust_rho(n))) - &
634  w(ixo^s,dust_rho(n)) * (w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_)))
635  elsewhere
636  w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + alpha(ixo^s, idir,n) * ( - &
637  w(ixo^s,dust_rho(n)) * (w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_)))
638  endwhere
639  else
640  w(ixo^s, gas_e_) = w(ixo^s, gas_e_) + vgas(ixo^s, idir) &
641  * tmp(ixo^s)
642  end if
643  end if
644  end if
645  end do
646  end do
647  end subroutine dust_terms
648 
649  !> Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
650  subroutine dust_implicit_update(dtfactor,qdt,qtC,psb,psa)
652  !use mod_ghostcells_update
653 
654  type(state), target :: psa(max_blocks)
655  type(state), target :: psb(max_blocks)
656  double precision, intent(in) :: qdt
657  double precision, intent(in) :: qtc
658  double precision, intent(in) :: dtfactor
659 
660  integer :: iigrid, igrid
661 
662  !call getbc(global_time,0.d0,psa,1,nw)
663  !$OMP PARALLEL DO PRIVATE(igrid)
664  do iigrid=1,igridstail; igrid=igrids(iigrid);
665  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
666  block=>psa(igrid)
667  call dust_advance_implicit_grid(ixg^ll, ixg^ll, psa(igrid)%w, psb(igrid)%w, psa(igrid)%x, dtfactor,qdt)
668  end do
669  !$OMP END PARALLEL DO
670 
671  end subroutine dust_implicit_update
672 
673  subroutine dust_advance_implicit_grid(ixI^L, ixO^L, w, wout, x, dtfactor,qdt)
675  integer, intent(in) :: ixI^L, ixO^L
676  double precision, intent(in) :: qdt
677  double precision, intent(in) :: dtfactor
678  double precision, intent(in) :: w(ixI^S,1:nw)
679  double precision, intent(in) :: x(ixI^S,1:ndim)
680  double precision, intent(out) :: wout(ixI^S,1:nw)
681 
682  integer :: n, m, idir
683  double precision :: alpha(ixI^S, 1:ndir, 1:dust_n_species)
684  double precision :: tmp(ixI^S),tmp2(ixI^S)
685  double precision :: tmp3(ixI^S)
686  double precision :: vgas(ixI^S, 1:ndir)
687 
688 
689  do idir = 1, ndir
690  vgas(ixo^s,idir)=w(ixo^s,gas_mom(idir))/w(ixo^s,gas_rho_)
691  end do
692  call get_alpha_dust(ixi^l, ixo^l, w, vgas, x, alpha)
693  !TODO this is still neeed?
694  wout(ixo^s,1:nw) = w(ixo^s,1:nw)
695 
696  do idir = 1, ndir
697  ! d1 from Eq 6
698  tmp2(ixo^s) = 0d0
699  do n = 1, dust_n_species
700  tmp2(ixo^s) = tmp2(ixo^s) + alpha(ixo^s, idir,n) * &
701  (w(ixo^s,gas_rho_) + w(ixo^s,dust_rho(n)))
702 
703  enddo
704  ! store D in tmp
705  tmp(ixo^s) = 1d0 + tmp2(ixo^s) * qdt
706  if(dust_implicit_second_order) then
707  ! d2 from Eq 6
708  tmp2(ixo^s) = 0d0
709  do n = 1, dust_n_species
710  do m = n+1, dust_n_species
711  tmp2(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) *&
712  (w(ixo^s,gas_rho_) + w(ixo^s,dust_rho(n))+w(ixo^s,dust_rho(m)))
713  enddo
714  enddo
715  ! multiplied at the end by rho_gas
716  tmp(ixo^s) = tmp(ixo^s) + w(ixo^s,gas_rho_)*tmp2(ixo^s) * (qdt**2)
717  endif
718 
719 
720 
721  do n = 1, dust_n_species
722  ! ni1 from eq 7
723  tmp2(ixo^s) = alpha(ixo^s, idir,n) * ( &
724  w(ixo^s,dust_rho(n)) * w(ixo^s, gas_mom(idir)) - &
725  w(ixo^s,gas_rho_) * w(ixo^s, dust_mom(idir, n))) * qdt
726 
727  if(dust_implicit_second_order) then
728  ! ni2 from eq 7
729  tmp3(ixo^s) = 0d0
730  do m = n+1, dust_n_species
731  tmp3(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) * &
732  ( w(ixo^s,dust_rho(n)) * (w(ixo^s, dust_mom(idir, n)) + w(ixo^s, gas_mom(idir))) - &
733  (w(ixo^s,gas_rho_) + w(ixo^s,dust_rho(m))) * w(ixo^s, dust_mom(idir, n)) )
734  enddo
735  ! tmp3 multiplied at the end by rho_gas
736  tmp2(ixo^s) = tmp2(ixo^s) + tmp3(ixo^s) * w(ixo^s,gas_rho_)* (qdt**2)
737  endif
738  tmp2(ixo^s) = tmp2(ixo^s)/tmp(ixo^s)
739  wout(ixo^s, dust_mom(idir,n)) = w(ixo^s, dust_mom(idir,n)) + tmp2(ixo^s)
740  enddo
741 
742  if (dust_backreaction) then
743  tmp2(ixo^s) = 0d0
744  !n1 from eq 8
745  do n = 1, dust_n_species
746  tmp2(ixo^s) = tmp2(ixo^s) + alpha(ixo^s, idir,n) * &
747  (w(ixo^s,gas_rho_) * w(ixo^s, dust_mom(idir,n)) - &
748  w(ixo^s,dust_rho(n)) * w(ixo^s, gas_mom(idir)))
749 
750  enddo
751  tmp2(ixo^s) = qdt * tmp2(ixo^s)
752  if(dust_implicit_second_order) then
753  !n2 from eq 8
754  tmp3(ixo^s) = 0d0
755  do n = 1, dust_n_species
756  do m = n+1, dust_n_species
757  tmp3(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) * &
758  (w(ixo^s,gas_rho_) * (w(ixo^s, dust_mom(idir, n)) + w(ixo^s, dust_mom(idir, m))) - &
759  (w(ixo^s,dust_rho(n)) + w(ixo^s,dust_rho(m)))* w(ixo^s, gas_mom(idir)))
760  enddo
761  enddo
762  ! tmp3 multiplied at the end by rho_gas
763  tmp2(ixo^s) = tmp2(ixo^s) + (qdt**2)*tmp3(ixo^s)* w(ixo^s,gas_rho_)
764  endif
765  ! store in tmp2 contribution to momentum
766  ! so that it is used when dust_backreaction_fh = .false.
767  tmp2(ixo^s) = tmp2(ixo^s) / tmp(ixo^s)
768  wout(ixo^s, gas_mom(idir)) = w(ixo^s, gas_mom(idir)) + tmp2(ixo^s)
769 
770  ! kinetic energy update
771  if (gas_e_ > 0) then
772  if(dust_backreaction_fh) then
773  ! add work done by coll terms + FrictionalHeating
774  tmp2(ixo^s) = 0d0
775  do n = 1, dust_n_species
776  ! 2*dust kinetic energy: dust rho can be 0
777  where(w(ixo^s,dust_rho(n)) > 0d0)
778  tmp3(ixo^s)= w(ixo^s, dust_mom(idir,n))**2/w(ixo^s,dust_rho(n))
779  elsewhere
780  tmp3(ixo^s) = 0d0
781  endwhere
782  tmp2(ixo^s) = tmp2(ixo^s) + alpha(ixo^s, idir,n) * &
783  (w(ixo^s,gas_rho_) * tmp3(ixo^s) - &
784  w(ixo^s,dust_rho(n)) * (w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_)))
785 
786  enddo
787  tmp2(ixo^s) = qdt * tmp2(ixo^s)
788  if(dust_implicit_second_order) then
789  tmp3(ixo^s) = 0d0
790  do n = 1, dust_n_species
791  do m = n+1, dust_n_species
792  tmp3(ixo^s) = tmp3(ixo^s) + alpha(ixo^s, idir,n) * alpha(ixo^s, idir,m) * &
793  (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))) - &
794  (w(ixo^s,dust_rho(n)) + w(ixo^s,dust_rho(m)))* w(ixo^s, gas_mom(idir))**2/w(ixo^s,gas_rho_))
795  enddo
796  enddo
797  ! tmp3 multiplied at the end by rho_gas
798  tmp2(ixo^s) = tmp2(ixo^s) + (qdt**2)*tmp3(ixo^s)* w(ixo^s,gas_rho_)
799  endif
800  wout(ixo^s, gas_e_) = wout(ixo^s, gas_e_) + 0.5d0 * tmp2(ixo^s) / tmp(ixo^s)
801  else
802  ! dust_backreaction_fh = .false.
803  ! add only work done by coll term by multiplyting the contribution in mom eq. by velocity
804  wout(ixo^s, gas_e_) = wout(ixo^s, gas_e_) + vgas(ixo^s, idir) * tmp2(ixo^s)
805  endif
806  end if
807  end if
808  end do !1..ndir
809 
810 
811  end subroutine dust_advance_implicit_grid
812 
813  ! copied from get_3d_dragforce subroutine
814  subroutine get_alpha_dust(ixI^L, ixO^L, w, vgas,x, alpha)
817  integer, intent(in) :: ixI^L, ixO^L
818  double precision, intent(in) :: x(ixI^S, 1:ndim)
819  double precision, intent(in) :: w(ixI^S, 1:nw)
820  double precision,intent(in) :: vgas(ixI^S, 1:ndir)
821  double precision, intent(out) :: &
822  alpha(ixI^S, 1:ndir, 1:dust_n_species)
823 
824  double precision :: ptherm(ixI^S)
825  double precision, dimension(ixI^S) :: vt2, deltav, fd, vdust
826  double precision :: alpha_T(ixI^S, 1:dust_n_species)
827  integer :: n, idir
828 
829  call phys_get_pthermal(w, x, ixi^l, ixo^l, ptherm)
830 
831  vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/w(ixo^s, gas_rho_)
832 
833  select case( trim(dust_method) )
834  case ('Kwok') ! assume sticking coefficient equals 0.25
835 
836  do idir = 1, ndir
837  do n = 1, dust_n_species
838  where(w(ixo^s, dust_rho(n)) > 0.0d0)
839 
840  ! 0.75 from sticking coefficient
841  fd(ixo^s) = 0.75d0 / (dust_density(n) * dust_size(n))
842 
843  ! 0.75 from spherical grainvolume
844  vdust(ixo^s) = w(ixo^s, dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
845  deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
846  fd(ixo^s) = fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
847  elsewhere
848  fd(ixo^s) = 0.0d0
849  end where
850  alpha(ixo^s, idir, n) = fd(ixo^s)
851  end do
852  end do
853 
854  case ('sticking') ! Calculate sticking coefficient based on the gas and dust temperatures
855 
856  call get_sticking(w, x, ixi^l, ixo^l, alpha_t, ptherm)
857 
858  do idir = 1, ndir
859  do n = 1, dust_n_species
860  where(w(ixo^s, dust_rho(n))>0.0d0)
861  ! sticking
862  fd(ixo^s) = (one-alpha_t(ixo^s,n)) / (dust_density(n)*dust_size(n))
863  ! 0.75 from spherical grainvolume
864  vdust(ixo^s) = w(ixo^s,dust_mom(idir, n)) / w(ixo^s, dust_rho(n))
865  deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
866  fd(ixo^s) = fd(ixo^s)*0.75d0*dsqrt(vt2(ixo^s) + deltav(ixo^s)**2)
867  else where
868  fd(ixo^s) = 0.0d0
869  end where
870  alpha(ixo^s, idir,n) = fd(ixo^s)
871  end do
872  end do
873 
874  case('linear') !linear with Deltav, for testing (see Laibe & Price 2011)
875  do idir = 1, ndir
876  do n = 1, dust_n_species
877  where(w(ixo^s, dust_rho(n))>0.0d0)
878  fd(ixo^s) = dust_k_lineardrag/(w(ixo^s,gas_rho_)*w(ixo^s, dust_rho(n)))
879  else where
880  fd(ixo^s) = 0.0d0
881  end where
882  alpha(ixo^s, idir,n) = fd(ixo^s)
883  end do
884  end do
885 
886  case('none')
887  alpha(ixo^s, :, :) = 0.0d0
888  case default
889  call mpistop( "=== This dust method has not been implemented===" )
890  end select
891 
892  end subroutine get_alpha_dust
893 
894 
895  !> Get dt related to dust and gas stopping time (Laibe 2011)
896  subroutine dust_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
899 
900  integer, intent(in) :: ixi^l, ixo^l
901  double precision, intent(in) :: dx^d, x(ixi^s, 1:ndim)
902  double precision, intent(in) :: w(ixi^s, 1:nw)
903  double precision, intent(inout) :: dtnew
904 
905  double precision :: ptherm(ixi^s), vgas(ixi^s, 1:ndir)
906  double precision, dimension(1:dust_n_species):: dtdust
907  double precision, dimension(ixI^S) :: vt2, deltav, tstop, vdust
908  double precision, dimension(ixI^S, 1:dust_n_species) :: alpha_t
909  integer :: n, idir
910 
911  if(dust_dtpar .le. 0) return
912 
913  call phys_get_pthermal(w, x, ixi^l, ixo^l, ptherm)
914  do idir = 1, ndir
915  vgas(ixo^s,idir)=w(ixo^s,gas_mom(idir))/w(ixo^s,gas_rho_)
916  end do
917 
918  select case( trim(dust_method) )
919 
920  case( 'Kwok' ) ! assume sticking coefficient equals 0.25
921  dtdust(:) = bigdouble
922 
923  vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/w(ixo^s, gas_rho_)
924 
925  do idir = 1, ndir
926  do n = 1, dust_n_species
927  where(w(ixo^s, dust_rho(n))>0.0d0)
928  vdust(ixo^s) = w(ixo^s,dust_mom(idir, n))/w(ixo^s, dust_rho(n))
929  deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
930  tstop(ixo^s) = 4.0d0*(dust_density(n)*dust_size(n))/ &
931  (3.0d0*(0.75d0)*dsqrt(vt2(ixo^s) + &
932  deltav(ixo^s)**2)*(w(ixo^s, dust_rho(n)) + &
933  w(ixo^s, gas_rho_)))
934  else where
935  tstop(ixo^s) = bigdouble
936  end where
937 
938  dtdust(n) = min(minval(tstop(ixo^s)), dtdust(n))
939  end do
940  end do
941 
942  dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
943 
944  case( 'sticking' ) ! Calculate sticking coefficient based on the gas temperature
945  dtdust(:) = bigdouble
946 
947  vt2(ixo^s) = gas_vtherm_factor*ptherm(ixo^s)/w(ixo^s, gas_rho_)
948 
949  ! Sticking coefficient
950  call get_sticking(w, x, ixi^l, ixo^l, alpha_t, ptherm)
951 
952  do idir = 1, ndir
953  do n = 1, dust_n_species
954  where(w(ixo^s, dust_rho(n))>0.0d0)
955  vdust(ixo^s) = w(ixo^s,dust_mom(idir, n))/w(ixo^s, dust_rho(n))
956  deltav(ixo^s) = vgas(ixo^s, idir)-vdust(ixo^s)
957  tstop(ixo^s) = 4.0d0*(dust_density(n)*dust_size(n))/ &
958  (3.0d0*(one-alpha_t(ixo^s,n))*dsqrt(vt2(ixo^s) + &
959  deltav(ixo^s)**2)*(w(ixo^s, dust_rho(n)) + &
960  w(ixo^s, gas_rho_)))
961  else where
962  tstop(ixo^s) = bigdouble
963  end where
964 
965  dtdust(n) = min(minval(tstop(ixo^s)), dtdust(n))
966  end do
967  end do
968 
969  dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
970 
971  case('linear') !linear with Deltav, for testing (see Laibe & Price 2011)
972  dtdust(:) = bigdouble
973 
974  do n = 1, dust_n_species
975  where(w(ixo^s, dust_rho(n))>0.0d0)
976  tstop(ixo^s) = (w(ixo^s, dust_rho(n))*w(ixo^s, gas_rho_))/ &
977  (dust_k_lineardrag*(w(ixo^s, dust_rho(n)) + w(ixo^s, gas_rho_)))
978  else where
979  tstop(ixo^s) = bigdouble
980  end where
981 
982  dtdust(n) = min(minval(tstop(ixo^s)), dtdust(n))
983  end do
984 
985  dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
986  case('usr')
987  dtdust(:) = bigdouble
988  call usr_dust_get_dt(w, ixi^l, ixo^l, dtdust, dx^d, x, dust_n_species)
989  dtnew = min(minval(dust_dtpar*dtdust(:)), dtnew)
990  case('none')
991  ! no dust timestep
992  case default
993  call mpistop( "=== This dust method has not been implemented===" )
994  end select
995 
996  if (dtnew < dtmin) then
997  write(unitterm,*)"-------------------------------------"
998  write(unitterm,*)"Warning: found DUST related time step too small! dtnew=", dtnew
999  write(unitterm,*)"on grid with index:", block%igrid," grid level=", block%level
1000  write(unitterm,*)"grid corners are=",{^d&rnode(rpxmin^d_, block%igrid), rnode(rpxmax^d_, block%igrid)}
1001  write(unitterm,*)" dtdust =", dtdust(:)
1002  write(unitterm,*)"on processor:", mype
1003  write(unitterm,*)"-------------------------------------"
1004  endif
1005 
1006  end subroutine dust_get_dt
1007 
1008  ! Note that cmax and cmin are assumed to be initialized
1009  subroutine dust_get_cmax(w, x, ixI^L, ixO^L, idim, cmax, cmin)
1011  use mod_variables
1012 
1013  integer, intent(in) :: ixi^l, ixo^l, idim
1014  double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
1015  double precision, intent(inout) :: cmax(ixi^s,1:number_species)
1016  double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
1017  double precision :: vdust(ixo^s)
1018  integer :: n
1019 
1020  do n = 1, dust_n_species
1021  vdust(ixo^s) = get_vdust(w, ixi^l, ixo^l, idim, n)
1022 
1023  if (present(cmin)) then
1024  cmin(ixo^s,1) = min(cmin(ixo^s,1), vdust(ixo^s))
1025  cmax(ixo^s,1) = max(cmax(ixo^s,1), vdust(ixo^s))
1026  else
1027  cmax(ixo^s,1) = max(cmax(ixo^s,1), abs(vdust(ixo^s)))
1028  end if
1029  end do
1030  end subroutine dust_get_cmax
1031 
1032  ! Note that cmax and cmin are assumed to be initialized
1033  subroutine dust_get_cmax_prim(w, x, ixI^L, ixO^L, idim, cmax, cmin)
1035 
1036  integer, intent(in) :: ixi^l, ixo^l, idim
1037  double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
1038  double precision, intent(inout) :: cmax(ixi^s)
1039  double precision, intent(inout), optional :: cmin(ixi^s)
1040  double precision :: vdust(ixo^s)
1041  integer :: n
1042 
1043  do n = 1, dust_n_species
1044  vdust(ixo^s) = get_vdust_prim(w, ixi^l, ixo^l, idim, n)
1045 
1046  if (present(cmin)) then
1047  cmin(ixo^s) = min(cmin(ixo^s), vdust(ixo^s))
1048  cmax(ixo^s) = max(cmax(ixo^s), vdust(ixo^s))
1049  else
1050  cmax(ixo^s) = max(cmax(ixo^s), abs(vdust(ixo^s)))
1051  end if
1052  end do
1053  end subroutine dust_get_cmax_prim
1054 
1055 end module mod_dust
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
Module for including dust species, which interact with the gas through a drag force.
Definition: mod_dust.t:3
double precision function, dimension(ixo^s) get_vdust_prim(w, ixIL, ixOL, idim, n)
Definition: mod_dust.t:312
double precision, public, protected dust_min_rho
Minimum dust density as used when dust_small_to_zero=T.
Definition: mod_dust.t:50
subroutine, public dust_check_w(ixIL, ixOL, w, flag)
Definition: mod_dust.t:194
subroutine dust_advance_implicit_grid(ixIL, ixOL, w, wout, x, dtfactor, qdt)
Definition: mod_dust.t:674
subroutine, public dust_get_flux(w, x, ixIL, ixOL, idim, f)
Definition: mod_dust.t:252
subroutine get_3d_dragforce(ixIL, ixOL, w, x, fdrag, ptherm, vgas)
Definition: mod_dust.t:352
subroutine, public dust_implicit_update(dtfactor, qdt, qtC, psb, psa)
Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
Definition: mod_dust.t:651
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:897
double precision, dimension(:), allocatable, public dust_size
Size of each dust species, dimensionless expression.
Definition: mod_dust.t:25
subroutine get_tdust(w, x, ixIL, ixOL, Td)
Returns dust temperature (in K), either as constant or based on equ. 5.41, 5.42 and 5....
Definition: mod_dust.t:468
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:59
subroutine dust_read_params(files)
Read dust_list module parameters from a file.
Definition: mod_dust.t:136
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:529
subroutine, public set_dusttozero(ixIL, ixOL, w, x)
Definition: mod_dust.t:327
subroutine, public dust_to_primitive(ixIL, ixOL, w, x)
Definition: mod_dust.t:228
subroutine get_sticking(w, x, ixIL, ixOL, alpha_T, ptherm)
Get sticking coefficient alpha_T (always between 0 and 1)
Definition: mod_dust.t:438
integer, dimension(:, :), allocatable, public, protected dust_mom
Indices of the dust momentum densities.
Definition: mod_dust.t:22
subroutine, public dust_get_cmax(w, x, ixIL, ixOL, idim, cmax, cmin)
Definition: mod_dust.t:1010
logical, public, protected dust_small_to_zero
Set small dust densities to zero to avoid numerical problems.
Definition: mod_dust.t:47
integer, public, protected dust_n_species
The number of dust species.
Definition: mod_dust.t:12
subroutine get_alpha_dust(ixIL, ixOL, w, vgas, x, alpha)
Definition: mod_dust.t:815
integer, dimension(:), allocatable, public, protected dust_rho
Indices of the dust densities.
Definition: mod_dust.t:19
double precision function, dimension(ixo^s) get_vdust(w, ixIL, ixOL, idim, n)
Definition: mod_dust.t:298
subroutine, public dust_get_flux_prim(w, x, ixIL, ixOL, idim, f)
Definition: mod_dust.t:275
subroutine, public dust_check_params()
Definition: mod_dust.t:154
subroutine, public dust_evaluate_implicit(qtC, psa)
inplace update of psa==>F_im(psa)
Definition: mod_dust.t:582
subroutine dust_terms(ixIL, ixOL, w, x)
Definition: mod_dust.t:604
subroutine, public dust_init(g_rho, g_mom, g_energy)
Definition: mod_dust.t:95
subroutine, public dust_get_cmax_prim(w, x, ixIL, ixOL, idim, cmax, cmin)
Definition: mod_dust.t:1034
double precision, dimension(:), allocatable, public dust_density
Internal density of each dust species, dimensionless expression.
Definition: mod_dust.t:28
subroutine, public dust_to_conserved(ixIL, ixOL, w, x)
Definition: mod_dust.t:208
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
integer coordinate
Definition: mod_geometry.t:7
integer, parameter spherical
Definition: mod_geometry.t:11
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.
integer, 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
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.
double precision, dimension(ndim) dxlevel
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.
Definition: mod_variables.t:23
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
Definition: mod_variables.t:71