aboutsummaryrefslogtreecommitdiff
path: root/src/linear_extrap_one_bndry.F
blob: f5ba75bbce5d033bd51d83f15cb7f83edeaf0379 (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
c this subroutine linearly extrapolates one Cactus variable
c on one boundary of the Cactus grid box
C $Header$

#include "cctk.h"
#include "cctk_Arguments.h"
#include "Exact.inc"

c #define-ing the symbol EXACT_NO_F90 will turn this file into a no-op
#ifndef EXACT_NO_F90
      subroutine Exact__linear_extrap_one_bndry(CCTK_ARGUMENTS,var)

      implicit none

      DECLARE_CCTK_ARGUMENTS

      integer nx,ny,nz

      CCTK_REAL var(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3))

C     Grid parameters.

      nx = cctk_lsh(1)
      ny = cctk_lsh(2)
      nz = cctk_lsh(3)

C     Linear extrapolation from the interiors to the boundaries.
C     Does not support octant or quadrant.

C     6 faces.

      if (cctk_bbox(1).eq.1 .and. nx.ge.4 ) then 
         var(1,:,:) = 2.d0 * var(2,:,:) - var(3,:,:)
      endif
      if (cctk_bbox(3).eq.1 .and. ny.ge.4 ) then 
         var(:,1,:) = 2.d0 * var(:,2,:) - var(:,3,:)
      endif
      if (cctk_bbox(5).eq.1 .and. nz.ge.4 ) then 
         var(:,:,1) = 2.d0 * var(:,:,2) - var(:,:,3)
      endif
      if (cctk_bbox(2).eq.1 .and. nx.ge.4 ) then 
         var(nx,:,:) = 2.d0 * var(nx-1,:,:) - var(nx-2,:,:)
      endif
      if (cctk_bbox(4).eq.1 .and. ny.ge.4 ) then 
         var(:,ny,:) = 2.d0 * var(:,ny-1,:) - var(:,ny-2,:)
      endif
      if (cctk_bbox(6).eq.1 .and. nz.ge.4 ) then 
         var(:,:,nz) = 2.d0 * var(:,:,nz-1) - var(:,:,nz-2)
      endif

C     12 edges.
C     4 round face x=min.
      if ( cctk_bbox(1).eq.1 .and. nx.ge.4 
     $     .and. cctk_bbox(3).eq.1 .and. ny.ge.4 ) then 
         var(1,1,:) = 2.d0 * var(2,2,:) - var(3,3,:)
      end if

      if ( cctk_bbox(1).eq.1 .and. nx.ge.4 
     $     .and. cctk_bbox(4).eq.1 .and. ny.ge.4 ) then 
         var(1,ny,:) = 2.d0 * var(2,ny-1,:) - var(2,ny-2,:)
      end if

      if ( cctk_bbox(1).eq.1 .and. nx.ge.4 
     $     .and. cctk_bbox(5).eq.1 .and. nz.ge.4 ) then 
         var(1,:,1) = 2.d0 * var(2,:,2) - var(3,:,3)
      end if

      if ( cctk_bbox(1).eq.1 .and. nx.ge.4 
     $     .and. cctk_bbox(6).eq.1 .and. nz.ge.4 ) then 
         var(1,:,nz) = 2.d0 * var(2,:,nz-1) - var(3,:,nz-2)
      end if

C     4 around face x=max.
      if ( cctk_bbox(2).eq.1 .and. nx.ge.4 
     $     .and. cctk_bbox(3).eq.1 .and. ny.ge.4 ) then 
         var(nx,1,:) = 2.d0 * var(nx-1,2,:) - var(nx-2,3,:)
      end if

      if ( cctk_bbox(2).eq.1 .and. nx.ge.4 
     $     .and. cctk_bbox(4).eq.1 .and. ny.ge.4 ) then 
         var(nx,ny,:) = 2.d0 * var(nx-1,ny-1,:) - var(nx-2,ny-2,:)
      end if

      if ( cctk_bbox(2).eq.1 .and. nx.ge.4 
     $     .and. cctk_bbox(5).eq.1 .and. nz.ge.4 ) then 
         var(nx,:,1) = 2.d0 * var(nx-1,:,2) - var(nx-2,:,3)
      end if

      if ( cctk_bbox(2).eq.1 .and. nx.ge.4 
     $     .and. cctk_bbox(6).eq.1 .and. nz.ge.4 ) then 
         var(nx,:,nz) = 2.d0 * var(nx-1,:,nz-1) - var(nx-2,:,nz-2)
      end if

