The Elk Code
mixerifc.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine mixerifc(mtype,n,v,dv,nwork,work)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: mtype,n
11 real(8), intent(inout) :: v(n)
12 real(8), intent(out) :: dv
13 integer, intent(inout) :: nwork
14 real(8), intent(inout) :: work(*)
15 select case(mtype)
16 case(0)
17 ! linear mixing
18  if (nwork < 1) then
19  nwork=n
20  return
21  end if
22  call mixlinear(iscl,amixpm(1),n,v,work,dv)
23 case(1)
24 ! adaptive linear mixing
25  if (nwork < 1) then
26  nwork=3*n
27  return
28  end if
29  call mixadapt(iscl,amixpm(1),amixpm(2),n,v,work,work(n+1),work(2*n+1),dv)
30 case(3)
31 ! Broyden mixing
32  if (nwork < 1) then
33  nwork=(4+2*mixsdb)*n+mixsdb**2
34  return
35  end if
36  call mixbroyden(iscl,n,mixsdb,broydpm(1),broydpm(2),v,work,work(2*n+1), &
37  work(4*n+1),work((4+mixsdb)*n+1),work((4+2*mixsdb)*n+1),dv)
38 case default
39  write(*,*)
40  write(*,'("Error(mixerifc): mtype not defined : ",I8)') mtype
41  write(*,*)
42  stop
43 end select
44 end subroutine
45 
46 subroutine getmixdata(mtype,mixdescr)
47 implicit none
48 ! arguments
49 integer, intent(in) :: mtype
50 character(*), intent(out) :: mixdescr
51 select case(mtype)
52 case(0)
53  mixdescr='Linear mixing'
54 case(1)
55  mixdescr='Adaptive linear mixing'
56 case(3)
57  mixdescr='Broyden mixing, J. Phys. A: Math. Gen. 17, L317 (1984)'
58 case default
59  write(*,*)
60  write(*,'("Error(getmixdata): mixtype not defined : ",I8)') mtype
61  write(*,*)
62  stop
63 end select
64 end subroutine
65 
real(8), dimension(2) amixpm
Definition: modmain.f90:702
integer iscl
Definition: modmain.f90:1051
pure subroutine mixadapt(iscl, beta0, betamax, n, nu, mu, beta, f, d)
Definition: mixadapt.f90:10
subroutine getmixdata(mtype, mixdescr)
Definition: mixerifc.f90:47
pure subroutine mixlinear(iscl, beta, n, nu, mu, d)
Definition: mixlinear.f90:7
real(8), dimension(2) broydpm
Definition: modmain.f90:706
integer mixsdb
Definition: modmain.f90:704
subroutine mixerifc(mtype, n, v, dv, nwork, work)
Definition: mixerifc.f90:7
subroutine mixbroyden(iscl, n, msd, alpha, w0, nu, mu, f, df, u, a, d)
Definition: mixbroyden.f90:7