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
|
! Initialisation of the level set function and various other things.
! $Header$
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
subroutine EHFinder_Init_F(CCTK_ARGUMENTS)
use EHFinder_mod
implicit none
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_FUNCTIONS
CCTK_INT :: i, j, k, status
CCTK_REAL, dimension(3) :: xp, xpt
CCTK_REAL, dimension(3,3) :: txyz
CCTK_REAL :: cosa, sina, cosb, sinb, cosc, sinc
CCTK_REAL :: last_time
CCTK_REAL :: theta, dtheta, thetamin, thetamax
CCTK_INT, dimension(1) :: lsh, lbnd
! Initialize the current_iteration variable.
current_iteration = last_iteration_number
! Initialize the last_time varible. Note the parameters should be
! chosen to be consistent with the run producing the numerical data.
! F.ex. if dt was 0.1 but the data was only stored every 4 iterations, the
! parameters should be chosen so that dt now is 0.4.
last_time = abs(cctk_delta_time) * last_iteration_number / &
saved_iteration_every
cctk_time = last_time
! If a sphere is requested...
if ( CCTK_EQUALS( initial_f, 'sphere' ) ) then
! Set up a sphere of radius initial_rad and translated
! (translate_x,translate_y,translate_z) away from the origin.
f = sqrt( ( x - translate_x )**2 + &
( y - translate_y )**2 + &
( z - translate_z )**2 ) - initial_rad
if ( evolve_generators .gt. 0 ) then
call CCTK_GrouplbndGN ( status, cctkGH, 1, lbnd, "ehfinder::generators" )
if ( status .lt. 0 ) then
call CCTK_WARN ( 0, "cannot get lower bounds for generator arrays" )
end if
call CCTK_GrouplshGN ( status, cctkGH, 1, lsh, "ehfinder::generators" )
if ( status .lt. 0 ) then
call CCTK_WARN ( 0, "cannot get local size for generator arrays" )
end if
if ( CCTK_EQUALS( generator_distribution, 'line' ) ) then
if ( CCTK_EQUALS( domain, 'full' ) ) then
thetamin = zero; thetamax = pi
else if ( CCTK_EQUALS( domain, 'bitant') ) then
if ( CCTK_EQUALS( bitant_plane, 'xy' ) ) then
thetamin = zero; thetamax = half * pi
else
thetamin = zero; thetamax = pi
end if
else if ( CCTK_EQUALS( domain, 'quadrant' ) ) then
if ( CCTK_EQUALS( quadrant_direction, 'x' ) .or. &
CCTK_EQUALS( quadrant_direction, 'y' ) ) then
thetamin = zero; thetamax = half * pi
else
thetamin = zero; thetamax = pi
end if
else if ( CCTK_EQUALS( domain, 'octant' ) ) then
thetamin = zero; thetamax = half * pi
end if
if ( number_of_generators .eq. 1 ) then
theta = half * ( thetamax - thetamin ) + thetamin
else
dtheta = ( thetamax - thetamin ) / ( number_of_generators - 1 )
end if
do i = 1, lsh(1)
theta = thetamin + dtheta * ( i + lbnd(1) - 1 )
xg(i) = initial_rad * sin(theta) + translate_x
yg(i) = translate_y
zg(i) = initial_rad * cos(theta) + translate_z
end do
end if
end if
end if
! If an ellipsoid is requested...
if ( CCTK_EQUALS( initial_f, 'ellipsoid' ) ) then
! Calculate sines and cosines of the rotation parameters.
cosa = cos(rotation_alpha)
sina = sin(rotation_alpha)
cosb = cos(rotation_beta)
sinb = sin(rotation_beta)
cosc = cos(rotation_gamma)
sinc = sin(rotation_gamma)
! Set up the rotation matrix. The order is alpha around the z-axis,
! beta around the y-axis and finally gamma around the x-axis.
txyz(1,1) = cosa * cosb
txyz(1,2) = sina * cosb
txyz(1,3) = -sinb
txyz(2,1) = cosa * sinb * sinc - sina * cosc
txyz(2,2) = sina * sinb * sinc + cosa * cosc
txyz(2,3) = cosb * sinc
txyz(3,1) = cosa * sinb * cosc + sina * sinc
txyz(3,2) = sina * sinb * cosc - cosa * sinc
txyz(3,3) = cosb * cosc
! print*,txyz
! Apply the rotations and translation for all points on the grid.
! Even though at first glance it looks like the translation is done
! first, the opposite is actually true.
do k = 1, nz
do j = 1, ny
do i = 1, nx
xp(1) = x(i,j,k) - translate_x
xp(2) = y(i,j,k) - translate_y
xp(3) = z(i,j,k) - translate_z
xpt = matmul ( txyz, xp )
f(i,j,k) = sqrt( xpt(1)**2 / initial_a**2 + &
xpt(2)**2 / initial_b**2 + &
xpt(3)**2 / initial_c**2) - 1.0
end do
end do
end do
end if
! if an ovaloid of Cassini is requested...
if ( CCTK_EQUALS( initial_f, 'cassini' ) ) then
f = (x**2+y**2+z**2)**2 + cas_a**4 - &
2*cas_a**2*(x**2 - (y**2+z**2)) - cas_b**4
end if
! Initialise the internal mask.
eh_mask = 0
return
end subroutine EHFinder_Init_F
subroutine EHFinder_Init(CCTK_ARGUMENTS)
use EHFinder_mod
implicit none
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_FUNCTIONS
! Set up the value used in interiour inactive cells.
ex_value = - ( one + shell_width ) * maxval(cctk_delta_space)
! Initialize nx, ny and nz based on the local gris size. They are defined
! in the module EHFinder_mod and will therefore be known from now on by all
! routine that uses this module.
nx = cctk_lsh(1)
ny = cctk_lsh(2)
nz = cctk_lsh(3)
! Find the maximal grid spacing.
delta = maxval ( cctk_delta_space )
! Initialise the ghostzone information. These variables are defined in
! the module EHFinder_mod.
ngx = cctk_nghostzones(1)
ngy = cctk_nghostzones(2)
ngz = cctk_nghostzones(3)
! Get handles for various reduction operations.
call CCTK_ReductionArrayHandle ( max_handle, 'maximum' )
if ( max_handle .lt. 0 ) then
call CCTK_WARN(0,'Could not obtain a handle for maximum reduction')
end if
call CCTK_ReductionArrayHandle ( min_handle, 'minimum' )
if ( min_handle .lt. 0 ) then
call CCTK_WARN(0,'Could not obtain a handle for minimum reduction')
end if
call CCTK_ReductionArrayHandle ( sum_handle, 'sum' )
if ( sum_handle .lt. 0 ) then
call CCTK_WARN(0,'Could not obtain a handle for sum reduction')
end if
! ! Register a coodinatsystem for the 2D surface grid arrays. This has no
! ! effect yet, but should do something useful, when TAGS can be used to
! ! assign coordinates to grid arrays.
! call CCTK_CoordRegisterSystem ( ierr, 2, 'cart2d' )
! if ( ierr .lt. 0 ) then
! call CCTK_WARN(1,'Could not register a 2D coordinate system as "cart2d"')
! end if
! call CCTK_CoordRegisterData ( ierr, 1, 'ehfinder::ctheta', &
! 'x', 'cart2d' )
! if ( ierr .lt. 0 ) then
! call CCTK_WARN(1,'Could not register ctheta as the 1st coordinate for "cart2d"')
! end if
! call CCTK_CoordRegisterData ( ierr, 2, 'ehfinder::cphi', &
! 'y', 'cart2d' )
! if ( ierr .lt. 0 ) then
! call CCTK_WARN(1,'Could not register cphi as the 2nd coordinate for "cart2d"')
! end if
!
! call CCTK_INFO('2d coordinate system registered')
end subroutine EHFinder_Init
|