xref: /petsc/src/mat/impls/aij/mpi/mpihashmat.h (revision 5ff6d247539c86491dc822dc7e845e819e6cc4a3)
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 
MatSetValues_MPI_Hash(Mat A,PetscInt m,const PetscInt * rows,PetscInt n,const PetscInt * cols,const PetscScalar * values,InsertMode addv)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;
129680b7ceSStefano Zampini   const PetscInt rStart         = A->rmap->rstart;
139680b7ceSStefano Zampini   const PetscInt rEnd           = A->rmap->rend;
149680b7ceSStefano Zampini   const PetscInt cStart         = A->cmap->rstart;
159680b7ceSStefano Zampini   const PetscInt cEnd           = A->cmap->rend;
1626cec326SBarry Smith #if defined(TYPE_SBAIJ)
179680b7ceSStefano Zampini   const PetscInt bs = A->rmap->bs;
1826cec326SBarry Smith #endif
199680b7ceSStefano Zampini   const PetscBool ignorezeroentries = ((Mat_SeqAIJ *)a->A->data)->ignorezeroentries;
2026cec326SBarry Smith 
2126cec326SBarry Smith   PetscFunctionBegin;
2226cec326SBarry Smith   for (PetscInt r = 0; r < m; ++r) {
2326cec326SBarry Smith     PetscScalar value;
2426cec326SBarry Smith     if (rows[r] < 0) continue;
2526cec326SBarry Smith     if (rows[r] < rStart || rows[r] >= rEnd) {
263ac1eb50SJunchao 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]);
273ac1eb50SJunchao Zhang       if (!a->donotstash) {
283ac1eb50SJunchao Zhang         A->assembled = PETSC_FALSE;
2926cec326SBarry Smith         if (a->roworiented) {
309680b7ceSStefano Zampini           PetscCall(MatStashValuesRow_Private(&A->stash, rows[r], n, cols, PetscSafePointerPlusOffset(values, r * n), (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
3126cec326SBarry Smith         } else {
329680b7ceSStefano Zampini           PetscCall(MatStashValuesCol_Private(&A->stash, rows[r], n, cols, PetscSafePointerPlusOffset(values, r), m, (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
3326cec326SBarry Smith         }
343ac1eb50SJunchao Zhang       }
3526cec326SBarry Smith     } else {
3626cec326SBarry Smith       for (PetscInt c = 0; c < n; ++c) {
3726cec326SBarry Smith #if defined(TYPE_SBAIJ)
3826cec326SBarry Smith         if (cols[c] / bs < rows[r] / bs) continue;
3926cec326SBarry Smith #else
4026cec326SBarry Smith         if (cols[c] < 0) continue;
4126cec326SBarry Smith #endif
42f4f49eeaSPierre Jolivet         value = values ? (a->roworiented ? values[r * n + c] : values[r + m * c]) : 0;
439680b7ceSStefano Zampini         if (ignorezeroentries && value == 0.0 && (addv == ADD_VALUES) && rows[r] != cols[c]) continue;
449680b7ceSStefano Zampini         if (cols[c] >= cStart && cols[c] < cEnd) {
459680b7ceSStefano Zampini           PetscCall(MatSetValue(a->A, rows[r] - rStart, cols[c] - cStart, value, addv));
469680b7ceSStefano Zampini         } else if (!ignorezeroentries || value != 0.0) {
479680b7ceSStefano Zampini           /* MPIAIJ never inserts in B if it is 0 */
489680b7ceSStefano Zampini           PetscCall(MatSetValue(a->B, rows[r] - rStart, cols[c], value, addv));
499680b7ceSStefano Zampini         }
5026cec326SBarry Smith       }
5126cec326SBarry Smith     }
5226cec326SBarry Smith   }
5326cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
5426cec326SBarry Smith }
5526cec326SBarry Smith 
MatAssemblyBegin_MPI_Hash(Mat A,PETSC_UNUSED MatAssemblyType type)56a57c2057SPierre Jolivet static PetscErrorCode MatAssemblyBegin_MPI_Hash(Mat A, PETSC_UNUSED MatAssemblyType type)
5726cec326SBarry Smith {
58cfe56a0cSJunchao Zhang   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
5926cec326SBarry Smith   PetscInt nstash, reallocs;
6026cec326SBarry Smith 
6126cec326SBarry Smith   PetscFunctionBegin;
62cfe56a0cSJunchao Zhang   if (a->donotstash || A->nooffprocentries) PetscFunctionReturn(PETSC_SUCCESS);
6326cec326SBarry Smith   PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range));
6426cec326SBarry Smith   PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs));
6526cec326SBarry Smith   PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
6626cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
6726cec326SBarry Smith }
6826cec326SBarry Smith 
MatFinishScatterAndSetValues_MPI_Hash(Mat A)69fe1fc275SAlexander static PetscErrorCode MatFinishScatterAndSetValues_MPI_Hash(Mat A)
7026cec326SBarry Smith {
7126cec326SBarry Smith   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
7226cec326SBarry Smith   PetscMPIInt  n;
7326cec326SBarry Smith   PetscScalar *val;
7426cec326SBarry Smith   PetscInt    *row, *col;
7526cec326SBarry Smith   PetscInt     j, ncols, flg, rstart;
7626cec326SBarry Smith 
7726cec326SBarry Smith   PetscFunctionBegin;
78cfe56a0cSJunchao Zhang   if (!a->donotstash && !A->nooffprocentries) {
7926cec326SBarry Smith     while (1) {
8026cec326SBarry Smith       PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg));
8126cec326SBarry Smith       if (!flg) break;
8226cec326SBarry Smith 
8326cec326SBarry Smith       for (PetscInt i = 0; i < n;) {
8426cec326SBarry Smith         /* Now identify the consecutive vals belonging to the same row */
8526cec326SBarry Smith         for (j = i, rstart = row[j]; j < n; j++) {
8626cec326SBarry Smith           if (row[j] != rstart) break;
8726cec326SBarry Smith         }
8826cec326SBarry Smith         if (j < n) ncols = j - i;
8926cec326SBarry Smith         else ncols = n - i;
9026cec326SBarry Smith         /* Now assemble all these values with a single function call */
9126cec326SBarry Smith         PetscCall(MatSetValues_MPI_Hash(A, 1, row + i, ncols, col + i, val + i, A->insertmode));
9226cec326SBarry Smith         i = j;
9326cec326SBarry Smith       }
9426cec326SBarry Smith     }
9526cec326SBarry Smith     PetscCall(MatStashScatterEnd_Private(&A->stash));
96cfe56a0cSJunchao Zhang   }
97fe1fc275SAlexander   PetscFunctionReturn(PETSC_SUCCESS);
98fe1fc275SAlexander }
99fe1fc275SAlexander 
MatAssemblyEnd_MPI_Hash(Mat A,MatAssemblyType type)100fe1fc275SAlexander static PetscErrorCode MatAssemblyEnd_MPI_Hash(Mat A, MatAssemblyType type)
101fe1fc275SAlexander {
102fe1fc275SAlexander   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
103fe1fc275SAlexander 
104fe1fc275SAlexander   PetscFunctionBegin;
105fe1fc275SAlexander   PetscCall(MatFinishScatterAndSetValues_MPI_Hash(A));
10626cec326SBarry Smith   if (type != MAT_FINAL_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
10726cec326SBarry Smith 
10826cec326SBarry Smith   A->insertmode = NOT_SET_VALUES; /* this was set by the previous calls to MatSetValues() */
10926cec326SBarry Smith 
110aea10558SJacob Faibussowitsch   A->ops[0]      = a->cops;
111ad79cf63SBarry Smith   A->hash_active = PETSC_FALSE;
11226cec326SBarry Smith 
11326cec326SBarry Smith   PetscCall(MatAssemblyBegin(a->A, MAT_FINAL_ASSEMBLY));
11426cec326SBarry Smith   PetscCall(MatAssemblyEnd(a->A, MAT_FINAL_ASSEMBLY));
11526cec326SBarry Smith   PetscCall(MatAssemblyBegin(a->B, MAT_FINAL_ASSEMBLY));
11626cec326SBarry Smith   PetscCall(MatAssemblyEnd(a->B, MAT_FINAL_ASSEMBLY));
11726cec326SBarry Smith   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
11826cec326SBarry Smith   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
11926cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
12026cec326SBarry Smith }
12126cec326SBarry Smith 
MatCopyHashToXAIJ_MPI_Hash(Mat A,Mat B)122fe1fc275SAlexander static PetscErrorCode MatCopyHashToXAIJ_MPI_Hash(Mat A, Mat B)
123fe1fc275SAlexander {
124fe1fc275SAlexander   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data, *b = (PetscConcat(Mat_MPI, TYPE) *)B->data;
125fe1fc275SAlexander 
126fe1fc275SAlexander   PetscFunctionBegin;
127fe1fc275SAlexander   /* Let's figure there's no harm done in doing the scatters for A now even if A != B */
128fe1fc275SAlexander   PetscCall(MatAssemblyBegin_MPI_Hash(A, /*unused*/ MAT_FINAL_ASSEMBLY));
129fe1fc275SAlexander   PetscCall(MatFinishScatterAndSetValues_MPI_Hash(A));
130fe1fc275SAlexander 
131fe1fc275SAlexander   PetscCall(MatCopyHashToXAIJ(a->A, b->A));
132fe1fc275SAlexander   PetscCall(MatCopyHashToXAIJ(a->B, b->B));
133fe1fc275SAlexander   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
134fe1fc275SAlexander   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
135fe1fc275SAlexander   PetscFunctionReturn(PETSC_SUCCESS);
136fe1fc275SAlexander }
137fe1fc275SAlexander 
MatDestroy_MPI_Hash(Mat A)13826cec326SBarry Smith static PetscErrorCode MatDestroy_MPI_Hash(Mat A)
13926cec326SBarry Smith {
14026cec326SBarry Smith   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
14126cec326SBarry Smith 
14226cec326SBarry Smith   PetscFunctionBegin;
14326cec326SBarry Smith   PetscCall(MatStashDestroy_Private(&A->stash));
14426cec326SBarry Smith   PetscCall(MatDestroy(&a->A));
14526cec326SBarry Smith   PetscCall(MatDestroy(&a->B));
14626cec326SBarry Smith   PetscCall((*a->cops.destroy)(A));
14726cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
14826cec326SBarry Smith }
14926cec326SBarry Smith 
MatZeroEntries_MPI_Hash(PETSC_UNUSED Mat A)150a57c2057SPierre Jolivet static PetscErrorCode MatZeroEntries_MPI_Hash(PETSC_UNUSED Mat A)
15126cec326SBarry Smith {
15226cec326SBarry Smith   PetscFunctionBegin;
15326cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
15426cec326SBarry Smith }
15526cec326SBarry Smith 
MatSetRandom_MPI_Hash(Mat A,PETSC_UNUSED PetscRandom r)156a57c2057SPierre Jolivet static PetscErrorCode MatSetRandom_MPI_Hash(Mat A, PETSC_UNUSED PetscRandom r)
15726cec326SBarry Smith {
15826cec326SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Must set preallocation first");
15926cec326SBarry Smith }
16026cec326SBarry Smith 
MatSetUp_MPI_Hash(Mat A)16126cec326SBarry Smith static PetscErrorCode MatSetUp_MPI_Hash(Mat A)
16226cec326SBarry Smith {
16326cec326SBarry Smith   PetscConcat(Mat_MPI, TYPE) *a = (PetscConcat(Mat_MPI, TYPE) *)A->data;
16426cec326SBarry Smith   PetscMPIInt size;
16526cec326SBarry Smith #if !defined(TYPE_AIJ)
16626cec326SBarry Smith   PetscInt bs;
16726cec326SBarry Smith #endif
16826cec326SBarry Smith 
16926cec326SBarry Smith   PetscFunctionBegin;
170*666f9b3fSPierre Jolivet   PetscCall(PetscInfo(A, "Using hash-based MatSetValues() for MATMPI" PetscStringize(TYPE) " because no preallocation provided\n"));
17126cec326SBarry Smith   PetscCall(PetscLayoutSetUp(A->rmap));
17226cec326SBarry Smith   PetscCall(PetscLayoutSetUp(A->cmap));
17326cec326SBarry Smith   if (A->rmap->bs < 1) A->rmap->bs = 1;
17426cec326SBarry Smith   if (A->cmap->bs < 1) A->cmap->bs = 1;
17526cec326SBarry Smith   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
17626cec326SBarry Smith 
17726cec326SBarry Smith #if !defined(TYPE_AIJ)
17826cec326SBarry Smith   PetscCall(MatGetBlockSize(A, &bs));
17926cec326SBarry Smith   /* these values are set in MatMPISBAIJSetPreallocation() */
18026cec326SBarry Smith   a->bs2 = bs * bs;
18126cec326SBarry Smith   a->mbs = A->rmap->n / bs;
18226cec326SBarry Smith   a->nbs = A->cmap->n / bs;
18326cec326SBarry Smith   a->Mbs = A->rmap->N / bs;
18426cec326SBarry Smith   a->Nbs = A->cmap->N / bs;
18526cec326SBarry Smith 
18626cec326SBarry Smith   for (PetscInt i = 0; i <= a->size; i++) a->rangebs[i] = A->rmap->range[i] / bs;
18726cec326SBarry Smith   a->rstartbs = A->rmap->rstart / bs;
18826cec326SBarry Smith   a->rendbs   = A->rmap->rend / bs;
18926cec326SBarry Smith   a->cstartbs = A->cmap->rstart / bs;
19026cec326SBarry Smith   a->cendbs   = A->cmap->rend / bs;
19126cec326SBarry Smith   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)A), A->rmap->bs, &A->bstash));
19226cec326SBarry Smith #endif
19326cec326SBarry Smith 
19426cec326SBarry Smith   PetscCall(MatCreate(PETSC_COMM_SELF, &a->A));
19526cec326SBarry Smith   PetscCall(MatSetSizes(a->A, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n));
19626cec326SBarry Smith   PetscCall(MatSetBlockSizesFromMats(a->A, A, A));
19726cec326SBarry Smith #if defined(SUB_TYPE_CUSPARSE)
19826cec326SBarry Smith   PetscCall(MatSetType(a->A, MATSEQAIJCUSPARSE));
19926cec326SBarry Smith #else
20026cec326SBarry Smith   PetscCall(MatSetType(a->A, PetscConcat(MATSEQ, TYPE)));
20126cec326SBarry Smith #endif
20226cec326SBarry Smith   PetscCall(MatSetUp(a->A));
20326cec326SBarry Smith 
20426cec326SBarry Smith   PetscCall(MatCreate(PETSC_COMM_SELF, &a->B));
20526cec326SBarry Smith   PetscCall(MatSetSizes(a->B, A->rmap->n, size > 1 ? A->cmap->N : 0, A->rmap->n, size > 1 ? A->cmap->N : 0));
20626cec326SBarry Smith   PetscCall(MatSetBlockSizesFromMats(a->B, A, A));
20726cec326SBarry Smith #if defined(TYPE_SBAIJ)
20826cec326SBarry Smith   PetscCall(MatSetType(a->B, MATSEQBAIJ));
20926cec326SBarry Smith #else
21026cec326SBarry Smith   #if defined(SUB_TYPE_CUSPARSE)
21126cec326SBarry Smith   PetscCall(MatSetType(a->B, MATSEQAIJCUSPARSE));
21226cec326SBarry Smith   #else
21326cec326SBarry Smith   PetscCall(MatSetType(a->B, PetscConcat(MATSEQ, TYPE)));
21426cec326SBarry Smith   #endif
21526cec326SBarry Smith #endif
21626cec326SBarry Smith   PetscCall(MatSetUp(a->B));
21726cec326SBarry Smith 
21826cec326SBarry Smith   /* keep a record of the operations so they can be reset when the hash handling is complete */
219aea10558SJacob Faibussowitsch   a->cops                  = A->ops[0];
22026cec326SBarry Smith   A->ops->assemblybegin    = MatAssemblyBegin_MPI_Hash;
22126cec326SBarry Smith   A->ops->assemblyend      = MatAssemblyEnd_MPI_Hash;
22226cec326SBarry Smith   A->ops->setvalues        = MatSetValues_MPI_Hash;
22326cec326SBarry Smith   A->ops->destroy          = MatDestroy_MPI_Hash;
22426cec326SBarry Smith   A->ops->zeroentries      = MatZeroEntries_MPI_Hash;
22526cec326SBarry Smith   A->ops->setrandom        = MatSetRandom_MPI_Hash;
226fe1fc275SAlexander   A->ops->copyhashtoxaij   = MatCopyHashToXAIJ_MPI_Hash;
22726cec326SBarry Smith   A->ops->setvaluesblocked = NULL;
22826cec326SBarry Smith 
22926cec326SBarry Smith   A->preallocated = PETSC_TRUE;
230ad79cf63SBarry Smith   A->hash_active  = PETSC_TRUE;
23126cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
23226cec326SBarry Smith }
233