MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_ard_phys.t
Go to the documentation of this file.
1 !> Module containing the physics routines for advection-reaction-diffusion equations
2 !>
3 !> This module can be seen as an extension of the reaction-diffusion (rd) module
4 !> and includes the same reaction systems and more: the Gray-Scott model, the
5 !> Schnakenberg model, the Brusselator model, the diffusive logistic equation,
6 !> an analytical testcase from "Numerical solution of time-dependent advection-
7 !> diffusion-reaction equations" by Hundsdorfer & Verwer, the Oregonator model,
8 !> the extended Brusselator model, the diffusive Lorenz system and the advection-
9 !> diffusion equation. See the documentation of the advection-reaction-diffusion
10 !> module for more information.
11 !>
12 !> An advection term can be aplied to these systems of the form:
13 !> nabla( (A1/adv_pow) * u^(adv_pow) ) (for the first unknown)
14 !> nabla( (A2/adv_pow) * v^(adv_pow) ) (for the second unknown, if applicable)
15 !> nabla( (A3/adv_pow) * w^(adv_pow) ) (for the third unknown, if applicable)
16 !>
17 !> IMEX methods are also supported. The implicit system is solved by a
18 !> multigrid solver coupled into MPI-AMRVAC.
19 !>
22  use mod_comm_lib, only: mpistop
23 
24  implicit none
25  private
26 
27  !> Indices of the unknowns
28  integer, protected, public :: u_ = 1
29  integer, protected, public :: v_ = 2 !< For 2 or more equations
30  integer, protected, public :: w_ = 3 !< For 3 or more equations
31 
32  !> Whether particles module is added
33  logical, public, protected :: ard_particles = .false.
34 
35  !> Parameter with which to multiply the reaction timestep restriction
36  double precision, public, protected :: dtreacpar = 0.5d0
37 
38  !> Name of the system to be solved
39  character(len=20), public, protected :: equation_name = "gray-scott"
40  integer :: number_of_species = 2
41  integer :: equation_type = 1
42  integer, parameter :: eq_gray_scott = 1 ! Gray-Scott model
43  integer, parameter :: eq_schnakenberg = 2 ! Schnakenberg model
44  integer, parameter :: eq_brusselator = 3 ! Brusselator model
45  integer, parameter :: eq_logistic = 4 ! Logistic equation
46  integer, parameter :: eq_analyt_hunds = 5
47  integer, parameter :: eq_belousov_fn = 6 ! Field-Noyes model, or Oregonator
48  integer, parameter :: eq_ext_brusselator = 7 ! Extended Brusselator
49  integer, parameter :: eq_lorenz = 8 ! Lorenz system
50  integer, parameter :: eq_no_reac = 9 ! Advection-diffusion equation
51 
52  !> Diffusion coefficient for first species (u)
53  double precision, public, protected :: d1 = 0.05d0
54  !> Diffusion coefficient for second species (v) (if applicable)
55  double precision, public, protected :: d2 = 1.0d0
56  !> Diffusion coefficient for third species (w) (if applicable)
57  double precision, public, protected :: d3 = 1.0d0
58 
59  !> Power of the unknown in the advection term (1 for linear)
60  integer, public, protected :: adv_pow = 1
61 
62  !> Advection coefficients for first species (u)
63  double precision, public, protected :: a1(^nd) = 0.0d0
64  !> Advection coefficients for second species (v) (if applicable)
65  double precision, public, protected :: a2(^nd) = 0.0d0
66  !> Advection coefficients for third species (w) (if applicable)
67  double precision, public, protected :: a3(^nd) = 0.0d0
68 
69  !> Parameters for Schnakenberg model
70  double precision, public, protected :: sb_alpha = 0.1305d0
71  double precision, public, protected :: sb_beta = 0.7695d0
72  double precision, public, protected :: sb_kappa = 100.0d0
73 
74  !> Feed rate for Gray-Scott model
75  double precision, public, protected :: gs_f = 0.046d0
76  !> Kill rate for Gray-Scott model
77  double precision, public, protected :: gs_k = 0.063d0
78 
79  !> Parameters for Brusselator model
80  double precision, public, protected :: br_a = 4.5d0
81  double precision, public, protected :: br_b = 8.0d0
82  double precision, public, protected :: br_c = 1.0d0
83  double precision, public, protected :: br_d = 1.0d0
84 
85  !> Parameter for logistic model (Fisher / KPP equation)
86  double precision, public, protected :: lg_lambda = 1.0d0
87 
88  !> Parameters for the Field-Noyes model of the Belousov-Zhabotinsky reaction
89  double precision, public, protected :: bzfn_epsilon = 1.0d0
90  double precision, public, protected :: bzfn_delta = 1.0d0
91  double precision, public, protected :: bzfn_lambda = 1.0d0
92  double precision, public, protected :: bzfn_mu = 1.0d0
93 
94  !> Parameter for Lorenz system (Rayleigh number)
95  double precision, public, protected :: lor_r = 28.0d0
96  !> Parameter for Lorenz system (Prandtl number)
97  double precision, public, protected :: lor_sigma = 10.0d0
98  !> Parameter for Lorenz system (aspect ratio of the convection rolls)
99  double precision, public, protected :: lor_b = 8.0d0 / 3.0d0
100 
101  !> Whether to handle the explicitly handled source term in split fashion
102  logical :: ard_source_split = .false.
103 
104  !> Boundary condition information for the multigrid method
105  type(mg_bc_t), public :: ard_mg_bc(3, mg_num_neighbors)
106 
107  ! Public methods
108  public :: ard_phys_init
109 
110 contains
111 
112  !> Read this module's parameters from a file
113  subroutine ard_params_read(files)
114  use mod_global_parameters, only: unitpar
115  character(len=*), intent(in) :: files(:)
116  integer :: n
117 
118  namelist /ard_list/ d1, d2, d3, adv_pow, a1, a2, a3, sb_alpha, sb_beta, sb_kappa, &
121  ard_source_split, dtreacpar
122 
123  do n = 1, size(files)
124  open(unitpar, file=trim(files(n)), status='old')
125  read(unitpar, ard_list, end=111)
126 111 close(unitpar)
127  end do
128 
129  !> Set the equation type and number of species
130  select case (equation_name)
131  case ("gray-scott")
132  equation_type = eq_gray_scott
133  number_of_species = 2
134  case ("schnakenberg")
135  equation_type = eq_schnakenberg
136  number_of_species = 2
137  case ("brusselator")
138  equation_type = eq_brusselator
139  number_of_species = 2
140  case ("logistic")
141  equation_type = eq_logistic
142  number_of_species = 1
143  case ("analyt_hunds")
144  equation_type = eq_analyt_hunds
145  number_of_species = 1
146  case ("belousov_fieldnoyes")
147  equation_type = eq_belousov_fn
148  number_of_species = 3
149  case ("ext_brusselator")
150  equation_type = eq_ext_brusselator
151  number_of_species = 3
152  case ("lorenz")
153  equation_type = eq_lorenz
154  number_of_species = 3
155  case ("no_reac")
156  equation_type = eq_no_reac
157  number_of_species = 1
158  case default
159  call mpistop("Unknown equation_name (valid: gray-scott, schnakenberg, ...)")
160  end select
161 
162  end subroutine ard_params_read
163 
164  !> Write this modules parameters to a snapshot
165  subroutine ard_write_info(fh)
167  integer, intent(in) :: fh
168  integer, parameter :: n_par = 0
169  integer, dimension(MPI_STATUS_SIZE) :: st
170  integer :: er
171  integer :: idim
172 
173  call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
174  end subroutine ard_write_info
175 
176  subroutine ard_phys_init()
178  use mod_physics
180  use mod_particles, only: particles_init
181 
182  call ard_params_read(par_files)
183 
184  physics_type = "ard"
185  phys_energy = .false.
186  ! Whether diagonal ghost cells are required for the physics
187  phys_req_diagonal = .false.
189 
190  allocate(start_indices(number_species),stop_indices(number_species))
191  ! set the index of the first flux variable for species 1
192  start_indices(1)=1
193  ! Use the first variable as a density
194  u_ = var_set_fluxvar("u", "u")
195  if (number_of_species >= 2) then
196  v_ = var_set_fluxvar("v", "v")
197  end if
198  if (number_of_species >= 3) then
199  w_ = var_set_fluxvar("w", "w")
200  end if
201 
202  ! set number of variables which need update ghostcells
203  nwgc=nwflux
204  ! set the index of the last flux variable for species 1
205  stop_indices(1)=nwflux
206 
207  ! Check whether custom flux types have been defined
208  if (.not. allocated(flux_type)) then
209  allocate(flux_type(ndir, nw))
211  else if (any(shape(flux_type) /= [ndir, nw])) then
212  call mpistop("phys_check error: flux_type has wrong shape")
213  end if
214 
227 
228  ! Initialize particles module
229  if (ard_particles) then
230  call particles_init()
231  phys_req_diagonal = .true.
232  end if
233 
234  end subroutine ard_phys_init
235 
236  subroutine ard_check_params
238  integer :: n, i, iw, species_list(number_of_species)
239 
240  if (use_imex_scheme) then
241  use_multigrid = .true.
242  select case(number_of_species)
243  case(1)
244  species_list = [u_]
245  if (d1 == 0.0d0) then
246  write(*, *) "Diffusion coefficient cannot be zero in IMEX scheme"
247  call mpistop("Zero diffusion in IMEX scheme")
248  end if
249  case(2)
250  species_list = [u_, v_]
251  if ((d1 == 0.0d0) .or. (d2 == 0.0d0)) then
252  write(*, *) "Diffusion coefficient cannot be zero in IMEX scheme"
253  call mpistop("Zero diffusion in IMEX scheme")
254  end if
255  case(3)
256  species_list = [u_, v_, w_]
257  if ((d1 == 0.0d0) .or. (d2 == 0.0d0) .or. (d3 == 0.0d0)) then
258  write(*, *) "Diffusion coefficient cannot be zero in IMEX scheme"
259  call mpistop("Zero diffusion in IMEX scheme")
260  end if
261  end select
262 
263  do i = 1, number_of_species
264  iw = species_list(i)
265 
266  ! Set boundary conditions for the multigrid solver
267  do n = 1, 2*ndim
268  select case (typeboundary(iw, n))
269  case (bc_symm)
270  ! d/dx u = 0
271  ard_mg_bc(i, n)%bc_type = mg_bc_neumann
272  ard_mg_bc(i, n)%bc_value = 0.0_dp
273  case (bc_asymm)
274  ! u = 0
275  ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
276  ard_mg_bc(i, n)%bc_value = 0.0_dp
277  case (bc_cont)
278  ! d/dx u = 0
279  ard_mg_bc(i, n)%bc_type = mg_bc_neumann
280  ard_mg_bc(i, n)%bc_value = 0.0_dp
281  case (bc_periodic)
282  ! Nothing to do here
283  case (bc_special)
284  if (.not. associated(ard_mg_bc(i, n)%boundary_cond)) then
285  write(*, "(A,I0,A,I0,A)") "typeboundary(", iw, ",", n, &
286  ") is 'special', but the corresponding method " // &
287  "ard_mg_bc(i, n)%boundary_cond is not set"
288  call mpistop("ard_mg_bc(i, n)%boundary_cond not set")
289  end if
290  case default
291  write(*,*) "ard_check_params warning: unknown boundary type"
292  ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
293  ard_mg_bc(i, n)%bc_value = 0.0_dp
294  end select
295  end do
296  end do
297  end if
298 
299  end subroutine ard_check_params
300 
301  subroutine ard_to_conserved(ixI^L, ixO^L, w, x)
303  integer, intent(in) :: ixI^L, ixO^L
304  double precision, intent(inout) :: w(ixI^S, nw)
305  double precision, intent(in) :: x(ixI^S, 1:^ND)
306 
307  ! Do nothing (primitive and conservative are equal for ard module)
308  end subroutine ard_to_conserved
309 
310  subroutine ard_to_primitive(ixI^L, ixO^L, w, x)
312  integer, intent(in) :: ixI^L, ixO^L
313  double precision, intent(inout) :: w(ixI^S, nw)
314  double precision, intent(in) :: x(ixI^S, 1:^ND)
315 
316  ! Do nothing (primitive and conservative are equal for ard module)
317  end subroutine ard_to_primitive
318 
319  subroutine ard_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
321  integer, intent(in) :: ixI^L, ixO^L, idim
322  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
323  double precision, intent(inout) :: cmax(ixI^S)
324 
325  cmax(ixo^s) = abs(a1(idim) * w(ixo^s,u_)**(adv_pow-1))
326  if (number_of_species >= 2) then
327  cmax(ixo^s) = max(cmax(ixo^s), abs(a2(idim) * w(ixo^s,v_)**(adv_pow-1)))
328  end if
329  if (number_of_species >= 3) then
330  cmax(ixo^s) = max(cmax(ixo^s), abs(a3(idim) * w(ixo^s,w_)**(adv_pow-1)))
331  end if
332 
333  end subroutine ard_get_cmax
334 
335  subroutine ard_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
337  use mod_variables
338  integer, intent(in) :: ixI^L, ixO^L, idim
339  double precision, intent(in) :: wLC(ixI^S, nw), wRC(ixI^S,nw)
340  double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S,nw)
341  double precision, intent(in) :: x(ixI^S, 1:^ND)
342  double precision, intent(in) :: Hspeed(ixI^S,1:number_species)
343  double precision, intent(inout) :: cmax(ixI^S,1:number_species)
344  double precision, intent(inout), optional :: cmin(ixI^S,1:number_species)
345 
346  double precision :: wmean(ixI^S,nw)
347 
348  ! Since the advection coefficient can depend on unknowns,
349  ! some average over the left and right state should be taken
350  wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
351 
352  if (present(cmin)) then
353  cmin(ixo^s,1) = min(a1(idim) * wmean(ixo^s,u_)**(adv_pow-1), zero)
354  cmax(ixo^s,1) = max(a1(idim) * wmean(ixo^s,u_)**(adv_pow-1), zero)
355  if (number_of_species >= 2) then
356  cmin(ixo^s,1) = min(cmin(ixo^s,1), a2(idim) * wmean(ixo^s,v_)**(adv_pow-1))
357  cmax(ixo^s,1) = max(cmax(ixo^s,1), a2(idim) * wmean(ixo^s,v_)**(adv_pow-1))
358  end if
359  if (number_of_species >= 3) then
360  cmin(ixo^s,1) = min(cmin(ixo^s,1), a3(idim) * wmean(ixo^s,w_)**(adv_pow-1))
361  cmax(ixo^s,1) = max(cmax(ixo^s,1), a3(idim) * wmean(ixo^s,w_)**(adv_pow-1))
362  end if
363  else
364  cmax(ixo^s,1) = maxval(abs(a1(idim) * wmean(ixo^s,u_)**(adv_pow-1)))
365  if (number_of_species >=2) then
366  cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(a2(idim) * wmean(ixo^s,v_)**(adv_pow-1))))
367  end if
368  if (number_of_species >=3) then
369  cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(a3(idim) * wmean(ixo^s,w_)**(adv_pow-1))))
370  end if
371  end if
372 
373  end subroutine ard_get_cbounds
374 
375  subroutine ard_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
377  integer, intent(in) :: ixI^L, ixO^L
378  double precision, intent(in) :: dx^D, x(ixI^S, 1:^ND)
379  double precision, intent(in) :: w(ixI^S, 1:nw)
380  double precision, intent(inout) :: dtnew
381  double precision :: maxrate
382  double precision :: maxD
383  double precision :: maxA
384 
385  dtnew = bigdouble
386 
387  ! dt < dx^2 / (2 * ndim * diffusion_coeff)
388  ! use dtdiffpar < 1 for explicit and > 1 for imex/split
389  maxd = d1
390  if (number_of_species >= 2) then
391  maxd = max(maxd, d2)
392  end if
393  if (number_of_species >= 3) then
394  maxd = max(maxd, d3)
395  end if
396  dtnew = min(dtnew, dtdiffpar * minval([ dx^d ])**2 / (2 * ndim * maxd))
397 
398  ! Estimate time step for reactions
399  select case (equation_type)
400  case (eq_gray_scott)
401  maxrate = max(maxval(w(ixo^s, v_))**2 + gs_f, &
402  maxval(w(ixo^s, v_) * w(ixo^s, u_)) - gs_f - gs_k)
403  case (eq_schnakenberg)
404  maxrate = max(maxval(abs(w(ixo^s, v_) * w(ixo^s, u_) - 1)), &
405  maxval(w(ixo^s, u_))**2)
406  case (eq_brusselator)
407  maxrate = max( maxval(w(ixo^s, u_)*w(ixo^s, v_) - (br_b+1)), &
408  maxval(w(ixo^s, u_)**2) )
409  case (eq_ext_brusselator)
410  maxrate = max( maxval(w(ixo^s, u_)*w(ixo^s, v_) - (br_b+1)) + br_c, &
411  maxval(w(ixo^s, u_)**2) )
412  maxrate = max(maxrate, br_d)
413  case (eq_logistic)
414  maxrate = lg_lambda*maxval(abs(1 - w(ixo^s, u_))) ! abs for safety, normally u < 1
415  case (eq_analyt_hunds)
416  maxrate = maxval(w(ixo^s, u_)*abs(1 - w(ixo^s, u_))) / d1
417  case (eq_belousov_fn)
418  maxrate = max(&
419  maxval(abs(1.0d0 - w(ixo^s, w_) - w(ixo^s, u_))) / bzfn_epsilon, &
420  maxval(bzfn_lambda + w(ixo^s, u_)) / bzfn_delta &
421  )
422  case (eq_lorenz)
423  ! det(J) = sigma(b(r-1) + x*(x*+y*))
424  maxrate = max(lor_sigma, 1.0d0, lor_b)
425  case (eq_no_reac)
426  ! No reaction term, so no influence on timestep
427  maxrate = zero
428  case default
429  maxrate = one
430  call mpistop("Unknown equation type")
431  end select
432 
433  dtnew = min(dtnew, dtreacpar / maxrate)
434 
435  end subroutine ard_get_dt
436 
437  ! Add the flux from the advection term
438  subroutine ard_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
440  integer, intent(in) :: ixI^L, ixO^L, idim
441  double precision, intent(in) :: wC(ixI^S, 1:nw)
442  double precision, intent(in) :: w(ixI^S, 1:nw)
443  double precision, intent(in) :: x(ixI^S, 1:^ND)
444  double precision, intent(out) :: f(ixI^S, nwflux)
445 
446  f(ixo^s, u_) = (a1(idim)/adv_pow) * w(ixo^s,u_)**adv_pow
447  if (number_of_species >=2) then
448  f(ixo^s, v_) = (a2(idim)/adv_pow) * w(ixo^s,v_)**adv_pow
449  end if
450  if (number_of_species >=3) then
451  f(ixo^s, w_) = (a3(idim)/adv_pow) * w(ixo^s,w_)**adv_pow
452  end if
453 
454  end subroutine ard_get_flux
455 
456  subroutine ard_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, w, x)
457 
458  ! Add geometrical source terms
459  ! There are no geometrical source terms
460 
462 
463  integer, intent(in) :: ixI^L, ixO^L
464  double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
465  double precision, intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
466 
467  end subroutine ard_add_source_geom
468 
469  ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
470  subroutine ard_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
472  integer, intent(in) :: ixI^L, ixO^L
473  double precision, intent(in) :: qdt, dtfactor
474  double precision, intent(in) :: wCT(ixI^S, 1:nw),wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
475  double precision, intent(inout) :: w(ixI^S, 1:nw)
476  double precision :: lpl_u(ixO^S), lpl_v(ixO^S), lpl_w(ixO^S)
477  logical, intent(in) :: qsourcesplit
478  logical, intent(inout) :: active
479 
480  ! here we add the reaction terms (always) and the diffusion if no imex is used
481  if (qsourcesplit .eqv. ard_source_split) then
482  if (.not.use_imex_scheme) then
483  call ard_laplacian(ixi^l, ixo^l, wct(ixi^s, u_), lpl_u)
484  if (number_of_species >= 2) then
485  call ard_laplacian(ixi^l, ixo^l, wct(ixi^s, v_), lpl_v)
486  end if
487  if (number_of_species >= 3) then
488  call ard_laplacian(ixi^l, ixo^l, wct(ixi^s, w_), lpl_w)
489  end if
490  else
491  ! for all IMEX scheme variants: only add the reactions
492  lpl_u = 0.0d0
493  lpl_v = 0.0d0
494  lpl_w = 0.0d0
495  end if
496 
497  select case (equation_type)
498  case (eq_gray_scott)
499  w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u - &
500  wct(ixo^s, u_) * wct(ixo^s, v_)**2 + &
501  gs_f * (1 - wct(ixo^s, u_)))
502  w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v + &
503  wct(ixo^s, u_) * wct(ixo^s, v_)**2 - &
504  (gs_f + gs_k) * wct(ixo^s, v_))
505  case (eq_schnakenberg)
506  w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
507  + sb_kappa * (sb_alpha - wct(ixo^s, u_) + &
508  wct(ixo^s, u_)**2 * wct(ixo^s, v_)))
509  w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
510  + sb_kappa * (sb_beta - wct(ixo^s, u_)**2 * wct(ixo^s, v_)))
511  case (eq_brusselator)
512  w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
513  + br_a - (br_b + 1) * wct(ixo^s, u_) &
514  + wct(ixo^s, u_)**2 * wct(ixo^s, v_))
515  w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
516  + br_b * wct(ixo^s, u_) - wct(ixo^s, u_)**2 * wct(ixo^s, v_))
517  case (eq_ext_brusselator)
518  w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
519  + br_a - (br_b + 1) * wct(ixo^s, u_) &
520  + wct(ixo^s, u_)**2 * wct(ixo^s, v_) &
521  - br_c * wct(ixo^s, u_) + br_d * w(ixo^s, w_))
522  w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
523  + br_b * wct(ixo^s, u_) &
524  - wct(ixo^s, u_)**2 * wct(ixo^s, v_))
525  w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
526  + br_c * wct(ixo^s, u_) - br_d * w(ixo^s, w_))
527  case (eq_logistic)
528  w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
529  + lg_lambda * w(ixo^s, u_) * (1 - w(ixo^s, u_)))
530  case (eq_analyt_hunds)
531  w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
532  + 1.0d0/d1 * w(ixo^s, u_)**2 * (1 - w(ixo^s, u_)))
533  case (eq_belousov_fn)
534  w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
535  + 1.0d0/bzfn_epsilon * (bzfn_lambda * wct(ixo^s, u_) &
536  - wct(ixo^s, u_)*wct(ixo^s, w_) + wct(ixo^s, u_) &
537  - wct(ixo^s, u_)**2))
538  w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
539  + wct(ixo^s, u_) - wct(ixo^s, v_))
540  w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
541  + 1.0d0/bzfn_delta * (-bzfn_lambda * wct(ixo^s, w_) &
542  - wct(ixo^s, u_)*wct(ixo^s, w_) + bzfn_mu * wct(ixo^s, v_)))
543  case (eq_lorenz)
544  ! xdot = sigma.(y-x)
545  w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
546  + lor_sigma * (wct(ixo^s, v_) - wct(ixo^s, u_)))
547  ! ydot = r.x - y - x.z
548  w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
549  + lor_r * wct(ixo^s, u_) - wct(ixo^s, v_) &
550  - wct(ixo^s, u_)*wct(ixo^s, w_))
551  ! zdot = x.y - b.z
552  w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
553  + wct(ixo^s, u_)*wct(ixo^s, v_) - lor_b * wct(ixo^s, w_))
554  case (eq_no_reac)
555  w(ixo^s, u_) = w(ixo^s, u_) + qdt * d1 * lpl_u
556  case default
557  call mpistop("Unknown equation type")
558  end select
559 
560  ! enforce getbc call after source addition
561  active = .true.
562  end if
563 
564  end subroutine ard_add_source
565 
566  !> Compute the Laplacian using a standard second order scheme. For now this
567  !> method only works in slab geometries. Requires one ghost cell only.
568  subroutine ard_laplacian(ixI^L,ixO^L,var,lpl)
570  integer, intent(in) :: ixI^L, ixO^L
571  double precision, intent(in) :: var(ixI^S)
572  double precision, intent(out) :: lpl(ixO^S)
573  integer :: idir, jxO^L, hxO^L
574  double precision :: h_inv2
575 
576  if (slab) then
577  lpl(ixo^s) = 0.0d0
578  do idir = 1, ndim
579  hxo^l=ixo^l-kr(idir,^d);
580  jxo^l=ixo^l+kr(idir,^d);
581  h_inv2 = 1/dxlevel(idir)**2
582  lpl(ixo^s) = lpl(ixo^s) + h_inv2 * &
583  (var(jxo^s) - 2 * var(ixo^s) + var(hxo^s))
584  end do
585  else
586  call mpistop("ard_laplacian not implemented in this geometry")
587  end if
588 
589  end subroutine ard_laplacian
590 
591  subroutine put_laplacians_onegrid(ixI^L,ixO^L,w)
593  integer, intent(in) :: ixI^L, ixO^L
594  double precision, intent(inout) :: w(ixI^S, 1:nw)
595 
596  double precision :: lpl_u(ixO^S), lpl_v(ixO^S), lpl_w(ixO^S)
597 
598  call ard_laplacian(ixi^l, ixo^l, w(ixi^s, u_), lpl_u)
599  if (number_of_species >= 2) then
600  call ard_laplacian(ixi^l, ixo^l, w(ixi^s, v_), lpl_v)
601  end if
602  if (number_of_species >= 3) then
603  call ard_laplacian(ixi^l, ixo^l, w(ixi^s, w_), lpl_w)
604  end if
605 
606  w(ixo^s,u_) = d1*lpl_u
607  if (number_of_species >= 2) then
608  w(ixo^s,v_) = d2*lpl_v
609  end if
610  if (number_of_species >= 3) then
611  w(ixo^s,w_) = d3*lpl_w
612  end if
613 
614  end subroutine put_laplacians_onegrid
615 
616  !> inplace update of psa==>F_im(psa)
617  subroutine ard_evaluate_implicit(qtC,psa)
619  type(state), target :: psa(max_blocks)
620  double precision, intent(in) :: qtC
621 
622  integer :: iigrid, igrid, level
623  integer :: ixO^L
624 
625  !ixO^L=ixG^LL^LSUB1;
626  ixo^l=ixm^ll;
627  !$OMP PARALLEL DO PRIVATE(igrid)
628  do iigrid=1,igridstail; igrid=igrids(iigrid);
629  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
630  call put_laplacians_onegrid(ixg^ll,ixo^l,psa(igrid)%w)
631  end do
632  !$OMP END PARALLEL DO
633 
634  end subroutine ard_evaluate_implicit
635 
636  !> Implicit solve of psa=psb+dtfactor*dt*F_im(psa)
637  subroutine ard_implicit_update(dtfactor,qdt,qtC,psa,psb)
639  use mod_forest
640 
641  type(state), target :: psa(max_blocks)
642  type(state), target :: psb(max_blocks)
643  double precision, intent(in) :: qdt
644  double precision, intent(in) :: qtC
645  double precision, intent(in) :: dtfactor
646 
647  integer :: n
648  double precision :: res, max_residual, lambda
649 
650  integer :: iw_to,iw_from
651  integer :: iigrid, igrid, id
652  integer :: nc, lvl
653  type(tree_node), pointer :: pnode
654  real(dp) :: fac
655 
656  ! Avoid setting a very restrictive limit to the residual when the time step
657  ! is small (as the operator is ~ 1/(D * qdt))
658  if (qdt < dtmin) then
659  if(mype==0)then
660  print *,'skipping implicit solve: dt too small!'
661  print *,'Currently at time=',global_time,' time step=',qdt,' dtmin=',dtmin
662  endif
663  return
664  endif
665  max_residual = 1d-7/qdt
666 
667  mg%operator_type = mg_helmholtz
668  call mg_set_methods(mg)
669 
670  if (.not. mg%is_allocated) call mpistop("multigrid tree not allocated yet")
671 
672  ! First handle the u variable ***************************************
673  lambda = 1/(dtfactor * qdt * d1)
674  call helmholtz_set_lambda(lambda)
675  mg%bc(:, mg_iphi) = ard_mg_bc(1, :)
676 
677  call mg_copy_to_tree(u_, mg_irhs, factor=-lambda, state_from=psb)
678  call mg_copy_to_tree(u_, mg_iphi, state_from=psb)
679 
680  call mg_fas_fmg(mg, .true., max_res=res)
681  do n = 1, 10
682  call mg_fas_vcycle(mg, max_res=res)
683  if (res < max_residual) exit
684  end do
685 
686  call mg_copy_from_tree_gc(mg_iphi, u_, state_to=psa)
687  ! Done with the u variable ***************************************
688 
689  ! Next handle the v variable ***************************************
690  if (number_of_species >= 2) then
691  lambda = 1/(dtfactor * qdt * d2)
692  call helmholtz_set_lambda(lambda)
693  mg%bc(:, mg_iphi) = ard_mg_bc(2, :)
694 
695  call mg_copy_to_tree(v_, mg_irhs, factor=-lambda, state_from=psb)
696  call mg_copy_to_tree(v_, mg_iphi, state_from=psb)
697 
698  call mg_fas_fmg(mg, .true., max_res=res)
699  do n = 1, 10
700  call mg_fas_vcycle(mg, max_res=res)
701  if (res < max_residual) exit
702  end do
703 
704  call mg_copy_from_tree_gc(mg_iphi, v_, state_to=psa)
705  end if
706  ! Done with the v variable ***************************************
707 
708  ! Next handle the w variable ***************************************
709  if (number_of_species >= 3) then
710  lambda = 1/(dtfactor * qdt * d3)
711  call helmholtz_set_lambda(lambda)
712 
713  call mg_copy_to_tree(w_, mg_irhs, factor=-lambda, state_from=psb)
714  call mg_copy_to_tree(w_, mg_iphi, state_from=psb)
715 
716  call mg_fas_fmg(mg, .true., max_res=res)
717  do n = 1, 10
718  call mg_fas_vcycle(mg, max_res=res)
719  if (res < max_residual) exit
720  end do
721 
722  call mg_copy_from_tree_gc(mg_iphi, w_, state_to=psa)
723  end if
724  ! Done with the w variable ***************************************
725 
726  end subroutine ard_implicit_update
727 
728 end module mod_ard_phys
Module containing the physics routines for advection-reaction-diffusion equations.
Definition: mod_ard_phys.t:20
subroutine ard_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, w, x)
Definition: mod_ard_phys.t:457
subroutine ard_get_flux(wC, w, x, ixIL, ixOL, idim, f)
Definition: mod_ard_phys.t:439
double precision, public, protected br_a
Parameters for Brusselator model.
Definition: mod_ard_phys.t:80
double precision, dimension(^nd), public, protected a3
Advection coefficients for third species (w) (if applicable)
Definition: mod_ard_phys.t:67
double precision, public, protected lg_lambda
Parameter for logistic model (Fisher / KPP equation)
Definition: mod_ard_phys.t:86
double precision, public, protected gs_k
Kill rate for Gray-Scott model.
Definition: mod_ard_phys.t:77
integer, public, protected v_
For 2 or more equations.
Definition: mod_ard_phys.t:29
subroutine ard_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
Definition: mod_ard_phys.t:376
double precision, public, protected sb_alpha
Parameters for Schnakenberg model.
Definition: mod_ard_phys.t:70
double precision, public, protected gs_f
Feed rate for Gray-Scott model.
Definition: mod_ard_phys.t:75
double precision, public, protected sb_beta
Definition: mod_ard_phys.t:71
double precision, public, protected br_b
Definition: mod_ard_phys.t:81
subroutine ard_get_cmax(w, x, ixIL, ixOL, idim, cmax)
Definition: mod_ard_phys.t:320
double precision, public, protected br_c
Definition: mod_ard_phys.t:82
subroutine, public ard_phys_init()
Definition: mod_ard_phys.t:177
subroutine ard_get_cbounds(wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
Definition: mod_ard_phys.t:336
double precision, public, protected sb_kappa
Definition: mod_ard_phys.t:72
character(len=20), public, protected equation_name
Name of the system to be solved.
Definition: mod_ard_phys.t:39
integer, public, protected u_
Indices of the unknowns.
Definition: mod_ard_phys.t:28
double precision, public, protected lor_b
Parameter for Lorenz system (aspect ratio of the convection rolls)
Definition: mod_ard_phys.t:99
double precision, public, protected dtreacpar
Parameter with which to multiply the reaction timestep restriction.
Definition: mod_ard_phys.t:36
double precision, public, protected bzfn_mu
Definition: mod_ard_phys.t:92
double precision, public, protected lor_sigma
Parameter for Lorenz system (Prandtl number)
Definition: mod_ard_phys.t:97
integer, public, protected adv_pow
Power of the unknown in the advection term (1 for linear)
Definition: mod_ard_phys.t:60
subroutine ard_laplacian(ixIL, ixOL, var, lpl)
Compute the Laplacian using a standard second order scheme. For now this method only works in slab ge...
Definition: mod_ard_phys.t:569
double precision, public, protected d1
Diffusion coefficient for first species (u)
Definition: mod_ard_phys.t:53
subroutine ard_add_source(qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
Definition: mod_ard_phys.t:471
double precision, public, protected d2
Diffusion coefficient for second species (v) (if applicable)
Definition: mod_ard_phys.t:55
type(mg_bc_t), dimension(3, mg_num_neighbors), public ard_mg_bc
Boundary condition information for the multigrid method.
Definition: mod_ard_phys.t:105
double precision, public, protected bzfn_epsilon
Parameters for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
Definition: mod_ard_phys.t:89
logical, public, protected ard_particles
Whether particles module is added.
Definition: mod_ard_phys.t:33
double precision, public, protected lor_r
Parameter for Lorenz system (Rayleigh number)
Definition: mod_ard_phys.t:95
double precision, dimension(^nd), public, protected a2
Advection coefficients for second species (v) (if applicable)
Definition: mod_ard_phys.t:65
subroutine put_laplacians_onegrid(ixIL, ixOL, w)
Definition: mod_ard_phys.t:592
double precision, public, protected d3
Diffusion coefficient for third species (w) (if applicable)
Definition: mod_ard_phys.t:57
double precision, public, protected bzfn_lambda
Definition: mod_ard_phys.t:91
subroutine ard_write_info(fh)
Write this modules parameters to a snapshot.
Definition: mod_ard_phys.t:166
integer, public, protected w_
For 3 or more equations.
Definition: mod_ard_phys.t:30
subroutine ard_implicit_update(dtfactor, qdt, qtC, psa, psb)
Implicit solve of psa=psb+dtfactor*dt*F_im(psa)
Definition: mod_ard_phys.t:638
subroutine ard_check_params
Definition: mod_ard_phys.t:237
double precision, dimension(^nd), public, protected a1
Advection coefficients for first species (u)
Definition: mod_ard_phys.t:63
subroutine ard_to_primitive(ixIL, ixOL, w, x)
Definition: mod_ard_phys.t:311
double precision, public, protected bzfn_delta
Definition: mod_ard_phys.t:90
subroutine ard_evaluate_implicit(qtC, psa)
inplace update of psa==>F_im(psa)
Definition: mod_ard_phys.t:618
subroutine ard_to_conserved(ixIL, ixOL, w, x)
Definition: mod_ard_phys.t:302
double precision, public, protected br_d
Definition: mod_ard_phys.t:83
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
Module with basic grid data structures.
Definition: mod_forest.t:2
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision dtdiffpar
For resistive MHD, the time step is also limited by the diffusion time: .
integer, parameter unitpar
file handle for IO
integer, parameter bc_asymm
double precision global_time
The global simulation time.
logical use_imex_scheme
whether IMEX in use or not
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, dimension(:, :), allocatable typeboundary
Array indicating the type of boundary condition per variable and per physical boundary.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical use_particles
Use particles module or not.
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
logical slab
Cartesian geometry or not.
integer, parameter bc_periodic
integer, parameter bc_special
boundary condition types
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, parameter bc_cont
integer, parameter bc_symm
logical use_multigrid
Use multigrid (only available in 2D and 3D)
double precision dtmin
Stop the simulation when the time step becomes smaller than this value.
double precision, dimension(ndim) dxlevel
Module to couple the octree-mg library to AMRVAC. This file uses the VACPP preprocessor,...
Module containing all the particle routines.
Definition: mod_particles.t:2
subroutine particles_init()
Initialize particle data and parameters.
Definition: mod_particles.t:15
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
procedure(sub_convert), pointer phys_to_primitive
Definition: mod_physics.t:59
procedure(sub_write_info), pointer phys_write_info
Definition: mod_physics.t:79
procedure(sub_get_flux), pointer phys_get_flux
Definition: mod_physics.t:67
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
Definition: mod_physics.t:85
logical phys_req_diagonal
Whether the physics routines require diagonal ghost cells, for example for computing a curl.
Definition: mod_physics.t:24
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition: mod_physics.t:66
procedure(sub_get_dt), pointer phys_get_dt
Definition: mod_physics.t:69
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition: mod_physics.t:70
procedure(sub_check_params), pointer phys_check_params
Definition: mod_physics.t:56
integer, parameter flux_default
Indicates a normal flux.
Definition: mod_physics.t:46
integer, dimension(:, :), allocatable flux_type
Array per direction per variable, which can be used to specify that certain fluxes have to be treated...
Definition: mod_physics.t:43
procedure(sub_convert), pointer phys_to_conserved
Definition: mod_physics.t:58
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition: mod_physics.t:16
procedure(sub_implicit_update), pointer phys_implicit_update
Definition: mod_physics.t:84
procedure(sub_add_source), pointer phys_add_source
Definition: mod_physics.t:71
procedure(sub_get_cmax), pointer phys_get_cmax
Definition: mod_physics.t:61
logical phys_energy
Solve energy equation or not.
Definition: mod_physics.t:27
integer nwflux
Number of flux variables.
Definition: mod_variables.t:8
The data structure that contains information about a tree node/grid block.
Definition: mod_forest.t:11