Actual source code: aijsell.c

  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: */

  8: #include <../src/mat/impls/aij/seq/aij.h>
  9: #include <../src/mat/impls/sell/seq/sell.h>

 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;

 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'. */
 22:   Mat            B        = *newmat;
 23:   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL*) A->spptr;

 26:   if (reuse == MAT_INITIAL_MATRIX) {
 27:     MatDuplicate(A,MAT_COPY_VALUES,&B);
 28:   }

 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;

 40:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijsell_seqaij_C",NULL);

 42:   if (reuse == MAT_INITIAL_MATRIX) aijsell = (Mat_SeqAIJSELL*)B->spptr;

 44:   /* Clean up the Mat_SeqAIJSELL data structure.
 45:    * Note that MatDestroy() simply returns if passed a NULL value, so it's OK to call even if the shadow matrix was never constructed. */
 46:   MatDestroy(&aijsell->S);
 47:   PetscFree(B->spptr);

 49:   /* Change the type of B to MATSEQAIJ. */
 50:   PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ);

 52:   *newmat = B;
 53:   return(0);
 54: }

 56: PetscErrorCode MatDestroy_SeqAIJSELL(Mat A)
 57: {
 59:   Mat_SeqAIJSELL  *aijsell = (Mat_SeqAIJSELL*) A->spptr;


 63:   /* If MatHeaderMerge() was used, then this SeqAIJSELL matrix will not have an
 64:    * spptr pointer. */
 65:   if (aijsell) {
 66:     /* Clean up everything in the Mat_SeqAIJSELL data structure, then free A->spptr. */
 67:     MatDestroy(&aijsell->S);
 68:     PetscFree(A->spptr);
 69:   }

 71:   /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ()
 72:    * to destroy everything that remains. */
 73:   PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ);
 74:   /* Note that I don't call MatSetType().  I believe this is because that
 75:    * is only to be called when *building* a matrix.  I could be wrong, but
 76:    * that is how things work for the SuperLU matrix class. */
 77:   MatDestroy_SeqAIJ(A);
 78:   return(0);
 79: }

 81: /* Build or update the shadow matrix if and only if needed.
 82:  * We track the ObjectState to determine when this needs to be done. */
 83: PETSC_INTERN PetscErrorCode MatSeqAIJSELL_build_shadow(Mat A)
 84: {
 85:   PetscErrorCode   ierr;
 86:   Mat_SeqAIJSELL   *aijsell = (Mat_SeqAIJSELL*) A->spptr;
 87:   PetscObjectState state;

 90:   PetscObjectStateGet((PetscObject)A,&state);
 91:   if (aijsell->S && aijsell->state == state) {
 92:     /* The existing shadow matrix is up-to-date, so simply exit. */
 93:     return(0);
 94:   }

 96:   PetscLogEventBegin(MAT_Convert,A,0,0,0);
 97:   if (aijsell->S) {
 98:     MatConvert_SeqAIJ_SeqSELL(A,MATSEQSELL,MAT_REUSE_MATRIX,&aijsell->S);
 99:   } else {
100:     MatConvert_SeqAIJ_SeqSELL(A,MATSEQSELL,MAT_INITIAL_MATRIX,&aijsell->S);
101:   }
102:   PetscLogEventEnd(MAT_Convert,A,0,0,0);

104:   /* Record the ObjectState so that we can tell when the shadow matrix needs updating */
105:   PetscObjectStateGet((PetscObject)A,&aijsell->state);

107:   return(0);
108: }

110: PetscErrorCode MatDuplicate_SeqAIJSELL(Mat A, MatDuplicateOption op, Mat *M)
111: {
113:   Mat_SeqAIJSELL *aijsell;
114:   Mat_SeqAIJSELL *aijsell_dest;

117:   MatDuplicate_SeqAIJ(A,op,M);
118:   aijsell      = (Mat_SeqAIJSELL*) A->spptr;
119:   aijsell_dest = (Mat_SeqAIJSELL*) (*M)->spptr;
120:   PetscArraycpy(aijsell_dest,aijsell,1);
121:   /* We don't duplicate the shadow matrix -- that will be constructed as needed. */
122:   aijsell_dest->S = NULL;
123:   if (aijsell->eager_shadow) {
124:     MatSeqAIJSELL_build_shadow(A);
125:   }
126:   return(0);
127: }

