aboutsummaryrefslogtreecommitdiff
path: root/src/CartGrid3D.F
blob: b04f5a73a1a6f83272aa7f42d46fff7c4431eca1 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
 /*@@
   @file      CartGrid3D.F
   @date      Thu Feb 18
   @author    Gabrielle Allen
   @desc 
   Set up coordinates for a 3D Cartesian grid 
   @enddesc 
 @@*/

/*#define CARTGRID3D_DEBUG*/

#include "cctk.h"
#include "declare_arguments.h"
#include "declare_parameters.h"


      subroutine CartGrid3D(CCTK_FARGUMENTS)

      implicit none

      DECLARE_CCTK_FARGUMENTS
      DECLARE_PARAMETERS

      integer :: CCTK_Equals
      integer :: iflag, iconv, i, j, k
      REAL    :: x0,y0,z0,dx,dy,dz

      iflag = 0

      iconv = 2**convlevel


c     --------------------------------------------------------------
c
c     BYRANGE
c
c     User gives: minimum and maximum values of coordinates and
c                 the number of gridpoints on the coarse grid
c
c     --------------------------------------------------------------

      if (CCTK_Equals(grid,'byrange') == 1) then

         iflag = iflag +1

         if (CCTK_Equals(symmetry,'octant') == 1) then

            iflag = iflag +1

c           Grid spacing on coarsest grid          
            coarse_dx = grid_xmax/(DBLE(global_sh(1))-1.5d0)
            coarse_dy = grid_ymax/(DBLE(global_sh(2))-1.5d0)
            coarse_dz = grid_zmax/(DBLE(global_sh(3))-1.5d0)

c           Grid spacing on this grid
            dx = -coarse_dx*0.5d0
            dy = -coarse_dy*0.5d0
            dz = -coarse_dz*0.5d0

c           Minimum coordinate values on this grid
            x0 = (0.5-dble(nghostzones))*dx
            y0 = (0.5-dble(nghostzones))*dy
            z0 = (0.5-dble(nghostzones))*dz

         else if (CCTK_Equals(symmetry,'quadrant') == 1) then

            iflag = iflag +1
            print *,"FIXME"

         else if (CCTK_Equals(symmetry,'bitant') == 1) then

            iflag = iflag +1
            print *,"FIXME"

         else if (CCTK_Equals(symmetry,'full') == 1) then

            iflag = iflag +1

c           Set minimum values of coordinates
            x0 = grid_xmin
            y0 = grid_ymin
            z0 = grid_zmin

c           dx,dy,dz on the coarsest grid of each GH
            coarse_dx = (grid_xmax-grid_xmin)/max(global_sh(1)-1,1)
            coarse_dy = (grid_ymax-grid_ymin)/max(global_sh(2)-1,1)
            coarse_dz = (grid_zmax-grid_zmin)/max(global_sh(3)-1,1)

         end if

c        dx,dy,dz on the grid we are on
         dx = coarse_dx/levfac
         dy = coarse_dy/levfac
         dz = coarse_dz/levfac

         write(*,'(1X,A,1X,3(A,G12.7,2X))') 
     &       'Grid by range: Setting ','dx=>',dx,'dy=>',dy,'dz=>',dz


c     -----------------------------------------------------------
c
c     BYSPACING
c
c     User gives: grid spacing on the coarsest GH and
c                 the number of gridpoints on the coarsest GH
c
c     -----------------------------------------------------------

      else if (CCTK_Equals(grid,'byspacing') == 1) then

c        dx, dy, dz on the grid we are on
         dx = grid_dx*iconv/levfac
         dy = grid_dy*iconv/levfac
         dz = grid_dz*iconv/levfac

         if (CCTK_Equals(symmetry,'bitant') == 1) then 

            iflag = iflag + 1
            x0 = (0.5 - global_sh(1)/2)*dx
            y0 = (0.5 - global_sh(2)/2)*dy
            z0 = (-dble(nghostzones)+0.5d0)*dz

         else if (CCTK_Equals(symmetry,'quadrant') == 1) then 

            iflag = iflag + 1

            x0 = (-dble(nghostzones)+0.5d0)*dx
            y0 = (-dble(nghostzones)+0.5d0)*dy
            z0 = (0.5D0 - global_sh(3)/2)*dz
	    
         else if (CCTK_Equals(symmetry,'octant') == 1) then 
         
            iflag = iflag + 1

            x0 = (-dble(nghostzones)+0.5d0)*dx
            y0 = (-dble(nghostzones)+0.5d0)*dy
            z0 = (-dble(nghostzones)+0.5d0)*dz

         else if (CCTK_Equals(symmetry,'full')==1) then 

            iflag = iflag + 1

            x0 = (0.5 - global_sh(1)/2)*dx
            y0 = (0.5 - global_sh(2)/2)*dy
            z0 = (0.5 - global_sh(3)/2)*dz

         end if


c     --------------------------------------------------------------
c
c     BOX (-0.5 to 0.5)
c
c     User gives: number of gridpoints on the coarse grid
c
c     --------------------------------------------------------------    

      elseif (CCTK_Equals(grid,'box')==1) then 

         iflag = iflag + 1
 
c        Coordinates are all -0.5 to 0.5
         x0 = -0.5
         y0 = -0.5
         z0 = -0.5

c        dx,dy,dz on the coarsest grid of each GH
         coarse_dx = 1.d0/max(global_sh(1)-1,1)
         coarse_dy = 1.d0/max(global_sh(1)-1,1)
         coarse_dz = 1.d0/max(global_sh(1)-1,1)

c        dx,dy,dz on the grid we are on
         dx = coarse_dx/levfac
         dy = coarse_dy/levfac
         dz = coarse_dz/levfac

c        Special cases
         if (global_sh(1) == 1) x0 = 0.0D0
         if (global_sh(2) == 1) y0 = 0.0D0
         if (global_sh(3) == 1) z0 = 0.0D0

         write(*,'(1X,A,1X,3(A,G12.7,2X))') 
     &       'Box grid: Setting ','dx=>',dx,'dy=>',dy,'dz=>',dz


      end if


c     No grid was set up
c     ------------------
      if (iflag.ne.1) then
        call CCTK_Warn(0,"No grid set up in CartGrid3D")
        call CCTK_Exit(GH)
      end if


c     Set spatial coordinates
c     -----------------------
      do i=1,sh(1)
         do j=1,sh(2)
            do k=1,sh(3)
               x(i,j,k) = dx*(i+lb(1)-1) + x0
               y(i,j,k) = dy*(j+lb(2)-1) + y0
               z(i,j,k) = dz*(k+lb(3)-1) + z0

#ifdef CARTGRID3D_DEBUG
               write(*,*) "Index: ",i,j,k
               write(*,*) "Coord: ",x(i,j,k),y(i,j,k),z(i,j,k)
#endif

            end do
         end do
      end do

      r = sqrt(x**2 + y**2 + z**2)

      delta_space(1) = dx
      delta_space(2) = dy
      delta_space(3) = dz

      return
      end subroutine CartGrid3D