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_boundary_adjust), pointer :: phys_boundary_adjust => null()
73 procedure(sub_write_info), pointer :: phys_write_info => null()
74 procedure(sub_small_values), pointer :: phys_handle_small_values => null()
75 procedure(sub_get_ct_velocity), pointer :: phys_get_ct_velocity => null()
76 procedure(sub_update_faces), pointer :: phys_update_faces => null()
77 procedure(sub_face_to_center), pointer :: phys_face_to_center => null()
78 procedure(sub_implicit_update), pointer :: phys_implicit_update => null()
79 procedure(sub_evaluate_implicit),pointer:: phys_evaluate_implicit => null()
80 procedure(sub_clean_divb), pointer :: phys_clean_divb => null()
81 ! set the equilibrium variables
82 procedure(sub_set_equi_vars), pointer :: phys_set_equi_vars => null()
83 ! subroutine with no parameters which creates EUV images
84 procedure(sub_check_params), pointer :: phys_te_images => null()
85 ! to update temperature variable with partial ionization
86 procedure(sub_update_temperature), pointer :: phys_update_temperature => null()
87 procedure(sub_get_auxiliary), pointer :: phys_get_auxiliary => null()
88 procedure(sub_get_auxiliary_prim), pointer :: phys_get_auxiliary_prim => null()
89
90 abstract interface
91
93 end subroutine sub_check_params
94
98 end subroutine sub_set_mg_bounds
99
100 subroutine sub_boundary_adjust(igrid,psb)
102 integer, intent(in) :: igrid
103 type(state), target :: psb(max_blocks)
104 end subroutine sub_boundary_adjust
105
106 subroutine sub_convert(ixI^L, ixO^L, w, x)
108 integer, intent(in) :: ixI^L, ixO^L
109 double precision, intent(inout) :: w(ixI^S, nw)
110 double precision, intent(in) :: x(ixI^S, 1:^ND)
111 end subroutine sub_convert
112
113 subroutine sub_modify_wlr(ixI^L, ixO^L, qt, wLC, wRC, wLp, wRp, s, idir)
115 integer, intent(in) :: ixI^L, ixO^L, idir
116 double precision, intent(in) :: qt
117 double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
118 double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
119 type(state) :: s
120 end subroutine sub_modify_wlr
121
122 subroutine sub_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
124 integer, intent(in) :: ixI^L, ixO^L, idim
125 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:^ND)
126 double precision, intent(inout) :: cmax(ixI^S)
127 end subroutine sub_get_cmax
128
129 subroutine sub_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
131 integer, intent(in) :: ixI^L, ixO^L
132 double precision, intent(inout) :: w(ixI^S, nw)
133 double precision, intent(in) :: x(ixI^S, 1:^ND)
134 double precision, intent(out) :: tco_local, Tmax_local
135 end subroutine sub_get_tcutoff
136
137 subroutine sub_trac_after_setdt(trac_alfa,tco,T_peak,T_bott)
138 double precision, intent(in) :: trac_alfa,tco,T_peak,T_bott
139 end subroutine sub_trac_after_setdt
140
141 subroutine sub_get_v(w,x,ixI^L,ixO^L,v)
143 integer, intent(in) :: ixI^L, ixO^L
144 double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
145 double precision, intent(out) :: v(ixI^S,1:ndir)
146 end subroutine sub_get_v
147
148 subroutine sub_get_rho(w,x,ixI^L,ixO^L,rho)
150 integer, intent(in) :: ixI^L, ixO^L
151 double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
152 double precision, intent(out) :: rho(ixI^S)
153 end subroutine sub_get_rho
154
155 subroutine sub_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
157 integer, intent(in) :: ixI^L, ixO^L, idim
158 double precision, intent(in) :: wprim(ixI^S, nw)
159 double precision, intent(in) :: x(ixI^S,1:ndim)
160 double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
161 end subroutine sub_get_h_speed
162
163 subroutine sub_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim, Hspeed, cmax, cmin)
165 use mod_variables
166 integer, intent(in) :: ixI^L, ixO^L, idim
167 double precision, intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
168 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
169 double precision, intent(in) :: x(ixI^S, 1:^ND)
170 double precision, intent(inout) :: cmax(ixI^S,1:number_species)
171 double precision, intent(inout), optional :: cmin(ixI^S,1:number_species)
172 double precision, intent(in) :: Hspeed(ixI^S,1:number_species)
173 end subroutine sub_get_cbounds
174
175 subroutine sub_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
177 integer, intent(in) :: ixI^L, ixO^L, idim
178 double precision, intent(in) :: wC(ixI^S, 1:nw)
179 double precision, intent(in) :: w(ixI^S, 1:nw)
180 double precision, intent(in) :: x(ixI^S, 1:^ND)
181 double precision, intent(out) :: f(ixI^S, nwflux)
182 end subroutine sub_get_flux
183
184 subroutine sub_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
186 integer, intent(in) :: ixI^L, ixO^L
187 double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
188 double precision, intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S, 1:nw), w(ixI^S, 1:nw)
189 end subroutine sub_add_source_geom
190
191 subroutine sub_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
192 qsourcesplit, active)
194 integer, intent(in) :: ixI^L, ixO^L
195 double precision, intent(in) :: qdt, dtfactor
196 double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
197 double precision, intent(inout) :: w(ixI^S, 1:nw)
198 logical, intent(in) :: qsourcesplit
199 logical, intent(inout) :: active !< Needs to be set to true when active
200 end subroutine sub_add_source
201
202 !> Add global source terms on complete domain (potentially implicit)
203 subroutine sub_global_source(qdt, qt, active)
205 double precision, intent(in) :: qdt !< Current time step
206 double precision, intent(in) :: qt !< Current time
207 logical, intent(inout) :: active !< Output if the source is active
208 end subroutine sub_global_source
209
210 !> clean initial divb
211 subroutine sub_clean_divb(qdt, qt, active)
213 double precision, intent(in) :: qdt !< Current time step
214 double precision, intent(in) :: qt !< Current time
215 logical, intent(inout) :: active !< Output if the source is active
216 end subroutine sub_clean_divb
217
218 !> set equilibrium variables other than b0 (e.g. p0 and rho0)
219 subroutine sub_set_equi_vars(igrid)
220 integer, intent(in) :: igrid
221 end subroutine sub_set_equi_vars
222
223
224 !> Add special advance in each advect step
225 subroutine sub_special_advance(qt, psa)
227 double precision, intent(in) :: qt !< Current time
228 type(state), target :: psa(max_blocks) !< Compute based on this state
229 end subroutine sub_special_advance
230
231 subroutine sub_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
233 integer, intent(in) :: ixI^L, ixO^L
234 double precision, intent(in) :: dx^D, x(ixI^S, 1:^ND)
235 double precision, intent(in) :: w(ixI^S, 1:nw)
236 double precision, intent(inout) :: dtnew
237 end subroutine sub_get_dt
238
239 subroutine sub_check_w(primitive, ixI^L, ixO^L, w, w_flag)
241 logical, intent(in) :: primitive
242 integer, intent(in) :: ixI^L, ixO^L
243 double precision, intent(in) :: w(ixI^S,1:nw)
244 logical, intent(inout) :: w_flag(ixI^S,1:nw)
245 end subroutine sub_check_w
246
247 subroutine sub_get_pthermal(w,x,ixI^L,ixO^L,pth)
249 integer, intent(in) :: ixI^L, ixO^L
250 double precision, intent(in) :: w(ixI^S,nw)
251 double precision, intent(in) :: x(ixI^S,1:ndim)
252 double precision, intent(out):: pth(ixI^S)
253 end subroutine sub_get_pthermal
254
255 subroutine sub_get_auxiliary(ixI^L,ixO^L,w,x)
257 integer, intent(in) :: ixI^L, ixO^L
258 double precision, intent(inout) :: w(ixI^S,nw)
259 double precision, intent(in) :: x(ixI^S,1:ndim)
260 end subroutine sub_get_auxiliary
261
262 subroutine sub_get_auxiliary_prim(ixI^L,ixO^L,w)
264 integer, intent(in) :: ixI^L, ixO^L
265 double precision, intent(inout) :: w(ixI^S,nw)
266 end subroutine sub_get_auxiliary_prim
267
268 subroutine sub_get_tgas(w,x,ixI^L,ixO^L,tgas)
270 integer, intent(in) :: ixI^L, ixO^L
271 double precision, intent(in) :: w(ixI^S,nw)
272 double precision, intent(in) :: x(ixI^S,1:ndim)
273 double precision, intent(out):: tgas(ixI^S)
274 end subroutine sub_get_tgas
275
276 subroutine sub_write_info(file_handle)
277 integer, intent(in) :: file_handle
278 end subroutine sub_write_info
279
280 subroutine sub_small_values(primitive, w, x, ixI^L, ixO^L, subname)
282 logical, intent(in) :: primitive
283 integer, intent(in) :: ixI^L,ixO^L
284 double precision, intent(inout) :: w(ixI^S,1:nw)
285 double precision, intent(in) :: x(ixI^S,1:ndim)
286 character(len=*), intent(in) :: subname
287 end subroutine sub_small_values
288
289 subroutine sub_get_var(ixI^L, ixO^L, w, out)
291 integer, intent(in) :: ixI^L, ixO^L
292 double precision, intent(in) :: w(ixI^S, nw)
293 double precision, intent(out) :: out(ixO^S)
294 end subroutine sub_get_var
295
296 subroutine sub_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
298
299 integer, intent(in) :: ixI^L, ixO^L, idim
300 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
301 double precision, intent(in) :: cmax(ixI^S)
302 double precision, intent(in), optional :: cmin(ixI^S)
303 type(ct_velocity), intent(inout):: vcts
304 end subroutine sub_get_ct_velocity
305
306 subroutine sub_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
308 integer, intent(in) :: ixI^L, ixO^L
309 double precision, intent(in) :: qt, qdt
310 ! cell-center primitive variables
311 double precision, intent(in) :: wprim(ixI^S,1:nw)
312 ! velocity structure
313 type(state) :: sCT, s
314 type(ct_velocity) :: vcts
315 double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
316 double precision, intent(inout) :: fE(ixI^S,sdim:3)
317 end subroutine sub_update_faces
318
319 subroutine sub_face_to_center(ixO^L,s)
321 integer, intent(in) :: ixO^L
322 type(state) :: s
323 end subroutine sub_face_to_center
324
325 subroutine sub_evaluate_implicit(qtC,psa)
327 type(state), target :: psa(max_blocks)
328 double precision, intent(in) :: qtC
329 end subroutine sub_evaluate_implicit
330
331 subroutine sub_implicit_update(dtfactor,qdt,qtC,psa,psb)
333 type(state), target :: psa(max_blocks)
334 type(state), target :: psb(max_blocks)
335 double precision, intent(in) :: qdt
336 double precision, intent(in) :: qtC
337 double precision, intent(in) :: dtfactor
338 end subroutine sub_implicit_update
339
340 subroutine sub_update_temperature(ixI^L,ixO^L,wCT,w,x)
342 integer, intent(in) :: ixI^L, ixO^L
343 double precision, intent(in) :: wCT(ixI^S,nw),x(ixI^S,1:ndim)
344 double precision, intent(inout) :: w(ixI^S,nw)
345 end subroutine sub_update_temperature
346
347 end interface
348
349contains
350
351 subroutine phys_check()
352 use mod_global_parameters, only: nw, ndir
353
356 use mod_comm_lib, only: mpistop
357
358 if (physics_type == "") call mpistop("Error: no physics module loaded")
359
360 call phys_hllc_check()
361 call phys_roe_check()
362
363 ! Checks whether the required physics methods have been defined
364 if (.not. associated(phys_check_params)) &
366
367 if (.not. associated(phys_to_conserved)) &
368 call mpistop("Error: phys_to_conserved not defined")
369
370 if (.not. associated(phys_to_primitive)) &
371 call mpistop("Error: phys_to_primitive not defined")
372
373 if (.not. associated(phys_modify_wlr)) &
375
376 if (.not. associated(phys_get_cmax)) &
377 call mpistop("Error: no phys_get_cmax not defined")
378
379 if (.not. associated(phys_get_h_speed)) &
381
382 if (.not. associated(phys_get_cbounds)) &
383 call mpistop("Error: no phys_get_cbounds not defined")
384
385 if (.not. associated(phys_get_flux)) &
386 call mpistop("Error: no phys_get_flux not defined")
387
388 if (.not. associated(phys_get_dt)) &
389 call mpistop("Error: no phys_get_dt not defined")
390
391 if (.not. associated(phys_add_source_geom)) &
393
394 if (.not. associated(phys_add_source)) &
396
397 if (.not. associated(phys_check_w)) &
399
400 if (.not. associated(phys_get_pthermal)) &
402 if (.not. associated(phys_get_auxiliary)) &
404 if (.not. associated(phys_get_auxiliary_prim)) &
406
407 if (.not. associated(phys_boundary_adjust)) &
409
410 if (.not. associated(phys_write_info)) &
412
413 if (.not. associated(phys_handle_small_values)) &
415
416 if (.not. associated(phys_get_ct_velocity)) &
418
419 if (.not. associated(phys_update_faces)) &
421
422 if (.not. associated(phys_face_to_center)) &
424
425 if (.not. associated(phys_evaluate_implicit)) &
427
428 if (.not. associated(phys_implicit_update)) &
430
431 end subroutine phys_check
432
434 end subroutine dummy_init_params
435
437 end subroutine dummy_check_params
438
439 subroutine dummy_modify_wlr(ixI^L, ixO^L, qt, wLC, wRC, wLp, wRp, s, idir)
441 integer, intent(in) :: ixI^L, ixO^L, idir
442 double precision, intent(in) :: qt
443 double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
444 double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
445 type(state) :: s
446 end subroutine dummy_modify_wlr
447
448 subroutine dummy_get_h_speed(wprim,x,ixI^L,ixO^L,idim,Hspeed)
450 integer, intent(in) :: ixI^L, ixO^L, idim
451 double precision, intent(in) :: wprim(ixI^S, nw)
452 double precision, intent(in) :: x(ixI^S,1:ndim)
453 double precision, intent(out) :: Hspeed(ixI^S,1:number_species)
454 end subroutine dummy_get_h_speed
455
456 subroutine dummy_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
458 integer, intent(in) :: ixI^L, ixO^L
459 double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:^ND)
460 double precision, intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S,1:nw),w(ixI^S, 1:nw)
461 end subroutine dummy_add_source_geom
462
463 subroutine dummy_add_source(qdt, dtfactor, ixI^L, ixO^L, wCT, wCTprim, w, x, &
464 qsourcesplit, active)
466 integer, intent(in) :: ixI^L, ixO^L
467 double precision, intent(in) :: qdt,dtfactor
468 double precision, intent(in) :: wCT(ixI^S, 1:nw), wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
469 double precision, intent(inout) :: w(ixI^S, 1:nw)
470 logical, intent(in) :: qsourcesplit
471 logical, intent(inout) :: active
472 ! Don't have to set active, since it starts as .false.
473 end subroutine dummy_add_source
474
475 subroutine dummy_check_w(primitive, ixI^L, ixO^L, w, w_flag)
477 logical, intent(in) :: primitive
478 integer, intent(in) :: ixI^L, ixO^L
479 double precision, intent(in) :: w(ixI^S,1:nw)
480 logical, intent(inout) :: w_flag(ixI^S,1:nw)
481
482 w_flag=.false. ! All okay
483 end subroutine dummy_check_w
484
485 subroutine dummy_get_pthermal(w, x, ixI^L, ixO^L, pth)
487 use mod_comm_lib, only: mpistop
488
489 integer, intent(in) :: ixI^L, ixO^L
490 double precision, intent(in) :: w(ixI^S, nw)
491 double precision, intent(in) :: x(ixI^S, 1:ndim)
492 double precision, intent(out):: pth(ixI^S)
493
494 call mpistop("No get_pthermal method specified")
495 end subroutine dummy_get_pthermal
496
497 subroutine dummy_get_auxiliary(ixI^L, ixO^L, w, x)
499 use mod_comm_lib, only: mpistop
500
501 integer, intent(in) :: ixI^L, ixO^L
502 double precision, intent(inout) :: w(ixI^S, nw)
503 double precision, intent(in) :: x(ixI^S, 1:ndim)
504
505 !call mpistop("No get_auxiliary method specified")
506 end subroutine dummy_get_auxiliary
507
508 subroutine dummy_get_auxiliary_prim(ixI^L, ixO^L, w)
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
515 call mpistop("No get_auxiliary_prim method specified")
516 end subroutine dummy_get_auxiliary_prim
517
518 subroutine dummy_boundary_adjust(igrid,psb)
520 integer, intent(in) :: igrid
521 type(state), target :: psb(max_blocks)
522 end subroutine dummy_boundary_adjust
523
524 subroutine dummy_write_info(fh)
526 integer, intent(in) :: fh !< File handle
527 integer, dimension(MPI_STATUS_SIZE) :: st
528 integer :: er
529
530 ! Number of physics parameters
531 integer, parameter :: n_par = 0
532
533 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
534 end subroutine dummy_write_info
535
536 subroutine dummy_small_values(primitive, w, x, ixI^L, ixO^L, subname)
538 logical, intent(in) :: primitive
539 integer, intent(in) :: ixI^L,ixO^L
540 double precision, intent(inout) :: w(ixI^S,1:nw)
541 double precision, intent(in) :: x(ixI^S,1:ndim)
542 character(len=*), intent(in) :: subname
543 end subroutine dummy_small_values
544
545 subroutine dummy_get_ct_velocity(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
547
548 integer, intent(in) :: ixI^L, ixO^L, idim
549 double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
550 double precision, intent(in) :: cmax(ixI^S)
551 double precision, intent(in), optional :: cmin(ixI^S)
552 type(ct_velocity), intent(inout):: vcts
553 end subroutine dummy_get_ct_velocity
554
555 subroutine dummy_update_faces(ixI^L,ixO^L,qt,qdt,wprim,fC,fE,sCT,s,vcts)
557 integer, intent(in) :: ixI^L, ixO^L
558 double precision, intent(in) :: qt, qdt
559 ! cell-center primitive variables
560 double precision, intent(in) :: wprim(ixI^S,1:nw)
561 type(state) :: sCT, s
562 type(ct_velocity) :: vcts
563 double precision, intent(in) :: fC(ixI^S,1:nwflux,1:ndim)
564 double precision, intent(inout) :: fE(ixI^S,sdim:3)
565 end subroutine dummy_update_faces
566
567 subroutine dummy_face_to_center(ixO^L,s)
569 integer, intent(in) :: ixO^L
570 type(state) :: s
571 end subroutine dummy_face_to_center
572
573 subroutine dummy_evaluate_implicit(qtC,psa)
575 type(state), target :: psa(max_blocks)
576 double precision, intent(in) :: qtC
577 integer :: iigrid, igrid
578
579 ! Just copy in nul state
580 !$OMP PARALLEL DO PRIVATE(igrid)
581 do iigrid=1,igridstail; igrid=igrids(iigrid);
582 psa(igrid)%w = 0.0d0*psa(igrid)%w
583 if(stagger_grid) psa(igrid)%ws = 0.0d0*psa(igrid)%ws
584 end do
585 !$OMP END PARALLEL DO
586
587 end subroutine dummy_evaluate_implicit
588
589 subroutine dummy_implicit_update(dtfactor,qdt,qtC,psa,psb)
591 type(state), target :: psa(max_blocks)
592 type(state), target :: psb(max_blocks)
593 double precision, intent(in) :: qdt
594 double precision, intent(in) :: qtC
595 double precision, intent(in) :: dtfactor
596 integer :: iigrid, igrid
597
598 ! Just copy in psb state when using the scheme without implicit part
599 !$OMP PARALLEL DO PRIVATE(igrid)
600 do iigrid=1,igridstail; igrid=igrids(iigrid);
601 psa(igrid)%w = psb(igrid)%w
602 if(stagger_grid) psa(igrid)%ws = psb(igrid)%ws
603 end do
604 !$OMP END PARALLEL DO
605
606 end subroutine dummy_implicit_update
607
608end 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:75
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:74
procedure(sub_write_info), pointer phys_write_info
Definition mod_physics.t:73
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:79
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:84
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:86
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:87
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:88
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:82
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:76
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:78
subroutine dummy_evaluate_implicit(qtc, psa)
subroutine dummy_check_params
procedure(sub_clean_divb), pointer phys_clean_divb
Definition mod_physics.t:80
procedure(sub_boundary_adjust), pointer phys_boundary_adjust
Definition mod_physics.t:72
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:77
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()
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.