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
|
c -*-Fortran-*-
c $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDSpaceTimeToy/src/InitialData.F77,v 1.5 2001/08/26 13:59:47 schnetter Exp $
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
subroutine IDSpaceTimeToy_InitialData (CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
CCTK_REAL pi
CCTK_REAL omega
CCTK_REAL dt
integer i,j,k
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)
phi(i,j,k) = amplitude
$ * cos((kx*cart3d_x(i,j,k) + ky*cart3d_y(i,j,k)
$ + kz*cart3d_z(i,j,k) + omega*cctk_time) * pi)
psi(i,j,k) = - amplitude
$ * sin((kx*cart3d_x(i,j,k) + ky*cart3d_y(i,j,k)
$ + kz*cart3d_z(i,j,k) + omega*cctk_time) * pi)
$ * pi * omega
phi_p(i,j,k) = amplitude
$ * cos((kx*cart3d_x(i,j,k) + ky*cart3d_y(i,j,k)
$ + kz*cart3d_z(i,j,k) + omega*(cctk_time-dt)) * pi)
psi_p(i,j,k) = - amplitude
$ * sin((kx*cart3d_x(i,j,k) + ky*cart3d_y(i,j,k)
$ + kz*cart3d_z(i,j,k) + omega*(cctk_time-dt)) * pi)
$ * pi * 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)
phi(i,j,k) = amplitude / spher3d_r(i,j,k)
$ * exp(- (spher3d_r(i,j,k) - radius - cctk_time)**2 / sigma**2)
psi(i,j,k) = phi(i,j,k)
$ * 2 * (spher3d_r(i,j,k) - radius - cctk_time) / sigma**2
phi_p(i,j,k) = amplitude / spher3d_r(i,j,k)
$ * exp(- (spher3d_r(i,j,k) - radius - (cctk_time-dt))**2 / sigma**2)
psi_p(i,j,k) = phi(i,j,k)
$ * 2 * (spher3d_r(i,j,k) - radius - (cctk_time-dt)) / sigma**2
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)
phi(i,j,k) = amplitude
$ * sin(kx * (cart3d_x(i,j,k) - 0.5d0) * pi)
$ * sin(ky * (cart3d_y(i,j,k) - 0.5d0) * pi)
$ * sin(kz * (cart3d_z(i,j,k) - 0.5d0) * pi)
$ * cos(omega * cctk_time * pi)
psi(i,j,k) = - amplitude
$ * sin(kx * (cart3d_x(i,j,k) - 0.5d0) * pi)
$ * sin(ky * (cart3d_y(i,j,k) - 0.5d0) * pi)
$ * sin(kz * (cart3d_z(i,j,k) - 0.5d0) * pi)
$ * sin(omega * cctk_time * pi)
$ * omega * pi
phi_p(i,j,k) = amplitude
$ * sin(kx * (cart3d_x(i,j,k) - 0.5d0) * pi)
$ * sin(ky * (cart3d_y(i,j,k) - 0.5d0) * pi)
$ * sin(kz * (cart3d_z(i,j,k) - 0.5d0) * pi)
$ * cos(omega * (cctk_time-dt) * pi)
psi_p(i,j,k) = - amplitude
$ * sin(kx * (cart3d_x(i,j,k) - 0.5d0) * pi)
$ * sin(ky * (cart3d_y(i,j,k) - 0.5d0) * pi)
$ * sin(kz * (cart3d_z(i,j,k) - 0.5d0) * pi)
$ * sin(omega * (cctk_time-dt) * pi)
$ * omega * pi
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)
phi(i,j,k) = 0
psi(i,j,k) = 0
phi_p(i,j,k) = 0
psi_p(i,j,k) = 0
end do
end do
end do
end if
if (hydrotoy_active.eq.1) then
do k=1,cctk_lsh(3)
do j=1,cctk_lsh(2)
do i=1,cctk_lsh(1)
psi(i,j,k) = psi(i,j,k) - u(i,j,k)
psi_p(i,j,k) = psi_p(i,j,k) - u_p(i,j,k)
end do
end do
end do
end if
end
|