TABLE OF CONTENTS


ABINIT/m_rec_tools [ Modules ]

[ Top ] [ Modules ]

NAME

  m_rec_tools

FUNCTION

  This module provides some functions more or less generic used
  in the Recursion Mathod

COPYRIGHT

 Copyright (C) 2002-2024 ABINIT group (MMancini)
 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.

NOTES

SOURCE

20 #if defined HAVE_CONFIG_H
21 #include "config.h"
22 #endif
23 
24 #include "abi_common.h"
25 
26 module m_rec_tools
27 
28  use defs_basis
29  use m_abicore
30 
31  use defs_rectypes, only : recparall_type
32 
33  implicit none
34 
35  private
36 
37  public ::            &
38    get_pt0_pt1,       &  !--To get pt0 pt1 from inf,sup
39    reshape_pot,       &  !--To rescale the potential between 2 grids
40    trottersum            !--To calculate the trotter sum
41 
42 CONTAINS  !===========================================================

m_rec_tools/get_pt0_pt1 [ Functions ]

[ Top ] [ m_rec_tools ] [ Functions ]

NAME

 get_pt0_pt1

FUNCTION

 utility function to get pt0 and pt1 for given inf and sup

INPUTS

 ngfft(3) = fine grid (corresponds to dtset%ngfft(1:3))
 inf = inferior point in-line coordinate on the coarse grid
 sup = superior point in-line coordinate on the coarse grid

OUTPUT

 recpar%pt0<type(vec_int)>=Intial point for this proc in x,y,z
 recpar%pt1<type(vec_int)>=Final point for this proc in x,y,z
 recpar%min_pt=inferior point in-line coordinate on the fine grid
 recpar%max_pt=superior point in-line coordinate on the dine grid

SIDE EFFECTS

SOURCE

 67 subroutine get_pt0_pt1(ngfft,gratio,inf,sup,recpar)
 68 
 69  implicit none
 70 
 71 !Arguments ------------------------------------
 72   integer,intent(in)        :: gratio
 73   integer,intent(in)        :: inf,sup
 74   type(recparall_type),intent(inout) :: recpar
 75   integer,intent(in)        :: ngfft(3)
 76 !Local ---------------------------
 77   integer :: x,y,z,count
 78   integer :: boz,boy,pt
 79 ! *********************************************************************
 80   count = 0
 81   do z = 0,ngfft(3)-1,gratio
 82     boz = z*ngfft(2)
 83     do y = 0,ngfft(2)-1,gratio
 84       boy = (boz+y)*ngfft(1)
 85       do x = 0,ngfft(1)-1,gratio
 86         pt = boy+x
 87         if(count >= inf) then
 88           if(count == inf) then
 89             recpar%pt0%x = x; recpar%pt0%y = y; recpar%pt0%z = z
 90             recpar%min_pt = pt
 91           end if
 92           recpar%pt1%x = x; recpar%pt1%y = y; recpar%pt1%z = z
 93           recpar%max_pt = pt
 94         endif
 95         count = count+1
 96         if(count == sup ) return
 97       end do
 98     end do
 99   end do
100 end subroutine get_pt0_pt1

m_rec_tools/reshape_pot [ Functions ]

[ Top ] [ m_rec_tools ] [ Functions ]

NAME

 reshape_pot

FUNCTION

 Reshape array on

INPUTS

 nfft=size of the big grid
 nfftrec=size of the cut grid
 trasl(3) center point
 ngfft(3) dimensions of the big grid
 ngfftrec(3) dimensions of the cut grid
 pot(0:nfft-1) 3d-array on the big grid

OUTPUT

 potloc is a cut of pot around trasl

SOURCE

