MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_ffhd_phys.t
Go to the documentation of this file.
1!> Frozen-field hydrodynamics module
3
4#include "amrvac.h"
5
6 use mod_global_parameters, only: std_len, const_c
10 use mod_physics
11 use mod_comm_lib, only: mpistop
12
13 implicit none
14 private
15
16 !> Whether an energy equation is used
17 logical, public, protected :: ffhd_energy = .true.
18
19 !> Whether thermal conduction is used
20 logical, public, protected :: ffhd_thermal_conduction = .false.
21 !> Whether hyperbolic type thermal conduction is used
22 logical, public, protected :: ffhd_hyperbolic_thermal_conduction = .false.
23 !> Wheterh saturation is considered for hyperbolic TC
24 logical, public, protected :: ffhd_htc_sat = .false.
25 !> type of fluid for thermal conduction
26 type(tc_fluid), public, allocatable :: tc_fl
27 !> type of fluid for thermal emission synthesis
28 type(te_fluid), public, allocatable :: te_fl_ffhd
29
30 !> Whether radiative cooling is added
31 logical, public, protected :: ffhd_radiative_cooling = .false.
32 !> type of fluid for radiative cooling
33 type(rc_fluid), public, allocatable :: rc_fl
34
35 !> Whether gravity is added
36 logical, public, protected :: ffhd_gravity = .false.
37
38 !> Whether TRAC method is used
39 logical, public, protected :: ffhd_trac = .false.
40
41 !> Which TRAC method is used
42 integer, public, protected :: ffhd_trac_type=1
43
44 !> Height of the mask used in the TRAC method
45 double precision, public, protected :: ffhd_trac_mask = 0.d0
46
47 !> Distance between two adjacent traced magnetic field lines (in finest cell size)
48 integer, public, protected :: ffhd_trac_finegrid=4
49
50 !> Whether plasma is partially ionized
51 logical, public, protected :: ffhd_partial_ionization = .false.
52
53 !> Index of the density (in the w array)
54 integer, public, protected :: rho_
55
56 !> Indices of the momentum density
57 integer, allocatable, public, protected :: mom(:)
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 integer, public, protected :: tweight_
71 integer, public, protected :: q_
72
73 !> The adiabatic index
74 double precision, public :: ffhd_gamma = 5.d0/3.0d0
75
76 !> The adiabatic constant
77 double precision, public :: ffhd_adiab = 1.0d0
78
79 !> The thermal conductivity kappa in hyperbolic thermal conduction
80 double precision, public :: hypertc_kappa
81
82 !> Helium abundance over Hydrogen
83 double precision, public, protected :: he_abundance=0.1d0
84 !> Ionization fraction of H
85 !> H_ion_fr = H+/(H+ + H)
86 double precision, public, protected :: h_ion_fr=1d0
87 !> Ionization fraction of He
88 !> He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
89 double precision, public, protected :: he_ion_fr=1d0
90 !> Ratio of number He2+ / number He+ + He2+
91 !> He_ion_fr2 = He2+/(He2+ + He+)
92 double precision, public, protected :: he_ion_fr2=1d0
93 ! used for eq of state when it is not defined by units,
94 ! the units do not contain terms related to ionization fraction
95 ! and it is p = RR * rho * T
96 double precision, public, protected :: rr=1d0
97 ! remove the below flag and assume default value = .false.
98 ! when eq state properly implemented everywhere
99 ! and not anymore through units
100 logical, public, protected :: eq_state_units = .true.
101
102 !> gamma minus one and its inverse
103 double precision :: gamma_1, inv_gamma_1
104
105 !define the function interface for the kinetic energy
106 abstract interface
107
108 function fun_kin_en(w, ixI^L, ixO^L, inv_rho) result(ke)
109 use mod_global_parameters, only: nw, ndim,block
110 integer, intent(in) :: ixi^l, ixo^l
111 double precision, intent(in) :: w(ixi^s, nw)
112 double precision :: ke(ixo^s)
113 double precision, intent(in), optional :: inv_rho(ixo^s)
114 end function fun_kin_en
115
116 end interface
117
118 procedure(sub_convert), pointer :: ffhd_to_primitive => null()
119 procedure(sub_convert), pointer :: ffhd_to_conserved => null()
120 procedure(sub_small_values), pointer :: ffhd_handle_small_values => null()
121 procedure(sub_get_pthermal), pointer :: ffhd_get_pthermal => null()
122 procedure(sub_get_pthermal), pointer :: ffhd_get_rfactor => null()
123 procedure(sub_get_pthermal), pointer :: ffhd_get_temperature => null()
124 procedure(sub_get_v), pointer :: ffhd_get_v => null()
125 procedure(fun_kin_en), pointer :: ffhd_kin_en => null()
126 ! Public methods
127 public :: ffhd_phys_init
128 public :: ffhd_kin_en
129 public :: ffhd_get_pthermal
130 public :: ffhd_get_temperature
131 public :: ffhd_get_v
132 public :: ffhd_get_rho
133 public :: ffhd_get_v_idim
134 public :: ffhd_to_conserved
135 public :: ffhd_to_primitive
136 public :: ffhd_get_csound2
137 public :: ffhd_e_to_ei
138 public :: ffhd_ei_to_e
139
140contains
141
142 subroutine ffhd_read_params(files)
144 character(len=*), intent(in) :: files(:)
145 integer :: n
146
147 namelist /ffhd_list/ ffhd_energy, ffhd_gamma, ffhd_adiab, &
151
152 do n = 1, size(files)
153 open(unitpar, file=trim(files(n)), status="old")
154 read(unitpar, ffhd_list, end=111)
155111 close(unitpar)
156 end do
157 end subroutine ffhd_read_params
158
159 !> Write this module's parameters to a snapsoht
160 subroutine ffhd_write_info(fh)
162 integer, intent(in) :: fh
163 integer, parameter :: n_par = 1
164 double precision :: values(n_par)
165 character(len=name_len) :: names(n_par)
166 integer, dimension(MPI_STATUS_SIZE) :: st
167 integer :: er
168
169 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
170
171 names(1) = "gamma"
172 values(1) = ffhd_gamma
173 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
174 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
175 end subroutine ffhd_write_info
176
177 subroutine ffhd_phys_init()
181 use mod_gravity, only: gravity_init
186 integer :: itr, idir
187
188 call ffhd_read_params(par_files)
189
190 if(.not. ffhd_energy) then
193 if(mype==0) write(*,*) 'WARNING: set ffhd_thermal_conduction=F when ffhd_energy=F'
194 end if
197 if(mype==0) write(*,*) 'WARNING: set ffhd_hyperbolic_thermal_conduction=F when ffhd_energy=F'
198 end if
201 if(mype==0) write(*,*) 'WARNING: set ffhd_radiative_cooling=F when ffhd_energy=F'
202 end if
203 if(ffhd_trac) then
204 ffhd_trac=.false.
205 if(mype==0) write(*,*) 'WARNING: set ffhd_trac=F when ffhd_energy=F'
206 end if
209 if(mype==0) write(*,*) 'WARNING: set ffhd_partial_ionization=F when ffhd_energy=F'
210 end if
211 end if
212 if(.not.eq_state_units) then
215 if(mype==0) write(*,*) 'WARNING: set ffhd_partial_ionization=F when eq_state_units=F'
216 end if
217 end if
218
221 if(mype==0) write(*,*) 'WARNING: turn off parabolic TC when using hyperbolic TC'
222 end if
223
224 physics_type = "ffhd"
225 phys_energy=ffhd_energy
226 phys_internal_e=.false.
229 phys_partial_ionization=ffhd_partial_ionization
230
231 phys_gamma = ffhd_gamma
232 phys_total_energy=ffhd_energy
234
235 {^ifoned
236 if(ffhd_trac .and. ffhd_trac_type .gt. 2) then
238 if(mype==0) write(*,*) 'WARNING: reset ffhd_trac_type=1 for 1D simulation'
239 end if
240 }
241 if(ffhd_trac .and. ffhd_trac_type .le. 4) then
242 ffhd_trac_mask=bigdouble
243 if(mype==0) write(*,*) 'WARNING: set ffhd_trac_mask==bigdouble for global TRAC method'
244 end if
246
247 allocate(start_indices(number_species),stop_indices(number_species))
248 start_indices(1)=1
249 ! Determine flux variables
250 rho_ = var_set_rho()
251
252 allocate(mom(1))
253 mom(:) = var_set_momentum(1)
254
255 ! Set index of energy variable
256 if(ffhd_energy) then
257 e_ = var_set_energy() ! energy density
258 p_ = e_ ! gas pressure
259 else
260 e_ = -1
261 p_ = -1
262 end if
263
265 q_ = var_set_q()
266 need_global_cmax=.true.
267 else
268 q_=-1
269 end if
270
271 ! set number of variables which need update ghostcells
272 nwgc=nwflux
273
274 ! set the index of the last flux variable for species 1
275 stop_indices(1)=nwflux
276
277 ! set temperature as an auxiliary variable to get ionization degree
279 te_ = var_set_auxvar('Te','Te')
280 else
281 te_ = -1
282 end if
283
284 ! set cutoff temperature when using the TRAC method, as well as an auxiliary weight
285 tweight_ = -1
286 if(ffhd_trac) then
287 tcoff_ = var_set_wextra()
288 iw_tcoff=tcoff_
289 if(ffhd_trac_type .ge. 3) then
290 tweight_ = var_set_wextra()
291 iw_tweight=tweight_
292 end if
293 else
294 tcoff_ = -1
295 end if
296
297 nvector = 0 ! No. vector vars
298
299 ! Check whether custom flux types have been defined
300 if(.not. allocated(flux_type)) then
301 allocate(flux_type(ndir, nwflux))
302 flux_type = flux_default
303 else if(any(shape(flux_type) /= [ndir, nwflux])) then
304 call mpistop("phys_check error: flux_type has wrong shape")
305 end if
306
307 phys_get_dt => ffhd_get_dt
308 phys_get_cmax => ffhd_get_cmax_origin
309 phys_get_tcutoff => ffhd_get_tcutoff
310 phys_get_cbounds => ffhd_get_cbounds
311 phys_to_primitive => ffhd_to_primitive_origin
312 ffhd_to_primitive => ffhd_to_primitive_origin
313 phys_to_conserved => ffhd_to_conserved_origin
314 ffhd_to_conserved => ffhd_to_conserved_origin
315 phys_get_flux => ffhd_get_flux
316 phys_get_v => ffhd_get_v_origin
317 ffhd_get_v => ffhd_get_v_origin
318 phys_get_rho => ffhd_get_rho
319 ffhd_kin_en => ffhd_kin_en_origin
320 phys_add_source_geom => ffhd_add_source_geom
321 phys_add_source => ffhd_add_source
322 phys_check_params => ffhd_check_params
323 phys_write_info => ffhd_write_info
324 phys_handle_small_values => ffhd_handle_small_values_origin
325 ffhd_handle_small_values => ffhd_handle_small_values_origin
326 phys_check_w => ffhd_check_w_origin
327
328 if(.not.ffhd_energy) then
329 phys_get_pthermal => ffhd_get_pthermal_iso
330 ffhd_get_pthermal => ffhd_get_pthermal_iso
331 else
332 phys_get_pthermal => ffhd_get_pthermal_origin
333 ffhd_get_pthermal => ffhd_get_pthermal_origin
334 end if
335
336 ! choose Rfactor in ideal gas law
338 ffhd_get_rfactor=>rfactor_from_temperature_ionization
339 phys_update_temperature => ffhd_update_temperature
340 else if(associated(usr_rfactor)) then
341 ffhd_get_rfactor=>usr_rfactor
342 else
343 ffhd_get_rfactor=>rfactor_from_constant_ionization
344 end if
345
347 ffhd_get_temperature => ffhd_get_temperature_from_te
348 else
349 ffhd_get_temperature => ffhd_get_temperature_from_etot
350 end if
351
352 ! derive units from basic units
353 call ffhd_physical_units()
354
356 if(si_unit)then
357 ! parallel conduction Spitzer
359 else
360 ! in cgs
362 endif
363 end if
364 if(.not. ffhd_energy .and. ffhd_thermal_conduction) then
365 call mpistop("thermal conduction needs ffhd_energy=T")
366 end if
368 call mpistop("hyperbolic thermal conduction needs ffhd_energy=T")
369 end if
370 if(.not. ffhd_energy .and. ffhd_radiative_cooling) then
371 call mpistop("radiative cooling needs ffhd_energy=T")
372 end if
373
374 ! initialize thermal conduction module
376 call sts_init()
378
379 allocate(tc_fl)
380 call tc_get_hd_params(tc_fl,tc_params_read_ffhd)
381 call add_sts_method(ffhd_get_tc_dt_ffhd,ffhd_sts_set_source_tc_ffhd,e_,1,e_,1,.false.)
382 tc_fl%get_temperature_from_conserved => ffhd_get_temperature_from_etot
383 tc_fl%get_temperature_from_eint => ffhd_get_temperature_from_eint
385 call set_error_handling_to_head(ffhd_tc_handle_small_e)
386 tc_fl%get_rho => ffhd_get_rho
387 tc_fl%e_ = e_
388 tc_fl%Tcoff_ = tcoff_
389 end if
390
391 ! Initialize radiative cooling module
394 allocate(rc_fl)
395 call radiative_cooling_init(rc_fl,rc_params_read)
396 rc_fl%get_rho => ffhd_get_rho
397 rc_fl%get_pthermal => ffhd_get_pthermal
398 rc_fl%get_var_Rfactor => ffhd_get_rfactor
399 rc_fl%e_ = e_
400 rc_fl%Tcoff_ = tcoff_
401 rc_fl%has_equi = .false.
402 rc_fl%subtract_equi = .false.
403 end if
404
405{^ifthreed
406 allocate(te_fl_ffhd)
407 te_fl_ffhd%get_rho=> ffhd_get_rho
408 te_fl_ffhd%get_pthermal=> ffhd_get_pthermal
409 te_fl_ffhd%get_var_Rfactor => ffhd_get_rfactor
410 phys_te_images => ffhd_te_images
411}
412
413 ! Initialize gravity module
414 if(ffhd_gravity) then
415 call gravity_init()
416 end if
417
418 ! initialize ionization degree table
420 end subroutine ffhd_phys_init
421
422{^ifthreed
423 subroutine ffhd_te_images
426
427 select case(convert_type)
428 case('EIvtiCCmpi','EIvtuCCmpi')
430 case('ESvtiCCmpi','ESvtuCCmpi')
432 case('SIvtiCCmpi','SIvtuCCmpi')
434 case('WIvtiCCmpi','WIvtuCCmpi')
436 case default
437 call mpistop("Error in synthesize emission: Unknown convert_type")
438 end select
439 end subroutine ffhd_te_images
440}
441
442 subroutine ffhd_sts_set_source_tc_ffhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
446 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
447 double precision, intent(in) :: x(ixi^s,1:ndim)
448 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
449 double precision, intent(in) :: my_dt
450 logical, intent(in) :: fix_conserve_at_step
451 call sts_set_source_tc_mhd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
452 end subroutine ffhd_sts_set_source_tc_ffhd
453
454 function ffhd_get_tc_dt_ffhd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
455 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
456 !where tc_k_para_i=tc_k_para*B_i**2/B**2
457 !and T=p/rho
460
461 integer, intent(in) :: ixi^l, ixo^l
462 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
463 double precision, intent(in) :: w(ixi^s,1:nw)
464 double precision :: dtnew
465
466 dtnew=get_tc_dt_mhd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
467 end function ffhd_get_tc_dt_ffhd
468
469 subroutine ffhd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
471
472 integer, intent(in) :: ixi^l,ixo^l
473 double precision, intent(inout) :: w(ixi^s,1:nw)
474 double precision, intent(in) :: x(ixi^s,1:ndim)
475 integer, intent(in) :: step
476 character(len=140) :: error_msg
477
478 write(error_msg,"(a,i3)") "Thermal conduction step ", step
479 call ffhd_handle_small_ei(w,x,ixi^l,ixo^l,e_,error_msg)
480 end subroutine ffhd_tc_handle_small_e
481
482 subroutine tc_params_read_ffhd(fl)
484 type(tc_fluid), intent(inout) :: fl
485 integer :: n
486 ! list parameters
487 logical :: tc_saturate=.false.
488 double precision :: tc_k_para=0d0
489 character(len=std_len) :: tc_slope_limiter="MC"
490
491 namelist /tc_list/ tc_saturate, tc_slope_limiter, tc_k_para
492
493 do n = 1, size(par_files)
494 open(unitpar, file=trim(par_files(n)), status="old")
495 read(unitpar, tc_list, end=111)
496111 close(unitpar)
497 end do
498
499 fl%tc_saturate = tc_saturate
500 fl%tc_k_para = tc_k_para
501 select case(tc_slope_limiter)
502 case ('no','none')
503 fl%tc_slope_limiter = 0
504 case ('MC')
505 ! montonized central limiter Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
506 fl%tc_slope_limiter = 1
507 case('minmod')
508 ! minmod limiter
509 fl%tc_slope_limiter = 2
510 case ('superbee')
511 ! Roes superbee limiter (eq.3.51i)
512 fl%tc_slope_limiter = 3
513 case ('koren')
514 ! Barry Koren Right variant
515 fl%tc_slope_limiter = 4
516 case default
517 call mpistop("Unknown tc_slope_limiter, choose MC, minmod")
518 end select
519 end subroutine tc_params_read_ffhd
520
521 subroutine rc_params_read(fl)
523 use mod_constants, only: bigdouble
524 type(rc_fluid), intent(inout) :: fl
525 integer :: n
526 integer :: ncool = 4000
527
528 !> Name of cooling curve
529 character(len=std_len) :: coolcurve='JCcorona'
530
531 !> Fixed temperature not lower than tlow
532 logical :: tfix=.false.
533
534 !> Lower limit of temperature
535 double precision :: tlow=bigdouble
536
537 !> Add cooling source in a split way (.true.) or un-split way (.false.)
538 logical :: rc_split=.false.
539 logical :: rad_cut=.false.
540 double precision :: rad_cut_hgt=0.5d0
541 double precision :: rad_cut_dey=0.15d0
542
543 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split, rad_cut, rad_cut_hgt, rad_cut_dey
544
545 do n = 1, size(par_files)
546 open(unitpar, file=trim(par_files(n)), status="old")
547 read(unitpar, rc_list, end=111)
548111 close(unitpar)
549 end do
550
551 fl%ncool=ncool
552 fl%coolcurve=coolcurve
553 fl%tlow=tlow
554 fl%Tfix=tfix
555 fl%rc_split=rc_split
556 fl%rad_cut=rad_cut
557 fl%rad_cut_hgt=rad_cut_hgt
558 fl%rad_cut_dey=rad_cut_dey
559 end subroutine rc_params_read
560
561 subroutine ffhd_check_params
565 use mod_geometry, only: coordinate
566
567
568 gamma_1=ffhd_gamma-1.d0
569 if (.not. ffhd_energy) then
570 if (ffhd_gamma <= 0.0d0) call mpistop ("Error: ffhd_gamma <= 0")
571 if (ffhd_adiab < 0.0d0) call mpistop ("Error: ffhd_adiab < 0")
573 else
574 if (ffhd_gamma <= 0.0d0 .or. ffhd_gamma == 1.0d0) &
575 call mpistop ("Error: ffhd_gamma <= 0 or ffhd_gamma == 1")
576 inv_gamma_1=1.d0/gamma_1
577 small_e = small_pressure * inv_gamma_1
578 end if
579
580 if (number_equi_vars > 0 .and. .not. associated(usr_set_equi_vars)) then
581 call mpistop("usr_set_equi_vars has to be implemented in the user file")
582 end if
583
584
585 if(mype==0)then
586 write(*,*)'====FFHD run with settings===================='
587 write(*,*)'Using mod_ffhd_phys with settings:'
588 write(*,*)'SI_unit=',si_unit
589 write(*,*)'Dimensionality :',ndim
590 write(*,*)'vector components:',ndir
591 write(*,*)'coordinate set to type,slab:',coordinate,slab
592 write(*,*)'number of variables nw=',nw
593 write(*,*)' start index iwstart=',iwstart
594 write(*,*)'number of vector variables=',nvector
595 write(*,*)'number of stagger variables nws=',nws
596 write(*,*)'number of variables with BCs=',nwgc
597 write(*,*)'number of vars with fluxes=',nwflux
598 write(*,*)'number of vars with flux + BC=',nwfluxbc
599 write(*,*)'number of auxiliary variables=',nwaux
600 write(*,*)'number of extra vars without flux=',nwextra
601 write(*,*)'number of extra vars for wextra=',nw_extra
602 write(*,*)'number of auxiliary I/O variables=',nwauxio
603 write(*,*)' ffhd_energy=',ffhd_energy
604 write(*,*)' ffhd_gravity=',ffhd_gravity
605 write(*,*)' ffhd_radiative_cooling=',ffhd_radiative_cooling
606 write(*,*)' ffhd_hyperbolic_thermal_conduction=',ffhd_hyperbolic_thermal_conduction
607 write(*,*)' ffhd_trac=',ffhd_trac
608 write(*,*)'number of ghostcells=',nghostcells
609 write(*,*)'number due to phys_wider_stencil=',phys_wider_stencil
610 write(*,*)'==========================================='
611 endif
612 end subroutine ffhd_check_params
613
614 subroutine ffhd_physical_units()
616 double precision :: mp,kb
617 double precision :: a,b
618
619 if(si_unit) then
620 mp=mp_si
621 kb=kb_si
622 else
623 mp=mp_cgs
624 kb=kb_cgs
625 end if
626 if(eq_state_units) then
627 a=1d0+4d0*he_abundance
629 b=1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0)
630 else
631 b=2d0+3d0*he_abundance
632 end if
633 rr=1d0
634 else
635 a=1d0
636 b=1d0
637 rr=(1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0))/(1d0+4d0*he_abundance)
638 end if
639 if(unit_density/=1.d0 .or. unit_numberdensity/=1.d0) then
640 if(unit_density/=1.d0) then
642 else if(unit_numberdensity/=1.d0) then
644 end if
645 if(unit_temperature/=1.d0) then
648 if(unit_length/=1.d0) then
650 else if(unit_time/=1.d0) then
652 end if
653 else if(unit_pressure/=1.d0) then
656 if(unit_length/=1.d0) then
658 else if(unit_time/=1.d0) then
660 end if
661 else if(unit_velocity/=1.d0) then
664 if(unit_length/=1.d0) then
666 else if(unit_time/=1.d0) then
668 end if
669 else if(unit_time/=1.d0) then
673 end if
674 else if(unit_temperature/=1.d0) then
675 ! units of temperature and velocity are dependent
676 if(unit_pressure/=1.d0) then
680 if(unit_length/=1.d0) then
682 else if(unit_time/=1.d0) then
684 end if
685 end if
686 else if(unit_pressure/=1.d0) then
687 if(unit_velocity/=1.d0) then
691 if(unit_length/=1.d0) then
693 else if(unit_time/=1.d0) then
695 end if
696 else if(unit_time/=0.d0) then
701 end if
702 end if
704 end subroutine ffhd_physical_units
705
706 subroutine ffhd_check_w_origin(primitive,ixI^L,ixO^L,w,flag)
708 logical, intent(in) :: primitive
709 integer, intent(in) :: ixi^l, ixo^l
710 double precision, intent(in) :: w(ixi^s,nw)
711 double precision :: tmp(ixi^s)
712 logical, intent(inout) :: flag(ixi^s,1:nw)
713
714 flag=.false.
715 where(w(ixo^s,rho_) < small_density) flag(ixo^s,rho_) = .true.
716
717 if(ffhd_energy) then
718 if(primitive) then
719 where(w(ixo^s,e_) < small_pressure) flag(ixo^s,e_) = .true.
720 else
721 tmp(ixo^s)=w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l)
722 where(tmp(ixo^s) < small_e) flag(ixo^s,e_) = .true.
723 end if
724 end if
725 end subroutine ffhd_check_w_origin
726
727 subroutine ffhd_to_conserved_origin(ixI^L,ixO^L,w,x)
729 integer, intent(in) :: ixi^l, ixo^l
730 double precision, intent(inout) :: w(ixi^s, nw)
731 double precision, intent(in) :: x(ixi^s, 1:ndim)
732
733 if(ffhd_energy) then
734 w(ixo^s,e_)=w(ixo^s,p_)*inv_gamma_1+half*w(ixo^s,mom(1))**2*w(ixo^s,rho_)
735 end if
736 w(ixo^s,mom(1))=w(ixo^s,rho_)*w(ixo^s,mom(1))
737 end subroutine ffhd_to_conserved_origin
738
739 subroutine ffhd_to_primitive_origin(ixI^L,ixO^L,w,x)
741 integer, intent(in) :: ixi^l, ixo^l
742 double precision, intent(inout) :: w(ixi^s, nw)
743 double precision, intent(in) :: x(ixi^s, 1:ndim)
744
745 if(fix_small_values) then
746 call ffhd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'ffhd_to_primitive_origin')
747 end if
748
749 w(ixo^s,mom(1)) = w(ixo^s,mom(1))/w(ixo^s,rho_)
750 if(ffhd_energy) then
751 w(ixo^s,p_)=gamma_1*(w(ixo^s,e_)-half*w(ixo^s,rho_)*w(ixo^s,mom(1))**2)
752 end if
753 end subroutine ffhd_to_primitive_origin
754
755 subroutine ffhd_ei_to_e(ixI^L,ixO^L,w,x)
757 integer, intent(in) :: ixi^l, ixo^l
758 double precision, intent(inout) :: w(ixi^s, nw)
759 double precision, intent(in) :: x(ixi^s, 1:ndim)
760
761 w(ixi^s,e_)=w(ixi^s,e_)+ffhd_kin_en(w,ixi^l,ixi^l)
762 end subroutine ffhd_ei_to_e
763
764 subroutine ffhd_e_to_ei(ixI^L,ixO^L,w,x)
766 integer, intent(in) :: ixi^l, ixo^l
767 double precision, intent(inout) :: w(ixi^s, nw)
768 double precision, intent(in) :: x(ixi^s, 1:ndim)
769
770 w(ixi^s,e_)=w(ixi^s,e_)-ffhd_kin_en(w,ixi^l,ixi^l)
771 if(fix_small_values) then
772 call ffhd_handle_small_ei(w,x,ixi^l,ixi^l,e_,'ffhd_e_to_ei')
773 end if
774 end subroutine ffhd_e_to_ei
775
776 subroutine ffhd_handle_small_values_origin(primitive, w, x, ixI^L, ixO^L, subname)
779 logical, intent(in) :: primitive
780 integer, intent(in) :: ixi^l,ixo^l
781 double precision, intent(inout) :: w(ixi^s,1:nw)
782 double precision, intent(in) :: x(ixi^s,1:ndim)
783 character(len=*), intent(in) :: subname
784
785 logical :: flag(ixi^s,1:nw)
786 double precision :: tmp2(ixi^s)
787
788 call phys_check_w(primitive, ixi^l, ixi^l, w, flag)
789
790 if(any(flag)) then
791 select case (small_values_method)
792 case ("replace")
793 where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
794 if(small_values_fix_iw(mom(1))) then
795 where(flag(ixo^s,rho_)) w(ixo^s, mom(1)) = 0.0d0
796 end if
797 if(ffhd_energy) then
798 if(primitive) then
799 where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
800 else
801 where(flag(ixo^s,e_))
802 w(ixo^s,e_) = small_e+ffhd_kin_en(w,ixi^l,ixo^l)
803 end where
804 end if
805 end if
806 case ("average")
807 call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
808 if(ffhd_energy) then
809 if(primitive) then
810 call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
811 else
812 w(ixi^s,e_)=w(ixi^s,e_)-ffhd_kin_en(w,ixi^l,ixi^l)
813 call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
814 w(ixi^s,e_)=w(ixi^s,e_)+ffhd_kin_en(w,ixi^l,ixi^l)
815 end if
816 end if
817 case default
818 if(.not.primitive) then
819 if(ffhd_energy) then
820 w(ixo^s,p_)=gamma_1*(w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l))
821 end if
822 w(ixo^s,mom(1))=w(ixo^s,mom(1))/w(ixo^s,rho_)
823 end if
824 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
825 end select
826 end if
827 end subroutine ffhd_handle_small_values_origin
828
829 subroutine ffhd_get_v_origin(w,x,ixI^L,ixO^L,v)
831 integer, intent(in) :: ixi^l, ixo^l
832 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
833 double precision, intent(out) :: v(ixi^s,ndir)
834 double precision :: rho(ixi^s)
835 integer :: idir
836
837 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
838 rho(ixo^s)=1.d0/rho(ixo^s)
839 do idir=1,ndir
840 v(ixo^s,ndir) = w(ixo^s,mom(1))*block%B0(ixo^s,idir,0)*rho(ixo^s)
841 end do
842 end subroutine ffhd_get_v_origin
843
844 subroutine ffhd_get_v_idim(w,x,ixI^L,ixO^L,idim,v)
846 integer, intent(in) :: ixi^l, ixo^l, idim
847 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
848 double precision, intent(out) :: v(ixi^s)
849 double precision :: rho(ixi^s)
850
851 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
852 v(ixo^s) = (w(ixo^s, mom(1))*block%B0(ixo^s,idim,0)) / rho(ixo^s)
853 end subroutine ffhd_get_v_idim
854
855 subroutine ffhd_get_cmax_origin(wprim,x,ixI^L,ixO^L,idim,cmax)
857 integer, intent(in) :: ixi^l, ixo^l, idim
858 ! w in primitive form
859 double precision, intent(in) :: wprim(ixi^s, nw), x(ixi^s,1:ndim)
860 double precision, intent(inout) :: cmax(ixi^s)
861
862 if(ffhd_energy) then
863 cmax(ixo^s)=dsqrt(ffhd_gamma*wprim(ixo^s,p_)/wprim(ixo^s,rho_))
864 else
865 cmax(ixo^s)=dsqrt(ffhd_gamma*ffhd_adiab*wprim(ixo^s,rho_)**gamma_1)
866 end if
867 cmax(ixo^s)=dabs(wprim(ixo^s,mom(1))*block%B0(ixo^s,idim,0))+cmax(ixo^s)
868
869 end subroutine ffhd_get_cmax_origin
870
871 subroutine ffhd_get_tcutoff(ixI^L,ixO^L,w,x,Tco_local,Tmax_local)
873 use mod_geometry
874 integer, intent(in) :: ixi^l,ixo^l
875 double precision, intent(in) :: x(ixi^s,1:ndim)
876 double precision, intent(inout) :: w(ixi^s,1:nw)
877 double precision, intent(out) :: tco_local,tmax_local
878 double precision, parameter :: trac_delta=0.25d0
879 double precision :: te(ixi^s),lts(ixi^s)
880 double precision, dimension(ixI^S,1:ndim) :: gradt
881 double precision :: bdir(ndim)
882 double precision :: ltrc,ltrp,altr
883 integer :: idims,ix^d,jxo^l,hxo^l,ixa^d,ixb^d
884 integer :: jxp^l,hxp^l,ixp^l,ixq^l
885 logical :: lrlt(ixi^s)
886
887 call ffhd_get_temperature(w,x,ixi^l,ixi^l,te)
888 tco_local=zero
889 tmax_local=maxval(te(ixo^s))
890
891 {^ifoned
892 select case(ffhd_trac_type)
893 case(0)
894 !> test case, fixed cutoff temperature
895 block%wextra(ixi^s,tcoff_)=2.5d5/unit_temperature
896 case(1)
897 hxo^l=ixo^l-1;
898 jxo^l=ixo^l+1;
899 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
900 lrlt=.false.
901 where(lts(ixo^s) > trac_delta)
902 lrlt(ixo^s)=.true.
903 end where
904 if(any(lrlt(ixo^s))) then
905 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
906 end if
907 case(2)
908 !> iijima et al. 2021, LTRAC method
909 ltrc=1.5d0
910 ltrp=4.d0
911 ixp^l=ixo^l^ladd1;
912 hxo^l=ixo^l-1;
913 jxo^l=ixo^l+1;
914 hxp^l=ixp^l-1;
915 jxp^l=ixp^l+1;
916 lts(ixp^s)=0.5d0*dabs(te(jxp^s)-te(hxp^s))/te(ixp^s)
917 lts(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
918 lts(ixo^s)=0.25d0*(lts(jxo^s)+two*lts(ixo^s)+lts(hxo^s))
919 block%wextra(ixo^s,tcoff_)=te(ixo^s)*lts(ixo^s)**0.4d0
920 case default
921 call mpistop("ffhd_trac_type not allowed for 1D simulation")
922 end select
923 }
924 {^nooned
925 select case(ffhd_trac_type)
926 case(0)
927 !> test case, fixed cutoff temperature
928 if(slab_uniform) then
929 !> assume cgs units
930 block%wextra(ixi^s,tcoff_)=max(min(3.d5/unit_temperature,6.d5/unit_temperature-3.d-4/unit_temperature*unit_length*x(ixi^s,ndim)),zero)
931 else
932 block%wextra(ixi^s,tcoff_)=2.5d5/unit_temperature
933 end if
934 case(1,4,6)
935 do idims=1,ndim
936 call gradient(te,ixi^l,ixo^l,idims,gradt(ixi^s,idims))
937 end do
938 if(ffhd_trac_type .gt. 1) then
939 ! B direction at cell center
940 bdir=zero
941 {do ixa^d=0,1\}
942 ixb^d=(ixomin^d+ixomax^d-1)/2+ixa^d;
943 bdir(1:ndim)=bdir(1:ndim)+block%B0(ixb^d,1:ndim,0)
944 {end do\}
945 if(sum(bdir(:)**2) .gt. zero) then
946 bdir(1:ndim)=bdir(1:ndim)/dsqrt(sum(bdir(:)**2))
947 end if
948 block%special_values(3:ndim+2)=bdir(1:ndim)
949 end if
950 {do ix^db=ixomin^db,ixomax^db\}
951 {^iftwod
952 ! temperature length scale inversed
953 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2))*&
954 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
955 }
956 {^ifthreed
957 ! temperature length scale inversed
958 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2),block%ds(ix^d,3))*&
959 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
960 }
961 if(lts(ix^d)>trac_delta) then
962 block%special_values(1)=max(block%special_values(1),te(ix^d))
963 end if
964 {end do\}
965 block%special_values(2)=tmax_local
966 case(2)
967 !> iijima et al. 2021, LTRAC method
968 ltrc=1.5d0
969 ltrp=4.d0
970 ixp^l=ixo^l^ladd2;
971 do idims=1,ndim
972 ixq^l=ixp^l;
973 hxp^l=ixp^l;
974 jxp^l=ixp^l;
975 select case(idims)
976 {case(^d)
977 ixqmin^d=ixqmin^d+1
978 ixqmax^d=ixqmax^d-1
979 hxpmax^d=ixpmin^d
980 jxpmin^d=ixpmax^d
981 \}
982 end select
983 call gradient(te,ixi^l,ixq^l,idims,gradt(ixi^s,idims))
984 call gradientf(te,x,ixi^l,hxp^l,idims,gradt(ixi^s,idims),nghostcells,.true.)
985 call gradientf(te,x,ixi^l,jxp^l,idims,gradt(ixi^s,idims),nghostcells,.false.)
986 end do
987 {do ix^db=ixpmin^db,ixpmax^db\}
988 ! temperature length scale inversed
989 lts(ix^d)=abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
990 ! fraction of cells size to temperature length scale
991 lts(ix^d)=min(^d&block%ds({ix^d},^d))*lts(ix^d)
992 lts(ix^d)=max(one,(exp(lts(ix^d))/ltrc)**ltrp)
993 {end do\}
994
995 ! need one ghost layer for thermal conductivity
996 ixp^l=ixo^l^ladd1;
997 {do ix^db=ixpmin^db,ixpmax^db\}
998 {^iftwod
999 altr=0.25d0*((lts(ix1-1,ix2)+two*lts(ix^d)+lts(ix1+1,ix2))*block%B0(ix^d,1,0)**2+&
1000 (lts(ix1,ix2-1)+two*lts(ix^d)+lts(ix1,ix2+1))*block%B0(ix^d,2,0)**2)
1001 block%wextra(ix^d,tcoff_)=te(ix^d)*altr**0.4d0
1002 }
1003 {^ifthreed
1004 altr=0.25d0*((lts(ix1-1,ix2,ix3)+two*lts(ix^d)+lts(ix1+1,ix2,ix3))*block%B0(ix^d,1,0)**2+&
1005 (lts(ix1,ix2-1,ix3)+two*lts(ix^d)+lts(ix1,ix2+1,ix3))*block%B0(ix^d,2,0)**2+&
1006 (lts(ix1,ix2,ix3-1)+two*lts(ix^d)+lts(ix1,ix2,ix3+1))*block%B0(ix^d,3,0)**2)
1007 block%wextra(ix^d,tcoff_)=te(ix^d)*altr**0.4d0
1008 }
1009 {end do\}
1010 case(3,5)
1011 !> do nothing here
1012 case default
1013 call mpistop("unknown ffhd_trac_type")
1014 end select
1015 }
1016 end subroutine ffhd_get_tcutoff
1017
1018 subroutine ffhd_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
1020 integer, intent(in) :: ixi^l, ixo^l, idim
1021 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
1022 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
1023 double precision, intent(in) :: x(ixi^s,1:ndim)
1024 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
1025 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
1026 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
1027 double precision :: wmean(ixi^s,nw)
1028 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1029
1030 select case (boundspeed)
1031 case (1)
1032 ! This implements formula (10.52) from "Riemann Solvers and Numerical
1033 ! Methods for Fluid Dynamics" by Toro.
1034 tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
1035 tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
1036 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1037 umean(ixo^s)=(wlp(ixo^s,mom(1))*tmp1(ixo^s)&
1038 +wrp(ixo^s,mom(1))*tmp2(ixo^s))*tmp3(ixo^s)
1039 umean(ixo^s)=umean(ixo^s)*block%B0(ixo^s,idim,idim)
1040 call ffhd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1041 call ffhd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1042 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1043 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1044 ((wrp(ixo^s,mom(1))-wlp(ixo^s,mom(1)))*block%B0(ixo^s,idim,idim))**2
1045 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1046 if(present(cmin)) then
1047 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1048 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1049 else
1050 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1051 end if
1052 case (2)
1053 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1054 tmp1(ixo^s)=wmean(ixo^s,mom(1))*block%B0(ixo^s,idim,idim)/wmean(ixo^s,rho_)
1055 call ffhd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
1056 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1057 if(present(cmin)) then
1058 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1059 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1060 else
1061 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1062 end if
1063 case (3)
1064 ! Miyoshi 2005 JCP 208, 315 equation (67)
1065 call ffhd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1066 call ffhd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1067 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1068 if(present(cmin)) then
1069 cmin(ixo^s,1)=min(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1070 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))-csoundl(ixo^s)
1071 cmax(ixo^s,1)=max(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1072 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1073 else
1074 cmax(ixo^s,1)=max(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1075 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1076 end if
1077 end select
1078 end subroutine ffhd_get_cbounds
1079
1080 subroutine ffhd_get_pthermal_iso(w,x,ixI^L,ixO^L,pth)
1082
1083 integer, intent(in) :: ixi^l, ixo^l
1084 double precision, intent(in) :: w(ixi^s,nw)
1085 double precision, intent(in) :: x(ixi^s,1:ndim)
1086 double precision, intent(out):: pth(ixi^s)
1087
1088 call ffhd_get_rho(w,x,ixi^l,ixo^l,pth)
1089 pth(ixo^s)=ffhd_adiab*pth(ixo^s)**ffhd_gamma
1090 end subroutine ffhd_get_pthermal_iso
1091
1092 subroutine ffhd_get_pthermal_origin(w,x,ixI^L,ixO^L,pth)
1095 integer, intent(in) :: ixi^l, ixo^l
1096 double precision, intent(in) :: w(ixi^s,nw)
1097 double precision, intent(in) :: x(ixi^s,1:ndim)
1098 double precision, intent(out):: pth(ixi^s)
1099 integer :: iw, ix^d
1100
1101 pth(ixo^s)=gamma_1*(w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l))
1102 if (fix_small_values) then
1103 {do ix^db= ixo^lim^db\}
1104 if(pth(ix^d)<small_pressure) then
1105 pth(ix^d)=small_pressure
1106 end if
1107 {end do^D&\}
1108 elseif(check_small_values) then
1109 {do ix^db= ixo^lim^db\}
1110 if(pth(ix^d)<small_pressure) then
1111 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1112 " encountered when call ffhd_get_pthermal"
1113 write(*,*) "Iteration: ", it, " Time: ", global_time
1114 write(*,*) "Location: ", x(ix^d,:)
1115 write(*,*) "Cell number: ", ix^d
1116 do iw=1,nw
1117 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1118 end do
1119 ! use erroneous arithmetic operation to crash the run
1120 if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1121 write(*,*) "Saving status at the previous time step"
1122 crash=.true.
1123 end if
1124 {end do^D&\}
1125 end if
1126 end subroutine ffhd_get_pthermal_origin
1127
1128 subroutine ffhd_get_temperature_from_te(w, x, ixI^L, ixO^L, res)
1130 integer, intent(in) :: ixi^l, ixo^l
1131 double precision, intent(in) :: w(ixi^s, 1:nw)
1132 double precision, intent(in) :: x(ixi^s, 1:ndim)
1133 double precision, intent(out):: res(ixi^s)
1134
1135 res(ixo^s) = w(ixo^s, te_)
1136 end subroutine ffhd_get_temperature_from_te
1137
1138 subroutine ffhd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1140 integer, intent(in) :: ixi^l, ixo^l
1141 double precision, intent(in) :: w(ixi^s, 1:nw)
1142 double precision, intent(in) :: x(ixi^s, 1:ndim)
1143 double precision, intent(out):: res(ixi^s)
1144 double precision :: r(ixi^s)
1145
1146 call ffhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1147 res(ixo^s) = gamma_1 * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1148 end subroutine ffhd_get_temperature_from_eint
1149
1150 subroutine ffhd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1152 integer, intent(in) :: ixi^l, ixo^l
1153 double precision, intent(in) :: w(ixi^s, 1:nw)
1154 double precision, intent(in) :: x(ixi^s, 1:ndim)
1155 double precision, intent(out):: res(ixi^s)
1156
1157 double precision :: r(ixi^s)
1158
1159 call ffhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1160 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,res)
1161 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1162 end subroutine ffhd_get_temperature_from_etot
1163
1164 subroutine ffhd_get_csound2(w,x,ixI^L,ixO^L,csound2)
1166 integer, intent(in) :: ixi^l, ixo^l
1167 double precision, intent(in) :: w(ixi^s,nw)
1168 double precision, intent(in) :: x(ixi^s,1:ndim)
1169 double precision, intent(out) :: csound2(ixi^s)
1170 double precision :: rho(ixi^s)
1171
1172 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
1173 if(ffhd_energy) then
1174 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,csound2)
1175 csound2(ixo^s)=ffhd_gamma*csound2(ixo^s)/rho(ixo^s)
1176 else
1177 csound2(ixo^s)=ffhd_gamma*ffhd_adiab*rho(ixo^s)**gamma_1
1178 end if
1179 end subroutine ffhd_get_csound2
1180
1181 subroutine ffhd_get_flux(wC,w,x,ixI^L,ixO^L,idim,f)
1183 use mod_geometry
1184 integer, intent(in) :: ixi^l, ixo^l, idim
1185 ! conservative w
1186 double precision, intent(in) :: wc(ixi^s,nw)
1187 ! primitive w
1188 double precision, intent(in) :: w(ixi^s,nw)
1189 double precision, intent(in) :: x(ixi^s,1:ndim)
1190 double precision,intent(out) :: f(ixi^s,nwflux)
1191 double precision :: ptotal(ixo^s)
1192
1193 f(ixo^s,rho_)=w(ixo^s,mom(1))*w(ixo^s,rho_)*block%B0(ixo^s,idim,idim)
1194
1195 if(ffhd_energy) then
1196 ptotal(ixo^s)=w(ixo^s,p_)
1197 else
1198 ptotal(ixo^s)=ffhd_adiab*w(ixo^s,rho_)**ffhd_gamma
1199 end if
1200
1201 ! Get flux of momentum
1202 f(ixo^s,mom(1))=(wc(ixo^s,mom(1))*w(ixo^s,mom(1))+ptotal(ixo^s))*block%B0(ixo^s,idim,idim)
1203
1204 ! Get flux of energy
1205 if(ffhd_energy) then
1206 f(ixo^s,e_)=w(ixo^s,mom(1))*(wc(ixo^s,e_)+ptotal(ixo^s))*block%B0(ixo^s,idim,idim)
1208 f(ixo^s,e_)=f(ixo^s,e_)+w(ixo^s,q_)*block%B0(ixo^s,idim,idim)
1209 f(ixo^s,q_)=zero
1210 end if
1211 end if
1212 end subroutine ffhd_get_flux
1213
1214 subroutine ffhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1218 integer, intent(in) :: ixi^l, ixo^l
1219 double precision, intent(in) :: qdt,dtfactor
1220 double precision, intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw), x(ixi^s,1:ndim)
1221 double precision, intent(inout) :: w(ixi^s,1:nw)
1222 logical, intent(in) :: qsourcesplit
1223 logical, intent(inout) :: active
1224
1225 if (.not. qsourcesplit) then
1226 active = .true.
1227 call add_punitb(qdt,ixi^l,ixo^l,wct,w,x,wctprim)
1229 call add_hypertc_source(qdt,ixi^l,ixo^l,wct,w,x,wctprim)
1230 end if
1231 end if
1232
1233 if(ffhd_radiative_cooling) then
1234 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
1235 w,x,qsourcesplit,active, rc_fl)
1236 end if
1237
1238 if(ffhd_gravity) then
1239 call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
1240 w,x,ffhd_energy,qsourcesplit,active)
1241 end if
1242
1243 ! update temperature from new pressure, density, and old ionization degree
1245 if(.not.qsourcesplit) then
1246 active = .true.
1247 call ffhd_update_temperature(ixi^l,ixo^l,wct,w,x)
1248 end if
1249 end if
1250 end subroutine ffhd_add_source
1251
1252 subroutine add_punitb(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1254 use mod_geometry
1255 integer, intent(in) :: ixi^l,ixo^l
1256 double precision, intent(in) :: qdt
1257 double precision, intent(in) :: wct(ixi^s,1:nw),x(ixi^s,1:ndim)
1258 double precision, intent(in) :: wctprim(ixi^s,1:nw)
1259 double precision, intent(inout) :: w(ixi^s,1:nw)
1260
1261 integer :: idims,hxo^l
1262 double precision :: divb(ixi^s)
1263 double precision :: rhovpar(ixi^s),gradrhov(ixi^s)
1264
1265 divb=zero
1266 if(slab_uniform) then
1267 do idims=1,ndim
1268 hxo^l=ixo^l-kr(idims,^d);
1269 divb(ixo^s)=divb(ixo^s)+(block%B0(ixo^s,idims,idims)-block%B0(hxo^s,idims,idims))/dxlevel(idims)
1270 end do
1271 else
1272 call divvector(block%B0(ixi^s,1:ndir,0),ixi^l,ixo^l,divb)
1273 end if
1274 w(ixo^s,mom(1))=w(ixo^s,mom(1))+qdt*wctprim(ixo^s,p_)*divb(ixo^s)
1275 end subroutine add_punitb
1276
1277 subroutine ffhd_get_rho(w,x,ixI^L,ixO^L,rho)
1279 integer, intent(in) :: ixi^l, ixo^l
1280 double precision, intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:ndim)
1281 double precision, intent(out) :: rho(ixi^s)
1282
1283 rho(ixo^s) = w(ixo^s,rho_)
1284 end subroutine ffhd_get_rho
1285
1286 subroutine ffhd_handle_small_ei(w, x, ixI^L, ixO^L, ie, subname)
1289 integer, intent(in) :: ixi^l,ixo^l, ie
1290 double precision, intent(inout) :: w(ixi^s,1:nw)
1291 double precision, intent(in) :: x(ixi^s,1:ndim)
1292 character(len=*), intent(in) :: subname
1293 integer :: idir
1294 logical :: flag(ixi^s,1:nw)
1295 double precision :: rho(ixi^s)
1296
1297 flag=.false.
1298 where(w(ixo^s,ie)<small_e) flag(ixo^s,ie)=.true.
1299 if(any(flag(ixo^s,ie))) then
1300 select case (small_values_method)
1301 case ("replace")
1302 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e
1303 case ("average")
1304 call small_values_average(ixi^l, ixo^l, w, x, flag, ie)
1305 case default
1306 w(ixo^s,e_)=w(ixo^s,e_)*gamma_1
1307 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
1308 w(ixo^s,mom(1)) = w(ixo^s,mom(1))/rho(ixo^s)
1309 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1310 end select
1311 end if
1312 end subroutine ffhd_handle_small_ei
1313
1314 subroutine ffhd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1317 integer, intent(in) :: ixi^l, ixo^l
1318 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1319 double precision, intent(inout) :: w(ixi^s,1:nw)
1320 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1321
1322 call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1323 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1324 w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1325 he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1326 end subroutine ffhd_update_temperature
1327
1328 subroutine ffhd_get_dt(wprim,ixI^L,ixO^L,dtnew,dx^D,x)
1330 use mod_usr_methods
1331 use mod_gravity, only: gravity_get_dt
1332 integer, intent(in) :: ixi^l, ixo^l
1333 double precision, intent(inout) :: dtnew
1334 double precision, intent(in) :: dx^d
1335 double precision, intent(in) :: wprim(ixi^s,1:nw)
1336 double precision, intent(in) :: x(ixi^s,1:ndim)
1337
1338 dtnew = bigdouble
1339
1340 if(ffhd_gravity) then
1341 call gravity_get_dt(wprim,ixi^l,ixo^l,dtnew,dx^d,x)
1342 end if
1343 end subroutine ffhd_get_dt
1344
1345 subroutine ffhd_add_source_geom(qdt,dtfactor,ixI^L,ixO^L,wCT,wprim,w,x)
1347 integer, intent(in) :: ixi^l, ixo^l
1348 double precision, intent(in) :: qdt, dtfactor,x(ixi^s,1:ndim)
1349 double precision, intent(inout) :: wct(ixi^s,1:nw), wprim(ixi^s,1:nw),w(ixi^s,1:nw)
1350
1351 ! no geometric source terms needed for ffhd, no divergences of tensors
1352 end subroutine ffhd_add_source_geom
1353
1354 function ffhd_kin_en_origin(w, ixI^L, ixO^L, inv_rho) result(ke)
1355 use mod_global_parameters, only: nw, ndim,block
1356 integer, intent(in) :: ixi^l, ixo^l
1357 double precision, intent(in) :: w(ixi^s, nw)
1358 double precision :: ke(ixo^s)
1359 double precision, intent(in), optional :: inv_rho(ixo^s)
1360
1361 if(present(inv_rho)) then
1362 ke(ixo^s)=0.5d0*w(ixo^s,mom(1))**2*inv_rho(ixo^s)
1363 else
1364 ke(ixo^s)=0.5d0*w(ixo^s,mom(1))**2/w(ixo^s,rho_)
1365 end if
1366 end function ffhd_kin_en_origin
1367
1368 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1371 integer, intent(in) :: ixi^l, ixo^l
1372 double precision, intent(in) :: w(ixi^s,1:nw)
1373 double precision, intent(in) :: x(ixi^s,1:ndim)
1374 double precision, intent(out):: rfactor(ixi^s)
1375 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1376
1377 call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1378 rfactor(ixo^s)=(1.d0+iz_h(ixo^s)+0.1d0*(1.d0+iz_he(ixo^s)*(1.d0+iz_he(ixo^s))))/(2.d0+3.d0*he_abundance)
1379 end subroutine rfactor_from_temperature_ionization
1380
1381 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1383 integer, intent(in) :: ixi^l, ixo^l
1384 double precision, intent(in) :: w(ixi^s,1:nw)
1385 double precision, intent(in) :: x(ixi^s,1:ndim)
1386 double precision, intent(out):: rfactor(ixi^s)
1387
1388 rfactor(ixo^s)=rr
1389 end subroutine rfactor_from_constant_ionization
1390
1391 subroutine add_hypertc_source(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1393
1394 integer, intent(in) :: ixi^l,ixo^l
1395 double precision, intent(in) :: qdt
1396 double precision, dimension(ixI^S,1:ndim), intent(in) :: x
1397 double precision, dimension(ixI^S,1:nw), intent(in) :: wct,wctprim
1398 double precision, dimension(ixI^S,1:nw), intent(inout) :: w
1399
1400 double precision, dimension(ixI^S) :: te
1401 double precision :: sigma_t5,sigma_t7,sigmat5_bgradt,f_sat,tau
1402 integer :: ix^d
1403
1404 te(ixi^s)=wctprim(ixi^s,p_)/wct(ixi^s,rho_)
1405 ! temperature on face T_(i+1/2)=(7(T_i+T_(i+1))-(T_(i-1)+T_(i+2)))/12
1406 ! T_(i+1/2)-T_(i-1/2)=(8(T_(i+1)-T_(i-1))-T_(i+2)+T_(i-2))/12
1407 {^iftwod
1408 do ix2=ixomin2,ixomax2
1409 do ix1=ixomin1,ixomax1
1410 if(ffhd_trac) then
1411 if(te(ix^d)<block%wextra(ix^d,tcoff_)) then
1412 sigma_t5=hypertc_kappa*sqrt(block%wextra(ix^d,tcoff_)**5)
1413 sigma_t7=sigma_t5*block%wextra(ix^d,tcoff_)
1414 else
1415 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1416 sigma_t7=sigma_t5*te(ix^d)
1417 end if
1418 else
1419 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1420 sigma_t7=sigma_t5*te(ix^d)
1421 end if
1422 sigmat5_bgradt=sigma_t5*(&
1423 block%B0(ix^d,1,0)*((8.d0*(te(ix1+1,ix2)-te(ix1-1,ix2))-te(ix1+2,ix2)+te(ix1-2,ix2))/12.d0)/block%ds(ix^d,1)&
1424 +block%B0(ix^d,2,0)*((8.d0*(te(ix1,ix2+1)-te(ix1,ix2-1))-te(ix1,ix2+2)+te(ix1,ix2-2))/12.d0)/block%ds(ix^d,2))
1425 if(ffhd_htc_sat) then
1426 ! 5 phi rho c^3, phi=0.3, c=sqrt(p/rho) isothermal sound speed
1427 f_sat=one/(one+dabs(sigmat5_bgradt)/(1.5d0*wct(ix^d,rho_)*(wctprim(ix^d,p_)/wct(ix^d,rho_))**1.5d0))
1428 tau=max(4.d0*dt, f_sat*sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1429 w(ix^d,q_)=w(ix^d,q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^d,q_))/tau
1430 else
1431 w(ix^d,q_)=w(ix^d,q_)-qdt*(sigmat5_bgradt+wct(ix^d,q_))/&
1432 max(4.d0*dt, sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1433 end if
1434 end do
1435 end do
1436 }
1437 {^ifthreed
1438 do ix3=ixomin3,ixomax3
1439 do ix2=ixomin2,ixomax2
1440 do ix1=ixomin1,ixomax1
1441 if(ffhd_trac) then
1442 if(te(ix^d)<block%wextra(ix^d,tcoff_)) then
1443 sigma_t5=hypertc_kappa*sqrt(block%wextra(ix^d,tcoff_)**5)
1444 sigma_t7=sigma_t5*block%wextra(ix^d,tcoff_)
1445 else
1446 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1447 sigma_t7=sigma_t5*te(ix^d)
1448 end if
1449 else
1450 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1451 sigma_t7=sigma_t5*te(ix^d)
1452 end if
1453 sigmat5_bgradt=sigma_t5*(&
1454 block%B0(ix^d,1,0)*((8.d0*(te(ix1+1,ix2,ix3)-te(ix1-1,ix2,ix3))-te(ix1+2,ix2,ix3)+te(ix1-2,ix2,ix3))/12.d0)/block%ds(ix^d,1)&
1455 +block%B0(ix^d,2,0)*((8.d0*(te(ix1,ix2+1,ix3)-te(ix1,ix2-1,ix3))-te(ix1,ix2+2,ix3)+te(ix1,ix2-2,ix3))/12.d0)/block%ds(ix^d,2)&
1456 +block%B0(ix^d,3,0)*((8.d0*(te(ix1,ix2,ix3+1)-te(ix1,ix2,ix3-1))-te(ix1,ix2,ix3+2)+te(ix1,ix2,ix3-2))/12.d0)/block%ds(ix^d,3))
1457 if(ffhd_htc_sat) then
1458 ! 5 phi rho c^3, phi=0.3, c=sqrt(p/rho) isothermal sound speed
1459 f_sat=one/(one+dabs(sigmat5_bgradt)/(1.5d0*wct(ix^d,rho_)*(wctprim(ix^d,p_)/wct(ix^d,rho_))**1.5d0))
1460 tau=max(4.d0*dt, f_sat*sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1461 w(ix^d,q_)=w(ix^d,q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^d,q_))/tau
1462 else
1463 w(ix^d,q_)=w(ix^d,q_)-qdt*(sigmat5_bgradt+wct(ix^d,q_))/&
1464 max(4.d0*dt, sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1465 end if
1466 end do
1467 end do
1468 end do
1469 }
1470 end subroutine add_hypertc_source
1471
1472end module mod_ffhd_phys
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.
subroutine add_convert_method(phys_convert_vars, nwc, dataset_names, file_suffix)
Definition mod_convert.t:59
Frozen-field hydrodynamics module.
integer, public, protected te_
Indices of temperature.
integer, public, protected ffhd_trac_type
Which TRAC method is used.
double precision, public hypertc_kappa
The thermal conductivity kappa in hyperbolic thermal conduction.
integer, public, protected e_
Index of the energy density (-1 if not present)
double precision, public ffhd_gamma
The adiabatic index.
logical, public, protected eq_state_units
double precision, public ffhd_adiab
The adiabatic constant.
procedure(sub_get_pthermal), pointer, public ffhd_get_temperature
logical, public, protected ffhd_hyperbolic_thermal_conduction
Whether hyperbolic type thermal conduction is used.
type(rc_fluid), allocatable, public rc_fl
type of fluid for radiative cooling
integer, public, protected rho_
Index of the density (in the w array)
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
procedure(sub_get_pthermal), pointer, public ffhd_get_pthermal
logical, public, protected ffhd_htc_sat
Wheterh saturation is considered for hyperbolic TC.
subroutine, public ffhd_get_rho(w, x, ixil, ixol, rho)
logical, public, protected ffhd_partial_ionization
Whether plasma is partially ionized.
procedure(sub_convert), pointer, public ffhd_to_conserved
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
double precision, public, protected h_ion_fr
Ionization fraction of H H_ion_fr = H+/(H+ + H)
procedure(fun_kin_en), pointer, public ffhd_kin_en
subroutine, public ffhd_phys_init()
procedure(sub_get_v), pointer, public ffhd_get_v
subroutine, public ffhd_get_csound2(w, x, ixil, ixol, csound2)
logical, public, protected ffhd_energy
Whether an energy equation is used.
double precision, public, protected rr
type(tc_fluid), allocatable, public tc_fl
type of fluid for thermal conduction
type(te_fluid), allocatable, public te_fl_ffhd
type of fluid for thermal emission synthesis
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
logical, public, protected ffhd_radiative_cooling
Whether radiative cooling is added.
subroutine, public ffhd_get_v_idim(w, x, ixil, ixol, idim, v)
integer, public, protected q_
procedure(sub_convert), pointer, public ffhd_to_primitive
integer, public, protected tweight_
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
logical, public, protected ffhd_trac
Whether TRAC method is used.
logical, public, protected ffhd_thermal_conduction
Whether thermal conduction is used.
subroutine, public ffhd_ei_to_e(ixil, ixol, w, x)
logical, public, protected ffhd_gravity
Whether gravity is added.
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
double precision, public, protected ffhd_trac_mask
Height of the mask used in the TRAC method.
double precision, public, protected he_ion_fr
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
integer, public, protected ffhd_trac_finegrid
Distance between two adjacent traced magnetic field lines (in finest cell size)
subroutine, public ffhd_e_to_ei(ixil, ixol, w, x)
Module for flux conservation near refinement boundaries.
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
subroutine divvector(qvec, ixil, ixol, divq, nth_in)
integer coordinate
Definition mod_geometry.t:7
subroutine gradient(q, ixil, ixol, idir, gradq, nth_in)
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.
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 unit_mass
Physical scaling factor for mass.
integer, dimension(3, 3) kr
Kronecker delta tensor.
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.
double precision cmax_global
global fastest wave speed needed in fd scheme and glm method
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
double precision dt
global time step
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision courantpar
The Courant (CFL) number used for the simulation.
logical slab
Cartesian geometry or not.
integer nwauxio
Number of auxiliary variables that are only included in output.
double precision unit_velocity
Physical scaling factor for velocity.
logical b0field
split magnetic field as background B0 field
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
integer nghostcells
Number of ghost cells surrounding a grid.
logical phys_trac
Use TRAC for MHD or 1D HD.
logical need_global_cmax
need global maximal wave speed
logical fix_small_values
fix small values with average or replace methods
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
integer, parameter unitconvert
integer number_equi_vars
number of equilibrium set variables, besides the mag field
Module for including gravity in (magneto)hydrodynamics simulations.
Definition mod_gravity.t:2
subroutine gravity_get_dt(wprim, 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
module ionization degree - get ionization degree for given temperature
subroutine ionization_degree_from_temperature(ixil, ixol, te, iz_h, iz_he)
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
subroutine radiative_cooling_init_params(phys_gamma, he_abund)
Radiative cooling initialization.
subroutine radiative_cooling_init(fl, read_params)
subroutine radiative_cooling_add_source(qdt, ixil, ixol, wct, wctprim, w, x, qsourcesplit, active, fl)
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 which can be used for multiple source terms in the governing equatio...
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_mhd(w, ixil, ixol, dxd, x, fl)
Get the explicit timestep for the TC (mhd implementation) Note: for multi-D MHD (1D MHD will use HD f...
subroutine tc_init_params(phys_gamma)
subroutine, public sts_set_source_tc_mhd(ixil, ixol, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
anisotropic thermal conduction with slope limited symmetric scheme Sharma 2007 Journal of Computation...
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_equi_vars), pointer usr_set_equi_vars