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

 24:   PetscFunctionBegin;
 25:   if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));

 27:   /* Reset the original function pointers. */
 28:   B->ops->duplicate        = MatDuplicate_SeqAIJ;
 29:   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJ;
 30:   B->ops->destroy          = MatDestroy_SeqAIJ;
 31:   B->ops->mult             = MatMult_SeqAIJ;
 32:   B->ops->multtranspose    = MatMultTranspose_SeqAIJ;
 33:   B->ops->multadd          = MatMultAdd_SeqAIJ;
 34:   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJ;
 35:   B->ops->sor              = MatSOR_SeqAIJ;

 37:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaijsell_seqaij_C", NULL));

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

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

 46:   /* Change the type of B to MATSEQAIJ. */
 47:   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ));

 49:   *newmat = B;
 50:   PetscFunctionReturn(PETSC_SUCCESS);
 51: }

 53: static PetscErrorCode MatDestroy_SeqAIJSELL(Mat A)
 54: {
 55:   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;

 57:   PetscFunctionBegin;
 58:   /* If MatHeaderMerge() was used, then this SeqAIJSELL matrix will not have an
 59:    * spptr pointer. */
 60:   if (aijsell) {
 61:     /* Clean up everything in the Mat_SeqAIJSELL data structure, then free A->spptr. */
 62:     PetscCall(MatDestroy(&aijsell->S));
 63:     PetscCall(PetscFree(A->spptr));
 64:   }

 66:   /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ()
 67:    * to destroy everything that remains. */
 68:   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ));
 69:   /* Note that I don't call MatSetType().  I believe this is because that
 70:    * is only to be called when *building* a matrix.  I could be wrong, but
 71:    * that is how things work for the SuperLU matrix class. */
 72:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaijsell_seqaij_C", NULL));
 73:   PetscCall(MatDestroy_SeqAIJ(A));
 74:   PetscFunctionReturn(PETSC_SUCCESS);
 75: }

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

 84:   PetscFunctionBegin;
 85:   PetscCall(PetscObjectStateGet((PetscObject)A, &state));
 86:   if (aijsell->S && aijsell->state == state) {
 87:     /* The existing shadow matrix is up-to-date, so simply exit. */
 88:     PetscFunctionReturn(PETSC_SUCCESS);
 89:   }

 91:   PetscCall(PetscLogEventBegin(MAT_Convert, A, 0, 0, 0));
 92:   if (aijsell->S) {
 93:     PetscCall(MatConvert_SeqAIJ_SeqSELL(A, MATSEQSELL, MAT_REUSE_MATRIX, &aijsell->S));
 94:   } else {
 95:     PetscCall(MatConvert_SeqAIJ_SeqSELL(A, MATSEQSELL, MAT_INITIAL_MATRIX, &aijsell->S));
 96:   }
 97:   PetscCall(PetscLogEventEnd(MAT_Convert, A, 0, 0, 0));

 99:   /* Record the ObjectState so that we can tell when the shadow matrix needs updating */
100:   PetscCall(PetscObjectStateGet((PetscObject)A, &aijsell->state));
101:   PetscFunctionReturn(PETSC_SUCCESS);
102: }

104: static PetscErrorCode MatDuplicate_SeqAIJSELL(Mat A, MatDuplicateOption op, Mat *M)
105: {
106:   Mat_SeqAIJSELL *aijsell;
107:   Mat_SeqAIJSELL *aijsell_dest;

109:   PetscFunctionBegin;
110:   PetscCall(MatDuplicate_SeqAIJ(A, op, M));
111:   aijsell      = (Mat_SeqAIJSELL *)A->spptr;
112:   aijsell_dest = (Mat_SeqAIJSELL *)(*M)->spptr;
113:   PetscCall(PetscArraycpy(aijsell_dest, aijsell, 1));
114:   /* We don't duplicate the shadow matrix -- that will be constructed as needed. */
115:   aijsell_dest->S = NULL;
116:   if (aijsell->eager_shadow) PetscCall(MatSeqAIJSELL_build_shadow(A));
117:   PetscFunctionReturn(PETSC_SUCCESS);
118: }

