Scrambler  1
PFFT.f90
Go to the documentation of this file.
00001 !#########################################################################
00002 !               
00003 !    Copyright (C) 2003-2012 Department of Physics and Astronomy,
00004 !                            University of Rochester,
00005 !                            Rochester, NY
00006 !
00007 !    PFFT.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 !#########################################################################
00023 MODULE PFFT
00024   USE globaldeclarations
00025   USE LayoutDeclarations
00026   USE LayoutComms
00027   USE TreeDeclarations
00028   USE DataDeclarations
00029   USE Fields
00030   IMPLICIT NONE
00031   INCLUDE 'fftw3.f'
00032   SAVE
00033   PRIVATE
00034 
00035   TYPE PFFTPlanDef
00036      TYPE(pLayoutDef), DIMENSION(:), ALLOCATABLE :: Layouts
00037      COMPLEX(8), DIMENSION(:,:,:,:), POINTER :: data !transforms done 'in place'     
00038      INTEGER, DIMENSION(3,2) :: mB !Global bounds
00039      INTEGER, DIMENSION(3,2) :: lmB !Local bounds
00040   END type PFFTPlanDef
00041 
00042   INTEGER, PUBLIC, PARAMETER :: FORWARD = FFTW_FORWARD, BACKWARD = FFTW_BACKWARD
00043   
00044   
00045   PUBLIC :: CreatePlan, DestroyPlan, PFFTPlanDef, LoadFieldIntoPFFT, UnLoadFieldFromPFFT, ExecutePlan, SpectralProlongation
00046 
00047 
00048 CONTAINS
00049 
00053   SUBROUTINE ExecutePlan(plan, direction)
00054     TYPE(PFFTPlanDef) :: plan
00055     INTEGER :: direction, p
00056     INTEGER :: i
00057     DO i=1,nDim
00058        CALL DoFFT(plan%data, i, direction)
00059        CALL TransferLayouts(plan%layouts(i)%p,plan%layouts(mod(i,nDim)+1)%p, plan%data)
00060     END DO
00061   END SUBROUTINE ExecutePlan
00062 
00063   SUBROUTINE DoFFT(data, dir, direction)
00064     COMPLEX(8), DIMENSION(:,:,:,:), POINTER :: data
00065     COMPLEX(8), DIMENSION(:,:,:), POINTER :: in, out
00066     INTEGER :: dir, mx(3), ip(3,2), i, j, k, direction
00067     INTEGER(8) :: p
00068     mx=1
00069     mx(dir)=size(data,dir)
00070     ALLOCATE(in(mx(1),mx(2),mx(3)), out(mx(1),mx(2),mx(3)))
00071     ip(dir,:)=(/lbound(data,dir),ubound(data,dir)/)    
00072     CALL dfftw_plan_dft_1d(p, mx(dir), in, out, direction, FFTW_ESTIMATE)
00073     DO i=lBound(data,mod(dir,3)+1), ubound(data,mod(dir,3)+1)
00074        DO j=lBound(data,mod(dir+1,3)+1), ubound(data,mod(dir+1,3)+1)
00075           DO k=1,size(data,4)
00076              ip(mod(dir,3)+1,:)=i
00077              ip(mod(dir+1,3)+1,:)=j
00078              in=data(ip(1,1):ip(1,2), ip(2,1):ip(2,2), ip(3,1):ip(3,2),k)
00079              CALL dfftw_execute(p)
00080              data(ip(1,1):ip(1,2), ip(2,1):ip(2,2), ip(3,1):ip(3,2),k)=out
00081           END DO
00082        END DO
00083     END DO
00084     CALL dfftw_destroy_plan(p)
00085     DEALLOCATE(in, out)
00086   END SUBROUTINE DoFFT
00087 
00088   SUBROUTINE CreatePlan(plan, level, mB, fields)
00089     INTEGER, DIMENSION(3,2) :: mB
00090     TYPE(PFFTPlanDef), POINTER :: plan
00091     INTEGER :: i, fields, level
00092     ALLOCATE(plan)
00093     plan%mB=mB
00094     ALLOCATE(plan%layouts(1:nDim))
00095     DO i=1,nDim
00096        mB(i,:)=plan%mB(i,1)  !Shrink the bounds along each direction so that the FFT can be done
00097        CALL CreateLayout(mB, plan%layouts(i)%p)
00098        plan%layouts(i)%p%mB(:,i,2)=plan%mB(i,2) !re-expand the layouts 
00099        plan%layouts(i)%p%level=level
00100        mB(i,:)=plan%mB(i,:)
00101     END DO
00102     plan%lmB=plan%layouts(1)%p%mB(MPI_ID,:,:)
00103     ALLOCATE(plan%data(plan%lmB(1,1):plan%lmB(1,2), plan%lmB(2,1):plan%lmB(2,2), plan%lmB(3,1):plan%lmB(3,2),fields))
00104   END SUBROUTINE CreatePlan
00105 
00106   SUBROUTINE DestroyPlan(plan)
00107     TYPE(PFFTPlanDef), POINTER :: plan
00108     INTEGER :: i
00109     DO i=1,nDim
00110        CALL DestroyLayout(plan%layouts(i)%p)
00111     END DO
00112     DEALLOCATE(plan%layouts)
00113     DEALLOCATE(plan%data)
00114     DEALLOCATE(plan)
00115     NULLIFY(plan)
00116   END SUBROUTINE DestroyPlan
00117 
00118   SUBROUTINE LoadFieldIntoPFFT(plan,  FieldID)
00119     TYPE(PFFTPlanDef), POINTER :: plan
00120     INTEGER :: FieldID(:)
00121     INTEGER, DIMENSION(3,2) :: mB
00122     CALL LoadFieldIntoLayout(plan%layouts(1)%p, plan%data, FieldID)
00123   END SUBROUTINE LoadFieldIntoPFFT
00124 
00125 
00126   SUBROUTINE UnLoadFieldFromPFFT(plan, FieldID, lPeriodic, rmbc)
00127     TYPE(PFFTPlanDef), POINTER :: plan
00128     INTEGER :: FieldID(:,:)
00129     INTEGER, DIMENSION(3,2) :: mB
00130     INTEGER :: rmbc
00131     LOGICAL, DIMENSION(:) :: lPeriodic
00132 !    mB=plan%layouts(1)%p%mB(MPI_ID,:,:)
00133 !    ALLOCATE(plan%data(mB(1,1):mB(1,2), mB(2,1):mB(2,2), mB(3,1):mB(3,2)))
00134     CALL UnLoadFieldFromLayout(plan%layouts(1)%p, plan%data, FieldID, lPeriodic, rmbc)
00135   END SUBROUTINE UnLoadFieldFromPFFT
00136   
00137 
00138 
00139   SUBROUTINE SpectralProlongation(plan, newplan)
00140      TYPE(PFFTPlanDef), POINTER :: plan, newplan     
00141      INTEGER :: i, j, k, n, l(3), mB(3,2), mC(3,2), mD(3,2), splitpoints(3,2), mA(3,2)
00142      REAL(KIND=qPREC) :: dphi(3), phase, kvec(3)
00143      !First do shifting
00144 
00145 
00146 
00147 
00148 !     mB=plan%lmB
00149 !     DO i=mB(1,1),mB(1,2)
00150 !        DO j=mB(2,1), mB(2,2)
00151 !           DO k=mB(3,1), mB(3,2)
00152 !              kvec=SpectraK((/i,j,k/),plan%mB)
00153 !              phase=sum(dphi * kvec)
00154 !              plan%data(i,j,k,:)=plan%data(i,j,k,:)*exp(cmplx(0d0,phase))
00155 !           END DO
00156 !        END DO
00157 !     END DO
00158 !     plan%data=plan%data/product(plan%mB(:,2)-plan%mB(:,1)+1)
00159 
00160      !consider an 8 point plan
00161      ! the wave vectors will be
00162      ! 0 1 2 3 4 -3 -2 -1
00163      ! and the split point is 5
00164      ! we want to divide the 5th entry in half and send 1:5 and 5:8
00165      
00166      !consider a 7 point plan
00167      ! 0 1 2 3 -3 -2 -1
00168      ! the split point is 4 
00169      ! and we want to send 1:4 and 5:7
00170 
00171      splitpoints(:,1)=(sum(plan%mB,2)+0)/2+1
00172      splitpoints(:,2)=(sum(plan%mB,2)+1)/2+1
00173     
00174      
00175      ! First cut duplicates reals in half.
00176      mA=plan%mB
00177      mB=plan%lmB
00178      DO i=1,nDim
00179         IF (splitpoints(i,1)==splitpoints(i,2) .AND. (splitpoints(i,1) >= mB(i,1) .AND. splitpoints(i,1) <= mB(i,2))) THEN
00180            mA(i,:)=splitpoints(i,1)
00181            plan%data(mA(1,1):mA(1,2), mA(2,1):mA(2,2), mA(3,1):mA(3,2),:)=half*&
00182                 plan%data(mA(1,1):mA(1,2), mA(2,1):mA(2,2), mA(3,1):mA(3,2),:)
00183            mA(i,:)=mB(i,:)
00184         END IF
00185      END DO
00186 
00187      ! Then stretch fourier components
00188 
00189      newplan%data=0d0
00190      l=0
00191      DO i=0,2**nDim-1
00192         mC=plan%mB
00193         mD=newplan%mB
00194         DO n=1,nDim
00195            l(n)=MOD(i/2**(n-1),2)
00196            IF (l(n) == 0) THEN
00197               mC(n,2)=splitpoints(n,1)
00198               mD(n,2)=mD(n,1)+(mC(n,2)-mC(n,1))
00199            ELSE
00200               mC(n,1)=splitpoints(n,2)+1
00201               mD(n,1)=mD(n,2)-(mC(n,2)-mC(n,1))
00202            END IF
00203         END DO
00204         CALL LayoutTransferC(mC, mD, plan%layouts(1)%p, newplan%layouts(1)%p, plan%data, newplan%data)
00205      END DO
00206 
00207      ! And finally adjust phases due to grid point locations.
00208 
00209      dphi=-Pi/REAL(newplan%mB(:,2)-newplan%mB(:,1)+1)
00210      mB=newplan%lmB
00211      DO i=mB(1,1),mB(1,2)
00212         DO j=mB(2,1), mB(2,2)
00213            DO k=mB(3,1), mB(3,2)
00214               kvec=SpectraK((/i,j,k/),newplan%mB)
00215               phase=sum(dphi * kvec)
00216               newplan%data(i,j,k,:)=newplan%data(i,j,k,:)*exp(cmplx(0d0,phase))
00217            END DO
00218         END DO
00219      END DO
00220 
00221      ! And divide final data by size of original FFT.
00222      newplan%data=newplan%data/product(plan%mB(:,2)-plan%mB(:,1)+1)
00223 
00224   END SUBROUTINE SpectralProlongation
00225 
00226 END MODULE PFFT
 All Classes Files Functions Variables