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: if (reuse == MAT_INITIAL_MATRIX) {
25: MatDuplicate(A,MAT_COPY_VALUES,&B);
26: }
28: /* Reset the original function pointers. */
29: B->ops->duplicate = MatDuplicate_SeqAIJ;
30: B->ops->assemblyend = MatAssemblyEnd_SeqAIJ;
31: B->ops->destroy = MatDestroy_SeqAIJ;
32: B->ops->mult = MatMult_SeqAIJ;
33: B->ops->multtranspose = MatMultTranspose_SeqAIJ;
34: B->ops->multadd = MatMultAdd_SeqAIJ;
35: B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJ;
36: B->ops->sor = MatSOR_SeqAIJ;
38: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijsell_seqaij_C",NULL);
40: if (reuse == MAT_INITIAL_MATRIX) aijsell = (Mat_SeqAIJSELL*)B->spptr;
42: /* Clean up the Mat_SeqAIJSELL data structure.
43: * Note that MatDestroy() simply returns if passed a NULL value, so it's OK to call even if the shadow matrix was never constructed. */
44: MatDestroy(&aijsell->S);
45: PetscFree(B->spptr);
47: /* Change the type of B to MATSEQAIJ. */
48: PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ);
50: *newmat = B;
51: return 0;
52: }
54: PetscErrorCode MatDestroy_SeqAIJSELL(Mat A)
55: {
56: Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL*) A->spptr;
59: /* If MatHeaderMerge() was used, then this SeqAIJSELL matrix will not have an
60: * spptr pointer. */
61: if (aijsell) {
62: /* Clean up everything in the Mat_SeqAIJSELL data structure, then free A->spptr. */
63: MatDestroy(&aijsell->S);
64: PetscFree(A->spptr);
65: }
67: /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ()
68: * to destroy everything that remains. */
69: PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ);
70: /* Note that I don't call MatSetType(). I believe this is because that
71: * is only to be called when *building* a matrix. I could be wrong, but
72: * that is how things work for the SuperLU matrix class. */
73: MatDestroy_SeqAIJ(A);
74: return 0;
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: PetscObjectStateGet((PetscObject)A,&state);
85: if (aijsell->S && aijsell->state == state) {
86: /* The existing shadow matrix is up-to-date, so simply exit. */
87: return 0;
88: }
90: PetscLogEventBegin(MAT_Convert,A,0,0,0);
91: if (aijsell->S) {
92: MatConvert_SeqAIJ_SeqSELL(A,MATSEQSELL,MAT_REUSE_MATRIX,&aijsell->S);
93: } else {
94: MatConvert_SeqAIJ_SeqSELL(A,MATSEQSELL,MAT_INITIAL_MATRIX,&aijsell->S);
95: }
96: PetscLogEventEnd(MAT_Convert,A,0,0,0);
98: /* Record the ObjectState so that we can tell when the shadow matrix needs updating */
99: PetscObjectStateGet((PetscObject)A,&aijsell->state);
101: return 0;
102: }
104: PetscErrorCode MatDuplicate_SeqAIJSELL(Mat A, MatDuplicateOption op, Mat *M)
105: {
106: Mat_SeqAIJSELL *aijsell;
107: Mat_SeqAIJSELL *aijsell_dest;
109: MatDuplicate_SeqAIJ(A,op,M);
110: aijsell = (Mat_SeqAIJSELL*) A->spptr;
111: aijsell_dest = (Mat_SeqAIJSELL*) (*M)->spptr;
112: PetscArraycpy(aijsell_dest,aijsell,1);
113: /* We don't duplicate the shadow matrix -- that will be constructed as needed. */
114: aijsell_dest->S = NULL;
115: if (aijsell->eager_shadow) {
116: MatSeqAIJSELL_build_shadow(A);
117: }
118: return 0;
119: }
121: PetscErrorCode MatAssemblyEnd_SeqAIJSELL(Mat A, MatAssemblyType mode)
122: {
123: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
124: Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL*)A->spptr;
126: if (mode == MAT_FLUSH_ASSEMBLY) return 0;
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: 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) {
144: MatSeqAIJSELL_build_shadow(A);
145: }
147: return 0;
148: }
150: PetscErrorCode MatMult_SeqAIJSELL(Mat A,Vec xx,Vec yy)
151: {
152: Mat_SeqAIJSELL *aijsell = (Mat_SeqAIJSELL*)A->spptr;
154: MatSeqAIJSELL_build_shadow(A);
155: MatMult_SeqSELL(aijsell->S,xx,yy);
156: return 0;
157: }
159: PetscErrorCode MatMultTranspose_SeqAIJSELL(Mat A,Vec xx,Vec yy)
160: {
161: Mat_SeqAIJSELL *aijsell=(Mat_SeqAIJSELL*)A->spptr;
163: MatSeqAIJSELL_build_shadow(A);
164: MatMultTranspose_SeqSELL(aijsell->S,xx,yy);
165: return 0;
166: }
168: PetscErrorCode MatMultAdd_SeqAIJSELL(Mat A,Vec xx,Vec yy,Vec zz)
169: {
170: Mat_SeqAIJSELL *aijsell=(Mat_SeqAIJSELL*)A->spptr;
172: MatSeqAIJSELL_build_shadow(A);
173: MatMultAdd_SeqSELL(aijsell->S,xx,yy,zz);
174: return 0;
175: }
177: PetscErrorCode MatMultTransposeAdd_SeqAIJSELL(Mat A,Vec xx,Vec yy,Vec zz)
178: {
179: Mat_SeqAIJSELL *aijsell=(Mat_SeqAIJSELL*)A->spptr;
181: MatSeqAIJSELL_build_shadow(A);
182: MatMultTransposeAdd_SeqSELL(aijsell->S,xx,yy,zz);
183: return 0;
184: }
186: PetscErrorCode MatSOR_SeqAIJSELL(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
187: {
188: Mat_SeqAIJSELL *aijsell=(Mat_SeqAIJSELL*)A->spptr;
190: MatSeqAIJSELL_build_shadow(A);
191: MatSOR_SeqSELL(aijsell->S,bb,omega,flag,fshift,its,lits,xx);
192: return 0;
193: }
195: /* MatConvert_SeqAIJ_SeqAIJSELL converts a SeqAIJ matrix into a
196: * SeqAIJSELL matrix. This routine is called by the MatCreate_SeqAIJSELL()
197: * routine, but can also be used to convert an assembled SeqAIJ matrix
198: * into a SeqAIJSELL one. */
199: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
200: {
202: Mat B = *newmat;
203: Mat_SeqAIJ *b;
204: Mat_SeqAIJSELL *aijsell;
205: PetscBool set;
206: PetscBool sametype;
208: if (reuse == MAT_INITIAL_MATRIX) {
209: MatDuplicate(A,MAT_COPY_VALUES,&B);
210: }
212: PetscObjectTypeCompare((PetscObject)A,type,&sametype);
213: if (sametype) return 0;
215: PetscNewLog(B,&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: 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) {
240: MatSeqAIJSELL_build_shadow(A);
241: }
243: B->ops->mult = MatMult_SeqAIJSELL;
244: B->ops->multtranspose = MatMultTranspose_SeqAIJSELL;
245: B->ops->multadd = MatMultAdd_SeqAIJSELL;
246: B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJSELL;
247: B->ops->sor = MatSOR_SeqAIJSELL;
249: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijsell_seqaij_C",MatConvert_SeqAIJSELL_SeqAIJ);
251: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJSELL);
252: *newmat = B;
253: return 0;
254: }
256: /*@C
257: MatCreateSeqAIJSELL - Creates a sparse matrix of type SEQAIJSELL.
258: This type inherits from AIJ and is largely identical, but keeps a "shadow"
259: copy of the matrix in SEQSELL format, which is used when this format
260: may be more suitable for a requested operation. Currently, SEQSELL format
261: is used for MatMult, MatMultTranspose, MatMultAdd, MatMultTransposeAdd,
262: and MatSOR operations.
263: Because SEQAIJSELL is a subtype of SEQAIJ, the option "-mat_seqaij_type seqaijsell" can be used to make
264: sequential AIJ matrices default to being instances of MATSEQAIJSELL.
266: Collective
268: Input Parameters:
269: + comm - MPI communicator, set to PETSC_COMM_SELF
270: . m - number of rows
271: . n - number of columns
272: . nz - number of nonzeros per row (same for all rows)
273: - nnz - array containing the number of nonzeros in the various rows
274: (possibly different for each row) or NULL
276: Output Parameter:
277: . A - the matrix
279: Options Database Keys:
280: . -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
282: Notes:
283: If nnz is given then nz is ignored
285: Level: intermediate
287: .seealso: MatCreate(), MatCreateMPIAIJSELL(), MatSetValues()
288: @*/
289: PetscErrorCode MatCreateSeqAIJSELL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
290: {
291: MatCreate(comm,A);
292: MatSetSizes(*A,m,n,m,n);
293: MatSetType(*A,MATSEQAIJSELL);
294: MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);
295: return 0;
296: }
298: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJSELL(Mat A)
299: {
300: MatSetType(A,MATSEQAIJ);
301: MatConvert_SeqAIJ_SeqAIJSELL(A,MATSEQAIJSELL,MAT_INPLACE_MATRIX,&A);
302: return 0;
303: }