aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
Diffstat (limited to 'src')
-rw-r--r--src/DistortedBHIVP.F2
-rw-r--r--src/Stab3d.F9
2 files changed, 8 insertions, 3 deletions
diff --git a/src/DistortedBHIVP.F b/src/DistortedBHIVP.F
index b76b184..4913781 100644
--- a/src/DistortedBHIVP.F
+++ b/src/DistortedBHIVP.F
@@ -620,7 +620,7 @@ c Compute eta,q,phi at the each points of cartesian grid
interp_coords(1) = CCTK_PointerTo(abseta)
interp_coords(2) = CCTK_PointerTo(q)
- interp_coords(2) = CCTK_PointerTo(phi)
+ interp_coords(3) = CCTK_PointerTo(phi)
in_array_dims(1) = ne
in_array_dims(2) = nq
diff --git a/src/Stab3d.F b/src/Stab3d.F
index f8ca5e0..aad61b6 100644
--- a/src/Stab3d.F
+++ b/src/Stab3d.F
@@ -97,12 +97,14 @@ c
real*8 x(im*jm*km),r(im*jm*km)
c Local variables
integer :: i,j,k,kk
- real*8 :: p(im*jm*km),Ap(im*jm*km),w(im*jm*km),As(im*jm*km)
+ real*8, allocatable :: p(:),Ap(:),w(:),As(:)
real*8 :: omega, chi,chi1,chi2, delta, deltap, pp
*
***********************************************************************
*
-
+c
+ allocate(p(im*jm*km),Ap(im*jm*km),w(im*jm*km),As(im*jm*km))
+
do i = 1,im*jm*km
p(i) = 0.
Ap(i) = 0.
@@ -237,6 +239,9 @@ c rnorm = sum(r**2)
rnorm=sqrt(rnorm)
if (rnorm .gt. tol) goto 1
+c
+ deallocate(p,Ap,w,As)
+c
return
end
c