aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorbentivegna <bentivegna@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-01-11 23:02:17 +0000
committerbentivegna <bentivegna@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-01-11 23:02:17 +0000
commit1a6a07d0ae8cecb53b0d0440a45bc682297bf6c1 (patch)
treed033a4a876fa6bbf9472e2c0883b3fa7cd1ecd32 /src
parent56522273ef0c1d3d2c4e6c1169a1c132ef5b9c34 (diff)
Deleted unused variables here and there.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@24 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src')
-rw-r--r--src/Whisky_Con2Prim.F908
-rw-r--r--src/Whisky_Eigenproblem.F9025
-rw-r--r--src/util/Riemann1d.f9029
3 files changed, 30 insertions, 32 deletions
diff --git a/src/Whisky_Con2Prim.F90 b/src/Whisky_Con2Prim.F90
index 53de9f5..6179048 100644
--- a/src/Whisky_Con2Prim.F90
+++ b/src/Whisky_Con2Prim.F90
@@ -866,10 +866,10 @@ subroutine Con2Prim_ptPolytype(handle, dens, sx, sy, sz, tau, rho, &
press, uxx, uxy, uxz, uyy, uyz, uzz, det, w_lorentz, x, &
y, z, r, whisky_rho_min
- CCTK_REAL s2, f, df, ftol, vlowx, vlowy, vlowz
- CCTK_INT count, i, handle, whisky_reflevel, whisky_C2P_failed
- CCTK_REAL udens, usx, usy, usz, rhoold, rhonew, epsold, epsnew, w2, &
- w2mhalf, enthalpy, denthalpy
+ CCTK_REAL s2, f, df, vlowx, vlowy, vlowz
+ CCTK_INT count, handle, whisky_reflevel, whisky_C2P_failed
+ CCTK_REAL udens, usx, usy, usz, rhoold, rhonew, epsold, epsnew, &
+ enthalpy, denthalpy
character(len=200) warnline
#ifdef _EOS_BASE_INC_
diff --git a/src/Whisky_Eigenproblem.F90 b/src/Whisky_Eigenproblem.F90
index 2abc3e7..e0a96ea 100644
--- a/src/Whisky_Eigenproblem.F90
+++ b/src/Whisky_Eigenproblem.F90
@@ -36,7 +36,7 @@
@@*/
subroutine eigenvalues(handle,rho,velx,vely,velz,eps,w_lorentz,&
- lam,gxx,gxy,gxz,gyy,gyz,gzz,u,det,alp,beta)
+ lam,gxx,gxy,gxz,gyy,gyz,gzz,u,alp,beta)
implicit none
@@ -44,7 +44,7 @@ subroutine eigenvalues(handle,rho,velx,vely,velz,eps,w_lorentz,&
CCTK_REAL rho,velx,vely,velz,eps,w_lorentz
CCTK_REAL lam(5)
- CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,det
+ CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz
CCTK_REAL alp,beta,u
CCTK_REAL cs2,one,two
@@ -114,7 +114,7 @@ end subroutine eigenvalues
@@*/
subroutine eigenproblem(handle,rho,velx,vely,velz,eps,&
- w_lorentz,gxx,gxy,gxz,gyy,gyz,gzz,u,det,&
+ w_lorentz,gxx,gxy,gxz,gyy,gyz,gzz,u,&
alp,beta,qdiff1,qdiff2,qdiff3,qdiff4,qdiff5,&
roeflux1,roeflux2,roeflux3,roeflux4,roeflux5)
@@ -126,7 +126,7 @@ subroutine eigenproblem(handle,rho,velx,vely,velz,eps,&
CCTK_REAL rho,velx,vely,velz,eps,w_lorentz
CCTK_REAL lam(5),p(5,5),q(5,5),dw(5),rflux(5)
- CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,u,det
+ CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,u
CCTK_REAL alp,beta,roeflux1,roeflux2,roeflux3,roeflux4,roeflux5
CCTK_REAL du(5),aa(5,5),qdiff1,qdiff2,qdiff3,qdiff4,qdiff5
@@ -484,7 +484,7 @@ end subroutine eigenproblem
@@*/
subroutine eigenproblem_leftright(handle,rho,velx,vely,velz,eps,&
- w_lorentz,gxx,gxy,gxz,gyy,gyz,gzz,u,det,&
+ w_lorentz,gxx,gxy,gxz,gyy,gyz,gzz,u,&
alp,beta,lambda,levec,revec)
USE Whisky_Scalars
@@ -495,7 +495,7 @@ subroutine eigenproblem_leftright(handle,rho,velx,vely,velz,eps,&
CCTK_REAL rho,velx,vely,velz,eps,w_lorentz
CCTK_REAL lambda(5),levec(5,5),revec(5,5)
- CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,u,det
+ CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,u
CCTK_REAL alp,beta
CCTK_REAL tmp1,tmp2
@@ -705,24 +705,23 @@ end subroutine eigenproblem_leftright
@@*/
subroutine eigenvalues_general(&
- rho,velx,vely,velz,eps,press,cs2,&
+ velx,vely,velz,cs2,&
lam,&
gxx,gxy,gxz,gyy,gyz,gzz,&
- u,det,alp,beta)
+ u,alp,beta)
implicit none
DECLARE_CCTK_PARAMETERS
- CCTK_REAL rho,velx,vely,velz,eps
+ CCTK_REAL velx,vely,velz
CCTK_REAL lam(5)
- CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,det
+ CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz
CCTK_REAL alp,beta,u
CCTK_REAL cs2,one,two
CCTK_REAL vlowx,vlowy,vlowz,v2,w
CCTK_REAL lam1,lam2,lam3,lamm,lamp,lamm_nobeta,lamp_nobeta
- CCTK_REAL press
character(len=256) NaN_WarnLine
@@ -793,7 +792,7 @@ subroutine eigenproblem_general(&
rho,velx,vely,velz,eps,&
press,cs2,dpdeps,&
gxx,gxy,gxz,gyy,gyz,gzz,&
- u,det,alp,beta,&
+ u,alp,beta,&
qdiff1,qdiff2,qdiff3,qdiff4,qdiff5,&
roeflux1,roeflux2,roeflux3,roeflux4,roeflux5)
@@ -805,7 +804,7 @@ subroutine eigenproblem_general(&
CCTK_REAL rho,velx,vely,velz,eps
CCTK_REAL lam(5),p(5,5),q(5,5),dw(5),rflux(5)
- CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,u,det
+ CCTK_REAL gxx,gxy,gxz,gyy,gyz,gzz,u
CCTK_REAL alp,beta,roeflux1,roeflux2,roeflux3,roeflux4,roeflux5
CCTK_REAL du(5),aa(5,5),qdiff1,qdiff2,qdiff3,qdiff4,qdiff5
diff --git a/src/util/Riemann1d.f90 b/src/util/Riemann1d.f90
index 5198f95..200b9a5 100644
--- a/src/util/Riemann1d.f90
+++ b/src/util/Riemann1d.f90
@@ -11,7 +11,7 @@ program Riemann1d
implicit none
- integer :: mn, i, j, n, iloop, step, test
+ integer :: mn, i, step, test
real(kind=8) :: tol, pmin, pmax, dvel1, dvel2, check, gamma, t, &
&gm1, mtol, xi, rad, vshockl, vshockr, csl, csls, csrs, &
&csr, vels, ul, uls, urs, ur, hl, hls, hrs, hr, &
@@ -155,7 +155,7 @@ program Riemann1d
guess = (0.5d0 - eps) * (left + right)
- call getp2(mtol, gamma, tol, left, right, starl, starr, guess, &
+ call getp2(gamma, tol, left, right, starl, starr, guess, &
&step)
else
@@ -238,7 +238,7 @@ program Riemann1d
end do
- call vshock(gamma, left, starl, starr, right, vshockl, vshockr)
+ call vshock(gamma, left, starl, starr, right)
ul = left(3) / gm1 / left(1)
uls = starl(3) / gm1 / starl(1)
@@ -382,7 +382,7 @@ contains
implicit none
real(kind=8), dimension(3) :: left, right, guess
- real(kind=8) :: dvel, pguess, gamma, gm1
+ real(kind=8) :: dvel, pguess, gamma
integer :: step
step = step + 1
@@ -413,7 +413,7 @@ contains
real(kind=8), dimension(3) :: known, guess
real(kind=8) :: gamma, gm1, cs, sign, a, b, c, h, e, ea, ua,&
- &ha, csa, k, sqgl1, v12, pp, n, d
+ &ha, csa, k, sqgl1, v12
gm1 = gamma - 1.d0
ua = known(3) / gm1 / known(1)
@@ -617,7 +617,6 @@ contains
real(kind=8), dimension(3) :: known
real(kind=8), dimension(4) :: out
integer, intent(in) :: sign
- integer :: step
real(kind=8) :: gm1, b, c, d, k, l, v, ua, ha, csa, cs2, &
&ocs2, fcs2, dfdcs2, tol
@@ -675,7 +674,7 @@ contains
implicit none
real(kind=8), dimension(3) :: left, right, guess
- real(kind=8) :: dvel, pguess, gamma, gm1, &
+ real(kind=8) :: dvel, pguess, gamma, &
&dfdp, dfdp1
integer :: step
@@ -708,10 +707,10 @@ contains
implicit none
real(kind=8), dimension(3) :: known, guess
- real(kind=8) :: gamma, gm1, cs, sign, a, b, c, h, u, j, ua, &
- &ha, csa, wa, vshock, wshock, k, sqgl1, dfdp, dhdp, dhordp, &
- &djdp, dwsojdp, dcsdp, dadp, d, e, dbdp, dcdp, drhodp, &
- &g, i, dvsdp, dwsdp, dv12dp, v12, dedp, ea
+ real(kind=8) :: gamma, gm1, cs, sign, a, b, c, h, ua, &
+ &ha, csa, wa, k, sqgl1, dfdp, dhdp, &
+ &dadp, e, dbdp, dcdp, &
+ &dv12dp, v12, dedp, ea
gm1 = gamma - 1.d0
ua = known(3) / gm1 / known(1)
@@ -778,13 +777,13 @@ contains
!!$Given an interval, computes the correct intermediate pressure and
!!$returns the intermediate state in star
- subroutine getp2(mtol, gamma, tol, left, right, starl, &
+ subroutine getp2(gamma, tol, left, right, starl, &
&starr, guess, step)
implicit none
real(kind=8), dimension(3) :: left, right, guess, starl, starr
- real(kind=8) :: tol, gamma, f, dfdp, mtol
+ real(kind=8) :: tol, gamma, f, dfdp
integer :: step
step = 0
@@ -831,12 +830,12 @@ contains
!!$Quick routine to find the shock velocities
- subroutine vshock(gamma, l, sl, sr, r, vl, vr)
+ subroutine vshock(gamma, l, sl, sr, r)
implicit none
real(kind=8), dimension(3) :: l, sl, sr, r
- real(kind=8) :: vl, vr, h, j, gamma, gm1, &
+ real(kind=8) :: h, j, gamma, gm1, &
&ha, wa2, a, b
gm1 = gamma - 1.d0