Skip to content

Commit 0017a2e

Browse files
committed
Bug fix: in "dReDistributre_A", dereference several uninitialized/unallocated arrays even though their sizes are zero.
Upgrade to v7.1.1
1 parent e3eba61 commit 0017a2e

File tree

6 files changed

+33
-29
lines changed

6 files changed

+33
-29
lines changed

CMakeLists.txt

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -12,7 +12,7 @@ cmake_minimum_required(VERSION 3.18.1 FATAL_ERROR)
1212
project(SuperLU_DIST C CXX)
1313
set(VERSION_MAJOR "7")
1414
set(VERSION_MINOR "1")
15-
set(VERSION_BugFix "0")
15+
set(VERSION_BugFix "1")
1616
set(PROJECT_VERSION ${VERSION_MAJOR}.${VERSION_MINOR}.${VERSION_BugFix})
1717

1818
list(APPEND CMAKE_MODULE_PATH "${PROJECT_SOURCE_DIR}/cmake")

README.md

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,4 @@
1-
# SuperLU_DIST (version 7.1.0) <img align=center width="55" alt="superlu" src="https://user-images.githubusercontent.com/11741943/103982988-5a9a9d00-5139-11eb-9ac4-a55e80a79f8d.png">
1+
# SuperLU_DIST (version 7.1.1) <img align=center width="55" alt="superlu" src="https://user-images.githubusercontent.com/11741943/103982988-5a9a9d00-5139-11eb-9ac4-a55e80a79f8d.png">
22

