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
|
/*@@
@file WaveToy.F77
@date
@author Tom Goodale
@desc
Evolution routines for the wave equation solver
@enddesc
@@*/
#include "cctk.h"
#include "cctk_parameters.h"
#include "cctk_arguments.h"
/*@@
@routine WaveToyF77_Boundaries
@date
@author Tom Goodale
@desc
Boundary conditions for the wave equation
@enddesc
@calls ApplyFlatBC,ApplyRadiativeBC
@calledby
@history
@endhistory
@@*/
subroutine WaveToyF77_Boundaries(CCTK_FARGUMENTS)
implicit none
DECLARE_CCTK_FARGUMENTS
DECLARE_CCTK_PARAMETERS
integer sw(3)
integer CCTK_Equals
c Set the stencil width
sw(1)=1
sw(2)=1
sw(3)=1
call ApplySymmetry(cctkGH,"wavetoy::scalarfields")
if (CCTK_EQUALS(bound,"flat")) then
call ApplyFlatBC(cctkGH,sw,"wavetoy::phi")
else if (CCTK_Equals(bound,"radiation").eq.1) then
call ApplyRadiativeBC(cctkGH,1.0,sw,"wavetoy::phi","wavetoy::phi_old")
end if
return
end
c --------------------------------------------------------------
/*@@
@routine WaveToyF77_Evolution
@date
@author Tom Goodale
@desc
Evolution for the wave equation
@enddesc
@calls CCTK_SyncGroup, wavetoy_boundaries
@calledby
@history
@endhistory
@@*/
subroutine WaveToyF77_evolution(CCTK_FARGUMENTS)
implicit none
DECLARE_CCTK_FARGUMENTS
DECLARE_CCTK_PARAMETERS
c Declare local variables
INTEGER i,j,k
INTEGER istart, jstart, kstart, iend, jend, kend
CCTK_REAL dx,dy,dz,dt
c Set up shorthands
c -----------------
dx = cctk_delta_space(1)
dy = cctk_delta_space(2)
dz = cctk_delta_space(3)
dt = cctk_delta_time
istart = 2
jstart = 2
kstart = 2
iend = cctk_lsh(1)-1
jend = cctk_lsh(2)-1
kend = cctk_lsh(3)-1
c Do the evolution
c ----------------
do k = kstart, kend
do j = jstart, jend
do i = istart, iend
phi_tmp(i,j,k) =
1 2.0*(1.0 - (dt**2)*(1.0/dx**2 +
2 1.0/dy**2 +1.0/dz**2))*phi(i,j,k) -
3 phi_old(i,j,k) + (dt**2) *
5 ((phi(i+1,j,k)+phi(i-1,j,k))/dx**2
6 +(phi(i,j+1,k)+phi(i,j-1,k))/dy**2
7 +(phi(i,j,k+1)+phi(i,j,k-1))/dz**2)
end do
end do
end do
c Update timeslices
c -----------------
do k = 1, cctk_lsh(3)
do j = 1, cctk_lsh(2)
do i = 1, cctk_lsh(1)
phi_old(i,j,k) = phi(i,j,k)
phi(i,j,k) = phi_tmp(i,j,k)
end do
end do
end do
c Apply boundary conditions
c -------------------------
call WaveToyF77_Boundaries(CCTK_FARGUMENTS)
c Synchronize
c -----------
call CCTK_SyncGroup(cctkGH,"wavetoy::scalarfields")
return
end
c --------------------------------------------------------------
|