00001
00002
00003 module radial_green_fn_grav_bhex_dd
00004 use phys_constant, only : long
00005 use grid_parameter, only : nrg, nlg, rgin, rgout
00006 use coordinate_grav_r, only : rg, rginv, hrg, hrginv
00007 use make_array_3d
00008 implicit none
00009 real(long), pointer :: hgfn_dd(:,:,:), gfnsf_dd(:,:,:)
00010 contains
00011 subroutine allocate_hgfn_bhex_dd
00012 implicit none
00013 call alloc_array3d(hgfn_dd,1,nrg,0,nlg,0,nrg)
00014 call alloc_array3d(gfnsf_dd,0,nlg,0,nrg,1,4)
00015 end subroutine allocate_hgfn_bhex_dd
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 subroutine calc_hgfn_bhex_dd
00027 implicit none
00028 integer :: ir, irr, il
00029 real(long) :: bhrad, bhradinv, bhsq
00030 real(long) :: radext, radextinv, radsq
00031 real(long) :: fac1, fac2, fac3, dil1, dil2
00032
00033 do ir = 0, nrg
00034 do il = 0, nlg
00035 do irr = 1, nrg
00036 hgfn_dd(irr,il,ir) = 0.d0
00037 end do
00038 end do
00039 end do
00040
00041
00042
00043 bhrad = rgin
00044 bhradinv = 1.0d0/rgin
00045 radext = rgout
00046 radextinv = 1.0d0/rgout
00047
00048 do ir = 0, nrg
00049 do irr = 1, nrg
00050 if (hrg(irr).lt.rg(ir)) then
00051 do il = 0, nlg
00052 fac1 = 1.0d0/(1.0d0 - (bhrad*radextinv)**(2*il+1))
00053 fac2 = (bhrad*radextinv)**il*radextinv
00054 hgfn_dd(irr,il,ir) = fac1*fac2 &
00055 & *((hrg(irr)*bhradinv)**il - (bhrad*hrginv(irr))**(il+1)) &
00056 & *((radext*rginv(ir))**(il+1) - (rg(ir)*radextinv)**il)
00057 end do
00058 end if
00059 if (hrg(irr).ge.rg(ir)) then
00060 do il = 0, nlg
00061 fac1 = 1.0d0/(1.0d0 - (bhrad*radextinv)**(2*il+1))
00062 fac2 = (bhrad*radextinv)**il*radextinv
00063 hgfn_dd(irr,il,ir) = fac1*fac2 &
00064 & *((rg(ir)*bhradinv)**il - (bhrad*rginv(ir))**(il+1)) &
00065 & *((radext*hrginv(irr))**(il+1) - (hrg(irr)*radextinv)**il)
00066 end do
00067 end if
00068 end do
00069 end do
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079 bhrad = rgin
00080 bhradinv = 1.0d0/rgin
00081 radext = rgout
00082 radextinv = 1.0d0/rgout
00083 bhsq = bhrad**2
00084 radsq = radext**2
00085
00086 do ir = 0, nrg
00087 il = 0
00088
00089 fac1 = 1.0d0/(1.0d0 - (bhrad*radextinv)**(2*il+1))
00090 fac2 = dble(2*il+1)*bhrad**(il-1)*radextinv**(il+1)
00091 fac3 = dble(2*il+1)*bhrad**il*radextinv**(il+2)
00092 gfnsf_dd(il,ir,1) = - 0.0d0
00093 gfnsf_dd(il,ir,2) = - fac1*fac2*(radext*rginv(ir) - 1.0d0)*bhsq
00094 gfnsf_dd(il,ir,3) = 0.0d0
00095 gfnsf_dd(il,ir,4) = - fac1*fac3*(1.0d0 - bhrad*rginv(ir))*radsq
00096
00097 do il = 1, nlg
00098
00099 fac1 = 1.0d0/(1.0d0 - (bhrad*radextinv)**(2*il+1))
00100 fac2 = dble(2*il+1)*bhrad**(il-1)*radextinv**(il+1)
00101 fac3 = dble(2*il+1)*bhrad**il*radextinv**(il+2)
00102 gfnsf_dd(il,ir,1) = - 0.0d0
00103 gfnsf_dd(il,ir,2) = - fac1*fac2 &
00104 & *((radext*rginv(ir))**(il+1)-(rg(ir)*radextinv)**il)*bhsq
00105 gfnsf_dd(il,ir,3) = 0.0d0
00106 gfnsf_dd(il,ir,4) = - fac1*fac3 &
00107 & *((rg(ir)*bhradinv)**il-(bhrad*rginv(ir))**(il+1))*radsq
00108
00109 end do
00110 end do
00111
00112 end subroutine calc_hgfn_bhex_dd
00113 end module radial_green_fn_grav_bhex_dd