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