120: static PetscErrorCode MatAssemblyEnd_SeqAIJSELL(Mat A, MatAssemblyType mode)
121: {
122:   Mat_SeqAIJ     *a       = (Mat_SeqAIJ *)A->data;
123:   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;

125:   PetscFunctionBegin;
126:   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);

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

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

139:   PetscCall(MatAssemblyEnd_SeqAIJ(A, mode));

141:   /* If the user has requested "eager" shadowing, create the SELL shadow matrix (if needed; the function checks).
142:    * (The default is to take a "lazy" approach, deferring this until something like MatMult() is called.) */
143:   if (aijsell->eager_shadow) PetscCall(MatSeqAIJSELL_build_shadow(A));
144:   PetscFunctionReturn(PETSC_SUCCESS);
145: }

147: static PetscErrorCode MatMult_SeqAIJSELL(Mat A, Vec xx, Vec yy)
148: {
149:   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;

151:   PetscFunctionBegin;
152:   PetscCall(MatSeqAIJSELL_build_shadow(A));
153:   PetscCall(MatMult_SeqSELL(aijsell->S, xx, yy));
154:   PetscFunctionReturn(PETSC_SUCCESS);
155: }

157: static PetscErrorCode MatMultTranspose_SeqAIJSELL(Mat A, Vec xx, Vec yy)
158: {
159:   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;

161:   PetscFunctionBegin;
162:   PetscCall(MatSeqAIJSELL_build_shadow(A));
163:   PetscCall(MatMultTranspose_SeqSELL(aijsell->S, xx, yy));
164:   PetscFunctionReturn(PETSC_SUCCESS);
165: }

167: static PetscErrorCode MatMultAdd_SeqAIJSELL(Mat A, Vec xx, Vec yy, Vec zz)
168: {
169:   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;

171:   PetscFunctionBegin;
172:   PetscCall(MatSeqAIJSELL_build_shadow(A));
173:   PetscCall(MatMultAdd_SeqSELL(aijsell->S, xx, yy, zz));
174:   PetscFunctionReturn(PETSC_SUCCESS);
175: }

177: static PetscErrorCode MatMultTransposeAdd_SeqAIJSELL(Mat A, Vec xx, Vec yy, Vec zz)
178: {
179:   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;

181:   PetscFunctionBegin;
182:   PetscCall(MatSeqAIJSELL_build_shadow(A));
183:   PetscCall(MatMultTransposeAdd_SeqSELL(aijsell->S, xx, yy, zz));
184:   PetscFunctionReturn(PETSC_SUCCESS);
185: }

187: static PetscErrorCode MatSOR_SeqAIJSELL(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
188: {
189:   Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL *)A->spptr;

191:   PetscFunctionBegin;
192:   PetscCall(MatSeqAIJSELL_build_shadow(A));
193:   PetscCall(MatSOR_SeqSELL(aijsell->S, bb, omega, flag, fshift, its, lits, xx));
194:   PetscFunctionReturn(PETSC_SUCCESS);
195: }

197: /* MatConvert_SeqAIJ_SeqAIJSELL converts a SeqAIJ matrix into a
198:  * SeqAIJSELL matrix.  This routine is called by the MatCreate_SeqAIJSELL()
199:  * routine, but can also be used to convert an assembled SeqAIJ matrix
200:  * into a SeqAIJSELL one. */
201: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat A, MatType type, MatReuse reuse, Mat *newmat)
202: {
203:   Mat             B = *newmat;
204:   Mat_SeqAIJ     *b;
205:   Mat_SeqAIJSELL *aijsell;
206:   PetscBool       set;
207:   PetscBool       sametype;

209:   PetscFunctionBegin;
210:   if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));

212:   PetscCall(PetscObjectTypeCompare((PetscObject)A, type, &sametype));
213:   if (sametype) PetscFunctionReturn(PETSC_SUCCESS);

