quantum-espresso/CPV/fft.f90

670 lines
22 KiB
Fortran
Raw Normal View History

!
! Copyright (C) 2002-2005 FPMD-CPV groups
! This file is distributed under the terms of the
! GNU General Public License. See the file `License'
! in the root directory of the present distribution,
! or http://www.gnu.org/copyleft/gpl.txt .
!
! ----------------------------------------------
! This Subroutines written by Carlo Cavazzoni
! Last modified December 2008
! ----------------------------------------------
#include "f_defs.h"
!-----------------------------------------------------------------------
subroutine invfft_x( grid_type, f, dfft, ia )
!-----------------------------------------------------------------------
! grid_type = 'Dense'
! inverse fourier transform of potentials and charge density
! on the dense grid . On output, f is overwritten
! grid_type = 'Smooth'
! inverse fourier transform of potentials and charge density
! on the smooth grid . On output, f is overwritten
! grid_type = 'Wave'
! inverse fourier transform of wave functions
! on the smooth grid . On output, f is overwritten
! grid_type = 'Box'
! not-so-parallel 3d fft for box grid, implemented only for sign=1
! G-space to R-space, output = \sum_G f(G)exp(+iG*R)
! The array f (overwritten on output) is NOT distributed:
! a copy is present on each processor.
! The fft along z is done on the entire grid.
! The fft along xy is done only on planes that have components on the
! dense grid for each processor. Note that the final array will no
! longer be the same on all processors.
!
!
USE kinds, ONLY: DP
use fft_base, only: dfftp, dffts, dfftb
use fft_scalar, only: cfft3d, cfft3ds, cft_b
use fft_parallel, only: tg_cft3s
use control_flags, only: use_task_groups
USE fft_types, only: fft_dlay_descriptor
IMPLICIT none
TYPE(fft_dlay_descriptor), INTENT(IN) :: dfft
INTEGER, OPTIONAL, INTENT(IN) :: ia
CHARACTER(LEN=*), INTENT(IN) :: grid_type
COMPLEX(DP) :: f(:)
!
INTEGER :: imin3, imax3, np3
IF( grid_type == 'Dense' ) THEN
IF( dfft%nr1 /= dfftp%nr1 .OR. dfft%nr2 /= dfftp%nr2 .OR. dfft%nr3 /= dfftp%nr3 .OR. &
dfft%nr1x /= dfftp%nr1x .OR. dfft%nr2x /= dfftp%nr2x .OR. dfft%nr3x /= dfftp%nr3x ) &
CALL errore( ' invfft ', ' inconsistent descriptor for Dense fft ' , 1 )
call start_clock( 'fft' )
ELSE IF( grid_type == 'Smooth' ) THEN
IF( dfft%nr1 /= dffts%nr1 .OR. dfft%nr2 /= dffts%nr2 .OR. dfft%nr3 /= dffts%nr3 .OR. &
dfft%nr1x /= dffts%nr1x .OR. dfft%nr2x /= dffts%nr2x .OR. dfft%nr3x /= dffts%nr3x ) &
CALL errore( ' invfft ', ' inconsistent descriptor for Smooth fft ' , 1 )
call start_clock( 'ffts' )
ELSE IF( grid_type == 'Wave' ) THEN
IF( dfft%nr1 /= dffts%nr1 .OR. dfft%nr2 /= dffts%nr2 .OR. dfft%nr3 /= dffts%nr3 .OR. &
dfft%nr1x /= dffts%nr1x .OR. dfft%nr2x /= dffts%nr2x .OR. dfft%nr3x /= dffts%nr3x ) &
CALL errore( ' invfft ', ' inconsistent descriptor for Wave fft ' , 1 )
call start_clock('fftw')
ELSE IF( grid_type == 'Box' ) THEN
IF( dfft%nr1 /= dfftb%nr1 .OR. dfft%nr2 /= dfftb%nr2 .OR. dfft%nr3 /= dfftb%nr3 .OR. &
dfft%nr1x /= dfftb%nr1x .OR. dfft%nr2x /= dfftb%nr2x .OR. dfft%nr3x /= dfftb%nr3x ) &
CALL errore( ' invfft ', ' inconsistent descriptor for Box fft ' , 1 )
call start_clock( 'fftb' )
ELSE
call errore( ' invfft ', ' unknown grid: '//grid_type , 1 )
END IF
#if defined __PARA && !defined __USE_3D_FFT
IF( grid_type == 'Box' ) THEN
imin3 = dfftb%imin3( ia )
imax3 = dfftb%imax3( ia )
np3 = dfftb%np3( ia ) ! imax3 - imin3 + 1
END IF
IF( grid_type == 'Dense' ) THEN
call tg_cft3s( f, dfftp, 1 )
ELSE IF( grid_type == 'Smooth' ) THEN
call tg_cft3s( f, dffts, 1 )
ELSE IF( grid_type == 'Wave' ) THEN
call tg_cft3s( f, dffts, 2, use_task_groups )
ELSE IF( grid_type == 'Box' .AND. np3 > 0 ) THEN
call cft_b( f, dfftb%nr1, dfftb%nr2, dfftb%nr3, dfftb%nr1x, dfftb%nr2x, dfftb%nr3x, imin3, imax3, 1 )
END IF
#else
# if defined __SCSL || __SX6 || __USE_3D_FFT
IF( grid_type == 'Dense' ) THEN
call cfft3d( f, dfftp%nr1, dfftp%nr2, dfftp%nr3, dfftp%nr1x, dfftp%nr2x, dfftp%nr3x, 1)
ELSE IF( grid_type == 'Smooth' ) THEN
call cfft3d( f, dffts%nr1, dffts%nr2, dffts%nr3, dffts%nr1x, dffts%nr2x, dffts%nr3x, 1)
ELSE IF( grid_type == 'Wave' ) THEN
call cfft3d( f, dffts%nr1, dffts%nr2, dffts%nr3, dffts%nr1x, dffts%nr2x, dffts%nr3x, 1)
ELSE IF( grid_type == 'Box' ) THEN
call cfft3d( f, dfftb%nr1, dfftb%nr2, dfftb%nr3, dfftb%nr1x, dfftb%nr2x, dfftb%nr3x, 1)
END IF
# elif defined __ESSL || __LINUX_ESSL || __FFTW || __FFTW3
IF( grid_type == 'Dense' ) THEN
call cfft3d( f, dfftp%nr1, dfftp%nr2, dfftp%nr3, dfftp%nr1x, dfftp%nr2x, dfftp%nr3x, 1)
ELSE IF( grid_type == 'Smooth' ) THEN
call cfft3d( f, dffts%nr1, dffts%nr2, dffts%nr3, dffts%nr1x, dffts%nr2x, dffts%nr3x, 1)
ELSE IF( grid_type == 'Wave' ) THEN
call cfft3ds( f, dffts%nr1, dffts%nr2, dffts%nr3, dffts%nr1x, dffts%nr2x, dffts%nr3x, 1, dffts%isind, dffts%iplw )
ELSE IF( grid_type == 'Box' ) THEN
call cfft3d( f, dfftb%nr1, dfftb%nr2, dfftb%nr3, dfftb%nr1x, dfftb%nr2x, dfftb%nr3x, 1)
END IF
# endif
#endif
IF( grid_type == 'Dense' ) THEN
call stop_clock( 'fft' )
ELSE IF( grid_type == 'Smooth' ) THEN
call stop_clock( 'ffts' )
ELSE IF( grid_type == 'Wave' ) THEN
call stop_clock('fftw')
ELSE IF( grid_type == 'Box' ) THEN
call stop_clock( 'fftb' )
END IF
!
return
end subroutine invfft_x
!-----------------------------------------------------------------------
subroutine fwfft_x( grid_type, f, dfft )
!-----------------------------------------------------------------------
! grid_type = 'Dense'
! forward fourier transform of potentials and charge density
! on the dense grid . On output, f is overwritten
! grid_type = 'Smooth'
! forward fourier transform of potentials and charge density
! on the smooth grid . On output, f is overwritten
! grid_type = 'Wave'
! forward fourier transform of wave functions
! on the smooth grid . On output, f is overwritten
!
USE kinds, ONLY: DP
use fft_base, only: dfftp, dffts
use fft_scalar, only: cfft3d, cfft3ds
use fft_parallel, only: tg_cft3s
use control_flags, only: use_task_groups
USE fft_types, only: fft_dlay_descriptor
implicit none
TYPE(fft_dlay_descriptor), INTENT(IN) :: dfft
CHARACTER(LEN=*), INTENT(IN) :: grid_type
COMPLEX(DP) :: f(:)
IF( grid_type == 'Dense' ) THEN
IF( dfft%nr1 /= dfftp%nr1 .OR. dfft%nr2 /= dfftp%nr2 .OR. dfft%nr3 /= dfftp%nr3 .OR. &
dfft%nr1x /= dfftp%nr1x .OR. dfft%nr2x /= dfftp%nr2x .OR. dfft%nr3x /= dfftp%nr3x ) &
CALL errore( ' fwfft ', ' inconsistent descriptor for Dense fft ' , 1 )
call start_clock( 'fft' )
ELSE IF( grid_type == 'Smooth' ) THEN
IF( dfft%nr1 /= dffts%nr1 .OR. dfft%nr2 /= dffts%nr2 .OR. dfft%nr3 /= dffts%nr3 .OR. &
dfft%nr1x /= dffts%nr1x .OR. dfft%nr2x /= dffts%nr2x .OR. dfft%nr3x /= dffts%nr3x ) &
CALL errore( ' fwfft ', ' inconsistent descriptor for Smooth fft ' , 1 )
call start_clock( 'ffts' )
ELSE IF( grid_type == 'Wave' ) THEN
IF( dfft%nr1 /= dffts%nr1 .OR. dfft%nr2 /= dffts%nr2 .OR. dfft%nr3 /= dffts%nr3 .OR. &
dfft%nr1x /= dffts%nr1x .OR. dfft%nr2x /= dffts%nr2x .OR. dfft%nr3x /= dffts%nr3x ) &
CALL errore( ' fwfft ', ' inconsistent descriptor for Wave fft ' , 1 )
call start_clock( 'fftw' )
ELSE
call errore( ' fwfft ', ' unknown grid: '//grid_type , 1 )
END IF
#if defined __PARA && !defined __USE_3D_FFT
IF( grid_type == 'Dense' ) THEN
call tg_cft3s(f,dfftp,-1)
ELSE IF( grid_type == 'Smooth' ) THEN
call tg_cft3s(f,dffts,-1)
ELSE IF( grid_type == 'Wave' ) THEN
call tg_cft3s(f,dffts,-2, use_task_groups )
END IF
#else
# if defined __SCSL || __SX6 || __USE_3D_FFT
IF( grid_type == 'Dense' ) THEN
call cfft3d( f, dfftp%nr1, dfftp%nr2, dfftp%nr3, dfftp%nr1x, dfftp%nr2x, dfftp%nr3x, -1)
ELSE IF( grid_type == 'Smooth' ) THEN
call cfft3d( f, dffts%nr1, dffts%nr2, dffts%nr3, dffts%nr1x, dffts%nr2x, dffts%nr3x, -1)
ELSE IF( grid_type == 'Wave' ) THEN
call cfft3d( f, dffts%nr1, dffts%nr2, dffts%nr3, dffts%nr1x, dffts%nr2x, dffts%nr3x, -1)
END IF
# elif defined __ESSL || __LINUX_ESSL || __FFTW || __FFTW3
IF( grid_type == 'Dense' ) THEN
call cfft3d( f, dfftp%nr1, dfftp%nr2, dfftp%nr3, dfftp%nr1x, dfftp%nr2x, dfftp%nr3x, -1)
ELSE IF( grid_type == 'Smooth' ) THEN
call cfft3d( f, dffts%nr1, dffts%nr2, dffts%nr3, dffts%nr1x, dffts%nr2x, dffts%nr3x, -1)
ELSE IF( grid_type == 'Wave' ) THEN
call cfft3ds( f, dffts%nr1, dffts%nr2, dffts%nr3, dffts%nr1x, dffts%nr2x, dffts%nr3x, -1, dffts%isind, dffts%iplw )
END IF
# endif
#endif
IF( grid_type == 'Dense' ) THEN
call stop_clock( 'fft' )
ELSE IF( grid_type == 'Smooth' ) THEN
call stop_clock( 'ffts' )
ELSE IF( grid_type == 'Wave' ) THEN
call stop_clock( 'fftw' )
END IF
return
end subroutine fwfft_x
!-----------------------------------------------------------------------
SUBROUTINE c2psi( psi, nnr, c, ca, ng, iflg )
!
use gvecs, only: nms, nps
use kinds, only: DP
implicit none
complex(DP) :: psi(*), c(*), ca(*)
integer, intent(in) :: nnr, ng, iflg
complex(DP), parameter :: ci=(0.0d0,1.0d0)
integer :: ig
psi( 1 : nnr ) = 0.0d0
!
! iflg "cases"
!
! 0 Do not use gamma symmetry
!
! 1 set psi using a wf with Gamma symmetry
! 2 set psi combining two wf with Gamma symmetry
!
SELECT CASE ( iflg )
!
! Case 0, 1 and 2 SMOOTH MESH
!
CASE ( 0 )
!
do ig = 1, ng
psi( nps( ig ) ) = c( ig )
end do
!
CASE ( 1 )
!
do ig = 1, ng
General cleanup of intrinsic functions: conversion to real => DBLE (including real part of a complex number) conversion to complex => CMPLX complex conjugate => CONJG imaginary part => AIMAG All functions are uppercase. CMPLX is preprocessed by f_defs.h and performs an explicit cast: #define CMPLX(a,b) cmplx(a,b,kind=DP) This implies that 1) f_defs.h must be included whenever a CMPLX is present, 2) CMPLX should stay in a single line, 3) DP must be defined. All occurrences of real, float, dreal, dfloat, dconjg, dimag, dcmplx removed - please do not reintroduce any of them. Tested only with ifc7 and g95 - beware unintended side effects Maybe not the best solution (explicit casts everywhere would be better) but it can be easily changed with a script if the need arises. The following code might be used to test for possible trouble: program test_intrinsic implicit none integer, parameter :: dp = selected_real_kind(14,200) real (kind=dp) :: a = 0.123456789012345_dp real (kind=dp) :: b = 0.987654321098765_dp complex (kind=dp) :: c = ( 0.123456789012345_dp, 0.987654321098765_dp) print *, ' A = ', a print *, ' DBLE(A)= ', DBLE(a) print *, ' C = ', c print *, 'CONJG(C)= ', CONJG(c) print *, 'DBLE(c),AIMAG(C) = ', DBLE(c), AIMAG(c) print *, 'CMPLX(A,B,kind=dp)= ', CMPLX( a, b, kind=dp) end program test_intrinsic Note that CMPLX and REAL without a cast yield single precision numbers on ifc7 and g95 !!! git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@2133 c92efa57-630b-4861-b058-cf58834340f0
2005-08-27 01:44:42 +08:00
psi( nms( ig ) ) = CONJG( c( ig ) )
psi( nps( ig ) ) = c( ig )
end do
!
CASE ( 2 )
!
do ig = 1, ng
General cleanup of intrinsic functions: conversion to real => DBLE (including real part of a complex number) conversion to complex => CMPLX complex conjugate => CONJG imaginary part => AIMAG All functions are uppercase. CMPLX is preprocessed by f_defs.h and performs an explicit cast: #define CMPLX(a,b) cmplx(a,b,kind=DP) This implies that 1) f_defs.h must be included whenever a CMPLX is present, 2) CMPLX should stay in a single line, 3) DP must be defined. All occurrences of real, float, dreal, dfloat, dconjg, dimag, dcmplx removed - please do not reintroduce any of them. Tested only with ifc7 and g95 - beware unintended side effects Maybe not the best solution (explicit casts everywhere would be better) but it can be easily changed with a script if the need arises. The following code might be used to test for possible trouble: program test_intrinsic implicit none integer, parameter :: dp = selected_real_kind(14,200) real (kind=dp) :: a = 0.123456789012345_dp real (kind=dp) :: b = 0.987654321098765_dp complex (kind=dp) :: c = ( 0.123456789012345_dp, 0.987654321098765_dp) print *, ' A = ', a print *, ' DBLE(A)= ', DBLE(a) print *, ' C = ', c print *, 'CONJG(C)= ', CONJG(c) print *, 'DBLE(c),AIMAG(C) = ', DBLE(c), AIMAG(c) print *, 'CMPLX(A,B,kind=dp)= ', CMPLX( a, b, kind=dp) end program test_intrinsic Note that CMPLX and REAL without a cast yield single precision numbers on ifc7 and g95 !!! git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@2133 c92efa57-630b-4861-b058-cf58834340f0
2005-08-27 01:44:42 +08:00
psi( nms( ig ) ) = CONJG( c( ig ) ) + ci * conjg( ca( ig ) )
psi( nps( ig ) ) = c( ig ) + ci * ca( ig )
end do
CASE DEFAULT
!
CALL errore(" c2psi "," wrong value for iflg ", ABS( iflg ) )
END SELECT
return
END SUBROUTINE c2psi
!
!
!
SUBROUTINE rho2psi( grid_type, psi, nnr, rho, ng )
!
use recvecs_indexes, only: nm, np
use gvecs, only: nms, nps
use kinds, only: DP
implicit none
complex(DP) :: psi(*), rho(*)
integer, intent(in) :: nnr, ng
character(len=*), intent(in) :: grid_type
integer :: ig
psi( 1 : nnr ) = 0.0d0
SELECT CASE ( grid_type )
!
! Case 0, 1 and 2 SMOOTH MESH
!
CASE ( 'Smooth' )
!
! without gamma sym
! do ig = 1, ng
! psi( nps( ig ) ) = rho( ig )
! end do
!
do ig = 1, ng
psi( nms( ig ) ) = CONJG( rho( ig ) )
psi( nps( ig ) ) = rho( ig )
end do
!
CASE ( 'Dense' )
!
! do ig = 1, ng
! psi( np( ig ) ) = rho( ig )
! end do
!
do ig = 1, ng
psi( nm( ig ) ) = CONJG( rho( ig ) )
psi( np( ig ) ) = rho( ig )
end do
!
CASE DEFAULT
!
CALL errore(" rho2psi "," wrong grid "//grid_type , 1 )
END SELECT
return
END SUBROUTINE rho2psi
!-----------------------------------------------------------------------
SUBROUTINE psi2c( psi, nnr, c, ca, ng, iflg )
use recvecs_indexes, only: nm, np
use gvecs, only: nms, nps
use kinds, only: DP
implicit none
complex(DP) :: psi(*), c(*), ca(*)
integer, intent(in) :: nnr, ng, iflg
complex(DP), parameter :: ci=(0.0d0,1.0d0)
integer :: ig
!
! iflg "cases"
!
! 0, 10 Do not use gamma symmetry
!
! 1, 11 set psi using a wf with Gamma symmetry
!
! 2, 12 set psi combining two wf with Gamma symmetry
!
SELECT CASE ( iflg )
!
! Case 0, 1 and 2 SMOOTH MESH
!
CASE ( 0 )
!
do ig = 1, ng
c( ig ) = psi( nps( ig ) )
end do
!
CASE ( 1 )
!
CALL errore(" psi2c "," wrong value for iflg ", 11 )
!
CASE ( 2 )
!
DO ig = 1, ng
ca(ig) = psi( nms( ig ) )
c (ig) = psi( nps( ig ) )
END DO
!
! Case 10, 11 and 12 DENSE MESH
!
CASE ( 10 )
!
do ig = 1, ng
c( ig ) = psi( np( ig ) )
end do
!
CASE ( 11 )
!
CALL errore(" psi2c "," wrong value for iflg ", 1 )
!
CASE ( 12 )
!
DO ig = 1, ng
ca(ig) = psi( nm( ig ) )
c (ig) = psi( np( ig ) )
END DO
CASE DEFAULT
!
CALL errore(" psi2c "," wrong value for iflg ", ABS( iflg ) )
END SELECT
return
END SUBROUTINE psi2c
!-----------------------------------------------------------------------
SUBROUTINE psi2rho( grid_type, psi, nnr, rho, ng )
use recvecs_indexes, only: nm, np
use gvecs, only: nms, nps
use kinds, only: DP
implicit none
complex(DP) :: psi(*), rho(*)
integer, intent(in) :: nnr, ng
character(len=*), intent(in) :: grid_type
integer :: ig
SELECT CASE ( grid_type )
!
CASE ( 'Smooth' )
!
do ig = 1, ng
rho( ig ) = psi( nps( ig ) )
end do
!
CASE ( 'Dense' )
!
do ig = 1, ng
rho( ig ) = psi( np( ig ) )
end do
!
CASE DEFAULT
!
CALL errore(" psi2rho "," wrong grid "//grid_type , 1 )
END SELECT
return
END SUBROUTINE psi2rho
!-----------------------------------------------------------------------
SUBROUTINE box2grid(irb,nfft,qv,vr)
!-----------------------------------------------------------------------
!
! add array qv(r) on box grid to array vr(r) on dense grid
! irb : position of the box in the dense grid
! nfft=1 add real part of qv(r) to real part of array vr(r)
! nfft=2 add imaginary part of qv(r) to real part of array vr(r)
!
USE grid_dimensions, ONLY: nr1, nr2, nr3, &
nr1x, nr2x, nnr => nnrx
USE smallbox_grid_dimensions, ONLY: nr1b, nr2b, nr3b, &
nr1bx, nr2bx, nnrb => nnrbx
USE fft_base, ONLY: dfftp
USE mp_global, ONLY: me_image
IMPLICIT NONE
INTEGER, INTENT(in):: nfft, irb(3)
REAL(8), INTENT(in):: qv(2,nnrb)
COMPLEX(8), INTENT(inout):: vr(nnr)
!
INTEGER ir1, ir2, ir3, ir, ibig1, ibig2, ibig3, ibig
INTEGER me
IF(nfft.LE.0.OR.nfft.GT.2) CALL errore('box2grid','wrong data',nfft)
me = me_image + 1
DO ir3=1,nr3b
ibig3=irb(3)+ir3-1
ibig3=1+MOD(ibig3-1,nr3)
IF(ibig3.LT.1.OR.ibig3.GT.nr3) &
& CALL errore('box2grid','ibig3 wrong',ibig3)
ibig3=ibig3-dfftp%ipp(me)
IF ( ibig3 .GT. 0 .AND. ibig3 .LE. ( dfftp%npp(me) ) ) THEN
DO ir2=1,nr2b
ibig2=irb(2)+ir2-1
ibig2=1+MOD(ibig2-1,nr2)
IF(ibig2.LT.1.OR.ibig2.GT.nr2) &
& CALL errore('box2grid','ibig2 wrong',ibig2)
DO ir1=1,nr1b
ibig1=irb(1)+ir1-1
ibig1=1+MOD(ibig1-1,nr1)
IF(ibig1.LT.1.OR.ibig1.GT.nr1) &
& CALL errore('box2grid','ibig1 wrong',ibig1)
ibig=ibig1+(ibig2-1)*nr1x+(ibig3-1)*nr1x*nr2x
ir=ir1+(ir2-1)*nr1bx+(ir3-1)*nr1bx*nr2bx
vr(ibig) = vr(ibig)+qv(nfft,ir)
END DO
END DO
END IF
END DO
!
RETURN
END SUBROUTINE box2grid
!-----------------------------------------------------------------------
SUBROUTINE box2grid2(irb,qv,v)
!-----------------------------------------------------------------------
!
! add array qv(r) on box grid to array v(r) on dense grid
! irb : position of the box in the dense grid
!
USE grid_dimensions, ONLY: nr1, nr2, nr3, &
nr1x, nr2x, nnr => nnrx
USE smallbox_grid_dimensions, ONLY: nr1b, nr2b, nr3b, &
nr1bx, nr2bx, nnrb => nnrbx
USE fft_base, ONLY: dfftp
USE mp_global, ONLY: me_image
!
IMPLICIT NONE
!
INTEGER, INTENT(in):: irb(3)
COMPLEX(8), INTENT(in):: qv(nnrb)
COMPLEX(8), INTENT(inout):: v(nnr)
!
INTEGER ir1, ir2, ir3, ir, ibig1, ibig2, ibig3, ibig
INTEGER me
me = me_image + 1
DO ir3=1,nr3b
ibig3=irb(3)+ir3-1
ibig3=1+MOD(ibig3-1,nr3)
IF(ibig3.LT.1.OR.ibig3.GT.nr3) &
& CALL errore('box2grid2','ibig3 wrong',ibig3)
ibig3=ibig3-dfftp%ipp(me)
IF (ibig3.GT.0.AND.ibig3.LE. dfftp%npp(me) ) THEN
DO ir2=1,nr2b
ibig2=irb(2)+ir2-1
ibig2=1+MOD(ibig2-1,nr2)
IF(ibig2.LT.1.OR.ibig2.GT.nr2) &
& CALL errore('box2grid2','ibig2 wrong',ibig2)
DO ir1=1,nr1b
ibig1=irb(1)+ir1-1
ibig1=1+MOD(ibig1-1,nr1)
IF(ibig1.LT.1.OR.ibig1.GT.nr1) &
& CALL errore('box2grid2','ibig1 wrong',ibig1)
ibig=ibig1+(ibig2-1)*nr1x+(ibig3-1)*nr1x*nr2x
ir=ir1+(ir2-1)*nr1bx+(ir3-1)*nr1bx*nr2bx
v(ibig) = v(ibig)+qv(ir)
END DO
END DO
END IF
END DO
RETURN
END SUBROUTINE box2grid2
!-----------------------------------------------------------------------
REAL(8) FUNCTION boxdotgrid(irb,nfft,qv,vr)
!-----------------------------------------------------------------------
!
! Calculate \sum_i qv(r_i)*vr(r_i) with r_i on box grid
! array qv(r) is defined on box grid, array vr(r)on dense grid
! irb : position of the box in the dense grid
! nfft=1 (2): use real (imaginary) part of qv(r)
! Parallel execution: remember to sum the contributions from other nodes
!
USE grid_dimensions, ONLY: nr1, nr2, nr3, &
nr1x, nr2x, nnr => nnrx
USE smallbox_grid_dimensions, ONLY: nr1b, nr2b, nr3b, &
nr1bx, nr2bx, nnrb => nnrbx
USE fft_base, ONLY: dfftp
USE mp_global, ONLY: me_image
IMPLICIT NONE
INTEGER, INTENT(in):: nfft, irb(3)
REAL(8), INTENT(in):: qv(2,nnrb), vr(nnr)
!
INTEGER ir1, ir2, ir3, ir, ibig1, ibig2, ibig3, ibig
INTEGER me
!
!
IF(nfft.LE.0.OR.nfft.GT.2) CALL errore('boxdotgrid','wrong data',nfft)
me = me_image + 1
boxdotgrid=0.d0
DO ir3=1,nr3b
ibig3=irb(3)+ir3-1
ibig3=1+MOD(ibig3-1,nr3)
ibig3=ibig3-dfftp%ipp(me)
IF (ibig3.GT.0.AND.ibig3.LE. dfftp%npp(me) ) THEN
DO ir2=1,nr2b
ibig2=irb(2)+ir2-1
ibig2=1+MOD(ibig2-1,nr2)
DO ir1=1,nr1b
ibig1=irb(1)+ir1-1
ibig1=1+MOD(ibig1-1,nr1)
ibig=ibig1 + (ibig2-1)*nr1x + (ibig3-1)*nr1x*nr2x
ir =ir1 + (ir2-1)*nr1bx + (ir3-1)*nr1bx*nr2bx
boxdotgrid = boxdotgrid + qv(nfft,ir)*vr(ibig)
END DO
END DO
ENDIF
END DO
RETURN
END FUNCTION boxdotgrid
!
!----------------------------------------------------------------------
subroutine parabox(nr3b,irb3,nr3,imin3,imax3)
!----------------------------------------------------------------------
!
! find if box grid planes in the z direction have component on the dense
! grid on this processor, and if, which range imin3-imax3
!
use mp_global, only: me_image
use fft_base, only: dfftp
! input
integer nr3b,irb3,nr3
! output
integer imin3,imax3
! local
integer ir3, ibig3, me
!
me = me_image + 1
imin3=nr3b
imax3=1
do ir3=1,nr3b
ibig3=1+mod(irb3+ir3-2,nr3)
if(ibig3.lt.1.or.ibig3.gt.nr3) &
& call errore('cfftpb','ibig3 wrong',ibig3)
ibig3=ibig3-dfftp%ipp(me)
if (ibig3.gt.0.and.ibig3.le.dfftp%npp(me)) then
imin3=min(imin3,ir3)
imax3=max(imax3,ir3)
end if
end do
!
return
end subroutine parabox