123 subroutine reshape_pot(trasl,nfft,nfftrec,ngfft,ngfftrec,pot,potloc)
124 
125  implicit none
126 
127  !Arguments ------------------------------------
128  integer,  intent(in) :: nfft,nfftrec
129  integer,  intent(in) :: trasl(3)
130  integer,  intent(in) :: ngfftrec(3),ngfft(3)
131  real(dp), intent(in) :: pot(0:nfft-1)
132  real(dp), intent(out):: potloc(0:nfftrec-1)
133  !Local ----------------------------------------
134  ! scalars
135  integer :: zz,yy,xx,parz,pary
136  integer :: modi,modj,modk
137  !character(len=500) :: msg
138  ! *********************************************************************
139  do zz = 0,ngfftrec(3)-1
140    modk = modulo(zz-trasl(3),ngfft(3))*ngfft(2)
141    parz = ngfftrec(2)*zz
142    do yy = 0,ngfftrec(2)-1
143      modj = (modulo(yy-trasl(2),ngfft(2))+modk)*ngfft(1)
144      pary = ngfftrec(1)*(yy+parz)
145      do xx = 0,ngfftrec(1)-1
146        modi = modulo(xx-trasl(1),ngfft(1))+modj
147        potloc(xx+pary) = pot(modi)
148      end do
149    end do
150  end do
151 
152 end subroutine reshape_pot

m_rec_tools/trottersum [ Functions ]

[ Top ] [ m_rec_tools ] [ Functions ]

NAME

 trottersum

FUNCTION

 Calculate the contribution to the partial fraction decomposition
 due to a recursion step.

INPUTS

  dim_trott=dimension of the partial fraction decomposition (PFD)
  pi_on_rtrotter=parameter pi/rtrotter
  an,bn2=recursion coefficients at irec
  exp1=numerical factor  (see density_rec)
  coeef_mu=numerical factor

OUTPUT

 SIZE EFFECTS
  D,N=denominator and numerator accumalator of PFD
  Dold,Nold=denominator and numerator of PFD (old values)
  facrec0=used to select irec=0
  error=estimated error of recursion at this step
  prod_b2=numerical factor

SOURCE

181 subroutine trottersum(dim_trott,error,&
182      &                prod_b2,pi_on_rtrotter,&
183      &                facrec0,coeef_mu,exp1,&
184      &                an,bn2,&
185      &                N,D,Nold,Dold)
186 
187  implicit none
188 
189  !Arguments ------------------------------------
190  !scalars
191  integer,  intent(in) :: dim_trott
192  real(dp), intent(in) :: an,bn2,exp1,pi_on_rtrotter
193  real(dp), intent(inout) :: error,prod_b2
194  complex(dp), intent(in) :: coeef_mu
195  complex(dp), intent(inout) :: facrec0
196  !arrays
197  complex(dpc),intent(inout) :: D(0:dim_trott),Dold(0:dim_trott)
198  complex(dpc),intent(inout) :: N(0:dim_trott),Nold(0:dim_trott)
199  !Local ----------------------------------------
200  ! scalars
201  integer :: itrot
202  real(dp) :: arg
203  complex(dpc) :: Dnew,Nnew,zj
204  !character(len=500) :: msg
205  ! *********************************************************************
206 
207  error = zero
208  prod_b2 = prod_b2 * exp1 * bn2
209 
210  do itrot=0,dim_trott
211    arg = pi_on_rtrotter*(real( itrot,dp) + half )
212    zj = cmplx(cos(arg),sin(arg),dp)*coeef_mu
213 
214    Nnew = zj*facrec0 + &
215         & (zj - cmplx(an  ,zero,dp))*N(itrot) - &
216         &       cmplx(bn2 ,zero,dp)*Nold(itrot)
217    Dnew = (zj - cmplx(an  ,zero,dp))*D(itrot) - &
218         &       cmplx(bn2 ,zero,dp)*Dold(itrot)
219    Nold(itrot) = N(itrot)
220    Dold(itrot) = D(itrot)
221    N(itrot) = Nnew
222    D(itrot) = Dnew
223 
224    !--Error estimator
225    error = error + abs(prod_b2/(D(itrot)*Dold(itrot)))
226  end do
227  facrec0 = czero
228 
229 end subroutine trottersum