The Elk Code
mixadapt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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 !BOP
7 ! !ROUTINE: mixadapt
8 ! !INTERFACE:
9 pure subroutine mixadapt(iscl,beta0,betamax,n,nu,mu,beta,f,d)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! iscl : self-consistent loop number (in,integer)
12 ! beta0 : mixing parameter (in,real)
13 ! betamax : maximum mixing parameter (in,real)
14 ! n : vector length (in,integer)
15 ! nu : current output vector as well as the next input vector of the
16 ! self-consistent loop (inout,real(n))
17 ! mu : used internally (inout,real(n))
18 ! beta : used internally (inout,real(n))
19 ! f : used internally (inout,real(n))
20 ! d : RMS difference between old and new output vectors (out,real)
21 ! !DESCRIPTION:
22 ! Given the input vector $\mu^i$ and output vector $\nu^i$ of the $i$th
23 ! self-consistent loop, this routine generates the next input vector to the
24 ! loop using an adaptive mixing scheme. The $j$th component of the output
25 ! vector is mixed with a fraction of the same component of the input vector:
26 ! $$ \mu^{i+1}_j=\beta^i_j\nu^i_j+(1-\beta^i_j)\mu^i_j, $$
27 ! where $\beta^{i+1}_j=\beta^i_j+\beta_0$ if $f^i_j\equiv\nu^i_j-\mu^i_j$ does
28 ! not change sign between loops. If $f^i_j$ does change sign, then
29 ! $\beta^{i+1}_j=(\beta^i_j+\beta_0)/2$. Note that the array {\tt nu} serves
30 ! for both input and output, and the arrays {\tt mu}, {\tt beta} and {\tt f}
31 ! are used internally and should not be changed between calls. The routine is
32 ! thread-safe so long as each thread has its own independent work arrays.
33 ! Complex arrays may be passed as real arrays with $n$ doubled.
34 !
35 ! !REVISION HISTORY:
36 ! Created March 2003 (JKD)
37 ! Modified, September 2008 (JKD)
38 ! Modified, August 2011 (JKD)
39 !EOP
40 !BOC
41 implicit none
42 ! arguments
43 integer, intent(in) :: iscl
44 real(8), intent(in) :: beta0,betamax
45 integer, intent(in) :: n
46 real(8), intent(inout) :: nu(n),mu(n),beta(n),f(n)
47 real(8), intent(out) :: d
48 ! local variables
49 integer i
50 real(8) t1
51 if (n < 1) return
52 ! initialise mixer
53 if (iscl < 1) then
54  mu(1:n)=nu(1:n)
55  f(1:n)=0.d0
56  beta(1:n)=beta0
57  d=1.d0
58  return
59 end if
60 d=0.d0
61 do i=1,n
62  t1=nu(i)-mu(i)
63  d=d+t1**2
64  if (t1*f(i) >= 0.d0) then
65  beta(i)=beta(i)+beta0
66  if (beta(i) > betamax) beta(i)=betamax
67  else
68  beta(i)=0.5d0*(beta(i)+beta0)
69  end if
70  f(i)=t1
71  nu(i)=beta(i)*nu(i)+(1.d0-beta(i))*mu(i)
72  mu(i)=nu(i)
73 end do
74 d=sqrt(d/dble(n))
75 end subroutine
76 !EOC
77 
pure subroutine mixadapt(iscl, beta0, betamax, n, nu, mu, beta, f, d)
Definition: mixadapt.f90:10