Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! multipole.f90 is part of AstroBEAR. 00008 ! 00009 ! AstroBEAR is free software: you can redistribute it and/or modify 00010 ! it under the terms of the GNU General Public License as published by 00011 ! the Free Software Foundation, either version 3 of the License, or 00012 ! (at your option) any later version. 00013 ! 00014 ! AstroBEAR is distributed in the hope that it will be useful, 00015 ! but WITHOUT ANY WARRANTY; without even the implied warranty of 00016 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00017 ! GNU General Public License for more details. 00018 ! 00019 ! You should have received a copy of the GNU General Public License 00020 ! along with AstroBEAR. If not, see <http://www.gnu.org/licenses/>. 00021 ! 00022 !######################################################################### 00025 00029 00030 00033 module multipole 00034 USE GlobalDeclarations 00035 USE PhysicsDeclarations 00036 USE TreeDeclarations 00037 USE DataDeclarations 00038 USE CommonFunctions 00039 USE Timing 00040 IMPLICIT NONE 00041 COMPLEX(8), DIMENSION(:), ALLOCATABLE :: moments 00042 INTEGER, DIMENSION(9,2) :: momentfactors 00043 REAL(KIND=qPREC) :: CM_00,CM_10,CM_11,CM_20,CM_21,CM_22 00044 LOGICAL :: lDipole, lQuadrupole,lReflect(3), lMultiPole(3) 00045 REAL(KIND=qPREC) :: massfact, costheta_fact, sintheta_fact, cosphi_fact, sinphi_fact, cos2phi_fact, sin2phi_fact, MultiPole_COM(3) 00046 REAL(KIND=qPREC) :: multipole_radius 00047 SAVE 00048 contains 00049 00051 SUBROUTINE InitMultiPoleMoments 00052 REAL(KIND=qPREC) :: dV 00053 lDipole=.true. 00054 lQuadrupole=.true. 00055 ! lDipole=.false. 00056 ! lQuadrupole=.false. 00057 dV=levels(0)%dx**nDim 00058 IF (nDim == 3) THEN 00059 IF (lDipole) THEN 00060 IF (lQuadrupole) THEN 00061 ALLOCATE(moments(9)) 00062 ELSE 00063 ALLOCATE(moments(4)) 00064 END IF 00065 ELSE 00066 ALLOCATE(moments(1)) 00067 END IF 00068 ELSE 00069 IF (lDipole) THEN 00070 write(*,*) 'WARNING*** Dipole expansion not supported in 2D' 00071 END IF 00072 lDipole=.false. 00073 lQuadrupole=.false. 00074 IF (lDipole) THEN 00075 IF (lQuadrupole) THEN 00076 ALLOCATE(moments(5)) 00077 ELSE 00078 ALLOCATE(moments(3)) 00079 END IF 00080 ELSE 00081 ALLOCATE(moments(1)) 00082 END IF 00083 END IF 00084 00085 00086 ! CM_00=-4d0*Pi*ScaleGrav*dV 00087 ! CM_10=-4d0*Pi*ScaleGrav*dV 00088 ! CM_11=-.5d0*4d0*Pi*ScaleGrav*dV 00089 ! CM_20=-.25d0*4d0*Pi*ScaleGrav*dV 00090 ! CM_21=-3d0/2d0*4d0*Pi*ScaleGrav*dV 00091 ! CM_22=-3d0/8d0*4d0*Pi*ScaleGrav*dV 00092 CM_00=-ScaleGrav*dV 00093 CM_10=-ScaleGrav*dV 00094 CM_11=-.5d0*ScaleGrav*dV 00095 CM_20=-.25d0*ScaleGrav*dV 00096 CM_21=-3d0/2d0*ScaleGrav*dV 00097 CM_22=-3d0/8d0*ScaleGrav*dV 00098 00099 IF (nDim == 3) THEN 00100 00101 !-----------------------------------------------------------------------! 00102 ! *** Symmetry of Spherical Harmonics *** ! 00103 !-----------------------------------------------------------------------! 00104 ! Term Symmetry Real/Imag! 00105 !-----------------------------------------------------------------------! 00106 ! Y_0 xyz / xyz ! 00107 !-----------------------------------------------------------------------! 00108 ! Y_1_-1 = sqrt(1-z^2)*(x-iy) yz / xz ! 00109 ! Y_1_0 = z xy / xy ! 00110 ! Y_1_1 = sqrt(1-z^2)*(x+iy) yz / xz ! 00111 !-----------------------------------------------------------------------! 00112 ! Y_2_-2 = (1-z^2)*(x-iy)^2 = (1-z^2)*(x^2+y^2-2*i*x*y) xyz / z ! 00113 ! Y_2_-1 = z*sqrt(1-z^2)*(x-iy) y / x ! 00114 ! Y_2_0 = (3*z^2-1) xyz / xyz ! 00115 ! Y_2_1 = z*sqrt(1-z^2)*(x-iy) y / x ! 00116 ! Y_2_-2 = (1-z^2)*(x+iy)^2 = (1-z^2)*(x^2+y^2+2*i*x*y) xyz / z ! 00117 !-----------------------------------------------------------------------! 00118 00119 00120 !----------------------------------------------------------| 00121 ! *** Multiplication factors given grid symmetry *** | 00122 !----------------------------------------------------------| 00123 ! symmetry | 0 | -1 0 1 | -2 -1 0 1 2 | 00124 ! none | 1/1 | 1/1 1/1 1/1 | 1/1 1/1 1/1 1/1 1/1 | 00125 ! x only | 2/2 | 0/2 2/2 0/2 ! 2/0 0/2 2/2 0/2 2/0 | 00126 ! y only | 2/2 | 2/0 2/2 2/0 | 2/0 0/2 2/2 0/2 2/0 | 00127 ! z only | 2/2 | 2/2 0/0 2/2 | 2/2 0/0 2/2 0/0 2/2 | 00128 ! xy | 4/4 | 0/0 4/4 0/0 | 4/0 0/4 4/4 0/4 4/0 | 00129 ! xz | 4/4 | 0/4 0/0 0/4 | 4/0 0/0 4/4 0/0 4/0 | 00130 ! yz | 4/4 | 4/0 0/0 4/0 | 4/0 0/0 4/4 0/0 4/0 | 00131 ! xyz | 8/8 | 0/0 0/0 0/0 | 8/0 0/0 8/8 0/0 8/0 | 00132 !----------------------------------------------------------| 00133 MomentFactors(:,:)=1 00134 IF (lReflect(1)) MomentFactors=MomentFactors*reshape((/2,2, 0,2, 2,2, 0,2, 2,0, 0,2, 2,2, 0,2, 2,0/),(/9,2/),(/0/),(/2,1/)) 00135 IF (lReflect(2)) MomentFactors=MomentFactors*reshape((/2,2, 2,0, 2,2, 2,0, 2,0, 0,2, 2,2, 0,2, 2,0/),(/9,2/),(/0/),(/2,1/)) 00136 IF (lReflect(3)) MomentFactors=MomentFactors*reshape((/2,2, 2,2, 0,0, 2,2, 2,2, 0,0, 2,2, 0,0, 2,2/),(/9,2/),(/0/),(/2,1/)) 00137 ELSE 00138 MomentFactors(:,:)=1 00139 IF (lReflect(1)) MomentFactors(1,:)=MomentFactors(1,:)*(/2,2/) 00140 IF (lReflect(2)) MomentFactors(1,:)=MomentFactors(1,:)*(/2,2/) 00141 END IF 00142 END SUBROUTINE InitMultiPoleMoments 00143 00144 00146 SUBROUTINE CalcMultiPoleMoments() 00147 REAL(KIND=qPREC) :: rmass(4) 00148 TYPE(NodeDef), POINTER :: node 00149 TYPE(NodeDefList), POINTER :: nodelist 00150 INTEGER :: iErr 00151 TYPE(InfoDef), POINTER :: info 00152 REAL(KIND=qPREC) :: x,y,z,s,r,costheta,sintheta,cosphi,sinphi,cos2phi,sin2phi,r2,r3, dx, dy, dz, rho, rrho, rrrho, pos(3), dm, theta 00153 INTEGER :: jmin, jmax, kmin, kmax, imin, imax, dim, i, j, k 00154 INTEGER, DIMENSION(3,2) :: ip 00155 nodelist=>Nodes(0)%p 00156 ! massfact=1 00157 rmass=0 00158 moments=0 00159 DO WHILE(ASSOCIATED(nodelist)) 00160 info=>nodelist%self%info 00161 DO dim=1,nDim 00162 IF (lMultiPole(dim)) THEN 00163 ip(:,1)=1 00164 ip(:,2)=Info%mX(:) 00165 DO i=1, Info%mX(dim) 00166 ip(dim,:)=i 00167 x=Info%xBounds(dim,1)+(REAL(i)-half)*levels(0)%dx 00168 dm=SUM(Info%q(ip(1,1):ip(1,2), ip(2,1):ip(2,2), ip(3,1):ip(3,2),1)) 00169 IF (dim == 1) rmass(4)=rmass(4)+dm 00170 rmass(dim) = rmass(dim) + dm*x 00171 END DO 00172 END IF 00173 END DO 00174 nodelist=>nodelist%next 00175 END DO 00176 CALL StartTimer(iBarrier, 0) 00177 CALL MPI_ALLREDUCE(MPI_IN_PLACE, rmass, 4, MPI_DOUBLE_PRECISION, MPI_SUM, levels(0)%MPI_COMM, iErr) 00178 CALL StopTimer(iBarrier, 0) 00179 multipole_com=rmass(1:3)/rmass(4) 00180 DO dim=1,nDim 00181 IF (lReflect(dim)) THEN !ANY(Gmthbc(dim,edge) == (/REFLECT_WALL, REFLECT_BPARALLEL, REFLECT_CYLINDRICAL/))) THEN 00182 multipole_com(dim)=GxBounds(dim,1) 00183 END IF 00184 END DO 00185 multipole_radius=huge(multipole_radius) 00186 DO dim=1,nDim 00187 IF (lMultiPole(dim)) THEN 00188 multipole_radius=min(multipole_radius,abs(multipole_com(dim)-GxBounds(dim,2))) 00189 END IF 00190 END DO 00191 ! write(*,*) multipole_com, multipole_radius 00192 SELECT CASE(nDim) 00193 CASE(2) 00194 z=GxBounds(3,1) 00195 nodelist=>Nodes(0)%p 00196 DO WHILE(ASSOCIATED(nodelist)) 00197 info=>nodelist%self%info 00198 00199 dx=multipole_radius 00200 imin=max(1,floor((multipole_com(1)-dx-Info%xBounds(1,1))/levels(0)%dx+half)) 00201 imax=min(Info%mX(1),ceiling((multipole_com(1)+dx-Info%xbounds(1,1))/levels(0)%dx+half)) 00202 DO i=imin,imax 00203 x=Info%Xbounds(1,1)+(REAL(i)-half)*levels(0)%dx - multipole_com(1) 00204 00205 dy=sqrt(multipole_radius**2-x**2) 00206 jmin=max(1,floor((multipole_com(2)-dy-Info%xbounds(2,1))/levels(0)%dx+half)) 00207 jmax=min(Info%mX(2),ceiling((multipole_com(2)+dy-Info%xbounds(2,1))/levels(0)%dx+half)) 00208 DO j=jmin,jmax 00209 y=Info%xbounds(2,1)+(REAL(j)-half)*levels(0)%dx - multipole_com(2) 00210 pos=(/x,y,z/) 00211 00212 r=sqrt(SUM(pos(1:nDim)**2)) 00213 IF (r < multipole_radius) THEN 00214 00215 rho=Info%q(i,j,1,1) 00216 moments(1)=moments(1)+rho 00217 ! write(*,*) moments(1), rho, massfact 00218 IF (lDipole) THEN 00219 moments(2)=moments(2)+rho*x 00220 moments(3)=moments(3)+rho*y 00221 IF (lQuadrupole) THEN 00222 rrho=rho/2d0*r**2 00223 theta=getphi(x,y) 00224 moments(4)=moments(4)+rrho*(cos(2d0*theta)) 00225 moments(5)=moments(5)+rrho*(sin(2d0*theta)) 00226 END IF 00227 END IF 00228 END IF 00229 END DO 00230 END DO 00231 nodelist=>nodelist%next 00232 END DO 00233 CASE(3) 00234 00235 00236 00237 00238 nodelist=>Nodes(0)%p 00239 DO WHILE(ASSOCIATED(nodelist)) 00240 info=>nodelist%self%info 00241 00242 dx=multipole_radius 00243 imin=max(1,floor((multipole_com(1)-dx-Info%xBounds(1,1))/levels(0)%dx+half)) 00244 imax=min(Info%mX(1),ceiling((multipole_com(1)+dx-Info%xbounds(1,1))/levels(0)%dx+half)) 00245 DO i=imin,imax 00246 x=Info%Xbounds(1,1)+(REAL(i)-half)*levels(0)%dx - multipole_com(1) 00247 00248 dy=sqrt(multipole_radius**2-x**2) 00249 jmin=max(1,floor((multipole_com(2)-dy-Info%xbounds(2,1))/levels(0)%dx+half)) 00250 jmax=min(Info%mX(2),ceiling((multipole_com(2)+dy-Info%xbounds(2,1))/levels(0)%dx+half)) 00251 DO j=jmin,jmax 00252 y=Info%xbounds(2,1)+(REAL(j)-half)*levels(0)%dx - multipole_com(2) 00253 00254 dz=sqrt(dy**2-y**2) 00255 kmin=max(1,floor((multipole_com(3)-dz-Info%xbounds(3,1))/levels(0)%dx+half)) 00256 kmax=min(Info%mX(3),ceiling((multipole_com(3)+dz-Info%xbounds(3,1))/levels(0)%dx+half)) 00257 DO k=kmin,kmax 00258 z=Info%xbounds(3,1)+(REAL(k)-half)*levels(0)%dx - multipole_com(3) 00259 00260 pos=(/x,y,z/) 00261 r=sqrt(SUM(pos(1:nDim)**2)) 00262 IF (r < multipole_radius) THEN 00263 00264 00265 rho=Info%q(i,j,k,1) 00266 moments(1)=moments(1)+rho 00267 00268 00269 IF (lDiPole) THEN 00270 s=sqrt(x**2+y**2) 00271 costheta=z/r 00272 sintheta=s/r 00273 cosphi=x/s 00274 sinphi=y/s 00275 IF (lQuadrupole) THEN 00276 cos2phi=(cosphi**2-sinphi**2) 00277 sin2phi=2d0*cosphi*sinphi 00278 END IF 00279 cosphi=cosphi 00280 sinphi=sinphi 00281 rrho=r*rho 00282 moments(2)=moments(2)+sintheta*cmplx(cosphi,+sinphi) * rrho 00283 moments(3)=moments(3)+costheta * rrho 00284 moments(4)=moments(4)+sintheta*cmplx(cosphi,-sinphi) * rrho 00285 00286 IF (lQuadrupole) THEN 00287 rrrho=r*rrho 00288 cos2phi=(cosphi**2-sinphi**2) * cos2phi_fact 00289 sin2phi=2d0*cosphi*sinphi * sin2phi_fact 00290 moments(5)=moments(5)+sintheta**2*cmplx(cos2phi,sin2phi)*rrrho 00291 moments(6)=moments(6)+sintheta*costheta*cmplx(cosphi,+sinphi) * rrrho 00292 moments(7)=moments(7)+(3d0*(z/r)**2-1d0) * rrrho 00293 moments(8)=moments(8)+sintheta*costheta*cmplx(cosphi,-sinphi) * rrrho 00294 moments(9)=moments(9)+sintheta**2*cmplx(cos2phi,-sin2phi) * rrrho 00295 END IF 00296 END IF 00297 END IF 00298 END DO 00299 END DO 00300 END DO 00301 nodelist=>nodelist%next 00302 END DO 00303 END SELECT 00304 DO i=1,size(moments) 00305 moments(i)=cmplx(real(moments(i))*momentfactors(i,1),aimag(moments(i))*momentfactors(i,2)) 00306 END DO 00307 CALL MPI_ALLREDUCE(MPI_IN_PLACE, moments, size(moments), MPI_DOUBLE_COMPLEX, MPI_SUM, levels(0)%MPI_COMM, iErr) 00308 ! write(*,*) 'moments=', moments*levels(0)%dx**nDim 00309 IF (nDim == 3) THEN 00310 moments(1)=moments(1)*CM_00 00311 !@ write(*,*) 'moments_now', moments(1), CM_00, ScaleGrav, massfact 00312 IF (lDipole) THEN 00313 moments(2)=moments(2)*CM_11 00314 moments(3)=moments(3)*CM_10 00315 moments(4)=-moments(4)*CM_11 00316 IF (lQuadrupole) THEN 00317 moments(5)=moments(5)*CM_22 00318 moments(6)=moments(6)*CM_21 00319 moments(7)=moments(7)*CM_20 00320 moments(8)=-moments(8)*CM_21 00321 moments(9)=moments(9)*CM_22 00322 END IF 00323 END IF 00324 ELSE 00325 ! write(*,*) 'moments=', moments*levels(0)%dx**nDim 00326 moments(:)=-moments(:)*2d0*CM_00 00327 00328 END IF 00329 ! moments(2:)=0d0 00330 END SUBROUTINE CalcMultiPoleMoments 00331 00332 00336 SUBROUTINE SetPhi(Info, ip) 00337 TYPE(InfoDef) :: Info 00338 INTEGER :: i,j,k 00339 INTEGER, DIMENSION(3,2) :: ip 00340 REAL(KIND=qPREC) :: x,y,z,s,r,costheta,sintheta,cosphi,sinphi,r2,r3,s2 00341 ! write(*,*) 'setting phi', ip 00342 DO i=ip(1,1),ip(1,2) 00343 x=Info%xbounds(1,1)+(REAL(i)-half)*levels(Info%level)%dx - multipole_com(1) 00344 00345 DO j=ip(2,1), ip(2,2) 00346 y=Info%xbounds(2,1)+(REAL(j)-half)*levels(Info%level)%dx - multipole_com(2) 00347 s2=x**2+y**2 00348 s=sqrt(s2) 00349 cosphi=x/s 00350 sinphi=y/s 00351 DO k=ip(3,1),ip(3,2) 00352 IF (nDim == 3) THEN 00353 z=Info%xbounds(3,1)+(REAL(k)-half)*levels(Info%level)%dx - multipole_com(3) 00354 r2=(s**2+z**2) 00355 r=sqrt(r2) 00356 costheta=z/r 00357 sintheta=s/r 00358 END IF 00359 ! write(*,*) r, moments(1) 00360 IF (nDim == 3) THEN 00361 Info%q(i,j,k,iPhiGas)=moments(1)/r 00362 00363 IF (lDiPole) THEN 00364 Info%q(i,j,k,iPhiGas)=Info%q(i,j,k,iPhiGas)+& 00365 abs(moments(2)*sintheta*cmplx(cosphi,-sinphi) + & 00366 moments(3)*costheta + & 00367 moments(4)*sintheta*cmplx(cosphi,+sinphi))/r2 00368 IF (lQuadrupole) THEN 00369 r3=r2*r 00370 Info%q(i,j,k,iPhiGas)=Info%q(i,j,k,iPhiGas)+ & 00371 abs(moments(5)*(sintheta*cmplx(cosphi,-sinphi))**2 + & 00372 moments(6)*sintheta*costheta*cmplx(cosphi,-sinphi) + & 00373 moments(7)*(3d0*costheta**2-1d0) + & 00374 moments(8)*sintheta*costheta*cmplx(cosphi,+sinphi) + & 00375 moments(9)*(sintheta*cmplx(cosphi,+sinphi))**2)/r3 00376 END IF 00377 END IF 00378 ELSE 00379 Info%q(i,j,k,iPhiGas)=moments(1)*log(s/R2DEff) 00380 ! write(*,*) moments(1), s, Info%q(i,j,k,iPhiGas) !REAL(moments(1)*log(s/R2Deff)) 00381 IF (lDiPole) THEN 00382 ! Info%q(i,j,k,iPhiGas)=Info%q(i,j,k,iPhiGas)-& 00383 ! (moments(2)*cosphi+moments(3)*sinphi)/s 00384 IF (lQuadrupole) THEN 00385 ! Info%q(i,j,k,iPhiGas)=Info%q(i,j,k,iPhiGas)-& 00386 ! (moments(4)*(cosphi**2-sinphi**2)+moments(5)*2d0*sinphi*cosphi)/s2 00387 END IF 00388 END IF 00389 END IF 00390 Info%q(i,j,k,iPhiDot)=0d0 00391 ! Info%q(i,j,k,iPhi)=Info%q(i,j,k,iPhiGas)+Info%q(i,j,k,iPhiSinks) 00392 END DO 00393 END DO 00394 END DO 00395 ! write(*,*) maxval(Info%q(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),iPhiGas)) 00396 ! write(*,*) ip 00397 ! pause 00398 END SUBROUTINE SetPhi 00399 END module multipole