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