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
66 allocate(start_indices(number_species),stop_indices(number_species))
75 stop_indices(1)=nwflux
82 call mpistop(
"phys_check error: flux_type has wrong shape")
98 subroutine rho_check_params
102 use mod_particles,
only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
110 write(*,*)
'====RHO run with settings===================='
111 write(*,*)
'Using mod_rho_phys with settings:'
112 write(*,*)
'Dimensionality :',
ndim
113 write(*,*)
'vector components:',
ndir
115 write(*,*)
'number of variables nw=',nw
116 write(*,*)
' start index iwstart=',iwstart
117 write(*,*)
'number of vector variables=',nvector
118 write(*,*)
'number of stagger variables nws=',nws
119 write(*,*)
'number of variables with BCs=',nwgc
120 write(*,*)
'number of vars with fluxes=',nwflux
121 write(*,*)
'number of vars with flux + BC=',nwfluxbc
122 write(*,*)
'number of auxiliary variables=',nwaux
123 write(*,*)
'number of extra vars without flux=',nwextra
124 write(*,*)
'number of extra vars for wextra=',nw_extra
125 write(*,*)
'number of auxiliary I/O variables=',
nwauxio
128 write(*,*)
'*****Using particles: npayload,ngridvars :', npayload,ngridvars
129 write(*,*)
'*****Using particles: nusrpayload :', nusrpayload
130 write(*,*)
'*****Using particles: num_particles :', num_particles
131 write(*,*)
'*****Using particles: physics_type_particles=',physics_type_particles
134 write(*,*)
'==========================================='
137 end subroutine rho_check_params
139 subroutine rho_to_conserved(ixI^L, ixO^L, w, x)
141 integer,
intent(in) :: ixi^
l, ixo^
l
142 double precision,
intent(inout) :: w(ixi^s, nw)
143 double precision,
intent(in) :: x(ixi^s, 1:^nd)
146 end subroutine rho_to_conserved
148 subroutine rho_to_primitive(ixI^L, ixO^L, w, x)
150 integer,
intent(in) :: ixi^
l, ixo^
l
151 double precision,
intent(inout) :: w(ixi^s, nw)
152 double precision,
intent(in) :: x(ixi^s, 1:^nd)
155 end subroutine rho_to_primitive
157 subroutine rho_get_v(w, x, ixI^L, ixO^L, idim, v, centered)
160 logical,
intent(in) :: centered
161 integer,
intent(in) :: ixi^
l, ixo^
l, idim
162 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
163 double precision,
intent(out) :: v(ixi^s)
165 double precision :: dtheta, dphi, halfdtheta, halfdphi, invdtheta, invdphi
167 double precision :: appcosphi(ixi^s), appsinphi(ixi^s), &
168 appcosthe(ixi^s), appsinthe(ixi^s)
174 call mpistop(
"advection in 1D cylindrical not available")
181 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2))+
rho_v(2)*dsin(x(ixo^s,2))
183 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,2))+
rho_v(2)*dcos(x(ixo^s,2))
187 dtheta=x(ixomin1,ixomin2+1,2)-x(ixomin1,ixomin2,2)
188 halfdtheta=0.5d0*dtheta
189 invdtheta=1.0d0/dtheta
192 v(ixo^s) =(
rho_v(1)*( dsin(x(ixo^s,2)+halfdtheta) &
193 -dsin(x(ixo^s,2)-halfdtheta)) &
194 +
rho_v(2)*(-dcos(x(ixo^s,2)+halfdtheta) &
195 +dcos(x(ixo^s,2)-halfdtheta)))*invdtheta
197 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,2)+halfdtheta) &
198 +
rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)
207 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,3))+
rho_v(2)*dsin(x(ixo^s,3))
211 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3))+
rho_v(2)*dcos(x(ixo^s,3))
215 dtheta=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
216 halfdtheta=0.5d0*dtheta
217 invdtheta=1.0d0/dtheta
220 v(ixo^s) =(
rho_v(1)*( dsin(x(ixo^s,3)+halfdtheta) &
221 -dsin(x(ixo^s,3)-halfdtheta)) &
222 +
rho_v(2)*(-dcos(x(ixo^s,3)+halfdtheta) &
223 +dcos(x(ixo^s,3)-halfdtheta)))*invdtheta
227 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)+halfdtheta) &
228 +
rho_v(2)*dcos(x(ixo^s,3)+halfdtheta)
234 call mpistop(
"advection in 1D spherical not available")
237 call mpistop(
"advection in 2D spherical not available")
244 v(ixo^s) =
rho_v(1)*dsin(x(ixo^s,2))*dcos(x(ixo^s,3)) &
245 +
rho_v(2)*dsin(x(ixo^s,2))*dsin(x(ixo^s,3)) &
246 +
rho_v(3)*dcos(x(ixo^s,2))
248 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2))*dcos(x(ixo^s,3)) &
249 +
rho_v(2)*dcos(x(ixo^s,2))*dsin(x(ixo^s,3)) &
250 -
rho_v(3)*dsin(x(ixo^s,2))
252 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)) &
253 +
rho_v(2)*dcos(x(ixo^s,3))
257 dtheta=x(ixomin1,ixomin2+1,ixomin3,2)-x(ixomin1,ixomin2,ixomin3,2)
258 dphi=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
259 halfdtheta=0.5d0*dtheta
261 invdtheta=1.0d0/dtheta
265 appcosphi(ixo^s)=( dsin(x(ixo^s,3)+halfdphi) &
266 -dsin(x(ixo^s,3)-halfdphi))*invdphi
267 appsinphi(ixo^s)=(-dcos(x(ixo^s,3)+halfdphi) &
268 +dcos(x(ixo^s,3)-halfdphi))*invdphi
269 appcosthe(ixo^s)=(dsin(x(ixo^s,2)+halfdtheta)**2 &
270 -dsin(x(ixo^s,2)-halfdtheta)**2) &
271 /(4.0d0*dabs(dsin(x(ixo^s,2)))*dsin(halfdtheta))
273 (-dsin(x(ixo^s,2)+halfdtheta)*dcos(x(ixo^s,2)+halfdtheta) &
274 +dsin(x(ixo^s,2)-halfdtheta)*dcos(x(ixo^s,2)-halfdtheta) &
275 +dtheta)/(4.0d0*dabs(dsin(x(ixo^s,2)))*dsin(halfdtheta))
276 v(ixo^s) =
rho_v(1)*appsinthe(ixo^s)*appcosphi(ixo^s) &
277 +
rho_v(2)*appsinthe(ixo^s)*appsinphi(ixo^s) &
278 +
rho_v(3)*appcosthe(ixo^s)
280 appcosphi(ixo^s)=( dsin(x(ixo^s,3)+halfdphi) &
281 -dsin(x(ixo^s,3)-halfdphi))*invdphi
282 appsinphi(ixo^s)=(-dcos(x(ixo^s,3)+halfdphi) &
283 +dcos(x(ixo^s,3)-halfdphi))*invdphi
284 v(ixo^s) =
rho_v(1)*dcos(x(ixo^s,2)+halfdtheta)*appcosphi(ixo^s) &
285 +
rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)*appsinphi(ixo^s) &
286 -
rho_v(3)*dsin(x(ixo^s,2)+halfdtheta)
288 v(ixo^s) =-
rho_v(1)*dsin(x(ixo^s,3)+halfdphi) &
289 +
rho_v(2)*dcos(x(ixo^s,3)+halfdphi)
294 v(ixo^s) =
rho_v(idim)
299 subroutine rho_get_v_idim(w,x,ixI^L,ixO^L,idim,v)
302 integer,
intent(in) :: ixi^
l, ixo^
l, idim
303 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:
ndim)
304 double precision,
intent(out) :: v(ixi^s)
306 v(ixo^s) =
rho_v(idim)
308 end subroutine rho_get_v_idim
310 subroutine rho_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
312 integer,
intent(in) :: ixi^
l, ixo^
l, idim
313 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
314 double precision,
intent(inout) :: cmax(ixi^s)
316 call rho_get_v(w, x, ixi^
l, ixo^
l, idim, cmax, .true.)
318 cmax(ixo^s) = abs(cmax(ixo^s))
320 end subroutine rho_get_cmax
322 subroutine rho_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
325 integer,
intent(in) :: ixi^
l, ixo^
l, idim
326 double precision,
intent(in) :: wlc(ixi^s,
nw), wrc(ixi^s,
nw)
327 double precision,
intent(in) :: wlp(ixi^s,
nw), wrp(ixi^s,
nw)
328 double precision,
intent(in) :: x(ixi^s, 1:^nd)
330 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
335 call rho_get_v(wlc, x, ixi^
l, ixo^
l, idim, cmax(ixi^s,1), .false.)
337 if (
present(cmin))
then
338 cmin(ixo^s,1) = min(cmax(ixo^s,1), zero)
339 cmax(ixo^s,1) = max(cmax(ixo^s,1), zero)
341 cmax(ixo^s,1) = maxval(abs(cmax(ixo^s,1)))
344 end subroutine rho_get_cbounds
346 subroutine rho_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
348 integer,
intent(in) :: ixi^
l, ixo^
l
349 double precision,
intent(in) ::
dx^
d, x(ixi^s, 1:^nd)
350 double precision,
intent(in) :: w(ixi^s, 1:nw)
351 double precision,
intent(inout) :: dtnew
354 end subroutine rho_get_dt
357 subroutine rho_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
359 integer,
intent(in) :: ixi^
l, ixo^
l, idim
360 double precision,
intent(in) :: wc(ixi^s, 1:nw)
361 double precision,
intent(in) :: w(ixi^s, 1:nw)
362 double precision,
intent(in) :: x(ixi^s, 1:^nd)
363 double precision,
intent(out) :: f(ixi^s, nwflux)
364 double precision :: v(ixi^s)
366 call rho_get_v(wc, x, ixi^
l, ixo^
l, idim, v, .false.)
368 f(ixo^s,
rho_) = w(ixo^s,
rho_) * v(ixo^s)
369 end subroutine rho_get_flux
371 subroutine rho_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT,wprim, w, x)
378 integer,
intent(in) :: ixi^
l, ixo^
l
379 double precision,
intent(in) :: qdt, dtfactor, x(ixi^s, 1:^nd)
380 double precision,
intent(inout) :: wct(ixi^s, 1:nw),wprim(ixi^s,1:nw), w(ixi^s, 1:nw)
382 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.
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.
logical slab
Cartesian geometry or not.
integer nwauxio
Number of auxiliary variables that are only included in output.
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
integer nghostcells
Number of ghost cells surrounding a grid.
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
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_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...