MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_rho_phys.t
Go to the documentation of this file.
1!> Module containing the physics routines for scalar advection
3
4 implicit none
5 private
6
7 integer, protected, public :: rho_ = 1
8
9 !> Whether particles module is added
10 logical, public, protected :: rho_particles = .false.
11
12 double precision, protected, public :: rho_v(^nd) = 1.0d0
13
14 ! Public methods
15 public :: rho_phys_init
16 public :: rho_get_v
17
18contains
19
20 subroutine rho_params_read(files)
22 character(len=*), intent(in) :: files(:)
23 integer :: n
24
25 namelist /rho_list/ rho_v, rho_particles
26
27 do n = 1, size(files)
28 open(unitpar, file=trim(files(n)), status='old')
29 read(unitpar, rho_list, end=111)
30111 close(unitpar)
31 end do
32
33 end subroutine rho_params_read
34
35 !> Write this module's parameters to a snapsoht
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
43 integer :: er
44 integer :: idim
45
46 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
47
48 do idim=1,ndim
49 write(names(idim),'(a,i1)') "v",idim
50 values(idim) = rho_v(idim)
51 end do
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
55
56 subroutine rho_phys_init()
58 use mod_physics
59
60 call rho_params_read(par_files)
61
62 physics_type = "rho"
63 phys_energy = .false.
65
66 allocate(start_indices(number_species),stop_indices(number_species))
67 ! set the index of the first flux variable for species 1
68 start_indices(1)=1
69 rho_ = var_set_rho()
70
71 ! set number of variables which need update ghostcells
72 nwgc=nwflux
73
74 ! set the index of the last flux variable for species 1
75 stop_indices(1)=nwflux
76
77 ! Check whether custom flux types have been defined
78 if (.not. allocated(flux_type)) then
79 allocate(flux_type(ndir, nw))
81 else if (any(shape(flux_type) /= [ndir, nw])) then
82 call mpistop("phys_check error: flux_type has wrong shape")
83 end if
84
85 phys_get_cmax => rho_get_cmax
86 phys_get_cbounds => rho_get_cbounds
87 phys_get_flux => rho_get_flux
88 phys_add_source_geom => rho_add_source_geom
89 phys_to_conserved => rho_to_conserved
90 phys_to_primitive => rho_to_primitive
91 phys_get_dt => rho_get_dt
92 phys_check_params => rho_check_params
93 phys_write_info => rho_write_info
94
95
96 end subroutine rho_phys_init
97
98 subroutine rho_check_params
100 use mod_geometry, only: coordinate
102 use mod_particles, only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
103
104 ! Initialize particles module here, so user extra variables are sample
105 if (rho_particles) then
106 call particles_init()
107 end if
108
109 if(mype==0)then
110 write(*,*)'====RHO run with settings===================='
111 write(*,*)'Using mod_rho_phys with settings:'
112 write(*,*)'Dimensionality :',ndim
113 write(*,*)'vector components:',ndir
114 write(*,*)'coordinate set to type,slab:',coordinate,slab
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
126 write(*,*)' rho_particles=',rho_particles
127 if(rho_particles) then
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
132 end if
133 write(*,*)'number of ghostcells=',nghostcells
134 write(*,*)'==========================================='
135 endif
136
137 end subroutine rho_check_params
138
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)
144
145 ! Do nothing (primitive and conservative are equal for rho module)
146 end subroutine rho_to_conserved
147
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)
153
154 ! Do nothing (primitive and conservative are equal for rho module)
155 end subroutine rho_to_primitive
156
157 subroutine rho_get_v(w, x, ixI^L, ixO^L, idim, v, centered)
159 use mod_geometry
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)
164
165 double precision :: dtheta, dphi, halfdtheta, halfdphi, invdtheta, invdphi
166 {^ifthreed
167 double precision :: appcosphi(ixi^s), appsinphi(ixi^s), &
168 appcosthe(ixi^s), appsinthe(ixi^s)
169 }
170
171 select case (coordinate)
172 case (cylindrical)
173 {^ifoned
174 call mpistop("advection in 1D cylindrical not available")
175 }
176 {^iftwod
177 ! advection in 2D cylindrical: polar grid: to v_R, v_varphi
178 if(centered)then
179 select case (idim)
180 case (1) ! radial velocity
181 v(ixo^s) = rho_v(1)*dcos(x(ixo^s,2))+rho_v(2)*dsin(x(ixo^s,2))
182 case (2) ! v_varphi
183 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,2))+rho_v(2)*dcos(x(ixo^s,2))
184 end select
185 else
186 ! assumed uniform in varphi=theta
187 dtheta=x(ixomin1,ixomin2+1,2)-x(ixomin1,ixomin2,2)
188 halfdtheta=0.5d0*dtheta
189 invdtheta=1.0d0/dtheta
190 select case (idim)
191 case (1) ! radial velocity
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
196 case (2) ! v_varphi
197 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,2)+halfdtheta) &
198 +rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)
199 end select
200 endif
201 }
202 {^ifthreed
203 ! advection in 3D cylindrical: convert to v_R, v_Z, v_varphi
204 if(centered)then
205 select case (idim)
206 case (1) ! v_R velocity
207 v(ixo^s) = rho_v(1)*dcos(x(ixo^s,3))+rho_v(2)*dsin(x(ixo^s,3))
208 case (2) ! v_Z velocity
209 v(ixo^s) = rho_v(3)
210 case (3) ! v_varphi velocity
211 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,3))+rho_v(2)*dcos(x(ixo^s,3))
212 end select
213 else
214 ! assumed uniform in varphi=theta
215 dtheta=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
216 halfdtheta=0.5d0*dtheta
217 invdtheta=1.0d0/dtheta
218 select case (idim)
219 case (1) ! radial velocity
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
224 case (2) ! v_Z velocity
225 v(ixo^s) = rho_v(3)
226 case (3) ! v_varphi
227 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,3)+halfdtheta) &
228 +rho_v(2)*dcos(x(ixo^s,3)+halfdtheta)
229 end select
230 endif
231 }
232 case (spherical)
233 {^ifoned
234 call mpistop("advection in 1D spherical not available")
235 }
236 {^iftwod
237 call mpistop("advection in 2D spherical not available")
238 }
239 {^ifthreed
240 ! advection in 3D spherical: convert to v_r, v_theta, v_phi
241 if(centered)then
242 select case (idim)
243 case (1) ! radial velocity
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))
247 case (2) ! theta velocity
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))
251 case (3) ! phi velocity
252 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,3)) &
253 +rho_v(2)*dcos(x(ixo^s,3))
254 end select
255 else
256 ! assumed uniform in theta and phi
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
260 halfdphi=0.5d0*dphi
261 invdtheta=1.0d0/dtheta
262 invdphi=1.0d0/dphi
263 select case (idim)
264 case (1) ! radial velocity
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))
272 appsinthe(ixo^s)= &
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)
279 case (2) ! theta velocity
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)
287 case (3) ! phi velocity
288 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,3)+halfdphi) &
289 +rho_v(2)*dcos(x(ixo^s,3)+halfdphi)
290 end select
291 endif
292 }
293 case default
294 v(ixo^s) = rho_v(idim)
295 end select
296 end subroutine rho_get_v
297
298 !> Calculate simple v component
299 subroutine rho_get_v_idim(w,x,ixI^L,ixO^L,idim,v)
301
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)
305
306 v(ixo^s) = rho_v(idim)
307
308 end subroutine rho_get_v_idim
309
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)
315
316 call rho_get_v(w, x, ixi^l, ixo^l, idim, cmax, .true.)
317
318 cmax(ixo^s) = abs(cmax(ixo^s))
319
320 end subroutine rho_get_cmax
321
322 subroutine rho_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
324 use mod_variables
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)
329 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
330 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
331 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
332
333 ! If get_v depends on w, the first argument should be some average over the
334 ! left and right state
335 call rho_get_v(wlc, x, ixi^l, ixo^l, idim, cmax(ixi^s,1), .false.)
336
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)
340 else
341 cmax(ixo^s,1) = maxval(abs(cmax(ixo^s,1)))
342 end if
343
344 end subroutine rho_get_cbounds
345
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
352
353 dtnew = bigdouble
354 end subroutine rho_get_dt
355
356 ! There is nothing to add to the transport flux in the transport equation
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)
365
366 call rho_get_v(wc, x, ixi^l, ixo^l, idim, v, .false.)
367
368 f(ixo^s, rho_) = w(ixo^s, rho_) * v(ixo^s)
369 end subroutine rho_get_flux
370
371 subroutine rho_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT,wprim, w, x)
372
373 ! Add geometrical source terms to w
374 ! There are no geometrical source terms in the transport equation
375
377
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)
381
382 end subroutine rho_add_source_geom
383
384end module mod_rho_phys
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
integer coordinate
Definition mod_geometry.t:7
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...
Definition mod_physics.t:4
procedure(sub_convert), pointer phys_to_primitive
Definition mod_physics.t:54
procedure(sub_write_info), pointer phys_write_info
Definition mod_physics.t:75
procedure(sub_get_flux), pointer phys_get_flux
Definition mod_physics.t:61
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition mod_physics.t:60
procedure(sub_get_dt), pointer phys_get_dt
Definition mod_physics.t:64
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition mod_physics.t:65
procedure(sub_check_params), pointer phys_check_params
Definition mod_physics.t:51
integer, parameter flux_default
Indicates a normal flux.
Definition mod_physics.t:24
integer, dimension(:, :), allocatable flux_type
Array per direction per variable, which can be used to specify that certain fluxes have to be treated...
Definition mod_physics.t:21
procedure(sub_convert), pointer phys_to_conserved
Definition mod_physics.t:53
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition mod_physics.t:49
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:56
logical phys_energy
Solve energy equation or not.
Definition mod_physics.t:37
Module containing the physics routines for scalar advection.
Definition mod_rho_phys.t:2
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_
Definition mod_rho_phys.t:7
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...