diff options
author | knarf <knarf@ac85fae7-cede-4708-beff-ae01c7fa1c26> | 2010-08-16 13:45:31 +0000 |
---|---|---|
committer | knarf <knarf@ac85fae7-cede-4708-beff-ae01c7fa1c26> | 2010-08-16 13:45:31 +0000 |
commit | c9f3deaaf486dcc02df03dd7b3e047524f7f5ca9 (patch) | |
tree | 68c0020a517f999fd6e4c41f04688e9c76da5d4c | |
parent | 22870a2eb56af28e5c5d9ef23e7b99fad0adb2ff (diff) |
add USE_EOS_OMNI checks
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/GRHydro_InitData/trunk@109 ac85fae7-cede-4708-beff-ae01c7fa1c26
-rw-r--r-- | src/GRHydro_C2P2C.F90 | 25 | ||||
-rw-r--r-- | src/GRHydro_Con2Prim.F90 | 25 | ||||
-rw-r--r-- | src/GRHydro_ReadConformalData.F90 | 30 | ||||
-rw-r--r-- | src/GRHydro_SimpleWave.F90 | 24 |
4 files changed, 98 insertions, 6 deletions
diff --git a/src/GRHydro_C2P2C.F90 b/src/GRHydro_C2P2C.F90 index 5808c1b..a2523a8 100644 --- a/src/GRHydro_C2P2C.F90 +++ b/src/GRHydro_C2P2C.F90 @@ -34,10 +34,12 @@ subroutine c2p2c(CCTK_ARGUMENTS) DECLARE_CCTK_ARGUMENTS DECLARE_CCTK_PARAMETERS +#if !USE_EOS_OMNI #ifdef _EOS_BASE_INC_ #undef _EOS_BASE_INC_ #endif #include "EOS_Base.inc" +#endif integer didit,i,j,k,nx,ny,nz CCTK_REAL det @@ -50,6 +52,21 @@ subroutine c2p2c(CCTK_ARGUMENTS) CCTK_INT C2P_failed logical epsnegative +#if USE_EOS_OMNI +! begin EOS Omni vars + integer :: n = 1 + integer :: poly_eoskey = 0 + integer :: keytemp = 0 + integer :: anyerr = 0 + integer :: keyerr(1) = 0 + real*8 :: xpress = 0.0d0 + real*8 :: xeps = 0.0d0 + real*8 :: xtemp = 0.0d0 + real*8 :: xye = 0.0d0 + poly_eoskey = GRHydro_poly_eoskey +! end EOS Omni vars +#endif + call CCTK_WARN(1,"This test works only with Ideal_Fluid EoS") nx = cctk_lsh(1) @@ -90,8 +107,16 @@ subroutine c2p2c(CCTK_ARGUMENTS) epsnegative = .false. GRHydro_rho_min = 1.e-10 +#if USE_EOS_OMNI + call EOS_Omni_press(poly_eoskey,keytemp,n,& + GRHydro_rho_min,1.0d0,xtemp,xye,pmin,keyerr,anyerr) + + call EOS_Omni_EpsFromPress(poly_eoskey,keytemp,n,& + GRHydro_rho_min,xeps,xtemp,xye,pmin,epsmin,keyerr,anyerr) +#else pmin = EOS_Pressure(GRHydro_eos_handle, GRHydro_rho_min, 1.0d0) epsmin = EOS_SpecificIntEnergy(GRHydro_eos_handle, GRHydro_rho_min, pmin) +#endif C2P_failed = 0 write(*,*) 'C2P2C test: initial values.' diff --git a/src/GRHydro_Con2Prim.F90 b/src/GRHydro_Con2Prim.F90 index 6c22850..5172d29 100644 --- a/src/GRHydro_Con2Prim.F90 +++ b/src/GRHydro_Con2Prim.F90 @@ -52,10 +52,27 @@ subroutine GRHydro_con2primtest(CCTK_ARGUMENTS) CCTK_INT C2P_failed logical epsnegative +#ifndef USE_EOS_OMNI #ifdef _EOS_BASE_INC_ #undef _EOS_BASE_INC_ #endif #include "EOS_Base.inc" +#endif + +#if USE_EOS_OMNI +! begin EOS Omni vars + integer :: n = 1 + integer :: poly_eoskey = 0 + integer :: keytemp = 0 + integer :: anyerr = 0 + integer :: keyerr(1) = 0 + real*8 :: xpress = 0.0d0 + real*8 :: xeps = 0.0d0 + real*8 :: xtemp = 0.0d0 + real*8 :: xye = 0.0d0 + poly_eoskey = GRHydro_poly_eoskey +! end EOS Omni vars +#endif call CCTK_WARN(1,"For this test, remember to use a polytropic EoS and to set eos_gamma = 2.0 and eos_k = 100.0") @@ -92,8 +109,16 @@ subroutine GRHydro_con2primtest(CCTK_ARGUMENTS) epsnegative = .false. +#if USE_EOS_OMNI + call EOS_Omni_press(poly_eoskey,keytemp,n,& + GRHydro_rho_min,1.0d0,xtemp,xye,pmin,keyerr,anyerr) + + call EOS_Omni_EpsFromPress(poly_eoskey,keytemp,n,& + GRHydro_rho_min,xeps,xtemp,xye,pmin,epsmin,keyerr,anyerr) +#else pmin = EOS_Pressure(GRHydro_eos_handle, GRHydro_rho_min, 1.0d0) epsmin = EOS_SpecificIntEnergy(GRHydro_eos_handle, GRHydro_rho_min, pmin) +#endif C2P_failed = 0 write(*,*) 'Con2Prim test: converting to primitive variables' diff --git a/src/GRHydro_ReadConformalData.F90 b/src/GRHydro_ReadConformalData.F90 index 5829d80..4513f47 100644 --- a/src/GRHydro_ReadConformalData.F90 +++ b/src/GRHydro_ReadConformalData.F90 @@ -51,9 +51,23 @@ subroutine GRHydro_ReadConformalData(CCTK_ARGUMENTS) CCTK_INT :: i,j,k,handle,ierr - CCTK_REAL :: eos_k, eos_gamma, rho_min, det - - +! CCTK_REAL :: eos_k, eos_gamma + CCTK_REAL :: rho_min, det + +#if USE_EOS_OMNI +! begin EOS Omni vars + integer :: n = 1 + integer :: poly_eoskey = 0 + integer :: keytemp = 0 + integer :: anyerr = 0 + integer :: keyerr(1) = 0 + real*8 :: xpress = 0.0d0 + real*8 :: xeps = 0.0d0 + real*8 :: xtemp = 0.0d0 + real*8 :: xye = 0.0d0 + poly_eoskey = GRHydro_poly_eoskey +! end EOS Omni vars +#endif ! only gxx has been read in; copy it into gyy and gzz as well gyy = gxx @@ -131,8 +145,8 @@ subroutine GRHydro_ReadConformalData(CCTK_ARGUMENTS) ! if (handle < 0) call CCTK_WARN(0,"For this hack you need to compile with EOS_Polytrope") - eos_k = EOS_Pressure(GRHydro_eos_handle, 1.0, 1.0) - eos_gamma = 1.0 + eos_k / EOS_SpecificIntEnergy(GRHydro_eos_handle,1.0,1.0) +! eos_k = EOS_Pressure(GRHydro_eos_handle, 1.0, 1.0) +! eos_gamma = 1.0 + eos_k / EOS_SpecificIntEnergy(GRHydro_eos_handle,1.0,1.0) ! press = eos_k * rho**eos_gamma @@ -140,7 +154,13 @@ subroutine GRHydro_ReadConformalData(CCTK_ARGUMENTS) do i=1,cctk_lsh(1) do j=1,cctk_lsh(2) do k=1,cctk_lsh(3) +#ifdef USE_EOS_OMNI + call EOS_Omni_EpsFromPress(poly_eoskey,keytemp,n,& + rho(i,j,k),xeps,xtemp,xye,press(i,j,k),eps(i,j,k),keyerr,anyerr) + +#else eps(i,j,k) = EOS_SpecificIntEnergy(GRHydro_eos_handle,rho(i,j,k),press(i,j,k)) +#endif end do end do end do diff --git a/src/GRHydro_SimpleWave.F90 b/src/GRHydro_SimpleWave.F90 index 829170f..d1f10d4 100644 --- a/src/GRHydro_SimpleWave.F90 +++ b/src/GRHydro_SimpleWave.F90 @@ -46,14 +46,31 @@ subroutine GRHydro_SimpleWave(CCTK_ARGUMENTS) DECLARE_CCTK_PARAMETERS DECLARE_CCTK_FUNCTIONS +#ifndef USE_EOS_OMNI #ifdef _EOS_BASE_INC_ #undef _EOS_BASE_INC_ #endif #include "EOS_Base.inc" - +#endif + CCTK_INT :: i, j, k, nx, ny, nz CCTK_REAL :: dr, k1, k2, k3, k4, in_data, old_data, source_data, new_data, c_0, det, pi +#if USE_EOS_OMNI +! begin EOS Omni vars + integer :: n = 1 + integer :: poly_eoskey = 0 + integer :: keytemp = 0 + integer :: anyerr = 0 + integer :: keyerr(1) = 0 + real*8 :: xpress = 0.0d0 + real*8 :: xeps = 0.0d0 + real*8 :: xtemp = 0.0d0 + real*8 :: xye = 0.0d0 + poly_eoskey = GRHydro_poly_eoskey +! end EOS Omni vars +#endif + call CCTK_INFO("Setting initial data for a simple wave as Anile Miller Motta") call CCTK_WARN(1, "The simple-wave initial-data routine works only for unigrid and on one node.") @@ -133,7 +150,12 @@ subroutine GRHydro_SimpleWave(CCTK_ARGUMENTS) eps(i,1,1) = rho_abs_min**(1.d0/3.d0) velx(i,1,1) = 0.d0 w_lorentz(i,1,1) = 1.d0 +#if USE_EOS_OMNI + call EOS_Omni_press(poly_eoskey,keytemp,n,& + rho(i,1,1),1.0d0,xtemp,xye,press(i,1,1),keyerr,anyerr) +#else press(i,1,1) = EOS_Pressure(GRHydro_polytrope_handle, rho(i,1,1), 1.0d0) +#endif ! polytrope only (initial data) end if |