00001 subroutine current_modify_MHD(cobj)
00002 use phys_constant, only : long, pi
00003 use grid_parameter
00004 use coordinate_grav_r, only : rg
00005 use trigonometry_grav_phi, only : sinphig, cosphig
00006 use def_metric, only : alph, psi
00007 use def_metric_hij, only : hxxu, hxyu, hxzu, hyyu, hyzu, hzzu
00008 use def_matter, only : rs
00009 use def_matter_parameter, only : radi
00010 use def_emfield, only : vaxd, vayd, vazd, jtu, jxu, jyu, jzu
00011 use def_emfield_derivatives, only : Lie_bAxd_grid, Lie_bAyd_grid, &
00012 & Lie_bAzd_grid
00013 use def_vector_phi, only : vec_phig
00014 use integrability_fnc_MHD
00015 use make_array_3d
00016 use interface_interpo_gr2fl
00017 use interface_grgrad_gridpoint_4th_type0
00018 use interface_calc_integrability_modify_At
00019 implicit none
00020 real(long), pointer :: vaphid(:,:,:), va_mod(:,:,:), alva_mod(:,:,:)
00021 real(long), pointer :: fxd_grid_mod(:,:,:), p6fxu(:,:,:)
00022 real(long), pointer :: fyd_grid_mod(:,:,:), p6fyu(:,:,:)
00023 real(long), pointer :: fzd_grid_mod(:,:,:), p6fzu(:,:,:)
00024 real(long) :: vecphi(3)
00025 real(long) :: dfxudx, dfxudy, dfxudz
00026 real(long) :: dfyudx, dfyudy, dfyudz
00027 real(long) :: dfzudx, dfzudy, dfzudz
00028 real(long) :: alps6inv, pi4invR2, jphiu, Aphi
00029 real(long) :: dalvadxgp, dalvadygp, dalvadzgp
00030 real(long) :: lie_bAx, lie_bAy, lie_bAz, ainvh
00031 integer :: irg, itg, ipg
00032 character(len=2), intent(in) :: cobj
00033
00034 call alloc_array3d(vaphid, 0, nrg, 0, ntg, 0, npg)
00035 call alloc_array3d(va_mod, 0, nrg, 0, ntg, 0, npg)
00036 call alloc_array3d(alva_mod, 0, nrg, 0, ntg, 0, npg)
00037 call alloc_array3d(fxd_grid_mod, 0, nrg, 0, ntg, 0, npg)
00038 call alloc_array3d(fyd_grid_mod, 0, nrg, 0, ntg, 0, npg)
00039 call alloc_array3d(fzd_grid_mod, 0, nrg, 0, ntg, 0, npg)
00040 call alloc_array3d(p6fxu, 0, nrg, 0, ntg, 0, npg)
00041 call alloc_array3d(p6fyu, 0, nrg, 0, ntg, 0, npg)
00042 call alloc_array3d(p6fzu, 0, nrg, 0, ntg, 0, npg)
00043
00044 vaphid(0:nrg,0:ntg,0:npg)=vayd(0:nrg,0:ntg,0:npg) &
00045 & *vec_phig(0:nrg,0:ntg,0:npg,2)
00046 call calc_integrability_modify_At(va_mod,'gr')
00047 alva_mod(0:nrg,0:ntg,0:npg) = alph(0:nrg,0:ntg,0:npg) &
00048 & *va_mod(0:nrg,0:ntg,0:npg)
00049
00050 do ipg = 0, npg
00051 do itg = 0, ntg
00052 do irg = 0, nrg
00053
00054 call grgrad_gridpoint_4th_type0(alva_mod, &
00055 & dalvadxgp,dalvadygp,dalvadzgp,irg,itg,ipg,cobj)
00056 ainvh = 1.0d0/alph(irg,itg,ipg)
00057 lie_bAx = Lie_bAxd_grid(irg,itg,ipg)
00058 lie_bAy = Lie_bAyd_grid(irg,itg,ipg)
00059 lie_bAz = Lie_bAzd_grid(irg,itg,ipg)
00060 fxd_grid_mod(irg,itg,ipg) = ainvh*(lie_bAx - dalvadxgp)
00061 fyd_grid_mod(irg,itg,ipg) = ainvh*(lie_bAy - dalvadygp)
00062 fzd_grid_mod(irg,itg,ipg) = ainvh*(lie_bAz - dalvadzgp)
00063
00064 end do
00065 end do
00066 end do
00067
00068 p6fxu(0:nrg,0:ntg,0:npg) = psi(0:nrg,0:ntg,0:npg)**10 &
00069 & *((1.0d0 + hxxu(0:nrg,0:ntg,0:npg))*fxd_grid_mod(0:nrg,0:ntg,0:npg) &
00070 & + hxyu(0:nrg,0:ntg,0:npg) *fyd_grid_mod(0:nrg,0:ntg,0:npg) &
00071 & + hxzu(0:nrg,0:ntg,0:npg) *fzd_grid_mod(0:nrg,0:ntg,0:npg))
00072 p6fyu(0:nrg,0:ntg,0:npg) = psi(0:nrg,0:ntg,0:npg)**10 &
00073 & *(hxyu(0:nrg,0:ntg,0:npg) *fxd_grid_mod(0:nrg,0:ntg,0:npg) &
00074 & + (1.0d0 + hyyu(0:nrg,0:ntg,0:npg))*fyd_grid_mod(0:nrg,0:ntg,0:npg) &
00075 & + hyzu(0:nrg,0:ntg,0:npg) *fzd_grid_mod(0:nrg,0:ntg,0:npg))
00076 p6fzu(0:nrg,0:ntg,0:npg) = psi(0:nrg,0:ntg,0:npg)**10 &
00077 & *(hxzu(0:nrg,0:ntg,0:npg) *fxd_grid_mod(0:nrg,0:ntg,0:npg) &
00078 & + hyzu(0:nrg,0:ntg,0:npg) *fyd_grid_mod(0:nrg,0:ntg,0:npg) &
00079 & + (1.0d0 + hzzu(0:nrg,0:ntg,0:npg))*fzd_grid_mod(0:nrg,0:ntg,0:npg))
00080
00081
00082 pi4invR2 = 1.0d0/(4.0d0*pi*radi**2)
00083 ipg = 0
00084 do itg = 0, ntg
00085 do irg = 0, nrg
00086
00087 vecphi(1)= vec_phig(irg,itg,ipg,1)
00088 vecphi(2)= vec_phig(irg,itg,ipg,2)
00089 vecphi(3)= vec_phig(irg,itg,ipg,3)
00090 alps6inv = 1.0d0/(alph(irg,itg,ipg)*psi(irg,itg,ipg)**6)
00091 Aphi = vaphid(irg,itg,ipg)
00092 call calc_integrability_fnc_MHD(Aphi)
00093 call grgrad_gridpoint_4th_type0(p6fxu,dfxudx,dfxudy,dfxudz,&
00094 & irg,itg,ipg,cobj)
00095 call grgrad_gridpoint_4th_type0(p6fyu,dfyudx,dfyudy,dfyudz,&
00096 & irg,itg,ipg,cobj)
00097 call grgrad_gridpoint_4th_type0(p6fzu,dfzudx,dfzudy,dfzudz,&
00098 & irg,itg,ipg,cobj)
00099
00100 jtu(irg,itg,ipg) = pi4invR2*alps6inv*(dfxudx + dfyudy +dfzudz)
00101
00102 if (irg.eq.0.or.itg.eq.0) then
00103 jphiu = 0.0d0
00104 else
00105 jphiu = jyu(irg,itg,ipg)/vecphi(2) - MHDfnc_dAt*jtu(irg,itg,ipg)
00106 end if
00107
00108 jyu(irg,itg,ipg) = jphiu*vecphi(2)
00109
00110 if (rg(irg).gt.rs(itg,ipg)) then
00111 jtu(irg,itg,ipg) = 0.0d0
00112 jyu(irg,itg,ipg) = 0.0d0
00113 end if
00114
00115 end do
00116 end do
00117
00118
00119
00120 do ipg = 0, npg
00121 do itg = 0, ntg
00122 do irg = 0, nrg
00123 jtu(irg,itg,ipg) = jtu(irg,itg,0)
00124 jxu(irg,itg,ipg) = cosphig(ipg)*jxu(irg,itg,0) &
00125 & - sinphig(ipg)*jyu(irg,itg,0)
00126 jyu(irg,itg,ipg) = sinphig(ipg)*jxu(irg,itg,0) &
00127 & + cosphig(ipg)*jyu(irg,itg,0)
00128 jzu(irg,itg,ipg) = jzu(irg,itg,0)
00129
00130 if (rg(irg).gt.rs(itg,ipg)) then
00131 jtu(irg,itg,ipg) = 0.0d0
00132 jxu(irg,itg,ipg) = 0.0d0
00133 jyu(irg,itg,ipg) = 0.0d0
00134 jzu(irg,itg,ipg) = 0.0d0
00135 end if
00136
00137 end do
00138 end do
00139 end do
00140
00141 deallocate(vaphid)
00142 deallocate(va_mod)
00143 deallocate(alva_mod)
00144 deallocate(fxd_grid_mod)
00145 deallocate(fyd_grid_mod)
00146 deallocate(fzd_grid_mod)
00147 deallocate(p6fxu)
00148 deallocate(p6fyu)
00149 deallocate(p6fzu)
00150
00151 end subroutine current_modify_MHD