aboutsummaryrefslogtreecommitdiff
path: root/src/Dissipation_4_2.F90
blob: d2c6d89f73cf7018fa04d42f7c1fc1b04630b7e4 (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
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
! $Header$

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

subroutine dissipation_4_2 (var, ni, nj, nk, bb, gsize, offset, delta, epsilon, rhs)

  implicit none

  DECLARE_CCTK_PARAMETERS

  integer ::  ni, nj, nk
  CCTK_REAL, dimension(ni,nj,nk), intent(in) :: var
  CCTK_REAL, dimension(ni,nj,nk), intent(inout) :: rhs
  CCTK_INT, dimension(6), intent(in) :: bb
  CCTK_INT, dimension(3), intent(in) :: gsize
  CCTK_INT, dimension(6), intent(in) :: offset
  CCTK_REAL, dimension(3), intent(in) :: delta
  CCTK_REAL, intent(in) :: epsilon 

  CCTK_REAL :: zero = 0.0
  integer, parameter :: wp = kind(zero)
  CCTK_REAL, dimension(6,4) :: a
  CCTK_REAL :: idel

  CCTK_INT :: il, ir, jl, jr, kl, kr, ol, or

  call set_coeff ( a )

  if ( scale_with_h > 0 ) then
    idel = epsilon / ( 16 * delta(1) )
  else
    idel = epsilon / 16
  end if

  if ( bb(1) == 0 ) then
    il = 1 + gsize(1)
  else
    ol = offset(1)
!$omp parallel workshare
    rhs(1+ol,:,:) = rhs(1+ol,:,:) + &
                 ( a(1,1) * var(1+ol,:,:) + a(2,1) * var(2+ol,:,:) + &
                   a(3,1) * var(3+ol,:,:) ) * idel
    rhs(2+ol,:,:) = rhs(2+ol,:,:) + &
                 ( a(1,2) * var(1+ol,:,:) + a(2,2) * var(2+ol,:,:) + &
                   a(3,2) * var(3+ol,:,:) + a(4,2) * var(4+ol,:,:) ) * idel
    rhs(3+ol,:,:) = rhs(3+ol,:,:) + &
                 ( a(1,3) * var(1+ol,:,:) + a(2,3) * var(2+ol,:,:) + &
                   a(3,3) * var(3+ol,:,:) + a(4,3) * var(4+ol,:,:) + &
                   a(5,3) * var(5+ol,:,:) ) * idel
    rhs(4+ol,:,:) = rhs(4+ol,:,:) + &
                 ( a(2,4) * var(2+ol,:,:) + a(3,4) * var(3+ol,:,:) + &
                   a(4,4) * var(4+ol,:,:) + a(5,4) * var(5+ol,:,:) + &
                   a(6,4) * var(6+ol,:,:) ) * idel
!$omp end parallel workshare

    il = 5 + ol
  end if
  if ( bb(2) == 0 ) then
    ir = ni - gsize(1)
  else
    or = ni - offset(2)
!$omp parallel workshare
    rhs(or-3,:,:) = rhs(or-3,:,:) + &
                 ( a(2,4) * var(or-1,:,:) + a(3,4) * var(or-2,:,:) + &
                   a(4,4) * var(or-3,:,:) + a(5,4) * var(or-4,:,:) + &
                   a(6,4) * var(or-5,:,:) ) * idel
    rhs(or-2,:,:) = rhs(or-2,:,:) + &
                 ( a(1,3) * var(or,:,:) + a(2,3) * var(or-1,:,:) + &
                   a(3,3) * var(or-2,:,:) + a(4,3) * var(or-3,:,:) + &
                   a(5,3) * var(or-4,:,:) ) * idel
    rhs(or-1,:,:) = rhs(or-1,:,:) + &
                 ( a(1,2) * var(or,:,:) + a(2,2) * var(or-1,:,:) + &
                   a(3,2) * var(or-2,:,:) + a(4,2) * var(or-3,:,:) ) * idel
    rhs(or,:,:)   = rhs(or,:,:) + &
                 ( a(1,1) * var(or,:,:) + a(2,1) * var(or-1,:,:) + &
                   a(3,1) * var(or-2,:,:) ) * idel
!$omp end parallel workshare

    ir = or - 4
  end if
!$omp parallel workshare
  rhs(il:ir,:,:) = rhs(il:ir,:,:) + &
                   ( -6.0_wp * var(il:ir,:,:) + &
                      4.0_wp * ( var(il-1:ir-1,:,:) + &
                                 var(il+1:ir+1,:,:) ) - &
                               ( var(il-2:ir-2,:,:) + &
                                 var(il+2:ir+2,:,:) ) ) * idel
!$omp end parallel workshare

  if ( zero_derivs_y == 0 ) then
    call set_coeff ( a )
  
    if ( scale_with_h > 0 ) then
      idel = epsilon / ( 16 * delta(2) )
    else
      idel = epsilon / 16
    end if
  
    if ( bb(3) == 0 ) then
      jl = 1 + gsize(2)
    else
      ol = offset(3)
!$omp parallel workshare
      rhs(:,1+ol,:) = rhs(:,1+ol,:) + &
                   ( a(1,1) * var(:,1+ol,:) + a(2,1) * var(:,2+ol,:) + &
                     a(3,1) * var(:,3+ol,:) ) * idel
      rhs(:,2+ol,:) = rhs(:,2+ol,:) + &
                   ( a(1,2) * var(:,1+ol,:) + a(2,2) * var(:,2+ol,:) + &
                     a(3,2) * var(:,3+ol,:) + a(4,2) * var(:,4+ol,:) ) * idel
      rhs(:,3+ol,:) = rhs(:,3+ol,:) + &
                   ( a(1,3) * var(:,1+ol,:) + a(2,3) * var(:,2+ol,:) + &
                     a(3,3) * var(:,3+ol,:) + a(4,3) * var(:,4+ol,:) + &
                     a(5,3) * var(:,5+ol,:) ) * idel
      rhs(:,4+ol,:) = rhs(:,4+ol,:) + &
                   ( a(2,4) * var(:,2+ol,:) + a(3,4) * var(:,3+ol,:) + &
                     a(4,4) * var(:,4+ol,:) + a(5,4) * var(:,5+ol,:) + &
                     a(6,4) * var(:,6+ol,:) ) * idel
!$omp end parallel workshare
  
      jl = 5 + ol
    end if
    if ( bb(4) == 0 ) then
      jr = nj - gsize(2)
    else
      or = nj - offset(4)
!$omp parallel workshare
      rhs(:,or-3,:) = rhs(:,or-3,:) + &
                   ( a(2,4) * var(:,or-1,:) + a(3,4) * var(:,or-2,:) + &
                     a(4,4) * var(:,or-3,:) + a(5,4) * var(:,or-4,:) + &
                     a(6,4) * var(:,or-5,:) ) * idel
      rhs(:,or-2,:) = rhs(:,or-2,:) + &
                   ( a(1,3) * var(:,or,:) + a(2,3) * var(:,or-1,:) + &
                     a(3,3) * var(:,or-2,:) + a(4,3) * var(:,or-3,:) + &
                     a(5,3) * var(:,or-4,:) ) * idel
      rhs(:,or-1,:) = rhs(:,or-1,:) + &
                   ( a(1,2) * var(:,or,:) + a(2,2) * var(:,or-1,:) + &
                     a(3,2) * var(:,or-2,:) + a(4,2) * var(:,or-3,:) ) * idel
      rhs(:,or,:)   = rhs(:,or,:) + &
                   ( a(1,1) * var(:,or,:) + a(2,1) * var(:,or-1,:) + &
                     a(3,1) * var(:,or-2,:) ) * idel
!$omp end parallel workshare
  
      jr = or - 4
    end if
!$omp parallel workshare
    rhs(:,jl:jr,:) = rhs(:,jl:jr,:) + &
                     ( -6.0_wp * var(:,jl:jr,:) + &
                        4.0_wp * ( var(:,jl-1:jr-1,:) + &
                                   var(:,jl+1:jr+1,:) ) - &
                                 ( var(:,jl-2:jr-2,:) + &
                                   var(:,jl+2:jr+2,:) ) ) * idel
!$omp end parallel workshare
  end if

  if ( zero_derivs_z == 0 ) then
    call set_coeff ( a )
  
    if ( scale_with_h > 0 ) then
      idel = epsilon / ( 16 * delta(3) )
    else
      idel = epsilon / 16
    end if
  
    if ( bb(5) == 0 ) then
      kl = 1 + gsize(3)
    else
      ol = offset(5)
!$omp parallel workshare
      rhs(:,:,1+ol) = rhs(:,:,1+ol) + &
                   ( a(1,1) * var(:,:,1+ol) + a(2,1) * var(:,:,2+ol) + &
                     a(3,1) * var(:,:,3+ol) ) * idel
      rhs(:,:,2+ol) = rhs(:,:,2+ol) + &
                   ( a(1,2) * var(:,:,1+ol) + a(2,2) * var(:,:,2+ol) + &
                     a(3,2) * var(:,:,3+ol) + a(4,2) * var(:,:,4+ol) ) * idel
      rhs(:,:,3+ol) = rhs(:,:,3+ol) + &
                   ( a(1,3) * var(:,:,1+ol) + a(2,3) * var(:,:,2+ol) + &
                     a(3,3) * var(:,:,3+ol) + a(4,3) * var(:,:,4+ol) + &
                     a(5,3) * var(:,:,5+ol) ) * idel
      rhs(:,:,4+ol) = rhs(:,:,4+ol) + &
                   ( a(2,4) * var(:,:,2+ol) + a(3,4) * var(:,:,3+ol) + &
                     a(4,4) * var(:,:,4+ol) + a(5,4) * var(:,:,5+ol) + &
                     a(6,4) * var(:,:,6+ol) ) * idel
!$omp end parallel workshare
  
      kl = 5 + ol
    end if
    if ( bb(6) == 0 ) then
      kr = nk - gsize(3)
    else
      or = nk - offset(6)
!$omp parallel workshare
      rhs(:,:,or-3) = rhs(:,:,or-3) + &
                   ( a(2,4) * var(:,:,or-1) + a(3,4) * var(:,:,or-2) + &
                     a(4,4) * var(:,:,or-3) + a(5,4) * var(:,:,or-4) + &
                     a(6,4) * var(:,:,or-5) ) * idel
      rhs(:,:,or-2) = rhs(:,:,or-2) + &
                   ( a(1,3) * var(:,:,or) + a(2,3) * var(:,:,or-1) + &
                     a(3,3) * var(:,:,or-2) + a(4,3) * var(:,:,or-3) + &
                     a(5,3) * var(:,:,or-4) ) * idel
      rhs(:,:,or-1) = rhs(:,:,or-1) + &
                   ( a(1,2) * var(:,:,or) + a(2,2) * var(:,:,or-1) + &
                     a(3,2) * var(:,:,or-2) + a(4,2) * var(:,:,or-3) ) * idel
      rhs(:,:,or)   = rhs(:,:,or) + &
                   ( a(1,1) * var(:,:,or) + a(2,1) * var(:,:,or-1) + &
                     a(3,1) * var(:,:,or-2) ) * idel
!$omp end parallel workshare
  
      kr = or - 4
    end if
!$omp parallel workshare
    rhs(:,:,kl:kr) = rhs(:,:,kl:kr) + &
                     ( -6.0_wp * var(:,:,kl:kr) + &
                        4.0_wp * ( var(:,:,kl-1:kr-1) + &
                                   var(:,:,kl+1:kr+1) ) - &
                                 ( var(:,:,kl-2:kr-2) + &
                                   var(:,:,kl+2:kr+2) ) ) * idel
!$omp end parallel workshare
  end if

contains

  subroutine set_coeff ( a )

    implicit none

    CCTK_REAL, dimension(6,4), intent(OUT) :: a
    CCTK_REAL :: zero = 0.0
    integer, parameter :: wp = kind(zero)

    a(1,1) = -2.8235294117647058823529411764705882352941176470588_wp
    a(2,1) = 5.6470588235294117647058823529411764705882352941176_wp
    a(3,1) = -2.8235294117647058823529411764705882352941176470588_wp
    a(4,1) = 0.0_wp
    a(5,1) = 0.0_wp
    a(6,1) = 0.0_wp
    a(1,2) = 1.6271186440677966101694915254237288135593220338983_wp
    a(2,2) = -4.0677966101694915254237288135593220338983050847458_wp
    a(3,2) = 3.2542372881355932203389830508474576271186440677966_wp
    a(4,2) = -0.81355932203389830508474576271186440677966101694915_wp
    a(5,2) = 0.0_wp
    a(6,2) = 0.0_wp
    a(1,3) = -1.1162790697674418604651162790697674418604651162791_wp
    a(2,3) = 4.4651162790697674418604651162790697674418604651163_wp
    a(3,3) = -6.6976744186046511627906976744186046511627906976744_wp
    a(4,3) = 4.4651162790697674418604651162790697674418604651163_wp
    a(5,3) = -1.1162790697674418604651162790697674418604651162791_wp
    a(6,3) = 0.0_wp
    a(1,4) = 0.0_wp
    a(2,4) = -0.97959183673469387755102040816326530612244897959184_wp
    a(3,4) = 3.9183673469387755102040816326530612244897959183673_wp
    a(4,4) = -5.8775510204081632653061224489795918367346938775510_wp
    a(5,4) = 3.9183673469387755102040816326530612244897959183673_wp
    a(6,4) = -0.97959183673469387755102040816326530612244897959184_wp

  end subroutine set_coeff

end subroutine dissipation_4_2