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
36 subroutine rho_write_info(fh)
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)
54 end subroutine rho_write_info
67 allocate(start_indices(number_species),stop_indices(number_species))
76 stop_indices(1)=nwflux
83 call mpistop(
"phys_check error: flux_type has wrong shape")
102 subroutine rho_to_conserved(ixI^L, ixO^L, w, x)
104 integer,
intent(in) :: ixi^
l, ixo^
l
105 double precision,
intent(inout) :: w(ixi^s, nw)
106 double precision,
intent(in) :: x(ixi^s, 1:^nd)
109 end subroutine rho_to_conserved
111 subroutine rho_to_primitive(ixI^L, ixO^L, w, x)
113 integer,
intent(in) :: ixi^
l, ixo^
l
114 double precision,
intent(inout) :: w(ixi^s, nw)
115 double precision,
intent(in) :: x(ixi^s, 1:^nd)
118 end subroutine rho_to_primitive
120 subroutine rho_get_v(w, x, ixI^L, ixO^L, idim, v, centered)
123 logical,
intent(in) :: centered
124 integer,
intent(in) :: ixi^
l, ixo^
l, idim
125 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
126 double precision,
intent(out) :: v(ixi^s)
128 double precision :: dtheta, dphi, halfdtheta, halfdphi, invdtheta, invdphi
130 double precision :: appcosphi(ixi^s), appsinphi(ixi^s), &
131 appcosthe(ixi^s), appsinthe(ixi^s)
137 call mpistop(
"advection in 1D cylindrical not available")
144 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2))+
rho_v(2)*dsin(x(ixo^s,2))
146 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,2))+
rho_v(2)*dcos(x(ixo^s,2))
150 dtheta=x(ixomin1,ixomin2+1,2)-x(ixomin1,ixomin2,2)
151 halfdtheta=0.5d0*dtheta
152 invdtheta=1.0d0/dtheta
155 v(ixo^s) =(
rho_v(1)*( dsin(x(ixo^s,2)+halfdtheta) &
156 -dsin(x(ixo^s,2)-halfdtheta)) &
157 +
rho_v(2)*(-dcos(x(ixo^s,2)+halfdtheta) &
158 +dcos(x(ixo^s,2)-halfdtheta)))*invdtheta
160 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,2)+halfdtheta) &
161 +
rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)
170 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,3))+
rho_v(2)*dsin(x(ixo^s,3))
174 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3))+
rho_v(2)*dcos(x(ixo^s,3))
178 dtheta=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
179 halfdtheta=0.5d0*dtheta
180 invdtheta=1.0d0/dtheta
183 v(ixo^s) =(
rho_v(1)*( dsin(x(ixo^s,3)+halfdtheta) &
184 -dsin(x(ixo^s,3)-halfdtheta)) &
185 +
rho_v(2)*(-dcos(x(ixo^s,3)+halfdtheta) &
186 +dcos(x(ixo^s,3)-halfdtheta)))*invdtheta
190 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)+halfdtheta) &
191 +
rho_v(2)*dcos(x(ixo^s,3)+halfdtheta)
197 call mpistop(
"advection in 1D spherical not available")
200 call mpistop(
"advection in 2D spherical not available")
207 v(ixo^s) =
rho_v(1)*dsin(x(ixo^s,2))*dcos(x(ixo^s,3)) &
208 +
rho_v(2)*dsin(x(ixo^s,2))*dsin(x(ixo^s,3)) &
209 +
rho_v(3)*dcos(x(ixo^s,2))
211 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2))*dcos(x(ixo^s,3)) &
212 +
rho_v(2)*dcos(x(ixo^s,2))*dsin(x(ixo^s,3)) &
213 -
rho_v(3)*dsin(x(ixo^s,2))
215 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)) &
216 +
rho_v(2)*dcos(x(ixo^s,3))
220 dtheta=x(ixomin1,ixomin2+1,ixomin3,2)-x(ixomin1,ixomin2,ixomin3,2)
221 dphi=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
222 halfdtheta=0.5d0*dtheta
224 invdtheta=1.0d0/dtheta
228 appcosphi(ixo^s)=( dsin(x(ixo^s,3)+halfdphi) &
229 -dsin(x(ixo^s,3)-halfdphi))*invdphi
230 appsinphi(ixo^s)=(-dcos(x(ixo^s,3)+halfdphi) &
231 +dcos(x(ixo^s,3)-halfdphi))*invdphi
232 appcosthe(ixo^s)=(dsin(x(ixo^s,2)+halfdtheta)**2 &
233 -dsin(x(ixo^s,2)-halfdtheta)**2) &
234 /(4.0d0*dabs(dsin(x(ixo^s,2)))*dsin(halfdtheta))
236 (-dsin(x(ixo^s,2)+halfdtheta)*dcos(x(ixo^s,2)+halfdtheta) &
237 +dsin(x(ixo^s,2)-halfdtheta)*dcos(x(ixo^s,2)-halfdtheta) &
238 +dtheta)/(4.0d0*dabs(dsin(x(ixo^s,2)))*dsin(halfdtheta))
239 v(ixo^s) =
rho_v(1)*appsinthe(ixo^s)*appcosphi(ixo^s) &
240 +
rho_v(2)*appsinthe(ixo^s)*appsinphi(ixo^s) &
241 +
rho_v(3)*appcosthe(ixo^s)
243 appcosphi(ixo^s)=( dsin(x(ixo^s,3)+halfdphi) &
244 -dsin(x(ixo^s,3)-halfdphi))*invdphi
245 appsinphi(ixo^s)=(-dcos(x(ixo^s,3)+halfdphi) &
246 +dcos(x(ixo^s,3)-halfdphi))*invdphi
247 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2)+halfdtheta)*appcosphi(ixo^s) &
248 +
rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)*appsinphi(ixo^s) &
249 -
rho_v(3)*dsin(x(ixo^s,2)+halfdtheta)
251 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)+halfdphi) &
252 +
rho_v(2)*dcos(x(ixo^s,3)+halfdphi)
257 v(ixo^s) =
rho_v(idim)
262 subroutine rho_get_v_idim(w,x,ixI^L,ixO^L,idim,v)
265 integer,
intent(in) :: ixi^
l, ixo^
l, idim
266 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:
ndim)
267 double precision,
intent(out) :: v(ixi^s)
269 v(ixo^s) =
rho_v(idim)
271 end subroutine rho_get_v_idim
273 subroutine rho_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
275 integer,
intent(in) :: ixi^
l, ixo^
l, idim
276 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
277 double precision,
intent(inout) :: cmax(ixi^s)
279 call rho_get_v(w, x, ixi^
l, ixo^
l, idim, cmax, .true.)
281 cmax(ixo^s) = abs(cmax(ixo^s))
283 end subroutine rho_get_cmax
285 subroutine rho_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
288 integer,
intent(in) :: ixi^
l, ixo^
l, idim
289 double precision,
intent(in) :: wlc(ixi^s,
nw), wrc(ixi^s,
nw)
290 double precision,
intent(in) :: wlp(ixi^s,
nw), wrp(ixi^s,
nw)
291 double precision,
intent(in) :: x(ixi^s, 1:^nd)
293 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
298 call rho_get_v(wlc, x, ixi^
l, ixo^
l, idim, cmax(ixi^s,1), .false.)
300 if (
present(cmin))
then
301 cmin(ixo^s,1) = min(cmax(ixo^s,1), zero)
302 cmax(ixo^s,1) = max(cmax(ixo^s,1), zero)
304 cmax(ixo^s,1) = maxval(abs(cmax(ixo^s,1)))
307 end subroutine rho_get_cbounds
309 subroutine rho_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
311 integer,
intent(in) :: ixi^
l, ixo^
l
312 double precision,
intent(in) ::
dx^
d, x(ixi^s, 1:^nd)
313 double precision,
intent(in) :: w(ixi^s, 1:nw)
314 double precision,
intent(inout) :: dtnew
317 end subroutine rho_get_dt
320 subroutine rho_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
322 integer,
intent(in) :: ixi^
l, ixo^
l, idim
323 double precision,
intent(in) :: wc(ixi^s, 1:nw)
324 double precision,
intent(in) :: w(ixi^s, 1:nw)
325 double precision,
intent(in) :: x(ixi^s, 1:^nd)
326 double precision,
intent(out) :: f(ixi^s, nwflux)
327 double precision :: v(ixi^s)
329 call rho_get_v(wc, x, ixi^
l, ixo^
l, idim, v, .false.)
331 f(ixo^s,
rho_) = w(ixo^s,
rho_) * v(ixo^s)
332 end subroutine rho_get_flux
334 subroutine rho_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT,wprim, w, x)
341 integer,
intent(in) :: ixi^
l, ixo^
l
342 double precision,
intent(in) :: qdt, dtfactor, x(ixi^s, 1:^nd)
343 double precision,
intent(inout) :: wct(ixi^s, 1:nw),wprim(ixi^s,1:nw), w(ixi^s, 1:nw)
345 end subroutine rho_add_source_geom
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.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision, dimension(:,:), allocatable dx
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_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, public rho_get_v(w, x, ixil, ixol, idim, v, centered)
double precision, dimension(^nd), public, protected rho_v
subroutine, public rho_phys_init()
logical, public, protected rho_particles
Whether particles module is added.
integer, public, protected rho_
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...