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_dmlprec_bld.f90

183 lines
7.0 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_dmlprec_bld.f90
!
! Subroutine: mld_dmlprec_bld
! Version: real
!
! This routine builds the base preconditioner corresponding to the current
! level of the multilevel preconditioner. The routine first builds the
! (coarse) matrix associated to the current level from the (fine) matrix
! associated to the previous level, then builds the related base preconditioner.
!
!
! Arguments:
! a - type(psb_dspmat_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_dbaseprc_type), input/output.
! The base preconditioner data structure containing the local
! part of the base preconditioner to be built.
! info - integer, output.
! Error code.
!
subroutine mld_dmlprec_bld(a,desc_a,p,info)
18 years ago
use psb_base_mod
use mld_inner_mod, mld_protect_name => mld_dmlprec_bld
18 years ago
implicit none
! Arguments
18 years ago
type(psb_dspmat_type), intent(in), target :: a
type(psb_desc_type), intent(in), target :: desc_a
type(mld_dbaseprc_type), intent(inout),target :: p
18 years ago
integer, intent(out) :: info
! Local variables
type(psb_desc_type) :: desc_ac
type(psb_dspmat_type) :: ac
character(len=20) :: name
integer :: ictxt, np, me, err_act
18 years ago
mld2p4: config/pac.m4 configure krylov/Makefile krylov/psb_prec_mod.F90 mlprec/Makefile mlprec/mld_basep_bld_mod.f90 mlprec/mld_caggrmap_bld.f90 mlprec/mld_caggrmat_asb.f90 mlprec/mld_caggrmat_raw_asb.F90 mlprec/mld_caggrmat_smth_asb.F90 mlprec/mld_cas_aply.f90 mlprec/mld_cas_bld.f90 mlprec/mld_cbaseprec_aply.f90 mlprec/mld_cbaseprec_bld.f90 mlprec/mld_cdiag_bld.f90 mlprec/mld_cfact_bld.f90 mlprec/mld_cilu0_fact.f90 mlprec/mld_cilu_bld.f90 mlprec/mld_ciluk_fact.f90 mlprec/mld_cilut_fact.f90 mlprec/mld_cmlprec_aply.f90 mlprec/mld_cmlprec_bld.f90 mlprec/mld_cprec_aply.f90 mlprec/mld_cprecbld.f90 mlprec/mld_cprecfree.f90 mlprec/mld_cprecinit.f90 mlprec/mld_cprecset.f90 mlprec/mld_cslu_bld.f90 mlprec/mld_cslu_interface.c mlprec/mld_cslud_bld.f90 mlprec/mld_cslud_interface.c mlprec/mld_csp_renum.f90 mlprec/mld_csub_aply.f90 mlprec/mld_csub_solve.f90 mlprec/mld_cumf_bld.f90 mlprec/mld_cumf_interface.c mlprec/mld_dmlprec_bld.f90 mlprec/mld_dprecset.f90 mlprec/mld_inner_mod.f90 mlprec/mld_prec_mod.f90 mlprec/mld_prec_type.f90 mlprec/mld_saggrmap_bld.f90 mlprec/mld_saggrmat_asb.f90 mlprec/mld_saggrmat_raw_asb.F90 mlprec/mld_saggrmat_smth_asb.F90 mlprec/mld_sas_aply.f90 mlprec/mld_sas_bld.f90 mlprec/mld_sbaseprec_aply.f90 mlprec/mld_sbaseprec_bld.f90 mlprec/mld_sdiag_bld.f90 mlprec/mld_sfact_bld.f90 mlprec/mld_silu0_fact.f90 mlprec/mld_silu_bld.f90 mlprec/mld_siluk_fact.f90 mlprec/mld_silut_fact.f90 mlprec/mld_smlprec_aply.f90 mlprec/mld_smlprec_bld.f90 mlprec/mld_sprec_aply.f90 mlprec/mld_sprecbld.f90 mlprec/mld_sprecfree.f90 mlprec/mld_sprecinit.f90 mlprec/mld_sprecset.f90 mlprec/mld_sslu_bld.f90 mlprec/mld_sslu_interface.c mlprec/mld_sslud_bld.f90 mlprec/mld_sslud_interface.c mlprec/mld_ssp_renum.f90 mlprec/mld_ssub_aply.f90 mlprec/mld_ssub_solve.f90 mlprec/mld_sumf_bld.f90 mlprec/mld_sumf_interface.c mlprec/mld_zilu0_fact.f90 mlprec/mld_zmlprec_aply.f90 mlprec/mld_zmlprec_bld.f90 mlprec/mld_zprecset.f90 test/fileread/Makefile test/fileread/cf_sample.f90 test/fileread/data_input.f90 test/fileread/df_sample.f90 test/fileread/runs/cfs.inp test/fileread/runs/dfs.inp test/fileread/runs/sfs.inp test/fileread/runs/zfs.inp test/fileread/sf_sample.f90 test/fileread/zf_sample.f90 test/pargen/Makefile test/pargen/runs/ppde.inp test/pargen/spde.f90 Merged single precision version.
17 years ago
name='mld_dmlprec_bld'
18 years ago
if (psb_get_errstatus().ne.0) return
call psb_erractionsave(err_act)
18 years ago
info = 0
ictxt = psb_cd_get_context(desc_a)
call psb_info(ictxt,me,np)
if (.not.allocated(p%iprcparm)) then
info = 2222
call psb_errpush(info,name)
goto 9999
endif
call mld_check_def(p%iprcparm(mld_ml_type_),'Multilevel type',&
& mld_mult_ml_,is_legal_ml_type)
call mld_check_def(p%iprcparm(mld_aggr_alg_),'Aggregation',&
& mld_dec_aggr_,is_legal_ml_aggr_alg)
call mld_check_def(p%iprcparm(mld_aggr_kind_),'Smoother',&
& mld_smooth_prol_,is_legal_ml_aggr_kind)
call mld_check_def(p%iprcparm(mld_coarse_mat_),'Coarse matrix',&
& mld_distr_mat_,is_legal_ml_coarse_mat)
17 years ago
call mld_check_def(p%iprcparm(mld_smoother_pos_),'smooth_pos',&
& mld_pre_smooth_,is_legal_ml_smooth_pos)
18 years ago
select case(p%iprcparm(mld_sub_solve_))
case(mld_ilu_n_,mld_milu_n_)
17 years ago
call mld_check_def(p%iprcparm(mld_sub_fillin_),'Level',0,is_legal_ml_lev)
case(mld_ilu_t_)
call mld_check_def(p%rprcparm(mld_fact_thrs_),'Eps',dzero,is_legal_fact_thrs)
18 years ago
end select
call mld_check_def(p%rprcparm(mld_aggr_damp_),'Omega',dzero,is_legal_omega)
call mld_check_def(p%rprcparm(mld_aggr_thresh_),'Aggr_Thresh',dzero,is_legal_aggr_thrs)
call mld_check_def(p%iprcparm(mld_smoother_sweeps_),'Jacobi sweeps',&
18 years ago
& 1,is_legal_jac_sweeps)
!
! Build a mapping between the row indices of the fine-level matrix
! and the row indices of the coarse-level matrix, according to a decoupled
! aggregation algorithm. This also defines a tentative prolongator from
! the coarse to the fine level.
!
call mld_aggrmap_bld(p%iprcparm(mld_aggr_alg_),p%rprcparm(mld_aggr_thresh_),&
& a,desc_a,p%nlaggr,p%mlia,info)
18 years ago
if(info /= 0) then
call psb_errpush(4010,name,a_err='mld_aggrmap_bld')
18 years ago
goto 9999
end if
!
! Build the coarse-level matrix from the fine level one, starting from
! the mapping defined by mld_aggrmap_bld and applying the aggregation
! algorithm specified by p%iprcparm(mld_aggr_kind_)
!
call mld_aggrmat_asb(a,desc_a,ac,desc_ac,p,info)
18 years ago
if(info /= 0) then
call psb_errpush(4010,name,a_err='mld_aggrmat_asb')
18 years ago
goto 9999
end if
!
! Build the 'base preconditioner' corresponding to the coarse level
!
call mld_baseprc_bld(ac,desc_ac,p,info)
if (info /= 0) then
call psb_errpush(4010,name,a_err='mld_baseprc_bld')
18 years ago
goto 9999
end if
!
! We have used a separate ac because
! 1. we want to reuse the same routines mld_ilu_bld, etc.,
! 2. we do NOT want to pass an argument twice to them (p%av(mld_ac_) and p),
! as this would violate the Fortran standard.
18 years ago
! Hence a separate AC and a TRANSFER function at the end.
!
call psb_sp_transfer(ac,p%av(mld_ac_),info)
p%base_a => p%av(mld_ac_)
if (info==0) call psb_cdtransfer(desc_ac,p%desc_ac,info)
18 years ago
p%map_desc = psb_inter_desc(psb_map_aggr_,desc_a,&
& p%desc_ac,p%av(mld_sm_pr_t_),p%av(mld_sm_pr_))
! The two matrices from p%av() have been copied, may free them.
if (info == 0) call psb_sp_free(p%av(mld_sm_pr_t_),info)
if (info == 0) call psb_sp_free(p%av(mld_sm_pr_),info)
18 years ago
if (info /= 0) then
call psb_errpush(4010,name,a_err='psb_cdtransfer')
18 years ago
goto 9999
end if
p%base_desc => p%desc_ac
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
end subroutine mld_dmlprec_bld