33
[![Build Status](https://travis-ci.org/xiaoyeli/superlu_dist.svg?branch=master)](https://travis-ci.org/xiaoyeli/superlu_dist)
44
[Nightly tests](http://my.cdash.org/index.php?project=superlu_dist)
@@ -25,7 +25,7 @@ acceleration capabilities.
2525
Table of Contents
2626
=================
2727

28-
* [SuperLU_DIST (version 7.1.0) <a href="https://user-images.githubusercontent.com/11741943/103982988-5a9a9d00-5139-11eb-9ac4-a55e80a79f8d.png" target="_blank" rel="nofollow"><img align="center" width="55" alt="superlu" src="https://user-images.githubusercontent.com/11741943/103982988-5a9a9d00-5139-11eb-9ac4-a55e80a79f8d.png" style="max-width:100%;"></a>](#superlu_dist-version-70---)
28+
* [SuperLU_DIST (version 7.1.1) <a href="https://user-images.githubusercontent.com/11741943/103982988-5a9a9d00-5139-11eb-9ac4-a55e80a79f8d.png" target="_blank" rel="nofollow"><img align="center" width="55" alt="superlu" src="https://user-images.githubusercontent.com/11741943/103982988-5a9a9d00-5139-11eb-9ac4-a55e80a79f8d.png" style="max-width:100%;"></a>](#superlu_dist-version-70---)
2929
* [Directory structure of the source code](#directory-structure-of-the-source-code)
3030
* [Installation](#installation)
3131
* [Installation option 1: Using CMake build system.](#installation-option-1-using-cmake-build-system)

SRC/pddistribute.c

Lines changed: 11 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -13,9 +13,10 @@ at the top-level directory.
1313
/*! @file
1414
* \brief Re-distribute A on the 2D process mesh.
1515
* <pre>
16-
* -- Distributed SuperLU routine (version 2.3) --
16+
* -- Distributed SuperLU routine (version 7.1.1) --
1717
* Lawrence Berkeley National Lab, Univ. of California Berkeley.
1818
* October 15, 2008
19+
* October 18, 2021, minor fix, v7.1.1
1920
* </pre>
2021
*/
2122

@@ -141,8 +142,8 @@ dReDistribute_A(SuperMatrix *A, dScalePermstruct_t *ScalePermstruct,
141142
ABORT("Malloc fails for ia[].");
142143
if ( !(aij = doubleMalloc_dist(k)) )
143144
ABORT("Malloc fails for aij[].");
145+
ja = ia + k;
144146
}
145-
ja = ia + k;
146147

147148
/* Allocate temporary storage for sending/receiving the A triplets. */
148149
if ( procs > 1 ) {
@@ -170,9 +171,9 @@ dReDistribute_A(SuperMatrix *A, dScalePermstruct_t *ScalePermstruct,
170171

171172
for (i = 0, j = 0, p = 0; p < procs; ++p) {
172173
if ( p != iam ) {
173-
ia_send[p] = &index[i];
174+
if (nnzToSend[p] > 0) ia_send[p] = &index[i];
174175
i += 2 * nnzToSend[p]; /* ia/ja indices alternate */
175-
aij_send[p] = &nzval[j];
176+
if (nnzToSend[p] > 0) aij_send[p] = &nzval[j];
176177
j += nnzToSend[p];
177178
}
178179
}
@@ -216,8 +217,8 @@ dReDistribute_A(SuperMatrix *A, dScalePermstruct_t *ScalePermstruct,
216217
NOTE: Can possibly use MPI_Alltoallv.
217218
------------------------------------------------------------*/
218219
for (p = 0; p < procs; ++p) {
219-
if ( p != iam && nnzToSend[p]>0 ) { // cause two of the tests to hang
220-
// if ( p != iam ) {
220+
if ( p != iam && nnzToSend[p] > 0 ) {
221+
//if ( p != iam ) {
221222
it = 2*nnzToSend[p];
222223
MPI_Isend( ia_send[p], it, mpi_int_t,
223224
p, iam, grid->comm, &send_req[p] );
@@ -228,8 +229,8 @@ dReDistribute_A(SuperMatrix *A, dScalePermstruct_t *ScalePermstruct,
228229
}
229230

230231
for (p = 0; p < procs; ++p) {
231-
if ( p != iam && nnzToRecv[p]>0 ) {
232-
//if ( p != iam ) {
232+
if ( p != iam && nnzToRecv[p] > 0 ) {
233+
//if ( p != iam ) {
233234
it = 2*nnzToRecv[p];
234235
MPI_Recv( itemp, it, mpi_int_t, p, p, grid->comm, &status );
235236
it = nnzToRecv[p];
@@ -248,8 +249,8 @@ dReDistribute_A(SuperMatrix *A, dScalePermstruct_t *ScalePermstruct,
248249
}
249250

250251
for (p = 0; p < procs; ++p) {
251-
if ( p != iam && nnzToSend[p] > 0 ) {
252-
//if ( p != iam ) {
252+
if ( p != iam && nnzToSend[p] > 0 ) { // cause two of the tests to hang
253+
//if ( p != iam ) {
253254
MPI_Wait( &send_req[p], &status);
254255
MPI_Wait( &send_req[procs+p], &status);
255256
}

SRC/psdistribute.c

Lines changed: 5 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -13,9 +13,10 @@ at the top-level directory.
1313
/*! @file
1414
* \brief Re-distribute A on the 2D process mesh.
1515
* <pre>
16-
* -- Distributed SuperLU routine (version 2.3) --
16+
* -- Distributed SuperLU routine (version 7.1.1) --
1717
* Lawrence Berkeley National Lab, Univ. of California Berkeley.
1818
* October 15, 2008
19+
* October 18, 2021, minor fix, v7.1.1
1920
* </pre>
2021
*/
2122

@@ -141,8 +142,8 @@ sReDistribute_A(SuperMatrix *A, sScalePermstruct_t *ScalePermstruct,
141142
ABORT("Malloc fails for ia[].");
142143
if ( !(aij = floatMalloc_dist(k)) )
143144
ABORT("Malloc fails for aij[].");
145+
ja = ia + k;
144146
}
145-
ja = ia + k;
146147

147148
/* Allocate temporary storage for sending/receiving the A triplets. */
148149
if ( procs > 1 ) {
@@ -170,9 +171,9 @@ sReDistribute_A(SuperMatrix *A, sScalePermstruct_t *ScalePermstruct,
170171

171172
for (i = 0, j = 0, p = 0; p < procs; ++p) {
172173
if ( p != iam ) {
173-
ia_send[p] = &index[i];
174+
if (nnzToSend[p] > 0) ia_send[p] = &index[i];
174175
i += 2 * nnzToSend[p]; /* ia/ja indices alternate */
175-
aij_send[p] = &nzval[j];
176+
if (nnzToSend[p] > 0) aij_send[p] = &nzval[j];
176177
j += nnzToSend[p];
177178
}
178179
}

SRC/pzdistribute.c

Lines changed: 11 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -12,9 +12,10 @@ at the top-level directory.
1212
/*! @file
1313
* \brief Re-distribute A on the 2D process mesh.
1414
* <pre>
15-
* -- Distributed SuperLU routine (version 2.3) --
15+
* -- Distributed SuperLU routine (version 7.1.1) --
1616
* Lawrence Berkeley National Lab, Univ. of California Berkeley.
1717
* October 15, 2008
18+
* October 18, 2021, minor fix, v7.1.1
1819
* </pre>
1920
*/
2021

@@ -140,8 +141,8 @@ zReDistribute_A(SuperMatrix *A, zScalePermstruct_t *ScalePermstruct,
140141
ABORT("Malloc fails for ia[].");
141142
if ( !(aij = doublecomplexMalloc_dist(k)) )
142143
ABORT("Malloc fails for aij[].");
144+
ja = ia + k;
143145
}
144-
ja = ia + k;
145146

146147
/* Allocate temporary storage for sending/receiving the A triplets. */
147148
if ( procs > 1 ) {
@@ -169,9 +170,9 @@ zReDistribute_A(SuperMatrix *A, zScalePermstruct_t *ScalePermstruct,
169170

170171
for (i = 0, j = 0, p = 0; p < procs; ++p) {
171172
if ( p != iam ) {
172-
ia_send[p] = &index[i];
173+
if (nnzToSend[p] > 0) ia_send[p] = &index[i];
173174
i += 2 * nnzToSend[p]; /* ia/ja indices alternate */
174-
aij_send[p] = &nzval[j];
175+
if (nnzToSend[p] > 0) aij_send[p] = &nzval[j];
175176
j += nnzToSend[p];
176177
}
177178
}
@@ -215,8 +216,8 @@ zReDistribute_A(SuperMatrix *A, zScalePermstruct_t *ScalePermstruct,
215216
NOTE: Can possibly use MPI_Alltoallv.
216217
------------------------------------------------------------*/
217218
for (p = 0; p < procs; ++p) {
218-
if ( p != iam && nnzToSend[p] > 0 ) {
219-
//if ( p != iam ) {
219+
if ( p != iam && nnzToSend[p] > 0 ) {
220+
//if ( p != iam ) {
220221
it = 2*nnzToSend[p];
221222
MPI_Isend( ia_send[p], it, mpi_int_t,
222223
p, iam, grid->comm, &send_req[p] );
@@ -227,8 +228,8 @@ zReDistribute_A(SuperMatrix *A, zScalePermstruct_t *ScalePermstruct,
227228
}
228229

229230
for (p = 0; p < procs; ++p) {
230-
if ( p != iam && nnzToRecv[p] > 0 ) {
231-
//if ( p != iam ) {
231+
if ( p != iam && nnzToRecv[p] > 0 ) {
232+
//if ( p != iam ) {
232233
it = 2*nnzToRecv[p];
233234
MPI_Recv( itemp, it, mpi_int_t, p, p, grid->comm, &status );
234235
it = nnzToRecv[p];
@@ -247,8 +248,8 @@ zReDistribute_A(SuperMatrix *A, zScalePermstruct_t *ScalePermstruct,
247248
}
248249

249250
for (p = 0; p < procs; ++p) {
250-
if ( p != iam && nnzToSend[p] > 0 ) {
251-
//if ( p != iam ) {
251+
if ( p != iam && nnzToSend[p] > 0 ) { // cause two of the tests to hang
252+
//if ( p != iam ) {
252253
MPI_Wait( &send_req[p], &status);
253254
MPI_Wait( &send_req[procs+p], &status);
254255
}

SRC/superlu_defs.h

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -25,6 +25,7 @@ at the top-level directory.
2525
* October 23, 2020 version 6.4.0
2626
* May 12, 2021 version 7.0.0
2727
* October 5, 2021 version 7.1.0
28+
* October 18, 2021 version 7.1.1
2829
* </pre>
2930
*/
3031

@@ -76,8 +77,8 @@ at the top-level directory.
7677
*/
7778
#define SUPERLU_DIST_MAJOR_VERSION 7
7879
#define SUPERLU_DIST_MINOR_VERSION 1
79-
#define SUPERLU_DIST_PATCH_VERSION 0
80-
#define SUPERLU_DIST_RELEASE_DATE "October 5, 2021"
80+
#define SUPERLU_DIST_PATCH_VERSION 1
81+
#define SUPERLU_DIST_RELEASE_DATE "October 18, 2021"
8182

8283
#include "superlu_dist_config.h"
8384
/* Define my integer size int_t */

0 commit comments

Comments
 (0)