Class | chrolog_class |
In: |
SolverClass/chrolog_class.F90
|
— Version
$Id: chrolog_class.F90,v 1.5 2011/05/21 10:18:11 ishikawa Exp $
Constant : | |||
CHRON_OP_CONVERGED = 3 : | integer, parameter
|
Constant : | |||
CHRON_OP_DO_MATVEC = 1 : | integer, parameter
|
Constant : | |||
CHRON_OP_MAXITER_REACHED = -9 : | integer, parameter
|
Constant : | |||
CHRON_OP_PRECOND = 4 : | integer, parameter
|
Constant : | |||
CHRON_OP_PRINT_STATUS = 2 : | integer, parameter
|
Constant : | |||
CHRON_OP_PROJECTION = 5 : | integer, parameter
|
Derived Type : | |||
store_vec(:) => NULL() : | complex(DP), pointer
| ||
src_vec(:) => NULL() : | complex(DP), pointer
| ||
dst_vec(:) => NULL() : | complex(DP), pointer
|
Subroutine : | |
this : | type(chrolog_alg), intent(inout) |
subroutine delete_chrolog(this) implicit none type(chrolog_alg), intent(inout) :: this integer :: i if (this%NMAX == 0) return if (associated(this%x)) then do i=1,SIZE(this%x(:)) if (associated(this%x(i)%v)) deallocate(this%x(i)%v) enddo deallocate(this%x) endif if (associated(this%xguess)) deallocate(this%xguess) if (associated(this%b)) deallocate(this%b) return end subroutine
Function : | |
ncron : | integer |
this : | type(chrolog_alg), intent(in) |
function get_ncron_chrolog(this) result(ncron) implicit none type(chrolog_alg), intent(in) :: this integer :: ncron ncron = this%NMAX return end function
Function : | |
res : | real(DP) |
this : | type(chrolog_alg), intent(in) |
function get_residual_norm_chrolog(this) result(res) implicit none type(chrolog_alg), intent(in) :: this real(DP) :: res res = this%res return end function
Function : | |
istat : | integer |
this : | type(chrolog_alg), intent(in) |
function get_status_chrolog(this) result(istat) implicit none type(chrolog_alg), intent(in) :: this integer :: istat istat = this%istat return end function
Function : | |
flag : | logical |
this : | type(chrolog_alg), intent(in) |
function has_history_chrolog(this) result(flag) implicit none type(chrolog_alg), intent(in) :: this logical :: flag if (this%n > 0) then flag = .true. else flag = .false. endif return end function
Subroutine : | |||
this : | type(chrolog_alg), intent(inout) | ||
NSIZE : | integer, intent(in)
| ||
NMAX : | integer, intent(in)
|
subroutine new_chrolog(this,NSIZE,NMAX) implicit none type(chrolog_alg), intent(inout) :: this integer, intent(in) :: NSIZE ! vector size integer, intent(in) :: NMAX ! dimension size integer :: i,nn,is this%NSIZE = NSIZE this%NMAX = NMAX this%n = 0 this%inow = 1 this%loop_count = 0 if (this%NSIZE == 0) then write(*,'("NSIZE is zero in new_chrolog.")') stop endif if (this%NMAX == 0) return if (associated(this%x)) then do i=1,SIZE(this%x(:)) if (associated(this%x(i)%v)) deallocate(this%x(i)%v) enddo deallocate(this%x) endif if (associated(this%xguess)) deallocate(this%xguess) if (associated(this%b)) deallocate(this%b) allocate(this%x(this%NMAX)) do i=1,this%NMAX allocate(this%x(i)%v(this%NSIZE)) enddo allocate(this%xguess(this%NSIZE)) allocate(this%b(this%NSIZE)) do i=1,this%NMAX !$OMP PARALLEL DO PRIVATE(is) do is=1,this%NSIZE this%x(i)%v(is) = Z0 enddo enddo this%src_vec => this%b this%dst_vec => this%xguess this%store_vec => this%x(this%inow)%v return end subroutine
Subroutine : | |
this : | type(chrolog_alg), intent(inout) |
Chronological initial solution guess (for full lattice solver)
Least squre fitting (Minimal residual guess) for Dy = D^(-1) y.
Minimize
Min_{c(i)}[ | D Dy - y |^2 ], with Dy = sum_{i=1,N} c(i) x(i),
where x(i)’s are previous solutions, c(i)’s are the linear fitting coefficients.
The fitting is done using QR factoryzation for (D x(i)).
this : chrolog_alg
subroutine solve_chrolog(this) ! ! Chronological initial solution guess (for full lattice solver) ! ! Least squre fitting (Minimal residual guess) for Dy = D^(-1) y. ! ! Minimize ! ! Min_{c(i)}[ | D Dy - y |^2 ], with Dy = sum_{i=1,N} c(i) x(i), ! ! where x(i)'s are previous solutions, c(i)'s are the linear fitting coefficients. ! ! The fitting is done using QR factoryzation for (D x(i)). ! ! this : chrolog_alg ! implicit none type(chrolog_alg), intent(inout) :: this complex(DP), allocatable :: R(:,:) complex(DP) :: ctmp real(DP) :: rtmp integer :: is,i,j,n,itop,nn select case(this%inext) case (0) if (this%nmax == 0) then this%istat = CHRON_OP_CONVERGED this%inext = 0 return endif n = this%n if (n == 0) then this%istat = CHRON_OP_CONVERGED this%inext = 0 return endif if (associated(this%vwork)) then do i=1,SIZE(this%vwork) deallocate(this%vwork(i)%v) enddo deallocate(this%vwork) endif nn = this%n+1 allocate(this%vwork(nn)) do i=1,nn allocate(this%vwork(i)%v(this%NSIZE)) enddo this%itop = mod(this%inow-1-1+this%nmax,this%nmax)+1 this%loop_count = 1 this%istat = CHRON_OP_NOP this%inext = this%inext + 1 return case (1) !------------------ ! dx(i) = D x(i) !------------------ this%src_vec => this%x(this%itop)%v this%dst_vec => this%vwork(this%loop_count)%v this%istat = CHRON_OP_DO_MATVEC this%inext = this%inext + 1 return case (2) this%loop_count = this%loop_count + 1 this%itop = mod(this%itop-1-1+this%nmax,this%nmax)+1 if (this%loop_count == this%n+1) then !-------------------------- ! set source vector !-------------------------- !$OMP PARALLEL DO PRIVATE(is) do is=1,this%NSIZE this%vwork(this%n+1)%v(is) = this%b(is) enddo this%istat = CHRON_OP_NOP this%inext = this%inext + 1 else this%istat = CHRON_OP_NOP this%inext = 1 endif return case (3) n = this%n nn = this%n+1 allocate(R(nn,nn)) R(:,:) = Z0 !------------------------------------------- ! QR decomposition for Q*R = dx, dx <= Q !------------------------------------------- do i=1,nn do j=1,i-1 ctmp = prod(this%vwork(j)%v,this%vwork(i)%v) !$OMP PARALLEL DO PRIVATE(is) do is=1,this%NSIZE this%vwork(i)%v(is) = this%vwork(i)%v(is) - this%vwork(j)%v(is)*ctmp enddo R(j,i) = ctmp enddo rtmp = abs(this%vwork(i)%v) R(i,i) = rtmp rtmp = 1.0_DP/rtmp !$OMP PARALLEL DO PRIVATE(is) do is=1,this%NSIZE this%vwork(i)%v(is) = this%vwork(i)%v(is)*rtmp enddo enddo !do j=1,nn !do i=1,nn ! ctmp = prod(this%vwork(i)%v,this%vwork(j)%v) ! write(*,'(2I3,2E24.15)')i,j,ctmp !enddo !enddo !------------------------------------- ! normalized residual |b-Ax|/|b| !------------------------------------- rtmp = real(dot_product(R(1:n,nn),R(1:n,nn)),kind=DP) this%res = sqrt(1.0_DP -rtmp/abs(this%b)**2) !------------------------------------- ! R(1:n,nn) <= R(1:n,1:n) \ R(1:n,nn) !------------------------------------- do i=n,1,-1 do j=i+1,n R(i,nn) = R(i,nn) - R(i,j)*R(j,nn) enddo rtmp = real(R(i,i),kind=KIND(rtmp)) R(i,nn) = R(i,nn) / rtmp enddo !------------------------------------- ! xguess = sum_{i=1,n}[x(i)*R(i,nn)] !------------------------------------- itop = mod(this%inow-1-1+this%nmax,this%nmax)+1 i=1 ctmp = R(i,nn) !$OMP PARALLEL DO PRIVATE(is) do is=1,this%NSIZE this%xguess(is) = this%x(itop)%v(is)*ctmp enddo itop = mod(itop-1-1+this%nmax,this%nmax)+1 do i=2,n ctmp = R(i,nn) !$OMP PARALLEL DO PRIVATE(is) do is=1,this%NSIZE this%xguess(is) = this%xguess(is) + this%x(itop)%v(is)*ctmp enddo itop = mod(itop-1-1+this%nmax,this%nmax)+1 enddo deallocate(R) if (associated(this%vwork)) then do i=1,SIZE(this%vwork) deallocate(this%vwork(i)%v) enddo deallocate(this%vwork) endif this%dst_vec => this%xguess this%src_vec => this%b this%istat = CHRON_OP_CONVERGED this%inext = 0 return end select return end subroutine
Subroutine : | |
this : | type(chrolog_alg), intent(inout) |
Store solution x for chronological guess.
subroutine store_chrolog(this) ! ! \Store solution x for chronological guess. ! implicit none type(chrolog_alg), intent(inout) :: this real(DP) :: rtmp integer :: is if (this%nmax==0) return #ifdef _DEBUG if (0 == nodeid) write(*,'("STORE START")') call print_stat_chrolog(this) #endif !------------------------------------- ! Normalize !------------------------------------- rtmp = 1.0_DP/abs(this%store_vec) !$OMP PARALLEL Do PRIVATE(is) do is=1,this%NSIZE this%store_vec(is) = this%store_vec(is)*rtmp enddo this%inow = mod(this%inow-1+1,this%nmax)+1 this%n = min(this%n + 1,this%nmax) this%store_vec => this%x(this%inow)%v this%inext = 0 #ifdef _DEBUG if (0 == nodeid ) write(*,'("STORE END")') call print_stat_chrolog(this) #endif return end subroutine