diff options
author | bmundim <bmundim@ac85fae7-cede-4708-beff-ae01c7fa1c26> | 2011-04-28 01:14:57 +0000 |
---|---|---|
committer | bmundim <bmundim@ac85fae7-cede-4708-beff-ae01c7fa1c26> | 2011-04-28 01:14:57 +0000 |
commit | d01d6ce889030ae0322409021d5c004b3fbcd649 (patch) | |
tree | cd4cae7d23bdb6d92615d6d7d0f1430bfc43f042 | |
parent | 58030c6ec544958da67c85831aa0e5024a5eff54 (diff) |
A few fixes on the interface calls and related bugs.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/GRHydro_InitData/trunk@127 ac85fae7-cede-4708-beff-ae01c7fa1c26
-rw-r--r-- | src/GRHydro_CylindricalExplosionM.F90 | 18 | ||||
-rw-r--r-- | src/GRHydro_MonopoleM.F90 | 39 | ||||
-rw-r--r-- | src/GRHydro_P2C2PM.F90 | 2 | ||||
-rw-r--r-- | src/GRHydro_ShockTubeM.F90 | 15 |
4 files changed, 45 insertions, 29 deletions
diff --git a/src/GRHydro_CylindricalExplosionM.F90 b/src/GRHydro_CylindricalExplosionM.F90 index d06f777..bc8ab61 100644 --- a/src/GRHydro_CylindricalExplosionM.F90 +++ b/src/GRHydro_CylindricalExplosionM.F90 @@ -22,6 +22,9 @@ #define Bvecx(i,j,k) Bvec(i,j,k,1) #define Bvecy(i,j,k) Bvec(i,j,k,2) #define Bvecz(i,j,k) Bvec(i,j,k,3) +#define Bconsx(i,j,k) Bcons(i,j,k,1) +#define Bconsy(i,j,k) Bcons(i,j,k,2) +#define Bconsz(i,j,k) Bcons(i,j,k,3) /*@@ @@ -117,25 +120,28 @@ subroutine GRHydro_cylindricalexplosionM(CCTK_ARGUMENTS) call Prim2ConPolyM(GRHydro_eos_handle,gxx(i,j,k),gxy(i,j,k),& gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),& det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),& - tau(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),rho(i,j,k),& + tau(i,j,k),Bconsx(i,j,k),Bconsy(i,j,k),Bconsz(i,j,k),rho(i,j,k),& velx(i,j,k),vely(i,j,k),velz(i,j,k),& - eps(i,j,k),press(i,j,k),w_lorentz(i,j,k)) + eps(i,j,k),press(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),& + w_lorentz(i,j,k)) else call Prim2ConGenM(GRHydro_eos_handle,gxx(i,j,k),gxy(i,j,k),& gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),& det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),& - tau(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),rho(i,j,k),& + tau(i,j,k),Bconsx(i,j,k),Bconsy(i,j,k),Bconsz(i,j,k),rho(i,j,k),& velx(i,j,k),vely(i,j,k),velz(i,j,k),& - eps(i,j,k),press(i,j,k),w_lorentz(i,j,k)) + eps(i,j,k),press(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),& + w_lorentz(i,j,k)) end if - enddo + + enddo enddo enddo densrhs = 0.d0 srhs = 0.d0 taurhs = 0.d0 - Bvecrhs = 0.d0 + Bconsrhs = 0.d0 return diff --git a/src/GRHydro_MonopoleM.F90 b/src/GRHydro_MonopoleM.F90 index a34a51c..6a6ea32 100644 --- a/src/GRHydro_MonopoleM.F90 +++ b/src/GRHydro_MonopoleM.F90 @@ -22,6 +22,9 @@ #define Bvecx(i,j,k) Bvec(i,j,k,1) #define Bvecy(i,j,k) Bvec(i,j,k,2) #define Bvecz(i,j,k) Bvec(i,j,k,3) +#define Bconsx(i,j,k) Bcons(i,j,k,1) +#define Bconsy(i,j,k) Bcons(i,j,k,2) +#define Bconsz(i,j,k) Bcons(i,j,k,3) /*@@ @routine GRHydro_MonopoleM @@ -77,22 +80,24 @@ subroutine GRHydro_MonopoleM(CCTK_ARGUMENTS) Bvecz(i,j,k)=0.0 det=SPATIAL_DETERMINANT(gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k)) - - if (CCTK_EQUALS(GRHydro_eos_type,"Polytype")) then - call Prim2ConPolyM(GRHydro_eos_handle,gxx(i,j,k),gxy(i,j,k),& - gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),& - det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),& - tau(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),rho(i,j,k),& - velx(i,j,k),vely(i,j,k),velz(i,j,k),& - eps(i,j,k),press(i,j,k),w_lorentz(i,j,k)) - else - call Prim2ConGenM(GRHydro_eos_handle,gxx(i,j,k),gxy(i,j,k),& - gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),& - det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),& - tau(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),rho(i,j,k),& - velx(i,j,k),vely(i,j,k),velz(i,j,k),& - eps(i,j,k),press(i,j,k),w_lorentz(i,j,k)) - end if + + if (CCTK_EQUALS(GRHydro_eos_type,"Polytype")) then + call Prim2ConPolyM(GRHydro_eos_handle,gxx(i,j,k),gxy(i,j,k),& + gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),& + det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),& + tau(i,j,k),Bconsx(i,j,k),Bconsy(i,j,k),Bconsz(i,j,k),rho(i,j,k),& + velx(i,j,k),vely(i,j,k),velz(i,j,k),& + eps(i,j,k),press(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),& + w_lorentz(i,j,k)) + else + call Prim2ConGenM(GRHydro_eos_handle,gxx(i,j,k),gxy(i,j,k),& + gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),& + det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),& + tau(i,j,k),Bconsx(i,j,k),Bconsy(i,j,k),Bconsz(i,j,k),rho(i,j,k),& + velx(i,j,k),vely(i,j,k),velz(i,j,k),& + eps(i,j,k),press(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),& + w_lorentz(i,j,k)) + end if enddo enddo @@ -101,7 +106,7 @@ subroutine GRHydro_MonopoleM(CCTK_ARGUMENTS) densrhs = 0.d0 srhs = 0.d0 taurhs = 0.d0 - Bvecrhs = 0.d0 + Bconsrhs = 0.d0 return diff --git a/src/GRHydro_P2C2PM.F90 b/src/GRHydro_P2C2PM.F90 index 15cacf0..1c1a9c7 100644 --- a/src/GRHydro_P2C2PM.F90 +++ b/src/GRHydro_P2C2PM.F90 @@ -131,7 +131,7 @@ subroutine p2c2pm(CCTK_ARGUMENTS) call Prim2ConGenM(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, & - bconsx_send,bcibsy_send,bconsz_send, & + bconsx_send,bconsy_send,bconsz_send, & rho_send(1),velx_send,vely_send,velz_send,eps_send(1),press_send(1), & bvcx_send,bvcy_send,bvcz_send,w_lorentz_send) diff --git a/src/GRHydro_ShockTubeM.F90 b/src/GRHydro_ShockTubeM.F90 index 72bac3e..a87c528 100644 --- a/src/GRHydro_ShockTubeM.F90 +++ b/src/GRHydro_ShockTubeM.F90 @@ -22,6 +22,9 @@ #define Bvecx(i,j,k) Bvec(i,j,k,1) #define Bvecy(i,j,k) Bvec(i,j,k,2) #define Bvecz(i,j,k) Bvec(i,j,k,3) +#define Bconsx(i,j,k) Bcons(i,j,k,1) +#define Bconsy(i,j,k) Bcons(i,j,k,2) +#define Bconsz(i,j,k) Bcons(i,j,k,3) #define OOSQRT2 (0.7071067811865475244008442) #define OOSQRT3 (0.5773502691896257645091489) @@ -611,16 +614,18 @@ subroutine GRHydro_shocktubeM(CCTK_ARGUMENTS) call Prim2ConPolyM(GRHydro_eos_handle,gxx(i,j,k),gxy(i,j,k),& gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),& det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),& - tau(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),rho(i,j,k),& + tau(i,j,k),Bconsx(i,j,k),Bconsy(i,j,k),Bconsz(i,j,k),rho(i,j,k),& velx(i,j,k),vely(i,j,k),velz(i,j,k),& - eps(i,j,k),press(i,j,k),w_lorentz(i,j,k)) + eps(i,j,k),press(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),& + w_lorentz(i,j,k)) else call Prim2ConGenM(GRHydro_eos_handle,gxx(i,j,k),gxy(i,j,k),& gxz(i,j,k),gyy(i,j,k),gyz(i,j,k),gzz(i,j,k),& det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),& - tau(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),rho(i,j,k),& + tau(i,j,k),Bconsx(i,j,k),Bconsy(i,j,k),Bconsz(i,j,k),rho(i,j,k),& velx(i,j,k),vely(i,j,k),velz(i,j,k),& - eps(i,j,k),press(i,j,k),w_lorentz(i,j,k)) + eps(i,j,k),press(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(i,j,k),& + w_lorentz(i,j,k)) end if enddo enddo @@ -629,7 +634,7 @@ subroutine GRHydro_shocktubeM(CCTK_ARGUMENTS) densrhs = 0.d0 srhs = 0.d0 taurhs = 0.d0 - Bvecrhs = 0.d0 + Bconsrhs = 0.d0 return |