aboutsummaryrefslogtreecommitdiff
path: root/src/covderivs2.F90
blob: 15e0143360938ce84300564ca389781c572e81dd (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
! $Header$

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

module covderivs2
  implicit none
  DECLARE_CCTK_PARAMETERS
  private
  
  public calc_2scalargrad
  public calc_2covectorgrad
  public calc_2vectorgrad
  public calc_2tensorgrad
  
  public calc_2scalargradgrad
  public calc_2covectorgradgrad
  public calc_2vectorgradgrad
  
  public calc_2longitudinal
  public calc_2gradlongitudinal
  
contains
  
  subroutine calc_2scalargrad (f, df, gamma, gf)
    CCTK_REAL, intent(in)  :: f
    CCTK_REAL, intent(in)  :: df(2)
    CCTK_REAL, intent(in)  :: gamma(2,2,2)
    CCTK_REAL, intent(out) :: gf(2)
    integer :: i
    ! f;i = f,i
    do i=1,2
       gf(i) = df(i)
    end do
  end subroutine calc_2scalargrad
  
  subroutine calc_2covectorgrad (f, df, gamma, gf)
    CCTK_REAL, intent(in)  :: f(2)
    CCTK_REAL, intent(in)  :: df(2,2)
    CCTK_REAL, intent(in)  :: gamma(2,2,2)
    CCTK_REAL, intent(out) :: gf(2,2)
    integer :: i,j,k
    ! f_i;j = f_i,j - Gamma^k_ij f_k
    do i=1,2
       do j=1,2
          gf(i,j) = df(i,j)
          do k=1,2
             gf(i,j) = gf(i,j) - gamma(k,i,j) * f(k)
          end do
       end do
    end do
  end subroutine calc_2covectorgrad
  
  subroutine calc_2vectorgrad (f, df, gamma, gf)
    CCTK_REAL, intent(in)  :: f(2)
    CCTK_REAL, intent(in)  :: df(2,2)
    CCTK_REAL, intent(in)  :: gamma(2,2,2)
    CCTK_REAL, intent(out) :: gf(2,2)
    integer :: i,j,k
    ! f^i;j = f^i,j + Gamma^i_kj f^k
    do i=1,2
       do j=1,2
          gf(i,j) = df(i,j)
          do k=1,2
             gf(i,j) = gf(i,j) + gamma(i,k,j) * f(k)
          end do
       end do
    end do
  end subroutine calc_2vectorgrad
  
  subroutine calc_2tensorgrad (f, df, gamma, gf)
    CCTK_REAL, intent(in)  :: f(2,2)
    CCTK_REAL, intent(in)  :: df(2,2,2)
    CCTK_REAL, intent(in)  :: gamma(2,2,2)
    CCTK_REAL, intent(out) :: gf(2,2,2)
    integer :: i,j,k,l
    ! f_ij;k = f_ij,k - Gamma^l_ik f_lj - Gamma^l_jk f_il
    do i=1,2
       do j=1,2
          do k=1,2
             gf(i,j,k) = df(i,j,k)
             do l=1,2
                gf(i,j,k) = gf(i,j,k) - gamma(l,i,k) * f(l,j) &
                     &                - gamma(l,j,k) * f(i,l)
             end do
          end do
       end do
    end do
  end subroutine calc_2tensorgrad
  
  
  
  subroutine calc_2scalargradgrad (f, df, ddf, gamma, ggf)
    CCTK_REAL, intent(in)  :: f
    CCTK_REAL, intent(in)  :: df(2)
    CCTK_REAL, intent(in)  :: ddf(2,2)
    CCTK_REAL, intent(in)  :: gamma(2,2,2)
    CCTK_REAL, intent(out) :: ggf(2,2)
    integer :: i,j,k
    ! f;ij = f,ij - Gamma^k_ij f,k
    do i=1,2
       do j=1,2
          ggf(i,j) = ddf(i,j)
          do k=1,2
             ggf(i,j) = ggf(i,j) - gamma(k,i,j) * df(k)
          end do
       end do
    end do
  end subroutine calc_2scalargradgrad
  
  subroutine calc_2covectorgradgrad (f, df, ddf, gamma, dgamma, ggf)
    CCTK_REAL, intent(in)  :: f(2)
    CCTK_REAL, intent(in)  :: df(2,2)
    CCTK_REAL, intent(in)  :: ddf(2,2,2)
    CCTK_REAL, intent(in)  :: gamma(2,2,2), dgamma(2,2,2,2)
    CCTK_REAL, intent(out) :: ggf(2,2,2)
    CCTK_REAL :: gf(2,2)
    integer :: i,j,k,l
    ! f_i;j  = f_i,j - Gamma^l_ij f_l
    ! f_i;jk = f_i,jk - Gamma^l_ij,k f_l - Gamma^l_ij f_l,k
    !          - Gamma^l_ik f_l;j - Gamma^l_jk f_i;l
    call calc_2covectorgrad (f, df, gamma, gf)
    do i=1,2
       do j=1,2
          do k=1,2
             ggf(i,j,k) = ddf(i,j,k)
             do l=1,2
                ggf(i,j,k) = ggf(i,j,k) &
                     - dgamma(l,i,j,k) * f(l) - gamma(l,i,j) * df(l,k) &
                     - gamma(l,i,k) * gf(l,j) - gamma(l,j,k) * gf(i,l)
             end do
          end do
       end do
    end do
  end subroutine calc_2covectorgradgrad
  
  subroutine calc_2vectorgradgrad (f, df, ddf, gamma, dgamma, ggf)
    CCTK_REAL, intent(in)  :: f(2)
    CCTK_REAL, intent(in)  :: df(2,2)
    CCTK_REAL, intent(in)  :: ddf(2,2,2)
    CCTK_REAL, intent(in)  :: gamma(2,2,2), dgamma(2,2,2,2)
    CCTK_REAL, intent(out) :: ggf(2,2,2)
    CCTK_REAL :: gf(2,2)
    integer :: i,j,k,l
    ! f^i;j  = f^i,j + Gamma^i_lj f^l
    ! f^i;jk = f^i,jk + Gamma^i_lj,k f^l + Gamma^i_lj f^l,k
    !          - Gamma^i_lk f^l;j - Gamma^l_jk f^i;l
    call calc_2vectorgrad (f, df, gamma, gf)
    do i=1,2
       do j=1,2
          do k=1,2
             ggf(i,j,k) = ddf(i,j,k)
             do l=1,2
                ggf(i,j,k) = ggf(i,j,k) &
                     + dgamma(i,l,j,k) * f(l) + gamma(i,l,j) * df(l,k) &
                     - gamma(i,l,k) * gf(l,j) - gamma(l,j,k) * gf(i,l)
             end do
          end do
       end do
    end do
  end subroutine calc_2vectorgradgrad
  
  
  
  subroutine calc_2longitudinal (gf, gg, gu, lf)
    CCTK_REAL, intent(in)  :: gf(2,2)
    CCTK_REAL, intent(in)  :: gg(2,2), gu(2,2)
    CCTK_REAL, intent(out) :: lf(2,2)
    integer :: i,j,k,l
    ! (Lf)_ij = D_i f_j + D_j f_i - 2/3 g_ij D^l f_l
    do i=1,2
       do j=1,2
          lf(i,j) = gf(i,j) + gf(j,i)
          do k=1,2
             do l=1,2
                lf(i,j) = lf(i,j) - (2.d0/3.d0) * gg(i,j) * gu(k,l) * gf(k,l)
             end do
          end do
       end do
    end do
  end subroutine calc_2longitudinal
  
  subroutine calc_2gradlongitudinal (ggf, gg, gu, glf)
    CCTK_REAL, intent(in)  :: ggf(2,2,2)
    CCTK_REAL, intent(in)  :: gg(2,2), gu(2,2)
    CCTK_REAL, intent(out) :: glf(2,2,2)
    integer :: i,j,k,l,m
    ! D_k (Lf)_ij = D_k (D_i f_j + D_j f_i - 2/3 g_ij D^l f_l)
    do i=1,2
       do j=1,2
          do k=1,2
             glf(i,j,k) = ggf(j,i,k) + ggf(i,j,k)
             do l=1,2
                do m=1,2
                   glf(i,j,k) = glf(i,j,k) &
                        - (2.d0/3.d0) * gg(i,j) * gu(l,m) * ggf(l,m,k)
                end do
             end do
          end do
       end do
    end do
  end subroutine calc_2gradlongitudinal
  
end module covderivs2