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