16 integer,
public,
protected ::
rho_
17 integer,
public,
protected ::
d_
20 integer,
allocatable,
public,
protected ::
mom(:)
23 integer,
allocatable,
public,
protected ::
tracer(:)
26 integer,
public,
protected ::
e_
28 integer,
public,
protected ::
p_
31 integer,
public,
protected ::
lfac_
34 integer,
public,
protected ::
xi_
52 double precision,
public ::
tolernr = 1.0d-9
53 double precision,
public ::
dmaxvel = 1.0d-7
72 subroutine srhd_read_params(files)
74 character(len=*),
intent(in) :: files(:)
82 open(
unitpar, file=trim(files(n)), status=
"old")
83 read(
unitpar, srhd_list,
end=111)
87 end subroutine srhd_read_params
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
99 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
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)
106 end subroutine srhd_write_info
115 physics_type =
"srhd"
117 phys_total_energy = .true.
121 phys_internal_e = .false.
122 phys_partial_ionization=.false.
128 allocate(start_indices(number_species),stop_indices(number_species))
137 mom(:) = var_set_momentum(
ndir)
140 e_ = var_set_energy()
144 call srhd_physical_units()
150 tracer(itr) = var_set_fluxvar(
"trc",
"trp", itr, need_bc=.false.)
155 xi_ = var_set_auxvar(
'xi',
'xi')
156 lfac_= var_set_auxvar(
'lfac',
'lfac')
162 stop_indices(1)=nwflux
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")
173 allocate(iw_vector(nvector))
174 iw_vector(1) =
mom(1) - 1
177 phys_add_source => srhd_add_source
178 phys_get_dt => srhd_get_dt
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
193 phys_get_v => srhd_get_v
194 phys_write_info => srhd_write_info
195 phys_handle_small_values => srhd_handle_small_values
204 use mod_particles,
only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
212 call mpistop (
"Error: srhd_gamma <= 0 or srhd_gamma == 1")
220 call srhd_get_smallvalues_eos
223 write(*,*)
'------------------------------------------------------------'
224 write(*,*)
'Using EOS set via srhd_eos=',
srhd_eos
225 write(*,*)
'Maximal lorentz factor (via dmaxvel) is=',
lfacmax
229 write(*,*)
'------------------------------------------------------------'
233 write(*,*)
'====SRHD run with settings===================='
234 write(*,*)
'Using mod_srhd_phys with settings:'
236 write(*,*)
'Dimensionality :',
ndim
237 write(*,*)
'vector components:',
ndir
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
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
259 write(*,*)
'number due to phys_wider_stencil=',phys_wider_stencil
260 write(*,*)
'==========================================='
266 subroutine srhd_physical_units
268 double precision :: mp,kb
281 call mpistop(
"Abort: must set positive values for unit length and numberdensity")
291 end subroutine srhd_physical_units
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)
309 where(w(ixo^s,
e_) <
small_e) flag(ixo^s,
e_) = .true.
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)
324 where(w(ixo^s,
lfac_) < one) flag(ixo^s,
lfac_) = .true.
326 if(any(flag(ixo^s,
xi_)))
then
327 write(*,*)
'auxiliary xi too low: abort'
328 call mpistop(
'auxiliary check failed')
330 if(any(flag(ixo^s,
lfac_)))
then
331 write(*,*)
'auxiliary lfac too low: abort'
332 call mpistop(
'auxiliary check failed')
335 end subroutine srhd_check_w_aux
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
346 rho(ixo^s) = sum(w(ixo^s,
mom(:))**2, dim=
ndim+1)
347 w(ixo^s,
lfac_) = dsqrt(1.0d0+rho(ixo^s))
349 rho(ixo^s)=w(ixo^s,
rho_)
350 pth(ixo^s)=w(ixo^s,
p_)
355 w(ixo^s,
xi_) = w(ixo^s,
lfac_)**2.0d0*rhoh(ixo^s)
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)
369 integer :: ix^
d,ierror,idir
370 integer :: flag_error(ixo^s)
371 double precision :: ssqr
374 {
do ix^db=ixomin^db,ixomax^db\}
378 ssqr= ssqr+w(ix^
d,
mom(idir))**2
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)
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)
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
401 {
do ix^db=ixomin^db,ixomax^db\}
405 ssqr= ssqr+w(ix^
d,
mom(idir))**2
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)
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)
424 w(ix^
d,
d_),ssqr,w(ix^
d,
e_),ierror)
425 flag_error(ix^
d) = ierror
430 if(any(flag_error(ixo^s)/=0))
then
431 print *,flag_error(ixo^s)
432 call mpistop(
'Problem when getting auxiliaries')
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)
446 double precision,
dimension(ixO^S) :: rhoh,rho,pth
450 rhoh(ixo^s) = sum(w(ixo^s,
mom(:))**2, dim=
ndim+1)
451 w(ixo^s,
lfac_) = dsqrt(1.0d0+rhoh(ixo^s))
453 rho(ixo^s)=w(ixo^s,
rho_)
454 pth(ixo^s)=w(ixo^s,
p_)
459 rhoh(ixo^s)= rhoh(ixo^s)*w(ixo^s,
lfac_)
461 w(ixo^s,
xi_) = w(ixo^s,
lfac_)*rhoh(ixo^s)
464 w(ixo^s,
d_)=w(ixo^s,
lfac_)*rho(ixo^s)
468 w(ixo^s,
mom(idir)) = rhoh(ixo^s)*w(ixo^s,
mom(idir))
472 w(ixo^s,
e_) = w(ixo^s,
xi_)-pth(ixo^s)-w(ixo^s,
d_)
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)
487 double precision,
dimension(ixO^S) :: rho,rhoh,e
488 double precision,
dimension(ixI^S) :: pth
489 character(len=30) :: subname_loc
495 rho(ixo^s) = w(ixo^s,
d_)/w(ixo^s,
lfac_)
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)
502 w(ixo^s,
rho_)=rho(ixo^s)
505 w(ixo^s,
mom(idir)) = w(ixo^s,
lfac_)*w(ixo^s,
mom(idir))&
508 w(ixo^s,
p_)=pth(ixo^s)
512 /(rho(ixo^s)*w(ixo^s,
lfac_))
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)
527 v(ixo^s,idir) = w(ixo^s,
mom(idir))/w(ixo^s,
xi_)
530 end subroutine srhd_get_v
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)
542 double precision :: rho(ixo^s)
545 call srhd_get_csound2_eos(ixi^
l,ixo^
l,rho,rhoh,csound2)
547 end subroutine srhd_get_csound2_rhoh
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)
559 double precision :: rho(ixo^s),rhoh(ixo^s)
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)
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)
582 double precision :: rho(ixo^s),rhoh(ixo^s),e(ixo^s)
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)
593 subroutine srhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
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
603 end subroutine srhd_add_source
606 subroutine srhd_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
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
616 end subroutine srhd_get_dt
620 subroutine srhd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
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)
627 double precision,
dimension(ixO^S) :: csound2,tmp1,tmp2,v2
628 double precision,
dimension(ixI^S) :: vidim, cmin
630 logical :: flag(ixi^s,1:nw)
635 tmp2=w(ixo^s,
xi_)/w(ixo^s,
lfac_)**2.0d0
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))
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)
653 cmax(ixo^s) = max(dabs(cmax(ixo^s)),dabs(cmin(ixo^s)))
655 end subroutine srhd_get_cmax
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)
667 double precision,
dimension(ixI^S):: tmp1,tmp2
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))
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))
680 cmin(ixo^s) = max(cmin(ixo^s), - 1.0d0)
681 cmin(ixo^s) = min(cmin(ixo^s), 1.0d0)
683 end subroutine srhd_get_cmax_loc
687 subroutine srhd_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
691 integer,
intent(in) :: ixi^
l, ixo^
l, idim
693 double precision,
intent(in) :: wlc(ixi^s,
nw), wrc(ixi^s,
nw)
695 double precision,
intent(in) :: wlp(ixi^s,
nw), wrp(ixi^s,
nw)
696 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
698 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
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
705 logical :: flag(ixi^s,1:
nw)
712 tmp2=wlp(ixo^s,
xi_)/wlp(ixo^s,
lfac_)**2.0d0
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)
722 tmp2=wrp(ixo^s,
xi_)/wrp(ixo^s,
lfac_)**2.0d0
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)
729 if(
present(cmin))
then
731 cmax(ixo^s,1)=max(cmaxl(ixo^s),cmaxr(ixo^s))
732 cmin(ixo^s,1)=min(cminl(ixo^s),cminr(ixo^s))
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))
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)
754 cmax(ixo^s,1)=max(cmaxl(ixo^s),dabs(cminl(ixo^s)))
762 tmp1=wmean(ixo^s,
rho_)
763 tmp2=wmean(ixo^s,
xi_)/wmean(ixo^s,
lfac_)**2.0d0
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)
773 cmax(ixo^s,1)=max(cmaxl(ixo^s),dabs(cminl(ixo^s)))
777 end subroutine srhd_get_cbounds
780 subroutine srhd_get_flux(wC,wP,x,ixI^L,ixO^L,idim,f)
782 integer,
intent(in) :: ixi^
l, ixo^
l, idim
784 double precision,
intent(in) :: wc(ixi^s,nw)
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)
790 double precision :: pth(ixi^s)
791 double precision :: v(ixi^s,1:
ndir)
794 pth(ixo^s)=wp(ixo^s,
p_)
796 v(ixo^s,idir) = wp(ixo^s,
mom(idir))/wp(ixo^s,
lfac_)
800 f(ixo^s,
d_)=v(ixo^s,idim)*wc(ixo^s,
rho_)
810 f(ixo^s,
mom(idir))= v(ixo^s,idim)*wc(ixo^s,
mom(idir))
812 f(ixo^s,
mom(idim))=pth(ixo^s)+f(ixo^s,
mom(idim))
816 f(ixo^s,
e_)=v(ixo^s,idim)*(wc(ixo^s,
e_) + pth(ixo^s))
818 end subroutine srhd_get_flux
821 subroutine srhd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
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)
829 double precision :: pth(ixi^s),
source(ixi^s), v(ixi^s,1:
ndir)
830 integer :: idir, h1x^
l{^nooned, h2x^
l}
832 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
843 source(ixo^s) =
source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
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_)
858 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
863 h1x^
l=ixo^
l-
kr(1,^
d); {^nooned h2x^
l=ixo^
l-
kr(2,^
d);}
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)
876 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s, 1)
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)
887 w(ixo^s,
mom(2)) = w(ixo^s,
mom(2)) + qdt *
source(ixo^s) / x(ixo^s, 1)
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)
898 end subroutine srhd_add_source_geom
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
911 logical :: flag(ixi^s,1:nw),flagall(ixi^s)
919 flagall(ixo^s)=(flag(ixo^s,
rho_).or.flag(ixo^s,
e_))
921 where(flagall(ixo^s))
933 where(flagall(ixo^s)) w(ixo^s,
e_) =
small_e
953 if(.not.primitive)
then
955 write(*,*)
"handle_small_values default: note reporting conservatives!"
961 end subroutine srhd_handle_small_values
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)
976 double precision,
dimension(ixO^S) :: pth,rho,e_th,e
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_)
983 e = e_th+dsqrt(e_th**2+rho**2)
985 (one-(rho(ixo^s)/e(ixo^s))**2)
989 e = e_th+dsqrt(e_th**2+w(ixo^s,
rho_)**2)
991 (one-(w(ixo^s,
rho_)/e(ixo^s))**2)
1000 subroutine srhd_get_smallvalues_eos
1004 double precision :: lsmalle,lsmallp,lsmallrho
1011 call mpistop(
"must set finite values small-density/pressure for small value treatments")
1021 gamma_1*lsmallrho*(lsmallrho/lsmalle))
1030 end subroutine srhd_get_smallvalues_eos
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)
1039 double precision,
dimension(ixO^S) :: e_th,e
1044 e = e_th+dsqrt(e_th**2+rho**2)
1053 {
do ix^db= ixo^lim^db\}
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)')
1064 {
do ix^db= ixo^lim^db\}
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)
1084 e = (rhoh+dsqrt(rhoh**2+(
srhd_gamma**2-one)*rho**2)) &
1086 p(ixo^s) = half*
gamma_1* (e-rho*(rho/e))
1092 {
do ix^db= ixo^lim^db\}
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)')
1105 {
do ix^db= ixo^lim^db\}
1113 end subroutine srhd_get_pressure_eos
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)
1123 double precision :: p(ixi^s)
1124 double precision :: e(ixo^s)
1127 call srhd_get_pressure_eos(ixi^
l,ixo^
l,rho,rhoh,p,e)
1130 +
gamma_1*(rho(ixo^s)/e(ixo^s))**2))&
1131 /(2.0d0*rhoh(ixo^s))
1133 csound2(ixo^s)=
srhd_gamma*p(ixo^s)/rhoh(ixo^s)
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)')
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
1153 if(csound2(ix^
d)<=0.0d0)
then
1159 end subroutine srhd_get_csound2_eos
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)
1169 double precision :: e(ixo^s)
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)')
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
1199 if(csound2(ix^
d)<=0.0d0)
then
1205 end subroutine srhd_get_csound2_prim_eos
1208 subroutine con2prim_eos(lfac,xi,myd,myssqr,mytau,ierror)
1211 double precision,
intent(in) :: myd, myssqr, mytau
1212 double precision,
intent(inout) :: lfac, xi
1213 integer,
intent(inout) :: ierror
1216 double precision:: f,df,lfacl
1220 d = myd;
ssqr = myssqr;
tau = mytau;
1226 call funcd_eos(xi,f,df,lfacl,
d,
ssqr,
tau,ierror)
1227 if (ierror == 0 .and. dabs(f/df)<
absaccnr)
then
1235 write(*,*)
'entering con2prim_eos with xi=',xi
1244 call con2primhydro_eos(lfac,xi,
d,
ssqr,
tau,ierror)
1246 end subroutine con2prim_eos
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
1254 double precision :: dlfac
1255 double precision :: vsqr,p,dpdxi
1261 mylfac = one/dsqrt(one-vsqr)
1262 dlfac = -mylfac**3*ssqr/(xi**3)
1264 call funcpressure_eos(xi,mylfac,d,dlfac,p,dpdxi)
1275 end subroutine funcd_eos
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
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
1306 pmin=dsqrt(sqrs)/(one-
dmaxvel)-tau-d
1307 plabs=max(
minp,pmin)
1326 pcurrent=half*(pcurrent+pprev)
1334 xicurrent=tau+d+pcurrent
1346 v2=sqrs/xicurrent**2
1348 if(lfac2inv>zero)
then
1349 lfac=one/dsqrt(lfac2inv)
1361 s2overcubeg2rh=sqrs/(xicurrent**3)
1363 call funcenthalpy_eos(pcurrent,lfac2inv,d,sqrs,xicurrent,&
1364 s2overcubeg2rh,h,dhdp)
1366 ff=-xicurrent*lfac2inv + h
1367 df=- two*sqrs/xicurrent**2 + dhdp - lfac2inv
1369 if (ff*df==zero)
then
1379 if (ff*df>zero)
then
1382 pnew=max(pnew,plabs)
1386 pnew=min(pnew,prabs)
1392 er=two*dabs(dp)/(pnew+pcurrent)
1398 if((dabs(oldff2-ff) < 1.0d-8 .or. niiter >=
maxitnr-
maxitnr/20).and.&
1399 ff * oldff1 < zero .and. dabs(ff)>
absaccnr)
then
1402 if(n2it<=3) pcurrent=half*(pnew+pcurrent)
1406 pcurrent=half*(pleft+pright)
1409 xicurrent=tau+d+pcurrent
1410 v2=sqrs/xicurrent**2
1413 if(lfac2inv>zero)
then
1414 lfac=one/dsqrt(lfac2inv)
1422 call bisection_enthalpy_eos(pnew,lfac2inv,d,xicurrent,h)
1423 nff=-xicurrent*lfac2inv + h
1426 if(ff * nff < zero)
then
1432 pcurrent=half*(pleft+pright)
1436 if(2.0d0*dabs(pleft-pright)/(pleft+pright)<
absaccnr &
1470 if(pcurrent<
minp)
then
1477 v2=sqrs/xicurrent**2
1479 if(lfac2inv>zero)
then
1480 lfac=one/dsqrt(lfac2inv)
1486 end subroutine con2primhydro_eos
1490 subroutine funcpressure_eos(xicurrent,lfac,d,dlfacdxi,p,dpdxi)
1492 double precision,
intent(in) :: xicurrent,lfac,d,dlfacdxi
1493 double precision,
intent(out) :: p,dpdxi
1495 double precision :: rho,h,e,dhdxi,rhotoe
1496 double precision :: dpdchi,dedxi
1499 h=xicurrent/(lfac**2)
1501 e = (h+dsqrt(h**2+(
srhd_gamma**2-one)*rho**2)) &
1505 p = half*
gamma_1*(e-rho*rhotoe)
1507 dhdxi = one/(lfac**2)-2.0d0*xicurrent/(lfac**2)*dlfacdxi/lfac
1509 dedxi=(dhdxi+(h*dhdxi-(
srhd_gamma**2-one)*rho**2*dlfacdxi/lfac)&
1514 dpdxi=half*
gamma_1*(2.0d0*rho*rhotoe*dlfacdxi/lfac+&
1515 (one+rhotoe**2)*dedxi)
1517 end subroutine funcpressure_eos
1521 subroutine funcenthalpy_eos(pcurrent,lfac2inv,d,sqrs,xicurrent,dv2d2p,h,dhdp)
1523 double precision,
intent(in) :: pcurrent,lfac2inv,d,sqrs,xicurrent,dv2d2p
1524 double precision,
intent(out):: h,dhdp
1527 double precision:: rho,e_th,e,de_thdp,dedp
1529 rho=d*dsqrt(lfac2inv)
1531 e = (e_th + dsqrt(e_th**2+rho**2))
1537 dedp = de_thdp * (one+e_th/dsqrt(e_th**2+rho**2))&
1538 + d**2*dv2d2p/dsqrt(e_th**2+rho**2)
1541 gamma_1*(rho*(rho/e))*(-2.0d0*dv2d2p/lfac2inv+dedp/e))
1542 end subroutine funcenthalpy_eos
1546 subroutine bisection_enthalpy_eos(pcurrent,lfac2inv,d,xicurrent,h)
1548 double precision,
intent(in) :: pcurrent,lfac2inv,d,xicurrent
1549 double precision,
intent(out):: h
1552 double precision:: rho,e_th,e
1554 rho=d*dsqrt(lfac2inv)
1556 e = (e_th + dsqrt(e_th**2+rho**2))
1561 end subroutine bisection_enthalpy_eos
1564 subroutine con2prim(lfac,xi,myd,myssqr,mytau,ierror)
1567 double precision,
intent(in) :: myd, myssqr, mytau
1568 double precision,
intent(inout) :: lfac, xi
1569 integer,
intent(inout) :: ierror
1572 double precision:: f,df,lfacl
1576 d = myd;
ssqr = myssqr;
tau = mytau;
1582 call funcd(xi,f,df,lfacl,
d,
ssqr,
tau,ierror)
1583 if (ierror == 0 .and. dabs(f/df)<
absaccnr)
then
1591 write(*,*)
'entering con2prim with xi=',xi
1600 call con2primhydro(lfac,xi,
d,
ssqr,
tau,ierror)
1602 end subroutine con2prim
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
1610 double precision :: dlfac
1611 double precision :: vsqr,p,dpdxi
1617 mylfac = one/dsqrt(one-vsqr)
1618 dlfac = -mylfac**3*ssqr/(xi**3)
1620 call funcpressure(xi,mylfac,d,dlfac,p,dpdxi)
1631 end subroutine funcd
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
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
1662 pmin=dsqrt(sqrs)/(one-
dmaxvel)-tau-d
1663 plabs=max(
minp,pmin)
1682 pcurrent=half*(pcurrent+pprev)
1690 xicurrent=tau+d+pcurrent
1702 v2=sqrs/xicurrent**2
1704 if(lfac2inv>zero)
then
1705 lfac=one/dsqrt(lfac2inv)
1717 s2overcubeg2rh=sqrs/(xicurrent**3)
1719 call funcenthalpy(pcurrent,lfac2inv,d,sqrs,xicurrent,&
1720 s2overcubeg2rh,h,dhdp)
1722 ff=-xicurrent*lfac2inv + h
1723 df=- two*sqrs/xicurrent**2 + dhdp - lfac2inv
1725 if (ff*df==zero)
then
1735 if (ff*df>zero)
then
1738 pnew=max(pnew,plabs)
1742 pnew=min(pnew,prabs)
1748 er=two*dabs(dp)/(pnew+pcurrent)
1754 if((dabs(oldff2-ff) < 1.0d-8 .or. niiter >=
maxitnr-
maxitnr/20).and.&
1755 ff * oldff1 < zero .and. dabs(ff)>
absaccnr)
then
1758 if(n2it<=3) pcurrent=half*(pnew+pcurrent)
1762 pcurrent=half*(pleft+pright)
1765 xicurrent=tau+d+pcurrent
1766 v2=sqrs/xicurrent**2
1769 if(lfac2inv>zero)
then
1770 lfac=one/dsqrt(lfac2inv)
1778 call bisection_enthalpy(pnew,lfac2inv,d,xicurrent,h)
1779 nff=-xicurrent*lfac2inv + h
1782 if(ff * nff < zero)
then
1788 pcurrent=half*(pleft+pright)
1792 if(2.0d0*dabs(pleft-pright)/(pleft+pright)<
absaccnr &
1826 if(pcurrent<
minp)
then
1833 v2=sqrs/xicurrent**2
1835 if(lfac2inv>zero)
then
1836 lfac=one/dsqrt(lfac2inv)
1842 end subroutine con2primhydro
1846 subroutine funcpressure(xicurrent,lfac,d,dlfacdxi,p,dpdxi)
1848 double precision,
intent(in) :: xicurrent,lfac,d,dlfacdxi
1849 double precision,
intent(out) :: p,dpdxi
1851 double precision :: rho,h,e,dhdxi,rhotoe
1852 double precision :: dpdchi,dedxi
1855 h=xicurrent/(lfac**2)
1860 dpdxi = dpdchi * one/lfac**2
1862 if (dlfacdxi /= 0.0d0) &
1863 dpdxi = dpdxi + dpdchi * ((d*lfac-2.0d0*xicurrent)/lfac**3) * dlfacdxi
1865 end subroutine funcpressure
1869 subroutine funcenthalpy(pcurrent,lfac2inv,d,sqrs,xicurrent,dv2d2p,h,dhdp)
1871 double precision,
intent(in) :: pcurrent,lfac2inv,d,sqrs,xicurrent,dv2d2p
1872 double precision,
intent(out):: h,dhdp
1875 double precision:: rho,e_th,e,de_thdp,dedp
1877 rho=d*dsqrt(lfac2inv)
1880 end subroutine funcenthalpy
1884 subroutine bisection_enthalpy(pcurrent,lfac2inv,d,xicurrent,h)
1886 double precision,
intent(in) :: pcurrent,lfac2inv,d,xicurrent
1887 double precision,
intent(out):: h
1890 double precision:: rho,e_th,e
1892 rho=d*dsqrt(lfac2inv)
1896 end subroutine bisection_enthalpy
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.
Module for physical and numeric constants.
Module with geometry-related routines (e.g., divergence, curl)
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 small_pressure
double precision unit_time
Physical scaling factor for time.
double precision unit_density
Physical scaling factor for density.
integer, parameter unitpar
file handle for IO
double precision unit_mass
Physical scaling factor for mass.
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision unit_numberdensity
Physical scaling factor for number density.
double precision unit_pressure
Physical scaling factor for pressure.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
logical use_particles
Use particles module or not.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
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
double precision small_density
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
Module containing all the particle routines.
subroutine particles_init()
Initialize particle data and parameters.
This module defines the procedures of a physics module. It contains function pointers for the various...
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.