00001 subroutine calc_quad_pole_peos
00002   use phys_constant, only  :   long, pi
00003   use grid_parameter, only :   nrf, ntf, npf
00004   use def_matter_parameter, only : radi, ome
00005   use make_array_3d
00006   use make_array_5d
00007   use def_quantities, only : restmass, Iij, Itf, dt1Itf, dt2Itf, dt3Itf, &
00008   &                          LGW, dJdt, hplus, hcross
00009   use interface_source_quad_pole_peos
00010   use interface_vol_int_fluid
00011   implicit none
00012   real(long)          :: fac2pi, fac4pi
00013   real(long)          :: volf
00014   real(long), pointer :: souf(:,:,:), souf5d(:,:,:,:,:)
00015   real(long)          :: theta, phi, rr
00016   integer :: ir, it, ip, i, j, k, l, iquad
00017 
00018   call alloc_array3d(souf, 0, nrf, 0, ntf, 0, npf)
00019   call alloc_array5d(souf5d, 0, nrf, 0, ntf, 0, npf, 1, 3, 1, 3)
00020 
00021   do iquad = 0, 4
00022     call source_quad_pole_peos(souf5d, iquad)
00023     do i = 1, 3
00024       do j = 1, 3
00025         souf(0:nrf,0:ntf,0:npf) = souf5d(0:nrf,0:ntf,0:npf,i,j)
00026         call vol_int_fluid(souf,volf)
00027         if (iquad == 0) Iij(i,j) = radi**5*volf
00028         if (iquad == 1) Itf(i,j) = radi**5*volf
00029         if (iquad == 2) dt1Itf(i,j) = ome*radi**4*volf
00030         if (iquad == 3) dt2Itf(i,j) = -ome**2*radi**3*volf
00031         if (iquad == 4) dt3Itf(i,j) = -ome**3*radi**2*volf
00032       end do
00033     end do
00034   end do
00035 
00036   LGW = 0.0d0
00037   dJdt(1:3) = 0.0d0
00038   do i = 1, 3
00039     do j = 1, 3
00040       LGW = LGW + dt3Itf(i,j)*dt3Itf(i,j)/5.0d0
00041     end do
00042     if (i.eq.1) then; j = 2; k = 3; end if
00043     if (i.eq.2) then; j = 3; k = 1; end if
00044     if (i.eq.3) then; j = 1; k = 2; end if
00045     do l = 1, 3
00046       dJdt(i) = dJdt(i) &
00047     & + (dt2Itf(j,l)*dt3Itf(l,k) - dt2Itf(k,l)*dt3Itf(l,j))*2.0d0/5.0d0
00048     end do
00049   end do
00050 
00051   theta = 0.0d0
00052   phi = 0.0d0
00053   rr = 1.0d0
00054   hplus = ((dt2Itf(1,1)-dt2Itf(2,2))*(cos(theta)**2+1.0d0)*cos(2*phi)/4.0d0 &
00055      &  -  (dt2Itf(1,1)+dt2Itf(2,2)-2.0d0*dt2Itf(3,3))*sin(theta)**2/4.0d0 &
00056      &  +   dt2Itf(1,2)*((cos(theta)**2+1.0d0)/2.0d0)*sin(2*phi) &
00057      &  -   dt2Itf(1,3)*sin(theta)*cos(theta)*cos(phi) &
00058      &  -   dt2Itf(2,3)*sin(theta)*cos(theta)*sin(phi) &
00059      &    )*2.0d0/rr
00060   hcross =( - (dt2Itf(1,1)-dt2Itf(2,2))*cos(theta)*sin(2*phi)/2.0d0 &
00061      &   +     dt2Itf(1,2)*cos(theta)*cos(2*phi) &
00062      &   +     dt2Itf(1,3)*sin(theta)*sin(phi) &
00063      &   -     dt2Itf(2,3)*sin(theta)*cos(phi) &
00064      &    )*2.0d0/rr
00065 
00066 
00067 
00068   deallocate(souf)
00069   deallocate(souf5d)
00070 end subroutine calc_quad_pole_peos