TABLE OF CONTENTS


ABINIT/tools_lotf [ Modules ]

[ Top ] [ Modules ]

NAME

 tools_lotf

FUNCTION

  Contains simple functions for LOTF

COPYRIGHT

 Copyright (C) 2005-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 .

SOURCE

15 #if defined HAVE_CONFIG_H
16 #include "config.h"
17 #endif
18 
19 #include "abi_common.h"
20 
21 module tools_lotf
22 
23  use defs_basis
24  use m_errors
25 
26  implicit none
27 
28  private
29 
30  public ::             &
31    pinterp,            &
32    pinterp_nolinear,   &
33    dlvsum,             &
34    icf
35 
36 
37 contains

tools_lotf/dlvsum [ Functions ]

[ Top ] [ tools_lotf ] [ Functions ]

NAME

 dlvsum

FUNCTION

INPUTS

SOURCE

115  subroutine dlvsum(local_pe,npes,sum_local_pe,ndimve)
116   !------------------------------------------------------------
117   ! (ADV) GLOBAL VECTOR SUM, REAL 28/5/94
118   implicit none
119 
120   !Arguments ------------------------
121   integer :: local_pe,npes,ndimve
122   real(dp) :: sum_local_pe(ndimve)
123 
124 ! *************************************************************************
125 
126   ! do nothing here !!! 
127   ! better remove dlvsum  everywhere in this version.
128   return
129  end subroutine dlvsum

tools_lotf/icf [ Functions ]

[ Top ] [ tools_lotf ] [ Functions ]

NAME

 icf

FUNCTION

INPUTS

SOURCE

141  FUNCTION icf(ix,iy,iz,ic1,ic2,ic3) 
142 
143   implicit none
144   !Arguments ------------------------
145   integer :: ix,iy,iz,ic1,ic2,ic3,icf  
146   icf = 1 &
147     + mod( ix - 1 + ic1, ic1 ) &
148     + mod( iy - 1 + ic2, ic2 ) * ic1 &
149     + mod( iz - 1 + ic3, ic3 ) * ic1 * ic2
150  end FUNCTION icf
151 
152 end module tools_lotf

tools_lotf/pinterp [ Functions ]

[ Top ] [ tools_lotf ] [ Functions ]

NAME

 pinterp

FUNCTION

INPUTS

SOURCE

48  subroutine  pinterp(a0,a1,ainterpoled,ndim,nitex,n)
49   ! interpolator : a0 for n=0
50   !                a1 for n=nx
51   !                a1 for nx = 0
52   implicit none
53 
54   !Arguments ------------------------
55   integer,intent(in) :: ndim,nitex,n 
56   real(dp),intent(in) :: a0(ndim),a1(ndim)
57   real(dp),intent(out) :: ainterpoled(ndim)
58   !Local --------------------------- 
59   character(len=500) :: message
60 
61 ! *************************************************************************
62 
63   !--Control
64   if(nitex < 0) then
65     write(message,'(a,i4,a)') '  LOTF: pinterp nitex =',nitex,'smaller than 0'
66     ABI_ERROR(message)
67 
68   elseif(nitex >= 1) then
69     ainterpoled = a0 + (n/real(nitex,dp))*(a1-a0) 
70   elseif(nitex == 0) then
71     ainterpoled = a1 
72   endif
73  end subroutine pinterp

tools_lotf/pinterp_nolinear [ Functions ]

[ Top ] [ tools_lotf ] [ Functions ]

NAME

 pinterp_nolinear

FUNCTION

INPUTS

SOURCE

 85  subroutine  pinterp_nolinear(a0,a1,ainterpoled,ndim,nitex,n)
 86   ! interpolator : a0 for n=0
 87   !                a1 for n=nx
 88   !                a1 for nx = 0
 89   implicit none
 90 
 91   !Arguments ------------------------
 92   integer,intent(in) :: ndim,nitex,n 
 93   real(dp),intent(in) :: a0(ndim),a1(ndim)
 94   real(dp),intent(out) :: ainterpoled(ndim)
 95   !Local --------------------------- 
 96   real(dp) :: lambda
 97 
 98 ! *************************************************************************
 99 
100   !--Control
101   lambda = n/real(nitex,dp)
102   ainterpoled = a0 - six*(a1-a0)*(lambda**3/three-lambda**2/two) 
103  end subroutine pinterp_nolinear