aboutsummaryrefslogtreecommitdiff
path: root/src/WaveToy.F77
blob: c6036a0a3818e5791d785c4fdc0b8b4c737ab813 (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
 /*@@
   @file      WaveToy.F77
   @date      
   @author    Tom Goodale
   @desc 
              Evolution routines for the wave equation solver
   @enddesc 
 @@*/
  
c Using Cactus infrastructure
#include "cctk.h" 

c Using Cactus parameters
#include "cctk_parameters.h"

c Using Cactus arguments lists
#include "cctk_arguments.h"


 /*@@
   @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

c     Declare variables in argument list
      DECLARE_CCTK_FARGUMENTS

      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_next(i,j,k) = 
     &		 2.0*(1.0 - (dt**2)*(1.0/dx**2 +   
     &              1.0/dy**2 +1.0/dz**2))*phi(i,j,k) -
     &              phi_old(i,j,k) + (dt**2) * 
     &              ((phi(i+1,j,k)+phi(i-1,j,k))/dx**2
     &              +(phi(i,j+1,k)+phi(i,j-1,k))/dy**2
     &              +(phi(i,j,k+1)+phi(i,j,k-1))/dz**2)

            end do
         end do
      end do
         
c     Synchronize
c     -----------
      call CCTK_SyncGroup(cctkGH,"wavetoyf77::scalartmps")

c     Apply boundary conditions
c     -------------------------
      call WaveToyF77_Boundaries(CCTK_PASS_FTOF)	

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_next(i,j,k)

            end do
         end do
      end do

      return
      end 


 /*@@
   @routine    WaveToyF77_Boundaries
   @date       
   @author     Tom Goodale
   @desc 
               Boundary conditions for the wave equation
   @enddesc 
   @calls      ApplySymmetry,ApplyFlatBC,ApplyRadiativeBC 
   @history 
 
   @endhistory 

@@*/

      subroutine WaveToyF77_Boundaries(CCTK_FARGUMENTS)

      implicit none

c     Declare arguement list
      DECLARE_CCTK_FARGUMENTS

c     Declare parameters
      DECLARE_CCTK_PARAMETERS

c     Local declarations
      CCTK_REAL zero
      integer ierr
      integer sw(3)

c     Cactus declarations
      integer CCTK_Equals

      zero = 0.0

c     Set the stencil width
c     ---------------------
      sw(1)=1
      sw(2)=1
      sw(3)=1

c     Apply the symmetry boundary conditions on any coordinate axes
c     -------------------------------------------------------------
      call ApplySymmetry(cctkGH,"wavetoyf77::scalartmps")

c     Apply the outer boundary conditions
c     -----------------------------------
      if (CCTK_EQUALS(bound,"flat")) then
         call ApplyFlatBC(ierr,cctkGH,sw,"wavetoyf77::phi_next")
      else if (CCTK_EQUALS(bound,"zero")) then
         call ApplyConstantBC(ierr,cctkGH,zero,sw,"wavetoyf77::phi_next")
      else if (CCTK_Equals(bound,"radiation").eq.1) then
         call ApplyRadiativeBC(ierr,cctkGH,zero,zero,sw,"wavetoyf77::phi_next",
     &	      "wavetoy::phi")
      end if

      if (ierr < 0)  then
        call CCTK_WARN(0,"Boundary conditions not applied - giving up!");
      end if

      return
      end