TABLE OF CONTENTS


ABINIT/gg1cc [ Functions ]

[ Top ] [ Functions ]

NAME

 gg1cc

FUNCTION

 gg1cc_xx=$(\frac{\sin(2\pi xx)}{(2\pi xx)(1-4xx^2)(1-xx^2)})^2$

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (XG, DCA, MM)
 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

  xx= abscisse to which gg1cc_xx is calculated

OUTPUT

  gg1cc_xx= gg1cc_x(xx)

PARENTS

      psp1cc

CHILDREN

SOURCE

29 #if defined HAVE_CONFIG_H
30 #include "config.h"
31 #endif
32 
33 #include "abi_common.h"
34 
35 
36 subroutine gg1cc(gg1cc_xx,xx)
37 
38  use defs_basis
39  use m_profiling_abi
40 
41 !This section has been created automatically by the script Abilint (TD).
42 !Do not modify the following lines by hand.
43 #undef ABI_FUNC
44 #define ABI_FUNC 'gg1cc'
45 !End of the abilint section
46 
47  implicit none
48 
49 !Arguments ------------------------------------
50 !scalars
51  real(dp),intent(in) :: xx
52  real(dp),intent(out) :: gg1cc_xx
53 
54 !Local variables -------------------------------------------
55 !The c s are coefficients for Taylor expansion of the analytic
56 !form near xx=0, 1/2, and 1.
57 !scalars
58  real(dp) :: c21=4.d0/9.d0,c22=-40.d0/27.d0,c23=20.d0/3.d0-16.d0*pi**2/27.d0
59  real(dp) :: c24=-4160.d0/243.d0+160.d0*pi**2/81.d0,c31=1.d0/36.d0
60  real(dp) :: c32=-25.d0/108.d0,c33=485.d0/432.d0-pi**2/27.d0
61  real(dp) :: c34=-4055.d0/972.d0+25.d0*pi**2/81.d0
62 
63 ! *************************************************************************
64 
65 !Cut off beyond 3/gcut=xcccrc
66  if (xx>3.0d0) then
67    gg1cc_xx=0.0d0
68 !  Take care of difficult limits near x=0, 1/2, and 1
69  else if (abs(xx)<=1.d-09) then
70    gg1cc_xx=1.d0
71  else if (abs(xx-0.5d0)<=1.d-04) then
72 !  (this limit and next are more troublesome for numerical cancellation)
73    gg1cc_xx=c21+(xx-0.5d0)*(c22+(xx-0.5d0)*(c23+(xx-0.5d0)*c24))
74  else if (abs(xx-1.d0)<=1.d-04) then
75    gg1cc_xx=c31+(xx-1.0d0)*(c32+(xx-1.0d0)*(c33+(xx-1.0d0)*c34))
76  else
77 !  The following is the square of the Fourier transform of a
78 !  function built out of two spherical bessel functions in G
79 !  space and cut off absolutely beyond gcut
80    gg1cc_xx=(sin(2.0d0*pi*xx)/( (2.0d0*pi*xx) * &
81 &   (1.d0-4.0d0*xx**2)*(1.d0-xx**2) )  )**2
82  end if
83 
84 end subroutine gg1cc