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