File: mixadapt.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (79 lines) | stat: -rw-r--r-- 2,608 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79

! Copyright (C) 2002-2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU Lesser General Public
! License. See the file COPYING for license details.

!BOP
! !ROUTINE: mixadapt
! !INTERFACE:
subroutine mixadapt(iscl,beta0,betamax,n,nu,mu,beta,f,d)
! !INPUT/OUTPUT PARAMETERS:
!   iscl    : self-consistent loop number (in,integer)
!   beta0   : mixing parameter (in,real)
!   betamax : maximum mixing parameter (in,real)
!   n       : vector length (in,integer)
!   nu      : current output vector as well as the next input vector of the
!             self-consistent loop (inout,real(n))
!   mu      : used internally (inout,real(n))
!   beta    : used internally (inout,real(n))
!   f       : used internally (inout,real(n))
!   d       : RMS difference between old and new output vectors (out,real)
! !DESCRIPTION:
!   Given the input vector $\mu^i$ and output vector $\nu^i$ of the $i$th
!   self-consistent loop, this routine generates the next input vector to the
!   loop using an adaptive mixing scheme. The $j$th component of the output
!   vector is mixed with a fraction of the same component of the input vector:
!   $$ \mu^{i+1}_j=\beta^i_j\nu^i_j+(1-\beta^i_j)\mu^i_j, $$
!   where $\beta^{i+1}_j=\beta^i_j+\beta_0$ if $f^i_j\equiv\nu^i_j-\mu^i_j$ does
!   not change sign between loops. If $f^i_j$ does change sign, then
!   $\beta^{i+1}_j=(\beta^i_j+\beta_0)/2$. Note that the array {\tt nu} serves
!   for both input and output, and the arrays {\tt mu}, {\tt beta} and {\tt f}
!   are used internally and should not be changed between calls. The routine is
!   thread-safe so long as each thread has its own independent work arrays.
!   Complex arrays may be passed as real arrays with $n$ doubled.
!
! !REVISION HISTORY:
!   Created March 2003 (JKD)
!   Modified, September 2008 (JKD)
!   Modified, August 2011 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: iscl
real(8), intent(in) :: beta0,betamax
integer, intent(in) :: n
real(8), intent(inout) :: nu(n),mu(n)
real(8), intent(inout) :: beta(n),f(n)
real(8), intent(out) :: d
! local variables
integer i
real(8) t1
if (n.le.0) return
! initialise mixer
if (iscl.le.0) then
  call dcopy(n,nu,1,mu,1)
  f(:)=0.d0
  beta(:)=beta0
  d=1.d0
  return
end if
d=0.d0
do i=1,n
  t1=nu(i)-mu(i)
  d=d+t1**2
  if (t1*f(i).ge.0.d0) then
    beta(i)=beta(i)+beta0
    if (beta(i).gt.betamax) beta(i)=betamax
  else
    beta(i)=(beta(i)+beta0)*0.5d0
  end if
  f(i)=t1
  nu(i)=beta(i)*nu(i)+(1.d0-beta(i))*mu(i)
  mu(i)=nu(i)
end do
d=sqrt(d/dble(n))
return
end subroutine
!EOC