xref: /petsc/src/mat/impls/aij/seq/aijsell/aijsell.c (revision e0f5bfbec699682fa3e8b8532b1176849ea4e12a)
1 /*
2   Defines basic operations for the MATSEQAIJSELL matrix class.
3   This class is derived from the MATAIJCLASS, but maintains a "shadow" copy
4   of the matrix stored in MATSEQSELL format, which is used as appropriate for
5   performing operations for which this format is more suitable.
6 */
7 
8 #include <../src/mat/impls/aij/seq/aij.h>
9 #include <../src/mat/impls/sell/seq/sell.h>
10 
11 typedef struct {
12   Mat              S; /* The SELL formatted "shadow" matrix. */
13   PetscBool        eager_shadow;
14   PetscObjectState state; /* State of the matrix when shadow matrix was last constructed. */
15 } Mat_SeqAIJSELL;
16 
17 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJSELL_SeqAIJ(Mat A, MatType type, MatReuse reuse, Mat *newmat) {
18   /* This routine is only called to convert a MATAIJSELL to its base PETSc type, */
19   /* so we will ignore 'MatType type'. */
20   Mat             B       = *newmat;
21   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;
22 
23   PetscFunctionBegin;
24   if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
25 
26   /* Reset the original function pointers. */
27   B->ops->duplicate        = MatDuplicate_SeqAIJ;
28   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJ;
29   B->ops->destroy          = MatDestroy_SeqAIJ;
30   B->ops->mult             = MatMult_SeqAIJ;
31   B->ops->multtranspose    = MatMultTranspose_SeqAIJ;
32   B->ops->multadd          = MatMultAdd_SeqAIJ;
33   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJ;
34   B->ops->sor              = MatSOR_SeqAIJ;
35 
36   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaijsell_seqaij_C", NULL));
37 
38   if (reuse == MAT_INITIAL_MATRIX) aijsell = (Mat_SeqAIJSELL *)B->spptr;
39 
40   /* Clean up the Mat_SeqAIJSELL data structure.
41    * Note that MatDestroy() simply returns if passed a NULL value, so it's OK to call even if the shadow matrix was never constructed. */
42   PetscCall(MatDestroy(&aijsell->S));
43   PetscCall(PetscFree(B->spptr));
44 
45   /* Change the type of B to MATSEQAIJ. */
46   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ));
47 
48   *newmat = B;
49   PetscFunctionReturn(0);
50 }
51 
52 PetscErrorCode MatDestroy_SeqAIJSELL(Mat A) {
53   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;
54 
55   PetscFunctionBegin;
56 
57   /* If MatHeaderMerge() was used, then this SeqAIJSELL matrix will not have an
58    * spptr pointer. */
59   if (aijsell) {
60     /* Clean up everything in the Mat_SeqAIJSELL data structure, then free A->spptr. */
61     PetscCall(MatDestroy(&aijsell->S));
62     PetscCall(PetscFree(A->spptr));
63   }
64 
65   /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ()
66    * to destroy everything that remains. */
67   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ));
68   /* Note that I don't call MatSetType().  I believe this is because that
69    * is only to be called when *building* a matrix.  I could be wrong, but
70    * that is how things work for the SuperLU matrix class. */
71   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaijsell_seqaij_C", NULL));
72   PetscCall(MatDestroy_SeqAIJ(A));
73   PetscFunctionReturn(0);
74 }
75 
76 /* Build or update the shadow matrix if and only if needed.
77  * We track the ObjectState to determine when this needs to be done. */
78 PETSC_INTERN PetscErrorCode MatSeqAIJSELL_build_shadow(Mat A) {
79   Mat_SeqAIJSELL  *aijsell = (Mat_SeqAIJSELL *)A->spptr;
80   PetscObjectState state;
81 
82   PetscFunctionBegin;
83   PetscCall(PetscObjectStateGet((PetscObject)A, &state));
84   if (aijsell->S && aijsell->state == state) {
85     /* The existing shadow matrix is up-to-date, so simply exit. */
86     PetscFunctionReturn(0);
87   }
88 
89   PetscCall(PetscLogEventBegin(MAT_Convert, A, 0, 0, 0));
90   if (aijsell->S) {
91     PetscCall(MatConvert_SeqAIJ_SeqSELL(A, MATSEQSELL, MAT_REUSE_MATRIX, &aijsell->S));
92   } else {
93     PetscCall(MatConvert_SeqAIJ_SeqSELL(A, MATSEQSELL, MAT_INITIAL_MATRIX, &aijsell->S));
94   }
95   PetscCall(PetscLogEventEnd(MAT_Convert, A, 0, 0, 0));
96 
97   /* Record the ObjectState so that we can tell when the shadow matrix needs updating */
98   PetscCall(PetscObjectStateGet((PetscObject)A, &aijsell->state));
99 
100   PetscFunctionReturn(0);
101 }
102 
103 PetscErrorCode MatDuplicate_SeqAIJSELL(Mat A, MatDuplicateOption op, Mat *M) {
104   Mat_SeqAIJSELL *aijsell;
105   Mat_SeqAIJSELL *aijsell_dest;
106 
107   PetscFunctionBegin;
108   PetscCall(MatDuplicate_SeqAIJ(A, op, M));
109   aijsell      = (Mat_SeqAIJSELL *)A->spptr;
110   aijsell_dest = (Mat_SeqAIJSELL *)(*M)->spptr;
111   PetscCall(PetscArraycpy(aijsell_dest, aijsell, 1));
112   /* We don't duplicate the shadow matrix -- that will be constructed as needed. */
113   aijsell_dest->S = NULL;
114   if (aijsell->eager_shadow) PetscCall(MatSeqAIJSELL_build_shadow(A));
115   PetscFunctionReturn(0);
116 }
117 
118 PetscErrorCode MatAssemblyEnd_SeqAIJSELL(Mat A, MatAssemblyType mode) {
119   Mat_SeqAIJ     *a       = (Mat_SeqAIJ *)A->data;
120   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;
121 
122   PetscFunctionBegin;
123   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
124 
125   /* I disable the use of the inode routines so that the AIJSELL ones will be
126    * used instead, but I wonder if it might make sense (and is feasible) to
127    * use some of them. */
128   a->inode.use = PETSC_FALSE;
129 
130   /* Since a MATSEQAIJSELL matrix is really just a MATSEQAIJ with some
131    * extra information and some different methods, call the AssemblyEnd
132    * routine for a MATSEQAIJ.
133    * I'm not sure if this is the best way to do this, but it avoids
134    * a lot of code duplication. */
135 
136   PetscCall(MatAssemblyEnd_SeqAIJ(A, mode));
137 
138   /* If the user has requested "eager" shadowing, create the SELL shadow matrix (if needed; the function checks).
139    * (The default is to take a "lazy" approach, deferring this until something like MatMult() is called.) */
140   if (aijsell->eager_shadow) PetscCall(MatSeqAIJSELL_build_shadow(A));
141 
142   PetscFunctionReturn(0);
143 }
144 
145 PetscErrorCode MatMult_SeqAIJSELL(Mat A, Vec xx, Vec yy) {
146   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;
147 
148   PetscFunctionBegin;
149   PetscCall(MatSeqAIJSELL_build_shadow(A));
150   PetscCall(MatMult_SeqSELL(aijsell->S, xx, yy));
151   PetscFunctionReturn(0);
152 }
153 
154 PetscErrorCode MatMultTranspose_SeqAIJSELL(Mat A, Vec xx, Vec yy) {
155   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;
156 
157   PetscFunctionBegin;
158   PetscCall(MatSeqAIJSELL_build_shadow(A));
159   PetscCall(MatMultTranspose_SeqSELL(aijsell->S, xx, yy));
160   PetscFunctionReturn(0);
161 }
162 
163 PetscErrorCode MatMultAdd_SeqAIJSELL(Mat A, Vec xx, Vec yy, Vec zz) {
164   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;
165 
166   PetscFunctionBegin;
167   PetscCall(MatSeqAIJSELL_build_shadow(A));
168   PetscCall(MatMultAdd_SeqSELL(aijsell->S, xx, yy, zz));
169   PetscFunctionReturn(0);
170 }
171 
172 PetscErrorCode MatMultTransposeAdd_SeqAIJSELL(Mat A, Vec xx, Vec yy, Vec zz) {
173   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;
174 
175   PetscFunctionBegin;
176   PetscCall(MatSeqAIJSELL_build_shadow(A));
177   PetscCall(MatMultTransposeAdd_SeqSELL(aijsell->S, xx, yy, zz));
178   PetscFunctionReturn(0);
179 }
180 
181 PetscErrorCode MatSOR_SeqAIJSELL(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) {
182   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;
183 
184   PetscFunctionBegin;
185   PetscCall(MatSeqAIJSELL_build_shadow(A));
186   PetscCall(MatSOR_SeqSELL(aijsell->S, bb, omega, flag, fshift, its, lits, xx));
187   PetscFunctionReturn(0);
188 }
189 
190 /* MatConvert_SeqAIJ_SeqAIJSELL converts a SeqAIJ matrix into a
191  * SeqAIJSELL matrix.  This routine is called by the MatCreate_SeqAIJSELL()
192  * routine, but can also be used to convert an assembled SeqAIJ matrix
193  * into a SeqAIJSELL one. */
194 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat A, MatType type, MatReuse reuse, Mat *newmat) {
195   Mat             B = *newmat;
196   Mat_SeqAIJ     *b;
197   Mat_SeqAIJSELL *aijsell;
198   PetscBool       set;
199   PetscBool       sametype;
200 
201   PetscFunctionBegin;
202   if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
203 
204   PetscCall(PetscObjectTypeCompare((PetscObject)A, type, &sametype));
205   if (sametype) PetscFunctionReturn(0);
206 
207   PetscCall(PetscNew(&aijsell));
208   b        = (Mat_SeqAIJ *)B->data;
209   B->spptr = (void *)aijsell;
210 
211   /* Disable use of the inode routines so that the AIJSELL ones will be used instead.
212    * This happens in MatAssemblyEnd_SeqAIJSELL as well, but the assembly end may not be called, so set it here, too.
213    * As noted elsewhere, I wonder if it might make sense and be feasible to use some of the inode routines. */
214   b->inode.use = PETSC_FALSE;
215 
216   /* Set function pointers for methods that we inherit from AIJ but override.
217    * We also parse some command line options below, since those determine some of the methods we point to. */
218   B->ops->duplicate   = MatDuplicate_SeqAIJSELL;
219   B->ops->assemblyend = MatAssemblyEnd_SeqAIJSELL;
220   B->ops->destroy     = MatDestroy_SeqAIJSELL;
221 
222   aijsell->S            = NULL;
223   aijsell->eager_shadow = PETSC_FALSE;
224 
225   /* Parse command line options. */
226   PetscOptionsBegin(PetscObjectComm((PetscObject)A), ((PetscObject)A)->prefix, "AIJSELL Options", "Mat");
227   PetscCall(PetscOptionsBool("-mat_aijsell_eager_shadow", "Eager Shadowing", "None", (PetscBool)aijsell->eager_shadow, (PetscBool *)&aijsell->eager_shadow, &set));
228   PetscOptionsEnd();
229 
230   /* If A has already been assembled and eager shadowing is specified, build the shadow matrix. */
231   if (A->assembled && aijsell->eager_shadow) PetscCall(MatSeqAIJSELL_build_shadow(A));
232 
233   B->ops->mult             = MatMult_SeqAIJSELL;
234   B->ops->multtranspose    = MatMultTranspose_SeqAIJSELL;
235   B->ops->multadd          = MatMultAdd_SeqAIJSELL;
236   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJSELL;
237   B->ops->sor              = MatSOR_SeqAIJSELL;
238 
239   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaijsell_seqaij_C", MatConvert_SeqAIJSELL_SeqAIJ));
240 
241   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJSELL));
242   *newmat = B;
243   PetscFunctionReturn(0);
244 }
245 
246 /*@C
247    MatCreateSeqAIJSELL - Creates a sparse matrix of type `MATSEQAIJSELL`.
248    This type inherits from AIJ and is largely identical, but keeps a "shadow"
249    copy of the matrix in `MATSEQSELL` format, which is used when this format
250    may be more suitable for a requested operation. Currently, `MATSEQSELL` format
251    is used for `MatMult()`, `MatMultTranspose()`, `MatMultAdd()`, `MatMultTransposeAdd()`,
252    and `MatSOR()` operations.
253 
254    Collective
255 
256    Input Parameters:
257 +  comm - MPI communicator, set to `PETSC_COMM_SELF`
258 .  m - number of rows
259 .  n - number of columns
260 .  nz - number of nonzeros per row (same for all rows)
261 -  nnz - array containing the number of nonzeros in the various rows
262          (possibly different for each row) or NULL
263 
264    Output Parameter:
265 .  A - the matrix
266 
267    Options Database Keys:
268 .  -mat_aijsell_eager_shadow - Construct shadow matrix upon matrix assembly; default is to take a "lazy" approach, performing this step the first time the matrix is applied
269 
270    Notes:
271    If nnz is given then nz is ignored
272 
273    Because `MATSEQAIJSELL is a subtype of `MATSEQAIJ`, the option "-mat_seqaij_type seqaijsell" can be used to make
274    sequential `MATSEAIJ` matrices default to being instances of `MATSEQAIJSELL`.
275 
276    Level: intermediate
277 
278 .seealso: `MatCreate()`, `MatCreateMPIAIJSELL()`, `MatSetValues()`
279 @*/
280 PetscErrorCode MatCreateSeqAIJSELL(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A) {
281   PetscFunctionBegin;
282   PetscCall(MatCreate(comm, A));
283   PetscCall(MatSetSizes(*A, m, n, m, n));
284   PetscCall(MatSetType(*A, MATSEQAIJSELL));
285   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*A, nz, nnz));
286   PetscFunctionReturn(0);
287 }
288 
289 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJSELL(Mat A) {
290   PetscFunctionBegin;
291   PetscCall(MatSetType(A, MATSEQAIJ));
292   PetscCall(MatConvert_SeqAIJ_SeqAIJSELL(A, MATSEQAIJSELL, MAT_INPLACE_MATRIX, &A));
293   PetscFunctionReturn(0);
294 }
295