TABLE OF CONTENTS


ABINIT/m_qplusg [ Modules ]

[ Top ] [ Modules ]

NAME

  m_qplusg

FUNCTION

COPYRIGHT

 Copyright (C) 1999-2024 ABINIT group (MG, FB, BG)
 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 .

SOURCE

16 #if defined HAVE_CONFIG_H
17 #include "config.h"
18 #endif
19 
20 #include "abi_common.h"
21 
22 module m_qplusg
23 
24  use defs_basis
25 
26  implicit none
27 
28  private

m_qplusg/cmod_qpg [ Functions ]

[ Top ] [ m_qplusg ] [ Functions ]

NAME

 cmod_qpg

FUNCTION

 Set up table of lengths |q+G| for the Coulomb potential.

INPUTS

 gvec(3,npwvec)=Reduced coordinates of the G vectors.
 gprimd(3,3)=Dimensional primitive translations for reciprocal space ($\textrm{bohr}^{-1}$)
 iq=Index specifying the q point.
 npwvec=Number of planewaves
 nq=Number of q points.
 q=Coordinates of q points.

OUTPUT

 qplusg(npwvec)=Norm of q+G vector

SOURCE

57 subroutine cmod_qpg(nq, iq, q, npwvec, gvec, gprimd, qplusg)
58 
59 !Arguments ------------------------------------
60 !scalars
61  integer,intent(in) :: iq,npwvec,nq
62 !arrays
63  integer,intent(in) :: gvec(3,npwvec)
64  real(dp),intent(in) :: gprimd(3,3),q(3,nq)
65  real(dp),intent(out) :: qplusg(npwvec)
66 
67 !Local variables ------------------------------
68 !scalars
69  integer :: ig,ii
70 !arrays
71  real(dp) :: gmet(3,3),gpq(3)
72 
73 !************************************************************************
74 
75  ! Compute reciprocal space metrics
76  do ii=1,3
77    gmet(ii,:)=gprimd(1,ii)*gprimd(1,:)+&
78               gprimd(2,ii)*gprimd(2,:)+&
79               gprimd(3,ii)*gprimd(3,:)
80  end do
81 
82  if (ALL(ABS(q(:,iq)) < tol3)) then
83    ! Treat q as if it were zero except when G=0
84    qplusg(1)=two_pi*SQRT(DOT_PRODUCT(q(:,iq),MATMUL(gmet,q(:,iq))))
85    do ig=2,npwvec
86      gpq(:)=gvec(:,ig)
87      qplusg(ig)=two_pi*SQRT(DOT_PRODUCT(gpq,MATMUL(gmet,gpq)))
88    end do
89  else
90    do ig=1,npwvec
91      gpq(:)=gvec(:,ig)+q(:,iq)
92      qplusg(ig)=two_pi*SQRT(DOT_PRODUCT(gpq,MATMUL(gmet,gpq)))
93    end do
94  end if
95 
96 end subroutine cmod_qpg