aboutsummaryrefslogtreecommitdiff
path: root/CarpetExtra/FOWaveToyF77/src/FOWaveToy.F77
blob: bc2a175499cd594aec82f92ce4902e0181975be5 (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
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
c     -*-Fortran-*-

 /*@@
   @file      FOWaveToy.F77
   @date      
   @author    Scott Hawley, based on WaveToy by Tom Goodale, Erik Schnetter
   @desc 
              Evolution routines for the wave equation solver
              written as a first order system.
              The field phi itself does not enter into the evolution,
              system but we carry it along using pi = partial phi / partial t.
   @enddesc 
 @@*/
  
#include "cctk.h" 
#include "cctk_Arguments.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"



 /*@@
   @routine    FOWaveToyF77_Evolution
   @date       
   @author     Scott Hawley, using code of Tom Goodale, Erik Schnetter
   @desc 
               Evolution for the wave equation as a 1st-order system
   @enddesc 
   @calls      CCTK_SyncGroup
   @calledby   
   @history 
 
   @endhistory 

@@*/

      subroutine FOWaveToyF77_Evolution (CCTK_ARGUMENTS)
      
      implicit none
      
c     Declare variables in argument list
      DECLARE_CCTK_ARGUMENTS
      DECLARE_CCTK_FUNCTIONS
      DECLARE_CCTK_PARAMETERS

      INTEGER   i,j,k
      INTEGER   istart(3), iend(3)
      INTEGER   d
      CCTK_REAL dx,dy,dz,dt
      CCTK_REAL dxi,dyi,dzi
      
c     call CCTK_INFO ("WaveToyF77_Evolution")
      
c     Set up shorthands
c     -----------------
      dx = 2*CCTK_DELTA_SPACE(1)
      dy = 2*CCTK_DELTA_SPACE(2)
      dz = 2*CCTK_DELTA_SPACE(3)
      dt = 2*CCTK_DELTA_TIME
      
      dxi = 1/dx
      dyi = 1/dy
      dzi = 1/dz
      
      do d = 1, 3
         if (cctk_bbox(2*d-1).eq.0) then
            istart(d) = 1+cctk_nghostzones(d)
         else
            istart(d) = 1+bound_width
         end if
         if (cctk_bbox(2*d).eq.0) then
            iend(d) = cctk_lsh(d)-cctk_nghostzones(d)
         else
            iend(d) = cctk_lsh(d)-bound_width
         end if
      end do
      
c     Do the evolution
c     ----------------
      do k = istart(3), iend(3)
         do j = istart(2), iend(2)
            do i = istart(1), iend(1)
               
               pi(i,j,k) = pi_p_p(i,j,k) + dt * ( 
     $                      (phix_p(i+1,j,k) - phix_p(i-1,j,k))*dxi
     $                    + (phiy_p(i,j+1,k) - phiy_p(i,j-1,k))*dyi
     $                    + (phiz_p(i,j,k+1) - phiz_p(i,j,k-1))*dzi
     $                                          )

               phix(i,j,k) = phix_p_p(i,j,k) + dt * dxi * ( 
     $                       pi_p(i+1,j,k) - pi_p(i-1,j,k)  )

               phiy(i,j,k) = phiy_p_p(i,j,k) + dt * dyi * ( 
     $                       pi_p(i,j+1,k) - pi_p(i,j-1,k)  )

               phiz(i,j,k) = phiz_p_p(i,j,k) + dt * dzi * ( 
     $                       pi_p(i,j,k+1) - pi_p(i,j,k-1)  )

               phi(i,j,k) = phi_p_p(i,j,k) + 2*dt * pi_p(i,j,i) 
               
            end do
         end do
      end do
      
      end



 /*@@
   @routine    FOWaveToyF77_Boundaries
   @date       
   @author     Scott Hawley stealing from Tom Goodale, Erik Schnetter
   @desc 
               Boundary conditions for the wave equation
   @enddesc 
   @history 
 
   @endhistory 

@@*/

      subroutine FOWaveToyF77_Boundaries (CCTK_ARGUMENTS)
      
      implicit none
      
      DECLARE_CCTK_ARGUMENTS
      DECLARE_CCTK_PARAMETERS
      DECLARE_CCTK_FUNCTIONS
      
c     Local declarations
      integer table
      data    table /-1/
      
      character fbound*100
      CCTK_INT  fboundlen
      
      integer i,j,k
      CCTK_REAL  spher3d_r
      
      integer ierr
      CCTK_REAL ri3
      
c     call CCTK_INFO ("FOWaveToyF77_Boundaries")
      
c     Apply the excision boundary condition
c     -------------------------------------
      if (CCTK_EQUALS(excision_bound, "none")) then
c     do nothing
      else if (CCTK_EQUALS(excision_bound, "1/r")) then
         do k=1,cctk_lsh(3)
            do j=1,cctk_lsh(2)
               do i=1,cctk_lsh(1)
                  spher3d_r = sqrt(x(i,j,k)**2 + y(i,j,k)**2 + z(i,j,k)**2)
                  if (spher3d_r .le. excision_radius) then
                     pi(i,j,k)  = 0.0
                     phi(i,j,k) = 1.0 / spher3d_r
                     ri3 = phi(i,j,k)**3
                     phix(i,j,k) = - x(i,j,k) * ri3
                     phiy(i,j,k) = - y(i,j,k) * ri3
                     phiz(i,j,k) = - z(i,j,k) * ri3
                  end if
               end do
            end do
         end do 
      else
         call CCTK_WARN (0, "internal error")
      end if

c     Apply the symmetry boundary conditions on any coordinate axes
c     -------------------------------------------------------------
      call CartSymGN(ierr,cctkGH,"fowavetoy::scalarevolve")
      if (ierr.ne.0) call CCTK_WARN (0, "internal error")
      call CartSymGN(ierr,cctkGH,"fowavetoy::scalarevolve_derivs")
      if (ierr.ne.0) call CCTK_WARN (0, "internal error")
      
c     Apply the outer boundary conditions
c     -----------------------------------
      
      if (table.eq.-1) then
         
         call Util_TableCreateFromString (table, "LIMIT=0.0 SPEED=1.0")
         if (table.lt.0) call CCTK_WARN (0, "internal error")
         
      end if
      
      call CCTK_FortranString (fboundlen, bound, fbound)
      if (fboundlen.lt.0) call CCTK_WARN (0, "internal error")
      
      ierr = Boundary_SelectGroupForBC (cctkGH, CCTK_ALL_FACES, bound_width,
     $     table, "fowavetoy::scalarevolve", fbound)
      if (ierr.ne.0) call CCTK_WARN (0, "internal error")
      ierr = Boundary_SelectGroupForBC (cctkGH, CCTK_ALL_FACES, bound_width,
     $     table, "fowavetoy::scalarevolve_derivs", fbound)
      if (ierr.ne.0) call CCTK_WARN (0, "internal error")
      
      end