aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@cct.lsu.edu>2006-12-12 16:02:00 +0000
committerErik Schnetter <schnetter@cct.lsu.edu>2006-12-12 16:02:00 +0000
commit1e08caa8a5e95a49273ce7fd90c2c32721e0fc87 (patch)
treeabdeef4a3e3708b63025ed53e8b5a9e64f8dccac /Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77
parent27164fb4a8f7108196cc4122ccbf2d238e68d09a (diff)
CarpetLib: Replace CHKIDX macros with calls to checkindex
darcs-hash:20061212160245-dae7b-19df81c29d911d9c77ae0aaa99ae999a0f6d27c9.gz
Diffstat (limited to 'Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77')
-rw-r--r--Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F7726
1 files changed, 9 insertions, 17 deletions
diff --git a/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77 b/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77
index 7ec6dc3a7..0bf91a371 100644
--- a/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77
+++ b/Carpet/CarpetLib/src/prolongate_3d_real8_2tl.F77
@@ -1,20 +1,10 @@
c -*-Fortran-*-
#include "cctk.h"
+#include "cctk_Parameters.h"
-#define CHKIDX(i,j,k, imax,jmax,kmax, where) \
- if ((i).lt.1 .or. (i).gt.(imax) \
- .or. (j).lt.1 .or. (j).gt.(jmax) \
- .or. (k).lt.1 .or. (k).gt.(kmax)) then &&\
- write (msg, '(a, " array index out of bounds: shape is (",i4,",",i4,",",i4,"), index is (",i4,",",i4,",",i4,")")') \
- (where), (imax), (jmax), (kmax), (i), (j), (k) &&\
- call CCTK_WARN (0, msg) &&\
- end if
-
-
-
subroutine prolongate_3d_real8_2tl (
$ src1, t1, src2, t2, srciext, srcjext, srckext,
$ dst, t, dstiext, dstjext, dstkext,
@@ -22,6 +12,8 @@ c -*-Fortran-*-
implicit none
+ DECLARE_CCTK_PARAMETERS
+
CCTK_REAL8 one
parameter (one = 1)
@@ -60,8 +52,6 @@ c bbox(:,3) is stride
CCTK_REAL8 res
integer d
- character msg*1000
-
do d=1,3
@@ -170,8 +160,9 @@ c Loop over fine region
fac = ifac(ii) * jfac(jj) * kfac(kk)
if (fac.ne.0) then
- CHKIDX (i0+ii, j0+jj, k0+kk, \
- srciext,srcjext,srckext, "source")
+ if (check_array_accesses.ne.0) then
+ call checkindex (i0+ii, j0+jj, k0+kk, 1,1,1, srciext,srcjext,srckext, "source")
+ end if
res = res
$ + fac * s1fac * src1(i0+ii, j0+jj, k0+kk)
$ + fac * s2fac * src2(i0+ii, j0+jj, k0+kk)
@@ -181,8 +172,9 @@ c Loop over fine region
end do
end do
- CHKIDX (dstioff+i+1, dstjoff+j+1, dstkoff+k+1, \
- dstiext,dstjext,dstkext, "destination")
+ if (check_array_accesses.ne.0) then
+ call checkindex (dstioff+i+1, dstjoff+j+1, dstkoff+k+1, 1,1,1, dstiext,dstjext,dstkext, "destination")
+ end if
dst (dstioff+i+1, dstjoff+j+1, dstkoff+k+1) = dstdiv * res
end do