MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
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
110contains
111
112 !> Read this module's parameters from a file
113 subroutine ard_params_read(files)
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)
126111 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
181
182 call ard_params_read(par_files)
183
184 physics_type = "ard"
185 phys_energy = .false.
187
188 allocate(start_indices(number_species),stop_indices(number_species))
189 ! set the index of the first flux variable for species 1
190 start_indices(1)=1
191 ! Use the first variable as a density
192 u_ = var_set_fluxvar("u", "u")
193 if (number_of_species >= 2) then
194 v_ = var_set_fluxvar("v", "v")
195 end if
196 if (number_of_species >= 3) then
197 w_ = var_set_fluxvar("w", "w")
198 end if
199
200 ! set number of variables which need update ghostcells
201 nwgc=nwflux
202 ! set the index of the last flux variable for species 1
203 stop_indices(1)=nwflux
204
205 ! Check whether custom flux types have been defined
206 if (.not. allocated(flux_type)) then
207 allocate(flux_type(ndir, nw))
209 else if (any(shape(flux_type) /= [ndir, nw])) then
210 call mpistop("phys_check error: flux_type has wrong shape")
211 end if
212
213 phys_get_cmax => ard_get_cmax
214 phys_get_cbounds => ard_get_cbounds
215 phys_get_flux => ard_get_flux
216 phys_to_conserved => ard_to_conserved
217 phys_to_primitive => ard_to_primitive
218 phys_add_source_geom => ard_add_source_geom
219 phys_add_source => ard_add_source
220 phys_get_dt => ard_get_dt
221 phys_write_info => ard_write_info
222 phys_check_params => ard_check_params
223 phys_implicit_update => ard_implicit_update
224 phys_evaluate_implicit => ard_evaluate_implicit
225
226 ! Initialize particles module
227 if (ard_particles) then
228 call particles_init()
229 end if
230
231 end subroutine ard_phys_init
232
233 subroutine ard_check_params
235 integer :: n, i, iw, species_list(number_of_species)
236
237 if (use_imex_scheme) then
238 use_multigrid = .true.
239 select case(number_of_species)
240 case(1)
241 species_list = [u_]
242 if (d1 == 0.0d0) then
243 write(*, *) "Diffusion coefficient cannot be zero in IMEX scheme"
244 call mpistop("Zero diffusion in IMEX scheme")
245 end if
246 case(2)
247 species_list = [u_, v_]
248 if ((d1 == 0.0d0) .or. (d2 == 0.0d0)) then
249 write(*, *) "Diffusion coefficient cannot be zero in IMEX scheme"
250 call mpistop("Zero diffusion in IMEX scheme")
251 end if
252 case(3)
253 species_list = [u_, v_, w_]
254 if ((d1 == 0.0d0) .or. (d2 == 0.0d0) .or. (d3 == 0.0d0)) then
255 write(*, *) "Diffusion coefficient cannot be zero in IMEX scheme"
256 call mpistop("Zero diffusion in IMEX scheme")
257 end if
258 end select
259
260 do i = 1, number_of_species
261 iw = species_list(i)
262
263 ! Set boundary conditions for the multigrid solver
264 do n = 1, 2*ndim
265 select case (typeboundary(iw, n))
266 case (bc_symm)
267 ! d/dx u = 0
268 ard_mg_bc(i, n)%bc_type = mg_bc_neumann
269 ard_mg_bc(i, n)%bc_value = 0.0_dp
270 case (bc_asymm)
271 ! u = 0
272 ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
273 ard_mg_bc(i, n)%bc_value = 0.0_dp
274 case (bc_cont)
275 ! d/dx u = 0
276 ard_mg_bc(i, n)%bc_type = mg_bc_neumann
277 ard_mg_bc(i, n)%bc_value = 0.0_dp
278 case (bc_periodic)
279 ! Nothing to do here
280 case (bc_special)
281 if (.not. associated(ard_mg_bc(i, n)%boundary_cond)) then
282 write(*, "(A,I0,A,I0,A)") "typeboundary(", iw, ",", n, &
283 ") is 'special', but the corresponding method " // &
284 "ard_mg_bc(i, n)%boundary_cond is not set"
285 call mpistop("ard_mg_bc(i, n)%boundary_cond not set")
286 end if
287 case default
288 write(*,*) "ard_check_params warning: unknown boundary type"
289 ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
290 ard_mg_bc(i, n)%bc_value = 0.0_dp
291 end select
292 end do
293 end do
294 end if
295
296 end subroutine ard_check_params
297
298 subroutine ard_to_conserved(ixI^L, ixO^L, w, x)
300 integer, intent(in) :: ixi^l, ixo^l
301 double precision, intent(inout) :: w(ixi^s, nw)
302 double precision, intent(in) :: x(ixi^s, 1:^nd)
303
304 ! Do nothing (primitive and conservative are equal for ard module)
305 end subroutine ard_to_conserved
306
307 subroutine ard_to_primitive(ixI^L, ixO^L, w, x)
309 integer, intent(in) :: ixi^l, ixo^l
310 double precision, intent(inout) :: w(ixi^s, nw)
311 double precision, intent(in) :: x(ixi^s, 1:^nd)
312
313 ! Do nothing (primitive and conservative are equal for ard module)
314 end subroutine ard_to_primitive
315
316 subroutine ard_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
318 integer, intent(in) :: ixi^l, ixo^l, idim
319 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
320 double precision, intent(inout) :: cmax(ixi^s)
321
322 cmax(ixo^s) = abs(a1(idim) * w(ixo^s,u_)**(adv_pow-1))
323 if (number_of_species >= 2) then
324 cmax(ixo^s) = max(cmax(ixo^s), abs(a2(idim) * w(ixo^s,v_)**(adv_pow-1)))
325 end if
326 if (number_of_species >= 3) then
327 cmax(ixo^s) = max(cmax(ixo^s), abs(a3(idim) * w(ixo^s,w_)**(adv_pow-1)))
328 end if
329
330 end subroutine ard_get_cmax
331
332 subroutine ard_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
334 use mod_variables
335 integer, intent(in) :: ixi^l, ixo^l, idim
336 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s,nw)
337 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s,nw)
338 double precision, intent(in) :: x(ixi^s, 1:^nd)
339 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
340 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
341 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
342
343 double precision :: wmean(ixi^s,nw)
344
345 ! Since the advection coefficient can depend on unknowns,
346 ! some average over the left and right state should be taken
347 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
348
349 if (present(cmin)) then
350 cmin(ixo^s,1) = min(a1(idim) * wmean(ixo^s,u_)**(adv_pow-1), zero)
351 cmax(ixo^s,1) = max(a1(idim) * wmean(ixo^s,u_)**(adv_pow-1), zero)
352 if (number_of_species >= 2) then
353 cmin(ixo^s,1) = min(cmin(ixo^s,1), a2(idim) * wmean(ixo^s,v_)**(adv_pow-1))
354 cmax(ixo^s,1) = max(cmax(ixo^s,1), a2(idim) * wmean(ixo^s,v_)**(adv_pow-1))
355 end if
356 if (number_of_species >= 3) then
357 cmin(ixo^s,1) = min(cmin(ixo^s,1), a3(idim) * wmean(ixo^s,w_)**(adv_pow-1))
358 cmax(ixo^s,1) = max(cmax(ixo^s,1), a3(idim) * wmean(ixo^s,w_)**(adv_pow-1))
359 end if
360 else
361 cmax(ixo^s,1) = maxval(abs(a1(idim) * wmean(ixo^s,u_)**(adv_pow-1)))
362 if (number_of_species >=2) then
363 cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(a2(idim) * wmean(ixo^s,v_)**(adv_pow-1))))
364 end if
365 if (number_of_species >=3) then
366 cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(a3(idim) * wmean(ixo^s,w_)**(adv_pow-1))))
367 end if
368 end if
369
370 end subroutine ard_get_cbounds
371
372 subroutine ard_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
374 integer, intent(in) :: ixi^l, ixo^l
375 double precision, intent(in) :: dx^d, x(ixi^s, 1:^nd)
376 double precision, intent(in) :: w(ixi^s, 1:nw)
377 double precision, intent(inout) :: dtnew
378 double precision :: maxrate
379 double precision :: maxd
380 double precision :: maxa
381
382 dtnew = bigdouble
383
384 ! dt < dx^2 / (2 * ndim * diffusion_coeff)
385 ! use dtdiffpar < 1 for explicit and > 1 for imex/split
386 maxd = d1
387 if (number_of_species >= 2) then
388 maxd = max(maxd, d2)
389 end if
390 if (number_of_species >= 3) then
391 maxd = max(maxd, d3)
392 end if
393 dtnew = min(dtnew, dtdiffpar * minval([ dx^d ])**2 / (2 * ndim * maxd))
394
395 ! Estimate time step for reactions
396 select case (equation_type)
397 case (eq_gray_scott)
398 maxrate = max(maxval(w(ixo^s, v_))**2 + gs_f, &
399 maxval(w(ixo^s, v_) * w(ixo^s, u_)) - gs_f - gs_k)
400 case (eq_schnakenberg)
401 maxrate = max(maxval(abs(w(ixo^s, v_) * w(ixo^s, u_) - 1)), &
402 maxval(w(ixo^s, u_))**2)
403 case (eq_brusselator)
404 maxrate = max( maxval(w(ixo^s, u_)*w(ixo^s, v_) - (br_b+1)), &
405 maxval(w(ixo^s, u_)**2) )
406 case (eq_ext_brusselator)
407 maxrate = max( maxval(w(ixo^s, u_)*w(ixo^s, v_) - (br_b+1)) + br_c, &
408 maxval(w(ixo^s, u_)**2) )
409 maxrate = max(maxrate, br_d)
410 case (eq_logistic)
411 maxrate = lg_lambda*maxval(abs(1 - w(ixo^s, u_))) ! abs for safety, normally u < 1
412 case (eq_analyt_hunds)
413 maxrate = maxval(w(ixo^s, u_)*abs(1 - w(ixo^s, u_))) / d1
414 case (eq_belousov_fn)
415 maxrate = max(&
416 maxval(abs(1.0d0 - w(ixo^s, w_) - w(ixo^s, u_))) / bzfn_epsilon, &
417 maxval(bzfn_lambda + w(ixo^s, u_)) / bzfn_delta &
418 )
419 case (eq_lorenz)
420 ! det(J) = sigma(b(r-1) + x*(x*+y*))
421 maxrate = max(lor_sigma, 1.0d0, lor_b)
422 case (eq_no_reac)
423 ! No reaction term, so no influence on timestep
424 maxrate = zero
425 case default
426 maxrate = one
427 call mpistop("Unknown equation type")
428 end select
429
430 dtnew = min(dtnew, dtreacpar / maxrate)
431
432 end subroutine ard_get_dt
433
434 ! Add the flux from the advection term
435 subroutine ard_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
437 integer, intent(in) :: ixi^l, ixo^l, idim
438 double precision, intent(in) :: wc(ixi^s, 1:nw)
439 double precision, intent(in) :: w(ixi^s, 1:nw)
440 double precision, intent(in) :: x(ixi^s, 1:^nd)
441 double precision, intent(out) :: f(ixi^s, nwflux)
442
443 f(ixo^s, u_) = (a1(idim)/adv_pow) * w(ixo^s,u_)**adv_pow
444 if (number_of_species >=2) then
445 f(ixo^s, v_) = (a2(idim)/adv_pow) * w(ixo^s,v_)**adv_pow
446 end if
447 if (number_of_species >=3) then
448 f(ixo^s, w_) = (a3(idim)/adv_pow) * w(ixo^s,w_)**adv_pow
449 end if
450
451 end subroutine ard_get_flux
452
453 subroutine ard_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim,w, x)
454
455 ! Add geometrical source terms
456 ! There are no geometrical source terms
457
459
460 integer, intent(in) :: ixi^l, ixo^l
461 double precision, intent(in) :: qdt, dtfactor, x(ixi^s, 1:^nd)
462 double precision, intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
463
464 end subroutine ard_add_source_geom
465
466 ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
467 subroutine ard_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
469 integer, intent(in) :: ixi^l, ixo^l
470 double precision, intent(in) :: qdt, dtfactor
471 double precision, intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:ndim)
472 double precision, intent(inout) :: w(ixi^s, 1:nw)
473 double precision :: lpl_u(ixo^s), lpl_v(ixo^s), lpl_w(ixo^s)
474 logical, intent(in) :: qsourcesplit
475 logical, intent(inout) :: active
476
477 ! here we add the reaction terms (always) and the diffusion if no imex is used
478 if (qsourcesplit .eqv. ard_source_split) then
479 if (.not.use_imex_scheme) then
480 call ard_laplacian(ixi^l, ixo^l, wct(ixi^s, u_), lpl_u)
481 if (number_of_species >= 2) then
482 call ard_laplacian(ixi^l, ixo^l, wct(ixi^s, v_), lpl_v)
483 end if
484 if (number_of_species >= 3) then
485 call ard_laplacian(ixi^l, ixo^l, wct(ixi^s, w_), lpl_w)
486 end if
487 else
488 ! for all IMEX scheme variants: only add the reactions
489 lpl_u = 0.0d0
490 lpl_v = 0.0d0
491 lpl_w = 0.0d0
492 end if
493
494 select case (equation_type)
495 case (eq_gray_scott)
496 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u - &
497 wct(ixo^s, u_) * wct(ixo^s, v_)**2 + &
498 gs_f * (1 - wct(ixo^s, u_)))
499 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v + &
500 wct(ixo^s, u_) * wct(ixo^s, v_)**2 - &
501 (gs_f + gs_k) * wct(ixo^s, v_))
502 case (eq_schnakenberg)
503 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
504 + sb_kappa * (sb_alpha - wct(ixo^s, u_) + &
505 wct(ixo^s, u_)**2 * wct(ixo^s, v_)))
506 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
507 + sb_kappa * (sb_beta - wct(ixo^s, u_)**2 * wct(ixo^s, v_)))
508 case (eq_brusselator)
509 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
510 + br_a - (br_b + 1) * wct(ixo^s, u_) &
511 + wct(ixo^s, u_)**2 * wct(ixo^s, v_))
512 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
513 + br_b * wct(ixo^s, u_) - wct(ixo^s, u_)**2 * wct(ixo^s, v_))
514 case (eq_ext_brusselator)
515 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
516 + br_a - (br_b + 1) * wct(ixo^s, u_) &
517 + wct(ixo^s, u_)**2 * wct(ixo^s, v_) &
518 - br_c * wct(ixo^s, u_) + br_d * w(ixo^s, w_))
519 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
520 + br_b * wct(ixo^s, u_) &
521 - wct(ixo^s, u_)**2 * wct(ixo^s, v_))
522 w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
523 + br_c * wct(ixo^s, u_) - br_d * w(ixo^s, w_))
524 case (eq_logistic)
525 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
526 + lg_lambda * w(ixo^s, u_) * (1 - w(ixo^s, u_)))
527 case (eq_analyt_hunds)
528 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
529 + 1.0d0/d1 * w(ixo^s, u_)**2 * (1 - w(ixo^s, u_)))
530 case (eq_belousov_fn)
531 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
532 + 1.0d0/bzfn_epsilon * (bzfn_lambda * wct(ixo^s, u_) &
533 - wct(ixo^s, u_)*wct(ixo^s, w_) + wct(ixo^s, u_) &
534 - wct(ixo^s, u_)**2))
535 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
536 + wct(ixo^s, u_) - wct(ixo^s, v_))
537 w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
538 + 1.0d0/bzfn_delta * (-bzfn_lambda * wct(ixo^s, w_) &
539 - wct(ixo^s, u_)*wct(ixo^s, w_) + bzfn_mu * wct(ixo^s, v_)))
540 case (eq_lorenz)
541 ! xdot = sigma.(y-x)
542 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
543 + lor_sigma * (wct(ixo^s, v_) - wct(ixo^s, u_)))
544 ! ydot = r.x - y - x.z
545 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
546 + lor_r * wct(ixo^s, u_) - wct(ixo^s, v_) &
547 - wct(ixo^s, u_)*wct(ixo^s, w_))
548 ! zdot = x.y - b.z
549 w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
550 + wct(ixo^s, u_)*wct(ixo^s, v_) - lor_b * wct(ixo^s, w_))
551 case (eq_no_reac)
552 w(ixo^s, u_) = w(ixo^s, u_) + qdt * d1 * lpl_u
553 case default
554 call mpistop("Unknown equation type")
555 end select
556
557 ! enforce getbc call after source addition
558 active = .true.
559 end if
560
561 end subroutine ard_add_source
562
563 !> Compute the Laplacian using a standard second order scheme. For now this
564 !> method only works in slab geometries. Requires one ghost cell only.
565 subroutine ard_laplacian(ixI^L,ixO^L,var,lpl)
567 integer, intent(in) :: ixi^l, ixo^l
568 double precision, intent(in) :: var(ixi^s)
569 double precision, intent(out) :: lpl(ixo^s)
570 integer :: idir, jxo^l, hxo^l
571 double precision :: h_inv2
572
573 if (slab) then
574 lpl(ixo^s) = 0.0d0
575 do idir = 1, ndim
576 hxo^l=ixo^l-kr(idir,^d);
577 jxo^l=ixo^l+kr(idir,^d);
578 h_inv2 = 1/dxlevel(idir)**2
579 lpl(ixo^s) = lpl(ixo^s) + h_inv2 * &
580 (var(jxo^s) - 2 * var(ixo^s) + var(hxo^s))
581 end do
582 else
583 call mpistop("ard_laplacian not implemented in this geometry")
584 end if
585
586 end subroutine ard_laplacian
587
588 subroutine put_laplacians_onegrid(ixI^L,ixO^L,w)
590 integer, intent(in) :: ixi^l, ixo^l
591 double precision, intent(inout) :: w(ixi^s, 1:nw)
592
593 double precision :: lpl_u(ixo^s), lpl_v(ixo^s), lpl_w(ixo^s)
594
595 call ard_laplacian(ixi^l, ixo^l, w(ixi^s, u_), lpl_u)
596 if (number_of_species >= 2) then
597 call ard_laplacian(ixi^l, ixo^l, w(ixi^s, v_), lpl_v)
598 end if
599 if (number_of_species >= 3) then
600 call ard_laplacian(ixi^l, ixo^l, w(ixi^s, w_), lpl_w)
601 end if
602
603 w(ixo^s,u_) = d1*lpl_u
604 if (number_of_species >= 2) then
605 w(ixo^s,v_) = d2*lpl_v
606 end if
607 if (number_of_species >= 3) then
608 w(ixo^s,w_) = d3*lpl_w
609 end if
610
611 end subroutine put_laplacians_onegrid
612
613 !> inplace update of psa==>F_im(psa)
614 subroutine ard_evaluate_implicit(qtC,psa)
616 type(state), target :: psa(max_blocks)
617 double precision, intent(in) :: qtc
618
619 integer :: iigrid, igrid, level
620 integer :: ixo^l
621
622 !ixO^L=ixG^LL^LSUB1;
623 ixo^l=ixm^ll;
624 !$OMP PARALLEL DO PRIVATE(igrid)
625 do iigrid=1,igridstail; igrid=igrids(iigrid);
626 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
627 call put_laplacians_onegrid(ixg^ll,ixo^l,psa(igrid)%w)
628 end do
629 !$OMP END PARALLEL DO
630
631 end subroutine ard_evaluate_implicit
632
633 !> Implicit solve of psa=psb+dtfactor*dt*F_im(psa)
634 subroutine ard_implicit_update(dtfactor,qdt,qtC,psa,psb)
636 use mod_forest
637
638 type(state), target :: psa(max_blocks)
639 type(state), target :: psb(max_blocks)
640 double precision, intent(in) :: qdt
641 double precision, intent(in) :: qtc
642 double precision, intent(in) :: dtfactor
643
644 integer :: n
645 double precision :: res, max_residual, lambda
646
647 integer :: iw_to,iw_from
648 integer :: iigrid, igrid, id
649 integer :: nc, lvl
650 type(tree_node), pointer :: pnode
651 real(dp) :: fac
652
653 ! Avoid setting a very restrictive limit to the residual when the time step
654 ! is small (as the operator is ~ 1/(D * qdt))
655 if (qdt < dtmin) then
656 if(mype==0)then
657 print *,'skipping implicit solve: dt too small!'
658 print *,'Currently at time=',global_time,' time step=',qdt,' dtmin=',dtmin
659 endif
660 return
661 endif
662 max_residual = 1d-7/qdt
663
664 mg%operator_type = mg_helmholtz
665 call mg_set_methods(mg)
666
667 if (.not. mg%is_allocated) call mpistop("multigrid tree not allocated yet")
668
669 ! First handle the u variable ***************************************
670 lambda = 1/(dtfactor * qdt * d1)
671 call helmholtz_set_lambda(lambda)
672 mg%bc(:, mg_iphi) = ard_mg_bc(1, :)
673
674 call mg_copy_to_tree(u_, mg_irhs, factor=-lambda, state_from=psb)
675 call mg_copy_to_tree(u_, mg_iphi, state_from=psb)
676
677 call mg_fas_fmg(mg, .true., max_res=res)
678 do n = 1, 10
679 call mg_fas_vcycle(mg, max_res=res)
680 if (res < max_residual) exit
681 end do
682
683 call mg_copy_from_tree_gc(mg_iphi, u_, state_to=psa)
684 ! Done with the u variable ***************************************
685
686 ! Next handle the v variable ***************************************
687 if (number_of_species >= 2) then
688 lambda = 1/(dtfactor * qdt * d2)
689 call helmholtz_set_lambda(lambda)
690 mg%bc(:, mg_iphi) = ard_mg_bc(2, :)
691
692 call mg_copy_to_tree(v_, mg_irhs, factor=-lambda, state_from=psb)
693 call mg_copy_to_tree(v_, mg_iphi, state_from=psb)
694
695 call mg_fas_fmg(mg, .true., max_res=res)
696 do n = 1, 10
697 call mg_fas_vcycle(mg, max_res=res)
698 if (res < max_residual) exit
699 end do
700
701 call mg_copy_from_tree_gc(mg_iphi, v_, state_to=psa)
702 end if
703 ! Done with the v variable ***************************************
704
705 ! Next handle the w variable ***************************************
706 if (number_of_species >= 3) then
707 lambda = 1/(dtfactor * qdt * d3)
708 call helmholtz_set_lambda(lambda)
709
710 call mg_copy_to_tree(w_, mg_irhs, factor=-lambda, state_from=psb)
711 call mg_copy_to_tree(w_, mg_iphi, state_from=psb)
712
713 call mg_fas_fmg(mg, .true., max_res=res)
714 do n = 1, 10
715 call mg_fas_vcycle(mg, max_res=res)
716 if (res < max_residual) exit
717 end do
718
719 call mg_copy_from_tree_gc(mg_iphi, w_, state_to=psa)
720 end if
721 ! Done with the w variable ***************************************
722
723 end subroutine ard_implicit_update
724
725end module mod_ard_phys
Module containing the physics routines for advection-reaction-diffusion equations.
double precision, public, protected br_a
Parameters for Brusselator model.
double precision, dimension(^nd), public, protected a3
Advection coefficients for third species (w) (if applicable)
double precision, public, protected lg_lambda
Parameter for logistic model (Fisher / KPP equation)
double precision, public, protected gs_k
Kill rate for Gray-Scott model.
integer, public, protected v_
For 2 or more equations.
double precision, public, protected sb_alpha
Parameters for Schnakenberg model.
double precision, public, protected gs_f
Feed rate for Gray-Scott model.
double precision, public, protected sb_beta
double precision, public, protected br_b
double precision, public, protected br_c
subroutine, public ard_phys_init()
double precision, public, protected sb_kappa
character(len=20), public, protected equation_name
Name of the system to be solved.
integer, public, protected u_
Indices of the unknowns.
double precision, public, protected lor_b
Parameter for Lorenz system (aspect ratio of the convection rolls)
double precision, public, protected dtreacpar
Parameter with which to multiply the reaction timestep restriction.
double precision, public, protected bzfn_mu
double precision, public, protected lor_sigma
Parameter for Lorenz system (Prandtl number)
integer, public, protected adv_pow
Power of the unknown in the advection term (1 for linear)
double precision, public, protected d1
Diffusion coefficient for first species (u)
double precision, public, protected d2
Diffusion coefficient for second species (v) (if applicable)
type(mg_bc_t), dimension(3, mg_num_neighbors), public ard_mg_bc
Boundary condition information for the multigrid method.
double precision, public, protected bzfn_epsilon
Parameters for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
logical, public, protected ard_particles
Whether particles module is added.
double precision, public, protected lor_r
Parameter for Lorenz system (Rayleigh number)
double precision, dimension(^nd), public, protected a2
Advection coefficients for second species (v) (if applicable)
double precision, public, protected d3
Diffusion coefficient for third species (w) (if applicable)
double precision, public, protected bzfn_lambda
integer, public, protected w_
For 3 or more equations.
double precision, dimension(^nd), public, protected a1
Advection coefficients for first species (u)
double precision, public, protected bzfn_delta
double precision, public, protected br_d
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
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
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.
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
logical slab
Cartesian geometry or not.
integer, parameter bc_periodic
integer, parameter bc_special
boundary condition types
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision, dimension(:,:), allocatable dx
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
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.
integer max_blocks
The maximum number of grid blocks in a processor.
Module to couple the octree-mg library to AMRVAC. This file uses the VACPP preprocessor,...
Module containing all the particle routines.
subroutine particles_init()
Initialize particle data and parameters.
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:55
procedure(sub_write_info), pointer phys_write_info
Definition mod_physics.t:77
procedure(sub_get_flux), pointer phys_get_flux
Definition mod_physics.t:64
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
Definition mod_physics.t:83
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition mod_physics.t:63
procedure(sub_get_dt), pointer phys_get_dt
Definition mod_physics.t:67
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition mod_physics.t:68
procedure(sub_check_params), pointer phys_check_params
Definition mod_physics.t:52
integer, parameter flux_default
Indicates a normal flux.
Definition mod_physics.t:24
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:21
procedure(sub_convert), pointer phys_to_conserved
Definition mod_physics.t:54
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition mod_physics.t:50
procedure(sub_implicit_update), pointer phys_implicit_update
Definition mod_physics.t:82
procedure(sub_add_source), pointer phys_add_source
Definition mod_physics.t:69
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:57
logical phys_energy
Solve energy equation or not.
Definition mod_physics.t:35
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
integer nwflux
Number of flux variables.
The data structure that contains information about a tree node/grid block.
Definition mod_forest.t:11