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
|
c -*-Fortran-*-
c $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/HydroToy/src/HydroToy.F77,v 1.3 2001/03/21 22:57:40 eschnett Exp $
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
subroutine HydroToy_EulerPredictor (CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_FUNCTIONS
DECLARE_CCTK_PARAMETERS
integer i,j,k
c Copy
do k=1,cctk_lsh(3)
do j=1,cctk_lsh(2)
do i=1,cctk_lsh(1)
u_i(i,j,k) = u(i,j,k)
vx_i(i,j,k) = vx(i,j,k)
vy_i(i,j,k) = vy(i,j,k)
vz_i(i,j,k) = vz(i,j,k)
end do
end do
end do
c Evolve
call HydroToy_Step (CCTK_PASS_FTOF)
c Apply boundaries
call HydroToy_Boundaries (CCTK_PASS_FTOF)
end
subroutine HydroToy_EulerCorrector (CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_FUNCTIONS
DECLARE_CCTK_PARAMETERS
CCTK_REAL two, half
parameter (two=2, half=1/two)
integer i,j,k
c Copy
do k=1,cctk_lsh(3)
do j=1,cctk_lsh(2)
do i=1,cctk_lsh(1)
u_i(i,j,k) = u_n(i,j,k)
vx_i(i,j,k) = vx_n(i,j,k)
vy_i(i,j,k) = vy_n(i,j,k)
vz_i(i,j,k) = vz_n(i,j,k)
end do
end do
end do
c Evolve
call HydroToy_Step (CCTK_PASS_FTOF)
c Average
do k=1,cctk_lsh(3)
do j=1,cctk_lsh(2)
do i=1,cctk_lsh(1)
u_n(i,j,k) = half * (u(i,j,k) + u_n(i,j,k))
vx_n(i,j,k) = half * (vx(i,j,k) + vx_n(i,j,k))
vy_n(i,j,k) = half * (vy(i,j,k) + vy_n(i,j,k))
vz_n(i,j,k) = half * (vz(i,j,k) + vz_n(i,j,k))
end do
end do
end do
c Apply boundaries
call HydroToy_Boundaries (CCTK_PASS_FTOF)
end
subroutine HydroToy_Step (CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_FUNCTIONS
DECLARE_CCTK_PARAMETERS
CCTK_REAL dx,dy,dz,dt
integer i,j,k
dx = CCTK_DELTA_SPACE(1)
dy = CCTK_DELTA_SPACE(2)
dz = CCTK_DELTA_SPACE(3)
dt = CCTK_DELTA_TIME
c Evolve
do k=1,cctk_lsh(3)
do j=1,cctk_lsh(2)
do i=1,cctk_lsh(1)
u_n(i,j,k) = u_i(i,j,k)
$ + dt * (vx_i(i+1,j,k) - vx_i(i-1,j,k)) / (2*dx)
$ + dt * (vy_i(i,j+1,k) - vy_i(i,j-1,k)) / (2*dy)
$ + dt * (vz_i(i,j,k+1) - vz_i(i,j,k-1)) / (2*dz)
vx_n(i,j,k) = vx_i(i,j,k)
$ + dt * (u_i(i+1,j,k) - u_i(i-1,j,k)) / (2*dx)
vy_n(i,j,k) = vy_i(i,j,k)
$ + dt * (u_i(i,j+1,k) - u_i(i,j-1,k)) / (2*dy)
vz_n(i,j,k) = vz_i(i,j,k)
$ + dt * (u_i(i,j,k+1) - u_i(i,j,k-1)) / (2*dz)
end do
end do
end do
end
subroutine HydroToy_Boundaries (CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_FUNCTIONS
DECLARE_CCTK_PARAMETERS
CCTK_REAL zero, one
parameter (zero=0, one=1)
CCTK_REAL finf, npow
parameter (finf=1, npow=1)
integer sw(3)
integer ierr
sw(1) = 1
sw(2) = 1
sw(3) = 1
c Apply boundary condition
if (CCTK_EQUALS(bound, "flat")) then
call BndFlatGN (ierr, cctkGH, sw, "hydrotoy::hydroevolve")
else if (CCTK_EQUALS(bound, "zero")) then
call BndScalarGN (ierr, cctkGH, sw, zero,
$ "hydrotoy::hydroevolve")
else if (CCTK_EQUALS(bound, "radiation")) then
call BndRadiativeGN (ierr, cctkGH, sw, zero, one,
$ "hydrotoy::hydroevolve", "hydrotoy::hydroevolve")
else if (CCTK_EQUALS(bound, "robin")) then
call BndRobinGN (ierr, cctkGH, sw, finf, npow,
$ "hydrotoy::hydroevolve")
else if (CCTK_EQUALS(bound, "none")) then
ierr = 0
else
call CCTK_WARN (0, "internal error")
end if
if (ierr .lt. 0) then
call CCTK_WARN (0, "Error while applying boundary condition")
end if
call Cart3dSymGN (ierr, cctkGH, "hydrotoy::hydroevolve")
if (ierr .lt. 0) then
call CCTK_WARN (0, "Error while applying boundary condition")
end if
end
|