28 integer,
protected,
public ::
u_ = 1
29 integer,
protected,
public ::
v_ = 2
30 integer,
protected,
public ::
w_ = 3
36 double precision,
public,
protected ::
dtreacpar = 0.5d0
40 integer :: number_of_species = 2
41 integer :: equation_type = 1
42 integer,
parameter :: eq_gray_scott = 1
43 integer,
parameter :: eq_schnakenberg = 2
44 integer,
parameter :: eq_brusselator = 3
45 integer,
parameter :: eq_logistic = 4
46 integer,
parameter :: eq_analyt_hunds = 5
47 integer,
parameter :: eq_belousov_fn = 6
48 integer,
parameter :: eq_ext_brusselator = 7
49 integer,
parameter :: eq_lorenz = 8
50 integer,
parameter :: eq_no_reac = 9
53 double precision,
public,
protected ::
d1 = 0.05d0
55 double precision,
public,
protected ::
d2 = 1.0d0
57 double precision,
public,
protected ::
d3 = 1.0d0
63 double precision,
public,
protected ::
a1(^nd) = 0.0d0
65 double precision,
public,
protected ::
a2(^nd) = 0.0d0
67 double precision,
public,
protected ::
a3(^nd) = 0.0d0
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
75 double precision,
public,
protected ::
gs_f = 0.046d0
77 double precision,
public,
protected ::
gs_k = 0.063d0
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
86 double precision,
public,
protected ::
lg_lambda = 1.0d0
90 double precision,
public,
protected ::
bzfn_delta = 1.0d0
92 double precision,
public,
protected ::
bzfn_mu = 1.0d0
95 double precision,
public,
protected ::
lor_r = 28.0d0
97 double precision,
public,
protected ::
lor_sigma = 10.0d0
99 double precision,
public,
protected ::
lor_b = 8.0d0 / 3.0d0
102 logical :: ard_source_split = .false.
113 subroutine ard_params_read(files)
115 character(len=*),
intent(in) :: files(:)
118 namelist /ard_list/
d1,
d2,
d3,
adv_pow,
a1,
a2,
a3,
sb_alpha,
sb_beta,
sb_kappa, &
123 do n = 1,
size(files)
124 open(
unitpar, file=trim(files(n)), status=
'old')
125 read(
unitpar, ard_list,
end=111)
132 equation_type = eq_gray_scott
133 number_of_species = 2
134 case (
"schnakenberg")
135 equation_type = eq_schnakenberg
136 number_of_species = 2
138 equation_type = eq_brusselator
139 number_of_species = 2
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
153 equation_type = eq_lorenz
154 number_of_species = 3
156 equation_type = eq_no_reac
157 number_of_species = 1
159 call mpistop(
"Unknown equation_name (valid: gray-scott, schnakenberg, ...)")
162 end subroutine ard_params_read
167 integer,
intent(in) :: fh
168 integer,
parameter :: n_par = 0
169 integer,
dimension(MPI_STATUS_SIZE) :: st
173 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
190 allocate(start_indices(number_species),stop_indices(number_species))
194 u_ = var_set_fluxvar(
"u",
"u")
195 if (number_of_species >= 2)
then
196 v_ = var_set_fluxvar(
"v",
"v")
198 if (number_of_species >= 3)
then
199 w_ = var_set_fluxvar(
"w",
"w")
205 stop_indices(1)=nwflux
212 call mpistop(
"phys_check error: flux_type has wrong shape")
238 integer :: n, i, iw, species_list(number_of_species)
242 select case(number_of_species)
245 if (
d1 == 0.0d0)
then
246 write(*, *)
"Diffusion coefficient cannot be zero in IMEX scheme"
247 call mpistop(
"Zero diffusion in IMEX scheme")
250 species_list = [
u_,
v_]
251 if ((
d1 == 0.0d0) .or. (
d2 == 0.0d0))
then
252 write(*, *)
"Diffusion coefficient cannot be zero in IMEX scheme"
253 call mpistop(
"Zero diffusion in IMEX scheme")
256 species_list = [
u_,
v_,
w_]
257 if ((
d1 == 0.0d0) .or. (
d2 == 0.0d0) .or. (
d3 == 0.0d0))
then
258 write(*, *)
"Diffusion coefficient cannot be zero in IMEX scheme"
259 call mpistop(
"Zero diffusion in IMEX scheme")
263 do i = 1, number_of_species
275 ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
284 if (.not.
associated(
ard_mg_bc(i, n)%boundary_cond))
then
285 write(*,
"(A,I0,A,I0,A)")
"typeboundary(", iw,
",", n, &
286 ") is 'special', but the corresponding method " // &
287 "ard_mg_bc(i, n)%boundary_cond is not set"
288 call mpistop(
"ard_mg_bc(i, n)%boundary_cond not set")
291 write(*,*)
"ard_check_params warning: unknown boundary type"
292 ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
303 integer,
intent(in) :: ixI^L, ixO^L
304 double precision,
intent(inout) :: w(ixI^S, nw)
305 double precision,
intent(in) :: x(ixI^S, 1:^ND)
312 integer,
intent(in) :: ixI^L, ixO^L
313 double precision,
intent(inout) :: w(ixI^S, nw)
314 double precision,
intent(in) :: x(ixI^S, 1:^ND)
321 integer,
intent(in) :: ixI^L, ixO^L, idim
322 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
323 double precision,
intent(inout) :: cmax(ixI^S)
325 cmax(ixo^s) = abs(
a1(idim) * w(ixo^s,
u_)**(
adv_pow-1))
326 if (number_of_species >= 2)
then
327 cmax(ixo^s) = max(cmax(ixo^s), abs(
a2(idim) * w(ixo^s,
v_)**(
adv_pow-1)))
329 if (number_of_species >= 3)
then
330 cmax(ixo^s) = max(cmax(ixo^s), abs(
a3(idim) * w(ixo^s,
w_)**(
adv_pow-1)))
335 subroutine ard_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
338 integer,
intent(in) :: ixI^L, ixO^L, idim
339 double precision,
intent(in) :: wLC(ixI^S, nw), wRC(ixI^S,nw)
340 double precision,
intent(in) :: wLp(ixI^S, nw), wRp(ixI^S,nw)
341 double precision,
intent(in) :: x(ixI^S, 1:^ND)
342 double precision,
intent(in) :: Hspeed(ixI^S,1:number_species)
343 double precision,
intent(inout) :: cmax(ixI^S,1:number_species)
344 double precision,
intent(inout),
optional :: cmin(ixI^S,1:number_species)
346 double precision :: wmean(ixI^S,nw)
352 if (
present(cmin))
then
353 cmin(ixo^s,1) = min(
a1(idim) * wmean(ixo^s,
u_)**(
adv_pow-1), zero)
354 cmax(ixo^s,1) = max(
a1(idim) * wmean(ixo^s,
u_)**(
adv_pow-1), zero)
355 if (number_of_species >= 2)
then
356 cmin(ixo^s,1) = min(cmin(ixo^s,1),
a2(idim) * wmean(ixo^s,
v_)**(
adv_pow-1))
357 cmax(ixo^s,1) = max(cmax(ixo^s,1),
a2(idim) * wmean(ixo^s,
v_)**(
adv_pow-1))
359 if (number_of_species >= 3)
then
360 cmin(ixo^s,1) = min(cmin(ixo^s,1),
a3(idim) * wmean(ixo^s,
w_)**(
adv_pow-1))
361 cmax(ixo^s,1) = max(cmax(ixo^s,1),
a3(idim) * wmean(ixo^s,
w_)**(
adv_pow-1))
364 cmax(ixo^s,1) = maxval(abs(
a1(idim) * wmean(ixo^s,
u_)**(
adv_pow-1)))
365 if (number_of_species >=2)
then
366 cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(
a2(idim) * wmean(ixo^s,
v_)**(
adv_pow-1))))
368 if (number_of_species >=3)
then
369 cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(
a3(idim) * wmean(ixo^s,
w_)**(
adv_pow-1))))
377 integer,
intent(in) :: ixI^L, ixO^L
378 double precision,
intent(in) :: dx^D, x(ixI^S, 1:^ND)
379 double precision,
intent(in) :: w(ixI^S, 1:nw)
380 double precision,
intent(inout) :: dtnew
381 double precision :: maxrate
382 double precision :: maxD
383 double precision :: maxA
390 if (number_of_species >= 2)
then
393 if (number_of_species >= 3)
then
396 dtnew = min(dtnew,
dtdiffpar * minval([ dx^d ])**2 / (2 *
ndim * maxd))
399 select case (equation_type)
401 maxrate = max(maxval(w(ixo^s,
v_))**2 +
gs_f, &
403 case (eq_schnakenberg)
404 maxrate = max(maxval(abs(w(ixo^s,
v_) * w(ixo^s,
u_) - 1)), &
405 maxval(w(ixo^s,
u_))**2)
406 case (eq_brusselator)
407 maxrate = max( maxval(w(ixo^s,
u_)*w(ixo^s,
v_) - (
br_b+1)), &
408 maxval(w(ixo^s,
u_)**2) )
409 case (eq_ext_brusselator)
410 maxrate = max( maxval(w(ixo^s,
u_)*w(ixo^s,
v_) - (
br_b+1)) +
br_c, &
411 maxval(w(ixo^s,
u_)**2) )
412 maxrate = max(maxrate,
br_d)
415 case (eq_analyt_hunds)
416 maxrate = maxval(w(ixo^s,
u_)*abs(1 - w(ixo^s,
u_))) /
d1
417 case (eq_belousov_fn)
430 call mpistop(
"Unknown equation type")
440 integer,
intent(in) :: ixI^L, ixO^L, idim
441 double precision,
intent(in) :: wC(ixI^S, 1:nw)
442 double precision,
intent(in) :: w(ixI^S, 1:nw)
443 double precision,
intent(in) :: x(ixI^S, 1:^ND)
444 double precision,
intent(out) :: f(ixI^S, nwflux)
447 if (number_of_species >=2)
then
450 if (number_of_species >=3)
then
463 integer,
intent(in) :: ixI^L, ixO^L
464 double precision,
intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
465 double precision,
intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
470 subroutine ard_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
472 integer,
intent(in) :: ixI^L, ixO^L
473 double precision,
intent(in) :: qdt, dtfactor
474 double precision,
intent(in) :: wCT(ixI^S, 1:nw),wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
475 double precision,
intent(inout) :: w(ixI^S, 1:nw)
476 double precision :: lpl_u(ixO^S), lpl_v(ixO^S), lpl_w(ixO^S)
477 logical,
intent(in) :: qsourcesplit
478 logical,
intent(inout) :: active
481 if (qsourcesplit .eqv. ard_source_split)
then
484 if (number_of_species >= 2)
then
487 if (number_of_species >= 3)
then
497 select case (equation_type)
499 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u - &
500 wct(ixo^s,
u_) * wct(ixo^s,
v_)**2 + &
501 gs_f * (1 - wct(ixo^s,
u_)))
502 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v + &
503 wct(ixo^s,
u_) * wct(ixo^s,
v_)**2 - &
505 case (eq_schnakenberg)
506 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
508 wct(ixo^s,
u_)**2 * wct(ixo^s,
v_)))
509 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
511 case (eq_brusselator)
512 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
514 + wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
515 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
516 +
br_b * wct(ixo^s,
u_) - wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
517 case (eq_ext_brusselator)
518 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
520 + wct(ixo^s,
u_)**2 * wct(ixo^s,
v_) &
522 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
524 - wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
525 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
528 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
530 case (eq_analyt_hunds)
531 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
532 + 1.0d0/
d1 * w(ixo^s,
u_)**2 * (1 - w(ixo^s,
u_)))
533 case (eq_belousov_fn)
534 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
536 - wct(ixo^s,
u_)*wct(ixo^s,
w_) + wct(ixo^s,
u_) &
537 - wct(ixo^s,
u_)**2))
538 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
539 + wct(ixo^s,
u_) - wct(ixo^s,
v_))
540 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
545 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
548 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
549 +
lor_r * wct(ixo^s,
u_) - wct(ixo^s,
v_) &
550 - wct(ixo^s,
u_)*wct(ixo^s,
w_))
552 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
553 + wct(ixo^s,
u_)*wct(ixo^s,
v_) -
lor_b * wct(ixo^s,
w_))
555 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt *
d1 * lpl_u
557 call mpistop(
"Unknown equation type")
570 integer,
intent(in) :: ixI^L, ixO^L
571 double precision,
intent(in) :: var(ixI^S)
572 double precision,
intent(out) :: lpl(ixO^S)
573 integer :: idir, jxO^L, hxO^L
574 double precision :: h_inv2
579 hxo^l=ixo^l-
kr(idir,^
d);
580 jxo^l=ixo^l+
kr(idir,^
d);
582 lpl(ixo^s) = lpl(ixo^s) + h_inv2 * &
583 (var(jxo^s) - 2 * var(ixo^s) + var(hxo^s))
586 call mpistop(
"ard_laplacian not implemented in this geometry")
593 integer,
intent(in) :: ixI^L, ixO^L
594 double precision,
intent(inout) :: w(ixI^S, 1:nw)
596 double precision :: lpl_u(ixO^S), lpl_v(ixO^S), lpl_w(ixO^S)
599 if (number_of_species >= 2)
then
602 if (number_of_species >= 3)
then
606 w(ixo^s,
u_) =
d1*lpl_u
607 if (number_of_species >= 2)
then
608 w(ixo^s,
v_) =
d2*lpl_v
610 if (number_of_species >= 3)
then
611 w(ixo^s,
w_) =
d3*lpl_w
619 type(state),
target :: psa(max_blocks)
620 double precision,
intent(in) :: qtC
622 integer :: iigrid, igrid, level
628 do iigrid=1,igridstail; igrid=igrids(iigrid);
641 type(state),
target :: psa(max_blocks)
642 type(state),
target :: psb(max_blocks)
643 double precision,
intent(in) :: qdt
644 double precision,
intent(in) :: qtC
645 double precision,
intent(in) :: dtfactor
648 double precision :: res, max_residual, lambda
650 integer :: iw_to,iw_from
651 integer :: iigrid, igrid, id
658 if (qdt <
dtmin)
then
660 print *,
'skipping implicit solve: dt too small!'
665 max_residual = 1
d-7/qdt
667 mg%operator_type = mg_helmholtz
668 call mg_set_methods(mg)
670 if (.not. mg%is_allocated)
call mpistop(
"multigrid tree not allocated yet")
673 lambda = 1/(dtfactor * qdt *
d1)
674 call helmholtz_set_lambda(lambda)
677 call mg_copy_to_tree(
u_, mg_irhs, factor=-lambda, state_from=psb)
678 call mg_copy_to_tree(
u_, mg_iphi, state_from=psb)
680 call mg_fas_fmg(mg, .true., max_res=res)
682 call mg_fas_vcycle(mg, max_res=res)
683 if (res < max_residual)
exit
686 call mg_copy_from_tree_gc(mg_iphi,
u_, state_to=psa)
690 if (number_of_species >= 2)
then
691 lambda = 1/(dtfactor * qdt *
d2)
692 call helmholtz_set_lambda(lambda)
695 call mg_copy_to_tree(
v_, mg_irhs, factor=-lambda, state_from=psb)
696 call mg_copy_to_tree(
v_, mg_iphi, state_from=psb)
698 call mg_fas_fmg(mg, .true., max_res=res)
700 call mg_fas_vcycle(mg, max_res=res)
701 if (res < max_residual)
exit
704 call mg_copy_from_tree_gc(mg_iphi,
v_, state_to=psa)
709 if (number_of_species >= 3)
then
710 lambda = 1/(dtfactor * qdt *
d3)
711 call helmholtz_set_lambda(lambda)
713 call mg_copy_to_tree(
w_, mg_irhs, factor=-lambda, state_from=psb)
714 call mg_copy_to_tree(
w_, mg_iphi, state_from=psb)
716 call mg_fas_fmg(mg, .true., max_res=res)
718 call mg_fas_vcycle(mg, max_res=res)
719 if (res < max_residual)
exit
722 call mg_copy_from_tree_gc(mg_iphi,
w_, state_to=psa)
Module containing the physics routines for advection-reaction-diffusion equations.
subroutine ard_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, w, x)
subroutine ard_get_flux(wC, w, x, ixIL, ixOL, idim, f)
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.
subroutine ard_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
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
subroutine ard_get_cmax(w, x, ixIL, ixOL, idim, cmax)
double precision, public, protected br_c
subroutine, public ard_phys_init()
subroutine ard_get_cbounds(wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
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)
subroutine ard_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 d1
Diffusion coefficient for first species (u)
subroutine ard_add_source(qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
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)
subroutine put_laplacians_onegrid(ixIL, ixOL, w)
double precision, public, protected d3
Diffusion coefficient for third species (w) (if applicable)
double precision, public, protected bzfn_lambda
subroutine ard_write_info(fh)
Write this modules parameters to a snapshot.
integer, public, protected w_
For 3 or more equations.
subroutine ard_implicit_update(dtfactor, qdt, qtC, psa, psb)
Implicit solve of psa=psb+dtfactor*dt*F_im(psa)
subroutine ard_check_params
double precision, dimension(^nd), public, protected a1
Advection coefficients for first species (u)
subroutine ard_to_primitive(ixIL, ixOL, w, x)
double precision, public, protected bzfn_delta
subroutine ard_evaluate_implicit(qtC, psa)
inplace update of psa==>F_im(psa)
subroutine ard_to_conserved(ixIL, ixOL, w, x)
double precision, public, protected br_d
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.
integer, 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
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, parameter bc_cont
integer, parameter bc_symm
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.
double precision, dimension(ndim) dxlevel
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_add_source_geom), pointer phys_add_source_geom
procedure(sub_check_params), pointer phys_check_params
integer, parameter flux_default
Indicates a normal flux.
integer, dimension(:, :), allocatable flux_type
Array per direction per variable, which can be used to specify that certain fluxes have to be treated...
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.
integer nwflux
Number of flux variables.
The data structure that contains information about a tree node/grid block.