21 integer,
protected,
public ::
u_ = 1
22 integer,
protected,
public ::
v_ = 2
23 integer,
protected,
public ::
w_ = 3
29 double precision,
public,
protected ::
dtreacpar = 0.5d0
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
41 integer,
parameter :: eq_ext_brusselator = 7
42 integer,
parameter :: eq_lorenz = 8
45 double precision,
public,
protected ::
d1 = 0.05d0
47 double precision,
public,
protected ::
d2 = 1.0d0
49 double precision,
public,
protected ::
d3 = 1.0d0
52 double precision,
public,
protected ::
sb_alpha = 0.1305d0
54 double precision,
public,
protected ::
sb_beta = 0.7695d0
56 double precision,
public,
protected ::
sb_kappa = 100.0d0
59 double precision,
public,
protected ::
gs_f = 0.046d0
61 double precision,
public,
protected ::
gs_k = 0.063d0
64 double precision,
public,
protected ::
br_a = 4.5d0
66 double precision,
public,
protected ::
br_b = 8.0d0
68 double precision,
public,
protected ::
br_c = 1.0d0
70 double precision,
public,
protected ::
br_d = 1.0d0
73 double precision,
public,
protected ::
lg_lambda = 1.0d0
78 double precision,
public,
protected ::
bzfn_delta = 1.0d0
82 double precision,
public,
protected ::
bzfn_mu = 1.0d0
85 double precision,
public,
protected ::
lor_r = 28.0d0
87 double precision,
public,
protected ::
lor_sigma = 10.0d0
89 double precision,
public,
protected ::
lor_b = 8.0d0 / 3.0d0
92 logical :: rd_source_split = .false.
95 type(mg_bc_t),
public ::
rd_mg_bc(3, mg_num_neighbors)
102 subroutine rd_params_read(files)
104 character(len=*),
intent(in) :: files(:)
112 do n = 1,
size(files)
113 open(
unitpar, file=trim(files(n)), status=
'old')
114 read(
unitpar, rd_list,
end=111)
120 equation_type = eq_gray_scott
121 number_of_species = 2
122 case (
"schnakenberg")
123 equation_type = eq_schnakenberg
124 number_of_species = 2
126 equation_type = eq_brusselator
127 number_of_species = 2
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
141 equation_type = eq_lorenz
142 number_of_species = 3
144 call mpistop(
"Unknown equation_name (valid: gray-scott, schnakenberg, ...)")
147 end subroutine rd_params_read
152 integer,
intent(in) :: fh
153 integer,
parameter :: n_par = 0
154 integer,
dimension(MPI_STATUS_SIZE) :: st
158 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
174 allocate(start_indices(number_species),stop_indices(number_species))
178 u_ = var_set_fluxvar(
"u",
"u")
179 if (number_of_species >= 2)
then
180 v_ = var_set_fluxvar(
"v",
"v")
182 if (number_of_species >= 3)
then
183 w_ = var_set_fluxvar(
"w",
"w")
190 stop_indices(1)=nwflux
217 integer :: n, i, iw, species_list(number_of_species)
221 call mpistop(
"mod_rd requires flux_scheme = source")
226 select case(number_of_species)
227 case(1); species_list = [
u_]
228 case(2); species_list = [
u_,
v_]
229 case(3); species_list = [
u_,
v_,
w_]
232 do i = 1, number_of_species
240 rd_mg_bc(i, n)%bc_type = mg_bc_neumann
244 rd_mg_bc(i, n)%bc_type = mg_bc_dirichlet
248 rd_mg_bc(i, n)%bc_type = mg_bc_neumann
253 if (.not.
associated(
rd_mg_bc(i, n)%boundary_cond))
then
254 write(*,
"(A,I0,A,I0,A)")
"typeboundary(", iw,
",", n, &
255 ") is 'special', but the corresponding method " // &
256 "rd_mg_bc(i, n)%boundary_cond is not set"
257 call mpistop(
"rd_mg_bc(i, n)%boundary_cond not set")
260 write(*,*)
"rd_check_params warning: unknown boundary type"
261 rd_mg_bc(i, n)%bc_type = mg_bc_dirichlet
272 integer,
intent(in) :: ixI^L, ixO^L
273 double precision,
intent(inout) :: w(ixI^S, nw)
274 double precision,
intent(in) :: x(ixI^S, 1:^ND)
281 integer,
intent(in) :: ixI^L, ixO^L
282 double precision,
intent(inout) :: w(ixI^S, nw)
283 double precision,
intent(in) :: x(ixI^S, 1:^ND)
290 integer,
intent(in) :: ixI^L, ixO^L, idim
291 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
292 double precision,
intent(inout) :: cmax(ixI^S)
297 subroutine rd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
300 integer,
intent(in) :: ixI^L, ixO^L, idim
301 double precision,
intent(in) :: wLC(ixI^S, nw), wRC(ixI^S,nw)
302 double precision,
intent(in) :: wLp(ixI^S, nw), wRp(ixI^S,nw)
303 double precision,
intent(in) :: x(ixI^S, 1:^ND)
304 double precision,
intent(inout) :: cmax(ixI^S,1:number_species)
305 double precision,
intent(inout),
optional :: cmin(ixI^S,1:number_species)
306 double precision,
intent(in) :: Hspeed(ixI^S,1:number_species)
308 if (
present(cmin))
then
309 cmin(ixo^s,1) = 0.0d0
310 cmax(ixo^s,1) = 0.0d0
312 cmax(ixo^s,1) = 0.0d0
319 integer,
intent(in) :: ixI^L, ixO^L
320 double precision,
intent(in) :: dx^D, x(ixI^S, 1:^ND)
321 double precision,
intent(in) :: w(ixI^S, 1:nw)
322 double precision,
intent(inout) :: dtnew
323 double precision :: maxrate
324 double precision :: maxD
329 if (number_of_species >= 2)
then
332 if (number_of_species >= 3)
then
338 select case (equation_type)
340 maxrate = max(maxval(w(ixo^s,
v_))**2 +
gs_f, &
342 case (eq_schnakenberg)
343 maxrate = max(maxval(abs(w(ixo^s,
v_) * w(ixo^s,
u_) - 1)), &
344 maxval(w(ixo^s,
u_))**2)
345 case (eq_brusselator)
346 maxrate = max( maxval(w(ixo^s,
u_)*w(ixo^s,
v_) - (
br_b+1)), &
347 maxval(w(ixo^s,
u_)**2) )
348 case (eq_ext_brusselator)
349 maxrate = max( maxval(w(ixo^s,
u_)*w(ixo^s,
v_) - (
br_b+1)) +
br_c, &
350 maxval(w(ixo^s,
u_)**2) )
351 maxrate = max(maxrate,
br_d)
354 case (eq_analyt_hunds)
355 maxrate = maxval(w(ixo^s,
u_)*abs(1 - w(ixo^s,
u_))) /
d1
356 case (eq_belousov_fn)
366 call mpistop(
"Unknown equation type")
376 integer,
intent(in) :: ixI^L, ixO^L, idim
377 double precision,
intent(in) :: wC(ixI^S, 1:nw)
378 double precision,
intent(in) :: w(ixI^S, 1:nw)
379 double precision,
intent(in) :: x(ixI^S, 1:^ND)
380 double precision,
intent(out) :: f(ixI^S, nwflux)
386 subroutine rd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
388 integer,
intent(in) :: ixI^L, ixO^L
389 double precision,
intent(in) :: qdt, dtfactor
390 double precision,
intent(in) :: wCT(ixI^S, 1:nw),wCTprim(ixI^S,1:nw),x(ixI^S, 1:ndim)
391 double precision,
intent(inout) :: w(ixI^S, 1:nw)
392 double precision :: lpl_u(ixO^S), lpl_v(ixO^S), lpl_w(ixO^S)
393 logical,
intent(in) :: qsourcesplit
394 logical,
intent(inout) :: active
397 if (qsourcesplit .eqv. rd_source_split)
then
400 if (number_of_species >= 2)
then
403 if (number_of_species >= 3)
then
413 select case (equation_type)
415 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u - &
416 wct(ixo^s,
u_) * wct(ixo^s,
v_)**2 + &
417 gs_f * (1 - wct(ixo^s,
u_)))
418 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v + &
419 wct(ixo^s,
u_) * wct(ixo^s,
v_)**2 - &
421 case (eq_schnakenberg)
422 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
424 wct(ixo^s,
u_)**2 * wct(ixo^s,
v_)))
425 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
427 case (eq_brusselator)
428 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
430 + wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
431 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
432 +
br_b * wct(ixo^s,
u_) - wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
433 case (eq_ext_brusselator)
434 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
436 + wct(ixo^s,
u_)**2 * wct(ixo^s,
v_) &
438 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
440 - wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
441 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
444 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
446 case (eq_analyt_hunds)
447 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
448 + 1.0d0/
d1 * w(ixo^s,
u_)**2 * (1 - w(ixo^s,
u_)))
449 case (eq_belousov_fn)
450 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
452 - wct(ixo^s,
u_)*wct(ixo^s,
w_) + wct(ixo^s,
u_) &
453 - wct(ixo^s,
u_)**2))
454 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
455 + wct(ixo^s,
u_) - wct(ixo^s,
v_))
456 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
461 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
464 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
465 +
lor_r * wct(ixo^s,
u_) - wct(ixo^s,
v_) &
466 - wct(ixo^s,
u_)*wct(ixo^s,
w_))
468 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
469 + wct(ixo^s,
u_)*wct(ixo^s,
v_) -
lor_b * wct(ixo^s,
w_))
471 call mpistop(
"Unknown equation type")
484 integer,
intent(in) :: ixI^L, ixO^L
485 double precision,
intent(in) :: var(ixI^S)
486 double precision,
intent(out) :: lpl(ixO^S)
487 integer :: idir, jxO^L, hxO^L
488 double precision :: h_inv2
493 hxo^l=ixo^l-
kr(idir,^
d);
494 jxo^l=ixo^l+
kr(idir,^
d);
496 lpl(ixo^s) = lpl(ixo^s) + h_inv2 * &
497 (var(jxo^s) - 2 * var(ixo^s) + var(hxo^s))
500 call mpistop(
"rd_laplacian not implemented in this geometry")
507 integer,
intent(in) :: ixI^L, ixO^L
508 double precision,
intent(inout) :: w(ixI^S, 1:nw)
510 double precision :: lpl_u(ixO^S), lpl_v(ixO^S), lpl_w(ixO^S)
513 if (number_of_species >= 2)
then
516 if (number_of_species >= 3)
then
520 w(ixo^s,
u_) =
d1*lpl_u
521 if (number_of_species >= 2)
then
522 w(ixo^s,
v_) =
d2*lpl_v
524 if (number_of_species >= 3)
then
525 w(ixo^s,
w_) =
d3*lpl_w
533 type(state),
target :: psa(max_blocks)
534 double precision,
intent(in) :: qtC
536 integer :: iigrid, igrid, level
542 do iigrid=1,igridstail; igrid=igrids(iigrid);
555 type(state),
target :: psa(max_blocks)
556 type(state),
target :: psb(max_blocks)
557 double precision,
intent(in) :: qdt
558 double precision,
intent(in) :: qtC
559 double precision,
intent(in) :: dtfactor
562 double precision :: res, max_residual, lambda
564 integer :: iw_to,iw_from
565 integer :: iigrid, igrid, id
572 if (qdt <
dtmin)
then
574 print *,
'skipping implicit solve: dt too small!'
579 max_residual = 1
d-7/qdt
581 mg%operator_type = mg_helmholtz
582 call mg_set_methods(mg)
584 if (.not. mg%is_allocated)
call mpistop(
"multigrid tree not allocated yet")
587 lambda = 1/(dtfactor * qdt *
d1)
588 call helmholtz_set_lambda(lambda)
591 call mg_copy_to_tree(
u_, mg_irhs, factor=-lambda, state_from=psb)
592 call mg_copy_to_tree(
u_, mg_iphi, state_from=psb)
594 call mg_fas_fmg(mg, .true., max_res=res)
596 call mg_fas_vcycle(mg, max_res=res)
597 if (res < max_residual)
exit
600 call mg_copy_from_tree_gc(mg_iphi,
u_, state_to=psa)
604 if (number_of_species >= 2)
then
605 lambda = 1/(dtfactor * qdt *
d2)
606 call helmholtz_set_lambda(lambda)
609 call mg_copy_to_tree(
v_, mg_irhs, factor=-lambda, state_from=psb)
610 call mg_copy_to_tree(
v_, mg_iphi, state_from=psb)
612 call mg_fas_fmg(mg, .true., max_res=res)
614 call mg_fas_vcycle(mg, max_res=res)
615 if (res < max_residual)
exit
618 call mg_copy_from_tree_gc(mg_iphi,
v_, state_to=psa)
623 if (number_of_species >= 3)
then
624 lambda = 1/(dtfactor * qdt *
d3)
625 call helmholtz_set_lambda(lambda)
627 call mg_copy_to_tree(
w_, mg_irhs, factor=-lambda, state_from=psb)
628 call mg_copy_to_tree(
w_, mg_iphi, state_from=psb)
630 call mg_fas_fmg(mg, .true., max_res=res)
632 call mg_fas_vcycle(mg, max_res=res)
633 if (res < max_residual)
exit
636 call mg_copy_from_tree_gc(mg_iphi,
w_, state_to=psa)
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with basic grid data structures.
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision dtdiffpar
For resistive MHD, the time step is also limited by the diffusion time: .
integer, parameter unitpar
file handle for IO
integer, parameter bc_asymm
double precision global_time
The global simulation time.
logical use_imex_scheme
whether IMEX in use or not
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, dimension(:, :), allocatable typeboundary
Array indicating the type of boundary condition per variable and per physical boundary.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical use_particles
Use particles module or not.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
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
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, parameter bc_cont
integer, parameter bc_symm
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, parameter fs_source
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...
procedure(sub_convert), pointer phys_to_primitive
procedure(sub_write_info), pointer phys_write_info
procedure(sub_get_flux), pointer phys_get_flux
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
logical phys_req_diagonal
Whether the physics routines require diagonal ghost cells, for example for computing a curl.
procedure(sub_get_cbounds), pointer phys_get_cbounds
procedure(sub_get_dt), pointer phys_get_dt
procedure(sub_check_params), pointer phys_check_params
procedure(sub_convert), pointer phys_to_conserved
character(len=name_len) physics_type
String describing the physics type of the simulation.
procedure(sub_implicit_update), pointer phys_implicit_update
procedure(sub_add_source), pointer phys_add_source
procedure(sub_get_cmax), pointer phys_get_cmax
logical phys_energy
Solve energy equation or not.
Reaction-diffusion module (physics routines)
subroutine rd_write_info(fh)
Write this module's parameters to a snapshot.
double precision, public, protected gs_f
Feed rate for Gray-Scott model.
double precision, public, protected lor_b
Parameter for Lorenz system (aspect ratio of the convection rolls)
integer, public, protected u_
double precision, public, protected bzfn_lambda
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
type(mg_bc_t), dimension(3, mg_num_neighbors), public rd_mg_bc
Boundary condition information for the multigrid method.
subroutine, public rd_phys_init()
subroutine rd_get_cmax(w, x, ixIL, ixOL, idim, cmax)
double precision, public, protected d3
Diffusion coefficient for third species (w) (if applicable)
logical, public, protected rd_particles
Whether particles module is added.
double precision, public, protected sb_kappa
Parameter for Schnakenberg model.
double precision, public, protected br_d
Parameter for extended Brusselator model.
double precision, public, protected gs_k
Kill rate for Gray-Scott model.
double precision, public, protected d2
Diffusion coefficient for second species (v) (if applicable)
subroutine rd_evaluate_implicit(qtC, psa)
inplace update of psa==>F_im(psa)
double precision, public, protected dtreacpar
Parameter with which to multiply the reaction timestep restriction.
integer, public, protected w_
For 3 or more equations.
subroutine rd_get_cbounds(wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
double precision, public, protected bzfn_epsilon
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
double precision, public, protected d1
Diffusion coefficient for first species (u)
double precision, public, protected lg_lambda
Parameter for logistic model (Fisher / KPP equation)
integer, public, protected v_
For 2 or more equations.
double precision, public, protected sb_beta
Parameter for Schnakenberg model.
subroutine rd_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
double precision, public, protected bzfn_mu
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
double precision, public, protected sb_alpha
Parameter for Schnakenberg model.
subroutine rd_check_params
subroutine rd_to_primitive(ixIL, ixOL, w, x)
subroutine rd_implicit_update(dtfactor, qdt, qtC, psa, psb)
Implicit solve of psa=psb+dtfactor*dt*F_im(psa)
double precision, public, protected bzfn_delta
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
double precision, public, protected lor_r
Parameter for Lorenz system (Rayleigh number)
subroutine rd_add_source(qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
subroutine put_laplacians_onegrid(ixIL, ixOL, w)
subroutine rd_to_conserved(ixIL, ixOL, w, x)
double precision, public, protected br_a
Parameter for Brusselator model.
subroutine rd_laplacian(ixIL, ixOL, var, lpl)
Compute the Laplacian using a standard second order scheme. For now this method only works in slab ge...
double precision, public, protected lor_sigma
Parameter for Lorenz system (Prandtl number)
subroutine rd_get_flux(wC, w, x, ixIL, ixOL, idim, f)
double precision, public, protected br_c
Parameter for extended Brusselator model.
character(len=20), public, protected equation_name
Name of the system to be solved.
double precision, public, protected br_b
Parameter for Brusselator model.
The data structure that contains information about a tree node/grid block.