129: PetscErrorCode MatAssemblyEnd_SeqAIJSELL(Mat A, MatAssemblyType mode)
130: {
131:   PetscErrorCode  ierr;
132:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
133:   Mat_SeqAIJSELL  *aijsell = (Mat_SeqAIJSELL*)A->spptr;

136:   if (mode == MAT_FLUSH_ASSEMBLY) return(0);

138:   /* I disable the use of the inode routines so that the AIJSELL ones will be
139:    * used instead, but I wonder if it might make sense (and is feasible) to
140:    * use some of them. */
141:   a->inode.use = PETSC_FALSE;

143:   /* Since a MATSEQAIJSELL matrix is really just a MATSEQAIJ with some
144:    * extra information and some different methods, call the AssemblyEnd
145:    * routine for a MATSEQAIJ.
146:    * I'm not sure if this is the best way to do this, but it avoids
147:    * a lot of code duplication. */

149:   MatAssemblyEnd_SeqAIJ(A, mode);

151:   /* If the user has requested "eager" shadowing, create the SELL shadow matrix (if needed; the function checks).
152:    * (The default is to take a "lazy" approach, deferring this until something like MatMult() is called.) */
153:   if (aijsell->eager_shadow) {
154:     MatSeqAIJSELL_build_shadow(A);
155:   }

157:   return(0);
158: }

160: PetscErrorCode MatMult_SeqAIJSELL(Mat A,Vec xx,Vec yy)
161: {
162:   Mat_SeqAIJSELL    *aijsell = (Mat_SeqAIJSELL*)A->spptr;
163:   PetscErrorCode    ierr;

166:   MatSeqAIJSELL_build_shadow(A);
167:   MatMult_SeqSELL(aijsell->S,xx,yy);
168:   return(0);
169: }

171: PetscErrorCode MatMultTranspose_SeqAIJSELL(Mat A,Vec xx,Vec yy)
172: {
173:   Mat_SeqAIJSELL    *aijsell=(Mat_SeqAIJSELL*)A->spptr;
174:   PetscErrorCode    ierr;

177:   MatSeqAIJSELL_build_shadow(A);
178:   MatMultTranspose_SeqSELL(aijsell->S,xx,yy);
179:   return(0);
180: }

182: PetscErrorCode MatMultAdd_SeqAIJSELL(Mat A,Vec xx,Vec yy,Vec zz)
183: {
184:   Mat_SeqAIJSELL    *aijsell=(Mat_SeqAIJSELL*)A->spptr;
185:   PetscErrorCode    ierr;

188:   MatSeqAIJSELL_build_shadow(A);
189:   MatMultAdd_SeqSELL(aijsell->S,xx,yy,zz);
190:   return(0);
191: }

193: PetscErrorCode MatMultTransposeAdd_SeqAIJSELL(Mat A,Vec xx,Vec yy,Vec zz)
194: {
195:   Mat_SeqAIJSELL    *aijsell=(Mat_SeqAIJSELL*)A->spptr;
196:   PetscErrorCode    ierr;

199:   MatSeqAIJSELL_build_shadow(A);
200:   MatMultTransposeAdd_SeqSELL(aijsell->S,xx,yy,zz);
201:   return(0);
202: }

