MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_rd_phys.t
Go to the documentation of this file.
1!> Reaction-diffusion module (physics routines)
2!>
3!> Multiple reaction-diffusion systems are included: the Gray-Scott model, the
4!> Schnakenberg model, the Brusselator model, the diffusive logistic equation,
5!> an analytical testcase from "Numerical solution of time-dependent advection-
6!> diffusion-reaction equations" by Hundsdorfer & Verwer, the Oregonator model,
7!> the extended Brusselator model, and the diffusive Lorenz system. See the
8!> documentation of the reaction-diffusion module for more information.
9!>
10!> IMEX methods are also supported. The implicit system is solved by a
11!> multigrid solver coupled into MPI-AMRVAC.
12!>
15 use mod_comm_lib, only: mpistop
16
17
18 implicit none
19 private
20
21 integer, protected, public :: u_ = 1
22 integer, protected, public :: v_ = 2 !< For 2 or more equations
23 integer, protected, public :: w_ = 3 !< For 3 or more equations
24
25 !> Whether particles module is added
26 logical, public, protected :: rd_particles = .false.
27
28 !> Parameter with which to multiply the reaction timestep restriction
29 double precision, public, protected :: dtreacpar = 0.5d0
30
31 !> Name of the system to be solved
32 character(len=20), public, protected :: equation_name = "gray-scott"
33 integer :: number_of_species = 2
34 integer :: equation_type = 1
35 integer, parameter :: eq_gray_scott = 1
36 integer, parameter :: eq_schnakenberg = 2
37 integer, parameter :: eq_brusselator = 3
38 integer, parameter :: eq_logistic = 4
39 integer, parameter :: eq_analyt_hunds = 5
40 integer, parameter :: eq_belousov_fn = 6 ! Field-Noyes model, or Oregonator
41 integer, parameter :: eq_ext_brusselator = 7
42 integer, parameter :: eq_lorenz = 8
43
44 !> Diffusion coefficient for first species (u)
45 double precision, public, protected :: d1 = 0.05d0
46 !> Diffusion coefficient for second species (v) (if applicable)
47 double precision, public, protected :: d2 = 1.0d0
48 !> Diffusion coefficient for third species (w) (if applicable)
49 double precision, public, protected :: d3 = 1.0d0
50
51 !> Parameter for Schnakenberg model
52 double precision, public, protected :: sb_alpha = 0.1305d0
53 !> Parameter for Schnakenberg model
54 double precision, public, protected :: sb_beta = 0.7695d0
55 !> Parameter for Schnakenberg model
56 double precision, public, protected :: sb_kappa = 100.0d0
57
58 !> Feed rate for Gray-Scott model
59 double precision, public, protected :: gs_f = 0.046d0
60 !> Kill rate for Gray-Scott model
61 double precision, public, protected :: gs_k = 0.063d0
62
63 !> Parameter for Brusselator model
64 double precision, public, protected :: br_a = 4.5d0
65 !> Parameter for Brusselator model
66 double precision, public, protected :: br_b = 8.0d0
67 !> Parameter for extended Brusselator model
68 double precision, public, protected :: br_c = 1.0d0
69 !> Parameter for extended Brusselator model
70 double precision, public, protected :: br_d = 1.0d0
71
72 !> Parameter for logistic model (Fisher / KPP equation)
73 double precision, public, protected :: lg_lambda = 1.0d0
74
75 !> Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction
76 double precision, public, protected :: bzfn_epsilon = 1.0d0
77 !> Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction
78 double precision, public, protected :: bzfn_delta = 1.0d0
79 !> Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction
80 double precision, public, protected :: bzfn_lambda = 1.0d0
81 !> Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction
82 double precision, public, protected :: bzfn_mu = 1.0d0
83
84 !> Parameter for Lorenz system (Rayleigh number)
85 double precision, public, protected :: lor_r = 28.0d0
86 !> Parameter for Lorenz system (Prandtl number)
87 double precision, public, protected :: lor_sigma = 10.0d0
88 !> Parameter for Lorenz system (aspect ratio of the convection rolls)
89 double precision, public, protected :: lor_b = 8.0d0 / 3.0d0
90
91 !> Whether to handle the explicitly handled source term in split fashion
92 logical :: rd_source_split = .false.
93
94 !> Boundary condition information for the multigrid method
95 type(mg_bc_t), public :: rd_mg_bc(3, mg_num_neighbors)
96
97 ! Public methods
98 public :: rd_phys_init
99
100contains
101
102 subroutine rd_params_read(files)
104 character(len=*), intent(in) :: files(:)
105 integer :: n
106
107 namelist /rd_list/ d1, d2, d3, sb_alpha, sb_beta, sb_kappa, gs_f, gs_k, &
110 equation_name, rd_particles, rd_source_split, dtreacpar
111
112 do n = 1, size(files)
113 open(unitpar, file=trim(files(n)), status='old')
114 read(unitpar, rd_list, end=111)
115111 close(unitpar)
116 end do
117
118 select case (equation_name)
119 case ("gray-scott")
120 equation_type = eq_gray_scott
121 number_of_species = 2
122 case ("schnakenberg")
123 equation_type = eq_schnakenberg
124 number_of_species = 2
125 case ("brusselator")
126 equation_type = eq_brusselator
127 number_of_species = 2
128 case ("logistic")
129 equation_type = eq_logistic
130 number_of_species = 1
131 case ("analyt_hunds")
132 equation_type = eq_analyt_hunds
133 number_of_species = 1
134 case ("belousov_fieldnoyes")
135 equation_type = eq_belousov_fn
136 number_of_species = 3
137 case ("ext_brusselator")
138 equation_type = eq_ext_brusselator
139 number_of_species = 3
140 case ("lorenz")
141 equation_type = eq_lorenz
142 number_of_species = 3
143 case default
144 call mpistop("Unknown equation_name (valid: gray-scott, schnakenberg, ...)")
145 end select
146
147 end subroutine rd_params_read
148
149 !> Write this module's parameters to a snapshot
150 subroutine rd_write_info(fh)
152 integer, intent(in) :: fh
153 integer, parameter :: n_par = 0
154 integer, dimension(MPI_STATUS_SIZE) :: st
155 integer :: er
156 integer :: idim
157
158 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
159 end subroutine rd_write_info
160
161 subroutine rd_phys_init()
163 use mod_physics
165
166 call rd_params_read(par_files)
167
168 physics_type = "rd"
169 phys_energy = .false.
171
172 allocate(start_indices(number_species),stop_indices(number_species))
173 ! set the index of the first flux variable for species 1
174 start_indices(1)=1
175 ! Use the first variable as a density
176 u_ = var_set_fluxvar("u", "u")
177 if (number_of_species >= 2) then
178 v_ = var_set_fluxvar("v", "v")
179 end if
180 if (number_of_species >= 3) then
181 w_ = var_set_fluxvar("w", "w")
182 end if
183
184 ! set number of variables which need update ghostcells
185 nwgc=nwflux
186
187 ! set the index of the last flux variable for species 1
188 stop_indices(1)=nwflux
189
190 ! Disable flux conservation near AMR boundaries, since we have no fluxes
191 fix_conserve_global = .false.
192
193 phys_get_cmax => rd_get_cmax
194 phys_get_cbounds => rd_get_cbounds
195 phys_get_flux => rd_get_flux
196 phys_to_conserved => rd_to_conserved
197 phys_to_primitive => rd_to_primitive
198 phys_add_source => rd_add_source
199 phys_get_dt => rd_get_dt
200 phys_write_info => rd_write_info
201 phys_check_params => rd_check_params
202 phys_implicit_update => rd_implicit_update
203 phys_evaluate_implicit => rd_evaluate_implicit
204
205
206 end subroutine rd_phys_init
207
208 subroutine rd_check_params
210 use mod_geometry, only: coordinate
212 use mod_particles, only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
213
214 integer :: n, i, iw, species_list(number_of_species)
215
216 if (any(flux_method /= fs_source)) then
217 ! there are no fluxes, only source terms in reaction-diffusion
218 call mpistop("mod_rd requires flux_scheme = source")
219 end if
220
221 ! Initialize particles module
222 if (rd_particles) then
223 call particles_init()
224 end if
225
226 if (use_imex_scheme) then
227 use_multigrid = .true.
228 select case(number_of_species)
229 case(1); species_list = [u_]
230 case(2); species_list = [u_, v_]
231 case(3); species_list = [u_, v_, w_]
232 end select
233
234 do i = 1, number_of_species
235 iw = species_list(i)
236
237 ! Set boundary conditions for the multigrid solver
238 do n = 1, 2*ndim
239 select case (typeboundary(iw, n))
240 case (bc_symm)
241 ! d/dx u = 0
242 rd_mg_bc(i, n)%bc_type = mg_bc_neumann
243 rd_mg_bc(i, n)%bc_value = 0.0_dp
244 case (bc_asymm)
245 ! u = 0
246 rd_mg_bc(i, n)%bc_type = mg_bc_dirichlet
247 rd_mg_bc(i, n)%bc_value = 0.0_dp
248 case (bc_cont)
249 ! d/dx u = 0
250 rd_mg_bc(i, n)%bc_type = mg_bc_neumann
251 rd_mg_bc(i, n)%bc_value = 0.0_dp
252 case (bc_periodic)
253 ! Nothing to do here
254 case (bc_special)
255 if (.not. associated(rd_mg_bc(i, n)%boundary_cond)) then
256 write(*, "(A,I0,A,I0,A)") "typeboundary(", iw, ",", n, &
257 ") is 'special', but the corresponding method " // &
258 "rd_mg_bc(i, n)%boundary_cond is not set"
259 call mpistop("rd_mg_bc(i, n)%boundary_cond not set")
260 end if
261 case default
262 write(*,*) "rd_check_params warning: unknown boundary type"
263 rd_mg_bc(i, n)%bc_type = mg_bc_dirichlet
264 rd_mg_bc(i, n)%bc_value = 0.0_dp
265 end select
266 end do
267 end do
268 end if
269
270
271 if(mype==0)then
272 write(*,*)'====RD run with settings===================='
273 write(*,*)'Using mod_rd_phys with settings:'
274 write(*,*)'Dimensionality :',ndim
275 write(*,*)'vector components:',ndir
276 write(*,*)'coordinate set to type,slab:',coordinate,slab
277 write(*,*)'number of variables nw=',nw
278 write(*,*)' start index iwstart=',iwstart
279 write(*,*)'number of vector variables=',nvector
280 write(*,*)'number of stagger variables nws=',nws
281 write(*,*)'number of variables with BCs=',nwgc
282 write(*,*)'number of vars with fluxes=',nwflux
283 write(*,*)'number of vars with flux + BC=',nwfluxbc
284 write(*,*)'number of auxiliary variables=',nwaux
285 write(*,*)'number of extra vars without flux=',nwextra
286 write(*,*)'number of extra vars for wextra=',nw_extra
287 write(*,*)'number of auxiliary I/O variables=',nwauxio
288 write(*,*)' rd_particles=',rd_particles
289 if(rd_particles) then
290 write(*,*) '*****Using particles: npayload,ngridvars :', npayload,ngridvars
291 write(*,*) '*****Using particles: nusrpayload :', nusrpayload
292 write(*,*) '*****Using particles: num_particles :', num_particles
293 write(*,*) '*****Using particles: physics_type_particles=',physics_type_particles
294 end if
295 write(*,*)'number of ghostcells=',nghostcells
296 write(*,*)'==========================================='
297 endif
298
299
300
301 end subroutine rd_check_params
302
303 subroutine rd_to_conserved(ixI^L, ixO^L, w, x)
305 integer, intent(in) :: ixi^l, ixo^l
306 double precision, intent(inout) :: w(ixi^s, nw)
307 double precision, intent(in) :: x(ixi^s, 1:^nd)
308
309 ! Do nothing (primitive and conservative are equal for rd module)
310 end subroutine rd_to_conserved
311
312 subroutine rd_to_primitive(ixI^L, ixO^L, w, x)
314 integer, intent(in) :: ixi^l, ixo^l
315 double precision, intent(inout) :: w(ixi^s, nw)
316 double precision, intent(in) :: x(ixi^s, 1:^nd)
317
318 ! Do nothing (primitive and conservative are equal for rd module)
319 end subroutine rd_to_primitive
320
321 subroutine rd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
323 integer, intent(in) :: ixi^l, ixo^l, idim
324 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
325 double precision, intent(inout) :: cmax(ixi^s)
326
327 cmax(ixo^s) = 0.0d0
328 end subroutine rd_get_cmax
329
330 subroutine rd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
332 use mod_variables
333 integer, intent(in) :: ixi^l, ixo^l, idim
334 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s,nw)
335 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s,nw)
336 double precision, intent(in) :: x(ixi^s, 1:^nd)
337 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
338 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
339 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
340
341 if (present(cmin)) then
342 cmin(ixo^s,1) = 0.0d0
343 cmax(ixo^s,1) = 0.0d0
344 else
345 cmax(ixo^s,1) = 0.0d0
346 end if
347
348 end subroutine rd_get_cbounds
349
350 subroutine rd_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
352 integer, intent(in) :: ixi^l, ixo^l
353 double precision, intent(in) :: dx^d, x(ixi^s, 1:^nd)
354 double precision, intent(in) :: w(ixi^s, 1:nw)
355 double precision, intent(inout) :: dtnew
356 double precision :: maxrate
357 double precision :: maxd
358
359 ! dt < dx^2 / (2 * ndim * diffusion_coeff)
360 ! use dtdiffpar < 1 for explicit and > 1 for imex/split
361 maxd = d1
362 if (number_of_species >= 2) then
363 maxd = max(maxd, d2)
364 end if
365 if (number_of_species >= 3) then
366 maxd = max(maxd, d3)
367 end if
368 dtnew = dtdiffpar * minval([ dx^d ])**2 / (2 * ndim * maxd)
369
370 ! Estimate time step for reactions
371 select case (equation_type)
372 case (eq_gray_scott)
373 maxrate = max(maxval(w(ixo^s, v_))**2 + gs_f, &
374 maxval(w(ixo^s, v_) * w(ixo^s, u_)) - gs_f - gs_k)
375 case (eq_schnakenberg)
376 maxrate = max(maxval(abs(w(ixo^s, v_) * w(ixo^s, u_) - 1)), &
377 maxval(w(ixo^s, u_))**2)
378 case (eq_brusselator)
379 maxrate = max( maxval(w(ixo^s, u_)*w(ixo^s, v_) - (br_b+1)), &
380 maxval(w(ixo^s, u_)**2) )
381 case (eq_ext_brusselator)
382 maxrate = max( maxval(w(ixo^s, u_)*w(ixo^s, v_) - (br_b+1)) + br_c, &
383 maxval(w(ixo^s, u_)**2) )
384 maxrate = max(maxrate, br_d)
385 case (eq_logistic)
386 maxrate = lg_lambda*maxval(abs(1 - w(ixo^s, u_))) ! abs for safety, normally u < 1
387 case (eq_analyt_hunds)
388 maxrate = maxval(w(ixo^s, u_)*abs(1 - w(ixo^s, u_))) / d1
389 case (eq_belousov_fn)
390 maxrate = max(&
391 maxval(abs(1.0d0 - w(ixo^s, w_) - w(ixo^s, u_))) / bzfn_epsilon, &
392 maxval(bzfn_lambda + w(ixo^s, u_)) / bzfn_delta &
393 )
394 case (eq_lorenz)
395 ! det(J) = sigma(b(r-1) + x*(x*+y*))
396 maxrate = max(lor_sigma, 1.0d0, lor_b)
397 case default
398 maxrate = one
399 call mpistop("Unknown equation type")
400 end select
401
402 dtnew = min(dtnew, dtreacpar / maxrate)
403
404 end subroutine rd_get_dt
405
406 ! There is nothing to add to the transport flux in the transport equation
407 subroutine rd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
409 integer, intent(in) :: ixi^l, ixo^l, idim
410 double precision, intent(in) :: wc(ixi^s, 1:nw)
411 double precision, intent(in) :: w(ixi^s, 1:nw)
412 double precision, intent(in) :: x(ixi^s, 1:^nd)
413 double precision, intent(out) :: f(ixi^s, nwflux)
414
415 f(ixo^s, :) = 0.0d0
416 end subroutine rd_get_flux
417
418 ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
419 subroutine rd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
421 integer, intent(in) :: ixi^l, ixo^l
422 double precision, intent(in) :: qdt, dtfactor
423 double precision, intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw),x(ixi^s, 1:ndim)
424 double precision, intent(inout) :: w(ixi^s, 1:nw)
425 double precision :: lpl_u(ixo^s), lpl_v(ixo^s), lpl_w(ixo^s)
426 logical, intent(in) :: qsourcesplit
427 logical, intent(inout) :: active
428
429 ! here we add the reaction terms (always) and the diffusion if no imex is used
430 if (qsourcesplit .eqv. rd_source_split) then
431 if (.not.use_imex_scheme) then
432 call rd_laplacian(ixi^l, ixo^l, wct(ixi^s, u_), lpl_u)
433 if (number_of_species >= 2) then
434 call rd_laplacian(ixi^l, ixo^l, wct(ixi^s, v_), lpl_v)
435 end if
436 if (number_of_species >= 3) then
437 call rd_laplacian(ixi^l, ixo^l, wct(ixi^s, w_), lpl_w)
438 end if
439 else
440 ! for all IMEX scheme variants: only add the reactions
441 lpl_u = 0.0d0
442 lpl_v = 0.0d0
443 lpl_w = 0.0d0
444 end if
445
446 select case (equation_type)
447 case (eq_gray_scott)
448 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u - &
449 wct(ixo^s, u_) * wct(ixo^s, v_)**2 + &
450 gs_f * (1 - wct(ixo^s, u_)))
451 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v + &
452 wct(ixo^s, u_) * wct(ixo^s, v_)**2 - &
453 (gs_f + gs_k) * wct(ixo^s, v_))
454 case (eq_schnakenberg)
455 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
456 + sb_kappa * (sb_alpha - wct(ixo^s, u_) + &
457 wct(ixo^s, u_)**2 * wct(ixo^s, v_)))
458 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
459 + sb_kappa * (sb_beta - wct(ixo^s, u_)**2 * wct(ixo^s, v_)))
460 case (eq_brusselator)
461 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
462 + br_a - (br_b + 1) * wct(ixo^s, u_) &
463 + wct(ixo^s, u_)**2 * wct(ixo^s, v_))
464 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
465 + br_b * wct(ixo^s, u_) - wct(ixo^s, u_)**2 * wct(ixo^s, v_))
466 case (eq_ext_brusselator)
467 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
468 + br_a - (br_b + 1) * wct(ixo^s, u_) &
469 + wct(ixo^s, u_)**2 * wct(ixo^s, v_) &
470 - br_c * wct(ixo^s, u_) + br_d * w(ixo^s, w_))
471 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
472 + br_b * wct(ixo^s, u_) &
473 - wct(ixo^s, u_)**2 * wct(ixo^s, v_))
474 w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
475 + br_c * wct(ixo^s, u_) - br_d * w(ixo^s, w_))
476 case (eq_logistic)
477 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
478 + lg_lambda * w(ixo^s, u_) * (1 - w(ixo^s, u_)))
479 case (eq_analyt_hunds)
480 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
481 + 1.0d0/d1 * w(ixo^s, u_)**2 * (1 - w(ixo^s, u_)))
482 case (eq_belousov_fn)
483 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
484 + 1.0d0/bzfn_epsilon * (bzfn_lambda * wct(ixo^s, u_) &
485 - wct(ixo^s, u_)*wct(ixo^s, w_) + wct(ixo^s, u_) &
486 - wct(ixo^s, u_)**2))
487 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
488 + wct(ixo^s, u_) - wct(ixo^s, v_))
489 w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
490 + 1.0d0/bzfn_delta * (-bzfn_lambda * wct(ixo^s, w_) &
491 - wct(ixo^s, u_)*wct(ixo^s, w_) + bzfn_mu * wct(ixo^s, v_)))
492 case (eq_lorenz)
493 ! xdot = sigma.(y-x)
494 w(ixo^s, u_) = w(ixo^s, u_) + qdt * (d1 * lpl_u &
495 + lor_sigma * (wct(ixo^s, v_) - wct(ixo^s, u_)))
496 ! ydot = r.x - y - x.z
497 w(ixo^s, v_) = w(ixo^s, v_) + qdt * (d2 * lpl_v &
498 + lor_r * wct(ixo^s, u_) - wct(ixo^s, v_) &
499 - wct(ixo^s, u_)*wct(ixo^s, w_))
500 ! zdot = x.y - b.z
501 w(ixo^s, w_) = w(ixo^s, w_) + qdt * (d3 * lpl_w &
502 + wct(ixo^s, u_)*wct(ixo^s, v_) - lor_b * wct(ixo^s, w_))
503 case default
504 call mpistop("Unknown equation type")
505 end select
506
507 ! enforce getbc call after source addition
508 active = .true.
509 end if
510
511 end subroutine rd_add_source
512
513 !> Compute the Laplacian using a standard second order scheme. For now this
514 !> method only works in slab geometries. Requires one ghost cell only.
515 subroutine rd_laplacian(ixI^L,ixO^L,var,lpl)
517 integer, intent(in) :: ixi^l, ixo^l
518 double precision, intent(in) :: var(ixi^s)
519 double precision, intent(out) :: lpl(ixo^s)
520 integer :: idir, jxo^l, hxo^l
521 double precision :: h_inv2
522
523 if (slab) then
524 lpl(ixo^s) = 0.0d0
525 do idir = 1, ndim
526 hxo^l=ixo^l-kr(idir,^d);
527 jxo^l=ixo^l+kr(idir,^d);
528 h_inv2 = 1/dxlevel(idir)**2
529 lpl(ixo^s) = lpl(ixo^s) + h_inv2 * &
530 (var(jxo^s) - 2 * var(ixo^s) + var(hxo^s))
531 end do
532 else
533 call mpistop("rd_laplacian not implemented in this geometry")
534 end if
535
536 end subroutine rd_laplacian
537
538 subroutine put_laplacians_onegrid(ixI^L,ixO^L,w)
540 integer, intent(in) :: ixi^l, ixo^l
541 double precision, intent(inout) :: w(ixi^s, 1:nw)
542
543 double precision :: lpl_u(ixo^s), lpl_v(ixo^s), lpl_w(ixo^s)
544
545 call rd_laplacian(ixi^l, ixo^l, w(ixi^s, u_), lpl_u)
546 if (number_of_species >= 2) then
547 call rd_laplacian(ixi^l, ixo^l, w(ixi^s, v_), lpl_v)
548 end if
549 if (number_of_species >= 3) then
550 call rd_laplacian(ixi^l, ixo^l, w(ixi^s, w_), lpl_w)
551 end if
552
553 w(ixo^s,u_) = d1*lpl_u
554 if (number_of_species >= 2) then
555 w(ixo^s,v_) = d2*lpl_v
556 end if
557 if (number_of_species >= 3) then
558 w(ixo^s,w_) = d3*lpl_w
559 end if
560
561 end subroutine put_laplacians_onegrid
562
563 !> inplace update of psa==>F_im(psa)
564 subroutine rd_evaluate_implicit(qtC,psa)
566 type(state), target :: psa(max_blocks)
567 double precision, intent(in) :: qtc
568
569 integer :: iigrid, igrid, level
570 integer :: ixo^l
571
572 !ixO^L=ixG^LL^LSUB1;
573 ixo^l=ixm^ll;
574 !$OMP PARALLEL DO PRIVATE(igrid)
575 do iigrid=1,igridstail; igrid=igrids(iigrid);
576 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
577 call put_laplacians_onegrid(ixg^ll,ixo^l,psa(igrid)%w)
578 end do
579 !$OMP END PARALLEL DO
580
581 end subroutine rd_evaluate_implicit
582
583 !> Implicit solve of psa=psb+dtfactor*dt*F_im(psa)
584 subroutine rd_implicit_update(dtfactor,qdt,qtC,psa,psb)
586 use mod_forest
587
588 type(state), target :: psa(max_blocks)
589 type(state), target :: psb(max_blocks)
590 double precision, intent(in) :: qdt
591 double precision, intent(in) :: qtc
592 double precision, intent(in) :: dtfactor
593
594 integer :: n
595 double precision :: res, max_residual, lambda
596
597 integer :: iw_to,iw_from
598 integer :: iigrid, igrid, id
599 integer :: nc, lvl
600 type(tree_node), pointer :: pnode
601 real(dp) :: fac
602
603 ! Avoid setting a very restrictive limit to the residual when the time step
604 ! is small (as the operator is ~ 1/(D * qdt))
605 if (qdt < dtmin) then
606 if(mype==0)then
607 print *,'skipping implicit solve: dt too small!'
608 print *,'Currently at time=',global_time,' time step=',qdt,' dtmin=',dtmin
609 endif
610 return
611 endif
612 max_residual = 1d-7/qdt
613
614 mg%operator_type = mg_helmholtz
615 call mg_set_methods(mg)
616
617 if (.not. mg%is_allocated) call mpistop("multigrid tree not allocated yet")
618
619 ! First handle the u variable ***************************************
620 lambda = 1/(dtfactor * qdt * d1)
621 call helmholtz_set_lambda(lambda)
622 mg%bc(:, mg_iphi) = rd_mg_bc(1, :)
623
624 call mg_copy_to_tree(u_, mg_irhs, factor=-lambda, state_from=psb)
625 call mg_copy_to_tree(u_, mg_iphi, state_from=psb)
626
627 call mg_fas_fmg(mg, .true., max_res=res)
628 do n = 1, 10
629 call mg_fas_vcycle(mg, max_res=res)
630 if (res < max_residual) exit
631 end do
632
633 if (res > max_residual) call mpistop("u_multigrid: no convergence")
634
635 call mg_copy_from_tree_gc(mg_iphi, u_, state_to=psa)
636 ! Done with the u variable ***************************************
637
638 ! Next handle the v variable ***************************************
639 if (number_of_species >= 2) then
640 lambda = 1/(dtfactor * qdt * d2)
641 call helmholtz_set_lambda(lambda)
642 mg%bc(:, mg_iphi) = rd_mg_bc(2, :)
643
644 call mg_copy_to_tree(v_, mg_irhs, factor=-lambda, state_from=psb)
645 call mg_copy_to_tree(v_, mg_iphi, state_from=psb)
646
647 call mg_fas_fmg(mg, .true., max_res=res)
648 do n = 1, 10
649 call mg_fas_vcycle(mg, max_res=res)
650 if (res < max_residual) exit
651 end do
652
653 if (res > max_residual) call mpistop("v_multigrid: no convergence")
654 call mg_copy_from_tree_gc(mg_iphi, v_, state_to=psa)
655 end if
656 ! Done with the v variable ***************************************
657
658 ! Next handle the w variable ***************************************
659 if (number_of_species >= 3) then
660 lambda = 1/(dtfactor * qdt * d3)
661 call helmholtz_set_lambda(lambda)
662
663 call mg_copy_to_tree(w_, mg_irhs, factor=-lambda, state_from=psb)
664 call mg_copy_to_tree(w_, mg_iphi, state_from=psb)
665
666 call mg_fas_fmg(mg, .true., max_res=res)
667 do n = 1, 10
668 call mg_fas_vcycle(mg, max_res=res)
669 if (res < max_residual) exit
670 end do
671
672 if (res > max_residual) call mpistop("w_multigrid: no convergence")
673 call mg_copy_from_tree_gc(mg_iphi, w_, state_to=psa)
674 end if
675 ! Done with the w variable ***************************************
676
677 end subroutine rd_implicit_update
678
679end module mod_rd_phys
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
integer, dimension(:), allocatable flux_method
Which flux scheme of spatial discretization to use (per grid level)
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.
logical fix_conserve_global
Whether to apply flux conservation at refinement boundaries.
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_check_params), pointer phys_check_params
Definition mod_physics.t:51
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
Reaction-diffusion module (physics routines)
Definition mod_rd_phys.t:13
double precision, public, protected gs_f
Feed rate for Gray-Scott model.
Definition mod_rd_phys.t:59
double precision, public, protected lor_b
Parameter for Lorenz system (aspect ratio of the convection rolls)
Definition mod_rd_phys.t:89
integer, public, protected u_
Definition mod_rd_phys.t:21
double precision, public, protected bzfn_lambda
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
Definition mod_rd_phys.t:80
type(mg_bc_t), dimension(3, mg_num_neighbors), public rd_mg_bc
Boundary condition information for the multigrid method.
Definition mod_rd_phys.t:95
subroutine, public rd_phys_init()
double precision, public, protected d3
Diffusion coefficient for third species (w) (if applicable)
Definition mod_rd_phys.t:49
logical, public, protected rd_particles
Whether particles module is added.
Definition mod_rd_phys.t:26
double precision, public, protected sb_kappa
Parameter for Schnakenberg model.
Definition mod_rd_phys.t:56
double precision, public, protected br_d
Parameter for extended Brusselator model.
Definition mod_rd_phys.t:70
double precision, public, protected gs_k
Kill rate for Gray-Scott model.
Definition mod_rd_phys.t:61
double precision, public, protected d2
Diffusion coefficient for second species (v) (if applicable)
Definition mod_rd_phys.t:47
double precision, public, protected dtreacpar
Parameter with which to multiply the reaction timestep restriction.
Definition mod_rd_phys.t:29
integer, public, protected w_
For 3 or more equations.
Definition mod_rd_phys.t:23
double precision, public, protected bzfn_epsilon
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
Definition mod_rd_phys.t:76
double precision, public, protected d1
Diffusion coefficient for first species (u)
Definition mod_rd_phys.t:45
double precision, public, protected lg_lambda
Parameter for logistic model (Fisher / KPP equation)
Definition mod_rd_phys.t:73
integer, public, protected v_
For 2 or more equations.
Definition mod_rd_phys.t:22
double precision, public, protected sb_beta
Parameter for Schnakenberg model.
Definition mod_rd_phys.t:54
double precision, public, protected bzfn_mu
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
Definition mod_rd_phys.t:82
double precision, public, protected sb_alpha
Parameter for Schnakenberg model.
Definition mod_rd_phys.t:52
double precision, public, protected bzfn_delta
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
Definition mod_rd_phys.t:78
double precision, public, protected lor_r
Parameter for Lorenz system (Rayleigh number)
Definition mod_rd_phys.t:85
double precision, public, protected br_a
Parameter for Brusselator model.
Definition mod_rd_phys.t:64
double precision, public, protected lor_sigma
Parameter for Lorenz system (Prandtl number)
Definition mod_rd_phys.t:87
double precision, public, protected br_c
Parameter for extended Brusselator model.
Definition mod_rd_phys.t:68
character(len=20), public, protected equation_name
Name of the system to be solved.
Definition mod_rd_phys.t:32
double precision, public, protected br_b
Parameter for Brusselator model.
Definition mod_rd_phys.t:66
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
The data structure that contains information about a tree node/grid block.
Definition mod_forest.t:11