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