00001 subroutine violation_gridpoint_MoC_CF_peos(MoC_vio)
00002 use phys_constant, only : long, pi
00003 use grid_parameter, only : nrg, ntg, npg
00004 use def_metric
00005 use def_metric_excurve_grid
00006 use def_bh_parameter, only : ome_bh
00007 use def_matter, only : emdg, omeg, vepxg, vepyg, vepzg, jomeg_int
00008 use def_matter_parameter, only : ome, ber, radi
00009 use def_vector_x
00010 use def_vector_phi
00011
00012 use make_array_3d
00013 use make_array_4d
00014 use interface_interpo_linear_type0
00015 use interface_grgrad_4th_gridpoint
00016 use interface_grgrad_gridpoint
00017 use interface_dadbscalar_3rd2nd_type0
00018
00019 implicit none
00020 real(long), pointer :: MoC_vio(:,:,:,:), ddivBeta(:,:,:,:), Lapla_beta(:,:,:,:)
00021 real(long), pointer :: fnc2(:,:,:), divBeta(:,:,:)
00022 real(long), pointer :: grad2x(:,:,:), grad2y(:,:,:), grad2z(:,:,:)
00023 real(long), pointer :: dbvxddx(:,:,:), dbvxddy(:,:,:), dbvxddz(:,:,:)
00024 real(long), pointer :: dbvyddx(:,:,:), dbvyddy(:,:,:), dbvyddz(:,:,:)
00025 real(long), pointer :: dbvzddx(:,:,:), dbvzddy(:,:,:), dbvzddz(:,:,:)
00026
00027 real(long) :: d2bxdxx, d2bxdxy, d2bxdxz
00028 real(long) :: d2bxdyx, d2bxdyy, d2bxdyz
00029 real(long) :: d2bxdzx, d2bxdzy, d2bxdzz
00030 real(long) :: d2bydxx, d2bydxy, d2bydxz
00031 real(long) :: d2bydyx, d2bydyy, d2bydyz
00032 real(long) :: d2bydzx, d2bydzy, d2bydzz
00033 real(long) :: d2bzdxx, d2bzdxy, d2bzdxz
00034 real(long) :: d2bzdyx, d2bzdyy, d2bzdyz
00035 real(long) :: d2bzdzx, d2bzdzy, d2bzdzz
00036
00037 real(long) :: vphig(3), san, san2, ovdgc(3), lam
00038 real(long) :: emdgc, rhogc, pregc, hhgc, utgc, oterm, zfac, rjj
00039 real(long) :: psigc, alpgc, fnc2gc, afnc2inv, dxafn, dyafn, dzafn
00040 real(long) :: bvxdgc, bvydgc, bvzdgc, tfkax, tfkay, tfkaz, ene
00041 real(long) :: omegc, jomeg_intgc
00042
00043 real(long) :: ovxu, ovyu, ovzu, vepxgc, vepygc, vepzgc
00044 real(long) :: dbxdx,dbxdy,dbxdz, dbydx,dbydy,dbydz, dbzdx,dbzdy,dbzdz
00045 integer :: ii, irg, itg, ipg, impt
00046 real(long) :: wxgc, wygc, wzgc, psigc4, psigcp, alpgc2, dvep2, wdvep, w2,
00047 wterm, uih2, hut
00048
00049 call alloc_array3d(fnc2 ,0,nrg,0,ntg,0,npg)
00050 call alloc_array3d(grad2x ,0,nrg,0,ntg,0,npg)
00051 call alloc_array3d(grad2y ,0,nrg,0,ntg,0,npg)
00052 call alloc_array3d(grad2z ,0,nrg,0,ntg,0,npg)
00053 call alloc_array3d(divBeta,0,nrg,0,ntg,0,npg)
00054
00055 call alloc_array3d(dbvxddx,0,nrg,0,ntg,0,npg)
00056 call alloc_array3d(dbvxddy,0,nrg,0,ntg,0,npg)
00057 call alloc_array3d(dbvxddz,0,nrg,0,ntg,0,npg)
00058
00059 call alloc_array3d(dbvyddx,0,nrg,0,ntg,0,npg)
00060 call alloc_array3d(dbvyddy,0,nrg,0,ntg,0,npg)
00061 call alloc_array3d(dbvyddz,0,nrg,0,ntg,0,npg)
00062
00063 call alloc_array3d(dbvzddx,0,nrg,0,ntg,0,npg)
00064 call alloc_array3d(dbvzddy,0,nrg,0,ntg,0,npg)
00065 call alloc_array3d(dbvzddz,0,nrg,0,ntg,0,npg)
00066
00067 call alloc_array4d(Lapla_beta,0,nrg,0,ntg,0,npg,1,3)
00068 call alloc_array4d(ddivBeta ,0,nrg,0,ntg,0,npg,1,3)
00069
00070 san = 1.0d0/3.0d0
00071 san2= 2.0d0/3.0d0
00072
00073
00074
00075
00076 do ipg = 0, npg
00077 do itg = 0, ntg
00078 do irg = 0, nrg
00079 psigc = psi(irg,itg,ipg)
00080 alpgc = alph(irg,itg,ipg)
00081 if (alpgc.le.0.0d0) alpgc = 1.0d-20
00082 fnc2(irg,itg,ipg) = psi(irg,itg,ipg)**6/alpgc
00083 call grgrad_4th_gridpoint(bvxd,dbxdx,dbxdy,dbxdz,irg,itg,ipg)
00084 call grgrad_4th_gridpoint(bvyd,dbydx,dbydy,dbydz,irg,itg,ipg)
00085 call grgrad_4th_gridpoint(bvzd,dbzdx,dbzdy,dbzdz,irg,itg,ipg)
00086
00087 divBeta(irg,itg,ipg) = dbxdx + dbydy + dbzdz
00088
00089 dbvxddx(irg,itg,ipg) = dbxdx
00090 dbvxddy(irg,itg,ipg) = dbxdy
00091 dbvxddz(irg,itg,ipg) = dbxdz
00092
00093 dbvyddx(irg,itg,ipg) = dbydx
00094 dbvyddy(irg,itg,ipg) = dbydy
00095 dbvyddz(irg,itg,ipg) = dbydz
00096
00097 dbvzddx(irg,itg,ipg) = dbzdx
00098 dbvzddy(irg,itg,ipg) = dbzdy
00099 dbvzddz(irg,itg,ipg) = dbzdz
00100 end do
00101 end do
00102 end do
00103
00104 call grgrad_gridpoint(fnc2,grad2x,grad2y,grad2z)
00105
00106 do ipg = 0, npg
00107 do itg = 0, ntg
00108 do irg = 0, nrg
00109 call grgrad_4th_gridpoint(divBeta,dxafn,dyafn,dzafn,irg,itg,ipg)
00110 ddivBeta(irg,itg,ipg,1) = dxafn
00111 ddivBeta(irg,itg,ipg,2) = dyafn
00112 ddivBeta(irg,itg,ipg,3) = dzafn
00113
00114 call grgrad_4th_gridpoint(dbvxddx,d2bxdxx,d2bxdxy,d2bxdxz,irg,itg,ipg)
00115 call grgrad_4th_gridpoint(dbvxddy,d2bxdyx,d2bxdyy,d2bxdyz,irg,itg,ipg)
00116 call grgrad_4th_gridpoint(dbvxddz,d2bxdzx,d2bxdzy,d2bxdzz,irg,itg,ipg)
00117 Lapla_beta(irg,itg,ipg,1) = d2bxdxx + d2bxdyy + d2bxdzz
00118
00119 call grgrad_4th_gridpoint(dbvyddx,d2bydxx,d2bydxy,d2bydxz,irg,itg,ipg)
00120 call grgrad_4th_gridpoint(dbvyddy,d2bydyx,d2bydyy,d2bydyz,irg,itg,ipg)
00121 call grgrad_4th_gridpoint(dbvyddz,d2bydzx,d2bydzy,d2bydzz,irg,itg,ipg)
00122 Lapla_beta(irg,itg,ipg,2) = d2bydxx + d2bydyy + d2bydzz
00123
00124 call grgrad_4th_gridpoint(dbvzddx,d2bzdxx,d2bzdxy,d2bzdxz,irg,itg,ipg)
00125 call grgrad_4th_gridpoint(dbvzddy,d2bzdyx,d2bzdyy,d2bzdyz,irg,itg,ipg)
00126 call grgrad_4th_gridpoint(dbvzddz,d2bzdzx,d2bzdzy,d2bzdzz,irg,itg,ipg)
00127 Lapla_beta(irg,itg,ipg,3) = d2bzdxx + d2bzdyy + d2bzdzz
00128 end do
00129 end do
00130 end do
00131
00132 do ii = 1, 3
00133 do ipg = 0, npg
00134 do itg = 0, ntg
00135 do irg = 0, nrg
00136 emdgc = emdg(irg,itg,ipg)
00137 omegc = omeg(irg,itg,ipg)
00138 jomeg_intgc = jomeg_int(irg,itg,ipg)
00139 psigc = psi(irg,itg,ipg)
00140 alpgc = alph(irg,itg,ipg)
00141 fnc2gc = fnc2(irg,itg,ipg)
00142 bvxdgc = bvxd(irg,itg,ipg)
00143 bvydgc = bvyd(irg,itg,ipg)
00144 bvzdgc = bvzd(irg,itg,ipg)
00145 afnc2inv = alpgc/fnc2gc
00146
00147 dxafn = grad2x(irg,itg,ipg)
00148 dyafn = grad2y(irg,itg,ipg)
00149 dzafn = grad2z(irg,itg,ipg)
00150 tfkax = tfkij_grid(irg,itg,ipg,ii,1)
00151 tfkay = tfkij_grid(irg,itg,ipg,ii,2)
00152 tfkaz = tfkij_grid(irg,itg,ipg,ii,3)
00153
00154 zfac = 1.0d0
00155 if (emdgc <= 1.0d-15) then
00156 emdgc = 1.0d-15
00157 zfac = 0.0d0
00158 end if
00159 call peos_q2hprho(emdgc, hhgc, pregc, rhogc, ene)
00160 utgc = hhgc/ber*exp(jomeg_intgc)
00161 vphig(1) = vec_phig(irg,itg,ipg,1)
00162 vphig(2) = vec_phig(irg,itg,ipg,2)
00163 vphig(3) = vec_phig(irg,itg,ipg,3)
00164 oterm = 0.0d0
00165 if (ii == 1) oterm = bvxdgc + omegc*vphig(1)
00166 if (ii == 2) oterm = bvydgc + omegc*vphig(2)
00167 if (ii == 3) oterm = bvzdgc + omegc*vphig(3)
00168
00169 rjj = hhgc*rhogc*alpgc*utgc**2*psigc**4*oterm
00170
00171 MoC_vio(irg,itg,ipg,ii) = - 2.0d0*afnc2inv &
00172 & *(tfkax*dxafn + tfkay*dyafn + tfkaz*dzafn) &
00173 & - san*ddivBeta(irg,itg,ipg,ii) &
00174 & + radi**2*16.0d0*pi*alpgc*rjj*zfac &
00175 & - Lapla_beta(irg,itg,ipg,ii)
00176 MoC_vio(irg,itg,ipg,ii) = MoC_vio(irg,itg,ipg,ii)/(2.0d0*alpgc*psigc**4)
00177
00178 end do
00179 end do
00180 end do
00181 end do
00182
00183 deallocate(fnc2)
00184 deallocate(grad2x)
00185 deallocate(grad2y)
00186 deallocate(grad2z)
00187 deallocate(divBeta)
00188
00189 deallocate(dbvxddx)
00190 deallocate(dbvxddy)
00191 deallocate(dbvxddz)
00192
00193 deallocate(dbvyddx)
00194 deallocate(dbvyddy)
00195 deallocate(dbvyddz)
00196
00197 deallocate(dbvzddx)
00198 deallocate(dbvzddy)
00199 deallocate(dbvzddz)
00200
00201 deallocate(Lapla_beta)
00202 deallocate(ddivBeta)
00203
00204 end subroutine violation_gridpoint_MoC_CF_peos