xref: /petsc/src/mat/impls/aij/seq/seqhashmat.h (revision 9f0612e409f6220a780be6348417bea34ef34962)
126cec326SBarry Smith /*
226cec326SBarry Smith    used by SEQAIJ, BAIJ and SBAIJ to reduce code duplication
326cec326SBarry Smith 
426cec326SBarry Smith      define TYPE to AIJ BAIJ or SBAIJ
526cec326SBarry Smith             TYPE_BS_ON for BAIJ and SBAIJ
626cec326SBarry Smith 
726cec326SBarry Smith */
826cec326SBarry Smith static PetscErrorCode MatAssemblyEnd_Seq_Hash(Mat A, MatAssemblyType type)
926cec326SBarry Smith {
1026cec326SBarry Smith   PetscConcat(Mat_Seq, TYPE) *a = (PetscConcat(Mat_Seq, TYPE) *)A->data;
1126cec326SBarry Smith   PetscHashIter  hi;
1226cec326SBarry Smith   PetscHashIJKey key;
1326cec326SBarry Smith   PetscScalar    value, *values;
1426cec326SBarry Smith   PetscInt       m, n, *cols, *rowstarts;
1526cec326SBarry Smith #if defined(TYPE_BS_ON)
1626cec326SBarry Smith   PetscInt bs;
1726cec326SBarry Smith #endif
1826cec326SBarry Smith 
1926cec326SBarry Smith   PetscFunctionBegin;
2026cec326SBarry Smith #if defined(TYPE_BS_ON)
2126cec326SBarry Smith   PetscCall(MatGetBlockSize(A, &bs));
2226cec326SBarry Smith   if (bs > 1) PetscCall(PetscHSetIJDestroy(&a->bht));
2326cec326SBarry Smith #endif
2426cec326SBarry Smith   A->preallocated = PETSC_FALSE; /* this was set to true for the MatSetValues_Hash() to work */
2526cec326SBarry Smith 
26aea10558SJacob Faibussowitsch   A->ops[0]      = a->cops;
27ad79cf63SBarry Smith   A->hash_active = PETSC_FALSE;
2826cec326SBarry Smith 
2926cec326SBarry Smith   /* move values from hash format to matrix type format */
3026cec326SBarry Smith   PetscCall(MatGetSize(A, &m, NULL));
3126cec326SBarry Smith #if defined(TYPE_BS_ON)
3226cec326SBarry Smith   if (bs > 1) PetscCall(PetscConcat(PetscConcat(MatSeq, TYPE), SetPreallocation)(A, bs, PETSC_DETERMINE, a->bdnz));
3326cec326SBarry Smith   else PetscCall(PetscConcat(PetscConcat(MatSeq, TYPE), SetPreallocation)(A, 1, PETSC_DETERMINE, a->dnz));
3426cec326SBarry Smith #else
3526cec326SBarry Smith   PetscCall(MatSeqAIJSetPreallocation(A, PETSC_DETERMINE, a->dnz));
3626cec326SBarry Smith #endif
3726cec326SBarry Smith   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
3826cec326SBarry Smith   PetscCall(PetscHMapIJVGetSize(a->ht, &n));
39*9f0612e4SBarry Smith   /* do not need PetscShmgetAllocateArray() since arrays are temporary */
4026cec326SBarry Smith   PetscCall(PetscMalloc3(n, &cols, m + 1, &rowstarts, n, &values));
4126cec326SBarry Smith   rowstarts[0] = 0;
4226cec326SBarry Smith   for (PetscInt i = 0; i < m; i++) rowstarts[i + 1] = rowstarts[i] + a->dnz[i];
4326cec326SBarry Smith 
4426cec326SBarry Smith   PetscHashIterBegin(a->ht, hi);
4526cec326SBarry Smith   while (!PetscHashIterAtEnd(a->ht, hi)) {
4626cec326SBarry Smith     PetscHashIterGetKey(a->ht, hi, key);
4726cec326SBarry Smith     PetscHashIterGetVal(a->ht, hi, value);
4826cec326SBarry Smith     cols[rowstarts[key.i]]     = key.j;
4926cec326SBarry Smith     values[rowstarts[key.i]++] = value;
5026cec326SBarry Smith     PetscHashIterNext(a->ht, hi);
5126cec326SBarry Smith   }
5226cec326SBarry Smith   PetscCall(PetscHMapIJVDestroy(&a->ht));
5326cec326SBarry Smith 
5426cec326SBarry Smith   for (PetscInt i = 0, start = 0; i < m; i++) {
558e3a54c0SPierre Jolivet     PetscCall(MatSetValues(A, 1, &i, a->dnz[i], PetscSafePointerPlusOffset(cols, start), PetscSafePointerPlusOffset(values, start), A->insertmode));
5626cec326SBarry Smith     start += a->dnz[i];
5726cec326SBarry Smith   }
5826cec326SBarry Smith   PetscCall(PetscFree3(cols, rowstarts, values));
5926cec326SBarry Smith   PetscCall(PetscFree(a->dnz));
6026cec326SBarry Smith #if defined(TYPE_BS_ON)
6126cec326SBarry Smith   if (bs > 1) PetscCall(PetscFree(a->bdnz));
6226cec326SBarry Smith #endif
6326cec326SBarry Smith   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
6426cec326SBarry Smith   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
6526cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
6626cec326SBarry Smith }
6726cec326SBarry Smith 
6826cec326SBarry Smith static PetscErrorCode MatDestroy_Seq_Hash(Mat A)
6926cec326SBarry Smith {
7026cec326SBarry Smith   PetscConcat(Mat_Seq, TYPE) *a = (PetscConcat(Mat_Seq, TYPE) *)A->data;
7126cec326SBarry Smith #if defined(TYPE_BS_ON)
7226cec326SBarry Smith   PetscInt bs;
7326cec326SBarry Smith #endif
7426cec326SBarry Smith 
7526cec326SBarry Smith   PetscFunctionBegin;
7626cec326SBarry Smith   PetscCall(PetscHMapIJVDestroy(&a->ht));
7726cec326SBarry Smith   PetscCall(PetscFree(a->dnz));
7826cec326SBarry Smith #if defined(TYPE_BS_ON)
7926cec326SBarry Smith   PetscCall(MatGetBlockSize(A, &bs));
8026cec326SBarry Smith   if (bs > 1) {
8126cec326SBarry Smith     PetscCall(PetscFree(a->bdnz));
8226cec326SBarry Smith     PetscCall(PetscHSetIJDestroy(&a->bht));
8326cec326SBarry Smith   }
8426cec326SBarry Smith #endif
8526cec326SBarry Smith   PetscCall((*a->cops.destroy)(A));
8626cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
8726cec326SBarry Smith }
8826cec326SBarry Smith 
8926cec326SBarry Smith static PetscErrorCode MatZeroEntries_Seq_Hash(Mat A)
9026cec326SBarry Smith {
9126cec326SBarry Smith   PetscFunctionBegin;
9226cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
9326cec326SBarry Smith }
9426cec326SBarry Smith 
9526cec326SBarry Smith static PetscErrorCode MatSetRandom_Seq_Hash(Mat A, PetscRandom r)
9626cec326SBarry Smith {
9726cec326SBarry Smith   PetscFunctionBegin;
9826cec326SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Must set preallocation first");
9926cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
10026cec326SBarry Smith }
10126cec326SBarry Smith 
10226cec326SBarry Smith static PetscErrorCode MatSetUp_Seq_Hash(Mat A)
10326cec326SBarry Smith {
10426cec326SBarry Smith   PetscConcat(Mat_Seq, TYPE) *a = (PetscConcat(Mat_Seq, TYPE) *)A->data;
10526cec326SBarry Smith   PetscInt m;
10626cec326SBarry Smith #if defined(TYPE_BS_ON)
10726cec326SBarry Smith   PetscInt bs;
10826cec326SBarry Smith #endif
10926cec326SBarry Smith 
11026cec326SBarry Smith   PetscFunctionBegin;
11126cec326SBarry Smith   PetscCall(PetscInfo(A, "Using hash-based MatSetValues() for MATSEQAIJ because no preallocation provided\n"));
11226cec326SBarry Smith   PetscCall(PetscLayoutSetUp(A->rmap));
11326cec326SBarry Smith   PetscCall(PetscLayoutSetUp(A->cmap));
11426cec326SBarry Smith   if (A->rmap->bs < 1) A->rmap->bs = 1;
11526cec326SBarry Smith   if (A->cmap->bs < 1) A->cmap->bs = 1;
11626cec326SBarry Smith 
11726cec326SBarry Smith   PetscCall(MatGetLocalSize(A, &m, NULL));
11826cec326SBarry Smith   PetscCall(PetscHMapIJVCreate(&a->ht));
11926cec326SBarry Smith   PetscCall(PetscCalloc1(m, &a->dnz));
12026cec326SBarry Smith #if defined(TYPE_BS_ON)
12126cec326SBarry Smith   PetscCall(MatGetBlockSize(A, &bs));
12226cec326SBarry Smith   if (bs > 1) {
12326cec326SBarry Smith     PetscCall(PetscCalloc1(m / bs, &a->bdnz));
12426cec326SBarry Smith     PetscCall(PetscHSetIJCreate(&a->bht));
12526cec326SBarry Smith   }
12626cec326SBarry Smith #endif
12726cec326SBarry Smith 
12826cec326SBarry Smith   /* keep a record of the operations so they can be reset when the hash handling is complete */
129aea10558SJacob Faibussowitsch   a->cops               = A->ops[0];
13026cec326SBarry Smith   A->ops->assemblybegin = NULL;
13126cec326SBarry Smith   A->ops->assemblyend   = MatAssemblyEnd_Seq_Hash;
13226cec326SBarry Smith   A->ops->destroy       = MatDestroy_Seq_Hash;
13326cec326SBarry Smith   A->ops->zeroentries   = MatZeroEntries_Seq_Hash;
13426cec326SBarry Smith   A->ops->setrandom     = MatSetRandom_Seq_Hash;
13526cec326SBarry Smith #if defined(TYPE_BS_ON)
13626cec326SBarry Smith   if (bs > 1) A->ops->setvalues = MatSetValues_Seq_Hash_BS;
13726cec326SBarry Smith   else
13826cec326SBarry Smith #endif
13926cec326SBarry Smith     A->ops->setvalues = MatSetValues_Seq_Hash;
14026cec326SBarry Smith   A->ops->setvaluesblocked = NULL;
14126cec326SBarry Smith 
14226cec326SBarry Smith   A->preallocated = PETSC_TRUE;
143ad79cf63SBarry Smith   A->hash_active  = PETSC_TRUE;
14426cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
14526cec326SBarry Smith }
146