aboutsummaryrefslogtreecommitdiff
path: root/src/tensor2.F90
blob: 277f81365ff6f030e5c02e3fb43d366ae5218e73 (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
! $Header$

#include "cctk.h"

module tensor2
  implicit none
  private
  
  public calc_2trace
  
  public calc_2det
  public calc_2detderiv
  public calc_2detdot
  
  public calc_2inv
  public calc_2invderiv
  public calc_2invdot
  
contains
  
  subroutine calc_2trace (kk, gu, trk)
    CCTK_REAL, intent(in)  :: kk(2,2)
    CCTK_REAL, intent(in)  :: gu(2,2)
    CCTK_REAL, intent(out) :: trk
    integer :: i,j
    trk = 0
    do i=1,2
       do j=1,2
          trk = trk + gu(i,j) * kk(i,j)
       end do
    end do
  end subroutine calc_2trace
  
  
  
  subroutine calc_2det (gg, dtg)
    CCTK_REAL, intent(in)  :: gg(2,2)
    CCTK_REAL, intent(out) :: dtg
    dtg =  gg(1,1) * gg(2,2) - gg(1,2)**2
  end subroutine calc_2det
  
  subroutine calc_2pdet (gg, pgg, pdtg)
    CCTK_REAL, intent(in)  :: gg(2,2), pgg(2,2)
    CCTK_REAL, intent(out) :: pdtg
    pdtg = pgg(1,1) * gg(2,2) + gg(1,1) * pgg(2,2) - 2 * gg(1,2) * pgg(1,2)
  end subroutine calc_2pdet
  
  subroutine calc_2detderiv (gg, dgg, ddtg)
    CCTK_REAL, intent(in)  :: gg(2,2), dgg(2,2,2)
    CCTK_REAL, intent(out) :: ddtg(2)
    integer :: i
    do i=1,2
       call calc_2pdet (gg, dgg(:,:,i), ddtg(i))
    end do
  end subroutine calc_2detderiv
  
  subroutine calc_2detdot (gg, gg_dot, dtg_dot)
    CCTK_REAL, intent(in)  :: gg(2,2), gg_dot(2,2)
    CCTK_REAL, intent(out) :: dtg_dot
    call calc_2pdet (gg, gg_dot, dtg_dot)
  end subroutine calc_2detdot
  
  
  
  subroutine calc_2inv (gg, dtg, gu)
    CCTK_REAL, intent(in)  :: gg(2,2), dtg
    CCTK_REAL, intent(out) :: gu(2,2)
    gu(1,1) = gg(2,2) / dtg
    gu(2,2) = gg(1,1) / dtg
    gu(1,2) = gg(1,2) / dtg
    gu(2,1) = gu(1,2)
  end subroutine calc_2inv
  
  subroutine calc_2pinv (gu, pgg, pgu)
    CCTK_REAL, intent(in)  :: gu(2,2), pgg(2,2)
    CCTK_REAL, intent(out) :: pgu(2,2)
    integer :: i,j,k,l
    do i=1,2
       do j=1,2
          pgu(i,j) = 0
          do k=1,2
             do l=1,2
                pgu(i,j) = pgu(i,j) - gu(i,k) * gu(j,l) * pgg(k,l)
             end do
          end do
       end do
    end do
  end subroutine calc_2pinv
  
  subroutine calc_2invderiv (gu, dgg, dgu)
    CCTK_REAL, intent(in)  :: gu(2,2), dgg(2,2,2)
    CCTK_REAL, intent(out) :: dgu(2,2,2)
    integer :: i
    do i=1,2
       call calc_2pinv (gu, dgg(:,:,i), dgu(:,:,i))
    end do
  end subroutine calc_2invderiv
  
  subroutine calc_2invdot (gu, gg_dot, gu_dot)
    CCTK_REAL, intent(in)  :: gu(2,2), gg_dot(2,2)
    CCTK_REAL, intent(out) :: gu_dot(2,2)
    call calc_2pinv (gu, gg_dot, gu_dot)
  end subroutine calc_2invdot
  
end module tensor2