MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_nonlinear_phys.t
Go to the documentation of this file.
1!> Module containing the physics routines for scalar nonlinear equation
3 use mod_comm_lib, only: mpistop
4
5 implicit none
6 private
7
8 !> index of the single scalar unknown
9 integer, protected, public :: rho_ = 1
10
11 !> switch between burgers (i.e. rho**2)
12 !> or nonconvex flux (i.e. rho**3)
13 integer, protected, public :: nonlinear_flux_type = 1
14
15 !> whether the KdV source term is added
16 logical, protected, public :: kdv_source_term = .false.
17
18 !> Whether particles module is added
19 logical, public, protected :: nonlinear_particles = .false.
20
21 ! Public methods
22 public :: nonlinear_phys_init
23 public :: nonlinear_get_v
24
25contains
26
27 !> Read this module's parameters from a file
28 subroutine nonlinear_params_read(files)
30 character(len=*), intent(in) :: files(:)
31 integer :: n
32
34
35 do n = 1, size(files)
36 open(unitpar, file=trim(files(n)), status='old')
37 read(unitpar, nonlinear_list, end=111)
38111 close(unitpar)
39 end do
40
41 end subroutine nonlinear_params_read
42
43 !> Write this module's parameters to a snapshot
44 subroutine nonlinear_write_info(fh)
45 ! for nonlinear scalar equation, nothing to write
46 ! note: this is info only stored at end of dat files,
47 ! is never read/used for restarts, only expects
48 ! an integer (number of parameters) and
49 ! corresponding double values and character names
50 ! and is meant for use in the python tools
52 integer, intent(in) :: fh
53 integer, dimension(MPI_STATUS_SIZE) :: st
54 integer :: er
55
56 ! Write zero parameters
57 call mpi_file_write(fh, 0, 1, mpi_integer, st, er)
58
59 end subroutine nonlinear_write_info
60
63 use mod_physics
64 use mod_kdv, only: kdv_init
65
66 call nonlinear_params_read(par_files)
67
68 physics_type = "nonlinear"
69 phys_energy = .false.
71
72 allocate(start_indices(number_species),stop_indices(number_species))
73 ! set the index of the first flux variable for species 1
74 start_indices(1)=1
75 rho_ = var_set_rho()
76
77 ! set number of variables which need update ghostcells
78 nwgc=nwflux
79
80 ! set the index of the last flux variable for species 1
81 stop_indices(1)=nwflux
82
83 ! Check whether custom flux types have been defined
84 if (.not. allocated(flux_type)) then
85 allocate(flux_type(ndir, nw))
87 else if (any(shape(flux_type) /= [ndir, nw])) then
88 call mpistop("phys_check error: flux_type has wrong shape")
89 end if
90
91 phys_get_cmax => nonlinear_get_cmax
92 phys_get_cbounds => nonlinear_get_cbounds
93 phys_get_flux => nonlinear_get_flux
94 phys_add_source_geom => nonlinear_add_source_geom
95 phys_add_source => nonlinear_add_source
96 phys_to_conserved => nonlinear_to_conserved
97 phys_to_primitive => nonlinear_to_primitive
98 phys_get_dt => nonlinear_get_dt
99 phys_check_params => nonlinear_check_params
100 phys_write_info => nonlinear_write_info
101
102 if (kdv_source_term) call kdv_init()
103
104
105 end subroutine nonlinear_phys_init
106
107 subroutine nonlinear_check_params
109 use mod_geometry, only: coordinate
111 use mod_particles, only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
112
113 ! Initialize particles module here, so user extra variables are sampled
114 if (nonlinear_particles) then
115 call particles_init()
116 end if
117
118 if(mype==0)then
119 write(*,*)'====NONLINEAR run with settings===================='
120 write(*,*)'Using mod_nonlinear_phys with settings:'
121 write(*,*)'Dimensionality :',ndim
122 write(*,*)'vector components:',ndir
123 write(*,*)'coordinate set to type,slab:',coordinate,slab
124 write(*,*)'number of variables nw=',nw
125 write(*,*)' start index iwstart=',iwstart
126 write(*,*)'number of vector variables=',nvector
127 write(*,*)'number of stagger variables nws=',nws
128 write(*,*)'number of variables with BCs=',nwgc
129 write(*,*)'number of vars with fluxes=',nwflux
130 write(*,*)'number of vars with flux + BC=',nwfluxbc
131 write(*,*)'number of auxiliary variables=',nwaux
132 write(*,*)'number of extra vars without flux=',nwextra
133 write(*,*)'number of extra vars for wextra=',nw_extra
134 write(*,*)'number of auxiliary I/O variables=',nwauxio
135 write(*,*)' nonlinear_particles=',nonlinear_particles
136 if(nonlinear_particles) then
137 write(*,*) '*****Using particles: npayload,ngridvars :', npayload,ngridvars
138 write(*,*) '*****Using particles: nusrpayload :', nusrpayload
139 write(*,*) '*****Using particles: num_particles :', num_particles
140 write(*,*) '*****Using particles: physics_type_particles=',physics_type_particles
141 end if
142 write(*,*)'number of ghostcells=',nghostcells
143 write(*,*)'==========================================='
144 endif
145
146 end subroutine nonlinear_check_params
147
148 subroutine nonlinear_to_conserved(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 nonlinear module)
155 end subroutine nonlinear_to_conserved
156
157 subroutine nonlinear_to_primitive(ixI^L, ixO^L, w, x)
159 integer, intent(in) :: ixi^l, ixo^l
160 double precision, intent(inout) :: w(ixi^s, nw)
161 double precision, intent(in) :: x(ixi^s, 1:^nd)
162
163 ! Do nothing (primitive and conservative are equal for nonlinear module)
164 end subroutine nonlinear_to_primitive
165
166 subroutine nonlinear_get_v(w, x, ixI^L, ixO^L, idim, v)
168 integer, intent(in) :: ixi^l, ixo^l, idim
169 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
170 double precision, intent(out) :: v(ixi^s)
171
172 select case(nonlinear_flux_type)
173 case(1)
174 v(ixo^s)=w(ixo^s,rho_)
175 case(2)
176 v(ixo^s)=3.0d0*w(ixo^s,rho_)**2
177 case default
178 call mpistop('Undefined fluxtype: set nonlinear_flux_type to 1 or 2')
179 end select
180
181 end subroutine nonlinear_get_v
182
183 subroutine nonlinear_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
185 integer, intent(in) :: ixi^l, ixo^l, idim
186 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:^nd)
187 double precision, intent(inout) :: cmax(ixi^s)
188
189 call nonlinear_get_v(w, x, ixi^l, ixo^l, idim, cmax)
190
191 cmax(ixo^s) = abs(cmax(ixo^s))
192
193 end subroutine nonlinear_get_cmax
194
195 subroutine nonlinear_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed, cmax, cmin)
197 use mod_variables
198 integer, intent(in) :: ixi^l, ixo^l, idim
199 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s,nw)
200 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s,nw)
201 double precision, intent(in) :: x(ixi^s, 1:^nd)
202 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
203 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
204 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
205 double precision :: wmean(ixi^s,nw)
206
207 ! since get_v depends on w, the first argument should be some average over the
208 ! left and right state
209 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
210 call nonlinear_get_v(wmean, x, ixi^l, ixo^l, idim, cmax(ixi^s,1))
211
212 if (present(cmin)) then
213 cmin(ixo^s,1) = min(cmax(ixo^s,1), zero)
214 cmax(ixo^s,1) = max(cmax(ixo^s,1), zero)
215 else
216 cmax(ixo^s,1) = maxval(abs(cmax(ixo^s,1)))
217 end if
218
219 end subroutine nonlinear_get_cbounds
220
221 subroutine nonlinear_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
223 use mod_kdv, only: kdv_get_dt
224
225 integer, intent(in) :: ixi^l, ixo^l
226 double precision, intent(in) :: dx^d, x(ixi^s, 1:^nd)
227 double precision, intent(in) :: w(ixi^s, 1:nw)
228 double precision, intent(inout) :: dtnew
229
230 dtnew = bigdouble
231
232 if(kdv_source_term) then
233 call kdv_get_dt(w, ixi^l, ixo^l, dtnew, dx^d, x)
234 endif
235 end subroutine nonlinear_get_dt
236
237 ! here we select the flux according to the nonlinear_flux_type parameter
238 subroutine nonlinear_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
240 integer, intent(in) :: ixi^l, ixo^l, idim
241 double precision, intent(in) :: wc(ixi^s, 1:nw)
242 double precision, intent(in) :: w(ixi^s, 1:nw)
243 double precision, intent(in) :: x(ixi^s, 1:^nd)
244 double precision, intent(out) :: f(ixi^s, nwflux)
245
246 select case(nonlinear_flux_type)
247 case(1)
248 f(ixo^s,rho_)=half*w(ixo^s,rho_)**2
249 case(2)
250 f(ixo^s,rho_)=w(ixo^s,rho_)**3
251 case default
252 call mpistop('Undefined fluxtype: set nonlinear_flux_type to 1 or 2')
253 end select
254
255 end subroutine nonlinear_get_flux
256
257 subroutine nonlinear_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT,wprim, w, x)
258
259 ! Add geometrical source terms to w
260 ! There are no geometrical source terms
261
263
264 integer, intent(in) :: ixi^l, ixo^l
265 double precision, intent(in) :: qdt, dtfactor, x(ixi^s, 1:^nd)
266 double precision, intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
267
268 end subroutine nonlinear_add_source_geom
269
270 subroutine nonlinear_add_source(qdt, dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
272 use mod_kdv, only: kdv_add_source
273
274 integer, intent(in) :: ixi^l, ixo^l
275 double precision, intent(in) :: qdt,dtfactor
276 double precision, intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:ndim)
277 double precision, intent(inout) :: w(ixi^s, 1:nw)
278 logical, intent(in) :: qsourcesplit
279 logical, intent(inout) :: active
280
281 if(kdv_source_term) then
282 call kdv_add_source(qdt,ixi^l,ixo^l,wct,w,x,qsourcesplit,active)
283 end if
284
285 end subroutine nonlinear_add_source
286
287end module mod_nonlinear_phys
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
integer coordinate
Definition mod_geometry.t:7
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 for including kdv source term in simulations adds over dimensions i.
Definition mod_kdv.t:3
subroutine kdv_add_source(qdt, ixil, ixol, wct, w, x, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
Definition mod_kdv.t:41
subroutine kdv_init()
Initialize the module.
Definition mod_kdv.t:33
subroutine kdv_get_dt(w, ixil, ixol, dtnew, dxd, x)
Definition mod_kdv.t:96
Module containing the physics routines for scalar nonlinear equation.
logical, public, protected nonlinear_particles
Whether particles module is added.
integer, public, protected nonlinear_flux_type
switch between burgers (i.e. rho**2) or nonconvex flux (i.e. rho**3)
logical, public, protected kdv_source_term
whether the KdV source term is added
integer, public, protected rho_
index of the single scalar unknown
subroutine, public nonlinear_phys_init()
subroutine, public nonlinear_get_v(w, x, ixil, ixol, idim, v)
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_add_source), pointer phys_add_source
Definition mod_physics.t:66
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
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
integer nwflux
Number of flux variables.