xref: /petsc/src/mat/impls/kaij/kaij.c (revision 910cf402deff9031709e244fa7c1e84ef2fdc3da)
149bd79ccSDebojyoti Ghosh 
249bd79ccSDebojyoti Ghosh /*
349bd79ccSDebojyoti Ghosh   Defines the basic matrix operations for the KAIJ  matrix storage format.
449bd79ccSDebojyoti Ghosh   This format is used to evaluate matrices of the form:
549bd79ccSDebojyoti Ghosh 
649bd79ccSDebojyoti Ghosh     [I \otimes S + A \otimes T]
749bd79ccSDebojyoti Ghosh 
849bd79ccSDebojyoti Ghosh   where
949bd79ccSDebojyoti Ghosh     S is a dense (p \times q) matrix
1049bd79ccSDebojyoti Ghosh     T is a dense (p \times q) matrix
1149bd79ccSDebojyoti Ghosh     A is an AIJ  (n \times n) matrix
1249bd79ccSDebojyoti Ghosh     I is the identity matrix
1349bd79ccSDebojyoti Ghosh 
1449bd79ccSDebojyoti Ghosh   The resulting matrix is (np \times nq)
1549bd79ccSDebojyoti Ghosh 
1649bd79ccSDebojyoti Ghosh   We provide:
1749bd79ccSDebojyoti Ghosh      MatMult()
1849bd79ccSDebojyoti Ghosh      MatMultAdd()
1949bd79ccSDebojyoti Ghosh      MatInvertBlockDiagonal()
2049bd79ccSDebojyoti Ghosh   and
2149bd79ccSDebojyoti Ghosh      MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*)
2249bd79ccSDebojyoti Ghosh 
2349bd79ccSDebojyoti Ghosh   This single directory handles both the sequential and parallel codes
2449bd79ccSDebojyoti Ghosh */
2549bd79ccSDebojyoti Ghosh 
2649bd79ccSDebojyoti Ghosh #include <../src/mat/impls/kaij/kaij.h> /*I "petscmat.h" I*/
2749bd79ccSDebojyoti Ghosh #include <../src/mat/utils/freespace.h>
2849bd79ccSDebojyoti Ghosh #include <petsc/private/vecimpl.h>
2949bd79ccSDebojyoti Ghosh 
3049bd79ccSDebojyoti Ghosh /*@C
3149bd79ccSDebojyoti Ghosh    MatKAIJGetAIJ - Get the AIJ matrix describing the blockwise action of the KAIJ matrix
3249bd79ccSDebojyoti Ghosh 
3349bd79ccSDebojyoti Ghosh    Not Collective, but if the KAIJ matrix is parallel, the AIJ matrix is also parallel
3449bd79ccSDebojyoti Ghosh 
3549bd79ccSDebojyoti Ghosh    Input Parameter:
3649bd79ccSDebojyoti Ghosh .  A - the KAIJ matrix
3749bd79ccSDebojyoti Ghosh 
3849bd79ccSDebojyoti Ghosh    Output Parameter:
3949bd79ccSDebojyoti Ghosh .  B - the AIJ matrix
4049bd79ccSDebojyoti Ghosh 
4149bd79ccSDebojyoti Ghosh    Level: advanced
4249bd79ccSDebojyoti Ghosh 
4349bd79ccSDebojyoti Ghosh    Notes: The reference count on the AIJ matrix is not increased so you should not destroy it.
4449bd79ccSDebojyoti Ghosh 
4549bd79ccSDebojyoti Ghosh .seealso: MatCreateKAIJ()
4649bd79ccSDebojyoti Ghosh @*/
4749bd79ccSDebojyoti Ghosh PetscErrorCode  MatKAIJGetAIJ(Mat A,Mat *B)
4849bd79ccSDebojyoti Ghosh {
4949bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
5049bd79ccSDebojyoti Ghosh   PetscBool      ismpikaij,isseqkaij;
5149bd79ccSDebojyoti Ghosh 
5249bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
5349bd79ccSDebojyoti Ghosh   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIKAIJ,&ismpikaij);CHKERRQ(ierr);
5449bd79ccSDebojyoti Ghosh   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQKAIJ,&isseqkaij);CHKERRQ(ierr);
5549bd79ccSDebojyoti Ghosh   if (ismpikaij) {
5649bd79ccSDebojyoti Ghosh     Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data;
5749bd79ccSDebojyoti Ghosh 
5849bd79ccSDebojyoti Ghosh     *B = b->A;
5949bd79ccSDebojyoti Ghosh   } else if (isseqkaij) {
6049bd79ccSDebojyoti Ghosh     Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data;
6149bd79ccSDebojyoti Ghosh 
6249bd79ccSDebojyoti Ghosh     *B = b->AIJ;
63b04351cbSRichard Tran Mills   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix passed in is not of type KAIJ");
6449bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
6549bd79ccSDebojyoti Ghosh }
6649bd79ccSDebojyoti Ghosh 
6749bd79ccSDebojyoti Ghosh /*@C
6849bd79ccSDebojyoti Ghosh    MatKAIJGetS - Get the S matrix describing the shift action of the KAIJ matrix
6949bd79ccSDebojyoti Ghosh 
700567c835SRichard Tran Mills    Not Collective; the entire S is stored and returned independently on all processes.
7149bd79ccSDebojyoti Ghosh 
7249bd79ccSDebojyoti Ghosh    Input Parameter:
7349bd79ccSDebojyoti Ghosh .  A - the KAIJ matrix
7449bd79ccSDebojyoti Ghosh 
75a5b5c723SRichard Tran Mills    Output Parameters:
76a5b5c723SRichard Tran Mills +  m - the number of rows in S
77a5b5c723SRichard Tran Mills .  n - the number of columns in S
78a5b5c723SRichard Tran Mills -  S - the S matrix, in form of a scalar array in column-major format
7949bd79ccSDebojyoti Ghosh 
80a5b5c723SRichard Tran Mills    Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired)
8131a97b9aSRichard Tran Mills 
8249bd79ccSDebojyoti Ghosh    Level: advanced
8349bd79ccSDebojyoti Ghosh 
8431a97b9aSRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes()
8549bd79ccSDebojyoti Ghosh @*/
86a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetS(Mat A,PetscInt *m,PetscInt *n,PetscScalar **S)
8749bd79ccSDebojyoti Ghosh {
8849bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data;
8949bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
90a5b5c723SRichard Tran Mills   if (m) *m = b->p;
91a5b5c723SRichard Tran Mills   if (n) *n = b->q;
92a5b5c723SRichard Tran Mills   if (S) *S = b->S;
93a5b5c723SRichard Tran Mills   PetscFunctionReturn(0);
94a5b5c723SRichard Tran Mills }
95a5b5c723SRichard Tran Mills 
96a5b5c723SRichard Tran Mills /*@C
97a5b5c723SRichard Tran Mills    MatKAIJGetSRead - Get a read-only pointer to the S matrix describing the shift action of the KAIJ matrix
98a5b5c723SRichard Tran Mills 
99a5b5c723SRichard Tran Mills    Not Collective; the entire S is stored and returned independently on all processes.
100a5b5c723SRichard Tran Mills 
101a5b5c723SRichard Tran Mills    Input Parameter:
102a5b5c723SRichard Tran Mills .  A - the KAIJ matrix
103a5b5c723SRichard Tran Mills 
104a5b5c723SRichard Tran Mills    Output Parameters:
105a5b5c723SRichard Tran Mills +  m - the number of rows in S
106a5b5c723SRichard Tran Mills .  n - the number of columns in S
107a5b5c723SRichard Tran Mills -  S - the S matrix, in form of a scalar array in column-major format
108a5b5c723SRichard Tran Mills 
109a5b5c723SRichard Tran Mills    Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired)
110a5b5c723SRichard Tran Mills 
111a5b5c723SRichard Tran Mills    Level: advanced
112a5b5c723SRichard Tran Mills 
113a5b5c723SRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes()
114a5b5c723SRichard Tran Mills @*/
115a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetSRead(Mat A,PetscInt *m,PetscInt *n,const PetscScalar **S)
116a5b5c723SRichard Tran Mills {
117a5b5c723SRichard Tran Mills   Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data;
118a5b5c723SRichard Tran Mills   PetscFunctionBegin;
119a5b5c723SRichard Tran Mills   if (m) *m = b->p;
120a5b5c723SRichard Tran Mills   if (n) *n = b->q;
121a5b5c723SRichard Tran Mills   if (S) *S = b->S;
122a5b5c723SRichard Tran Mills   PetscFunctionReturn(0);
123a5b5c723SRichard Tran Mills }
124a5b5c723SRichard Tran Mills 
125a5b5c723SRichard Tran Mills /*@C
126a5b5c723SRichard Tran Mills   MatKAIJRestoreS - Restore array obtained with MatKAIJGetS()
127a5b5c723SRichard Tran Mills 
128a5b5c723SRichard Tran Mills   Not collective
129a5b5c723SRichard Tran Mills 
130a5b5c723SRichard Tran Mills   Input Parameter:
131a5b5c723SRichard Tran Mills . A - the KAIJ matrix
132a5b5c723SRichard Tran Mills 
133a5b5c723SRichard Tran Mills   Output Parameter:
134a5b5c723SRichard Tran Mills . S - location of pointer to array obtained with MatKAIJGetS()
135a5b5c723SRichard Tran Mills 
136a5b5c723SRichard Tran Mills   Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored.
137a5b5c723SRichard Tran Mills   If NULL is passed, it will not attempt to zero the array pointer.
138a5b5c723SRichard Tran Mills 
139a5b5c723SRichard Tran Mills   Level: advanced
140a5b5c723SRichard Tran Mills .seealso: MatKAIJGetS(), MatKAIJGetSRead(), MatKAIJRestoreSRead()
141a5b5c723SRichard Tran Mills @*/
142a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreS(Mat A,PetscScalar **S)
143a5b5c723SRichard Tran Mills {
14466f58c76SRichard Tran Mills   PetscErrorCode ierr;
14566f58c76SRichard Tran Mills 
146a5b5c723SRichard Tran Mills   PetscFunctionBegin;
147a5b5c723SRichard Tran Mills   if (S) *S = NULL;
14866f58c76SRichard Tran Mills   ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr);
149a5b5c723SRichard Tran Mills   PetscFunctionReturn(0);
150a5b5c723SRichard Tran Mills }
151a5b5c723SRichard Tran Mills 
152a5b5c723SRichard Tran Mills /*@C
153a5b5c723SRichard Tran Mills   MatKAIJRestoreSRead - Restore array obtained with MatKAIJGetSRead()
154a5b5c723SRichard Tran Mills 
155a5b5c723SRichard Tran Mills   Not collective
156a5b5c723SRichard Tran Mills 
157a5b5c723SRichard Tran Mills   Input Parameter:
158a5b5c723SRichard Tran Mills . A - the KAIJ matrix
159a5b5c723SRichard Tran Mills 
160a5b5c723SRichard Tran Mills   Output Parameter:
161a5b5c723SRichard Tran Mills . S - location of pointer to array obtained with MatKAIJGetS()
162a5b5c723SRichard Tran Mills 
163a5b5c723SRichard Tran Mills   Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored.
164a5b5c723SRichard Tran Mills   If NULL is passed, it will not attempt to zero the array pointer.
165a5b5c723SRichard Tran Mills 
166a5b5c723SRichard Tran Mills   Level: advanced
167a5b5c723SRichard Tran Mills .seealso: MatKAIJGetS(), MatKAIJGetSRead(), MatKAIJRestoreSRead()
168a5b5c723SRichard Tran Mills @*/
169a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreSRead(Mat A,const PetscScalar **S)
170a5b5c723SRichard Tran Mills {
171a5b5c723SRichard Tran Mills   PetscFunctionBegin;
172a5b5c723SRichard Tran Mills   if (S) *S = NULL;
17349bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
17449bd79ccSDebojyoti Ghosh }
17549bd79ccSDebojyoti Ghosh 
17649bd79ccSDebojyoti Ghosh /*@C
17731a97b9aSRichard Tran Mills    MatKAIJGetT - Get the transformation matrix T associated with the KAIJ matrix
17849bd79ccSDebojyoti Ghosh 
1790567c835SRichard Tran Mills    Not Collective; the entire T is stored and returned independently on all processes
18049bd79ccSDebojyoti Ghosh 
18149bd79ccSDebojyoti Ghosh    Input Parameter:
18249bd79ccSDebojyoti Ghosh .  A - the KAIJ matrix
18349bd79ccSDebojyoti Ghosh 
18449bd79ccSDebojyoti Ghosh    Output Parameter:
185a5b5c723SRichard Tran Mills +  m - the number of rows in T
186a5b5c723SRichard Tran Mills .  n - the number of columns in T
187a5b5c723SRichard Tran Mills -  T - the T matrix, in form of a scalar array in column-major format
18849bd79ccSDebojyoti Ghosh 
189a5b5c723SRichard Tran Mills    Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired)
19031a97b9aSRichard Tran Mills 
19149bd79ccSDebojyoti Ghosh    Level: advanced
19249bd79ccSDebojyoti Ghosh 
19331a97b9aSRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes()
19449bd79ccSDebojyoti Ghosh @*/
195a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetT(Mat A,PetscInt *m,PetscInt *n,PetscScalar **T)
19649bd79ccSDebojyoti Ghosh {
19749bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data;
19849bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
199a5b5c723SRichard Tran Mills   if (m) *m = b->p;
200a5b5c723SRichard Tran Mills   if (n) *n = b->q;
201a5b5c723SRichard Tran Mills   if (T) *T = b->T;
202a5b5c723SRichard Tran Mills   PetscFunctionReturn(0);
203a5b5c723SRichard Tran Mills }
204a5b5c723SRichard Tran Mills 
205a5b5c723SRichard Tran Mills /*@C
206a5b5c723SRichard Tran Mills    MatKAIJGetTRead - Get a read-only pointer to the transformation matrix T associated with the KAIJ matrix
207a5b5c723SRichard Tran Mills 
208a5b5c723SRichard Tran Mills    Not Collective; the entire T is stored and returned independently on all processes
209a5b5c723SRichard Tran Mills 
210a5b5c723SRichard Tran Mills    Input Parameter:
211a5b5c723SRichard Tran Mills .  A - the KAIJ matrix
212a5b5c723SRichard Tran Mills 
213a5b5c723SRichard Tran Mills    Output Parameter:
214a5b5c723SRichard Tran Mills +  m - the number of rows in T
215a5b5c723SRichard Tran Mills .  n - the number of columns in T
216a5b5c723SRichard Tran Mills -  T - the T matrix, in form of a scalar array in column-major format
217a5b5c723SRichard Tran Mills 
218a5b5c723SRichard Tran Mills    Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired)
219a5b5c723SRichard Tran Mills 
220a5b5c723SRichard Tran Mills    Level: advanced
221a5b5c723SRichard Tran Mills 
222a5b5c723SRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes()
223a5b5c723SRichard Tran Mills @*/
224a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetTRead(Mat A,PetscInt *m,PetscInt *n,const PetscScalar **T)
225a5b5c723SRichard Tran Mills {
226a5b5c723SRichard Tran Mills   Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data;
227a5b5c723SRichard Tran Mills   PetscFunctionBegin;
228a5b5c723SRichard Tran Mills   if (m) *m = b->p;
229a5b5c723SRichard Tran Mills   if (n) *n = b->q;
230a5b5c723SRichard Tran Mills   if (T) *T = b->T;
231a5b5c723SRichard Tran Mills   PetscFunctionReturn(0);
232a5b5c723SRichard Tran Mills }
233a5b5c723SRichard Tran Mills 
234a5b5c723SRichard Tran Mills /*@C
235a5b5c723SRichard Tran Mills   MatKAIJRestoreT - Restore array obtained with MatKAIJGetT()
236a5b5c723SRichard Tran Mills 
237a5b5c723SRichard Tran Mills   Not collective
238a5b5c723SRichard Tran Mills 
239a5b5c723SRichard Tran Mills   Input Parameter:
240a5b5c723SRichard Tran Mills . A - the KAIJ matrix
241a5b5c723SRichard Tran Mills 
242a5b5c723SRichard Tran Mills   Output Parameter:
243a5b5c723SRichard Tran Mills . T - location of pointer to array obtained with MatKAIJGetS()
244a5b5c723SRichard Tran Mills 
245a5b5c723SRichard Tran Mills   Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored.
246a5b5c723SRichard Tran Mills   If NULL is passed, it will not attempt to zero the array pointer.
247a5b5c723SRichard Tran Mills 
248a5b5c723SRichard Tran Mills   Level: advanced
249a5b5c723SRichard Tran Mills .seealso: MatKAIJGetT(), MatKAIJGetTRead(), MatKAIJRestoreTRead()
250a5b5c723SRichard Tran Mills @*/
251a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreT(Mat A,PetscScalar **T)
252a5b5c723SRichard Tran Mills {
25366f58c76SRichard Tran Mills   PetscErrorCode ierr;
25466f58c76SRichard Tran Mills 
255a5b5c723SRichard Tran Mills   PetscFunctionBegin;
256a5b5c723SRichard Tran Mills   if (T) *T = NULL;
25766f58c76SRichard Tran Mills   ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr);
258a5b5c723SRichard Tran Mills   PetscFunctionReturn(0);
259a5b5c723SRichard Tran Mills }
260a5b5c723SRichard Tran Mills 
261a5b5c723SRichard Tran Mills /*@C
262a5b5c723SRichard Tran Mills   MatKAIJRestoreTRead - Restore array obtained with MatKAIJGetTRead()
263a5b5c723SRichard Tran Mills 
264a5b5c723SRichard Tran Mills   Not collective
265a5b5c723SRichard Tran Mills 
266a5b5c723SRichard Tran Mills   Input Parameter:
267a5b5c723SRichard Tran Mills . A - the KAIJ matrix
268a5b5c723SRichard Tran Mills 
269a5b5c723SRichard Tran Mills   Output Parameter:
270a5b5c723SRichard Tran Mills . T - location of pointer to array obtained with MatKAIJGetS()
271a5b5c723SRichard Tran Mills 
272a5b5c723SRichard Tran Mills   Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored.
273a5b5c723SRichard Tran Mills   If NULL is passed, it will not attempt to zero the array pointer.
274a5b5c723SRichard Tran Mills 
275a5b5c723SRichard Tran Mills   Level: advanced
276a5b5c723SRichard Tran Mills .seealso: MatKAIJGetT(), MatKAIJGetTRead(), MatKAIJRestoreTRead()
277a5b5c723SRichard Tran Mills @*/
278a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreTRead(Mat A,const PetscScalar **T)
279a5b5c723SRichard Tran Mills {
280a5b5c723SRichard Tran Mills   PetscFunctionBegin;
281a5b5c723SRichard Tran Mills   if (T) *T = NULL;
28249bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
28349bd79ccSDebojyoti Ghosh }
28449bd79ccSDebojyoti Ghosh 
2850567c835SRichard Tran Mills /*@
2860567c835SRichard Tran Mills    MatKAIJSetAIJ - Set the AIJ matrix describing the blockwise action of the KAIJ matrix
2870567c835SRichard Tran Mills 
2880567c835SRichard Tran Mills    Logically Collective; if the AIJ matrix is parallel, the KAIJ matrix is also parallel
2890567c835SRichard Tran Mills 
2900567c835SRichard Tran Mills    Input Parameters:
2910567c835SRichard Tran Mills +  A - the KAIJ matrix
2920567c835SRichard Tran Mills -  B - the AIJ matrix
2930567c835SRichard Tran Mills 
29415b9d025SRichard Tran Mills    Notes:
29515b9d025SRichard Tran Mills    This function increases the reference count on the AIJ matrix, so the user is free to destroy the matrix if it is not needed.
29615b9d025SRichard Tran Mills    Changes to the entries of the AIJ matrix will immediately affect the KAIJ matrix.
29715b9d025SRichard Tran Mills 
2980567c835SRichard Tran Mills    Level: advanced
2990567c835SRichard Tran Mills 
3000567c835SRichard Tran Mills .seealso: MatKAIJGetAIJ(), MatKAIJSetS(), MatKAIJSetT()
3010567c835SRichard Tran Mills @*/
3020567c835SRichard Tran Mills PetscErrorCode MatKAIJSetAIJ(Mat A,Mat B)
3030567c835SRichard Tran Mills {
3040567c835SRichard Tran Mills   PetscErrorCode ierr;
3050567c835SRichard Tran Mills   PetscMPIInt    size;
3060567c835SRichard Tran Mills 
3070567c835SRichard Tran Mills   PetscFunctionBegin;
3080567c835SRichard Tran Mills   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
3090567c835SRichard Tran Mills   if (size == 1) {
3100567c835SRichard Tran Mills     Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data;
3110567c835SRichard Tran Mills     a->AIJ = B;
3120567c835SRichard Tran Mills   } else {
3130567c835SRichard Tran Mills     Mat_MPIKAIJ *a = (Mat_MPIKAIJ*)A->data;
3140567c835SRichard Tran Mills     a->A = B;
3150567c835SRichard Tran Mills   }
31615b9d025SRichard Tran Mills   ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
3170567c835SRichard Tran Mills   PetscFunctionReturn(0);
3180567c835SRichard Tran Mills }
3190567c835SRichard Tran Mills 
3200567c835SRichard Tran Mills /*@C
3210567c835SRichard Tran Mills    MatKAIJSetS - Set the S matrix describing the shift action of the KAIJ matrix
3220567c835SRichard Tran Mills 
3230567c835SRichard Tran Mills    Logically Collective; the entire S is stored independently on all processes.
3240567c835SRichard Tran Mills 
3250567c835SRichard Tran Mills    Input Parameters:
3260567c835SRichard Tran Mills +  A - the KAIJ matrix
3270567c835SRichard Tran Mills .  p - the number of rows in S
3280567c835SRichard Tran Mills .  q - the number of columns in S
3290567c835SRichard Tran Mills -  S - the S matrix, in form of a scalar array in column-major format
3300567c835SRichard Tran Mills 
3310567c835SRichard Tran Mills    Notes: The dimensions p and q must match those of the transformation matrix T associated with the KAIJ matrix.
33288f48298SRichard Tran Mills    The S matrix is copied, so the user can destroy this array.
3330567c835SRichard Tran Mills 
3340567c835SRichard Tran Mills    Level: Advanced
3350567c835SRichard Tran Mills 
3360567c835SRichard Tran Mills .seealso: MatKAIJGetS(), MatKAIJSetT(), MatKAIJSetAIJ()
3370567c835SRichard Tran Mills @*/
3380567c835SRichard Tran Mills PetscErrorCode MatKAIJSetS(Mat A,PetscInt p,PetscInt q,const PetscScalar S[])
3390567c835SRichard Tran Mills {
3400567c835SRichard Tran Mills   PetscErrorCode ierr;
3410567c835SRichard Tran Mills   Mat_SeqKAIJ    *a = (Mat_SeqKAIJ*)A->data;
3420567c835SRichard Tran Mills 
3430567c835SRichard Tran Mills   PetscFunctionBegin;
3440567c835SRichard Tran Mills   ierr = PetscFree(a->S);CHKERRQ(ierr);
3450567c835SRichard Tran Mills   if (S) {
346a84f8069SRichard Tran Mills     ierr = PetscMalloc1(p*q*sizeof(PetscScalar),&a->S);CHKERRQ(ierr);
3470567c835SRichard Tran Mills     ierr = PetscMemcpy(a->S,S,p*q*sizeof(PetscScalar));CHKERRQ(ierr);
3480567c835SRichard Tran Mills   } else  a->S = NULL;
3490567c835SRichard Tran Mills 
3500567c835SRichard Tran Mills   a->p = p;
3510567c835SRichard Tran Mills   a->q = q;
3520567c835SRichard Tran Mills   PetscFunctionReturn(0);
3530567c835SRichard Tran Mills }
3540567c835SRichard Tran Mills 
3550567c835SRichard Tran Mills /*@C
356*910cf402Sprj-    MatKAIJGetScaledIdentity - Check if both S and T are scaled identities.
357*910cf402Sprj- 
358*910cf402Sprj-    Logically Collective.
359*910cf402Sprj- 
360*910cf402Sprj-    Input Parameter:
361*910cf402Sprj- .  A - the KAIJ matrix
362*910cf402Sprj- 
363*910cf402Sprj-   Output Parameter:
364*910cf402Sprj- .  identity - the Boolean value
365*910cf402Sprj- 
366*910cf402Sprj-    Level: Advanced
367*910cf402Sprj- 
368*910cf402Sprj- .seealso: MatKAIJGetS(), MatKAIJGetT()
369*910cf402Sprj- @*/
370*910cf402Sprj- PetscErrorCode MatKAIJGetScaledIdentity(Mat A,PetscBool* identity)
371*910cf402Sprj- {
372*910cf402Sprj-   Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data;
373*910cf402Sprj-   PetscInt    i,j;
374*910cf402Sprj- 
375*910cf402Sprj-   PetscFunctionBegin;
376*910cf402Sprj-   if (a->p != a->q) {
377*910cf402Sprj-     *identity = PETSC_FALSE;
378*910cf402Sprj-     PetscFunctionReturn(0);
379*910cf402Sprj-   } else *identity = PETSC_TRUE;
380*910cf402Sprj-   if (!a->isTI || a->S) {
381*910cf402Sprj-     for (i=0; i<a->p && *identity; i++) {
382*910cf402Sprj-       for (j=0; j<a->p && *identity; j++) {
383*910cf402Sprj-         if (i != j) {
384*910cf402Sprj-           if(a->S && PetscAbsScalar(a->S[i+j*a->p]) > PETSC_SMALL) *identity = PETSC_FALSE;
385*910cf402Sprj-           if(a->T && PetscAbsScalar(a->T[i+j*a->p]) > PETSC_SMALL) *identity = PETSC_FALSE;
386*910cf402Sprj-         } else {
387*910cf402Sprj-           if(a->S && PetscAbsScalar(a->S[i*(a->p+1)]-a->S[0]) > PETSC_SMALL) *identity = PETSC_FALSE;
388*910cf402Sprj-           if(a->T && PetscAbsScalar(a->T[i*(a->p+1)]-a->T[0]) > PETSC_SMALL) *identity = PETSC_FALSE;
389*910cf402Sprj-         }
390*910cf402Sprj-       }
391*910cf402Sprj-     }
392*910cf402Sprj-   }
393*910cf402Sprj-   PetscFunctionReturn(0);
394*910cf402Sprj- }
395*910cf402Sprj- 
396*910cf402Sprj- /*@C
3970567c835SRichard Tran Mills    MatKAIJSetT - Set the transformation matrix T associated with the KAIJ matrix
3980567c835SRichard Tran Mills 
3990567c835SRichard Tran Mills    Logically Collective; the entire T is stored independently on all processes.
4000567c835SRichard Tran Mills 
4010567c835SRichard Tran Mills    Input Parameters:
4020567c835SRichard Tran Mills +  A - the KAIJ matrix
4030567c835SRichard Tran Mills .  p - the number of rows in S
4040567c835SRichard Tran Mills .  q - the number of columns in S
4050567c835SRichard Tran Mills -  T - the T matrix, in form of a scalar array in column-major format
4060567c835SRichard Tran Mills 
4070567c835SRichard Tran Mills    Notes: The dimensions p and q must match those of the shift matrix S associated with the KAIJ matrix.
40888f48298SRichard Tran Mills    The T matrix is copied, so the user can destroy this array.
4090567c835SRichard Tran Mills 
4100567c835SRichard Tran Mills    Level: Advanced
4110567c835SRichard Tran Mills 
4120567c835SRichard Tran Mills .seealso: MatKAIJGetT(), MatKAIJSetS(), MatKAIJSetAIJ()
4130567c835SRichard Tran Mills @*/
4140567c835SRichard Tran Mills PetscErrorCode MatKAIJSetT(Mat A,PetscInt p,PetscInt q,const PetscScalar T[])
4150567c835SRichard Tran Mills {
4160567c835SRichard Tran Mills   PetscErrorCode ierr;
4170567c835SRichard Tran Mills   PetscInt       i,j;
4180567c835SRichard Tran Mills   Mat_SeqKAIJ    *a = (Mat_SeqKAIJ*)A->data;
4190567c835SRichard Tran Mills   PetscBool      isTI = PETSC_FALSE;
4200567c835SRichard Tran Mills 
4210567c835SRichard Tran Mills   PetscFunctionBegin;
4220567c835SRichard Tran Mills   /* check if T is an identity matrix */
4230567c835SRichard Tran Mills   if (T && (p == q)) {
4240567c835SRichard Tran Mills     isTI = PETSC_TRUE;
4250567c835SRichard Tran Mills     for (i=0; i<p; i++) {
4260567c835SRichard Tran Mills       for (j=0; j<q; j++) {
4270567c835SRichard Tran Mills         if (i == j) {
4280567c835SRichard Tran Mills           /* diagonal term must be 1 */
4290567c835SRichard Tran Mills           if (T[i+j*p] != 1.0) isTI = PETSC_FALSE;
4300567c835SRichard Tran Mills         } else {
4310567c835SRichard Tran Mills           /* off-diagonal term must be 0 */
4320567c835SRichard Tran Mills           if (T[i+j*p] != 0.0) isTI = PETSC_FALSE;
4330567c835SRichard Tran Mills         }
4340567c835SRichard Tran Mills       }
4350567c835SRichard Tran Mills     }
4360567c835SRichard Tran Mills   }
4370567c835SRichard Tran Mills   a->isTI = isTI;
4380567c835SRichard Tran Mills 
4390567c835SRichard Tran Mills   ierr = PetscFree(a->T);CHKERRQ(ierr);
4400567c835SRichard Tran Mills   if (T && (!isTI)) {
441a84f8069SRichard Tran Mills     ierr = PetscMalloc1(p*q*sizeof(PetscScalar),&a->T);CHKERRQ(ierr);
4420567c835SRichard Tran Mills     ierr = PetscMemcpy(a->T,T,p*q*sizeof(PetscScalar));CHKERRQ(ierr);
44350d19d74SRichard Tran Mills   } else a->T = NULL;
4440567c835SRichard Tran Mills 
4450567c835SRichard Tran Mills   a->p = p;
4460567c835SRichard Tran Mills   a->q = q;
4470567c835SRichard Tran Mills   PetscFunctionReturn(0);
4480567c835SRichard Tran Mills }
4490567c835SRichard Tran Mills 
45049bd79ccSDebojyoti Ghosh PetscErrorCode MatDestroy_SeqKAIJ(Mat A)
45149bd79ccSDebojyoti Ghosh {
45249bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
45349bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ    *b = (Mat_SeqKAIJ*)A->data;
45449bd79ccSDebojyoti Ghosh 
45549bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
45649bd79ccSDebojyoti Ghosh   ierr = MatDestroy(&b->AIJ);CHKERRQ(ierr);
457a84f8069SRichard Tran Mills   ierr = PetscFree(b->S);CHKERRQ(ierr);
458a84f8069SRichard Tran Mills   ierr = PetscFree(b->T);CHKERRQ(ierr);
459a84f8069SRichard Tran Mills   ierr = PetscFree(b->ibdiag);CHKERRQ(ierr);
46049bd79ccSDebojyoti Ghosh   ierr = PetscFree5(b->sor.w,b->sor.y,b->sor.work,b->sor.t,b->sor.arr);CHKERRQ(ierr);
46149bd79ccSDebojyoti Ghosh   ierr = PetscFree(A->data);CHKERRQ(ierr);
46249bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
46349bd79ccSDebojyoti Ghosh }
46449bd79ccSDebojyoti Ghosh 
46549bd79ccSDebojyoti Ghosh PetscErrorCode MatSetUp_KAIJ(Mat A)
46649bd79ccSDebojyoti Ghosh {
4670567c835SRichard Tran Mills   PetscErrorCode ierr;
4680567c835SRichard Tran Mills   PetscInt       n;
4690567c835SRichard Tran Mills   PetscMPIInt    size;
4700567c835SRichard Tran Mills   Mat_SeqKAIJ    *seqkaij = (Mat_SeqKAIJ*)A->data;
4710567c835SRichard Tran Mills 
47249bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
4730567c835SRichard Tran Mills   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
4740567c835SRichard Tran Mills   if (size == 1) {
4750567c835SRichard Tran Mills     ierr = MatSetSizes(A,seqkaij->p*seqkaij->AIJ->rmap->n,seqkaij->q*seqkaij->AIJ->cmap->n,seqkaij->p*seqkaij->AIJ->rmap->N,seqkaij->q*seqkaij->AIJ->cmap->N);CHKERRQ(ierr);
4760567c835SRichard Tran Mills     ierr = PetscLayoutSetBlockSize(A->rmap,seqkaij->p);CHKERRQ(ierr);
4770567c835SRichard Tran Mills     ierr = PetscLayoutSetBlockSize(A->cmap,seqkaij->q);CHKERRQ(ierr);
4780567c835SRichard Tran Mills     ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
4790567c835SRichard Tran Mills     ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
4800567c835SRichard Tran Mills   } else {
4810567c835SRichard Tran Mills     Mat_MPIKAIJ *a;
4820567c835SRichard Tran Mills     Mat_MPIAIJ  *mpiaij;
4830567c835SRichard Tran Mills     IS          from,to;
4840567c835SRichard Tran Mills     Vec         gvec;
4850567c835SRichard Tran Mills     PetscScalar *T;
4860567c835SRichard Tran Mills     PetscInt    i,j;
4870567c835SRichard Tran Mills 
4880567c835SRichard Tran Mills     a = (Mat_MPIKAIJ*)A->data;
489d3f912faSRichard Tran Mills     mpiaij = (Mat_MPIAIJ*)a->A->data;
4900567c835SRichard Tran Mills     ierr = MatSetSizes(A,a->p*a->A->rmap->n,a->q*a->A->cmap->n,a->p*a->A->rmap->N,a->q*a->A->cmap->N);CHKERRQ(ierr);
4910567c835SRichard Tran Mills     ierr = PetscLayoutSetBlockSize(A->rmap,seqkaij->p);CHKERRQ(ierr);
4920567c835SRichard Tran Mills     ierr = PetscLayoutSetBlockSize(A->cmap,seqkaij->q);CHKERRQ(ierr);
4930567c835SRichard Tran Mills     ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
4940567c835SRichard Tran Mills     ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
4950567c835SRichard Tran Mills 
4960567c835SRichard Tran Mills     if (a->isTI) {
4970567c835SRichard Tran Mills       /* If the transformation matrix associated with the parallel matrix A is the identity matrix, then a->T will be NULL.
4980567c835SRichard Tran Mills        * In this case, if we pass a->T directly to the MatCreateKAIJ() calls to create the sequential submatrices, the routine will
4990567c835SRichard Tran Mills        * not be able to tell that transformation matrix should be set to the identity; thus we create a temporary identity matrix
5000567c835SRichard Tran Mills        * to pass in. */
501a84f8069SRichard Tran Mills       ierr = PetscMalloc1(a->p*a->q*sizeof(PetscScalar),&T);CHKERRQ(ierr);
5020567c835SRichard Tran Mills       for (i=0; i<a->p; i++) {
5030567c835SRichard Tran Mills         for (j=0; j<a->q; j++) {
5040567c835SRichard Tran Mills           if (i==j) T[i+j*a->p] = 1.0;
5050567c835SRichard Tran Mills           else      T[i+j*a->p] = 0.0;
5060567c835SRichard Tran Mills         }
5070567c835SRichard Tran Mills       }
5080567c835SRichard Tran Mills     } else T = a->T;
5090567c835SRichard Tran Mills     ierr = MatCreateKAIJ(mpiaij->A,a->p,a->q,a->S,T,&a->AIJ);CHKERRQ(ierr);
5100567c835SRichard Tran Mills     ierr = MatCreateKAIJ(mpiaij->B,a->p,a->q,NULL,T,&a->OAIJ);CHKERRQ(ierr);
511c138d2acSRichard Tran Mills     if (a->isTI) {
5120567c835SRichard Tran Mills       ierr = PetscFree(T);CHKERRQ(ierr);
513c138d2acSRichard Tran Mills     }
5140567c835SRichard Tran Mills 
5150567c835SRichard Tran Mills     ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr);
5160567c835SRichard Tran Mills     ierr = VecCreate(PETSC_COMM_SELF,&a->w);CHKERRQ(ierr);
5170567c835SRichard Tran Mills     ierr = VecSetSizes(a->w,n*a->q,n*a->q);CHKERRQ(ierr);
5180567c835SRichard Tran Mills     ierr = VecSetBlockSize(a->w,a->q);CHKERRQ(ierr);
5190567c835SRichard Tran Mills     ierr = VecSetType(a->w,VECSEQ);CHKERRQ(ierr);
5200567c835SRichard Tran Mills 
5210567c835SRichard Tran Mills     /* create two temporary Index sets for build scatter gather */
5220567c835SRichard Tran Mills     ierr = ISCreateBlock(PetscObjectComm((PetscObject)a->A),a->q,n,mpiaij->garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
5230567c835SRichard Tran Mills     ierr = ISCreateStride(PETSC_COMM_SELF,n*a->q,0,1,&to);CHKERRQ(ierr);
5240567c835SRichard Tran Mills 
5250567c835SRichard Tran Mills     /* create temporary global vector to generate scatter context */
5260567c835SRichard Tran Mills     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)a->A),a->q,a->q*a->A->cmap->n,a->q*a->A->cmap->N,NULL,&gvec);CHKERRQ(ierr);
5270567c835SRichard Tran Mills 
5280567c835SRichard Tran Mills     /* generate the scatter context */
5294589b4e5SRichard Tran Mills     ierr = VecScatterCreate(gvec,from,a->w,to,&a->ctx);CHKERRQ(ierr);
5300567c835SRichard Tran Mills 
5310567c835SRichard Tran Mills     ierr = ISDestroy(&from);CHKERRQ(ierr);
5320567c835SRichard Tran Mills     ierr = ISDestroy(&to);CHKERRQ(ierr);
5330567c835SRichard Tran Mills     ierr = VecDestroy(&gvec);CHKERRQ(ierr);
5340567c835SRichard Tran Mills   }
5350567c835SRichard Tran Mills 
5360567c835SRichard Tran Mills   A->assembled = PETSC_TRUE;
53749bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
53849bd79ccSDebojyoti Ghosh }
53949bd79ccSDebojyoti Ghosh 
540e6985dafSRichard Tran Mills PetscErrorCode MatView_KAIJ(Mat A,PetscViewer viewer)
54149bd79ccSDebojyoti Ghosh {
542e6985dafSRichard Tran Mills   PetscViewerFormat format;
543e6985dafSRichard Tran Mills   Mat_SeqKAIJ       *a = (Mat_SeqKAIJ*)A->data;
54449bd79ccSDebojyoti Ghosh   Mat               B;
545e6985dafSRichard Tran Mills   PetscInt          i;
546e6985dafSRichard Tran Mills   PetscErrorCode    ierr;
547e6985dafSRichard Tran Mills   PetscBool         ismpikaij;
54849bd79ccSDebojyoti Ghosh 
54949bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
550e6985dafSRichard Tran Mills   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIKAIJ,&ismpikaij);CHKERRQ(ierr);
551e6985dafSRichard Tran Mills   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
552e6985dafSRichard Tran Mills   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_IMPL) {
553e6985dafSRichard Tran Mills     ierr = PetscViewerASCIIPrintf(viewer,"S and T have %D rows and %D columns\n",a->p,a->q);CHKERRQ(ierr);
554e6985dafSRichard Tran Mills 
555e6985dafSRichard Tran Mills     /* Print appropriate details for S. */
556e6985dafSRichard Tran Mills     if (!a->S) {
5572ae760e3SRichard Tran Mills       ierr = PetscViewerASCIIPrintf(viewer,"S is NULL\n");CHKERRQ(ierr);
558e6985dafSRichard Tran Mills     } else if (format == PETSC_VIEWER_ASCII_IMPL) {
559e6985dafSRichard Tran Mills       ierr = PetscViewerASCIIPrintf(viewer,"Entries of S are ");CHKERRQ(ierr);
560e6985dafSRichard Tran Mills       for (i=0; i<(a->p * a->q); i++) {
561e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
562e6985dafSRichard Tran Mills         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e ",(double)PetscRealPart(a->S[i]),(double)PetscImaginaryPart(a->S[i]));CHKERRQ(ierr);
563e6985dafSRichard Tran Mills #else
564e6985dafSRichard Tran Mills         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(a->S[i]));CHKERRQ(ierr);
565e6985dafSRichard Tran Mills #endif
566e6985dafSRichard Tran Mills       }
567e6985dafSRichard Tran Mills       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
56849bd79ccSDebojyoti Ghosh     }
56949bd79ccSDebojyoti Ghosh 
570e6985dafSRichard Tran Mills     /* Print appropriate details for T. */
571e6985dafSRichard Tran Mills     if (a->isTI) {
5722ae760e3SRichard Tran Mills       ierr = PetscViewerASCIIPrintf(viewer,"T is the identity matrix\n");CHKERRQ(ierr);
573e6985dafSRichard Tran Mills     } else if (!a->T) {
5742ae760e3SRichard Tran Mills       ierr = PetscViewerASCIIPrintf(viewer,"T is NULL\n");CHKERRQ(ierr);
575e6985dafSRichard Tran Mills     } else if (format == PETSC_VIEWER_ASCII_IMPL) {
576e6985dafSRichard Tran Mills       ierr = PetscViewerASCIIPrintf(viewer,"Entries of T are ");CHKERRQ(ierr);
577e6985dafSRichard Tran Mills       for (i=0; i<(a->p * a->q); i++) {
578e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
579e6985dafSRichard Tran Mills         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e ",(double)PetscRealPart(a->T[i]),(double)PetscImaginaryPart(a->T[i]));CHKERRQ(ierr);
580e6985dafSRichard Tran Mills #else
581e6985dafSRichard Tran Mills         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(a->T[i]));CHKERRQ(ierr);
582e6985dafSRichard Tran Mills #endif
583e6985dafSRichard Tran Mills       }
584e6985dafSRichard Tran Mills       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
585e6985dafSRichard Tran Mills     }
58649bd79ccSDebojyoti Ghosh 
587e6985dafSRichard Tran Mills     /* Now print details for the AIJ matrix, using the AIJ viewer. */
588e6985dafSRichard Tran Mills     ierr = PetscViewerASCIIPrintf(viewer,"Now viewing the associated AIJ matrix:\n");CHKERRQ(ierr);
589e6985dafSRichard Tran Mills     if (ismpikaij) {
590e6985dafSRichard Tran Mills       Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data;
591e6985dafSRichard Tran Mills       ierr = MatView(b->A,viewer);CHKERRQ(ierr);
592e6985dafSRichard Tran Mills     } else {
593e6985dafSRichard Tran Mills       ierr = MatView(a->AIJ,viewer);CHKERRQ(ierr);
594e6985dafSRichard Tran Mills     }
595e6985dafSRichard Tran Mills 
596e6985dafSRichard Tran Mills   } else {
597e6985dafSRichard Tran Mills     /* For all other matrix viewer output formats, simply convert to an AIJ matrix and call MatView() on that. */
598e6985dafSRichard Tran Mills     if (ismpikaij) {
59949bd79ccSDebojyoti Ghosh       ierr = MatConvert(A,MATMPIAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
600e6985dafSRichard Tran Mills     } else {
601e6985dafSRichard Tran Mills       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
602e6985dafSRichard Tran Mills     }
60349bd79ccSDebojyoti Ghosh     ierr = MatView(B,viewer);CHKERRQ(ierr);
60449bd79ccSDebojyoti Ghosh     ierr = MatDestroy(&B);CHKERRQ(ierr);
605e6985dafSRichard Tran Mills   }
60649bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
60749bd79ccSDebojyoti Ghosh }
60849bd79ccSDebojyoti Ghosh 
60949bd79ccSDebojyoti Ghosh PetscErrorCode MatDestroy_MPIKAIJ(Mat A)
61049bd79ccSDebojyoti Ghosh {
61149bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
61249bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ    *b = (Mat_MPIKAIJ*)A->data;
61349bd79ccSDebojyoti Ghosh 
61449bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
61549bd79ccSDebojyoti Ghosh   ierr = MatDestroy(&b->AIJ);CHKERRQ(ierr);
61649bd79ccSDebojyoti Ghosh   ierr = MatDestroy(&b->OAIJ);CHKERRQ(ierr);
61749bd79ccSDebojyoti Ghosh   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
61849bd79ccSDebojyoti Ghosh   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
61949bd79ccSDebojyoti Ghosh   ierr = VecDestroy(&b->w);CHKERRQ(ierr);
620a84f8069SRichard Tran Mills   ierr = PetscFree(b->S);CHKERRQ(ierr);
621a84f8069SRichard Tran Mills   ierr = PetscFree(b->T);CHKERRQ(ierr);
622a84f8069SRichard Tran Mills   ierr = PetscFree(b->ibdiag);CHKERRQ(ierr);
62349bd79ccSDebojyoti Ghosh   ierr = PetscFree(A->data);CHKERRQ(ierr);
62449bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
62549bd79ccSDebojyoti Ghosh }
62649bd79ccSDebojyoti Ghosh 
62749bd79ccSDebojyoti Ghosh /* --------------------------------------------------------------------------------------*/
62849bd79ccSDebojyoti Ghosh 
62949bd79ccSDebojyoti Ghosh /* zz = yy + Axx */
630836168d5SRichard Tran Mills PetscErrorCode MatMultAdd_SeqKAIJ(Mat A,Vec xx,Vec yy,Vec zz)
63149bd79ccSDebojyoti Ghosh {
63249bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ       *b = (Mat_SeqKAIJ*)A->data;
63349bd79ccSDebojyoti Ghosh   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)b->AIJ->data;
63449bd79ccSDebojyoti Ghosh   const PetscScalar *s = b->S, *t = b->T;
63549bd79ccSDebojyoti Ghosh   const PetscScalar *x,*v,*bx;
63649bd79ccSDebojyoti Ghosh   PetscScalar       *y,*sums;
63749bd79ccSDebojyoti Ghosh   PetscErrorCode    ierr;
63849bd79ccSDebojyoti Ghosh   const PetscInt    m = b->AIJ->rmap->n,*idx,*ii;
63949bd79ccSDebojyoti Ghosh   PetscInt          n,i,jrow,j,l,p=b->p,q=b->q,k;
64049bd79ccSDebojyoti Ghosh 
64149bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
64249bd79ccSDebojyoti Ghosh   if (!yy) {
64349bd79ccSDebojyoti Ghosh     ierr = VecSet(zz,0.0);CHKERRQ(ierr);
64449bd79ccSDebojyoti Ghosh   } else {
64549bd79ccSDebojyoti Ghosh     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
64649bd79ccSDebojyoti Ghosh   }
64749bd79ccSDebojyoti Ghosh   if ((!s) && (!t) && (!b->isTI)) PetscFunctionReturn(0);
64849bd79ccSDebojyoti Ghosh 
64949bd79ccSDebojyoti Ghosh   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
65049bd79ccSDebojyoti Ghosh   ierr = VecGetArray(zz,&y);CHKERRQ(ierr);
65149bd79ccSDebojyoti Ghosh   idx  = a->j;
65249bd79ccSDebojyoti Ghosh   v    = a->a;
65349bd79ccSDebojyoti Ghosh   ii   = a->i;
65449bd79ccSDebojyoti Ghosh 
65549bd79ccSDebojyoti Ghosh   if (b->isTI) {
65649bd79ccSDebojyoti Ghosh     for (i=0; i<m; i++) {
65749bd79ccSDebojyoti Ghosh       jrow = ii[i];
65849bd79ccSDebojyoti Ghosh       n    = ii[i+1] - jrow;
65949bd79ccSDebojyoti Ghosh       sums = y + p*i;
66049bd79ccSDebojyoti Ghosh       for (j=0; j<n; j++) {
66149bd79ccSDebojyoti Ghosh         for (k=0; k<p; k++) {
66249bd79ccSDebojyoti Ghosh           sums[k] += v[jrow+j]*x[q*idx[jrow+j]+k];
66349bd79ccSDebojyoti Ghosh         }
66449bd79ccSDebojyoti Ghosh       }
66549bd79ccSDebojyoti Ghosh     }
666234d9204SRichard Tran Mills     ierr = PetscLogFlops((a->nz)*3*p);CHKERRQ(ierr);
66749bd79ccSDebojyoti Ghosh   } else if (t) {
66849bd79ccSDebojyoti Ghosh     for (i=0; i<m; i++) {
66949bd79ccSDebojyoti Ghosh       jrow = ii[i];
67049bd79ccSDebojyoti Ghosh       n    = ii[i+1] - jrow;
67149bd79ccSDebojyoti Ghosh       sums = y + p*i;
67249bd79ccSDebojyoti Ghosh       for (j=0; j<n; j++) {
67349bd79ccSDebojyoti Ghosh         for (k=0; k<p; k++) {
67449bd79ccSDebojyoti Ghosh           for (l=0; l<q; l++) {
67549bd79ccSDebojyoti Ghosh             sums[k] += v[jrow+j]*t[k+l*p]*x[q*idx[jrow+j]+l];
67649bd79ccSDebojyoti Ghosh           }
67749bd79ccSDebojyoti Ghosh         }
67849bd79ccSDebojyoti Ghosh       }
67949bd79ccSDebojyoti Ghosh     }
680234d9204SRichard Tran Mills     /* The flop count below assumes that v[jrow+j] is hoisted out (which an optimizing compiler is likely to do),
681234d9204SRichard Tran Mills      * and also that T part is hoisted outside this loop (in exchange for temporary storage) as (A \otimes I) (I \otimes T),
682234d9204SRichard Tran Mills      * so that this multiply doesn't have to be redone for each matrix entry, but just once per column. The latter
683234d9204SRichard Tran Mills      * transformation is much less likely to be applied, but we nonetheless count the minimum flops required. */
684234d9204SRichard Tran Mills     ierr = PetscLogFlops((2.0*p*q-p)*m+2*p*a->nz);CHKERRQ(ierr);
68549bd79ccSDebojyoti Ghosh   }
68649bd79ccSDebojyoti Ghosh   if (s) {
68749bd79ccSDebojyoti Ghosh     for (i=0; i<m; i++) {
68849bd79ccSDebojyoti Ghosh       sums = y + p*i;
68949bd79ccSDebojyoti Ghosh       bx   = x + q*i;
69049bd79ccSDebojyoti Ghosh       if (i < b->AIJ->cmap->n) {
69149bd79ccSDebojyoti Ghosh         for (j=0; j<q; j++) {
69249bd79ccSDebojyoti Ghosh           for (k=0; k<p; k++) {
69349bd79ccSDebojyoti Ghosh             sums[k] += s[k+j*p]*bx[j];
69449bd79ccSDebojyoti Ghosh           }
69549bd79ccSDebojyoti Ghosh         }
69649bd79ccSDebojyoti Ghosh       }
69749bd79ccSDebojyoti Ghosh     }
698234d9204SRichard Tran Mills     ierr = PetscLogFlops(m*2*p*q);CHKERRQ(ierr);
69949bd79ccSDebojyoti Ghosh   }
70049bd79ccSDebojyoti Ghosh 
70149bd79ccSDebojyoti Ghosh   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
70249bd79ccSDebojyoti Ghosh   ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr);
70349bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
70449bd79ccSDebojyoti Ghosh }
70549bd79ccSDebojyoti Ghosh 
706bb6fb833SRichard Tran Mills PetscErrorCode MatMult_SeqKAIJ(Mat A,Vec xx,Vec yy)
70749bd79ccSDebojyoti Ghosh {
70849bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
70949bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
710836168d5SRichard Tran Mills   ierr = MatMultAdd_SeqKAIJ(A,xx,PETSC_NULL,yy);CHKERRQ(ierr);
71149bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
71249bd79ccSDebojyoti Ghosh }
71349bd79ccSDebojyoti Ghosh 
71449bd79ccSDebojyoti Ghosh #include <petsc/private/kernels/blockinvert.h>
71549bd79ccSDebojyoti Ghosh 
716bb6fb833SRichard Tran Mills PetscErrorCode MatInvertBlockDiagonal_SeqKAIJ(Mat A,const PetscScalar **values)
71749bd79ccSDebojyoti Ghosh {
71849bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ       *b  = (Mat_SeqKAIJ*)A->data;
71949bd79ccSDebojyoti Ghosh   Mat_SeqAIJ        *a  = (Mat_SeqAIJ*)b->AIJ->data;
72049bd79ccSDebojyoti Ghosh   const PetscScalar *S  = b->S;
72149bd79ccSDebojyoti Ghosh   const PetscScalar *T  = b->T;
72249bd79ccSDebojyoti Ghosh   const PetscScalar *v  = a->a;
72349bd79ccSDebojyoti Ghosh   const PetscInt     p  = b->p, q = b->q, m = b->AIJ->rmap->n, *idx = a->j, *ii = a->i;
72449bd79ccSDebojyoti Ghosh   PetscErrorCode    ierr;
72549bd79ccSDebojyoti Ghosh   PetscInt          i,j,*v_pivots,dof,dof2;
72649bd79ccSDebojyoti Ghosh   PetscScalar       *diag,aval,*v_work;
72749bd79ccSDebojyoti Ghosh 
72849bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
72949bd79ccSDebojyoti Ghosh   if (p != q) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATKAIJ: Block size must be square to calculate inverse.");
73031a97b9aSRichard Tran Mills   if ((!S) && (!T) && (!b->isTI)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATKAIJ: Cannot invert a zero matrix.");
73149bd79ccSDebojyoti Ghosh 
73249bd79ccSDebojyoti Ghosh   dof  = p;
73349bd79ccSDebojyoti Ghosh   dof2 = dof*dof;
73449bd79ccSDebojyoti Ghosh 
73549bd79ccSDebojyoti Ghosh   if (b->ibdiagvalid) {
73649bd79ccSDebojyoti Ghosh     if (values) *values = b->ibdiag;
73749bd79ccSDebojyoti Ghosh     PetscFunctionReturn(0);
73849bd79ccSDebojyoti Ghosh   }
73949bd79ccSDebojyoti Ghosh   if (!b->ibdiag) {
740a84f8069SRichard Tran Mills     ierr = PetscMalloc1(dof2*m*sizeof(PetscScalar),&b->ibdiag);CHKERRQ(ierr);
74149bd79ccSDebojyoti Ghosh     ierr = PetscLogObjectMemory((PetscObject)A,dof2*m*sizeof(PetscScalar));CHKERRQ(ierr);
74249bd79ccSDebojyoti Ghosh   }
74349bd79ccSDebojyoti Ghosh   if (values) *values = b->ibdiag;
74449bd79ccSDebojyoti Ghosh   diag = b->ibdiag;
74549bd79ccSDebojyoti Ghosh 
74649bd79ccSDebojyoti Ghosh   ierr = PetscMalloc2(dof,&v_work,dof,&v_pivots);CHKERRQ(ierr);
74749bd79ccSDebojyoti Ghosh   for (i=0; i<m; i++) {
74849bd79ccSDebojyoti Ghosh     if (S) {
74949bd79ccSDebojyoti Ghosh       ierr = PetscMemcpy(diag,S,dof2*sizeof(PetscScalar));CHKERRQ(ierr);
75049bd79ccSDebojyoti Ghosh     } else {
75149bd79ccSDebojyoti Ghosh       ierr = PetscMemzero(diag,dof2*sizeof(PetscScalar));CHKERRQ(ierr);
75249bd79ccSDebojyoti Ghosh     }
75349bd79ccSDebojyoti Ghosh     if (b->isTI) {
75449bd79ccSDebojyoti Ghosh       aval = 0;
75549bd79ccSDebojyoti Ghosh       for (j=ii[i]; j<ii[i+1]; j++) if (idx[j] == i) aval = v[j];
75649bd79ccSDebojyoti Ghosh       for (j=0; j<dof; j++) diag[j+dof*j] += aval;
75749bd79ccSDebojyoti Ghosh     } else if (T) {
75849bd79ccSDebojyoti Ghosh       aval = 0;
75949bd79ccSDebojyoti Ghosh       for (j=ii[i]; j<ii[i+1]; j++) if (idx[j] == i) aval = v[j];
76049bd79ccSDebojyoti Ghosh       for (j=0; j<dof2; j++) diag[j] += aval*T[j];
76149bd79ccSDebojyoti Ghosh     }
76249bd79ccSDebojyoti Ghosh     ierr = PetscKernel_A_gets_inverse_A(dof,diag,v_pivots,v_work,PETSC_FALSE,NULL);CHKERRQ(ierr);
76349bd79ccSDebojyoti Ghosh     diag += dof2;
76449bd79ccSDebojyoti Ghosh   }
76549bd79ccSDebojyoti Ghosh   ierr = PetscFree2(v_work,v_pivots);CHKERRQ(ierr);
76649bd79ccSDebojyoti Ghosh 
76749bd79ccSDebojyoti Ghosh   b->ibdiagvalid = PETSC_TRUE;
76849bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
76949bd79ccSDebojyoti Ghosh }
77049bd79ccSDebojyoti Ghosh 
77149bd79ccSDebojyoti Ghosh static PetscErrorCode MatGetDiagonalBlock_MPIKAIJ(Mat A,Mat *B)
77249bd79ccSDebojyoti Ghosh {
77349bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ *kaij = (Mat_MPIKAIJ*) A->data;
77449bd79ccSDebojyoti Ghosh 
77549bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
77649bd79ccSDebojyoti Ghosh   *B = kaij->AIJ;
77749bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
77849bd79ccSDebojyoti Ghosh }
77949bd79ccSDebojyoti Ghosh 
78049bd79ccSDebojyoti Ghosh PetscErrorCode MatSOR_SeqKAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
78149bd79ccSDebojyoti Ghosh {
78249bd79ccSDebojyoti Ghosh   PetscErrorCode    ierr;
78349bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ       *kaij = (Mat_SeqKAIJ*) A->data;
78449bd79ccSDebojyoti Ghosh   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)kaij->AIJ->data;
78549bd79ccSDebojyoti Ghosh   const PetscScalar *aa = a->a, *T = kaij->T, *v;
78649bd79ccSDebojyoti Ghosh   const PetscInt    m  = kaij->AIJ->rmap->n, *ai=a->i, *aj=a->j, p = kaij->p, q = kaij->q, *diag, *vi;
78749bd79ccSDebojyoti Ghosh   const PetscScalar *b, *xb, *idiag;
78849bd79ccSDebojyoti Ghosh   PetscScalar       *x, *work, *workt, *w, *y, *arr, *t, *arrt;
78949bd79ccSDebojyoti Ghosh   PetscInt          i, j, k, i2, bs, bs2, nz;
79049bd79ccSDebojyoti Ghosh 
79149bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
79249bd79ccSDebojyoti Ghosh   its = its*lits;
79349bd79ccSDebojyoti Ghosh   if (flag & SOR_EISENSTAT) SETERRQ (PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
79449bd79ccSDebojyoti Ghosh   if (its <= 0)             SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
7956a375485SRichard Tran Mills   if (fshift)               SETERRQ (PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for diagonal shift");
7966a375485SRichard Tran Mills   if ((flag & SOR_APPLY_UPPER) || (flag & SOR_APPLY_LOWER)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for applying upper or lower triangular parts");
7976a375485SRichard Tran Mills   if (p != q) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSOR for KAIJ: No support for non-square dense blocks");
79849bd79ccSDebojyoti Ghosh   else        {bs = p; bs2 = bs*bs; }
79949bd79ccSDebojyoti Ghosh 
80049bd79ccSDebojyoti Ghosh   if (!m) PetscFunctionReturn(0);
80149bd79ccSDebojyoti Ghosh 
802bb6fb833SRichard Tran Mills   if (!kaij->ibdiagvalid) { ierr = MatInvertBlockDiagonal_SeqKAIJ(A,NULL);CHKERRQ(ierr); }
80349bd79ccSDebojyoti Ghosh   idiag = kaij->ibdiag;
80449bd79ccSDebojyoti Ghosh   diag  = a->diag;
80549bd79ccSDebojyoti Ghosh 
80649bd79ccSDebojyoti Ghosh   if (!kaij->sor.setup) {
80749bd79ccSDebojyoti Ghosh     ierr = PetscMalloc5(bs,&kaij->sor.w,bs,&kaij->sor.y,m*bs,&kaij->sor.work,m*bs,&kaij->sor.t,m*bs2,&kaij->sor.arr);CHKERRQ(ierr);
80849bd79ccSDebojyoti Ghosh     kaij->sor.setup = PETSC_TRUE;
80949bd79ccSDebojyoti Ghosh   }
81049bd79ccSDebojyoti Ghosh   y     = kaij->sor.y;
81149bd79ccSDebojyoti Ghosh   w     = kaij->sor.w;
81249bd79ccSDebojyoti Ghosh   work  = kaij->sor.work;
81349bd79ccSDebojyoti Ghosh   t     = kaij->sor.t;
81449bd79ccSDebojyoti Ghosh   arr   = kaij->sor.arr;
81549bd79ccSDebojyoti Ghosh 
81649bd79ccSDebojyoti Ghosh   ierr = VecGetArray(xx,&x);    CHKERRQ(ierr);
81749bd79ccSDebojyoti Ghosh   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
81849bd79ccSDebojyoti Ghosh 
81949bd79ccSDebojyoti Ghosh   if (flag & SOR_ZERO_INITIAL_GUESS) {
82049bd79ccSDebojyoti Ghosh     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
82149bd79ccSDebojyoti Ghosh       PetscKernel_w_gets_Ar_times_v(bs,bs,b,idiag,x);                            /* x[0:bs] <- D^{-1} b[0:bs] */
82249bd79ccSDebojyoti Ghosh       ierr   =  PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr);
82349bd79ccSDebojyoti Ghosh       i2     =  bs;
82449bd79ccSDebojyoti Ghosh       idiag  += bs2;
82549bd79ccSDebojyoti Ghosh       for (i=1; i<m; i++) {
82649bd79ccSDebojyoti Ghosh         v  = aa + ai[i];
82749bd79ccSDebojyoti Ghosh         vi = aj + ai[i];
82849bd79ccSDebojyoti Ghosh         nz = diag[i] - ai[i];
82949bd79ccSDebojyoti Ghosh 
83049bd79ccSDebojyoti Ghosh         if (T) {                /* b - T (Arow * x) */
8312ae760e3SRichard Tran Mills           ierr = PetscMemzero(w,bs*sizeof(PetscScalar));CHKERRQ(ierr);
83249bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
83349bd79ccSDebojyoti Ghosh             for (k=0; k<bs; k++) w[k] -= v[j] * x[vi[j]*bs+k];
83449bd79ccSDebojyoti Ghosh           }
83549bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs,w,T,&t[i2]);
83649bd79ccSDebojyoti Ghosh           for (k=0; k<bs; k++) t[i2+k] += b[i2+k];
83749bd79ccSDebojyoti Ghosh         } else if (kaij->isTI) {
8389eb573c1SRichard Tran Mills           ierr = PetscMemcpy(t+i2,b+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr);
83949bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
84049bd79ccSDebojyoti Ghosh             for (k=0; k<bs; k++) t[i2+k] -= v[j] * x[vi[j]*bs+k];
84149bd79ccSDebojyoti Ghosh           }
84249bd79ccSDebojyoti Ghosh         } else {
8439eb573c1SRichard Tran Mills           ierr = PetscMemcpy(t+i2,b+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr);
84449bd79ccSDebojyoti Ghosh         }
84549bd79ccSDebojyoti Ghosh 
84649bd79ccSDebojyoti Ghosh         PetscKernel_w_gets_Ar_times_v(bs,bs,t+i2,idiag,y);
84749bd79ccSDebojyoti Ghosh         for (j=0; j<bs; j++) x[i2+j] = omega * y[j];
84849bd79ccSDebojyoti Ghosh 
84949bd79ccSDebojyoti Ghosh         idiag += bs2;
85049bd79ccSDebojyoti Ghosh         i2    += bs;
85149bd79ccSDebojyoti Ghosh       }
85249bd79ccSDebojyoti Ghosh       /* for logging purposes assume number of nonzero in lower half is 1/2 of total */
85349bd79ccSDebojyoti Ghosh       ierr = PetscLogFlops(1.0*bs2*a->nz);CHKERRQ(ierr);
85449bd79ccSDebojyoti Ghosh       xb = t;
85549bd79ccSDebojyoti Ghosh     } else xb = b;
85649bd79ccSDebojyoti Ghosh     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
85749bd79ccSDebojyoti Ghosh       idiag = kaij->ibdiag+bs2*(m-1);
85849bd79ccSDebojyoti Ghosh       i2    = bs * (m-1);
85949bd79ccSDebojyoti Ghosh       ierr  = PetscMemcpy(w,xb+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr);
86049bd79ccSDebojyoti Ghosh       PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,x+i2);
86149bd79ccSDebojyoti Ghosh       i2    -= bs;
86249bd79ccSDebojyoti Ghosh       idiag -= bs2;
86349bd79ccSDebojyoti Ghosh       for (i=m-2; i>=0; i--) {
86449bd79ccSDebojyoti Ghosh         v  = aa + diag[i] + 1 ;
86549bd79ccSDebojyoti Ghosh         vi = aj + diag[i] + 1;
86649bd79ccSDebojyoti Ghosh         nz = ai[i+1] - diag[i] - 1;
86749bd79ccSDebojyoti Ghosh 
86849bd79ccSDebojyoti Ghosh         if (T) {                /* FIXME: This branch untested */
86949bd79ccSDebojyoti Ghosh           ierr = PetscMemcpy(w,xb+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr);
87049bd79ccSDebojyoti Ghosh           /* copy all rows of x that are needed into contiguous space */
87149bd79ccSDebojyoti Ghosh           workt = work;
87249bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
87349bd79ccSDebojyoti Ghosh             ierr   = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr);
87449bd79ccSDebojyoti Ghosh             workt += bs;
87549bd79ccSDebojyoti Ghosh           }
87649bd79ccSDebojyoti Ghosh           arrt = arr;
87749bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
87849bd79ccSDebojyoti Ghosh             ierr  = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
87949bd79ccSDebojyoti Ghosh             for (k=0; k<bs2; k++) arrt[k] *= v[j];
88049bd79ccSDebojyoti Ghosh             arrt += bs2;
88149bd79ccSDebojyoti Ghosh           }
88249bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
88349bd79ccSDebojyoti Ghosh         } else if (kaij->isTI) {
8849eb573c1SRichard Tran Mills           ierr = PetscMemcpy(w,t+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr);
88549bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
88649bd79ccSDebojyoti Ghosh             for (k=0; k<bs; k++) w[k] -= v[j] * x[vi[j]*bs+k];
88749bd79ccSDebojyoti Ghosh           }
88849bd79ccSDebojyoti Ghosh         }
88949bd79ccSDebojyoti Ghosh 
89049bd79ccSDebojyoti Ghosh         PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y); /* RHS incorrect for omega != 1.0 */
89149bd79ccSDebojyoti Ghosh         for (j=0; j<bs; j++) x[i2+j] = (1.0-omega) * x[i2+j] + omega * y[j];
89249bd79ccSDebojyoti Ghosh 
89349bd79ccSDebojyoti Ghosh         idiag -= bs2;
89449bd79ccSDebojyoti Ghosh         i2    -= bs;
89549bd79ccSDebojyoti Ghosh       }
89649bd79ccSDebojyoti Ghosh       ierr = PetscLogFlops(1.0*bs2*(a->nz));CHKERRQ(ierr);
89749bd79ccSDebojyoti Ghosh     }
89849bd79ccSDebojyoti Ghosh     its--;
89949bd79ccSDebojyoti Ghosh   }
90049bd79ccSDebojyoti Ghosh   while (its--) {               /* FIXME: This branch not updated */
90149bd79ccSDebojyoti Ghosh     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
90249bd79ccSDebojyoti Ghosh       i2     =  0;
90349bd79ccSDebojyoti Ghosh       idiag  = kaij->ibdiag;
90449bd79ccSDebojyoti Ghosh       for (i=0; i<m; i++) {
90549bd79ccSDebojyoti Ghosh         ierr = PetscMemcpy(w,b+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr);
90649bd79ccSDebojyoti Ghosh 
90749bd79ccSDebojyoti Ghosh         v  = aa + ai[i];
90849bd79ccSDebojyoti Ghosh         vi = aj + ai[i];
90949bd79ccSDebojyoti Ghosh         nz = diag[i] - ai[i];
91049bd79ccSDebojyoti Ghosh         workt = work;
91149bd79ccSDebojyoti Ghosh         for (j=0; j<nz; j++) {
91249bd79ccSDebojyoti Ghosh           ierr   = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr);
91349bd79ccSDebojyoti Ghosh           workt += bs;
91449bd79ccSDebojyoti Ghosh         }
91549bd79ccSDebojyoti Ghosh         arrt = arr;
91649bd79ccSDebojyoti Ghosh         if (T) {
91749bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
91849bd79ccSDebojyoti Ghosh             ierr  = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
91949bd79ccSDebojyoti Ghosh             for (k=0; k<bs2; k++) arrt[k] *= v[j];
92049bd79ccSDebojyoti Ghosh             arrt += bs2;
92149bd79ccSDebojyoti Ghosh           }
92249bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
92349bd79ccSDebojyoti Ghosh         } else if (kaij->isTI) {
92449bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
92549bd79ccSDebojyoti Ghosh             ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
92649bd79ccSDebojyoti Ghosh             for (k=0; k<bs; k++) arrt[k+bs*k] = v[j];
92749bd79ccSDebojyoti Ghosh             arrt += bs2;
92849bd79ccSDebojyoti Ghosh           }
92949bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
93049bd79ccSDebojyoti Ghosh         }
93149bd79ccSDebojyoti Ghosh         ierr = PetscMemcpy(t+i2,w,bs*sizeof(PetscScalar));CHKERRQ(ierr);
93249bd79ccSDebojyoti Ghosh 
93349bd79ccSDebojyoti Ghosh         v  = aa + diag[i] + 1;
93449bd79ccSDebojyoti Ghosh         vi = aj + diag[i] + 1;
93549bd79ccSDebojyoti Ghosh         nz = ai[i+1] - diag[i] - 1;
93649bd79ccSDebojyoti Ghosh         workt = work;
93749bd79ccSDebojyoti Ghosh         for (j=0; j<nz; j++) {
93849bd79ccSDebojyoti Ghosh           ierr   = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr);
93949bd79ccSDebojyoti Ghosh           workt += bs;
94049bd79ccSDebojyoti Ghosh         }
94149bd79ccSDebojyoti Ghosh         arrt = arr;
94249bd79ccSDebojyoti Ghosh         if (T) {
94349bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
94449bd79ccSDebojyoti Ghosh             ierr  = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
94549bd79ccSDebojyoti Ghosh             for (k=0; k<bs2; k++) arrt[k] *= v[j];
94649bd79ccSDebojyoti Ghosh             arrt += bs2;
94749bd79ccSDebojyoti Ghosh           }
94849bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
94949bd79ccSDebojyoti Ghosh         } else if (kaij->isTI) {
95049bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
95149bd79ccSDebojyoti Ghosh             ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
95249bd79ccSDebojyoti Ghosh             for (k=0; k<bs; k++) arrt[k+bs*k] = v[j];
95349bd79ccSDebojyoti Ghosh             arrt += bs2;
95449bd79ccSDebojyoti Ghosh           }
95549bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
95649bd79ccSDebojyoti Ghosh         }
95749bd79ccSDebojyoti Ghosh 
95849bd79ccSDebojyoti Ghosh         PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y);
95949bd79ccSDebojyoti Ghosh         for (j=0; j<bs; j++) *(x+i2+j) = (1.0-omega) * *(x+i2+j) + omega * *(y+j);
96049bd79ccSDebojyoti Ghosh 
96149bd79ccSDebojyoti Ghosh         idiag += bs2;
96249bd79ccSDebojyoti Ghosh         i2    += bs;
96349bd79ccSDebojyoti Ghosh       }
96449bd79ccSDebojyoti Ghosh       xb = t;
96549bd79ccSDebojyoti Ghosh     }
96649bd79ccSDebojyoti Ghosh     else xb = b;
96749bd79ccSDebojyoti Ghosh     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
96849bd79ccSDebojyoti Ghosh       idiag = kaij->ibdiag+bs2*(m-1);
96949bd79ccSDebojyoti Ghosh       i2    = bs * (m-1);
97049bd79ccSDebojyoti Ghosh       if (xb == b) {
97149bd79ccSDebojyoti Ghosh         for (i=m-1; i>=0; i--) {
97249bd79ccSDebojyoti Ghosh           ierr = PetscMemcpy(w,b+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr);
97349bd79ccSDebojyoti Ghosh 
97449bd79ccSDebojyoti Ghosh           v  = aa + ai[i];
97549bd79ccSDebojyoti Ghosh           vi = aj + ai[i];
97649bd79ccSDebojyoti Ghosh           nz = diag[i] - ai[i];
97749bd79ccSDebojyoti Ghosh           workt = work;
97849bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
97949bd79ccSDebojyoti Ghosh             ierr   = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr);
98049bd79ccSDebojyoti Ghosh             workt += bs;
98149bd79ccSDebojyoti Ghosh           }
98249bd79ccSDebojyoti Ghosh           arrt = arr;
98349bd79ccSDebojyoti Ghosh           if (T) {
98449bd79ccSDebojyoti Ghosh             for (j=0; j<nz; j++) {
98549bd79ccSDebojyoti Ghosh               ierr  = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
98649bd79ccSDebojyoti Ghosh               for (k=0; k<bs2; k++) arrt[k] *= v[j];
98749bd79ccSDebojyoti Ghosh               arrt += bs2;
98849bd79ccSDebojyoti Ghosh             }
98949bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
99049bd79ccSDebojyoti Ghosh           } else if (kaij->isTI) {
99149bd79ccSDebojyoti Ghosh             for (j=0; j<nz; j++) {
99249bd79ccSDebojyoti Ghosh               ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
99349bd79ccSDebojyoti Ghosh               for (k=0; k<bs; k++) arrt[k+bs*k] = v[j];
99449bd79ccSDebojyoti Ghosh               arrt += bs2;
99549bd79ccSDebojyoti Ghosh             }
99649bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
99749bd79ccSDebojyoti Ghosh           }
99849bd79ccSDebojyoti Ghosh 
99949bd79ccSDebojyoti Ghosh           v  = aa + diag[i] + 1;
100049bd79ccSDebojyoti Ghosh           vi = aj + diag[i] + 1;
100149bd79ccSDebojyoti Ghosh           nz = ai[i+1] - diag[i] - 1;
100249bd79ccSDebojyoti Ghosh           workt = work;
100349bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
100449bd79ccSDebojyoti Ghosh             ierr   = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr);
100549bd79ccSDebojyoti Ghosh             workt += bs;
100649bd79ccSDebojyoti Ghosh           }
100749bd79ccSDebojyoti Ghosh           arrt = arr;
100849bd79ccSDebojyoti Ghosh           if (T) {
100949bd79ccSDebojyoti Ghosh             for (j=0; j<nz; j++) {
101049bd79ccSDebojyoti Ghosh               ierr  = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
101149bd79ccSDebojyoti Ghosh               for (k=0; k<bs2; k++) arrt[k] *= v[j];
101249bd79ccSDebojyoti Ghosh               arrt += bs2;
101349bd79ccSDebojyoti Ghosh             }
101449bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
101549bd79ccSDebojyoti Ghosh           } else if (kaij->isTI) {
101649bd79ccSDebojyoti Ghosh             for (j=0; j<nz; j++) {
101749bd79ccSDebojyoti Ghosh               ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
101849bd79ccSDebojyoti Ghosh               for (k=0; k<bs; k++) arrt[k+bs*k] = v[j];
101949bd79ccSDebojyoti Ghosh               arrt += bs2;
102049bd79ccSDebojyoti Ghosh             }
102149bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
102249bd79ccSDebojyoti Ghosh           }
102349bd79ccSDebojyoti Ghosh 
102449bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y);
102549bd79ccSDebojyoti Ghosh           for (j=0; j<bs; j++) *(x+i2+j) = (1.0-omega) * *(x+i2+j) + omega * *(y+j);
102649bd79ccSDebojyoti Ghosh         }
102749bd79ccSDebojyoti Ghosh       } else {
102849bd79ccSDebojyoti Ghosh         for (i=m-1; i>=0; i--) {
102949bd79ccSDebojyoti Ghosh           ierr = PetscMemcpy(w,xb+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr);
103049bd79ccSDebojyoti Ghosh           v  = aa + diag[i] + 1;
103149bd79ccSDebojyoti Ghosh           vi = aj + diag[i] + 1;
103249bd79ccSDebojyoti Ghosh           nz = ai[i+1] - diag[i] - 1;
103349bd79ccSDebojyoti Ghosh           workt = work;
103449bd79ccSDebojyoti Ghosh           for (j=0; j<nz; j++) {
103549bd79ccSDebojyoti Ghosh             ierr   = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr);
103649bd79ccSDebojyoti Ghosh             workt += bs;
103749bd79ccSDebojyoti Ghosh           }
103849bd79ccSDebojyoti Ghosh           arrt = arr;
103949bd79ccSDebojyoti Ghosh           if (T) {
104049bd79ccSDebojyoti Ghosh             for (j=0; j<nz; j++) {
104149bd79ccSDebojyoti Ghosh               ierr  = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
104249bd79ccSDebojyoti Ghosh               for (k=0; k<bs2; k++) arrt[k] *= v[j];
104349bd79ccSDebojyoti Ghosh               arrt += bs2;
104449bd79ccSDebojyoti Ghosh             }
104549bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
104649bd79ccSDebojyoti Ghosh           } else if (kaij->isTI) {
104749bd79ccSDebojyoti Ghosh             for (j=0; j<nz; j++) {
104849bd79ccSDebojyoti Ghosh               ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr);
104949bd79ccSDebojyoti Ghosh               for (k=0; k<bs; k++) arrt[k+bs*k] = v[j];
105049bd79ccSDebojyoti Ghosh               arrt += bs2;
105149bd79ccSDebojyoti Ghosh             }
105249bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work);
105349bd79ccSDebojyoti Ghosh           }
105449bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y);
105549bd79ccSDebojyoti Ghosh           for (j=0; j<bs; j++) *(x+i2+j) = (1.0-omega) * *(x+i2+j) + omega * *(y+j);
105649bd79ccSDebojyoti Ghosh         }
105749bd79ccSDebojyoti Ghosh       }
105849bd79ccSDebojyoti Ghosh       ierr = PetscLogFlops(1.0*bs2*(a->nz));CHKERRQ(ierr);
105949bd79ccSDebojyoti Ghosh     }
106049bd79ccSDebojyoti Ghosh   }
106149bd79ccSDebojyoti Ghosh 
106249bd79ccSDebojyoti Ghosh   ierr = VecRestoreArray(xx,&x);    CHKERRQ(ierr);
106349bd79ccSDebojyoti Ghosh   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
106449bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
106549bd79ccSDebojyoti Ghosh }
106649bd79ccSDebojyoti Ghosh 
106749bd79ccSDebojyoti Ghosh /*===================================================================================*/
106849bd79ccSDebojyoti Ghosh 
1069836168d5SRichard Tran Mills PetscErrorCode MatMultAdd_MPIKAIJ(Mat A,Vec xx,Vec yy,Vec zz)
107049bd79ccSDebojyoti Ghosh {
107149bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ    *b = (Mat_MPIKAIJ*)A->data;
107249bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
107349bd79ccSDebojyoti Ghosh 
107449bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
107549bd79ccSDebojyoti Ghosh   if (!yy) {
107649bd79ccSDebojyoti Ghosh     ierr = VecSet(zz,0.0);CHKERRQ(ierr);
107749bd79ccSDebojyoti Ghosh   } else {
107849bd79ccSDebojyoti Ghosh     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
107949bd79ccSDebojyoti Ghosh   }
108049bd79ccSDebojyoti Ghosh   /* start the scatter */
108149bd79ccSDebojyoti Ghosh   ierr = VecScatterBegin(b->ctx,xx,b->w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
108249bd79ccSDebojyoti Ghosh   ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,zz,zz);CHKERRQ(ierr);
108349bd79ccSDebojyoti Ghosh   ierr = VecScatterEnd(b->ctx,xx,b->w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
108449bd79ccSDebojyoti Ghosh   ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,zz,zz);CHKERRQ(ierr);
108549bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
108649bd79ccSDebojyoti Ghosh }
108749bd79ccSDebojyoti Ghosh 
1088bb6fb833SRichard Tran Mills PetscErrorCode MatMult_MPIKAIJ(Mat A,Vec xx,Vec yy)
108949bd79ccSDebojyoti Ghosh {
109049bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
109149bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
1092836168d5SRichard Tran Mills   ierr = MatMultAdd_MPIKAIJ(A,xx,PETSC_NULL,yy);CHKERRQ(ierr);
109349bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
109449bd79ccSDebojyoti Ghosh }
109549bd79ccSDebojyoti Ghosh 
1096bb6fb833SRichard Tran Mills PetscErrorCode MatInvertBlockDiagonal_MPIKAIJ(Mat A,const PetscScalar **values)
109749bd79ccSDebojyoti Ghosh {
109849bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ     *b = (Mat_MPIKAIJ*)A->data;
109949bd79ccSDebojyoti Ghosh   PetscErrorCode  ierr;
110049bd79ccSDebojyoti Ghosh 
110149bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
110249bd79ccSDebojyoti Ghosh   ierr = (*b->AIJ->ops->invertblockdiagonal)(b->AIJ,values);CHKERRQ(ierr);
110349bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
110449bd79ccSDebojyoti Ghosh }
110549bd79ccSDebojyoti Ghosh 
110649bd79ccSDebojyoti Ghosh /* ----------------------------------------------------------------*/
110749bd79ccSDebojyoti Ghosh 
110849bd79ccSDebojyoti Ghosh PetscErrorCode MatGetRow_SeqKAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **values)
110949bd79ccSDebojyoti Ghosh {
111049bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ     *b   = (Mat_SeqKAIJ*) A->data;
11111ca5ffdbSRichard Tran Mills   PetscErrorCode  diag = PETSC_FALSE;
11121ca5ffdbSRichard Tran Mills   PetscErrorCode  ierr;
111349bd79ccSDebojyoti Ghosh   PetscInt        nzaij,nz,*colsaij,*idx,i,j,p=b->p,q=b->q,r=row/p,s=row%p,c;
111449bd79ccSDebojyoti Ghosh   PetscScalar     *vaij,*v,*S=b->S,*T=b->T;
111549bd79ccSDebojyoti Ghosh 
111649bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
111749bd79ccSDebojyoti Ghosh   if (b->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
111849bd79ccSDebojyoti Ghosh   b->getrowactive = PETSC_TRUE;
111949bd79ccSDebojyoti Ghosh   if (row < 0 || row >= A->rmap->n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D out of range",row);
112049bd79ccSDebojyoti Ghosh 
112149bd79ccSDebojyoti Ghosh   if ((!S) && (!T) && (!b->isTI)) {
112249bd79ccSDebojyoti Ghosh     if (ncols)    *ncols  = 0;
112349bd79ccSDebojyoti Ghosh     if (cols)     *cols   = NULL;
112449bd79ccSDebojyoti Ghosh     if (values)   *values = NULL;
112549bd79ccSDebojyoti Ghosh     PetscFunctionReturn(0);
112649bd79ccSDebojyoti Ghosh   }
112749bd79ccSDebojyoti Ghosh 
112849bd79ccSDebojyoti Ghosh   if (T || b->isTI) {
112949bd79ccSDebojyoti Ghosh     ierr  = MatGetRow_SeqAIJ(b->AIJ,r,&nzaij,&colsaij,&vaij);CHKERRQ(ierr);
113049bd79ccSDebojyoti Ghosh     c     = nzaij;
113149bd79ccSDebojyoti Ghosh     for (i=0; i<nzaij; i++) {
113249bd79ccSDebojyoti Ghosh       /* check if this row contains a diagonal entry */
113349bd79ccSDebojyoti Ghosh       if (colsaij[i] == r) {
113449bd79ccSDebojyoti Ghosh         diag = PETSC_TRUE;
113549bd79ccSDebojyoti Ghosh         c = i;
113649bd79ccSDebojyoti Ghosh       }
113749bd79ccSDebojyoti Ghosh     }
113849bd79ccSDebojyoti Ghosh   } else nzaij = c = 0;
113949bd79ccSDebojyoti Ghosh 
114049bd79ccSDebojyoti Ghosh   /* calculate size of row */
114149bd79ccSDebojyoti Ghosh   nz = 0;
114249bd79ccSDebojyoti Ghosh   if (S)            nz += q;
114349bd79ccSDebojyoti Ghosh   if (T || b->isTI) nz += (diag && S ? (nzaij-1)*q : nzaij*q);
114449bd79ccSDebojyoti Ghosh 
114549bd79ccSDebojyoti Ghosh   if (cols || values) {
114649bd79ccSDebojyoti Ghosh     ierr = PetscMalloc2(nz,&idx,nz,&v);CHKERRQ(ierr);
114738822f9dSRichard Tran Mills     for (i=0; i<q; i++) {
114838822f9dSRichard Tran Mills       /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */
114938822f9dSRichard Tran Mills       v[i] = 0.0;
115038822f9dSRichard Tran Mills     }
115149bd79ccSDebojyoti Ghosh     if (b->isTI) {
115249bd79ccSDebojyoti Ghosh       for (i=0; i<nzaij; i++) {
115349bd79ccSDebojyoti Ghosh         for (j=0; j<q; j++) {
115449bd79ccSDebojyoti Ghosh           idx[i*q+j] = colsaij[i]*q+j;
115549bd79ccSDebojyoti Ghosh           v[i*q+j]   = (j==s ? vaij[i] : 0);
115649bd79ccSDebojyoti Ghosh         }
115749bd79ccSDebojyoti Ghosh       }
115849bd79ccSDebojyoti Ghosh     } else if (T) {
115949bd79ccSDebojyoti Ghosh       for (i=0; i<nzaij; i++) {
116049bd79ccSDebojyoti Ghosh         for (j=0; j<q; j++) {
116149bd79ccSDebojyoti Ghosh           idx[i*q+j] = colsaij[i]*q+j;
116249bd79ccSDebojyoti Ghosh           v[i*q+j]   = vaij[i]*T[s+j*p];
116349bd79ccSDebojyoti Ghosh         }
116449bd79ccSDebojyoti Ghosh       }
116549bd79ccSDebojyoti Ghosh     }
116649bd79ccSDebojyoti Ghosh     if (S) {
116749bd79ccSDebojyoti Ghosh       for (j=0; j<q; j++) {
116849bd79ccSDebojyoti Ghosh         idx[c*q+j] = r*q+j;
116949bd79ccSDebojyoti Ghosh         v[c*q+j]  += S[s+j*p];
117049bd79ccSDebojyoti Ghosh       }
117149bd79ccSDebojyoti Ghosh     }
117249bd79ccSDebojyoti Ghosh   }
117349bd79ccSDebojyoti Ghosh 
117449bd79ccSDebojyoti Ghosh   if (ncols)    *ncols  = nz;
117549bd79ccSDebojyoti Ghosh   if (cols)     *cols   = idx;
117649bd79ccSDebojyoti Ghosh   if (values)   *values = v;
117749bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
117849bd79ccSDebojyoti Ghosh }
117949bd79ccSDebojyoti Ghosh 
118049bd79ccSDebojyoti Ghosh PetscErrorCode MatRestoreRow_SeqKAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
118149bd79ccSDebojyoti Ghosh {
118249bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
118349bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
118449bd79ccSDebojyoti Ghosh   ierr = PetscFree2(*idx,*v);CHKERRQ(ierr);
118549bd79ccSDebojyoti Ghosh   ((Mat_SeqKAIJ*)A->data)->getrowactive = PETSC_FALSE;
118649bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
118749bd79ccSDebojyoti Ghosh }
118849bd79ccSDebojyoti Ghosh 
118949bd79ccSDebojyoti Ghosh PetscErrorCode MatGetRow_MPIKAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **values)
119049bd79ccSDebojyoti Ghosh {
119149bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ     *b      = (Mat_MPIKAIJ*) A->data;
119249bd79ccSDebojyoti Ghosh   Mat             MatAIJ  = ((Mat_SeqKAIJ*)b->AIJ->data)->AIJ;
119349bd79ccSDebojyoti Ghosh   Mat             MatOAIJ = ((Mat_SeqKAIJ*)b->OAIJ->data)->AIJ;
119449bd79ccSDebojyoti Ghosh   Mat             AIJ     = b->A;
1195fc64b2cfSRichard Tran Mills   PetscBool       diag    = PETSC_FALSE;
119649bd79ccSDebojyoti Ghosh   PetscErrorCode  ierr;
119749bd79ccSDebojyoti Ghosh   const PetscInt  rstart=A->rmap->rstart,rend=A->rmap->rend,p=b->p,q=b->q,*garray;
119849bd79ccSDebojyoti Ghosh   PetscInt        nz,*idx,ncolsaij,ncolsoaij,*colsaij,*colsoaij,r,s,c,i,j,lrow;
119949bd79ccSDebojyoti Ghosh   PetscScalar     *v,*vals,*ovals,*S=b->S,*T=b->T;
120049bd79ccSDebojyoti Ghosh 
120149bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
120249bd79ccSDebojyoti Ghosh   if (b->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
120349bd79ccSDebojyoti Ghosh   b->getrowactive = PETSC_TRUE;
120449bd79ccSDebojyoti Ghosh   if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only local rows");
120549bd79ccSDebojyoti Ghosh   lrow = row - rstart;
120649bd79ccSDebojyoti Ghosh 
120749bd79ccSDebojyoti Ghosh   if ((!S) && (!T) && (!b->isTI)) {
120849bd79ccSDebojyoti Ghosh     if (ncols)    *ncols  = 0;
120949bd79ccSDebojyoti Ghosh     if (cols)     *cols   = NULL;
121049bd79ccSDebojyoti Ghosh     if (values)   *values = NULL;
121149bd79ccSDebojyoti Ghosh     PetscFunctionReturn(0);
121249bd79ccSDebojyoti Ghosh   }
121349bd79ccSDebojyoti Ghosh 
121449bd79ccSDebojyoti Ghosh   r = lrow/p;
121549bd79ccSDebojyoti Ghosh   s = lrow%p;
121649bd79ccSDebojyoti Ghosh 
121749bd79ccSDebojyoti Ghosh   if (T || b->isTI) {
12182ae760e3SRichard Tran Mills     ierr = MatMPIAIJGetSeqAIJ(AIJ,NULL,NULL,&garray);CHKERRQ(ierr);
121949bd79ccSDebojyoti Ghosh     ierr = MatGetRow_SeqAIJ(MatAIJ,lrow/p,&ncolsaij,&colsaij,&vals);CHKERRQ(ierr);
122049bd79ccSDebojyoti Ghosh     ierr = MatGetRow_SeqAIJ(MatOAIJ,lrow/p,&ncolsoaij,&colsoaij,&ovals);CHKERRQ(ierr);
122149bd79ccSDebojyoti Ghosh 
122249bd79ccSDebojyoti Ghosh     c     = ncolsaij + ncolsoaij;
122349bd79ccSDebojyoti Ghosh     for (i=0; i<ncolsaij; i++) {
122449bd79ccSDebojyoti Ghosh       /* check if this row contains a diagonal entry */
122549bd79ccSDebojyoti Ghosh       if (colsaij[i] == r) {
122649bd79ccSDebojyoti Ghosh         diag = PETSC_TRUE;
122749bd79ccSDebojyoti Ghosh         c = i;
122849bd79ccSDebojyoti Ghosh       }
122949bd79ccSDebojyoti Ghosh     }
123049bd79ccSDebojyoti Ghosh   } else c = 0;
123149bd79ccSDebojyoti Ghosh 
123249bd79ccSDebojyoti Ghosh   /* calculate size of row */
123349bd79ccSDebojyoti Ghosh   nz = 0;
123449bd79ccSDebojyoti Ghosh   if (S)            nz += q;
123549bd79ccSDebojyoti Ghosh   if (T || b->isTI) nz += (diag && S ? (ncolsaij+ncolsoaij-1)*q : (ncolsaij+ncolsoaij)*q);
123649bd79ccSDebojyoti Ghosh 
123749bd79ccSDebojyoti Ghosh   if (cols || values) {
123849bd79ccSDebojyoti Ghosh     ierr = PetscMalloc2(nz,&idx,nz,&v);CHKERRQ(ierr);
1239a437a796SRichard Tran Mills     for (i=0; i<q; i++) {
1240a437a796SRichard Tran Mills       /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */
1241a437a796SRichard Tran Mills       v[i] = 0.0;
1242a437a796SRichard Tran Mills     }
124349bd79ccSDebojyoti Ghosh     if (b->isTI) {
124449bd79ccSDebojyoti Ghosh       for (i=0; i<ncolsaij; i++) {
124549bd79ccSDebojyoti Ghosh         for (j=0; j<q; j++) {
124649bd79ccSDebojyoti Ghosh           idx[i*q+j] = (colsaij[i]+rstart/p)*q+j;
124749bd79ccSDebojyoti Ghosh           v[i*q+j]   = (j==s ? vals[i] : 0.0);
124849bd79ccSDebojyoti Ghosh         }
124949bd79ccSDebojyoti Ghosh       }
125049bd79ccSDebojyoti Ghosh       for (i=0; i<ncolsoaij; i++) {
125149bd79ccSDebojyoti Ghosh         for (j=0; j<q; j++) {
125249bd79ccSDebojyoti Ghosh           idx[(i+ncolsaij)*q+j] = garray[colsoaij[i]]*q+j;
125349bd79ccSDebojyoti Ghosh           v[(i+ncolsaij)*q+j]   = (j==s ? ovals[i]: 0.0);
125449bd79ccSDebojyoti Ghosh         }
125549bd79ccSDebojyoti Ghosh       }
125649bd79ccSDebojyoti Ghosh     } else if (T) {
125749bd79ccSDebojyoti Ghosh       for (i=0; i<ncolsaij; i++) {
125849bd79ccSDebojyoti Ghosh         for (j=0; j<q; j++) {
125949bd79ccSDebojyoti Ghosh           idx[i*q+j] = (colsaij[i]+rstart/p)*q+j;
126049bd79ccSDebojyoti Ghosh           v[i*q+j]   = vals[i]*T[s+j*p];
126149bd79ccSDebojyoti Ghosh         }
126249bd79ccSDebojyoti Ghosh       }
126349bd79ccSDebojyoti Ghosh       for (i=0; i<ncolsoaij; i++) {
126449bd79ccSDebojyoti Ghosh         for (j=0; j<q; j++) {
126549bd79ccSDebojyoti Ghosh           idx[(i+ncolsaij)*q+j] = garray[colsoaij[i]]*q+j;
126649bd79ccSDebojyoti Ghosh           v[(i+ncolsaij)*q+j]   = ovals[i]*T[s+j*p];
126749bd79ccSDebojyoti Ghosh         }
126849bd79ccSDebojyoti Ghosh       }
126949bd79ccSDebojyoti Ghosh     }
127049bd79ccSDebojyoti Ghosh     if (S) {
127149bd79ccSDebojyoti Ghosh       for (j=0; j<q; j++) {
127249bd79ccSDebojyoti Ghosh         idx[c*q+j] = (r+rstart/p)*q+j;
127349bd79ccSDebojyoti Ghosh         v[c*q+j]  += S[s+j*p];
127449bd79ccSDebojyoti Ghosh       }
127549bd79ccSDebojyoti Ghosh     }
127649bd79ccSDebojyoti Ghosh   }
127749bd79ccSDebojyoti Ghosh 
127849bd79ccSDebojyoti Ghosh   if (ncols)  *ncols  = nz;
127949bd79ccSDebojyoti Ghosh   if (cols)   *cols   = idx;
128049bd79ccSDebojyoti Ghosh   if (values) *values = v;
128149bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
128249bd79ccSDebojyoti Ghosh }
128349bd79ccSDebojyoti Ghosh 
128449bd79ccSDebojyoti Ghosh PetscErrorCode MatRestoreRow_MPIKAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
128549bd79ccSDebojyoti Ghosh {
128649bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
128749bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
128849bd79ccSDebojyoti Ghosh   ierr = PetscFree2(*idx,*v);CHKERRQ(ierr);
128949bd79ccSDebojyoti Ghosh   ((Mat_SeqKAIJ*)A->data)->getrowactive = PETSC_FALSE;
129049bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
129149bd79ccSDebojyoti Ghosh }
129249bd79ccSDebojyoti Ghosh 
129349bd79ccSDebojyoti Ghosh PetscErrorCode  MatCreateSubMatrix_KAIJ(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat)
129449bd79ccSDebojyoti Ghosh {
129549bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
129649bd79ccSDebojyoti Ghosh   Mat            A;
129749bd79ccSDebojyoti Ghosh 
129849bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
129949bd79ccSDebojyoti Ghosh   ierr = MatConvert(mat,MATAIJ,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
130049bd79ccSDebojyoti Ghosh   ierr = MatCreateSubMatrix(A,isrow,iscol,cll,newmat);CHKERRQ(ierr);
130149bd79ccSDebojyoti Ghosh   ierr = MatDestroy(&A);CHKERRQ(ierr);
130249bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
130349bd79ccSDebojyoti Ghosh }
130449bd79ccSDebojyoti Ghosh 
130549bd79ccSDebojyoti Ghosh /* ---------------------------------------------------------------------------------- */
130649bd79ccSDebojyoti Ghosh /*@C
130749bd79ccSDebojyoti Ghosh   MatCreateKAIJ - Creates a matrix type to be used for matrices of the following form:
130849bd79ccSDebojyoti Ghosh 
130949bd79ccSDebojyoti Ghosh     [I \otimes S + A \otimes T]
131049bd79ccSDebojyoti Ghosh 
131149bd79ccSDebojyoti Ghosh   where
131249bd79ccSDebojyoti Ghosh     S is a dense (p \times q) matrix
131349bd79ccSDebojyoti Ghosh     T is a dense (p \times q) matrix
131449bd79ccSDebojyoti Ghosh     A is an AIJ  (n \times n) matrix
131549bd79ccSDebojyoti Ghosh     I is the identity matrix
131649bd79ccSDebojyoti Ghosh   The resulting matrix is (np \times nq)
131749bd79ccSDebojyoti Ghosh 
1318d3b90ce1SRichard Tran Mills   S and T are always stored independently on all processes as PetscScalar arrays in column-major format.
131949bd79ccSDebojyoti Ghosh 
132049bd79ccSDebojyoti Ghosh   Collective
132149bd79ccSDebojyoti Ghosh 
132249bd79ccSDebojyoti Ghosh   Input Parameters:
132349bd79ccSDebojyoti Ghosh + A - the AIJ matrix
132449bd79ccSDebojyoti Ghosh . p - number of rows in S and T
1325d3b90ce1SRichard Tran Mills . q - number of columns in S and T
1326d3b90ce1SRichard Tran Mills . S - the S matrix (can be PETSC_NULL), stored as a PetscScalar array (column-major)
1327d3b90ce1SRichard Tran Mills - T - the T matrix (can be PETSC_NULL), stored as a PetscScalar array (column-major)
132849bd79ccSDebojyoti Ghosh 
132949bd79ccSDebojyoti Ghosh   Output Parameter:
133049bd79ccSDebojyoti Ghosh . kaij - the new KAIJ matrix
133149bd79ccSDebojyoti Ghosh 
1332d3b90ce1SRichard Tran Mills   Notes:
1333d3b90ce1SRichard Tran Mills   This function increases the reference count on the AIJ matrix, so the user is free to destroy the matrix if it is not needed.
1334d3b90ce1SRichard Tran Mills   Changes to the entries of the AIJ matrix will immediately affect the KAIJ matrix.
133549bd79ccSDebojyoti Ghosh 
133649bd79ccSDebojyoti Ghosh   Level: advanced
133749bd79ccSDebojyoti Ghosh 
13380567c835SRichard Tran Mills .seealso: MatKAIJSetAIJ(), MatKAIJSetS(), MatKAIJSetT(), MatKAIJGetAIJ(), MatKAIJGetS(), MatKAIJGetT(), MATKAIJ
133949bd79ccSDebojyoti Ghosh @*/
134049bd79ccSDebojyoti Ghosh PetscErrorCode  MatCreateKAIJ(Mat A,PetscInt p,PetscInt q,const PetscScalar S[],const PetscScalar T[],Mat *kaij)
134149bd79ccSDebojyoti Ghosh {
134249bd79ccSDebojyoti Ghosh   PetscErrorCode ierr;
134349bd79ccSDebojyoti Ghosh   PetscMPIInt    size;
134449bd79ccSDebojyoti Ghosh 
134549bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
13460567c835SRichard Tran Mills   ierr = MatCreate(PetscObjectComm((PetscObject)A),kaij);CHKERRQ(ierr);
134749bd79ccSDebojyoti Ghosh   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
134849bd79ccSDebojyoti Ghosh   if (size == 1) {
13490567c835SRichard Tran Mills     ierr = MatSetType(*kaij,MATSEQKAIJ);CHKERRQ(ierr);
135049bd79ccSDebojyoti Ghosh   } else {
13510567c835SRichard Tran Mills     ierr = MatSetType(*kaij,MATMPIKAIJ);CHKERRQ(ierr);
135249bd79ccSDebojyoti Ghosh   }
13530567c835SRichard Tran Mills   ierr = MatKAIJSetAIJ(*kaij,A);CHKERRQ(ierr);
13540567c835SRichard Tran Mills   ierr = MatKAIJSetS(*kaij,p,q,S);CHKERRQ(ierr);
13550567c835SRichard Tran Mills   ierr = MatKAIJSetT(*kaij,p,q,T);CHKERRQ(ierr);
13562ae760e3SRichard Tran Mills   ierr = MatSetUp(*kaij);CHKERRQ(ierr);
13570567c835SRichard Tran Mills   PetscFunctionReturn(0);
13580567c835SRichard Tran Mills }
135949bd79ccSDebojyoti Ghosh 
13600567c835SRichard Tran Mills /*MC
13615881e567SRichard Tran Mills   MATKAIJ - MATKAIJ = "kaij" - A matrix type to be used to evaluate matrices of form
13625881e567SRichard Tran Mills     [I \otimes S + A \otimes T],
13630567c835SRichard Tran Mills   where
13645881e567SRichard Tran Mills     S is a dense (p \times q) matrix,
13655881e567SRichard Tran Mills     T is a dense (p \times q) matrix,
13665881e567SRichard Tran Mills     A is an AIJ  (n \times n) matrix,
13675881e567SRichard Tran Mills     and I is the identity matrix.
13685881e567SRichard Tran Mills   The resulting matrix is (np \times nq).
13690567c835SRichard Tran Mills 
1370d3b90ce1SRichard Tran Mills   S and T are always stored independently on all processes as PetscScalar arrays in column-major format.
13710567c835SRichard Tran Mills 
13725881e567SRichard Tran Mills   Notes:
13735881e567SRichard Tran Mills   A linear system with multiple right-hand sides, AX = B, can be expressed in the KAIJ-friendly form of (A \otimes I) x = b,
13745881e567SRichard Tran Mills   where x and b are column vectors containing the row-major representations of X and B.
13755881e567SRichard Tran Mills 
13760567c835SRichard Tran Mills   Level: advanced
13770567c835SRichard Tran Mills 
13780567c835SRichard Tran Mills .seealso: MatKAIJSetAIJ(), MatKAIJSetS(), MatKAIJSetT(), MatKAIJGetAIJ(), MatKAIJGetS(), MatKAIJGetT(), MatCreateKAIJ()
13790567c835SRichard Tran Mills M*/
13800567c835SRichard Tran Mills 
13810567c835SRichard Tran Mills PETSC_EXTERN PetscErrorCode MatCreate_KAIJ(Mat A)
13820567c835SRichard Tran Mills {
13830567c835SRichard Tran Mills   PetscErrorCode ierr;
13840567c835SRichard Tran Mills   Mat_MPIKAIJ    *b;
13850567c835SRichard Tran Mills   PetscMPIInt    size;
13860567c835SRichard Tran Mills 
13870567c835SRichard Tran Mills   PetscFunctionBegin;
13880567c835SRichard Tran Mills   ierr     = PetscNewLog(A,&b);CHKERRQ(ierr);
13890567c835SRichard Tran Mills   A->data  = (void*)b;
13900567c835SRichard Tran Mills 
13910567c835SRichard Tran Mills   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
13920567c835SRichard Tran Mills 
13930567c835SRichard Tran Mills   A->ops->setup = MatSetUp_KAIJ;
13940567c835SRichard Tran Mills 
13950567c835SRichard Tran Mills   b->w    = 0;
13960567c835SRichard Tran Mills   ierr    = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
13970567c835SRichard Tran Mills   if (size == 1) {
13980567c835SRichard Tran Mills     ierr = PetscObjectChangeTypeName((PetscObject)A,MATSEQKAIJ);CHKERRQ(ierr);
13990567c835SRichard Tran Mills     A->ops->setup               = MatSetUp_KAIJ;
14000567c835SRichard Tran Mills     A->ops->destroy             = MatDestroy_SeqKAIJ;
1401e6985dafSRichard Tran Mills     A->ops->view                = MatView_KAIJ;
1402bb6fb833SRichard Tran Mills     A->ops->mult                = MatMult_SeqKAIJ;
1403bb6fb833SRichard Tran Mills     A->ops->multadd             = MatMultAdd_SeqKAIJ;
1404bb6fb833SRichard Tran Mills     A->ops->invertblockdiagonal = MatInvertBlockDiagonal_SeqKAIJ;
14050567c835SRichard Tran Mills     A->ops->getrow              = MatGetRow_SeqKAIJ;
14060567c835SRichard Tran Mills     A->ops->restorerow          = MatRestoreRow_SeqKAIJ;
14070567c835SRichard Tran Mills     A->ops->sor                 = MatSOR_SeqKAIJ;
14080567c835SRichard Tran Mills   } else {
14090567c835SRichard Tran Mills     ierr = PetscObjectChangeTypeName((PetscObject)A,MATMPIKAIJ);CHKERRQ(ierr);
14100567c835SRichard Tran Mills     A->ops->setup               = MatSetUp_KAIJ;
14110567c835SRichard Tran Mills     A->ops->destroy             = MatDestroy_MPIKAIJ;
1412e6985dafSRichard Tran Mills     A->ops->view                = MatView_KAIJ;
1413bb6fb833SRichard Tran Mills     A->ops->mult                = MatMult_MPIKAIJ;
1414bb6fb833SRichard Tran Mills     A->ops->multadd             = MatMultAdd_MPIKAIJ;
1415bb6fb833SRichard Tran Mills     A->ops->invertblockdiagonal = MatInvertBlockDiagonal_MPIKAIJ;
14160567c835SRichard Tran Mills     A->ops->getrow              = MatGetRow_MPIKAIJ;
14170567c835SRichard Tran Mills     A->ops->restorerow          = MatRestoreRow_MPIKAIJ;
14180567c835SRichard Tran Mills     ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPIKAIJ);CHKERRQ(ierr);
14190567c835SRichard Tran Mills   }
14200567c835SRichard Tran Mills   A->ops->createsubmatrix = MatCreateSubMatrix_KAIJ;
142149bd79ccSDebojyoti Ghosh   PetscFunctionReturn(0);
142249bd79ccSDebojyoti Ghosh }
1423