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

subroutine set_coeff_6_5 ( 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(3), save :: a
  CCTK_REAL, dimension(10,7), save :: q 

  CCTK_INT :: i, il, ir

  logical, save :: first = .true.

  if ( first ) then
    a(1) = 3.0_wp/4.0_wp; a(2) = -3.0_wp/20.0_wp; a(3) = 1.0_wp/60.0_wp

    q(1,1) = -137.0_wp/60.0_wp; q(2,1) = 5.0_wp; q(3,1) = -5.0_wp;
    q(4,1) = 10.0_wp/3.0_wp; q(5,1) = -5.0_wp/4.0_wp; q(6,1) = 1.0_wp/5.0_wp;
    q(7,1) = zero; q(8,1) = zero; q(9,1) = zero; q(10,1) = zero;

    q(1,2) = -0.3209098489859668118347_wp
    q(2,2) = -0.359441670715467075894_wp
    q(3,2) = 0.195756852998105503889_wp
    q(4,2) = 1.394685510250317033169_wp
    q(5,2) = -1.448965775497476572821_wp
    q(6,2) = 0.6519476244467816674835_wp
    q(7,2) = -0.1115052611983591304249_wp
    q(8,2) = -0.00156743129793461356831_wp
    q(9,2) = zero; q(10,2) = zero

    q(1,3) = 0.0980796259621024870281_wp
    q(2,3) = -0.784752101846956387607_wp
    q(3,3) = 0.364206221878666751300_wp
    q(4,3) = 0.102097753636344358379_wp
    q(5,3) = 0.377167650934576412896_wp
    q(6,3) = -0.173241400242683302020_wp
    q(7,3) = 0.0062120424243610783650_wp
    q(8,3) = 0.01153111791917461507898_wp
    q(9,3) = -0.00130091066558601341956_wp
    q(10,3) = zero

    q(1,4) = 0.0771947815148234674634_wp
    q(2,4) = -0.4124719432622107102831_wp
    q(3,4) = 0.654798717867622763906_wp
    q(4,4) = -1.873474657556531216053_wp
    q(5,4) = 2.161961221935152129759_wp
    q(6,4) = -0.729147814207596325277_wp
    q(7,4) = 0.1292509053479867634674_wp
    q(8,4) = -0.01092898521375837831427_wp
    q(9,4) = 0.00316983426828354261009_wp
    q(10,4) = -0.000352060693772037278516_wp

    q(1,5) = -0.01135836398682271847110_wp
    q(2,5) = 0.0114974985688122284189_wp
    q(3,5) = 0.228765565489867047381_wp
    q(4,5) = -1.179122428698534075609_wp
    q(5,5) = 0.774618306056840861215_wp
    q(6,5) = 0.016612403904600041087_wp
    q(7,5) = 0.2399304273448163041220_wp
    q(8,5) = -0.0959180672407791626759_wp
    q(9,5) = 0.01612460763754670343830_wp
    q(10,5) = -0.001149949076347228906108_wp

    q(1,6) = -0.01912046962023373131259_wp
    q(2,6) = 0.1569245474161812198493_wp
    q(3,6) = -0.567250839387340928718_wp
    q(4,6) = 1.230411310316952671353_wp
    q(5,6) = -2.048795717924254691221_wp
    q(6,6) = 0.982715833107347313532_wp
    q(7,6) = 0.2876584784911885215053_wp
    q(8,6) = -0.0207657038198929292487_wp
    q(9,6) = -0.00335290957957120811713_wp
    q(10,6) = 0.001575470999623762377395_wp

    q(1,7) = -0.0076072322555736921424_wp 
    q(2,7) = 0.034503414345314189270_wp
    q(3,7) = -0.044377233426986476366_wp
    q(4,7) = -0.049408595803628143797_wp
    q(5,7) = 0.304693537653175477283_wp
    q(6,7) = -0.935863321236636595622_wp
    q(7,7) = 0.111375967229143365899_wp 
    q(8,7) = 0.7149322477536284813995_wp
    q(9,7) = -0.1444768161656941519630_wp
    q(10,7) = 0.01622803190725754603783_wp

    first = .false.
  end if

  dd = zero
  imin = 0
  imax = -1

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

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

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

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

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

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

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

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

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

    dd(nsize-9:nsize,nsize-4) = -q(10:1:-1,5)
    imin(nsize-4) = nsize-9; imax(nsize-4) = nsize

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

    dd(nsize-8:nsize,nsize-2) = -q(9:1:-1,3)
    imin(nsize-2) = nsize-8; imax(nsize-2) = nsize

    dd(nsize-7:nsize,nsize-1) = -q(8:1:-1,2)
    imin(nsize-1) = nsize-7; imax(nsize-1) = nsize

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

    ir = nsize - 7
  end if
  do i = il, ir
    dd(i-3:i-1,i) = -a(3:1:-1); dd(i+1:i+3,i) = a(1:3)
    imin(i) = i-3; imax(i) = i+3
  end do
end subroutine set_coeff_6_5