TABLE OF CONTENTS


ABINIT/dotprodm_v [ Functions ]

[ Top ] [ Functions ]

NAME

 dotprodm_v

FUNCTION

 For two sets of potentials,
 compute dot product of each pair of two potentials (integral over FFT grid), to obtain
 a series of square residual-like quantity (so the sum of product of values
 is NOT divided by the number of FFT points, and NOT multiplied by the primitive cell volume).
 Take into account the spin components of the potentials (nspden),
 and sum over them.
 Need the index of the first pair of potentials to be treated, in each array
 of potentials, and the number of potentials to be treated.
 Might be used to compute just one square of norm, in
 a big array, such as to avoid copying a potential from a big array
 to a temporary place.

COPYRIGHT

 Copyright (C) 1999-2018 ABINIT group (XG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  cplex=if 1, real space functions on FFT grid are REAL, if 2, COMPLEX
  cpldot=if 1, the dot array is real, if 2, the dot array is complex
  index1=index of the first potential to be treated in the potarr1 array
  index2=index of the first potential to be treated in the potarr2 array
  mpicomm=the mpi communicator used for the summation
  mpi_summarize=set it to .true. if parallelisation is done over FFT
  mult1=number of potentials to be treated in the first set
  mult2=number of potentials to be treated in the second set
  nfft= (effective) number of FFT grid points (for this processor)
  npot1= third dimension of the potarr1 array
  npot2= third dimension of the potarr2 array
  nspden=number of spin-density components
  opt_storage: 0, if potentials are stored as V^up-up, V^dn-dn, Re[V^up-dn], Im[V^up-dn]
               1, if potentials are stored as V, B_x, B_y, Bz  (B=magn. field)
  potarr1(cplex*nfft,nspden,npot)=first array of real space potentials on FFT grid
    (if cplex=2 and cpldot=2, potarr1 is the array that will be complex conjugated)
  potarr2(cplex*nfft,nspden,npot)=second array of real space potentials on FFT grid

OUTPUT

  dot(cpldot,mult1,mult2)= series of values of the dot product

SIDE EFFECTS

NOTES

  Concerning storage when nspden=4:
   cplex=1:
     opt_storage=0: V are stored as : V^11, V^22, Re[V^12], Im[V^12] (complex, hermitian)
     opt_storage=1: V are stored as : V, B_x, B_y, B_z               (real)
   cplex=2:
     opt_storage=0: V are stored as : V^11, V^22, V^12, i.V^21 (complex)
     opt_storage=1: V are stored as : V, B_x, B_y, B_z         (complex)

PARENTS

      scfopt

CHILDREN

      timab,xmpi_sum

SOURCE

 68 #if defined HAVE_CONFIG_H
 69 #include "config.h"
 70 #endif
 71 
 72 #include "abi_common.h"
 73 
 74 
 75 subroutine dotprodm_v(cplex,cpldot,dot,index1,index2,mpicomm,mpi_summarize,&
 76 &   mult1,mult2,nfft,npot1,npot2,nspden,opt_storage,potarr1,potarr2)
 77 
 78  use defs_basis
 79  use m_profiling_abi
 80  use m_errors
 81  use m_xmpi
 82 
 83 !This section has been created automatically by the script Abilint (TD).
 84 !Do not modify the following lines by hand.
 85 #undef ABI_FUNC
 86 #define ABI_FUNC 'dotprodm_v'
 87  use interfaces_18_timing
 88 !End of the abilint section
 89 
 90  implicit none
 91 
 92 !Arguments ------------------------------------
 93 !scalars
 94  integer,intent(in) :: cpldot,cplex,index1,index2,mult1,mult2,nfft,npot1,npot2
 95  integer,intent(in) :: nspden,opt_storage,mpicomm
 96  logical, intent(in) :: mpi_summarize
 97 !arrays
 98  real(dp),intent(in) :: potarr1(cplex*nfft,nspden,npot1)
 99  real(dp),intent(in) :: potarr2(cplex*nfft,nspden,npot2)
100  real(dp),intent(out) :: dot(cpldot,mult1,mult2)
101 
102 !Local variables-------------------------------
103 !scalars
104  integer :: i1,i2,ierr,ifft,ispden
105  real(dp) :: ai,ar
106 !arrays
107  real(dp) :: tsec(2)
108 
109 ! *************************************************************************
110 
111 !Real or complex inputs are coded
112  DBG_CHECK(ANY(cplex==(/1,2/)),"Wrong cplex")
113 
114 !Real or complex outputs are coded
115  DBG_CHECK(ANY(cpldot==(/1,2/)),"Wrong cpldot")
116  DBG_CHECK(ANY(nspden==(/1,2,4/)),"Wrong nspden")
117  DBG_CHECK( npot1-index1-mult1 >= -1,"npot1-index1-mult1")
118  DBG_CHECK( npot2-index2-mult2 >= -1,"npot2-index2-mult2")
119 
120  if(cplex==1 .or. cpldot==1)then
121 
122    do i1=1,mult1
123      do i2=1,mult2
124        ar=zero
125        do ispden=1,min(nspden,2)
126 !$OMP PARALLEL DO PRIVATE(ifft) SHARED(cplex,i1,i2,index1,index2,ispden,nfft,potarr1,potarr2) REDUCTION(+:ar)
127          do ifft=1,cplex*nfft
128            ar=ar + potarr1(ifft,ispden,index1+i1-1)*potarr2(ifft,ispden,index2+i2-1)
129          end do
130        end do
131        dot(1,i1,i2)=ar
132        if (nspden==4) then
133          ar=zero
134          do ispden=3,4
135 !$OMP PARALLEL DO PRIVATE(ifft) SHARED(cplex,i1,i2,index1,index2,ispden,nfft,potarr1,potarr2) REDUCTION(+:ar)
136            do ifft=1,cplex*nfft
137              ar=ar + potarr1(ifft,ispden,index1+i1-1)*potarr2(ifft,ispden,index2+i2-1)
138            end do
139          end do
140          if (opt_storage==0) then
141            if (cplex==1) then
142              dot(1,i1,i2)=dot(1,i1,i2)+two*ar
143            else
144              dot(1,i1,i2)=dot(1,i1,i2)+ar
145            end if
146          else
147            dot(1,i1,i2)=half*(dot(1,i1,i2)+ar)
148          end if
149        end if
150      end do
151    end do
152 
153  else ! if (cplex==2 .and. cpldot==2)
154 
155    do i1=1,mult1
156      do i2=1,mult2
157        ar=zero ; ai=zero
158        do ispden=1,min(nspden,2)
159 !$OMP PARALLEL DO PRIVATE(ifft) SHARED(cplex,i1,i2,index1,index2,ispden,nfft,potarr1,potarr2) REDUCTION(+:ar,ai)
160          do ifft=1,nfft
161            ar=ar + potarr1(2*ifft-1,ispden,index1+i1-1)*potarr2(2*ifft-1,ispden,index2+i2-1) &
162 &           + potarr1(2*ifft  ,ispden,index1+i1-1)*potarr2(2*ifft  ,ispden,index2+i2-1)
163            ai=ai + potarr1(2*ifft-1,ispden,index1+i1-1)*potarr2(2*ifft  ,ispden,index2+i2-1) &
164 &           - potarr1(2*ifft  ,ispden,index1+i1-1)*potarr2(2*ifft-1,ispden,index2+i2-1)
165          end do
166        end do
167        dot(1,i1,i2)=ar ; dot(2,i1,i2)=ai
168        if (nspden==4) then
169          ar=zero
170          do ispden=3,4
171 !$OMP PARALLEL DO PRIVATE(ifft) SHARED(cplex,i1,i2,index1,index2,ispden,nfft,potarr1,potarr2) REDUCTION(+:ar,ai)
172            do ifft=1,nfft
173              ar=ar + potarr1(2*ifft-1,ispden,index1+i1-1)*potarr2(2*ifft-1,ispden,index2+i2-1) &
174 &             + potarr1(2*ifft  ,ispden,index1+i1-1)*potarr2(2*ifft  ,ispden,index2+i2-1)
175              ai=ai + potarr1(2*ifft-1,ispden,index1+i1-1)*potarr2(2*ifft  ,ispden,index2+i2-1) &
176 &             - potarr1(2*ifft  ,ispden,index1+i1-1)*potarr2(2*ifft-1,ispden,index2+i2-1)
177            end do
178          end do
179          if (opt_storage==0) then
180            dot(1,i1,i2)=dot(1,i1,i2)+ar
181            dot(2,i1,i2)=dot(2,i1,i2)+ai
182          else
183            dot(1,i1,i2)=half*(dot(1,i1,i2)+ar)
184            dot(2,i1,i2)=half*(dot(2,i1,i2)+ai)
185          end if
186        end if
187      end do
188    end do
189  end if
190 
191 !XG030513 : MPIWF reduction (addition) on dot is needed here
192  if (mpi_summarize) then
193    call timab(48,1,tsec)
194    call xmpi_sum(dot,mpicomm ,ierr)
195    call timab(48,2,tsec)
196  end if
197 
198  if(cpldot==2 .and. cplex==1)dot(2,:,:)=zero
199 
200 end subroutine dotprodm_v