1*26cec326SBarry Smith /* 2*26cec326SBarry Smith used by MPIAIJ, BAIJ and SBAIJ to reduce code duplication 3*26cec326SBarry Smith 4*26cec326SBarry Smith define TYPE to AIJ BAIJ or SBAIJ 5*26cec326SBarry Smith TYPE_SBAIJ for SBAIJ matrix 6*26cec326SBarry Smith 7*26cec326SBarry Smith */ 8*26cec326SBarry Smith 9*26cec326SBarry Smith static PetscErrorCode MatSetValues_MPI_Hash(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode addv) 10*26cec326SBarry Smith { 11*26cec326SBarry Smith PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data; 12*26cec326SBarry Smith PetscInt rStart, rEnd, cStart, cEnd; 13*26cec326SBarry Smith #if defined(TYPE_SBAIJ) 14*26cec326SBarry Smith PetscInt bs; 15*26cec326SBarry Smith #endif 16*26cec326SBarry Smith 17*26cec326SBarry Smith PetscFunctionBegin; 18*26cec326SBarry Smith PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 19*26cec326SBarry Smith PetscCall(MatGetOwnershipRangeColumn(A, &cStart, &cEnd)); 20*26cec326SBarry Smith #if defined(TYPE_SBAIJ) 21*26cec326SBarry Smith PetscCall(MatGetBlockSize(A, &bs)); 22*26cec326SBarry Smith #endif 23*26cec326SBarry Smith for (PetscInt r = 0; r < m; ++r) { 24*26cec326SBarry Smith PetscScalar value; 25*26cec326SBarry Smith if (rows[r] < 0) continue; 26*26cec326SBarry Smith if (rows[r] < rStart || rows[r] >= rEnd) { 27*26cec326SBarry Smith if (a->roworiented) { 28*26cec326SBarry Smith PetscCall(MatStashValuesRow_Private(&A->stash, rows[r], n, cols, values + r * n, PETSC_FALSE)); 29*26cec326SBarry Smith } else { 30*26cec326SBarry Smith PetscCall(MatStashValuesCol_Private(&A->stash, rows[r], n, cols, values + r, m, PETSC_FALSE)); 31*26cec326SBarry Smith } 32*26cec326SBarry Smith } else { 33*26cec326SBarry Smith for (PetscInt c = 0; c < n; ++c) { 34*26cec326SBarry Smith #if defined(TYPE_SBAIJ) 35*26cec326SBarry Smith if (cols[c] / bs < rows[r] / bs) continue; 36*26cec326SBarry Smith #else 37*26cec326SBarry Smith if (cols[c] < 0) continue; 38*26cec326SBarry Smith #endif 39*26cec326SBarry Smith value = values ? ((a->roworiented) ? values[r * n + c] : values[r + m * c]) : 0; 40*26cec326SBarry Smith if (cols[c] >= cStart && cols[c] < cEnd) PetscCall(MatSetValue(a->A, rows[r] - rStart, cols[c] - cStart, value, addv)); 41*26cec326SBarry Smith else PetscCall(MatSetValue(a->B, rows[r] - rStart, cols[c], value, addv)); 42*26cec326SBarry Smith } 43*26cec326SBarry Smith } 44*26cec326SBarry Smith } 45*26cec326SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 46*26cec326SBarry Smith } 47*26cec326SBarry Smith 48*26cec326SBarry Smith static PetscErrorCode MatAssemblyBegin_MPI_Hash(Mat A, MatAssemblyType type) 49*26cec326SBarry Smith { 50*26cec326SBarry Smith PetscInt nstash, reallocs; 51*26cec326SBarry Smith 52*26cec326SBarry Smith PetscFunctionBegin; 53*26cec326SBarry Smith PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range)); 54*26cec326SBarry Smith PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs)); 55*26cec326SBarry Smith PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs)); 56*26cec326SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 57*26cec326SBarry Smith } 58*26cec326SBarry Smith 59*26cec326SBarry Smith static PetscErrorCode MatAssemblyEnd_MPI_Hash(Mat A, MatAssemblyType type) 60*26cec326SBarry Smith { 61*26cec326SBarry Smith PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data; 62*26cec326SBarry Smith PetscMPIInt n; 63*26cec326SBarry Smith PetscScalar *val; 64*26cec326SBarry Smith PetscInt *row, *col; 65*26cec326SBarry Smith PetscInt j, ncols, flg, rstart; 66*26cec326SBarry Smith 67*26cec326SBarry Smith PetscFunctionBegin; 68*26cec326SBarry Smith while (1) { 69*26cec326SBarry Smith PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg)); 70*26cec326SBarry Smith if (!flg) break; 71*26cec326SBarry Smith 72*26cec326SBarry Smith for (PetscInt i = 0; i < n;) { 73*26cec326SBarry Smith /* Now identify the consecutive vals belonging to the same row */ 74*26cec326SBarry Smith for (j = i, rstart = row[j]; j < n; j++) { 75*26cec326SBarry Smith if (row[j] != rstart) break; 76*26cec326SBarry Smith } 77*26cec326SBarry Smith if (j < n) ncols = j - i; 78*26cec326SBarry Smith else ncols = n - i; 79*26cec326SBarry Smith /* Now assemble all these values with a single function call */ 80*26cec326SBarry Smith PetscCall(MatSetValues_MPI_Hash(A, 1, row + i, ncols, col + i, val + i, A->insertmode)); 81*26cec326SBarry Smith i = j; 82*26cec326SBarry Smith } 83*26cec326SBarry Smith } 84*26cec326SBarry Smith PetscCall(MatStashScatterEnd_Private(&A->stash)); 85*26cec326SBarry Smith if (type != MAT_FINAL_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS); 86*26cec326SBarry Smith 87*26cec326SBarry Smith A->insertmode = NOT_SET_VALUES; /* this was set by the previous calls to MatSetValues() */ 88*26cec326SBarry Smith 89*26cec326SBarry Smith PetscCall(PetscMemcpy(&A->ops, &a->cops, sizeof(*(A->ops)))); 90*26cec326SBarry Smith 91*26cec326SBarry Smith PetscCall(MatAssemblyBegin(a->A, MAT_FINAL_ASSEMBLY)); 92*26cec326SBarry Smith PetscCall(MatAssemblyEnd(a->A, MAT_FINAL_ASSEMBLY)); 93*26cec326SBarry Smith PetscCall(MatAssemblyBegin(a->B, MAT_FINAL_ASSEMBLY)); 94*26cec326SBarry Smith PetscCall(MatAssemblyEnd(a->B, MAT_FINAL_ASSEMBLY)); 95*26cec326SBarry Smith PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 96*26cec326SBarry Smith PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 97*26cec326SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 98*26cec326SBarry Smith } 99*26cec326SBarry Smith 100*26cec326SBarry Smith static PetscErrorCode MatDestroy_MPI_Hash(Mat A) 101*26cec326SBarry Smith { 102*26cec326SBarry Smith PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data; 103*26cec326SBarry Smith 104*26cec326SBarry Smith PetscFunctionBegin; 105*26cec326SBarry Smith PetscCall(MatStashDestroy_Private(&A->stash)); 106*26cec326SBarry Smith PetscCall(MatDestroy(&a->A)); 107*26cec326SBarry Smith PetscCall(MatDestroy(&a->B)); 108*26cec326SBarry Smith PetscCall((*a->cops.destroy)(A)); 109*26cec326SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 110*26cec326SBarry Smith } 111*26cec326SBarry Smith 112*26cec326SBarry Smith static PetscErrorCode MatZeroEntries_MPI_Hash(Mat A) 113*26cec326SBarry Smith { 114*26cec326SBarry Smith PetscFunctionBegin; 115*26cec326SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 116*26cec326SBarry Smith } 117*26cec326SBarry Smith 118*26cec326SBarry Smith static PetscErrorCode MatSetRandom_MPI_Hash(Mat A, PetscRandom r) 119*26cec326SBarry Smith { 120*26cec326SBarry Smith PetscFunctionBegin; 121*26cec326SBarry Smith SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Must set preallocation first"); 122*26cec326SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 123*26cec326SBarry Smith } 124*26cec326SBarry Smith 125*26cec326SBarry Smith static PetscErrorCode MatSetUp_MPI_Hash(Mat A) 126*26cec326SBarry Smith { 127*26cec326SBarry Smith PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data; 128*26cec326SBarry Smith PetscMPIInt size; 129*26cec326SBarry Smith #if !defined(TYPE_AIJ) 130*26cec326SBarry Smith PetscInt bs; 131*26cec326SBarry Smith #endif 132*26cec326SBarry Smith 133*26cec326SBarry Smith PetscFunctionBegin; 134*26cec326SBarry Smith PetscCall(PetscInfo(A, "Using hash-based MatSetValues() for MATMPISBAIJ because no preallocation provided\n")); 135*26cec326SBarry Smith PetscCall(PetscLayoutSetUp(A->rmap)); 136*26cec326SBarry Smith PetscCall(PetscLayoutSetUp(A->cmap)); 137*26cec326SBarry Smith if (A->rmap->bs < 1) A->rmap->bs = 1; 138*26cec326SBarry Smith if (A->cmap->bs < 1) A->cmap->bs = 1; 139*26cec326SBarry Smith PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 140*26cec326SBarry Smith 141*26cec326SBarry Smith #if !defined(TYPE_AIJ) 142*26cec326SBarry Smith PetscCall(MatGetBlockSize(A, &bs)); 143*26cec326SBarry Smith /* these values are set in MatMPISBAIJSetPreallocation() */ 144*26cec326SBarry Smith a->bs2 = bs * bs; 145*26cec326SBarry Smith a->mbs = A->rmap->n / bs; 146*26cec326SBarry Smith a->nbs = A->cmap->n / bs; 147*26cec326SBarry Smith a->Mbs = A->rmap->N / bs; 148*26cec326SBarry Smith a->Nbs = A->cmap->N / bs; 149*26cec326SBarry Smith 150*26cec326SBarry Smith for (PetscInt i = 0; i <= a->size; i++) a->rangebs[i] = A->rmap->range[i] / bs; 151*26cec326SBarry Smith a->rstartbs = A->rmap->rstart / bs; 152*26cec326SBarry Smith a->rendbs = A->rmap->rend / bs; 153*26cec326SBarry Smith a->cstartbs = A->cmap->rstart / bs; 154*26cec326SBarry Smith a->cendbs = A->cmap->rend / bs; 155*26cec326SBarry Smith PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)A), A->rmap->bs, &A->bstash)); 156*26cec326SBarry Smith #endif 157*26cec326SBarry Smith 158*26cec326SBarry Smith PetscCall(MatCreate(PETSC_COMM_SELF, &a->A)); 159*26cec326SBarry Smith PetscCall(MatSetSizes(a->A, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n)); 160*26cec326SBarry Smith PetscCall(MatSetBlockSizesFromMats(a->A, A, A)); 161*26cec326SBarry Smith #if defined(SUB_TYPE_CUSPARSE) 162*26cec326SBarry Smith PetscCall(MatSetType(a->A, MATSEQAIJCUSPARSE)); 163*26cec326SBarry Smith #else 164*26cec326SBarry Smith PetscCall(MatSetType(a->A, PetscConcat(MATSEQ, TYPE))); 165*26cec326SBarry Smith #endif 166*26cec326SBarry Smith PetscCall(MatSetUp(a->A)); 167*26cec326SBarry Smith 168*26cec326SBarry Smith PetscCall(MatCreate(PETSC_COMM_SELF, &a->B)); 169*26cec326SBarry Smith PetscCall(MatSetSizes(a->B, A->rmap->n, size > 1 ? A->cmap->N : 0, A->rmap->n, size > 1 ? A->cmap->N : 0)); 170*26cec326SBarry Smith PetscCall(MatSetBlockSizesFromMats(a->B, A, A)); 171*26cec326SBarry Smith #if defined(TYPE_SBAIJ) 172*26cec326SBarry Smith PetscCall(MatSetType(a->B, MATSEQBAIJ)); 173*26cec326SBarry Smith #else 174*26cec326SBarry Smith #if defined(SUB_TYPE_CUSPARSE) 175*26cec326SBarry Smith PetscCall(MatSetType(a->B, MATSEQAIJCUSPARSE)); 176*26cec326SBarry Smith #else 177*26cec326SBarry Smith PetscCall(MatSetType(a->B, PetscConcat(MATSEQ, TYPE))); 178*26cec326SBarry Smith #endif 179*26cec326SBarry Smith #endif 180*26cec326SBarry Smith PetscCall(MatSetUp(a->B)); 181*26cec326SBarry Smith 182*26cec326SBarry Smith /* keep a record of the operations so they can be reset when the hash handling is complete */ 183*26cec326SBarry Smith PetscCall(PetscMemcpy(&a->cops, &A->ops, sizeof(*(A->ops)))); 184*26cec326SBarry Smith 185*26cec326SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_MPI_Hash; 186*26cec326SBarry Smith A->ops->assemblyend = MatAssemblyEnd_MPI_Hash; 187*26cec326SBarry Smith A->ops->setvalues = MatSetValues_MPI_Hash; 188*26cec326SBarry Smith A->ops->destroy = MatDestroy_MPI_Hash; 189*26cec326SBarry Smith A->ops->zeroentries = MatZeroEntries_MPI_Hash; 190*26cec326SBarry Smith A->ops->setrandom = MatSetRandom_MPI_Hash; 191*26cec326SBarry Smith A->ops->setvaluesblocked = NULL; 192*26cec326SBarry Smith 193*26cec326SBarry Smith A->preallocated = PETSC_TRUE; 194*26cec326SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 195*26cec326SBarry Smith } 196