8 integer,
protected,
public ::
rho_ = 1
27 subroutine nonlinear_params_read(files)
29 character(len=*),
intent(in) :: files(:)
35 open(
unitpar, file=trim(files(n)), status=
'old')
36 read(
unitpar, nonlinear_list,
end=111)
40 end subroutine nonlinear_params_read
51 integer,
intent(in) :: fh
52 integer,
dimension(MPI_STATUS_SIZE) :: st
56 call mpi_file_write(fh, 0, 1, mpi_integer, st, er)
74 allocate(start_indices(number_species),stop_indices(number_species))
83 stop_indices(1)=nwflux
93 call mpistop(
"phys_check error: flux_type has wrong shape")
118 integer,
intent(in) :: ixI^L, ixO^L
119 double precision,
intent(inout) :: w(ixI^S, nw)
120 double precision,
intent(in) :: x(ixI^S, 1:^ND)
127 integer,
intent(in) :: ixI^L, ixO^L
128 double precision,
intent(inout) :: w(ixI^S, nw)
129 double precision,
intent(in) :: x(ixI^S, 1:^ND)
136 integer,
intent(in) :: ixi^
l, ixo^
l, idim
137 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
138 double precision,
intent(out) :: v(ixi^s)
142 v(ixo^s)=w(ixo^s,
rho_)
144 v(ixo^s)=3.0d0*w(ixo^s,
rho_)**2
146 call mpistop(
'Undefined fluxtype: set nonlinear_flux_type to 1 or 2')
153 integer,
intent(in) :: ixI^L, ixO^L, idim
154 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
155 double precision,
intent(inout) :: cmax(ixI^S)
159 cmax(ixo^s) = abs(cmax(ixo^s))
163 subroutine nonlinear_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
166 integer,
intent(in) :: ixI^L, ixO^L, idim
167 double precision,
intent(in) :: wLC(ixI^S, nw), wRC(ixI^S,nw)
168 double precision,
intent(in) :: wLp(ixI^S, nw), wRp(ixI^S,nw)
169 double precision,
intent(in) :: x(ixI^S, 1:^ND)
170 double precision,
intent(inout) :: cmax(ixI^S,1:number_species)
171 double precision,
intent(inout),
optional :: cmin(ixI^S,1:number_species)
172 double precision,
intent(in) :: Hspeed(ixI^S,1:number_species)
173 double precision :: wmean(ixI^S,nw)
180 if (
present(cmin))
then
181 cmin(ixo^s,1) = min(cmax(ixo^s,1), zero)
182 cmax(ixo^s,1) = max(cmax(ixo^s,1), zero)
184 cmax(ixo^s,1) = maxval(abs(cmax(ixo^s,1)))
193 integer,
intent(in) :: ixI^L, ixO^L
194 double precision,
intent(in) :: dx^D, x(ixI^S, 1:^ND)
195 double precision,
intent(in) :: w(ixI^S, 1:nw)
196 double precision,
intent(inout) :: dtnew
201 call kdv_get_dt(w, ixi^l, ixo^l, dtnew, dx^d, x)
208 integer,
intent(in) :: ixI^L, ixO^L, idim
209 double precision,
intent(in) :: wC(ixI^S, 1:nw)
210 double precision,
intent(in) :: w(ixI^S, 1:nw)
211 double precision,
intent(in) :: x(ixI^S, 1:^ND)
212 double precision,
intent(out) :: f(ixI^S, nwflux)
220 call mpistop(
'Undefined fluxtype: set nonlinear_flux_type to 1 or 2')
232 integer,
intent(in) :: ixI^L, ixO^L
233 double precision,
intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
234 double precision,
intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S,1:nw),w(ixI^S, 1:nw)
242 integer,
intent(in) :: ixI^L, ixO^L
243 double precision,
intent(in) :: qdt,dtfactor
244 double precision,
intent(in) :: wCT(ixI^S, 1:nw),wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
245 double precision,
intent(inout) :: w(ixI^S, 1:nw)
246 logical,
intent(in) :: qsourcesplit
247 logical,
intent(inout) :: active
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, parameter unitpar
file handle for IO
logical use_particles
Use particles module or not.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer ndir
Number of spatial dimensions (components) for vector variables.
Module for including kdv source term in simulations adds over dimensions i.
subroutine kdv_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine kdv_init()
Initialize the module.
subroutine kdv_add_source(qdt, ixIL, ixOL, wCT, w, x, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
Module containing the physics routines for scalar nonlinear equation.
subroutine nonlinear_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, wprim, w, x)
subroutine nonlinear_get_flux(wC, w, x, ixIL, ixOL, idim, f)
logical, public, protected nonlinear_particles
Whether particles module is added.
subroutine nonlinear_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
integer, public, protected nonlinear_flux_type
switch between burgers (i.e. rho**2) or nonconvex flux (i.e. rho**3)
logical, public, protected kdv_source_term
whether the KdV source term is added
integer, public, protected rho_
index of the single scalar unknown
subroutine, public nonlinear_phys_init()
subroutine, public nonlinear_get_v(w, x, ixIL, ixOL, idim, v)
subroutine nonlinear_get_cmax(w, x, ixIL, ixOL, idim, cmax)
subroutine nonlinear_get_cbounds(wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
subroutine nonlinear_add_source(qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
subroutine nonlinear_to_conserved(ixIL, ixOL, w, x)
subroutine nonlinear_write_info(fh)
Write this module's parameters to a snapshot.
subroutine nonlinear_to_primitive(ixIL, ixOL, w, x)
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
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
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_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.