MPI-AMRVAC 3.2
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
180
181 call ard_params_read(par_files)
182
183 physics_type = "ard"
184 phys_energy = .false.
186
187 allocate(start_indices(number_species),stop_indices(number_species))
188 ! set the index of the first flux variable for species 1
189 start_indices(1)=1
190 ! Use the first variable as a density
191 u_ = var_set_fluxvar("u", "u")
192 if (number_of_species >= 2) then
193 v_ = var_set_fluxvar("v", "v")
194 end if
195 if (number_of_species >= 3) then
196 w_ = var_set_fluxvar("w", "w")
197 end if
198
199 ! set number of variables which need update ghostcells
200 nwgc=nwflux
201 ! set the index of the last flux variable for species 1
202 stop_indices(1)=nwflux
203
204 ! Check whether custom flux types have been defined
205 if (.not. allocated(flux_type)) then
206 allocate(flux_type(ndir, nw))
208 else if (any(shape(flux_type) /= [ndir, nw])) then
209 call mpistop("phys_check error: flux_type has wrong shape")
210 end if
211
212 phys_get_cmax => ard_get_cmax
213 phys_get_cbounds => ard_get_cbounds
214 phys_get_flux => ard_get_flux
215 phys_to_conserved => ard_to_conserved
216 phys_to_primitive => ard_to_primitive
217 phys_add_source_geom => ard_add_source_geom
218 phys_add_source => ard_add_source
219 phys_get_dt => ard_get_dt
220 phys_write_info => ard_write_info
221 phys_check_params => ard_check_params
222 phys_implicit_update => ard_implicit_update
223 phys_evaluate_implicit => ard_evaluate_implicit
224
225
226 end subroutine ard_phys_init
227
228 subroutine ard_check_params
231 use mod_geometry, only: coordinate
233 use mod_particles, only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
234 integer :: n, i, iw, species_list(number_of_species)
235
236 ! Initialize particles module
237 if (ard_particles) then
238 call particles_init()
239 end if
240
241 if (use_imex_scheme) then
242 use_multigrid = .true.
243 select case(number_of_species)
244 case(1)
245 species_list = [u_]
246 if (d1 == 0.0d0) then
247 write(*, *) "Diffusion coefficient cannot be zero in IMEX scheme"
248 call mpistop("Zero diffusion in IMEX scheme")
249 end if
250 case(2)
251 species_list = [u_, v_]
252 if ((d1 == 0.0d0) .or. (d2 == 0.0d0)) then
253 write(*, *) "Diffusion coefficient cannot be zero in IMEX scheme"
254 call mpistop("Zero diffusion in IMEX scheme")
255 end if
256 case(3)
257 species_list = [u_, v_, w_]
258 if ((d1 == 0.0d0) .or. (d2 == 0.0d0) .or. (d3 == 0.0d0)) then
259 write(*, *) "Diffusion coefficient cannot be zero in IMEX scheme"
260 call mpistop("Zero diffusion in IMEX scheme")
261 end if
262 end select
263
264 do i = 1, number_of_species
265 iw = species_list(i)
266
267 ! Set boundary conditions for the multigrid solver
268 do n = 1, 2*ndim
269 select case (typeboundary(iw, n))
270 case (bc_symm)
271 ! d/dx u = 0
272 ard_mg_bc(i, n)%bc_type = mg_bc_neumann
273 ard_mg_bc(i, n)%bc_value = 0.0_dp
274 case (bc_asymm)
275 ! u = 0
276 ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
277 ard_mg_bc(i, n)%bc_value = 0.0_dp
278 case (bc_cont)
279 ! d/dx u = 0
280 ard_mg_bc(i, n)%bc_type = mg_bc_neumann
281 ard_mg_bc(i, n)%bc_value = 0.0_dp
282 case (bc_periodic)
283 ! Nothing to do here
284 case (bc_special)
285 if (.not. associated(usr_special_mg_bc)) then
286 call mpistop("special BC for MG not defined")
287 else
288 call usr_special_mg_bc(n)
289 endif
290 !!if (.not. associated(ard_mg_bc(i, n)%boundary_cond)) then
291 !! write(*, "(A,I0,A,I0,A)") "typeboundary(", iw, ",", n, &
292 !! ") is 'special', but the corresponding method " // &
293 !! "ard_mg_bc(i, n)%boundary_cond is not set"
294 !! call mpistop("ard_mg_bc(i, n)%boundary_cond not set")
295 !!end if
296 case default
297 write(*,*) "ard_check_params warning: unknown boundary type"
298 ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
299 ard_mg_bc(i, n)%bc_value = 0.0_dp
300 end select
301 end do
302 end do
303 end if
304
305
306 if(mype==0)then
307 write(*,*)'====ARD run with settings===================='
308 write(*,*)'Using mod_ard_phys with settings:'
309 write(*,*)'Dimensionality :',ndim
310 write(*,*)'vector components:',ndir
311 write(*,*)'coordinate set to type,slab:',coordinate,slab
312 write(*,*)'number of variables nw=',nw
313 write(*,*)' start index iwstart=',iwstart
314 write(*,*)'number of vector variables=',nvector
315 write(*,*)'number of stagger variables nws=',nws
316 write(*,*)'number of variables with BCs=',nwgc
317 write(*,*)'number of vars with fluxes=',nwflux
318 write(*,*)'number of vars with flux + BC=',nwfluxbc
319 write(*,*)'number of auxiliary variables=',nwaux
320 write(*,*)'number of extra vars without flux=',nwextra
321 write(*,*)'number of extra vars for wextra=',nw_extra
322 write(*,*)'number of auxiliary I/O variables=',nwauxio
323 write(*,*)' ard_particles=',ard_particles
324 if(ard_particles) then
325 write(*,*) '*****Using particles: npayload,ngridvars :', npayload,ngridvars
326 write(*,*) '*****Using particles: nusrpayload :', nusrpayload
327 write(*,*) '*****Using particles: num_particles :', num_particles
328 write(*,*) '*****Using particles: physics_type_particles=',physics_type_particles
329 end if
330 write(*,*)'number of ghostcells=',nghostcells
331 write(*,*)'==========================================='
332 endif
333
334
335 end subroutine ard_check_params
336
337 subroutine ard_to_conserved(ixI^L, ixO^L, w, x)
339 integer, intent(in) :: ixi^l, ixo^l
340 double precision, intent(inout) :: w(ixi^s, nw)
341 double precision, intent(in) :: x(ixi^s, 1:^nd)
342
343 ! Do nothing (primitive and conservative are equal for ard module)
344 end subroutine ard_to_conserved
345
346 subroutine ard_to_primitive(ixI^L, ixO^L, w, x)
348 integer, intent(in) :: ixi^l, ixo^l
349 double precision, intent(inout) :: w(ixi^s, nw)
350 double precision, intent(in) :: x(ixi^s, 1:^nd)
351
352 ! Do nothing (primitive and conservative are equal for ard module)
353 end subroutine ard_to_primitive
354
355 subroutine ard_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
357 integer, intent(in) :: ixi^l, ixo^l, idim
358 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
359 double precision, intent(inout) :: cmax(ixi^s)
360
361 cmax(ixo^s) = abs(a1(idim) * w(ixo^s,u_)**(adv_pow-1))
362 if (number_of_species >= 2) then
363 cmax(ixo^s) = max(cmax(ixo^s), abs(a2(idim) * w(ixo^s,v_)**(adv_pow-1)))
364 end if
365 if (number_of_species >= 3) then
366 cmax(ixo^s) = max(cmax(ixo^s), abs(a3(idim) * w(ixo^s,w_)**(adv_pow-1)))
367 end if
368
369 end subroutine ard_get_cmax
370
371 subroutine ard_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
373 use mod_variables
374 integer, intent(in) :: ixi^l, ixo^l, idim
375 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s,nw)
376 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s,nw)
377 double precision, intent(in) :: x(ixi^s, 1:^nd)
378 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
379 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
380 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
381
382 double precision :: wmean(ixi^s,nw)
383
384 ! Since the advection coefficient can depend on unknowns,
385 ! some average over the left and right state should be taken
386 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
387
388 if (present(cmin)) then
389 cmin(ixo^s,1) = min(a1(idim) * wmean(ixo^s,u_)**(adv_pow-1), zero)
390 cmax(ixo^s,1) = max(a1(idim) * wmean(ixo^s,u_)**(adv_pow-1), zero)
391 if (number_of_species >= 2) then
392 cmin(ixo^s,1) = min(cmin(ixo^s,1), a2(idim) * wmean(ixo^s,v_)**(adv_pow-1))
393 cmax(ixo^s,1) = max(cmax(ixo^s,1), a2(idim) * wmean(ixo^s,v_)**(adv_pow-1))
394 end if
395 if (number_of_species >= 3) then
396 cmin(ixo^s,1) = min(cmin(ixo^s,1), a3(idim) * wmean(ixo^s,w_)**(adv_pow-1))
397 cmax(ixo^s,1) = max(cmax(ixo^s,1), a3(idim) * wmean(ixo^s,w_)**(adv_pow-1))
398 end if
399 else
400 cmax(ixo^s,1) = maxval(abs(a1(idim) * wmean(ixo^s,u_)**(adv_pow-1)))
401 if (number_of_species >=2) then
402 cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(a2(idim) * wmean(ixo^s,v_)**(adv_pow-1))))
403 end if
404 if (number_of_species >=3) then
405 cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(a3(idim) * wmean(ixo^s,w_)**(adv_pow-1))))
406 end if
407 end if
408
409 end subroutine ard_get_cbounds
410
411 subroutine ard_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
413 integer, intent(in) :: ixi^l, ixo^l
414 double precision, intent(in) :: dx^d, x(ixi^s, 1:^nd)
415 double precision, intent(in) :: w(ixi^s, 1:nw)
416 double precision, intent(inout) :: dtnew
417 double precision :: maxrate
418 double precision :: maxd
419 double precision :: maxa
420
421 dtnew = bigdouble
422
423 ! dt < dx^2 / (2 * ndim * diffusion_coeff)
424 ! use dtdiffpar < 1 for explicit and > 1 for imex/split
425 maxd = d1
426 if (number_of_species >= 2) then
427 maxd = max(maxd, d2)
428 end if
429 if (number_of_species >= 3) then
430 maxd = max(maxd, d3)
431 end if
432 dtnew = min(dtnew, dtdiffpar * minval([ dx^d ])**2 / (2 * ndim * maxd))
433
434 ! Estimate time step for reactions
435 select case (equation_type)
436 case (eq_gray_scott)
437 maxrate = max(maxval(w(ixo^s, v_))**2 + gs_f, &
438 maxval(w(ixo^s, v_) * w(ixo^s, u_)) - gs_f - gs_k)
439 case (eq_schnakenberg)
440 maxrate = max(maxval(abs(w(ixo^s, v_) * w(ixo^s, u_) - 1)), &
441 maxval(w(ixo^s, u_))**2)
442 case (eq_brusselator)
443 maxrate = max( maxval(w(ixo^s, u_)*w(ixo^s, v_) - (br_b+1)), &
444 maxval(w(ixo^s, u_)**2) )
445 case (eq_ext_brusselator)
446 maxrate = max( maxval(w(ixo^s, u_)*w(ixo^s, v_) - (br_b+1)) + br_c, &
447 maxval(w(ixo^s, u_)**2) )
448 maxrate = max(maxrate, br_d)
449 case (eq_logistic)
450 maxrate = lg_lambda*maxval(abs(1 - w(ixo^s, u_))) ! abs for safety, normally u < 1
451 case (eq_analyt_hunds)
452 maxrate = maxval(w(ixo^s, u_)*abs(1 - w(ixo^s, u_))) / d1
453 case (eq_belousov_fn)
454 maxrate = max(&
455 maxval(abs(1.0d0 - w(ixo^s, w_) - w(ixo^s, u_))) / bzfn_epsilon, &
456 maxval(bzfn_lambda + w(ixo^s, u_)) / bzfn_delta &
457 )
458 case (eq_lorenz)
459 ! det(J) = sigma(b(r-1) + x*(x*+y*))
460 maxrate = max(lor_sigma, 1.0d0, lor_b)
461 case (eq_no_reac)
462 ! No reaction term, so no influence on timestep
463 maxrate = zero
464 case default
465 maxrate = one
466 call mpistop("Unknown equation type")
467 end select
468
469 dtnew = min(dtnew, dtreacpar / maxrate)
470
471 end subroutine ard_get_dt
472
473 ! Add the flux from the advection term
474 subroutine ard_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
476 integer, intent(in) :: ixi^l, ixo^l, idim
477 double precision, intent(in) :: wc(ixi^s, 1:nw)
478 double precision, intent(in) :: w(ixi^s, 1:nw)
479 double precision, intent(in) :: x(ixi^s, 1:^nd)
480 double precision, intent(out) :: f(ixi^s, nwflux)
481
482 f(ixo^s, u_) = (a1(idim)/adv_pow) * w(ixo^s,u_)**adv_pow
483 if (number_of_species >=2) then
484 f(ixo^s, v_) = (a2(idim)/adv_pow) * w(ixo^s,v_)**adv_pow
485 end if
486 if (number_of_species >=3) then
487 f(ixo^s, w_) = (a3(idim)/adv_pow) * w(ixo^s,w_)**adv_pow
488 end if
489
490 end subroutine ard_get_flux
491
492 subroutine ard_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim,w, x)
493
494 ! Add geometrical source terms
495 ! There are no geometrical source terms
496
498
499 integer, intent(in) :: ixi^l, ixo^l
500 double precision, intent(in) :: qdt, dtfactor, x(ixi^s, 1:^nd)
501 double precision, intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
502
503 end subroutine ard_add_source_geom
504
505 ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
506 subroutine ard_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
508 integer, intent(in) :: ixi^l, ixo^l
509 double precision, intent(in) :: qdt, dtfactor
510 double precision, intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:ndim)
511 double precision, intent(inout) :: w(ixi^s, 1:nw)
512 double precision :: lpl_u(ixo^s), lpl_v(ixo^s), lpl_w(ixo^s)
513 logical, intent(in) :: qsourcesplit
514 logical, intent(inout) :: active
515
516 ! here we add the reaction terms (always) and the diffusion if no imex is used
517 if (qsourcesplit .eqv. ard_source_split) then
518 if (.not.use_imex_scheme) then
519 call ard_laplacian(ixi^l, ixo^l, wct(ixi^s, u_), lpl_u)
520 if (number_of_species >= 2) then
521 call ard_laplacian(ixi^l, ixo^l, wct(ixi^s, v_), lpl_v)
522 end if
523 if (number_of_species >= 3) then
524 call ard_laplacian(ixi^l, ixo^l, wct(ixi^s, w_), lpl_w)
525 end if
526 else
527 ! for all IMEX scheme variants: only add the reactions
528 lpl_u = 0.0d0
529 lpl_v = 0.0d0
530 lpl_w = 0.0d0
531 end if
532
533 select case (equation_type)
534 case (eq_gray_scott)
535 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u - &
536 wct(ixo^s, u_) * wct(ixo^s, v_)**2 + &
537 gs_f * (1 - wct(ixo^s, u_)))
538 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v + &
539 wct(ixo^s, u_) * wct(ixo^s, v_)**2 - &
540 (gs_f + gs_k) * wct(ixo^s, v_))
541 case (eq_schnakenberg)
542 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
543 + sb_kappa * (sb_alpha - wct(ixo^s, u_) + &
544 wct(ixo^s, u_)**2 * wct(ixo^s, v_)))
545 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
546 + sb_kappa * (sb_beta - wct(ixo^s, u_)**2 * wct(ixo^s, v_)))
547 case (eq_brusselator)
548 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
549 + br_a - (br_b + 1) * wct(ixo^s, u_) &
550 + wct(ixo^s, u_)**2 * wct(ixo^s, v_))
551 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
552 + br_b * wct(ixo^s, u_) - wct(ixo^s, u_)**2 * wct(ixo^s, v_))
553 case (eq_ext_brusselator)
554 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
555 + br_a - (br_b + 1) * wct(ixo^s, u_) &
556 + wct(ixo^s, u_)**2 * wct(ixo^s, v_) &
557 - br_c * wct(ixo^s, u_) + br_d * w(ixo^s, w_))
558 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
559 + br_b * wct(ixo^s, u_) &
560 - wct(ixo^s, u_)**2 * wct(ixo^s, v_))
561 w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
562 + br_c * wct(ixo^s, u_) - br_d * w(ixo^s, w_))
563 case (eq_logistic)
564 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
565 + lg_lambda * w(ixo^s, u_) * (1 - w(ixo^s, u_)))
566 case (eq_analyt_hunds)
567 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
568 + 1.0d0/d1 * w(ixo^s, u_)**2 * (1 - w(ixo^s, u_)))
569 case (eq_belousov_fn)
570 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
571 + 1.0d0/bzfn_epsilon * (bzfn_lambda * wct(ixo^s, u_) &
572 - wct(ixo^s, u_)*wct(ixo^s, w_) + wct(ixo^s, u_) &
573 - wct(ixo^s, u_)**2))
574 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
575 + wct(ixo^s, u_) - wct(ixo^s, v_))
576 w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
577 + 1.0d0/bzfn_delta * (-bzfn_lambda * wct(ixo^s, w_) &
578 - wct(ixo^s, u_)*wct(ixo^s, w_) + bzfn_mu * wct(ixo^s, v_)))
579 case (eq_lorenz)
580 ! xdot = sigma.(y-x)
581 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
582 + lor_sigma * (wct(ixo^s, v_) - wct(ixo^s, u_)))
583 ! ydot = r.x - y - x.z
584 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
585 + lor_r * wct(ixo^s, u_) - wct(ixo^s, v_) &
586 - wct(ixo^s, u_)*wct(ixo^s, w_))
587 ! zdot = x.y - b.z
588 w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
589 + wct(ixo^s, u_)*wct(ixo^s, v_) - lor_b * wct(ixo^s, w_))
590 case (eq_no_reac)
591 w(ixo^s, u_) = w(ixo^s, u_) + qdt * d1 * lpl_u
592 case default
593 call mpistop("Unknown equation type")
594 end select
595
596 ! enforce getbc call after source addition
597 active = .true.
598 end if
599
600 end subroutine ard_add_source
601
602 !> Compute the Laplacian using a standard second order scheme. For now this
603 !> method only works in slab geometries. Requires one ghost cell only.
604 subroutine ard_laplacian(ixI^L,ixO^L,var,lpl)
606 integer, intent(in) :: ixi^l, ixo^l
607 double precision, intent(in) :: var(ixi^s)
608 double precision, intent(out) :: lpl(ixo^s)
609 integer :: idir, jxo^l, hxo^l
610 double precision :: h_inv2
611
612 if (slab) then
613 lpl(ixo^s) = 0.0d0
614 do idir = 1, ndim
615 hxo^l=ixo^l-kr(idir,^d);
616 jxo^l=ixo^l+kr(idir,^d);
617 h_inv2 = 1/dxlevel(idir)**2
618 lpl(ixo^s) = lpl(ixo^s) + h_inv2 * &
619 (var(jxo^s) - 2 * var(ixo^s) + var(hxo^s))
620 end do
621 else
622 call mpistop("ard_laplacian not implemented in this geometry")
623 end if
624
625 end subroutine ard_laplacian
626
627 subroutine put_laplacians_onegrid(ixI^L,ixO^L,w)
629 integer, intent(in) :: ixi^l, ixo^l
630 double precision, intent(inout) :: w(ixi^s, 1:nw)
631
632 double precision :: lpl_u(ixo^s), lpl_v(ixo^s), lpl_w(ixo^s)
633
634 call ard_laplacian(ixi^l, ixo^l, w(ixi^s, u_), lpl_u)
635 if (number_of_species >= 2) then
636 call ard_laplacian(ixi^l, ixo^l, w(ixi^s, v_), lpl_v)
637 end if
638 if (number_of_species >= 3) then
639 call ard_laplacian(ixi^l, ixo^l, w(ixi^s, w_), lpl_w)
640 end if
641
642 w(ixo^s,u_) = d1*lpl_u
643 if (number_of_species >= 2) then
644 w(ixo^s,v_) = d2*lpl_v
645 end if
646 if (number_of_species >= 3) then
647 w(ixo^s,w_) = d3*lpl_w
648 end if
649
650 end subroutine put_laplacians_onegrid
651
652 !> inplace update of psa==>F_im(psa)
653 subroutine ard_evaluate_implicit(qtC,psa)
655 type(state), target :: psa(max_blocks)
656 double precision, intent(in) :: qtc
657
658 integer :: iigrid, igrid, level
659 integer :: ixo^l
660
661 !ixO^L=ixG^LL^LSUB1;
662 ixo^l=ixm^ll;
663 !$OMP PARALLEL DO PRIVATE(igrid)
664 do iigrid=1,igridstail; igrid=igrids(iigrid);
665 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
666 call put_laplacians_onegrid(ixg^ll,ixo^l,psa(igrid)%w)
667 end do
668 !$OMP END PARALLEL DO
669
670 end subroutine ard_evaluate_implicit
671
672 !> Implicit solve of psa=psb+dtfactor*dt*F_im(psa)
673 subroutine ard_implicit_update(dtfactor,qdt,qtC,psa,psb)
675 use mod_forest
676
677 type(state), target :: psa(max_blocks)
678 type(state), target :: psb(max_blocks)
679 double precision, intent(in) :: qdt
680 double precision, intent(in) :: qtc
681 double precision, intent(in) :: dtfactor
682
683 integer :: n
684 double precision :: res, max_residual, lambda
685
686 integer :: iw_to,iw_from
687 integer :: iigrid, igrid, id
688 integer :: nc, lvl
689 type(tree_node), pointer :: pnode
690 real(dp) :: fac
691
692 ! Avoid setting a very restrictive limit to the residual when the time step
693 ! is small (as the operator is ~ 1/(D * qdt))
694 if (qdt < dtmin) then
695 if(mype==0)then
696 print *,'skipping implicit solve: dt too small!'
697 print *,'Currently at time=',global_time,' time step=',qdt,' dtmin=',dtmin
698 endif
699 return
700 endif
701 max_residual = 1d-7/qdt
702
703 mg%operator_type = mg_helmholtz
704 call mg_set_methods(mg)
705
706 if (.not. mg%is_allocated) call mpistop("multigrid tree not allocated yet")
707
708 ! First handle the u variable ***************************************
709 lambda = 1/(dtfactor * qdt * d1)
710 call helmholtz_set_lambda(lambda)
711 mg%bc(:, mg_iphi) = ard_mg_bc(1, :)
712
713 call mg_copy_to_tree(u_, mg_irhs, factor=-lambda, state_from=psb)
714 call mg_copy_to_tree(u_, 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, u_, state_to=psa)
723 ! Done with the u variable ***************************************
724
725 ! Next handle the v variable ***************************************
726 if (number_of_species >= 2) then
727 lambda = 1/(dtfactor * qdt * d2)
728 call helmholtz_set_lambda(lambda)
729 mg%bc(:, mg_iphi) = ard_mg_bc(2, :)
730
731 call mg_copy_to_tree(v_, mg_irhs, factor=-lambda, state_from=psb)
732 call mg_copy_to_tree(v_, mg_iphi, state_from=psb)
733
734 call mg_fas_fmg(mg, .true., max_res=res)
735 do n = 1, 10
736 call mg_fas_vcycle(mg, max_res=res)
737 if (res < max_residual) exit
738 end do
739
740 call mg_copy_from_tree_gc(mg_iphi, v_, state_to=psa)
741 end if
742 ! Done with the v variable ***************************************
743
744 ! Next handle the w variable ***************************************
745 if (number_of_species >= 3) then
746 lambda = 1/(dtfactor * qdt * d3)
747 call helmholtz_set_lambda(lambda)
748
749 call mg_copy_to_tree(w_, mg_irhs, factor=-lambda, state_from=psb)
750 call mg_copy_to_tree(w_, mg_iphi, state_from=psb)
751
752 call mg_fas_fmg(mg, .true., max_res=res)
753 do n = 1, 10
754 call mg_fas_vcycle(mg, max_res=res)
755 if (res < max_residual) exit
756 end do
757
758 call mg_copy_from_tree_gc(mg_iphi, w_, state_to=psa)
759 end if
760 ! Done with the w variable ***************************************
761
762 end subroutine ard_implicit_update
763
764end 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
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
integer coordinate
Definition mod_geometry.t:7
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
integer nwauxio
Number of auxiliary variables that are only included in output.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
integer nghostcells
Number of ghost cells surrounding a grid.
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:54
procedure(sub_write_info), pointer phys_write_info
Definition mod_physics.t:75
procedure(sub_get_flux), pointer phys_get_flux
Definition mod_physics.t:61
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
Definition mod_physics.t:81
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition mod_physics.t:60
procedure(sub_get_dt), pointer phys_get_dt
Definition mod_physics.t:64
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition mod_physics.t:65
procedure(sub_check_params), pointer phys_check_params
Definition mod_physics.t:51
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:53
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition mod_physics.t:49
procedure(sub_implicit_update), pointer phys_implicit_update
Definition mod_physics.t:80
procedure(sub_add_source), pointer phys_add_source
Definition mod_physics.t:66
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:56
logical phys_energy
Solve energy equation or not.
Definition mod_physics.t:37
Module with all the methods that users can customize in AMRVAC.
procedure(special_mg_bc), pointer usr_special_mg_bc
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