00001 subroutine sourceterm_trfreeG_WL(souten)
00002 use grid_parameter, only : long, nrg, ntg, npg
00003 use phys_constant, only : pi
00004 use coordinate_grav_r, only : rg, hrg
00005 use def_matter_parameter, only : ome, ber, radi
00006 use def_metric, only : psi, alph, alps, tfkij, &
00007 & bvxu, bvyu, bvzu, bvxd, bvyd, bvzd, alps2
00008 use def_metric_excurve_grid, only : tfkij_grid
00009 use def_metric_rotshift, only : ovxu, ovyu, ovzu, &
00010 & ovxd, ovyd, ovzd
00011 use def_cristoffel, only : cri
00012 use def_ricci_tensor, only : rabnl
00013 use def_shift_derivatives, only : cdbvxd, cdbvyd, cdbvzd, cdivbv
00014 use def_Lie_derivatives, only : rlpxx, rlpxy, rlpxz, rlpyy, rlpyz, rlpzz, &
00015 & rlbxx, rlbxy, rlbxz, rlbyy, rlbyz, rlbzz
00016 use def_Lie_derivatives_grid, only : rlbxx_grid, rlbxy_grid, rlbxz_grid, &
00017 & rlbyy_grid, rlbyz_grid, rlbzz_grid, &
00018 & rlpxx_grid, rlpxy_grid, rlpxz_grid, &
00019 & rlpyy_grid, rlpyz_grid, rlpzz_grid
00020 use def_formulation, only : swlp, swls, iswl
00021 use def_cutsw, only : cutfac
00022 use def_dvphi, only : dphiu
00023 use def_metric_hij, only : hxxu, hxyu, hxzu, hyyu, hyzu, hzzu, &
00024 & hxxd, hxyd, hxzd, hyyd, hyzd, hzzd
00025 use def_formulation, only : chgra, chope
00026 use make_array_3d
00027 use make_array_4d
00028 use make_array_5d
00029 use interface_interpo_linear_type0
00030 use interface_grgrad_midpoint
00031 use interface_grgrad1g_midpoint
00032 use interface_grdphi_midpoint_type0
00033 use interface_dadbscalar_type0
00034 use interface_dadbscalar_type3
00035 use interface_grd2phi_midpoint_type0
00036 implicit none
00037
00038 real(8), pointer :: souten(:,:,:,:)
00039 real(8), pointer :: dfdx(:,:,:), dfdy(:,:,:), dfdz(:,:,:)
00040 real(8), pointer :: fnc0(:,:,:)
00041 real(8), pointer :: grada(:,:,:,:), gradp(:,:,:,:)
00042 real(8), pointer :: gradap(:,:,:,:), gradap2(:,:,:,:)
00043 real(8), pointer :: gradk(:,:,:,:,:), gradb(:,:,:,:,:)
00044 real(8) :: dabalps2(1:3,1:3), daalps2(3),
00045 dapsi(3), daalph(3), daalps(3), cdbv(3,3)
00046 real(8) :: dbv(1:3,1:3), dov(1:3,1:3)=0.0d0
00047 real(8) :: daij(1:6,1:3)
00048 real(8) :: gamu(1:3,1:3), gamd(1:3,1:3), aij(1:3,1:3)
00049 real(8) :: rlieaij(3,3), rnlab(3,3), rpsiab(3,3),
00050 aacabc(3,3), sab(3,3), ovd(3), ovu(3), rpddal(3,3)
00051 real(8) :: heli(3,3), dphabd(3,3), dphabu(3,3), d2phab(3,3),
00052 rlbgab(3,3), rlpgab(3,3), dliepg(6,3), dliebg(6,3), &
00053 rlplbg(3,3), rlblpg(3,3), rlblbg(3,3), palidp(3,3)
00054 real(8) :: ica(6), icb(6)
00055 real(8) :: tm1a(3,3),tm2a(3,3),tm3a(3,3),tm4a(3,3),tm5a(3,3)
00056
00057 real(8) :: alpgc, pregc, psigc,
00058 hhxxu, hhxyu, hhxzu, hhyxu, hhyyu, hhyzu, &
00059 hhzxu, hhzyu, hhzzu, &
00060 bvxdc, bvxgc, bvydc, bvygc, bvzdc, bvzgc, &
00061 san, san2, cutoff, &
00062 alpgcinv, ap2gc, ap2gcinv, c1ab, c2ab, c3ab, cdalps2, cdivbc, &
00063 d2phxx, d2phxy, d2phxz, d2phyx, d2phyy, d2phyz, &
00064 d2phzx, d2phzy, d2phzz, &
00065 dapda, dhdh, &
00066 dhxxddx, dhxxddy, dhxxddz, dhxxudx, dhxxudy, dhxxudz, &
00067 dhxyddx, dhxyddy, dhxyddz, dhxyudx, dhxyudy, dhxyudz, &
00068 dhxzddx, dhxzddy, dhxzddz, dhxzudx, dhxzudy, dhxzudz, &
00069 dhyxddx, dhyxddy, dhyxddz, dhyxudx, dhyxudy, dhyxudz, &
00070 dhyyddx, dhyyddy, dhyyddz, dhyyudx, dhyyudy, dhyyudz, &
00071 dhyzddx, dhyzddy, dhyzddz, dhyzudx, dhyzudy, dhyzudz, &
00072 dhzxddx, dhzxddy, dhzxddz, dhzxudx, dhzxudy, dhzxudz, &
00073 dhzyddx, dhzyddy, dhzyddz, dhzyudx, dhzyudy, dhzyudz, &
00074 dhzzddx, dhzzddy, dhzzddz, dhzzudx, dhzzudy, dhzzudz, &
00075 dpdbgxx, dpdbgxy, dpdbgxz, dpdbgyx, dpdbgyy, dpdbgyz, &
00076 dpdbgzx, dpdbgzy, dpdbgzz, &
00077 dphdph, dphxxd, dphxyd, dphxzd, dphxxu, dphxyu, dphxzu, &
00078 dphyxd, dphyyd, dphyzd, dphyxu, dphyyu, dphyzu, &
00079 dphzxd, dphzyd, dphzzd, dphzxu, dphzyu, dphzzu, &
00080 gadadp, gamdhdh, hhxxd, hhxyd, hhxzd, &
00081 ogdphdph, ovpadpsi, ovxgc, ovygc, ovzgc, &
00082 ps43oal, ps4oal, ps4oal2, psigc4, psiinv, psiinv2, &
00083 rlielna, rlielnp, roku, sgam, sw1, sw2, sw3, &
00084 term1, term2, term3, term4, term5, tfaa, &
00085 tfaij, tfd2ph, tfheli, tflieaij, tfliep4aij, tfrnlab, tfrpddal, &
00086 trsab, traa, traij, trd2ph, trheli, trlieaij, trrnlab, trrpddal, &
00087 hhyxd, hhyyd, hhyzd, hhzxd, hhzyd, hhzzd, &
00088 ene, grad(1:3)
00089 integer :: ipg, irg, itg, iph, ii, ia, ib, ic
00090
00091 call alloc_array3d(dfdx,1,nrg,1,ntg,1,npg)
00092 call alloc_array3d(dfdy,1,nrg,1,ntg,1,npg)
00093 call alloc_array3d(dfdz,1,nrg,1,ntg,1,npg)
00094 call alloc_array3d(fnc0,0,nrg,0,ntg,0,npg)
00095 call alloc_array4d(grada,1,nrg,1,ntg,1,npg,1,3)
00096 call alloc_array4d(gradp,1,nrg,1,ntg,1,npg,1,3)
00097 call alloc_array4d(gradap,1,nrg,1,ntg,1,npg,1,3)
00098 call alloc_array4d(gradap2,1,nrg,1,ntg,1,npg,1,3)
00099 call alloc_array5d(gradk,1,nrg,1,ntg,1,npg,1,6,1,3)
00100 call alloc_array5d(gradb,1,nrg,1,ntg,1,npg,1,3,1,3)
00101
00102
00103
00104 do ii = 1, 3
00105
00106 if (ii == 1) call grgrad_midpoint(bvxu,dfdx,dfdy,dfdz)
00107 if (ii == 2) call grgrad_midpoint(bvyu,dfdx,dfdy,dfdz)
00108 if (ii == 3) call grgrad_midpoint(bvzu,dfdx,dfdy,dfdz)
00109 gradb(1:nrg,1:ntg,1:npg,ii,1) = dfdx(1:nrg,1:ntg,1:npg)
00110 gradb(1:nrg,1:ntg,1:npg,ii,2) = dfdy(1:nrg,1:ntg,1:npg)
00111 gradb(1:nrg,1:ntg,1:npg,ii,3) = dfdz(1:nrg,1:ntg,1:npg)
00112 end do
00113
00114 do ic = 1, 6
00115 ia = 1 + ic/4 + ic/6
00116 ib = ic - (ic/4)*2 - ic/6
00117 fnc0(0:nrg,0:ntg,0:npg) = tfkij_grid(0:nrg,0:ntg,0:npg,ia,ib)
00118 call grgrad_midpoint(fnc0,dfdx,dfdy,dfdz)
00119 gradk(1:nrg,1:ntg,1:npg,ic,1) = dfdx(1:nrg,1:ntg,1:npg)
00120 gradk(1:nrg,1:ntg,1:npg,ic,2) = dfdy(1:nrg,1:ntg,1:npg)
00121 gradk(1:nrg,1:ntg,1:npg,ic,3) = dfdz(1:nrg,1:ntg,1:npg)
00122 end do
00123
00124 alps2(0:nrg,0:ntg,0:npg) = &
00125 & alph(0:nrg,0:ntg,0:npg)*psi(0:nrg,0:ntg,0:npg)**2
00126
00127 call grgrad_midpoint(psi,dfdx,dfdy,dfdz)
00128 gradp(1:nrg,1:ntg,1:npg,1) = dfdx(1:nrg,1:ntg,1:npg)
00129 gradp(1:nrg,1:ntg,1:npg,2) = dfdy(1:nrg,1:ntg,1:npg)
00130 gradp(1:nrg,1:ntg,1:npg,3) = dfdz(1:nrg,1:ntg,1:npg)
00131 call grgrad_midpoint(alph,dfdx,dfdy,dfdz)
00132 grada(1:nrg,1:ntg,1:npg,1) = dfdx(1:nrg,1:ntg,1:npg)
00133 grada(1:nrg,1:ntg,1:npg,2) = dfdy(1:nrg,1:ntg,1:npg)
00134 grada(1:nrg,1:ntg,1:npg,3) = dfdz(1:nrg,1:ntg,1:npg)
00135 call grgrad_midpoint(alps,dfdx,dfdy,dfdz)
00136 gradap(1:nrg,1:ntg,1:npg,1) = dfdx(1:nrg,1:ntg,1:npg)
00137 gradap(1:nrg,1:ntg,1:npg,2) = dfdy(1:nrg,1:ntg,1:npg)
00138 gradap(1:nrg,1:ntg,1:npg,3) = dfdz(1:nrg,1:ntg,1:npg)
00139 call grgrad_midpoint(alps2,dfdx,dfdy,dfdz)
00140 gradap2(1:nrg,1:ntg,1:npg,1) = dfdx(1:nrg,1:ntg,1:npg)
00141 gradap2(1:nrg,1:ntg,1:npg,2) = dfdy(1:nrg,1:ntg,1:npg)
00142 gradap2(1:nrg,1:ntg,1:npg,3) = dfdz(1:nrg,1:ntg,1:npg)
00143
00144 san = 1.0d0/3.0d0
00145 san2= 2.0d0/3.0d0
00146 roku= 1.0d0/6.0d0
00147
00148 do ipg = 1, npg
00149 do itg = 1, ntg
00150 do irg = 1, nrg
00151
00152 cutoff = 1.0d0
00153 if (chgra == 'C'.or.chgra == 'H'.or.chgra == 'W') then
00154 if (rg(irg) > cutfac*pi/ome) cutoff = 0.0d0
00155 end if
00156
00157 call interpo_linear_type0(psigc,psi,irg,itg,ipg)
00158 call interpo_linear_type0(alpgc,alph,irg,itg,ipg)
00159 call interpo_linear_type0(ap2gc,alps2,irg,itg,ipg)
00160 psigc4 = psigc**4
00161 psiinv = 1.0d0/psigc
00162 psiinv2 = psiinv**2
00163 alpgcinv = 1.0d0/alpgc
00164 ap2gcinv = 1.0d0/ap2gc
00165 ps43oal = 4.0d0*psigc**3*alpgcinv
00166 ps4oal = psigc4*alpgcinv
00167 ps4oal2= psigc4*alpgcinv**2
00168
00169 call interpo_linear_type0(hhxxu,hxxu,irg,itg,ipg)
00170 call interpo_linear_type0(hhxyu,hxyu,irg,itg,ipg)
00171 call interpo_linear_type0(hhxzu,hxzu,irg,itg,ipg)
00172 call interpo_linear_type0(hhyyu,hyyu,irg,itg,ipg)
00173 call interpo_linear_type0(hhyzu,hyzu,irg,itg,ipg)
00174 call interpo_linear_type0(hhzzu,hzzu,irg,itg,ipg)
00175 hhyxu = hhxyu
00176 hhzxu = hhxzu
00177 hhzyu = hhyzu
00178 gamu(1,1) = hhxxu + 1.0d0
00179 gamu(1,2) = hhxyu
00180 gamu(1,3) = hhxzu
00181 gamu(2,2) = hhyyu + 1.0d0
00182 gamu(2,3) = hhyzu
00183 gamu(3,3) = hhzzu + 1.0d0
00184 gamu(2,1) = gamu(1,2)
00185 gamu(3,1) = gamu(1,3)
00186 gamu(3,2) = gamu(2,3)
00187
00188 call interpo_linear_type0(hhxxd,hxxd,irg,itg,ipg)
00189 call interpo_linear_type0(hhxyd,hxyd,irg,itg,ipg)
00190 call interpo_linear_type0(hhxzd,hxzd,irg,itg,ipg)
00191 call interpo_linear_type0(hhyyd,hyyd,irg,itg,ipg)
00192 call interpo_linear_type0(hhyzd,hyzd,irg,itg,ipg)
00193 call interpo_linear_type0(hhzzd,hzzd,irg,itg,ipg)
00194 hhyxd = hhxyd
00195 hhzxd = hhxzd
00196 hhzyd = hhyzd
00197 gamd(1,1) = hhxxd + 1.0d0
00198 gamd(1,2) = hhxyd
00199 gamd(1,3) = hhxzd
00200 gamd(2,2) = hhyyd + 1.0d0
00201 gamd(2,3) = hhyzd
00202 gamd(3,3) = hhzzd + 1.0d0
00203 gamd(2,1) = gamd(1,2)
00204 gamd(3,1) = gamd(1,3)
00205 gamd(3,2) = gamd(2,3)
00206
00207 call interpo_linear_type0(ovxgc,ovxd,irg,itg,ipg)
00208 call interpo_linear_type0(ovygc,ovyd,irg,itg,ipg)
00209 call interpo_linear_type0(ovzgc,ovzd,irg,itg,ipg)
00210 ovd(1) = ovxgc
00211 ovd(2) = ovygc
00212 ovd(3) = ovzgc
00213 call interpo_linear_type0(ovxgc,ovxu,irg,itg,ipg)
00214 call interpo_linear_type0(ovygc,ovyu,irg,itg,ipg)
00215 call interpo_linear_type0(ovzgc,ovzu,irg,itg,ipg)
00216 ovu(1) = ovxgc
00217 ovu(2) = ovygc
00218 ovu(3) = ovzgc
00219
00220
00221 call interpo_linear_type0(bvxgc,bvxu,irg,itg,ipg)
00222 call interpo_linear_type0(bvygc,bvyu,irg,itg,ipg)
00223 call interpo_linear_type0(bvzgc,bvzu,irg,itg,ipg)
00224
00225 dapsi(1:3) = gradp(irg,itg,ipg,1:3)
00226 daalph(1:3) = grada(irg,itg,ipg,1:3)
00227 daalps(1:3) = gradap(irg,itg,ipg,1:3)
00228 daalps2(1:3) = gradap2(irg,itg,ipg,1:3)
00229 cdbv(1,1:3) = cdbvxd(irg,itg,ipg,1:3)
00230 cdbv(2,1:3) = cdbvyd(irg,itg,ipg,1:3)
00231 cdbv(3,1:3) = cdbvzd(irg,itg,ipg,1:3)
00232 do ia = 1, 3
00233 aij(ia,1:3) = tfkij(irg,itg,ipg,ia,1:3)
00234 dbv(ia,1:3) = gradb(irg,itg,ipg,ia,1:3)
00235 end do
00236 daij(1:6,1:3) = gradk(irg,itg,ipg,1:6,1:3)
00237
00238
00239
00240 cdivbc = cdivbv(irg,itg,ipg)
00241
00242
00243
00244
00245 call dadbscalar_type3(alps2,dabalps2,irg,itg,ipg)
00246
00247 gadadp=gamu(1,1)*daalph(1)*dapsi(1)+gamu(1,2)*daalph(1)*dapsi(2) &
00248 & +gamu(1,3)*daalph(1)*dapsi(3)+gamu(2,1)*daalph(2)*dapsi(1) &
00249 & +gamu(2,2)*daalph(2)*dapsi(2)+gamu(2,3)*daalph(2)*dapsi(3) &
00250 & +gamu(3,1)*daalph(3)*dapsi(1)+gamu(3,2)*daalph(3)*dapsi(2) &
00251 & +gamu(3,3)*daalph(3)*dapsi(3)
00252
00253
00254
00255 sw1 = swlp(iswl)
00256 sw2 = 1.0d0 - sw1
00257 sw3 = swls(iswl)
00258
00259
00260
00261 ovpadpsi=ps43oal*(ovxgc*dapsi(1)+ovygc*dapsi(2)+ovzgc*dapsi(3))
00262 dov(1:3,1:3) = dbv(1:3,1:3)
00263 dov(1,2) = dbv(1,2) + ome*dphiu(1,2)
00264 dov(2,1) = dbv(2,1) + ome*dphiu(2,1)
00265
00266
00267
00268 if (chgra == 'h'.or.chgra == 'c'.or.chgra == 'C' &
00269 &.or.chgra == 'H'.or.chgra == 'W') then
00270
00271 rlbgab(1,1) = rlbxx(irg,itg,ipg)
00272 rlbgab(1,2) = rlbxy(irg,itg,ipg)
00273 rlbgab(1,3) = rlbxz(irg,itg,ipg)
00274 rlbgab(2,1) = rlbxy(irg,itg,ipg)
00275 rlbgab(2,2) = rlbyy(irg,itg,ipg)
00276 rlbgab(2,3) = rlbyz(irg,itg,ipg)
00277 rlbgab(3,1) = rlbxz(irg,itg,ipg)
00278 rlbgab(3,2) = rlbyz(irg,itg,ipg)
00279 rlbgab(3,3) = rlbzz(irg,itg,ipg)
00280 rlpgab(1,1) = rlpxx(irg,itg,ipg)
00281 rlpgab(1,2) = rlpxy(irg,itg,ipg)
00282 rlpgab(1,3) = rlpxz(irg,itg,ipg)
00283 rlpgab(2,1) = rlpxy(irg,itg,ipg)
00284 rlpgab(2,2) = rlpyy(irg,itg,ipg)
00285 rlpgab(2,3) = rlpyz(irg,itg,ipg)
00286 rlpgab(3,1) = rlpxz(irg,itg,ipg)
00287 rlpgab(3,2) = rlpyz(irg,itg,ipg)
00288 rlpgab(3,3) = rlpzz(irg,itg,ipg)
00289
00290 call grdphi_midpoint_type0(hxxd,dphxxd,irg,itg,ipg)
00291 call grdphi_midpoint_type0(hxyd,dphxyd,irg,itg,ipg)
00292 call grdphi_midpoint_type0(hxzd,dphxzd,irg,itg,ipg)
00293 call grdphi_midpoint_type0(hyyd,dphyyd,irg,itg,ipg)
00294 call grdphi_midpoint_type0(hyzd,dphyzd,irg,itg,ipg)
00295 call grdphi_midpoint_type0(hzzd,dphzzd,irg,itg,ipg)
00296 dphyxd = dphxyd
00297 dphzxd = dphxzd
00298 dphzyd = dphyzd
00299 dphabd(1,1) = dphxxd
00300 dphabd(1,2) = dphxyd
00301 dphabd(1,3) = dphxzd
00302 dphabd(2,1) = dphyxd
00303 dphabd(2,2) = dphyyd
00304 dphabd(2,3) = dphyzd
00305 dphabd(3,1) = dphzxd
00306 dphabd(3,2) = dphzyd
00307 dphabd(3,3) = dphzzd
00308 call grdphi_midpoint_type0(hxxu,dphxxu,irg,itg,ipg)
00309 call grdphi_midpoint_type0(hxyu,dphxyu,irg,itg,ipg)
00310 call grdphi_midpoint_type0(hxzu,dphxzu,irg,itg,ipg)
00311 call grdphi_midpoint_type0(hyyu,dphyyu,irg,itg,ipg)
00312 call grdphi_midpoint_type0(hyzu,dphyzu,irg,itg,ipg)
00313 call grdphi_midpoint_type0(hzzu,dphzzu,irg,itg,ipg)
00314 dphyxu = dphxyu
00315 dphzxu = dphxzu
00316 dphzyu = dphyzu
00317 dphabu(1,1) = dphxxu
00318 dphabu(1,2) = dphxyu
00319 dphabu(1,3) = dphxzu
00320 dphabu(2,1) = dphyxu
00321 dphabu(2,2) = dphyyu
00322 dphabu(2,3) = dphyzu
00323 dphabu(3,1) = dphzxu
00324 dphabu(3,2) = dphzyu
00325 dphabu(3,3) = dphzzu
00326 call grd2phi_midpoint_type0(hxxd,d2phxx,irg,itg,ipg)
00327 call grd2phi_midpoint_type0(hxyd,d2phxy,irg,itg,ipg)
00328 call grd2phi_midpoint_type0(hxzd,d2phxz,irg,itg,ipg)
00329 call grd2phi_midpoint_type0(hyyd,d2phyy,irg,itg,ipg)
00330 call grd2phi_midpoint_type0(hyzd,d2phyz,irg,itg,ipg)
00331 call grd2phi_midpoint_type0(hzzd,d2phzz,irg,itg,ipg)
00332 d2phyx = d2phxy
00333 d2phzx = d2phxz
00334 d2phzy = d2phyz
00335 d2phab(1,1) = d2phxx
00336 d2phab(1,2) = d2phxy
00337 d2phab(1,3) = d2phxz
00338 d2phab(2,1) = d2phyx
00339 d2phab(2,2) = d2phyy
00340 d2phab(2,3) = d2phyz
00341 d2phab(3,1) = d2phzx
00342 d2phab(3,2) = d2phzy
00343 d2phab(3,3) = d2phzz
00344
00345
00346
00347
00348 call grdphi_midpoint_type0(rlbxx_grid,dpdbgxx,irg,itg,ipg)
00349 call grdphi_midpoint_type0(rlbxy_grid,dpdbgxy,irg,itg,ipg)
00350 call grdphi_midpoint_type0(rlbxz_grid,dpdbgxz,irg,itg,ipg)
00351 call grdphi_midpoint_type0(rlbyy_grid,dpdbgyy,irg,itg,ipg)
00352 call grdphi_midpoint_type0(rlbyz_grid,dpdbgyz,irg,itg,ipg)
00353 call grdphi_midpoint_type0(rlbzz_grid,dpdbgzz,irg,itg,ipg)
00354 dpdbgyx = dpdbgxy
00355 dpdbgzx = dpdbgxz
00356 dpdbgzy = dpdbgyz
00357 rlplbg(1,1)=dpdbgxx+rlbgab(1,2)*dphiu(2,1)+rlbgab(2,1)*dphiu(2,1)
00358 rlplbg(1,2)=dpdbgxy+rlbgab(1,1)*dphiu(1,2)+rlbgab(2,2)*dphiu(2,1)
00359 rlplbg(1,3)=dpdbgxz+rlbgab(2,3)*dphiu(2,1)
00360 rlplbg(2,1)=dpdbgyx+rlbgab(2,2)*dphiu(2,1)+rlbgab(1,1)*dphiu(1,2)
00361 rlplbg(2,2)=dpdbgyy+rlbgab(2,1)*dphiu(1,2)+rlbgab(1,2)*dphiu(1,2)
00362 rlplbg(2,3)=dpdbgyz+rlbgab(1,3)*dphiu(1,2)
00363 rlplbg(3,1) = dpdbgzx
00364 rlplbg(3,2) = dpdbgzy
00365 rlplbg(3,3) = dpdbgzz
00366
00367 rlielnp = 8.0d0*psiinv &
00368 & *(ovu(1)* dapsi(1) + ovu(2)* dapsi(2) + ovu(3)* dapsi(3))
00369 rlielna = alpgcinv &
00370 & *(ovu(1)*daalph(1) + ovu(2)*daalph(2) + ovu(3)*daalph(3))
00371
00372 call grgrad1g_midpoint(rlbxx_grid,grad,irg,itg,ipg)
00373 dliebg(1,1:3) = grad(1:3)
00374 call grgrad1g_midpoint(rlbxy_grid,grad,irg,itg,ipg)
00375 dliebg(2,1:3) = grad(1:3)
00376 call grgrad1g_midpoint(rlbxz_grid,grad,irg,itg,ipg)
00377 dliebg(3,1:3) = grad(1:3)
00378 call grgrad1g_midpoint(rlbyy_grid,grad,irg,itg,ipg)
00379 dliebg(4,1:3) = grad(1:3)
00380 call grgrad1g_midpoint(rlbyz_grid,grad,irg,itg,ipg)
00381 dliebg(5,1:3) = grad(1:3)
00382 call grgrad1g_midpoint(rlbzz_grid,grad,irg,itg,ipg)
00383 dliebg(6,1:3) = grad(1:3)
00384
00385 call grgrad1g_midpoint(rlpxx_grid,grad,irg,itg,ipg)
00386 dliepg(1,1:3) = grad(1:3)
00387 call grgrad1g_midpoint(rlpxy_grid,grad,irg,itg,ipg)
00388 dliepg(2,1:3) = grad(1:3)
00389 call grgrad1g_midpoint(rlpxz_grid,grad,irg,itg,ipg)
00390 dliepg(3,1:3) = grad(1:3)
00391 call grgrad1g_midpoint(rlpyy_grid,grad,irg,itg,ipg)
00392 dliepg(4,1:3) = grad(1:3)
00393 call grgrad1g_midpoint(rlpyz_grid,grad,irg,itg,ipg)
00394 dliepg(5,1:3) = grad(1:3)
00395 call grgrad1g_midpoint(rlpzz_grid,grad,irg,itg,ipg)
00396 dliepg(6,1:3) = grad(1:3)
00397
00398 end if
00399
00400 do ic = 1, 6
00401
00402 ia = 1 + ic/4 + ic/6
00403 ib = ic - (ic/4)*2 - ic/6
00404
00405
00406
00407
00408 rlieaij(ia,ib) = &
00409 & ovxgc*daij(ic,1) + ovygc*daij(ic,2) + ovzgc*daij(ic,3) &
00410 & + aij(ia,1)*dov(1,ib)+aij(ia,2)*dov(2,ib)+aij(ia,3)*dov(3,ib) &
00411 & + aij(1,ib)*dov(1,ia)+aij(2,ib)*dov(2,ia)+aij(3,ib)*dov(3,ia)
00412
00413
00414
00415
00416
00417
00418 rnlab(ia,ib) = rabnl(irg,itg,ipg,ic)
00419
00420
00421
00422 c1ab = cri(irg,itg,ipg,1,ic)
00423 c2ab = cri(irg,itg,ipg,2,ic)
00424 c3ab = cri(irg,itg,ipg,3,ic)
00425
00426 cdalps2 = c1ab*daalps2(1) + c2ab*daalps2(2) + c3ab*daalps2(3)
00427 dapda = 4.0d0*(daalps(ia)*dapsi(ib)+daalps(ib)*dapsi(ia))
00428 rpddal(ia,ib) = ap2gcinv*(- dabalps2(ia,ib) + cdalps2 + dapda)
00429
00430
00431
00432 aacabc(ia,ib) = gamu(1,1)*aij(ia,1)*aij(ib,1) &
00433 & + gamu(1,2)*aij(ia,1)*aij(ib,2) &
00434 & + gamu(1,3)*aij(ia,1)*aij(ib,3) &
00435 & + gamu(2,1)*aij(ia,2)*aij(ib,1) &
00436 & + gamu(2,2)*aij(ia,2)*aij(ib,2) &
00437 & + gamu(2,3)*aij(ia,2)*aij(ib,3) &
00438 & + gamu(3,1)*aij(ia,3)*aij(ib,1) &
00439 & + gamu(3,2)*aij(ia,3)*aij(ib,2) &
00440 & + gamu(3,3)*aij(ia,3)*aij(ib,3)
00441
00442
00443
00444 heli(ia,ib) = 0.0d0
00445 if (chgra == 'h'.or.chgra == 'H'.or.chgra == 'W') then
00446
00447 rlblpg(ia,ib) = &
00448 & bvxgc*dliepg(ic,1)+bvygc*dliepg(ic,2)+bvzgc*dliepg(ic,3) &
00449 & + rlpgab(ia,1)*dov(1,ib) + rlpgab(ia,2)*dov(2,ib) &
00450 & + rlpgab(ia,3)*dov(3,ib) + rlpgab(1,ib)*dov(1,ia) &
00451 & + rlpgab(2,ib)*dov(2,ia) + rlpgab(3,ib)*dov(3,ia)
00452 rlblbg(ia,ib) = &
00453 & bvxgc*dliebg(ic,1)+bvygc*dliebg(ic,2)+bvzgc*dliebg(ic,3) &
00454 & + rlbgab(ia,1)*dov(1,ib) + rlbgab(ia,2)*dov(2,ib) &
00455 & + rlbgab(ia,3)*dov(3,ib) + rlbgab(1,ib)*dov(1,ia) &
00456 & + rlbgab(2,ib)*dov(2,ia) + rlbgab(3,ib)*dov(3,ia)
00457
00458 palidp(ia,ib) = (dphabd(ia,1) + rlpgab(ia,1))*dphiu(1,ib) &
00459 & + (dphabd(ia,2) + rlpgab(ia,2))*dphiu(2,ib) &
00460 & + (dphabd(ia,3) + rlpgab(ia,3))*dphiu(3,ib) &
00461 & + (dphabd(1,ib) + rlpgab(1,ib))*dphiu(1,ia) &
00462 & + (dphabd(2,ib) + rlpgab(2,ib))*dphiu(2,ia) &
00463 & + (dphabd(3,ib) + rlpgab(3,ib))*dphiu(3,ia)
00464
00465 term1 = 0.5d0*(ps4oal2 - 1.0d0)*ome**2*d2phab(ia,ib)
00466 term2 = 0.5d0*ps4oal2*ome**2*palidp(ia,ib)
00467 term3 = 0.5d0*ps4oal2*ome*(rlplbg(ia,ib) + rlblpg(ia,ib))
00468 term4 = 0.5d0*ps4oal2*rlblbg(ia,ib)
00469 term5 = ps4oal*aij(ia,ib)*(rlielnp - rlielna)
00470
00471
00472 if (irg == 2.and.itg == 1.and.ipg == 0) then
00473 write(6,'(1p,6e14.6)')term1,term2,term3,term4,term5
00474 end if
00475
00476 tm1a(ia,ib) = term1
00477 tm2a(ia,ib) = term2
00478 tm3a(ia,ib) = term3
00479 tm4a(ia,ib) = term4
00480 tm5a(ia,ib) = term5
00481
00482
00483 if (chgra == 'H') then
00484 term2 = term2*cutoff
00485 term3 = term3*cutoff
00486 end if
00487
00488 heli(ia,ib) = term1 + term2 + term3 + term4 + term5
00489
00490 end if
00491
00492 end do
00493
00494
00495
00496 call grgrad1g_midpoint(hxxd,grad,irg,itg,ipg)
00497 dhxxddx = grad(1)
00498 dhxxddy = grad(2)
00499 dhxxddz = grad(3)
00500 call grgrad1g_midpoint(hxyd,grad,irg,itg,ipg)
00501 dhxyddx = grad(1)
00502 dhxyddy = grad(2)
00503 dhxyddz = grad(3)
00504 call grgrad1g_midpoint(hxzd,grad,irg,itg,ipg)
00505 dhxzddx = grad(1)
00506 dhxzddy = grad(2)
00507 dhxzddz = grad(3)
00508 call grgrad1g_midpoint(hyyd,grad,irg,itg,ipg)
00509 dhyyddx = grad(1)
00510 dhyyddy = grad(2)
00511 dhyyddz = grad(3)
00512 call grgrad1g_midpoint(hyzd,grad,irg,itg,ipg)
00513 dhyzddx = grad(1)
00514 dhyzddy = grad(2)
00515 dhyzddz = grad(3)
00516 call grgrad1g_midpoint(hzzd,grad,irg,itg,ipg)
00517 dhzzddx = grad(1)
00518 dhzzddy = grad(2)
00519 dhzzddz = grad(3)
00520
00521 dhyxddx = dhxyddx
00522 dhyxddy = dhxyddy
00523 dhyxddz = dhxyddz
00524 dhzxddx = dhxzddx
00525 dhzxddy = dhxzddy
00526 dhzxddz = dhxzddz
00527 dhzyddx = dhyzddx
00528 dhzyddy = dhyzddy
00529 dhzyddz = dhyzddz
00530
00531 call grgrad1g_midpoint(hxxu,grad,irg,itg,ipg)
00532 dhxxudx = grad(1)
00533 dhxxudy = grad(2)
00534 dhxxudz = grad(3)
00535 call grgrad1g_midpoint(hxyu,grad,irg,itg,ipg)
00536 dhxyudx = grad(1)
00537 dhxyudy = grad(2)
00538 dhxyudz = grad(3)
00539 call grgrad1g_midpoint(hxzu,grad,irg,itg,ipg)
00540 dhxzudx = grad(1)
00541 dhxzudy = grad(2)
00542 dhxzudz = grad(3)
00543 call grgrad1g_midpoint(hyyu,grad,irg,itg,ipg)
00544 dhyyudx = grad(1)
00545 dhyyudy = grad(2)
00546 dhyyudz = grad(3)
00547 call grgrad1g_midpoint(hyzu,grad,irg,itg,ipg)
00548 dhyzudx = grad(1)
00549 dhyzudy = grad(2)
00550 dhyzudz = grad(3)
00551 call grgrad1g_midpoint(hzzu,grad,irg,itg,ipg)
00552 dhzzudx = grad(1)
00553 dhzzudy = grad(2)
00554 dhzzudz = grad(3)
00555
00556 dhyxudx = dhxyudx
00557 dhyxudy = dhxyudy
00558 dhyxudz = dhxyudz
00559 dhzxudx = dhxzudx
00560 dhzxudy = dhxzudy
00561 dhzxudz = dhxzudz
00562 dhzyudx = dhyzudx
00563 dhzyudy = dhyzudy
00564 dhzyudz = dhyzudz
00565
00566
00567
00568 dhdh = dhxxudx*dhxxddx + dhxxudy*dhxxddy + dhxxudz*dhxxddz &
00569 & + dhxyudx*dhxyddx + dhxyudy*dhxyddy + dhxyudz*dhxyddz &
00570 & + dhxzudx*dhxzddx + dhxzudy*dhxzddy + dhxzudz*dhxzddz &
00571 & + dhyxudx*dhyxddx + dhyxudy*dhyxddy + dhyxudz*dhyxddz &
00572 & + dhyyudx*dhyyddx + dhyyudy*dhyyddy + dhyyudz*dhyyddz &
00573 & + dhyzudx*dhyzddx + dhyzudy*dhyzddy + dhyzudz*dhyzddz &
00574 & + dhzxudx*dhzxddx + dhzxudy*dhzxddy + dhzxudz*dhzxddz &
00575 & + dhzyudx*dhzyddx + dhzyudy*dhzyddy + dhzyudz*dhzyddz &
00576 & + dhzzudx*dhzzddx + dhzzudy*dhzzddy + dhzzudz*dhzzddz
00577
00578 dphdph = dphxxu*dphxxd + dphxyu*dphxyd + dphxzu*dphxzd &
00579 & + dphyxu*dphyxd + dphyyu*dphyyd + dphyzu*dphyzd &
00580 & + dphzxu*dphzxd + dphzyu*dphzyd + dphzzu*dphzzd
00581
00582 heli(2,1) = heli(1,2)
00583 heli(3,1) = heli(1,3)
00584 heli(3,2) = heli(2,3)
00585 rlieaij(2,1) = rlieaij(1,2)
00586 rlieaij(3,1) = rlieaij(1,3)
00587 rlieaij(3,2) = rlieaij(2,3)
00588 rnlab(2,1) = rnlab(1,2)
00589 rnlab(3,1) = rnlab(1,3)
00590 rnlab(3,2) = rnlab(2,3)
00591 rpddal(2,1) = rpddal(1,2)
00592 rpddal(3,1) = rpddal(1,3)
00593 rpddal(3,2) = rpddal(2,3)
00594 aacabc(2,1) = aacabc(1,2)
00595 aacabc(3,1) = aacabc(1,3)
00596 aacabc(3,2) = aacabc(2,3)
00597 sab(2,1) = sab(1,2)
00598 sab(3,1) = sab(1,3)
00599 sab(3,2) = sab(2,3)
00600
00601 trheli = 0.0d0
00602 traij = 0.0d0
00603 trlieaij = 0.0d0
00604 trrnlab = 0.0d0
00605 trrpddal = 0.0d0
00606 traa = 0.0d0
00607 trsab = 0.0d0
00608 trd2ph = 0.0d0
00609 do ib = 1, 3
00610 do ia = 1, 3
00611 trheli = trheli + gamu(ia,ib)* heli(ia,ib)
00612 traij = traij + gamu(ia,ib)* aij(ia,ib)
00613 trlieaij = trlieaij + gamu(ia,ib)*rlieaij(ia,ib)
00614 trrnlab = trrnlab + gamu(ia,ib)* rnlab(ia,ib)
00615 trrpddal = trrpddal + gamu(ia,ib)* rpddal(ia,ib)
00616 traa = traa + gamu(ia,ib)* aacabc(ia,ib)
00617 trsab = trsab + gamu(ia,ib)* sab(ia,ib)
00618 trd2ph = trd2ph + gamu(ia,ib)* d2phab(ia,ib)
00619 end do
00620 end do
00621
00622 do ic = 1, 6
00623
00624 ia = 1 + ic/4 + ic/6
00625 ib = ic - (ic/4)*2 - ic/6
00626
00627 sgam = san*gamd(ia,ib)
00628 tfaij = aij(ia,ib) - sgam*traij
00629 tflieaij = rlieaij(ia,ib) - sgam*trlieaij
00630 tfrnlab = rnlab(ia,ib) - sgam*trrnlab
00631 tfrpddal = rpddal(ia,ib) - sgam*trrpddal
00632 tfaa = aacabc(ia,ib) - sgam*traa
00633
00634 tfd2ph = d2phab(ia,ib) - sgam*trd2ph
00635
00636 gamdhdh = gamd(ia,ib)*dhdh
00637 ogdphdph = 0.0d0
00638 if (chope == 'H') ogdphdph = gamd(ia,ib)*ome**2*dphdph
00639
00640 tfliep4aij = 0.0d0
00641 tfheli = 0.0d0
00642 if (chgra == 'w') tfliep4aij = tfaij*ovpadpsi + ps4oal*tflieaij
00643
00644
00645 if (chgra == 'c') tfliep4aij = tfaij*ovpadpsi + ps4oal*tflieaij
00646 if (chgra == 'C') tfliep4aij = tfaij*ovpadpsi + ps4oal*tflieaij &
00647 & - 0.5d0*ome**2*tfd2ph*cutoff
00648 if (chgra == 'h') tfheli = heli(ia,ib) - sgam*trheli
00649 if (chgra == 'H') tfheli = heli(ia,ib) - sgam*trheli
00650 if (chgra == 'W') then
00651 tfheli = heli(ia,ib) - sgam*trheli
00652 tfheli = tfheli*cutoff
00653 tfliep4aij = tfaij*ovpadpsi + ps4oal*tflieaij
00654 tfliep4aij = tfliep4aij*(1.0d0 - cutoff)
00655 end if
00656
00657 souten(irg,itg,ipg,ic) = &
00658 & 2.0d0*( tfheli + tfliep4aij &
00659 & + tfrnlab + tfrpddal - 2.0d0*psigc4*tfaa &
00660 & - roku*gamdhdh + roku*ogdphdph )
00661
00662
00663
00664
00665
00666
00667
00668
00669
00670
00671 end do
00672 end do
00673 end do
00674 end do
00675
00676
00677 deallocate(dfdx)
00678 deallocate(dfdy)
00679 deallocate(dfdz)
00680 deallocate(fnc0)
00681 deallocate(grada)
00682 deallocate(gradp)
00683 deallocate(gradap)
00684 deallocate(gradap2)
00685 deallocate(gradk)
00686 deallocate(gradb)
00687
00688 end subroutine sourceterm_trfreeG_WL