00001 subroutine sourceterm_MoC_CF_with_divshift_peos(souvec)
00002 use phys_constant, only : long, pi
00003 use grid_parameter, only : nrg, ntg, npg
00004 use def_metric
00005 use def_bh_parameter, only : ome_bh
00006 use def_matter, only : emdg, omeg, jomeg_int
00007 use def_matter_parameter, only : ome, ber, radi
00008 use def_vector_x
00009 use def_vector_phi
00010 use make_array_3d
00011 use make_array_4d
00012 use interface_interpo_linear_type0
00013 use interface_grgrad_4th_gridpoint
00014 use interface_grgrad_midpoint_r3rd_nsbh
00015 use interface_grgrad_midpoint_r3rd_type0_ns
00016 use interface_dadbscalar_3rd2nd_type0
00017
00018 implicit none
00019 real(long), pointer :: souvec(:,:,:,:), ddivBeta(:,:,:,:)
00020 real(long), pointer :: fnc2(:,:,:), divBeta(:,:,:), gradtrK(:,:,:,:)
00021 real(long), pointer :: grad2x(:,:,:), grad2y(:,:,:), grad2z(:,:,:)
00022 real(long) :: vphig(3), san, san2, san4
00023 real(long) :: emdgc, rhogc, pregc, hhgc, utgc, oterm, zfac, rjj
00024 real(long) :: psigc, alpgc, fnc2gc, afnc2inv, dxafn, dyafn, dzafn
00025 real(long) :: bvxdgc, bvydgc, bvzdgc, tfkax, tfkay, tfkaz, dtrK, ene
00026 real(long) :: omegc, jomeg_intgc
00027
00028 real(long) :: ovxu, ovyu, ovzu
00029 real(long) :: dbxdx,dbxdy,dbxdz, dbydx,dbydy,dbydz, dbzdx,dbzdy,dbzdz
00030 real(long) :: d2bvx(3,3), d2bvy(3,3), d2bvz(3,3)
00031 integer :: ii, irg, itg, ipg, impt
00032
00033 call alloc_array3d(fnc2,0,nrg,0,ntg,0,npg)
00034 call alloc_array3d(grad2x,1,nrg,1,ntg,1,npg)
00035 call alloc_array3d(grad2y,1,nrg,1,ntg,1,npg)
00036 call alloc_array3d(grad2z,1,nrg,1,ntg,1,npg)
00037 call alloc_array3d(divBeta,0,nrg,0,ntg,0,npg)
00038 call alloc_array4d(ddivBeta,1,nrg,1,ntg,1,npg,1,3)
00039 call alloc_array4d(gradtrK,1,nrg,1,ntg,1,npg,1,3)
00040
00041 san = 1.0d0/3.0d0
00042 san2= 2.0d0/3.0d0
00043 san4= 4.0d0/3.0d0
00044
00045
00046
00047
00048 call grgrad_midpoint(trk_grid,grad2x,grad2y,grad2z)
00049 gradtrK(1:nrg,1:ntg,1:npg,1) = grad2x(1:nrg,1:ntg,1:npg)
00050 gradtrK(1:nrg,1:ntg,1:npg,2) = grad2y(1:nrg,1:ntg,1:npg)
00051 gradtrK(1:nrg,1:ntg,1:npg,3) = grad2z(1:nrg,1:ntg,1:npg)
00052
00053 do ipg = 0, npg
00054 do itg = 0, ntg
00055 do irg = 0, nrg
00056
00057
00058 psigc = psi(irg,itg,ipg)
00059 alpgc = alph(irg,itg,ipg)
00060 if (alpgc.le.0.0d0) alpgc = 1.0d-20
00061 fnc2(irg,itg,ipg) = psi(irg,itg,ipg)**6/alpgc
00062
00063
00064
00065 call grgrad_4th_gridpoint(bvxd,dbxdx,dbxdy,dbxdz,irg,itg,ipg)
00066 call grgrad_4th_gridpoint(bvyd,dbydx,dbydy,dbydz,irg,itg,ipg)
00067 call grgrad_4th_gridpoint(bvzd,dbzdx,dbzdy,dbzdz,irg,itg,ipg)
00068 divBeta(irg,itg,ipg) = dbxdx + dbydy + dbzdz
00069
00070
00071
00072
00073
00074
00075 end do
00076 end do
00077 end do
00078
00079 call grgrad_midpoint_r3rd_nsbh(fnc2,grad2x,grad2y,grad2z,'ns')
00080
00081 do ipg = 1, npg
00082 do itg = 1, ntg
00083 do irg = 1, nrg
00084
00085 call grgrad_midpoint_r3rd_type0_ns(divBeta,dxafn,dyafn,dzafn,irg,itg,ipg)
00086 ddivBeta(irg,itg,ipg,1) = dxafn
00087 ddivBeta(irg,itg,ipg,2) = dyafn
00088 ddivBeta(irg,itg,ipg,3) = dzafn
00089
00090
00091
00092
00093
00094
00095 end do
00096 end do
00097 end do
00098
00099 do ii = 1, 3
00100 do ipg = 1, npg
00101 do itg = 1, ntg
00102 do irg = 1, nrg
00103 call interpo_linear_type0(emdgc,emdg,irg,itg,ipg)
00104 call interpo_linear_type0(omegc,omeg,irg,itg,ipg)
00105 call interpo_linear_type0(jomeg_intgc,jomeg_int,irg,itg,ipg)
00106 call interpo_linear_type0(psigc,psi,irg,itg,ipg)
00107 call interpo_linear_type0(alpgc,alph,irg,itg,ipg)
00108 call interpo_linear_type0(fnc2gc,fnc2,irg,itg,ipg)
00109 call interpo_linear_type0(bvxdgc,bvxd,irg,itg,ipg)
00110 call interpo_linear_type0(bvydgc,bvyd,irg,itg,ipg)
00111 call interpo_linear_type0(bvzdgc,bvzd,irg,itg,ipg)
00112 afnc2inv = alpgc/fnc2gc
00113
00114 dxafn = grad2x(irg,itg,ipg)
00115 dyafn = grad2y(irg,itg,ipg)
00116 dzafn = grad2z(irg,itg,ipg)
00117 tfkax = tfkij(irg,itg,ipg,ii,1)
00118 tfkay = tfkij(irg,itg,ipg,ii,2)
00119 tfkaz = tfkij(irg,itg,ipg,ii,3)
00120 dtrK = gradtrK(irg,itg,ipg,ii)
00121
00122 zfac = 1.0d0
00123 if (emdgc <= 1.0d-15) then
00124 emdgc = 1.0d-15
00125 zfac = 0.0d0
00126 end if
00127 call peos_q2hprho(emdgc, hhgc, pregc, rhogc, ene)
00128 utgc = hhgc/ber*exp(jomeg_intgc)
00129 vphig(1) = hvec_phig(irg,itg,ipg,1)
00130 vphig(2) = hvec_phig(irg,itg,ipg,2)
00131 vphig(3) = hvec_phig(irg,itg,ipg,3)
00132 oterm = 0.0d0
00133 if (ii == 1) oterm = bvxdgc + ome*vphig(1)
00134 if (ii == 2) oterm = bvydgc + ome*vphig(2)
00135 if (ii == 3) oterm = bvzdgc + ome*vphig(3)
00136
00137
00138
00139
00140 rjj = hhgc*rhogc*alpgc*utgc**2*psigc**4*oterm
00141
00142 souvec(irg,itg,ipg,ii) = - 2.0d0*afnc2inv &
00143 & *(tfkax*dxafn + tfkay*dyafn + tfkaz*dzafn) &
00144 & + san4*alpgc*dtrK &
00145 & - san*ddivBeta(irg,itg,ipg,ii) &
00146 & + radi**2*16.0d0*pi*alpgc*rjj*zfac
00147
00148 end do
00149 end do
00150 end do
00151 end do
00152
00153 deallocate(fnc2)
00154 deallocate(grad2x)
00155 deallocate(grad2y)
00156 deallocate(grad2z)
00157 deallocate(divBeta)
00158 deallocate(ddivBeta)
00159 deallocate(gradtrK)
00160
00161 end subroutine sourceterm_MoC_CF_with_divshift_peos