xref: /petsc/src/mat/impls/aij/mpi/mpihashmat.h (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
126cec326SBarry Smith /*
226cec326SBarry Smith    used by MPIAIJ, BAIJ and SBAIJ to reduce code duplication
326cec326SBarry Smith 
426cec326SBarry Smith      define TYPE to AIJ BAIJ or SBAIJ
526cec326SBarry Smith             TYPE_SBAIJ for SBAIJ matrix
626cec326SBarry Smith 
726cec326SBarry Smith */
826cec326SBarry Smith 
926cec326SBarry Smith static PetscErrorCode MatSetValues_MPI_Hash(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode addv)
1026cec326SBarry Smith {
1126cec326SBarry Smith   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
1226cec326SBarry Smith   PetscInt rStart, rEnd, cStart, cEnd;
1326cec326SBarry Smith #if defined(TYPE_SBAIJ)
1426cec326SBarry Smith   PetscInt bs;
1526cec326SBarry Smith #endif
1626cec326SBarry Smith 
1726cec326SBarry Smith   PetscFunctionBegin;
1826cec326SBarry Smith   PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd));
1926cec326SBarry Smith   PetscCall(MatGetOwnershipRangeColumn(A, &cStart, &cEnd));
2026cec326SBarry Smith #if defined(TYPE_SBAIJ)
2126cec326SBarry Smith   PetscCall(MatGetBlockSize(A, &bs));
2226cec326SBarry Smith #endif
2326cec326SBarry Smith   for (PetscInt r = 0; r < m; ++r) {
2426cec326SBarry Smith     PetscScalar value;
2526cec326SBarry Smith     if (rows[r] < 0) continue;
2626cec326SBarry Smith     if (rows[r] < rStart || rows[r] >= rEnd) {
273ac1eb50SJunchao Zhang       PetscCheck(!A->nooffprocentries, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Setting off process row %" PetscInt_FMT " even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set", rows[r]);
283ac1eb50SJunchao Zhang       if (!a->donotstash) {
293ac1eb50SJunchao Zhang         A->assembled = PETSC_FALSE;
3026cec326SBarry Smith         if (a->roworiented) {
3126cec326SBarry Smith           PetscCall(MatStashValuesRow_Private(&A->stash, rows[r], n, cols, values + r * n, PETSC_FALSE));
3226cec326SBarry Smith         } else {
3326cec326SBarry Smith           PetscCall(MatStashValuesCol_Private(&A->stash, rows[r], n, cols, values + r, m, PETSC_FALSE));
3426cec326SBarry Smith         }
353ac1eb50SJunchao Zhang       }
3626cec326SBarry Smith     } else {
3726cec326SBarry Smith       for (PetscInt c = 0; c < n; ++c) {
3826cec326SBarry Smith #if defined(TYPE_SBAIJ)
3926cec326SBarry Smith         if (cols[c] / bs < rows[r] / bs) continue;
4026cec326SBarry Smith #else
4126cec326SBarry Smith         if (cols[c] < 0) continue;
4226cec326SBarry Smith #endif
43*f4f49eeaSPierre Jolivet         value = values ? (a->roworiented ? values[r * n + c] : values[r + m * c]) : 0;
4426cec326SBarry Smith         if (cols[c] >= cStart && cols[c] < cEnd) PetscCall(MatSetValue(a->A, rows[r] - rStart, cols[c] - cStart, value, addv));
4526cec326SBarry Smith         else PetscCall(MatSetValue(a->B, rows[r] - rStart, cols[c], value, addv));
4626cec326SBarry Smith       }
4726cec326SBarry Smith     }
4826cec326SBarry Smith   }
4926cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
5026cec326SBarry Smith }
5126cec326SBarry Smith 
52a57c2057SPierre Jolivet static PetscErrorCode MatAssemblyBegin_MPI_Hash(Mat A, PETSC_UNUSED MatAssemblyType type)
5326cec326SBarry Smith {
54cfe56a0cSJunchao Zhang   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
5526cec326SBarry Smith   PetscInt nstash, reallocs;
5626cec326SBarry Smith 
5726cec326SBarry Smith   PetscFunctionBegin;
58cfe56a0cSJunchao Zhang   if (a->donotstash || A->nooffprocentries) PetscFunctionReturn(PETSC_SUCCESS);
5926cec326SBarry Smith   PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range));
6026cec326SBarry Smith   PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs));
6126cec326SBarry Smith   PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
6226cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
6326cec326SBarry Smith }
6426cec326SBarry Smith 
6526cec326SBarry Smith static PetscErrorCode MatAssemblyEnd_MPI_Hash(Mat A, MatAssemblyType type)
6626cec326SBarry Smith {
6726cec326SBarry Smith   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
6826cec326SBarry Smith   PetscMPIInt  n;
6926cec326SBarry Smith   PetscScalar *val;
7026cec326SBarry Smith   PetscInt    *row, *col;
7126cec326SBarry Smith   PetscInt     j, ncols, flg, rstart;
7226cec326SBarry Smith 
7326cec326SBarry Smith   PetscFunctionBegin;
74cfe56a0cSJunchao Zhang   if (!a->donotstash && !A->nooffprocentries) {
7526cec326SBarry Smith     while (1) {
7626cec326SBarry Smith       PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg));
7726cec326SBarry Smith       if (!flg) break;
7826cec326SBarry Smith 
7926cec326SBarry Smith       for (PetscInt i = 0; i < n;) {
8026cec326SBarry Smith         /* Now identify the consecutive vals belonging to the same row */
8126cec326SBarry Smith         for (j = i, rstart = row[j]; j < n; j++) {
8226cec326SBarry Smith           if (row[j] != rstart) break;
8326cec326SBarry Smith         }
8426cec326SBarry Smith         if (j < n) ncols = j - i;
8526cec326SBarry Smith         else ncols = n - i;
8626cec326SBarry Smith         /* Now assemble all these values with a single function call */
8726cec326SBarry Smith         PetscCall(MatSetValues_MPI_Hash(A, 1, row + i, ncols, col + i, val + i, A->insertmode));
8826cec326SBarry Smith         i = j;
8926cec326SBarry Smith       }
9026cec326SBarry Smith     }
9126cec326SBarry Smith     PetscCall(MatStashScatterEnd_Private(&A->stash));
92cfe56a0cSJunchao Zhang   }
9326cec326SBarry Smith   if (type != MAT_FINAL_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
9426cec326SBarry Smith 
9526cec326SBarry Smith   A->insertmode = NOT_SET_VALUES; /* this was set by the previous calls to MatSetValues() */
9626cec326SBarry Smith 
97aea10558SJacob Faibussowitsch   A->ops[0]      = a->cops;
98ad79cf63SBarry Smith   A->hash_active = PETSC_FALSE;
9926cec326SBarry Smith 
10026cec326SBarry Smith   PetscCall(MatAssemblyBegin(a->A, MAT_FINAL_ASSEMBLY));
10126cec326SBarry Smith   PetscCall(MatAssemblyEnd(a->A, MAT_FINAL_ASSEMBLY));
10226cec326SBarry Smith   PetscCall(MatAssemblyBegin(a->B, MAT_FINAL_ASSEMBLY));
10326cec326SBarry Smith   PetscCall(MatAssemblyEnd(a->B, MAT_FINAL_ASSEMBLY));
10426cec326SBarry Smith   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
10526cec326SBarry Smith   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
10626cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
10726cec326SBarry Smith }
10826cec326SBarry Smith 
10926cec326SBarry Smith static PetscErrorCode MatDestroy_MPI_Hash(Mat A)
11026cec326SBarry Smith {
11126cec326SBarry Smith   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
11226cec326SBarry Smith 
11326cec326SBarry Smith   PetscFunctionBegin;
11426cec326SBarry Smith   PetscCall(MatStashDestroy_Private(&A->stash));
11526cec326SBarry Smith   PetscCall(MatDestroy(&a->A));
11626cec326SBarry Smith   PetscCall(MatDestroy(&a->B));
11726cec326SBarry Smith   PetscCall((*a->cops.destroy)(A));
11826cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
11926cec326SBarry Smith }
12026cec326SBarry Smith 
121a57c2057SPierre Jolivet static PetscErrorCode MatZeroEntries_MPI_Hash(PETSC_UNUSED Mat A)
12226cec326SBarry Smith {
12326cec326SBarry Smith   PetscFunctionBegin;
12426cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
12526cec326SBarry Smith }
12626cec326SBarry Smith 
127a57c2057SPierre Jolivet static PetscErrorCode MatSetRandom_MPI_Hash(Mat A, PETSC_UNUSED PetscRandom r)
12826cec326SBarry Smith {
12926cec326SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Must set preallocation first");
13026cec326SBarry Smith }
13126cec326SBarry Smith 
13226cec326SBarry Smith static PetscErrorCode MatSetUp_MPI_Hash(Mat A)
13326cec326SBarry Smith {
13426cec326SBarry Smith   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
13526cec326SBarry Smith   PetscMPIInt size;
13626cec326SBarry Smith #if !defined(TYPE_AIJ)
13726cec326SBarry Smith   PetscInt bs;
13826cec326SBarry Smith #endif
13926cec326SBarry Smith 
14026cec326SBarry Smith   PetscFunctionBegin;
14126cec326SBarry Smith   PetscCall(PetscInfo(A, "Using hash-based MatSetValues() for MATMPISBAIJ because no preallocation provided\n"));
14226cec326SBarry Smith   PetscCall(PetscLayoutSetUp(A->rmap));
14326cec326SBarry Smith   PetscCall(PetscLayoutSetUp(A->cmap));
14426cec326SBarry Smith   if (A->rmap->bs < 1) A->rmap->bs = 1;
14526cec326SBarry Smith   if (A->cmap->bs < 1) A->cmap->bs = 1;
14626cec326SBarry Smith   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
14726cec326SBarry Smith 
14826cec326SBarry Smith #if !defined(TYPE_AIJ)
14926cec326SBarry Smith   PetscCall(MatGetBlockSize(A, &bs));
15026cec326SBarry Smith   /* these values are set in MatMPISBAIJSetPreallocation() */
15126cec326SBarry Smith   a->bs2 = bs * bs;
15226cec326SBarry Smith   a->mbs = A->rmap->n / bs;
15326cec326SBarry Smith   a->nbs = A->cmap->n / bs;
15426cec326SBarry Smith   a->Mbs = A->rmap->N / bs;
15526cec326SBarry Smith   a->Nbs = A->cmap->N / bs;
15626cec326SBarry Smith 
15726cec326SBarry Smith   for (PetscInt i = 0; i <= a->size; i++) a->rangebs[i] = A->rmap->range[i] / bs;
15826cec326SBarry Smith   a->rstartbs = A->rmap->rstart / bs;
15926cec326SBarry Smith   a->rendbs   = A->rmap->rend / bs;
16026cec326SBarry Smith   a->cstartbs = A->cmap->rstart / bs;
16126cec326SBarry Smith   a->cendbs   = A->cmap->rend / bs;
16226cec326SBarry Smith   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)A), A->rmap->bs, &A->bstash));
16326cec326SBarry Smith #endif
16426cec326SBarry Smith 
16526cec326SBarry Smith   PetscCall(MatCreate(PETSC_COMM_SELF, &a->A));
16626cec326SBarry Smith   PetscCall(MatSetSizes(a->A, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n));
16726cec326SBarry Smith   PetscCall(MatSetBlockSizesFromMats(a->A, A, A));
16826cec326SBarry Smith #if defined(SUB_TYPE_CUSPARSE)
16926cec326SBarry Smith   PetscCall(MatSetType(a->A, MATSEQAIJCUSPARSE));
17026cec326SBarry Smith #else
17126cec326SBarry Smith   PetscCall(MatSetType(a->A, PetscConcat(MATSEQ, TYPE)));
17226cec326SBarry Smith #endif
17326cec326SBarry Smith   PetscCall(MatSetUp(a->A));
17426cec326SBarry Smith 
17526cec326SBarry Smith   PetscCall(MatCreate(PETSC_COMM_SELF, &a->B));
17626cec326SBarry Smith   PetscCall(MatSetSizes(a->B, A->rmap->n, size > 1 ? A->cmap->N : 0, A->rmap->n, size > 1 ? A->cmap->N : 0));
17726cec326SBarry Smith   PetscCall(MatSetBlockSizesFromMats(a->B, A, A));
17826cec326SBarry Smith #if defined(TYPE_SBAIJ)
17926cec326SBarry Smith   PetscCall(MatSetType(a->B, MATSEQBAIJ));
18026cec326SBarry Smith #else
18126cec326SBarry Smith   #if defined(SUB_TYPE_CUSPARSE)
18226cec326SBarry Smith   PetscCall(MatSetType(a->B, MATSEQAIJCUSPARSE));
18326cec326SBarry Smith   #else
18426cec326SBarry Smith   PetscCall(MatSetType(a->B, PetscConcat(MATSEQ, TYPE)));
18526cec326SBarry Smith   #endif
18626cec326SBarry Smith #endif
18726cec326SBarry Smith   PetscCall(MatSetUp(a->B));
18826cec326SBarry Smith 
18926cec326SBarry Smith   /* keep a record of the operations so they can be reset when the hash handling is complete */
190aea10558SJacob Faibussowitsch   a->cops                  = A->ops[0];
19126cec326SBarry Smith   A->ops->assemblybegin    = MatAssemblyBegin_MPI_Hash;
19226cec326SBarry Smith   A->ops->assemblyend      = MatAssemblyEnd_MPI_Hash;
19326cec326SBarry Smith   A->ops->setvalues        = MatSetValues_MPI_Hash;
19426cec326SBarry Smith   A->ops->destroy          = MatDestroy_MPI_Hash;
19526cec326SBarry Smith   A->ops->zeroentries      = MatZeroEntries_MPI_Hash;
19626cec326SBarry Smith   A->ops->setrandom        = MatSetRandom_MPI_Hash;
19726cec326SBarry Smith   A->ops->setvaluesblocked = NULL;
19826cec326SBarry Smith 
19926cec326SBarry Smith   A->preallocated = PETSC_TRUE;
200ad79cf63SBarry Smith   A->hash_active  = PETSC_TRUE;
20126cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
20226cec326SBarry Smith }
203