aboutsummaryrefslogtreecommitdiff
path: root/src/Coefficients_4_3.F90
blob: 89c824ddf7cb98622d3a6186f8ac3eca1f98662a (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
#include "cctk.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"

subroutine set_coeff_4_3 ( nsize, loc_order, bb, gsize, imin, imax, dd )

  implicit none

  DECLARE_CCTK_FUNCTIONS
  DECLARE_CCTK_PARAMETERS

  CCTK_REAL, parameter :: zero = 0.0
  integer, parameter :: wp = kind(zero)
  CCTK_INT, intent(IN) :: nsize, loc_order
  CCTK_INT, dimension(2), intent(IN) :: bb
  CCTK_INT, intent(IN) :: gsize
  CCTK_INT, dimension(nsize), intent(OUT) :: imin, imax
  CCTK_REAL, dimension(nsize,nsize), intent(OUT) :: dd

  CCTK_REAL, dimension(2), save :: a
  CCTK_REAL, dimension(7,5), save :: q 
  CCTK_REAL :: f0, f1, f2, f3

  CCTK_INT :: i, il, ir

  logical, save :: first = .true.

  if ( first ) then
    a(1) = 2.0_wp/3.0_wp; a(2) = -1.0_wp/12.0_wp

    f0 = sqrt(26116897.0_wp)
    f1 = -56764003702447356523.0_wp + 8154993476273221.0_wp * f0
    f2 = -55804550303.0_wp + 9650225.0_wp * f0
    f3 = 3262210757.0_wp + 271861.0_wp * f0

    q(1,1) = -11._wp/6.0_wp; q(2,1) = 3.0_wp
    q(3,1) = -3.0_wp/2.0_wp; q(4,1) = 1.0_wp/3.0_wp
    q(5,1) = zero; q(6,1) = zero; q(7,1) = zero

    q(1,2) =  -24.0_wp * ( -779042810827742869.0_wp + &
                            104535124033147.0_wp * f0 ) / f1
    q(2,2) = -( -176530817412806109689.0_wp + &
                 29768274816875927.0_wp * f0 ) / ( 6.0_wp * f1 )
    q(3,2) = 343.0_wp * ( -171079116122226871.0_wp + &
                           27975630462649.0_wp * f0 ) / f1
    q(4,2) = -3.0_wp * ( -7475554291248533227.0_wp + &
                          1648464218793925.0_wp * f0 ) / ( 2.0_wp * f1 )
    q(5,2) = ( -2383792768180030915.0_wp + &
                1179620587812973.0_wp * f0 ) / ( 3.0_wp * f1 )
    q(6,2) = -1232.0_wp * ( -115724529581315.0_wp + 37280576429.0_wp * f0 ) / f1
    q(7,2) = zero

    q(1,3) = -12.0_wp * ( -380966843.0_wp + 86315.0_wp * f0 ) / f2
    q(2,3) = ( 5024933015.0_wp + 2010631.0_wp * f0 ) / ( 3.0_wp * f2 )
    q(3,3) = -231.0_wp * ( -431968921.0_wp + 86711.0_wp * f0 ) / ( 2.0_wp * f2 )
    q(4,3) = ( -65931742559.0_wp + 12256337.0_wp * f0 ) / f2
    q(5,3) = -( -50597298167.0_wp + 9716873.0_wp * f0 ) / ( 6.0_wp * f2 )
    q(6,3) = -88.0_wp * ( -15453061.0_wp + 2911.0_wp * f0 ) / f2
    q(7,3) = zero

    q(1,4) = 48.0_wp * ( -56020909845192541.0_wp + &
                          9790180507043.0_wp * f0 ) / f1
    q(2,4) = ( -9918249049237586011.0_wp + &
                1463702013196501.0_wp * f0 ) / ( 6.0_wp * f1)
    q(3,4) = -13.0_wp * ( -4130451756851441723.0_wp + &
                           664278707201077.0_wp * f0 ) / f1
    q(4,4) = 3.0_wp * ( -26937108467782666617.0_wp + &
                         5169063172799767.0_wp * f0 ) / ( 2.0_wp * f1 )
    q(5,4) = -( 6548308508012371315.0_wp + &
                3968886380989379.0_wp * f0 ) / ( 3.0_wp * f1 )
    q(6,4) = 88.0_wp * ( -91337851897923397.0_wp + &
                          19696768305507.0_wp * f0 ) / f1
    q(7,4) = 242.0_wp * ( -120683.0_wp + 15.0_wp * f0 ) / f3

    q(1,5) = 264.0_wp * ( -120683.0_wp + 15.0_wp * f0 ) / f3
    q(2,5) = ( -43118111.0_wp + 23357.0_wp * f0 ) / ( 3.0_wp * f3 )
    q(3,5) = -47.0_wp * ( -28770085.0_wp + 2259.0_wp * f0 ) / ( 2.0_wp * f3 ) 
    q(4,5) = -3.0_wp * ( 1003619433.0_wp + 11777.0_wp * f0 ) / f3
    q(5,5) = -11.0_wp * ( -384168269.0_wp + 65747.0_wp * f0 ) / ( 6.0_wp * f3 )
    q(6,5) = 22.0_wp * ( 87290207.0_wp + 10221.0_wp * f0 ) / f3
    q(7,5) = -66.0_wp * ( 3692405.0_wp + 419.0_wp * f0 ) / f3
    first = .false.
  end if

  dd = zero
  imin = 0
  imax = -1

  if ( bb(1) == 0 ) then
    il = 1 + gsize
  else
    dd(1:4,1) = q(1:4,1)
    imin(1) = 1; imax(1) = 4

    dd(1:6,2) = q(1:6,2)
    imin(2) = 1; imax(2) = 6

    dd(1:6,3) = q(1:6,3)
    imin(3) = 1; imax(3) = 6

    dd(1:7,4) = q(1:7,4)
    imin(4) = 1; imax(4) = 7

    dd(1:7,5) = q(1:7,5)
    imin(5) = 1; imax(5) = 7

    il = 6
  end if
  if ( bb(2) == 0 ) then
    ir = nsize - gsize
  else
    dd(nsize-6:nsize,nsize-4) = -q(7:1:-1,5)
    imin(nsize-4) = nsize-6; imax(nsize-4) = nsize

    dd(nsize-6:nsize,nsize-3) = -q(7:1:-1,4)
    imin(nsize-3) = nsize-6; imax(nsize-3) = nsize
   
    dd(nsize-5:nsize,nsize-2) = -q(6:1:-1,3)
    imin(nsize-2) = nsize-5; imax(nsize-2) = nsize

    dd(nsize-5:nsize,nsize-1) = -q(6:1:-1,2)
    imin(nsize-1) = nsize-5; imax(nsize-1) = nsize

    dd(nsize-3:nsize,nsize) = -q(4:1:-1,1)
    imin(nsize) = nsize-3; imax(nsize) = nsize

    ir = nsize - 5
  end if
  do i = il, ir
    dd(i-2:i-1,i) = -a(2:1:-1); dd(i+1:i+2,i) = a(1:2)
    imin(i) = i-2; imax(i) = i+2
  end do
end subroutine set_coeff_4_3