TABLE OF CONTENTS


ABINIT/bracketing [ Functions ]

[ Top ] [ Functions ]

NAME

 bracketing

FUNCTION

 bracket a minimun of a function f

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, 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 .
 For the initials of contributors, see ~ABINIT/Infos/contributors .

INPUTS

 dp_dum_vdp: the function of which the mimimum should be bracketted

OUTPUT

 b= last member of the bracketing triplet a < x < b
 fa,fx,fb= value of the function at dp_dum_vdp(v(:)+y*grad(:))

SIDE EFFECTS

 v: the initial vector for the function (return unchanged)
 grad: the direction on which the bracketting is to be performed (return unchanged)
 a,x: two members of the bracketing triplet (see b)

PARENTS

      linmin

CHILDREN

SOURCE

 37 #if defined HAVE_CONFIG_H
 38 #include "config.h"
 39 #endif
 40 
 41 #include "abi_common.h"
 42 
 43 
 44 subroutine bracketing (nv1,nv2,dp_dum_v2dp,v,grad,a,x,b,fa,fx,fb)
 45 
 46  use defs_basis
 47  use m_profiling_abi
 48 
 49 !This section has been created automatically by the script Abilint (TD).
 50 !Do not modify the following lines by hand.
 51 #undef ABI_FUNC
 52 #define ABI_FUNC 'bracketing'
 53 !End of the abilint section
 54 
 55  implicit none
 56 
 57 !Arguments ------------------------------------
 58  include "dummy_functions.inc"
 59 !scalars
 60  integer,intent(in) :: nv1,nv2
 61  real(dp),intent(inout) :: a,x
 62  real(dp),intent(out) :: b,fa,fb,fx
 63 !arrays
 64  real(dp),intent(inout) :: grad(nv1,nv2),v(nv1,nv2)
 65 
 66 !Local variables-------------------------------
 67 !scalars
 68  real(dp),parameter :: maglimit=10000.0_dp
 69  real(dp) :: c,fu,q,r,u,ulim
 70 
 71 ! *************************************************************************
 72 
 73  fa=dp_dum_v2dp(nv1,nv2,v(:,:)+(a*grad(:,:)))
 74  fx=dp_dum_v2dp(nv1,nv2,(x*grad(:,:))+v(:,:))
 75  if(fx > fa) then
 76   c=a
 77   a=x
 78   x=c
 79   c=fa
 80   fa=fx
 81   fx=c
 82  end if
 83  b=x+gold*(x-a)
 84  fb=dp_dum_v2dp(nv1,nv2,(b*grad(:,:))+v(:,:))
 85  do
 86   if (fx <= fb) return
 87   r=(x-a)*(fx-fb)
 88   q=(x-b)*(fx-fa)
 89   u=x-((x-b)*q-(x-a)*r)/(two*sign(max(abs(q-r),smallest_real),q-r))
 90   ulim=x+maglimit*(b-x)
 91   if((x-u)*(u-b) > zero) then
 92    fu=dp_dum_v2dp(nv1,nv2,(u*grad(:,:))+v(:,:))
 93    if(fu < fb) then
 94     a=x
 95     fa=fx
 96     x=u
 97     fx=fu
 98     return
 99    else if (fx < fu) then
100     b=u
101     fb=fu
102     return
103    end if
104    u=b+gold*(b-x)
105    fu=dp_dum_v2dp(nv1,nv2,(u*grad(:,:))+v(:,:))
106   else if((b-u)*(u-ulim) > zero) then
107    fu=dp_dum_v2dp(nv1,nv2,u*grad(:,:)+v(:,:))
108    if(fu<fb) then
109     x=b
110     b=u
111     u=b+gold*(b-x)
112     fx=fb
113     fb=fu
114     fu=dp_dum_v2dp(nv1,nv2,(u*grad(:,:))+v(:,:))
115    end if
116   else if((u-ulim)*(ulim-b) >= zero) then
117    u=ulim
118    fu=dp_dum_v2dp(nv1,nv2,(u*grad(:,:))+v(:,:))
119   else
120    u=b+gold*(b-x)
121    fu=dp_dum_v2dp(nv1,nv2,(u*grad(:,:))+v(:,:))
122   end if
123   a=x
124   x=b
125   b=u
126   fa=fx
127   fx=fb
128   fb=fu
129  end do
130 
131 end subroutine bracketing