TABLE OF CONTENTS


ABINIT/spline_paw_fncs [ Functions ]

[ Top ] [ Functions ]

NAME

 spline_paw_fncs

FUNCTION

 Compute radial PAW functions and their derivatives on a set of points in the PAW sphere.

COPYRIGHT

 Copyright (C) 2005-2018 ABINIT group (JJ,MT)
 This file is distributed under the terms of the
 GNU General Public License, see ~ABINIT/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

 integer :: nnl : number of nl PAW basis functions in set
 integer :: npts : number of points to perform fits on
 real(dp) :: points(npts) : SORTED vector of points to perform fits on
 type(pawrad_type) :: pawrad : paw radial mesh data
 type(pawtab_type) :: pawtab : paw wavefunctions around each type of atom

OUTPUT

 real(dp) :: phi(npts,nnl), dphi(npts,nnl), tphi(npts,nnl), dtphi(npts,nnl) : PAW functions
             phi, tphi and their radial derivatives evaluated at the input points by spline fits

NOTES

 The PAW basis functions are defined by $<r|\phi_i>=(u_i(r)/r)S_{lm}(\hat{r})$, evaluated on a radial
 grid. This subroutine computes $u(r)$ and $d u(r)/dr$ on the set of points provided on input. Typically
 the input points will be the fine grid points in the PAW sphere. They are presumed to be sorted
 already on input.

PARENTS

CHILDREN

      nderiv_gen,spline,splint

SOURCE

 39 #if defined HAVE_CONFIG_H
 40 #include "config.h"
 41 #endif
 42 
 43 #include "abi_common.h"
 44 
 45 subroutine spline_paw_fncs(dphi,dtphi,nnl,npts,pawrad,pawtab,points,phi,tphi)
 46 
 47  use m_profiling_abi
 48 
 49  use defs_basis
 50  use m_errors
 51  use m_splines
 52  use m_pawrad, only : pawrad_type, nderiv_gen
 53  use m_pawtab, only : pawtab_type
 54 
 55 !This section has been created automatically by the script Abilint (TD).
 56 !Do not modify the following lines by hand.
 57 #undef ABI_FUNC
 58 #define ABI_FUNC 'spline_paw_fncs'
 59 !End of the abilint section
 60 
 61  implicit none
 62 
 63 !Arguments ------------------------------------
 64 !scalars
 65  integer,intent(in) :: nnl,npts
 66 !arrays
 67  real(dp),intent(in) :: points(npts)
 68  real(dp),intent(out) :: phi(npts,nnl),dphi(npts,nnl),tphi(npts,nnl),dtphi(npts,nnl)
 69  type(pawrad_type),intent(in) :: pawrad
 70  type(pawtab_type),intent(in) :: pawtab
 71 
 72 !Local variables-------------------------------
 73 !scalars
 74  integer :: inl
 75  real(dp) :: ybcbeg, ybcend
 76 !arrays
 77  real(dp),allocatable :: der(:),diag(:),ypp(:)
 78 
 79 ! ************************************************************************
 80 
 81  DBG_ENTER("COLL")
 82 
 83  ABI_ALLOCATE(der,(pawtab%mesh_size))
 84  ABI_ALLOCATE(ypp,(pawtab%mesh_size))
 85  ABI_ALLOCATE(diag,(pawtab%mesh_size))
 86  do inl = 1, nnl
 87 
 88 !  spline phi onto points
 89    ypp(:) = zero; diag(:) = zero; ybcbeg = zero; ybcend = zero
 90    call spline(pawrad%rad,pawtab%phi(:,inl),pawtab%mesh_size,ybcbeg,ybcend,ypp)
 91    call splint(pawtab%mesh_size,pawrad%rad,pawtab%phi(:,inl),ypp,npts,points,phi(:,inl))
 92 
 93 !  next spline d phi/dr onto points
 94 !  need derivative of phi with respect to radius
 95    der(:) = zero
 96    call nderiv_gen(der,pawtab%phi(:,inl),pawrad)
 97    ypp(:) = zero; diag(:) = zero; ybcbeg = zero; ybcend = zero
 98    call spline(pawrad%rad,der,pawtab%mesh_size,ybcbeg,ybcend,ypp)
 99    call splint(pawtab%mesh_size,pawrad%rad,der,ypp,npts,points,dphi(:,inl))
100 
101 !  next splint tphi onto points
102    ypp(:) = zero; diag(:) = zero;
103    call spline(pawrad%rad,pawtab%tphi(:,inl),pawtab%mesh_size,ybcbeg,ybcend,ypp)
104    call splint(pawtab%mesh_size,pawrad%rad,pawtab%tphi(:,inl),ypp,npts,points,tphi(:,inl))
105 
106 !  finally spline d tphi/dr onto points
107 !  need derivative of tphi with respect to radius
108    der(:) = zero
109    call nderiv_gen(der,pawtab%tphi(:,inl),pawrad)
110    ypp(:) = zero; diag(:) = zero; ybcbeg = zero; ybcend = zero
111    call spline(pawrad%rad,der,pawtab%mesh_size,ybcbeg,ybcend,ypp)
112    call splint(pawtab%mesh_size,pawrad%rad,der,ypp,npts,points,dtphi(:,inl))
113 
114  end do ! end loop over nnl basis functions
115 
116  ABI_DEALLOCATE(der)
117  ABI_DEALLOCATE(ypp)
118  ABI_DEALLOCATE(diag)
119 
120  DBG_EXIT("COLL")
121 
122  end subroutine spline_paw_fncs