quack/src/utils/DIIS_extrapolation.f90

56 lines
1.1 KiB
Fortran
Raw Normal View History

2019-03-19 11:21:34 +01:00
subroutine DIIS_extrapolation(rcond,n_err,n_e,n_diis,error,e,error_in,e_inout)
2019-03-19 10:13:33 +01:00
! Perform DIIS extrapolation
implicit none
include 'parameters.h'
! Input variables
integer,intent(in) :: n_err,n_e
double precision,intent(in) :: error_in(n_err),error(n_err,n_diis),e(n_e,n_diis)
! Local variables
double precision,allocatable :: A(:,:),b(:),w(:)
! Output variables
2019-03-19 11:21:34 +01:00
double precision,intent(out) :: rcond
2019-03-19 10:13:33 +01:00
integer,intent(inout) :: n_diis
double precision,intent(inout):: e_inout(n_e)
! Memory allocaiton
allocate(A(n_diis+1,n_diis+1),b(n_diis+1),w(n_diis+1))
! Update DIIS "history"
call prepend(n_err,n_diis,error,error_in)
call prepend(n_e,n_diis,e,e_inout)
! Build A matrix
A(1:n_diis,1:n_diis) = matmul(transpose(error),error)
A(1:n_diis,n_diis+1) = -1d0
A(n_diis+1,1:n_diis) = -1d0
A(n_diis+1,n_diis+1) = +0d0
! Build x matrix
b(1:n_diis) = +0d0
b(n_diis+1) = -1d0
! Solve linear system
2022-11-29 16:01:02 +01:00
! call easy_linear_solve(n_diis+1,A,b,w)
call linear_solve(n_diis+1,A,b,w,rcond)
2019-03-19 10:13:33 +01:00
! Extrapolate
2022-11-29 14:47:51 +01:00
e_inout(:) = matmul(w(1:n_diis),transpose(e(:,1:n_diis)))
2019-03-19 10:13:33 +01:00
end subroutine DIIS_extrapolation