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
86 call mpistop(
"phys_check error: flux_type has wrong shape")
105 subroutine rho_to_conserved(ixI^L, ixO^L, w, x)
107 integer,
intent(in) :: ixi^
l, ixo^
l
108 double precision,
intent(inout) :: w(ixi^s, nw)
109 double precision,
intent(in) :: x(ixi^s, 1:^nd)
112 end subroutine rho_to_conserved
114 subroutine rho_to_primitive(ixI^L, ixO^L, w, x)
116 integer,
intent(in) :: ixi^
l, ixo^
l
117 double precision,
intent(inout) :: w(ixi^s, nw)
118 double precision,
intent(in) :: x(ixi^s, 1:^nd)
121 end subroutine rho_to_primitive
123 subroutine rho_get_v(w, x, ixI^L, ixO^L, idim, v, centered)
126 logical,
intent(in) :: centered
127 integer,
intent(in) :: ixi^
l, ixo^
l, idim
128 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
129 double precision,
intent(out) :: v(ixi^s)
131 double precision :: dtheta, dphi, halfdtheta, halfdphi, invdtheta, invdphi
133 double precision :: appcosphi(ixi^s), appsinphi(ixi^s), &
134 appcosthe(ixi^s), appsinthe(ixi^s)
140 call mpistop(
"advection in 1D cylindrical not available")
147 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2))+
rho_v(2)*dsin(x(ixo^s,2))
149 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,2))+
rho_v(2)*dcos(x(ixo^s,2))
153 dtheta=x(ixomin1,ixomin2+1,2)-x(ixomin1,ixomin2,2)
154 halfdtheta=0.5d0*dtheta
155 invdtheta=1.0d0/dtheta
158 v(ixo^s) =(
rho_v(1)*( dsin(x(ixo^s,2)+halfdtheta) &
159 -dsin(x(ixo^s,2)-halfdtheta)) &
160 +
rho_v(2)*(-dcos(x(ixo^s,2)+halfdtheta) &
161 +dcos(x(ixo^s,2)-halfdtheta)))*invdtheta
163 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,2)+halfdtheta) &
164 +
rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)
173 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,3))+
rho_v(2)*dsin(x(ixo^s,3))
177 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3))+
rho_v(2)*dcos(x(ixo^s,3))
181 dtheta=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
182 halfdtheta=0.5d0*dtheta
183 invdtheta=1.0d0/dtheta
186 v(ixo^s) =(
rho_v(1)*( dsin(x(ixo^s,3)+halfdtheta) &
187 -dsin(x(ixo^s,3)-halfdtheta)) &
188 +
rho_v(2)*(-dcos(x(ixo^s,3)+halfdtheta) &
189 +dcos(x(ixo^s,3)-halfdtheta)))*invdtheta
193 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)+halfdtheta) &
194 +
rho_v(2)*dcos(x(ixo^s,3)+halfdtheta)
200 call mpistop(
"advection in 1D spherical not available")
203 call mpistop(
"advection in 2D spherical not available")
210 v(ixo^s) =
rho_v(1)*dsin(x(ixo^s,2))*dcos(x(ixo^s,3)) &
211 +
rho_v(2)*dsin(x(ixo^s,2))*dsin(x(ixo^s,3)) &
212 +
rho_v(3)*dcos(x(ixo^s,2))
214 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2))*dcos(x(ixo^s,3)) &
215 +
rho_v(2)*dcos(x(ixo^s,2))*dsin(x(ixo^s,3)) &
216 -
rho_v(3)*dsin(x(ixo^s,2))
218 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)) &
219 +
rho_v(2)*dcos(x(ixo^s,3))
223 dtheta=x(ixomin1,ixomin2+1,ixomin3,2)-x(ixomin1,ixomin2,ixomin3,2)
224 dphi=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
225 halfdtheta=0.5d0*dtheta
227 invdtheta=1.0d0/dtheta
231 appcosphi(ixo^s)=( dsin(x(ixo^s,3)+halfdphi) &
232 -dsin(x(ixo^s,3)-halfdphi))*invdphi
233 appsinphi(ixo^s)=(-dcos(x(ixo^s,3)+halfdphi) &
234 +dcos(x(ixo^s,3)-halfdphi))*invdphi
235 appcosthe(ixo^s)=(dsin(x(ixo^s,2)+halfdtheta)**2 &
236 -dsin(x(ixo^s,2)-halfdtheta)**2) &
237 /(4.0d0*dabs(dsin(x(ixo^s,2)))*dsin(halfdtheta))
239 (-dsin(x(ixo^s,2)+halfdtheta)*dcos(x(ixo^s,2)+halfdtheta) &
240 +dsin(x(ixo^s,2)-halfdtheta)*dcos(x(ixo^s,2)-halfdtheta) &
241 +dtheta)/(4.0d0*dabs(dsin(x(ixo^s,2)))*dsin(halfdtheta))
242 v(ixo^s) =
rho_v(1)*appsinthe(ixo^s)*appcosphi(ixo^s) &
243 +
rho_v(2)*appsinthe(ixo^s)*appsinphi(ixo^s) &
244 +
rho_v(3)*appcosthe(ixo^s)
246 appcosphi(ixo^s)=( dsin(x(ixo^s,3)+halfdphi) &
247 -dsin(x(ixo^s,3)-halfdphi))*invdphi
248 appsinphi(ixo^s)=(-dcos(x(ixo^s,3)+halfdphi) &
249 +dcos(x(ixo^s,3)-halfdphi))*invdphi
250 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2)+halfdtheta)*appcosphi(ixo^s) &
251 +
rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)*appsinphi(ixo^s) &
252 -
rho_v(3)*dsin(x(ixo^s,2)+halfdtheta)
254 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)+halfdphi) &
255 +
rho_v(2)*dcos(x(ixo^s,3)+halfdphi)
260 v(ixo^s) =
rho_v(idim)
265 subroutine rho_get_v_idim(w,x,ixI^L,ixO^L,idim,v)
268 integer,
intent(in) :: ixi^
l, ixo^
l, idim
269 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:
ndim)
270 double precision,
intent(out) :: v(ixi^s)
272 v(ixo^s) =
rho_v(idim)
274 end subroutine rho_get_v_idim
276 subroutine rho_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
278 integer,
intent(in) :: ixi^
l, ixo^
l, idim
279 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
280 double precision,
intent(inout) :: cmax(ixi^s)
282 call rho_get_v(w, x, ixi^
l, ixo^
l, idim, cmax, .true.)
284 cmax(ixo^s) = abs(cmax(ixo^s))
286 end subroutine rho_get_cmax
288 subroutine rho_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
291 integer,
intent(in) :: ixi^
l, ixo^
l, idim
292 double precision,
intent(in) :: wlc(ixi^s,
nw), wrc(ixi^s,
nw)
293 double precision,
intent(in) :: wlp(ixi^s,
nw), wrp(ixi^s,
nw)
294 double precision,
intent(in) :: x(ixi^s, 1:^nd)
296 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
301 call rho_get_v(wlc, x, ixi^
l, ixo^
l, idim, cmax(ixi^s,1), .false.)
303 if (
present(cmin))
then
304 cmin(ixo^s,1) = min(cmax(ixo^s,1), zero)
305 cmax(ixo^s,1) = max(cmax(ixo^s,1), zero)
307 cmax(ixo^s,1) = maxval(abs(cmax(ixo^s,1)))
310 end subroutine rho_get_cbounds
312 subroutine rho_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
314 integer,
intent(in) :: ixi^
l, ixo^
l
315 double precision,
intent(in) ::
dx^
d, x(ixi^s, 1:^nd)
316 double precision,
intent(in) :: w(ixi^s, 1:nw)
317 double precision,
intent(inout) :: dtnew
320 end subroutine rho_get_dt
323 subroutine rho_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
325 integer,
intent(in) :: ixi^
l, ixo^
l, idim
326 double precision,
intent(in) :: wc(ixi^s, 1:nw)
327 double precision,
intent(in) :: w(ixi^s, 1:nw)
328 double precision,
intent(in) :: x(ixi^s, 1:^nd)
329 double precision,
intent(out) :: f(ixi^s, nwflux)
330 double precision :: v(ixi^s)
332 call rho_get_v(wc, x, ixi^
l, ixo^
l, idim, v, .false.)
334 f(ixo^s,
rho_) = w(ixo^s,
rho_) * v(ixo^s)
335 end subroutine rho_get_flux
337 subroutine rho_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT,wprim, w, x)
344 integer,
intent(in) :: ixi^
l, ixo^
l
345 double precision,
intent(in) :: qdt, dtfactor, x(ixi^s, 1:^nd)
346 double precision,
intent(inout) :: wct(ixi^s, 1:nw),wprim(ixi^s,1:nw), w(ixi^s, 1:nw)
348 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...