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.
psblas3/base/internals/psi_compute_size.f90

138 lines
4.3 KiB
Fortran

!!$
!!$ Parallel Sparse BLAS version 2.2
!!$ (C) Copyright 2006/2007/2008
!!$ Salvatore Filippone University of Rome Tor Vergata
!!$ Alfredo Buttari University of Rome Tor Vergata
!!$
!!$ 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 PSBLAS 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 PSBLAS 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.
!!$
!!$
!
! Compute maximum data exchange size; small utility for assembly of descriptors.
!
subroutine psi_compute_size(desc_data, index_in, dl_lda, info)
20 years ago
use psi_mod, psb_protect_name => psi_compute_size
20 years ago
use psb_const_mod
use psb_descriptor_type
20 years ago
use psb_error_mod
use psb_penv_mod
20 years ago
implicit none
! ....scalars parameters....
integer :: info, dl_lda
! .....array parameters....
integer :: desc_data(:), index_in(:)
! ....local scalars....
integer :: i,np,me,proc, max_index
integer :: ictxt, err_act
20 years ago
! ...local array...
integer :: int_err(5)
integer, allocatable :: counter_recv(:), counter_dl(:)
20 years ago
! ...parameters
integer :: debug_level, debug_unit
20 years ago
character(len=20) :: name
name='psi_compute_size'
call psb_get_erraction(err_act)
debug_unit = psb_get_debug_unit()
debug_level = psb_get_debug_level()
20 years ago
info = psb_success_
ictxt = desc_data(psb_ctxt_)
20 years ago
call psb_info(ictxt,me,np)
if (np == -1) then
psblas3: base/comm/psb_cgather.f90 base/comm/psb_chalo.f90 base/comm/psb_covrl.f90 base/comm/psb_cscatter.F90 base/comm/psb_dgather.f90 base/comm/psb_dhalo.f90 base/comm/psb_dovrl.f90 base/comm/psb_dscatter.F90 base/comm/psb_igather.f90 base/comm/psb_ihalo.f90 base/comm/psb_iovrl.f90 base/comm/psb_iscatter.F90 base/comm/psb_sgather.f90 base/comm/psb_shalo.f90 base/comm/psb_sovrl.f90 base/comm/psb_sscatter.F90 base/comm/psb_zgather.f90 base/comm/psb_zhalo.f90 base/comm/psb_zovrl.f90 base/comm/psb_zscatter.F90 base/internals/psi_bld_g2lmap.f90 base/internals/psi_bld_tmphalo.f90 base/internals/psi_bld_tmpovrl.f90 base/internals/psi_compute_size.f90 base/internals/psi_crea_index.f90 base/internals/psi_cswapdata.F90 base/internals/psi_cswaptran.F90 base/internals/psi_desc_index.F90 base/internals/psi_dswapdata.F90 base/internals/psi_dswaptran.F90 base/internals/psi_fnd_owner.F90 base/internals/psi_iswapdata.F90 base/internals/psi_iswaptran.F90 base/internals/psi_ldsc_pre_halo.f90 base/internals/psi_sswapdata.F90 base/internals/psi_sswaptran.F90 base/internals/psi_zswapdata.F90 base/internals/psi_zswaptran.F90 base/modules/psb_const_mod.F90 base/modules/psb_desc_type.f90 base/modules/psb_error_mod.F90 base/psblas/psb_camax.f90 base/psblas/psb_casum.f90 base/psblas/psb_caxpby.f90 base/psblas/psb_cdot.f90 base/psblas/psb_cnrm2.f90 base/psblas/psb_cnrmi.f90 base/psblas/psb_cspmm.f90 base/psblas/psb_cspsm.f90 base/psblas/psb_damax.f90 base/psblas/psb_dasum.f90 base/psblas/psb_daxpby.f90 base/psblas/psb_ddot.f90 base/psblas/psb_dnrm2.f90 base/psblas/psb_dnrmi.f90 base/psblas/psb_dspmm.f90 base/psblas/psb_dspsm.f90 base/psblas/psb_samax.f90 base/psblas/psb_sasum.f90 base/psblas/psb_saxpby.f90 base/psblas/psb_sdot.f90 base/psblas/psb_snrm2.f90 base/psblas/psb_snrmi.f90 base/psblas/psb_sspmm.f90 base/psblas/psb_sspsm.f90 base/psblas/psb_sxdot.f90 base/psblas/psb_zamax.f90 base/psblas/psb_zasum.f90 base/psblas/psb_zaxpby.f90 base/psblas/psb_zdot.f90 base/psblas/psb_znrm2.f90 base/psblas/psb_znrmi.f90 base/psblas/psb_zspmm.f90 base/psblas/psb_zspsm.f90 base/serial/psi_impl.f90 base/tools/psb_callc.f90 base/tools/psb_casb.f90 base/tools/psb_cdcpy.f90 base/tools/psb_cdren.f90 base/tools/psb_cfree.f90 base/tools/psb_cins.f90 base/tools/psb_cspalloc.f90 base/tools/psb_cspasb.f90 base/tools/psb_dallc.f90 base/tools/psb_dasb.f90 base/tools/psb_dfree.f90 base/tools/psb_dins.f90 base/tools/psb_dspalloc.f90 base/tools/psb_dspasb.f90 base/tools/psb_ialloc.f90 base/tools/psb_iasb.f90 base/tools/psb_icdasb.F90 base/tools/psb_ifree.f90 base/tools/psb_iins.f90 base/tools/psb_sallc.f90 base/tools/psb_sasb.f90 base/tools/psb_sfree.f90 base/tools/psb_sins.f90 base/tools/psb_sspalloc.f90 base/tools/psb_sspasb.f90 base/tools/psb_zallc.f90 base/tools/psb_zasb.f90 base/tools/psb_zfree.f90 base/tools/psb_zins.f90 base/tools/psb_zspalloc.f90 base/tools/psb_zspasb.f90 test/serial/Makefile test/serial/d_coo_matgen.f03 Fix blacs_error into context_error
15 years ago
info = psb_err_context_error_
call psb_errpush(info,name)
goto 9999
20 years ago
endif
allocate(counter_dl(0:np-1),counter_recv(0:np-1),stat=info)
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_,name,a_err='Allocate')
goto 9999
end if
20 years ago
! ..initialize counters...
do i=0,np-1
counter_recv(i)=0
counter_dl(i)=0
20 years ago
enddo
! ....verify local correctness of halo_in....
i=1
do while (index_in(i) /= -1)
proc=index_in(i)
if ((proc > np-1).or.(proc < 0)) then
info = psb_err_invalid_pid_arg_
int_err(1) = 11
int_err(2) = proc
call psb_errpush(info,name,i_err=int_err)
goto 9999
endif
counter_dl(proc)=1
! ..update no of elements to receive from proc proc..
counter_recv(proc)=counter_recv(proc)+&
& index_in(i+1)
i=i+index_in(i+1)+2
20 years ago
enddo
! ...computing max_halo: max halo points to be received from
! same processor
max_index=0
dl_lda=0
do i=0,np-1
if (counter_recv(i) > max_index) max_index = counter_recv(i)
if (counter_dl(i) == 1) dl_lda = dl_lda+1
20 years ago
enddo
! computing max global value of dl_lda
call psb_amx(ictxt, dl_lda)
20 years ago
if (debug_level>=psb_debug_inner_) then
write(debug_unit,*) me,' ',trim(name),': ',dl_lda
20 years ago
endif
call psb_erractionrestore(err_act)
return
9999 continue
call psb_erractionrestore(err_act)
if (err_act == psb_act_abort_) then
call psb_error(ictxt)
return
20 years ago
end if
return
end subroutine psi_compute_size