MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
Loading...
Searching...
No Matches
mod_physics.t
Go to the documentation of this file.
1!> This module defines the procedures of a physics module. It contains function
2!> pointers for the various supported routines. An actual physics module has to
3!> set these pointers to its implementation of these routines.
5 use mod_global_parameters, only: name_len, max_nw
8
9 implicit none
10 public
11
12
13 double precision :: phys_gamma=5.d0/3.d0
14
15 !> To use wider stencils in flux calculations. A value of 1 will extend it by
16 !> one cell in both directions, in any dimension
17 integer :: phys_wider_stencil = 0
18
19 !> Array per direction per variable, which can be used to specify that certain
20 !> fluxes have to be treated differently
21 integer, allocatable :: flux_type(:, :)
22
23 !> Indicates a normal flux
24 integer, parameter :: flux_default = 0
25 !> Indicates the flux should be treated with tvdlf
26 integer, parameter :: flux_tvdlf = 1
27 !> Indicates dissipation should be omitted
28 integer, parameter :: flux_no_dissipation = 2
29 !> Indicates the flux should be specially treated
30 integer, parameter :: flux_special = 3
31 !> Indicates the flux should be treated with hll
32 integer, parameter :: flux_hll = 4
33
34 !> Solve energy equation or not
35 logical :: phys_energy=.false.
36
37 !> Solve total energy equation or not
38 logical :: phys_total_energy=.false.
39
40 !> Solve internal energy instead of total energy
41 logical :: phys_internal_e=.false.
42
43 !> Solve partially ionized one-fluid plasma
44 logical :: phys_partial_ionization=.false.
45
46 !> if equilibrium pressure is splitted
47 logical :: phys_equi_pe=.false.
48
49 !> String describing the physics type of the simulation
50 character(len=name_len) :: physics_type = ""
51
52 procedure(sub_check_params), pointer :: phys_check_params => null()
53 procedure(sub_set_mg_bounds), pointer :: phys_set_mg_bounds => null()
54 procedure(sub_convert), pointer :: phys_to_conserved => null()
55 procedure(sub_convert), pointer :: phys_to_primitive => null()
56 procedure(sub_modify_wlr), pointer :: phys_modify_wlr => null()
57 procedure(sub_get_cmax), pointer :: phys_get_cmax => null()
58 procedure(sub_get_a2max), pointer :: phys_get_a2max => null()
59 procedure(sub_get_tcutoff), pointer :: phys_get_tcutoff => null()
60 procedure(sub_trac_after_setdt), pointer:: phys_trac_after_setdt => null()
61 procedure(sub_get_h_speed), pointer :: phys_get_h_speed => null()
62 procedure(sub_get_cbounds), pointer :: phys_get_cbounds => null()
63 procedure(sub_get_flux), pointer :: phys_get_flux => null()
64 procedure(sub_get_v), pointer :: phys_get_v => null()
65 procedure(sub_get_rho), pointer :: phys_get_rho => null()
66 procedure(sub_get_dt), pointer :: phys_get_dt => null()
67 procedure(sub_add_source_geom), pointer :: phys_add_source_geom => null()
68 procedure(sub_add_source), pointer :: phys_add_source => null()
69 procedure(sub_global_source), pointer :: phys_global_source_after => null()
70 procedure(sub_special_advance), pointer :: phys_special_advance => null()
71 procedure(sub_check_w), pointer :: phys_check_w => null()
72 procedure(sub_get_pthermal), pointer :: phys_get_pthermal => null()
73 procedure(sub_get_tgas), pointer :: phys_get_tgas => null()
74 procedure(sub_get_trad), pointer :: phys_get_trad => null()
75 procedure(sub_boundary_adjust), pointer :: phys_boundary_adjust => null()
76 procedure(sub_write_info), pointer :: phys_write_info => null()
77 procedure(sub_small_values), pointer :: phys_handle_small_values => null()
78 procedure(sub_get_ct_velocity), pointer :: phys_get_ct_velocity => null()
79 procedure(sub_update_faces), pointer :: phys_update_faces => null()
80 procedure(sub_face_to_center), pointer :: phys_face_to_center => null()
81 procedure(sub_implicit_update), pointer :: phys_implicit_update => null()
82 procedure(sub_evaluate_implicit),pointer:: phys_evaluate_implicit => null()
83 procedure(sub_clean_divb), pointer :: phys_clean_divb => null()
84 ! set the equilibrium variables
85 procedure(sub_set_equi_vars), pointer :: phys_set_equi_vars => null()
86 ! subroutine with no parameters which creates EUV images
87 procedure(sub_check_params), pointer :: phys_te_images => null()
88 ! to update temperature variable with partial ionization
89 procedure(sub_update_temperature), pointer :: phys_update_temperature => null()
90 procedure(sub_get_auxiliary), pointer :: phys_get_auxiliary => null()
91 procedure(sub_get_auxiliary_prim), pointer :: phys_get_auxiliary_prim => null()
92
93 abstract interface
94
96 end subroutine sub_check_params
97
101 end subroutine sub_set_mg_bounds
102
103 subroutine sub_boundary_adjust(igrid,psb)
105 integer, intent(in) :: igrid
106 type(state), target :: psb(max_blocks)
107 end subroutine sub_boundary_adjust
108
109 subroutine sub_convert(ixI^L, ixO^L, w, x)
111 integer, intent(in) :: ixI^L, ixO^L
112 double precision, intent(inout) :: w(ixI^S, nw)
113 double precision, intent(in) :: x(ixI^S, 1:^ND)
114 end subroutine sub_convert
115
116 subroutine sub_modify_wlr(ixI^L, ixO^L, qt, wLC, wRC, wLp, wRp, s, idir)
118 integer, intent(in) :: ixI^L, ixO^L, idir
119 double precision, intent(in) :: qt
120 double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
121 double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
122 type(state) :: s
123 end subroutine sub_modify_wlr
124
125 subroutine sub_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
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(inout) :: cmax(ixI^S)
130 end subroutine sub_get_cmax
131
132 subroutine sub_get_a2max(w, x, ixI^L, ixO^L, a2max)
134 integer, intent(in) :: ixI^L, ixO^L
135 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
136 double precision, intent(inout) :: a2max(ndim)
137 end subroutine sub_get_a2max
138
139 subroutine sub_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
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 double precision, intent(out) :: tco_local, Tmax_local
145 end subroutine sub_get_tcutoff
146
147 subroutine sub_trac_after_setdt(trac_alfa,tco,T_peak,T_bott)
148 double precision, intent(in) :: trac_alfa,tco,T_peak,T_bott
149 end subroutine sub_trac_after_setdt
150
151 subroutine sub_get_v(w,x,ixI^L,ixO^L,v)
153 integer, intent(in) :: ixI^L, ixO^L
154 double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
155 double precision, intent(out) :: v(ixI^S,1:ndir)
156 end subroutine sub_get_v
157
158 subroutine sub_get_rho(w,x,ixI^L,ixO^L,rho)
160 integer, intent(in) :: ixI^L, ixO^L
161 double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
162 double precision, intent(out) :: rho(ixI^S)
163 end subroutine sub_get_rho
164
165 subroutine sub_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
167 integer, intent(in) :: ixI^L, ixO^L, idim
168 double precision, intent(in) :: wprim(ixI^S, nw)
169 double precision, intent(in) :: x(ixI^S,1:ndim)
170 double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
171 end subroutine sub_get_h_speed
172
173 subroutine sub_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim, Hspeed, cmax, cmin)
175 use mod_variables
176 integer, intent(in) :: ixI^L, ixO^L, idim
177 double precision, intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
178 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
179 double precision, intent(in) :: x(ixI^S, 1:^ND)
180 double precision, intent(inout) :: cmax(ixI^S,1:number_species)
181 double precision, intent(inout), optional :: cmin(ixI^S,1:number_species)
182 double precision, intent(in) :: Hspeed(ixI^S,1:number_species)
183 end subroutine sub_get_cbounds
184
185 subroutine sub_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
187 integer, intent(in) :: ixI^L, ixO^L, idim
188 double precision, intent(in) :: wC(ixI^S, 1:nw)
189 double precision, intent(in) :: w(ixI^S, 1:nw)
190 double precision, intent(in) :: x(ixI^S, 1:^ND)
191 double precision, intent(out) :: f(ixI^S, nwflux)
192 end subroutine sub_get_flux
193
194 subroutine sub_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
196 integer, intent(in) :: ixI^L, ixO^L
197 double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
198 double precision, intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S, 1:nw), w(ixI^S, 1:nw)
199 end subroutine sub_add_source_geom
200
201 subroutine sub_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
202 qsourcesplit, active)
204 integer, intent(in) :: ixI^L, ixO^L
205 double precision, intent(in) :: qdt, dtfactor
206 double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
207 double precision, intent(inout) :: w(ixI^S, 1:nw)
208 logical, intent(in) :: qsourcesplit
209 logical, intent(inout) :: active !< Needs to be set to true when active
210 end subroutine sub_add_source
211
212 !> Add global source terms on complete domain (potentially implicit)
213 subroutine sub_global_source(qdt, qt, active)
215 double precision, intent(in) :: qdt !< Current time step
216 double precision, intent(in) :: qt !< Current time
217 logical, intent(inout) :: active !< Output if the source is active
218 end subroutine sub_global_source
219
220 !> clean initial divb
221 subroutine sub_clean_divb(qdt, qt, active)
223 double precision, intent(in) :: qdt !< Current time step
224 double precision, intent(in) :: qt !< Current time
225 logical, intent(inout) :: active !< Output if the source is active
226 end subroutine sub_clean_divb
227
228 !> set equilibrium variables other than b0 (e.g. p0 and rho0)
229 subroutine sub_set_equi_vars(igrid)
230 integer, intent(in) :: igrid
231 end subroutine sub_set_equi_vars
232
233
234 !> Add special advance in each advect step
235 subroutine sub_special_advance(qt, psa)
237 double precision, intent(in) :: qt !< Current time
238 type(state), target :: psa(max_blocks) !< Compute based on this state
239 end subroutine sub_special_advance
240
241 subroutine sub_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
243 integer, intent(in) :: ixI^L, ixO^L
244 double precision, intent(in) :: dx^D, x(ixI^S, 1:^ND)
245 double precision, intent(in) :: w(ixI^S, 1:nw)
246 double precision, intent(inout) :: dtnew
247 end subroutine sub_get_dt
248
249 subroutine sub_check_w(primitive, ixI^L, ixO^L, w, w_flag)
251 logical, intent(in) :: primitive
252 integer, intent(in) :: ixI^L, ixO^L
253 double precision, intent(in) :: w(ixI^S,1:nw)
254 logical, intent(inout) :: w_flag(ixI^S,1:nw)
255 end subroutine sub_check_w
256
257 subroutine sub_get_pthermal(w,x,ixI^L,ixO^L,pth)
259 integer, intent(in) :: ixI^L, ixO^L
260 double precision, intent(in) :: w(ixI^S,nw)
261 double precision, intent(in) :: x(ixI^S,1:ndim)
262 double precision, intent(out):: pth(ixI^S)
263 end subroutine sub_get_pthermal
264
265 subroutine sub_get_auxiliary(ixI^L,ixO^L,w,x)
267 integer, intent(in) :: ixI^L, ixO^L
268 double precision, intent(inout) :: w(ixI^S,nw)
269 double precision, intent(in) :: x(ixI^S,1:ndim)
270 end subroutine sub_get_auxiliary
271
272 subroutine sub_get_auxiliary_prim(ixI^L,ixO^L,w)
274 integer, intent(in) :: ixI^L, ixO^L
275 double precision, intent(inout) :: w(ixI^S,nw)
276 end subroutine sub_get_auxiliary_prim
277
278 subroutine sub_get_tgas(w,x,ixI^L,ixO^L,tgas)
280 integer, intent(in) :: ixI^L, ixO^L
281 double precision, intent(in) :: w(ixI^S,nw)
282 double precision, intent(in) :: x(ixI^S,1:ndim)
283 double precision, intent(out):: tgas(ixI^S)
284 end subroutine sub_get_tgas
285
286 subroutine sub_get_trad(w,x,ixI^L,ixO^L,trad)
288 integer, intent(in) :: ixI^L, ixO^L
289 double precision, intent(in) :: w(ixI^S,nw)
290 double precision, intent(in) :: x(ixI^S,1:ndim)
291 double precision, intent(out):: trad(ixI^S)
292 end subroutine sub_get_trad
293
294 subroutine sub_write_info(file_handle)
295 integer, intent(in) :: file_handle
296 end subroutine sub_write_info
297
298 subroutine sub_small_values(primitive, w, x, ixI^L, ixO^L, subname)
300 logical, intent(in) :: primitive
301 integer, intent(in) :: ixI^L,ixO^L
302 double precision, intent(inout) :: w(ixI^S,1:nw)
303 double precision, intent(in) :: x(ixI^S,1:ndim)
304 character(len=*), intent(in) :: subname
305 end subroutine sub_small_values
306
307 subroutine sub_get_var(ixI^L, ixO^L, w, out)
309 integer, intent(in) :: ixI^L, ixO^L
310 double precision, intent(in) :: w(ixI^S, nw)
311 double precision, intent(out) :: out(ixO^S)
312 end subroutine sub_get_var
313
314 subroutine sub_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
316
317 integer, intent(in) :: ixI^L, ixO^L, idim
318 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
319 double precision, intent(in) :: cmax(ixI^S)
320 double precision, intent(in), optional :: cmin(ixI^S)
321 type(ct_velocity), intent(inout):: vcts
322 end subroutine sub_get_ct_velocity
323
324 subroutine sub_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
326 integer, intent(in) :: ixI^L, ixO^L
327 double precision, intent(in) :: qt, qdt
328 ! cell-center primitive variables
329 double precision, intent(in) :: wprim(ixI^S,1:nw)
330 ! velocity structure
331 type(state) :: sCT, s
332 type(ct_velocity) :: vcts
333 double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
334 double precision, intent(inout) :: fE(ixI^S,sdim:3)
335 end subroutine sub_update_faces
336
337 subroutine sub_face_to_center(ixO^L,s)
339 integer, intent(in) :: ixO^L
340 type(state) :: s
341 end subroutine sub_face_to_center
342
343 subroutine sub_evaluate_implicit(qtC,psa)
345 type(state), target :: psa(max_blocks)
346 double precision, intent(in) :: qtC
347 end subroutine sub_evaluate_implicit
348
349 subroutine sub_implicit_update(dtfactor,qdt,qtC,psa,psb)
351 type(state), target :: psa(max_blocks)
352 type(state), target :: psb(max_blocks)
353 double precision, intent(in) :: qdt
354 double precision, intent(in) :: qtC
355 double precision, intent(in) :: dtfactor
356 end subroutine sub_implicit_update
357
358 subroutine sub_update_temperature(ixI^L,ixO^L,wCT,w,x)
360 integer, intent(in) :: ixI^L, ixO^L
361 double precision, intent(in) :: wCT(ixI^S,nw),x(ixI^S,1:ndim)
362 double precision, intent(inout) :: w(ixI^S,nw)
363 end subroutine sub_update_temperature
364
365 end interface
366
367contains
368
369 subroutine phys_check()
370 use mod_global_parameters, only: nw, ndir
371
374 use mod_comm_lib, only: mpistop
375
376 if (physics_type == "") call mpistop("Error: no physics module loaded")
377
378 call phys_hllc_check()
379 call phys_roe_check()
380
381 ! Checks whether the required physics methods have been defined
382 if (.not. associated(phys_check_params)) &
384
385 if (.not. associated(phys_to_conserved)) &
386 call mpistop("Error: phys_to_conserved not defined")
387
388 if (.not. associated(phys_to_primitive)) &
389 call mpistop("Error: phys_to_primitive not defined")
390
391 if (.not. associated(phys_modify_wlr)) &
393
394 if (.not. associated(phys_get_cmax)) &
395 call mpistop("Error: no phys_get_cmax not defined")
396
397 if (.not. associated(phys_get_a2max)) &
399
400 if (.not. associated(phys_get_h_speed)) &
402
403 if (.not. associated(phys_get_cbounds)) &
404 call mpistop("Error: no phys_get_cbounds not defined")
405
406 if (.not. associated(phys_get_flux)) &
407 call mpistop("Error: no phys_get_flux not defined")
408
409 if (.not. associated(phys_get_dt)) &
410 call mpistop("Error: no phys_get_dt not defined")
411
412 if (.not. associated(phys_add_source_geom)) &
414
415 if (.not. associated(phys_add_source)) &
417
418 if (.not. associated(phys_check_w)) &
420
421 if (.not. associated(phys_get_pthermal)) &
423 if (.not. associated(phys_get_auxiliary)) &
425 if (.not. associated(phys_get_auxiliary_prim)) &
427
428 if (.not. associated(phys_boundary_adjust)) &
430
431 if (.not. associated(phys_write_info)) &
433
434 if (.not. associated(phys_handle_small_values)) &
436
437 if (.not. associated(phys_get_ct_velocity)) &
439
440 if (.not. associated(phys_update_faces)) &
442
443 if (.not. associated(phys_face_to_center)) &
445
446 if (.not. associated(phys_evaluate_implicit)) &
448
449 if (.not. associated(phys_implicit_update)) &
451
452 end subroutine phys_check
453
455 end subroutine dummy_init_params
456
458 end subroutine dummy_check_params
459
460 subroutine dummy_modify_wlr(ixI^L, ixO^L, qt, wLC, wRC, wLp, wRp, s, idir)
462 integer, intent(in) :: ixI^L, ixO^L, idir
463 double precision, intent(in) :: qt
464 double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
465 double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
466 type(state) :: s
467 end subroutine dummy_modify_wlr
468
469 subroutine dummy_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
471 integer, intent(in) :: ixI^L, ixO^L, idim
472 double precision, intent(in) :: wprim(ixI^S, nw)
473 double precision, intent(in) :: x(ixI^S,1:ndim)
474 double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
475 end subroutine dummy_get_h_speed
476
477 subroutine dummy_get_a2max(w, x, ixI^L, ixO^L, a2max)
479 use mod_comm_lib, only: mpistop
480 integer, intent(in) :: ixI^L, ixO^L
481 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
482 double precision, intent(inout) :: a2max(ndim)
483 call mpistop("Error: entered dummy_get_a2max")
484 end subroutine dummy_get_a2max
485
486 subroutine dummy_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
488 integer, intent(in) :: ixI^L, ixO^L
489 double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
490 double precision, intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S,1:nw),w(ixI^S, 1:nw)
491 end subroutine dummy_add_source_geom
492
493 subroutine dummy_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
494 qsourcesplit, active)
496 integer, intent(in) :: ixI^L, ixO^L
497 double precision, intent(in) :: qdt,dtfactor
498 double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
499 double precision, intent(inout) :: w(ixI^S, 1:nw)
500 logical, intent(in) :: qsourcesplit
501 logical, intent(inout) :: active
502 ! Don't have to set active, since it starts as .false.
503 end subroutine dummy_add_source
504
505 subroutine dummy_check_w(primitive, ixI^L, ixO^L, w, w_flag)
507 logical, intent(in) :: primitive
508 integer, intent(in) :: ixI^L, ixO^L
509 double precision, intent(in) :: w(ixI^S,1:nw)
510 logical, intent(inout) :: w_flag(ixI^S,1:nw)
511
512 w_flag=.false. ! All okay
513 end subroutine dummy_check_w
514
515 subroutine dummy_get_pthermal(w, x, ixI^L, ixO^L, pth)
517 use mod_comm_lib, only: mpistop
518
519 integer, intent(in) :: ixI^L, ixO^L
520 double precision, intent(in) :: w(ixI^S, nw)
521 double precision, intent(in) :: x(ixI^S, 1:ndim)
522 double precision, intent(out):: pth(ixI^S)
523
524 call mpistop("No get_pthermal method specified")
525 end subroutine dummy_get_pthermal
526
527 subroutine dummy_get_auxiliary(ixI^L, ixO^L, w, x)
529 use mod_comm_lib, only: mpistop
530
531 integer, intent(in) :: ixI^L, ixO^L
532 double precision, intent(inout) :: w(ixI^S, nw)
533 double precision, intent(in) :: x(ixI^S, 1:ndim)
534
535 !call mpistop("No get_auxiliary method specified")
536 end subroutine dummy_get_auxiliary
537
538 subroutine dummy_get_auxiliary_prim(ixI^L, ixO^L, w)
540 use mod_comm_lib, only: mpistop
541
542 integer, intent(in) :: ixI^L, ixO^L
543 double precision, intent(inout) :: w(ixI^S, nw)
544
545 call mpistop("No get_auxiliary_prim method specified")
546 end subroutine dummy_get_auxiliary_prim
547
548 subroutine dummy_boundary_adjust(igrid,psb)
550 integer, intent(in) :: igrid
551 type(state), target :: psb(max_blocks)
552 end subroutine dummy_boundary_adjust
553
554 subroutine dummy_write_info(fh)
556 integer, intent(in) :: fh !< File handle
557 integer, dimension(MPI_STATUS_SIZE) :: st
558 integer :: er
559
560 ! Number of physics parameters
561 integer, parameter :: n_par = 0
562
563 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
564 end subroutine dummy_write_info
565
566 subroutine dummy_small_values(primitive, w, x, ixI^L, ixO^L, subname)
568 logical, intent(in) :: primitive
569 integer, intent(in) :: ixI^L,ixO^L
570 double precision, intent(inout) :: w(ixI^S,1:nw)
571 double precision, intent(in) :: x(ixI^S,1:ndim)
572 character(len=*), intent(in) :: subname
573 end subroutine dummy_small_values
574
575 subroutine dummy_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
577
578 integer, intent(in) :: ixI^L, ixO^L, idim
579 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
580 double precision, intent(in) :: cmax(ixI^S)
581 double precision, intent(in), optional :: cmin(ixI^S)
582 type(ct_velocity), intent(inout):: vcts
583 end subroutine dummy_get_ct_velocity
584
585 subroutine dummy_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
587 integer, intent(in) :: ixI^L, ixO^L
588 double precision, intent(in) :: qt, qdt
589 ! cell-center primitive variables
590 double precision, intent(in) :: wprim(ixI^S,1:nw)
591 type(state) :: sCT, s
592 type(ct_velocity) :: vcts
593 double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
594 double precision, intent(inout) :: fE(ixI^S,sdim:3)
595 end subroutine dummy_update_faces
596
597 subroutine dummy_face_to_center(ixO^L,s)
599 integer, intent(in) :: ixO^L
600 type(state) :: s
601 end subroutine dummy_face_to_center
602
603 subroutine dummy_evaluate_implicit(qtC,psa)
605 type(state), target :: psa(max_blocks)
606 double precision, intent(in) :: qtC
607 integer :: iigrid, igrid
608
609 ! Just copy in nul state
610 !$OMP PARALLEL DO PRIVATE(igrid)
611 do iigrid=1,igridstail; igrid=igrids(iigrid);
612 psa(igrid)%w = 0.0d0*psa(igrid)%w
613 if(stagger_grid) psa(igrid)%ws = 0.0d0*psa(igrid)%ws
614 end do
615 !$OMP END PARALLEL DO
616
617 end subroutine dummy_evaluate_implicit
618
619 subroutine dummy_implicit_update(dtfactor,qdt,qtC,psa,psb)
621 type(state), target :: psa(max_blocks)
622 type(state), target :: psb(max_blocks)
623 double precision, intent(in) :: qdt
624 double precision, intent(in) :: qtC
625 double precision, intent(in) :: dtfactor
626 integer :: iigrid, igrid
627
628 ! Just copy in psb state when using the scheme without implicit part
629 !$OMP PARALLEL DO PRIVATE(igrid)
630 do iigrid=1,igridstail; igrid=igrids(iigrid);
631 psa(igrid)%w = psb(igrid)%w
632 if(stagger_grid) psa(igrid)%ws = psb(igrid)%ws
633 end do
634 !$OMP END PARALLEL DO
635
636 end subroutine dummy_implicit_update
637
638end module mod_physics
Add global source terms on complete domain (potentially implicit)
set equilibrium variables other than b0 (e.g. p0 and rho0)
Add special advance in each advect step.
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
This module contains definitions of global parameters and variables and some generic functions/subrou...
logical stagger_grid
True for using stagger grid.
integer ndir
Number of spatial dimensions (components) for vector variables.
subroutine phys_hllc_check
subroutine phys_roe_check()
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
procedure(sub_get_a2max), pointer phys_get_a2max
Definition mod_physics.t:58
procedure(sub_get_ct_velocity), pointer phys_get_ct_velocity
Definition mod_physics.t:78
subroutine dummy_implicit_update(dtfactor, qdt, qtc, psa, psb)
subroutine dummy_boundary_adjust(igrid, psb)
procedure(sub_convert), pointer phys_to_primitive
Definition mod_physics.t:55
procedure(sub_get_tgas), pointer phys_get_tgas
Definition mod_physics.t:73
integer, parameter flux_tvdlf
Indicates the flux should be treated with tvdlf.
Definition mod_physics.t:26
procedure(sub_small_values), pointer phys_handle_small_values
Definition mod_physics.t:77
procedure(sub_write_info), pointer phys_write_info
Definition mod_physics.t:76
procedure(sub_get_flux), pointer phys_get_flux
Definition mod_physics.t:63
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
Definition mod_physics.t:82
subroutine dummy_get_auxiliary(ixil, ixol, w, x)
subroutine dummy_face_to_center(ixol, s)
subroutine dummy_modify_wlr(ixil, ixol, qt, wlc, wrc, wlp, wrp, s, idir)
subroutine dummy_get_a2max(w, x, ixil, ixol, a2max)
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition mod_physics.t:62
procedure(sub_check_params), pointer phys_te_images
Definition mod_physics.t:87
procedure(sub_get_dt), pointer phys_get_dt
Definition mod_physics.t:66
subroutine dummy_update_faces(ixil, ixol, qt, qdt, wprim, fc, fe, sct, s, vcts)
logical phys_partial_ionization
Solve partially ionized one-fluid plasma.
Definition mod_physics.t:44
logical phys_total_energy
Solve total energy equation or not.
Definition mod_physics.t:38
integer, parameter flux_hll
Indicates the flux should be treated with hll.
Definition mod_physics.t:32
procedure(sub_get_pthermal), pointer phys_get_pthermal
Definition mod_physics.t:72
procedure(sub_update_temperature), pointer phys_update_temperature
Definition mod_physics.t:89
procedure(sub_get_tcutoff), pointer phys_get_tcutoff
Definition mod_physics.t:59
subroutine dummy_get_pthermal(w, x, ixil, ixol, pth)
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition mod_physics.t:67
procedure(sub_check_params), pointer phys_check_params
Definition mod_physics.t:52
integer, parameter flux_default
Indicates a normal flux.
Definition mod_physics.t:24
procedure(sub_get_auxiliary), pointer phys_get_auxiliary
Definition mod_physics.t:90
subroutine dummy_get_h_speed(wprim, x, ixil, ixol, idim, hspeed)
subroutine dummy_small_values(primitive, w, x, ixil, ixol, subname)
procedure(sub_get_auxiliary_prim), pointer phys_get_auxiliary_prim
Definition mod_physics.t:91
integer, parameter flux_special
Indicates the flux should be specially treated.
Definition mod_physics.t:30
procedure(sub_trac_after_setdt), pointer phys_trac_after_setdt
Definition mod_physics.t:60
integer, parameter flux_no_dissipation
Indicates dissipation should be omitted.
Definition mod_physics.t:28
subroutine dummy_add_source_geom(qdt, dtfactor, ixil, ixol, wct, wprim, w, x)
procedure(sub_set_equi_vars), pointer phys_set_equi_vars
Definition mod_physics.t:85
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
integer phys_wider_stencil
To use wider stencils in flux calculations. A value of 1 will extend it by one cell in both direction...
Definition mod_physics.t:17
subroutine dummy_init_params
procedure(sub_update_faces), pointer phys_update_faces
Definition mod_physics.t:79
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_check_w), pointer phys_check_w
Definition mod_physics.t:71
subroutine dummy_get_ct_velocity(vcts, wlp, wrp, ixil, ixol, idim, cmax, cmin)
procedure(sub_set_mg_bounds), pointer phys_set_mg_bounds
Definition mod_physics.t:53
logical phys_equi_pe
if equilibrium pressure is splitted
Definition mod_physics.t:47
procedure(sub_implicit_update), pointer phys_implicit_update
Definition mod_physics.t:81
procedure(sub_get_trad), pointer phys_get_trad
Definition mod_physics.t:74
subroutine dummy_evaluate_implicit(qtc, psa)
subroutine dummy_check_params
procedure(sub_clean_divb), pointer phys_clean_divb
Definition mod_physics.t:83
procedure(sub_boundary_adjust), pointer phys_boundary_adjust
Definition mod_physics.t:75
subroutine dummy_add_source(qdt, dtfactor, ixil, ixol, wct, wctprim, w, x, qsourcesplit, active)
double precision phys_gamma
Definition mod_physics.t:13
procedure(sub_global_source), pointer phys_global_source_after
Definition mod_physics.t:69
procedure(sub_add_source), pointer phys_add_source
Definition mod_physics.t:68
procedure(sub_face_to_center), pointer phys_face_to_center
Definition mod_physics.t:80
procedure(sub_modify_wlr), pointer phys_modify_wlr
Definition mod_physics.t:56
subroutine dummy_write_info(fh)
procedure(sub_get_h_speed), pointer phys_get_h_speed
Definition mod_physics.t:61
subroutine dummy_get_auxiliary_prim(ixil, ixol, w)
logical phys_internal_e
Solve internal energy instead of total energy.
Definition mod_physics.t:41
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:57
procedure(sub_get_rho), pointer phys_get_rho
Definition mod_physics.t:65
procedure(sub_get_v), pointer phys_get_v
Definition mod_physics.t:64
subroutine dummy_check_w(primitive, ixil, ixol, w, w_flag)
subroutine phys_check()
logical phys_energy
Solve energy equation or not.
Definition mod_physics.t:35
procedure(sub_special_advance), pointer phys_special_advance
Definition mod_physics.t:70
Module with all the methods that users can customize in AMRVAC.