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.
399 lines
15 KiB
Fortran
399 lines
15 KiB
Fortran
8 years ago
|
!
|
||
|
!
|
||
7 years ago
|
! MLD2P4 version 2.2
|
||
8 years ago
|
! MultiLevel Domain Decomposition Parallel Preconditioners Package
|
||
8 years ago
|
! based on PSBLAS (Parallel Sparse BLAS version 3.5)
|
||
8 years ago
|
!
|
||
7 years ago
|
! (C) Copyright 2008-2018
|
||
8 years ago
|
!
|
||
7 years ago
|
! Salvatore Filippone
|
||
|
! Pasqua D'Ambra
|
||
|
! Daniela di Serafino
|
||
8 years ago
|
!
|
||
|
! 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 MLD2P4 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 MLD2P4 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.
|
||
|
!
|
||
|
!
|
||
14 years ago
|
!
|
||
4 years ago
|
! File: amg_s_diag_solver_mod.f90
|
||
14 years ago
|
!
|
||
4 years ago
|
! Module: amg_s_diag_solver_mod
|
||
14 years ago
|
!
|
||
8 years ago
|
! This module defines:
|
||
4 years ago
|
! - the amg_s_diag_solver_type data structure containing the
|
||
8 years ago
|
! simple diagonal solver. This extracts the main diagonal of a matrix
|
||
|
! and precomputes its inverse. Combined with a Jacobi "smoother" generates
|
||
|
! what are commonly known as the classic Jacobi iterations
|
||
14 years ago
|
!
|
||
4 years ago
|
module amg_s_diag_solver
|
||
14 years ago
|
|
||
4 years ago
|
use amg_s_base_solver_mod
|
||
14 years ago
|
|
||
4 years ago
|
type, extends(amg_s_base_solver_type) :: amg_s_diag_solver_type
|
||
13 years ago
|
type(psb_s_vect_type), allocatable :: dv
|
||
|
real(psb_spk_), allocatable :: d(:)
|
||
14 years ago
|
contains
|
||
4 years ago
|
procedure, pass(sv) :: dump => amg_s_diag_solver_dmp
|
||
|
procedure, pass(sv) :: build => amg_s_diag_solver_bld
|
||
|
procedure, pass(sv) :: cnv => amg_s_diag_solver_cnv
|
||
|
procedure, pass(sv) :: clone => amg_s_diag_solver_clone
|
||
|
procedure, pass(sv) :: clear_data => amg_s_diag_solver_clear_data
|
||
|
procedure, pass(sv) :: apply_v => amg_s_diag_solver_apply_vect
|
||
|
procedure, pass(sv) :: apply_a => amg_s_diag_solver_apply
|
||
13 years ago
|
procedure, pass(sv) :: free => s_diag_solver_free
|
||
|
procedure, pass(sv) :: descr => s_diag_solver_descr
|
||
|
procedure, pass(sv) :: sizeof => s_diag_solver_sizeof
|
||
|
procedure, pass(sv) :: get_nzeros => s_diag_solver_get_nzeros
|
||
12 years ago
|
procedure, nopass :: get_fmt => s_diag_solver_get_fmt
|
||
8 years ago
|
procedure, nopass :: get_id => s_diag_solver_get_id
|
||
4 years ago
|
end type amg_s_diag_solver_type
|
||
14 years ago
|
|
||
|
|
||
12 years ago
|
private :: s_diag_solver_free, s_diag_solver_descr, &
|
||
12 years ago
|
& s_diag_solver_sizeof, s_diag_solver_get_nzeros, &
|
||
8 years ago
|
& s_diag_solver_get_fmt, s_diag_solver_get_id
|
||
14 years ago
|
|
||
|
|
||
13 years ago
|
interface
|
||
4 years ago
|
subroutine amg_s_diag_solver_apply_vect(alpha,sv,x,beta,y,desc_data,&
|
||
7 years ago
|
& trans,work,wv,info,init,initu)
|
||
13 years ago
|
import :: psb_desc_type, psb_sspmat_type, psb_s_base_sparse_mat, &
|
||
12 years ago
|
& psb_s_vect_type, psb_s_base_vect_type, psb_spk_, &
|
||
4 years ago
|
& amg_s_diag_solver_type, psb_ipk_
|
||
12 years ago
|
type(psb_desc_type), intent(in) :: desc_data
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(inout) :: sv
|
||
13 years ago
|
type(psb_s_vect_type), intent(inout) :: x
|
||
|
type(psb_s_vect_type), intent(inout) :: y
|
||
12 years ago
|
real(psb_spk_),intent(in) :: alpha,beta
|
||
|
character(len=1),intent(in) :: trans
|
||
|
real(psb_spk_),target, intent(inout) :: work(:)
|
||
7 years ago
|
type(psb_s_vect_type),intent(inout) :: wv(:)
|
||
12 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
9 years ago
|
character, intent(in), optional :: init
|
||
|
type(psb_s_vect_type),intent(inout), optional :: initu
|
||
4 years ago
|
end subroutine amg_s_diag_solver_apply_vect
|
||
13 years ago
|
end interface
|
||
13 years ago
|
|
||
13 years ago
|
interface
|
||
4 years ago
|
subroutine amg_s_diag_solver_apply(alpha,sv,x,beta,y,desc_data,&
|
||
9 years ago
|
& trans,work,info,init,initu)
|
||
13 years ago
|
import :: psb_desc_type, psb_sspmat_type, psb_s_base_sparse_mat, &
|
||
12 years ago
|
& psb_s_vect_type, psb_s_base_vect_type, psb_spk_, &
|
||
4 years ago
|
& amg_s_diag_solver_type, psb_ipk_
|
||
12 years ago
|
type(psb_desc_type), intent(in) :: desc_data
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(inout) :: sv
|
||
13 years ago
|
real(psb_spk_), intent(inout) :: x(:)
|
||
|
real(psb_spk_), intent(inout) :: y(:)
|
||
|
real(psb_spk_),intent(in) :: alpha,beta
|
||
12 years ago
|
character(len=1),intent(in) :: trans
|
||
13 years ago
|
real(psb_spk_),target, intent(inout) :: work(:)
|
||
12 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
9 years ago
|
character, intent(in), optional :: init
|
||
|
real(psb_spk_),intent(inout), optional :: initu(:)
|
||
4 years ago
|
end subroutine amg_s_diag_solver_apply
|
||
13 years ago
|
end interface
|
||
13 years ago
|
|
||
13 years ago
|
interface
|
||
4 years ago
|
subroutine amg_s_diag_solver_bld(a,desc_a,sv,info,b,amold,vmold,imold)
|
||
13 years ago
|
import :: psb_desc_type, psb_sspmat_type, psb_s_base_sparse_mat, &
|
||
12 years ago
|
& psb_s_vect_type, psb_s_base_vect_type, psb_spk_, &
|
||
4 years ago
|
& amg_s_diag_solver_type, psb_ipk_, psb_i_base_vect_type
|
||
13 years ago
|
type(psb_sspmat_type), intent(in), target :: a
|
||
6 years ago
|
Type(psb_desc_type), Intent(inout) :: desc_a
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(inout) :: sv
|
||
12 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
13 years ago
|
type(psb_sspmat_type), intent(in), target, optional :: b
|
||
|
class(psb_s_base_sparse_mat), intent(in), optional :: amold
|
||
|
class(psb_s_base_vect_type), intent(in), optional :: vmold
|
||
11 years ago
|
class(psb_i_base_vect_type), intent(in), optional :: imold
|
||
4 years ago
|
end subroutine amg_s_diag_solver_bld
|
||
13 years ago
|
end interface
|
||
11 years ago
|
|
||
|
interface
|
||
4 years ago
|
subroutine amg_s_diag_solver_cnv(sv,info,amold,vmold,imold)
|
||
11 years ago
|
import :: psb_s_base_sparse_mat, psb_s_base_vect_type, psb_spk_, &
|
||
4 years ago
|
& amg_s_diag_solver_type, psb_ipk_, psb_i_base_vect_type
|
||
|
class(amg_s_diag_solver_type), intent(inout) :: sv
|
||
11 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
|
class(psb_s_base_sparse_mat), intent(in), optional :: amold
|
||
|
class(psb_s_base_vect_type), intent(in), optional :: vmold
|
||
|
class(psb_i_base_vect_type), intent(in), optional :: imold
|
||
4 years ago
|
end subroutine amg_s_diag_solver_cnv
|
||
11 years ago
|
end interface
|
||
9 years ago
|
|
||
|
interface
|
||
4 years ago
|
subroutine amg_s_diag_solver_dmp(sv,desc,level,info,prefix,head,solver,global_num)
|
||
|
import :: psb_desc_type, amg_s_diag_solver_type, psb_s_vect_type, psb_spk_, &
|
||
9 years ago
|
& psb_sspmat_type, psb_s_base_sparse_mat, psb_s_base_vect_type, &
|
||
|
& psb_ipk_
|
||
|
implicit none
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(in) :: sv
|
||
5 years ago
|
type(psb_desc_type), intent(in) :: desc
|
||
9 years ago
|
integer(psb_ipk_), intent(in) :: level
|
||
|
integer(psb_ipk_), intent(out) :: info
|
||
|
character(len=*), intent(in), optional :: prefix, head
|
||
5 years ago
|
logical, optional, intent(in) :: solver, global_num
|
||
4 years ago
|
end subroutine amg_s_diag_solver_dmp
|
||
9 years ago
|
end interface
|
||
|
|
||
12 years ago
|
interface
|
||
4 years ago
|
subroutine amg_s_diag_solver_clone(sv,svout,info)
|
||
12 years ago
|
import :: psb_desc_type, psb_sspmat_type, psb_s_base_sparse_mat, &
|
||
|
& psb_s_vect_type, psb_s_base_vect_type, psb_spk_, &
|
||
4 years ago
|
& amg_s_base_solver_type, amg_s_diag_solver_type, psb_ipk_
|
||
12 years ago
|
Implicit None
|
||
|
|
||
|
! Arguments
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(inout) :: sv
|
||
|
class(amg_s_base_solver_type), allocatable, intent(inout) :: svout
|
||
12 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
4 years ago
|
end subroutine amg_s_diag_solver_clone
|
||
12 years ago
|
end interface
|
||
5 years ago
|
|
||
|
interface
|
||
4 years ago
|
subroutine amg_s_diag_solver_clear_data(sv,info)
|
||
5 years ago
|
import :: psb_desc_type, psb_sspmat_type, psb_s_base_sparse_mat, &
|
||
|
& psb_s_vect_type, psb_s_base_vect_type, psb_spk_, &
|
||
4 years ago
|
& amg_s_diag_solver_type, psb_ipk_
|
||
5 years ago
|
Implicit None
|
||
|
|
||
|
! Arguments
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(inout) :: sv
|
||
5 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
4 years ago
|
end subroutine amg_s_diag_solver_clear_data
|
||
5 years ago
|
end interface
|
||
13 years ago
|
|
||
|
|
||
14 years ago
|
contains
|
||
|
|
||
|
subroutine s_diag_solver_free(sv,info)
|
||
|
|
||
|
Implicit None
|
||
|
|
||
|
! Arguments
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(inout) :: sv
|
||
12 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
|
integer(psb_ipk_) :: err_act
|
||
|
character(len=20) :: name='s_diag_solver_free'
|
||
14 years ago
|
|
||
|
call psb_erractionsave(err_act)
|
||
|
info = psb_success_
|
||
13 years ago
|
|
||
9 years ago
|
if (allocated(sv%dv)) call sv%dv%free(info)
|
||
|
|
||
14 years ago
|
if (allocated(sv%d)) then
|
||
|
deallocate(sv%d,stat=info)
|
||
|
if (info /= psb_success_) then
|
||
|
info = psb_err_alloc_dealloc_
|
||
|
call psb_errpush(info,name)
|
||
|
goto 9999
|
||
|
end if
|
||
|
end if
|
||
|
|
||
|
call psb_erractionrestore(err_act)
|
||
|
return
|
||
|
|
||
10 years ago
|
9999 call psb_error_handler(err_act)
|
||
14 years ago
|
return
|
||
10 years ago
|
|
||
14 years ago
|
end subroutine s_diag_solver_free
|
||
|
|
||
14 years ago
|
subroutine s_diag_solver_descr(sv,info,iout,coarse)
|
||
14 years ago
|
|
||
|
Implicit None
|
||
|
|
||
|
! Arguments
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(in) :: sv
|
||
12 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
|
integer(psb_ipk_), intent(in), optional :: iout
|
||
|
logical, intent(in), optional :: coarse
|
||
14 years ago
|
|
||
|
! Local variables
|
||
12 years ago
|
integer(psb_ipk_) :: err_act
|
||
4 years ago
|
character(len=20), parameter :: name='amg_s_diag_solver_descr'
|
||
12 years ago
|
integer(psb_ipk_) :: iout_
|
||
14 years ago
|
|
||
|
info = psb_success_
|
||
|
if (present(iout)) then
|
||
|
iout_ = iout
|
||
|
else
|
||
8 years ago
|
iout_ = psb_out_unit
|
||
14 years ago
|
endif
|
||
|
|
||
|
write(iout_,*) ' Diagonal local solver '
|
||
|
|
||
|
return
|
||
|
|
||
|
end subroutine s_diag_solver_descr
|
||
|
|
||
|
function s_diag_solver_sizeof(sv) result(val)
|
||
|
implicit none
|
||
|
! Arguments
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(in) :: sv
|
||
7 years ago
|
integer(psb_epk_) :: val
|
||
12 years ago
|
integer(psb_ipk_) :: i
|
||
14 years ago
|
|
||
|
val = 0
|
||
13 years ago
|
if (allocated(sv%dv)) val = val + sv%dv%sizeof()
|
||
14 years ago
|
|
||
|
return
|
||
|
end function s_diag_solver_sizeof
|
||
|
|
||
13 years ago
|
function s_diag_solver_get_nzeros(sv) result(val)
|
||
|
implicit none
|
||
|
! Arguments
|
||
4 years ago
|
class(amg_s_diag_solver_type), intent(in) :: sv
|
||
7 years ago
|
integer(psb_epk_) :: val
|
||
12 years ago
|
integer(psb_ipk_) :: i
|
||
13 years ago
|
|
||
|
val = 0
|
||
|
if (allocated(sv%dv)) val = val + sv%dv%get_nrows()
|
||
|
|
||
|
return
|
||
|
end function s_diag_solver_get_nzeros
|
||
|
|
||
12 years ago
|
function s_diag_solver_get_fmt() result(val)
|
||
|
implicit none
|
||
|
character(len=32) :: val
|
||
|
|
||
|
val = "Diag solver"
|
||
|
end function s_diag_solver_get_fmt
|
||
|
|
||
8 years ago
|
function s_diag_solver_get_id() result(val)
|
||
|
implicit none
|
||
|
integer(psb_ipk_) :: val
|
||
|
|
||
4 years ago
|
val = amg_diag_scale_
|
||
8 years ago
|
end function s_diag_solver_get_id
|
||
12 years ago
|
|
||
4 years ago
|
end module amg_s_diag_solver
|
||
6 years ago
|
|
||
5 years ago
|
!
|
||
4 years ago
|
! Module: amg_s_l1_diag_solver_mod
|
||
5 years ago
|
!
|
||
|
! This module defines:
|
||
4 years ago
|
! - the amg_s_l1_diag_solver_type data structure containing the
|
||
5 years ago
|
! L1 diagonal solver.
|
||
|
! The solver is defined as a diagonal containing in each element the
|
||
|
! inverse of the sum of the absolute values of the matrix entries
|
||
|
! along the corresponding row.
|
||
|
! Combined with a Jacobi "smoother" generates
|
||
|
! what are commonly known as the L1-Jacobi iterations
|
||
|
!
|
||
6 years ago
|
|
||
4 years ago
|
module amg_s_l1_diag_solver
|
||
6 years ago
|
|
||
4 years ago
|
use amg_s_diag_solver
|
||
6 years ago
|
|
||
4 years ago
|
type, extends(amg_s_diag_solver_type) :: amg_s_l1_diag_solver_type
|
||
6 years ago
|
contains
|
||
4 years ago
|
procedure, pass(sv) :: dump => amg_s_l1_diag_solver_dmp
|
||
|
procedure, pass(sv) :: build => amg_s_l1_diag_solver_bld
|
||
6 years ago
|
procedure, pass(sv) :: descr => s_l1_diag_solver_descr
|
||
|
procedure, nopass :: get_fmt => s_l1_diag_solver_get_fmt
|
||
|
procedure, nopass :: get_id => s_l1_diag_solver_get_id
|
||
4 years ago
|
end type amg_s_l1_diag_solver_type
|
||
6 years ago
|
|
||
|
|
||
|
private :: s_l1_diag_solver_descr, &
|
||
|
& s_l1_diag_solver_get_fmt, s_l1_diag_solver_get_id
|
||
|
|
||
|
interface
|
||
4 years ago
|
subroutine amg_s_l1_diag_solver_bld(a,desc_a,sv,info,b,amold,vmold,imold)
|
||
6 years ago
|
import :: psb_desc_type, psb_sspmat_type, psb_s_base_sparse_mat, &
|
||
|
& psb_s_vect_type, psb_s_base_vect_type, psb_spk_, &
|
||
4 years ago
|
& amg_s_l1_diag_solver_type, psb_ipk_, psb_i_base_vect_type
|
||
6 years ago
|
type(psb_sspmat_type), intent(in), target :: a
|
||
|
Type(psb_desc_type), Intent(inout) :: desc_a
|
||
4 years ago
|
class(amg_s_l1_diag_solver_type), intent(inout) :: sv
|
||
6 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
|
type(psb_sspmat_type), intent(in), target, optional :: b
|
||
|
class(psb_s_base_sparse_mat), intent(in), optional :: amold
|
||
|
class(psb_s_base_vect_type), intent(in), optional :: vmold
|
||
|
class(psb_i_base_vect_type), intent(in), optional :: imold
|
||
4 years ago
|
end subroutine amg_s_l1_diag_solver_bld
|
||
6 years ago
|
end interface
|
||
|
|
||
|
interface
|
||
4 years ago
|
subroutine amg_s_l1_diag_solver_dmp(sv,desc,level,info,prefix,head,solver,global_num)
|
||
|
import :: psb_desc_type, amg_s_l1_diag_solver_type, psb_s_vect_type, psb_spk_, &
|
||
6 years ago
|
& psb_sspmat_type, psb_s_base_sparse_mat, psb_s_base_vect_type, &
|
||
|
& psb_ipk_
|
||
|
implicit none
|
||
4 years ago
|
class(amg_s_l1_diag_solver_type), intent(in) :: sv
|
||
5 years ago
|
type(psb_desc_type), intent(in) :: desc
|
||
6 years ago
|
integer(psb_ipk_), intent(in) :: level
|
||
|
integer(psb_ipk_), intent(out) :: info
|
||
|
character(len=*), intent(in), optional :: prefix, head
|
||
5 years ago
|
logical, optional, intent(in) :: solver, global_num
|
||
4 years ago
|
end subroutine amg_s_l1_diag_solver_dmp
|
||
6 years ago
|
end interface
|
||
|
|
||
|
contains
|
||
|
|
||
|
subroutine s_l1_diag_solver_descr(sv,info,iout,coarse)
|
||
|
|
||
|
Implicit None
|
||
|
|
||
|
! Arguments
|
||
4 years ago
|
class(amg_s_l1_diag_solver_type), intent(in) :: sv
|
||
6 years ago
|
integer(psb_ipk_), intent(out) :: info
|
||
|
integer(psb_ipk_), intent(in), optional :: iout
|
||
|
logical, intent(in), optional :: coarse
|
||
|
|
||
|
! Local variables
|
||
|
integer(psb_ipk_) :: err_act
|
||
4 years ago
|
character(len=20), parameter :: name='amg_s_l1_diag_solver_descr'
|
||
6 years ago
|
integer(psb_ipk_) :: iout_
|
||
|
|
||
|
info = psb_success_
|
||
|
if (present(iout)) then
|
||
|
iout_ = iout
|
||
|
else
|
||
|
iout_ = psb_out_unit
|
||
|
endif
|
||
|
|
||
|
write(iout_,*) ' L1 Diagonal solver '
|
||
|
|
||
|
return
|
||
|
|
||
|
end subroutine s_l1_diag_solver_descr
|
||
|
|
||
|
function s_l1_diag_solver_get_fmt() result(val)
|
||
|
implicit none
|
||
|
character(len=32) :: val
|
||
|
|
||
|
val = "L1 Diag solver"
|
||
|
end function s_l1_diag_solver_get_fmt
|
||
|
|
||
|
function s_l1_diag_solver_get_id() result(val)
|
||
|
implicit none
|
||
|
integer(psb_ipk_) :: val
|
||
|
|
||
4 years ago
|
val = amg_l1_diag_scale_
|
||
6 years ago
|
end function s_l1_diag_solver_get_id
|
||
|
|
||
4 years ago
|
end module amg_s_l1_diag_solver
|
||
6 years ago
|
|