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