aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_C2P2C.F90
diff options
context:
space:
mode:
Diffstat (limited to 'src/GRHydro_C2P2C.F90')
-rw-r--r--src/GRHydro_C2P2C.F9059
1 files changed, 27 insertions, 32 deletions
diff --git a/src/GRHydro_C2P2C.F90 b/src/GRHydro_C2P2C.F90
index a1cc36d..731dbb8 100644
--- a/src/GRHydro_C2P2C.F90
+++ b/src/GRHydro_C2P2C.F90
@@ -10,6 +10,7 @@
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
+#include "cctk_Functions.h"
/*@@
@routine c2p2c
@@ -33,34 +34,26 @@ subroutine c2p2c(CCTK_ARGUMENTS)
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
+ DECLARE_CCTK_FUNCTIONS
- integer didit,i,j,k,nx,ny,nz
- CCTK_REAL det
- CCTK_REAL uxx,uxy,uxz,uyy,uyz,uzz
- CCTK_REAL gxx_send,gxy_send,gxz_send,gyy_send,gyz_send,gzz_send
- CCTK_REAL dens_send,sx_send,sy_send,sz_send,tau_send
- CCTK_REAL rho_send,velx_send,vely_send,velz_send,eps_send
- CCTK_REAL press_send,w_lorentz_send,x_send,y_send,z_send,r_send
- CCTK_REAL pmin, epsmin
- CCTK_INT C2P_failed
- logical epsnegative
+ CCTK_REAL :: det
+ CCTK_REAL :: uxx,uxy,uxz,uyy,uyz,uzz
+ CCTK_REAL :: gxx_send,gxy_send,gxz_send,gyy_send,gyz_send,gzz_send
+ CCTK_REAL :: dens_send,sx_send,sy_send,sz_send,tau_send
+ CCTK_REAL :: rho_send,velx_send,vely_send,velz_send,eps_send
+ CCTK_REAL :: press_send,w_lorentz_send,x_send,y_send,z_send,r_send
+ CCTK_REAL :: C2P_failed
+ CCTK_INT :: epsnegative
! begin EOS Omni vars
- integer :: n = 1
- 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
+ CCTK_REAL :: pmin(1), epsmin(1), epsval(1)
+ CCTK_INT :: n,keytemp,anyerr,keyerr(1)
+ CCTK_REAL :: xpress(1),xtemp(1),xye(1),xeps(1),xrho(1)
+ n=1;keytemp=0;anyerr=0;keyerr(1)=0
+ xpress(1)=0.0d0;xtemp(1)=0.0d0;xye(1)=0.0d0;xeps(1)=0.0d0
! end EOS Omni vars
call CCTK_WARN(1,"This test works only with Ideal_Fluid EoS")
-
- nx = cctk_lsh(1)
- ny = cctk_lsh(2)
- nz = cctk_lsh(3)
x_send = 0.0d0
y_send = 0.0d0
@@ -93,17 +86,18 @@ subroutine c2p2c(CCTK_ARGUMENTS)
press_send = 6.666666666666667d-7
w_lorentz_send = 1.0d0
- epsnegative = .false.
+ epsnegative = 0
- GRHydro_rho_min = 1.e-10
+ xrho = 1.0d-10
+ epsval = 1.0d0
call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
- GRHydro_rho_min,1.0d0,xtemp,xye,pmin,keyerr,anyerr)
+ xrho,epsval,xtemp,xye,pmin,keyerr,anyerr)
call EOS_Omni_EpsFromPress(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
- GRHydro_rho_min,xeps,xtemp,xye,pmin,epsmin,keyerr,anyerr)
+ xrho,xeps,xtemp,xye,pmin,epsmin,keyerr,anyerr)
- C2P_failed = 0
+ C2P_failed = 0.0d0
write(*,*) 'C2P2C test: initial values.'
write(*,*) ' conservative variables: '
@@ -116,11 +110,11 @@ subroutine c2p2c(CCTK_ARGUMENTS)
write(*,*) ' W : ',w_lorentz_send
write(*,*) 'C2P2C test: getting the associated primitive variables.'
- call Con2Prim_pt(GRHydro_eos_handle,dens_send,sx_send,sy_send,sz_send, &
+ call Con2PrimGen(GRHydro_eos_handle,dens_send,sx_send,sy_send,sz_send, &
tau_send,rho_send,velx_send,vely_send,velz_send, &
eps_send,press_send,w_lorentz_send, &
uxx,uxy,uxz,uyy,uyz,uzz,det,x_send,y_send,z_send,r_send,&
- epsnegative,GRHydro_rho_min, pmin, epsmin, GRHydro_init_data_reflevel,C2P_failed)
+ epsnegative,xrho(1),pmin(1),epsmin(1),GRHydro_init_data_reflevel,C2P_failed)
write(*,*) 'C2P2C test: the primitive variables are'
write(*,*) ' primitive variables: '
@@ -133,9 +127,10 @@ subroutine c2p2c(CCTK_ARGUMENTS)
write(*,*) ' W : ',w_lorentz_send
write(*,*) 'C2P2C test: converting back to conserved variables.'
- call prim2con(GRHydro_eos_handle,gxx_send, gxy_send, gxz_send, gyy_send, gyz_send, gzz_send, det, &
- dens_send, sx_send, sy_send, sz_send, tau_send, rho_send, &
- velx_send, vely_send, velz_send, eps_send, press_send, w_lorentz_send)
+ call Prim2ConGen(GRHydro_eos_handle, &
+ gxx_send,gxy_send,gxz_send,gyy_send,gyz_send,gzz_send,det, &
+ dens_send,sx_send,sy_send,sz_send,tau_send,rho_send, &
+ velx_send,vely_send,velz_send,eps_send,press_send,w_lorentz_send)
write(*,*) 'C2P2C test: the conserved variables are'
write(*,*) ' conservative variables: '