7 integer,
protected,
public ::
rho_ = 1
12 double precision,
protected,
public ::
rho_v(^nd) = 1.0d0
20 subroutine rho_params_read(files)
22 character(len=*),
intent(in) :: files(:)
28 open(
unitpar, file=trim(files(n)), status=
'old')
29 read(
unitpar, rho_list,
end=111)
33 end subroutine rho_params_read
38 integer,
intent(in) :: fh
39 integer,
parameter :: n_par = ^nd
40 double precision :: values(n_par)
41 character(len=name_len) :: names(n_par)
42 integer,
dimension(MPI_STATUS_SIZE) :: st
46 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
49 write(names(idim),
'(a,i1)')
"v",idim
50 values(idim) =
rho_v(idim)
52 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
53 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
68 allocate(start_indices(number_species),stop_indices(number_species))
77 stop_indices(1)=nwflux
84 call mpistop(
"phys_check error: flux_type has wrong shape")
106 integer,
intent(in) :: ixI^L, ixO^L
107 double precision,
intent(inout) :: w(ixI^S, nw)
108 double precision,
intent(in) :: x(ixI^S, 1:^ND)
115 integer,
intent(in) :: ixI^L, ixO^L
116 double precision,
intent(inout) :: w(ixI^S, nw)
117 double precision,
intent(in) :: x(ixI^S, 1:^ND)
122 subroutine rho_get_v(w, x, ixI^L, ixO^L, idim, v, centered)
125 logical,
intent(in) :: centered
126 integer,
intent(in) :: ixi^
l, ixo^
l, idim
127 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
128 double precision,
intent(out) :: v(ixi^s)
130 double precision :: dtheta, dphi, halfdtheta, halfdphi, invdtheta, invdphi
132 double precision :: appcosphi(ixi^s), appsinphi(ixi^s), &
133 appcosthe(ixi^s), appsinthe(ixi^s)
139 call mpistop(
"advection in 1D cylindrical not available")
146 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2))+
rho_v(2)*dsin(x(ixo^s,2))
148 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,2))+
rho_v(2)*dcos(x(ixo^s,2))
152 dtheta=x(ixomin1,ixomin2+1,2)-x(ixomin1,ixomin2,2)
153 halfdtheta=0.5d0*dtheta
154 invdtheta=1.0d0/dtheta
157 v(ixo^s) =(
rho_v(1)*( dsin(x(ixo^s,2)+halfdtheta) &
158 -dsin(x(ixo^s,2)-halfdtheta)) &
159 +
rho_v(2)*(-dcos(x(ixo^s,2)+halfdtheta) &
160 +dcos(x(ixo^s,2)-halfdtheta)))*invdtheta
162 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,2)+halfdtheta) &
163 +
rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)
172 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,3))+
rho_v(2)*dsin(x(ixo^s,3))
176 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3))+
rho_v(2)*dcos(x(ixo^s,3))
180 dtheta=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
181 halfdtheta=0.5d0*dtheta
182 invdtheta=1.0d0/dtheta
185 v(ixo^s) =(
rho_v(1)*( dsin(x(ixo^s,3)+halfdtheta) &
186 -dsin(x(ixo^s,3)-halfdtheta)) &
187 +
rho_v(2)*(-dcos(x(ixo^s,3)+halfdtheta) &
188 +dcos(x(ixo^s,3)-halfdtheta)))*invdtheta
192 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)+halfdtheta) &
193 +
rho_v(2)*dcos(x(ixo^s,3)+halfdtheta)
199 call mpistop(
"advection in 1D spherical not available")
202 call mpistop(
"advection in 2D spherical not available")
209 v(ixo^s) =
rho_v(1)*dsin(x(ixo^s,2))*dcos(x(ixo^s,3)) &
210 +
rho_v(2)*dsin(x(ixo^s,2))*dsin(x(ixo^s,3)) &
211 +
rho_v(3)*dcos(x(ixo^s,2))
213 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2))*dcos(x(ixo^s,3)) &
214 +
rho_v(2)*dcos(x(ixo^s,2))*dsin(x(ixo^s,3)) &
215 -
rho_v(3)*dsin(x(ixo^s,2))
217 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)) &
218 +
rho_v(2)*dcos(x(ixo^s,3))
222 dtheta=x(ixomin1,ixomin2+1,ixomin3,2)-x(ixomin1,ixomin2,ixomin3,2)
223 dphi=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
224 halfdtheta=0.5d0*dtheta
226 invdtheta=1.0d0/dtheta
230 appcosphi(ixo^s)=( dsin(x(ixo^s,3)+halfdphi) &
231 -dsin(x(ixo^s,3)-halfdphi))*invdphi
232 appsinphi(ixo^s)=(-dcos(x(ixo^s,3)+halfdphi) &
233 +dcos(x(ixo^s,3)-halfdphi))*invdphi
234 appcosthe(ixo^s)=(dsin(x(ixo^s,2)+halfdtheta)**2 &
235 -dsin(x(ixo^s,2)-halfdtheta)**2) &
236 /(4.0d0*dabs(dsin(x(ixo^s,2)))*dsin(halfdtheta))
238 (-dsin(x(ixo^s,2)+halfdtheta)*dcos(x(ixo^s,2)+halfdtheta) &
239 +dsin(x(ixo^s,2)-halfdtheta)*dcos(x(ixo^s,2)-halfdtheta) &
240 +dtheta)/(4.0d0*dabs(dsin(x(ixo^s,2)))*dsin(halfdtheta))
241 v(ixo^s) =
rho_v(1)*appsinthe(ixo^s)*appcosphi(ixo^s) &
242 +
rho_v(2)*appsinthe(ixo^s)*appsinphi(ixo^s) &
243 +
rho_v(3)*appcosthe(ixo^s)
245 appcosphi(ixo^s)=( dsin(x(ixo^s,3)+halfdphi) &
246 -dsin(x(ixo^s,3)-halfdphi))*invdphi
247 appsinphi(ixo^s)=(-dcos(x(ixo^s,3)+halfdphi) &
248 +dcos(x(ixo^s,3)-halfdphi))*invdphi
249 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2)+halfdtheta)*appcosphi(ixo^s) &
250 +
rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)*appsinphi(ixo^s) &
251 -
rho_v(3)*dsin(x(ixo^s,2)+halfdtheta)
253 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)+halfdphi) &
254 +
rho_v(2)*dcos(x(ixo^s,3)+halfdphi)
259 v(ixo^s) =
rho_v(idim)
267 integer,
intent(in) :: ixI^L, ixO^L, idim
268 double precision,
intent(in) :: w(ixI^S,nw), x(ixI^S,1:ndim)
269 double precision,
intent(out) :: v(ixI^S)
271 v(ixo^s) =
rho_v(idim)
277 integer,
intent(in) :: ixI^L, ixO^L, idim
278 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
279 double precision,
intent(inout) :: cmax(ixI^S)
281 call rho_get_v(w, x, ixi^l, ixo^l, idim, cmax, .true.)
283 cmax(ixo^s) = abs(cmax(ixo^s))
287 subroutine rho_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
290 integer,
intent(in) :: ixI^L, ixO^L, idim
291 double precision,
intent(in) :: wLC(ixI^S, nw), wRC(ixI^S,nw)
292 double precision,
intent(in) :: wLp(ixI^S, nw), wRp(ixI^S,nw)
293 double precision,
intent(in) :: x(ixI^S, 1:^ND)
294 double precision,
intent(inout) :: cmax(ixI^S,1:number_species)
295 double precision,
intent(inout),
optional :: cmin(ixI^S,1:number_species)
296 double precision,
intent(in) :: Hspeed(ixI^S,1:number_species)
300 call rho_get_v(wlc, x, ixi^l, ixo^l, idim, cmax(ixi^s,1), .false.)
302 if (
present(cmin))
then
303 cmin(ixo^s,1) = min(cmax(ixo^s,1), zero)
304 cmax(ixo^s,1) = max(cmax(ixo^s,1), zero)
306 cmax(ixo^s,1) = maxval(abs(cmax(ixo^s,1)))
313 integer,
intent(in) :: ixI^L, ixO^L
314 double precision,
intent(in) :: dx^D, x(ixI^S, 1:^ND)
315 double precision,
intent(in) :: w(ixI^S, 1:nw)
316 double precision,
intent(inout) :: dtnew
324 integer,
intent(in) :: ixI^L, ixO^L, idim
325 double precision,
intent(in) :: wC(ixI^S, 1:nw)
326 double precision,
intent(in) :: w(ixI^S, 1:nw)
327 double precision,
intent(in) :: x(ixI^S, 1:^ND)
328 double precision,
intent(out) :: f(ixI^S, nwflux)
329 double precision :: v(ixI^S)
331 call rho_get_v(wc, x, ixi^l, ixo^l, idim, v, .false.)
333 f(ixo^s,
rho_) = w(ixo^s,
rho_) * v(ixo^s)
343 integer,
intent(in) :: ixI^L, ixO^L
344 double precision,
intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
345 double precision,
intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
Module with geometry-related routines (e.g., divergence, curl)
integer, parameter spherical
integer, parameter cylindrical
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, parameter unitpar
file handle for IO
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 ndir
Number of spatial dimensions (components) for vector variables.
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_get_cmax), pointer phys_get_cmax
logical phys_energy
Solve energy equation or not.
Module containing the physics routines for scalar advection.
subroutine rho_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, w, x)
subroutine, public rho_get_v(w, x, ixIL, ixOL, idim, v, centered)
subroutine rho_to_conserved(ixIL, ixOL, w, x)
subroutine rho_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine rho_to_primitive(ixIL, ixOL, w, x)
double precision, dimension(^nd), public, protected rho_v
subroutine rho_get_v_idim(w, x, ixIL, ixOL, idim, v)
Calculate simple v component.
subroutine rho_write_info(fh)
Write this module's parameters to a snapsoht.
subroutine rho_get_flux(wC, w, x, ixIL, ixOL, idim, f)
subroutine, public rho_phys_init()
subroutine rho_get_cbounds(wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
subroutine rho_get_cmax(w, x, ixIL, ixOL, idim, cmax)
logical, public, protected rho_particles
Whether particles module is added.
integer, public, protected rho_