diff options
-rw-r--r-- | src/GRHydro_Prim2Con.F90 | 9 | ||||
-rw-r--r-- | src/GRHydro_UpdateMask.F90 | 12 |
2 files changed, 13 insertions, 8 deletions
diff --git a/src/GRHydro_Prim2Con.F90 b/src/GRHydro_Prim2Con.F90 index bbbab06..805e486 100644 --- a/src/GRHydro_Prim2Con.F90 +++ b/src/GRHydro_Prim2Con.F90 @@ -305,8 +305,8 @@ subroutine prim2con_hot(handle, GRHydro_reflevel, cctk_iteration, ii, jj, kk, & ! This is a way of recovering even on finer refinement levels: ! Use the average temperature at the interface instead of the ! reconstructed specific internal energy. - !$OMP CRITICAL if(GRHydro_eos_hot_prim2con_warn.ne.0) then + !$OMP CRITICAL call CCTK_WARN(1,"EOS error in prim2con_hot: NOW using averaged temp!") write(warnline,"(i8,3i5,1P10E15.6)") cctk_iteration,ii,jj,kk,x,y,z,r call CCTK_WARN(1,warnline) @@ -316,8 +316,8 @@ subroutine prim2con_hot(handle, GRHydro_reflevel, cctk_iteration, ii, jj, kk, & call CCTK_WARN(1,warnline) write(warnline,"(A10,i5)") "reflevel: ", GRHydro_reflevel call CCTK_WARN(1,warnline) - endif !$OMP END CRITICAL + endif keytemp=1 temp = temp0 call EOS_Omni_press(handle,keytemp,GRHydro_eos_rf_prec,n,& @@ -441,8 +441,9 @@ subroutine Primitive2ConservativeCells(CCTK_ARGUMENTS) det = SPATIAL_DETERMINANT(g11(i,j,k),g12(i,j,k),g13(i,j,k), \ g22(i,j,k),g23(i,j,k),g33(i,j,k)) - call prim2con_hot(GRHydro_eos_handle,GRHydro_reflevel,i,j,k,& - x(i,j,k),y(i,j,k),z(i,j,k),g11(i,j,k),& + call prim2con_hot(GRHydro_eos_handle,GRHydro_reflevel,cctk_iteration,& + i,j,k,& + x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k),g11(i,j,k),& g12(i,j,k),g13(i,j,k),& g22(i,j,k),g23(i,j,k),g33(i,j,k),& det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),& diff --git a/src/GRHydro_UpdateMask.F90 b/src/GRHydro_UpdateMask.F90 index 02934d7..ae13d52 100644 --- a/src/GRHydro_UpdateMask.F90 +++ b/src/GRHydro_UpdateMask.F90 @@ -303,7 +303,8 @@ subroutine GRHydro_AtmosphereReset(CCTK_ARGUMENTS) press(i,j,k),keyerr,anyerr) call prim2con_hot(GRHydro_eos_handle, GRHydro_reflevel,& - i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11(i,j,k),g12(i,j,k),& + cctk_iteration,i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k),& + g11(i,j,k),g12(i,j,k),& g13(i,j,k),g22(i,j,k),g23(i,j,k),g33(i,j,k), & det,dens(i,j,k),scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), & tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), & @@ -465,7 +466,8 @@ subroutine GRHydro_InitialAtmosphereReset(CCTK_ARGUMENTS) press(i,j,k),keyerr,anyerr) call prim2con_hot(GRHydro_eos_handle, GRHydro_reflevel,& - i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11(i,j,k),g12(i,j,k),& + cctk_iteration,i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k),& + g11(i,j,k),g12(i,j,k),& g13(i,j,k),g22(i,j,k),g23(i,j,k),g33(i,j,k), & det,dens(i,j,k),scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), & tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), & @@ -505,7 +507,8 @@ subroutine GRHydro_InitialAtmosphereReset(CCTK_ARGUMENTS) press_p(i,j,k),keyerr,anyerr) call prim2con_hot(GRHydro_eos_handle, GRHydro_reflevel,& - i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11_p(i,j,k),g12_p(i,j,k),& + cctk_iteration,i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k),& + g11_p(i,j,k),g12_p(i,j,k),& g13_p(i,j,k),g22_p(i,j,k),g23_p(i,j,k),g33_p(i,j,k), & det,dens_p(i,j,k),scon_p(i,j,k,1), scon_p(i,j,k,2), scon_p(i,j,k,3), & tau_p(i,j,k), rho_p(i,j,k), velx_p(i,j,k), vely_p(i,j,k), & @@ -545,7 +548,8 @@ subroutine GRHydro_InitialAtmosphereReset(CCTK_ARGUMENTS) rho_p_p(i,j,k),eps_p_p(i,j,k),temperature_p_p(i,j,k),y_e_p_p(i,j,k),& press_p_p(i,j,k),keyerr,anyerr) call prim2con_hot(GRHydro_eos_handle, GRHydro_reflevel,& - i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),g11_p_p(i,j,k),g12_p_p(i,j,k),& + cctk_iteration,i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),r(i,j,k),& + g11_p_p(i,j,k),g12_p_p(i,j,k),& g13_p_p(i,j,k),g22_p_p(i,j,k),g23_p_p(i,j,k),g33_p_p(i,j,k), & det,dens_p_p(i,j,k),scon_p_p(i,j,k,1), scon_p_p(i,j,k,2), scon_p_p(i,j,k,3), & tau_p_p(i,j,k), rho_p_p(i,j,k), velx_p_p(i,j,k), vely_p_p(i,j,k), & |