00001 subroutine sourceterm_MWspatial_WL(souvec)
00002 
00003   use phys_constant, only : pi
00004   use def_matter_parameter, only : ome
00005   use def_metric, only : psi, alph, tfkij
00006   use def_emfield, only : vaxd, vayd, vazd, vaxu, vayu, vazu
00007   use grid_parameter, only : nrg, ntg, npg
00008   use coordinate_grav_r, only : rg
00009   use def_cristoffel_grid, only : cri_grid
00010   use def_gamma_crist, only : gmcrix, gmcriy, gmcriz
00011   use def_ricci_tensor, only : rab
00012   use def_faraday_tensor, only : fxd, fyd, fzd, fijd
00013   use def_emfield_derivatives, only : cdvaxd, cdvayd, cdvazd
00014   use def_formulation, only : swflu
00015   use def_cristoffel, only : cri
00016   use def_cutsw, only : cutfac
00017   use def_metric_hij, only : hxxu, hxyu, hxzu, hyyu, hyzu, hzzu
00018   use def_metric_rotshift, only : ovxd, ovyd, ovzd
00019   use def_formulation, only : chgra
00020   use interface_interpo_linear_type0
00021   use interface_grgrad_midpoint
00022   use interface_grgrad1g_midpoint
00023   use interface_dadbscalar_type0
00024   use make_array_3d
00025   implicit none
00026 
00027   real(8), pointer :: souvec(:,:,:,:)
00028   real(8), pointer :: grad2x(:,:,:), grad2y(:,:,:), grad2z(:,:,:)
00029   real(8), pointer :: fnc2(:,:,:)
00030   real(8) :: grad(1:3), dabfnc(3,3), Fdps2al(3,3), AabFc(3,3), 
00031             fijdgc(3,3), tfkijgc(3,3), hiju(3,3), gmiju(3,3)
00032 
00033   real(8), pointer :: cri1va(:,:,:), cri2va(:,:,:), cri3va(:,:,:), 
00034                      cri4va(:,:,:), cri5va(:,:,:), cri6va(:,:,:)
00035   real(8) :: gmxxu, gmxyu, gmxzu, gmyyu, gmyzu, gmzzu, 
00036             gmyxu, gmzxu, gmzyu, &
00037             hhxxu, hhxyu, hhxzu, hhyxu, hhyyu, hhyzu, &
00038             hhzxu, hhzyu, hhzzu, &
00039             vaxdc, vaxgc, vaydc, vaygc, vazdc, vazgc, &
00040             c11, c12, c13, c14, c15, c16, &
00041             c21, c22, c23, c24, c25, c26, &
00042             c31, c32, c33, c34, c35, c36, &
00043             gc1, gc2, gc3, gclp, gcrdb1, gcrdb2, gdcb, &
00044             gmxvax, gmxvay, gmxvaz, gmyvax, gmyvay, gmyvaz, &
00045             gmzvax, gmzvay, gmzvaz, &
00046             hddb, rax, ray, raz, &
00047             san, san2, zfac, &
00048             cutoff, p2al, p2alinv, divlp, &
00049             dxc1b, dxc2b, dxc3b, dxc4b, dxc5b, dxc6b, &
00050             dyc1b, dyc2b, dyc3b, dyc4b, dyc5b, dyc6b, &
00051             dzc1b, dzc2b, dzc3b, dzc4b, dzc5b, dzc6b, &
00052             dvaxdx, dvaxdy, dvaxdz, dvaydx, dvaydy, dvaydz, &
00053             dvazdx, dvazdy, dvazdz, &
00054             dxafn, dyafn, dzafn, &
00055             vaxdgc, vaydgc, vazdgc, fxdgc, fydgc, fzdgc
00056   real(8) :: psigc, ricciva, Fap2dp2a, p4AF
00057   integer :: ipg, irg, itg, ic1, ic2, ic3, ii, ic0(1:3)
00058 
00059   call alloc_array3d(fnc2,0,nrg,0,ntg,0,npg)
00060   call alloc_array3d(grad2x,0,nrg,0,ntg,0,npg)
00061   call alloc_array3d(grad2y,0,nrg,0,ntg,0,npg)
00062   call alloc_array3d(grad2z,0,nrg,0,ntg,0,npg)
00063   call alloc_array3d(cri1va,0,nrg,0,ntg,0,npg)
00064   call alloc_array3d(cri2va,0,nrg,0,ntg,0,npg)
00065   call alloc_array3d(cri3va,0,nrg,0,ntg,0,npg)
00066   call alloc_array3d(cri4va,0,nrg,0,ntg,0,npg)
00067   call alloc_array3d(cri5va,0,nrg,0,ntg,0,npg)
00068   call alloc_array3d(cri6va,0,nrg,0,ntg,0,npg)
00069 
00070 
00071 
00072   do ipg = 0, npg
00073     do itg = 0, ntg
00074       do irg = 0, nrg
00075 
00076         vaxdc=vaxd(irg,itg,ipg)
00077         vaydc=vayd(irg,itg,ipg)
00078         vazdc=vazd(irg,itg,ipg)
00079         c11 = cri_grid(irg,itg,ipg,1,1)
00080         c12 = cri_grid(irg,itg,ipg,1,2)
00081         c13 = cri_grid(irg,itg,ipg,1,3)
00082         c14 = cri_grid(irg,itg,ipg,1,4)
00083         c15 = cri_grid(irg,itg,ipg,1,5)
00084         c16 = cri_grid(irg,itg,ipg,1,6)
00085         c21 = cri_grid(irg,itg,ipg,2,1)
00086         c22 = cri_grid(irg,itg,ipg,2,2)
00087         c23 = cri_grid(irg,itg,ipg,2,3)
00088         c24 = cri_grid(irg,itg,ipg,2,4)
00089         c25 = cri_grid(irg,itg,ipg,2,5)
00090         c26 = cri_grid(irg,itg,ipg,2,6)
00091         c31 = cri_grid(irg,itg,ipg,3,1)
00092         c32 = cri_grid(irg,itg,ipg,3,2)
00093         c33 = cri_grid(irg,itg,ipg,3,3)
00094         c34 = cri_grid(irg,itg,ipg,3,4)
00095         c35 = cri_grid(irg,itg,ipg,3,5)
00096         c36 = cri_grid(irg,itg,ipg,3,6)
00097         cri1va(irg,itg,ipg) = c11*vaxdc + c21*vaydc + c31*vazdc
00098         cri2va(irg,itg,ipg) = c12*vaxdc + c22*vaydc + c32*vazdc
00099         cri3va(irg,itg,ipg) = c13*vaxdc + c23*vaydc + c33*vazdc
00100         cri4va(irg,itg,ipg) = c14*vaxdc + c24*vaydc + c34*vazdc
00101         cri5va(irg,itg,ipg) = c15*vaxdc + c25*vaydc + c35*vazdc
00102         cri6va(irg,itg,ipg) = c16*vaxdc + c26*vaydc + c36*vazdc
00103 
00104       end do
00105     end do
00106   end do
00107 
00108   fnc2(0:nrg,0:ntg,0:npg) = psi(0:nrg,0:ntg,0:npg)**2/alph(0:nrg,0:ntg,0:npg)
00109   call grgrad_midpoint(fnc2,grad2x,grad2y,grad2z)
00110 
00111   san = 1.0d0/3.0d0
00112   san2= 2.0d0/3.0d0
00113 
00114   do ii = 1, 3
00115 
00116     do ipg = 1, npg
00117       do itg = 1, ntg
00118         do irg = 1, nrg
00119 
00120           cutoff = 1.0d0
00121           if (chgra == 'c'.or.chgra == 'C'.or.chgra == 'W') then
00122             if (rg(irg) > cutfac*pi/ome) cutoff = 0.0d0
00123           end if
00124 
00125           call interpo_linear_type0(vaxgc,vaxu,irg,itg,ipg)
00126           call interpo_linear_type0(vaygc,vayu,irg,itg,ipg)
00127           call interpo_linear_type0(vazgc,vazu,irg,itg,ipg)
00128           call interpo_linear_type0(psigc,psi,irg,itg,ipg)
00129           call interpo_linear_type0(p2al,fnc2,irg,itg,ipg)
00130           p2alinv = 1.0d0/p2al
00131 
00132           call interpo_linear_type0(hhxxu,hxxu,irg,itg,ipg)
00133           call interpo_linear_type0(hhxyu,hxyu,irg,itg,ipg)
00134           call interpo_linear_type0(hhxzu,hxzu,irg,itg,ipg)
00135           call interpo_linear_type0(hhyyu,hyyu,irg,itg,ipg)
00136           call interpo_linear_type0(hhyzu,hyzu,irg,itg,ipg)
00137           call interpo_linear_type0(hhzzu,hzzu,irg,itg,ipg)
00138           hhyxu = hhxyu
00139           hhzxu = hhxzu
00140           hhzyu = hhyzu
00141           gmxxu = hhxxu + 1.0d0
00142           gmxyu = hhxyu
00143           gmxzu = hhxzu
00144           gmyyu = hhyyu + 1.0d0
00145           gmyzu = hhyzu
00146           gmzzu = hhzzu + 1.0d0
00147           gmyxu = gmxyu
00148           gmzxu = gmxzu
00149           gmzyu = gmyzu
00150           hiju(1,1) = hhxxu ; hiju(1,2) = hhxyu ; hiju(1,3) = hhxzu
00151           hiju(2,1) = hhxyu ; hiju(2,2) = hhyyu ; hiju(2,3) = hhyzu
00152           hiju(3,1) = hhxzu ; hiju(3,2) = hhyzu ; hiju(3,3) = hhzzu
00153           gmiju(1,1) = gmxxu; gmiju(1,2) = gmxyu; gmiju(1,3) = gmxzu
00154           gmiju(2,1) = gmxyu; gmiju(2,2) = gmyyu; gmiju(2,3) = gmyzu
00155           gmiju(3,1) = gmxzu; gmiju(3,2) = gmyzu; gmiju(3,3) = gmzzu
00156 
00157           dxafn = grad2x(irg,itg,ipg)
00158           dyafn = grad2y(irg,itg,ipg)
00159           dzafn = grad2z(irg,itg,ipg)
00160           tfkijgc(1,1) = tfkij(irg,itg,ipg,1,1)
00161           tfkijgc(1,2) = tfkij(irg,itg,ipg,1,2)
00162           tfkijgc(1,3) = tfkij(irg,itg,ipg,1,3)
00163           tfkijgc(2,1) = tfkij(irg,itg,ipg,2,1)
00164           tfkijgc(2,2) = tfkij(irg,itg,ipg,2,2)
00165           tfkijgc(2,3) = tfkij(irg,itg,ipg,2,3)
00166           tfkijgc(3,1) = tfkij(irg,itg,ipg,3,1)
00167           tfkijgc(3,2) = tfkij(irg,itg,ipg,3,2)
00168           tfkijgc(3,3) = tfkij(irg,itg,ipg,3,3)
00169 
00170           dvaxdx = cdvaxd(irg,itg,ipg,1)
00171           dvaydx = cdvayd(irg,itg,ipg,1)
00172           dvazdx = cdvazd(irg,itg,ipg,1)
00173           dvaxdy = cdvaxd(irg,itg,ipg,2)
00174           dvaydy = cdvayd(irg,itg,ipg,2)
00175           dvazdy = cdvazd(irg,itg,ipg,2)
00176           dvaxdz = cdvaxd(irg,itg,ipg,3)
00177           dvaydz = cdvayd(irg,itg,ipg,3)
00178           dvazdz = cdvazd(irg,itg,ipg,3)
00179 
00180           fxdgc = fxd(irg,itg,ipg)  
00181           fydgc = fyd(irg,itg,ipg)  
00182           fzdgc = fzd(irg,itg,ipg)  
00183           fijdgc(1,1) = 0.0d0 ; fijdgc(2,2) = 0.0d0 ; fijdgc(3,3) = 0.0d0
00184           fijdgc(1,2) = fijd(irg,itg,ipg,1) ; fijdgc(2,1) = - fijdgc(1,2)
00185           fijdgc(1,3) = fijd(irg,itg,ipg,2) ; fijdgc(3,1) = - fijdgc(1,3)
00186           fijdgc(2,3) = fijd(irg,itg,ipg,3) ; fijdgc(3,2) = - fijdgc(2,3)
00187 
00188           hddb = 0.0d0
00189           gdcb = 0.0d0
00190           gcrdb1 = 0.0d0
00191           gcrdb2 = 0.0d0
00192           rax = 0.0d0
00193           ray = 0.0d0
00194           raz = 0.0d0
00195 
00196           c11 = cri(irg,itg,ipg,1,1) ; c12 = cri(irg,itg,ipg,1,2)
00197           c13 = cri(irg,itg,ipg,1,3) ; c14 = cri(irg,itg,ipg,1,4)
00198           c15 = cri(irg,itg,ipg,1,5) ; c16 = cri(irg,itg,ipg,1,6)
00199           c21 = cri(irg,itg,ipg,2,1) ; c22 = cri(irg,itg,ipg,2,2)
00200           c23 = cri(irg,itg,ipg,2,3) ; c24 = cri(irg,itg,ipg,2,4)
00201           c25 = cri(irg,itg,ipg,2,5) ; c26 = cri(irg,itg,ipg,2,6)
00202           c31 = cri(irg,itg,ipg,3,1) ; c32 = cri(irg,itg,ipg,3,2)
00203           c33 = cri(irg,itg,ipg,3,3) ; c34 = cri(irg,itg,ipg,3,4)
00204           c35 = cri(irg,itg,ipg,3,5) ; c36 = cri(irg,itg,ipg,3,6)
00205 
00206           gc1 = gmcrix(irg,itg,ipg)
00207           gc2 = gmcriy(irg,itg,ipg)
00208           gc3 = gmcriz(irg,itg,ipg)
00209 
00210           ic1 = 1*(ii-2)*(ii-3)/2 - 2*(ii-1)*(ii-3) + 3*(ii-1)*(ii-2)/2
00211           ic2 = 2*(ii-2)*(ii-3)/2 - 4*(ii-1)*(ii-3) + 5*(ii-1)*(ii-2)/2
00212           ic3 = 3*(ii-2)*(ii-3)/2 - 5*(ii-1)*(ii-3) + 6*(ii-1)*(ii-2)/2
00213 
00214           rax = rab(irg,itg,ipg,ic1)
00215           ray = rab(irg,itg,ipg,ic2)
00216           raz = rab(irg,itg,ipg,ic3)
00217 
00218           ic0(1) = ic1
00219           ic0(2) = ic2
00220           ic0(3) = ic3
00221 
00222 
00223 
00224           if(ii == 1) call dadbscalar_type0(vaxd,dabfnc,irg,itg,ipg)
00225           if(ii == 2) call dadbscalar_type0(vayd,dabfnc,irg,itg,ipg)
00226           if(ii == 3) call dadbscalar_type0(vazd,dabfnc,irg,itg,ipg)
00227           hddb = hhxxu*dabfnc(1,1) + hhxyu*dabfnc(1,2) + hhxzu*dabfnc(1,3) &
00228           &    + hhyxu*dabfnc(2,1) + hhyyu*dabfnc(2,2) + hhyzu*dabfnc(2,3) &
00229           &    + hhzxu*dabfnc(3,1) + hhzyu*dabfnc(3,2) + hhzzu*dabfnc(3,3)
00230 
00231 
00232 
00233           if (ii == 1) then
00234             call grgrad1g_midpoint(cri1va,grad,irg,itg,ipg)
00235             dxc1b = grad(1)
00236             dyc1b = grad(2)
00237             dzc1b = grad(3)
00238             call grgrad1g_midpoint(cri2va,grad,irg,itg,ipg)
00239             dxc2b = grad(1)
00240             dyc2b = grad(2)
00241             dzc2b = grad(3)
00242             call grgrad1g_midpoint(cri3va,grad,irg,itg,ipg)
00243             dxc3b = grad(1)
00244             dyc3b = grad(2)
00245             dzc3b = grad(3)
00246             gdcb = gmxxu*dxc1b + gmxyu*dxc2b + gmxzu*dxc3b &
00247             &    + gmyxu*dyc1b + gmyyu*dyc2b + gmyzu*dyc3b &
00248             &    + gmzxu*dzc1b + gmzyu*dzc2b + gmzzu*dzc3b
00249           end if
00250           if (ii == 2) then
00251             call grgrad1g_midpoint(cri2va,grad,irg,itg,ipg)
00252             dxc2b = grad(1)
00253             dyc2b = grad(2)
00254             dzc2b = grad(3)
00255             call grgrad1g_midpoint(cri4va,grad,irg,itg,ipg)
00256             dxc4b = grad(1)
00257             dyc4b = grad(2)
00258             dzc4b = grad(3)
00259             call grgrad1g_midpoint(cri5va,grad,irg,itg,ipg)
00260             dxc5b = grad(1)
00261             dyc5b = grad(2)
00262             dzc5b = grad(3)
00263             gdcb = gmxxu*dxc2b + gmxyu*dxc4b + gmxzu*dxc5b &
00264             &    + gmyxu*dyc2b + gmyyu*dyc4b + gmyzu*dyc5b &
00265             &    + gmzxu*dzc2b + gmzyu*dzc4b + gmzzu*dzc5b
00266           end if
00267           if (ii == 3) then
00268             call grgrad1g_midpoint(cri3va,grad,irg,itg,ipg)
00269             dxc3b = grad(1)
00270             dyc3b = grad(2)
00271             dzc3b = grad(3)
00272             call grgrad1g_midpoint(cri5va,grad,irg,itg,ipg)
00273             dxc5b = grad(1)
00274             dyc5b = grad(2)
00275             dzc5b = grad(3)
00276             call grgrad1g_midpoint(cri6va,grad,irg,itg,ipg)
00277             dxc6b = grad(1)
00278             dyc6b = grad(2)
00279             dzc6b = grad(3)
00280             gdcb = gmxxu*dxc3b + gmxyu*dxc5b + gmxzu*dxc6b &
00281             &    + gmyxu*dyc3b + gmyyu*dyc5b + gmyzu*dyc6b &
00282             &    + gmzxu*dzc3b + gmzyu*dzc5b + gmzzu*dzc6b
00283           end if
00284 
00285 
00286 
00287           if (ii == 1) gcrdb1 = gc1*dvaxdx + gc2*dvaxdy + gc3*dvaxdz
00288           if (ii == 2) gcrdb1 = gc1*dvaydx + gc2*dvaydy + gc3*dvaydz
00289           if (ii == 3) gcrdb1 = gc1*dvazdx + gc2*dvazdy + gc3*dvazdz
00290 
00291           gmxvax = gmxxu*dvaxdx + gmxyu*dvaxdy + gmxzu*dvaxdz
00292           gmxvay = gmxxu*dvaydx + gmxyu*dvaydy + gmxzu*dvaydz
00293           gmxvaz = gmxxu*dvazdx + gmxyu*dvazdy + gmxzu*dvazdz
00294           gmyvax = gmyxu*dvaxdx + gmyyu*dvaxdy + gmyzu*dvaxdz
00295           gmyvay = gmyxu*dvaydx + gmyyu*dvaydy + gmyzu*dvaydz
00296           gmyvaz = gmyxu*dvazdx + gmyyu*dvazdy + gmyzu*dvazdz
00297           gmzvax = gmzxu*dvaxdx + gmzyu*dvaxdy + gmzzu*dvaxdz
00298           gmzvay = gmzxu*dvaydx + gmzyu*dvaydy + gmzzu*dvaydz
00299           gmzvaz = gmzxu*dvazdx + gmzyu*dvazdy + gmzzu*dvazdz
00300 
00301           if (ii == 1) gcrdb2 = c11*gmxvax + c21*gmxvay + c31*gmxvaz &
00302           &                   + c12*gmyvax + c22*gmyvay + c32*gmyvaz &
00303           &                   + c13*gmzvax + c23*gmzvay + c33*gmzvaz
00304 
00305           if (ii == 2) gcrdb2 = c12*gmxvax + c22*gmxvay + c32*gmxvaz &
00306           &                   + c14*gmyvax + c24*gmyvay + c34*gmyvaz &
00307           &                   + c15*gmzvax + c25*gmzvay + c35*gmzvaz
00308 
00309           if (ii == 3) gcrdb2 = c13*gmxvax + c23*gmxvay + c33*gmxvaz &
00310           &                   + c15*gmyvax + c25*gmyvay + c35*gmyvaz &
00311           &                   + c16*gmzvax + c26*gmzvay + c36*gmzvaz
00312 
00313 
00314           ricciva = rax*vaxgc + ray*vaygc + raz*vazgc
00315 
00316           Fdps2al(1,1) = fijdgc(ii,1)*dxafn
00317           Fdps2al(1,2) = fijdgc(ii,1)*dyafn
00318           Fdps2al(1,3) = fijdgc(ii,1)*dzafn
00319           Fdps2al(2,1) = fijdgc(ii,2)*dxafn
00320           Fdps2al(2,2) = fijdgc(ii,2)*dyafn
00321           Fdps2al(2,3) = fijdgc(ii,2)*dzafn
00322           Fdps2al(3,1) = fijdgc(ii,3)*dxafn
00323           Fdps2al(3,2) = fijdgc(ii,3)*dyafn
00324           Fdps2al(3,3) = fijdgc(ii,3)*dzafn
00325           call compute_trace(hiju,Fdps2al,Fap2dp2a)
00326           Fap2dp2a = p2alinv*Fap2dp2a
00327 
00328           AabFc(1,1) = tfkijgc(ii,1)*fxdgc
00329           AabFc(1,2) = tfkijgc(ii,1)*fydgc
00330           AabFc(1,3) = tfkijgc(ii,1)*fzdgc
00331           AabFc(2,1) = tfkijgc(ii,2)*fxdgc
00332           AabFc(2,2) = tfkijgc(ii,2)*fydgc
00333           AabFc(2,3) = tfkijgc(ii,2)*fzdgc
00334           AabFc(3,1) = tfkijgc(ii,3)*fxdgc
00335           AabFc(3,2) = tfkijgc(ii,3)*fydgc
00336           AabFc(3,3) = tfkijgc(ii,3)*fzdgc
00337           call compute_trace(hiju,AabFc,p4AF)
00338           p4AF = 2.0d0*psigc**4*p4AF
00339 
00340           souvec(irg,itg,ipg,ii) = - hddb + gdcb + gcrdb1 + gcrdb2 &
00341           &                      + ricciva + Fap2dp2a - p4AF
00342 
00343         end do
00344       end do
00345     end do
00346   end do
00347 
00348   deallocate(fnc2)
00349   deallocate(grad2x)
00350   deallocate(grad2y)
00351   deallocate(grad2z)
00352   deallocate(cri1va)
00353   deallocate(cri2va)
00354   deallocate(cri3va)
00355   deallocate(cri4va)
00356   deallocate(cri5va)
00357   deallocate(cri6va)
00358 
00359 end subroutine sourceterm_MWspatial_WL