204: PetscErrorCode MatSOR_SeqAIJSELL(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
205: {
206:   Mat_SeqAIJSELL    *aijsell=(Mat_SeqAIJSELL*)A->spptr;
207:   PetscErrorCode    ierr;

210:   MatSeqAIJSELL_build_shadow(A);
211:   MatSOR_SeqSELL(aijsell->S,bb,omega,flag,fshift,its,lits,xx);
212:   return(0);
213: }

215: /* MatConvert_SeqAIJ_SeqAIJSELL converts a SeqAIJ matrix into a
216:  * SeqAIJSELL matrix.  This routine is called by the MatCreate_SeqAIJSELL()
217:  * routine, but can also be used to convert an assembled SeqAIJ matrix
218:  * into a SeqAIJSELL one. */
219: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
220: {
222:   Mat            B = *newmat;
223:   Mat_SeqAIJ     *b;
224:   Mat_SeqAIJSELL *aijsell;
225:   PetscBool      set;
226:   PetscBool      sametype;

229:   if (reuse == MAT_INITIAL_MATRIX) {
230:     MatDuplicate(A,MAT_COPY_VALUES,&B);
231:   }

233:   PetscObjectTypeCompare((PetscObject)A,type,&sametype);
234:   if (sametype) return(0);

236:   PetscNewLog(B,&aijsell);
237:   b        = (Mat_SeqAIJ*) B->data;
238:   B->spptr = (void*) aijsell;

240:   /* Disable use of the inode routines so that the AIJSELL ones will be used instead.
241:    * This happens in MatAssemblyEnd_SeqAIJSELL as well, but the assembly end may not be called, so set it here, too.
242:    * As noted elsewhere, I wonder if it might make sense and be feasible to use some of the inode routines. */
243:   b->inode.use = PETSC_FALSE;

245:   /* Set function pointers for methods that we inherit from AIJ but override.
246:    * We also parse some command line options below, since those determine some of the methods we point to. */
247:   B->ops->duplicate        = MatDuplicate_SeqAIJSELL;
248:   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJSELL;
249:   B->ops->destroy          = MatDestroy_SeqAIJSELL;

251:   aijsell->S = NULL;
252:   aijsell->eager_shadow = PETSC_FALSE;

254:   /* Parse command line options. */
255:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"AIJSELL Options","Mat");
256:   PetscOptionsBool("-mat_aijsell_eager_shadow","Eager Shadowing","None",(PetscBool)aijsell->eager_shadow,(PetscBool*)&aijsell->eager_shadow,&set);
257:   PetscOptionsEnd();

259:   /* If A has already been assembled and eager shadowing is specified, build the shadow matrix. */
260:   if (A->assembled && aijsell->eager_shadow) {
261:     MatSeqAIJSELL_build_shadow(A);
262:   }

264:   B->ops->mult             = MatMult_SeqAIJSELL;
265:   B->ops->multtranspose    = MatMultTranspose_SeqAIJSELL;
266:   B->ops->multadd          = MatMultAdd_SeqAIJSELL;
267:   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJSELL;
268:   B->ops->sor              = MatSOR_SeqAIJSELL;

270:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijsell_seqaij_C",MatConvert_SeqAIJSELL_SeqAIJ);

272:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJSELL);
273:   *newmat = B;
274:   return(0);
275: }

277: /*@C
278:    MatCreateSeqAIJSELL - Creates a sparse matrix of type SEQAIJSELL.
279:    This type inherits from AIJ and is largely identical, but keeps a "shadow"
280:    copy of the matrix in SEQSELL format, which is used when this format
281:    may be more suitable for a requested operation. Currently, SEQSELL format
282:    is used for MatMult, MatMultTranspose, MatMultAdd, MatMultTransposeAdd,
283:    and MatSOR operations.
284:    Because SEQAIJSELL is a subtype of SEQAIJ, the option "-mat_seqaij_type seqaijsell" can be used to make
285:    sequential AIJ matrices default to being instances of MATSEQAIJSELL.

287:    Collective

289:    Input Parameters:
290: +  comm - MPI communicator, set to PETSC_COMM_SELF
291: .  m - number of rows
292: .  n - number of columns
293: .  nz - number of nonzeros per row (same for all rows)
294: -  nnz - array containing the number of nonzeros in the various rows
295:          (possibly different for each row) or NULL

297:    Output Parameter:
298: .  A - the matrix

300:    Options Database Keys:
301: .  -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

303:    Notes:
304:    If nnz is given then nz is ignored

306:    Level: intermediate

308: .seealso: MatCreate(), MatCreateMPIAIJSELL(), MatSetValues()
309: @*/
310: PetscErrorCode  MatCreateSeqAIJSELL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
311: {

315:   MatCreate(comm,A);
316:   MatSetSizes(*A,m,n,m,n);
317:   MatSetType(*A,MATSEQAIJSELL);
318:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);
319:   return(0);
320: }

322: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJSELL(Mat A)
323: {

327:   MatSetType(A,MATSEQAIJ);
328:   MatConvert_SeqAIJ_SeqAIJSELL(A,MATSEQAIJSELL,MAT_INPLACE_MATRIX,&A);
329:   return(0);
330: }