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

      implicit none

c     Declare variables in argument list
      DECLARE_CCTK_ARGUMENTS

      INTEGER   i,j,k,ierr
      INTEGER   istart, jstart, kstart, iend, jend, kend
      CCTK_REAL dx,dy,dz,dt
      CCTK_REAL dx2,dy2,dz2,dt2
      CCTK_REAL dx2i,dy2i,dz2i

      CCTK_REAL factor

c     Set up shorthands
c     -----------------
      dx = CCTK_DELTA_SPACE(1)
      dy = CCTK_DELTA_SPACE(2)
      dz = CCTK_DELTA_SPACE(3)
      dt = CCTK_DELTA_TIME

      dx2 = dx*dx
      dy2 = dy*dy
      dz2 = dz*dz
      dt2 = dt*dt

      dx2i = 1.0/dx2
      dy2i = 1.0/dy2
      dz2i = 1.0/dz2

      istart = 2
      jstart = 2
      kstart = 2
      
      iend = cctk_lsh(1)-1
      jend = cctk_lsh(2)-1
      kend = cctk_lsh(3)-1

      factor = 2*(1 - (dt2)*(dx2i + dy2i + dz2i))

c     Do the evolution
c     ----------------
      do k = kstart, kend
         do j = jstart, jend
            do i = istart, iend
               
               phi(i,j,k) = factor*phi_p(i,j,k) -
     &              phi_p_p(i,j,k) + (dt2) * 
     &              ((phi_p(i+1,j,k)+phi_p(i-1,j,k))*dx2i
     &              +(phi_p(i,j+1,k)+phi_p(i,j-1,k))*dy2i
     &              +(phi_p(i,j,k+1)+phi_p(i,j,k-1))*dz2i)

            end do
         end do
      end do
         
      return
      end 


 /*@@
   @routine    WaveToyF77_Boundaries
   @date       
   @author     Tom Goodale
   @desc 
               Boundary conditions for the wave equation
   @enddesc 
   @history 
 
   @endhistory 

@@*/

      subroutine WaveToyF77_Boundaries(CCTK_ARGUMENTS)

      implicit none

      DECLARE_CCTK_ARGUMENTS
      DECLARE_CCTK_PARAMETERS
      DECLARE_CCTK_FUNCTIONS
      CCTK_INT Boundary_SelectVarForBC

c     Local declarations
      CCTK_REAL zero,one
      CCTK_REAL finf
      integer npow
      integer ierr
      integer sw(3)

      ierr = -1
      zero = 0.0
      one = 1.0

      npow = 1
      finf = 1.0d0

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 CartSymGN(ierr,cctkGH,"wavetoy::scalarevolve")

c     Apply the outer boundary conditions
c     -----------------------------------
c     Note: In each of the following calls to Boundary_SelectVarForBC,
c     default arguments are used, so an invalid table handle of -1 can
c     be passed
      if (CCTK_EQUALS(bound,"flat")) then
         ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, -1,
     $        "wavetoy::phi", "Flat");
      else if (CCTK_EQUALS(bound,"static")) then
         ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, -1,
     $        "wavetoy::phi", "Static");
      else if (CCTK_EQUALS(bound,"radiation")) then
         print *, "selecting for rediation cond"
         ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, -1,
     $        "wavetoy::phi", "Radiation");
      else if (CCTK_EQUALS(bound,"robin")) then
         ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, -1,
     $        "wavetoy::phi", "Robin");
      else if (CCTK_EQUALS(bound,"zero")) then
         ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, -1,
     $        "wavetoy::phi", "Scalar");
      else if (CCTK_EQUALS(bound,"none")) then
         ierr = Boundary_SelectVarForBC(cctkGH, CCTK_ALL_FACES, -1,
     $        "wavetoy::phi", "None");
      else
   	 call CCTK_WARN(0,"Unrecognized boundary condition")
      end if
      if (ierr < 0)  then
        call CCTK_WARN(0,"Boundary conditions not applied - giving up!")
      end if

      return
      end