aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorallen <allen@c78560ca-4b45-4335-b268-5f3340f3cb52>1999-09-14 11:36:40 +0000
committerallen <allen@c78560ca-4b45-4335-b268-5f3340f3cb52>1999-09-14 11:36:40 +0000
commit04555d7fe6cfcfbdda5cb31f4c0e95afcba4f3e2 (patch)
tree79aa95a473a8977a5efb7337c129262b0b746a5a /src
parentc393c416f9bac40811e07f2ed95999969aa0af73 (diff)
Now calculated the base grid spacings
Changed the CCTK_INFO lines a tad git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/CartGrid3D/trunk@53 c78560ca-4b45-4335-b268-5f3340f3cb52
Diffstat (limited to 'src')
-rw-r--r--src/CartGrid3D.F7770
1 files changed, 36 insertions, 34 deletions
diff --git a/src/CartGrid3D.F77 b/src/CartGrid3D.F77
index a3e3830..d580037 100644
--- a/src/CartGrid3D.F77
+++ b/src/CartGrid3D.F77
@@ -26,8 +26,8 @@
CCTK_REAL x_origin,y_origin,z_origin
CCTK_REAL this_dx,this_dy,this_dz
CCTK_REAL xmax1,ymax1,zmax1
- CCTK_REAL lower,upper,olower,oupper
-c Joan was here
+ CCTK_REAL lowerx,upperx,lowery,uppery,lowerz,upperz
+ CCTK_REAL olower,oupper
character*80 infoline
iconv = 2**(cctk_convlevel-1)
@@ -68,9 +68,9 @@ c Grid spacing on coarsest grid
end if
c Grid spacing on this grid
- this_dx = coarse_dx/cctk_levfac(1)
- this_dy = coarse_dy/cctk_levfac(2)
- this_dz = coarse_dz/cctk_levfac(3)
+ this_dx = coarse_dx
+ this_dy = coarse_dy
+ this_dz = coarse_dz
c Minimum coordinate values on this grid
if (no_origin.eq.1) then
@@ -106,14 +106,9 @@ c dx,dy,dz on the coarsest grid of each GH
end if
c dx,dy,dz on the grid we are on
- this_dx = coarse_dx*iconv/cctk_levfac(1)
- this_dy = coarse_dy*iconv/cctk_levfac(2)
- this_dz = coarse_dz*iconv/cctk_levfac(3)
-
- write(infoline,'(1X,A,1X,3(A,G12.7,2X))')
- & 'Grid by range: Setting ','dx=>',this_dx,'dy=>',
- & this_dy,'dz=>',this_dz
- call CCTK_INFO(infoline)
+ this_dx = coarse_dx*iconv
+ this_dy = coarse_dy*iconv
+ this_dz = coarse_dz*iconv
c -----------------------------------------------------------
c
@@ -139,9 +134,9 @@ c Dx, Dy, Dx on the coarsest grid
end if
c dx, dy, dz on the grid we are on
- this_dx = coarse_dx*iconv/cctk_levfac(1)
- this_dy = coarse_dy*iconv/cctk_levfac(2)
- this_dz = coarse_dz*iconv/cctk_levfac(3)
+ this_dx = coarse_dx*iconv
+ this_dy = coarse_dy*iconv
+ this_dz = coarse_dz*iconv
if (CCTK_Equals(domain,'bitant') .eq. 1) then
@@ -215,20 +210,15 @@ c dx,dy,dz on the coarsest grid of each GH
coarse_dz = 1.d0/max(cctk_gsh(3)-1,1)
c dx,dy,dz on the grid we are on
- this_dx = coarse_dx*iconv/cctk_levfac(1)
- this_dy = coarse_dy*iconv/cctk_levfac(2)
- this_dz = coarse_dz*iconv/cctk_levfac(3)
+ this_dx = coarse_dx*iconv
+ this_dy = coarse_dy*iconv
+ this_dz = coarse_dz*iconv
c Special cases
if (cctk_gsh(1) .eq. 1) x_origin = 0.0D0
if (cctk_gsh(2) .eq. 1) y_origin = 0.0D0
if (cctk_gsh(3) .eq. 1) z_origin = 0.0D0
- write(infoline,'(1X,A,1X,3(A,G9.4,1X))')
- & 'Box grid: Setting ','dx=>',this_dx,'dy=>',this_dy,
- & 'dz=>',this_dz
- CALL CCTK_INFO(infoline)
-
end if
@@ -253,15 +243,28 @@ c -----------------------
cctk_origin_space(2) = y_origin
cctk_origin_space(3) = z_origin
- lower = x_origin
- upper = x_origin+this_dx*(cctk_gsh(1)-1)
- call CCTK_RegisterCoordRange(cctkGH,lower,upper,"x")
- lower = y_origin
- upper = y_origin+this_dy*(cctk_gsh(2)-1)
- call CCTK_RegisterCoordRange(cctkGH,lower,upper,"y")
- lower = z_origin
- upper = z_origin+this_dz*(cctk_gsh(3)-1)
- call CCTK_RegisterCoordRange(cctkGH,lower,upper,"z")
+ lowerx = x_origin
+ upperx = x_origin+this_dx*(cctk_gsh(1)-1)
+ call CCTK_RegisterCoordRange(cctkGH,lowerx,upperx,"x")
+ lowery = y_origin
+ uppery = y_origin+this_dy*(cctk_gsh(2)-1)
+ call CCTK_RegisterCoordRange(cctkGH,lowery,uppery,"y")
+ lowerz = z_origin
+ upperz = z_origin+this_dz*(cctk_gsh(3)-1)
+ call CCTK_RegisterCoordRange(cctkGH,lowerz,upperz,"z")
+
+ write(infoline,'(1X,3(A,G12.7,2X))')
+ & 'dx=>',cctk_delta_space(1),
+ & 'dy=>',cctk_delta_space(2),
+ & 'dz=>',cctk_delta_space(3)
+ call CCTK_INFO(infoline)
+
+ write(infoline,'(1X,3(A,F6.3,A,F6.3,A,2X))')
+ & 'x=>[',lowerx,',',upperx,']',
+ & 'y=>[',lowery,',',uppery,']',
+ & 'z=>[',lowerz,',',upperz,']'
+ call CCTK_INFO(infoline)
+
#ifdef CCTK_DEBUG
write(*,*)
@@ -272,7 +275,6 @@ c -----------------------
write(*,101) cctk_delta_space(1),cctk_delta_space(2),cctk_delta_space(3)
101 format(2X,"Dx,Dy,Dz on this grid ",f6.3,f6.3,f6.3)
write(*,*) " Convergence level = ",cctk_convlevel
- write(*,*) " Grid level = ",cctk_levfac(1),cctk_levfac(2),cctk_levfac(3)
write(*,102) x_origin,y_origin,z_origin
102 format(2X,"Minimum global coords ",f6.3,f6.3,f6.3)
write(*,103) x_origin+this_dx*(cctk_gsh(1)-1),