MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
|
Reaction-diffusion module (physics routines) More...
Functions/Subroutines | |
subroutine, public | rd_phys_init () |
Variables | |
integer, public, protected | u_ = 1 |
integer, public, protected | v_ = 2 |
For 2 or more equations. | |
integer, public, protected | w_ = 3 |
For 3 or more equations. | |
logical, public, protected | rd_particles = .false. |
Whether particles module is added. | |
double precision, public, protected | dtreacpar = 0.5d0 |
Parameter with which to multiply the reaction timestep restriction. | |
character(len=20), public, protected | equation_name = "gray-scott" |
Name of the system to be solved. | |
double precision, public, protected | d1 = 0.05d0 |
Diffusion coefficient for first species (u) | |
double precision, public, protected | d2 = 1.0d0 |
Diffusion coefficient for second species (v) (if applicable) | |
double precision, public, protected | d3 = 1.0d0 |
Diffusion coefficient for third species (w) (if applicable) | |
double precision, public, protected | sb_alpha = 0.1305d0 |
Parameter for Schnakenberg model. | |
double precision, public, protected | sb_beta = 0.7695d0 |
Parameter for Schnakenberg model. | |
double precision, public, protected | sb_kappa = 100.0d0 |
Parameter for Schnakenberg model. | |
double precision, public, protected | gs_f = 0.046d0 |
Feed rate for Gray-Scott model. | |
double precision, public, protected | gs_k = 0.063d0 |
Kill rate for Gray-Scott model. | |
double precision, public, protected | br_a = 4.5d0 |
Parameter for Brusselator model. | |
double precision, public, protected | br_b = 8.0d0 |
Parameter for Brusselator model. | |
double precision, public, protected | br_c = 1.0d0 |
Parameter for extended Brusselator model. | |
double precision, public, protected | br_d = 1.0d0 |
Parameter for extended Brusselator model. | |
double precision, public, protected | lg_lambda = 1.0d0 |
Parameter for logistic model (Fisher / KPP equation) | |
double precision, public, protected | bzfn_epsilon = 1.0d0 |
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction. | |
double precision, public, protected | bzfn_delta = 1.0d0 |
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction. | |
double precision, public, protected | bzfn_lambda = 1.0d0 |
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction. | |
double precision, public, protected | bzfn_mu = 1.0d0 |
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction. | |
double precision, public, protected | lor_r = 28.0d0 |
Parameter for Lorenz system (Rayleigh number) | |
double precision, public, protected | lor_sigma = 10.0d0 |
Parameter for Lorenz system (Prandtl number) | |
double precision, public, protected | lor_b = 8.0d0 / 3.0d0 |
Parameter for Lorenz system (aspect ratio of the convection rolls) | |
type(mg_bc_t), dimension(3, mg_num_neighbors), public | rd_mg_bc |
Boundary condition information for the multigrid method. | |
Reaction-diffusion module (physics routines)
Multiple reaction-diffusion systems are included: the Gray-Scott model, the Schnakenberg model, the Brusselator model, the diffusive logistic equation, an analytical testcase from "Numerical solution of time-dependent advection- diffusion-reaction equations" by Hundsdorfer & Verwer, the Oregonator model, the extended Brusselator model, and the diffusive Lorenz system. See the documentation of the reaction-diffusion module for more information.
IMEX methods are also supported. The implicit system is solved by a multigrid solver coupled into MPI-AMRVAC.
subroutine, public mod_rd_phys::rd_phys_init |
double precision, public, protected mod_rd_phys::br_a = 4.5d0 |
Parameter for Brusselator model.
Definition at line 64 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::br_b = 8.0d0 |
Parameter for Brusselator model.
Definition at line 66 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::br_c = 1.0d0 |
Parameter for extended Brusselator model.
Definition at line 68 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::br_d = 1.0d0 |
Parameter for extended Brusselator model.
Definition at line 70 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::bzfn_delta = 1.0d0 |
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
Definition at line 78 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::bzfn_epsilon = 1.0d0 |
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
Definition at line 76 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::bzfn_lambda = 1.0d0 |
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
Definition at line 80 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::bzfn_mu = 1.0d0 |
Parameter for the Field-Noyes model of the Belousov-Zhabotinsky reaction.
Definition at line 82 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::d1 = 0.05d0 |
Diffusion coefficient for first species (u)
Definition at line 45 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::d2 = 1.0d0 |
Diffusion coefficient for second species (v) (if applicable)
Definition at line 47 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::d3 = 1.0d0 |
Diffusion coefficient for third species (w) (if applicable)
Definition at line 49 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::dtreacpar = 0.5d0 |
Parameter with which to multiply the reaction timestep restriction.
Definition at line 29 of file mod_rd_phys.t.
character(len=20), public, protected mod_rd_phys::equation_name = "gray-scott" |
Name of the system to be solved.
Definition at line 32 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::gs_f = 0.046d0 |
Feed rate for Gray-Scott model.
Definition at line 59 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::gs_k = 0.063d0 |
Kill rate for Gray-Scott model.
Definition at line 61 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::lg_lambda = 1.0d0 |
Parameter for logistic model (Fisher / KPP equation)
Definition at line 73 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::lor_b = 8.0d0 / 3.0d0 |
Parameter for Lorenz system (aspect ratio of the convection rolls)
Definition at line 89 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::lor_r = 28.0d0 |
Parameter for Lorenz system (Rayleigh number)
Definition at line 85 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::lor_sigma = 10.0d0 |
Parameter for Lorenz system (Prandtl number)
Definition at line 87 of file mod_rd_phys.t.
type(mg_bc_t), dimension(3, mg_num_neighbors), public mod_rd_phys::rd_mg_bc |
Boundary condition information for the multigrid method.
Definition at line 95 of file mod_rd_phys.t.
logical, public, protected mod_rd_phys::rd_particles = .false. |
Whether particles module is added.
Definition at line 26 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::sb_alpha = 0.1305d0 |
Parameter for Schnakenberg model.
Definition at line 52 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::sb_beta = 0.7695d0 |
Parameter for Schnakenberg model.
Definition at line 54 of file mod_rd_phys.t.
double precision, public, protected mod_rd_phys::sb_kappa = 100.0d0 |
Parameter for Schnakenberg model.
Definition at line 56 of file mod_rd_phys.t.
integer, public, protected mod_rd_phys::u_ = 1 |
Definition at line 21 of file mod_rd_phys.t.
integer, public, protected mod_rd_phys::v_ = 2 |
For 2 or more equations.
Definition at line 22 of file mod_rd_phys.t.
integer, public, protected mod_rd_phys::w_ = 3 |
For 3 or more equations.
Definition at line 23 of file mod_rd_phys.t.