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