xref: /petsc/src/mat/impls/aij/seq/seqhashmat.h (revision 65179781167c53baf4bbcd2ba89ddd86e7bb5b3d)
1 /*
2    used by SEQAIJ, BAIJ and SBAIJ to reduce code duplication
3 
4      define TYPE to AIJ BAIJ or SBAIJ
5             TYPE_BS_ON for BAIJ and SBAIJ
6 
7 */
8 static PetscErrorCode MatAssemblyEnd_Seq_Hash(Mat A, MatAssemblyType type)
9 {
10   PetscConcat(Mat_Seq, TYPE) *a = (PetscConcat(Mat_Seq, TYPE) *)A->data;
11   PetscHashIter  hi;
12   PetscHashIJKey key;
13   PetscScalar    value, *values;
14   PetscInt       m, n, *cols, *rowstarts;
15 #if defined(TYPE_BS_ON)
16   PetscInt bs;
17 #endif
18 
19   PetscFunctionBegin;
20 #if defined(TYPE_BS_ON)
21   PetscCall(MatGetBlockSize(A, &bs));
22   if (bs > 1) PetscCall(PetscHSetIJDestroy(&a->bht));
23 #endif
24   A->preallocated = PETSC_FALSE; /* this was set to true for the MatSetValues_Hash() to work */
25 
26   PetscCall(PetscMemcpy(&A->ops, &a->cops, sizeof(*(A->ops))));
27 
28   /* move values from hash format to matrix type format */
29   PetscCall(MatGetSize(A, &m, NULL));
30 #if defined(TYPE_BS_ON)
31   if (bs > 1) PetscCall(PetscConcat(PetscConcat(MatSeq, TYPE), SetPreallocation)(A, bs, PETSC_DETERMINE, a->bdnz));
32   else PetscCall(PetscConcat(PetscConcat(MatSeq, TYPE), SetPreallocation)(A, 1, PETSC_DETERMINE, a->dnz));
33 #else
34   PetscCall(MatSeqAIJSetPreallocation(A, PETSC_DETERMINE, a->dnz));
35 #endif
36   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
37   PetscCall(PetscHMapIJVGetSize(a->ht, &n));
38   PetscCall(PetscMalloc3(n, &cols, m + 1, &rowstarts, n, &values));
39   rowstarts[0] = 0;
40   for (PetscInt i = 0; i < m; i++) rowstarts[i + 1] = rowstarts[i] + a->dnz[i];
41 
42   PetscHashIterBegin(a->ht, hi);
43   while (!PetscHashIterAtEnd(a->ht, hi)) {
44     PetscHashIterGetKey(a->ht, hi, key);
45     PetscHashIterGetVal(a->ht, hi, value);
46     cols[rowstarts[key.i]]     = key.j;
47     values[rowstarts[key.i]++] = value;
48     PetscHashIterNext(a->ht, hi);
49   }
50   PetscCall(PetscHMapIJVDestroy(&a->ht));
51 
52   for (PetscInt i = 0, start = 0; i < m; i++) {
53     PetscCall(MatSetValues(A, 1, &i, a->dnz[i], &cols[start], &values[start], A->insertmode));
54     start += a->dnz[i];
55   }
56   PetscCall(PetscFree3(cols, rowstarts, values));
57   PetscCall(PetscFree(a->dnz));
58 #if defined(TYPE_BS_ON)
59   if (bs > 1) PetscCall(PetscFree(a->bdnz));
60 #endif
61   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
62   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
63   PetscFunctionReturn(PETSC_SUCCESS);
64 }
65 
66 static PetscErrorCode MatDestroy_Seq_Hash(Mat A)
67 {
68   PetscConcat(Mat_Seq, TYPE) *a = (PetscConcat(Mat_Seq, TYPE) *)A->data;
69 #if defined(TYPE_BS_ON)
70   PetscInt bs;
71 #endif
72 
73   PetscFunctionBegin;
74   PetscCall(PetscHMapIJVDestroy(&a->ht));
75   PetscCall(PetscFree(a->dnz));
76 #if defined(TYPE_BS_ON)
77   PetscCall(MatGetBlockSize(A, &bs));
78   if (bs > 1) {
79     PetscCall(PetscFree(a->bdnz));
80     PetscCall(PetscHSetIJDestroy(&a->bht));
81   }
82 #endif
83   PetscCall((*a->cops.destroy)(A));
84   PetscFunctionReturn(PETSC_SUCCESS);
85 }
86 
87 static PetscErrorCode MatZeroEntries_Seq_Hash(Mat A)
88 {
89   PetscFunctionBegin;
90   PetscFunctionReturn(PETSC_SUCCESS);
91 }
92 
93 static PetscErrorCode MatSetRandom_Seq_Hash(Mat A, PetscRandom r)
94 {
95   PetscFunctionBegin;
96   SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Must set preallocation first");
97   PetscFunctionReturn(PETSC_SUCCESS);
98 }
99 
100 static PetscErrorCode MatSetUp_Seq_Hash(Mat A)
101 {
102   PetscConcat(Mat_Seq, TYPE) *a = (PetscConcat(Mat_Seq, TYPE) *)A->data;
103   PetscInt m;
104 #if defined(TYPE_BS_ON)
105   PetscInt bs;
106 #endif
107 
108   PetscFunctionBegin;
109   PetscCall(PetscInfo(A, "Using hash-based MatSetValues() for MATSEQAIJ because no preallocation provided\n"));
110   PetscCall(PetscLayoutSetUp(A->rmap));
111   PetscCall(PetscLayoutSetUp(A->cmap));
112   if (A->rmap->bs < 1) A->rmap->bs = 1;
113   if (A->cmap->bs < 1) A->cmap->bs = 1;
114 
115   PetscCall(MatGetLocalSize(A, &m, NULL));
116   PetscCall(PetscHMapIJVCreate(&a->ht));
117   PetscCall(PetscCalloc1(m, &a->dnz));
118 #if defined(TYPE_BS_ON)
119   PetscCall(MatGetBlockSize(A, &bs));
120   if (bs > 1) {
121     PetscCall(PetscCalloc1(m / bs, &a->bdnz));
122     PetscCall(PetscHSetIJCreate(&a->bht));
123   }
124 #endif
125 
126   /* keep a record of the operations so they can be reset when the hash handling is complete */
127   PetscCall(PetscMemcpy(&a->cops, &A->ops, sizeof(*(A->ops))));
128 
129   A->ops->assemblybegin = NULL;
130   A->ops->assemblyend   = MatAssemblyEnd_Seq_Hash;
131   A->ops->destroy       = MatDestroy_Seq_Hash;
132   A->ops->zeroentries   = MatZeroEntries_Seq_Hash;
133   A->ops->setrandom     = MatSetRandom_Seq_Hash;
134 #if defined(TYPE_BS_ON)
135   if (bs > 1) A->ops->setvalues = MatSetValues_Seq_Hash_BS;
136   else
137 #endif
138     A->ops->setvalues = MatSetValues_Seq_Hash;
139   A->ops->setvaluesblocked = NULL;
140 
141   A->preallocated = PETSC_TRUE;
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144