aboutsummaryrefslogtreecommitdiff
path: root/src/Dissipation_4_3_min_err_coeff.F90
blob: d7d011f26e8994389dfd653ab667d0b1b984ac05 (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
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
! $Header$

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


subroutine dissipation_4_3_opt (var, lsh, gsh, lbnd, bb, gsize, delta, epsdisl, dfl, rhs)

  implicit none

  DECLARE_CCTK_PARAMETERS

  CCTK_INT, dimension(3), intent(in) ::  lsh, gsh, lbnd
  CCTK_INT :: ni, nj, nk
  CCTK_REAL, dimension(lsh(1),lsh(2),lsh(3)), intent(in) :: var
  CCTK_REAL, dimension(lsh(1),lsh(2),lsh(3)), intent(inout) :: rhs
  CCTK_INT, dimension(6), intent(in) :: bb
  CCTK_INT, dimension(3), intent(in) :: gsize
  CCTK_REAL, dimension(3), intent(in) :: delta
  CCTK_REAL, intent(in) :: epsdisl, dfl

  CCTK_REAL :: zero = 0.0
  integer, parameter :: wp = kind(zero)
  integer :: i, j, k
  CCTK_REAL, dimension(:,:), allocatable :: a, d, b, h
  CCTK_REAL :: idel

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

  ni = lsh(1); nj = lsh(2); nk = lsh(3)
  allocate ( a(ni,ni), d(ni,ni), b(ni,ni), h(ni,ni) )
  a = zero; d = zero; b = zero; h = zero;

  call set_dmatrix ( d, bb(1:2) )

  call set_bmatrix ( b, bb(1:2), lsh(1), gsh(1), lbnd(1), delta(1), dfl )
  
  call set_hmatrix ( h, bb(1:2) )

  a = - transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) )

  idel = epsdisl /  16.0_wp

  if ( bb(1) == 0 ) then
    il = 1 + gsize(1)
  else
    rhs(1,:,:) = rhs(1,:,:) + &
                 ( a(1,1) * var(1,:,:) + a(2,1) * var(2,:,:) + &
                   a(3,1) * var(3,:,:) ) * idel
    rhs(2,:,:) = rhs(2,:,:) + &
                 ( a(1,2) * var(1,:,:) + a(2,2) * var(2,:,:) + &
                   a(3,2) * var(3,:,:) + a(4,2) * var(4,:,:) + &
                   a(5,2) * var(5,:,:) + a(6,2) * var(6,:,:) + &
                   a(7,2) * var(7,:,:) ) * idel
    rhs(3,:,:) = rhs(3,:,:) + &
                 ( a(1,3) * var(1,:,:) + a(2,3) * var(2,:,:) + &
                   a(3,3) * var(3,:,:) + a(4,3) * var(4,:,:) + &
                   a(5,3) * var(5,:,:) + a(6,3) * var(6,:,:) + &
                   a(7,3) * var(7,:,:) ) * idel
    rhs(4,:,:) = rhs(4,:,:) + &
                 ( a(1,4) * var(1,:,:) + a(2,4) * var(2,:,:) + &
                   a(3,4) * var(3,:,:) + a(4,4) * var(4,:,:) + &
                   a(5,4) * var(5,:,:) + a(6,4) * var(6,:,:) + &
                   a(7,4) * var(7,:,:) ) * idel
    rhs(5,:,:) = rhs(5,:,:) + &
                 ( a(1,5) * var(1,:,:) + a(2,5) * var(2,:,:) + &
                   a(3,5) * var(3,:,:) + a(4,5) * var(4,:,:) + &
                   a(5,5) * var(5,:,:) + a(6,5) * var(6,:,:) + &
                   a(7,5) * var(7,:,:) ) * idel
    il = 6
  end if
  if ( bb(2) == 0 ) then
    ir = ni - gsize(1)
  else
    rhs(ni-4,:,:) = rhs(ni-4,:,:) + &
                    ( a(ni-6,ni-4) * var(ni-6,:,:) + &
                      a(ni-5,ni-4) * var(ni-5,:,:) + &
                      a(ni-4,ni-4) * var(ni-4,:,:) + &
                      a(ni-3,ni-4) * var(ni-3,:,:) + &
                      a(ni-2,ni-4) * var(ni-2,:,:) + &
                      a(ni-1,ni-4) * var(ni-1,:,:) + &
                      a(ni,ni-4) * var(ni,:,:) ) * idel
    rhs(ni-3,:,:) = rhs(ni-3,:,:) + &
                    ( a(ni-6,ni-3) * var(ni-6,:,:) + &
                      a(ni-5,ni-3) * var(ni-5,:,:) + &
                      a(ni-4,ni-3) * var(ni-4,:,:) + &
                      a(ni-3,ni-3) * var(ni-3,:,:) + &
                      a(ni-2,ni-3) * var(ni-2,:,:) + &
                      a(ni-1,ni-3) * var(ni-1,:,:) + &
                      a(ni,ni-3) * var(ni,:,:) ) * idel
    rhs(ni-2,:,:) = rhs(ni-2,:,:) + &
                    ( a(ni-6,ni-2) * var(ni-6,:,:) + &
                      a(ni-5,ni-2) * var(ni-5,:,:) + &
                      a(ni-4,ni-2) * var(ni-4,:,:) + &
                      a(ni-3,ni-2) * var(ni-3,:,:) + &
                      a(ni-2,ni-2) * var(ni-2,:,:) + &
                      a(ni-1,ni-2) * var(ni-1,:,:) + &
                      a(ni,ni-2) * var(ni,:,:) ) * idel
    rhs(ni-1,:,:) = rhs(ni-1,:,:) + &
                    ( a(ni-6,ni-1) * var(ni-6,:,:) + &
                      a(ni-5,ni-1) * var(ni-5,:,:) + &
                      a(ni-4,ni-1) * var(ni-4,:,:) + &
                      a(ni-3,ni-1) * var(ni-3,:,:) + &
                      a(ni-2,ni-1) * var(ni-2,:,:) + &
                      a(ni-1,ni-1) * var(ni-1,:,:) + &
                      a(ni,ni-1) * var(ni,:,:) ) * idel
    rhs(ni,:,:)   = rhs(ni,:,:) + &
                    ( a(ni-2,ni) * var(ni-2,:,:) + &
                      a(ni-1,ni) * var(ni-1,:,:) + &
                      a(ni,ni) * var(ni,:,:) ) * idel
    ir = ni - 5
  end if
  do i = il, ir
    rhs(i,:,:) = rhs(i,:,:) + &
                   ( a(i-2,i) * var(i-2,:,:) + &
                     a(i-1,i) * var(i-1,:,:) + &
                     a(i,i) * var(i,:,:) + &
                     a(i+1,i) * var(i+1,:,:) + &
                     a(i+2,i) * var(i+2,:,:) ) * idel
  end do

  deallocate ( a, d, b, h )

  allocate ( a(nj,nj), d(nj,nj), b(nj,nj), h(nj,nj) )
  a = zero; d = zero; b = zero; h = zero;

  call set_dmatrix ( d, bb(3:4) )

  call set_bmatrix ( b, bb(3:4), lsh(2), gsh(2), lbnd(2), delta(2), dfl )
  
  call set_hmatrix ( h, bb(3:4) )

  a = - transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) )

  idel = epsdisl /  16.0_wp

  if ( bb(3) == 0 ) then
    jl = 1 + gsize(2)
  else
    rhs(:,1,:) = rhs(:,1,:) + &
                 ( a(1,1) * var(:,1,:) + a(2,1) * var(:,2,:) + &
                   a(3,1) * var(:,3,:) ) * idel
    rhs(:,2,:) = rhs(:,2,:) + &
                 ( a(1,2) * var(:,1,:) + a(2,2) * var(:,2,:) + &
                   a(3,2) * var(:,3,:) + a(4,2) * var(:,4,:) + &
                   a(5,2) * var(:,5,:) + a(6,2) * var(:,6,:) + &
                   a(7,2) * var(:,7,:) ) * idel
    rhs(:,3,:) = rhs(:,3,:) + &
                 ( a(1,3) * var(:,1,:) + a(2,3) * var(:,2,:) + &
                   a(3,3) * var(:,3,:) + a(4,3) * var(:,4,:) + &
                   a(5,3) * var(:,5,:) + a(6,3) * var(:,6,:) + &
                   a(7,3) * var(:,7,:) ) * idel
    rhs(:,4,:) = rhs(:,4,:) + &
                 ( a(1,4) * var(:,1,:) + a(2,4) * var(:,2,:) + &
                   a(3,4) * var(:,3,:) + a(4,4) * var(:,4,:) + &
                   a(5,4) * var(:,5,:) + a(6,4) * var(:,6,:) + &
                   a(7,4) * var(:,7,:) ) * idel
    rhs(:,5,:) = rhs(:,5,:) + &
                 ( a(1,5) * var(:,1,:) + a(2,5) * var(:,2,:) + &
                   a(3,5) * var(:,3,:) + a(4,5) * var(:,4,:) + &
                   a(5,5) * var(:,5,:) + a(6,5) * var(:,6,:) + &
                   a(7,5) * var(:,7,:) ) * idel
    jl = 6
  end if
  if ( bb(4) == 0 ) then
    jr = nj - gsize(2)
  else
    rhs(:,nj-4,:) = rhs(:,nj-4,:) + &
                    ( a(nj-6,nj-4) * var(:,nj-6,:) + &
                      a(nj-5,nj-4) * var(:,nj-5,:) + &
                      a(nj-4,nj-4) * var(:,nj-4,:) + &
                      a(nj-3,nj-4) * var(:,nj-3,:) + &
                      a(nj-2,nj-4) * var(:,nj-2,:) + &
                      a(nj-1,nj-4) * var(:,nj-1,:) + &
                      a(nj,nj-4) * var(:,nj,:) ) * idel
    rhs(:,nj-3,:) = rhs(:,nj-3,:) + &
                    ( a(nj-6,nj-3) * var(:,nj-6,:) + &
                      a(nj-5,nj-3) * var(:,nj-5,:) + &
                      a(nj-4,nj-3) * var(:,nj-4,:) + &
                      a(nj-3,nj-3) * var(:,nj-3,:) + &
                      a(nj-2,nj-3) * var(:,nj-2,:) + &
                      a(nj-1,nj-3) * var(:,nj-1,:) + &
                      a(nj,nj-3) * var(:,nj,:) ) * idel
    rhs(:,nj-2,:) = rhs(:,nj-2,:) + &
                    ( a(nj-6,nj-2) * var(:,nj-6,:) + &
                      a(nj-5,nj-2) * var(:,nj-5,:) + &
                      a(nj-4,nj-2) * var(:,nj-4,:) + &
                      a(nj-3,nj-2) * var(:,nj-3,:) + &
                      a(nj-2,nj-2) * var(:,nj-2,:) + &
                      a(nj-1,nj-2) * var(:,nj-1,:) + &
                      a(nj,nj-2) * var(:,nj,:) ) * idel
    rhs(:,nj-1,:) = rhs(:,nj-1,:) + &
                    ( a(nj-6,nj-1) * var(:,nj-6,:) + &
                      a(nj-5,nj-1) * var(:,nj-5,:) + &
                      a(nj-4,nj-1) * var(:,nj-4,:) + &
                      a(nj-3,nj-1) * var(:,nj-3,:) + &
                      a(nj-2,nj-1) * var(:,nj-2,:) + &
                      a(nj-1,nj-1) * var(:,nj-1,:) + &
                      a(nj,nj-1) * var(:,nj,:) ) * idel
    rhs(:,nj,:)   = rhs(:,nj,:) + &
                    ( a(nj-2,nj) * var(:,nj-2,:) + &
                      a(nj-1,nj) * var(:,nj-1,:) + &
                      a(nj,nj) * var(:,nj,:) ) * idel
    jr = nj - 5
  end if

  do j = jl, jr
    rhs(:,j,:) = rhs(:,j,:) + &
                   ( a(j-2,j) * var(:,j-2,:) + &
                     a(j-1,j) * var(:,j-1,:) + &
                     a(j,j) * var(:,j,:) + &
                     a(j+1,j) * var(:,j+1,:) + &
                     a(j+2,j) * var(:,j+2,:) ) * idel
  end do

  deallocate ( a, d, b, h )

  allocate ( a(nk,nk), d(nk,nk), b(nk,nk), h(nk,nk) )
  a = zero; d = zero; b = zero; h = zero;

  call set_dmatrix ( d, bb(5:6) )

  call set_bmatrix ( b, bb(5:6), lsh(3), gsh(3), lbnd(3), delta(3), dfl )

  call set_hmatrix ( h, bb(5:6) )

  a = - transpose ( matmul ( h, matmul ( transpose(d), matmul ( b, d ) ) ) )

    idel = epsdisl /  16.0_wp


  if ( bb(5) == 0 ) then
    kl = 1 + gsize(3)
  else
    rhs(:,:,1) = rhs(:,:,1) + &
                 ( a(1,1) * var(:,:,1) + a(2,1) * var(:,:,2) + &
                   a(3,1) * var(:,:,3) ) * idel
    rhs(:,:,2) = rhs(:,:,2) + &
                 ( a(1,2) * var(:,:,1) + a(2,2) * var(:,:,2) + &
                   a(3,2) * var(:,:,3) + a(4,2) * var(:,:,4) + &
                   a(5,2) * var(:,:,5) + a(6,2) * var(:,:,6) + &
                   a(7,2) * var(:,:,7) ) * idel
    rhs(:,:,3) = rhs(:,:,3) + &
                 ( a(1,3) * var(:,:,1) + a(2,3) * var(:,:,2) + &
                   a(3,3) * var(:,:,3) + a(4,3) * var(:,:,4) + &
                   a(5,3) * var(:,:,5) + a(6,3) * var(:,:,6) + &
                   a(7,3) * var(:,:,7) ) * idel
    rhs(:,:,4) = rhs(:,:,4) + &
                 ( a(1,4) * var(:,:,1) + a(2,4) * var(:,:,2) + &
                   a(3,4) * var(:,:,3) + a(4,4) * var(:,:,4) + &
                   a(5,4) * var(:,:,5) + a(6,4) * var(:,:,6) + &
                   a(7,4) * var(:,:,7) ) * idel
    rhs(:,:,5) = rhs(:,:,5) + &
                 ( a(1,5) * var(:,:,1) + a(2,5) * var(:,:,2) + &
                   a(3,5) * var(:,:,3) + a(4,5) * var(:,:,4) + &
                   a(5,5) * var(:,:,5) + a(6,5) * var(:,:,6) + &
                   a(7,5) * var(:,:,7) ) * idel
    kl = 6
  end if
  if ( bb(6) == 0 ) then
    kr = nk - gsize(3)
  else
    rhs(:,:,nk-4) = rhs(:,:,nk-4) + &
                    ( a(nk-6,nk-4) * var(:,:,nk-6) + &
                      a(nk-5,nk-4) * var(:,:,nk-5) + &
                      a(nk-4,nk-4) * var(:,:,nk-4) + &
                      a(nk-3,nk-4) * var(:,:,nk-3) + &
                      a(nk-2,nk-4) * var(:,:,nk-2) + &
                      a(nk-1,nk-4) * var(:,:,nk-1) + &
                      a(nk,nk-4) * var(:,:,nk) ) * idel
    rhs(:,:,nk-3) = rhs(:,:,nk-3) + &
                    ( a(nk-6,nk-3) * var(:,:,nk-6) + &
                      a(nk-5,nk-3) * var(:,:,nk-5) + &
                      a(nk-4,nk-3) * var(:,:,nk-4) + &
                      a(nk-3,nk-3) * var(:,:,nk-3) + &
                      a(nk-2,nk-3) * var(:,:,nk-2) + &
                      a(nk-1,nk-3) * var(:,:,nk-1) + &
                      a(nk,nk-3) * var(:,:,nk) ) * idel
    rhs(:,:,nk-2) = rhs(:,:,nk-2) + &
                    ( a(nk-6,nk-2) * var(:,:,nk-6) + &
                      a(nk-5,nk-2) * var(:,:,nk-5) + &
                      a(nk-4,nk-2) * var(:,:,nk-4) + &
                      a(nk-3,nk-2) * var(:,:,nk-3) + &
                      a(nk-2,nk-2) * var(:,:,nk-2) + &
                      a(nk-1,nk-2) * var(:,:,nk-1) + &
                      a(nk,nk-2) * var(:,:,nk) ) * idel
    rhs(:,:,nk-1) = rhs(:,:,nk-1) + &
                    ( a(nk-6,nk-1) * var(:,:,nk-6) + &
                      a(nk-5,nk-1) * var(:,:,nk-5) + &
                      a(nk-4,nk-1) * var(:,:,nk-4) + &
                      a(nk-3,nk-1) * var(:,:,nk-3) + &
                      a(nk-2,nk-1) * var(:,:,nk-2) + &
                      a(nk-1,nk-1) * var(:,:,nk-1) + &
                      a(nk,nk-1) * var(:,:,nk) ) * idel
    rhs(:,:,nk)   = rhs(:,:,nk) + &
                    ( a(nk-2,nk) * var(:,:,nk-2) + &
                      a(nk-1,nk) * var(:,:,nk-1) + &
                      a(nk,nk) * var(:,:,nk) ) * idel
    kr = nk - 5
  end if

  do k = kl, kr
    rhs(:,:,k) = rhs(:,:,k) + &
                   ( a(k-2,k) * var(:,:,k-2) + &
                     a(k-1,k) * var(:,:,k-1) + &
                     a(k,k) * var(:,:,k) + &
                     a(k+1,k) * var(:,:,k+1) + &
                     a(k+2,k) * var(:,:,k+2) ) * idel
  end do

  deallocate ( a, d, b, h )

