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