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
165 subroutine ard_write_info(fh)
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)
174 end subroutine ard_write_info
187 allocate(start_indices(number_species),stop_indices(number_species))
191 u_ = var_set_fluxvar(
"u",
"u")
192 if (number_of_species >= 2)
then
193 v_ = var_set_fluxvar(
"v",
"v")
195 if (number_of_species >= 3)
then
196 w_ = var_set_fluxvar(
"w",
"w")
202 stop_indices(1)=nwflux
209 call mpistop(
"phys_check error: flux_type has wrong shape")
228 subroutine ard_check_params
233 use mod_particles,
only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
234 integer :: n, i, iw, species_list(number_of_species)
243 select case(number_of_species)
246 if (
d1 == 0.0d0)
then
247 write(*, *)
"Diffusion coefficient cannot be zero in IMEX scheme"
248 call mpistop(
"Zero diffusion in IMEX scheme")
251 species_list = [
u_,
v_]
252 if ((
d1 == 0.0d0) .or. (
d2 == 0.0d0))
then
253 write(*, *)
"Diffusion coefficient cannot be zero in IMEX scheme"
254 call mpistop(
"Zero diffusion in IMEX scheme")
257 species_list = [
u_,
v_,
w_]
258 if ((
d1 == 0.0d0) .or. (
d2 == 0.0d0) .or. (
d3 == 0.0d0))
then
259 write(*, *)
"Diffusion coefficient cannot be zero in IMEX scheme"
260 call mpistop(
"Zero diffusion in IMEX scheme")
264 do i = 1, number_of_species
276 ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
286 call mpistop(
"special BC for MG not defined")
297 write(*,*)
"ard_check_params warning: unknown boundary type"
298 ard_mg_bc(i, n)%bc_type = mg_bc_dirichlet
307 write(*,*)
'====ARD run with settings===================='
308 write(*,*)
'Using mod_ard_phys with settings:'
309 write(*,*)
'Dimensionality :',
ndim
310 write(*,*)
'vector components:',
ndir
312 write(*,*)
'number of variables nw=',nw
313 write(*,*)
' start index iwstart=',iwstart
314 write(*,*)
'number of vector variables=',nvector
315 write(*,*)
'number of stagger variables nws=',nws
316 write(*,*)
'number of variables with BCs=',nwgc
317 write(*,*)
'number of vars with fluxes=',nwflux
318 write(*,*)
'number of vars with flux + BC=',nwfluxbc
319 write(*,*)
'number of auxiliary variables=',nwaux
320 write(*,*)
'number of extra vars without flux=',nwextra
321 write(*,*)
'number of extra vars for wextra=',nw_extra
322 write(*,*)
'number of auxiliary I/O variables=',
nwauxio
325 write(*,*)
'*****Using particles: npayload,ngridvars :', npayload,ngridvars
326 write(*,*)
'*****Using particles: nusrpayload :', nusrpayload
327 write(*,*)
'*****Using particles: num_particles :', num_particles
328 write(*,*)
'*****Using particles: physics_type_particles=',physics_type_particles
331 write(*,*)
'==========================================='
335 end subroutine ard_check_params
337 subroutine ard_to_conserved(ixI^L, ixO^L, w, x)
339 integer,
intent(in) :: ixi^
l, ixo^
l
340 double precision,
intent(inout) :: w(ixi^s, nw)
341 double precision,
intent(in) :: x(ixi^s, 1:^nd)
344 end subroutine ard_to_conserved
346 subroutine ard_to_primitive(ixI^L, ixO^L, w, x)
348 integer,
intent(in) :: ixi^
l, ixo^
l
349 double precision,
intent(inout) :: w(ixi^s, nw)
350 double precision,
intent(in) :: x(ixi^s, 1:^nd)
353 end subroutine ard_to_primitive
355 subroutine ard_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
357 integer,
intent(in) :: ixi^
l, ixo^
l, idim
358 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
359 double precision,
intent(inout) :: cmax(ixi^s)
361 cmax(ixo^s) = abs(
a1(idim) * w(ixo^s,
u_)**(
adv_pow-1))
362 if (number_of_species >= 2)
then
363 cmax(ixo^s) = max(cmax(ixo^s), abs(
a2(idim) * w(ixo^s,
v_)**(
adv_pow-1)))
365 if (number_of_species >= 3)
then
366 cmax(ixo^s) = max(cmax(ixo^s), abs(
a3(idim) * w(ixo^s,
w_)**(
adv_pow-1)))
369 end subroutine ard_get_cmax
371 subroutine ard_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
374 integer,
intent(in) :: ixi^
l, ixo^
l, idim
375 double precision,
intent(in) :: wlc(ixi^s,
nw), wrc(ixi^s,
nw)
376 double precision,
intent(in) :: wlp(ixi^s,
nw), wrp(ixi^s,
nw)
377 double precision,
intent(in) :: x(ixi^s, 1:^nd)
380 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
382 double precision :: wmean(ixi^s,
nw)
388 if (
present(cmin))
then
389 cmin(ixo^s,1) = min(
a1(idim) * wmean(ixo^s,
u_)**(
adv_pow-1), zero)
390 cmax(ixo^s,1) = max(
a1(idim) * wmean(ixo^s,
u_)**(
adv_pow-1), zero)
391 if (number_of_species >= 2)
then
392 cmin(ixo^s,1) = min(cmin(ixo^s,1),
a2(idim) * wmean(ixo^s,
v_)**(
adv_pow-1))
393 cmax(ixo^s,1) = max(cmax(ixo^s,1),
a2(idim) * wmean(ixo^s,
v_)**(
adv_pow-1))
395 if (number_of_species >= 3)
then
396 cmin(ixo^s,1) = min(cmin(ixo^s,1),
a3(idim) * wmean(ixo^s,
w_)**(
adv_pow-1))
397 cmax(ixo^s,1) = max(cmax(ixo^s,1),
a3(idim) * wmean(ixo^s,
w_)**(
adv_pow-1))
400 cmax(ixo^s,1) = maxval(abs(
a1(idim) * wmean(ixo^s,
u_)**(
adv_pow-1)))
401 if (number_of_species >=2)
then
402 cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(
a2(idim) * wmean(ixo^s,
v_)**(
adv_pow-1))))
404 if (number_of_species >=3)
then
405 cmax(ixo^s,1) = max(cmax(ixo^s,1), maxval(abs(
a3(idim) * wmean(ixo^s,
w_)**(
adv_pow-1))))
409 end subroutine ard_get_cbounds
411 subroutine ard_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
413 integer,
intent(in) :: ixi^
l, ixo^
l
414 double precision,
intent(in) ::
dx^
d, x(ixi^s, 1:^nd)
415 double precision,
intent(in) :: w(ixi^s, 1:nw)
416 double precision,
intent(inout) :: dtnew
417 double precision :: maxrate
418 double precision :: maxd
419 double precision :: maxa
426 if (number_of_species >= 2)
then
429 if (number_of_species >= 3)
then
435 select case (equation_type)
437 maxrate = max(maxval(w(ixo^s,
v_))**2 +
gs_f, &
439 case (eq_schnakenberg)
440 maxrate = max(maxval(abs(w(ixo^s,
v_) * w(ixo^s,
u_) - 1)), &
441 maxval(w(ixo^s,
u_))**2)
442 case (eq_brusselator)
443 maxrate = max( maxval(w(ixo^s,
u_)*w(ixo^s,
v_) - (
br_b+1)), &
444 maxval(w(ixo^s,
u_)**2) )
445 case (eq_ext_brusselator)
446 maxrate = max( maxval(w(ixo^s,
u_)*w(ixo^s,
v_) - (
br_b+1)) +
br_c, &
447 maxval(w(ixo^s,
u_)**2) )
448 maxrate = max(maxrate,
br_d)
451 case (eq_analyt_hunds)
452 maxrate = maxval(w(ixo^s,
u_)*abs(1 - w(ixo^s,
u_))) /
d1
453 case (eq_belousov_fn)
466 call mpistop(
"Unknown equation type")
471 end subroutine ard_get_dt
474 subroutine ard_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
476 integer,
intent(in) :: ixi^
l, ixo^
l, idim
477 double precision,
intent(in) :: wc(ixi^s, 1:nw)
478 double precision,
intent(in) :: w(ixi^s, 1:nw)
479 double precision,
intent(in) :: x(ixi^s, 1:^nd)
480 double precision,
intent(out) :: f(ixi^s, nwflux)
483 if (number_of_species >=2)
then
486 if (number_of_species >=3)
then
490 end subroutine ard_get_flux
492 subroutine ard_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim,w, x)
499 integer,
intent(in) :: ixi^
l, ixo^
l
500 double precision,
intent(in) :: qdt, dtfactor, x(ixi^s, 1:^nd)
501 double precision,
intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
503 end subroutine ard_add_source_geom
506 subroutine ard_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
508 integer,
intent(in) :: ixi^
l, ixo^
l
509 double precision,
intent(in) :: qdt, dtfactor
510 double precision,
intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:
ndim)
511 double precision,
intent(inout) :: w(ixi^s, 1:nw)
512 double precision :: lpl_u(ixo^s), lpl_v(ixo^s), lpl_w(ixo^s)
513 logical,
intent(in) :: qsourcesplit
514 logical,
intent(inout) :: active
517 if (qsourcesplit .eqv. ard_source_split)
then
519 call ard_laplacian(ixi^
l, ixo^
l, wct(ixi^s,
u_), lpl_u)
520 if (number_of_species >= 2)
then
521 call ard_laplacian(ixi^
l, ixo^
l, wct(ixi^s,
v_), lpl_v)
523 if (number_of_species >= 3)
then
524 call ard_laplacian(ixi^
l, ixo^
l, wct(ixi^s,
w_), lpl_w)
533 select case (equation_type)
535 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u - &
536 wct(ixo^s,
u_) * wct(ixo^s,
v_)**2 + &
537 gs_f * (1 - wct(ixo^s,
u_)))
538 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v + &
539 wct(ixo^s,
u_) * wct(ixo^s,
v_)**2 - &
541 case (eq_schnakenberg)
542 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
544 wct(ixo^s,
u_)**2 * wct(ixo^s,
v_)))
545 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
547 case (eq_brusselator)
548 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
550 + wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
551 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
552 +
br_b * wct(ixo^s,
u_) - wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
553 case (eq_ext_brusselator)
554 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
556 + wct(ixo^s,
u_)**2 * wct(ixo^s,
v_) &
558 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
560 - wct(ixo^s,
u_)**2 * wct(ixo^s,
v_))
561 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
564 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
566 case (eq_analyt_hunds)
567 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
568 + 1.0d0/
d1 * w(ixo^s,
u_)**2 * (1 - w(ixo^s,
u_)))
569 case (eq_belousov_fn)
570 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
572 - wct(ixo^s,
u_)*wct(ixo^s,
w_) + wct(ixo^s,
u_) &
573 - wct(ixo^s,
u_)**2))
574 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
575 + wct(ixo^s,
u_) - wct(ixo^s,
v_))
576 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
581 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt * (
d1 * lpl_u &
584 w(ixo^s,
v_) = w(ixo^s,
v_) + qdt * (
d2 * lpl_v &
585 +
lor_r * wct(ixo^s,
u_) - wct(ixo^s,
v_) &
586 - wct(ixo^s,
u_)*wct(ixo^s,
w_))
588 w(ixo^s,
w_) = w(ixo^s,
w_) + qdt * (
d3 * lpl_w &
589 + wct(ixo^s,
u_)*wct(ixo^s,
v_) -
lor_b * wct(ixo^s,
w_))
591 w(ixo^s,
u_) = w(ixo^s,
u_) + qdt *
d1 * lpl_u
593 call mpistop(
"Unknown equation type")
600 end subroutine ard_add_source
604 subroutine ard_laplacian(ixI^L,ixO^L,var,lpl)
606 integer,
intent(in) :: ixi^
l, ixo^
l
607 double precision,
intent(in) :: var(ixi^s)
608 double precision,
intent(out) :: lpl(ixo^s)
609 integer :: idir, jxo^
l, hxo^
l
610 double precision :: h_inv2
615 hxo^
l=ixo^
l-
kr(idir,^
d);
616 jxo^
l=ixo^
l+
kr(idir,^
d);
618 lpl(ixo^s) = lpl(ixo^s) + h_inv2 * &
619 (var(jxo^s) - 2 * var(ixo^s) + var(hxo^s))
622 call mpistop(
"ard_laplacian not implemented in this geometry")
625 end subroutine ard_laplacian
627 subroutine put_laplacians_onegrid(ixI^L,ixO^L,w)
629 integer,
intent(in) :: ixi^
l, ixo^
l
630 double precision,
intent(inout) :: w(ixi^s, 1:nw)
632 double precision :: lpl_u(ixo^s), lpl_v(ixo^s), lpl_w(ixo^s)
634 call ard_laplacian(ixi^
l, ixo^
l, w(ixi^s,
u_), lpl_u)
635 if (number_of_species >= 2)
then
636 call ard_laplacian(ixi^
l, ixo^
l, w(ixi^s,
v_), lpl_v)
638 if (number_of_species >= 3)
then
639 call ard_laplacian(ixi^
l, ixo^
l, w(ixi^s,
w_), lpl_w)
642 w(ixo^s,
u_) =
d1*lpl_u
643 if (number_of_species >= 2)
then
644 w(ixo^s,
v_) =
d2*lpl_v
646 if (number_of_species >= 3)
then
647 w(ixo^s,
w_) =
d3*lpl_w
650 end subroutine put_laplacians_onegrid
653 subroutine ard_evaluate_implicit(qtC,psa)
656 double precision,
intent(in) :: qtc
658 integer :: iigrid, igrid, level
664 do iigrid=1,igridstail; igrid=igrids(iigrid);
666 call put_laplacians_onegrid(ixg^
ll,ixo^
l,psa(igrid)%w)
670 end subroutine ard_evaluate_implicit
673 subroutine ard_implicit_update(dtfactor,qdt,qtC,psa,psb)
679 double precision,
intent(in) :: qdt
680 double precision,
intent(in) :: qtc
681 double precision,
intent(in) :: dtfactor
684 double precision :: res, max_residual, lambda
686 integer :: iw_to,iw_from
687 integer :: iigrid, igrid, id
694 if (qdt <
dtmin)
then
696 print *,
'skipping implicit solve: dt too small!'
701 max_residual = 1
d-7/qdt
703 mg%operator_type = mg_helmholtz
704 call mg_set_methods(mg)
706 if (.not. mg%is_allocated)
call mpistop(
"multigrid tree not allocated yet")
709 lambda = 1/(dtfactor * qdt *
d1)
710 call helmholtz_set_lambda(lambda)
713 call mg_copy_to_tree(
u_, mg_irhs, factor=-lambda, state_from=psb)
714 call mg_copy_to_tree(
u_, 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,
u_, state_to=psa)
726 if (number_of_species >= 2)
then
727 lambda = 1/(dtfactor * qdt *
d2)
728 call helmholtz_set_lambda(lambda)
731 call mg_copy_to_tree(
v_, mg_irhs, factor=-lambda, state_from=psb)
732 call mg_copy_to_tree(
v_, mg_iphi, state_from=psb)
734 call mg_fas_fmg(mg, .true., max_res=res)
736 call mg_fas_vcycle(mg, max_res=res)
737 if (res < max_residual)
exit
740 call mg_copy_from_tree_gc(mg_iphi,
v_, state_to=psa)
745 if (number_of_species >= 3)
then
746 lambda = 1/(dtfactor * qdt *
d3)
747 call helmholtz_set_lambda(lambda)
749 call mg_copy_to_tree(
w_, mg_irhs, factor=-lambda, state_from=psb)
750 call mg_copy_to_tree(
w_, mg_iphi, state_from=psb)
752 call mg_fas_fmg(mg, .true., max_res=res)
754 call mg_fas_vcycle(mg, max_res=res)
755 if (res < max_residual)
exit
758 call mg_copy_from_tree_gc(mg_iphi,
w_, state_to=psa)
762 end subroutine ard_implicit_update
Module containing the physics routines for advection-reaction-diffusion equations.
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.
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
double precision, public, protected br_c
subroutine, public ard_phys_init()
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)
double precision, public, protected d1
Diffusion coefficient for first species (u)
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)
double precision, public, protected d3
Diffusion coefficient for third species (w) (if applicable)
double precision, public, protected bzfn_lambda
integer, public, protected w_
For 3 or more equations.
double precision, dimension(^nd), public, protected a1
Advection coefficients for first species (u)
double precision, public, protected bzfn_delta
double precision, public, protected br_d
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
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
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...
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_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.
Module with all the methods that users can customize in AMRVAC.
procedure(special_mg_bc), pointer usr_special_mg_bc
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
integer nwflux
Number of flux variables.
The data structure that contains information about a tree node/grid block.