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/mlprec/mld_zprecbld.f90

235 lines
7.4 KiB
Fortran

18 years ago
!!$
!!$
!!$ MLD2P4 version 1.0
!!$ MultiLevel Domain Decomposition Parallel Preconditioners Package
!!$ based on PSBLAS (Parallel Sparse BLAS version 2.2)
!!$
!!$ (C) Copyright 2008
!!$
!!$ Salvatore Filippone University of Rome Tor Vergata
!!$ Alfredo Buttari University of Rome Tor Vergata
!!$ Pasqua D'Ambra ICAR-CNR, Naples
!!$ Daniela di Serafino Second University of Naples
18 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
18 years ago
!!$ 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
18 years ago
!!$ 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: mld_zprecbld.f90
!
! Subroutine: mld_zprecbld
! Version: complex
! Contains: subroutine init_baseprc_av
!
! This routine builds the preconditioner according to the requirements made by
! the user trough the subroutines mld_precinit and mld_precset.
!
! A multilevel preconditioner is regarded as an array of 'base preconditioners',
! each representing the part of the preconditioner associated to a certain level.
! The levels are numbered in increasing order starting from the finest one, i.e.
! level 1 is the finest level.
!
!
! Arguments:
! a - type(psb_zspmat_type).
! The sparse matrix structure containing the local part of the
! matrix to be preconditioned.
! desc_a - type(psb_desc_type), input.
! The communication descriptor of a.
! p - type(mld_zprec_type), input/output.
! The preconditioner data structure containing the local part
! of the preconditioner to be built.
! info - integer, output.
! Error code.
!
subroutine mld_zprecbld(a,desc_a,p,info)
18 years ago
use psb_base_mod
use mld_inner_mod
use mld_prec_mod, mld_protect_name => mld_zprecbld
18 years ago
Implicit None
! Arguments
18 years ago
type(psb_zspmat_type), target :: a
type(psb_desc_type), intent(in), target :: desc_a
type(mld_zprec_type),intent(inout) :: p
18 years ago
integer, intent(out) :: info
!!$ character, intent(in), optional :: upd
18 years ago
! Local Variables
Integer :: err,i,k,ictxt, me,np, err_act, iszv
18 years ago
integer :: int_err(5)
character :: upd_
integer :: debug_level, debug_unit
character(len=20) :: name, ch_err
18 years ago
if(psb_get_errstatus().ne.0) return
info=0
err=0
call psb_erractionsave(err_act)
debug_unit = psb_get_debug_unit()
debug_level = psb_get_debug_level()
18 years ago
name = 'mld_zprecbld'
18 years ago
info = 0
int_err(1) = 0
ictxt = psb_cd_get_context(desc_a)
call psb_info(ictxt, me, np)
if (debug_level >= psb_debug_outer_) &
& write(debug_unit,*) me,' ',trim(name),&
& 'Entering ',desc_a%matrix_data(:)
!
! For the time being we are commenting out the UPDATE argument
! we plan to resurrect it later.
!!$ if (present(upd)) then
!!$ if (debug_level >= psb_debug_outer_) &
!!$ & write(debug_unit,*) me,' ',trim(name),'UPD ', upd
!!$
!!$ if ((toupper(upd).eq.'F').or.(toupper(upd).eq.'T')) then
!!$ upd_=toupper(upd)
!!$ else
!!$ upd_='F'
!!$ endif
!!$ else
!!$ upd_='F'
!!$ endif
upd_ = 'F'
18 years ago
if (.not.allocated(p%baseprecv)) then
!! Error: should have called mld_dprecinit
info=3111
call psb_errpush(info,name)
goto 9999
18 years ago
end if
18 years ago
!
! Check to ensure all procs have the same
!
iszv = size(p%baseprecv)
call psb_bcast(ictxt,iszv)
if (iszv /= size(p%baseprecv)) then
info=4001
call psb_errpush(info,name,a_err='Inconsistent size of baseprecv')
goto 9999
end if
if (iszv >= 1) then
!
! Allocate and build the fine level preconditioner
!
18 years ago
call init_baseprc_av(p%baseprecv(1),info)
if (info == 0) call mld_baseprc_bld(a,desc_a,p%baseprecv(1),info,upd_)
18 years ago
if (info /= 0) then
call psb_errpush(4001,name,a_err='Base level precbuild.')
18 years ago
goto 9999
end if
18 years ago
else
info=4010
ch_err='size bpv'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
18 years ago
endif
if (iszv > 1) then
18 years ago
!
! Build the base preconditioners corresponding to the remaining
! levels
!
do i=2, iszv
!
! Allocate the av component of the preconditioner data type
! at level i
!
if (i<iszv) then
!
! A replicated matrix only makes sense at the coarsest level
!
call mld_check_def(p%baseprecv(i)%iprcparm(mld_coarse_mat_),'Coarse matrix',&
& mld_distr_mat_,is_distr_ml_coarse_mat)
end if
call init_baseprc_av(p%baseprecv(i),info)
if (debug_level >= psb_debug_outer_) &
& write(debug_unit,*) me,' ',trim(name),&
& 'Calling mlprcbld at level ',i
!
! Build the base preconditioner corresponding to level i
!
if (info == 0) call mld_mlprec_bld(p%baseprecv(i-1)%base_a,&
& p%baseprecv(i-1)%base_desc, p%baseprecv(i),info)
18 years ago
if (info /= 0) then
call psb_errpush(4001,name,a_err='Init & build upper level preconditioner')
18 years ago
goto 9999
endif
if (debug_level >= psb_debug_outer_) &
& write(debug_unit,*) me,' ',trim(name),&
& 'Return from ',i,' call to mlprcbld ',info
18 years ago
end do
endif
call psb_erractionrestore(err_act)
return
9999 continue
call psb_erractionrestore(err_act)
if (err_act.eq.psb_act_abort_) then
18 years ago
call psb_error()
return
end if
return
contains
subroutine init_baseprc_av(p,info)
type(mld_zbaseprc_type), intent(inout) :: p
18 years ago
integer :: info
if (allocated(p%av)) then
if (size(p%av) /= mld_max_avsz_) then
deallocate(p%av,stat=info)
if (info /= 0) return
endif
end if
if (.not.(allocated(p%av))) then
allocate(p%av(mld_max_avsz_),stat=info)
if (info /= 0) return
18 years ago
end if
do k=1,size(p%av)
call psb_nullify_sp(p%av(k))
end do
18 years ago
end subroutine init_baseprc_av
end subroutine mld_zprecbld
18 years ago