Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! poisson.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 00028 00032 00036 00039 MODULE Poisson 00040 USE EllipticDeclarations 00041 USE GlobalDeclarations 00042 USE TreeDeclarations 00043 USE DataDeclarations 00044 USE PhysicsDeclarations 00045 USE MultiPole 00046 USE EllipticComms 00047 IMPLICIT NONE 00048 INTEGER, DIMENSION(:,:), ALLOCATABLE :: offsets 00049 INTEGER, DIMENSION(:), ALLOCATABLE :: StencilValues 00050 INTEGER, DIMENSION(:), ALLOCATABLE :: iEntries 00051 INTEGER, DIMENSION(3,2) :: Poisson_mthbc 00052 REAL(KIND=qPREC) :: PoissonScale 00053 INTEGER,PARAMETER,PUBLIC :: ZEROSLOPE=0,REFLECTING=1,PERIODIC=2,MULTIPOLE_EXPANSION=3,USERSPECIFIED=4 00054 LOGICAL, DIMENSION(3) :: lPoissonPeriodic 00055 INTEGER :: UpdateFreq=1 00056 SAVE 00057 ! Initializes local variables and creates stencil needed for poisson equation 00058 ! INCLUDE 'mpif.h' 00059 INCLUDE 'fftw3.f' 00060 CONTAINS 00061 00064 00066 SUBROUTINE Poisson_Init(EllipticObject) 00067 ! EllipticObject%stencil - handle for hyper stencil created here 00068 ! EllipticObject%Interface - (Currently only supports StructInterface==1) 00069 ! EllipticObject%Solver - (StuctPCG=1, StructGMRES=2) 00070 ! EllipticObject%tolerance - (Solver tolerance) 00071 ! EllipticObject%MaxIters - (Maximum # of iterations) 00072 ! EllipticObject%hverbosity - (hypre verbosity) 00073 00074 TYPE(EllipticObjectDef) :: EllipticObject 00075 INTEGER :: iErr 00076 INTEGER :: i,j,m 00077 INTEGER :: solver ! type of HYPRE solver 00078 INTEGER :: printLevel ! verbosity 00079 INTEGER :: maxIters ! maximum number of iterations 00080 INTEGER :: hVerbosity ! how verbose user-side HYPRE should be 00081 REAL :: tolerance ! convergence tolerance 00082 INTEGER, DIMENSION(3,2) :: mthbc 00083 NAMELIST /PoissonData/ solver,tolerance,printlevel,MaxIters,hverbosity, mthbc 00084 00085 ! Read in scheme value. 00086 READ(GLOBAL_DATA_HANDLE,NML=PoissonData,IOStat=iErr) 00087 00088 IF(iErr/=0) THEN 00089 PRINT *, "PoissonInit() error: unable to read ", GLOBAL_DATA_FILE, "." 00090 STOP 00091 END IF 00092 PoissonScale=4d0*Pi*ScaleGrav 00093 EllipticObject%Interface=StructInterface 00094 EllipticObject%Solver=solver 00095 EllipticObject%tolerance=tolerance 00096 EllipticObject%PrintLevel=printlevel 00097 EllipticObject%MaxIters=MaxIters 00098 EllipticObject%hverbosity=hverbosity 00099 Poisson_mthbc=mthbc 00100 CALL Checkmthbc(mthbc) 00101 ! lPoissonPeriodic=ANY(mthbc(1:nDim,:)==2,2) 00102 lEllipticPeriodic(1:nDim)=lEllipticPeriodic(1:nDim) .OR. lPoissonPeriodic(1:nDim) !ANY(mthbc(1:nDim,:)==2,2) 00103 lGravityPeriodic=ANY(lPoissonPeriodic) 00104 ! write(*,*) 'lEllipticPeriodic=', lEllipticPeriodic 00105 ! write(*,*) mthbc, lReflect 00106 ! Sets up stencil offsets for poisson solve 00107 ALLOCATE(offsets(nDim,0:2*nDim)) 00108 ALLOCATE(StencilValues(0:2*nDim)) 00109 ALLOCATE(iEntries(1:2*nDim+1)) 00110 iEntries=(/(i,i=0,2*nDim)/) 00111 00112 CALL C_StructStencilCreate(ndim, 2*nDim+1, EllipticObject%stencil, iErr) ;CALL CheckErr('C_StructStencilCreate',iErr) 00113 offsets(:,0)=0 00114 stencilvalues(0)=-2*nDim 00115 m=0 00116 DO i=1,nDim 00117 DO j=-1,1,2 00118 m=m+1 00119 offsets(:,m)=0 00120 offsets(i,m)=j 00121 stencilvalues(m)=1 00122 END DO 00123 END DO 00124 DO m=0,2*nDim 00125 CALL C_StructStencilSetElement(EllipticObject%stencil, m,offsets(:,m), iErr) 00126 CALL CheckErr('C_StructStencilSetElement',iErr) 00127 END DO 00128 CALL InitMultiPoleMoments 00129 END SUBROUTINE Poisson_Init 00131 00134 00138 SUBROUTINE PoissonPreSolveComm(n,done_solving) 00139 INTEGER :: n 00140 LOGICAL :: done_solving 00141 ! IF (n == 0 .AND. ALL(lPoissonPeriodic(1:nDim))) CALL UpdateMeanDensity(n) 00142 done_solving=.false. 00143 END SUBROUTINE PoissonPreSolveComm 00144 00148 SUBROUTINE PoissonBetweenSolveComm(n,done_solving) 00149 INTEGER :: n 00150 LOGICAL :: done_solving 00151 TYPE(NodeDef), POINTER :: node, overlap 00152 TYPE(NodeDefList), POINTER :: nodelist, overlaplist 00153 done_solving=.true. 00154 END SUBROUTINE PoissonBetweenSolveComm 00155 00156 00159 SUBROUTINE PoissonPostSolveComm(n) 00160 INTEGER :: n 00161 TYPE(NodeDef), POINTER :: node, overlap 00162 TYPE(NodeDefList), POINTER :: nodelist, overlaplist 00163 00164 ! Need to ghost variables q(:,:,:,iPhi) and q(:,:,:,iPhiDot) into regions to calculate source terms. 00165 ! CALL EllipticTransfer(n,(/iPhiGas,iPhi,iPhiDot/),levels(n)%gmbc(levels(n)%step+1)+1) 00166 END SUBROUTINE PoissonPostSolveComm 00167 00169 00172 00176 SUBROUTINE Poisson_Setup(n,EllipticLevelObject) 00177 INTEGER :: n 00178 TYPE(EllipticLevelObjectDef) :: EllipticLevelObject 00179 TYPE(NodeDef), POINTER :: node 00180 TYPE(NodeDefList), POINTER :: nodelist 00181 nodelist=>nodes(n)%p 00182 DO WHILE (associated(nodelist)) 00183 node=>nodelist%self 00184 CALL PoissonMatrixSetBoxValues(n,node%info,EllipticLevelObject) 00185 nodelist=>nodelist%next 00186 END DO 00187 END SUBROUTINE Poisson_Setup 00188 00192 SUBROUTINE Poisson_LoadVectors(n,EllipticLevelObject) 00193 INTEGER :: n 00194 TYPE(EllipticLevelObjectDef) :: EllipticLevelObject 00195 TYPE(NodeDef), POINTER :: node 00196 TYPE(NodeDefList), POINTER :: nodelist 00197 nodelist=>nodes(n)%p 00198 DO WHILE (associated(nodelist)) 00199 node=>nodelist%self 00200 CALL PoissonVectorSetBoxValues(n,node%info,EllipticLevelObject) 00201 nodelist=>nodelist%next 00202 END DO 00203 END SUBROUTINE Poisson_LoadVectors 00204 00208 SUBROUTINE Poisson_UnLoadVectors(n,EllipticLevelObject) 00209 INTEGER :: n 00210 TYPE(EllipticLevelObjectDef) :: EllipticLevelObject 00211 TYPE(NodeDef), POINTER :: node 00212 TYPE(NodeDefList), POINTER :: nodelist 00213 nodelist=>nodes(n)%p 00214 DO WHILE (associated(nodelist)) 00215 node=>nodelist%self 00216 CALL PoissonVectorGetBoxValues(n,node%info,EllipticLevelObject) 00217 nodelist=>nodelist%next 00218 END DO 00219 ! write(*,*) 'solved for phigas and updated phidot and phi on level', n 00220 END SUBROUTINE Poisson_UnLoadVectors 00221 00222 00225 SUBROUTINE PoissonSetBC(n) 00226 INTEGER :: n, rmbc 00227 TYPE(NodeDef), POINTER :: node 00228 TYPE(NodeDefList), POINTER :: nodelist 00229 rmbc=levels(n)%egmbc(levels(n)%step) 00230 IF (n == 0 .AND. ANY(lMultiPole)) CALL CalcMultiPoleMoments() 00231 nodelist=>Nodes(n)%p 00232 DO WHILE (associated(nodelist)) 00233 node=>nodelist%self 00234 IF (ANY(lReflect(1:nDim))) CALL ApplyPoissonPhysicalBC(node%info, rmbc) 00235 IF (ANY(lMultiPole(1:nDim))) CALL ApplyAnalyticBC(node%info, rmbc) 00236 nodelist=>nodelist%next 00237 END DO 00238 END SUBROUTINE PoissonSetBC 00239 00242 SUBROUTINE PoissonPostElliptic(n) 00243 INTEGER :: n 00244 TYPE(NodeDef), POINTER :: node 00245 TYPE(NodeDefList), POINTER :: nodelist 00246 CALL EllipticTransfer(n, (/iPhiGas, iPhiDot/), levels(n)%egmbc(levels(n)%step)) 00247 00248 nodelist=>Nodes(n)%p 00249 DO WHILE (associated(nodelist)) 00250 node=>nodelist%self 00251 CALL PoissonGetMaxSpeed(node%info) 00252 nodelist=>nodelist%next 00253 END DO 00254 00255 ! write(*,*) 'updating mass fluxes with new potential on level', n 00256 00257 nodelist=>Nodes(n)%p 00258 DO WHILE (associated(nodelist)) 00259 node=>nodelist%self 00260 CALL SelfGravCorrection(node%info) 00261 nodelist=>nodelist%next 00262 END DO 00263 00264 00265 END SUBROUTINE PoissonPostElliptic 00266 00269 SUBROUTINE UpdatePhiAtBoundaries(n) 00270 INTEGER :: n 00271 TYPE(NodeDef), POINTER :: node 00272 TYPE(NodeDefList), POINTER :: nodelist 00273 nodelist=>Nodes(n)%p 00274 DO WHILE (ASSOCIATED(nodelist)) 00275 CALL ApplyAnalyticBC(nodelist%self%info, 1) 00276 nodelist=>nodelist%next 00277 END DO 00278 END SUBROUTINE UpdatePhiAtBoundaries 00279 00280 00282 SUBROUTINE SelfGravCorrection(info) 00283 TYPE(InfoDef) :: Info 00284 INTEGER, DIMENSION(3,2) :: mb ! bounds to update 00285 INTEGER :: i,j,k,l,m,mbc,gpmb(3,3,2) ! size of ghost regions 00286 REAL(KIND=qPREC) :: dx, iPhiSign(2),dtdx 00287 INTEGER, DIMENSION(2) :: iPhiIndex(2) 00288 REAL(KIND=qPREC), DIMENSION(:,:,:,:), ALLOCATABLE :: f2x_, f2y_, f2z_, phi 00289 REAL(KIND=qPREC), DIMENSION(:,:,:), ALLOCATABLE :: gradphix_, gradphiy_, gradphiz_ 00290 IF (levels(Info%level)%dt == 0) RETURN 00291 ! RETURN 00292 dx=levels(Info%level)%dx 00293 mb=1 00294 mbc=levels(Info%level)%ambc(levels(Info%level)%step) 00295 00296 !Expand bounds to include row of ghost zones for phi 00297 mb(1:nDim,1)=1-mbc-1 00298 mB(1:nDim,2)=Info%mx(1:nDim)+mbc+1 00299 00300 allocate(phi(mb(1,1):mb(1,2), mb(2,1):mb(2,2), mb(3,1):mb(3,2),2)) 00301 phi(:,:,:,2)=Info%q(mb(1,1):mb(1,2), mb(2,1):mb(2,2), mb(3,1):mb(3,2), iPhiGas) 00302 phi(:,:,:,1)=phi(:,:,:,2)-Info%q(mb(1,1):mb(1,2), mb(2,1):mb(2,2), mb(3,1):mb(3,2), iPhiDot)*levels(Info%level)%dt 00303 ! IF (Info%level == 0) THEN 00304 ! write(*,*) mb 00305 ! write(*,*) minval(phi), maxval(phi) 00306 ! END IF 00307 !Set bounds for grad phi 00308 gpmb=1 00309 DO i=1,nDim 00310 gpmb(i,1:nDim,1)=mb(1:nDim,1) 00311 gpmb(i,1:nDim,2)=mb(1:nDim,2) 00312 gpmb(i,i,1)=mb(i,1)+1 00313 END DO 00314 ! write(*,*) 'gpmb=', gpmb 00315 !Readjust bounds for region of q to update 00316 mb(1:nDim,1)=mB(1:nDim,1)+1 00317 mB(1:nDim,2)=mB(1:nDim,2)-1 00318 ! write(*,*) 'reset mb to', mB 00319 ALLOCATE(gradphix_(gpmb(1,1,1):gpmb(1,1,2), gpmb(1,2,1):gpmb(1,2,2), gpmb(1,3,1):gpmb(1,3,2))) 00320 ALLOCATE(f2x_(mb(1,1):mb(1,2)+1,mb(2,1):mb(2,2),mb(3,1):mb(3,2),1:nDim)) 00321 f2x_=0d0 00322 IF (nDim >= 2) THEN 00323 ALLOCATE(gradphiy_(gpmb(2,1,1):gpmb(2,1,2), gpmb(2,2,1):gpmb(2,2,2), gpmb(2,3,1):gpmb(2,3,2))) 00324 ALLOCATE(f2y_(mb(1,1):mb(1,2),mb(2,1):mb(2,2)+1,mb(3,1):mb(3,2),1:nDim)) 00325 f2y_=0d0 00326 IF (nDim >= 3) THEN 00327 ALLOCATE(gradphiz_(gpmb(3,1,1):gpmb(3,1,2), gpmb(3,2,1):gpmb(3,2,2), gpmb(3,3,1):gpmb(3,3,2))) 00328 ALLOCATE(f2z_(mb(1,1):mb(1,2),mb(2,1):mb(2,2),mb(3,1):mb(3,2)+1,1:nDim)) 00329 f2z_=0d0 00330 END IF 00331 END IF 00332 00333 iPhiSign=(/-1d0,1d0/) 00334 00335 DO l=1,2 00336 FORALL (i=gpmb(1,1,1):gpmb(1,1,2), j=gpmb(1,2,1):gpmb(1,2,2), k=gpmb(1,3,1):gpmb(1,3,2)) 00337 gradphix_(i,j,k)=(phi(i,j,k,l)-phi(i-1,j,k,l))/dx 00338 END FORALL 00339 IF (nDim >= 2) THEN 00340 FORALL (i=gpmb(2,1,1):gpmb(2,1,2), j=gpmb(2,2,1):gpmb(2,2,2), k=gpmb(2,3,1):gpmb(2,3,2)) 00341 gradphiy_(i,j,k)=(phi(i,j,k,l)-phi(i,j-1,k,l))/dx 00342 END FORALL 00343 IF (nDim >= 3) THEN 00344 FORALL (i=gpmb(3,1,1):gpmb(3,1,2), j=gpmb(3,2,1):gpmb(3,2,2), k=gpmb(3,3,1):gpmb(3,3,2)) 00345 gradphiz_(i,j,k)=(phi(i,j,k,l)-phi(i,j,k-1,l))/dx 00346 END FORALL 00347 END IF 00348 END IF 00349 dtdx=half*levels(Info%level)%dt/dx*iPhiSign(l) 00350 IF (nDim == 1) THEN 00351 mB(1,2)=mb(1,2)+1 00352 FORALL (i=mb(1,1):mb(1,2), j=mb(2,1):mb(2,2), k=mb(3,1):mb(3,2)) 00353 f2x_(i,j,k,1) = f2x_(i,j,k,1) + & 00354 dtdx*(.125d0/Pi/ScaleGrav*gradphix_(i,j,k)**2+half*mean_density*(phi(i,j,k,l)+phi(i-1,j,k,l))) 00355 END FORALL 00356 mB(1,2)=mb(1,2)-1 00357 ELSEIF (nDim == 2) THEN 00358 mB(1,2)=mb(1,2)+1 00359 FORALL (i=mb(1,1):mb(1,2), j=mb(2,1):mb(2,2), k=mb(3,1):mb(3,2)) 00360 f2x_(i,j,k,1) = f2x_(i,j,k,1) + dtdx*(.125d0/Pi/ScaleGrav*& 00361 (gradphix_(i,j,k)**2-(.25d0*sum(gradphiy_(i-1:i,j:j+1,k)))**2)+& 00362 half*mean_density*SUM(phi(i-1:i,j,k,l))) 00363 f2x_(i,j,k,2) = f2x_(i,j,k,2) + dtdx*(.25d0/Pi/ScaleGrav*& 00364 (gradphix_(i,j,k)*.25d0*sum(gradphiy_(i-1:i,j:j+1,k)))) 00365 END FORALL 00366 mB(1,2)=mb(1,2)-1 00367 mB(2,2)=mb(2,2)+1 00368 FORALL (i=mb(1,1):mb(1,2), j=mb(2,1):mb(2,2), k=mb(3,1):mb(3,2)) 00369 f2y_(i,j,k,2) = f2y_(i,j,k,2) + dtdx*(.125d0/Pi/ScaleGrav*& 00370 (gradphiy_(i,j,k)**2-(.25d0*sum(gradphix_(i:i+1,j-1:j,k)))**2)+& 00371 half*mean_density*SUM(phi(i,j-1:j,k,l))) 00372 f2y_(i,j,k,1) = f2y_(i,j,k,1) + dtdx*(.25d0/Pi/ScaleGrav*& 00373 (gradphiy_(i,j,k)*.25d0*sum(gradphix_(i:i+1,j-1:j,k)))) 00374 END FORALL 00375 mB(2,2)=mb(2,2)-1 00376 ELSE! nDim == 3 00377 mB(1,2)=mb(1,2)+1 00378 FORALL (i=mb(1,1):mb(1,2), j=mb(2,1):mb(2,2), k=mb(3,1):mb(3,2)) 00379 f2x_(i,j,k,1) = f2x_(i,j,k,1) + dtdx*(.125d0/Pi/ScaleGrav*& 00380 (gradphix_(i,j,k)**2-(.25d0*sum(gradphiy_(i-1:i,j:j+1,k)))**2-(.25d0*sum(gradphiz_(i-1:i,j,k:k+1)))**2)+& 00381 half*mean_density*SUM(phi(i-1:i,j,k,l))) 00382 f2x_(i,j,k,2) = f2x_(i,j,k,2) + dtdx*(.25d0/Pi/ScaleGrav*& 00383 (gradphix_(i,j,k)*.25d0*sum(gradphiy_(i-1:i,j:j+1,k)))) 00384 f2x_(i,j,k,3) = f2x_(i,j,k,3) + dtdx*(.25d0/Pi/ScaleGrav*& 00385 (gradphix_(i,j,k)*.25d0*sum(gradphiz_(i-1:i,j,k:k+1)))) 00386 END FORALL 00387 mB(1,2)=mb(1,2)-1 00388 mB(2,2)=mb(2,2)+1 00389 FORALL (i=mb(1,1):mb(1,2), j=mb(2,1):mb(2,2), k=mb(3,1):mb(3,2)) 00390 f2y_(i,j,k,2) = f2y_(i,j,k,2) + dtdx*(.125d0/Pi/ScaleGrav*& 00391 (gradphiy_(i,j,k)**2-(.25d0*sum(gradphix_(i:i+1,j-1:j,k)))**2-(.25d0*sum(gradphiz_(i,j-1:j,k:k+1)))**2)+& 00392 half*mean_density*SUM(phi(i,j-1:j,k,l))) 00393 f2y_(i,j,k,1) = f2y_(i,j,k,1) + dtdx*(.25d0/Pi/ScaleGrav*& 00394 (gradphiy_(i,j,k)*.25d0*sum(gradphix_(i:i+1,j-1:j,k)))) 00395 f2y_(i,j,k,3) = f2y_(i,j,k,3) + dtdx*(.25d0/Pi/ScaleGrav*& 00396 (gradphiy_(i,j,k)*.25d0*sum(gradphiz_(i,j-1:j,k:k+1)))) 00397 END FORALL 00398 mB(2,2)=mb(2,2)-1 00399 mB(3,2)=mb(3,2)+1 00400 FORALL (i=mb(1,1):mb(1,2), j=mb(2,1):mb(2,2), k=mb(3,1):mb(3,2)) 00401 f2z_(i,j,k,3) = f2z_(i,j,k,3) + dtdx*(.125d0/Pi/ScaleGrav*& 00402 (gradphiz_(i,j,k)**2-(.25d0*sum(gradphix_(i:i+1,j,k-1:k)))**2-(.25d0*sum(gradphiy_(i,j:j+1,k-1:k)))**2)+& 00403 half*mean_density*SUM(phi(i,j,k-1:k,l))) 00404 f2z_(i,j,k,1) = f2z_(i,j,k,1) + dtdx*(.25d0/Pi/ScaleGrav*& 00405 (gradphiz_(i,j,k)*.25d0*sum(gradphix_(i:i+1,j,k-1:k)))) 00406 f2z_(i,j,k,2) = f2z_(i,j,k,2) + dtdx*(.25d0/Pi/ScaleGrav*& 00407 (gradphiz_(i,j,k)*.25d0*sum(gradphiy_(i,j:j+1,k-1:k)))) 00408 END FORALL 00409 mB(3,2)=mb(3,2)-1 00410 END IF 00411 END DO 00412 FORALL(i=mB(1,1):mB(1,2),j=mb(2,1):mB(2,2),k=mB(3,1):mB(3,2),m=1:nDim) 00413 Info%q(i,j,k,imom(m))=Info%q(i,j,k,imom(m))+f2x_(i,j,k,m)-f2x_(i+1,j,k,m) 00414 END FORALL 00415 mB(1,2)=mB(1,2)+1 00416 CALL storefixupfluxes(Info,mB,1,f2x_,imom(1:nDim)) 00417 mb(1,2)=mb(1,2)-1 00418 DEALLOCATE(f2x_, gradphix_) 00419 IF (nDim >= 2) THEN 00420 FORALL(i=mB(1,1):mB(1,2),j=mb(2,1):mB(2,2),k=mB(3,1):mB(3,2),m=1:nDim) 00421 Info%q(i,j,k,imom(m))=Info%q(i,j,k,imom(m))+f2y_(i,j,k,m)-f2y_(i,j+1,k,m) 00422 END FORALL 00423 mB(2,2)=mB(2,2)+1 00424 CALL storefixupfluxes(Info,mB,2,f2y_,imom(1:nDim)) 00425 mB(2,2)=mB(2,2)-1 00426 DEALLOCATE(f2y_, gradphiy_) 00427 IF (nDim >= 3) THEN 00428 FORALL(i=mB(1,1):mB(1,2),j=mb(2,1):mB(2,2),k=mB(3,1):mB(3,2),m=1:nDim) 00429 Info%q(i,j,k,imom(m))=Info%q(i,j,k,imom(m))+f2z_(i,j,k,m)-f2z_(i,j,k+1,m) 00430 END FORALL 00431 mB(3,2)=mB(3,2)+1 00432 CALL storefixupfluxes(Info,mB,3,f2z_,imom(1:nDim)) 00433 mb(3,2)=mb(3,2)-1 00434 DEALLOCATE(f2z_, gradphiz_) 00435 END IF 00436 END IF 00437 00438 IF (iE /= 0) THEN 00439 FORALL(i=mb(1,1):mb(1,2), j=mb(2,1):mb(2,2), k=mb(3,1):mb(3,2)) 00440 Info%q(i,j,k,iE)=Info%q(i,j,k,iE)-.25d0*( & 00441 Info%MassFlux(i,j,k,1)*(phi(i,j,k,2)-phi(i,j,k,1)-phi(i-1,j,k,2)+phi(i-1,j,k,1)) + & 00442 Info%MassFlux(i+1,j,k,1)*(phi(i+1,j,k,2)-phi(i+1,j,k,1)-phi(i,j,k,2)+phi(i,j,k,1))) 00443 END FORALL 00444 IF (nDim >= 2) THEN 00445 FORALL(i=mb(1,1):mb(1,2), j=mb(2,1):mb(2,2), k=mb(3,1):mb(3,2)) 00446 Info%q(i,j,k,iE)=Info%q(i,j,k,iE)-.25d0*( & 00447 Info%MassFlux(i,j,k,2)*(phi(i,j,k,2)-phi(i,j,k,1)-phi(i,j-1,k,2)+phi(i,j-1,k,1)) + & 00448 Info%MassFlux(i,j+1,k,2)*(phi(i,j+1,k,2)-phi(i,j+1,k,1)-phi(i,j,k,2)+phi(i,j,k,1))) 00449 END FORALL 00450 IF (nDim == 3) THEN 00451 FORALL(i=mb(1,1):mb(1,2), j=mb(2,1):mb(2,2), k=mb(3,1):mb(3,2)) 00452 Info%q(i,j,k,iE)=Info%q(i,j,k,iE)-.25d0*( & 00453 Info%MassFlux(i,j,k,3)*(phi(i,j,k,2)-phi(i,j,k,1)-phi(i,j,k-1,2)+phi(i,j,k-1,1)) + & 00454 Info%MassFlux(i,j,k+1,3)*(phi(i,j,k+1,2)-phi(i,j,k+1,1)-phi(i,j,k,2)+phi(i,j,k,1))) 00455 END FORALL 00456 END IF 00457 END IF 00458 END IF 00459 DEALLOCATE(phi) 00460 END SUBROUTINE SelfGravCorrection 00462 00463 00466 00471 SUBROUTINE PoissonMatrixSetBoxValues(n,Info,EllipticLevelObject) 00472 INTEGER :: n ! Level 00473 TYPE(InfoDef) :: Info 00474 TYPE(EllipticLevelObjectDef) :: EllipticLevelObject 00475 00476 INTEGER :: nvalues,nCells 00477 REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: matrixvalues 00478 INTEGER :: i,j,k,m,edge,p,temp,ip(3),l 00479 LOGICAL, DIMENSION(3,2) :: internal_to_grid, internal_to_domain 00480 INTEGER :: ierr 00481 nvalues=size(offsets,2) 00482 nCells=product(Info%mX(1:nDim)) 00483 ALLOCATE(matrixvalues(nCells*nvalues)) 00484 00485 matrixvalues=RESHAPE(SPREAD(StencilValues,2,nCells),shape(matrixvalues)) !default values 00486 m=1-nvalues 00487 DO k=Info%mGlobal(3,1),Info%mGlobal(3,2) 00488 internal_to_grid(3,1) = (k > Info%mGlobal(3,1)) 00489 internal_to_grid(3,2) = (k < Info%mGlobal(3,2)) 00490 IF (.NOT. internal_to_grid(3,1)) internal_to_domain(3,1) = lPoissonPeriodic(3) .OR. k > 1 00491 IF (.NOT. internal_to_grid(3,2)) internal_to_domain(3,2) = lPoissonPeriodic(3) .OR. k < levels(n)%mX(3) 00492 00493 DO j=Info%mGlobal(2,1),Info%mGlobal(2,2) 00494 internal_to_grid(2,1) = (j > Info%mGlobal(2,1)) 00495 internal_to_grid(2,2) = (j < Info%mGlobal(2,2)) 00496 IF (.NOT. internal_to_grid(2,1)) internal_to_domain(2,1) = lPoissonPeriodic(2) .OR. j > 1 00497 IF (.NOT. internal_to_grid(2,2)) internal_to_domain(2,2) = lPoissonPeriodic(2) .OR. j < levels(n)%mX(2) 00498 00499 DO i=Info%mGlobal(1,1), Info%mGlobal(1,2) 00500 internal_to_grid(1,1) = (i > Info%mGlobal(1,1)) 00501 internal_to_grid(1,2) = (i < Info%mGlobal(1,2)) 00502 IF (.NOT. internal_to_grid(1,1)) internal_to_domain(1,1) = lPoissonPeriodic(1) .OR. i > 1 00503 IF (.NOT. internal_to_grid(1,2)) internal_to_domain(1,2) = lPoissonPeriodic(1) .OR. i < levels(n)%mX(1) 00504 m=m+nvalues 00505 IF (ALL(internal_to_grid)) CYCLE 00506 ip=(/i,j,k/)-Info%mGlobal(1:3,1)+1 00507 p=0 00508 DO l=1,nDim 00509 temp=ip(l) 00510 DO edge=1,2 00511 p=p+1 00512 IF (internal_to_grid(l,edge)) CYCLE 00513 IF (internal_to_domain(l,edge)) THEN !zero out connections to ancestors 00514 ip(l)=temp+(-1)**edge 00515 IF (isAncestor(Info%ChildMask(ip(1),ip(2),ip(3)))) matrixvalues(m+p) = 0 !zero out matrix connections to other levels 00516 00517 ELSE !physical boundary condition exists 00518 matrixvalues(m+p) = 0 00519 SELECT CASE (Poisson_mthbc(l,edge)) 00520 CASE(MULTIPOLE_EXPANSION) 00521 !Do nothing - handled by source terms 00522 CASE(REFLECTING) 00523 !matrixvalues(m+p-(-1)**edge)=matrixvalues(m+p-(-1)**edge)+StencilValues(p) 00524 matrixvalues(m)=matrixvalues(m)+StencilValues(p) 00525 CASE(ZEROSLOPE) 00526 matrixvalues(m)=matrixvalues(m)-StencilValues(p) 00527 END SELECT 00528 END IF 00529 END DO 00530 ip(l)=temp 00531 END DO 00532 END DO 00533 END DO 00534 END DO 00535 CALL C_StructMatrixSetBoxValues(EllipticLevelObject%Matrix, Info%mGlobal(:,1), Info%mGlobal(:,2), nvalues, iEntries, matrixValues, iErr) 00536 CALL CheckErr('C_StructMatrixSetBoxValues',iErr) 00537 DEALLOCATE(matrixValues) 00538 END SUBROUTINE PoissonMatrixSetBoxValues 00539 00545 SUBROUTINE PoissonVectorSetBoxValues(n,Info,EllipticLevelObject) 00546 USE EOS 00547 INTEGER :: n 00548 TYPE(InfoDef) :: Info 00549 TYPE(EllipticLevelObjectDef) :: EllipticLevelObject 00550 00551 INTEGER ::nCells 00552 REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: vectorvalues 00553 INTEGER :: i,j,k,m,edge,p,temp,mb(3,2),ip(3),l,mbc 00554 LOGICAL, DIMENSION(3,2) :: internal_to_grid, internal_to_domain 00555 INTEGER :: iErr 00556 REAL(KIND=qPREC), DIMENSION(3) :: pos 00557 REAL(KIND=qPREC), DIMENSION(:,:,:), ALLOCATABLE :: oldphi 00558 mbc=levels(Info%level)%egmbc(levels(Info%level)%step) 00559 mb(nDim+1:3,:)=1 00560 mb(1:nDim,1)=1-mbc 00561 mb(1:nDim,2)=Info%mX(1:nDim)+mbc 00562 CALL Protectq(Info, mb, 'setboxvalues') 00563 00564 nCells=product(Info%mX(1:nDim)) 00565 !Update phiGas with predictor and store old phigas in phidot 00566 IF (levels(n)%dt > 0) THEN 00567 ALLOCATE(OldPhi(Info%mX(1),Info%mX(2),Info%mX(3))) 00568 Oldphi(1:Info%mX(1),1:Info%mX(2), 1:Info%mX(3))= & 00569 Info%q(1:Info%mX(1),1:Info%MX(2),1:Info%mX(3),iPhiGas) 00570 Info%q(mb(1,1):mb(1,2),mb(2,1):mb(2,2),mb(3,1):mb(3,2),iPhiGas)=& 00571 Info%q(mb(1,1):mb(1,2),mb(2,1):mb(2,2),mb(3,1):mb(3,2),iPhiGas)+& 00572 levels(n)%dt*Info%q(mb(1,1):mb(1,2),mb(2,1):mb(2,2),mb(3,1):mb(3,2),iPhiDot) 00573 Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),iPhiDot)=Oldphi 00574 DEALLOCATE(oldphi) 00575 END IF 00576 00577 ALLOCATE(vectorvalues(nCells)) 00578 00579 IF (ALL(lPoissonPeriodic(1:nDim))) THEN 00580 vectorvalues=reshape(Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),irho)-mean_density,shape(vectorvalues))*levels(n)%dx**2*PoissonScale 00581 00582 ELSEIF (ALL(lMultipole(1:nDim))) THEN 00583 m=0 00584 DO k=1,Info%mX(3) 00585 DO j=1,Info%mX(2) 00586 DO i=1,Info%mX(1) 00587 m=m+1 00588 pos=Info%xbounds(:,1)+(REAL((/i,j,k/))-half)*levels(Info%level)%dX - Multipole_COM(:) 00589 IF (SUM(pos(1:nDim)**2) < multipole_radius**2) THEN 00590 vectorvalues(m)=Info%q(i,j,k,1)*levels(n)%dx**2*PoissonScale 00591 ELSE 00592 vectorvalues(m)=0 00593 END IF 00594 END DO 00595 END DO 00596 END DO 00597 ELSE 00598 vectorvalues=reshape(Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),irho),shape(vectorvalues))*levels(n)%dx**2*PoissonScale 00599 END IF 00600 00601 m=0 00602 DO k=Info%mGlobal(3,1),Info%mGlobal(3,2) 00603 internal_to_grid(3,1) = (k > Info%mGlobal(3,1)) 00604 internal_to_grid(3,2) = (k < Info%mGlobal(3,2)) 00605 IF (.NOT. internal_to_grid(3,1)) internal_to_domain(3,1) = lPoissonPeriodic(3) .OR. k > 1 00606 IF (.NOT. internal_to_grid(3,2)) internal_to_domain(3,2) = lPoissonPeriodic(3) .OR. k < levels(n)%mX(3) 00607 00608 DO j=Info%mGlobal(2,1),Info%mGlobal(2,2) 00609 internal_to_grid(2,1) = (j > Info%mGlobal(2,1)) 00610 internal_to_grid(2,2) = (j < Info%mGlobal(2,2)) 00611 IF (.NOT. internal_to_grid(2,1)) internal_to_domain(2,1) = lPoissonPeriodic(2) .OR. j > 1 00612 IF (.NOT. internal_to_grid(2,2)) internal_to_domain(2,2) = lPoissonPeriodic(2) .OR. j < levels(n)%mX(2) 00613 00614 DO i=Info%mGlobal(1,1), Info%mGlobal(1,2) 00615 internal_to_grid(1,1) = (i > Info%mGlobal(1,1)) 00616 internal_to_grid(1,2) = (i < Info%mGlobal(1,2)) 00617 IF (.NOT. internal_to_grid(1,1)) internal_to_domain(1,1) = lPoissonPeriodic(1) .OR. i > 1 00618 IF (.NOT. internal_to_grid(1,2)) internal_to_domain(1,2) = lPoissonPeriodic(1) .OR. i < levels(n)%mX(1) 00619 m=m+1 00620 IF (ALL(internal_to_grid)) CYCLE 00621 ip=(/i,j,k/)-Info%mGlobal(:,1)+1 00622 p=0 00623 DO l=1,nDim 00624 temp=ip(l) 00625 DO edge=1,2 00626 p=p+1 00627 IF (internal_to_grid(l,edge)) CYCLE 00628 ip(l)=temp+(-1)**edge 00629 IF ((internal_to_domain(l,edge) .AND. isAncestor(Info%childmask(ip(1),ip(2),ip(3)))) .OR. ((.NOT. internal_to_domain(l,edge)) .AND. (Poisson_mthbc(l,edge) == USERSPECIFIED .OR. Poisson_mthbc(l,edge) == MULTIPOLE_EXPANSION))) THEN !zero out connections to ancestors 00630 vectorvalues(m) = vectorvalues(m)-StencilValues(p)*Info%q(ip(1),ip(2),ip(3),iPhiGas) !add boundary term 00631 END IF 00632 END DO 00633 ip(l)=temp 00634 END DO 00635 END DO 00636 END DO 00637 END DO 00638 CALL C_StructVectorSetBoxValues(EllipticLevelObject%VariableVector, Info%mGlobal(:,1), Info%mGlobal(:,2),vectorValues, iErr) 00639 CALL CheckErr('C_StructVectorSetBoxValues',iErr) 00640 vectorvalues=reshape(Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),iPhiGas),shape(vectorvalues)) 00641 CALL C_StructVectorSetBoxValues(EllipticLevelObject%SolutionVector, Info%mGlobal(:,1), Info%mGlobal(:,2), vectorvalues, iErr) 00642 CALL CheckErr('C_StructVectorSetBoxValues',iErr) 00643 00644 DEALLOCATE(vectorValues) 00645 00646 END SUBROUTINE PoissonVectorSetBoxValues 00647 00653 SUBROUTINE PoissonVectorGetBoxValues(n,Info,EllipticLevelObject) 00654 INTEGER :: n 00655 TYPE(InfoDef) :: Info 00656 TYPE(EllipticLevelObjectDef) :: EllipticLevelObject 00657 INTEGER ::nCells, mb(3,2), mbc 00658 REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: vectorvalues 00659 INTEGER :: iErr 00660 ! RETURN 00661 nCells=product(Info%mX(1:nDim)) 00662 ALLOCATE(vectorvalues(nCells)) 00663 00664 CALL C_StructVectorGetBoxValues(EllipticLevelObject%SolutionVector, Info%mGlobal(:,1), Info%mGlobal(:,2),vectorValues, iErr) 00665 CALL CheckErr('C_StructVectorGetBoxValues',iErr) 00666 00667 mbc=levels(Info%level)%egmbc(levels(Info%level)%step) 00668 mb(nDim+1:3,:)=1 00669 mb(1:nDim,1)=1-mbc 00670 mb(1:nDim,2)=Info%mX(1:nDim)+mbc 00671 !Update iPhiGas 00672 Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),iPhiGas)=reshape(vectorvalues,Info%mX) 00673 IF (levels(n)%dt > 0) THEN !predictive phidot 00674 !Update iPhiDot using old Phi stored in iPhiDot and new phi 00675 Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),iPhiDot)=& 00676 (Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),iPhiGas) - & 00677 Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),iPhiDot))/levels(n)%dt 00678 ELSE 00679 ! Info%q(:,:,:,iPhiDot)=0 00680 END IF 00681 DEALLOCATE(vectorValues) 00682 00683 END SUBROUTINE PoissonVectorGetBoxValues 00684 00685 00688 SUBROUTINE PoissonGetMaxSpeed(Info) 00689 USE ParticleDeclarations 00690 USE ParticleControl 00691 TYPE(InfoDef) :: Info 00692 INTEGER :: level 00693 INTEGER, DIMENSION(3,2) :: ip,ir,il 00694 INTEGER :: m 00695 level=Info%level 00696 ip(nDim+1:3,:)=1 00697 ip(1:nDim,1)=1!-levels(level)%gmbc(levels(level)%step) 00698 ip(1:nDim,2)=Info%mX(1:nDim)!+levels(level)%gmbc(levels(level)%step) 00699 ir=ip 00700 il=ip 00701 DO m=1,nDim 00702 ir(m,1)=ip(m,1)+1 00703 il(m,2)=ip(m,2)-1 00704 elliptic_maxspeed(level)=max(elliptic_maxspeed(level),sqrt(nDim*& 00705 maxval(abs(Info%q(ir(1,1):ir(1,2),ir(2,1):ir(2,2),ir(3,1):ir(3,2),iPhiGas)-& 00706 Info%q(il(1,1):il(1,2),il(2,1):il(2,2),il(3,1):il(3,2),iPhiGas))))) 00707 ir(m,:)=ip(m,:) 00708 il(m,:)=ip(m,:) 00709 END DO 00710 !~ write(*,*) 'elliptic_maxspeed=', elliptic_maxspeed(level) 00711 END SUBROUTINE PoissonGetMaxSpeed 00712 00713 00716 SUBROUTINE ApplyAnalyticBC(Info, rmbc) 00717 ! Interface declarations 00718 TYPE (InfoDef) :: Info 00719 ! Internal declarations 00720 INTEGER rmbc, level,dim,edge,start,ncells,i,dir 00721 INTEGER, DIMENSION(3,2):: lGmGlobal, ip 00722 level=Info%level 00723 ! rmbc=levels(level)%egmbc(1)!levels(level)%step) 00724 lGmGlobal(:,1)=GmGlobal(:,1) 00725 lGmGlobal(:,2)=GmGlobal(:,2)*PRODUCT(levels(0:level-1)%CoarsenRatio) 00726 ip=1 00727 DO dim=1,nDim 00728 DO edge=1,2 00729 IF (Poisson_mthbc(dim,edge) == MULTIPOLE_EXPANSION) THEN 00730 IF (edge == 1) THEN 00731 start=lGmGlobal(dim,1)-Info%mGlobal(dim,1) !first cell on left boundary 00732 nCells=start-(1-rmbc)+1 00733 ELSE 00734 start=(lGmGlobal(dim,2)+1)-(Info%mGlobal(dim,1)-1) !first ghost cell on right boundary 00735 nCells=Info%mx(dim)+rmbc-start+1 00736 END IF 00737 DO i=1,nDim 00738 IF (i==dim) CYCLE 00739 ip(i,1)=1-rmbc 00740 ip(i,2)=Info%mX(i)+rmbc 00741 END DO 00742 dir=(-1)**edge !(edge*2-3) !direction of edge (1=>-1, 2=>1) 00743 ip(dim,edge)=start+dir*(nCells-1) 00744 ip(dim,3-edge)=start 00745 IF (nCells > 0) CALL SetPhi(Info, ip) 00746 END IF 00747 END DO 00748 END DO 00749 END SUBROUTINE ApplyAnalyticBC 00750 00753 SUBROUTINE ApplyPoissonPhysicalBC(Info,rmbc) 00754 ! Interface declarations 00755 TYPE (InfoDef) :: Info 00756 ! Internal declarations 00757 INTEGER rmbc, level,dim,edge,start,ncells 00758 INTEGER, DIMENSION(3,2):: lGmGlobal 00759 level=Info%level 00760 lGmGlobal(:,1)=GmGlobal(:,1) 00761 lGmGlobal(:,2)=GmGlobal(:,2)*PRODUCT(levels(0:level-1)%CoarsenRatio) 00762 DO dim=1,nDim 00763 DO edge=1,2 00764 IF (Poisson_mthbc(dim,edge)==REFLECTING) THEN 00765 IF (edge == 1) THEN 00766 start=lGmGlobal(dim,1)-Info%mGlobal(dim,1) !first cell on left boundary 00767 nCells=start-(1-rmbc)+1 00768 ELSE 00769 start=(lGmGlobal(dim,2)+1)-(Info%mGlobal(dim,1)-1) !first ghost cell on right boundary 00770 nCells=Info%mx(dim)+rmbc-start+1 00771 END IF 00772 IF (nCells > 0) CALL MirrorPhiGhost(dim,edge,start,nCells,Info) 00773 END IF 00774 END DO 00775 END DO 00776 END SUBROUTINE ApplyPoissonPhysicalBC 00777 00779 00782 00783 ! Consistency Check on Poisson boundary conditions 00784 SUBROUTINE Checkmthbc(mthbc) 00785 INTEGER, DIMENSION(:,:) :: mthbc 00786 INTEGER :: dim,edge 00787 lMultiPole=.false. 00788 lReflect=.false. 00789 00790 IF (ANY(mthbc(1:nDim,:)==MULTIPOLE_EXPANSION)) THEN 00791 IF (ANY(ALL(mthbc(1:nDim,:)/=MULTIPOLE_EXPANSION, 2))) THEN 00792 PRINT*, 'must have at least one boundary in each dimension being multipole if any are multipole' 00793 STOP 00794 END IF 00795 END IF 00796 DO dim=1,nDim 00797 DO edge=1,2 00798 lReflect(dim) = lReflect(dim) .OR. mthbc(dim,edge) == REFLECTING !ANY(Gmthbc(dim,edge) == (/REFLECT_WALL, REFLECT_BPARALLEL, REFLECT_CYLINDRICAL/)) 00799 lMultiPole(dim) = lMultiPole(dim) .OR. mthbc(dim,edge) == MULTIPOLE_EXPANSION 00800 END DO 00801 END DO 00802 lPoissonPeriodic(1:nDim)=ANY(mthbc(1:nDim,:)==PERIODIC,2) 00803 IF (ALL(lPoissonPeriodic(1:nDim))) lNeedMeanDensity=.true. 00804 IF (ANY(lPoissonPeriodic(1:nDim) .AND. (mthbc(1:nDim,1) /= mthbc(1:nDim,2)))) THEN 00805 write(*,*) 'Boundaries must be periodic on both sides of any given direction' 00806 STOP 00807 END IF 00808 ! write(*,*) lMultiPole 00809 END SUBROUTINE Checkmthbc 00810 00811 SUBROUTINE MirrorPhiGhost(dim, edge, start, nCells, Info) 00812 INTEGER :: dim, start, nCells, edge,dir, level, rmbc,i,j 00813 TYPE(InfoDef) :: Info 00814 INTEGER, DIMENSION(3,2) :: ip,iq,ir,is 00815 LOGICAL :: lFirst 00816 INTEGER, DIMENSION(2) :: AuxParFields 00817 REAL(KIND=qPREC) :: aux_psign, aux_nsign 00818 INTEGER, DIMENSION(4) :: ReflectVars 00819 INTEGER :: nReflect, MirrorFields(3), nmirror 00820 IF (nCells == 0) RETURN 00821 nMirror=2 00822 MirrorFields(1:2)=(/iPhiGas, iPhiDot/) 00823 level=Info%level 00824 dir=(-1)**edge !(edge*2-3) !direction of edge (1=>-1, 2=>1) 00825 ip(nDim+1:3,:)=1 00826 iq(nDim+1:3,:)=1 00827 ! Stretch bounds by nCells 00828 DO i=1,nDim 00829 IF (i==dim) CYCLE 00830 ip(i,1)=1-nCells 00831 ip(i,2)=Info%mX(i)+nCells 00832 iq(i,:)=ip(i,:) 00833 END DO 00834 DO j=1,nCells 00835 iq(dim,:)=start+dir*(j-1) 00836 ip(dim,:)=start-dir*(j) 00837 Info%q(iq(1,1):iq(1,2),iq(2,1):iq(2,2),iq(3,1):iq(3,2),MirrorFields(1:nMirror))=& 00838 Info%q(ip(1,1):ip(1,2),ip(2,1):ip(2,2),ip(3,1):ip(3,2),MirrorFields(1:nMirror)) 00839 END DO 00840 END SUBROUTINE MirrorPhiGhost 00841 00842 00843 SUBROUTINE SelfGravCorrections(n) 00844 INTEGER :: n 00845 TYPE(NodeDef), POINTER :: node 00846 TYPE(NodeDefList), POINTER :: nodelist 00847 TYPE(InfoDef), POINTER :: info 00848 nodelist=>Nodes(n)%p 00849 DO WHILE(ASSOCIATED(nodelist)) 00850 info=>nodelist%self%info 00851 CALL SelfGravCorrection(info) 00852 nodelist=>nodelist%next 00853 END DO 00854 END SUBROUTINE SelfGravCorrections 00855 00857 00858 00859 SUBROUTINE FFTSOLVE(Info) 00860 TYPE(InfoDef) :: Info 00861 COMPLEX(8), dimension(:,:,:), allocatable :: in ! maps to transform of density 00862 COMPLEX(8), dimension(:,:,:), allocatable :: phi ! 00863 INTEGER :: i,j,k,frame,ikx,iky,ikz 00864 INTEGER(8) :: p1, p2 00865 REAL, DIMENSION(:), ALLOCATABLE :: kx,ky,kz 00866 REAL :: C1,dk(3) 00867 ALLOCATE(in(Info%mX(1),Info%mX(2),Info%mX(3))) 00868 ALLOCATE(phi(Info%mx(1),Info%mX(2),Info%mX(3))) 00869 C1=4d0*Pi*ScaleGrav/(PRODUCT(Info%mX(:))) !extra scaling factors (scale_grav is just G_ scaled to computational units) 00870 dk(:)=2d0*Pi/(GXBounds(1:3,2)-GxBounds(1:3,1))/(sqrt(C1)) 00871 CALL dfftw_plan_dft_3d(p1, Info%mx(1), Info%mx(2), Info%mx(3), in, phi, FFTW_FORWARD, FFTW_ESTIMATE) 00872 CALL dfftw_plan_dft_3d(p2, Info%mx(1), Info%mx(2), Info%mx(3), phi, in, FFTW_BACKWARD, FFTW_ESTIMATE) 00873 00874 in=Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),1) 00875 CALL dfftw_execute(p1) 00876 00877 ALLOCATE (kx(Info%mX(1)),ky(Info%mX(2)),kz(Info%mX(3))) 00878 DO ikx = 0, Info%mX(1)-1 00879 kx(ikx+1)=dk(1)*REAL((mod((ikx+Info%mX(1)/2),Info%mX(1))-Info%mX(1)/2),8) 00880 END DO 00881 DO iky = 0, Info%mX(2)-1 00882 ky(iky+1)=dk(2)*REAL((mod((iky+Info%mX(2)/2),Info%mX(2))-Info%mX(2)/2),8) 00883 END DO 00884 DO ikz = 0, Info%mX(3)-1 00885 kz(ikz+1)=dk(3)*REAL((mod((ikz+Info%mX(3)/2),Info%mX(3))-Info%mX(3)/2),8) 00886 END DO 00887 DO i=1,Info%mX(1) 00888 DO j=1,Info%mX(2) 00889 DO k=1,Info%mX(3) 00890 phi(i,j,k)=-1d0*phi(i,j,k)/REAL(kx(i)**2+ky(j)**2+kz(k)**2) 00891 END DO 00892 END DO 00893 END DO 00894 phi(1,1,1)=0d0 !periodic bc's 00895 CALL dfftw_execute(p2) 00896 Info%q(1:Info%mX(1),1:Info%mX(2),1:Info%mX(3),iPhiGas)=REAL(in) 00897 END SUBROUTINE FFTSOLVE 00898 00899 END MODULE Poisson