contains

  subroutine set_dmatrix ( d, bb )

    implicit none

    CCTK_REAL, dimension(:,:), intent(out) :: d
    CCTK_INT, dimension(2), intent(in) :: bb
    CCTK_INT :: n
    CCTK_REAL, dimension(3), save :: ac = (/ 1.0_wp, -2.0_wp, 1.0_wp /)
    CCTK_INT :: i

    n = size(d,1)
    if ( bb(1) == 0 ) then
      d(1,1:2) = ac(2:3)
    else
      d(1,1:3) = ac
    end if
    if ( bb(2) == 0 ) then
      d(n,n-1:n) = ac(1:2)
    else
      d(n,n-2:n) = ac
    end if
    do i = 2, n-1
      d(i,i-1:i+1) = ac
    end do
  end subroutine set_dmatrix

  subroutine set_bmatrix ( b, bb, lsh, gsh, lbnd, h, dfl )

    implicit none

    CCTK_REAL, dimension(:,:), intent(out) :: b
    CCTK_INT, dimension(2), intent(in) :: bb
    CCTK_INT, intent(in) :: lsh, gsh, lbnd
    CCTK_REAL, intent(in) :: h, dfl
    CCTK_INT :: n, nwt
    CCTK_INT :: i
    CCTK_REAL :: xp, xmax

    CCTK_REAL :: zero = 0.0
    integer, parameter :: wp = kind(zero)

    n = size(b,1)
    xp = ( gsh - 1 ) * dfl - 1
    do i = 1, n
      b(i,i) = bf ( real(i+lbnd-1,wp), xp, real(gsh-1,wp), h )
    end do

  end subroutine set_bmatrix
    
  subroutine set_hmatrix ( h, bb )

    implicit none

    CCTK_REAL, dimension(:,:), intent(out) :: h
    CCTK_INT, dimension(2), intent(in) :: bb
    CCTK_INT :: n, i

    CCTK_REAL :: zero = 0.0
    integer, parameter :: wp = kind(zero)

    n = size(h,1)
    do i = 1, n
      h(i,i) = 1.0_wp
    end do
    if ( bb(1) /= 0 ) then
      h(1:5,1:5) = sigma ( )
    end if
    if ( bb(2) /= 0 ) then
      h(n:n-4:-1,n:n-4:-1) = sigma ( )
    end if

  end subroutine set_hmatrix
    
  function bf ( x, xp, xmax, h )

    implicit none

    CCTK_REAL :: bf
    CCTK_REAL, intent(in) :: x, xp, xmax, h

    CCTK_REAL :: zero = 0.0
    integer, parameter :: wp = kind(zero)

    if ( x < xp ) then
      bf = ( 1.0_wp - h ) / xp * x + h
    else if ( x > xmax - xp ) then
      bf = ( h - 1.0_wp ) / xp * ( x + xp - xmax ) + 1
    else
      bf = 1.0_wp
    end if

  end function bf

  function sigma ( )

    implicit none

    CCTK_REAL, dimension(5,5) :: sigma

    CCTK_REAL :: zero = 0.0
    integer, parameter :: wp = kind(zero)

    sigma(1,1) = 4.186595370392226897362216859769846226369_wp
    sigma(2,1) = 0_wp
    sigma(3,1) = 0_wp
    sigma(4,1) = 0_wp
    sigma(5,1) = 0_wp
    sigma(1,2) = 0_wp
    sigma(2,2) = 0.6725191921225620731888714836983116420871_wp
    sigma(3,2) = 0.3613418181134949259370502966736306984367_wp
    sigma(4,2) = -0.2021316117293899791481674539631879662707_wp
    sigma(5,2) = 0.03455320708729270824077678274955265350304_wp
    sigma(1,3) = 0_wp
    sigma(2,3) = 0.3613418181134949259370502966736306984367_wp
    sigma(3,3) = 0.7206133711630147057720442098623847362950_wp
    sigma(4,3) = 0.1376472340546569368321616389764958792591_wp
    sigma(5,3) = -0.04136405531324488624637892257286207044784_wp
    sigma(1,4) = 0_wp
    sigma(2,4) = -0.2021316117293899791481674539631879662707_wp
    sigma(3,4) = 0.1376472340546569368321616389764958792591_wp
    sigma(4,4) = 0.9578653607931026822074133441449909689509_wp
    sigma(5,4) = 0.02069353627247161734563597102894256809696_wp
    sigma(1,5) = 0_wp
    sigma(2,5) = 0.03455320708729270824077678274955265350304_wp
    sigma(3,5) = -0.04136405531324488624637892257286207044784_wp
    sigma(4,5) = 0.02069353627247161734563597102894256809696_wp
    sigma(5,5) = 0.9908272703370861473007798925906968380654_wp

  end function sigma

  subroutine set_coeff ( a, dx )

    implicit none

    CCTK_REAL, dimension(7,5), intent(OUT) :: a
    CCTK_REAL, intent(in) :: dx
    CCTK_REAL :: zero = 0.0
    integer, parameter :: wp = kind(zero)
    CCTK_REAL, dimension(4) :: b

    b(1) = dx
    b(2) = ( 1 + 3*b(1) ) / 4
    b(3) = ( 1 + b(1) ) / 2
    b(4) = ( 3 + b(1) )/4

    a(1,1) = -2.001729494283060564327579028280434295231079652820_wp*b(1)-2.001729494283060564327579028280434295231079652820_wp*b(2)
    a(2,1) = 4.0034589885661211286551580565608685904621593056401_wp*b(1)+4.0034589885661211286551580565608685904621593056401_wp*b(2)
    a(3,1) = -2.001729494283060564327579028280434295231079652820_wp*b(1)-2.001729494283060564327579028280434295231079652820_wp*b(2)
    a(4,1) = 0.0_wp
    a(5,1) = 0.0_wp
    a(6,1) = 0.0_wp
    a(7,1) = 0.0_wp
    a(1,2) = 33.349009885955840248954469595024528895341083236112_wp*b(1)+33.349009885955840248954469595024528895341083236112_wp*b(2)
    a(2,2) = -66.698019771911680497908939190049057790682166472223_wp*b(1)-66.698019771911680497908939190049057790682166472223_wp*b(2)-14.236744515267393430983451831106076148131079582685_wp*b(3)
    a(3,2) = 33.349009885955840248954469595024528895341083236112_wp*b(1)+33.349009885955840248954469595024528895341083236112_wp*b(2)+28.473489030534786861966903662212152296262159165371_wp*b(3)-5.0440331538485788446043928053469767557227163958272_wp*b(4)
    a(4,2) = -1.3951356569032137488502820136180622005437534045293_wp-14.236744515267393430983451831106076148131079582685_wp*b(3)+10.088066307697157689208785610693953511445432791654_wp*b(4)
    a(5,2) = 5.9576955605539351261459121347762848951403545560302_wp-5.0440331538485788446043928053469767557227163958272_wp*b(4)
    a(6,2) = -7.7299841503982290057409782286983831886494488984724_wp
    a(7,2) = 3.1674242467475076284453481075401604940528477469716_wp
    a(1,3) = -3.9903407976170909924931618645027940552738274990905_wp*b(1)-3.9903407976170909924931618645027940552738274990905_wp*b(2)
    a(2,3) = 7.9806815952341819849863237290055881105476549981809_wp*b(1)+7.9806815952341819849863237290055881105476549981809_wp*b(2)+2.0281925060753259562659269077917853627407035645458_wp*b(3)
    a(3,3) = -3.9903407976170909924931618645027940552738274990905_wp*b(1)-3.9903407976170909924931618645027940552738274990905_wp*b(2)-4.0563850121506519125318538155835707254814071290916_wp*b(3)+0.85749863949531149345224149788443686605331530370993_wp*b(4)
    a(4,3) = -0.17805296968581960088505371131995522597503148483988_wp+2.0281925060753259562659269077917853627407035645458_wp*b(3)-1.7149972789906229869044829957688737321066306074199_wp*b(4)
    a(5,3) = -0.029159711483604197600932897250942820943416950674515_wp+0.85749863949531149345224149788443686605331530370993_wp*b(4)
    a(6,3) = 0.59247833202466719785702692846175131981192835586867_wp
    a(7,3) = -0.38526565085524339937104031989085327289347992035428_wp
    a(1,4) = -10.828009944579909415707962754026179346885875504376_wp*b(1)-10.828009944579909415707962754026179346885875504376_wp*b(2)
    a(2,4) = 21.656019889159818831415925508052358693771751008752_wp*b(1)+21.656019889159818831415925508052358693771751008752_wp*b(2)+4.8030111910833762941985823652201696605573253856589_wp*b(3)
    a(3,4) = -10.828009944579909415707962754026179346885875504376_wp*b(1)-10.828009944579909415707962754026179346885875504376_wp*b(2)-9.6060223821667525883971647304403393211146507713177_wp*b(3)+1.8413742894937777907392755981492717093567250921496_wp*b(4)
    a(4,4) = 0.52169344312294565855548171683437104934836522915934_wp+4.8030111910833762941985823652201696605573253856589_wp*b(3)-3.6827485789875555814785511962985434187134501842992_wp*b(4)
    a(5,4) = -2.3211687019578831529367545613086339586352127503243_wp+1.8413742894937777907392755981492717093567250921496_wp*b(4)
    a(6,4) = 3.0772570745469293302070639721141547692253298131705_wp
    a(7,4) = -1.2777818157119918358257911276398918599384822920056_wp
    a(1,5) = -6.7201141443502586562617365349711742609991754142974_wp*b(1)-6.7201141443502586562617365349711742609991754142974_wp*b(2)
    a(2,5) = 13.440228288700517312523473069942348521998350828595_wp*b(1)+13.440228288700517312523473069942348521998350828595_wp*b(2)+2.6601737327460025913616376196819751799013252956745_wp*b(3)
    a(3,5) = -6.7201141443502586562617365349711742609991754142974_wp*b(1)-6.7201141443502586562617365349711742609991754142974_wp*b(2)-5.3203474654920051827232752393639503598026505913491_wp*b(3)+0.68697583481584777698277591803888262651185956794776_wp*b(4)
    a(4,5) = 1.6888624757937931547697409070602037810047678994127_wp+2.6601737327460025913616376196819751799013252956745_wp*b(3)-1.3739516696316955539655518360777652530237191358955_wp*b(4)
    a(5,5) = -4.8610470973404788048372478314704553824811608945344_wp+0.68697583481584777698277591803888262651185956794776_wp*b(4)
    a(6,5) = 4.6555067672995781453652729417602994219480180908309_wp
    a(7,5) = -1.4833221457528924952977660173500478204716250957091_wp

  end subroutine set_coeff
end subroutine dissipation_4_3_opt