MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_srhd_phys.t
Go to the documentation of this file.
1 !> Special Relativistic Hydrodynamics (with EOS) physics module
3  use mod_physics
4  use mod_constants
5  use mod_comm_lib, only: mpistop
6  implicit none
7  private
8 
9  !> Whether particles module is added
10  logical, public, protected :: srhd_particles = .false.
11 
12  !> Number of tracer species
13  integer, public,protected :: srhd_n_tracer = 0
14 
15  !> Index of the density (in the w array)
16  integer, public,protected :: rho_
17  integer, public,protected :: d_
18 
19  !> Indices of the momentum density
20  integer, allocatable, public, protected :: mom(:)
21 
22  !> Indices of the tracers
23  integer, allocatable, public, protected :: tracer(:)
24 
25  !> Index of the energy density
26  integer, public,protected :: e_
27  !> Index of the gas pressure should equal e_
28  integer, public,protected :: p_
29 
30  !> Index of the Lorentz factor
31  integer, public,protected :: lfac_
32 
33  !> Index of the inertia
34  integer, public,protected :: xi_
35 
36  !> Whether synge eos is used
37  logical, public :: srhd_eos = .false.
38 
39  !> The adiabatic index and derived values
40  double precision, public :: srhd_gamma = 5.d0/3.0d0
41  double precision, public :: gamma_1,inv_gamma_1,gamma_to_gamma_1
42 
43  !> The smallest allowed energy
44  double precision, public :: small_e
45  !> The smallest allowed inertia
46  double precision, public :: small_xi
47 
48  !> Allows overruling default corner filling (for debug mode, otherwise corner primitives fail)
49  logical, public, protected :: srhd_force_diagonal = .false.
50 
51  !> Helium abundance over Hydrogen
52  double precision, public, protected :: he_abundance=0.0d0
53 
54  !> parameters for NR in con2prim
55  integer, public :: maxitnr = 100
56  double precision, public :: absaccnr = 1.0d-8
57  double precision, public :: tolernr = 1.0d-9
58  double precision, public :: dmaxvel = 1.0d-7
59  double precision, public :: lfacmax
60  double precision, public :: minp, minrho, smalltau, smallxi
61 
62  ! Public methods
63  public :: srhd_phys_init
64  public :: srhd_get_pthermal
65  public :: srhd_get_auxiliary
66  public :: srhd_get_auxiliary_prim
67  public :: srhd_get_csound2
68  public :: srhd_to_conserved
69  public :: srhd_to_primitive
70  public :: srhd_check_params
71  public :: srhd_check_w
72  public :: srhd_get_geff_eos
73  public :: srhd_get_enthalpy_eos
74 contains
75 
76  !> Read this module's parameters from a file
77  subroutine srhd_read_params(files)
79  character(len=*), intent(in) :: files(:)
80  integer :: n
81 
82  namelist /srhd_list/ srhd_n_tracer, srhd_eos, srhd_gamma, &
85 
86  do n = 1, size(files)
87  open(unitpar, file=trim(files(n)), status="old")
88  read(unitpar, srhd_list, end=111)
89 111 close(unitpar)
90  end do
91 
92  end subroutine srhd_read_params
93 
94  !> Write this module's parameters to a snapshot
95  subroutine srhd_write_info(fh)
97  integer, intent(in) :: fh
98  integer, parameter :: n_par = 1
99  double precision :: values(n_par)
100  character(len=name_len) :: names(n_par)
101  integer, dimension(MPI_STATUS_SIZE) :: st
102  integer :: er
103 
104  call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
105 
106  names(1) = "gamma"
107  values(1) = srhd_gamma
108  call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
109  call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
110 
111  end subroutine srhd_write_info
112 
113  !> Initialize the module
114  subroutine srhd_phys_init()
116  use mod_particles, only: particles_init
117  integer :: itr,idir
118 
119  call srhd_read_params(par_files)
120 
121  physics_type = "srhd"
122  phys_energy = .true.
123  phys_total_energy = .true.
124  phys_gamma = srhd_gamma
125 
126  ! unused physics options
127  phys_internal_e = .false.
128  phys_partial_ionization=.false.
129  phys_trac=.false.
130 
132 
133  ! note: number_species is 1 for srhd
134  allocate(start_indices(number_species),stop_indices(number_species))
135  ! set the index of the first flux variable for species 1
136  start_indices(1)=1
137 
138  ! Determine flux variables
139  rho_ = var_set_rho()
140  d_=rho_
141 
142  allocate(mom(ndir))
143  mom(:) = var_set_momentum(ndir)
144 
145  ! Set index of energy variable
146  e_ = var_set_energy()
147  p_ = e_
148 
149  ! Whether diagonal ghost cells are required for the physics
150  phys_req_diagonal = .false.
151 
152  ! derive units from basic units
153  call srhd_physical_units()
154 
155  if (srhd_force_diagonal) then
156  ! ensure corners are filled, otherwise divide by zero when getting primitives
157  ! --> only for debug purposes
158  phys_req_diagonal = .true.
159  endif
160 
161  allocate(tracer(srhd_n_tracer))
162 
163  ! Set starting index of tracers
164  do itr = 1, srhd_n_tracer
165  tracer(itr) = var_set_fluxvar("trc", "trp", itr, need_bc=.false.)
166  end do
167 
168  ! Set index for auxiliary variables
169  ! MUST be after the possible tracers (which have fluxes)
170  xi_ = var_set_auxvar('xi','xi')
171  lfac_= var_set_auxvar('lfac','lfac')
172 
173  ! set number of variables which need update ghostcells
174  nwgc=nwflux+nwaux
175 
176  ! set the index of the last flux variable for species 1
177  stop_indices(1)=nwflux
178 
179  ! Number of variables need reconstruction in w
180  nw_recon=nwflux
181 
182  ! Check whether custom flux types have been defined
183  if (.not. allocated(flux_type)) then
184  allocate(flux_type(ndir, nw))
185  flux_type = flux_default
186  else if (any(shape(flux_type) /= [ndir, nw])) then
187  call mpistop("phys_check error: flux_type has wrong shape")
188  end if
189 
190  nvector = 1 ! No. vector vars
191  allocate(iw_vector(nvector))
192  iw_vector(1) = mom(1) - 1
193 
194  ! dummy for now, no extra source terms precoded
195  phys_add_source => srhd_add_source
196  phys_get_dt => srhd_get_dt
197  ! copied in from HD/MHD, for certain limiters
198  phys_get_a2max => srhd_get_a2max
199 
200  ! actual srhd routines
201  phys_check_params => srhd_check_params
202  phys_check_w => srhd_check_w
203  phys_get_cmax => srhd_get_cmax
204  phys_get_cbounds => srhd_get_cbounds
205  phys_get_flux => srhd_get_flux
206  phys_add_source_geom => srhd_add_source_geom
207  phys_to_conserved => srhd_to_conserved
208  phys_to_primitive => srhd_to_primitive
209  phys_get_pthermal => srhd_get_pthermal
210  phys_get_auxiliary => srhd_get_auxiliary
211  phys_get_auxiliary_prim => srhd_get_auxiliary_prim
212  phys_get_pthermal => srhd_get_pthermal
213  phys_get_v => srhd_get_v
214  phys_write_info => srhd_write_info
215  phys_handle_small_values => srhd_handle_small_values
216 
217  ! Initialize particles module
218  if (srhd_particles) then
219  call particles_init()
220  phys_req_diagonal = .true.
221  end if
222 
223  end subroutine srhd_phys_init
224 
225  subroutine srhd_check_params
227 
228  if (srhd_gamma <= 0.0d0 .or. srhd_gamma == 1.0d0) &
229  call mpistop ("Error: srhd_gamma <= 0 or srhd_gamma == 1")
230  ! additional useful values
231  gamma_1=srhd_gamma-1.0d0
232  inv_gamma_1=1.0d0/gamma_1
234 
235  ! the following sets small_e and small_xi from small_density/small_pressure
236  ! according to the srhd_eos used
238 
239  if(mype==0)then
240  write(*,*)'------------------------------------------------------------'
241  write(*,*)'Using EOS set via srhd_eos=',srhd_eos
242  write(*,*)'Maximal lorentz factor (via dmaxvel) is=',lfacmax
243  write(*,*)'Use fixes set through check/fix small values:', check_small_values,fix_small_values
244  write(*,*)'Controlled with small pressure/density:', small_pressure,small_density
245  write(*,*)'Derived small values: xi and e ',small_xi,small_e
246  write(*,*)'------------------------------------------------------------'
247  endif
248 
249  end subroutine srhd_check_params
250 
253  double precision :: mp,kB
254 
255  ! Derive scaling units
256  if(si_unit) then
257  mp=mp_si
258  kb=kb_si
259  unit_velocity=c_si
260  else
261  mp=mp_cgs
262  kb=kb_cgs
263  unit_velocity=const_c
264  end if
265  if(unit_numberdensity*unit_length<=0.0d0)then
266  call mpistop("Abort: must set positive values for unit length and numberdensity")
267  endif
268  ! we assume user sets: unit_numberdensity, unit_length, He_abundance
269  ! then together with light speed c, all units fixed
275 
276  end subroutine srhd_physical_units
277 
278  !> Returns logical argument flag T where values are not ok
279  subroutine srhd_check_w(primitive, ixI^L, ixO^L, w, flag)
281  logical, intent(in) :: primitive
282  integer, intent(in) :: ixi^l, ixo^l
283  double precision, intent(in) :: w(ixi^s, nw)
284  logical, intent(inout) :: flag(ixi^s,1:nw)
285 
286  flag=.false.
287 
288  ! NOTE: we should not check or use nwaux variables here
289  if(primitive) then
290  where(w(ixo^s,rho_) < small_density) flag(ixo^s,rho_) = .true.
291  where(w(ixo^s,p_) < small_pressure) flag(ixo^s,p_) = .true.
292  else
293  where(w(ixo^s,d_) < small_density) flag(ixo^s,d_) = .true.
294  where(w(ixo^s,e_) < small_e) flag(ixo^s,e_) = .true.
295  endif
296 
297  end subroutine srhd_check_w
298 
299  !> Returns logical argument flag T where auxiliary values are not ok
300  subroutine srhd_check_w_aux(ixI^L, ixO^L, w, flag)
302  integer, intent(in) :: ixI^L, ixO^L
303  double precision, intent(in) :: w(ixI^S, nw)
304  logical, intent(inout) :: flag(ixI^S,1:nw)
305 
306  flag=.false.
307 
308  where(w(ixo^s,xi_) < small_xi) flag(ixo^s,xi_) = .true.
309  where(w(ixo^s,lfac_) < one) flag(ixo^s,lfac_) = .true.
310 
311  if(any(flag(ixo^s,xi_)))then
312  write(*,*)'auxiliary xi too low: abort'
313  call mpistop('auxiliary check failed')
314  endif
315  if(any(flag(ixo^s,lfac_)))then
316  write(*,*)'auxiliary lfac too low: abort'
317  call mpistop('auxiliary check failed')
318  endif
319 
320  end subroutine srhd_check_w_aux
321 
322  !> Set auxiliary variables lfac and xi from a primitive state
323  !> only used when handle_small_values average on primitives
324  subroutine srhd_get_auxiliary_prim(ixI^L,ixO^L,w)
326  integer, intent(in) :: ixi^l, ixo^l
327  double precision, intent(inout) :: w(ixi^s, nw)
328  double precision, dimension(ixO^S) :: rho,rhoh,pth
329 
330  ! assume four-velocity in momentum vector (i.e. lfac*v)
331  rho(ixo^s) = sum(w(ixo^s, mom(:))**2, dim=ndim+1)
332  w(ixo^s,lfac_) = dsqrt(1.0d0+rho(ixo^s))
333 
334  rho(ixo^s)=w(ixo^s,rho_)
335  pth(ixo^s)=w(ixo^s,p_)
336  ! compute rho*h (enthalpy h) from density-pressure
337  call srhd_get_enthalpy_eos(ixo^l,rho,pth,rhoh)
338 
339  ! fill auxiliary variable xi= lfac^2 rhoh
340  w(ixo^s,xi_) = w(ixo^s,lfac_)**2.0d0*rhoh(ixo^s)
341 
342  end subroutine srhd_get_auxiliary_prim
343 
344  !> Compute auxiliary variables lfac and xi from a conservative state
345  !> using srhd_con2prim to calculate enthalpy and lorentz factor
346  subroutine srhd_get_auxiliary(ixI^L,ixO^L,w,x)
348  implicit none
349 
350  integer, intent(in) :: ixi^l,ixo^l
351  double precision, intent(inout) :: w(ixi^s, nw)
352  double precision, intent(in) :: x(ixi^s, 1:ndim)
353 
354  integer :: ix^d,ierror,idir
355  integer :: flag_error(ixo^s)
356  double precision :: ssqr
357 
358  if(srhd_eos)then
359  {do ix^db=ixomin^db,ixomax^db\}
360  ierror=0
361  ssqr=0.0d0
362  do idir=1,ndir
363  ssqr= ssqr+w(ix^d,mom(idir))**2
364  enddo
365  if(w(ix^d,d_)<small_density)then
366  print *,'entering con2prim with', w(ix^d,lfac_),w(ix^d,xi_), &
367  w(ix^d,d_),ssqr,w(ix^d,e_)
368  print *,'in position',ix^d,x(ix^d,1:ndim)
370  if(check_small_values) call mpistop('small density on entry con2prim')
372  endif
373  if(w(ix^d,e_)<small_e)then
374  print *,'entering con2prim with', w(ix^d,lfac_),w(ix^d,xi_), &
375  w(ix^d,d_),ssqr,w(ix^d,e_)
376  print *,'in position',ix^d,x(ix^d,1:ndim)
378  if(check_small_values) call mpistop('small energy on entry con2prim')
379  if(fix_small_values) w(ix^d,e_)=small_e
380  endif
381  call con2prim_eos(w(ix^d,lfac_),w(ix^d,xi_), &
382  w(ix^d,d_),ssqr,w(ix^d,e_),ierror)
383  flag_error(ix^d) = ierror
384  {enddo^d&\}
385  else
386  {do ix^db=ixomin^db,ixomax^db\}
387  ierror=0
388  ssqr=0.0d0
389  do idir=1,ndir
390  ssqr= ssqr+w(ix^d,mom(idir))**2
391  enddo
392  if(w(ix^d,d_)<small_density)then
393  print *,'entering con2prim with', w(ix^d,lfac_),w(ix^d,xi_), &
394  w(ix^d,d_),ssqr,w(ix^d,e_)
395  print *,'in position',ix^d,x(ix^d,1:ndim)
397  if(check_small_values) call mpistop('small density on entry con2prim')
399  endif
400  if(w(ix^d,e_)<small_e)then
401  print *,'entering con2prim with', w(ix^d,lfac_),w(ix^d,xi_), &
402  w(ix^d,d_),ssqr,w(ix^d,e_)
403  print *,'in position',ix^d,x(ix^d,1:ndim)
405  if(check_small_values) call mpistop('small energy on entry con2prim')
406  if(fix_small_values) w(ix^d,e_)=small_e
407  endif
408  call con2prim(w(ix^d,lfac_),w(ix^d,xi_), &
409  w(ix^d,d_),ssqr,w(ix^d,e_),ierror)
410  flag_error(ix^d) = ierror
411  {enddo^d&\}
412  endif
413 
414  if(check_small_values)then
415  if(any(flag_error(ixo^s)/=0))then
416  print *,flag_error(ixo^s)
417  call mpistop('Problem when getting auxiliaries')
418  ! call srhd_handle_small_values(.false.,w,x,ixI^L,ixO^L,'srhd_get_auxiliary')
419  end if
420  end if
421 
422  end subroutine srhd_get_auxiliary
423 
424  !> Transform primitive variables into conservative ones
425  subroutine srhd_to_conserved(ixI^L, ixO^L, w, x)
427  integer, intent(in) :: ixi^l, ixo^l
428  double precision, intent(inout) :: w(ixi^s, nw)
429  double precision, intent(in) :: x(ixi^s, 1:ndim)
430  integer :: idir,itr
431  double precision, dimension(ixO^S) :: rhoh,rho,pth
432 
433  ! assume four-velocity in momentum vector (i.e. lfac*v)
434  ! use rhoh slot for temporary array
435  rhoh(ixo^s) = sum(w(ixo^s, mom(:))**2, dim=ndim+1)
436  w(ixo^s,lfac_) = dsqrt(1.0d0+rhoh(ixo^s))
437 
438  rho(ixo^s)=w(ixo^s,rho_)
439  pth(ixo^s)=w(ixo^s,p_)
440  ! compute rho*h (enthalpy h) from density-pressure
441  call srhd_get_enthalpy_eos(ixo^l,rho,pth,rhoh)
442 
443  ! compute rhoh*lfac (recycle rhoh)
444  rhoh(ixo^s)= rhoh(ixo^s)*w(ixo^s,lfac_)
445  ! fill auxiliary variable xi= lfac^2 rhoh
446  w(ixo^s,xi_) = w(ixo^s,lfac_)*rhoh(ixo^s)
447 
448  ! set conservative density: d = lfac * rho
449  w(ixo^s,d_)=w(ixo^s,lfac_)*rho(ixo^s)
450 
451  ! Convert four-velocity (lfac*v) to momentum (xi*v=[rho*h*lfac^2]*v)
452  do idir = 1, ndir
453  w(ixo^s, mom(idir)) = rhoh(ixo^s)*w(ixo^s, mom(idir))
454  end do
455 
456  ! set tau = xi-p-d energy variable
457  w(ixo^s,e_) = w(ixo^s,xi_)-pth(ixo^s)-w(ixo^s,d_)
458 
459  do itr=1,srhd_n_tracer
460  w(ixo^s,tracer(itr)) = w(ixo^s,d_)*w(ixo^s,tracer(itr))
461  end do
462 
463  end subroutine srhd_to_conserved
464 
465  !> Transform conservative variables into primitive ones
466  subroutine srhd_to_primitive(ixI^L, ixO^L, w, x)
468  integer, intent(in) :: ixi^l, ixo^l
469  double precision, intent(inout) :: w(ixi^s, nw)
470  double precision, intent(in) :: x(ixi^s, 1:ndim)
471  integer :: idir,itr
472  double precision, dimension(ixO^S) :: rho,rhoh,e
473  double precision, dimension(ixI^S) :: pth
474  character(len=30) :: subname_loc
475 
476  ! get auxiliary variables lfac and xi from conserved set
477  call srhd_get_auxiliary(ixi^l,ixo^l,w,x)
478 
479  ! from d to rho (d=rho*lfac)
480  rho(ixo^s) = w(ixo^s,d_)/w(ixo^s,lfac_)
481 
482  ! compute pressure
483  ! deduce rho*h from xi/lfac^2
484  rhoh(ixo^s) = w(ixo^s,xi_)/w(ixo^s,lfac_)**2.0d0
485  call srhd_get_pressure_eos(ixi^l,ixo^l,rho,rhoh,pth,e)
486 
487  w(ixo^s,rho_)=rho(ixo^s)
488  ! from xi*v to U=lfac*v (four-velocity)
489  do idir=1,ndir
490  w(ixo^s,mom(idir)) = w(ixo^s,lfac_)*w(ixo^s,mom(idir))&
491  /w(ixo^s,xi_)
492  end do
493  w(ixo^s,p_)=pth(ixo^s)
494 
495  do itr=1,srhd_n_tracer
496  w(ixo^s,tracer(itr)) = w(ixo^s,tracer(itr)) &
497  /(rho(ixo^s)*w(ixo^s,lfac_))
498  end do
499 
500  end subroutine srhd_to_primitive
501 
502  !> Calculate v vector from conservatives
503  subroutine srhd_get_v(w,x,ixI^L,ixO^L,v)
505  integer, intent(in) :: ixI^L, ixO^L
506  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
507  double precision, intent(out) :: v(ixI^S,1:ndir)
508  integer :: idir
509 
510  ! get v from xi*v
511  do idir=1,ndir
512  v(ixo^s,idir) = w(ixo^s, mom(idir))/w(ixo^s,xi_)
513  end do
514 
515  end subroutine srhd_get_v
516 
517  !> Calculate the square of the thermal sound speed csound2 within ixO^L.
518  !> here computed from conservative set WITH ADDED rho*h
519  !> local version: does not do con2prim
520  subroutine srhd_get_csound2_rhoh(w,x,ixI^L,ixO^L,rhoh,csound2)
522  integer, intent(in) :: ixI^L, ixO^L
523  double precision, intent(in) :: w(ixI^S,nw),rhoh(ixO^S)
524  double precision, intent(in) :: x(ixI^S,1:ndim)
525  double precision, intent(out) :: csound2(ixO^S)
526 
527  double precision :: rho(ixO^S)
528 
529  rho=w(ixo^s,d_)/w(ixo^s,lfac_)
530  call srhd_get_csound2_eos(ixi^l,ixo^l,rho,rhoh,csound2)
531 
532  end subroutine srhd_get_csound2_rhoh
533 
534  !> Calculate the square of the thermal sound speed csound2 within ixO^L.
535  !> here computed from conservative set and uses con2prim!!!
536  !> public version!
537  subroutine srhd_get_csound2(w,x,ixI^L,ixO^L,csound2)
539  integer, intent(in) :: ixi^l, ixo^l
540  double precision, intent(inout) :: w(ixi^s,nw)
541  double precision, intent(in) :: x(ixi^s,1:ndim)
542  double precision, intent(out) :: csound2(ixo^s)
543 
544  double precision :: rho(ixo^s),rhoh(ixo^s)
545 
546  ! get auxiliary variables lfac and xi from conserved set
547  call srhd_get_auxiliary(ixi^l,ixo^l,w,x)
548  ! quantify rho, rho*h
549  rho = w(ixo^s,d_)/w(ixo^s,lfac_)
550  rhoh = w(ixo^s,xi_)/w(ixo^s,lfac_)**2.0d0
551  call srhd_get_csound2_eos(ixi^l,ixo^l,rho,rhoh,csound2)
552 
553  end subroutine srhd_get_csound2
554 
555  !> Calculate thermal pressure p within ixO^L
556  !> must follow after update conservative with auxiliaries
557  subroutine srhd_get_pthermal(w, x, ixI^L, ixO^L, pth)
560 
561  integer, intent(in) :: ixi^l, ixo^l
562  double precision, intent(in) :: w(ixi^s, nw)
563  double precision, intent(in) :: x(ixi^s, 1:ndim)
564  double precision, intent(out) :: pth(ixi^s)
565 
566  integer :: iw, ix^d
567  double precision :: rho(ixo^s),rhoh(ixo^s),e(ixo^s)
568 
569  ! quantify rho, rho*h, tau and get pthermal
570  rho = w(ixo^s,d_)/w(ixo^s,lfac_)
571  rhoh = w(ixo^s,xi_)/w(ixo^s,lfac_)**2.0d0
572  call srhd_get_pressure_eos(ixi^l,ixo^l,rho,rhoh,pth,e)
573 
574  end subroutine srhd_get_pthermal
575 
576  !> dummy addsource subroutine
577  ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
578  subroutine srhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
580 
581  integer, intent(in) :: ixI^L, ixO^L
582  double precision, intent(in) :: qdt,dtfactor
583  double precision, intent(in) :: wCT(ixI^S, 1:nw),wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
584  double precision, intent(inout) :: w(ixI^S, 1:nw)
585  logical, intent(in) :: qsourcesplit
586  logical, intent(inout) :: active
587 
588  end subroutine srhd_add_source
589 
590  !> dummy get_dt subroutine
591  subroutine srhd_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
593 
594  integer, intent(in) :: ixI^L, ixO^L
595  double precision, intent(in) :: dx^D, x(ixI^S, 1:^ND)
596  double precision, intent(in) :: w(ixI^S, 1:nw)
597  double precision, intent(inout) :: dtnew
598 
599  dtnew = bigdouble
600 
601  end subroutine srhd_get_dt
602 
603  !> Calculate cmax_idim within ixO^L
604  !> used especially for setdt CFL limit
605  subroutine srhd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
607 
608  integer, intent(in) :: ixI^L, ixO^L, idim
609  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:ndim)
610  double precision, intent(inout) :: cmax(ixI^S)
611 
612  double precision :: wc(ixI^S,nw)
613  double precision, dimension(ixO^S) :: csound2,tmp1,tmp2,v2
614  double precision, dimension(ixI^S) :: vidim, cmin
615 
616  logical :: flag(ixI^S,1:nw)
617 
618  !!call srhd_check_w_aux(ixI^L, ixO^L, w, flag)
619 
620  ! input w is in primitive form TODO use it
621  wc=w
622  call srhd_to_conserved(ixi^l, ixo^l, wc, x)
623  ! auxiliaries are filled here
624  tmp1(ixo^s)=wc(ixo^s,xi_)/wc(ixo^s,lfac_)**2.0d0
625  v2(ixo^s)=1.0d0-1.0d0/wc(ixo^s,lfac_)**2
626  call srhd_get_csound2_rhoh(wc,x,ixi^l,ixo^l,tmp1,csound2)
627  vidim(ixo^s) = wc(ixo^s, mom(idim))/wc(ixo^s, xi_)
628  tmp2(ixo^s)=vidim(ixo^s)**2.0d0
629  tmp1(ixo^s)=1.0d0-v2(ixo^s)*csound2(ixo^s) &
630  -tmp2(ixo^s)*(1.0d0-csound2(ixo^s))
631  tmp2(ixo^s)=dsqrt(csound2(ixo^s)*(one-v2(ixo^s))*tmp1(ixo^s))
632  tmp1(ixo^s)=vidim(ixo^s)*(one-csound2(ixo^s))
633  cmax(ixo^s)=(tmp1(ixo^s)+tmp2(ixo^s))/(one-v2(ixo^s)*csound2(ixo^s))
634  cmin(ixo^s)=(tmp1(ixo^s)-tmp2(ixo^s))/(one-v2(ixo^s)*csound2(ixo^s))
635  ! Limit by speed of light
636  cmin(ixo^s) = max(cmin(ixo^s), - 1.0d0)
637  cmin(ixo^s) = min(cmin(ixo^s), 1.0d0)
638  cmax(ixo^s) = max(cmax(ixo^s), - 1.0d0)
639  cmax(ixo^s) = min(cmax(ixo^s), 1.0d0)
640  ! now take extremal value only for dt limit
641  cmax(ixo^s) = max(dabs(cmax(ixo^s)),dabs(cmin(ixo^s)))
642 
643  end subroutine srhd_get_cmax
644 
645  subroutine srhd_get_a2max(w,x,ixI^L,ixO^L,a2max)
647 
648  integer, intent(in) :: ixI^L, ixO^L
649  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
650  double precision, intent(inout) :: a2max(ndim)
651  double precision :: a2(ixI^S,ndim,nw)
652  integer :: gxO^L,hxO^L,jxO^L,kxO^L,i
653 
654  a2=zero
655  do i = 1,ndim
656  !> 4th order
657  hxo^l=ixo^l-kr(i,^d);
658  gxo^l=hxo^l-kr(i,^d);
659  jxo^l=ixo^l+kr(i,^d);
660  kxo^l=jxo^l+kr(i,^d);
661  a2(ixo^s,i,1:nwflux)=dabs(-w(kxo^s,1:nwflux)+16.d0*w(jxo^s,1:nwflux)&
662  -30.d0*w(ixo^s,1:nwflux)+16.d0*w(hxo^s,1:nwflux)-w(gxo^s,1:nwflux))
663  a2max(i)=maxval(a2(ixo^s,i,1:nwflux))/12.d0/dxlevel(i)**2
664  end do
665 
666  end subroutine srhd_get_a2max
667 
668  !> local version for recycling code when computing cmax-cmin
669  subroutine srhd_get_cmax_loc(ixI^L,ixO^L,vidim,csound2,v2,cmax,cmin)
671  integer, intent(in) :: ixI^L, ixO^L
672  double precision, intent(in) :: vidim(ixI^S)
673  double precision, intent(in), dimension(ixO^S) :: csound2
674  double precision, intent(in) :: v2(ixI^S)
675  double precision, intent(out) :: cmax(ixI^S)
676  double precision, intent(out) :: cmin(ixI^S)
677 
678  double precision, dimension(ixI^S):: tmp1,tmp2
679 
680  tmp2(ixo^s)=vidim(ixo^s)**2.0d0
681  tmp1(ixo^s)=1.0d0-v2(ixo^s)*csound2(ixo^s) &
682  -tmp2(ixo^s)*(1.0d0-csound2(ixo^s))
683  tmp2(ixo^s)=dsqrt(csound2(ixo^s)*(one-v2(ixo^s))*tmp1(ixo^s))
684  tmp1(ixo^s)=vidim(ixo^s)*(one-csound2(ixo^s))
685  cmax(ixo^s)=(tmp1(ixo^s)+tmp2(ixo^s))/(one-v2(ixo^s)*csound2(ixo^s))
686  ! Limit by speed of light
687  cmax(ixo^s) = max(cmax(ixo^s), - 1.0d0)
688  cmax(ixo^s) = min(cmax(ixo^s), 1.0d0)
689  cmin(ixo^s)=(tmp1(ixo^s)-tmp2(ixo^s))/(one-v2(ixo^s)*csound2(ixo^s))
690  ! Limit by speed of light
691  cmin(ixo^s) = max(cmin(ixo^s), - 1.0d0)
692  cmin(ixo^s) = min(cmin(ixo^s), 1.0d0)
693 
694  end subroutine srhd_get_cmax_loc
695 
696  !> Estimating bounds for the minimum and maximum signal velocities
697  !> here we will not use Hspeed at all (one species only)
698  subroutine srhd_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
700  use mod_variables
701 
702  integer, intent(in) :: ixI^L, ixO^L, idim
703  ! conservative left and right status
704  double precision, intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
705  ! primitive left and right status
706  double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
707  double precision, intent(in) :: x(ixI^S, 1:ndim)
708  double precision, intent(inout) :: cmax(ixI^S,1:number_species)
709  double precision, intent(inout), optional :: cmin(ixI^S,1:number_species)
710  double precision, intent(in) :: Hspeed(ixI^S,1:number_species)
711 
712  double precision :: wmean(ixI^S,nw)
713  double precision, dimension(ixO^S) :: csound2,tmp1,tmp2,tmp3
714  double precision, dimension(ixI^S) :: vidim,cmaxL,cmaxR,cminL,cminR,v2
715 
716  logical :: flag(ixI^S,1:nw)
717 
718  select case(boundspeed)
719  case(1) ! we do left-right first and take maximals
720  !!call srhd_check_w(.true.,ixI^L, ixO^L, wLp, flag)
721  !!call srhd_check_w_aux(ixI^L, ixO^L, wLp, flag)
722  tmp1=wlp(ixo^s,rho_)
723  tmp2=wlp(ixo^s,xi_)/wlp(ixo^s,lfac_)**2.0d0
724  tmp3=wlp(ixo^s,p_)
725  call srhd_get_csound2_prim_eos(ixo^l,tmp1,tmp2,tmp3,csound2)
726  vidim(ixo^s) = wlp(ixo^s, mom(idim))/wlp(ixo^s, lfac_)
727  v2(ixo^s) = 1.0d0-1.0d0/wlp(ixo^s,lfac_)**2
728  call srhd_get_cmax_loc(ixi^l,ixo^l,vidim,csound2,v2,cmaxl,cminl)
729 
730  !!call srhd_check_w(.true.,ixI^L, ixO^L, wRp, flag)
731  !!call srhd_check_w_aux(ixI^L, ixO^L, wRp, flag)
732  tmp1=wrp(ixo^s,rho_)
733  tmp2=wrp(ixo^s,xi_)/wrp(ixo^s,lfac_)**2.0d0
734  tmp3=wrp(ixo^s,p_)
735  call srhd_get_csound2_prim_eos(ixo^l,tmp1,tmp2,tmp3,csound2)
736  vidim(ixo^s) = wrp(ixo^s, mom(idim))/wrp(ixo^s, lfac_)
737  v2(ixo^s) = 1.0d0-1.0d0/wrp(ixo^s,lfac_)**2
738  call srhd_get_cmax_loc(ixi^l,ixo^l,vidim,csound2,v2,cmaxr,cminr)
739 
740  if(present(cmin))then
741  ! for HLL
742  cmax(ixo^s,1)=max(cmaxl(ixo^s),cmaxr(ixo^s))
743  cmin(ixo^s,1)=min(cminl(ixo^s),cminr(ixo^s))
744  else
745  ! for TVDLF
746  cmaxl(ixo^s)=max(cmaxl(ixo^s),dabs(cminl(ixo^s)))
747  cmaxr(ixo^s)=max(cmaxr(ixo^s),dabs(cminr(ixo^s)))
748  cmax(ixo^s,1)=max(cmaxl(ixo^s),cmaxr(ixo^s))
749  endif
750  case(2) ! this is cmaxmean from conservatives
751  ! here we do arithmetic mean of conservative vars
752  wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
753  ! get auxiliary variables
754  call srhd_get_auxiliary(ixi^l,ixo^l,wmean,x)
755  ! here tmp1 is rhoh
756  tmp1=wmean(ixo^s,xi_)/wmean(ixo^s,lfac_)**2.0d0
757  call srhd_get_csound2_rhoh(wmean,x,ixi^l,ixo^l,tmp1,csound2)
758  vidim(ixo^s) = wmean(ixo^s, mom(idim))/wmean(ixo^s, xi_)
759  v2(ixo^s)=1.0d0-1.0d0/wmean(ixo^s,lfac_)**2
760  call srhd_get_cmax_loc(ixi^l,ixo^l,vidim,csound2,v2,cmaxl,cminl)
761  if(present(cmin)) then
762  cmax(ixo^s,1)=cmaxl(ixo^s)
763  cmin(ixo^s,1)=cminl(ixo^s)
764  else
765  cmax(ixo^s,1)=max(cmaxl(ixo^s),dabs(cminl(ixo^s)))
766  endif
767  case(3) ! this is cmaxmean from primitives
768  ! here we do arithmetic mean of primitive vars
769  wmean(ixo^s,1:nwflux)=0.5d0*(wlp(ixo^s,1:nwflux)+wrp(ixo^s,1:nwflux))
770  ! get auxiliary variables for wmean (primitive array)
771  call srhd_get_auxiliary_prim(ixi^l,ixo^l,wmean)
772  ! here tmp1 is rhoh
773  tmp1=wmean(ixo^s,rho_)
774  tmp2=wmean(ixo^s,xi_)/wmean(ixo^s,lfac_)**2.0d0
775  tmp3=wmean(ixo^s,p_)
776  call srhd_get_csound2_prim_eos(ixo^l,tmp1,tmp2,tmp3,csound2)
777  vidim(ixo^s) = wmean(ixo^s, mom(idim))/wmean(ixo^s, lfac_)
778  v2(ixo^s) = 1.0d0-1.0d0/wmean(ixo^s,lfac_)**2
779  call srhd_get_cmax_loc(ixi^l,ixo^l,vidim,csound2,v2,cmaxl,cminl)
780  if(present(cmin)) then
781  cmax(ixo^s,1)=cmaxl(ixo^s)
782  cmin(ixo^s,1)=cminl(ixo^s)
783  else
784  cmax(ixo^s,1)=max(cmaxl(ixo^s),dabs(cminl(ixo^s)))
785  endif
786  end select
787 
788  end subroutine srhd_get_cbounds
789 
790  !> Calculate fluxes within ixO^L.
791  subroutine srhd_get_flux(wC,wP,x,ixI^L,ixO^L,idim,f)
793  integer, intent(in) :: ixI^L, ixO^L, idim
794  ! conservative w
795  double precision, intent(in) :: wC(ixI^S,nw)
796  ! primitive w
797  double precision, intent(in) :: wP(ixI^S,nw)
798  double precision, intent(in) :: x(ixI^S,1:ndim)
799  double precision,intent(out) :: f(ixI^S,nwflux)
800 
801  double precision :: pth(ixI^S)
802  double precision :: v(ixI^S,1:ndir)
803  integer :: iw,idir
804 
805  pth(ixo^s)=wp(ixo^s,p_)
806  do idir=1,ndir
807  v(ixo^s,idir) = wp(ixo^s, mom(idir))/wp(ixo^s,lfac_)
808  end do
809 
810  ! Get flux of density d, namely D*v
811  f(ixo^s,d_)=v(ixo^s,idim)*wc(ixo^s,rho_)
812 
813  ! Get flux of tracer
814  do iw=1,srhd_n_tracer
815  f(ixo^s,tracer(iw))=v(ixo^s,idim)*wc(ixo^s,tracer(iw))
816  end do
817 
818  ! Get flux of momentum
819  ! f_i[m_k]=v_i*m_k [+pth if i==k]
820  do idir=1,ndir
821  f(ixo^s,mom(idir))= v(ixo^s,idim)*wc(ixo^s,mom(idir))
822  end do
823  f(ixo^s,mom(idim))=pth(ixo^s)+f(ixo^s,mom(idim))
824 
825  ! Get flux of energy
826  ! f_i[e]=v_i*e+v_i*pth
827  f(ixo^s,e_)=v(ixo^s,idim)*(wc(ixo^s,e_) + pth(ixo^s))
828 
829  end subroutine srhd_get_flux
830 
831  !> Add geometrical source terms to w
832  subroutine srhd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
835  use mod_geometry
836  integer, intent(in) :: ixI^L, ixO^L
837  double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:ndim)
838  double precision, intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S, 1:nw), w(ixI^S, 1:nw)
839 
840  double precision :: pth(ixI^S), source(ixI^S), v(ixI^S,1:ndir)
841  integer :: idir, h1x^L{^NOONED, h2x^L}
842  integer :: mr_,mphi_ ! Polar var. names
843  double precision :: exp_factor(ixI^S), del_exp_factor(ixI^S), exp_factor_primitive(ixI^S)
844 
845  select case (coordinate)
846 
847  case(cartesian_expansion)
848  !the user provides the functions of exp_factor and del_exp_factor
849  if(associated(usr_set_surface)) &
850  call usr_set_surface(ixi^l,x,block%dx,exp_factor,del_exp_factor,exp_factor_primitive)
851  ! get auxiliary variables lfac and xi from conserved set
852  call srhd_get_auxiliary(ixi^l,ixo^l,wct,x)
853  call srhd_get_pthermal(wct, x, ixi^l, ixo^l, source)
854  source(ixo^s) = source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
855  w(ixo^s,mom(1)) = w(ixo^s,mom(1)) + qdt*source(ixo^s)
856 
857  case (cylindrical)
858  mr_ = mom(r_)
859  ! get auxiliary variables lfac and xi from conserved set
860  call srhd_get_auxiliary(ixi^l,ixo^l,wct,x)
861  call srhd_get_pthermal(wct, x, ixi^l, ixo^l, source)
862  if (phi_ > 0) then
863  mphi_ = mom(phi_)
864  source(ixo^s) = source(ixo^s) + wct(ixo^s, mphi_)*wprim(ixo^s,mom(phi_))
865  w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
866  source(ixo^s) = -wct(ixo^s, mphi_) * wprim(ixo^s,mom(r_))
867  w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt * source(ixo^s) / x(ixo^s, r_)
868  else
869  w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
870  end if
871 
872  case (spherical)
873  mr_ = mom(r_)
874  h1x^l=ixo^l-kr(1,^d); {^nooned h2x^l=ixo^l-kr(2,^d);}
875  ! s[mr]=((stheta*vtheta+sphi*vphi)+2*p)/r
876  ! get auxiliary variables lfac and xi from conserved set
877  call srhd_get_auxiliary(ixi^l,ixo^l,wct,x)
878  call srhd_get_pthermal(wct, x, ixi^l, ixo^l, pth)
879  source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
880  *(block%surfaceC(ixo^s, 1) - block%surfaceC(h1x^s, 1)) &
881  /block%dvolume(ixo^s)
882  if (ndir > 1) then
883  do idir = 2, ndir
884  source(ixo^s) = source(ixo^s) + wct(ixo^s, mom(idir))*wprim(ixo^s,mom(idir))
885  end do
886  end if
887  w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, 1)
888 
889  {^nooned
890  ! s[mtheta]=-(stheta*vr)/r+cot(theta)*(sphi*vphi+p)/r
891  source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
892  * (block%surfaceC(ixo^s, 2) - block%surfaceC(h2x^s, 2)) &
893  / block%dvolume(ixo^s)
894  if (ndir == 3) then
895  source(ixo^s) = source(ixo^s) + (wct(ixo^s, mom(3))*v(ixo^s,ndir)) / dtan(x(ixo^s, 2))
896  end if
897  source(ixo^s) = source(ixo^s) - wct(ixo^s, mom(2)) * wprim(ixo^s, mom(1))
898  w(ixo^s, mom(2)) = w(ixo^s, mom(2)) + qdt * source(ixo^s) / x(ixo^s, 1)
899 
900  if (ndir == 3) then
901  ! s[mphi]=-(sphi*vr)/r-cot(theta)*(sphi*vtheta)/r
902  source(ixo^s) = -(wct(ixo^s, mom(3)) * wprim(ixo^s, mom(1))) &
903  - (wct(ixo^s, mom(3)) * wprim(ixo^s, mom(2))) / dtan(x(ixo^s, 2))
904  w(ixo^s, mom(3)) = w(ixo^s, mom(3)) + qdt * source(ixo^s) / x(ixo^s, 1)
905  end if
906  }
907  end select
908 
909  end subroutine srhd_add_source_geom
910 
911  !> handles bootstrapping
912  subroutine srhd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
914  use mod_small_values
915  logical, intent(in) :: primitive
916  integer, intent(in) :: ixI^L,ixO^L
917  double precision, intent(inout) :: w(ixI^S,1:nw)
918  double precision, intent(in) :: x(ixI^S,1:ndim)
919  character(len=*), intent(in) :: subname
920 
921  integer :: n,idir
922  logical :: flag(ixI^S,1:nw),flagall(ixI^S)
923 
924  call srhd_check_w(primitive, ixi^l, ixo^l, w, flag)
925 
926  if (any(flag)) then
927  select case (small_values_method)
928  case ("replace")
929  ! any faulty cell is replaced by physical lower limit
930  flagall(ixo^s)=(flag(ixo^s,rho_).or.flag(ixo^s,e_))
931 
932  where(flagall(ixo^s))
933  ! D or rho: no difference primitive-conservative
934  w(ixo^s,rho_) = small_density
935  !w(ixO^S,lfac_)= 1.0d0
936  !w(ixO^S,xi_) = small_xi
937  endwhere
938  !do idir = 1, ndir
939  ! where(flagall(ixO^S)) w(ixO^S, mom(idir)) = 0.0d0
940  !end do
941  if(primitive) then
942  where(flagall(ixo^s)) w(ixo^s, p_) = small_pressure
943  else
944  where(flagall(ixo^s)) w(ixo^s, e_) = small_e
945  endif
946 
947  case ("average")
948  ! note: in small_values_average we use
949  ! small_values_fix_iw(1:nw) and small_values_daverage
950  ! when fails, may use small_pressure/small_density
951  if(primitive)then
952  ! averaging for all primitive fields (p, lfac*v, tau))
953  call small_values_average(ixi^l, ixo^l, w, x, flag)
954  ! update the auxiliaries from primitives
955  call srhd_get_auxiliary_prim(ixi^l,ixo^l,w)
956  else
957  ! do averaging of density d
958  call small_values_average(ixi^l, ixo^l, w, x, flag, d_)
959  ! do averaging of energy tau
960  call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
961  ! and now hope for the best....
962  endif
963  case default
964  if(.not.primitive) then
965  ! note that we throw error here, which assumes w is primitive
966  write(*,*) "handle_small_values default: note reporting conservatives!"
967  end if
968  call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
969  end select
970  end if
971 
972  end subroutine srhd_handle_small_values
973 
974  !> calculate effective gamma
975  subroutine srhd_get_geff_eos(w,ixI^L,ixO^L,varconserve,Geff)
976  use mod_global_parameters, only: nw
977  !================== IMPORTANT ==================!
978  !This subroutine is used with conserved variables in w when varconserve=T
979  !This subroutine is used with primitive variables in w when varconserve=F
980  ! both cases assume updated auxiliary variables xi_ en lfac_
981  !===============================================!
982  integer, intent(in) :: ixi^l,ixo^l
983  double precision, intent(in) :: w(ixi^s, 1:nw)
984  logical, intent(in) :: varconserve
985  double precision, intent(out) :: geff(ixi^s)
986 
987  double precision, dimension(ixO^S) :: pth,rho,e_th,e
988 
989  if (srhd_eos) then
990  if (varconserve) then
991  pth(ixo^s)=w(ixo^s,xi_)-w(ixo^s,e_)-w(ixo^s,d_)
992  rho(ixo^s)=w(ixo^s,d_)/w(ixo^s,lfac_)
993  e_th = pth*inv_gamma_1
994  e = e_th+dsqrt(e_th**2+rho**2)
995  geff(ixo^s) = srhd_gamma-half*gamma_1 * &
996  (one-(rho(ixo^s)/e(ixo^s))**2)
997  else
998  ! primitives available
999  e_th = w(ixo^s,p_)*inv_gamma_1
1000  e = e_th+dsqrt(e_th**2+w(ixo^s,rho_)**2)
1001  geff(ixo^s) = srhd_gamma-half*gamma_1 * &
1002  (one-(w(ixo^s,rho_)/e(ixo^s))**2)
1003  end if
1004  else
1005  geff(ixo^s) = srhd_gamma
1006  endif
1007 
1008  end subroutine srhd_get_geff_eos
1009 
1010  !> Compute the small value limits
1013  implicit none
1014  ! local small values
1015  double precision :: LsmallE,Lsmallp,Lsmallrho
1016 
1017  ! the maximal allowed Lorentz factor
1018  lfacmax=one/dsqrt(one-(one-dmaxvel)**2)
1021  if(small_density*small_pressure<=0.0d0)then
1022  call mpistop("must set finite values small-density/pressure for small value treatments")
1023  endif
1024  if(srhd_eos)then
1025  lsmallp=(one+10.d0*small_pressure)*small_pressure
1026  lsmallrho=(one+10.d0*small_density)*small_density
1027  !!Lsmallp=small_pressure
1028  !!Lsmallrho=small_density
1029  lsmalle=lsmallp*inv_gamma_1+&
1030  dsqrt((lsmallp*inv_gamma_1)**2+lsmallrho**2)
1031  small_xi=half*((srhd_gamma+one)*lsmalle-&
1032  gamma_1*lsmallrho*(lsmallrho/lsmalle))
1033  small_e=small_xi-lsmallp-lsmallrho
1034  else
1037  endif
1040 
1041  end subroutine srhd_get_smallvalues_eos
1042 
1043  !> Compute the enthalpy rho*h from rho and pressure p
1044  subroutine srhd_get_enthalpy_eos(ixO^L,rho,p,rhoh)
1046  integer, intent(in) :: ixo^l
1047  double precision, intent(in) :: rho(ixo^s),p(ixo^s)
1048  double precision, intent(out) :: rhoh(ixo^s)
1049 
1050  double precision, dimension(ixO^S) :: e_th,e
1051  integer :: ix^d
1052 
1053  if(srhd_eos) then
1054  e_th = p*inv_gamma_1
1055  e = e_th+dsqrt(e_th**2+rho**2)
1056  ! writing rho/E on purpose, for numerics
1057  rhoh = half*((srhd_gamma+one)*e &
1058  - gamma_1*rho*(rho/e))
1059  else
1060  rhoh = rho+gamma_to_gamma_1*p
1061  end if
1062 
1063  if (check_small_values) then
1064  {do ix^db= ixo^lim^db\}
1065  if(rhoh(ix^d)<small_xi) then
1066  write(*,*) "local pressure and density",p(ix^d),rho(ix^d)
1067  write(*,*) "Error: small value of enthalpy rho*h=",rhoh(ix^d),&
1068  " encountered when call srhd_get_enthalpy_eos"
1069  call mpistop('enthalpy below small_xi: stop (may need to turn on fixes)')
1070  end if
1071  {enddo^d&\}
1072  end if
1073 
1074  if (fix_small_values) then
1075  {do ix^db= ixo^lim^db\}
1076  if(rhoh(ix^d)<small_xi) then
1077  rhoh(ix^d)=small_xi
1078  endif
1079  {enddo^d&\}
1080  endif
1081 
1082  end subroutine srhd_get_enthalpy_eos
1083 
1084  !> Calculate thermal pressure p from density rho and enthalpy rho*h
1085  !> will provide p (and E if srhd_eos)
1086  subroutine srhd_get_pressure_eos(ixI^L,ixO^L,rho,rhoh,p,E)
1088  integer, intent(in) :: ixI^L, ixO^L
1089  double precision, intent(in) :: rho(ixO^S),rhoh(ixO^S)
1090  double precision, intent(out) :: p(ixI^S)
1091  double precision, intent(out) :: E(ixO^S)
1092  integer :: ix^D
1093 
1094  if(srhd_eos) then
1095  e = (rhoh+dsqrt(rhoh**2+(srhd_gamma**2-one)*rho**2)) &
1096  /(srhd_gamma+one)
1097  p(ixo^s) = half*gamma_1* (e-rho*(rho/e))
1098  else
1099  p(ixo^s) = (rhoh-rho)/gamma_to_gamma_1
1100  end if
1101 
1102  if (check_small_values) then
1103  {do ix^db= ixo^lim^db\}
1104  if(p(ix^d)<small_pressure) then
1105  write(*,*) "local enthalpy rho*h and density rho",rhoh(ix^d),rho(ix^d)
1106  if(srhd_eos) write(*,*) 'E, rho^2/E, difference', &
1107  e(ix^d),rho(ix^d)**2/e(ix^d),e(ix^d)-rho(ix^d)**2/e(ix^d)
1108  write(*,*) "Error: small value of gas pressure",p(ix^d),&
1109  " encountered when call srhd_get_pressure_eos"
1110  call mpistop('pressure below small_pressure: stop (may need to turn on fixes)')
1111  end if
1112  {enddo^d&\}
1113  end if
1114 
1115  if (fix_small_values) then
1116  {do ix^db= ixo^lim^db\}
1117  if(p(ix^d)<small_pressure) then
1118  p(ix^d)=small_pressure
1119  if(srhd_eos)e(ix^d)=max(small_e,e(ix^d))
1120  endif
1121  {enddo^d&\}
1122  endif
1123 
1124  end subroutine srhd_get_pressure_eos
1125 
1126  !> Calculate the square of the thermal sound speed csound2 within ixO^L.
1127  !> available rho - rho*h
1128  subroutine srhd_get_csound2_eos(ixI^L,ixO^L,rho,rhoh,csound2)
1130  integer, intent(in) :: ixI^L, ixO^L
1131  double precision, intent(in) :: rho(ixO^S),rhoh(ixO^S)
1132  double precision, intent(out) :: csound2(ixO^S)
1133 
1134  double precision :: p(ixI^S)
1135  double precision :: E(ixO^S)
1136  integer :: ix^D
1137 
1138  call srhd_get_pressure_eos(ixi^l,ixo^l,rho,rhoh,p,e)
1139  if(srhd_eos) then
1140  csound2(ixo^s)=(p(ixo^s)*((srhd_gamma+one)&
1141  +gamma_1*(rho(ixo^s)/e(ixo^s))**2))&
1142  /(2.0d0*rhoh(ixo^s))
1143  else
1144  csound2(ixo^s)=srhd_gamma*p(ixo^s)/rhoh(ixo^s)
1145  end if
1146 
1147  if (check_small_values) then
1148  {do ix^db= ixo^lim^db\}
1149  if(csound2(ix^d)>=1.0d0.or.csound2(ix^d)<=0.0d0) then
1150  write(*,*) "sound speed error with p - rho - rhoh",p(ix^d),rhoh(ix^d),rho(ix^d)
1151  if(srhd_eos) write(*,*) 'and E', e(ix^d)
1152  write(*,*) "Error: value of csound2",csound2(ix^d),&
1153  " encountered when call srhd_get_csound2_eos"
1154  call mpistop('sound speed stop (may need to turn on fixes)')
1155  end if
1156  {enddo^d&\}
1157  end if
1158 
1159  if (fix_small_values) then
1160  {do ix^db= ixo^lim^db\}
1161  if(csound2(ix^d)>=1.0d0) then
1162  csound2(ix^d)=1.0d0-1.0d0/lfacmax**2
1163  endif
1164  if(csound2(ix^d)<=0.0d0) then
1165  csound2(ix^d)=srhd_gamma*small_pressure/small_xi
1166  endif
1167  {enddo^d&\}
1168  endif
1169 
1170  end subroutine srhd_get_csound2_eos
1171 
1172  !> Calculate the square of the thermal sound speed csound2 within ixO^L.
1173  !> available rho - rho*h - p
1174  subroutine srhd_get_csound2_prim_eos(ixO^L,rho,rhoh,p,csound2)
1176  integer, intent(in) :: ixO^L
1177  double precision, intent(in) :: rho(ixO^S),rhoh(ixO^S),p(ixO^S)
1178  double precision, intent(out) :: csound2(ixO^S)
1179 
1180  double precision :: E(ixO^S)
1181  integer :: ix^D
1182 
1183  if(srhd_eos) then
1184  e = (rhoh+dsqrt(rhoh**2+(srhd_gamma**2-one)&
1185  *rho**2))/(srhd_gamma+one)
1186  csound2(ixo^s)=(p*((srhd_gamma+one)&
1187  +gamma_1*(rho/e)**2))&
1188  /(2.0d0*rhoh)
1189  else
1190  csound2(ixo^s)=srhd_gamma*p/rhoh
1191  end if
1192 
1193  if (check_small_values) then
1194  {do ix^db= ixo^lim^db\}
1195  if(csound2(ix^d)>=1.0d0.or.csound2(ix^d)<=0.0d0) then
1196  write(*,*) "sound speed error with p - rho - rhoh",p(ix^d),rhoh(ix^d),rho(ix^d)
1197  if(srhd_eos) write(*,*) 'and E', e(ix^d)
1198  write(*,*) "Error: value of csound2",csound2(ix^d),&
1199  " encountered when call srhd_get_csound2_prim_eos"
1200  call mpistop('sound speed stop (may need to turn on fixes)')
1201  end if
1202  {enddo^d&\}
1203  end if
1204 
1205  if (fix_small_values) then
1206  {do ix^db= ixo^lim^db\}
1207  if(csound2(ix^d)>=1.0d0) then
1208  csound2(ix^d)=1.0d0-1.0d0/lfacmax**2
1209  endif
1210  if(csound2(ix^d)<=0.0d0) then
1211  csound2(ix^d)=srhd_gamma*small_pressure/small_xi
1212  endif
1213  {enddo^d&\}
1214  endif
1215 
1216  end subroutine srhd_get_csound2_prim_eos
1217 
1218  !> con2prim: (D,S**2,tau) --> compute auxiliaries lfac and xi
1219  subroutine con2prim_eos(lfac,xi,myd,myssqr,mytau,ierror)
1220  use mod_con2prim_vars
1221 
1222  double precision, intent(in) :: myd, myssqr, mytau
1223  double precision, intent(inout) :: lfac, xi
1224  integer, intent(inout) :: ierror
1225 
1226  ! .. local ..
1227  double precision:: f,df,lfacl
1228  !------------------------------------------------------------------
1229 
1230  ! Save the input-state in mod_con2prim_vars
1231  d = myd; ssqr = myssqr; tau = mytau;
1232 
1233  ierror=0
1234 
1235  ! Check if guess is close enough: gives f,df,lfacl
1236  if(xi>smallxi)then
1237  call funcd_eos(xi,f,df,lfacl,d,ssqr,tau,ierror)
1238  if (ierror == 0 .and. dabs(f/df)<absaccnr) then
1239  xi = xi - f/df
1240  lfac = lfacl
1241  return
1242  else
1243  ierror = 0
1244  end if
1245  else
1246  write(*,*)'entering con2prim_eos with xi=',xi
1247  end if
1248 
1249  ! ierror=1 : must have D>=minrho, tau>=smalltau
1250  !if(d<minrho .or. tau<smalltau)then
1251  ! ierror=1
1252  ! return
1253  !endif
1254 
1255  call con2primhydro_eos(lfac,xi,d,ssqr,tau,ierror)
1256 
1257  end subroutine con2prim_eos
1258 
1259  subroutine funcd_eos(xi,f,df,mylfac,d,ssqr,tau,ierror)
1260  double precision, intent(in) :: xi,d,ssqr,tau
1261  double precision, intent(out) :: f,df,mylfac
1262  integer, intent(inout) :: ierror
1263 
1264  ! .. local ..
1265  double precision :: dlfac
1266  double precision :: vsqr,p,dpdxi
1267  !-----------------------------------------------------------------
1268 
1269  vsqr = ssqr/xi**2
1270 
1271  if (vsqr<one) then
1272  mylfac = one/dsqrt(one-vsqr)
1273  dlfac = -mylfac**3*ssqr/(xi**3)
1274  !===== Pressure, calculate using EOS =====!
1275  call funcpressure_eos(xi,mylfac,d,dlfac,p,dpdxi)
1276  !=========================================!
1277  f = xi-tau-d-p
1278  df = one-dpdxi
1279  else
1280  ! print *,'Erroneous input to funcd since vsqr=',vsqr,' >=1'
1281  ! print *,'input values d, ssqr, tau:',d,ssqr,tau
1282  ierror =6
1283  return
1284  end if
1285 
1286  end subroutine funcd_eos
1287 
1288  !> SRHD iteration solves for p via NR, and then gives xi as output
1289  subroutine con2primhydro_eos(lfac,xi,d,sqrs,tau,ierror)
1290  double precision, intent(out) :: xi,lfac
1291  double precision, intent(in) :: d,sqrs,tau
1292  integer,intent(inout) :: ierror
1293 
1294  ! .. local ..
1295  integer :: ni,niiter,nit,n2it,ni3
1296  double precision :: pcurrent,pnew
1297  double precision :: er,er1,ff,df,dp,v2
1298  double precision :: pmin,lfac2inv,pLabs,pRabs,pprev
1299  double precision :: s2overcubeG2rh
1300  double precision :: xicurrent,h,dhdp
1301  double precision :: oldff1,oldff2,Nff
1302  double precision :: pleft,pright
1303  !---------------------------------------------------------------------
1304 
1305  ierror=0
1306  ! ierror=0 : ok
1307  ! we already checked D>=minrho, tau>=smalltau (ierror=1)
1308  !
1309  ! ierror<>0
1310  !
1311  ! ierror=2 : maxitnr reached without convergence
1312  ! ierror=3 : final pressure value < smallp or xi<smallxi during iteration
1313  ! ierror=4 : final v^2=1 hence problem as lfac=1/0
1314  ! ierror=5 : nonmonotonic function f (as df=0)
1315 
1316  ! left and right brackets for p-range
1317  pmin=dsqrt(sqrs)/(one-dmaxvel)-tau-d
1318  plabs=max(minp,pmin)
1319  prabs=1.0d99
1320  ! start value from input
1321  pcurrent=plabs
1322 
1323  er1=one
1324  pprev=pcurrent
1325 
1326  ! Fudge Parameters
1327  oldff1=1.0d7 ! High number
1328  oldff2=1.0d9 ! High number bigger then oldff1
1329  n2it = 0
1330  nit = 0
1331 
1332  loopnr: do ni=1,maxitnr
1333  nit = nit + 1
1334  !=== Relax NR iteration accuracy=======!
1335  if(nit>maxitnr/4)then
1336  ! mix pressure value for convergence
1337  pcurrent=half*(pcurrent+pprev)
1338  ! relax accuracy requirement
1339  er1=10.0d0*er1
1340  nit = nit - maxitnr/10
1341  endif
1342  !=======================================!
1343 
1344  niiter=ni
1345  xicurrent=tau+d+pcurrent
1346 
1347  if(xicurrent<smallxi) then
1348  ! print *,'stop: too small xi iterate:',xicurrent
1349  ! print *,'for pressure iterate p',pcurrent
1350  ! print *,'pressure bracket pLabs pRabs',pLabs,pRabs
1351  ! print *,'iteration number:',ni
1352  ! print *,'values for d,s,tau,s2:',d,sqrs,tau,sqrs
1353  ierror=3
1354  return
1355  endif
1356 
1357  v2=sqrs/xicurrent**2
1358  lfac2inv=one - v2
1359  if(lfac2inv>zero) then
1360  lfac=one/dsqrt(lfac2inv)
1361  else
1362  ! print *,'stop: negative or zero factor 1-v2:',lfac2inv
1363  ! print *,'for pressure iterate p',pcurrent
1364  ! print *,'absolute pressure bracket pLabs pRabs',pLabs,pRabs
1365  ! print *,'iteration number:',ni
1366  ! print *,'values for d,s,tau,s2:',d,sqrs,tau,sqrs
1367  ! print *,'values for v2,xi:',v2,xicurrent
1368  ierror=4
1369  return
1370  endif
1371 
1372  s2overcubeg2rh=sqrs/(xicurrent**3)
1373  !== calculation done using the EOS ==!
1374  call funcenthalpy_eos(pcurrent,lfac2inv,d,sqrs,xicurrent,&
1375  s2overcubeg2rh,h,dhdp)
1376  !=======================================!
1377  ff=-xicurrent*lfac2inv + h
1378  df=- two*sqrs/xicurrent**2 + dhdp - lfac2inv
1379 
1380  if (ff*df==zero) then
1381  if (ff==zero) then
1382  exit ! zero found
1383  else
1384  ! print *,'stop: df becomes zero, non-monotonic f(p)!'
1385  ierror=5
1386  return
1387  endif
1388  else
1389  pnew=pcurrent-ff/df
1390  if (ff*df>zero) then
1391  ! pressure iterate has decreased
1392  ! restrict to left
1393  pnew=max(pnew,plabs)
1394  else ! ff*df<0
1395  ! pressure iterate has increased
1396  ! restrict to right
1397  pnew=min(pnew,prabs)
1398  endif
1399  endif
1400 
1401  !===============================================!
1402  dp=pcurrent-pnew
1403  er=two*dabs(dp)/(pnew+pcurrent)
1404  if(((er<tolernr*er1).or.(dabs(dp)<absaccnr))) exit loopnr
1405  !===============================================!
1406 
1407  ! For very small values of pressure, NR algorithm is not efficient to
1408  ! find root, use Euler algorithm to find precise value of pressure
1409  if((dabs(oldff2-ff) < 1.0d-8 .or. niiter >= maxitnr-maxitnr/20).and.&
1410  ff * oldff1 < zero .and. dabs(ff)>absaccnr)then
1411 
1412  n2it=n2it+1
1413  if(n2it<=3) pcurrent=half*(pnew+pcurrent)
1414  if(n2it>3)then
1415  pright =pcurrent
1416  pleft=pprev
1417  pcurrent=half*(pleft+pright)
1418  dicho: do ni3=1,maxitnr
1419  !===================!
1420  xicurrent=tau+d+pcurrent
1421  v2=sqrs/xicurrent**2
1422  lfac2inv=one - v2
1423 
1424  if(lfac2inv>zero)then
1425  lfac=one/dsqrt(lfac2inv)
1426  else
1427  ierror=4
1428  return
1429  endif
1430  !===================!
1431 
1432  !== calculation done using the EOS ==!
1433  call bisection_enthalpy_eos(pnew,lfac2inv,d,xicurrent,h)
1434  nff=-xicurrent*lfac2inv + h
1435  !=======================================!
1436  !==== Iterate ====!
1437  if(ff * nff < zero)then
1438  pleft=pcurrent
1439  else
1440  pright=pcurrent
1441  endif
1442 
1443  pcurrent=half*(pleft+pright)
1444  !==================!
1445 
1446  !=== The iteration converged ===!
1447  if(2.0d0*dabs(pleft-pright)/(pleft+pright)< absaccnr &
1448  .or. dabs(ff)<absaccnr)then
1449  pnew=pcurrent
1450  exit loopnr
1451  endif
1452  !==============================!
1453 
1454  !==============================!
1455 
1456  !=== conserve the last value of Nff ===!
1457  ff=nff
1458  !======================================!
1459  enddo dicho
1460  endif
1461 
1462  else
1463  !====== There is no problems, continue the NR iteration ======!
1464  pprev=pcurrent
1465  pcurrent=pnew
1466  !=============================================================!
1467  endif
1468 
1469 
1470  !=== keep the values of the 2 last ff ===!
1471  oldff2=oldff1
1472  oldff1=ff
1473  !========================================!
1474  enddo loopnr
1475 
1476  if(niiter==maxitnr)then
1477  ierror=2
1478  return
1479  endif
1480 
1481  if(pcurrent<minp) then
1482  ierror=3
1483  return
1484  endif
1485 
1486  !------------------------------!
1487  xi=tau+d+pcurrent
1488  v2=sqrs/xicurrent**2
1489  lfac2inv=one - v2
1490  if(lfac2inv>zero) then
1491  lfac=one/dsqrt(lfac2inv)
1492  else
1493  ierror=4
1494  return
1495  endif
1496 
1497  end subroutine con2primhydro_eos
1498 
1499  !> pointwise evaluations used in con2prim
1500  !> compute pointwise value for pressure p and dpdxi
1501  subroutine funcpressure_eos(xicurrent,lfac,d,dlfacdxi,p,dpdxi)
1502 
1503  double precision, intent(in) :: xicurrent,lfac,d,dlfacdxi
1504  double precision, intent(out) :: p,dpdxi
1505  ! .. local ..
1506  double precision :: rho,h,E,dhdxi,rhotoE
1507  double precision :: dpdchi,dEdxi
1508 
1509  ! rhoh here called h
1510  h=xicurrent/(lfac**2)
1511  rho=d/lfac
1512  e = (h+dsqrt(h**2+(srhd_gamma**2-one)*rho**2)) &
1513  /(srhd_gamma+one)
1514  ! output pressure
1515  rhotoe = rho/e
1516  p = half*gamma_1*(e-rho*rhotoe)
1517 
1518  dhdxi = one/(lfac**2)-2.0d0*xicurrent/(lfac**2)*dlfacdxi/lfac
1519 
1520  dedxi=(dhdxi+(h*dhdxi-(srhd_gamma**2-one)*rho**2*dlfacdxi/lfac)&
1521  /dsqrt(h**2+(srhd_gamma**2-one)*rho**2))&
1522  /(srhd_gamma+one)
1523 
1524  ! output pressure derivative to xi
1525  dpdxi=half*gamma_1*(2.0d0*rho*rhotoe*dlfacdxi/lfac+&
1526  (one+rhotoe**2)*dedxi)
1527 
1528  end subroutine funcpressure_eos
1529 
1530  !> pointwise evaluations used in con2prim
1531  !> returns enthalpy rho*h (h) and derivative d(rho*h)/dp (dhdp)
1532  subroutine funcenthalpy_eos(pcurrent,lfac2inv,d,sqrs,xicurrent,dv2d2p,h,dhdp)
1533 
1534  double precision, intent(in) :: pcurrent,lfac2inv,d,sqrs,xicurrent,dv2d2p
1535  double precision, intent(out):: h,dhdp
1536 
1537  ! local
1538  double precision:: rho,E_th,E,dE_thdp,dEdp
1539 
1540  rho=d*dsqrt(lfac2inv)
1541  e_th = pcurrent*inv_gamma_1
1542  e = (e_th + dsqrt(e_th**2+rho**2))
1543  !== Enthalpy ==!
1544  h = half*((srhd_gamma+one)*e-gamma_1*rho*(rho/e))
1545  !=== Derivative of thermal energy ===!
1546  de_thdp = one*inv_gamma_1
1547  !=== Derivative of internal energy ===!
1548  dedp = de_thdp * (one+e_th/dsqrt(e_th**2+rho**2))&
1549  + d**2*dv2d2p/dsqrt(e_th**2+rho**2)
1550  !====== Derivative of Enthalpy ======!
1551  dhdp = half*((srhd_gamma+one)*dedp + &
1552  gamma_1*(rho*(rho/e))*(-2.0d0*dv2d2p/lfac2inv+dedp/e))
1553  end subroutine funcenthalpy_eos
1554 
1555  !> pointwise evaluations used in con2prim
1556  !> returns enthalpy rho*h (h)
1557  subroutine bisection_enthalpy_eos(pcurrent,lfac2inv,d,xicurrent,h)
1558 
1559  double precision, intent(in) :: pcurrent,lfac2inv,d,xicurrent
1560  double precision, intent(out):: h
1561 
1562  ! local
1563  double precision:: rho,E_th,E
1564 
1565  rho=d*dsqrt(lfac2inv)
1566  e_th = pcurrent*inv_gamma_1
1567  e = (e_th + dsqrt(e_th**2+rho**2))
1568  !== Enthalpy ==!
1569  h = half*((srhd_gamma+one)*e-gamma_1*rho*(rho/e))
1570 
1571  return
1572  end subroutine bisection_enthalpy_eos
1573 
1574  !> con2prim: (D,S**2,tau) --> compute auxiliaries lfac and xi
1575  subroutine con2prim(lfac,xi,myd,myssqr,mytau,ierror)
1576  use mod_con2prim_vars
1577 
1578  double precision, intent(in) :: myd, myssqr, mytau
1579  double precision, intent(inout) :: lfac, xi
1580  integer, intent(inout) :: ierror
1581 
1582  ! .. local ..
1583  double precision:: f,df,lfacl
1584  !------------------------------------------------------------------
1585 
1586  ! Save the input-state in mod_con2prim_vars
1587  d = myd; ssqr = myssqr; tau = mytau;
1588 
1589  ierror=0
1590 
1591  ! Check if guess is close enough: gives f,df,lfacl
1592  if(xi>smallxi)then
1593  call funcd(xi,f,df,lfacl,d,ssqr,tau,ierror)
1594  if (ierror == 0 .and. dabs(f/df)<absaccnr) then
1595  xi = xi - f/df
1596  lfac = lfacl
1597  return
1598  else
1599  ierror = 0
1600  end if
1601  else
1602  write(*,*) 'entering con2prim with xi=',xi
1603  end if
1604 
1605  ! ierror=1 : must have D>=minrho, tau>=smalltau
1606  !if(d<minrho .or. tau<smalltau)then
1607  ! ierror=1
1608  ! return
1609  !endif
1610 
1611  call con2primhydro(lfac,xi,d,ssqr,tau,ierror)
1612 
1613  end subroutine con2prim
1614 
1615  subroutine funcd(xi,f,df,mylfac,d,ssqr,tau,ierror)
1616  double precision, intent(in) :: xi,d,ssqr,tau
1617  double precision, intent(out) :: f,df,mylfac
1618  integer, intent(inout) :: ierror
1619 
1620  ! .. local ..
1621  double precision :: dlfac
1622  double precision :: vsqr,p,dpdxi
1623  !-----------------------------------------------------------------
1624 
1625  vsqr = ssqr/xi**2
1626 
1627  if (vsqr<one) then
1628  mylfac = one/dsqrt(one-vsqr)
1629  dlfac = -mylfac**3*ssqr/(xi**3)
1630  !===== Pressure, calculate using EOS =====!
1631  call funcpressure(xi,mylfac,d,dlfac,p,dpdxi)
1632  !=========================================!
1633  f = xi-tau-d-p
1634  df = one-dpdxi
1635  else
1636  ! print *,'Erroneous input to funcd since vsqr=',vsqr,' >=1'
1637  ! print *,'input values d, ssqr, tau:',d,ssqr,tau
1638  ierror =6
1639  return
1640  end if
1641 
1642  end subroutine funcd
1643 
1644  !> SRHD iteration solves for p via NR, and then gives xi as output
1645  subroutine con2primhydro(lfac,xi,d,sqrs,tau,ierror)
1646  double precision, intent(out) :: xi,lfac
1647  double precision, intent(in) :: d,sqrs,tau
1648  integer,intent(inout) :: ierror
1649 
1650  ! .. local ..
1651  integer :: ni,niiter,nit,n2it,ni3
1652  double precision :: pcurrent,pnew
1653  double precision :: er,er1,ff,df,dp,v2
1654  double precision :: pmin,lfac2inv,pLabs,pRabs,pprev
1655  double precision :: s2overcubeG2rh
1656  double precision :: xicurrent,h,dhdp
1657  double precision :: oldff1,oldff2,Nff
1658  double precision :: pleft,pright
1659  !---------------------------------------------------------------------
1660 
1661  ierror=0
1662  ! ierror=0 : ok
1663  ! we already checked D>=minrho, tau>=smalltau (ierror=1)
1664  !
1665  ! ierror<>0
1666  !
1667  ! ierror=2 : maxitnr reached without convergence
1668  ! ierror=3 : final pressure value < smallp or xi<smallxi during iteration
1669  ! ierror=4 : final v^2=1 hence problem as lfac=1/0
1670  ! ierror=5 : nonmonotonic function f (as df=0)
1671 
1672  ! left and right brackets for p-range
1673  pmin=dsqrt(sqrs)/(one-dmaxvel)-tau-d
1674  plabs=max(minp,pmin)
1675  prabs=1.0d99
1676  ! start value from input
1677  pcurrent=plabs
1678 
1679  er1=one
1680  pprev=pcurrent
1681 
1682  ! Fudge Parameters
1683  oldff1=1.0d7 ! High number
1684  oldff2=1.0d9 ! High number bigger then oldff1
1685  n2it = 0
1686  nit = 0
1687 
1688  loopnr: do ni=1,maxitnr
1689  nit = nit + 1
1690  !=== Relax NR iteration accuracy=======!
1691  if(nit>maxitnr/4)then
1692  ! mix pressure value for convergence
1693  pcurrent=half*(pcurrent+pprev)
1694  ! relax accuracy requirement
1695  er1=10.0d0*er1
1696  nit = nit - maxitnr/10
1697  endif
1698  !=======================================!
1699 
1700  niiter=ni
1701  xicurrent=tau+d+pcurrent
1702 
1703  if(xicurrent<smallxi) then
1704  ! print *,'stop: too small xi iterate:',xicurrent
1705  ! print *,'for pressure iterate p',pcurrent
1706  ! print *,'pressure bracket pLabs pRabs',pLabs,pRabs
1707  ! print *,'iteration number:',ni
1708  ! print *,'values for d,s,tau,s2:',d,sqrs,tau,sqrs
1709  ierror=3
1710  return
1711  endif
1712 
1713  v2=sqrs/xicurrent**2
1714  lfac2inv=one - v2
1715  if(lfac2inv>zero) then
1716  lfac=one/dsqrt(lfac2inv)
1717  else
1718  ! print *,'stop: negative or zero factor 1-v2:',lfac2inv
1719  ! print *,'for pressure iterate p',pcurrent
1720  ! print *,'absolute pressure bracket pLabs pRabs',pLabs,pRabs
1721  ! print *,'iteration number:',ni
1722  ! print *,'values for d,s,tau,s2:',d,sqrs,tau,sqrs
1723  ! print *,'values for v2,xi:',v2,xicurrent
1724  ierror=4
1725  return
1726  endif
1727 
1728  s2overcubeg2rh=sqrs/(xicurrent**3)
1729  !== calculation done using the EOS ==!
1730  call funcenthalpy(pcurrent,lfac2inv,d,sqrs,xicurrent,&
1731  s2overcubeg2rh,h,dhdp)
1732  !=======================================!
1733  ff=-xicurrent*lfac2inv + h
1734  df=- two*sqrs/xicurrent**2 + dhdp - lfac2inv
1735 
1736  if (ff*df==zero) then
1737  if (ff==zero) then
1738  exit ! zero found
1739  else
1740  ! print *,'stop: df becomes zero, non-monotonic f(p)!'
1741  ierror=5
1742  return
1743  endif
1744  else
1745  pnew=pcurrent-ff/df
1746  if (ff*df>zero) then
1747  ! pressure iterate has decreased
1748  ! restrict to left
1749  pnew=max(pnew,plabs)
1750  else ! ff*df<0
1751  ! pressure iterate has increased
1752  ! restrict to right
1753  pnew=min(pnew,prabs)
1754  endif
1755  endif
1756 
1757  !===============================================!
1758  dp=pcurrent-pnew
1759  er=two*dabs(dp)/(pnew+pcurrent)
1760  if(((er<tolernr*er1).or.(dabs(dp)<absaccnr))) exit loopnr
1761  !===============================================!
1762 
1763  ! For very small values of pressure, NR algorithm is not efficient to
1764  ! find root, use Euler algorithm to find precise value of pressure
1765  if((dabs(oldff2-ff) < 1.0d-8 .or. niiter >= maxitnr-maxitnr/20).and.&
1766  ff * oldff1 < zero .and. dabs(ff)>absaccnr)then
1767 
1768  n2it=n2it+1
1769  if(n2it<=3) pcurrent=half*(pnew+pcurrent)
1770  if(n2it>3)then
1771  pright =pcurrent
1772  pleft=pprev
1773  pcurrent=half*(pleft+pright)
1774  dicho: do ni3=1,maxitnr
1775  !===================!
1776  xicurrent=tau+d+pcurrent
1777  v2=sqrs/xicurrent**2
1778  lfac2inv=one - v2
1779 
1780  if(lfac2inv>zero)then
1781  lfac=one/dsqrt(lfac2inv)
1782  else
1783  ierror=4
1784  return
1785  endif
1786  !===================!
1787 
1788  !== calculation done using the EOS ==!
1789  call bisection_enthalpy(pnew,lfac2inv,d,xicurrent,h)
1790  nff=-xicurrent*lfac2inv + h
1791  !=======================================!
1792  !==== Iterate ====!
1793  if(ff * nff < zero)then
1794  pleft=pcurrent
1795  else
1796  pright=pcurrent
1797  endif
1798 
1799  pcurrent=half*(pleft+pright)
1800  !==================!
1801 
1802  !=== The iteration converged ===!
1803  if(2.0d0*dabs(pleft-pright)/(pleft+pright)< absaccnr &
1804  .or. dabs(ff)<absaccnr)then
1805  pnew=pcurrent
1806  exit loopnr
1807  endif
1808  !==============================!
1809 
1810  !==============================!
1811 
1812  !=== conserve the last value of Nff ===!
1813  ff=nff
1814  !======================================!
1815  enddo dicho
1816  endif
1817 
1818  else
1819  !====== There is no problems, continue the NR iteration ======!
1820  pprev=pcurrent
1821  pcurrent=pnew
1822  !=============================================================!
1823  endif
1824 
1825 
1826  !=== keep the values of the 2 last ff ===!
1827  oldff2=oldff1
1828  oldff1=ff
1829  !========================================!
1830  enddo loopnr
1831 
1832  if(niiter==maxitnr)then
1833  ierror=2
1834  return
1835  endif
1836 
1837  if(pcurrent<minp) then
1838  ierror=3
1839  return
1840  endif
1841 
1842  !------------------------------!
1843  xi=tau+d+pcurrent
1844  v2=sqrs/xicurrent**2
1845  lfac2inv=one - v2
1846  if(lfac2inv>zero) then
1847  lfac=one/dsqrt(lfac2inv)
1848  else
1849  ierror=4
1850  return
1851  endif
1852 
1853  end subroutine con2primhydro
1854 
1855  !> pointwise evaluations used in con2prim
1856  !> compute pointwise value for pressure p and dpdxi
1857  subroutine funcpressure(xicurrent,lfac,d,dlfacdxi,p,dpdxi)
1858 
1859  double precision, intent(in) :: xicurrent,lfac,d,dlfacdxi
1860  double precision, intent(out) :: p,dpdxi
1861  ! .. local ..
1862  double precision :: rho,h,E,dhdxi,rhotoE
1863  double precision :: dpdchi,dEdxi
1864 
1865  ! rhoh here called h
1866  h=xicurrent/(lfac**2)
1867  rho=d/lfac
1868  ! output pressure
1869  p = (h - rho)/gamma_to_gamma_1
1870  dpdchi = one/gamma_to_gamma_1
1871  dpdxi = dpdchi * one/lfac**2
1872  ! zero case dlfacdxi implies zero velocity (ssqr=0)
1873  if (dlfacdxi /= 0.0d0) &
1874  dpdxi = dpdxi + dpdchi * ((d*lfac-2.0d0*xicurrent)/lfac**3) * dlfacdxi
1875 
1876  end subroutine funcpressure
1877 
1878  !> pointwise evaluations used in con2prim
1879  !> returns enthalpy rho*h (h) and derivative d(rho*h)/dp (dhdp)
1880  subroutine funcenthalpy(pcurrent,lfac2inv,d,sqrs,xicurrent,dv2d2p,h,dhdp)
1881 
1882  double precision, intent(in) :: pcurrent,lfac2inv,d,sqrs,xicurrent,dv2d2p
1883  double precision, intent(out):: h,dhdp
1884 
1885  ! local
1886  double precision:: rho,E_th,E,dE_thdp,dEdp
1887 
1888  rho=d*dsqrt(lfac2inv)
1889  h = rho + gamma_to_gamma_1 * pcurrent
1890  dhdp = gamma_to_gamma_1 + d/dsqrt(lfac2inv)*sqrs/xicurrent**3
1891  end subroutine funcenthalpy
1892 
1893  !> pointwise evaluations used in con2prim
1894  !> returns enthalpy rho*h (h)
1895  subroutine bisection_enthalpy(pcurrent,lfac2inv,d,xicurrent,h)
1896 
1897  double precision, intent(in) :: pcurrent,lfac2inv,d,xicurrent
1898  double precision, intent(out):: h
1899 
1900  ! local
1901  double precision:: rho,E_th,E
1902 
1903  rho=d*dsqrt(lfac2inv)
1904  h = rho + gamma_to_gamma_1 * pcurrent
1905 
1906  return
1907  end subroutine bisection_enthalpy
1908 
1909 end module mod_srhd_phys
Calculate w(iw)=w(iw)+qdt*SOURCE[wCT,qtC,x] within ixO for all indices iw=iwmin......
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
double precision ssqr
double precision d
double precision tau
Module for physical and numeric constants.
Definition: mod_constants.t:2
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
integer coordinate
Definition: mod_geometry.t:7
integer, parameter spherical
Definition: mod_geometry.t:11
integer, parameter cylindrical
Definition: mod_geometry.t:10
integer, parameter cartesian_expansion
Definition: mod_geometry.t:12
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 small_pressure
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.
double precision unit_pressure
Physical scaling factor for pressure.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
logical use_particles
Use particles module or not.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision unit_velocity
Physical scaling factor for velocity.
double precision unit_temperature
Physical scaling factor for temperature.
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
logical phys_trac
Use TRAC for MHD or 1D HD.
logical fix_small_values
fix small values with average or replace methods
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
double precision small_density
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
Module containing all the particle routines.
Definition: mod_particles.t:2
subroutine particles_init()
Initialize particle data and parameters.
Definition: mod_particles.t:15
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
Module for handling problematic values in simulations, such as negative pressures.
logical, public trace_small_values
trace small values in the source file using traceback flag of compiler
subroutine, public small_values_average(ixIL, ixOL, w, x, w_flag, windex)
subroutine, public small_values_error(wprim, x, ixIL, ixOL, w_flag, subname)
character(len=20), public small_values_method
How to handle small values.
Special Relativistic Hydrodynamics (with EOS) physics module.
Definition: mod_srhd_phys.t:2
subroutine srhd_add_source(qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
dummy addsource subroutine
subroutine srhd_get_flux(wC, wP, x, ixIL, ixOL, idim, f)
Calculate fluxes within ixO^L.
logical, public, protected srhd_force_diagonal
Allows overruling default corner filling (for debug mode, otherwise corner primitives fail)
Definition: mod_srhd_phys.t:49
subroutine srhd_physical_units
double precision, public tolernr
Definition: mod_srhd_phys.t:57
subroutine srhd_get_cmax_loc(ixIL, ixOL, vidim, csound2, v2, cmax, cmin)
local version for recycling code when computing cmax-cmin
integer, public, protected srhd_n_tracer
Number of tracer species.
Definition: mod_srhd_phys.t:13
double precision, public inv_gamma_1
Definition: mod_srhd_phys.t:41
subroutine srhd_write_info(fh)
Write this module's parameters to a snapshot.
Definition: mod_srhd_phys.t:96
subroutine srhd_get_cmax(w, x, ixIL, ixOL, idim, cmax)
Calculate cmax_idim within ixO^L used especially for setdt CFL limit.
integer, dimension(:), allocatable, public, protected tracer
Indices of the tracers.
Definition: mod_srhd_phys.t:23
subroutine, public srhd_get_auxiliary(ixIL, ixOL, w, x)
Compute auxiliary variables lfac and xi from a conservative state using srhd_con2prim to calculate en...
subroutine, public srhd_get_pthermal(w, x, ixIL, ixOL, pth)
Calculate thermal pressure p within ixO^L must follow after update conservative with auxiliaries.
double precision, public absaccnr
Definition: mod_srhd_phys.t:56
double precision, public small_e
The smallest allowed energy.
Definition: mod_srhd_phys.t:44
integer, public, protected lfac_
Index of the Lorentz factor.
Definition: mod_srhd_phys.t:31
subroutine funcd_eos(xi, f, df, mylfac, d, ssqr, tau, ierror)
double precision, public gamma_to_gamma_1
Definition: mod_srhd_phys.t:41
subroutine funcd(xi, f, df, mylfac, d, ssqr, tau, ierror)
subroutine srhd_get_csound2_eos(ixIL, ixOL, rho, rhoh, csound2)
Calculate the square of the thermal sound speed csound2 within ixO^L. available rho - rho*h.
subroutine, public srhd_get_auxiliary_prim(ixIL, ixOL, w)
Set auxiliary variables lfac and xi from a primitive state only used when handle_small_values average...
integer, public maxitnr
parameters for NR in con2prim
Definition: mod_srhd_phys.t:55
subroutine, public srhd_phys_init()
Initialize the module.
double precision, public gamma_1
Definition: mod_srhd_phys.t:41
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
Definition: mod_srhd_phys.t:52
subroutine, public srhd_get_csound2(w, x, ixIL, ixOL, csound2)
Calculate the square of the thermal sound speed csound2 within ixO^L. here computed from conservative...
double precision, public small_xi
The smallest allowed inertia.
Definition: mod_srhd_phys.t:46
subroutine, public srhd_to_conserved(ixIL, ixOL, w, x)
Transform primitive variables into conservative ones.
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
Definition: mod_srhd_phys.t:20
subroutine srhd_get_csound2_rhoh(w, x, ixIL, ixOL, rhoh, csound2)
Calculate the square of the thermal sound speed csound2 within ixO^L. here computed from conservative...
subroutine srhd_get_v(w, x, ixIL, ixOL, v)
Calculate v vector from conservatives.
subroutine, public srhd_get_enthalpy_eos(ixOL, rho, p, rhoh)
Compute the enthalpy rho*h from rho and pressure p.
double precision, public minp
Definition: mod_srhd_phys.t:60
subroutine, public srhd_to_primitive(ixIL, ixOL, w, x)
Transform conservative variables into primitive ones.
logical, public srhd_eos
Whether synge eos is used.
Definition: mod_srhd_phys.t:37
double precision, public lfacmax
Definition: mod_srhd_phys.t:59
integer, public, protected d_
Definition: mod_srhd_phys.t:17
subroutine srhd_handle_small_values(primitive, w, x, ixIL, ixOL, subname)
handles bootstrapping
subroutine srhd_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
dummy get_dt subroutine
subroutine srhd_get_smallvalues_eos
Compute the small value limits.
subroutine srhd_get_pressure_eos(ixIL, ixOL, rho, rhoh, p, E)
Calculate thermal pressure p from density rho and enthalpy rho*h will provide p (and E if srhd_eos)
subroutine, public srhd_check_params
subroutine srhd_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, wprim, w, x)
Add geometrical source terms to w.
integer, public, protected p_
Index of the gas pressure should equal e_.
Definition: mod_srhd_phys.t:28
subroutine con2prim_eos(lfac, xi, myd, myssqr, mytau, ierror)
con2prim: (D,S**2,tau) --> compute auxiliaries lfac and xi
subroutine srhd_get_csound2_prim_eos(ixOL, rho, rhoh, p, csound2)
Calculate the square of the thermal sound speed csound2 within ixO^L. available rho - rho*h - p.
subroutine srhd_get_cbounds(wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
Estimating bounds for the minimum and maximum signal velocities here we will not use Hspeed at all (o...
logical, public, protected srhd_particles
Whether particles module is added.
Definition: mod_srhd_phys.t:10
double precision, public smalltau
Definition: mod_srhd_phys.t:60
double precision, public dmaxvel
Definition: mod_srhd_phys.t:58
integer, public, protected e_
Index of the energy density.
Definition: mod_srhd_phys.t:26
subroutine srhd_check_w_aux(ixIL, ixOL, w, flag)
Returns logical argument flag T where auxiliary values are not ok.
double precision, public minrho
Definition: mod_srhd_phys.t:60
subroutine, public srhd_check_w(primitive, ixIL, ixOL, w, flag)
Returns logical argument flag T where values are not ok.
subroutine srhd_get_a2max(w, x, ixIL, ixOL, a2max)
double precision, public smallxi
Definition: mod_srhd_phys.t:60
integer, public, protected xi_
Index of the inertia.
Definition: mod_srhd_phys.t:34
double precision, public srhd_gamma
The adiabatic index and derived values.
Definition: mod_srhd_phys.t:40
integer, public, protected rho_
Index of the density (in the w array)
Definition: mod_srhd_phys.t:16
subroutine, public srhd_get_geff_eos(w, ixIL, ixOL, varconserve, Geff)
calculate effective gamma
Module with all the methods that users can customize in AMRVAC.
procedure(set_surface), pointer usr_set_surface
integer nwflux
Number of flux variables.
Definition: mod_variables.t:8