MPI-AMRVAC 3.1
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
60
61 call rho_params_read(par_files)
62
63 physics_type = "rho"
64 phys_energy = .false.
66
67 allocate(start_indices(number_species),stop_indices(number_species))
68 ! set the index of the first flux variable for species 1
69 start_indices(1)=1
70 rho_ = var_set_rho()
71
72 ! set number of variables which need update ghostcells
73 nwgc=nwflux
74
75 ! set the index of the last flux variable for species 1
76 stop_indices(1)=nwflux
77
78 ! Check whether custom flux types have been defined
79 if (.not. allocated(flux_type)) then
80 allocate(flux_type(ndir, nw))
82 else if (any(shape(flux_type) /= [ndir, nw])) then
83 call mpistop("phys_check error: flux_type has wrong shape")
84 end if
85
86 phys_get_cmax => rho_get_cmax
87 phys_get_cbounds => rho_get_cbounds
88 phys_get_flux => rho_get_flux
89 phys_add_source_geom => rho_add_source_geom
90 phys_to_conserved => rho_to_conserved
91 phys_to_primitive => rho_to_primitive
92 phys_get_dt => rho_get_dt
93 phys_write_info => rho_write_info
94
95 ! Initialize particles module
96 if (rho_particles) then
97 call particles_init()
98 end if
99
100 end subroutine rho_phys_init
101
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)
107
108 ! Do nothing (primitive and conservative are equal for rho module)
109 end subroutine rho_to_conserved
110
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)
116
117 ! Do nothing (primitive and conservative are equal for rho module)
118 end subroutine rho_to_primitive
119
120 subroutine rho_get_v(w, x, ixI^L, ixO^L, idim, v, centered)
122 use mod_geometry
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)
127
128 double precision :: dtheta, dphi, halfdtheta, halfdphi, invdtheta, invdphi
129 {^ifthreed
130 double precision :: appcosphi(ixi^s), appsinphi(ixi^s), &
131 appcosthe(ixi^s), appsinthe(ixi^s)
132 }
133
134 select case (coordinate)
135 case (cylindrical)
136 {^ifoned
137 call mpistop("advection in 1D cylindrical not available")
138 }
139 {^iftwod
140 ! advection in 2D cylindrical: polar grid: to v_R, v_varphi
141 if(centered)then
142 select case (idim)
143 case (1) ! radial velocity
144 v(ixo^s) = rho_v(1)*dcos(x(ixo^s,2))+rho_v(2)*dsin(x(ixo^s,2))
145 case (2) ! v_varphi
146 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,2))+rho_v(2)*dcos(x(ixo^s,2))
147 end select
148 else
149 ! assumed uniform in varphi=theta
150 dtheta=x(ixomin1,ixomin2+1,2)-x(ixomin1,ixomin2,2)
151 halfdtheta=0.5d0*dtheta
152 invdtheta=1.0d0/dtheta
153 select case (idim)
154 case (1) ! radial velocity
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
159 case (2) ! v_varphi
160 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,2)+halfdtheta) &
161 +rho_v(2)*dcos(x(ixo^s,2)+halfdtheta)
162 end select
163 endif
164 }
165 {^ifthreed
166 ! advection in 3D cylindrical: convert to v_R, v_Z, v_varphi
167 if(centered)then
168 select case (idim)
169 case (1) ! v_R velocity
170 v(ixo^s) = rho_v(1)*dcos(x(ixo^s,3))+rho_v(2)*dsin(x(ixo^s,3))
171 case (2) ! v_Z velocity
172 v(ixo^s) = rho_v(3)
173 case (3) ! v_varphi velocity
174 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,3))+rho_v(2)*dcos(x(ixo^s,3))
175 end select
176 else
177 ! assumed uniform in varphi=theta
178 dtheta=x(ixomin1,ixomin2,ixomin3+1,3)-x(ixomin1,ixomin2,ixomin3,3)
179 halfdtheta=0.5d0*dtheta
180 invdtheta=1.0d0/dtheta
181 select case (idim)
182 case (1) ! radial velocity
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
187 case (2) ! v_Z velocity
188 v(ixo^s) = rho_v(3)
189 case (3) ! v_varphi
190 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,3)+halfdtheta) &
191 +rho_v(2)*dcos(x(ixo^s,3)+halfdtheta)
192 end select
193 endif
194 }
195 case (spherical)
196 {^ifoned
197 call mpistop("advection in 1D spherical not available")
198 }
199 {^iftwod
200 call mpistop("advection in 2D spherical not available")
201 }
202 {^ifthreed
203 ! advection in 3D spherical: convert to v_r, v_theta, v_phi
204 if(centered)then
205 select case (idim)
206 case (1) ! radial velocity
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))
210 case (2) ! theta velocity
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))
214 case (3) ! phi velocity
215 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,3)) &
216 +rho_v(2)*dcos(x(ixo^s,3))
217 end select
218 else
219 ! assumed uniform in theta and phi
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
223 halfdphi=0.5d0*dphi
224 invdtheta=1.0d0/dtheta
225 invdphi=1.0d0/dphi
226 select case (idim)
227 case (1) ! radial velocity
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))
235 appsinthe(ixo^s)= &
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)
242 case (2) ! theta velocity
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)
250 case (3) ! phi velocity
251 v(ixo^s) =-rho_v(1)*dsin(x(ixo^s,3)+halfdphi) &
252 +rho_v(2)*dcos(x(ixo^s,3)+halfdphi)
253 end select
254 endif
255 }
256 case default
257 v(ixo^s) = rho_v(idim)
258 end select
259 end subroutine rho_get_v
260
261 !> Calculate simple v component
262 subroutine rho_get_v_idim(w,x,ixI^L,ixO^L,idim,v)
264
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)
268
269 v(ixo^s) = rho_v(idim)
270
271 end subroutine rho_get_v_idim
272
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)
278
279 call rho_get_v(w, x, ixi^l, ixo^l, idim, cmax, .true.)
280
281 cmax(ixo^s) = abs(cmax(ixo^s))
282
283 end subroutine rho_get_cmax
284
285 subroutine rho_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
287 use mod_variables
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)
292 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
293 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
294 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
295
296 ! If get_v depends on w, the first argument should be some average over the
297 ! left and right state
298 call rho_get_v(wlc, x, ixi^l, ixo^l, idim, cmax(ixi^s,1), .false.)
299
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)
303 else
304 cmax(ixo^s,1) = maxval(abs(cmax(ixo^s,1)))
305 end if
306
307 end subroutine rho_get_cbounds
308
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
315
316 dtnew = bigdouble
317 end subroutine rho_get_dt
318
319 ! There is nothing to add to the transport flux in the transport equation
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)
328
329 call rho_get_v(wc, x, ixi^l, ixo^l, idim, v, .false.)
330
331 f(ixo^s, rho_) = w(ixo^s, rho_) * v(ixo^s)
332 end subroutine rho_get_flux
333
334 subroutine rho_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT,wprim, w, x)
335
336 ! Add geometrical source terms to w
337 ! There are no geometrical source terms in the transport equation
338
340
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)
344
345 end subroutine rho_add_source_geom
346
347end 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.
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...
Definition mod_physics.t:4
procedure(sub_convert), pointer phys_to_primitive
Definition mod_physics.t:55
procedure(sub_write_info), pointer phys_write_info
Definition mod_physics.t:77
procedure(sub_get_flux), pointer phys_get_flux
Definition mod_physics.t:64
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition mod_physics.t:63
procedure(sub_get_dt), pointer phys_get_dt
Definition mod_physics.t:67
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition mod_physics.t:68
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:54
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition mod_physics.t:50
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:57
logical phys_energy
Solve energy equation or not.
Definition mod_physics.t:35
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...