You cannot select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
amg4psblas/amgprec/amg_d_sludist_solver.F90

494 lines
15 KiB
Fortran

!
!
! AMG4PSBLAS version 1.0
! Algebraic Multigrid Package
! based on PSBLAS (Parallel Sparse BLAS version 3.7)
!
! (C) Copyright 2021
!
! Salvatore Filippone
! Pasqua D'Ambra
! Fabio Durastante
!
! Redistribution and use in source and binary forms, with or without
! modification, are permitted provided that the following conditions
! are met:
! 1. Redistributions of source code must retain the above copyright
! notice, this list of conditions and the following disclaimer.
! 2. Redistributions in binary form must reproduce the above copyright
! notice, this list of conditions, and the following disclaimer in the
! documentation and/or other materials provided with the distribution.
! 3. The name of the AMG4PSBLAS group or the names of its contributors may
! not be used to endorse or promote products derived from this
! software without specific written permission.
!
! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
! ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
! TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
! PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE AMG4PSBLAS GROUP OR ITS CONTRIBUTORS
! BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
! CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
! SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
! INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
! CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
! ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
! POSSIBILITY OF SUCH DAMAGE.
!
!
!
! File: amg_d_sludist_solver_mod.f90
!
! Module: amg_d_sludist_solver_mod
!
! This module defines:
! - the amg_d_sludist_solver_type data structure containing the ingredients
! to interface with the SuperLU_Dist package.
! 1. The factorization is distributed (and thus exact)
!
!
!
module amg_d_sludist_solver
use iso_c_binding
use amg_d_base_solver_mod
#if (!defined(HAVE_SLUDIST_)) || defined(IPK8)
type, extends(amg_d_base_solver_type) :: amg_d_sludist_solver_type
end type amg_d_sludist_solver_type
#else
type, extends(amg_d_base_solver_type) :: amg_d_sludist_solver_type
type(c_ptr) :: lufactors=c_null_ptr
integer(c_long_long) :: symbsize=0, numsize=0
contains
procedure, pass(sv) :: build => d_sludist_solver_bld
procedure, pass(sv) :: apply_a => d_sludist_solver_apply
procedure, pass(sv) :: apply_v => d_sludist_solver_apply_vect
procedure, pass(sv) :: free => d_sludist_solver_free
procedure, pass(sv) :: clear_data => d_sludist_solver_clear_data
procedure, pass(sv) :: descr => d_sludist_solver_descr
procedure, pass(sv) :: sizeof => d_sludist_solver_sizeof
procedure, nopass :: get_fmt => d_sludist_solver_get_fmt
procedure, nopass :: get_id => d_sludist_solver_get_id
procedure, pass(sv) :: is_global => d_sludist_solver_is_global
final :: d_sludist_solver_finalize
end type amg_d_sludist_solver_type
private :: d_sludist_solver_bld, d_sludist_solver_apply, &
& d_sludist_solver_free, d_sludist_solver_descr, &
& d_sludist_solver_sizeof, d_sludist_solver_apply_vect, &
& d_sludist_solver_get_fmt, d_sludist_solver_get_id, &
& d_sludist_solver_is_global, d_sludist_solver_clear_data
private :: d_sludist_solver_finalize
interface
function amg_dsludist_fact(n,nl,nnz,ifrst, &
& values,rowptr,colind,lufactors,npr,npc) &
& bind(c,name='amg_dsludist_fact') result(info)
use iso_c_binding
integer(c_int), value :: n,nl,nnz,ifrst,npr,npc
integer(c_int) :: info
integer(c_int) :: rowptr(*),colind(*)
real(c_double) :: values(*)
type(c_ptr) :: lufactors
end function amg_dsludist_fact
end interface
interface
function amg_dsludist_solve(itrans,n,nrhs, b, ldb, lufactors)&
& bind(c,name='amg_dsludist_solve') result(info)
use iso_c_binding
integer(c_int) :: info
integer(c_int), value :: itrans,n,nrhs,ldb
real(c_double) :: b(ldb,*)
type(c_ptr), value :: lufactors
end function amg_dsludist_solve
end interface
interface
function amg_dsludist_free(lufactors)&
& bind(c,name='amg_dsludist_free') result(info)
use iso_c_binding
integer(c_int) :: info
type(c_ptr), value :: lufactors
end function amg_dsludist_free
end interface
contains
subroutine d_sludist_solver_apply(alpha,sv,x,beta,y,desc_data,&
& trans,work,info,init,initu)
use psb_base_mod
implicit none
type(psb_desc_type), intent(in) :: desc_data
class(amg_d_sludist_solver_type), intent(inout) :: sv
real(psb_dpk_),intent(inout) :: x(:)
real(psb_dpk_),intent(inout) :: y(:)
real(psb_dpk_),intent(in) :: alpha,beta
character(len=1),intent(in) :: trans
real(psb_dpk_),target, intent(inout) :: work(:)
integer, intent(out) :: info
character, intent(in), optional :: init
real(psb_dpk_),intent(inout), optional :: initu(:)
integer :: n_row,n_col
real(psb_dpk_), pointer :: ww(:)
type(psb_ctxt_type) :: ctxt
integer :: np,me,i, err_act
character :: trans_
character(len=20) :: name='d_sludist_solver_apply'
call psb_erractionsave(err_act)
info = psb_success_
trans_ = psb_toupper(trans)
select case(trans_)
case('N')
case('T','C')
case default
call psb_errpush(psb_err_iarg_invalid_i_,name)
goto 9999
end select
!
! For non-iterative solvers, init and initu are ignored.
!
n_row = desc_data%get_local_rows()
n_col = desc_data%get_local_cols()
if (n_col <= size(work)) then
ww => work(1:n_col)
else
allocate(ww(n_col),stat=info)
if (info /= psb_success_) then
info=psb_err_alloc_request_
call psb_errpush(info,name,i_err=(/n_col/),&
& a_err='real(psb_dpk_)')
goto 9999
end if
endif
if (info == psb_success_)&
& call psb_geaxpby(done,x,dzero,ww,desc_data,info)
select case(trans_)
case('N')
info = amg_dsludist_solve(0,n_row,1,ww,n_row,sv%lufactors)
case('T')
info = amg_dsludist_solve(1,n_row,1,ww,n_row,sv%lufactors)
case('C')
info = amg_dsludist_solve(2,n_row,1,ww,n_row,sv%lufactors)
case default
call psb_errpush(psb_err_internal_error_,&
& name,a_err='Invalid TRANS in subsolve')
goto 9999
end select
if (info == psb_success_)&
& call psb_geaxpby(alpha,ww,beta,y,desc_data,info)
if (info /= psb_success_) then
call psb_errpush(psb_err_internal_error_,&
& name,a_err='Error in subsolve')
goto 9999
endif
if (n_col > size(work)) then
deallocate(ww)
endif
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(err_act)
return
end subroutine d_sludist_solver_apply
subroutine d_sludist_solver_apply_vect(alpha,sv,x,beta,y,desc_data,&
& trans,work,wv,info,init,initu)
use psb_base_mod
implicit none
type(psb_desc_type), intent(in) :: desc_data
class(amg_d_sludist_solver_type), intent(inout) :: sv
type(psb_d_vect_type),intent(inout) :: x
type(psb_d_vect_type),intent(inout) :: y
real(psb_dpk_),intent(in) :: alpha,beta
character(len=1),intent(in) :: trans
real(psb_dpk_),target, intent(inout) :: work(:)
type(psb_d_vect_type),intent(inout) :: wv(:)
integer, intent(out) :: info
character, intent(in), optional :: init
type(psb_d_vect_type),intent(inout), optional :: initu
integer :: err_act
character(len=20) :: name='d_sludist_solver_apply_vect'
call psb_erractionsave(err_act)
info = psb_success_
!
! For non-iterative solvers, init and initu are ignored.
!
call x%v%sync()
call y%v%sync()
call sv%apply(alpha,x%v%v,beta,y%v%v,desc_data,trans,work,info)
call y%v%set_host()
if (info /= 0) goto 9999
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(err_act)
return
end subroutine d_sludist_solver_apply_vect
subroutine d_sludist_solver_bld(a,desc_a,sv,info,b,amold,vmold,imold)
use psb_base_mod
Implicit None
! Arguments
type(psb_dspmat_type), intent(in), target :: a
Type(psb_desc_type), Intent(inout) :: desc_a
class(amg_d_sludist_solver_type), intent(inout) :: sv
integer, intent(out) :: info
type(psb_dspmat_type), intent(in), target, optional :: b
class(psb_d_base_sparse_mat), intent(in), optional :: amold
class(psb_d_base_vect_type), intent(in), optional :: vmold
class(psb_i_base_vect_type), intent(in), optional :: imold
! Local variables
type(psb_dspmat_type) :: atmp
type(psb_d_csr_sparse_mat) :: acsr
type(psb_ctxt_type) :: ctxt
integer(psb_lpk_), allocatable :: gia(:), gja(:)
integer(psb_lpk_) :: lfrst
integer(psb_ipk_) :: n_row,n_col, nrow_a, nztota, nglob, nzt, npr, npc
integer(psb_ipk_) :: ifrst, ibcheck
integer(psb_ipk_) :: np,me,i, err_act, debug_unit, debug_level
character(len=20) :: name='d_sludist_solver_bld', ch_err
info=psb_success_
call psb_erractionsave(err_act)
debug_unit = psb_get_debug_unit()
debug_level = psb_get_debug_level()
ctxt = desc_a%get_context()
call psb_info(ctxt, me, np)
npr = np
npc = 1
if (debug_level >= psb_debug_outer_) &
& write(debug_unit,*) me,' ',trim(name),' start'
n_row = desc_a%get_local_rows()
n_col = desc_a%get_local_cols()
nglob = desc_a%get_global_rows()
!
! Strategy here is as follows: because a call to SLUDIST
! as a gobal solver is mostly done at the coarsest level,
! even if we start from a problem requiring 8 bytes, chances
! are that the global size will be suitable for 4 bytes
! anyway, so we hope for the best, and throw an error
! if something goes wrong.
!
if (nglob > huge(1_psb_ipk_)) then
write(0,*) me,' ',trim(name),': Error: overflow of local indices '
info=psb_err_internal_error_
call psb_errpush(info,name)
goto 9999
end if
call a%cscnv(atmp,info,type='csr')
! This in case we are dealing with AS
call psb_rwextd(n_row,atmp,info,b=b)
call atmp%mv_to(acsr)
nrow_a = acsr%get_nrows()
nztota = acsr%get_nzeros()
call psb_loc_to_glob(ione,lfrst,desc_a,info)
! Fix the entries to call C-base SuperLU
call psb_realloc(nztota,gja,info)
call psb_loc_to_glob(acsr%ja(1:nztota),gja(1:nztota), desc_a, info, iact='I')
acsr%ja(1:nztota) = gja(1:nztota)
acsr%ja(:) = acsr%ja(:) - 1
acsr%irp(:) = acsr%irp(:) - 1
ifrst = lfrst - 1
info = amg_dsludist_fact(nglob,nrow_a,nztota,ifrst,&
& acsr%val,acsr%irp,acsr%ja,sv%lufactors,&
& npr,npc)
if (info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='amg_dsludist_fact'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
call acsr%free()
if (debug_level >= psb_debug_outer_) &
& write(debug_unit,*) me,' ',trim(name),' end'
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(err_act)
return
end subroutine d_sludist_solver_bld
subroutine d_sludist_solver_free(sv,info)
Implicit None
! Arguments
class(amg_d_sludist_solver_type), intent(inout) :: sv
integer, intent(out) :: info
Integer :: err_act
character(len=20) :: name='d_sludist_solver_free'
call psb_erractionsave(err_act)
info = 0
call sv%clear_data(info)
if (info /= psb_success_) goto 9999
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(err_act)
return
end subroutine d_sludist_solver_free
subroutine d_sludist_solver_clear_data(sv,info)
Implicit None
! Arguments
class(amg_d_sludist_solver_type), intent(inout) :: sv
integer, intent(out) :: info
Integer :: err_act
character(len=20) :: name='d_sludist_solver_clear_data'
call psb_erractionsave(err_act)
info = psb_success_
if (c_associated(sv%lufactors)) info = amg_dsludist_free(sv%lufactors)
sv%lufactors = c_null_ptr
if (info /= psb_success_) goto 9999
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(err_act)
return
end subroutine d_sludist_solver_clear_data
!
function d_sludist_solver_is_global(sv) result(val)
implicit none
class(amg_d_sludist_solver_type), intent(in) :: sv
logical :: val
val = .true.
end function d_sludist_solver_is_global
subroutine d_sludist_solver_finalize(sv)
Implicit None
! Arguments
type(amg_d_sludist_solver_type), intent(inout) :: sv
integer :: info
Integer :: err_act
character(len=20) :: name='d_sludist_solver_finalize'
call sv%free(info)
return
end subroutine d_sludist_solver_finalize
subroutine d_sludist_solver_descr(sv,info,iout,coarse,prefix)
Implicit None
! Arguments
class(amg_d_sludist_solver_type), intent(in) :: sv
integer, intent(out) :: info
integer, intent(in), optional :: iout
logical, intent(in), optional :: coarse
character(len=*), intent(in), optional :: prefix
! Local variables
integer :: err_act
type(psb_ctxt_type) :: ctxt
integer :: me, np
character(len=20), parameter :: name='amg_d_sludist_solver_descr'
integer :: iout_
character(1024) :: prefix_
call psb_erractionsave(err_act)
info = psb_success_
if (present(iout)) then
iout_ = iout
else
iout_ = psb_out_unit
endif
if (present(prefix)) then
prefix_ = prefix
else
prefix_ = ""
end if
write(iout_,*) trim(prefix_), ' SuperLU_Dist Sparse Factorization Solver. '
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(err_act)
return
end subroutine d_sludist_solver_descr
function d_sludist_solver_sizeof(sv) result(val)
implicit none
! Arguments
class(amg_d_sludist_solver_type), intent(in) :: sv
integer(psb_epk_) :: val
integer :: i
val = 2*psb_sizeof_ip + psb_sizeof_dp
val = val + sv%symbsize
val = val + sv%numsize
return
end function d_sludist_solver_sizeof
function d_sludist_solver_get_fmt() result(val)
implicit none
character(len=32) :: val
val = "SuperLU_Dist solver"
end function d_sludist_solver_get_fmt
function d_sludist_solver_get_id() result(val)
implicit none
integer(psb_ipk_) :: val
val = amg_sludist_
end function d_sludist_solver_get_id
#endif
end module amg_d_sludist_solver