@@ -153,6 +153,11 @@ module scalapackfx_module
153
153
module procedure scalafx_creatematrix_complex, scalafx_creatematrix_dcomplex
154
154
end interface scalafx_creatematrix
155
155
156
+ !> Maps global position in a distributed matrix to local one.
157
+ interface scalafx_infog2l
158
+ module procedure scalafx_infog2l_single, scalafx_infog2l_array
159
+ end interface scalafx_infog2l
160
+
156
161
!************************************************************************
157
162
!*** ppotrf
158
163
!************************************************************************
@@ -1965,7 +1970,8 @@ contains
1965
1970
!! \param rsrc Row of the process owning the local matrix.
1966
1971
!! \param csrc Column of the process owning the local matrix.
1967
1972
!!
1968
- subroutine scalafx_infog2l (mygrid , desc , grow , gcol , lrow , lcol , rsrc , csrc )
1973
+ subroutine scalafx_infog2l_single (mygrid , desc , grow , gcol ,&
1974
+ & lrow , lcol , rsrc , csrc )
1969
1975
type(blacsgrid), intent (in ) :: mygrid
1970
1976
integer , intent (in ) :: desc(DLEN_)
1971
1977
integer , intent (in ) :: grow, gcol
@@ -1975,8 +1981,92 @@ contains
1975
1981
call infog2l(grow, gcol, desc, mygrid%nrow, mygrid%ncol, mygrid%myrow,&
1976
1982
& mygrid%mycol, lrow, lcol, rsrc, csrc)
1977
1983
1978
- end subroutine scalafx_infog2l
1984
+ end subroutine scalafx_infog2l_single
1985
+
1986
+ !> Maps global positions in a distributed matrix to local one.
1987
+ !!
1988
+ !! \param mygrid BLACS descriptor.
1989
+ !! \param desc Descriptor of the distributed matrix.
1990
+ !! \param grow Global row indices.
1991
+ !! \param gcol Global column indices.
1992
+ !! \param lrow Local row indices on output.
1993
+ !! \param lcol Local column indices on output.
1994
+ !! \param rsrc Rows of the process owning the local matrix.
1995
+ !! \param csrc Columns of the process owning the local matrix.
1996
+ !! \param calcAllIndices Whether to calculate all lrow and lcol,
1997
+ !! even if the current process does not own them. (default: true)
1998
+ !!
1999
+ subroutine scalafx_infog2l_array (mygrid , desc , grow , gcol ,&
2000
+ & lrow , lcol , rsrc , csrc , calcAllIndices )
2001
+ type(blacsgrid), intent (in ) :: mygrid
2002
+ integer , intent (in ) :: desc(DLEN_)
2003
+ integer , intent (in ) :: grow(:), gcol(:)
2004
+ integer , intent (out ) :: lrow(:), rsrc(:)
2005
+ integer , intent (out ) :: lcol(:), csrc(:)
2006
+ logical , intent (in ), optional :: calcAllIndices
2007
+
2008
+ call scalapackfx_infog2l_helper(grow, desc(MB_), desc(RSRC_),&
2009
+ & mygrid%myrow, mygrid%nrow, lrow, rsrc, calcAllIndices)
2010
+
2011
+ call scalapackfx_infog2l_helper(gcol, desc(NB_), desc(CSRC_),&
2012
+ & mygrid%mycol, mygrid%ncol, lcol, csrc, calcAllIndices)
2013
+
2014
+ end subroutine scalafx_infog2l_array
2015
+
2016
+ !> Helper routine for scalafx_infog2l_array.
2017
+ !!
2018
+ !! \param globalInd Global row/ column indices.
2019
+ !! \param descB Either desc(MB_) or desc(NB_).
2020
+ !! \param descSRC Either desc(RSRC_) or desc(CSRC_).
2021
+ !! \param myPos Row/ column of the current process.
2022
+ !! \param nPos Number of rows/ columns.
2023
+ !! \param localInd Local row/ column indices on output.
2024
+ !! \param localPos Rows/ columns of the process owning the local matrix.
2025
+ !! \param calcAllIndices Whether to calculate all local indices,
2026
+ !! even if the current process does not own them. (default: true)
2027
+ !!
2028
+ subroutine scalapackfx_infog2l_helper (globalInd , descB , descSRC ,&
2029
+ & myPos , nPos , localInd , localPos , calcAllIndices )
2030
+ integer , intent (in ) :: globalInd(:)
2031
+ integer , intent (in ) :: descB, descSRC
2032
+ integer , intent (in ) :: myPos, nPos
2033
+ integer , intent (out ) :: localInd(:), localPos(:)
2034
+ logical , intent (in ), optional :: calcAllIndices
2035
+
2036
+ real (dp) :: inv
2037
+ integer , dimension (size (globalInd)) :: blk
2038
+ integer :: check, i
2039
+ logical :: calcAllIndices_
2040
+
2041
+ ! Note that we explicitly multiply with a double here instead of
2042
+ ! dividing by an integer to enhance performance.
2043
+ inv = 1.0_dp / real (descB, kind= dp)
2044
+ blk = (globalInd - 1 ) * inv
2045
+
2046
+ check = modulo(myPos - descSRC, nPos)
2047
+
2048
+ localPos = mod (blk + descSRC, nPos)
2049
+
2050
+ calcAllIndices_ = .true.
2051
+ if (present (calcAllIndices)) then
2052
+ calcAllIndices_ = calcAllIndices
2053
+ end if
2054
+
2055
+ do i = 1 , size (globalInd)
2056
+ if (calcAllIndices_ .or. myPos == localPos(i)) then
2057
+ localInd(i) = (blk(i) / nPos + 1 ) * descB + 1
2058
+ if (check >= mod (blk(i), nPos)) then
2059
+ if (myPos == localPos(i)) then
2060
+ localInd(i) = localInd(i) + mod (globalInd(i) - 1 , descB)
2061
+ end if
2062
+ localInd(i) = localInd(i) - descB
2063
+ end if
2064
+ else
2065
+ localInd(i) = - 1
2066
+ end if
2067
+ end do
1979
2068
2069
+ end subroutine scalapackfx_infog2l_helper
1980
2070
1981
2071
!> Maps local row or column index onto global matrix position.
1982
2072
!!
0 commit comments