The Elk Code
mixbroyden.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2011 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 subroutine mixbroyden(iscl,n,msd,alpha,w0,nu,mu,f,df,u,a,d)
7 use modomp
8 implicit none
9 ! arguments
10 integer, intent(in) :: iscl,n,msd
11 real(8), intent(in) :: alpha,w0
12 real(8), intent(inout) :: nu(n),mu(n,2)
13 real(8), intent(inout) :: f(n,2),df(n,msd)
14 real(8), intent(inout) :: u(n,msd),a(msd,msd)
15 real(8), intent(out) :: d
16 ! local variables
17 integer jc,kp,kc,k,m
18 integer info,nthd
19 real(8) t1
20 ! automatic arrays
21 integer ipiv(msd)
22 real(8) c(msd),beta(msd,msd),gamma(msd),work(msd)
23 ! external functions
24 real(8), external :: ddot
25 if (n < 1) then
26  write(*,*)
27  write(*,'("Error(mixbroyden): n < 1 : ",I8)') n
28  write(*,*)
29  stop
30 end if
31 if (msd < 2) then
32  write(*,*)
33  write(*,'("Error(mixbroyden): msd < 2 : ",I8)') msd
34  write(*,*)
35  stop
36 end if
37 ! initialise mixer
38 if (iscl < 1) then
39  mu(1:n,1)=nu(1:n)
40  mu(1:n,2)=nu(1:n)
41  f(1:n,1)=0.d0
42  df(1:n,1)=0.d0
43  u(1:n,1)=0.d0
44  a(1:msd,1:msd)=0.d0
45  d=1.d0
46  return
47 end if
48 ! current subspace dimension
49 m=min(iscl+1,msd)
50 ! current index modulo m
51 jc=mod(iscl,m)+1
52 ! previous index modulo 2
53 kp=mod(iscl-1,2)+1
54 ! current index modulo 2
55 kc=mod(iscl,2)+1
56 f(1:n,kc)=nu(1:n)-mu(1:n,kp)
57 d=sum(f(1:n,kc)**2)
58 d=sqrt(d/dble(n))
59 df(1:n,jc)=f(1:n,kc)-f(1:n,kp)
60 t1=norm2(df(1:n,jc))
61 if (t1 > 1.d-8) t1=1.d0/t1
62 df(1:n,jc)=t1*df(1:n,jc)
63 u(1:n,jc)=alpha*df(1:n,jc)+t1*(mu(1:n,kp)-mu(1:n,kc))
64 call holdthd(2*m,nthd)
65 !$OMP PARALLEL DEFAULT(SHARED) &
66 !$OMP NUM_THREADS(nthd)
67 !$OMP DO SCHEDULE(DYNAMIC)
68 do k=1,m
69  c(k)=ddot(n,df(1,k),1,f(1,kc),1)
70 end do
71 !$OMP END DO NOWAIT
72 !$OMP DO SCHEDULE(DYNAMIC)
73 do k=1,m
74  a(k,jc)=ddot(n,df(1,k),1,df(1,jc),1)
75  a(jc,k)=a(k,jc)
76 end do
77 !$OMP END DO
78 !$OMP END PARALLEL
79 call freethd(nthd)
80 beta(1:msd,1:msd)=a(1:msd,1:msd)
81 do k=1,m
82  beta(k,k)=beta(k,k)+w0**2
83 end do
84 ! invert beta
85 call dgetrf(m,m,beta,msd,ipiv,info)
86 if (info == 0) call dgetri(m,beta,msd,ipiv,work,m,info)
87 if (info /= 0) then
88  write(*,*)
89  write(*,'("Error(mixbroyden): could not invert matrix")')
90  write(*,*)
91  stop
92 end if
93 do k=1,m
94  gamma(k)=sum(c(1:m)*beta(1:m,k))
95 end do
96 nu(1:n)=mu(1:n,kp)+alpha*f(1:n,kc)
97 do k=1,m
98  t1=-gamma(k)
99  nu(1:n)=nu(1:n)+t1*u(1:n,k)
100 end do
101 mu(1:n,kc)=nu(1:n)
102 end subroutine
103 
Definition: modomp.f90:6
subroutine mixbroyden(iscl, n, msd, alpha, w0, nu, mu, f, df, u, a, d)
Definition: mixbroyden.f90:7
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78