MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_hd_phys.t
Go to the documentation of this file.
1!> Hydrodynamics physics module
7 use mod_comm_lib, only: mpistop
8 implicit none
9 private
10
11 !> Whether an energy equation is used
12 logical, public, protected :: hd_energy = .true.
13
14 !> Whether thermal conduction is added
15 logical, public, protected :: hd_thermal_conduction = .false.
16 type(tc_fluid), allocatable, public :: tc_fl
17 type(te_fluid), allocatable, public :: te_fl_hd
18
19 !> Whether radiative cooling is added
20 logical, public, protected :: hd_radiative_cooling = .false.
21 type(rc_fluid), allocatable, public :: rc_fl
22
23 !> Whether dust is added
24 logical, public, protected :: hd_dust = .false.
25
26 !> Whether viscosity is added
27 logical, public, protected :: hd_viscosity = .false.
28
29 !> Whether gravity is added
30 logical, public, protected :: hd_gravity = .false.
31
32 !> Whether particles module is added
33 logical, public, protected :: hd_particles = .false.
34
35 !> Whether rotating frame is activated
36 logical, public, protected :: hd_rotating_frame = .false.
37
38 !> Whether CAK radiation line force is activated
39 logical, public, protected :: hd_cak_force = .false.
40
41 !> Number of tracer species
42 integer, public, protected :: hd_n_tracer = 0
43
44 !> Whether plasma is partially ionized
45 logical, public, protected :: hd_partial_ionization = .false.
46
47 !> Index of the density (in the w array)
48 integer, public, protected :: rho_
49
50 !> Indices of the momentum density
51 integer, allocatable, public, protected :: mom(:)
52
53 !> Indices of the tracers
54 integer, allocatable, public, protected :: tracer(:)
55
56 !> Index of the energy density (-1 if not present)
57 integer, public, protected :: e_
58
59 !> Index of the gas pressure (-1 if not present) should equal e_
60 integer, public, protected :: p_
61
62 !> Indices of temperature
63 integer, public, protected :: te_
64
65 !> Index of the cutoff temperature for the TRAC method
66 integer, public, protected :: tcoff_
67
68 !> The adiabatic index
69 double precision, public :: hd_gamma = 5.d0/3.0d0
70
71 !> The adiabatic constant
72 double precision, public :: hd_adiab = 1.0d0
73
74 !> The small_est allowed energy
75 double precision, protected :: small_e
76
77 !> Whether TRAC method is used
78 logical, public, protected :: hd_trac = .false.
79 integer, public, protected :: hd_trac_type = 1
80
81 !> Helium abundance over Hydrogen
82 double precision, public, protected :: he_abundance=0.1d0
83 !> Ionization fraction of H
84 !> H_ion_fr = H+/(H+ + H)
85 double precision, public, protected :: h_ion_fr=1d0
86 !> Ionization fraction of He
87 !> He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
88 double precision, public, protected :: he_ion_fr=1d0
89 !> Ratio of number He2+ / number He+ + He2+
90 !> He_ion_fr2 = He2+/(He2+ + He+)
91 double precision, public, protected :: he_ion_fr2=1d0
92 ! used for eq of state when it is not defined by units,
93 ! the units do not contain terms related to ionization fraction
94 ! and it is p = RR * rho * T
95 double precision, public, protected :: rr=1d0
96 ! remove the below flag and assume default value = .false.
97 ! when eq state properly implemented everywhere
98 ! and not anymore through units
99 logical, public, protected :: eq_state_units = .true.
100
101 procedure(sub_get_pthermal), pointer :: hd_get_rfactor => null()
102 ! Public methods
103 public :: hd_phys_init
104 public :: hd_kin_en
105 public :: hd_get_pthermal
106 public :: hd_get_csound2
107 public :: hd_to_conserved
108 public :: hd_to_primitive
109 public :: hd_check_params
110 public :: hd_check_w
111
112contains
113
114 !> Read this module's parameters from a file
115 subroutine hd_read_params(files)
117 character(len=*), intent(in) :: files(:)
118 integer :: n
119
120 namelist /hd_list/ hd_energy, hd_n_tracer, hd_gamma, hd_adiab, &
125
126 do n = 1, size(files)
127 open(unitpar, file=trim(files(n)), status="old")
128 read(unitpar, hd_list, end=111)
129111 close(unitpar)
130 end do
131
132 end subroutine hd_read_params
133
134 !> Write this module's parameters to a snapsoht
135 subroutine hd_write_info(fh)
137 integer, intent(in) :: fh
138 integer, parameter :: n_par = 1
139 double precision :: values(n_par)
140 character(len=name_len) :: names(n_par)
141 integer, dimension(MPI_STATUS_SIZE) :: st
142 integer :: er
143
144 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
145
146 names(1) = "gamma"
147 values(1) = hd_gamma
148 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
149 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
150 end subroutine hd_write_info
151
152 !> Initialize the module
153 subroutine hd_phys_init()
157 use mod_dust, only: dust_init
159 use mod_gravity, only: gravity_init
162 use mod_cak_force, only: cak_init
167
168 integer :: itr, idir
169
170 call hd_read_params(par_files)
171
172 physics_type = "hd"
173 phys_energy = hd_energy
174 phys_total_energy = hd_energy
175 phys_internal_e = .false.
176 phys_gamma = hd_gamma
177 phys_partial_ionization=hd_partial_ionization
178
180 if(phys_trac) then
181 if(ndim .eq. 1) then
182 if(hd_trac_type .gt. 2) then
184 if(mype==0) write(*,*) 'WARNING: set hd_trac_type=1'
185 end if
187 else
188 phys_trac=.false.
189 if(mype==0) write(*,*) 'WARNING: set hd_trac=F when ndim>=2'
190 end if
191 end if
192
193 ! set default gamma for polytropic/isothermal process
194 if(.not.hd_energy) then
195 if(hd_thermal_conduction) then
197 if(mype==0) write(*,*) 'WARNING: set hd_thermal_conduction=F when hd_energy=F'
198 end if
199 if(hd_radiative_cooling) then
201 if(mype==0) write(*,*) 'WARNING: set hd_radiative_cooling=F when hd_energy=F'
202 end if
203 end if
204 if(.not.eq_state_units) then
205 if(hd_partial_ionization) then
207 if(mype==0) write(*,*) 'WARNING: set hd_partial_ionization=F when eq_state_units=F'
208 end if
209 end if
211
212 allocate(start_indices(number_species),stop_indices(number_species))
213
214 ! set the index of the first flux variable for species 1
215 start_indices(1)=1
216 ! Determine flux variables
217 rho_ = var_set_rho()
218
219 allocate(mom(ndir))
220 mom(:) = var_set_momentum(ndir)
221
222 ! Set index of energy variable
223 if (hd_energy) then
224 e_ = var_set_energy()
225 p_ = e_
226 else
227 e_ = -1
228 p_ = -1
229 end if
230
231 phys_get_dt => hd_get_dt
232 phys_get_cmax => hd_get_cmax
233 phys_get_a2max => hd_get_a2max
234 phys_get_tcutoff => hd_get_tcutoff
235 phys_get_cbounds => hd_get_cbounds
236 phys_get_flux => hd_get_flux
237 phys_add_source_geom => hd_add_source_geom
238 phys_add_source => hd_add_source
239 phys_to_conserved => hd_to_conserved
240 phys_to_primitive => hd_to_primitive
241 phys_check_params => hd_check_params
242 phys_check_w => hd_check_w
243 phys_get_pthermal => hd_get_pthermal
244 phys_get_v => hd_get_v
245 phys_get_rho => hd_get_rho
246 phys_write_info => hd_write_info
247 phys_handle_small_values => hd_handle_small_values
248
249 ! derive units from basic units
250 call hd_physical_units()
251
252 if (hd_dust) then
253 call dust_init(rho_, mom(:), e_)
254 endif
255
256 allocate(tracer(hd_n_tracer))
257
258 ! Set starting index of tracers
259 do itr = 1, hd_n_tracer
260 tracer(itr) = var_set_fluxvar("trc", "trp", itr, need_bc=.false.)
261 end do
262
263 ! set number of variables which need update ghostcells
264 nwgc=nwflux+nwaux
265
266 ! set the index of the last flux variable for species 1
267 stop_indices(1)=nwflux
268
269 ! set temperature as an auxiliary variable to get ionization degree
270 if(hd_partial_ionization) then
271 te_ = var_set_auxvar('Te','Te')
272 else
273 te_ = -1
274 end if
275
276 if(hd_trac) then
277 tcoff_ = var_set_wextra()
278 iw_tcoff=tcoff_
279 else
280 tcoff_ = -1
281 end if
282
283 ! choose Rfactor in ideal gas law
284 if(hd_partial_ionization) then
285 hd_get_rfactor=>rfactor_from_temperature_ionization
286 phys_update_temperature => hd_update_temperature
287 else if(associated(usr_rfactor)) then
288 hd_get_rfactor=>usr_rfactor
289 else
290 hd_get_rfactor=>rfactor_from_constant_ionization
291 end if
292
293 ! initialize thermal conduction module
294 if (hd_thermal_conduction) then
295 if (.not. hd_energy) &
296 call mpistop("thermal conduction needs hd_energy=T")
297
298 call sts_init()
300
301 allocate(tc_fl)
302 call tc_get_hd_params(tc_fl,tc_params_read_hd)
303 call add_sts_method(hd_get_tc_dt_hd,hd_sts_set_source_tc_hd,e_,1,e_,1,.false.)
304 call set_conversion_methods_to_head(hd_e_to_ei, hd_ei_to_e)
305 call set_error_handling_to_head(hd_tc_handle_small_e)
306 tc_fl%get_temperature_from_conserved => hd_get_temperature_from_etot
307 tc_fl%get_temperature_from_eint => hd_get_temperature_from_eint
308 tc_fl%get_rho => hd_get_rho
309 tc_fl%e_ = e_
310 end if
311
312 ! Initialize radiative cooling module
313 if (hd_radiative_cooling) then
314 if (.not. hd_energy) &
315 call mpistop("radiative cooling needs hd_energy=T")
317 allocate(rc_fl)
318 call radiative_cooling_init(rc_fl,rc_params_read)
319 rc_fl%get_rho => hd_get_rho
320 rc_fl%get_pthermal => hd_get_pthermal
321 rc_fl%get_var_Rfactor => hd_get_rfactor
322 rc_fl%e_ = e_
323 rc_fl%Tcoff_ = tcoff_
324 end if
325 allocate(te_fl_hd)
326 te_fl_hd%get_rho=> hd_get_rho
327 te_fl_hd%get_pthermal=> hd_get_pthermal
328 te_fl_hd%get_var_Rfactor => hd_get_rfactor
329{^ifthreed
330 phys_te_images => hd_te_images
331}
332 ! Initialize viscosity module
333 if (hd_viscosity) call viscosity_init(phys_wider_stencil)
334
335 ! Initialize gravity module
336 if (hd_gravity) call gravity_init()
337
338 ! Initialize rotating_frame module
340
341 ! Initialize CAK radiation force module
343
344 ! Initialize particles module
345 if (hd_particles) then
346 call particles_init()
347 end if
348
349 ! Check whether custom flux types have been defined
350 if (.not. allocated(flux_type)) then
351 allocate(flux_type(ndir, nw))
352 flux_type = flux_default
353 else if (any(shape(flux_type) /= [ndir, nw])) then
354 call mpistop("phys_check error: flux_type has wrong shape")
355 end if
356
357 nvector = 1 ! No. vector vars
358 allocate(iw_vector(nvector))
359 iw_vector(1) = mom(1) - 1
360 ! initialize ionization degree table
362
363 end subroutine hd_phys_init
364
365{^ifthreed
366 subroutine hd_te_images
369 select case(convert_type)
370 case('EIvtiCCmpi','EIvtuCCmpi')
372 case('ESvtiCCmpi','ESvtuCCmpi')
374 case('SIvtiCCmpi','SIvtuCCmpi')
376 case('WIvtiCCmpi','WIvtuCCmpi')
378 case default
379 call mpistop("Error in synthesize emission: Unknown convert_type")
380 end select
381 end subroutine hd_te_images
382}
383!!start th cond
384 ! wrappers for STS functions in thermal_conductivity module
385 ! which take as argument the tc_fluid (defined in the physics module)
386 subroutine hd_sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
390 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
391 double precision, intent(in) :: x(ixi^s,1:ndim)
392 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
393 double precision, intent(in) :: my_dt
394 logical, intent(in) :: fix_conserve_at_step
395 call sts_set_source_tc_hd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
396 end subroutine hd_sts_set_source_tc_hd
397
398 function hd_get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
399 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
400 !where tc_k_para_i=tc_k_para*B_i**2/B**2
401 !and T=p/rho
404
405 integer, intent(in) :: ixi^l, ixo^l
406 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
407 double precision, intent(in) :: w(ixi^s,1:nw)
408 double precision :: dtnew
409
410 dtnew=get_tc_dt_hd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
411 end function hd_get_tc_dt_hd
412
413 subroutine hd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
414 ! move this in a different routine as in mhd if needed in more places
417
418 integer, intent(in) :: ixi^l,ixo^l
419 double precision, intent(inout) :: w(ixi^s,1:nw)
420 double precision, intent(in) :: x(ixi^s,1:ndim)
421 integer, intent(in) :: step
422
423 integer :: idir
424 logical :: flag(ixi^s,1:nw)
425 character(len=140) :: error_msg
426
427 flag=.false.
428 where(w(ixo^s,e_)<small_e) flag(ixo^s,e_)=.true.
429 if(any(flag(ixo^s,e_))) then
430 select case (small_values_method)
431 case ("replace")
432 where(flag(ixo^s,e_)) w(ixo^s,e_)=small_e
433 case ("average")
434 call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
435 case default
436 ! small values error shows primitive variables
437 w(ixo^s,e_)=w(ixo^s,e_)*(hd_gamma - 1.0d0)
438 do idir = 1, ndir
439 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,rho_)
440 end do
441 write(error_msg,"(a,i3)") "Thermal conduction step ", step
442 call small_values_error(w, x, ixi^l, ixo^l, flag, error_msg)
443 end select
444 end if
445 end subroutine hd_tc_handle_small_e
446
447 ! fill in tc_fluid fields from namelist
448 subroutine tc_params_read_hd(fl)
450 type(tc_fluid), intent(inout) :: fl
451 integer :: n
452 logical :: tc_saturate=.false.
453 double precision :: tc_k_para=0d0
454
455 namelist /tc_list/ tc_saturate, tc_k_para
456
457 do n = 1, size(par_files)
458 open(unitpar, file=trim(par_files(n)), status="old")
459 read(unitpar, tc_list, end=111)
460111 close(unitpar)
461 end do
462 fl%tc_saturate = tc_saturate
463 fl%tc_k_para = tc_k_para
464
465 end subroutine tc_params_read_hd
466
467 subroutine hd_get_rho(w,x,ixI^L,ixO^L,rho)
469 integer, intent(in) :: ixi^l, ixo^l
470 double precision, intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:ndim)
471 double precision, intent(out) :: rho(ixi^s)
472
473 rho(ixo^s) = w(ixo^s,rho_)
474
475 end subroutine hd_get_rho
476
477!!end th cond
478!!rad cool
479 subroutine rc_params_read(fl)
481 use mod_constants, only: bigdouble
482 use mod_basic_types, only: std_len
483 type(rc_fluid), intent(inout) :: fl
484 integer :: n
485 ! list parameters
486 integer :: ncool = 4000
487 double precision :: cfrac=0.1d0
488
489 !> Name of cooling curve
490 character(len=std_len) :: coolcurve='JCcorona'
491
492 !> Name of cooling method
493 character(len=std_len) :: coolmethod='exact'
494
495 !> Fixed temperature not lower than tlow
496 logical :: tfix=.false.
497
498 !> Lower limit of temperature
499 double precision :: tlow=bigdouble
500
501 !> Add cooling source in a split way (.true.) or un-split way (.false.)
502 logical :: rc_split=.false.
503
504
505 namelist /rc_list/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split
506
507 do n = 1, size(par_files)
508 open(unitpar, file=trim(par_files(n)), status="old")
509 read(unitpar, rc_list, end=111)
510111 close(unitpar)
511 end do
512
513 fl%ncool=ncool
514 fl%coolcurve=coolcurve
515 fl%coolmethod=coolmethod
516 fl%tlow=tlow
517 fl%Tfix=tfix
518 fl%rc_split=rc_split
519 fl%cfrac=cfrac
520 end subroutine rc_params_read
521!! end rad cool
522
526
527 if (.not. hd_energy) then
528 if (hd_gamma <= 0.0d0) call mpistop ("Error: hd_gamma <= 0")
529 if (hd_adiab < 0.0d0) call mpistop ("Error: hd_adiab < 0")
531 else
532 if (hd_gamma <= 0.0d0 .or. hd_gamma == 1.0d0) &
533 call mpistop ("Error: hd_gamma <= 0 or hd_gamma == 1.0")
534 small_e = small_pressure/(hd_gamma - 1.0d0)
535 end if
536
537 if (hd_dust) call dust_check_params()
538 if(use_imex_scheme) then
539 ! implicit dust update
540 phys_implicit_update => dust_implicit_update
541 phys_evaluate_implicit => dust_evaluate_implicit
542 endif
543
544 end subroutine hd_check_params
545
546 subroutine hd_physical_units
548 double precision :: mp,kb
549 double precision :: a,b
550 ! Derive scaling units
551 if(si_unit) then
552 mp=mp_si
553 kb=kb_si
554 else
555 mp=mp_cgs
556 kb=kb_cgs
557 end if
558 if(eq_state_units) then
559 a = 1d0 + 4d0 * he_abundance
560 if(hd_partial_ionization) then
561 b = 2.d0+3.d0*he_abundance
562 else
563 b = 1d0 + h_ion_fr + he_abundance*(he_ion_fr*(he_ion_fr2 + 1d0)+1d0)
564 end if
565 rr = 1d0
566 else
567 a = 1d0
568 b = 1d0
569 rr = (1d0 + h_ion_fr + he_abundance*(he_ion_fr*(he_ion_fr2 + 1d0)+1d0))/(1d0 + 4d0 * he_abundance)
570 end if
571 if(unit_density/=1.d0) then
573 else
574 ! unit of numberdensity is independent by default
576 end if
577 if(unit_velocity/=1.d0) then
580 else if(unit_pressure/=1.d0) then
583 else
584 ! unit of temperature is independent by default
587 end if
588 if(unit_time/=1.d0) then
590 else
591 ! unit of length is independent by default
593 end if
595
596 end subroutine hd_physical_units
597
598 !> Returns logical argument flag where values are ok
599 subroutine hd_check_w(primitive, ixI^L, ixO^L, w, flag)
601 use mod_dust, only: dust_check_w
602
603 logical, intent(in) :: primitive
604 integer, intent(in) :: ixi^l, ixo^l
605 double precision, intent(in) :: w(ixi^s, nw)
606 logical, intent(inout) :: flag(ixi^s,1:nw)
607 double precision :: tmp(ixi^s)
608
609 flag=.false.
610
611 if (hd_energy) then
612 if (primitive) then
613 where(w(ixo^s, e_) < small_pressure) flag(ixo^s,e_) = .true.
614 else
615 tmp(ixo^s) = (hd_gamma - 1.0d0)*(w(ixo^s, e_) - &
616 hd_kin_en(w, ixi^l, ixo^l))
617 where(tmp(ixo^s) < small_pressure) flag(ixo^s,e_) = .true.
618 endif
619 end if
620
621 where(w(ixo^s, rho_) < small_density) flag(ixo^s,rho_) = .true.
622
623 if(hd_dust) call dust_check_w(ixi^l,ixo^l,w,flag)
624
625 end subroutine hd_check_w
626
627 !> Transform primitive variables into conservative ones
628 subroutine hd_to_conserved(ixI^L, ixO^L, w, x)
630 use mod_dust, only: dust_to_conserved
631 integer, intent(in) :: ixi^l, ixo^l
632 double precision, intent(inout) :: w(ixi^s, nw)
633 double precision, intent(in) :: x(ixi^s, 1:ndim)
634 double precision :: invgam
635 integer :: idir, itr
636
637 !!if (fix_small_values) then
638 !! call hd_handle_small_values(.true., w, x, ixI^L, ixO^L, 'hd_to_conserved')
639 !!end if
640
641 if (hd_energy) then
642 invgam = 1.d0/(hd_gamma - 1.0d0)
643 ! Calculate total energy from pressure and kinetic energy
644 w(ixo^s, e_) = w(ixo^s, e_) * invgam + &
645 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * w(ixo^s, rho_)
646 end if
647
648 ! Convert velocity to momentum
649 do idir = 1, ndir
650 w(ixo^s, mom(idir)) = w(ixo^s, rho_) * w(ixo^s, mom(idir))
651 end do
652
653 if (hd_dust) then
654 call dust_to_conserved(ixi^l, ixo^l, w, x)
655 end if
656
657 end subroutine hd_to_conserved
658
659 !> Transform conservative variables into primitive ones
660 subroutine hd_to_primitive(ixI^L, ixO^L, w, x)
662 use mod_dust, only: dust_to_primitive
663 integer, intent(in) :: ixi^l, ixo^l
664 double precision, intent(inout) :: w(ixi^s, nw)
665 double precision, intent(in) :: x(ixi^s, 1:ndim)
666 integer :: itr, idir
667 double precision :: inv_rho(ixo^s)
668
669 if (fix_small_values) then
670 call hd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'hd_to_primitive')
671 end if
672
673 inv_rho = 1.0d0 / w(ixo^s, rho_)
674
675 if (hd_energy) then
676 ! Compute pressure
677 w(ixo^s, e_) = (hd_gamma - 1.0d0) * (w(ixo^s, e_) - &
678 hd_kin_en(w, ixi^l, ixo^l, inv_rho))
679 end if
680
681 ! Convert momentum to velocity
682 do idir = 1, ndir
683 w(ixo^s, mom(idir)) = w(ixo^s, mom(idir)) * inv_rho
684 end do
685
686 ! Convert dust momentum to dust velocity
687 if (hd_dust) then
688 call dust_to_primitive(ixi^l, ixo^l, w, x)
689 end if
690
691 end subroutine hd_to_primitive
692
693 !> Transform internal energy to total energy
694 subroutine hd_ei_to_e(ixI^L,ixO^L,w,x)
696 integer, intent(in) :: ixi^l, ixo^l
697 double precision, intent(inout) :: w(ixi^s, nw)
698 double precision, intent(in) :: x(ixi^s, 1:ndim)
699
700 ! Calculate total energy from internal and kinetic energy
701 w(ixo^s,e_)=w(ixo^s,e_)&
702 +hd_kin_en(w,ixi^l,ixo^l)
703
704 end subroutine hd_ei_to_e
705
706 !> Transform total energy to internal energy
707 subroutine hd_e_to_ei(ixI^L,ixO^L,w,x)
709 integer, intent(in) :: ixi^l, ixo^l
710 double precision, intent(inout) :: w(ixi^s, nw)
711 double precision, intent(in) :: x(ixi^s, 1:ndim)
712
713 ! Calculate ei = e - ek
714 w(ixo^s,e_)=w(ixo^s,e_)&
715 -hd_kin_en(w,ixi^l,ixo^l)
716
717 end subroutine hd_e_to_ei
718
719 subroutine e_to_rhos(ixI^L, ixO^L, w, x)
721
722 integer, intent(in) :: ixi^l, ixo^l
723 double precision :: w(ixi^s, nw)
724 double precision, intent(in) :: x(ixi^s, 1:ndim)
725
726 if (hd_energy) then
727 w(ixo^s, e_) = (hd_gamma - 1.0d0) * w(ixo^s, rho_)**(1.0d0 - hd_gamma) * &
728 (w(ixo^s, e_) - hd_kin_en(w, ixi^l, ixo^l))
729 else
730 call mpistop("energy from entropy can not be used with -eos = iso !")
731 end if
732 end subroutine e_to_rhos
733
734 subroutine rhos_to_e(ixI^L, ixO^L, w, x)
736
737 integer, intent(in) :: ixi^l, ixo^l
738 double precision :: w(ixi^s, nw)
739 double precision, intent(in) :: x(ixi^s, 1:ndim)
740
741 if (hd_energy) then
742 w(ixo^s, e_) = w(ixo^s, rho_)**(hd_gamma - 1.0d0) * w(ixo^s, e_) &
743 / (hd_gamma - 1.0d0) + hd_kin_en(w, ixi^l, ixo^l)
744 else
745 call mpistop("entropy from energy can not be used with -eos = iso !")
746 end if
747 end subroutine rhos_to_e
748
749 !> Calculate v_i = m_i / rho within ixO^L
750 subroutine hd_get_v_idim(w, x, ixI^L, ixO^L, idim, v)
752 integer, intent(in) :: ixi^l, ixo^l, idim
753 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
754 double precision, intent(out) :: v(ixi^s)
755
756 v(ixo^s) = w(ixo^s, mom(idim)) / w(ixo^s, rho_)
757 end subroutine hd_get_v_idim
758
759 !> Calculate velocity vector v_i = m_i / rho within ixO^L
760 subroutine hd_get_v(w,x,ixI^L,ixO^L,v)
762
763 integer, intent(in) :: ixi^l, ixo^l
764 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:^nd)
765 double precision, intent(out) :: v(ixi^s,1:ndir)
766
767 integer :: idir
768
769 do idir=1,ndir
770 v(ixo^s,idir) = w(ixo^s, mom(idir)) / w(ixo^s, rho_)
771 end do
772
773 end subroutine hd_get_v
774
775 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
776 subroutine hd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
780
781 integer, intent(in) :: ixi^l, ixo^l, idim
782 ! w in primitive form
783 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
784 double precision, intent(inout) :: cmax(ixi^s)
785
786 if(hd_energy) then
787 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(hd_gamma*w(ixo^s,p_)/w(ixo^s,rho_))
788 else
789 if (.not. associated(usr_set_pthermal)) then
790 cmax(ixo^s) = hd_adiab * w(ixo^s, rho_)**hd_gamma
791 else
792 call usr_set_pthermal(w,x,ixi^l,ixo^l,cmax)
793 end if
794 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(hd_gamma*cmax(ixo^s)/w(ixo^s,rho_))
795 end if
796
797 if (hd_dust) then
798 call dust_get_cmax_prim(w, x, ixi^l, ixo^l, idim, cmax)
799 end if
800 end subroutine hd_get_cmax
801
802 subroutine hd_get_a2max(w,x,ixI^L,ixO^L,a2max)
804
805 integer, intent(in) :: ixi^l, ixo^l
806 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
807 double precision, intent(inout) :: a2max(ndim)
808 double precision :: a2(ixi^s,ndim,nw)
809 integer :: gxo^l,hxo^l,jxo^l,kxo^l,i,j
810
811 a2=zero
812 do i = 1,ndim
813 !> 4th order
814 hxo^l=ixo^l-kr(i,^d);
815 gxo^l=hxo^l-kr(i,^d);
816 jxo^l=ixo^l+kr(i,^d);
817 kxo^l=jxo^l+kr(i,^d);
818 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
819 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
820 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/dxlevel(i)**2
821 end do
822 end subroutine hd_get_a2max
823
824 !> get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
825 subroutine hd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
827 integer, intent(in) :: ixi^l,ixo^l
828 double precision, intent(in) :: x(ixi^s,1:ndim)
829 ! in primitive form
830 double precision, intent(inout) :: w(ixi^s,1:nw)
831 double precision, intent(out) :: tco_local, tmax_local
832
833 double precision, parameter :: trac_delta=0.25d0
834 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
835 double precision :: ltr(ixi^s),ltrc,ltrp,tcoff(ixi^s)
836 integer :: jxo^l,hxo^l
837 integer :: jxp^l,hxp^l,ixp^l
838 logical :: lrlt(ixi^s)
839
840 {^ifoned
841 call hd_get_rfactor(w,x,ixi^l,ixi^l,te)
842 te(ixi^s)=w(ixi^s,p_)/(te(ixi^s)*w(ixi^s,rho_))
843
844 tco_local=zero
845 tmax_local=maxval(te(ixo^s))
846 select case(hd_trac_type)
847 case(0)
848 w(ixi^s,tcoff_)=3.d5/unit_temperature
849 case(1)
850 hxo^l=ixo^l-1;
851 jxo^l=ixo^l+1;
852 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
853 lrlt=.false.
854 where(lts(ixo^s) > trac_delta)
855 lrlt(ixo^s)=.true.
856 end where
857 if(any(lrlt(ixo^s))) then
858 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
859 end if
860 case(2)
861 !> iijima et al. 2021, LTRAC method
862 ltrc=1.5d0
863 ltrp=2.5d0
864 ixp^l=ixo^l^ladd1;
865 hxo^l=ixo^l-1;
866 jxo^l=ixo^l+1;
867 hxp^l=ixp^l-1;
868 jxp^l=ixp^l+1;
869 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
870 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
871 w(ixo^s,tcoff_)=te(ixo^s)*&
872 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
873 case default
874 call mpistop("mhd_trac_type not allowed for 1D simulation")
875 end select
876 }
877 end subroutine hd_get_tcutoff
878
879 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
880 subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
882 use mod_dust, only: dust_get_cmax
883 use mod_variables
884
885 integer, intent(in) :: ixi^l, ixo^l, idim
886 ! conservative left and right status
887 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
888 ! primitive left and right status
889 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
890 double precision, intent(in) :: x(ixi^s, 1:ndim)
891 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
892 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
893 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
894
895 double precision :: wmean(ixi^s,nw)
896 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
897 integer :: ix^d
898
899 select case(boundspeed)
900 case (1)
901 ! This implements formula (10.52) from "Riemann Solvers and Numerical
902 ! Methods for Fluid Dynamics" by Toro.
903
904 tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
905 tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
906 tmp3(ixo^s)=1.d0/(dsqrt(wlp(ixo^s,rho_))+dsqrt(wrp(ixo^s,rho_)))
907 umean(ixo^s)=(wlp(ixo^s,mom(idim))*tmp1(ixo^s)+wrp(ixo^s,mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
908
909 if(hd_energy) then
910 csoundl(ixo^s)=hd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
911 csoundr(ixo^s)=hd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
912 else
913 call hd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
914 call hd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
915 end if
916
917 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
918 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
919 (wrp(ixo^s,mom(idim))-wlp(ixo^s,mom(idim)))**2
920
921 dmean(ixo^s)=dsqrt(dmean(ixo^s))
922 if(present(cmin)) then
923 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
924 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
925 if(h_correction) then
926 {do ix^db=ixomin^db,ixomax^db\}
927 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
928 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
929 {end do\}
930 end if
931 else
932 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
933 end if
934
935 if (hd_dust) then
936 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
937 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
938 end if
939
940 case (2)
941 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
942 tmp1(ixo^s)=wmean(ixo^s,mom(idim))/wmean(ixo^s,rho_)
943 call hd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
944 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
945
946 if(present(cmin)) then
947 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
948 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
949 if(h_correction) then
950 {do ix^db=ixomin^db,ixomax^db\}
951 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
952 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
953 {end do\}
954 end if
955 else
956 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
957 end if
958
959 if (hd_dust) then
960 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
961 end if
962 case (3)
963 ! Miyoshi 2005 JCP 208, 315 equation (67)
964 if(hd_energy) then
965 csoundl(ixo^s)=hd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
966 csoundr(ixo^s)=hd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
967 else
968 call hd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
969 call hd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
970 end if
971 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
972 if(present(cmin)) then
973 cmin(ixo^s,1)=min(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))-csoundl(ixo^s)
974 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
975 if(h_correction) then
976 {do ix^db=ixomin^db,ixomax^db\}
977 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
978 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
979 {end do\}
980 end if
981 else
982 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
983 end if
984 if (hd_dust) then
985 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
986 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
987 end if
988 end select
989
990 end subroutine hd_get_cbounds
991
992 !> Calculate the square of the thermal sound speed csound2 within ixO^L.
993 !> csound2=gamma*p/rho
994 subroutine hd_get_csound2(w,x,ixI^L,ixO^L,csound2)
996 integer, intent(in) :: ixi^l, ixo^l
997 double precision, intent(in) :: w(ixi^s,nw)
998 double precision, intent(in) :: x(ixi^s,1:ndim)
999 double precision, intent(out) :: csound2(ixi^s)
1000
1001 call hd_get_pthermal(w,x,ixi^l,ixo^l,csound2)
1002 csound2(ixo^s)=hd_gamma*csound2(ixo^s)/w(ixo^s,rho_)
1003
1004 end subroutine hd_get_csound2
1005
1006 !> Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L
1007 subroutine hd_get_pthermal(w, x, ixI^L, ixO^L, pth)
1011
1012 integer, intent(in) :: ixi^l, ixo^l
1013 double precision, intent(in) :: w(ixi^s, 1:nw)
1014 double precision, intent(in) :: x(ixi^s, 1:ndim)
1015 double precision, intent(out):: pth(ixi^s)
1016 integer :: iw, ix^d
1017
1018 if (hd_energy) then
1019 pth(ixo^s) = (hd_gamma - 1.0d0) * (w(ixo^s, e_) - &
1020 hd_kin_en(w, ixi^l, ixo^l))
1021 else
1022 if (.not. associated(usr_set_pthermal)) then
1023 pth(ixo^s) = hd_adiab * w(ixo^s, rho_)**hd_gamma
1024 else
1025 call usr_set_pthermal(w,x,ixi^l,ixo^l,pth)
1026 end if
1027 end if
1028
1029 if (fix_small_values) then
1030 {do ix^db= ixo^lim^db\}
1031 if(pth(ix^d)<small_pressure) then
1032 pth(ix^d)=small_pressure
1033 endif
1034 {enddo^d&\}
1035 else if (check_small_values) then
1036 {do ix^db= ixo^lim^db\}
1037 if(pth(ix^d)<small_pressure) then
1038 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1039 " encountered when call hd_get_pthermal"
1040 write(*,*) "Iteration: ", it, " Time: ", global_time
1041 write(*,*) "Location: ", x(ix^d,:)
1042 write(*,*) "Cell number: ", ix^d
1043 do iw=1,nw
1044 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1045 end do
1046 ! use erroneous arithmetic operation to crash the run
1047 if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1048 write(*,*) "Saving status at the previous time step"
1049 crash=.true.
1050 end if
1051 {enddo^d&\}
1052 end if
1053
1054 end subroutine hd_get_pthermal
1055
1056 !> Calculate temperature=p/rho when in e_ the total energy is stored
1057 subroutine hd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1059 integer, intent(in) :: ixi^l, ixo^l
1060 double precision, intent(in) :: w(ixi^s, 1:nw)
1061 double precision, intent(in) :: x(ixi^s, 1:ndim)
1062 double precision, intent(out):: res(ixi^s)
1063
1064 double precision :: r(ixi^s)
1065
1066 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1067 call hd_get_pthermal(w, x, ixi^l, ixo^l, res)
1068 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1069 end subroutine hd_get_temperature_from_etot
1070
1071 !> Calculate temperature=p/rho when in e_ the internal energy is stored
1072 subroutine hd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1074 integer, intent(in) :: ixi^l, ixo^l
1075 double precision, intent(in) :: w(ixi^s, 1:nw)
1076 double precision, intent(in) :: x(ixi^s, 1:ndim)
1077 double precision, intent(out):: res(ixi^s)
1078
1079 double precision :: r(ixi^s)
1080
1081 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1082 res(ixo^s) = (hd_gamma - 1.0d0) * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1083 end subroutine hd_get_temperature_from_eint
1084
1085 ! Calculate flux f_idim[iw]
1086 subroutine hd_get_flux_cons(w, x, ixI^L, ixO^L, idim, f)
1088 use mod_dust, only: dust_get_flux
1089
1090 integer, intent(in) :: ixi^l, ixo^l, idim
1091 double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
1092 double precision, intent(out) :: f(ixi^s, nwflux)
1093 double precision :: pth(ixi^s), v(ixi^s),frame_vel(ixi^s)
1094 integer :: idir, itr
1095
1096 call hd_get_pthermal(w, x, ixi^l, ixo^l, pth)
1097 call hd_get_v_idim(w, x, ixi^l, ixo^l, idim, v)
1098
1099 f(ixo^s, rho_) = v(ixo^s) * w(ixo^s, rho_)
1100
1101 ! Momentum flux is v_i*m_i, +p in direction idim
1102 do idir = 1, ndir
1103 f(ixo^s, mom(idir)) = v(ixo^s) * w(ixo^s, mom(idir))
1104 end do
1105
1106 f(ixo^s, mom(idim)) = f(ixo^s, mom(idim)) + pth(ixo^s)
1107
1108 if(hd_energy) then
1109 ! Energy flux is v_i*(e + p)
1110 f(ixo^s, e_) = v(ixo^s) * (w(ixo^s, e_) + pth(ixo^s))
1111 end if
1112
1113 do itr = 1, hd_n_tracer
1114 f(ixo^s, tracer(itr)) = v(ixo^s) * w(ixo^s, tracer(itr))
1115 end do
1116
1117 ! Dust fluxes
1118 if (hd_dust) then
1119 call dust_get_flux(w, x, ixi^l, ixo^l, idim, f)
1120 end if
1121
1122 end subroutine hd_get_flux_cons
1123
1124 ! Calculate flux f_idim[iw]
1125 subroutine hd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1127 use mod_dust, only: dust_get_flux_prim
1128 use mod_viscosity, only: visc_get_flux_prim ! viscInDiv
1129
1130 integer, intent(in) :: ixi^l, ixo^l, idim
1131 ! conservative w
1132 double precision, intent(in) :: wc(ixi^s, 1:nw)
1133 ! primitive w
1134 double precision, intent(in) :: w(ixi^s, 1:nw)
1135 double precision, intent(in) :: x(ixi^s, 1:ndim)
1136 double precision, intent(out) :: f(ixi^s, nwflux)
1137 double precision :: pth(ixi^s),frame_vel(ixi^s)
1138 integer :: idir, itr
1139
1140 if (hd_energy) then
1141 pth(ixo^s) = w(ixo^s,p_)
1142 else
1143 call hd_get_pthermal(wc, x, ixi^l, ixo^l, pth)
1144 end if
1145
1146 f(ixo^s, rho_) = w(ixo^s,mom(idim)) * w(ixo^s, rho_)
1147
1148 ! Momentum flux is v_i*m_i, +p in direction idim
1149 do idir = 1, ndir
1150 f(ixo^s, mom(idir)) = w(ixo^s,mom(idim)) * wc(ixo^s, mom(idir))
1151 end do
1152
1153 f(ixo^s, mom(idim)) = f(ixo^s, mom(idim)) + pth(ixo^s)
1154
1155 if(hd_energy) then
1156 ! Energy flux is v_i*(e + p)
1157 f(ixo^s, e_) = w(ixo^s,mom(idim)) * (wc(ixo^s, e_) + w(ixo^s,p_))
1158 end if
1159
1160 do itr = 1, hd_n_tracer
1161 f(ixo^s, tracer(itr)) = w(ixo^s,mom(idim)) * w(ixo^s, tracer(itr))
1162 end do
1163
1164 ! Dust fluxes
1165 if (hd_dust) then
1166 call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1167 end if
1168
1169 ! Viscosity fluxes - viscInDiv
1170 if (hd_viscosity) then
1171 call visc_get_flux_prim(w, x, ixi^l, ixo^l, idim, f, hd_energy)
1172 endif
1173
1174 end subroutine hd_get_flux
1175
1176 !> Add geometrical source terms to w
1177 !>
1178 !> Notice that the expressions of the geometrical terms depend only on ndir,
1179 !> not ndim. Eg, they are the same in 2.5D and in 3D, for any geometry.
1180 !>
1181 !> Ileyk : to do :
1182 !> - address the source term for the dust in case (coordinate == spherical)
1183 subroutine hd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
1186 use mod_viscosity, only: visc_add_source_geom ! viscInDiv
1189 use mod_geometry
1190 integer, intent(in) :: ixi^l, ixo^l
1191 double precision, intent(in) :: qdt, dtfactor, x(ixi^s, 1:ndim)
1192 double precision, intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
1193 ! to change and to set as a parameter in the parfile once the possibility to
1194 ! solve the equations in an angular momentum conserving form has been
1195 ! implemented (change tvdlf.t eg)
1196 double precision :: pth(ixi^s), source(ixi^s), minrho
1197 integer :: iw,idir, h1x^l{^nooned, h2x^l}
1198 integer :: mr_,mphi_ ! Polar var. names
1199 integer :: irho, ifluid, n_fluids
1200 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
1201
1202 if (hd_dust) then
1203 n_fluids = 1 + dust_n_species
1204 else
1205 n_fluids = 1
1206 end if
1207
1208 select case (coordinate)
1209
1211 !the user provides the functions of exp_factor and del_exp_factor
1212 if(associated(usr_set_surface)) call usr_set_surface(ixi^l,x,block%dx,exp_factor,del_exp_factor,exp_factor_primitive)
1213 if(hd_energy) then
1214 source(ixo^s)=wprim(ixo^s, p_)
1215 else
1216 if(.not. associated(usr_set_pthermal)) then
1217 source(ixo^s)=hd_adiab * wprim(ixo^s, rho_)**hd_gamma
1218 else
1219 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1220 end if
1221 end if
1222 source(ixo^s) = source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1223 w(ixo^s,mom(1)) = w(ixo^s,mom(1)) + qdt*source(ixo^s)
1224
1225 case (cylindrical)
1226 do ifluid = 0, n_fluids-1
1227 ! s[mr]=(pthermal+mphi**2/rho)/radius
1228 if (ifluid == 0) then
1229 ! gas
1230 irho = rho_
1231 mr_ = mom(r_)
1232 if(phi_>0) mphi_ = mom(phi_)
1233 if(hd_energy) then
1234 source(ixo^s)=wprim(ixo^s, p_)
1235 else
1236 if(.not. associated(usr_set_pthermal)) then
1237 source(ixo^s)=hd_adiab * wprim(ixo^s, rho_)**hd_gamma
1238 else
1239 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1240 end if
1241 end if
1242 minrho = 0.0d0
1243 else
1244 ! dust : no pressure
1245 irho = dust_rho(ifluid)
1246 mr_ = dust_mom(r_, ifluid)
1247 if(phi_>0) mphi_ = dust_mom(phi_, ifluid)
1248 source(ixi^s) = zero
1249 minrho = 0.0d0
1250 end if
1251 if(phi_ > 0) then
1252 where (wct(ixo^s, irho) > minrho)
1253 source(ixo^s) = source(ixo^s) + wct(ixo^s,mphi_)*wprim(ixo^s,mphi_)
1254 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt*source(ixo^s)/x(ixo^s,r_)
1255 end where
1256 ! s[mphi]=(-mphi*vr)/radius
1257 where (wct(ixo^s, irho) > minrho)
1258 source(ixo^s) = -wct(ixo^s, mphi_) * wprim(ixo^s, mr_)
1259 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt * source(ixo^s) / x(ixo^s, r_)
1260 end where
1261 else
1262 ! s[mr]=2pthermal/radius
1263 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
1264 end if
1265 end do
1266 case (spherical)
1267 if (hd_dust) then
1268 call mpistop("Dust geom source terms not implemented yet with spherical geometries")
1269 end if
1270 mr_ = mom(r_)
1271 if(phi_>0) mphi_ = mom(phi_)
1272 h1x^l=ixo^l-kr(1,^d); {^nooned h2x^l=ixo^l-kr(2,^d);}
1273 if(hd_energy) then
1274 pth(ixo^s)=wprim(ixo^s, p_)
1275 else
1276 if(.not. associated(usr_set_pthermal)) then
1277 pth(ixo^s)=hd_adiab * wprim(ixo^s, rho_)**hd_gamma
1278 else
1279 call usr_set_pthermal(wct,x,ixi^l,ixo^l,pth)
1280 end if
1281 end if
1282 ! s[mr]=((vtheta**2+vphi**2)*rho+2*p)/r
1283 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1284 *(block%surfaceC(ixo^s, 1) - block%surfaceC(h1x^s, 1)) &
1285 /block%dvolume(ixo^s)
1286 do idir = 2, ndir
1287 source(ixo^s) = source(ixo^s) + wprim(ixo^s, mom(idir))**2 * wprim(ixo^s, rho_)
1288 end do
1289 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, 1)
1290
1291 {^nooned
1292 ! s[mtheta]=-(vr*vtheta*rho)/r+cot(theta)*(vphi**2*rho+p)/r
1293 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1294 * (block%surfaceC(ixo^s, 2) - block%surfaceC(h2x^s, 2)) &
1295 / block%dvolume(ixo^s)
1296 if (ndir == 3) then
1297 source(ixo^s) = source(ixo^s) + (wprim(ixo^s, mom(3))**2 * wprim(ixo^s, rho_)) / tan(x(ixo^s, 2))
1298 end if
1299 source(ixo^s) = source(ixo^s) - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)
1300 w(ixo^s, mom(2)) = w(ixo^s, mom(2)) + qdt * source(ixo^s) / x(ixo^s, 1)
1301
1302 if (ndir == 3) then
1303 ! s[mphi]=-(vphi*vr/rho)/r-cot(theta)*(vtheta*vphi/rho)/r
1304 source(ixo^s) = -(wprim(ixo^s, mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)&
1305 - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mom(3))) * wprim(ixo^s, rho_) / tan(x(ixo^s, 2))
1306 w(ixo^s, mom(3)) = w(ixo^s, mom(3)) + qdt * source(ixo^s) / x(ixo^s, 1)
1307 end if
1308 }
1309 end select
1310
1311 if (hd_viscosity) call visc_add_source_geom(qdt,ixi^l,ixo^l,wprim,w,x)
1312
1313 if (hd_rotating_frame) then
1314 if (hd_dust) then
1315 call mpistop("Rotating frame not implemented yet with dust")
1316 else
1317 call rotating_frame_add_source(qdt,dtfactor,ixi^l,ixo^l,wprim,w,x)
1318 end if
1319 end if
1320
1321 end subroutine hd_add_source_geom
1322
1323 ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
1324 subroutine hd_add_source(qdt,dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1329 use mod_usr_methods, only: usr_gravity
1331 use mod_cak_force, only: cak_add_source
1332
1333 integer, intent(in) :: ixi^l, ixo^l
1334 double precision, intent(in) :: qdt, dtfactor
1335 double precision, intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:ndim)
1336 double precision, intent(inout) :: w(ixi^s, 1:nw)
1337 logical, intent(in) :: qsourcesplit
1338 logical, intent(inout) :: active
1339
1340 double precision :: gravity_field(ixi^s, 1:ndim)
1341 integer :: idust, idim
1342
1343 if(hd_dust .and. .not. use_imex_scheme) then
1344 call dust_add_source(qdt,ixi^l,ixo^l,wct,w,x,qsourcesplit,active)
1345 end if
1346
1347 if(hd_radiative_cooling) then
1348 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1349 qsourcesplit,active, rc_fl)
1350 end if
1351
1352 if(hd_viscosity) then
1353 call viscosity_add_source(qdt,ixi^l,ixo^l,wct,w,x,&
1354 hd_energy,qsourcesplit,active)
1355 end if
1356
1357 if (hd_gravity) then
1358 call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1359 hd_energy,.false.,qsourcesplit,active)
1360
1361 if (hd_dust .and. qsourcesplit .eqv. grav_split) then
1362 active = .true.
1363
1364 call usr_gravity(ixi^l, ixo^l, wct, x, gravity_field)
1365 do idust = 1, dust_n_species
1366 do idim = 1, ndim
1367 w(ixo^s, dust_mom(idim, idust)) = w(ixo^s, dust_mom(idim, idust)) &
1368 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s, dust_rho(idust))
1369 end do
1370 end do
1371 end if
1372 end if
1373
1374 if (hd_cak_force) then
1375 call cak_add_source(qdt,ixi^l,ixo^l,wct,w,x,hd_energy,qsourcesplit,active)
1376 end if
1377
1378 if(hd_partial_ionization) then
1379 if(.not.qsourcesplit) then
1380 active = .true.
1381 call hd_update_temperature(ixi^l,ixo^l,wct,w,x)
1382 end if
1383 end if
1384
1385 end subroutine hd_add_source
1386
1387 subroutine hd_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
1389 use mod_dust, only: dust_get_dt
1392 use mod_gravity, only: gravity_get_dt
1393 use mod_cak_force, only: cak_get_dt
1394
1395 integer, intent(in) :: ixi^l, ixo^l
1396 double precision, intent(in) :: dx^d, x(ixi^s, 1:^nd)
1397 double precision, intent(in) :: w(ixi^s, 1:nw)
1398 double precision, intent(inout) :: dtnew
1399
1400 dtnew = bigdouble
1401
1402 if(hd_dust) then
1403 call dust_get_dt(w, ixi^l, ixo^l, dtnew, dx^d, x)
1404 end if
1405
1406 if(hd_radiative_cooling) then
1407 call cooling_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x,rc_fl)
1408 end if
1409
1410 if(hd_viscosity) then
1411 call viscosity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1412 end if
1413
1414 if(hd_gravity) then
1415 call gravity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1416 end if
1417
1418 if (hd_cak_force) then
1419 call cak_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1420 end if
1421
1422 end subroutine hd_get_dt
1423
1424 function hd_kin_en(w, ixI^L, ixO^L, inv_rho) result(ke)
1425 use mod_global_parameters, only: nw, ndim
1426 integer, intent(in) :: ixi^l, ixo^l
1427 double precision, intent(in) :: w(ixi^s, nw)
1428 double precision :: ke(ixo^s)
1429 double precision, intent(in), optional :: inv_rho(ixo^s)
1430
1431 if (present(inv_rho)) then
1432 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * inv_rho
1433 else
1434 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) / w(ixo^s, rho_)
1435 end if
1436 end function hd_kin_en
1437
1438 function hd_inv_rho(w, ixI^L, ixO^L) result(inv_rho)
1439 use mod_global_parameters, only: nw, ndim
1440 integer, intent(in) :: ixi^l, ixo^l
1441 double precision, intent(in) :: w(ixi^s, nw)
1442 double precision :: inv_rho(ixo^s)
1443
1444 ! Can make this more robust
1445 inv_rho = 1.0d0 / w(ixo^s, rho_)
1446 end function hd_inv_rho
1447
1448 subroutine hd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1449 ! handles hydro (density,pressure,velocity) bootstrapping
1450 ! any negative dust density is flagged as well (and throws an error)
1451 ! small_values_method=replace also for dust
1455 logical, intent(in) :: primitive
1456 integer, intent(in) :: ixi^l,ixo^l
1457 double precision, intent(inout) :: w(ixi^s,1:nw)
1458 double precision, intent(in) :: x(ixi^s,1:ndim)
1459 character(len=*), intent(in) :: subname
1460
1461 integer :: n,idir
1462 logical :: flag(ixi^s,1:nw)
1463
1464 call hd_check_w(primitive, ixi^l, ixo^l, w, flag)
1465
1466 if (any(flag)) then
1467 select case (small_values_method)
1468 case ("replace")
1469 where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
1470 do idir = 1, ndir
1471 if(small_values_fix_iw(mom(idir))) then
1472 where(flag(ixo^s,rho_)) w(ixo^s, mom(idir)) = 0.0d0
1473 end if
1474 end do
1475 if(hd_energy)then
1476 if(small_values_fix_iw(e_)) then
1477 if(primitive) then
1478 where(flag(ixo^s,rho_)) w(ixo^s, p_) = small_pressure
1479 else
1480 where(flag(ixo^s,rho_)) w(ixo^s, e_) = small_e + hd_kin_en(w,ixi^l,ixo^l)
1481 endif
1482 end if
1483 endif
1484
1485 if(hd_energy) then
1486 if(primitive) then
1487 where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
1488 else
1489 where(flag(ixo^s,e_))
1490 ! Add kinetic energy
1491 w(ixo^s,e_) = small_e + hd_kin_en(w,ixi^l,ixo^l)
1492 end where
1493 end if
1494 end if
1495
1496 if(hd_dust)then
1497 do n=1,dust_n_species
1498 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1499 do idir = 1, ndir
1500 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1501 enddo
1502 enddo
1503 endif
1504 case ("average")
1505 if(primitive)then
1506 ! averaging for all primitive fields, including dust
1507 call small_values_average(ixi^l, ixo^l, w, x, flag)
1508 else
1509 ! do averaging of density
1510 call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
1511 if(hd_energy) then
1512 ! do averaging of pressure
1513 w(ixi^s,p_)=(hd_gamma-1.d0)*(w(ixi^s,e_) &
1514 -0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_))
1515 call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
1516 w(ixi^s,e_)=w(ixi^s,p_)/(hd_gamma-1.d0) &
1517 +0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_)
1518 end if
1519 if(hd_dust)then
1520 do n=1,dust_n_species
1521 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1522 do idir = 1, ndir
1523 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1524 enddo
1525 enddo
1526 endif
1527 endif
1528 case default
1529 if(.not.primitive) then
1530 !convert w to primitive
1531 ! Calculate pressure = (gamma-1) * (e-ek)
1532 if(hd_energy) then
1533 w(ixo^s,p_)=(hd_gamma-1.d0)*(w(ixo^s,e_)-hd_kin_en(w,ixi^l,ixo^l))
1534 end if
1535 ! Convert gas momentum to velocity
1536 do idir = 1, ndir
1537 w(ixo^s, mom(idir)) = w(ixo^s, mom(idir))/w(ixo^s,rho_)
1538 end do
1539 end if
1540 ! NOTE: dust entries may still have conserved values here
1541 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1542 end select
1543 end if
1544 end subroutine hd_handle_small_values
1545
1546 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1549 integer, intent(in) :: ixi^l, ixo^l
1550 double precision, intent(in) :: w(ixi^s,1:nw)
1551 double precision, intent(in) :: x(ixi^s,1:ndim)
1552 double precision, intent(out):: rfactor(ixi^s)
1553
1554 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1555
1556 call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1557 ! assume the first and second ionization of Helium have the same degree
1558 rfactor(ixo^s)=(1.d0+iz_h(ixo^s)+0.1d0*(1.d0+iz_he(ixo^s)*(1.d0+iz_he(ixo^s))))/2.3d0
1559
1560 end subroutine rfactor_from_temperature_ionization
1561
1562 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1564 integer, intent(in) :: ixi^l, ixo^l
1565 double precision, intent(in) :: w(ixi^s,1:nw)
1566 double precision, intent(in) :: x(ixi^s,1:ndim)
1567 double precision, intent(out):: rfactor(ixi^s)
1568
1569 rfactor(ixo^s)=rr
1570
1571 end subroutine rfactor_from_constant_ionization
1572
1573 subroutine hd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1576
1577 integer, intent(in) :: ixi^l, ixo^l
1578 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1579 double precision, intent(inout) :: w(ixi^s,1:nw)
1580
1581 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1582
1583 call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1584
1585 call hd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1586
1587 w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1588 he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1589
1590 end subroutine hd_update_temperature
1591
1592end module mod_hd_phys
Calculate w(iw)=w(iw)+qdt*SOURCE[wCT,qtC,x] within ixO for all indices iw=iwmin......
Module with basic data types used in amrvac.
integer, parameter std_len
Default length for strings.
Module to include CAK radiation line force in (magneto)hydrodynamic models Computes both the force fr...
subroutine cak_get_dt(w, ixil, ixol, dtnew, dxd, x)
Check time step for total radiation contribution.
subroutine cak_init(phys_gamma)
Initialize the module.
subroutine cak_add_source(qdt, ixil, ixol, wct, w, x, energy, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter bigdouble
A very large real number.
Module for including dust species, which interact with the gas through a drag force.
Definition mod_dust.t:3
subroutine, public dust_add_source(qdt, ixil, ixol, wct, w, x, qsourcesplit, active)
w[iw]= w[iw]+qdt*S[wCT, x] where S is the source based on wCT within ixO
Definition mod_dust.t:530
subroutine, public dust_evaluate_implicit(qtc, psa)
inplace update of psa==>F_im(psa)
Definition mod_dust.t:583
subroutine, public dust_to_primitive(ixil, ixol, w, x)
Definition mod_dust.t:228
subroutine, public dust_get_dt(w, ixil, ixol, dtnew, dxd, x)
Get dt related to dust and gas stopping time (Laibe 2011)
Definition mod_dust.t:898
subroutine, public dust_get_flux(w, x, ixil, ixol, idim, f)
Definition mod_dust.t:252
integer, dimension(:, :), allocatable, public, protected dust_mom
Indices of the dust momentum densities.
Definition mod_dust.t:47
subroutine, public dust_to_conserved(ixil, ixol, w, x)
Definition mod_dust.t:208
integer, public, protected dust_n_species
The number of dust species.
Definition mod_dust.t:37
subroutine, public dust_get_flux_prim(w, x, ixil, ixol, idim, f)
Definition mod_dust.t:275
subroutine, public dust_check_w(ixil, ixol, w, flag)
Definition mod_dust.t:194
integer, dimension(:), allocatable, public, protected dust_rho
Indices of the dust densities.
Definition mod_dust.t:44
subroutine, public dust_get_cmax(w, x, ixil, ixol, idim, cmax, cmin)
Definition mod_dust.t:1011
subroutine, public dust_check_params()
Definition mod_dust.t:154
subroutine, public dust_get_cmax_prim(w, x, ixil, ixol, idim, cmax, cmin)
Definition mod_dust.t:1035
subroutine, public dust_init(g_rho, g_mom, g_energy)
Definition mod_dust.t:95
subroutine, public dust_implicit_update(dtfactor, qdt, qtc, psb, psa)
Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
Definition mod_dust.t:652
Module for flux conservation near refinement boundaries.
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
integer coordinate
Definition mod_geometry.t:7
integer, parameter spherical
integer, parameter cylindrical
integer, parameter cartesian_expansion
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
logical h_correction
If true, do H-correction to fix the carbuncle problem at grid-aligned shocks.
double precision unit_time
Physical scaling factor for time.
double precision unit_density
Physical scaling factor for density.
integer, parameter unitpar
file handle for IO
double precision global_time
The global simulation time.
double precision unit_mass
Physical scaling factor for mass.
logical use_imex_scheme
whether IMEX in use or not
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer it
Number of time steps taken.
double precision unit_numberdensity
Physical scaling factor for number density.
character(len=std_len) convert_type
Which format to use when converting.
double precision unit_pressure
Physical scaling factor for pressure.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
logical use_particles
Use particles module or not.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision unit_velocity
Physical scaling factor for velocity.
double precision unit_temperature
Physical scaling factor for temperature.
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
double precision, dimension(:,:), allocatable dx
logical phys_trac
Use TRAC for MHD or 1D HD.
logical fix_small_values
fix small values with average or replace methods
logical crash
Save a snapshot before crash a run met unphysical values.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
integer, parameter unitconvert
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
Module for including gravity in (magneto)hydrodynamics simulations.
Definition mod_gravity.t:2
logical grav_split
source split or not
Definition mod_gravity.t:6
subroutine gravity_get_dt(w, ixil, ixol, dtnew, dxd, x)
Definition mod_gravity.t:87
subroutine gravity_init()
Initialize the module.
Definition mod_gravity.t:26
subroutine gravity_add_source(qdt, ixil, ixol, wct, wctprim, w, x, energy, rhov, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
Definition mod_gravity.t:43
Hydrodynamics physics module.
Definition mod_hd_phys.t:2
subroutine, public hd_check_params
logical, public, protected hd_energy
Whether an energy equation is used.
Definition mod_hd_phys.t:12
logical, public, protected hd_dust
Whether dust is added.
Definition mod_hd_phys.t:24
integer, public, protected e_
Index of the energy density (-1 if not present)
Definition mod_hd_phys.t:57
logical, public, protected hd_radiative_cooling
Whether radiative cooling is added.
Definition mod_hd_phys.t:20
double precision, public, protected rr
Definition mod_hd_phys.t:95
double precision, public hd_gamma
The adiabatic index.
Definition mod_hd_phys.t:69
integer, public, protected hd_trac_type
Definition mod_hd_phys.t:79
logical, public, protected hd_particles
Whether particles module is added.
Definition mod_hd_phys.t:33
type(tc_fluid), allocatable, public tc_fl
Definition mod_hd_phys.t:16
subroutine, public hd_check_w(primitive, ixil, ixol, w, flag)
Returns logical argument flag where values are ok.
logical, public, protected hd_viscosity
Whether viscosity is added.
Definition mod_hd_phys.t:27
subroutine, public hd_get_csound2(w, x, ixil, ixol, csound2)
Calculate the square of the thermal sound speed csound2 within ixO^L. csound2=gamma*p/rho.
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
Definition mod_hd_phys.t:66
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
Definition mod_hd_phys.t:91
integer, public, protected te_
Indices of temperature.
Definition mod_hd_phys.t:63
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
Definition mod_hd_phys.t:51
double precision, public, protected h_ion_fr
Ionization fraction of H H_ion_fr = H+/(H+ + H)
Definition mod_hd_phys.t:85
double precision function, dimension(ixo^s), public hd_kin_en(w, ixil, ixol, inv_rho)
subroutine, public hd_to_conserved(ixil, ixol, w, x)
Transform primitive variables into conservative ones.
logical, public, protected hd_cak_force
Whether CAK radiation line force is activated.
Definition mod_hd_phys.t:39
subroutine, public hd_phys_init()
Initialize the module.
integer, dimension(:), allocatable, public, protected tracer
Indices of the tracers.
Definition mod_hd_phys.t:54
logical, public, protected hd_thermal_conduction
Whether thermal conduction is added.
Definition mod_hd_phys.t:15
logical, public, protected eq_state_units
Definition mod_hd_phys.t:99
double precision, public hd_adiab
The adiabatic constant.
Definition mod_hd_phys.t:72
subroutine, public hd_to_primitive(ixil, ixol, w, x)
Transform conservative variables into primitive ones.
integer, public, protected rho_
Index of the density (in the w array)
Definition mod_hd_phys.t:48
logical, public, protected hd_partial_ionization
Whether plasma is partially ionized.
Definition mod_hd_phys.t:45
double precision, public, protected he_ion_fr
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
Definition mod_hd_phys.t:88
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
Definition mod_hd_phys.t:82
logical, public, protected hd_gravity
Whether gravity is added.
Definition mod_hd_phys.t:30
type(rc_fluid), allocatable, public rc_fl
Definition mod_hd_phys.t:21
logical, public, protected hd_trac
Whether TRAC method is used.
Definition mod_hd_phys.t:78
integer, public, protected hd_n_tracer
Number of tracer species.
Definition mod_hd_phys.t:42
type(te_fluid), allocatable, public te_fl_hd
Definition mod_hd_phys.t:17
logical, public, protected hd_rotating_frame
Whether rotating frame is activated.
Definition mod_hd_phys.t:36
subroutine, public hd_get_pthermal(w, x, ixil, ixol, pth)
Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L.
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
Definition mod_hd_phys.t:60
module ionization degree - get ionization degree for given temperature
subroutine ionization_degree_from_temperature(ixil, ixol, te, iz_h, iz_he)
Module containing all the particle routines.
subroutine particles_init()
Initialize particle data and parameters.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
module radiative cooling – add optically thin radiative cooling for HD and MHD
subroutine radiative_cooling_init_params(phys_gamma, he_abund)
Radiative cooling initialization.
subroutine cooling_get_dt(w, ixil, ixol, dtnew, dxd, x, fl)
subroutine radiative_cooling_init(fl, read_params)
subroutine radiative_cooling_add_source(qdt, ixil, ixol, wct, wctprim, w, x, qsourcesplit, active, fl)
Module for including rotating frame in (magneto)hydrodynamics simulations The rotation vector is assu...
subroutine rotating_frame_add_source(qdt, dtfactor, ixil, ixol, wct, w, x)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
subroutine rotating_frame_init()
Initialize the module.
Module for handling problematic values in simulations, such as negative pressures.
subroutine, public small_values_average(ixil, ixol, w, x, w_flag, windex)
logical, public trace_small_values
trace small values in the source file using traceback flag of compiler
subroutine, public small_values_error(wprim, x, ixil, ixol, w_flag, subname)
logical, dimension(:), allocatable, public small_values_fix_iw
Whether to apply small value fixes to certain variables.
character(len=20), public small_values_method
How to handle small values.
Generic supertimestepping method 1) in amrvac.par in sts_list set the following parameters which have...
subroutine, public add_sts_method(sts_getdt, sts_set_sources, startvar, nflux, startwbc, nwbc, evolve_b)
subroutine which added programatically a term to be calculated using STS Params: sts_getdt function c...
subroutine, public set_conversion_methods_to_head(sts_before_first_cycle, sts_after_last_cycle)
Set the hooks called before the first cycle and after the last cycle in the STS update This method sh...
subroutine, public set_error_handling_to_head(sts_error_handling)
Set the hook of error handling in the STS update. This method is called before updating the BC....
subroutine, public sts_init()
Initialize sts module.
Thermal conduction for HD and MHD or RHD and RMHD or twofl (plasma-neutral) module Adaptation of mod_...
subroutine, public tc_get_hd_params(fl, read_hd_params)
Init TC coefficients: HD case.
double precision function, public get_tc_dt_hd(w, ixil, ixol, dxd, x, fl)
Get the explicit timestep for the TC (hd implementation)
subroutine tc_init_params(phys_gamma)
subroutine, public sts_set_source_tc_hd(ixil, ixol, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
subroutine get_euv_image(qunit, fl)
subroutine get_sxr_image(qunit, fl)
subroutine get_euv_spectrum(qunit, fl)
subroutine get_whitelight_image(qunit, fl)
Module with all the methods that users can customize in AMRVAC.
procedure(rfactor), pointer usr_rfactor
procedure(set_surface), pointer usr_set_surface
procedure(phys_gravity), pointer usr_gravity
procedure(hd_pthermal), pointer usr_set_pthermal
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
The module add viscous source terms and check time step.
subroutine, public visc_get_flux_prim(w, x, ixil, ixol, idim, f, energy)
subroutine viscosity_add_source(qdt, ixil, ixol, wct, w, x, energy, qsourcesplit, active)
subroutine viscosity_init(phys_wider_stencil)
Initialize the module.
subroutine viscosity_get_dt(w, ixil, ixol, dtnew, dxd, x)
subroutine visc_add_source_geom(qdt, ixil, ixol, wct, w, x)