C     Remaining 2 in y=min.
      if ( cctk_bbox(3).eq.1 .and. ny.ge.4
     $     .and. cctk_bbox(5).eq.1 .and. nz.ge.4 ) then
         var(:,1,1) = 2.d0 * var(:,2,2) - var(:,3,3)
      end if

      if ( cctk_bbox(3).eq.1 .and. ny.ge.4
     $     .and. cctk_bbox(6).eq.1 .and. nz.ge.4 ) then 
         var(:,1,nz) = 2.d0 * var(:,2,nz-1) - var(:,2,nz-2)
      end if

C     Remaining 2 in y=ymax.      
      if ( cctk_bbox(4).eq.1 .and. ny.ge.4
     $     .and. cctk_bbox(5).eq.1 .and. nz.ge.4 ) then
         var(:,ny,1) = 2.d0 * var(:,ny-1,2) - var(:,ny-2,3)  
      end if

      if ( cctk_bbox(4).eq.1 .and. ny.ge.4
     $     .and. cctk_bbox(6).eq.1 .and. nz.ge.4 ) then 
         var(:,ny,nz) = 2.d0 * var(:,ny-1,nz-1) - var(:,ny-2,nz-2)
      end if

C     8 corners.

      if (nx.ge.4 .and. ny.ge.4 .and. nz.ge.4) then
         if (cctk_bbox(1).eq.1 .and. cctk_bbox(3).eq.1 .and. cctk_bbox(5).eq.1) then
            var(1,1,1) = 2.d0*var(2,2,2) - var(3,3,3)
         end if
         if (cctk_bbox(1).eq.1 .and. cctk_bbox(3).eq.1 .and. cctk_bbox(6).eq.1) then
            var(1,1,nz) = 2.d0*var(2,2,nz-1) - var(3,3,nz-2)
         end if
         if (cctk_bbox(1).eq.1 .and. cctk_bbox(4).eq.1 .and. cctk_bbox(5).eq.1) then
            var(1,ny,1) = 2.d0*var(2,ny-1,2) - var(3,ny-2,3)
         end if
         if (cctk_bbox(1).eq.1 .and. cctk_bbox(4).eq.1 .and. cctk_bbox(6).eq.1) then
            var(1,ny,nz) = 2.d0*var(2,ny-1,nz-1) - var(3,ny-2,nz-2)
         end if
        if (cctk_bbox(2).eq.1 .and. cctk_bbox(3).eq.1 .and. cctk_bbox(5).eq.1) then
            var(nx,1,1) = 2.d0*var(nx-1,2,2) - var(nx-2,3,3)
         end if
         if (cctk_bbox(2).eq.1 .and. cctk_bbox(3).eq.1 .and. cctk_bbox(6).eq.1) then
            var(nx,1,nz) = 2.d0*var(nx-1,2,nz-1) - var(nx-2,3,nz-2)
         end if
         if (cctk_bbox(2).eq.1 .and. cctk_bbox(4).eq.1 .and. cctk_bbox(5).eq.1) then
            var(nx,ny,1) = 2.d0*var(nx-1,ny-1,2) - var(nx-2,ny-2,3)
         end if
         if (cctk_bbox(2).eq.1 .and. cctk_bbox(4).eq.1 .and. cctk_bbox(6).eq.1) then
            var(nx,ny,nz) = 2.d0*var(nx-1,ny-1,nz-1) - var(nx-2,ny-2,nz-2)
         end if
      end if


      return
      end
#endif