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