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