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