aboutsummaryrefslogtreecommitdiff
path: root/src/SetSym.F
blob: a7c331b6d74ed86838bad738079047a20abbd85c (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
197
 /*@@
   @file    SetSym.F
   @date    March 2000
   @author  Miguel Alcubierre
   @desc
            Sets the symmetries for the AHF grid functions.
   @enddesc
 @@*/

#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"

 /*@@
   @routine    AHFinder_SetSym
   @date       March 2000
   @author     Miguel Alcubierre
   @desc 
               Sets the symmetries for the AHF grid functions.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

      subroutine AHFinder_SetSym(CCTK_ARGUMENTS)

      implicit none
      
      DECLARE_CCTK_ARGUMENTS 
      DECLARE_CCTK_PARAMETERS

      integer sym(3)
      CCTK_INT CCTK_Equals
      integer, PARAMETER :: one = 1
      INTEGER :: ierr

c     +,+,+

      sym(1)=+1
      sym(2)=+1
      sym(3)=+1

      call SetCartSymVN(ierr,cctkGH,sym,'ahfinder::ahfgrid')
      call SetCartSymVN(ierr,cctkGH,sym,'ahfinder::ahf_exp')
      call SetCartSymVN(ierr,cctkGH,sym,'ahfinder::ahfgradn')
      call SetCartSymVN(ierr,cctkGH,sym,'ahfinder::ahfgauss')
      call SetCartSymVN(ierr,cctkGH,sym,'ahfinder::ahmask')

c     -,+,+

      sym(1)=-1
      sym(2)=+1
      sym(3)=+1

      call SetCartSymVN(ierr,cctkGH,sym,'ahfinder::ahfgradx')

c     +,-,+

      sym(1)=+1
      sym(2)=-1
      sym(3)=+1

      call SetCartSymVN(ierr,cctkGH,sym,'ahfinder::ahfgrady')

c     +,+,-

      sym(1)=+1
      sym(2)=+1
      sym(3)=-1

      call SetCartSymVN(ierr,cctkGH,sym,'ahfinder::ahfgradz')

c     End.

      end subroutine AHFinder_SetSym


 /*@@
   @file    SetSym.F
   @date    Jan 2003
   @author  Denis Pollney
   @desc
            Sets the reflection symmetry flags.
   @enddesc
 @@*/

      subroutine AHFinder_SetReflections(CCTK_ARGUMENTS)
      use AHFinder_dat

      implicit none

      DECLARE_CCTK_ARGUMENTS
      DECLARE_CCTK_PARAMETERS
      DECLARE_CCTK_FUNCTIONS

      CCTK_REAL, parameter :: zero = 0.0D0
      CCTK_REAL, parameter :: half = 0.5D0
      CCTK_REAL, parameter :: one  = 1.0D0

!     Find out reflection symmetries depending on
!     the values of the parameters.

      if ((CCTK_Equals(ahf_octant,"yes").eq.1).or.
     .    (CCTK_Equals(ahf_octant,"high").eq.1)) then

         refx = .true.
         refy = .true.
         refz = .true.

      else

         if (ahf_refx.ne.0) then
            refx = .true.
         else
            refx = .false.
         end if

         if (ahf_refy.ne.0) then
            refy = .true.
         else
            refy = .false.
         end if

         if (ahf_refz.ne.0) then
            refz = .true.
         else
            refz = .false.
         end if

      end if

      if (ahf_phi.eq.0) then
         refx = .true.
         refy = .true.
      end if

!     Force grid symmetries.  If the grid has some
!     symmetry (octant, quadrant, bitant), and the
!     horizon is centered on one of the symmetry
!     planes, then the grid symmetries are forced
!     on it regardless of the values of the symmetry
!     parameters.

      if (CCTK_Equals(domain,"octant").eq.1) then

         if (xc.eq.zero) refx = .true.
         if (yc.eq.zero) refy = .true.
         if (zc.eq.zero) refz = .true.

      else if (CCTK_Equals(domain,"quadrant").eq.1) then

         if (xc.eq.zero) refx = .true.
         if (yc.eq.zero) refy = .true.

      else if (CCTK_Equals(domain,"bitant").eq.1) then

         if (CCTK_Equals(bitant_plane,"xy").eq.1) then
            if (zc.eq.zero) refz = .true.
         else if (CCTK_Equals(bitant_plane,"xz").eq.1) then
            if (yc.eq.zero) refy = .true.
         else
            if (xc.eq.zero) refx = .true.
         end if

      end if

!     Check if we are using cartoon.  If we are, then the
!     corresponding symmetries are forced on the horizon.

      if (ahf_cartoon.ne.0) then
         cartoon = .true.
      else
         cartoon = .false.
      end if

      if (cartoon) then
         nonaxi = .false.
         refx = .true.
         refy = .true.
      end if

!     If desired, write values of symmetry flags to screen.

      if (veryver) then
         write(*,*)
         write(*,*) 'Symmetries used for horizon:'
         write(*,*) 'refx = ',refx
         write(*,*) 'refy = ',refy
         write(*,*) 'refz = ',refz
         write(*,*)
      end if

      end subroutine AHFinder_SetReflections