TABLE OF CONTENTS


ABINIT/sqnormm_v [ Functions ]

[ Top ] [ Functions ]

NAME

 sqnormm_v

FUNCTION

 For a series of potentials,
 compute square of the norm (integral over FFT grid), to obtain
 a 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 density and potentials (nspden), and sum over them.
 Need the index of the first potential to be treated, in the provided 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 function on FFT grid is REAL, if 2, COMPLEX
  index=index of the first potential to be treated
  mpicomm=the mpi communicator used for the summation
  mpi_summarize=set it to .true. if parallelisation is done over FFT
  mult=number of potentials to be treated
  nfft= (effective) number of FFT grid points (for this processor)
  npot= third dimension of the potarr array
  nspden=number of spin-density components
  opt_storage: 0, if potential is stored as V^up-up, V^dn-dn, Re[V^up-dn], Im[V^up-dn]
               1, if potential is stored as V, B_x, B_y, Bz  (B=magn. field)
  potarr(cplex*nfft,nspden,npot)=array of real space potentials on FFT grid

OUTPUT

  norm2(mult)= value of the square of the norm of the different potentials

PARENTS

      scfcge,scfopt

CHILDREN

      timab,xmpi_sum

SOURCE

 48 #if defined HAVE_CONFIG_H
 49 #include "config.h"
 50 #endif
 51 
 52 #include "abi_common.h"
 53 
 54 
 55 subroutine sqnormm_v(cplex,index,mpicomm, mpi_summarize,mult,nfft,norm2,npot,nspden,opt_storage,potarr)
 56 
 57  use defs_basis
 58  use m_errors
 59  use m_profiling_abi
 60  use m_xmpi
 61 
 62 !This section has been created automatically by the script Abilint (TD).
 63 !Do not modify the following lines by hand.
 64 #undef ABI_FUNC
 65 #define ABI_FUNC 'sqnormm_v'
 66  use interfaces_18_timing
 67 !End of the abilint section
 68 
 69  implicit none
 70 
 71 !Arguments ------------------------------------
 72 !scalars
 73  integer,intent(in) :: cplex,index,mult,nfft,npot,nspden,opt_storage,mpicomm
 74  logical, intent(in) :: mpi_summarize
 75 !arrays
 76  real(dp),intent(in) :: potarr(cplex*nfft,nspden,npot)
 77  real(dp),intent(out) :: norm2(mult)
 78 
 79 !Local variables-------------------------------
 80 !scalars
 81  integer :: ierr,ifft,ii,ispden
 82  real(dp) :: ar
 83 !arrays
 84  real(dp) :: tsec(2)
 85 
 86 ! *************************************************************************
 87 
 88 !Real or complex inputs are coded
 89  DBG_CHECK(ANY(cplex==(/1,2/)),"Wrong cplex")
 90  DBG_CHECK(ANY(nspden==(/1,2,4/)),"Wrong nspden")
 91 
 92  DBG_CHECK(index>=1,"wrong index")
 93  DBG_CHECK(mult>=1,"wrong mult")
 94  DBG_CHECK(npot>=1,"wrong npot")
 95 
 96  DBG_CHECK(npot-index-mult>=-1,'npot-index-mult')
 97 
 98  do ii=1,mult
 99    ar=zero
100    do ispden=1,min(nspden,2)
101 !$OMP PARALLEL DO PRIVATE(ifft) SHARED(cplex,ii,index,ispden,nfft,potarr) REDUCTION(+:ar)
102      do ifft=1,cplex*nfft
103        ar=ar + potarr(ifft,ispden,index+ii-1)**2
104      end do
105    end do
106    norm2(ii)=ar
107    if (nspden==4) then
108      ar=zero
109      do ispden=3,4
110 !$OMP PARALLEL DO PRIVATE(ifft) SHARED(cplex,ii,index,ispden,nfft,potarr) REDUCTION(+:ar)
111        do ifft=1,cplex*nfft
112          ar=ar + potarr(ifft,ispden,index+ii-1)**2
113        end do
114      end do
115      if (opt_storage==0) then
116        if (cplex==1) then
117          norm2(ii)=norm2(ii)+two*ar
118        else
119          norm2(ii)=norm2(ii)+ar
120        end if
121      else
122        norm2(ii)=half*(norm2(ii)+ar)
123      end if
124    end if
125  end do
126 
127 !XG030513 : MPIWF reduction (addition) on norm2 is needed here
128  if (mpi_summarize) then
129    call timab(48,1,tsec)
130    call xmpi_sum(norm2,mpicomm ,ierr)
131    call timab(48,2,tsec)
132  end if
133 
134 end subroutine sqnormm_v