MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
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_cs2max), pointer :: phys_get_cs2max => null()
64  procedure(sub_get_tcutoff), pointer :: phys_get_tcutoff => null()
65  procedure(sub_trac_after_setdt), pointer:: phys_trac_after_setdt => null()
66  procedure(sub_get_h_speed), pointer :: phys_get_h_speed => null()
67  procedure(sub_get_cbounds), pointer :: phys_get_cbounds => null()
68  procedure(sub_get_flux), pointer :: phys_get_flux => null()
69  procedure(sub_get_v), pointer :: phys_get_v => null()
70  procedure(sub_get_rho), pointer :: phys_get_rho => null()
71  procedure(sub_get_dt), pointer :: phys_get_dt => null()
72  procedure(sub_add_source_geom), pointer :: phys_add_source_geom => null()
73  procedure(sub_add_source), pointer :: phys_add_source => null()
74  procedure(sub_global_source), pointer :: phys_global_source_after => null()
75  procedure(sub_special_advance), pointer :: phys_special_advance => null()
76  procedure(sub_check_w), pointer :: phys_check_w => null()
77  procedure(sub_get_pthermal), pointer :: phys_get_pthermal => null()
78  procedure(sub_get_tgas), pointer :: phys_get_tgas => null()
79  procedure(sub_get_trad), pointer :: phys_get_trad => null()
80  procedure(sub_boundary_adjust), pointer :: phys_boundary_adjust => null()
81  procedure(sub_write_info), pointer :: phys_write_info => null()
82  procedure(sub_small_values), pointer :: phys_handle_small_values => null()
83  procedure(sub_get_ct_velocity), pointer :: phys_get_ct_velocity => null()
84  procedure(sub_update_faces), pointer :: phys_update_faces => null()
85  procedure(sub_face_to_center), pointer :: phys_face_to_center => null()
86  procedure(sub_implicit_update), pointer :: phys_implicit_update => null()
87  procedure(sub_evaluate_implicit),pointer:: phys_evaluate_implicit => null()
88  procedure(sub_clean_divb), pointer :: phys_clean_divb => null()
89  ! set the equilibrium variables
90  procedure(sub_set_equi_vars), pointer :: phys_set_equi_vars => null()
91  ! subroutine with no parameters which creates EUV images
92  procedure(sub_check_params), pointer :: phys_te_images => null()
93  ! to update temperature variable with partial ionization
94  procedure(sub_update_temperature), pointer :: phys_update_temperature => null()
95  procedure(sub_get_auxiliary), pointer :: phys_get_auxiliary => null()
96  procedure(sub_get_auxiliary_prim), pointer :: phys_get_auxiliary_prim => null()
97 
98  abstract interface
99 
100  subroutine sub_check_params
101  end subroutine sub_check_params
102 
103  subroutine sub_set_mg_bounds
105  use mod_usr_methods
106  end subroutine sub_set_mg_bounds
107 
108  subroutine sub_boundary_adjust(igrid,psb)
110  integer, intent(in) :: igrid
111  type(state), target :: psb(max_blocks)
112  end subroutine sub_boundary_adjust
113 
114  subroutine sub_convert(ixI^L, ixO^L, w, x)
116  integer, intent(in) :: ixI^L, ixO^L
117  double precision, intent(inout) :: w(ixI^S, nw)
118  double precision, intent(in) :: x(ixI^S, 1:^ND)
119  end subroutine sub_convert
120 
121  subroutine sub_modify_wlr(ixI^L, ixO^L, qt, wLC, wRC, wLp, wRp, s, idir)
123  integer, intent(in) :: ixI^L, ixO^L, idir
124  double precision, intent(in) :: qt
125  double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
126  double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
127  type(state) :: s
128  end subroutine sub_modify_wlr
129 
130  subroutine sub_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
132  integer, intent(in) :: ixI^L, ixO^L, idim
133  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
134  double precision, intent(inout) :: cmax(ixI^S)
135  end subroutine sub_get_cmax
136 
137  subroutine sub_get_a2max(w, x, ixI^L, ixO^L, a2max)
139  integer, intent(in) :: ixI^L, ixO^L
140  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
141  double precision, intent(inout) :: a2max(ndim)
142  end subroutine sub_get_a2max
143 
144  subroutine sub_get_cs2max(w, x, ixI^L, ixO^L, cs2max)
146  integer, intent(in) :: ixI^L, ixO^L
147  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
148  double precision, intent(inout) :: cs2max
149  end subroutine sub_get_cs2max
150 
151  subroutine sub_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
153  integer, intent(in) :: ixI^L, ixO^L
154  double precision, intent(inout) :: w(ixI^S, nw)
155  double precision, intent(in) :: x(ixI^S, 1:^ND)
156  double precision, intent(out) :: tco_local, Tmax_local
157  end subroutine sub_get_tcutoff
158 
159  subroutine sub_trac_after_setdt(trac_alfa,tco,T_peak,T_bott)
160  double precision, intent(in) :: trac_alfa,tco,T_peak,T_bott
161  end subroutine sub_trac_after_setdt
162 
163  subroutine sub_get_v(w,x,ixI^L,ixO^L,v)
165  integer, intent(in) :: ixI^L, ixO^L
166  double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
167  double precision, intent(out) :: v(ixI^S,1:ndir)
168  end subroutine sub_get_v
169 
170  subroutine sub_get_rho(w,x,ixI^L,ixO^L,rho)
172  integer, intent(in) :: ixI^L, ixO^L
173  double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
174  double precision, intent(out) :: rho(ixI^S)
175  end subroutine sub_get_rho
176 
177  subroutine sub_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
179  integer, intent(in) :: ixI^L, ixO^L, idim
180  double precision, intent(in) :: wprim(ixI^S, nw)
181  double precision, intent(in) :: x(ixI^S,1:ndim)
182  double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
183  end subroutine sub_get_h_speed
184 
185  subroutine sub_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim, Hspeed, cmax, cmin)
187  use mod_variables
188  integer, intent(in) :: ixI^L, ixO^L, idim
189  double precision, intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
190  double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
191  double precision, intent(in) :: x(ixI^S, 1:^ND)
192  double precision, intent(inout) :: cmax(ixI^S,1:number_species)
193  double precision, intent(inout), optional :: cmin(ixI^S,1:number_species)
194  double precision, intent(in) :: Hspeed(ixI^S,1:number_species)
195  end subroutine sub_get_cbounds
196 
197  subroutine sub_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
199  integer, intent(in) :: ixI^L, ixO^L, idim
200  double precision, intent(in) :: wC(ixI^S, 1:nw)
201  double precision, intent(in) :: w(ixI^S, 1:nw)
202  double precision, intent(in) :: x(ixI^S, 1:^ND)
203  double precision, intent(out) :: f(ixI^S, nwflux)
204  end subroutine sub_get_flux
205 
206  subroutine sub_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, w, x)
208  integer, intent(in) :: ixI^L, ixO^L
209  double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
210  double precision, intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
211  end subroutine sub_add_source_geom
212 
213  subroutine sub_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
214  qsourcesplit, active)
216  integer, intent(in) :: ixI^L, ixO^L
217  double precision, intent(in) :: qdt, dtfactor
218  double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
219  double precision, intent(inout) :: w(ixI^S, 1:nw)
220  logical, intent(in) :: qsourcesplit
221  logical, intent(inout) :: active !< Needs to be set to true when active
222  end subroutine sub_add_source
223 
224  !> Add global source terms on complete domain (potentially implicit)
225  subroutine sub_global_source(qdt, qt, active)
227  double precision, intent(in) :: qdt !< Current time step
228  double precision, intent(in) :: qt !< Current time
229  logical, intent(inout) :: active !< Output if the source is active
230  end subroutine sub_global_source
231 
232  !> clean initial divb
233  subroutine sub_clean_divb(qdt, qt, active)
235  double precision, intent(in) :: qdt !< Current time step
236  double precision, intent(in) :: qt !< Current time
237  logical, intent(inout) :: active !< Output if the source is active
238  end subroutine sub_clean_divb
239 
240  !> set equilibrium variables other than b0 (e.g. p0 and rho0)
241  subroutine sub_set_equi_vars(igrid)
242  integer, intent(in) :: igrid
243  end subroutine sub_set_equi_vars
244 
245 
246  !> Add special advance in each advect step
247  subroutine sub_special_advance(qt, psa)
249  double precision, intent(in) :: qt !< Current time
250  type(state), target :: psa(max_blocks) !< Compute based on this state
251  end subroutine sub_special_advance
252 
253  subroutine sub_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
255  integer, intent(in) :: ixI^L, ixO^L
256  double precision, intent(in) :: dx^D, x(ixI^S, 1:^ND)
257  double precision, intent(in) :: w(ixI^S, 1:nw)
258  double precision, intent(inout) :: dtnew
259  end subroutine sub_get_dt
260 
261  subroutine sub_check_w(primitive, ixI^L, ixO^L, w, w_flag)
263  logical, intent(in) :: primitive
264  integer, intent(in) :: ixI^L, ixO^L
265  double precision, intent(in) :: w(ixI^S,1:nw)
266  logical, intent(inout) :: w_flag(ixI^S,1:nw)
267  end subroutine sub_check_w
268 
269  subroutine sub_get_pthermal(w,x,ixI^L,ixO^L,pth)
271  integer, intent(in) :: ixI^L, ixO^L
272  double precision, intent(in) :: w(ixI^S,nw)
273  double precision, intent(in) :: x(ixI^S,1:ndim)
274  double precision, intent(out):: pth(ixI^S)
275  end subroutine sub_get_pthermal
276 
277  subroutine sub_get_auxiliary(ixI^L,ixO^L,w,x)
279  integer, intent(in) :: ixI^L, ixO^L
280  double precision, intent(inout) :: w(ixI^S,nw)
281  double precision, intent(in) :: x(ixI^S,1:ndim)
282  end subroutine sub_get_auxiliary
283 
284  subroutine sub_get_auxiliary_prim(ixI^L,ixO^L,w)
286  integer, intent(in) :: ixI^L, ixO^L
287  double precision, intent(inout) :: w(ixI^S,nw)
288  end subroutine sub_get_auxiliary_prim
289 
290  subroutine sub_get_tgas(w,x,ixI^L,ixO^L,tgas)
292  integer, intent(in) :: ixI^L, ixO^L
293  double precision, intent(in) :: w(ixI^S,nw)
294  double precision, intent(in) :: x(ixI^S,1:ndim)
295  double precision, intent(out):: tgas(ixI^S)
296  end subroutine sub_get_tgas
297 
298  subroutine sub_get_trad(w,x,ixI^L,ixO^L,trad)
300  integer, intent(in) :: ixI^L, ixO^L
301  double precision, intent(in) :: w(ixI^S,nw)
302  double precision, intent(in) :: x(ixI^S,1:ndim)
303  double precision, intent(out):: trad(ixI^S)
304  end subroutine sub_get_trad
305 
306  subroutine sub_write_info(file_handle)
307  integer, intent(in) :: file_handle
308  end subroutine sub_write_info
309 
310  subroutine sub_small_values(primitive, w, x, ixI^L, ixO^L, subname)
312  logical, intent(in) :: primitive
313  integer, intent(in) :: ixI^L,ixO^L
314  double precision, intent(inout) :: w(ixI^S,1:nw)
315  double precision, intent(in) :: x(ixI^S,1:ndim)
316  character(len=*), intent(in) :: subname
317  end subroutine sub_small_values
318 
319  subroutine sub_get_var(ixI^L, ixO^L, w, out)
321  integer, intent(in) :: ixI^L, ixO^L
322  double precision, intent(in) :: w(ixI^S, nw)
323  double precision, intent(out) :: out(ixO^S)
324  end subroutine sub_get_var
325 
326  subroutine sub_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
328 
329  integer, intent(in) :: ixI^L, ixO^L, idim
330  double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
331  double precision, intent(in) :: cmax(ixI^S)
332  double precision, intent(in), optional :: cmin(ixI^S)
333  type(ct_velocity), intent(inout):: vcts
334  end subroutine sub_get_ct_velocity
335 
336  subroutine sub_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
338  integer, intent(in) :: ixI^L, ixO^L
339  double precision, intent(in) :: qt, qdt
340  ! cell-center primitive variables
341  double precision, intent(in) :: wprim(ixI^S,1:nw)
342  ! velocity structure
343  type(state) :: sCT, s
344  type(ct_velocity) :: vcts
345  double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
346  double precision, intent(inout) :: fE(ixI^S,sdim:3)
347  end subroutine sub_update_faces
348 
349  subroutine sub_face_to_center(ixO^L,s)
351  integer, intent(in) :: ixO^L
352  type(state) :: s
353  end subroutine sub_face_to_center
354 
355  subroutine sub_evaluate_implicit(qtC,psa)
357  type(state), target :: psa(max_blocks)
358  double precision, intent(in) :: qtC
359  end subroutine sub_evaluate_implicit
360 
361  subroutine sub_implicit_update(dtfactor,qdt,qtC,psa,psb)
363  type(state), target :: psa(max_blocks)
364  type(state), target :: psb(max_blocks)
365  double precision, intent(in) :: qdt
366  double precision, intent(in) :: qtC
367  double precision, intent(in) :: dtfactor
368  end subroutine sub_implicit_update
369 
370  subroutine sub_update_temperature(ixI^L,ixO^L,wCT,w,x)
372  integer, intent(in) :: ixI^L, ixO^L
373  double precision, intent(in) :: wCT(ixI^S,nw),x(ixI^S,1:ndim)
374  double precision, intent(inout) :: w(ixI^S,nw)
375  end subroutine sub_update_temperature
376 
377  end interface
378 
379 contains
380 
381  subroutine phys_check()
382  use mod_global_parameters, only: nw, ndir
383 
386  use mod_comm_lib, only: mpistop
387 
388  if (physics_type == "") call mpistop("Error: no physics module loaded")
389 
390  call phys_hllc_check()
391  call phys_roe_check()
392 
393  ! Checks whether the required physics methods have been defined
394  if (.not. associated(phys_check_params)) &
396 
397  if (.not. associated(phys_to_conserved)) &
398  call mpistop("Error: phys_to_conserved not defined")
399 
400  if (.not. associated(phys_to_primitive)) &
401  call mpistop("Error: phys_to_primitive not defined")
402 
403  if (.not. associated(phys_modify_wlr)) &
405 
406  if (.not. associated(phys_get_cmax)) &
407  call mpistop("Error: no phys_get_cmax not defined")
408 
409  if (.not. associated(phys_get_a2max)) &
411 
412  if (.not. associated(phys_get_cs2max)) &
414 
415  if (.not. associated(phys_get_h_speed)) &
417 
418  if (.not. associated(phys_get_cbounds)) &
419  call mpistop("Error: no phys_get_cbounds not defined")
420 
421  if (.not. associated(phys_get_flux)) &
422  call mpistop("Error: no phys_get_flux not defined")
423 
424  if (.not. associated(phys_get_dt)) &
425  call mpistop("Error: no phys_get_dt not defined")
426 
427  if (.not. associated(phys_add_source_geom)) &
429 
430  if (.not. associated(phys_add_source)) &
432 
433  if (.not. associated(phys_check_w)) &
435 
436  if (.not. associated(phys_get_pthermal)) &
438  if (.not. associated(phys_get_auxiliary)) &
440  if (.not. associated(phys_get_auxiliary_prim)) &
442 
443  if (.not. associated(phys_boundary_adjust)) &
445 
446  if (.not. associated(phys_write_info)) &
448 
449  if (.not. associated(phys_handle_small_values)) &
451 
452  if (.not. associated(phys_get_ct_velocity)) &
454 
455  if (.not. associated(phys_update_faces)) &
457 
458  if (.not. associated(phys_face_to_center)) &
460 
461  if (.not. associated(phys_evaluate_implicit)) &
463 
464  if (.not. associated(phys_implicit_update)) &
466 
467  end subroutine phys_check
468 
469  subroutine dummy_init_params
470  end subroutine dummy_init_params
471 
473  end subroutine dummy_check_params
474 
475  subroutine dummy_modify_wlr(ixI^L, ixO^L, qt, wLC, wRC, wLp, wRp, s, idir)
477  integer, intent(in) :: ixI^L, ixO^L, idir
478  double precision, intent(in) :: qt
479  double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
480  double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
481  type(state) :: s
482  end subroutine dummy_modify_wlr
483 
484  subroutine dummy_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
486  integer, intent(in) :: ixI^L, ixO^L, idim
487  double precision, intent(in) :: wprim(ixI^S, nw)
488  double precision, intent(in) :: x(ixI^S,1:ndim)
489  double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
490  end subroutine dummy_get_h_speed
491 
492  subroutine dummy_get_a2max(w, x, ixI^L, ixO^L, a2max)
494  use mod_comm_lib, only: mpistop
495  integer, intent(in) :: ixI^L, ixO^L
496  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
497  double precision, intent(inout) :: a2max(ndim)
498  call mpistop("Error: entered dummy_get_a2max")
499  end subroutine dummy_get_a2max
500 
501  subroutine dummy_get_cs2max(w, x, ixI^L, ixO^L, cs2max)
503  use mod_comm_lib, only: mpistop
504  integer, intent(in) :: ixI^L, ixO^L
505  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
506  double precision, intent(inout) :: cs2max
507  call mpistop("Error: entered dummy_get_cs2max")
508  end subroutine dummy_get_cs2max
509 
510  subroutine dummy_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, w, x)
512  integer, intent(in) :: ixI^L, ixO^L
513  double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
514  double precision, intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
515  end subroutine dummy_add_source_geom
516 
517  subroutine dummy_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
518  qsourcesplit, active)
520  integer, intent(in) :: ixI^L, ixO^L
521  double precision, intent(in) :: qdt,dtfactor
522  double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
523  double precision, intent(inout) :: w(ixI^S, 1:nw)
524  logical, intent(in) :: qsourcesplit
525  logical, intent(inout) :: active
526  ! Don't have to set active, since it starts as .false.
527  end subroutine dummy_add_source
528 
529  subroutine dummy_check_w(primitive, ixI^L, ixO^L, w, w_flag)
531  logical, intent(in) :: primitive
532  integer, intent(in) :: ixI^L, ixO^L
533  double precision, intent(in) :: w(ixI^S,1:nw)
534  logical, intent(inout) :: w_flag(ixI^S,1:nw)
535 
536  w_flag=.false. ! All okay
537  end subroutine dummy_check_w
538 
539  subroutine dummy_get_pthermal(w, x, ixI^L, ixO^L, pth)
541  use mod_comm_lib, only: mpistop
542 
543  integer, intent(in) :: ixI^L, ixO^L
544  double precision, intent(in) :: w(ixI^S, nw)
545  double precision, intent(in) :: x(ixI^S, 1:ndim)
546  double precision, intent(out):: pth(ixI^S)
547 
548  call mpistop("No get_pthermal method specified")
549  end subroutine dummy_get_pthermal
550 
551  subroutine dummy_get_auxiliary(ixI^L, ixO^L, w, x)
553  use mod_comm_lib, only: mpistop
554 
555  integer, intent(in) :: ixI^L, ixO^L
556  double precision, intent(inout) :: w(ixI^S, nw)
557  double precision, intent(in) :: x(ixI^S, 1:ndim)
558 
559  !call mpistop("No get_auxiliary method specified")
560  end subroutine dummy_get_auxiliary
561 
562  subroutine dummy_get_auxiliary_prim(ixI^L, ixO^L, w)
564  use mod_comm_lib, only: mpistop
565 
566  integer, intent(in) :: ixI^L, ixO^L
567  double precision, intent(inout) :: w(ixI^S, nw)
568 
569  call mpistop("No get_auxiliary_prim method specified")
570  end subroutine dummy_get_auxiliary_prim
571 
572  subroutine dummy_boundary_adjust(igrid,psb)
574  integer, intent(in) :: igrid
575  type(state), target :: psb(max_blocks)
576  end subroutine dummy_boundary_adjust
577 
578  subroutine dummy_write_info(fh)
580  integer, intent(in) :: fh !< File handle
581  integer, dimension(MPI_STATUS_SIZE) :: st
582  integer :: er
583 
584  ! Number of physics parameters
585  integer, parameter :: n_par = 0
586 
587  call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
588  end subroutine dummy_write_info
589 
590  subroutine dummy_small_values(primitive, w, x, ixI^L, ixO^L, subname)
592  logical, intent(in) :: primitive
593  integer, intent(in) :: ixI^L,ixO^L
594  double precision, intent(inout) :: w(ixI^S,1:nw)
595  double precision, intent(in) :: x(ixI^S,1:ndim)
596  character(len=*), intent(in) :: subname
597  end subroutine dummy_small_values
598 
599  subroutine dummy_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
601 
602  integer, intent(in) :: ixI^L, ixO^L, idim
603  double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
604  double precision, intent(in) :: cmax(ixI^S)
605  double precision, intent(in), optional :: cmin(ixI^S)
606  type(ct_velocity), intent(inout):: vcts
607  end subroutine dummy_get_ct_velocity
608 
609  subroutine dummy_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
611  integer, intent(in) :: ixI^L, ixO^L
612  double precision, intent(in) :: qt, qdt
613  ! cell-center primitive variables
614  double precision, intent(in) :: wprim(ixI^S,1:nw)
615  type(state) :: sCT, s
616  type(ct_velocity) :: vcts
617  double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
618  double precision, intent(inout) :: fE(ixI^S,sdim:3)
619  end subroutine dummy_update_faces
620 
621  subroutine dummy_face_to_center(ixO^L,s)
623  integer, intent(in) :: ixO^L
624  type(state) :: s
625  end subroutine dummy_face_to_center
626 
627  subroutine dummy_evaluate_implicit(qtC,psa)
629  type(state), target :: psa(max_blocks)
630  double precision, intent(in) :: qtC
631  integer :: iigrid, igrid
632 
633  ! Just copy in nul state
634  !$OMP PARALLEL DO PRIVATE(igrid)
635  do iigrid=1,igridstail; igrid=igrids(iigrid);
636  psa(igrid)%w = 0.0d0*psa(igrid)%w
637  if(stagger_grid) psa(igrid)%ws = 0.0d0*psa(igrid)%ws
638  end do
639  !$OMP END PARALLEL DO
640 
641  end subroutine dummy_evaluate_implicit
642 
643  subroutine dummy_implicit_update(dtfactor,qdt,qtC,psa,psb)
645  type(state), target :: psa(max_blocks)
646  type(state), target :: psb(max_blocks)
647  double precision, intent(in) :: qdt
648  double precision, intent(in) :: qtC
649  double precision, intent(in) :: dtfactor
650  integer :: iigrid, igrid
651 
652  ! Just copy in psb state when using the scheme without implicit part
653  !$OMP PARALLEL DO PRIVATE(igrid)
654  do iigrid=1,igridstail; igrid=igrids(iigrid);
655  psa(igrid)%w = psb(igrid)%w
656  if(stagger_grid) psa(igrid)%ws = psb(igrid)%ws
657  end do
658  !$OMP END PARALLEL DO
659 
660  end subroutine dummy_implicit_update
661 
662 end module mod_physics
clean initial divb
Definition: mod_physics.t:233
Add global source terms on complete domain (potentially implicit)
Definition: mod_physics.t:225
set equilibrium variables other than b0 (e.g. p0 and rho0)
Definition: mod_physics.t:241
Add special advance in each advect step.
Definition: mod_physics.t:247
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:83
subroutine dummy_update_faces(ixIL, ixOL, qt, qdt, wprim, fC, fE, sCT, s, vcts)
Definition: mod_physics.t:610
subroutine dummy_boundary_adjust(igrid, psb)
Definition: mod_physics.t:573
subroutine dummy_get_ct_velocity(vcts, wLp, wRp, ixIL, ixOL, idim, cmax, cmin)
Definition: mod_physics.t:600
procedure(sub_convert), pointer phys_to_primitive
Definition: mod_physics.t:59
procedure(sub_get_tgas), pointer phys_get_tgas
Definition: mod_physics.t:78
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:82
procedure(sub_write_info), pointer phys_write_info
Definition: mod_physics.t:81
procedure(sub_get_flux), pointer phys_get_flux
Definition: mod_physics.t:68
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
Definition: mod_physics.t:87
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:67
procedure(sub_check_params), pointer phys_te_images
Definition: mod_physics.t:92
subroutine dummy_face_to_center(ixOL, s)
Definition: mod_physics.t:622
procedure(sub_get_dt), pointer phys_get_dt
Definition: mod_physics.t:71
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:77
subroutine dummy_get_h_speed(wprim, x, ixIL, ixOL, idim, Hspeed)
Definition: mod_physics.t:485
procedure(sub_update_temperature), pointer phys_update_temperature
Definition: mod_physics.t:94
procedure(sub_get_tcutoff), pointer phys_get_tcutoff
Definition: mod_physics.t:64
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition: mod_physics.t:72
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_cs2max), pointer phys_get_cs2max
Definition: mod_physics.t:63
procedure(sub_get_auxiliary), pointer phys_get_auxiliary
Definition: mod_physics.t:95
procedure(sub_get_auxiliary_prim), pointer phys_get_auxiliary_prim
Definition: mod_physics.t:96
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:65
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:90
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:470
procedure(sub_update_faces), pointer phys_update_faces
Definition: mod_physics.t:84
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:76
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:628
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:86
procedure(sub_get_trad), pointer phys_get_trad
Definition: mod_physics.t:79
subroutine dummy_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, w, x)
Definition: mod_physics.t:511
subroutine dummy_check_params
Definition: mod_physics.t:473
subroutine dummy_get_auxiliary(ixIL, ixOL, w, x)
Definition: mod_physics.t:552
procedure(sub_clean_divb), pointer phys_clean_divb
Definition: mod_physics.t:88
subroutine dummy_implicit_update(dtfactor, qdt, qtC, psa, psb)
Definition: mod_physics.t:644
procedure(sub_boundary_adjust), pointer phys_boundary_adjust
Definition: mod_physics.t:80
double precision phys_gamma
Definition: mod_physics.t:13
procedure(sub_global_source), pointer phys_global_source_after
Definition: mod_physics.t:74
procedure(sub_add_source), pointer phys_add_source
Definition: mod_physics.t:73
procedure(sub_face_to_center), pointer phys_face_to_center
Definition: mod_physics.t:85
subroutine dummy_get_pthermal(w, x, ixIL, ixOL, pth)
Definition: mod_physics.t:540
procedure(sub_modify_wlr), pointer phys_modify_wlr
Definition: mod_physics.t:60
subroutine dummy_write_info(fh)
Definition: mod_physics.t:579
procedure(sub_get_h_speed), pointer phys_get_h_speed
Definition: mod_physics.t:66
logical phys_internal_e
Solve internal energy instead of total energy.
Definition: mod_physics.t:33
subroutine dummy_get_cs2max(w, x, ixIL, ixOL, cs2max)
Definition: mod_physics.t:502
procedure(sub_get_cmax), pointer phys_get_cmax
Definition: mod_physics.t:61
procedure(sub_get_rho), pointer phys_get_rho
Definition: mod_physics.t:70
procedure(sub_get_v), pointer phys_get_v
Definition: mod_physics.t:69
subroutine dummy_add_source(qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
Definition: mod_physics.t:519
subroutine dummy_get_auxiliary_prim(ixIL, ixOL, w)
Definition: mod_physics.t:563
subroutine dummy_get_a2max(w, x, ixIL, ixOL, a2max)
Definition: mod_physics.t:493
subroutine phys_check()
Definition: mod_physics.t:382
subroutine dummy_small_values(primitive, w, x, ixIL, ixOL, subname)
Definition: mod_physics.t:591
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:530
procedure(sub_special_advance), pointer phys_special_advance
Definition: mod_physics.t:75
subroutine dummy_modify_wlr(ixIL, ixOL, qt, wLC, wRC, wLp, wRp, s, idir)
Definition: mod_physics.t:476
Module with all the methods that users can customize in AMRVAC.