215:   PetscCall(PetscNew(&aijsell));
216:   b        = (Mat_SeqAIJ *)B->data;
217:   B->spptr = (void *)aijsell;

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

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

230:   aijsell->S            = NULL;
231:   aijsell->eager_shadow = PETSC_FALSE;

233:   /* Parse command line options. */
234:   PetscOptionsBegin(PetscObjectComm((PetscObject)A), ((PetscObject)A)->prefix, "AIJSELL Options", "Mat");
235:   PetscCall(PetscOptionsBool("-mat_aijsell_eager_shadow", "Eager Shadowing", "None", (PetscBool)aijsell->eager_shadow, (PetscBool *)&aijsell->eager_shadow, &set));
236:   PetscOptionsEnd();

238:   /* If A has already been assembled and eager shadowing is specified, build the shadow matrix. */
239:   if (A->assembled && aijsell->eager_shadow) PetscCall(MatSeqAIJSELL_build_shadow(A));

241:   B->ops->mult             = MatMult_SeqAIJSELL;
242:   B->ops->multtranspose    = MatMultTranspose_SeqAIJSELL;
243:   B->ops->multadd          = MatMultAdd_SeqAIJSELL;
244:   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJSELL;
245:   B->ops->sor              = MatSOR_SeqAIJSELL;

247:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaijsell_seqaij_C", MatConvert_SeqAIJSELL_SeqAIJ));

249:   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJSELL));
250:   *newmat = B;
251:   PetscFunctionReturn(PETSC_SUCCESS);
252: }

254: /*@C
255:   MatCreateSeqAIJSELL - Creates a sparse matrix of type `MATSEQAIJSELL`.

257:   Collective

259:   Input Parameters:
260: + comm - MPI communicator, set to `PETSC_COMM_SELF`
261: . m    - number of rows
262: . n    - number of columns
263: . nz   - number of nonzeros per row (same for all rows)
264: - nnz  - array containing the number of nonzeros in the various rows
265:          (possibly different for each row) or `NULL`

267:   Output Parameter:
268: . A - the matrix

270:   Options Database Keys:
271: . -mat_aijsell_eager_shadow - Construct shadow matrix upon matrix assembly; default is to take a "lazy" approach,
272:                                performing this step the first time the matrix is applied

274:   Level: intermediate

276:   Notes:
277:   This type inherits from AIJ and is largely identical, but keeps a "shadow" copy of the matrix
278:   in `MATSEQSELL` format, which is used when this format may be more suitable for a requested
279:   operation. Currently, `MATSEQSELL` format is used for `MatMult()`, `MatMultTranspose()`,
280:   `MatMultAdd()`, `MatMultTransposeAdd()`, and `MatSOR()` operations.

282:   If `nnz` is given then `nz` is ignored

284:   Because `MATSEQAIJSELL` is a subtype of `MATSEQAIJ`, the option `-mat_seqaij_type seqaijsell` can be used to make
285:   sequential `MATSEQAIJ` matrices default to being instances of `MATSEQAIJSELL`.

287: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateMPIAIJSELL()`, `MatSetValues()`
288: @*/
289: PetscErrorCode MatCreateSeqAIJSELL(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A)
290: {
291:   PetscFunctionBegin;
292:   PetscCall(MatCreate(comm, A));
293:   PetscCall(MatSetSizes(*A, m, n, m, n));
294:   PetscCall(MatSetType(*A, MATSEQAIJSELL));
295:   PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*A, nz, nnz));
296:   PetscFunctionReturn(PETSC_SUCCESS);
297: }

299: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJSELL(Mat A)
300: {
301:   PetscFunctionBegin;
302:   PetscCall(MatSetType(A, MATSEQAIJ));
303:   PetscCall(MatConvert_SeqAIJ_SeqAIJSELL(A, MATSEQAIJSELL, MAT_INPLACE_MATRIX, &A));
304:   PetscFunctionReturn(PETSC_SUCCESS);
305: }