TABLE OF CONTENTS
ABINIT/inspln [ Functions ]
NAME
inspln
FUNCTION
This procedure gives the values of the spline coefficients (second derivatives) in the 1D grid with periodic boundary conditions at rsid - the values of the unknown functions specified in the vector valf of direction idir
COPYRIGHT
Copyright (C) 2002-2017 ABINIT group (PCasek,FF,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
idir= direction following which the derivatives are evaluated snn, tnn=remaining bi-dimensional coordinates of the line along which the derivative is to be computed
OUTPUT
(see side effects)
SIDE EFFECTS
This routine works on the data contained in the aimfields module WARNING This file does not follow the ABINIT coding rules (yet)
PARENTS
initaim
CHILDREN
SOURCE
40 #if defined HAVE_CONFIG_H 41 #include "config.h" 42 #endif 43 44 #include "abi_common.h" 45 46 47 subroutine inspln(idir,snn,tnn) 48 49 use m_profiling_abi 50 51 use defs_basis 52 use defs_aimfields 53 54 !This section has been created automatically by the script Abilint (TD). 55 !Do not modify the following lines by hand. 56 #undef ABI_FUNC 57 #define ABI_FUNC 'inspln' 58 !End of the abilint section 59 60 implicit none 61 62 !Arguments ------------------------------------ 63 !scalars 64 integer,intent(in) :: idir,snn,tnn 65 66 !Local variables------------------------------- 67 !scalars 68 integer :: dim,ii 69 real(dp) :: ss 70 !arrays 71 real(dp) :: rsid(ngfft(idir)),valf(ngfft(idir)) 72 real(dp),pointer :: ptc(:),ptd(:),ptp(:) 73 74 ! ************************************************************************* 75 76 !POINTER INITIALIZATION 77 78 if (idir==1) then 79 valf(:)=dvl(:,snn,tnn) 80 elseif (idir==2) then 81 valf(:)=dvl(tnn,:,snn) 82 else 83 valf(:)=dvl(snn,tnn,:) 84 end if 85 86 nullify(ptd,ptc,ptp) 87 if(idir==1) then 88 ptd=>dig1;ptc=>cdig1;ptp=>llg1 89 elseif (idir==2) then 90 ptd=>dig2;ptc=>cdig2;ptp=>llg2 91 else 92 ptd=>dig3;ptc=>cdig3;ptp=>llg3 93 end if 94 95 dim=ngfft(idir) 96 97 !FIRST CYCLE OF RECURRENCE 98 99 rsid(1)=valf(2)+valf(dim)-2.*valf(1) 100 rsid(1)=rsid(1)/ptd(1) 101 do ii=2,dim-1 102 rsid(ii)=valf(ii+1)+valf(ii-1)-2.*valf(ii) 103 rsid(ii)=(rsid(ii)-ptc(ii-1)*rsid(ii-1))/ptd(ii) 104 end do 105 ss=0._dp 106 do ii=1,dim-1 107 ss=ss+rsid(ii)*ptp(ii) 108 end do 109 rsid(dim)=valf(1)+valf(dim-1)-2.*valf(dim) 110 rsid(dim)=(rsid(dim)-ss)/ptd(dim) 111 112 !SECOND CYCLE WITH TRANSPOSED MATRIX 113 114 rsid(dim)=rsid(dim)/ptd(dim) 115 rsid(dim-1)=(rsid(dim-1)-ptc(dim-1)*rsid(dim))/ptd(dim-1) 116 do ii=dim-2,1,-1 117 rsid(ii)=(rsid(ii)-ptc(ii)*rsid(ii+1)-ptp(ii)*rsid(dim))/ptd(ii) 118 end do 119 120 if (idir==1) then 121 ddx(:,snn,tnn)=rsid(:) 122 elseif (idir==2) then 123 ddy(tnn,:,snn)=rsid(:) 124 else 125 ddz(snn,tnn,:)=rsid(:) 126 end if 127 128 end subroutine inspln