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 !> extra ghost cells in the transverse dimensions for fluxes of CT
35
36 !> Solve energy equation or not
37 logical :: phys_energy=.false.
38
39 !> Solve total energy equation or not
40 logical :: phys_total_energy=.false.
41
42 !> Solve internal energy instead of total energy
43 logical :: phys_internal_e=.false.
44
45 !> Solve partially ionized one-fluid plasma
46 logical :: phys_partial_ionization=.false.
47
48 !> if equilibrium pressure is splitted
49 logical :: phys_equi_pe=.false.
50
51 !> String describing the physics type of the simulation
52 character(len=name_len) :: physics_type = ""
53
54 procedure(sub_check_params), pointer :: phys_check_params => null()
55 procedure(sub_set_mg_bounds), pointer :: phys_set_mg_bounds => null()
56 procedure(sub_convert), pointer :: phys_to_conserved => null()
57 procedure(sub_convert), pointer :: phys_to_primitive => null()
58 procedure(sub_modify_wlr), pointer :: phys_modify_wlr => null()
59 procedure(sub_get_cmax), pointer :: phys_get_cmax => null()
60 procedure(sub_get_a2max), pointer :: phys_get_a2max => null()
61 procedure(sub_get_tcutoff), pointer :: phys_get_tcutoff => null()
62 procedure(sub_trac_after_setdt), pointer:: phys_trac_after_setdt => null()
63 procedure(sub_get_h_speed), pointer :: phys_get_h_speed => null()
64 procedure(sub_get_cbounds), pointer :: phys_get_cbounds => null()
65 procedure(sub_get_flux), pointer :: phys_get_flux => null()
66 procedure(sub_get_v), pointer :: phys_get_v => null()
67 procedure(sub_get_rho), pointer :: phys_get_rho => null()
68 procedure(sub_get_dt), pointer :: phys_get_dt => null()
69 procedure(sub_add_source_geom), pointer :: phys_add_source_geom => null()
70 procedure(sub_add_source), pointer :: phys_add_source => null()
71 procedure(sub_global_source), pointer :: phys_global_source_after => null()
72 procedure(sub_special_advance), pointer :: phys_special_advance => null()
73 procedure(sub_check_w), pointer :: phys_check_w => null()
74 procedure(sub_get_pthermal), pointer :: phys_get_pthermal => null()
75 procedure(sub_get_tgas), pointer :: phys_get_tgas => null()
76 procedure(sub_get_trad), pointer :: phys_get_trad => null()
77 procedure(sub_boundary_adjust), pointer :: phys_boundary_adjust => null()
78 procedure(sub_write_info), pointer :: phys_write_info => null()
79 procedure(sub_small_values), pointer :: phys_handle_small_values => null()
80 procedure(sub_get_ct_velocity), pointer :: phys_get_ct_velocity => null()
81 procedure(sub_update_faces), pointer :: phys_update_faces => null()
82 procedure(sub_face_to_center), pointer :: phys_face_to_center => null()
83 procedure(sub_implicit_update), pointer :: phys_implicit_update => null()
84 procedure(sub_evaluate_implicit),pointer:: phys_evaluate_implicit => null()
85 procedure(sub_clean_divb), pointer :: phys_clean_divb => null()
86 ! set the equilibrium variables
87 procedure(sub_set_equi_vars), pointer :: phys_set_equi_vars => null()
88 ! subroutine with no parameters which creates EUV images
89 procedure(sub_check_params), pointer :: phys_te_images => null()
90 ! to update temperature variable with partial ionization
91 procedure(sub_update_temperature), pointer :: phys_update_temperature => null()
92 procedure(sub_get_auxiliary), pointer :: phys_get_auxiliary => null()
93 procedure(sub_get_auxiliary_prim), pointer :: phys_get_auxiliary_prim => null()
94
95 abstract interface
96
98 end subroutine sub_check_params
99
103 end subroutine sub_set_mg_bounds
104
105 subroutine sub_boundary_adjust(igrid,psb)
107 integer, intent(in) :: igrid
108 type(state), target :: psb(max_blocks)
109 end subroutine sub_boundary_adjust
110
111 subroutine sub_convert(ixI^L, ixO^L, w, x)
113 integer, intent(in) :: ixI^L, ixO^L
114 double precision, intent(inout) :: w(ixI^S, nw)
115 double precision, intent(in) :: x(ixI^S, 1:^ND)
116 end subroutine sub_convert
117
118 subroutine sub_modify_wlr(ixI^L, ixO^L, qt, wLC, wRC, wLp, wRp, s, idir)
120 integer, intent(in) :: ixI^L, ixO^L, idir
121 double precision, intent(in) :: qt
122 double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
123 double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
124 type(state) :: s
125 end subroutine sub_modify_wlr
126
127 subroutine sub_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
129 integer, intent(in) :: ixI^L, ixO^L, idim
130 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
131 double precision, intent(inout) :: cmax(ixI^S)
132 end subroutine sub_get_cmax
133
134 subroutine sub_get_a2max(w, x, ixI^L, ixO^L, a2max)
136 integer, intent(in) :: ixI^L, ixO^L
137 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
138 double precision, intent(inout) :: a2max(ndim)
139 end subroutine sub_get_a2max
140
141 subroutine sub_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
143 integer, intent(in) :: ixI^L, ixO^L
144 double precision, intent(inout) :: w(ixI^S, nw)
145 double precision, intent(in) :: x(ixI^S, 1:^ND)
146 double precision, intent(out) :: tco_local, Tmax_local
147 end subroutine sub_get_tcutoff
148
149 subroutine sub_trac_after_setdt(trac_alfa,tco,T_peak,T_bott)
150 double precision, intent(in) :: trac_alfa,tco,T_peak,T_bott
151 end subroutine sub_trac_after_setdt
152
153 subroutine sub_get_v(w,x,ixI^L,ixO^L,v)
155 integer, intent(in) :: ixI^L, ixO^L
156 double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
157 double precision, intent(out) :: v(ixI^S,1:ndir)
158 end subroutine sub_get_v
159
160 subroutine sub_get_rho(w,x,ixI^L,ixO^L,rho)
162 integer, intent(in) :: ixI^L, ixO^L
163 double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
164 double precision, intent(out) :: rho(ixI^S)
165 end subroutine sub_get_rho
166
167 subroutine sub_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
169 integer, intent(in) :: ixI^L, ixO^L, idim
170 double precision, intent(in) :: wprim(ixI^S, nw)
171 double precision, intent(in) :: x(ixI^S,1:ndim)
172 double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
173 end subroutine sub_get_h_speed
174
175 subroutine sub_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim, Hspeed, cmax, cmin)
177 use mod_variables
178 integer, intent(in) :: ixI^L, ixO^L, idim
179 double precision, intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
180 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
181 double precision, intent(in) :: x(ixI^S, 1:^ND)
182 double precision, intent(inout) :: cmax(ixI^S,1:number_species)
183 double precision, intent(inout), optional :: cmin(ixI^S,1:number_species)
184 double precision, intent(in) :: Hspeed(ixI^S,1:number_species)
185 end subroutine sub_get_cbounds
186
187 subroutine sub_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
189 integer, intent(in) :: ixI^L, ixO^L, idim
190 double precision, intent(in) :: wC(ixI^S, 1:nw)
191 double precision, intent(in) :: w(ixI^S, 1:nw)
192 double precision, intent(in) :: x(ixI^S, 1:^ND)
193 double precision, intent(out) :: f(ixI^S, nwflux)
194 end subroutine sub_get_flux
195
196 subroutine sub_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
198 integer, intent(in) :: ixI^L, ixO^L
199 double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
200 double precision, intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S, 1:nw), w(ixI^S, 1:nw)
201 end subroutine sub_add_source_geom
202
203 subroutine sub_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
204 qsourcesplit, active)
206 integer, intent(in) :: ixI^L, ixO^L
207 double precision, intent(in) :: qdt, dtfactor
208 double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
209 double precision, intent(inout) :: w(ixI^S, 1:nw)
210 logical, intent(in) :: qsourcesplit
211 logical, intent(inout) :: active !< Needs to be set to true when active
212 end subroutine sub_add_source
213
214 !> Add global source terms on complete domain (potentially implicit)
215 subroutine sub_global_source(qdt, qt, active)
217 double precision, intent(in) :: qdt !< Current time step
218 double precision, intent(in) :: qt !< Current time
219 logical, intent(inout) :: active !< Output if the source is active
220 end subroutine sub_global_source
221
222 !> clean initial divb
223 subroutine sub_clean_divb(qdt, qt, active)
225 double precision, intent(in) :: qdt !< Current time step
226 double precision, intent(in) :: qt !< Current time
227 logical, intent(inout) :: active !< Output if the source is active
228 end subroutine sub_clean_divb
229
230 !> set equilibrium variables other than b0 (e.g. p0 and rho0)
231 subroutine sub_set_equi_vars(igrid)
232 integer, intent(in) :: igrid
233 end subroutine sub_set_equi_vars
234
235
236 !> Add special advance in each advect step
237 subroutine sub_special_advance(qt, psa)
239 double precision, intent(in) :: qt !< Current time
240 type(state), target :: psa(max_blocks) !< Compute based on this state
241 end subroutine sub_special_advance
242
243 subroutine sub_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
245 integer, intent(in) :: ixI^L, ixO^L
246 double precision, intent(in) :: dx^D, x(ixI^S, 1:^ND)
247 double precision, intent(in) :: w(ixI^S, 1:nw)
248 double precision, intent(inout) :: dtnew
249 end subroutine sub_get_dt
250
251 subroutine sub_check_w(primitive, ixI^L, ixO^L, w, w_flag)
253 logical, intent(in) :: primitive
254 integer, intent(in) :: ixI^L, ixO^L
255 double precision, intent(in) :: w(ixI^S,1:nw)
256 logical, intent(inout) :: w_flag(ixI^S,1:nw)
257 end subroutine sub_check_w
258
259 subroutine sub_get_pthermal(w,x,ixI^L,ixO^L,pth)
261 integer, intent(in) :: ixI^L, ixO^L
262 double precision, intent(in) :: w(ixI^S,nw)
263 double precision, intent(in) :: x(ixI^S,1:ndim)
264 double precision, intent(out):: pth(ixI^S)
265 end subroutine sub_get_pthermal
266
267 subroutine sub_get_auxiliary(ixI^L,ixO^L,w,x)
269 integer, intent(in) :: ixI^L, ixO^L
270 double precision, intent(inout) :: w(ixI^S,nw)
271 double precision, intent(in) :: x(ixI^S,1:ndim)
272 end subroutine sub_get_auxiliary
273
274 subroutine sub_get_auxiliary_prim(ixI^L,ixO^L,w)
276 integer, intent(in) :: ixI^L, ixO^L
277 double precision, intent(inout) :: w(ixI^S,nw)
278 end subroutine sub_get_auxiliary_prim
279
280 subroutine sub_get_tgas(w,x,ixI^L,ixO^L,tgas)
282 integer, intent(in) :: ixI^L, ixO^L
283 double precision, intent(in) :: w(ixI^S,nw)
284 double precision, intent(in) :: x(ixI^S,1:ndim)
285 double precision, intent(out):: tgas(ixI^S)
286 end subroutine sub_get_tgas
287
288 subroutine sub_get_trad(w,x,ixI^L,ixO^L,trad)
290 integer, intent(in) :: ixI^L, ixO^L
291 double precision, intent(in) :: w(ixI^S,nw)
292 double precision, intent(in) :: x(ixI^S,1:ndim)
293 double precision, intent(out):: trad(ixI^S)
294 end subroutine sub_get_trad
295
296 subroutine sub_write_info(file_handle)
297 integer, intent(in) :: file_handle
298 end subroutine sub_write_info
299
300 subroutine sub_small_values(primitive, w, x, ixI^L, ixO^L, subname)
302 logical, intent(in) :: primitive
303 integer, intent(in) :: ixI^L,ixO^L
304 double precision, intent(inout) :: w(ixI^S,1:nw)
305 double precision, intent(in) :: x(ixI^S,1:ndim)
306 character(len=*), intent(in) :: subname
307 end subroutine sub_small_values
308
309 subroutine sub_get_var(ixI^L, ixO^L, w, out)
311 integer, intent(in) :: ixI^L, ixO^L
312 double precision, intent(in) :: w(ixI^S, nw)
313 double precision, intent(out) :: out(ixO^S)
314 end subroutine sub_get_var
315
316 subroutine sub_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
318
319 integer, intent(in) :: ixI^L, ixO^L, idim
320 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
321 double precision, intent(in) :: cmax(ixI^S)
322 double precision, intent(in), optional :: cmin(ixI^S)
323 type(ct_velocity), intent(inout):: vcts
324 end subroutine sub_get_ct_velocity
325
326 subroutine sub_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
328 integer, intent(in) :: ixI^L, ixO^L
329 double precision, intent(in) :: qt, qdt
330 ! cell-center primitive variables
331 double precision, intent(in) :: wprim(ixI^S,1:nw)
332 ! velocity structure
333 type(state) :: sCT, s
334 type(ct_velocity) :: vcts
335 double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
336 double precision, intent(inout) :: fE(ixI^S,sdim:3)
337 end subroutine sub_update_faces
338
339 subroutine sub_face_to_center(ixO^L,s)
341 integer, intent(in) :: ixO^L
342 type(state) :: s
343 end subroutine sub_face_to_center
344
345 subroutine sub_evaluate_implicit(qtC,psa)
347 type(state), target :: psa(max_blocks)
348 double precision, intent(in) :: qtC
349 end subroutine sub_evaluate_implicit
350
351 subroutine sub_implicit_update(dtfactor,qdt,qtC,psa,psb)
353 type(state), target :: psa(max_blocks)
354 type(state), target :: psb(max_blocks)
355 double precision, intent(in) :: qdt
356 double precision, intent(in) :: qtC
357 double precision, intent(in) :: dtfactor
358 end subroutine sub_implicit_update
359
360 subroutine sub_update_temperature(ixI^L,ixO^L,wCT,w,x)
362 integer, intent(in) :: ixI^L, ixO^L
363 double precision, intent(in) :: wCT(ixI^S,nw),x(ixI^S,1:ndim)
364 double precision, intent(inout) :: w(ixI^S,nw)
365 end subroutine sub_update_temperature
366
367 end interface
368
369contains
370
371 subroutine phys_check()
372 use mod_global_parameters, only: nw, ndir
373
376 use mod_comm_lib, only: mpistop
377
378 if (physics_type == "") call mpistop("Error: no physics module loaded")
379
380 call phys_hllc_check()
381 call phys_roe_check()
382
383 ! Checks whether the required physics methods have been defined
384 if (.not. associated(phys_check_params)) &
386
387 if (.not. associated(phys_to_conserved)) &
388 call mpistop("Error: phys_to_conserved not defined")
389
390 if (.not. associated(phys_to_primitive)) &
391 call mpistop("Error: phys_to_primitive not defined")
392
393 if (.not. associated(phys_modify_wlr)) &
395
396 if (.not. associated(phys_get_cmax)) &
397 call mpistop("Error: no phys_get_cmax not defined")
398
399 if (.not. associated(phys_get_a2max)) &
401
402 if (.not. associated(phys_get_h_speed)) &
404
405 if (.not. associated(phys_get_cbounds)) &
406 call mpistop("Error: no phys_get_cbounds not defined")
407
408 if (.not. associated(phys_get_flux)) &
409 call mpistop("Error: no phys_get_flux not defined")
410
411 if (.not. associated(phys_get_dt)) &
412 call mpistop("Error: no phys_get_dt not defined")
413
414 if (.not. associated(phys_add_source_geom)) &
416
417 if (.not. associated(phys_add_source)) &
419
420 if (.not. associated(phys_check_w)) &
422
423 if (.not. associated(phys_get_pthermal)) &
425 if (.not. associated(phys_get_auxiliary)) &
427 if (.not. associated(phys_get_auxiliary_prim)) &
429
430 if (.not. associated(phys_boundary_adjust)) &
432
433 if (.not. associated(phys_write_info)) &
435
436 if (.not. associated(phys_handle_small_values)) &
438
439 if (.not. associated(phys_get_ct_velocity)) &
441
442 if (.not. associated(phys_update_faces)) &
444
445 if (.not. associated(phys_face_to_center)) &
447
448 if (.not. associated(phys_evaluate_implicit)) &
450
451 if (.not. associated(phys_implicit_update)) &
453
454 end subroutine phys_check
455
457 end subroutine dummy_init_params
458
460 end subroutine dummy_check_params
461
462 subroutine dummy_modify_wlr(ixI^L, ixO^L, qt, wLC, wRC, wLp, wRp, s, idir)
464 integer, intent(in) :: ixI^L, ixO^L, idir
465 double precision, intent(in) :: qt
466 double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
467 double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
468 type(state) :: s
469 end subroutine dummy_modify_wlr
470
471 subroutine dummy_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
473 integer, intent(in) :: ixI^L, ixO^L, idim
474 double precision, intent(in) :: wprim(ixI^S, nw)
475 double precision, intent(in) :: x(ixI^S,1:ndim)
476 double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
477 end subroutine dummy_get_h_speed
478
479 subroutine dummy_get_a2max(w, x, ixI^L, ixO^L, a2max)
481 use mod_comm_lib, only: mpistop
482 integer, intent(in) :: ixI^L, ixO^L
483 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
484 double precision, intent(inout) :: a2max(ndim)
485 call mpistop("Error: entered dummy_get_a2max")
486 end subroutine dummy_get_a2max
487
488 subroutine dummy_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
490 integer, intent(in) :: ixI^L, ixO^L
491 double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
492 double precision, intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S,1:nw),w(ixI^S, 1:nw)
493 end subroutine dummy_add_source_geom
494
495 subroutine dummy_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
496 qsourcesplit, active)
498 integer, intent(in) :: ixI^L, ixO^L
499 double precision, intent(in) :: qdt,dtfactor
500 double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
501 double precision, intent(inout) :: w(ixI^S, 1:nw)
502 logical, intent(in) :: qsourcesplit
503 logical, intent(inout) :: active
504 ! Don't have to set active, since it starts as .false.
505 end subroutine dummy_add_source
506
507 subroutine dummy_check_w(primitive, ixI^L, ixO^L, w, w_flag)
509 logical, intent(in) :: primitive
510 integer, intent(in) :: ixI^L, ixO^L
511 double precision, intent(in) :: w(ixI^S,1:nw)
512 logical, intent(inout) :: w_flag(ixI^S,1:nw)
513
514 w_flag=.false. ! All okay
515 end subroutine dummy_check_w
516
517 subroutine dummy_get_pthermal(w, x, ixI^L, ixO^L, pth)
519 use mod_comm_lib, only: mpistop
520
521 integer, intent(in) :: ixI^L, ixO^L
522 double precision, intent(in) :: w(ixI^S, nw)
523 double precision, intent(in) :: x(ixI^S, 1:ndim)
524 double precision, intent(out):: pth(ixI^S)
525
526 call mpistop("No get_pthermal method specified")
527 end subroutine dummy_get_pthermal
528
529 subroutine dummy_get_auxiliary(ixI^L, ixO^L, w, x)
531 use mod_comm_lib, only: mpistop
532
533 integer, intent(in) :: ixI^L, ixO^L
534 double precision, intent(inout) :: w(ixI^S, nw)
535 double precision, intent(in) :: x(ixI^S, 1:ndim)
536
537 !call mpistop("No get_auxiliary method specified")
538 end subroutine dummy_get_auxiliary
539
540 subroutine dummy_get_auxiliary_prim(ixI^L, ixO^L, w)
542 use mod_comm_lib, only: mpistop
543
544 integer, intent(in) :: ixI^L, ixO^L
545 double precision, intent(inout) :: w(ixI^S, nw)
546
547 call mpistop("No get_auxiliary_prim method specified")
548 end subroutine dummy_get_auxiliary_prim
549
550 subroutine dummy_boundary_adjust(igrid,psb)
552 integer, intent(in) :: igrid
553 type(state), target :: psb(max_blocks)
554 end subroutine dummy_boundary_adjust
555
556 subroutine dummy_write_info(fh)
558 integer, intent(in) :: fh !< File handle
559 integer, dimension(MPI_STATUS_SIZE) :: st
560 integer :: er
561
562 ! Number of physics parameters
563 integer, parameter :: n_par = 0
564
565 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
566 end subroutine dummy_write_info
567
568 subroutine dummy_small_values(primitive, w, x, ixI^L, ixO^L, subname)
570 logical, intent(in) :: primitive
571 integer, intent(in) :: ixI^L,ixO^L
572 double precision, intent(inout) :: w(ixI^S,1:nw)
573 double precision, intent(in) :: x(ixI^S,1:ndim)
574 character(len=*), intent(in) :: subname
575 end subroutine dummy_small_values
576
577 subroutine dummy_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
579
580 integer, intent(in) :: ixI^L, ixO^L, idim
581 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
582 double precision, intent(in) :: cmax(ixI^S)
583 double precision, intent(in), optional :: cmin(ixI^S)
584 type(ct_velocity), intent(inout):: vcts
585 end subroutine dummy_get_ct_velocity
586
587 subroutine dummy_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
589 integer, intent(in) :: ixI^L, ixO^L
590 double precision, intent(in) :: qt, qdt
591 ! cell-center primitive variables
592 double precision, intent(in) :: wprim(ixI^S,1:nw)
593 type(state) :: sCT, s
594 type(ct_velocity) :: vcts
595 double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
596 double precision, intent(inout) :: fE(ixI^S,sdim:3)
597 end subroutine dummy_update_faces
598
599 subroutine dummy_face_to_center(ixO^L,s)
601 integer, intent(in) :: ixO^L
602 type(state) :: s
603 end subroutine dummy_face_to_center
604
605 subroutine dummy_evaluate_implicit(qtC,psa)
607 type(state), target :: psa(max_blocks)
608 double precision, intent(in) :: qtC
609 integer :: iigrid, igrid
610
611 ! Just copy in nul state
612 !$OMP PARALLEL DO PRIVATE(igrid)
613 do iigrid=1,igridstail; igrid=igrids(iigrid);
614 psa(igrid)%w = 0.0d0*psa(igrid)%w
615 if(stagger_grid) psa(igrid)%ws = 0.0d0*psa(igrid)%ws
616 end do
617 !$OMP END PARALLEL DO
618
619 end subroutine dummy_evaluate_implicit
620
621 subroutine dummy_implicit_update(dtfactor,qdt,qtC,psa,psb)
623 type(state), target :: psa(max_blocks)
624 type(state), target :: psb(max_blocks)
625 double precision, intent(in) :: qdt
626 double precision, intent(in) :: qtC
627 double precision, intent(in) :: dtfactor
628 integer :: iigrid, igrid
629
630 ! Just copy in psb state when using the scheme without implicit part
631 !$OMP PARALLEL DO PRIVATE(igrid)
632 do iigrid=1,igridstail; igrid=igrids(iigrid);
633 psa(igrid)%w = psb(igrid)%w
634 if(stagger_grid) psa(igrid)%ws = psb(igrid)%ws
635 end do
636 !$OMP END PARALLEL DO
637
638 end subroutine dummy_implicit_update
639
640end 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:60
procedure(sub_get_ct_velocity), pointer phys_get_ct_velocity
Definition mod_physics.t:80
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:57
procedure(sub_get_tgas), pointer phys_get_tgas
Definition mod_physics.t:75
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:79
procedure(sub_write_info), pointer phys_write_info
Definition mod_physics.t:78
procedure(sub_get_flux), pointer phys_get_flux
Definition mod_physics.t:65
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
Definition mod_physics.t:84
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:64
procedure(sub_check_params), pointer phys_te_images
Definition mod_physics.t:89
procedure(sub_get_dt), pointer phys_get_dt
Definition mod_physics.t:68
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:46
logical phys_total_energy
Solve total energy equation or not.
Definition mod_physics.t:40
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:74
procedure(sub_update_temperature), pointer phys_update_temperature
Definition mod_physics.t:91
procedure(sub_get_tcutoff), pointer phys_get_tcutoff
Definition mod_physics.t:61
subroutine dummy_get_pthermal(w, x, ixil, ixol, pth)
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition mod_physics.t:69
procedure(sub_check_params), pointer phys_check_params
Definition mod_physics.t:54
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:92
integer transverse_ghost_cells
extra ghost cells in the transverse dimensions for fluxes of CT
Definition mod_physics.t:34
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:93
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:62
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:87
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:81
procedure(sub_convert), pointer phys_to_conserved
Definition mod_physics.t:56
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition mod_physics.t:52
procedure(sub_check_w), pointer phys_check_w
Definition mod_physics.t:73
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:55
logical phys_equi_pe
if equilibrium pressure is splitted
Definition mod_physics.t:49
procedure(sub_implicit_update), pointer phys_implicit_update
Definition mod_physics.t:83
procedure(sub_get_trad), pointer phys_get_trad
Definition mod_physics.t:76
subroutine dummy_evaluate_implicit(qtc, psa)
subroutine dummy_check_params
procedure(sub_clean_divb), pointer phys_clean_divb
Definition mod_physics.t:85
procedure(sub_boundary_adjust), pointer phys_boundary_adjust
Definition mod_physics.t:77
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:71
procedure(sub_add_source), pointer phys_add_source
Definition mod_physics.t:70
procedure(sub_face_to_center), pointer phys_face_to_center
Definition mod_physics.t:82
procedure(sub_modify_wlr), pointer phys_modify_wlr
Definition mod_physics.t:58
subroutine dummy_write_info(fh)
procedure(sub_get_h_speed), pointer phys_get_h_speed
Definition mod_physics.t:63
subroutine dummy_get_auxiliary_prim(ixil, ixol, w)
logical phys_internal_e
Solve internal energy instead of total energy.
Definition mod_physics.t:43
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:59
procedure(sub_get_rho), pointer phys_get_rho
Definition mod_physics.t:67
procedure(sub_get_v), pointer phys_get_v
Definition mod_physics.t:66
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:37
procedure(sub_special_advance), pointer phys_special_advance
Definition mod_physics.t:72
Module with all the methods that users can customize in AMRVAC.