MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
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 !> String describing the physics type of the simulation
49 character(len=name_len) :: physics_type = ""
50
51 procedure(sub_check_params), pointer :: phys_check_params => null()
52 procedure(sub_set_mg_bounds), pointer :: phys_set_mg_bounds => null()
53 procedure(sub_convert), pointer :: phys_to_conserved => null()
54 procedure(sub_convert), pointer :: phys_to_primitive => null()
55 procedure(sub_modify_wlr), pointer :: phys_modify_wlr => null()
56 procedure(sub_get_cmax), pointer :: phys_get_cmax => null()
57 procedure(sub_get_cs2), pointer :: phys_get_cs2 => 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_cs2(w, x, ixI^L, ixO^L, cs2)
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) :: cs2(ixI^S)
137 end subroutine sub_get_cs2
138
139 subroutine sub_get_a2max(w, x, ixI^L, ixO^L, a2max)
141 integer, intent(in) :: ixI^L, ixO^L
142 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
143 double precision, intent(inout) :: a2max(ndim)
144 end subroutine sub_get_a2max
145
146 subroutine sub_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
148 integer, intent(in) :: ixI^L, ixO^L
149 double precision, intent(inout) :: w(ixI^S, nw)
150 double precision, intent(in) :: x(ixI^S, 1:^ND)
151 double precision, intent(out) :: tco_local, Tmax_local
152 end subroutine sub_get_tcutoff
153
154 subroutine sub_trac_after_setdt(trac_alfa,tco,T_peak,T_bott)
155 double precision, intent(in) :: trac_alfa,tco,T_peak,T_bott
156 end subroutine sub_trac_after_setdt
157
158 subroutine sub_get_v(w,x,ixI^L,ixO^L,v)
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) :: v(ixI^S,1:ndir)
163 end subroutine sub_get_v
164
165 subroutine sub_get_rho(w,x,ixI^L,ixO^L,rho)
167 integer, intent(in) :: ixI^L, ixO^L
168 double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
169 double precision, intent(out) :: rho(ixI^S)
170 end subroutine sub_get_rho
171
172 subroutine sub_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
174 integer, intent(in) :: ixI^L, ixO^L, idim
175 double precision, intent(in) :: wprim(ixI^S, nw)
176 double precision, intent(in) :: x(ixI^S,1:ndim)
177 double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
178 end subroutine sub_get_h_speed
179
180 subroutine sub_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim, Hspeed, cmax, cmin)
182 use mod_variables
183 integer, intent(in) :: ixI^L, ixO^L, idim
184 double precision, intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
185 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
186 double precision, intent(in) :: x(ixI^S, 1:^ND)
187 double precision, intent(inout) :: cmax(ixI^S,1:number_species)
188 double precision, intent(inout), optional :: cmin(ixI^S,1:number_species)
189 double precision, intent(in) :: Hspeed(ixI^S,1:number_species)
190 end subroutine sub_get_cbounds
191
192 subroutine sub_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
194 integer, intent(in) :: ixI^L, ixO^L, idim
195 double precision, intent(in) :: wC(ixI^S, 1:nw)
196 double precision, intent(in) :: w(ixI^S, 1:nw)
197 double precision, intent(in) :: x(ixI^S, 1:^ND)
198 double precision, intent(out) :: f(ixI^S, nwflux)
199 end subroutine sub_get_flux
200
201 subroutine sub_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
203 integer, intent(in) :: ixI^L, ixO^L
204 double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
205 double precision, intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S, 1:nw), w(ixI^S, 1:nw)
206 end subroutine sub_add_source_geom
207
208 subroutine sub_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
209 qsourcesplit, active)
211 integer, intent(in) :: ixI^L, ixO^L
212 double precision, intent(in) :: qdt, dtfactor
213 double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
214 double precision, intent(inout) :: w(ixI^S, 1:nw)
215 logical, intent(in) :: qsourcesplit
216 logical, intent(inout) :: active !< Needs to be set to true when active
217 end subroutine sub_add_source
218
219 !> Add global source terms on complete domain (potentially implicit)
220 subroutine sub_global_source(qdt, qt, active)
222 double precision, intent(in) :: qdt !< Current time step
223 double precision, intent(in) :: qt !< Current time
224 logical, intent(inout) :: active !< Output if the source is active
225 end subroutine sub_global_source
226
227 !> clean initial divb
228 subroutine sub_clean_divb(qdt, qt, active)
230 double precision, intent(in) :: qdt !< Current time step
231 double precision, intent(in) :: qt !< Current time
232 logical, intent(inout) :: active !< Output if the source is active
233 end subroutine sub_clean_divb
234
235 !> set equilibrium variables other than b0 (e.g. p0 and rho0)
236 subroutine sub_set_equi_vars(igrid)
237 integer, intent(in) :: igrid
238 end subroutine sub_set_equi_vars
239
240
241 !> Add special advance in each advect step
242 subroutine sub_special_advance(qt, psa)
244 double precision, intent(in) :: qt !< Current time
245 type(state), target :: psa(max_blocks) !< Compute based on this state
246 end subroutine sub_special_advance
247
248 subroutine sub_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
250 integer, intent(in) :: ixI^L, ixO^L
251 double precision, intent(in) :: dx^D, x(ixI^S, 1:^ND)
252 double precision, intent(in) :: w(ixI^S, 1:nw)
253 double precision, intent(inout) :: dtnew
254 end subroutine sub_get_dt
255
256 subroutine sub_check_w(primitive, ixI^L, ixO^L, w, w_flag)
258 logical, intent(in) :: primitive
259 integer, intent(in) :: ixI^L, ixO^L
260 double precision, intent(in) :: w(ixI^S,1:nw)
261 logical, intent(inout) :: w_flag(ixI^S,1:nw)
262 end subroutine sub_check_w
263
264 subroutine sub_get_pthermal(w,x,ixI^L,ixO^L,pth)
266 integer, intent(in) :: ixI^L, ixO^L
267 double precision, intent(in) :: w(ixI^S,nw)
268 double precision, intent(in) :: x(ixI^S,1:ndim)
269 double precision, intent(out):: pth(ixI^S)
270 end subroutine sub_get_pthermal
271
272 subroutine sub_get_auxiliary(ixI^L,ixO^L,w,x)
274 integer, intent(in) :: ixI^L, ixO^L
275 double precision, intent(inout) :: w(ixI^S,nw)
276 double precision, intent(in) :: x(ixI^S,1:ndim)
277 end subroutine sub_get_auxiliary
278
279 subroutine sub_get_auxiliary_prim(ixI^L,ixO^L,w)
281 integer, intent(in) :: ixI^L, ixO^L
282 double precision, intent(inout) :: w(ixI^S,nw)
283 end subroutine sub_get_auxiliary_prim
284
285 subroutine sub_get_tgas(w,x,ixI^L,ixO^L,tgas)
287 integer, intent(in) :: ixI^L, ixO^L
288 double precision, intent(in) :: w(ixI^S,nw)
289 double precision, intent(in) :: x(ixI^S,1:ndim)
290 double precision, intent(out):: tgas(ixI^S)
291 end subroutine sub_get_tgas
292
293 subroutine sub_get_trad(w,x,ixI^L,ixO^L,trad)
295 integer, intent(in) :: ixI^L, ixO^L
296 double precision, intent(in) :: w(ixI^S,nw)
297 double precision, intent(in) :: x(ixI^S,1:ndim)
298 double precision, intent(out):: trad(ixI^S)
299 end subroutine sub_get_trad
300
301 subroutine sub_write_info(file_handle)
302 integer, intent(in) :: file_handle
303 end subroutine sub_write_info
304
305 subroutine sub_small_values(primitive, w, x, ixI^L, ixO^L, subname)
307 logical, intent(in) :: primitive
308 integer, intent(in) :: ixI^L,ixO^L
309 double precision, intent(inout) :: w(ixI^S,1:nw)
310 double precision, intent(in) :: x(ixI^S,1:ndim)
311 character(len=*), intent(in) :: subname
312 end subroutine sub_small_values
313
314 subroutine sub_get_var(ixI^L, ixO^L, w, out)
316 integer, intent(in) :: ixI^L, ixO^L
317 double precision, intent(in) :: w(ixI^S, nw)
318 double precision, intent(out) :: out(ixO^S)
319 end subroutine sub_get_var
320
321 subroutine sub_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
323
324 integer, intent(in) :: ixI^L, ixO^L, idim
325 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
326 double precision, intent(in) :: cmax(ixI^S)
327 double precision, intent(in), optional :: cmin(ixI^S)
328 type(ct_velocity), intent(inout):: vcts
329 end subroutine sub_get_ct_velocity
330
331 subroutine sub_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
333 integer, intent(in) :: ixI^L, ixO^L
334 double precision, intent(in) :: qt, qdt
335 ! cell-center primitive variables
336 double precision, intent(in) :: wprim(ixI^S,1:nw)
337 ! velocity structure
338 type(state) :: sCT, s
339 type(ct_velocity) :: vcts
340 double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
341 double precision, intent(inout) :: fE(ixI^S,sdim:3)
342 end subroutine sub_update_faces
343
344 subroutine sub_face_to_center(ixO^L,s)
346 integer, intent(in) :: ixO^L
347 type(state) :: s
348 end subroutine sub_face_to_center
349
350 subroutine sub_evaluate_implicit(qtC,psa)
352 type(state), target :: psa(max_blocks)
353 double precision, intent(in) :: qtC
354 end subroutine sub_evaluate_implicit
355
356 subroutine sub_implicit_update(dtfactor,qdt,qtC,psa,psb)
358 type(state), target :: psa(max_blocks)
359 type(state), target :: psb(max_blocks)
360 double precision, intent(in) :: qdt
361 double precision, intent(in) :: qtC
362 double precision, intent(in) :: dtfactor
363 end subroutine sub_implicit_update
364
365 subroutine sub_update_temperature(ixI^L,ixO^L,wCT,w,x)
367 integer, intent(in) :: ixI^L, ixO^L
368 double precision, intent(in) :: wCT(ixI^S,nw),x(ixI^S,1:ndim)
369 double precision, intent(inout) :: w(ixI^S,nw)
370 end subroutine sub_update_temperature
371
372 end interface
373
374contains
375
376 subroutine phys_check()
377 use mod_global_parameters, only: nw, ndir
378
381 use mod_comm_lib, only: mpistop
382
383 if (physics_type == "") call mpistop("Error: no physics module loaded")
384
385 call phys_hllc_check()
386 call phys_roe_check()
387
388 ! Checks whether the required physics methods have been defined
389 if (.not. associated(phys_check_params)) &
391
392 if (.not. associated(phys_to_conserved)) &
393 call mpistop("Error: phys_to_conserved not defined")
394
395 if (.not. associated(phys_to_primitive)) &
396 call mpistop("Error: phys_to_primitive not defined")
397
398 if (.not. associated(phys_modify_wlr)) &
400
401 if (.not. associated(phys_get_cmax)) &
402 call mpistop("Error: no phys_get_cmax not defined")
403
404 if (.not. associated(phys_get_a2max)) &
406
407 if (.not. associated(phys_get_h_speed)) &
409
410 if (.not. associated(phys_get_cbounds)) &
411 call mpistop("Error: no phys_get_cbounds not defined")
412
413 if (.not. associated(phys_get_flux)) &
414 call mpistop("Error: no phys_get_flux not defined")
415
416 if (.not. associated(phys_get_dt)) &
417 call mpistop("Error: no phys_get_dt not defined")
418
419 if (.not. associated(phys_add_source_geom)) &
421
422 if (.not. associated(phys_add_source)) &
424
425 if (.not. associated(phys_check_w)) &
427
428 if (.not. associated(phys_get_pthermal)) &
430 if (.not. associated(phys_get_auxiliary)) &
432 if (.not. associated(phys_get_auxiliary_prim)) &
434
435 if (.not. associated(phys_boundary_adjust)) &
437
438 if (.not. associated(phys_write_info)) &
440
441 if (.not. associated(phys_handle_small_values)) &
443
444 if (.not. associated(phys_get_ct_velocity)) &
446
447 if (.not. associated(phys_update_faces)) &
449
450 if (.not. associated(phys_face_to_center)) &
452
453 if (.not. associated(phys_evaluate_implicit)) &
455
456 if (.not. associated(phys_implicit_update)) &
458
459 end subroutine phys_check
460
462 end subroutine dummy_init_params
463
465 end subroutine dummy_check_params
466
467 subroutine dummy_modify_wlr(ixI^L, ixO^L, qt, wLC, wRC, wLp, wRp, s, idir)
469 integer, intent(in) :: ixI^L, ixO^L, idir
470 double precision, intent(in) :: qt
471 double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
472 double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
473 type(state) :: s
474 end subroutine dummy_modify_wlr
475
476 subroutine dummy_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
478 integer, intent(in) :: ixI^L, ixO^L, idim
479 double precision, intent(in) :: wprim(ixI^S, nw)
480 double precision, intent(in) :: x(ixI^S,1:ndim)
481 double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
482 end subroutine dummy_get_h_speed
483
484 subroutine dummy_get_a2max(w, x, ixI^L, ixO^L, a2max)
486 use mod_comm_lib, only: mpistop
487 integer, intent(in) :: ixI^L, ixO^L
488 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
489 double precision, intent(inout) :: a2max(ndim)
490 call mpistop("Error: entered dummy_get_a2max")
491 end subroutine dummy_get_a2max
492
493 subroutine dummy_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
495 integer, intent(in) :: ixI^L, ixO^L
496 double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
497 double precision, intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S,1:nw),w(ixI^S, 1:nw)
498 end subroutine dummy_add_source_geom
499
500 subroutine dummy_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
501 qsourcesplit, active)
503 integer, intent(in) :: ixI^L, ixO^L
504 double precision, intent(in) :: qdt,dtfactor
505 double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
506 double precision, intent(inout) :: w(ixI^S, 1:nw)
507 logical, intent(in) :: qsourcesplit
508 logical, intent(inout) :: active
509 ! Don't have to set active, since it starts as .false.
510 end subroutine dummy_add_source
511
512 subroutine dummy_check_w(primitive, ixI^L, ixO^L, w, w_flag)
514 logical, intent(in) :: primitive
515 integer, intent(in) :: ixI^L, ixO^L
516 double precision, intent(in) :: w(ixI^S,1:nw)
517 logical, intent(inout) :: w_flag(ixI^S,1:nw)
518
519 w_flag=.false. ! All okay
520 end subroutine dummy_check_w
521
522 subroutine dummy_get_pthermal(w, x, ixI^L, ixO^L, pth)
524 use mod_comm_lib, only: mpistop
525
526 integer, intent(in) :: ixI^L, ixO^L
527 double precision, intent(in) :: w(ixI^S, nw)
528 double precision, intent(in) :: x(ixI^S, 1:ndim)
529 double precision, intent(out):: pth(ixI^S)
530
531 call mpistop("No get_pthermal method specified")
532 end subroutine dummy_get_pthermal
533
534 subroutine dummy_get_auxiliary(ixI^L, ixO^L, w, x)
536 use mod_comm_lib, only: mpistop
537
538 integer, intent(in) :: ixI^L, ixO^L
539 double precision, intent(inout) :: w(ixI^S, nw)
540 double precision, intent(in) :: x(ixI^S, 1:ndim)
541
542 !call mpistop("No get_auxiliary method specified")
543 end subroutine dummy_get_auxiliary
544
545 subroutine dummy_get_auxiliary_prim(ixI^L, ixO^L, w)
547 use mod_comm_lib, only: mpistop
548
549 integer, intent(in) :: ixI^L, ixO^L
550 double precision, intent(inout) :: w(ixI^S, nw)
551
552 call mpistop("No get_auxiliary_prim method specified")
553 end subroutine dummy_get_auxiliary_prim
554
555 subroutine dummy_boundary_adjust(igrid,psb)
557 integer, intent(in) :: igrid
558 type(state), target :: psb(max_blocks)
559 end subroutine dummy_boundary_adjust
560
561 subroutine dummy_write_info(fh)
563 integer, intent(in) :: fh !< File handle
564 integer, dimension(MPI_STATUS_SIZE) :: st
565 integer :: er
566
567 ! Number of physics parameters
568 integer, parameter :: n_par = 0
569
570 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
571 end subroutine dummy_write_info
572
573 subroutine dummy_small_values(primitive, w, x, ixI^L, ixO^L, subname)
575 logical, intent(in) :: primitive
576 integer, intent(in) :: ixI^L,ixO^L
577 double precision, intent(inout) :: w(ixI^S,1:nw)
578 double precision, intent(in) :: x(ixI^S,1:ndim)
579 character(len=*), intent(in) :: subname
580 end subroutine dummy_small_values
581
582 subroutine dummy_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
584
585 integer, intent(in) :: ixI^L, ixO^L, idim
586 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
587 double precision, intent(in) :: cmax(ixI^S)
588 double precision, intent(in), optional :: cmin(ixI^S)
589 type(ct_velocity), intent(inout):: vcts
590 end subroutine dummy_get_ct_velocity
591
592 subroutine dummy_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
594 integer, intent(in) :: ixI^L, ixO^L
595 double precision, intent(in) :: qt, qdt
596 ! cell-center primitive variables
597 double precision, intent(in) :: wprim(ixI^S,1:nw)
598 type(state) :: sCT, s
599 type(ct_velocity) :: vcts
600 double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
601 double precision, intent(inout) :: fE(ixI^S,sdim:3)
602 end subroutine dummy_update_faces
603
604 subroutine dummy_face_to_center(ixO^L,s)
606 integer, intent(in) :: ixO^L
607 type(state) :: s
608 end subroutine dummy_face_to_center
609
610 subroutine dummy_evaluate_implicit(qtC,psa)
612 type(state), target :: psa(max_blocks)
613 double precision, intent(in) :: qtC
614 integer :: iigrid, igrid
615
616 ! Just copy in nul state
617 !$OMP PARALLEL DO PRIVATE(igrid)
618 do iigrid=1,igridstail; igrid=igrids(iigrid);
619 psa(igrid)%w = 0.0d0*psa(igrid)%w
620 if(stagger_grid) psa(igrid)%ws = 0.0d0*psa(igrid)%ws
621 end do
622 !$OMP END PARALLEL DO
623
624 end subroutine dummy_evaluate_implicit
625
626 subroutine dummy_implicit_update(dtfactor,qdt,qtC,psa,psb)
628 type(state), target :: psa(max_blocks)
629 type(state), target :: psb(max_blocks)
630 double precision, intent(in) :: qdt
631 double precision, intent(in) :: qtC
632 double precision, intent(in) :: dtfactor
633 integer :: iigrid, igrid
634
635 ! Just copy in psb state when using the scheme without implicit part
636 !$OMP PARALLEL DO PRIVATE(igrid)
637 do iigrid=1,igridstail; igrid=igrids(iigrid);
638 psa(igrid)%w = psb(igrid)%w
639 if(stagger_grid) psa(igrid)%ws = psb(igrid)%ws
640 end do
641 !$OMP END PARALLEL DO
642
643 end subroutine dummy_implicit_update
644
645end 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:54
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: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: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:51
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
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: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:53
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition mod_physics.t:49
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:52
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_get_cs2), pointer phys_get_cs2
Definition mod_physics.t:57
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:55
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:43
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:56
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:37
procedure(sub_special_advance), pointer phys_special_advance
Definition mod_physics.t:70
Module with all the methods that users can customize in AMRVAC.