aboutsummaryrefslogtreecommitdiff
path: root/CarpetExtra/IDHydroToy/src/InitialData.F77
blob: 1a10654651da2f7ac1e84484bfa651d22911ec67 (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
226
227
228
229
230
231
c     -*-Fortran-*-
c     $Header:$
      
#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"
      
      subroutine IDHydroToy_InitialData (CCTK_ARGUMENTS)
      
      implicit none
      
      DECLARE_CCTK_ARGUMENTS
      DECLARE_CCTK_PARAMETERS
      DECLARE_CCTK_FUNCTIONS
      
      CCTK_REAL pi
      CCTK_REAL omega
      CCTK_REAL dt
      CCTK_REAL x,y,z, r
      integer i,j,k
      
      CCTK_REAL vr
      
      external erf
      real*8 erf
      
      pi = 4*atan(1.d0)
      
      omega = sqrt(kx**2+ky**2+kz**2)
      
      dt = CCTK_DELTA_TIME
      
      if (CCTK_EQUALS(initial_data,"plane")) then 
         
         do k=1,cctk_lsh(3)
            do j=1,cctk_lsh(2)
               do i=1,cctk_lsh(1)
                  
                  x = cart3d_x(i,j,k)
                  y = cart3d_y(i,j,k)
                  z = cart3d_z(i,j,k)
                  
                  u(i,j,k) = amplitude
     $                 * cos((kx*x + ky*y + kz*z + omega*cctk_time) * pi)
                  vx(i,j,k) = u(i,j,k) * kx / omega
                  vy(i,j,k) = u(i,j,k) * ky / omega
                  vz(i,j,k) = u(i,j,k) * kz / omega
                  
                  u_p(i,j,k) = amplitude
     $                 * cos((kx*x + ky*y + kz*z + omega*(cctk_time-dt)) * pi)
                  vx_p(i,j,k) = u_p(i,j,k) * kx / omega
                  vy_p(i,j,k) = u_p(i,j,k) * ky / omega
                  vz_p(i,j,k) = u_p(i,j,k) * kz / omega
                  
                  u_p_p(i,j,k) = amplitude
     $                 * cos((kx*x + ky*y + kz*z + omega*(cctk_time-2*dt)) * pi)
                  vx_p_p(i,j,k) = u_p_p(i,j,k) * kx / omega
                  vy_p_p(i,j,k) = u_p_p(i,j,k) * ky / omega
                  vz_p_p(i,j,k) = u_p_p(i,j,k) * kz / omega
                  
               end do
            end do
         end do
         
      else if (CCTK_EQUALS(initial_data,"gaussian")) then 
         
         do k=1,cctk_lsh(3)
            do j=1,cctk_lsh(2)
               do i=1,cctk_lsh(1)
                  
                  x = cart3d_x(i,j,k)
                  y = cart3d_y(i,j,k)
                  z = cart3d_z(i,j,k)
                  r = spher3d_r(i,j,k)
                  
                  u(i,j,k) = amplitude
     $                 * exp(- (r - radius)**2 / sigma**2)
                  
                  vr = - 2*amplitude * (r - radius) / sigma**2
     $                 * exp(- (r - radius)**2 / sigma**2)
                  vx(i,j,k) = vr * x/r
                  vy(i,j,k) = vr * y/r
                  vz(i,j,k) = vr * z/r
                  
                  u_p(i,j,k) = amplitude/2 * (r - dt) / r
     $                 * exp(- (r - radius - dt)**2 / sigma**2)
     $                 + amplitude/2 * (r + dt) / r
     $                 * exp(- (r - radius + dt)**2 / sigma**2)
                  
                  vr = - amplitude/2 * (-dt / r**2 + (r - dt) * (r - radius - dt) / (r * sigma**2))
     $                 * exp(- (r - radius - dt)**2 / sigma**2)
     $                 - amplitude/2 * ( dt / r**2 + (r + dt) * (r - radius + dt) / (r * sigma**2))
     $                 * exp(- (r - radius - dt)**2 / sigma**2)
                  vx_p(i,j,k) = vr * x/r
                  vy_p(i,j,k) = vr * y/r
                  vz_p(i,j,k) = vr * z/r
                  
                  u_p_p(i,j,k) = amplitude/2 * (r - 2*dt) / r
     $                 * exp(- (r - radius - 2*dt)**2 / sigma**2)
     $                 + amplitude/2 * (r + 2*dt) / r
     $                 * exp(- (r - radius + 2*dt)**2 / sigma**2)
                  
                  vr = - amplitude/2 * (-2*dt / r**2 + (r - 2*dt) * (r - radius - 2*dt) / (r * sigma**2))
     $                 * exp(- (r - radius - 2*dt)**2 / sigma**2)
     $                 - amplitude/2 * ( 2*dt / r**2 + (r + 2*dt) * (r - radius + 2*dt) / (r * sigma**2))
     $                 * exp(- (r - radius - 2*dt)**2 / sigma**2)
                  vx_p_p(i,j,k) = vr * x/r
                  vy_p_p(i,j,k) = vr * y/r
                  vz_p_p(i,j,k) = vr * z/r
                  
               end do
            end do
         end do
         
      else if (CCTK_EQUALS(initial_data, "box")) then
         
c     Use kx,ky,kz as number of modes in each direction.
         
         do k=1,cctk_lsh(3)
            do j=1,cctk_lsh(2)
               do i=1,cctk_lsh(1)
                  
                  x = cart3d_x(i,j,k)
                  y = cart3d_y(i,j,k)
                  z = cart3d_z(i,j,k)
                  
                  u(i,j,k) = amplitude
     $                 * sin(kx * (x - 0.5d0) * pi)
     $                 * sin(ky * (y - 0.5d0) * pi)
     $                 * sin(kz * (z - 0.5d0) * pi)
     $                 * cos(omega * cctk_time * pi)
                  vx(i,j,k) = amplitude
     $                 * cos(kx * (x - 0.5d0) * pi)
     $                 * sin(ky * (y - 0.5d0) * pi)
     $                 * sin(kz * (z - 0.5d0) * pi)
     $                 * sin(omega * cctk_time * pi)
     $                 * kx / omega
                  vy(i,j,k) = amplitude
     $                 * sin(kx * (x - 0.5d0) * pi)
     $                 * cos(ky * (y - 0.5d0) * pi)
     $                 * sin(kz * (z - 0.5d0) * pi)
     $                 * sin(omega * cctk_time * pi)
     $                 * ky / omega
                  vz(i,j,k) = amplitude
     $                 * sin(kx * (x - 0.5d0) * pi)
     $                 * sin(ky * (y - 0.5d0) * pi)
     $                 * cos(kz * (z - 0.5d0) * pi)
     $                 * sin(omega * cctk_time * pi)
     $                 * kz / omega
                  
                  u_p(i,j,k) = amplitude
     $                 * sin(kx * (x - 0.5d0) * pi)
     $                 * sin(ky * (y - 0.5d0) * pi)
     $                 * sin(kz * (z - 0.5d0) * pi)
     $                 * cos(omega * (cctk_time - dt) * pi)
                  vx_p(i,j,k) = amplitude
     $                 * cos(kx * (x - 0.5d0) * pi)
     $                 * sin(ky * (y - 0.5d0) * pi)
     $                 * sin(kz * (z - 0.5d0) * pi)
     $                 * sin(omega * (cctk_time - dt) * pi)
     $                 * kx / omega
                  vy_p(i,j,k) = amplitude
     $                 * sin(kx * (x - 0.5d0) * pi)
     $                 * cos(ky * (y - 0.5d0) * pi)
     $                 * sin(kz * (z - 0.5d0) * pi)
     $                 * sin(omega * (cctk_time - dt) * pi)
     $                 * ky / omega
                  vz_p(i,j,k) = amplitude
     $                 * sin(kx * (x - 0.5d0) * pi)
     $                 * sin(ky * (y - 0.5d0) * pi)
     $                 * cos(kz * (z - 0.5d0) * pi)
     $                 * sin(omega * (cctk_time - dt) * pi)
     $                 * kz / omega
                  
                  u_p_p(i,j,k) = amplitude
     $                 * sin(kx * (x - 0.5d0) * pi)
     $                 * sin(ky * (y - 0.5d0) * pi)
     $                 * sin(kz * (z - 0.5d0) * pi)
     $                 * cos(omega * (cctk_time - 2*dt) * pi)
                  vx_p_p(i,j,k) = amplitude
     $                 * cos(kx * (x - 0.5d0) * pi)
     $                 * sin(ky * (y - 0.5d0) * pi)
     $                 * sin(kz * (z - 0.5d0) * pi)
     $                 * sin(omega * (cctk_time - 2*dt) * pi)
     $                 * kx / omega
                  vy_p_p(i,j,k) = amplitude
     $                 * sin(kx * (x - 0.5d0) * pi)
     $                 * cos(ky * (y - 0.5d0) * pi)
     $                 * sin(kz * (z - 0.5d0) * pi)
     $                 * sin(omega * (cctk_time - 2*dt) * pi)
     $                 * ky / omega
                  vz_p_p(i,j,k) = amplitude
     $                 * sin(kx * (x - 0.5d0) * pi)
     $                 * sin(ky * (y - 0.5d0) * pi)
     $                 * cos(kz * (z - 0.5d0) * pi)
     $                 * sin(omega * (cctk_time - 2*dt) * pi)
     $                 * kz / omega
                  
               end do
            end do
         end do
         
      else
         
         do k=1,cctk_lsh(3)
            do j=1,cctk_lsh(2)
               do i=1,cctk_lsh(1)
                  
                  u(i,j,k)  = 0
                  vx(i,j,k) = 0
                  vy(i,j,k) = 0
                  vz(i,j,k) = 0
                  
                  u_p(i,j,k)  = 0
                  vx_p(i,j,k) = 0
                  vy_p(i,j,k) = 0
                  vz_p(i,j,k) = 0
                  
                  u_p_p(i,j,k)  = 0
                  vx_p_p(i,j,k) = 0
                  vy_p_p(i,j,k) = 0
                  vz_p_p(i,j,k) = 0
                  
               end do
            end do
         end do 
         
      end if
      
      end