Class | solver_class |
In: |
SolverClass/solver_class.F90
OBSOLETES/PREPROC/sample.F90 OBSOLETES/PREPROC/solver_class.F90 |
This module contains solver algorithm with reverse communication interface (RCI). BiCGStab and BiCGstab(L=2) are included. This module uses some Fortran 2003 features: type extension and class polymorphism.
To solve, A x = b , program flow is as follows.
integer, parameter :: N = 1000 integer :: istat type(bicgstab_alg) :: solver complex(8) :: A(N,N) complex(8) :: x(N),b(N) ..... set A(:,:) something set b(:) something ..... ! set matrix size to 1000, maximul iteration counts to 100 ! set solver tolerance to 10^-12 for |Ax-b|/|b| < tol call new(solver,NSIZE=N,max_iter=100,guess=0,tol=1.0d-12) ! set source vector b(n) to solver solver%src_vec(1:N) = b(:) ! set solution vector x(n) = 0 (guess = 0 case) solver%dst_vec(1:N) = 0.0d0 ! iteration loop with RCI do call solve(solver) istat = get_status(solver) select case(istat) case(OP_NOP) cycle ! DO NOTHING case(OP_DO_MATVEC) ! Multiply matrix A on source vector and store the result in destination vector solver%dst_vec(:) = MATMUL(A(:,:),solver%src_vec(:)) case(OP_PRINT_STATUS) ! Print out iteration status write(*,'("ITER=",I4," RES=",E24.15)') get_current_iteration(solver), get_residual_norm(solver) case(OP_CONVERGED) ! Solver converged write(*,'("ITER=",I4," RES=",E24.15)') get_current_iteration(solver), get_residual_norm(solver) write(*,'("SOLVER CONVERGED!")') exit ! exit do loop case(OP_MAXITER_REACHED) ! Solver reaches maximum iteration counts. write(*,'("ITER=",I4," RES=",E24.15)') get_current_iteration(solver), get_residual_norm(solver) write(*,'("SOLVER DID NOT CONVERGE!")') stop ! stop program end select enddo ! get solution to x(:) x(:) = solver%dst_vec(:) stop
$Id: solver_class.F90,v 1.22 2010/08/27 09:08:03 ishikawa Exp ishikawa $
Constant : | |||
OP_DO_MATVEC = 1 : | integer, parameter
|
Constant : | |||
OP_DO_MATVEC = 1 : | integer, parameter
|
Constant : | |||
OP_DO_MATVEC = 1 : | integer, parameter
|
Constant : | |||
OP_MAXITER_REACHED = -9 : | integer, parameter
|
Constant : | |||
OP_MAXITER_REACHED = -9 : | integer, parameter
|
Constant : | |||
OP_MAXITER_REACHED = -9 : | integer, parameter
|
Constant : | |||
OP_PRECOND = 4 : | integer, parameter
|
Constant : | |||
OP_PRECOND = 4 : | integer, parameter
|
Constant : | |||
OP_PRECOND = 4 : | integer, parameter
|
Constant : | |||
OP_PRINT_STATUS = 2 : | integer, parameter
|
Constant : | |||
OP_PRINT_STATUS = 2 : | integer, parameter
|
Constant : | |||
OP_PRINT_STATUS = 2 : | integer, parameter
|
Constant : | |||
OP_PROJECTION = 5 : | integer, parameter
|
Constant : | |||
OP_PROJECTION = 5 : | integer, parameter
|
Constant : | |||
OP_PROJECTION = 5 : | integer, parameter
|
Subroutine : | |
this : | class(solver_algorithm), intent(inout) |
delete solver algorithm
subroutine delete_alg(this) ! ! delete solver algorithm ! implicit none class(solver_algorithm), intent(inout) :: this integer :: i if (this%debug) then write(*,fmt='("I am ",A,": delete")')TRIM(this%name) endif do i=1,this%NVWORK deallocate(this%vwork(i)%v) this%vwork(i)%v => NULL() enddo deallocate(this%vwork) this%vwork => NULL() deallocate(this%rwork) this%rwork => NULL() deallocate(this%zwork) this%zwork => NULL() this%src_vec => NULL() this%dst_vec => NULL() return end subroutine
Subroutine : | |
this : | class(solver_algorithm), intent(inout) |
delete solver algorithm
subroutine delete_alg(this) ! ! delete solver algorithm ! implicit none class(solver_algorithm), intent(inout) :: this integer :: i if (this%debug) then write(*,fmt='("I am ",A,": delete")')TRIM(this%name) endif do i=1,this%NVWORK deallocate(this%vwork(i)%v) this%vwork(i)%v => NULL() enddo deallocate(this%vwork) this%vwork => NULL() deallocate(this%rwork) this%rwork => NULL() deallocate(this%zwork) this%zwork => NULL() this%src_vec => NULL() this%dst_vec => NULL() return end subroutine
Subroutine : | |
this : | class(solver_algorithm), intent(inout) |
delete solver algorithm
subroutine delete_alg(this) ! ! delete solver algorithm ! implicit none class(solver_algorithm), intent(inout) :: this integer :: i if (this%debug) then write(*,fmt='("I am ",A,": delete")')TRIM(this%name) endif do i=1,this%NVWORK deallocate(this%vwork(i)%v) this%vwork(i)%v => NULL() enddo deallocate(this%vwork) this%vwork => NULL() deallocate(this%rwork) this%rwork => NULL() deallocate(this%zwork) this%zwork => NULL() this%src_vec => NULL() this%dst_vec => NULL() return end subroutine
Function : | |
iter : | integer |
this : | class(solver_algorithm), intent(in) |
return current iteration count of the solver
function get_current_iteration_alg(this) result(iter) ! ! return current iteration count of the solver ! implicit none class(solver_algorithm), intent(in) :: this integer :: iter iter = this%iter if (this%debug) then write(*,fmt='("I am ",A,": get_current_iteration")')TRIM(this%name) endif return end function
Function : | |
iter : | integer |
this : | class(solver_algorithm), intent(in) |
return current iteration count of the solver
function get_current_iteration_alg(this) result(iter) ! ! return current iteration count of the solver ! implicit none class(solver_algorithm), intent(in) :: this integer :: iter iter = this%iter if (this%debug) then write(*,fmt='("I am ",A,": get_current_iteration")')TRIM(this%name) endif return end function
Function : | |
iter : | integer |
this : | class(solver_algorithm), intent(in) |
return current iteration count of the solver
function get_current_iteration_alg(this) result(iter) ! ! return current iteration count of the solver ! implicit none class(solver_algorithm), intent(in) :: this integer :: iter iter = this%iter if (this%debug) then write(*,fmt='("I am ",A,": get_current_iteration")')TRIM(this%name) endif return end function
Function : | |
name : | character(ALG_NAMELEN) |
this : | class(solver_algorithm), intent(in) |
return solver name in string
function get_name_alg(this) result(name) ! ! return solver name in string ! implicit none class(solver_algorithm), intent(in) :: this character(ALG_NAMELEN) :: name name = this%name if (this%debug) then write(*,fmt='("I am ",A,": get_name")')TRIM(this%name) endif return end function
Function : | |
name : | character(ALG_NAMELEN) |
this : | class(solver_algorithm), intent(in) |
return solver name in string
function get_name_alg(this) result(name) ! ! return solver name in string ! implicit none class(solver_algorithm), intent(in) :: this character(ALG_NAMELEN) :: name name = this%name if (this%debug) then write(*,fmt='("I am ",A,": get_name")')TRIM(this%name) endif return end function
Function : | |
name : | character(ALG_NAMELEN) |
this : | class(solver_algorithm), intent(in) |
return solver name in string
function get_name_alg(this) result(name) ! ! return solver name in string ! implicit none class(solver_algorithm), intent(in) :: this character(ALG_NAMELEN) :: name name = this%name if (this%debug) then write(*,fmt='("I am ",A,": get_name")')TRIM(this%name) endif return end function
Function : | |
res : | real(DP) |
this : | class(solver_algorithm), intent(in) |
return current residual norm
function get_residual_norm_alg(this) result(res) ! ! return current residual norm ! implicit none class(solver_algorithm), intent(in) :: this real(DP) :: res res = this%res if (this%debug) then write(*,fmt='("I am ",A,": get_residual_norm")')TRIM(this%name) endif return end function
Function : | |
res : | real(DP) |
this : | class(solver_algorithm), intent(in) |
return current residual norm
function get_residual_norm_alg(this) result(res) ! ! return current residual norm ! implicit none class(solver_algorithm), intent(in) :: this real(DP) :: res res = this%res if (this%debug) then write(*,fmt='("I am ",A,": get_residual_norm")')TRIM(this%name) endif return end function
Function : | |
res : | real(DP) |
this : | class(solver_algorithm), intent(in) |
return current residual norm
function get_residual_norm_alg(this) result(res) ! ! return current residual norm ! implicit none class(solver_algorithm), intent(in) :: this real(DP) :: res res = this%res if (this%debug) then write(*,fmt='("I am ",A,": get_residual_norm")')TRIM(this%name) endif return end function
Function : | |
istat : | integer |
this : | class(solver_algorithm), intent(in) |
return RCI status
function get_status_alg(this) result(istat) ! ! return RCI status ! implicit none class(solver_algorithm), intent(in) :: this integer :: istat istat = this%istat if (this%debug) then write(*,fmt='("I am ",A,": get_status")')TRIM(this%name) endif return end function
Function : | |
istat : | integer |
this : | class(solver_algorithm), intent(in) |
return RCI status
function get_status_alg(this) result(istat) ! ! return RCI status ! implicit none class(solver_algorithm), intent(in) :: this integer :: istat istat = this%istat if (this%debug) then write(*,fmt='("I am ",A,": get_status")')TRIM(this%name) endif return end function
Function : | |
istat : | integer |
this : | class(solver_algorithm), intent(in) |
return RCI status
function get_status_alg(this) result(istat) ! ! return RCI status ! implicit none class(solver_algorithm), intent(in) :: this integer :: istat istat = this%istat if (this%debug) then write(*,fmt='("I am ",A,": get_status")')TRIM(this%name) endif return end function
Subroutine : | |||
this : | type(bicgstab_alg), intent(inout) | ||
NSIZE : | integer, intent(in)
| ||
max_iter : | integer, optional, intent(in)
| ||
guess : | integer, optional, intent(in)
| ||
mode : | integer, optional, intent(in)
| ||
tol : | real(DP), optional, intent(in)
| ||
debug : | logical, optional, intent(in)
|
Initialize BiCGstab algorithm instance.
User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x as descrived above after calling this method,
subroutine new_bicgstab(this,NSIZE,max_iter,guess,mode,tol,debug) ! ! Initialize BiCGstab algorithm instance. ! ! ! User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x ! as descrived above after calling this method, ! implicit none type(bicgstab_alg), intent(inout) :: this integer, intent(in) :: NSIZE ! set vector length / matrix size integer, optional, intent(in) :: max_iter ! set maximum iteration counts integer, optional, intent(in) :: guess ! set switch for initial guess. ! guess = 0, zero initial vector is assumed. ! guess = 1, user should supply inital guess in this%dst_vec(:). integer, optional, intent(in) :: mode ! set mode for preconditioner ! mode = 0, normal ! mode = 1, flexible preconditioner ! mode = 2, deflation is used real(DP), optional, intent(in) :: tol ! set stopping condition logical, optional, intent(in) :: debug ! debug=.true. debug mode character(CHARLEN) :: str this%NSIZE = NSIZE this%NVWORK = 8 this%NRWORK = 1 this%NZWORK = 5 this%name = "BiCGstab" if (present(tol)) then this%tol = tol else this%tol = EPSILON(1.0_DP)*2.0_DP endif if (present(max_iter)) then this%max_iter = max_iter endif if (present(guess)) then this%guess = guess endif if (present(mode)) then if (mode == MODE_PROJECTION) then write(*,'("PROJECTION(DEFLATION) IS NOT IMPREMENTED IN BICGSTAB ALG.")') stop endif this%mode = mode endif if (present(debug)) then this%debug = debug endif call new_alg(this) ! this initialize algorithm independent part of BiCGStab. return end subroutine
Subroutine : | |||
this : | type(bicgstab_alg), intent(inout) | ||
NSIZE : | integer, intent(in)
| ||
max_iter : | integer, optional, intent(in)
| ||
guess : | integer, optional, intent(in)
| ||
mode : | integer, optional, intent(in)
| ||
tol : | real(DP), optional, intent(in)
| ||
debug : | logical, optional, intent(in)
|
Initialize BiCGstab algorithm instance.
User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x as descrived above after calling this method,
subroutine new_bicgstab(this,NSIZE,max_iter,guess,mode,tol,debug) ! ! Initialize BiCGstab algorithm instance. ! ! ! User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x ! as descrived above after calling this method, ! implicit none type(bicgstab_alg), intent(inout) :: this integer, intent(in) :: NSIZE ! set vector length / matrix size integer, optional, intent(in) :: max_iter ! set maximum iteration counts integer, optional, intent(in) :: guess ! set switch for initial guess. ! guess = 0, zero initial vector is assumed. ! guess = 1, user should supply inital guess in this%dst_vec(:). integer, optional, intent(in) :: mode ! set mode for preconditioner ! mode = 0, normal ! mode = 1, flexible preconditioner ! mode = 2, deflation is used real(DP), optional, intent(in) :: tol ! set stopping condition logical, optional, intent(in) :: debug ! debug=.true. debug mode character(CHARLEN) :: str this%NSIZE = NSIZE this%NVWORK = 8 this%NRWORK = 1 this%NZWORK = 5 this%name = "BiCGstab" if (present(tol)) then this%tol = tol else this%tol = EPSILON(1.0_DP)*2.0_DP endif if (present(max_iter)) then this%max_iter = max_iter endif if (present(guess)) then this%guess = guess endif if (present(mode)) then if (mode == MODE_PROJECTION) then write(*,'("PROJECTION(DEFLATION) IS NOT IMPREMENTED IN BICGSTAB ALG.")') stop endif this%mode = mode endif if (present(debug)) then this%debug = debug endif call new_alg(this) ! this initialize algorithm independent part of BiCGStab. return end subroutine
Subroutine : | |||
this : | type(bicgstab_alg), intent(inout) | ||
NSIZE : | integer, intent(in)
| ||
max_iter : | integer, optional, intent(in)
| ||
guess : | integer, optional, intent(in)
| ||
mode : | integer, optional, intent(in)
| ||
tol : | real(DP), optional, intent(in)
| ||
debug : | logical, optional, intent(in)
|
Initialize BiCGstab algorithm instance.
User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x as descrived above after calling this method,
subroutine new_bicgstab(this,NSIZE,max_iter,guess,mode,tol,debug) ! ! Initialize BiCGstab algorithm instance. ! ! ! User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x ! as descrived above after calling this method, ! implicit none type(bicgstab_alg), intent(inout) :: this integer, intent(in) :: NSIZE ! set vector length / matrix size integer, optional, intent(in) :: max_iter ! set maximum iteration counts integer, optional, intent(in) :: guess ! set switch for initial guess. ! guess = 0, zero initial vector is assumed. ! guess = 1, user should supply inital guess in this%dst_vec(:). integer, optional, intent(in) :: mode ! set mode for preconditioner ! mode = 0, normal ! mode = 1, flexible preconditioner ! mode = 2, deflation is used real(DP), optional, intent(in) :: tol ! set stopping condition logical, optional, intent(in) :: debug ! debug=.true. debug mode character(CHARLEN) :: str this%NSIZE = NSIZE this%NVWORK = 8 this%NRWORK = 1 this%NZWORK = 5 this%name = "BiCGstab" if (present(tol)) then this%tol = tol else this%tol = EPSILON(1.0_DP)*2.0_DP endif if (present(max_iter)) then this%max_iter = max_iter endif if (present(guess)) then this%guess = guess endif if (present(mode)) then if (mode == MODE_PROJECTION) then write(*,'("PROJECTION(DEFLATION) IS NOT IMPREMENTED IN BICGSTAB ALG.")') stop endif this%mode = mode endif if (present(debug)) then this%debug = debug endif call new_alg(this) ! this initialize algorithm independent part of BiCGStab. return end subroutine
Subroutine : | |||
this : | type(cg_alg), intent(inout) | ||
NSIZE : | integer, intent(in)
| ||
max_iter : | integer, optional, intent(in)
| ||
guess : | integer, optional, intent(in)
| ||
mode : | integer, optional, intent(in)
| ||
tol : | real(DP), optional, intent(in)
| ||
debug : | logical, optional, intent(in)
|
Initialize CG algorithm instance.
User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x as descrived above after calling this method,
subroutine new_cg(this,NSIZE,max_iter,guess,mode,tol,debug) ! ! Initialize CG algorithm instance. ! ! ! User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x ! as descrived above after calling this method, ! implicit none type(cg_alg), intent(inout) :: this integer, intent(in) :: NSIZE ! set vector lenght / matrix size integer, optional, intent(in) :: max_iter ! set maximum iteration counts integer, optional, intent(in) :: guess ! set switch for initial guess. ! guess = 1, user should supply inital guess in this%dst_vec(:). ! guess = 0, user should set zero in this%dst_vec(:). integer, optional, intent(in) :: mode ! set mode for solver type ! mode = 0, normal solver ! mode = 1, flexible preconditioner ! mode = 2, deflation real(DP), optional, intent(in) :: tol ! set stopping condition logical, optional, intent(in) :: debug ! if true debug mode character(CHARLEN) :: str this%NSIZE = NSIZE this%NVWORK = 6 this%NRWORK = 5 this%NZWORK = 0 this%name = "CG" if (present(tol)) then this%tol = tol else this%tol = EPSILON(1.0_DP)*2.0_DP endif if (present(max_iter)) then this%max_iter = max_iter endif if (present(guess)) then this%guess = guess endif if (present(mode)) then if (mode == MODE_PRECOND) then write(*,'("PRECONDITIONER IS NOT IMPREMENTED IN CG ALG.")') stop endif this%mode = mode endif if (present(debug)) then this%debug = debug endif call new_alg(this) ! this initialize algorithm independent part of BiCGStab. return end subroutine
Subroutine : | |||
this : | type(cg_alg), intent(inout) | ||
NSIZE : | integer, intent(in)
| ||
max_iter : | integer, optional, intent(in)
| ||
guess : | integer, optional, intent(in)
| ||
mode : | integer, optional, intent(in)
| ||
tol : | real(DP), optional, intent(in)
| ||
debug : | logical, optional, intent(in)
|
Initialize CG algorithm instance.
User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x as descrived above after calling this method,
subroutine new_cg(this,NSIZE,max_iter,guess,mode,tol,debug) ! ! Initialize CG algorithm instance. ! ! ! User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x ! as descrived above after calling this method, ! implicit none type(cg_alg), intent(inout) :: this integer, intent(in) :: NSIZE ! set vector lenght / matrix size integer, optional, intent(in) :: max_iter ! set maximum iteration counts integer, optional, intent(in) :: guess ! set switch for initial guess. ! guess = 1, user should supply inital guess in this%dst_vec(:). ! guess = 0, user should set zero in this%dst_vec(:). integer, optional, intent(in) :: mode ! set mode for solver type ! mode = 0, normal solver ! mode = 1, flexible preconditioner ! mode = 2, deflation real(DP), optional, intent(in) :: tol ! set stopping condition logical, optional, intent(in) :: debug ! if true debug mode character(CHARLEN) :: str this%NSIZE = NSIZE this%NVWORK = 6 this%NRWORK = 5 this%NZWORK = 0 this%name = "CG" if (present(tol)) then this%tol = tol else this%tol = EPSILON(1.0_DP)*2.0_DP endif if (present(max_iter)) then this%max_iter = max_iter endif if (present(guess)) then this%guess = guess endif if (present(mode)) then if (mode == MODE_PRECOND) then write(*,'("PRECONDITIONER IS NOT IMPREMENTED IN CG ALG.")') stop endif this%mode = mode endif if (present(debug)) then this%debug = debug endif call new_alg(this) ! this initialize algorithm independent part of BiCGStab. return end subroutine
Subroutine : | |||
this : | type(cg_alg), intent(inout) | ||
NSIZE : | integer, intent(in)
| ||
max_iter : | integer, optional, intent(in)
| ||
guess : | integer, optional, intent(in)
| ||
mode : | integer, optional, intent(in)
| ||
tol : | real(DP), optional, intent(in)
| ||
debug : | logical, optional, intent(in)
|
Initialize CG algorithm instance.
User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x as descrived above after calling this method,
subroutine new_cg(this,NSIZE,max_iter,guess,mode,tol,debug) ! ! Initialize CG algorithm instance. ! ! ! User should set source vector b (rhs of Ax = b) in this%src_srcv(:), and set x ! as descrived above after calling this method, ! implicit none type(cg_alg), intent(inout) :: this integer, intent(in) :: NSIZE ! set vector lenght / matrix size integer, optional, intent(in) :: max_iter ! set maximum iteration counts integer, optional, intent(in) :: guess ! set switch for initial guess. ! guess = 1, user should supply inital guess in this%dst_vec(:). ! guess = 0, user should set zero in this%dst_vec(:). integer, optional, intent(in) :: mode ! set mode for solver type ! mode = 0, normal solver ! mode = 1, flexible preconditioner ! mode = 2, deflation real(DP), optional, intent(in) :: tol ! set stopping condition logical, optional, intent(in) :: debug ! if true debug mode character(CHARLEN) :: str this%NSIZE = NSIZE this%NVWORK = 6 this%NRWORK = 5 this%NZWORK = 0 this%name = "CG" if (present(tol)) then this%tol = tol else this%tol = EPSILON(1.0_DP)*2.0_DP endif if (present(max_iter)) then this%max_iter = max_iter endif if (present(guess)) then this%guess = guess endif if (present(mode)) then if (mode == MODE_PRECOND) then write(*,'("PRECONDITIONER IS NOT IMPREMENTED IN CG ALG.")') stop endif this%mode = mode endif if (present(debug)) then this%debug = debug endif call new_alg(this) ! this initialize algorithm independent part of BiCGStab. return end subroutine
Subroutine : | |
this : | type(bicgstab_alg), intent(inout) |
BiCGStab Solver with reverse communication interface (RCI)
subroutine solve_bicgstab(this) ! ! BiCGStab Solver with reverse communication interface (RCI) ! implicit none type(bicgstab_alg), intent(inout) :: this ! ! index to real scalar variables ! integer, parameter :: idxrwork_rsrc = 1 ! ! index to complex scalar variables ! integer, parameter :: idxzwork_rho = 1 integer, parameter :: idxzwork_rho0 = 2 integer, parameter :: idxzwork_bet = 3 integer, parameter :: idxzwork_omg = 4 integer, parameter :: idxzwork_alp = 5 ! ! index to complex vector variables ! integer, parameter :: idx_b = 1 integer, parameter :: idx_x = 2 integer, parameter :: idx_r = 3 integer, parameter :: idx_rt= 4 integer, parameter :: idx_p = 5 integer, parameter :: idx_q = 6 integer, parameter :: idx_v = 7 integer, parameter :: idx_t = 8 integer :: isite complex(DP) :: rho,rho0,bet,omg,alp complex(DP) :: ctmp0 real(DP) :: rtmp0,rtmp1,rsrc if (this%debug) write(*,*)"INEXT=",this%inext select case (this%inext) case (0) this%iter=0 if (this%guess == GUESS_USE) then !----------------- ! p = A x !----------------- this%istat = OP_DO_MATVEC this%src_vec => this%vwork(idx_x)%v this%dst_vec => this%vwork(idx_p)%v else this%istat = OP_NOP this%src_vec => NULL() this%dst_vec => NULL() endif this%inext=this%inext+1 return case (1) if (this%guess == GUESS_USE) then !----------------------- ! r = b - p = b - A x !----------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) -this%vwork(idx_p)%v(isite) enddo else !----------------------- ! r = b ! x = 0 !----------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) this%vwork(idx_x)%v(isite) = Z0 enddo endif !--------------- ! rt = r ! p = r ! rho0 = <rt|r> ! rsrc = |b| !--------------- rtmp0 = 0.0_DP rtmp1 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0,rtmp1) do isite=1,this%NSIZE this%vwork(idx_rt)%v(isite) = this%vwork(idx_r)%v(isite) this%vwork(idx_p )%v(isite) = this%vwork(idx_r)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_b)%v(isite))**2 +aimag(this%vwork(idx_b)%v(isite))**2 rtmp1 = rtmp1 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) call comlib_sumcast(rtmp1) #endif rsrc = sqrt(rtmp0) rho0 = rtmp1 this%rwork(idxrwork_rsrc) = rsrc this%zwork(idxzwork_rho0) = rho0 this%res = sqrt(rtmp1)/rsrc this%istat = OP_PRINT_STATUS this%inext = this%inext+1 return case (2) !--------------------------------------- ! start main loop !--------------------------------------- if (this%mode == MODE_PRECOND) then !------------------- ! v = M p !------------------- this%src_vec => this%vwork(idx_p)%v this%dst_vec => this%vwork(idx_v)%v this%istat = OP_PRECOND else this%src_vec => NULL() this%dst_vec => NULL() this%istat = OP_NOP endif this%inext = this%inext+1 return case (3) if (this%mode == MODE_PRECOND) then !------------------- ! q = A v !------------------- this%src_vec => this%vwork(idx_v)%v this%dst_vec => this%vwork(idx_q)%v this%istat = OP_DO_MATVEC else !------------------- ! q = A p !------------------- this%src_vec => this%vwork(idx_p)%v this%dst_vec => this%vwork(idx_q)%v this%istat = OP_DO_MATVEC endif this%inext = this%inext+1 this%iter = this%iter+1 return case (4) !-------------------------------- ! alp = rho0 / <rt|q> !-------------------------------- ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_rt)%v(isite)) *this%vwork(idx_q )%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) #endif rho0 = this%zwork(idxzwork_rho0) if (this%debug) then write(*,'("rho0=",2E24.15," <rt|q>=",2E24.15)')rho0,ctmp0 endif alp = rho0/ctmp0 if (this%debug) then write(*,'("alp=",2E24.15)')alp endif this%zwork(idxzwork_alp) = alp if (this%mode == MODE_PRECOND) then !---------------- ! r = r - alp q ! x = x + alp v !---------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -alp*this%vwork(idx_q)%v(isite) this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +alp*this%vwork(idx_v)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 + aimag(this%vwork(idx_r)%v(isite))**2 enddo else !---------------- ! x = x + alp p ! r = r - alp q !---------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +alp*this%vwork(idx_p)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -alp*this%vwork(idx_q)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 + aimag(this%vwork(idx_r)%v(isite))**2 enddo endif #ifndef _singlePU call comlib_sumcast(rtmp0) #endif rsrc = this%rwork(idxrwork_rsrc) this%res = sqrt(rtmp0)/rsrc if (this%debug) write(*,'("res=",E24.15," tol=",E24.15)')this%res,this%tol if ( this%res <= this%tol ) then this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_CONVERGED this%inext = 0 else if ( this%iter <= this%max_iter ) then this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_PRINT_STATUS this%inext = this%inext+1 else this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_MAXITER_REACHED this%inext = 0 endif return case(5) if (this%mode == MODE_PRECOND) then !------------ ! v = M r !------------ this%src_vec => this%vwork(idx_r)%v this%dst_vec => this%vwork(idx_v)%v this%istat = OP_PRECOND else this%src_vec => NULL() this%dst_vec => NULL() this%istat = OP_NOP endif this%inext = this%inext+1 return case(6) if (this%mode == MODE_PRECOND) then !------------ ! t = A v !------------ this%src_vec => this%vwork(idx_v)%v this%dst_vec => this%vwork(idx_t)%v this%istat = OP_DO_MATVEC else !------------ ! t = A r !------------ this%src_vec => this%vwork(idx_r)%v this%dst_vec => this%vwork(idx_t)%v this%istat = OP_DO_MATVEC endif this%inext = this%inext+1 this%iter = this%iter+1 return case(7) !------------------------ ! omg = <t|r>/<t|t> !------------------------ rtmp0 = 0.0_DP ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0,ctmp0) do isite=1,this%NSIZE rtmp0 = rtmp0 + real(this%vwork(idx_t)%v(isite))**2 +aimag(this%vwork(idx_t)%v(isite))**2 ctmp0 = ctmp0 +conjg(this%vwork(idx_t)%v(isite)) *this%vwork(idx_r)%v(isite) enddo #ifndef _singlePU call comlib_sumcast(rtmp0) call comlib_sumcast(ctmp0) #endif if (this%debug) then write(*,'("<t|r>=",2E24.15," <t|t>=",E24.15)')ctmp0,rtmp0 endif omg = ctmp0/rtmp0 this%zwork(idxzwork_omg) = omg if (this%mode == MODE_PRECOND) then !---------------------------- ! x = x + omg v ! r = r - omg t ! res = |r|/rsrc !---------------------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +omg*this%vwork(idx_v)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -omg*this%vwork(idx_t)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo else !---------------------------- ! x = x + omg r ! r = r - omg t ! res = |r|/rsrc !---------------------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +omg*this%vwork(idx_r)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -omg*this%vwork(idx_t)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo endif #ifndef _singlePU call comlib_sumcast(rtmp0) #endif rsrc = this%rwork(idxrwork_rsrc) this%res=sqrt(rtmp0)/rsrc if (this%debug) write(*,'("res=",E24.15," tol=",E24.15)')this%res,this%tol if ( this%res <= this%tol ) then this%istat = OP_CONVERGED this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%inext = 0 else if ( this%iter <= this%max_iter ) then this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_PRINT_STATUS this%inext = this%inext+1 else this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_MAXITER_REACHED this%inext = 0 endif return case(8) !------------------------------- ! rho = <rt|r> ! bet = (alp/omg)*(rho/rho0) ! rho0 = rho !------------------------------- ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_rt)%v(isite)) *this%vwork(idx_r )%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) #endif rho = ctmp0 alp = this%zwork(idxzwork_alp) omg = this%zwork(idxzwork_omg) rho0 = this%zwork(idxzwork_rho0) if (this%debug) then write(*,'("alp=",2E24.15," omg=",2E24.15," rho=",2E24.15," rho0=",2E24.15)')alp,omg,rho,rho0 endif bet = (alp/omg)*(rho/rho0) rho0 = rho this%zwork(idxzwork_rho0) = rho0 !------------------------------- ! p = r + bet (p - omg q) !------------------------------- !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) +bet*( this%vwork(idx_p)%v(isite) -omg*this%vwork(idx_q)%v(isite)) enddo this%istat = OP_NOP this%inext = 2 return end select return end subroutine
Subroutine : | |
this : | type(bicgstab_alg), intent(inout) |
BiCGStab Solver with reverse communication interface (RCI)
subroutine solve_bicgstab(this) ! ! BiCGStab Solver with reverse communication interface (RCI) ! implicit none type(bicgstab_alg), intent(inout) :: this ! ! index to real scalar variables ! integer, parameter :: idxrwork_rsrc = 1 ! ! index to complex scalar variables ! integer, parameter :: idxzwork_rho = 1 integer, parameter :: idxzwork_rho0 = 2 integer, parameter :: idxzwork_bet = 3 integer, parameter :: idxzwork_omg = 4 integer, parameter :: idxzwork_alp = 5 ! ! index to complex vector variables ! integer, parameter :: idx_b = 1 integer, parameter :: idx_x = 2 integer, parameter :: idx_r = 3 integer, parameter :: idx_rt= 4 integer, parameter :: idx_p = 5 integer, parameter :: idx_q = 6 integer, parameter :: idx_v = 7 integer, parameter :: idx_t = 8 integer :: isite complex(DP) :: rho,rho0,bet,omg,alp complex(DP) :: ctmp0 real(DP) :: rtmp0,rtmp1,rsrc if (this%debug) write(*,*)"INEXT=",this%inext select case (this%inext) case (0) this%iter=0 if (this%guess == GUESS_USE) then !----------------- ! p = A x !----------------- this%istat = OP_DO_MATVEC this%src_vec => this%vwork(idx_x)%v this%dst_vec => this%vwork(idx_p)%v else this%istat = OP_NOP this%src_vec => NULL() this%dst_vec => NULL() endif this%inext=this%inext+1 return case (1) if (this%guess == GUESS_USE) then !----------------------- ! r = b - p = b - A x !----------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) -this%vwork(idx_p)%v(isite) enddo else !----------------------- ! r = b ! x = 0 !----------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) this%vwork(idx_x)%v(isite) = Z0 enddo endif !--------------- ! rt = r ! p = r ! rho0 = <rt|r> ! rsrc = |b| !--------------- rtmp0 = 0.0_DP rtmp1 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0,rtmp1) do isite=1,this%NSIZE this%vwork(idx_rt)%v(isite) = this%vwork(idx_r)%v(isite) this%vwork(idx_p )%v(isite) = this%vwork(idx_r)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_b)%v(isite))**2 +aimag(this%vwork(idx_b)%v(isite))**2 rtmp1 = rtmp1 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) call comlib_sumcast(rtmp1) #endif rsrc = sqrt(rtmp0) rho0 = rtmp1 this%rwork(idxrwork_rsrc) = rsrc this%zwork(idxzwork_rho0) = rho0 this%res = sqrt(rtmp1)/rsrc this%istat = OP_PRINT_STATUS this%inext = this%inext+1 return case (2) !--------------------------------------- ! start main loop !--------------------------------------- if (this%mode == MODE_PRECOND) then !------------------- ! v = M p !------------------- this%src_vec => this%vwork(idx_p)%v this%dst_vec => this%vwork(idx_v)%v this%istat = OP_PRECOND else this%src_vec => NULL() this%dst_vec => NULL() this%istat = OP_NOP endif this%inext = this%inext+1 return case (3) if (this%mode == MODE_PRECOND) then !------------------- ! q = A v !------------------- this%src_vec => this%vwork(idx_v)%v this%dst_vec => this%vwork(idx_q)%v this%istat = OP_DO_MATVEC else !------------------- ! q = A p !------------------- this%src_vec => this%vwork(idx_p)%v this%dst_vec => this%vwork(idx_q)%v this%istat = OP_DO_MATVEC endif this%inext = this%inext+1 this%iter = this%iter+1 return case (4) !-------------------------------- ! alp = rho0 / <rt|q> !-------------------------------- ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_rt)%v(isite)) *this%vwork(idx_q )%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) #endif rho0 = this%zwork(idxzwork_rho0) if (this%debug) then write(*,'("rho0=",2E24.15," <rt|q>=",2E24.15)')rho0,ctmp0 endif alp = rho0/ctmp0 if (this%debug) then write(*,'("alp=",2E24.15)')alp endif this%zwork(idxzwork_alp) = alp if (this%mode == MODE_PRECOND) then !---------------- ! r = r - alp q ! x = x + alp v !---------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -alp*this%vwork(idx_q)%v(isite) this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +alp*this%vwork(idx_v)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 + aimag(this%vwork(idx_r)%v(isite))**2 enddo else !---------------- ! x = x + alp p ! r = r - alp q !---------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +alp*this%vwork(idx_p)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -alp*this%vwork(idx_q)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 + aimag(this%vwork(idx_r)%v(isite))**2 enddo endif #ifndef _singlePU call comlib_sumcast(rtmp0) #endif rsrc = this%rwork(idxrwork_rsrc) this%res = sqrt(rtmp0)/rsrc if (this%debug) write(*,'("res=",E24.15," tol=",E24.15)')this%res,this%tol if ( this%res <= this%tol ) then this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_CONVERGED this%inext = 0 else if ( this%iter <= this%max_iter ) then this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_PRINT_STATUS this%inext = this%inext+1 else this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_MAXITER_REACHED this%inext = 0 endif return case(5) if (this%mode == MODE_PRECOND) then !------------ ! v = M r !------------ this%src_vec => this%vwork(idx_r)%v this%dst_vec => this%vwork(idx_v)%v this%istat = OP_PRECOND else this%src_vec => NULL() this%dst_vec => NULL() this%istat = OP_NOP endif this%inext = this%inext+1 return case(6) if (this%mode == MODE_PRECOND) then !------------ ! t = A v !------------ this%src_vec => this%vwork(idx_v)%v this%dst_vec => this%vwork(idx_t)%v this%istat = OP_DO_MATVEC else !------------ ! t = A r !------------ this%src_vec => this%vwork(idx_r)%v this%dst_vec => this%vwork(idx_t)%v this%istat = OP_DO_MATVEC endif this%inext = this%inext+1 this%iter = this%iter+1 return case(7) !------------------------ ! omg = <t|r>/<t|t> !------------------------ rtmp0 = 0.0_DP ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0,ctmp0) do isite=1,this%NSIZE rtmp0 = rtmp0 + real(this%vwork(idx_t)%v(isite))**2 +aimag(this%vwork(idx_t)%v(isite))**2 ctmp0 = ctmp0 +conjg(this%vwork(idx_t)%v(isite)) *this%vwork(idx_r)%v(isite) enddo #ifndef _singlePU call comlib_sumcast(rtmp0) call comlib_sumcast(ctmp0) #endif if (this%debug) then write(*,'("<t|r>=",2E24.15," <t|t>=",E24.15)')ctmp0,rtmp0 endif omg = ctmp0/rtmp0 this%zwork(idxzwork_omg) = omg if (this%mode == MODE_PRECOND) then !---------------------------- ! x = x + omg v ! r = r - omg t ! res = |r|/rsrc !---------------------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +omg*this%vwork(idx_v)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -omg*this%vwork(idx_t)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo else !---------------------------- ! x = x + omg r ! r = r - omg t ! res = |r|/rsrc !---------------------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +omg*this%vwork(idx_r)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -omg*this%vwork(idx_t)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo endif #ifndef _singlePU call comlib_sumcast(rtmp0) #endif rsrc = this%rwork(idxrwork_rsrc) this%res=sqrt(rtmp0)/rsrc if (this%debug) write(*,'("res=",E24.15," tol=",E24.15)')this%res,this%tol if ( this%res <= this%tol ) then this%istat = OP_CONVERGED this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%inext = 0 else if ( this%iter <= this%max_iter ) then this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_PRINT_STATUS this%inext = this%inext+1 else this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_MAXITER_REACHED this%inext = 0 endif return case(8) !------------------------------- ! rho = <rt|r> ! bet = (alp/omg)*(rho/rho0) ! rho0 = rho !------------------------------- ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_rt)%v(isite)) *this%vwork(idx_r )%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) #endif rho = ctmp0 alp = this%zwork(idxzwork_alp) omg = this%zwork(idxzwork_omg) rho0 = this%zwork(idxzwork_rho0) if (this%debug) then write(*,'("alp=",2E24.15," omg=",2E24.15," rho=",2E24.15," rho0=",2E24.15)')alp,omg,rho,rho0 endif bet = (alp/omg)*(rho/rho0) rho0 = rho this%zwork(idxzwork_rho0) = rho0 !------------------------------- ! p = r + bet (p - omg q) !------------------------------- !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) +bet*( this%vwork(idx_p)%v(isite) -omg*this%vwork(idx_q)%v(isite)) enddo this%istat = OP_NOP this%inext = 2 return end select return end subroutine
Subroutine : | |
this : | type(bicgstab_alg), intent(inout) |
BiCGStab Solver with reverse communication interface (RCI)
subroutine solve_bicgstab(this) ! ! BiCGStab Solver with reverse communication interface (RCI) ! implicit none type(bicgstab_alg), intent(inout) :: this ! ! index to real scalar variables ! integer, parameter :: idxrwork_rsrc = 1 ! ! index to complex scalar variables ! integer, parameter :: idxzwork_rho = 1 integer, parameter :: idxzwork_rho0 = 2 integer, parameter :: idxzwork_bet = 3 integer, parameter :: idxzwork_omg = 4 integer, parameter :: idxzwork_alp = 5 ! ! index to complex vector variables ! integer, parameter :: idx_b = 1 integer, parameter :: idx_x = 2 integer, parameter :: idx_r = 3 integer, parameter :: idx_rt= 4 integer, parameter :: idx_p = 5 integer, parameter :: idx_q = 6 integer, parameter :: idx_v = 7 integer, parameter :: idx_t = 8 integer :: isite complex(DP) :: rho,rho0,bet,omg,alp complex(DP) :: ctmp0 real(DP) :: rtmp0,rtmp1,rsrc if (this%debug) write(*,*)"INEXT=",this%inext select case (this%inext) case (0) this%iter=0 if (this%guess == GUESS_USE) then !----------------- ! p = A x !----------------- this%istat = OP_DO_MATVEC this%src_vec => this%vwork(idx_x)%v this%dst_vec => this%vwork(idx_p)%v else this%istat = OP_NOP this%src_vec => NULL() this%dst_vec => NULL() endif this%inext=this%inext+1 return case (1) if (this%guess == GUESS_USE) then !----------------------- ! r = b - p = b - A x !----------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) -this%vwork(idx_p)%v(isite) enddo else !----------------------- ! r = b ! x = 0 !----------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) this%vwork(idx_x)%v(isite) = Z0 enddo endif !--------------- ! rt = r ! p = r ! rho0 = <rt|r> ! rsrc = |b| !--------------- rtmp0 = 0.0_DP rtmp1 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0,rtmp1) do isite=1,this%NSIZE this%vwork(idx_rt)%v(isite) = this%vwork(idx_r)%v(isite) this%vwork(idx_p )%v(isite) = this%vwork(idx_r)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_b)%v(isite))**2 +aimag(this%vwork(idx_b)%v(isite))**2 rtmp1 = rtmp1 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) call comlib_sumcast(rtmp1) #endif rsrc = sqrt(rtmp0) rho0 = rtmp1 this%rwork(idxrwork_rsrc) = rsrc this%zwork(idxzwork_rho0) = rho0 this%res = sqrt(rtmp1)/rsrc this%istat = OP_PRINT_STATUS this%inext = this%inext+1 return case (2) !--------------------------------------- ! start main loop !--------------------------------------- if (this%mode == MODE_PRECOND) then !------------------- ! v = M p !------------------- this%src_vec => this%vwork(idx_p)%v this%dst_vec => this%vwork(idx_v)%v this%istat = OP_PRECOND else this%src_vec => NULL() this%dst_vec => NULL() this%istat = OP_NOP endif this%inext = this%inext+1 return case (3) if (this%mode == MODE_PRECOND) then !------------------- ! q = A v !------------------- this%src_vec => this%vwork(idx_v)%v this%dst_vec => this%vwork(idx_q)%v this%istat = OP_DO_MATVEC else !------------------- ! q = A p !------------------- this%src_vec => this%vwork(idx_p)%v this%dst_vec => this%vwork(idx_q)%v this%istat = OP_DO_MATVEC endif this%inext = this%inext+1 this%iter = this%iter+1 return case (4) !-------------------------------- ! alp = rho0 / <rt|q> !-------------------------------- ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_rt)%v(isite)) *this%vwork(idx_q )%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) #endif rho0 = this%zwork(idxzwork_rho0) if (this%debug) then write(*,'("rho0=",2E24.15," <rt|q>=",2E24.15)')rho0,ctmp0 endif alp = rho0/ctmp0 if (this%debug) then write(*,'("alp=",2E24.15)')alp endif this%zwork(idxzwork_alp) = alp if (this%mode == MODE_PRECOND) then !---------------- ! r = r - alp q ! x = x + alp v !---------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -alp*this%vwork(idx_q)%v(isite) this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +alp*this%vwork(idx_v)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 + aimag(this%vwork(idx_r)%v(isite))**2 enddo else !---------------- ! x = x + alp p ! r = r - alp q !---------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +alp*this%vwork(idx_p)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -alp*this%vwork(idx_q)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 + aimag(this%vwork(idx_r)%v(isite))**2 enddo endif #ifndef _singlePU call comlib_sumcast(rtmp0) #endif rsrc = this%rwork(idxrwork_rsrc) this%res = sqrt(rtmp0)/rsrc if (this%debug) write(*,'("res=",E24.15," tol=",E24.15)')this%res,this%tol if ( this%res <= this%tol ) then this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_CONVERGED this%inext = 0 else if ( this%iter <= this%max_iter ) then this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_PRINT_STATUS this%inext = this%inext+1 else this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_MAXITER_REACHED this%inext = 0 endif return case(5) if (this%mode == MODE_PRECOND) then !------------ ! v = M r !------------ this%src_vec => this%vwork(idx_r)%v this%dst_vec => this%vwork(idx_v)%v this%istat = OP_PRECOND else this%src_vec => NULL() this%dst_vec => NULL() this%istat = OP_NOP endif this%inext = this%inext+1 return case(6) if (this%mode == MODE_PRECOND) then !------------ ! t = A v !------------ this%src_vec => this%vwork(idx_v)%v this%dst_vec => this%vwork(idx_t)%v this%istat = OP_DO_MATVEC else !------------ ! t = A r !------------ this%src_vec => this%vwork(idx_r)%v this%dst_vec => this%vwork(idx_t)%v this%istat = OP_DO_MATVEC endif this%inext = this%inext+1 this%iter = this%iter+1 return case(7) !------------------------ ! omg = <t|r>/<t|t> !------------------------ rtmp0 = 0.0_DP ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0,ctmp0) do isite=1,this%NSIZE rtmp0 = rtmp0 + real(this%vwork(idx_t)%v(isite))**2 +aimag(this%vwork(idx_t)%v(isite))**2 ctmp0 = ctmp0 +conjg(this%vwork(idx_t)%v(isite)) *this%vwork(idx_r)%v(isite) enddo #ifndef _singlePU call comlib_sumcast(rtmp0) call comlib_sumcast(ctmp0) #endif if (this%debug) then write(*,'("<t|r>=",2E24.15," <t|t>=",E24.15)')ctmp0,rtmp0 endif omg = ctmp0/rtmp0 this%zwork(idxzwork_omg) = omg if (this%mode == MODE_PRECOND) then !---------------------------- ! x = x + omg v ! r = r - omg t ! res = |r|/rsrc !---------------------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +omg*this%vwork(idx_v)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -omg*this%vwork(idx_t)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo else !---------------------------- ! x = x + omg r ! r = r - omg t ! res = |r|/rsrc !---------------------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +omg*this%vwork(idx_r)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -omg*this%vwork(idx_t)%v(isite) rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo endif #ifndef _singlePU call comlib_sumcast(rtmp0) #endif rsrc = this%rwork(idxrwork_rsrc) this%res=sqrt(rtmp0)/rsrc if (this%debug) write(*,'("res=",E24.15," tol=",E24.15)')this%res,this%tol if ( this%res <= this%tol ) then this%istat = OP_CONVERGED this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%inext = 0 else if ( this%iter <= this%max_iter ) then this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_PRINT_STATUS this%inext = this%inext+1 else this%src_vec => NULL() this%dst_vec => this%vwork(idx_x)%v this%istat = OP_MAXITER_REACHED this%inext = 0 endif return case(8) !------------------------------- ! rho = <rt|r> ! bet = (alp/omg)*(rho/rho0) ! rho0 = rho !------------------------------- ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_rt)%v(isite)) *this%vwork(idx_r )%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) #endif rho = ctmp0 alp = this%zwork(idxzwork_alp) omg = this%zwork(idxzwork_omg) rho0 = this%zwork(idxzwork_rho0) if (this%debug) then write(*,'("alp=",2E24.15," omg=",2E24.15," rho=",2E24.15," rho0=",2E24.15)')alp,omg,rho,rho0 endif bet = (alp/omg)*(rho/rho0) rho0 = rho this%zwork(idxzwork_rho0) = rho0 !------------------------------- ! p = r + bet (p - omg q) !------------------------------- !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) +bet*( this%vwork(idx_p)%v(isite) -omg*this%vwork(idx_q)%v(isite)) enddo this%istat = OP_NOP this%inext = 2 return end select return end subroutine
Subroutine : | |
this : | type(cg_alg), intent(inout) |
CG Solver with reverse communication interface (RCI)
subroutine solve_cg(this) ! ! CG Solver with reverse communication interface (RCI) ! implicit none type(cg_alg), intent(inout) :: this ! ! index to scalar variables ! integer, parameter :: idxrwork_rsrc =1 integer, parameter :: idxrwork_rho =2 integer, parameter :: idxrwork_rho0 =3 integer, parameter :: idxrwork_res0 =4 integer, parameter :: idxrwork_flag =5 ! ! index to vector variables ! integer, parameter :: idx_b =1 integer, parameter :: idx_x =2 integer, parameter :: idx_r =3 integer, parameter :: idx_p =4 integer, parameter :: idx_q =5 integer, parameter :: idx_v =6 integer :: isite real(DP) :: rho,rho0,beta complex(DP) :: ctmp0,alpha,cbeta,ctmp1 real(DP) :: rtmp0,rtmp1,rsrc select case (this%inext) case (0) this%iter=0 if (this%guess == GUESS_USE) then !----------------- ! p = A x !----------------- this%istat = OP_DO_MATVEC this%src_vec => this%vwork(idx_x)%v this%dst_vec => this%vwork(idx_p)%v else this%istat = OP_NOP this%src_vec => NULL() this%dst_vec => NULL() endif this%inext=this%inext+1 return case (1) if (this%guess == GUESS_USE) then !---------------------- ! r = b - p = b - A x !---------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) -this%vwork(idx_p)%v(isite) enddo else !---------------------- ! r = b, x = 0 !---------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) this%vwork(idx_x)%v(isite) = Z0 enddo endif call projection() this%inext = this%inext+1 return case (2) !----------------- ! p = r ! rsrc = |b| ! rho0 = |r|^2 !----------------- rtmp0 = 0.0_DP rtmp1 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0,rtmp1) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) rtmp0 = rtmp0+ real(this%vwork(idx_b)%v(isite))**2 +aimag(this%vwork(idx_b)%v(isite))**2 rtmp1 = rtmp1+ real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) call comlib_sumcast(rtmp1) #endif rsrc=sqrt(rtmp0) rho0=rtmp1 this%res = sqrt(rtmp1)/rsrc this%istat = OP_PRINT_STATUS this%inext = this%inext+1 this%rwork(idxrwork_rsrc) = rsrc this%rwork(idxrwork_rho0) = rho0 this%rwork(idxrwork_res0) = this%res this%rwork(idxrwork_flag) = 0 return case (3) !----------------------- ! start main loop !----------------------- this%iter=this%iter+1 !------------- ! q = A p !------------- this%src_vec => this%vwork(idx_p)%v this%dst_vec => this%vwork(idx_q)%v this%istat = OP_DO_MATVEC this%inext = this%inext+1 return case (4) !--------------------- ! alpha = rho0/<p|q> ! alpha = rho0/<r|q> !--------------------- ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_p)%v(isite)) *this%vwork(idx_q)%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) #endif rho0 = this%rwork(idxrwork_rho0) alpha=rho0/ctmp0 !------------------- ! x = x + alpha p ! r = r - alpha q ! rho = |r|^2 ! err = |r|/|b| !------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +alpha*this%vwork(idx_p)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -alpha*this%vwork(idx_q)%v(isite) enddo call projection() this%inext = this%inext+1 return case (5) rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) #endif rho = rtmp0 rsrc = this%rwork(idxrwork_rsrc) this%res = sqrt(rtmp0)/rsrc this%rwork(idxrwork_rho) = rho if ( this%res <= this%tol ) then this%istat = OP_CONVERGED this%dst_vec => this%vwork(idx_x)%v this%inext = 0 else if ( this%iter <= this%max_iter ) then this%istat = OP_PRINT_STATUS this%inext = this%inext+1 else this%istat = OP_MAXITER_REACHED this%inext = 0 endif return case (6) !----------------- ! true residual ! r = A x !----------------- if ( this%res/this%rwork(idxrwork_res0) < 1.0e-5_DP) then this%rwork(idxrwork_res0) = this%res this%rwork(idxrwork_flag) = 1 this%istat = OP_DO_MATVEC this%src_vec => this%vwork(idx_x)%v this%dst_vec => this%vwork(idx_r)%v this%inext = this%inext+1 else this%rwork(idxrwork_flag) = 0 this%istat = OP_NOP this%inext = this%inext+2 endif return case (7) !----------------- ! true residual ! r = b - A x ! = b - r !----------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) -this%vwork(idx_r)%v(isite) rtmp0 = rtmp0+ real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) #endif this%rwork(idxrwork_rho) = rtmp0 rsrc = this%rwork(idxrwork_rsrc) this%res = sqrt(rtmp0)/rsrc call projection() this%inext = this%inext+1 return case (8) rho = this%rwork(idxrwork_rho) rho0 = this%rwork(idxrwork_rho0) !------------------- ! p = r + beta p !------------------- if (this%rwork(idxrwork_flag) == 0) then beta = rho/rho0 rho0 = rho !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) +beta*this%vwork(idx_p)%v(isite) enddo else ctmp0 = Z0 ctmp1 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0,ctmp1) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_r)%v(isite)) *this%vwork(idx_q)%v(isite) ctmp1 = ctmp1 + conjg(this%vwork(idx_p)%v(isite)) *this%vwork(idx_q)%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) call comlib_sumcast(ctmp1) #endif cbeta =-conjg(ctmp0/ctmp1) rho0 = rho !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) +cbeta*this%vwork(idx_p)%v(isite) enddo endif this%istat = OP_NOP this%inext = 3 ! back to iteration loop this%rwork(idxrwork_rho) = rho this%rwork(idxrwork_rho0) = rho0 return end select return contains subroutine projection() if (this%mode == MODE_PROJECTION) then !---------------------- ! r <= r - A*V*E\V'*r ! x <= x + V*E\V'*r !---------------------- this%istat = OP_PROJECTION this%src_vec => this%vwork(idx_r)%v this%dst_vec => this%vwork(idx_x)%v else this%istat = OP_NOP endif return end subroutine end subroutine
Subroutine : | |
this : | type(cg_alg), intent(inout) |
CG Solver with reverse communication interface (RCI)
subroutine solve_cg(this) ! ! CG Solver with reverse communication interface (RCI) ! implicit none type(cg_alg), intent(inout) :: this ! ! index to scalar variables ! integer, parameter :: idxrwork_rsrc =1 integer, parameter :: idxrwork_rho =2 integer, parameter :: idxrwork_rho0 =3 integer, parameter :: idxrwork_res0 =4 integer, parameter :: idxrwork_flag =5 ! ! index to vector variables ! integer, parameter :: idx_b =1 integer, parameter :: idx_x =2 integer, parameter :: idx_r =3 integer, parameter :: idx_p =4 integer, parameter :: idx_q =5 integer, parameter :: idx_v =6 integer :: isite real(DP) :: rho,rho0,beta complex(DP) :: ctmp0,alpha,cbeta,ctmp1 real(DP) :: rtmp0,rtmp1,rsrc select case (this%inext) case (0) this%iter=0 if (this%guess == GUESS_USE) then !----------------- ! p = A x !----------------- this%istat = OP_DO_MATVEC this%src_vec => this%vwork(idx_x)%v this%dst_vec => this%vwork(idx_p)%v else this%istat = OP_NOP this%src_vec => NULL() this%dst_vec => NULL() endif this%inext=this%inext+1 return case (1) if (this%guess == GUESS_USE) then !---------------------- ! r = b - p = b - A x !---------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) -this%vwork(idx_p)%v(isite) enddo else !---------------------- ! r = b, x = 0 !---------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) this%vwork(idx_x)%v(isite) = Z0 enddo endif call projection() this%inext = this%inext+1 return case (2) !----------------- ! p = r ! rsrc = |b| ! rho0 = |r|^2 !----------------- rtmp0 = 0.0_DP rtmp1 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0,rtmp1) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) rtmp0 = rtmp0+ real(this%vwork(idx_b)%v(isite))**2 +aimag(this%vwork(idx_b)%v(isite))**2 rtmp1 = rtmp1+ real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) call comlib_sumcast(rtmp1) #endif rsrc=sqrt(rtmp0) rho0=rtmp1 this%res = sqrt(rtmp1)/rsrc this%istat = OP_PRINT_STATUS this%inext = this%inext+1 this%rwork(idxrwork_rsrc) = rsrc this%rwork(idxrwork_rho0) = rho0 this%rwork(idxrwork_res0) = this%res this%rwork(idxrwork_flag) = 0 return case (3) !----------------------- ! start main loop !----------------------- this%iter=this%iter+1 !------------- ! q = A p !------------- this%src_vec => this%vwork(idx_p)%v this%dst_vec => this%vwork(idx_q)%v this%istat = OP_DO_MATVEC this%inext = this%inext+1 return case (4) !--------------------- ! alpha = rho0/<p|q> ! alpha = rho0/<r|q> !--------------------- ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_p)%v(isite)) *this%vwork(idx_q)%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) #endif rho0 = this%rwork(idxrwork_rho0) alpha=rho0/ctmp0 !------------------- ! x = x + alpha p ! r = r - alpha q ! rho = |r|^2 ! err = |r|/|b| !------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +alpha*this%vwork(idx_p)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -alpha*this%vwork(idx_q)%v(isite) enddo call projection() this%inext = this%inext+1 return case (5) rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) #endif rho = rtmp0 rsrc = this%rwork(idxrwork_rsrc) this%res = sqrt(rtmp0)/rsrc this%rwork(idxrwork_rho) = rho if ( this%res <= this%tol ) then this%istat = OP_CONVERGED this%dst_vec => this%vwork(idx_x)%v this%inext = 0 else if ( this%iter <= this%max_iter ) then this%istat = OP_PRINT_STATUS this%inext = this%inext+1 else this%istat = OP_MAXITER_REACHED this%inext = 0 endif return case (6) !----------------- ! true residual ! r = A x !----------------- if ( this%res/this%rwork(idxrwork_res0) < 1.0e-5_DP) then this%rwork(idxrwork_res0) = this%res this%rwork(idxrwork_flag) = 1 this%istat = OP_DO_MATVEC this%src_vec => this%vwork(idx_x)%v this%dst_vec => this%vwork(idx_r)%v this%inext = this%inext+1 else this%rwork(idxrwork_flag) = 0 this%istat = OP_NOP this%inext = this%inext+2 endif return case (7) !----------------- ! true residual ! r = b - A x ! = b - r !----------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) -this%vwork(idx_r)%v(isite) rtmp0 = rtmp0+ real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) #endif this%rwork(idxrwork_rho) = rtmp0 rsrc = this%rwork(idxrwork_rsrc) this%res = sqrt(rtmp0)/rsrc call projection() this%inext = this%inext+1 return case (8) rho = this%rwork(idxrwork_rho) rho0 = this%rwork(idxrwork_rho0) !------------------- ! p = r + beta p !------------------- if (this%rwork(idxrwork_flag) == 0) then beta = rho/rho0 rho0 = rho !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) +beta*this%vwork(idx_p)%v(isite) enddo else ctmp0 = Z0 ctmp1 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0,ctmp1) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_r)%v(isite)) *this%vwork(idx_q)%v(isite) ctmp1 = ctmp1 + conjg(this%vwork(idx_p)%v(isite)) *this%vwork(idx_q)%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) call comlib_sumcast(ctmp1) #endif cbeta =-conjg(ctmp0/ctmp1) rho0 = rho !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) +cbeta*this%vwork(idx_p)%v(isite) enddo endif this%istat = OP_NOP this%inext = 3 ! back to iteration loop this%rwork(idxrwork_rho) = rho this%rwork(idxrwork_rho0) = rho0 return end select return contains subroutine projection() if (this%mode == MODE_PROJECTION) then !---------------------- ! r <= r - A*V*E\V'*r ! x <= x + V*E\V'*r !---------------------- this%istat = OP_PROJECTION this%src_vec => this%vwork(idx_r)%v this%dst_vec => this%vwork(idx_x)%v else this%istat = OP_NOP endif return end subroutine end subroutine
Subroutine : | |
this : | type(cg_alg), intent(inout) |
CG Solver with reverse communication interface (RCI)
subroutine solve_cg(this) ! ! CG Solver with reverse communication interface (RCI) ! implicit none type(cg_alg), intent(inout) :: this ! ! index to scalar variables ! integer, parameter :: idxrwork_rsrc =1 integer, parameter :: idxrwork_rho =2 integer, parameter :: idxrwork_rho0 =3 integer, parameter :: idxrwork_res0 =4 integer, parameter :: idxrwork_flag =5 ! ! index to vector variables ! integer, parameter :: idx_b =1 integer, parameter :: idx_x =2 integer, parameter :: idx_r =3 integer, parameter :: idx_p =4 integer, parameter :: idx_q =5 integer, parameter :: idx_v =6 integer :: isite real(DP) :: rho,rho0,beta complex(DP) :: ctmp0,alpha,cbeta,ctmp1 real(DP) :: rtmp0,rtmp1,rsrc select case (this%inext) case (0) this%iter=0 if (this%guess == GUESS_USE) then !----------------- ! p = A x !----------------- this%istat = OP_DO_MATVEC this%src_vec => this%vwork(idx_x)%v this%dst_vec => this%vwork(idx_p)%v else this%istat = OP_NOP this%src_vec => NULL() this%dst_vec => NULL() endif this%inext=this%inext+1 return case (1) if (this%guess == GUESS_USE) then !---------------------- ! r = b - p = b - A x !---------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) -this%vwork(idx_p)%v(isite) enddo else !---------------------- ! r = b, x = 0 !---------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) this%vwork(idx_x)%v(isite) = Z0 enddo endif call projection() this%inext = this%inext+1 return case (2) !----------------- ! p = r ! rsrc = |b| ! rho0 = |r|^2 !----------------- rtmp0 = 0.0_DP rtmp1 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0,rtmp1) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) rtmp0 = rtmp0+ real(this%vwork(idx_b)%v(isite))**2 +aimag(this%vwork(idx_b)%v(isite))**2 rtmp1 = rtmp1+ real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) call comlib_sumcast(rtmp1) #endif rsrc=sqrt(rtmp0) rho0=rtmp1 this%res = sqrt(rtmp1)/rsrc this%istat = OP_PRINT_STATUS this%inext = this%inext+1 this%rwork(idxrwork_rsrc) = rsrc this%rwork(idxrwork_rho0) = rho0 this%rwork(idxrwork_res0) = this%res this%rwork(idxrwork_flag) = 0 return case (3) !----------------------- ! start main loop !----------------------- this%iter=this%iter+1 !------------- ! q = A p !------------- this%src_vec => this%vwork(idx_p)%v this%dst_vec => this%vwork(idx_q)%v this%istat = OP_DO_MATVEC this%inext = this%inext+1 return case (4) !--------------------- ! alpha = rho0/<p|q> ! alpha = rho0/<r|q> !--------------------- ctmp0 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_p)%v(isite)) *this%vwork(idx_q)%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) #endif rho0 = this%rwork(idxrwork_rho0) alpha=rho0/ctmp0 !------------------- ! x = x + alpha p ! r = r - alpha q ! rho = |r|^2 ! err = |r|/|b| !------------------- !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_x)%v(isite) = this%vwork(idx_x)%v(isite) +alpha*this%vwork(idx_p)%v(isite) this%vwork(idx_r)%v(isite) = this%vwork(idx_r)%v(isite) -alpha*this%vwork(idx_q)%v(isite) enddo call projection() this%inext = this%inext+1 return case (5) rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE rtmp0 = rtmp0 + real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) #endif rho = rtmp0 rsrc = this%rwork(idxrwork_rsrc) this%res = sqrt(rtmp0)/rsrc this%rwork(idxrwork_rho) = rho if ( this%res <= this%tol ) then this%istat = OP_CONVERGED this%dst_vec => this%vwork(idx_x)%v this%inext = 0 else if ( this%iter <= this%max_iter ) then this%istat = OP_PRINT_STATUS this%inext = this%inext+1 else this%istat = OP_MAXITER_REACHED this%inext = 0 endif return case (6) !----------------- ! true residual ! r = A x !----------------- if ( this%res/this%rwork(idxrwork_res0) < 1.0e-5_DP) then this%rwork(idxrwork_res0) = this%res this%rwork(idxrwork_flag) = 1 this%istat = OP_DO_MATVEC this%src_vec => this%vwork(idx_x)%v this%dst_vec => this%vwork(idx_r)%v this%inext = this%inext+1 else this%rwork(idxrwork_flag) = 0 this%istat = OP_NOP this%inext = this%inext+2 endif return case (7) !----------------- ! true residual ! r = b - A x ! = b - r !----------------- rtmp0 = 0.0_DP !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:rtmp0) do isite=1,this%NSIZE this%vwork(idx_r)%v(isite) = this%vwork(idx_b)%v(isite) -this%vwork(idx_r)%v(isite) rtmp0 = rtmp0+ real(this%vwork(idx_r)%v(isite))**2 +aimag(this%vwork(idx_r)%v(isite))**2 enddo #ifndef _singlePU call comlib_sumcast(rtmp0) #endif this%rwork(idxrwork_rho) = rtmp0 rsrc = this%rwork(idxrwork_rsrc) this%res = sqrt(rtmp0)/rsrc call projection() this%inext = this%inext+1 return case (8) rho = this%rwork(idxrwork_rho) rho0 = this%rwork(idxrwork_rho0) !------------------- ! p = r + beta p !------------------- if (this%rwork(idxrwork_flag) == 0) then beta = rho/rho0 rho0 = rho !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) +beta*this%vwork(idx_p)%v(isite) enddo else ctmp0 = Z0 ctmp1 = Z0 !$OMP PARALLEL DO PRIVATE(isite) REDUCTION(+:ctmp0,ctmp1) do isite=1,this%NSIZE ctmp0 = ctmp0 + conjg(this%vwork(idx_r)%v(isite)) *this%vwork(idx_q)%v(isite) ctmp1 = ctmp1 + conjg(this%vwork(idx_p)%v(isite)) *this%vwork(idx_q)%v(isite) enddo #ifndef _singlePU call comlib_sumcast(ctmp0) call comlib_sumcast(ctmp1) #endif cbeta =-conjg(ctmp0/ctmp1) rho0 = rho !$OMP PARALLEL DO PRIVATE(isite) do isite=1,this%NSIZE this%vwork(idx_p)%v(isite) = this%vwork(idx_r)%v(isite) +cbeta*this%vwork(idx_p)%v(isite) enddo endif this%istat = OP_NOP this%inext = 3 ! back to iteration loop this%rwork(idxrwork_rho) = rho this%rwork(idxrwork_rho0) = rho0 return end select return contains subroutine projection() if (this%mode == MODE_PROJECTION) then !---------------------- ! r <= r - A*V*E\V'*r ! x <= x + V*E\V'*r !---------------------- this%istat = OP_PROJECTION this%src_vec => this%vwork(idx_r)%v this%dst_vec => this%vwork(idx_x)%v else this%istat = OP_NOP endif return end subroutine end subroutine
Derived Type : | |||
src_vec(:) => NULL() : | complex(DP), public, pointer
| ||
dst_vec(:) => NULL() : | complex(DP), public, pointer
|
base type for solver algorithm