00001
00002
00003 module legendre_fn_grav
00004 use phys_constant, only : long
00005 implicit none
00006 real(long), pointer :: plmg(:,:,:), hplmg(:,:,:)
00007 real(long), pointer :: dtplmg(:,:,:), hdtplmg(:,:,:)
00008 real(long), pointer :: yplmg(:,:,:), hyplmg(:,:,:)
00009 real(long), pointer :: dtyplmg(:,:,:), hdtyplmg(:,:,:)
00010 real(long), pointer :: facnmg(:,:), epsig(:)
00011 contains
00012 subroutine allocate_legendre
00013 use phys_constant, only : long
00014 use grid_parameter, only : ntg, nlg
00015 use make_array_3d
00016 use make_array_2d
00017 use make_array_1d
00018 implicit none
00019
00020
00021
00022 call alloc_array3d(plmg, 0, nlg, 0, nlg, 0, ntg)
00023 call alloc_array3d(hplmg, 0, nlg, 0, nlg, 1, ntg)
00024 call alloc_array3d(dtplmg, 0, nlg, 0, nlg, 0, ntg)
00025 call alloc_array3d(hdtplmg, 0, nlg, 0, nlg, 1, ntg)
00026 call alloc_array3d(yplmg, 0, nlg, 0, nlg, 0, ntg)
00027 call alloc_array3d(hyplmg, 0, nlg, 0, nlg, 1, ntg)
00028 call alloc_array3d(dtyplmg, 0, nlg, 0, nlg, 0, ntg)
00029 call alloc_array3d(hdtyplmg, 0, nlg, 0, nlg, 1, ntg)
00030 call alloc_array2d(facnmg, 0, nlg, 0, nlg)
00031 call alloc_array1d(epsig, 0, nlg)
00032
00033 end subroutine allocate_legendre
00034
00035 SUBROUTINE legendre
00036 use phys_constant, only : nnlg, long, pi
00037 use grid_parameter, only : ntg, nlg
00038 use coordinate_grav_theta, only : thg, hthg
00039 use make_array_3d
00040 use make_array_2d
00041 use make_array_1d
00042 implicit none
00043 integer :: it, mm, nn
00044 real(long) :: fmm, fnn, fnmfm, yplm_fac, theta
00045 real(long) :: pnag(0:nnlg,0:nnlg), dtpnag(0:nnlg,0:nnlg)
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057 do it = 0, ntg
00058 theta = thg(it)
00059 call legendre_theta(pnag,dtpnag,theta)
00060 plmg(0:nlg,0:nlg,it) = pnag(0:nlg,0:nlg)
00061 dtplmg(0:nlg,0:nlg,it) = dtpnag(0:nlg,0:nlg)
00062 end do
00063 do it = 1, ntg
00064 theta = hthg(it)
00065 call legendre_theta(pnag,dtpnag,theta)
00066 hplmg(0:nlg,0:nlg,it) = pnag(0:nlg,0:nlg)
00067 hdtplmg(0:nlg,0:nlg,it) = dtpnag(0:nlg,0:nlg)
00068 end do
00069
00070
00071
00072
00073 do nn = 0, nlg
00074 do mm = 0, nlg
00075 facnmg(nn,mm) = 0.0d+0
00076 end do
00077 end do
00078
00079 facnmg(0,0) = 1.0d0
00080 do nn = 1, nlg
00081 fnn = real(nn)
00082 facnmg(nn,0) = 1.0d0
00083 do mm = 1, nn
00084 fmm = real(mm)
00085 fnmfm= fnn-fmm + 1.0d0
00086 facnmg(nn,mm) = facnmg(nn,mm-1)/(fnn+fmm)/fnmfm
00087 end do
00088 end do
00089
00090
00091 do mm = 0, nlg
00092 epsig(mm) = 2.0d+0
00093 if (mm == 0) epsig(mm) = 1.0d+0
00094 end do
00095
00096 do nn = 0, nlg
00097 do mm = 0, nn
00098 do it = 0, ntg
00099 yplm_fac = (epsig(mm)*(2.0d0*dble(nn)+1.0d0)/(4.0d0*pi)* &
00100 & facnmg(nn,mm))**0.5d0
00101 yplmg(nn,mm,it) = yplm_fac* plmg(nn,mm,it)
00102 dtyplmg(nn,mm,it) = yplm_fac* dtplmg(nn,mm,it)
00103 if (it.eq.0) cycle
00104 hyplmg(nn,mm,it) = yplm_fac* hplmg(nn,mm,it)
00105 hdtyplmg(nn,mm,it) = yplm_fac*hdtplmg(nn,mm,it)
00106 end do
00107 end do
00108 end do
00109
00110 end subroutine legendre
00111 SUBROUTINE legendre_theta(pnag,dtpnag,theta)
00112 use phys_constant, only : nnlg, long
00113 use grid_parameter, only : nlg
00114 implicit none
00115 integer :: it, mm, nn, kk
00116 real(long) :: fmm, fkk, cc, ss, q1, q2
00117 real(long) :: pnag(0:nnlg,0:nnlg), dtpnag(0:nnlg,0:nnlg)
00118 real(long) :: theta, fac(0:nnlg)
00119
00120 cc = cos(theta)
00121 ss = sin(theta)
00122
00123 pnag(0:nlg,0:nlg) = 0.0d0
00124 dtpnag(0:nlg,0:nlg) = 0.0d0
00125
00126 fac(0) = 1.0d0
00127 do mm = 1, nlg
00128 fmm = real(mm)
00129 fac(mm) = (2.0d0*fmm-1.0d0) * fac(mm-1)
00130 end do
00131
00132 pnag(0,0) = 1.d0
00133 dtpnag(0,0) = 0.d0
00134 do mm = 1, nlg
00135 pnag(mm,mm) = fac(mm) * (-ss)**mm
00136 end do
00137
00138 mm = 1
00139 dtpnag(mm,mm) = fac(mm) * dble(mm)*(-cc)
00140 do mm = 2, nlg
00141 dtpnag(mm,mm) = fac(mm) * dble(mm)*(-ss)**(mm-1)*(-cc)
00142 end do
00143
00144 do mm = 0, nlg-1
00145 fmm = real(mm)
00146 pnag(mm+1,mm) = (2.0d0*fmm + 1.0d0)* cc* pnag(mm,mm)
00147 dtpnag(mm+1,mm) = (2.0d0*fmm + 1.0d0)*(cc*dtpnag(mm,mm) &
00148 & -ss* pnag(mm,mm))
00149 end do
00150
00151 do mm = 0, nlg-2
00152 fmm = real(mm)
00153 do kk = 2, nlg-mm
00154 fkk = real(kk)
00155 q1 = ( 2.0d0 * fmm + 2.0d0 * fkk - 1.0d0 ) / fkk
00156 q2 = ( 2.0d0 * fmm + fkk - 1.0d0 ) / fkk
00157 pnag(mm+kk,mm) = q1* cc* pnag(mm+kk-1,mm) - q2* pnag(mm+kk-2,mm)
00158 dtpnag(mm+kk,mm) = q1*(cc*dtpnag(mm+kk-1,mm) - ss * pnag(mm+kk-1,mm)) &
00159 & - q2*dtpnag(mm+kk-2,mm)
00160 end do
00161 end do
00162 end subroutine legendre_theta
00163 end module legendre_fn_grav