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
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
|
/*@@
@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 "cctk_arguments.h"
#include "cctk_parameters.h"
subroutine CartGrid3D(CCTK_FARGUMENTS)
implicit none
DECLARE_CCTK_FARGUMENTS
DECLARE_CCTK_PARAMETERS
integer :: CCTK_Equals
integer :: iconv, i, j, k
CCTK_REAL :: x_origin,y_origin,z_origin
CCTK_REAL :: this_dx,this_dy,this_dz
CCTK_REAL :: xmax1,ymax1,zmax1
iconv = 2**(cctk_convlevel-1)
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(type,'byrange') == 1) then
if (xyzmax > 0) then
xmax1 = xyzmax
ymax1 = xyzmax
zmax1 = xyzmax
else
xmax1 = xmax
ymax1 = ymax
zmax1 = zmax
end if
if (CCTK_Equals(symmetry,'octant') == 1) then
c Grid spacing on coarsest grid
if (no_origin==1) then
coarse_dx = xmax/(DBLE(cctk_gsh(1))-1.5d0)
coarse_dy = ymax/(DBLE(cctk_gsh(2))-1.5d0)
coarse_dz = zmax/(DBLE(cctk_gsh(3))-1.5d0)
else
coarse_dx = xmax/(DBLE(cctk_gsh(1))-1.0d0)
coarse_dy = ymax/(DBLE(cctk_gsh(2))-1.0d0)
coarse_dz = zmax/(DBLE(cctk_gsh(3))-1.0d0)
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)
c Minimum coordinate values on this grid
if (no_origin==1) then
x_origin = (0.5-dble(cctk_nghostzones(1)))*this_dx
y_origin = (0.5-dble(cctk_nghostzones(2)))*this_dy
z_origin = (0.5-dble(cctk_nghostzones(3)))*this_dz
else
x_origin = (-dble(cctk_nghostzones(1)))*this_dx
y_origin = (-dble(cctk_nghostzones(2)))*this_dy
z_origin = (-dble(cctk_nghostzones(3)))*this_dz
end if
else if (CCTK_Equals(symmetry,'quadrant') == 1) then
print *,"FIXME"
else if (CCTK_Equals(symmetry,'bitant') == 1) then
print *,"FIXME"
else if (CCTK_Equals(symmetry,'full') == 1) then
c Set minimum values of coordinates
x_origin = xmin
y_origin = ymin
z_origin = zmin
c dx,dy,dz on the coarsest grid of each GH
coarse_dx = (xmax-xmin)/max(cctk_gsh(1)-1,1)
coarse_dy = (ymax-ymin)/max(cctk_gsh(2)-1,1)
coarse_dz = (zmax-zmin)/max(cctk_gsh(3)-1,1)
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(*,'(1X,A,1X,3(A,G12.7,2X))')
& 'Grid by range: Setting ','dx=>',this_dx,'dy=>',
& this_dy,'dz=>',this_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(type,'byspacing') == 1) then
c Dx, Dy, Dx on the coarsest grid
if (dxyz > 0) then
coarse_dx = dxyz
coarse_dy = dxyz
coarse_dz = dxyz
else
coarse_dx = dx
coarse_dy = dy
coarse_dz = dz
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)
if (CCTK_Equals(symmetry,'bitant') == 1) then
if (no_origin==1) then
x_origin = (0.5 - cctk_gsh(1)/2)*this_dx
y_origin = (0.5 - cctk_gsh(2)/2)*this_dy
z_origin = (-dble(cctk_nghostzones(3))+0.5d0)*this_dz
else
x_origin = (- cctk_gsh(1)/2)*this_dx
y_origin = (- cctk_gsh(2)/2)*this_dy
z_origin = (-dble(cctk_nghostzones(3)))*this_dz
end if
else if (CCTK_Equals(symmetry,'quadrant') == 1) then
if (no_origin==1) then
x_origin = (-dble(cctk_nghostzones(1))+0.5d0)*this_dx
y_origin = (-dble(cctk_nghostzones(2))+0.5d0)*this_dy
z_origin = (0.5D0 - cctk_gsh(3)/2)*this_dz
else
x_origin = (-dble(cctk_nghostzones(1)))*this_dx
y_origin = (-dble(cctk_nghostzones(2)))*this_dy
z_origin = (- cctk_gsh(3)/2)*this_dz
end if
else if (CCTK_Equals(symmetry,'octant') == 1) then
if (no_origin==1) then
x_origin = (-dble(cctk_nghostzones(1))+0.5d0)*this_dx
y_origin = (-dble(cctk_nghostzones(2))+0.5d0)*this_dy
z_origin = (-dble(cctk_nghostzones(3))+0.5d0)*this_dz
else
x_origin = (-dble(cctk_nghostzones(1)))*this_dx
y_origin = (-dble(cctk_nghostzones(2)))*this_dy
z_origin = (-dble(cctk_nghostzones(3)))*this_dz
end if
else if (CCTK_Equals(symmetry,'full')==1) then
if (no_origin==1) then
x_origin = (0.5 - cctk_gsh(1)/2)*this_dx
y_origin = (0.5 - cctk_gsh(2)/2)*this_dy
z_origin = (0.5 - cctk_gsh(3)/2)*this_dz
else
x_origin = (- cctk_gsh(1)/2)*this_dx
y_origin = (- cctk_gsh(2)/2)*this_dy
z_origin = (- cctk_gsh(3)/2)*this_dz
end if
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(type,'box') == 1) then
c Coordinates are all -0.5 to 0.5
x_origin = -0.5
y_origin = -0.5
z_origin = -0.5
c dx,dy,dz on the coarsest grid of each GH
coarse_dx = 1.d0/max(cctk_gsh(1)-1,1)
coarse_dy = 1.d0/max(cctk_gsh(1)-1,1)
coarse_dz = 1.d0/max(cctk_gsh(1)-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)
c Special cases
if (cctk_gsh(1) == 1) x_origin = 0.0D0
if (cctk_gsh(2) == 1) y_origin = 0.0D0
if (cctk_gsh(3) == 1) z_origin = 0.0D0
write(*,'(1X,A,1X,3(A,G12.7,2X))')
& 'Box grid: Setting ','dx=>',this_dx,'dy=>',this_dy,
& 'dz=>',this_dz
end if
c Set spatial coordinates
c -----------------------
do i=1,cctk_lsh(1)
do j=1,cctk_lsh(2)
do k=1,cctk_lsh(3)
x(i,j,k) = this_dx*(i+cctk_lbnd(1)-1) + x_origin
y(i,j,k) = this_dy*(j+cctk_lbnd(2)-1) + y_origin
z(i,j,k) = this_dz*(k+cctk_lbnd(3)-1) + z_origin
end do
end do
end do
r = sqrt(x**2 + y**2 + z**2)
cctk_delta_space(1) = this_dx
cctk_delta_space(2) = this_dy
cctk_delta_space(3) = this_dz
cctk_origin_space(1) = x_origin
cctk_origin_space(2) = y_origin
cctk_origin_space(3) = z_origin
#ifdef CARTGRID3D_DEBUG
write(*,*)
write(*,*) " CartGrid3D"
write(*,*) " ----------"
write(*,100) coarse_dx,coarse_dy,coarse_dz
100 format(2X,"Dx,Dy,Dz on coarse grid ",f6.3,f6.3,f6.3)
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),
& y_origin+this_dy*(cctk_gsh(2)-1),z_origin+this_dz*(cctk_gsh(3)-1)
103 format(2X,"Maximum global coords ",f6.3,f6.3,f6.3)
write(*,104) x(1,1,1),y(1,1,1),z(1,1,1)
104 format(2X,"Minimum local coords ",f6.3,f6.3,f6.3)
write(*,105) x(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),
& y(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),z(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3))
105 format(2X,"Maximum local coords ",f6.3,f6.3,f6.3)
#endif
end subroutine CartGrid3D
|