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
150 subroutine rd_write_info(fh)
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)
159 end subroutine rd_write_info
172 allocate(start_indices(number_species),stop_indices(number_species))
176 u_ = var_set_fluxvar(
"u",
"u")
177 if (number_of_species >= 2)
then
178 v_ = var_set_fluxvar(
"v",
"v")
180 if (number_of_species >= 3)
then
181 w_ = var_set_fluxvar(
"w",
"w")
188 stop_indices(1)=nwflux
208 subroutine rd_check_params
212 use mod_particles,
only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
214 integer :: n, i, iw, species_list(number_of_species)
218 call mpistop(
"mod_rd requires flux_scheme = source")
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_]
234 do i = 1, number_of_species
242 rd_mg_bc(i, n)%bc_type = mg_bc_neumann
246 rd_mg_bc(i, n)%bc_type = mg_bc_dirichlet
250 rd_mg_bc(i, n)%bc_type = mg_bc_neumann
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")
262 write(*,*)
"rd_check_params warning: unknown boundary type"
263 rd_mg_bc(i, n)%bc_type = mg_bc_dirichlet
272 write(*,*)
'====RD run with settings===================='
273 write(*,*)
'Using mod_rd_phys with settings:'
274 write(*,*)
'Dimensionality :',
ndim
275 write(*,*)
'vector components:',
ndir
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
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
296 write(*,*)
'==========================================='
301 end subroutine rd_check_params
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)
310 end subroutine rd_to_conserved
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)
319 end subroutine rd_to_primitive
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)
328 end subroutine rd_get_cmax
330 subroutine rd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
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)
338 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
341 if (
present(cmin))
then
342 cmin(ixo^s,1) = 0.0d0
343 cmax(ixo^s,1) = 0.0d0
345 cmax(ixo^s,1) = 0.0d0
348 end subroutine rd_get_cbounds
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
362 if (number_of_species >= 2)
then
365 if (number_of_species >= 3)
then
371 select case (equation_type)
373 maxrate = max(maxval(w(ixo^s,
v_))**2 +
gs_f, &
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)
387 case (eq_analyt_hunds)
388 maxrate = maxval(w(ixo^s,
u_)*abs(1 - w(ixo^s,
u_))) /
d1
389 case (eq_belousov_fn)
399 call mpistop(
"Unknown equation type")
404 end subroutine rd_get_dt
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)
416 end subroutine rd_get_flux
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
430 if (qsourcesplit .eqv. rd_source_split)
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)
436 if (number_of_species >= 3)
then
437 call rd_laplacian(ixi^
l, ixo^
l, wct(ixi^s,
w_), lpl_w)
446 select case (equation_type)
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 - &
454 case (eq_schnakenberg)
455 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
457 wct(ixo^s,
u_)**2 * wct(ixo^s,
v_)))
458 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
460 case (eq_brusselator)
461 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_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 &
469 + wct(ixo^s,
u_)**2 * wct(ixo^s,
v_) &
471 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
473 - wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
474 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
477 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_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 &
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 &
494 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
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_))
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_))
504 call mpistop(
"Unknown equation type")
511 end subroutine rd_add_source
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
526 hxo^
l=ixo^
l-
kr(idir,^
d);
527 jxo^
l=ixo^
l+
kr(idir,^
d);
529 lpl(ixo^s) = lpl(ixo^s) + h_inv2 * &
530 (var(jxo^s) - 2 * var(ixo^s) + var(hxo^s))
533 call mpistop(
"rd_laplacian not implemented in this geometry")
536 end subroutine rd_laplacian
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)
543 double precision :: lpl_u(ixo^s), lpl_v(ixo^s), lpl_w(ixo^s)
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)
549 if (number_of_species >= 3)
then
550 call rd_laplacian(ixi^
l, ixo^
l, w(ixi^s,
w_), lpl_w)
553 w(ixo^s,
u_) =
d1*lpl_u
554 if (number_of_species >= 2)
then
555 w(ixo^s,
v_) =
d2*lpl_v
557 if (number_of_species >= 3)
then
558 w(ixo^s,
w_) =
d3*lpl_w
561 end subroutine put_laplacians_onegrid
564 subroutine rd_evaluate_implicit(qtC,psa)
567 double precision,
intent(in) :: qtc
569 integer :: iigrid, igrid, level
575 do iigrid=1,igridstail; igrid=igrids(iigrid);
577 call put_laplacians_onegrid(ixg^
ll,ixo^
l,psa(igrid)%w)
581 end subroutine rd_evaluate_implicit
584 subroutine rd_implicit_update(dtfactor,qdt,qtC,psa,psb)
590 double precision,
intent(in) :: qdt
591 double precision,
intent(in) :: qtc
592 double precision,
intent(in) :: dtfactor
595 double precision :: res, max_residual, lambda
597 integer :: iw_to,iw_from
598 integer :: iigrid, igrid, id
605 if (qdt <
dtmin)
then
607 print *,
'skipping implicit solve: dt too small!'
612 max_residual = 1
d-7/qdt
614 mg%operator_type = mg_helmholtz
615 call mg_set_methods(mg)
617 if (.not. mg%is_allocated)
call mpistop(
"multigrid tree not allocated yet")
620 lambda = 1/(dtfactor * qdt *
d1)
621 call helmholtz_set_lambda(lambda)
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)
627 call mg_fas_fmg(mg, .true., max_res=res)
629 call mg_fas_vcycle(mg, max_res=res)
630 if (res < max_residual)
exit
633 if (res > max_residual)
call mpistop(
"u_multigrid: no convergence")
635 call mg_copy_from_tree_gc(mg_iphi,
u_, state_to=psa)
639 if (number_of_species >= 2)
then
640 lambda = 1/(dtfactor * qdt *
d2)
641 call helmholtz_set_lambda(lambda)
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)
647 call mg_fas_fmg(mg, .true., max_res=res)
649 call mg_fas_vcycle(mg, max_res=res)
650 if (res < max_residual)
exit
653 if (res > max_residual)
call mpistop(
"v_multigrid: no convergence")
654 call mg_copy_from_tree_gc(mg_iphi,
v_, state_to=psa)
659 if (number_of_species >= 3)
then
660 lambda = 1/(dtfactor * qdt *
d3)
661 call helmholtz_set_lambda(lambda)
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)
666 call mg_fas_fmg(mg, .true., max_res=res)
668 call mg_fas_vcycle(mg, max_res=res)
669 if (res < max_residual)
exit
672 if (res > max_residual)
call mpistop(
"w_multigrid: no convergence")
673 call mg_copy_from_tree_gc(mg_iphi,
w_, state_to=psa)
677 end subroutine rd_implicit_update
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with basic grid data structures.
Module with geometry-related routines (e.g., divergence, curl)
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 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.
integer, parameter bc_cont
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
integer nghostcells
Number of ghost cells surrounding a grid.
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 max_blocks
The maximum number of grid blocks in a processor.
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
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)
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()
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)
double precision, public, protected dtreacpar
Parameter with which to multiply the reaction timestep restriction.
integer, public, protected w_
For 3 or more equations.
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.
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.
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)
double precision, public, protected br_a
Parameter for Brusselator model.
double precision, public, protected lor_sigma
Parameter for Lorenz system (Prandtl number)
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.
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.