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