Skip to content
GitLab
Projects
Groups
Snippets
Help
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
fleur
Project overview
Project overview
Details
Activity
Releases
Repository
Repository
Files
Commits
Branches
Tags
Contributors
Graph
Compare
Issues
52
Issues
52
List
Boards
Labels
Service Desk
Milestones
Operations
Operations
Incidents
Packages & Registries
Packages & Registries
Container Registry
Analytics
Analytics
Repository
Value Stream
Wiki
Wiki
Members
Members
Collapse sidebar
Close sidebar
Activity
Graph
Create a new issue
Commits
Issue Boards
Open sidebar
fleur
fleur
Commits
7af771fa
Commit
7af771fa
authored
Oct 04, 2018
by
Daniel Wortmann
Browse files
Options
Browse Files
Download
Email Patches
Plain Diff
Unifying use of blacs processor grid
parent
34d390b2
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
with
40 additions
and
46 deletions
+40
-46
diagonalization/elpa.F90
diagonalization/elpa.F90
+12
-14
types/types_mpimat.F90
types/types_mpimat.F90
+28
-32
No files found.
diagonalization/elpa.F90
View file @
7af771fa
...
...
@@ -47,7 +47,7 @@ CONTAINS
INCLUDE
'mpif.h'
!... Local variables
!
INTEGER
::
num
,
num2
,
myrow
,
mycol
INTEGER
::
num
,
num2
INTEGER
::
nb
,
myid
,
np
INTEGER
::
n_col
,
n_row
LOGICAL
::
ok
...
...
@@ -77,18 +77,16 @@ CONTAINS
CALL
MPI_BARRIER
(
hmat
%
blacsdata
%
mpi_com
,
err
)
CALL
MPI_COMM_RANK
(
hmat
%
blacsdata
%
mpi_com
,
myid
,
err
)
CALL
MPI_COMM_SIZE
(
hmat
%
blacsdata
%
mpi_com
,
np
,
err
)
myrow
=
myid
/
hmat
%
blacsdata
%
npcol
mycol
=
myid
-
(
myid
/
hmat
%
blacsdata
%
npcol
)
*
hmat
%
blacsdata
%
npcol
!Create communicators for ELPA
#if defined (CPP_ELPA_201705003)
mpi_comm_rows
=
-1
mpi_comm_cols
=
-1
#elif defined (CPP_ELPA_201605004) || defined (CPP_ELPA_201605003)||defined(CPP_ELPA_NEW)
err
=
get_elpa_row_col_comms
(
hmat
%
blacsdata
%
mpi_com
,
myrow
,
mycol
,
mpi_comm_rows
,
mpi_comm_cols
)
err
=
get_elpa_row_col_comms
(
hmat
%
blacsdata
%
mpi_com
,
hmat
%
blacsdata
%
myrow
,
hmat
%
blacsdata
%
mycol
,
mpi_comm_rows
,
mpi_comm_cols
)
#else
CALL
get_elpa_row_col_comms
(
hmat
%
blacsdata
%
mpi_com
,
myrow
,
mycol
,
mpi_comm_rows
,
mpi_comm_cols
)
CALL
get_elpa_row_col_comms
(
hmat
%
blacsdata
%
mpi_com
,
hmat
%
blacsdata
%
myrow
,
hmat
%
blacsdata
%
mycol
,
mpi_comm_rows
,
mpi_comm_cols
)
#endif
!print *,"creating ELPA comms -- done"
...
...
@@ -119,8 +117,8 @@ CONTAINS
CALL
elpa_obj
%
set
(
"local_ncols"
,
hmat
%
matsize2
,
err
)
CALL
elpa_obj
%
set
(
"nblk"
,
nb
,
err
)
CALL
elpa_obj
%
set
(
"mpi_comm_parent"
,
hmat
%
blacsdata
%
mpi_com
,
err
)
CALL
elpa_obj
%
set
(
"process_row"
,
myrow
,
err
)
CALL
elpa_obj
%
set
(
"process_col"
,
mycol
,
err
)
CALL
elpa_obj
%
set
(
"process_row"
,
hmat
%
blacsdata
%
myrow
,
err
)
CALL
elpa_obj
%
set
(
"process_col"
,
hmat
%
blacsdata
%
mycol
,
err
)
#ifdef CPP_ELPA2
CALL
elpa_obj
%
set
(
"solver"
,
ELPA_SOLVER_2STAGE
)
#else
...
...
@@ -188,8 +186,8 @@ CONTAINS
DO
i
=
1
,
hmat
%
matsize2
! Get global column corresponding to i and number of local rows up to
! and including the diagonal, these are unchanged in H
n_col
=
indxl2g
(
i
,
nb
,
mycol
,
0
,
hmat
%
blacsdata
%
npcol
)
n_row
=
numroc
(
n_col
,
nb
,
myrow
,
0
,
hmat
%
blacsdata
%
nprow
)
n_col
=
indxl2g
(
i
,
nb
,
hmat
%
blacsdata
%
mycol
,
0
,
hmat
%
blacsdata
%
npcol
)
n_row
=
numroc
(
n_col
,
nb
,
hmat
%
blacsdata
%
myrow
,
0
,
hmat
%
blacsdata
%
nprow
)
IF
(
hmat
%
l_real
)
THEN
hmat
%
data_r
(
n_row
+1
:
hmat
%
matsize1
,
i
)
=
0.d0
ELSE
...
...
@@ -210,8 +208,8 @@ CONTAINS
DO
i
=
1
,
hmat
%
matsize2
! Get global column corresponding to i and number of local rows up to
! and including the diagonal, these are unchanged in H
n_col
=
indxl2g
(
i
,
nb
,
mycol
,
0
,
hmat
%
blacsdata
%
npcol
)
n_row
=
numroc
(
n_col
,
nb
,
myrow
,
0
,
hmat
%
blacsdata
%
nprow
)
n_col
=
indxl2g
(
i
,
nb
,
hmat
%
blacsdata
%
mycol
,
0
,
hmat
%
blacsdata
%
npcol
)
n_row
=
numroc
(
n_col
,
nb
,
hmat
%
blacsdata
%
myrow
,
0
,
hmat
%
blacsdata
%
nprow
)
IF
(
hmat
%
l_real
)
THEN
hmat
%
data_r
(
n_row
+1
:
hmat
%
matsize1
,
i
)
=
ev_dist
%
data_r
(
n_row
+1
:
ev_dist
%
matsize1
,
i
)
ELSE
...
...
@@ -317,8 +315,8 @@ CONTAINS
DO
i
=
1
,
hmat
%
matsize2
! Get global column corresponding to i and number of local rows up to
! and including the diagonal, these are unchanged in A
n_col
=
indxl2g
(
i
,
nb
,
mycol
,
0
,
hmat
%
blacsdata
%
npcol
)
n_row
=
numroc
(
n_col
,
nb
,
myrow
,
0
,
hmat
%
blacsdata
%
nprow
)
n_col
=
indxl2g
(
i
,
nb
,
hmat
%
blacsdata
%
mycol
,
0
,
hmat
%
blacsdata
%
npcol
)
n_row
=
numroc
(
n_col
,
nb
,
hmat
%
blacsdata
%
myrow
,
0
,
hmat
%
blacsdata
%
nprow
)
IF
(
hmat
%
l_real
)
THEN
hmat
%
data_r
(
n_row
+1
:
hmat
%
matsize1
,
i
)
=
ev_dist
%
data_r
(
n_row
+1
:
ev_dist
%
matsize1
,
i
)
ELSE
...
...
types/types_mpimat.F90
View file @
7af771fa
...
...
@@ -28,6 +28,7 @@ MODULE m_types_mpimat
!> 7,8: row/colum of grid for first row/colum of matrix
!> 9: leading dimension of local matrix
INTEGER
::
npcol
,
nprow
!> the number of columns/rows in the processor grid
INTEGER
::
mycol
,
myrow
END
TYPE
t_blacsdata
...
...
@@ -91,7 +92,7 @@ CONTAINS
SUBROUTINE
generate_full_matrix
(
mat
)
CLASS
(
t_mpimat
),
INTENT
(
INOUT
)
::
mat
INTEGER
::
i
,
j
,
i_glob
,
j_glob
,
npcol
,
nprow
,
myid
,
err
,
myrow
,
mycol
,
np
INTEGER
::
i
,
j
,
i_glob
,
j_glob
,
myid
,
err
,
np
COMPLEX
,
ALLOCATABLE
::
tmp_c
(:,:)
REAL
,
ALLOCATABLE
::
tmp_r
(:,:)
#ifdef CPP_SCALAPACK
...
...
@@ -100,16 +101,14 @@ CONTAINS
CALL
MPI_COMM_RANK
(
mat
%
blacsdata
%
mpi_com
,
myid
,
err
)
CALL
MPI_COMM_SIZE
(
mat
%
blacsdata
%
mpi_com
,
np
,
err
)
CALL
blacs_gridinfo
(
mat
%
blacsdata
%
blacs_desc
(
2
),
nprow
,
npcol
,
myrow
,
mycol
)
!Set lower part of matrix to zero
DO
i
=
1
,
mat
%
matsize1
DO
j
=
1
,
mat
%
matsize2
! Get global column corresponding to i and number of local rows up to
! and including the diagonal, these are unchanged in A
i_glob
=
indxl2g
(
i
,
mat
%
blacsdata
%
blacs_desc
(
5
),
m
yrow
,
0
,
nprow
)
j_glob
=
indxl2g
(
j
,
mat
%
blacsdata
%
blacs_desc
(
6
),
m
ycol
,
0
,
npcol
)
i_glob
=
indxl2g
(
i
,
mat
%
blacsdata
%
blacs_desc
(
5
),
m
at
%
blacsdata
%
myrow
,
0
,
mat
%
blacsdata
%
nprow
)
j_glob
=
indxl2g
(
j
,
mat
%
blacsdata
%
blacs_desc
(
6
),
m
at
%
blacsdata
%
mycol
,
0
,
mat
%
blacsdata
%
npcol
)
IF
(
i_glob
>
j_glob
)
THEN
IF
(
mat
%
l_real
)
THEN
...
...
@@ -299,9 +298,7 @@ CONTAINS
mat
%
blacsdata
%
no_use
=
1
mat
%
blacsdata
%
mpi_com
=
mpi_subcom
CALL
priv_create_blacsgrid
(
mat
%
blacsdata
%
mpi_com
,
l_2d
,
matsize1
,
matsize2
,
nbx
,
nby
,&
mat
%
blacsdata
%
blacs_desc
,&
mat
%
matsize1
,
mat
%
matsize2
,&
mat
%
blacsdata
%
npcol
,
mat
%
blacsdata
%
nprow
)
mat
%
blacsdata
,
mat
%
matsize1
,
mat
%
matsize2
)
CALL
mat
%
alloc
(
l_real
)
!Attention,sizes determined in call to priv_create_blacsgrid
!check if this matrix is actually distributed over MPI_COMM_SELF
IF
(
mpi_subcom
==
MPI_COMM_SELF
)
THEN
...
...
@@ -333,21 +330,20 @@ CONTAINS
END
SUBROUTINE
mpimat_init_template
SUBROUTINE
priv_create_blacsgrid
(
mpi_subcom
,
l_2d
,
m1
,
m2
,
nbc
,
nbr
,
sc_desc
,
local_size1
,
local_size2
,
nprow
,
npcol
)
SUBROUTINE
priv_create_blacsgrid
(
mpi_subcom
,
l_2d
,
m1
,
m2
,
nbc
,
nbr
,
blacsdata
,
local_size1
,
local_size2
)
IMPLICIT
NONE
INTEGER
,
INTENT
(
IN
)
::
mpi_subcom
INTEGER
,
INTENT
(
IN
)
::
m1
,
m2
INTEGER
,
INTENT
(
INOUT
)::
nbc
,
nbr
LOGICAL
,
INTENT
(
IN
)
::
l_2d
INTEGER
,
INTENT
(
OUT
)::
sc_desc
(:)
type
(
t_blacsdata
),
INTENT
(
OUT
)::
blacsdata
INTEGER
,
INTENT
(
OUT
)::
local_size1
,
local_size2
INTEGER
,
INTENT
(
OUT
)::
npcol
,
nprow
#ifdef CPP_MPI
INCLUDE
'mpif.h'
INTEGER
::
myrowssca
,
mycolssca
,
myrow
,
mycol
INTEGER
::
myrowssca
,
mycolssca
INTEGER
::
iamblacs
,
npblacs
,
np
,
myid
INTEGER
::
nprow2
,
npcol2
,
myrowblacs
,
mycolblacs
INTEGER
::
nprow2
,
npcol2
INTEGER
::
k
,
i
,
j
,
ictextblacs
INTEGER
::
ierr
...
...
@@ -366,34 +362,34 @@ CONTAINS
IF
(
l_2d
)
THEN
distloop
:
DO
j
=
INT
(
SQRT
(
REAL
(
np
))),
1
,
-1
IF
(
(
np
/
j
)
*
j
==
np
)
THEN
npcol
=
np
/
j
nprow
=
j
blacsdata
%
npcol
=
np
/
j
blacsdata
%
nprow
=
j
EXIT
distloop
ENDIF
ENDDO
distloop
ELSE
nbc
=
1
nbr
=
MAX
(
m1
,
m2
)
npcol
=
np
nprow
=
1
blacsdata
%
npcol
=
np
blacsdata
%
nprow
=
1
ENDIF
ALLOCATE
(
iblacsnums
(
np
),
ihelp
(
np
),
iusermap
(
nprow
,
npcol
))
ALLOCATE
(
iblacsnums
(
np
),
ihelp
(
np
),
iusermap
(
blacsdata
%
nprow
,
blacsdata
%
npcol
))
! An
nprow*
npcol processor grid will be created
! An
blacsdata%nprow*blacsdata%
npcol processor grid will be created
! Row and column index myrow, mycol of this processor in the grid
! and distribution of A and B in ScaLAPACK
! The local processor will get myrowssca rows and mycolssca columns
! of A and B
!
myrow
=
myid
/
npcol
! my row number in the BLACS nprow*
npcol grid
mycol
=
myid
-
(
myid
/
npcol
)
*
npcol
! my column number in the BLACS nprow*
npcol grid
myrow
=
myid
/
blacsdata
%
npcol
! my row number in the BLACS blacsdata%nprow*blacsdata%
npcol grid
mycol
=
myid
-
(
myid
/
blacsdata
%
npcol
)
*
blacsdata
%
npcol
! my column number in the BLACS blacsdata%nprow*blacsdata%
npcol grid
!
! Now allocate Asca to put the elements of Achi or receivebuffer to
!
myrowssca
=
(
m1
-1
)/(
nbr
*
nprow
)
*
nbr
+
MIN
(
MAX
(
m1
-
(
m1
-1
)/(
nbr
*
nprow
)
*
nbr
*
nprow
-
nbr
*
myrow
,
0
),
nbr
)
myrowssca
=
(
m1
-1
)/(
nbr
*
blacsdata
%
nprow
)
*
nbr
+
MIN
(
MAX
(
m1
-
(
m1
-1
)/(
nbr
*
blacsdata
%
nprow
)
*
nbr
*
blacsdata
%
nprow
-
nbr
*
myrow
,
0
),
nbr
)
! Number of rows the local process gets in ScaLAPACK distribution
mycolssca
=
(
m2
-1
)/(
nbc
*
npcol
)
*
nbc
+
MIN
(
MAX
(
m2
-
(
m2
-1
)/(
nbc
*
npcol
)
*
nbc
*
npcol
-
nbc
*
mycol
,
0
),
nbc
)
mycolssca
=
(
m2
-1
)/(
nbc
*
blacsdata
%
npcol
)
*
nbc
+
MIN
(
MAX
(
m2
-
(
m2
-1
)/(
nbc
*
blacsdata
%
npcol
)
*
nbc
*
blacsdata
%
npcol
-
nbc
*
mycol
,
0
),
nbc
)
!Get BLACS ranks for all MPI ranks
...
...
@@ -407,8 +403,8 @@ CONTAINS
! iblacsnums(i) is the BLACS-process number of MPI-process i-1
k
=
1
DO
i
=
1
,
nprow
DO
j
=
1
,
npcol
DO
i
=
1
,
blacsdata
%
nprow
DO
j
=
1
,
blacsdata
%
npcol
iusermap
(
i
,
j
)
=
iblacsnums
(
k
)
k
=
k
+
1
ENDDO
...
...
@@ -416,23 +412,23 @@ CONTAINS
!Get the Blacs default context
CALL
BLACS_GET
(
0
,
0
,
ictextblacs
)
! Create the Grid
CALL
BLACS_GRIDMAP
(
ictextblacs
,
iusermap
,
nprow
,
nprow
,
npcol
)
CALL
BLACS_GRIDMAP
(
ictextblacs
,
iusermap
,
size
(
iusermap
,
1
),
blacsdata
%
nprow
,
blacsdata
%
npcol
)
! Now control, whether the BLACS grid is the one we wanted
CALL
BLACS_GRIDINFO
(
ictextblacs
,
nprow2
,
npcol2
,
myrowblacs
,
mycolblacs
)
IF
(
nprow2
/
=
nprow
)
THEN
CALL
BLACS_GRIDINFO
(
ictextblacs
,
nprow2
,
npcol2
,
blacsdata
%
myrow
,
blacsdata
%
mycol
)
IF
(
nprow2
/
=
blacsdata
%
nprow
)
THEN
WRITE
(
6
,
*
)
'Wrong number of rows in BLACS grid'
WRITE
(
6
,
*
)
'nprow='
,
nprow
,
' nprow2='
,
nprow2
WRITE
(
6
,
*
)
'nprow='
,
blacsdata
%
nprow
,
' nprow2='
,
nprow2
call
judft_error
(
'Wrong number of rows in BLACS grid'
)
ENDIF
IF
(
npcol2
/
=
npcol
)
THEN
IF
(
npcol2
/
=
blacsdata
%
npcol
)
THEN
WRITE
(
6
,
*
)
'Wrong number of columns in BLACS grid'
WRITE
(
6
,
*
)
'npcol='
,
npcol
,
' npcol2='
,
npcol2
WRITE
(
6
,
*
)
'npcol='
,
blacsdata
%
npcol
,
' npcol2='
,
npcol2
call
judft_error
(
'Wrong number of columns in BLACS grid'
)
ENDIF
!Create the descriptors
CALL
descinit
(
sc
_desc
,
m1
,
m2
,
nbr
,
nbc
,
0
,
0
,
ictextblacs
,
myrowssca
,
ierr
)
CALL
descinit
(
blacsdata
%
blacs
_desc
,
m1
,
m2
,
nbr
,
nbc
,
0
,
0
,
ictextblacs
,
myrowssca
,
ierr
)
IF
(
ierr
/
=
0
)
call
judft_error
(
'Creation of BLACS descriptor failed'
)
local_size1
=
myrowssca
local_size2
=
mycolssca
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment