xref: /petsc/src/mat/impls/kaij/kaij.c (revision ac530a7e429a3ef5a9263376acf6071236a5db98)
149bd79ccSDebojyoti Ghosh /*
249bd79ccSDebojyoti Ghosh   Defines the basic matrix operations for the KAIJ  matrix storage format.
349bd79ccSDebojyoti Ghosh   This format is used to evaluate matrices of the form:
449bd79ccSDebojyoti Ghosh 
549bd79ccSDebojyoti Ghosh     [I \otimes S + A \otimes T]
649bd79ccSDebojyoti Ghosh 
749bd79ccSDebojyoti Ghosh   where
849bd79ccSDebojyoti Ghosh     S is a dense (p \times q) matrix
949bd79ccSDebojyoti Ghosh     T is a dense (p \times q) matrix
1049bd79ccSDebojyoti Ghosh     A is an AIJ  (n \times n) matrix
1149bd79ccSDebojyoti Ghosh     I is the identity matrix
1249bd79ccSDebojyoti Ghosh 
1349bd79ccSDebojyoti Ghosh   The resulting matrix is (np \times nq)
1449bd79ccSDebojyoti Ghosh 
1549bd79ccSDebojyoti Ghosh   We provide:
1649bd79ccSDebojyoti Ghosh      MatMult()
1749bd79ccSDebojyoti Ghosh      MatMultAdd()
1849bd79ccSDebojyoti Ghosh      MatInvertBlockDiagonal()
1949bd79ccSDebojyoti Ghosh   and
2049bd79ccSDebojyoti Ghosh      MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*)
2149bd79ccSDebojyoti Ghosh 
2249bd79ccSDebojyoti Ghosh   This single directory handles both the sequential and parallel codes
2349bd79ccSDebojyoti Ghosh */
2449bd79ccSDebojyoti Ghosh 
2549bd79ccSDebojyoti Ghosh #include <../src/mat/impls/kaij/kaij.h> /*I "petscmat.h" I*/
2649bd79ccSDebojyoti Ghosh #include <../src/mat/utils/freespace.h>
2749bd79ccSDebojyoti Ghosh #include <petsc/private/vecimpl.h>
2849bd79ccSDebojyoti Ghosh 
29cc4c1da9SBarry Smith /*@
3011a5261eSBarry Smith   MatKAIJGetAIJ - Get the `MATAIJ` matrix describing the blockwise action of the `MATKAIJ` matrix
3149bd79ccSDebojyoti Ghosh 
3211a5261eSBarry Smith   Not Collective, but if the `MATKAIJ` matrix is parallel, the `MATAIJ` matrix is also parallel
3349bd79ccSDebojyoti Ghosh 
3449bd79ccSDebojyoti Ghosh   Input Parameter:
3511a5261eSBarry Smith . A - the `MATKAIJ` matrix
3649bd79ccSDebojyoti Ghosh 
3749bd79ccSDebojyoti Ghosh   Output Parameter:
3811a5261eSBarry Smith . B - the `MATAIJ` matrix
3949bd79ccSDebojyoti Ghosh 
4049bd79ccSDebojyoti Ghosh   Level: advanced
4149bd79ccSDebojyoti Ghosh 
4211a5261eSBarry Smith   Note:
4311a5261eSBarry Smith   The reference count on the `MATAIJ` matrix is not increased so you should not destroy it.
4449bd79ccSDebojyoti Ghosh 
451cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateKAIJ()`, `MATKAIJ`, `MATAIJ`
4649bd79ccSDebojyoti Ghosh @*/
47d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetAIJ(Mat A, Mat *B)
48d71ae5a4SJacob Faibussowitsch {
4949bd79ccSDebojyoti Ghosh   PetscBool ismpikaij, isseqkaij;
5049bd79ccSDebojyoti Ghosh 
5149bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIKAIJ, &ismpikaij));
539566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQKAIJ, &isseqkaij));
5449bd79ccSDebojyoti Ghosh   if (ismpikaij) {
5549bd79ccSDebojyoti Ghosh     Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
5649bd79ccSDebojyoti Ghosh 
5749bd79ccSDebojyoti Ghosh     *B = b->A;
5849bd79ccSDebojyoti Ghosh   } else if (isseqkaij) {
5949bd79ccSDebojyoti Ghosh     Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data;
6049bd79ccSDebojyoti Ghosh 
6149bd79ccSDebojyoti Ghosh     *B = b->AIJ;
62b04351cbSRichard Tran Mills   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix passed in is not of type KAIJ");
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6449bd79ccSDebojyoti Ghosh }
6549bd79ccSDebojyoti Ghosh 
6649bd79ccSDebojyoti Ghosh /*@C
672ef1f0ffSBarry Smith   MatKAIJGetS - Get the `S` matrix describing the shift action of the `MATKAIJ` matrix
6849bd79ccSDebojyoti Ghosh 
692ef1f0ffSBarry Smith   Not Collective; the entire `S` is stored and returned independently on all processes.
7049bd79ccSDebojyoti Ghosh 
7149bd79ccSDebojyoti Ghosh   Input Parameter:
7211a5261eSBarry Smith . A - the `MATKAIJ` matrix
7349bd79ccSDebojyoti Ghosh 
74a5b5c723SRichard Tran Mills   Output Parameters:
752ef1f0ffSBarry Smith + m - the number of rows in `S`
762ef1f0ffSBarry Smith . n - the number of columns in `S`
77a5b5c723SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format
7849bd79ccSDebojyoti Ghosh 
7949bd79ccSDebojyoti Ghosh   Level: advanced
8049bd79ccSDebojyoti Ghosh 
812ef1f0ffSBarry Smith   Note:
822ef1f0ffSBarry Smith   All output parameters are optional (pass `NULL` if not desired)
832ef1f0ffSBarry Smith 
841cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()`
8549bd79ccSDebojyoti Ghosh @*/
865d83a8b1SBarry Smith PetscErrorCode MatKAIJGetS(Mat A, PetscInt *m, PetscInt *n, PetscScalar *S[])
87d71ae5a4SJacob Faibussowitsch {
8849bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data;
894d86920dSPierre Jolivet 
9049bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
91a5b5c723SRichard Tran Mills   if (m) *m = b->p;
92a5b5c723SRichard Tran Mills   if (n) *n = b->q;
93a5b5c723SRichard Tran Mills   if (S) *S = b->S;
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95a5b5c723SRichard Tran Mills }
96a5b5c723SRichard Tran Mills 
97a5b5c723SRichard Tran Mills /*@C
982ef1f0ffSBarry Smith   MatKAIJGetSRead - Get a read-only pointer to the `S` matrix describing the shift action of the `MATKAIJ` matrix
99a5b5c723SRichard Tran Mills 
1002ef1f0ffSBarry Smith   Not Collective; the entire `S` is stored and returned independently on all processes.
101a5b5c723SRichard Tran Mills 
102a5b5c723SRichard Tran Mills   Input Parameter:
10311a5261eSBarry Smith . A - the `MATKAIJ` matrix
104a5b5c723SRichard Tran Mills 
105a5b5c723SRichard Tran Mills   Output Parameters:
1062ef1f0ffSBarry Smith + m - the number of rows in `S`
1072ef1f0ffSBarry Smith . n - the number of columns in `S`
108a5b5c723SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format
109a5b5c723SRichard Tran Mills 
110a5b5c723SRichard Tran Mills   Level: advanced
111a5b5c723SRichard Tran Mills 
1122ef1f0ffSBarry Smith   Note:
1132ef1f0ffSBarry Smith   All output parameters are optional (pass `NULL` if not desired)
1142ef1f0ffSBarry Smith 
1151cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()`
116a5b5c723SRichard Tran Mills @*/
1175d83a8b1SBarry Smith PetscErrorCode MatKAIJGetSRead(Mat A, PetscInt *m, PetscInt *n, const PetscScalar *S[])
118d71ae5a4SJacob Faibussowitsch {
119a5b5c723SRichard Tran Mills   Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data;
1204d86920dSPierre Jolivet 
121a5b5c723SRichard Tran Mills   PetscFunctionBegin;
122a5b5c723SRichard Tran Mills   if (m) *m = b->p;
123a5b5c723SRichard Tran Mills   if (n) *n = b->q;
124a5b5c723SRichard Tran Mills   if (S) *S = b->S;
1253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126a5b5c723SRichard Tran Mills }
127a5b5c723SRichard Tran Mills 
128a5b5c723SRichard Tran Mills /*@C
12911a5261eSBarry Smith   MatKAIJRestoreS - Restore array obtained with `MatKAIJGetS()`
130a5b5c723SRichard Tran Mills 
1312ef1f0ffSBarry Smith   Not Collective
132a5b5c723SRichard Tran Mills 
1332ef1f0ffSBarry Smith   Input Parameters:
1342ef1f0ffSBarry Smith + A - the `MATKAIJ` matrix
1352ef1f0ffSBarry Smith - S - location of pointer to array obtained with `MatKAIJGetS()`
136a5b5c723SRichard Tran Mills 
137a5b5c723SRichard Tran Mills   Level: advanced
13811a5261eSBarry Smith 
1392ef1f0ffSBarry Smith   Note:
1402ef1f0ffSBarry Smith   This routine zeros the array pointer to prevent accidental reuse after it has been restored.
1412ef1f0ffSBarry Smith   If `NULL` is passed, it will not attempt to zero the array pointer.
1422ef1f0ffSBarry Smith 
1431cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJGetSRead()`, `MatKAIJRestoreSRead()`
144a5b5c723SRichard Tran Mills @*/
1455d83a8b1SBarry Smith PetscErrorCode MatKAIJRestoreS(Mat A, PetscScalar *S[])
146d71ae5a4SJacob Faibussowitsch {
147a5b5c723SRichard Tran Mills   PetscFunctionBegin;
148a5b5c723SRichard Tran Mills   if (S) *S = NULL;
1499566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)A));
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151a5b5c723SRichard Tran Mills }
152a5b5c723SRichard Tran Mills 
153a5b5c723SRichard Tran Mills /*@C
15411a5261eSBarry Smith   MatKAIJRestoreSRead - Restore array obtained with `MatKAIJGetSRead()`
155a5b5c723SRichard Tran Mills 
1562ef1f0ffSBarry Smith   Not Collective
157a5b5c723SRichard Tran Mills 
1582ef1f0ffSBarry Smith   Input Parameters:
1592ef1f0ffSBarry Smith + A - the `MATKAIJ` matrix
1602ef1f0ffSBarry Smith - S - location of pointer to array obtained with `MatKAIJGetS()`
161a5b5c723SRichard Tran Mills 
162a5b5c723SRichard Tran Mills   Level: advanced
16311a5261eSBarry Smith 
1642ef1f0ffSBarry Smith   Note:
1652ef1f0ffSBarry Smith   This routine zeros the array pointer to prevent accidental reuse after it has been restored.
1662ef1f0ffSBarry Smith   If `NULL` is passed, it will not attempt to zero the array pointer.
1672ef1f0ffSBarry Smith 
16842747ad1SJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJGetSRead()`
169a5b5c723SRichard Tran Mills @*/
1705d83a8b1SBarry Smith PetscErrorCode MatKAIJRestoreSRead(Mat A, const PetscScalar *S[])
171d71ae5a4SJacob Faibussowitsch {
172a5b5c723SRichard Tran Mills   PetscFunctionBegin;
173a5b5c723SRichard Tran Mills   if (S) *S = NULL;
1743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17549bd79ccSDebojyoti Ghosh }
17649bd79ccSDebojyoti Ghosh 
17749bd79ccSDebojyoti Ghosh /*@C
1782ef1f0ffSBarry Smith   MatKAIJGetT - Get the transformation matrix `T` associated with the `MATKAIJ` matrix
17949bd79ccSDebojyoti Ghosh 
1802ef1f0ffSBarry Smith   Not Collective; the entire `T` is stored and returned independently on all processes
18149bd79ccSDebojyoti Ghosh 
18249bd79ccSDebojyoti Ghosh   Input Parameter:
18311a5261eSBarry Smith . A - the `MATKAIJ` matrix
18449bd79ccSDebojyoti Ghosh 
185d8d19677SJose E. Roman   Output Parameters:
1862ef1f0ffSBarry Smith + m - the number of rows in `T`
1872ef1f0ffSBarry Smith . n - the number of columns in `T`
188a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format
18949bd79ccSDebojyoti Ghosh 
19049bd79ccSDebojyoti Ghosh   Level: advanced
19149bd79ccSDebojyoti Ghosh 
1922ef1f0ffSBarry Smith   Note:
1932ef1f0ffSBarry Smith   All output parameters are optional (pass `NULL` if not desired)
1942ef1f0ffSBarry Smith 
1951cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()`
19649bd79ccSDebojyoti Ghosh @*/
1975d83a8b1SBarry Smith PetscErrorCode MatKAIJGetT(Mat A, PetscInt *m, PetscInt *n, PetscScalar *T[])
198d71ae5a4SJacob Faibussowitsch {
19949bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data;
2004d86920dSPierre Jolivet 
20149bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
202a5b5c723SRichard Tran Mills   if (m) *m = b->p;
203a5b5c723SRichard Tran Mills   if (n) *n = b->q;
204a5b5c723SRichard Tran Mills   if (T) *T = b->T;
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
206a5b5c723SRichard Tran Mills }
207a5b5c723SRichard Tran Mills 
208a5b5c723SRichard Tran Mills /*@C
2092ef1f0ffSBarry Smith   MatKAIJGetTRead - Get a read-only pointer to the transformation matrix `T` associated with the `MATKAIJ` matrix
210a5b5c723SRichard Tran Mills 
2112ef1f0ffSBarry Smith   Not Collective; the entire `T` is stored and returned independently on all processes
212a5b5c723SRichard Tran Mills 
213a5b5c723SRichard Tran Mills   Input Parameter:
21411a5261eSBarry Smith . A - the `MATKAIJ` matrix
215a5b5c723SRichard Tran Mills 
216d8d19677SJose E. Roman   Output Parameters:
2172ef1f0ffSBarry Smith + m - the number of rows in `T`
2182ef1f0ffSBarry Smith . n - the number of columns in `T`
219a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format
220a5b5c723SRichard Tran Mills 
221a5b5c723SRichard Tran Mills   Level: advanced
222a5b5c723SRichard Tran Mills 
2232ef1f0ffSBarry Smith   Note:
2242ef1f0ffSBarry Smith   All output parameters are optional (pass `NULL` if not desired)
2252ef1f0ffSBarry Smith 
2261cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()`
227a5b5c723SRichard Tran Mills @*/
2285d83a8b1SBarry Smith PetscErrorCode MatKAIJGetTRead(Mat A, PetscInt *m, PetscInt *n, const PetscScalar *T[])
229d71ae5a4SJacob Faibussowitsch {
230a5b5c723SRichard Tran Mills   Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data;
2314d86920dSPierre Jolivet 
232a5b5c723SRichard Tran Mills   PetscFunctionBegin;
233a5b5c723SRichard Tran Mills   if (m) *m = b->p;
234a5b5c723SRichard Tran Mills   if (n) *n = b->q;
235a5b5c723SRichard Tran Mills   if (T) *T = b->T;
2363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
237a5b5c723SRichard Tran Mills }
238a5b5c723SRichard Tran Mills 
239a5b5c723SRichard Tran Mills /*@C
24011a5261eSBarry Smith   MatKAIJRestoreT - Restore array obtained with `MatKAIJGetT()`
241a5b5c723SRichard Tran Mills 
2422ef1f0ffSBarry Smith   Not Collective
243a5b5c723SRichard Tran Mills 
2442ef1f0ffSBarry Smith   Input Parameters:
2452ef1f0ffSBarry Smith + A - the `MATKAIJ` matrix
2462ef1f0ffSBarry Smith - T - location of pointer to array obtained with `MatKAIJGetS()`
247a5b5c723SRichard Tran Mills 
248a5b5c723SRichard Tran Mills   Level: advanced
24911a5261eSBarry Smith 
2502ef1f0ffSBarry Smith   Note:
2512ef1f0ffSBarry Smith   This routine zeros the array pointer to prevent accidental reuse after it has been restored.
2522ef1f0ffSBarry Smith   If `NULL` is passed, it will not attempt to zero the array pointer.
2532ef1f0ffSBarry Smith 
2541cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetT()`, `MatKAIJGetTRead()`, `MatKAIJRestoreTRead()`
255a5b5c723SRichard Tran Mills @*/
2565d83a8b1SBarry Smith PetscErrorCode MatKAIJRestoreT(Mat A, PetscScalar *T[])
257d71ae5a4SJacob Faibussowitsch {
258a5b5c723SRichard Tran Mills   PetscFunctionBegin;
259a5b5c723SRichard Tran Mills   if (T) *T = NULL;
2609566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)A));
2613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
262a5b5c723SRichard Tran Mills }
263a5b5c723SRichard Tran Mills 
264a5b5c723SRichard Tran Mills /*@C
26511a5261eSBarry Smith   MatKAIJRestoreTRead - Restore array obtained with `MatKAIJGetTRead()`
266a5b5c723SRichard Tran Mills 
2672ef1f0ffSBarry Smith   Not Collective
268a5b5c723SRichard Tran Mills 
2692ef1f0ffSBarry Smith   Input Parameters:
2702ef1f0ffSBarry Smith + A - the `MATKAIJ` matrix
2712ef1f0ffSBarry Smith - T - location of pointer to array obtained with `MatKAIJGetS()`
272a5b5c723SRichard Tran Mills 
273a5b5c723SRichard Tran Mills   Level: advanced
27411a5261eSBarry Smith 
2752ef1f0ffSBarry Smith   Note:
2762ef1f0ffSBarry Smith   This routine zeros the array pointer to prevent accidental reuse after it has been restored.
2772ef1f0ffSBarry Smith   If `NULL` is passed, it will not attempt to zero the array pointer.
2782ef1f0ffSBarry Smith 
27942747ad1SJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetT()`, `MatKAIJGetTRead()`
280a5b5c723SRichard Tran Mills @*/
2815d83a8b1SBarry Smith PetscErrorCode MatKAIJRestoreTRead(Mat A, const PetscScalar *T[])
282d71ae5a4SJacob Faibussowitsch {
283a5b5c723SRichard Tran Mills   PetscFunctionBegin;
284a5b5c723SRichard Tran Mills   if (T) *T = NULL;
2853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28649bd79ccSDebojyoti Ghosh }
28749bd79ccSDebojyoti Ghosh 
2880567c835SRichard Tran Mills /*@
28911a5261eSBarry Smith   MatKAIJSetAIJ - Set the `MATAIJ` matrix describing the blockwise action of the `MATKAIJ` matrix
2900567c835SRichard Tran Mills 
29111a5261eSBarry Smith   Logically Collective; if the `MATAIJ` matrix is parallel, the `MATKAIJ` matrix is also parallel
2920567c835SRichard Tran Mills 
2930567c835SRichard Tran Mills   Input Parameters:
29411a5261eSBarry Smith + A - the `MATKAIJ` matrix
29511a5261eSBarry Smith - B - the `MATAIJ` matrix
2960567c835SRichard Tran Mills 
2972ef1f0ffSBarry Smith   Level: advanced
2982ef1f0ffSBarry Smith 
29915b9d025SRichard Tran Mills   Notes:
30011a5261eSBarry Smith   This function increases the reference count on the `MATAIJ` matrix, so the user is free to destroy the matrix if it is not needed.
30111a5261eSBarry Smith 
30211a5261eSBarry Smith   Changes to the entries of the `MATAIJ` matrix will immediately affect the `MATKAIJ` matrix.
30315b9d025SRichard Tran Mills 
3041cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()`
3050567c835SRichard Tran Mills @*/
306d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJSetAIJ(Mat A, Mat B)
307d71ae5a4SJacob Faibussowitsch {
3080567c835SRichard Tran Mills   PetscMPIInt size;
30906d61a37SPierre Jolivet   PetscBool   flg;
3100567c835SRichard Tran Mills 
3110567c835SRichard Tran Mills   PetscFunctionBegin;
3129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
3130567c835SRichard Tran Mills   if (size == 1) {
3149566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &flg));
31528b400f6SJacob Faibussowitsch     PetscCheck(flg, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "MatKAIJSetAIJ() with MATSEQKAIJ does not support %s as the AIJ mat", ((PetscObject)B)->type_name);
3160567c835SRichard Tran Mills     Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data;
3170567c835SRichard Tran Mills     a->AIJ         = B;
3180567c835SRichard Tran Mills   } else {
3190567c835SRichard Tran Mills     Mat_MPIKAIJ *a = (Mat_MPIKAIJ *)A->data;
3200567c835SRichard Tran Mills     a->A           = B;
3210567c835SRichard Tran Mills   }
3229566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)B));
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3240567c835SRichard Tran Mills }
3250567c835SRichard Tran Mills 
326cc4c1da9SBarry Smith /*@
3272ef1f0ffSBarry Smith   MatKAIJSetS - Set the `S` matrix describing the shift action of the `MATKAIJ` matrix
3280567c835SRichard Tran Mills 
3292ef1f0ffSBarry Smith   Logically Collective; the entire `S` is stored independently on all processes.
3300567c835SRichard Tran Mills 
3310567c835SRichard Tran Mills   Input Parameters:
33211a5261eSBarry Smith + A - the `MATKAIJ` matrix
3332ef1f0ffSBarry Smith . p - the number of rows in `S`
3342ef1f0ffSBarry Smith . q - the number of columns in `S`
3350567c835SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format
3360567c835SRichard Tran Mills 
33711a5261eSBarry Smith   Level: advanced
3380567c835SRichard Tran Mills 
3392ef1f0ffSBarry Smith   Notes:
3402ef1f0ffSBarry Smith   The dimensions `p` and `q` must match those of the transformation matrix `T` associated with the `MATKAIJ` matrix.
3412ef1f0ffSBarry Smith 
3422ef1f0ffSBarry Smith   The `S` matrix is copied, so the user can destroy this array.
3432ef1f0ffSBarry Smith 
3441cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJSetT()`, `MatKAIJSetAIJ()`
3450567c835SRichard Tran Mills @*/
346d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJSetS(Mat A, PetscInt p, PetscInt q, const PetscScalar S[])
347d71ae5a4SJacob Faibussowitsch {
3480567c835SRichard Tran Mills   Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data;
3490567c835SRichard Tran Mills 
3500567c835SRichard Tran Mills   PetscFunctionBegin;
3519566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->S));
3520567c835SRichard Tran Mills   if (S) {
3539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(p * q, &a->S));
3549566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(a->S, S, p * q * sizeof(PetscScalar)));
3550567c835SRichard Tran Mills   } else a->S = NULL;
3560567c835SRichard Tran Mills 
3570567c835SRichard Tran Mills   a->p = p;
3580567c835SRichard Tran Mills   a->q = q;
3593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3600567c835SRichard Tran Mills }
3610567c835SRichard Tran Mills 
362cc4c1da9SBarry Smith /*@
3632ef1f0ffSBarry Smith   MatKAIJGetScaledIdentity - Check if both `S` and `T` are scaled identities.
364910cf402Sprj- 
365910cf402Sprj-   Logically Collective.
366910cf402Sprj- 
367910cf402Sprj-   Input Parameter:
36811a5261eSBarry Smith . A - the `MATKAIJ` matrix
369910cf402Sprj- 
370910cf402Sprj-   Output Parameter:
371910cf402Sprj- . identity - the Boolean value
372910cf402Sprj- 
373fe59aa6dSJacob Faibussowitsch   Level: advanced
374910cf402Sprj- 
3751cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJGetT()`
376910cf402Sprj- @*/
377d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetScaledIdentity(Mat A, PetscBool *identity)
378d71ae5a4SJacob Faibussowitsch {
379910cf402Sprj-   Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data;
380910cf402Sprj-   PetscInt     i, j;
381910cf402Sprj- 
382910cf402Sprj-   PetscFunctionBegin;
383910cf402Sprj-   if (a->p != a->q) {
384910cf402Sprj-     *identity = PETSC_FALSE;
3853ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
386910cf402Sprj-   } else *identity = PETSC_TRUE;
387910cf402Sprj-   if (!a->isTI || a->S) {
388910cf402Sprj-     for (i = 0; i < a->p && *identity; i++) {
389910cf402Sprj-       for (j = 0; j < a->p && *identity; j++) {
390910cf402Sprj-         if (i != j) {
391910cf402Sprj-           if (a->S && PetscAbsScalar(a->S[i + j * a->p]) > PETSC_SMALL) *identity = PETSC_FALSE;
392910cf402Sprj-           if (a->T && PetscAbsScalar(a->T[i + j * a->p]) > PETSC_SMALL) *identity = PETSC_FALSE;
393910cf402Sprj-         } else {
394910cf402Sprj-           if (a->S && PetscAbsScalar(a->S[i * (a->p + 1)] - a->S[0]) > PETSC_SMALL) *identity = PETSC_FALSE;
395910cf402Sprj-           if (a->T && PetscAbsScalar(a->T[i * (a->p + 1)] - a->T[0]) > PETSC_SMALL) *identity = PETSC_FALSE;
396910cf402Sprj-         }
397910cf402Sprj-       }
398910cf402Sprj-     }
399910cf402Sprj-   }
4003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
401910cf402Sprj- }
402910cf402Sprj- 
403cc4c1da9SBarry Smith /*@
4042ef1f0ffSBarry Smith   MatKAIJSetT - Set the transformation matrix `T` associated with the `MATKAIJ` matrix
4050567c835SRichard Tran Mills 
4062ef1f0ffSBarry Smith   Logically Collective; the entire `T` is stored independently on all processes.
4070567c835SRichard Tran Mills 
4080567c835SRichard Tran Mills   Input Parameters:
40911a5261eSBarry Smith + A - the `MATKAIJ` matrix
4102ef1f0ffSBarry Smith . p - the number of rows in `S`
4112ef1f0ffSBarry Smith . q - the number of columns in `S`
4122ef1f0ffSBarry Smith - T - the `T` matrix, in form of a scalar array in column-major format
4130567c835SRichard Tran Mills 
414fe59aa6dSJacob Faibussowitsch   Level: advanced
4150567c835SRichard Tran Mills 
4162ef1f0ffSBarry Smith   Notes:
4172ef1f0ffSBarry Smith   The dimensions `p` and `q` must match those of the shift matrix `S` associated with the `MATKAIJ` matrix.
4182ef1f0ffSBarry Smith 
4192ef1f0ffSBarry Smith   The `T` matrix is copied, so the user can destroy this array.
4202ef1f0ffSBarry Smith 
4211cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetT()`, `MatKAIJSetS()`, `MatKAIJSetAIJ()`
4220567c835SRichard Tran Mills @*/
423d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJSetT(Mat A, PetscInt p, PetscInt q, const PetscScalar T[])
424d71ae5a4SJacob Faibussowitsch {
4250567c835SRichard Tran Mills   PetscInt     i, j;
4260567c835SRichard Tran Mills   Mat_SeqKAIJ *a    = (Mat_SeqKAIJ *)A->data;
4270567c835SRichard Tran Mills   PetscBool    isTI = PETSC_FALSE;
4280567c835SRichard Tran Mills 
4290567c835SRichard Tran Mills   PetscFunctionBegin;
4300567c835SRichard Tran Mills   /* check if T is an identity matrix */
4310567c835SRichard Tran Mills   if (T && (p == q)) {
4320567c835SRichard Tran Mills     isTI = PETSC_TRUE;
4330567c835SRichard Tran Mills     for (i = 0; i < p; i++) {
4340567c835SRichard Tran Mills       for (j = 0; j < q; j++) {
4350567c835SRichard Tran Mills         if (i == j) {
4360567c835SRichard Tran Mills           /* diagonal term must be 1 */
4370567c835SRichard Tran Mills           if (T[i + j * p] != 1.0) isTI = PETSC_FALSE;
4380567c835SRichard Tran Mills         } else {
4390567c835SRichard Tran Mills           /* off-diagonal term must be 0 */
4400567c835SRichard Tran Mills           if (T[i + j * p] != 0.0) isTI = PETSC_FALSE;
4410567c835SRichard Tran Mills         }
4420567c835SRichard Tran Mills       }
4430567c835SRichard Tran Mills     }
4440567c835SRichard Tran Mills   }
4450567c835SRichard Tran Mills   a->isTI = isTI;
4460567c835SRichard Tran Mills 
4479566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->T));
4480567c835SRichard Tran Mills   if (T && (!isTI)) {
4499566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(p * q, &a->T));
4509566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(a->T, T, p * q * sizeof(PetscScalar)));
45150d19d74SRichard Tran Mills   } else a->T = NULL;
4520567c835SRichard Tran Mills 
4530567c835SRichard Tran Mills   a->p = p;
4540567c835SRichard Tran Mills   a->q = q;
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4560567c835SRichard Tran Mills }
4570567c835SRichard Tran Mills 
458ff6a9541SJacob Faibussowitsch static PetscErrorCode MatDestroy_SeqKAIJ(Mat A)
459d71ae5a4SJacob Faibussowitsch {
46049bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data;
46149bd79ccSDebojyoti Ghosh 
46249bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
4639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&b->AIJ));
4649566063dSJacob Faibussowitsch   PetscCall(PetscFree(b->S));
4659566063dSJacob Faibussowitsch   PetscCall(PetscFree(b->T));
4669566063dSJacob Faibussowitsch   PetscCall(PetscFree(b->ibdiag));
4679566063dSJacob Faibussowitsch   PetscCall(PetscFree5(b->sor.w, b->sor.y, b->sor.work, b->sor.t, b->sor.arr));
4689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqkaij_seqaij_C", NULL));
4699566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47149bd79ccSDebojyoti Ghosh }
47249bd79ccSDebojyoti Ghosh 
473ff6a9541SJacob Faibussowitsch static PetscErrorCode MatKAIJ_build_AIJ_OAIJ(Mat A)
474d71ae5a4SJacob Faibussowitsch {
475e0e5a793SRichard Tran Mills   Mat_MPIKAIJ     *a;
476e0e5a793SRichard Tran Mills   Mat_MPIAIJ      *mpiaij;
477e0e5a793SRichard Tran Mills   PetscScalar     *T;
478e0e5a793SRichard Tran Mills   PetscInt         i, j;
479e0e5a793SRichard Tran Mills   PetscObjectState state;
480e0e5a793SRichard Tran Mills 
481e0e5a793SRichard Tran Mills   PetscFunctionBegin;
482e0e5a793SRichard Tran Mills   a      = (Mat_MPIKAIJ *)A->data;
483e0e5a793SRichard Tran Mills   mpiaij = (Mat_MPIAIJ *)a->A->data;
484e0e5a793SRichard Tran Mills 
4859566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)a->A, &state));
486e0e5a793SRichard Tran Mills   if (state == a->state) {
487e0e5a793SRichard Tran Mills     /* The existing AIJ and KAIJ members are up-to-date, so simply exit. */
4883ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
489e0e5a793SRichard Tran Mills   } else {
4909566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&a->AIJ));
4919566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&a->OAIJ));
492e0e5a793SRichard Tran Mills     if (a->isTI) {
493e0e5a793SRichard Tran Mills       /* If the transformation matrix associated with the parallel matrix A is the identity matrix, then a->T will be NULL.
494e0e5a793SRichard Tran Mills        * In this case, if we pass a->T directly to the MatCreateKAIJ() calls to create the sequential submatrices, the routine will
495e0e5a793SRichard Tran Mills        * not be able to tell that transformation matrix should be set to the identity; thus we create a temporary identity matrix
496e0e5a793SRichard Tran Mills        * to pass in. */
4979566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(a->p * a->q, &T));
498e0e5a793SRichard Tran Mills       for (i = 0; i < a->p; i++) {
499e0e5a793SRichard Tran Mills         for (j = 0; j < a->q; j++) {
500e0e5a793SRichard Tran Mills           if (i == j) T[i + j * a->p] = 1.0;
501e0e5a793SRichard Tran Mills           else T[i + j * a->p] = 0.0;
502e0e5a793SRichard Tran Mills         }
503e0e5a793SRichard Tran Mills       }
504e0e5a793SRichard Tran Mills     } else T = a->T;
5059566063dSJacob Faibussowitsch     PetscCall(MatCreateKAIJ(mpiaij->A, a->p, a->q, a->S, T, &a->AIJ));
5069566063dSJacob Faibussowitsch     PetscCall(MatCreateKAIJ(mpiaij->B, a->p, a->q, NULL, T, &a->OAIJ));
5071baa6e33SBarry Smith     if (a->isTI) PetscCall(PetscFree(T));
508e0e5a793SRichard Tran Mills     a->state = state;
509e0e5a793SRichard Tran Mills   }
5103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
511e0e5a793SRichard Tran Mills }
512e0e5a793SRichard Tran Mills 
513ff6a9541SJacob Faibussowitsch static PetscErrorCode MatSetUp_KAIJ(Mat A)
514d71ae5a4SJacob Faibussowitsch {
5150567c835SRichard Tran Mills   PetscInt     n;
5160567c835SRichard Tran Mills   PetscMPIInt  size;
5170567c835SRichard Tran Mills   Mat_SeqKAIJ *seqkaij = (Mat_SeqKAIJ *)A->data;
5180567c835SRichard Tran Mills 
51949bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
5209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
5210567c835SRichard Tran Mills   if (size == 1) {
5229566063dSJacob Faibussowitsch     PetscCall(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));
5239566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(A->rmap, seqkaij->p));
5249566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(A->cmap, seqkaij->q));
5259566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(A->rmap));
5269566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(A->cmap));
5270567c835SRichard Tran Mills   } else {
5280567c835SRichard Tran Mills     Mat_MPIKAIJ *a;
5290567c835SRichard Tran Mills     Mat_MPIAIJ  *mpiaij;
5300567c835SRichard Tran Mills     IS           from, to;
5310567c835SRichard Tran Mills     Vec          gvec;
5320567c835SRichard Tran Mills 
5330567c835SRichard Tran Mills     a      = (Mat_MPIKAIJ *)A->data;
534d3f912faSRichard Tran Mills     mpiaij = (Mat_MPIAIJ *)a->A->data;
5359566063dSJacob Faibussowitsch     PetscCall(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));
5369566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(A->rmap, seqkaij->p));
5379566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(A->cmap, seqkaij->q));
5389566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(A->rmap));
5399566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(A->cmap));
5400567c835SRichard Tran Mills 
5419566063dSJacob Faibussowitsch     PetscCall(MatKAIJ_build_AIJ_OAIJ(A));
5420567c835SRichard Tran Mills 
5439566063dSJacob Faibussowitsch     PetscCall(VecGetSize(mpiaij->lvec, &n));
5449566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &a->w));
5459566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(a->w, n * a->q, n * a->q));
5469566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(a->w, a->q));
5479566063dSJacob Faibussowitsch     PetscCall(VecSetType(a->w, VECSEQ));
5480567c835SRichard Tran Mills 
5490567c835SRichard Tran Mills     /* create two temporary Index sets for build scatter gather */
5509566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)a->A), a->q, n, mpiaij->garray, PETSC_COPY_VALUES, &from));
5519566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, n * a->q, 0, 1, &to));
5520567c835SRichard Tran Mills 
5530567c835SRichard Tran Mills     /* create temporary global vector to generate scatter context */
5549566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)a->A), a->q, a->q * a->A->cmap->n, a->q * a->A->cmap->N, NULL, &gvec));
5550567c835SRichard Tran Mills 
5560567c835SRichard Tran Mills     /* generate the scatter context */
5579566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(gvec, from, a->w, to, &a->ctx));
5580567c835SRichard Tran Mills 
5599566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&from));
5609566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&to));
5619566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&gvec));
5620567c835SRichard Tran Mills   }
5630567c835SRichard Tran Mills 
5640567c835SRichard Tran Mills   A->assembled = PETSC_TRUE;
5653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56649bd79ccSDebojyoti Ghosh }
56749bd79ccSDebojyoti Ghosh 
568ff6a9541SJacob Faibussowitsch static PetscErrorCode MatView_KAIJ(Mat A, PetscViewer viewer)
569d71ae5a4SJacob Faibussowitsch {
570e6985dafSRichard Tran Mills   PetscViewerFormat format;
571e6985dafSRichard Tran Mills   Mat_SeqKAIJ      *a = (Mat_SeqKAIJ *)A->data;
57249bd79ccSDebojyoti Ghosh   Mat               B;
573e6985dafSRichard Tran Mills   PetscInt          i;
574e6985dafSRichard Tran Mills   PetscBool         ismpikaij;
57549bd79ccSDebojyoti Ghosh 
57649bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
5779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIKAIJ, &ismpikaij));
5789566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
579e6985dafSRichard Tran Mills   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_IMPL) {
5809566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "S and T have %" PetscInt_FMT " rows and %" PetscInt_FMT " columns\n", a->p, a->q));
581e6985dafSRichard Tran Mills 
582e6985dafSRichard Tran Mills     /* Print appropriate details for S. */
583e6985dafSRichard Tran Mills     if (!a->S) {
5849566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "S is NULL\n"));
585e6985dafSRichard Tran Mills     } else if (format == PETSC_VIEWER_ASCII_IMPL) {
5869566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Entries of S are "));
587e6985dafSRichard Tran Mills       for (i = 0; i < (a->p * a->q); i++) {
588e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
5899566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e %18.16e ", (double)PetscRealPart(a->S[i]), (double)PetscImaginaryPart(a->S[i])));
590e6985dafSRichard Tran Mills #else
5919566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e ", (double)PetscRealPart(a->S[i])));
592e6985dafSRichard Tran Mills #endif
593e6985dafSRichard Tran Mills       }
5949566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
59549bd79ccSDebojyoti Ghosh     }
59649bd79ccSDebojyoti Ghosh 
597e6985dafSRichard Tran Mills     /* Print appropriate details for T. */
598e6985dafSRichard Tran Mills     if (a->isTI) {
5999566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "T is the identity matrix\n"));
600e6985dafSRichard Tran Mills     } else if (!a->T) {
6019566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "T is NULL\n"));
602e6985dafSRichard Tran Mills     } else if (format == PETSC_VIEWER_ASCII_IMPL) {
6039566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Entries of T are "));
604e6985dafSRichard Tran Mills       for (i = 0; i < (a->p * a->q); i++) {
605e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
6069566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e %18.16e ", (double)PetscRealPart(a->T[i]), (double)PetscImaginaryPart(a->T[i])));
607e6985dafSRichard Tran Mills #else
6089566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e ", (double)PetscRealPart(a->T[i])));
609e6985dafSRichard Tran Mills #endif
610e6985dafSRichard Tran Mills       }
6119566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
612e6985dafSRichard Tran Mills     }
61349bd79ccSDebojyoti Ghosh 
614e6985dafSRichard Tran Mills     /* Now print details for the AIJ matrix, using the AIJ viewer. */
6159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Now viewing the associated AIJ matrix:\n"));
616e6985dafSRichard Tran Mills     if (ismpikaij) {
617e6985dafSRichard Tran Mills       Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
6189566063dSJacob Faibussowitsch       PetscCall(MatView(b->A, viewer));
619e6985dafSRichard Tran Mills     } else {
6209566063dSJacob Faibussowitsch       PetscCall(MatView(a->AIJ, viewer));
621e6985dafSRichard Tran Mills     }
622e6985dafSRichard Tran Mills 
623e6985dafSRichard Tran Mills   } else {
624e6985dafSRichard Tran Mills     /* For all other matrix viewer output formats, simply convert to an AIJ matrix and call MatView() on that. */
6259566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
6269566063dSJacob Faibussowitsch     PetscCall(MatView(B, viewer));
6279566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
628e6985dafSRichard Tran Mills   }
6293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63049bd79ccSDebojyoti Ghosh }
63149bd79ccSDebojyoti Ghosh 
632ff6a9541SJacob Faibussowitsch static PetscErrorCode MatDestroy_MPIKAIJ(Mat A)
633d71ae5a4SJacob Faibussowitsch {
63449bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
63549bd79ccSDebojyoti Ghosh 
63649bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
6379566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&b->AIJ));
6389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&b->OAIJ));
6399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&b->A));
6409566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&b->ctx));
6419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->w));
6429566063dSJacob Faibussowitsch   PetscCall(PetscFree(b->S));
6439566063dSJacob Faibussowitsch   PetscCall(PetscFree(b->T));
6449566063dSJacob Faibussowitsch   PetscCall(PetscFree(b->ibdiag));
6459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetDiagonalBlock_C", NULL));
6469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_mpikaij_mpiaij_C", NULL));
6479566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
6483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64949bd79ccSDebojyoti Ghosh }
65049bd79ccSDebojyoti Ghosh 
65149bd79ccSDebojyoti Ghosh /* zz = yy + Axx */
652ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMultAdd_SeqKAIJ(Mat A, Vec xx, Vec yy, Vec zz)
653d71ae5a4SJacob Faibussowitsch {
65449bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ       *b = (Mat_SeqKAIJ *)A->data;
65549bd79ccSDebojyoti Ghosh   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)b->AIJ->data;
65649bd79ccSDebojyoti Ghosh   const PetscScalar *s = b->S, *t = b->T;
65749bd79ccSDebojyoti Ghosh   const PetscScalar *x, *v, *bx;
65849bd79ccSDebojyoti Ghosh   PetscScalar       *y, *sums;
65949bd79ccSDebojyoti Ghosh   const PetscInt     m = b->AIJ->rmap->n, *idx, *ii;
66049bd79ccSDebojyoti Ghosh   PetscInt           n, i, jrow, j, l, p = b->p, q = b->q, k;
66149bd79ccSDebojyoti Ghosh 
66249bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
66349bd79ccSDebojyoti Ghosh   if (!yy) {
6649566063dSJacob Faibussowitsch     PetscCall(VecSet(zz, 0.0));
66549bd79ccSDebojyoti Ghosh   } else {
6669566063dSJacob Faibussowitsch     PetscCall(VecCopy(yy, zz));
66749bd79ccSDebojyoti Ghosh   }
6683ba16761SJacob Faibussowitsch   if ((!s) && (!t) && (!b->isTI)) PetscFunctionReturn(PETSC_SUCCESS);
66949bd79ccSDebojyoti Ghosh 
6709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx, &x));
6719566063dSJacob Faibussowitsch   PetscCall(VecGetArray(zz, &y));
67249bd79ccSDebojyoti Ghosh   idx = a->j;
67349bd79ccSDebojyoti Ghosh   v   = a->a;
67449bd79ccSDebojyoti Ghosh   ii  = a->i;
67549bd79ccSDebojyoti Ghosh 
67649bd79ccSDebojyoti Ghosh   if (b->isTI) {
67749bd79ccSDebojyoti Ghosh     for (i = 0; i < m; i++) {
67849bd79ccSDebojyoti Ghosh       jrow = ii[i];
67949bd79ccSDebojyoti Ghosh       n    = ii[i + 1] - jrow;
68049bd79ccSDebojyoti Ghosh       sums = y + p * i;
68149bd79ccSDebojyoti Ghosh       for (j = 0; j < n; j++) {
682ad540459SPierre Jolivet         for (k = 0; k < p; k++) sums[k] += v[jrow + j] * x[q * idx[jrow + j] + k];
68349bd79ccSDebojyoti Ghosh       }
68449bd79ccSDebojyoti Ghosh     }
6859566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(3.0 * (a->nz) * p));
68649bd79ccSDebojyoti Ghosh   } else if (t) {
68749bd79ccSDebojyoti Ghosh     for (i = 0; i < m; i++) {
68849bd79ccSDebojyoti Ghosh       jrow = ii[i];
68949bd79ccSDebojyoti Ghosh       n    = ii[i + 1] - jrow;
69049bd79ccSDebojyoti Ghosh       sums = y + p * i;
69149bd79ccSDebojyoti Ghosh       for (j = 0; j < n; j++) {
69249bd79ccSDebojyoti Ghosh         for (k = 0; k < p; k++) {
693ad540459SPierre Jolivet           for (l = 0; l < q; l++) sums[k] += v[jrow + j] * t[k + l * p] * x[q * idx[jrow + j] + l];
69449bd79ccSDebojyoti Ghosh         }
69549bd79ccSDebojyoti Ghosh       }
69649bd79ccSDebojyoti Ghosh     }
697234d9204SRichard Tran Mills     /* The flop count below assumes that v[jrow+j] is hoisted out (which an optimizing compiler is likely to do),
698234d9204SRichard Tran Mills      * and also that T part is hoisted outside this loop (in exchange for temporary storage) as (A \otimes I) (I \otimes T),
699234d9204SRichard Tran Mills      * so that this multiply doesn't have to be redone for each matrix entry, but just once per column. The latter
700234d9204SRichard Tran Mills      * transformation is much less likely to be applied, but we nonetheless count the minimum flops required. */
7019566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops((2.0 * p * q - p) * m + 2.0 * p * a->nz));
70249bd79ccSDebojyoti Ghosh   }
70349bd79ccSDebojyoti Ghosh   if (s) {
70449bd79ccSDebojyoti Ghosh     for (i = 0; i < m; i++) {
70549bd79ccSDebojyoti Ghosh       sums = y + p * i;
70649bd79ccSDebojyoti Ghosh       bx   = x + q * i;
70749bd79ccSDebojyoti Ghosh       if (i < b->AIJ->cmap->n) {
70849bd79ccSDebojyoti Ghosh         for (j = 0; j < q; j++) {
709ad540459SPierre Jolivet           for (k = 0; k < p; k++) sums[k] += s[k + j * p] * bx[j];
71049bd79ccSDebojyoti Ghosh         }
71149bd79ccSDebojyoti Ghosh       }
71249bd79ccSDebojyoti Ghosh     }
7139566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(2.0 * m * p * q));
71449bd79ccSDebojyoti Ghosh   }
71549bd79ccSDebojyoti Ghosh 
7169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx, &x));
7179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(zz, &y));
7183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71949bd79ccSDebojyoti Ghosh }
72049bd79ccSDebojyoti Ghosh 
721ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMult_SeqKAIJ(Mat A, Vec xx, Vec yy)
722d71ae5a4SJacob Faibussowitsch {
72349bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
724f3fa974cSJacob Faibussowitsch   PetscCall(MatMultAdd_SeqKAIJ(A, xx, NULL, yy));
7253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72649bd79ccSDebojyoti Ghosh }
72749bd79ccSDebojyoti Ghosh 
72849bd79ccSDebojyoti Ghosh #include <petsc/private/kernels/blockinvert.h>
72949bd79ccSDebojyoti Ghosh 
730ff6a9541SJacob Faibussowitsch static PetscErrorCode MatInvertBlockDiagonal_SeqKAIJ(Mat A, const PetscScalar **values)
731d71ae5a4SJacob Faibussowitsch {
73249bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ       *b = (Mat_SeqKAIJ *)A->data;
73349bd79ccSDebojyoti Ghosh   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)b->AIJ->data;
73449bd79ccSDebojyoti Ghosh   const PetscScalar *S = b->S;
73549bd79ccSDebojyoti Ghosh   const PetscScalar *T = b->T;
73649bd79ccSDebojyoti Ghosh   const PetscScalar *v = a->a;
73749bd79ccSDebojyoti Ghosh   const PetscInt     p = b->p, q = b->q, m = b->AIJ->rmap->n, *idx = a->j, *ii = a->i;
73849bd79ccSDebojyoti Ghosh   PetscInt           i, j, *v_pivots, dof, dof2;
73949bd79ccSDebojyoti Ghosh   PetscScalar       *diag, aval, *v_work;
74049bd79ccSDebojyoti Ghosh 
74149bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
74208401ef6SPierre Jolivet   PetscCheck(p == q, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MATKAIJ: Block size must be square to calculate inverse.");
743aed4548fSBarry Smith   PetscCheck(S || T || b->isTI, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MATKAIJ: Cannot invert a zero matrix.");
74449bd79ccSDebojyoti Ghosh 
74549bd79ccSDebojyoti Ghosh   dof  = p;
74649bd79ccSDebojyoti Ghosh   dof2 = dof * dof;
74749bd79ccSDebojyoti Ghosh 
74849bd79ccSDebojyoti Ghosh   if (b->ibdiagvalid) {
74949bd79ccSDebojyoti Ghosh     if (values) *values = b->ibdiag;
7503ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
75149bd79ccSDebojyoti Ghosh   }
752aa624791SPierre Jolivet   if (!b->ibdiag) PetscCall(PetscMalloc1(dof2 * m, &b->ibdiag));
75349bd79ccSDebojyoti Ghosh   if (values) *values = b->ibdiag;
75449bd79ccSDebojyoti Ghosh   diag = b->ibdiag;
75549bd79ccSDebojyoti Ghosh 
7569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dof, &v_work, dof, &v_pivots));
75749bd79ccSDebojyoti Ghosh   for (i = 0; i < m; i++) {
75849bd79ccSDebojyoti Ghosh     if (S) {
7599566063dSJacob Faibussowitsch       PetscCall(PetscMemcpy(diag, S, dof2 * sizeof(PetscScalar)));
76049bd79ccSDebojyoti Ghosh     } else {
7619566063dSJacob Faibussowitsch       PetscCall(PetscMemzero(diag, dof2 * sizeof(PetscScalar)));
76249bd79ccSDebojyoti Ghosh     }
76349bd79ccSDebojyoti Ghosh     if (b->isTI) {
76449bd79ccSDebojyoti Ghosh       aval = 0;
7659371c9d4SSatish Balay       for (j = ii[i]; j < ii[i + 1]; j++)
7669371c9d4SSatish Balay         if (idx[j] == i) aval = v[j];
76749bd79ccSDebojyoti Ghosh       for (j = 0; j < dof; j++) diag[j + dof * j] += aval;
76849bd79ccSDebojyoti Ghosh     } else if (T) {
76949bd79ccSDebojyoti Ghosh       aval = 0;
7709371c9d4SSatish Balay       for (j = ii[i]; j < ii[i + 1]; j++)
7719371c9d4SSatish Balay         if (idx[j] == i) aval = v[j];
77249bd79ccSDebojyoti Ghosh       for (j = 0; j < dof2; j++) diag[j] += aval * T[j];
77349bd79ccSDebojyoti Ghosh     }
7749566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A(dof, diag, v_pivots, v_work, PETSC_FALSE, NULL));
77549bd79ccSDebojyoti Ghosh     diag += dof2;
77649bd79ccSDebojyoti Ghosh   }
7779566063dSJacob Faibussowitsch   PetscCall(PetscFree2(v_work, v_pivots));
77849bd79ccSDebojyoti Ghosh 
77949bd79ccSDebojyoti Ghosh   b->ibdiagvalid = PETSC_TRUE;
7803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78149bd79ccSDebojyoti Ghosh }
78249bd79ccSDebojyoti Ghosh 
783d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_MPIKAIJ(Mat A, Mat *B)
784d71ae5a4SJacob Faibussowitsch {
78549bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ *kaij = (Mat_MPIKAIJ *)A->data;
78649bd79ccSDebojyoti Ghosh 
78749bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
78849bd79ccSDebojyoti Ghosh   *B = kaij->AIJ;
7893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79049bd79ccSDebojyoti Ghosh }
79149bd79ccSDebojyoti Ghosh 
792d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatConvert_KAIJ_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
793d71ae5a4SJacob Faibussowitsch {
794fac53fa1SPierre Jolivet   Mat_SeqKAIJ   *a = (Mat_SeqKAIJ *)A->data;
795fac53fa1SPierre Jolivet   Mat            AIJ, OAIJ, B;
796fac53fa1SPierre Jolivet   PetscInt      *d_nnz, *o_nnz = NULL, nz, i, j, m, d;
797fac53fa1SPierre Jolivet   const PetscInt p = a->p, q = a->q;
798fac53fa1SPierre Jolivet   PetscBool      ismpikaij, missing;
799fac53fa1SPierre Jolivet 
800fac53fa1SPierre Jolivet   PetscFunctionBegin;
801fac53fa1SPierre Jolivet   if (reuse != MAT_REUSE_MATRIX) {
8029566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIKAIJ, &ismpikaij));
803fac53fa1SPierre Jolivet     if (ismpikaij) {
804fac53fa1SPierre Jolivet       Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
805fac53fa1SPierre Jolivet       AIJ            = ((Mat_SeqKAIJ *)b->AIJ->data)->AIJ;
806fac53fa1SPierre Jolivet       OAIJ           = ((Mat_SeqKAIJ *)b->OAIJ->data)->AIJ;
807fac53fa1SPierre Jolivet     } else {
808fac53fa1SPierre Jolivet       AIJ  = a->AIJ;
809fac53fa1SPierre Jolivet       OAIJ = NULL;
810fac53fa1SPierre Jolivet     }
8119566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
8129566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
8139566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATAIJ));
8149566063dSJacob Faibussowitsch     PetscCall(MatGetSize(AIJ, &m, NULL));
8159566063dSJacob Faibussowitsch     PetscCall(MatMissingDiagonal(AIJ, &missing, &d)); /* assumption that all successive rows will have a missing diagonal */
816fac53fa1SPierre Jolivet     if (!missing || !a->S) d = m;
8179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m * p, &d_nnz));
818fac53fa1SPierre Jolivet     for (i = 0; i < m; ++i) {
8199566063dSJacob Faibussowitsch       PetscCall(MatGetRow_SeqAIJ(AIJ, i, &nz, NULL, NULL));
820fac53fa1SPierre Jolivet       for (j = 0; j < p; ++j) d_nnz[i * p + j] = nz * q + (i >= d) * q;
8219566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow_SeqAIJ(AIJ, i, &nz, NULL, NULL));
822fac53fa1SPierre Jolivet     }
823fac53fa1SPierre Jolivet     if (OAIJ) {
8249566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m * p, &o_nnz));
825fac53fa1SPierre Jolivet       for (i = 0; i < m; ++i) {
8269566063dSJacob Faibussowitsch         PetscCall(MatGetRow_SeqAIJ(OAIJ, i, &nz, NULL, NULL));
827fac53fa1SPierre Jolivet         for (j = 0; j < p; ++j) o_nnz[i * p + j] = nz * q;
8289566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow_SeqAIJ(OAIJ, i, &nz, NULL, NULL));
829fac53fa1SPierre Jolivet       }
8309566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJSetPreallocation(B, 0, d_nnz, 0, o_nnz));
831fac53fa1SPierre Jolivet     } else {
8329566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(B, 0, d_nnz));
833fac53fa1SPierre Jolivet     }
8349566063dSJacob Faibussowitsch     PetscCall(PetscFree(d_nnz));
8359566063dSJacob Faibussowitsch     PetscCall(PetscFree(o_nnz));
836fac53fa1SPierre Jolivet   } else B = *newmat;
8379566063dSJacob Faibussowitsch   PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &B));
838*ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
839*ac530a7eSPierre Jolivet   else *newmat = B;
8403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
841fac53fa1SPierre Jolivet }
842fac53fa1SPierre Jolivet 
843ff6a9541SJacob Faibussowitsch static PetscErrorCode MatSOR_SeqKAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
844d71ae5a4SJacob Faibussowitsch {
84549bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ       *kaij = (Mat_SeqKAIJ *)A->data;
84649bd79ccSDebojyoti Ghosh   Mat_SeqAIJ        *a    = (Mat_SeqAIJ *)kaij->AIJ->data;
84749bd79ccSDebojyoti Ghosh   const PetscScalar *aa = a->a, *T = kaij->T, *v;
84849bd79ccSDebojyoti Ghosh   const PetscInt     m = kaij->AIJ->rmap->n, *ai = a->i, *aj = a->j, p = kaij->p, q = kaij->q, *diag, *vi;
84949bd79ccSDebojyoti Ghosh   const PetscScalar *b, *xb, *idiag;
85049bd79ccSDebojyoti Ghosh   PetscScalar       *x, *work, *workt, *w, *y, *arr, *t, *arrt;
85149bd79ccSDebojyoti Ghosh   PetscInt           i, j, k, i2, bs, bs2, nz;
85249bd79ccSDebojyoti Ghosh 
85349bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
85449bd79ccSDebojyoti Ghosh   its = its * lits;
855aed4548fSBarry Smith   PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat");
85608401ef6SPierre Jolivet   PetscCheck(its > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive", its, lits);
85728b400f6SJacob Faibussowitsch   PetscCheck(!fshift, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for diagonal shift");
85808401ef6SPierre Jolivet   PetscCheck(!(flag & SOR_APPLY_UPPER) && !(flag & SOR_APPLY_LOWER), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for applying upper or lower triangular parts");
85908401ef6SPierre Jolivet   PetscCheck(p == q, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatSOR for KAIJ: No support for non-square dense blocks");
860f7d195e4SLawrence Mitchell   bs  = p;
861f7d195e4SLawrence Mitchell   bs2 = bs * bs;
86249bd79ccSDebojyoti Ghosh 
8633ba16761SJacob Faibussowitsch   if (!m) PetscFunctionReturn(PETSC_SUCCESS);
86449bd79ccSDebojyoti Ghosh 
8659566063dSJacob Faibussowitsch   if (!kaij->ibdiagvalid) PetscCall(MatInvertBlockDiagonal_SeqKAIJ(A, NULL));
86649bd79ccSDebojyoti Ghosh   idiag = kaij->ibdiag;
86749bd79ccSDebojyoti Ghosh   diag  = a->diag;
86849bd79ccSDebojyoti Ghosh 
86949bd79ccSDebojyoti Ghosh   if (!kaij->sor.setup) {
8709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc5(bs, &kaij->sor.w, bs, &kaij->sor.y, m * bs, &kaij->sor.work, m * bs, &kaij->sor.t, m * bs2, &kaij->sor.arr));
87149bd79ccSDebojyoti Ghosh     kaij->sor.setup = PETSC_TRUE;
87249bd79ccSDebojyoti Ghosh   }
87349bd79ccSDebojyoti Ghosh   y    = kaij->sor.y;
87449bd79ccSDebojyoti Ghosh   w    = kaij->sor.w;
87549bd79ccSDebojyoti Ghosh   work = kaij->sor.work;
87649bd79ccSDebojyoti Ghosh   t    = kaij->sor.t;
87749bd79ccSDebojyoti Ghosh   arr  = kaij->sor.arr;
87849bd79ccSDebojyoti Ghosh 
879d0609cedSBarry Smith   PetscCall(VecGetArray(xx, &x));
8809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
88149bd79ccSDebojyoti Ghosh 
88249bd79ccSDebojyoti Ghosh   if (flag & SOR_ZERO_INITIAL_GUESS) {
88349bd79ccSDebojyoti Ghosh     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
88449bd79ccSDebojyoti Ghosh       PetscKernel_w_gets_Ar_times_v(bs, bs, b, idiag, x); /* x[0:bs] <- D^{-1} b[0:bs] */
8859566063dSJacob Faibussowitsch       PetscCall(PetscMemcpy(t, b, bs * sizeof(PetscScalar)));
88649bd79ccSDebojyoti Ghosh       i2 = bs;
88749bd79ccSDebojyoti Ghosh       idiag += bs2;
88849bd79ccSDebojyoti Ghosh       for (i = 1; i < m; i++) {
88949bd79ccSDebojyoti Ghosh         v  = aa + ai[i];
89049bd79ccSDebojyoti Ghosh         vi = aj + ai[i];
89149bd79ccSDebojyoti Ghosh         nz = diag[i] - ai[i];
89249bd79ccSDebojyoti Ghosh 
89349bd79ccSDebojyoti Ghosh         if (T) { /* b - T (Arow * x) */
8949566063dSJacob Faibussowitsch           PetscCall(PetscMemzero(w, bs * sizeof(PetscScalar)));
89549bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
89649bd79ccSDebojyoti Ghosh             for (k = 0; k < bs; k++) w[k] -= v[j] * x[vi[j] * bs + k];
89749bd79ccSDebojyoti Ghosh           }
89849bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs, w, T, &t[i2]);
89949bd79ccSDebojyoti Ghosh           for (k = 0; k < bs; k++) t[i2 + k] += b[i2 + k];
90049bd79ccSDebojyoti Ghosh         } else if (kaij->isTI) {
9019566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(t + i2, b + i2, bs * sizeof(PetscScalar)));
90249bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
90349bd79ccSDebojyoti Ghosh             for (k = 0; k < bs; k++) t[i2 + k] -= v[j] * x[vi[j] * bs + k];
90449bd79ccSDebojyoti Ghosh           }
90549bd79ccSDebojyoti Ghosh         } else {
9069566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(t + i2, b + i2, bs * sizeof(PetscScalar)));
90749bd79ccSDebojyoti Ghosh         }
90849bd79ccSDebojyoti Ghosh 
90949bd79ccSDebojyoti Ghosh         PetscKernel_w_gets_Ar_times_v(bs, bs, t + i2, idiag, y);
91049bd79ccSDebojyoti Ghosh         for (j = 0; j < bs; j++) x[i2 + j] = omega * y[j];
91149bd79ccSDebojyoti Ghosh 
91249bd79ccSDebojyoti Ghosh         idiag += bs2;
91349bd79ccSDebojyoti Ghosh         i2 += bs;
91449bd79ccSDebojyoti Ghosh       }
91549bd79ccSDebojyoti Ghosh       /* for logging purposes assume number of nonzero in lower half is 1/2 of total */
9169566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(1.0 * bs2 * a->nz));
91749bd79ccSDebojyoti Ghosh       xb = t;
91849bd79ccSDebojyoti Ghosh     } else xb = b;
91949bd79ccSDebojyoti Ghosh     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
92049bd79ccSDebojyoti Ghosh       idiag = kaij->ibdiag + bs2 * (m - 1);
92149bd79ccSDebojyoti Ghosh       i2    = bs * (m - 1);
9229566063dSJacob Faibussowitsch       PetscCall(PetscMemcpy(w, xb + i2, bs * sizeof(PetscScalar)));
92349bd79ccSDebojyoti Ghosh       PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, x + i2);
92449bd79ccSDebojyoti Ghosh       i2 -= bs;
92549bd79ccSDebojyoti Ghosh       idiag -= bs2;
92649bd79ccSDebojyoti Ghosh       for (i = m - 2; i >= 0; i--) {
92749bd79ccSDebojyoti Ghosh         v  = aa + diag[i] + 1;
92849bd79ccSDebojyoti Ghosh         vi = aj + diag[i] + 1;
92949bd79ccSDebojyoti Ghosh         nz = ai[i + 1] - diag[i] - 1;
93049bd79ccSDebojyoti Ghosh 
93149bd79ccSDebojyoti Ghosh         if (T) { /* FIXME: This branch untested */
9329566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(w, xb + i2, bs * sizeof(PetscScalar)));
93349bd79ccSDebojyoti Ghosh           /* copy all rows of x that are needed into contiguous space */
93449bd79ccSDebojyoti Ghosh           workt = work;
93549bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
9369566063dSJacob Faibussowitsch             PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar)));
93749bd79ccSDebojyoti Ghosh             workt += bs;
93849bd79ccSDebojyoti Ghosh           }
93949bd79ccSDebojyoti Ghosh           arrt = arr;
94049bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
9419566063dSJacob Faibussowitsch             PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar)));
94249bd79ccSDebojyoti Ghosh             for (k = 0; k < bs2; k++) arrt[k] *= v[j];
94349bd79ccSDebojyoti Ghosh             arrt += bs2;
94449bd79ccSDebojyoti Ghosh           }
94549bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
94649bd79ccSDebojyoti Ghosh         } else if (kaij->isTI) {
9479566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(w, t + i2, bs * sizeof(PetscScalar)));
94849bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
94949bd79ccSDebojyoti Ghosh             for (k = 0; k < bs; k++) w[k] -= v[j] * x[vi[j] * bs + k];
95049bd79ccSDebojyoti Ghosh           }
95149bd79ccSDebojyoti Ghosh         }
95249bd79ccSDebojyoti Ghosh 
95349bd79ccSDebojyoti Ghosh         PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); /* RHS incorrect for omega != 1.0 */
95449bd79ccSDebojyoti Ghosh         for (j = 0; j < bs; j++) x[i2 + j] = (1.0 - omega) * x[i2 + j] + omega * y[j];
95549bd79ccSDebojyoti Ghosh 
95649bd79ccSDebojyoti Ghosh         idiag -= bs2;
95749bd79ccSDebojyoti Ghosh         i2 -= bs;
95849bd79ccSDebojyoti Ghosh       }
9599566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz)));
96049bd79ccSDebojyoti Ghosh     }
96149bd79ccSDebojyoti Ghosh     its--;
96249bd79ccSDebojyoti Ghosh   }
96349bd79ccSDebojyoti Ghosh   while (its--) { /* FIXME: This branch not updated */
96449bd79ccSDebojyoti Ghosh     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
96549bd79ccSDebojyoti Ghosh       i2    = 0;
96649bd79ccSDebojyoti Ghosh       idiag = kaij->ibdiag;
96749bd79ccSDebojyoti Ghosh       for (i = 0; i < m; i++) {
9689566063dSJacob Faibussowitsch         PetscCall(PetscMemcpy(w, b + i2, bs * sizeof(PetscScalar)));
96949bd79ccSDebojyoti Ghosh 
97049bd79ccSDebojyoti Ghosh         v     = aa + ai[i];
97149bd79ccSDebojyoti Ghosh         vi    = aj + ai[i];
97249bd79ccSDebojyoti Ghosh         nz    = diag[i] - ai[i];
97349bd79ccSDebojyoti Ghosh         workt = work;
97449bd79ccSDebojyoti Ghosh         for (j = 0; j < nz; j++) {
9759566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar)));
97649bd79ccSDebojyoti Ghosh           workt += bs;
97749bd79ccSDebojyoti Ghosh         }
97849bd79ccSDebojyoti Ghosh         arrt = arr;
97949bd79ccSDebojyoti Ghosh         if (T) {
98049bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
9819566063dSJacob Faibussowitsch             PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar)));
98249bd79ccSDebojyoti Ghosh             for (k = 0; k < bs2; k++) arrt[k] *= v[j];
98349bd79ccSDebojyoti Ghosh             arrt += bs2;
98449bd79ccSDebojyoti Ghosh           }
98549bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
98649bd79ccSDebojyoti Ghosh         } else if (kaij->isTI) {
98749bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
9889566063dSJacob Faibussowitsch             PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar)));
98949bd79ccSDebojyoti Ghosh             for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
99049bd79ccSDebojyoti Ghosh             arrt += bs2;
99149bd79ccSDebojyoti Ghosh           }
99249bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
99349bd79ccSDebojyoti Ghosh         }
9949566063dSJacob Faibussowitsch         PetscCall(PetscMemcpy(t + i2, w, bs * sizeof(PetscScalar)));
99549bd79ccSDebojyoti Ghosh 
99649bd79ccSDebojyoti Ghosh         v     = aa + diag[i] + 1;
99749bd79ccSDebojyoti Ghosh         vi    = aj + diag[i] + 1;
99849bd79ccSDebojyoti Ghosh         nz    = ai[i + 1] - diag[i] - 1;
99949bd79ccSDebojyoti Ghosh         workt = work;
100049bd79ccSDebojyoti Ghosh         for (j = 0; j < nz; j++) {
10019566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar)));
100249bd79ccSDebojyoti Ghosh           workt += bs;
100349bd79ccSDebojyoti Ghosh         }
100449bd79ccSDebojyoti Ghosh         arrt = arr;
100549bd79ccSDebojyoti Ghosh         if (T) {
100649bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
10079566063dSJacob Faibussowitsch             PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar)));
100849bd79ccSDebojyoti Ghosh             for (k = 0; k < bs2; k++) arrt[k] *= v[j];
100949bd79ccSDebojyoti Ghosh             arrt += bs2;
101049bd79ccSDebojyoti Ghosh           }
101149bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
101249bd79ccSDebojyoti Ghosh         } else if (kaij->isTI) {
101349bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
10149566063dSJacob Faibussowitsch             PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar)));
101549bd79ccSDebojyoti Ghosh             for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
101649bd79ccSDebojyoti Ghosh             arrt += bs2;
101749bd79ccSDebojyoti Ghosh           }
101849bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
101949bd79ccSDebojyoti Ghosh         }
102049bd79ccSDebojyoti Ghosh 
102149bd79ccSDebojyoti Ghosh         PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y);
102249bd79ccSDebojyoti Ghosh         for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j);
102349bd79ccSDebojyoti Ghosh 
102449bd79ccSDebojyoti Ghosh         idiag += bs2;
102549bd79ccSDebojyoti Ghosh         i2 += bs;
102649bd79ccSDebojyoti Ghosh       }
102749bd79ccSDebojyoti Ghosh       xb = t;
10289371c9d4SSatish Balay     } else xb = b;
102949bd79ccSDebojyoti Ghosh     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
103049bd79ccSDebojyoti Ghosh       idiag = kaij->ibdiag + bs2 * (m - 1);
103149bd79ccSDebojyoti Ghosh       i2    = bs * (m - 1);
103249bd79ccSDebojyoti Ghosh       if (xb == b) {
103349bd79ccSDebojyoti Ghosh         for (i = m - 1; i >= 0; i--) {
10349566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(w, b + i2, bs * sizeof(PetscScalar)));
103549bd79ccSDebojyoti Ghosh 
103649bd79ccSDebojyoti Ghosh           v     = aa + ai[i];
103749bd79ccSDebojyoti Ghosh           vi    = aj + ai[i];
103849bd79ccSDebojyoti Ghosh           nz    = diag[i] - ai[i];
103949bd79ccSDebojyoti Ghosh           workt = work;
104049bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
10419566063dSJacob Faibussowitsch             PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar)));
104249bd79ccSDebojyoti Ghosh             workt += bs;
104349bd79ccSDebojyoti Ghosh           }
104449bd79ccSDebojyoti Ghosh           arrt = arr;
104549bd79ccSDebojyoti Ghosh           if (T) {
104649bd79ccSDebojyoti Ghosh             for (j = 0; j < nz; j++) {
10479566063dSJacob Faibussowitsch               PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar)));
104849bd79ccSDebojyoti Ghosh               for (k = 0; k < bs2; k++) arrt[k] *= v[j];
104949bd79ccSDebojyoti Ghosh               arrt += bs2;
105049bd79ccSDebojyoti Ghosh             }
105149bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
105249bd79ccSDebojyoti Ghosh           } else if (kaij->isTI) {
105349bd79ccSDebojyoti Ghosh             for (j = 0; j < nz; j++) {
10549566063dSJacob Faibussowitsch               PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar)));
105549bd79ccSDebojyoti Ghosh               for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
105649bd79ccSDebojyoti Ghosh               arrt += bs2;
105749bd79ccSDebojyoti Ghosh             }
105849bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
105949bd79ccSDebojyoti Ghosh           }
106049bd79ccSDebojyoti Ghosh 
106149bd79ccSDebojyoti Ghosh           v     = aa + diag[i] + 1;
106249bd79ccSDebojyoti Ghosh           vi    = aj + diag[i] + 1;
106349bd79ccSDebojyoti Ghosh           nz    = ai[i + 1] - diag[i] - 1;
106449bd79ccSDebojyoti Ghosh           workt = work;
106549bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
10669566063dSJacob Faibussowitsch             PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar)));
106749bd79ccSDebojyoti Ghosh             workt += bs;
106849bd79ccSDebojyoti Ghosh           }
106949bd79ccSDebojyoti Ghosh           arrt = arr;
107049bd79ccSDebojyoti Ghosh           if (T) {
107149bd79ccSDebojyoti Ghosh             for (j = 0; j < nz; j++) {
10729566063dSJacob Faibussowitsch               PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar)));
107349bd79ccSDebojyoti Ghosh               for (k = 0; k < bs2; k++) arrt[k] *= v[j];
107449bd79ccSDebojyoti Ghosh               arrt += bs2;
107549bd79ccSDebojyoti Ghosh             }
107649bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
107749bd79ccSDebojyoti Ghosh           } else if (kaij->isTI) {
107849bd79ccSDebojyoti Ghosh             for (j = 0; j < nz; j++) {
10799566063dSJacob Faibussowitsch               PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar)));
108049bd79ccSDebojyoti Ghosh               for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
108149bd79ccSDebojyoti Ghosh               arrt += bs2;
108249bd79ccSDebojyoti Ghosh             }
108349bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
108449bd79ccSDebojyoti Ghosh           }
108549bd79ccSDebojyoti Ghosh 
108649bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y);
108749bd79ccSDebojyoti Ghosh           for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j);
108849bd79ccSDebojyoti Ghosh         }
108949bd79ccSDebojyoti Ghosh       } else {
109049bd79ccSDebojyoti Ghosh         for (i = m - 1; i >= 0; i--) {
10919566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(w, xb + i2, bs * sizeof(PetscScalar)));
109249bd79ccSDebojyoti Ghosh           v     = aa + diag[i] + 1;
109349bd79ccSDebojyoti Ghosh           vi    = aj + diag[i] + 1;
109449bd79ccSDebojyoti Ghosh           nz    = ai[i + 1] - diag[i] - 1;
109549bd79ccSDebojyoti Ghosh           workt = work;
109649bd79ccSDebojyoti Ghosh           for (j = 0; j < nz; j++) {
10979566063dSJacob Faibussowitsch             PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar)));
109849bd79ccSDebojyoti Ghosh             workt += bs;
109949bd79ccSDebojyoti Ghosh           }
110049bd79ccSDebojyoti Ghosh           arrt = arr;
110149bd79ccSDebojyoti Ghosh           if (T) {
110249bd79ccSDebojyoti Ghosh             for (j = 0; j < nz; j++) {
11039566063dSJacob Faibussowitsch               PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar)));
110449bd79ccSDebojyoti Ghosh               for (k = 0; k < bs2; k++) arrt[k] *= v[j];
110549bd79ccSDebojyoti Ghosh               arrt += bs2;
110649bd79ccSDebojyoti Ghosh             }
110749bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
110849bd79ccSDebojyoti Ghosh           } else if (kaij->isTI) {
110949bd79ccSDebojyoti Ghosh             for (j = 0; j < nz; j++) {
11109566063dSJacob Faibussowitsch               PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar)));
111149bd79ccSDebojyoti Ghosh               for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
111249bd79ccSDebojyoti Ghosh               arrt += bs2;
111349bd79ccSDebojyoti Ghosh             }
111449bd79ccSDebojyoti Ghosh             PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
111549bd79ccSDebojyoti Ghosh           }
111649bd79ccSDebojyoti Ghosh           PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y);
111749bd79ccSDebojyoti Ghosh           for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j);
111849bd79ccSDebojyoti Ghosh         }
111949bd79ccSDebojyoti Ghosh       }
11209566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz)));
112149bd79ccSDebojyoti Ghosh     }
112249bd79ccSDebojyoti Ghosh   }
112349bd79ccSDebojyoti Ghosh 
1124d0609cedSBarry Smith   PetscCall(VecRestoreArray(xx, &x));
11259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
11263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112749bd79ccSDebojyoti Ghosh }
112849bd79ccSDebojyoti Ghosh 
112949bd79ccSDebojyoti Ghosh /*===================================================================================*/
113049bd79ccSDebojyoti Ghosh 
1131ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMultAdd_MPIKAIJ(Mat A, Vec xx, Vec yy, Vec zz)
1132d71ae5a4SJacob Faibussowitsch {
113349bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
113449bd79ccSDebojyoti Ghosh 
113549bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
113649bd79ccSDebojyoti Ghosh   if (!yy) {
11379566063dSJacob Faibussowitsch     PetscCall(VecSet(zz, 0.0));
113849bd79ccSDebojyoti Ghosh   } else {
11399566063dSJacob Faibussowitsch     PetscCall(VecCopy(yy, zz));
114049bd79ccSDebojyoti Ghosh   }
11419566063dSJacob Faibussowitsch   PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ and b->OAIJ are up to date. */
114249bd79ccSDebojyoti Ghosh   /* start the scatter */
11439566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD));
11449566063dSJacob Faibussowitsch   PetscCall((*b->AIJ->ops->multadd)(b->AIJ, xx, zz, zz));
11459566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD));
11469566063dSJacob Faibussowitsch   PetscCall((*b->OAIJ->ops->multadd)(b->OAIJ, b->w, zz, zz));
11473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114849bd79ccSDebojyoti Ghosh }
114949bd79ccSDebojyoti Ghosh 
1150ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMult_MPIKAIJ(Mat A, Vec xx, Vec yy)
1151d71ae5a4SJacob Faibussowitsch {
115249bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
1153f3fa974cSJacob Faibussowitsch   PetscCall(MatMultAdd_MPIKAIJ(A, xx, NULL, yy));
11543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115549bd79ccSDebojyoti Ghosh }
115649bd79ccSDebojyoti Ghosh 
1157ff6a9541SJacob Faibussowitsch static PetscErrorCode MatInvertBlockDiagonal_MPIKAIJ(Mat A, const PetscScalar **values)
1158d71ae5a4SJacob Faibussowitsch {
115949bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
116049bd79ccSDebojyoti Ghosh 
116149bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
11629566063dSJacob Faibussowitsch   PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ is up to date. */
11639566063dSJacob Faibussowitsch   PetscCall((*b->AIJ->ops->invertblockdiagonal)(b->AIJ, values));
11643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
116549bd79ccSDebojyoti Ghosh }
116649bd79ccSDebojyoti Ghosh 
1167ff6a9541SJacob Faibussowitsch static PetscErrorCode MatGetRow_SeqKAIJ(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **values)
1168d71ae5a4SJacob Faibussowitsch {
116949bd79ccSDebojyoti Ghosh   Mat_SeqKAIJ *b    = (Mat_SeqKAIJ *)A->data;
11703ba16761SJacob Faibussowitsch   PetscBool    diag = PETSC_FALSE;
117149bd79ccSDebojyoti Ghosh   PetscInt     nzaij, nz, *colsaij, *idx, i, j, p = b->p, q = b->q, r = row / p, s = row % p, c;
117249bd79ccSDebojyoti Ghosh   PetscScalar *vaij, *v, *S = b->S, *T = b->T;
117349bd79ccSDebojyoti Ghosh 
117449bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
117528b400f6SJacob Faibussowitsch   PetscCheck(!b->getrowactive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Already active");
117649bd79ccSDebojyoti Ghosh   b->getrowactive = PETSC_TRUE;
1177aed4548fSBarry Smith   PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " out of range", row);
117849bd79ccSDebojyoti Ghosh 
117949bd79ccSDebojyoti Ghosh   if ((!S) && (!T) && (!b->isTI)) {
118049bd79ccSDebojyoti Ghosh     if (ncols) *ncols = 0;
118149bd79ccSDebojyoti Ghosh     if (cols) *cols = NULL;
118249bd79ccSDebojyoti Ghosh     if (values) *values = NULL;
11833ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
118449bd79ccSDebojyoti Ghosh   }
118549bd79ccSDebojyoti Ghosh 
118649bd79ccSDebojyoti Ghosh   if (T || b->isTI) {
11879566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqAIJ(b->AIJ, r, &nzaij, &colsaij, &vaij));
118849bd79ccSDebojyoti Ghosh     c = nzaij;
118949bd79ccSDebojyoti Ghosh     for (i = 0; i < nzaij; i++) {
119049bd79ccSDebojyoti Ghosh       /* check if this row contains a diagonal entry */
119149bd79ccSDebojyoti Ghosh       if (colsaij[i] == r) {
119249bd79ccSDebojyoti Ghosh         diag = PETSC_TRUE;
119349bd79ccSDebojyoti Ghosh         c    = i;
119449bd79ccSDebojyoti Ghosh       }
119549bd79ccSDebojyoti Ghosh     }
119649bd79ccSDebojyoti Ghosh   } else nzaij = c = 0;
119749bd79ccSDebojyoti Ghosh 
119849bd79ccSDebojyoti Ghosh   /* calculate size of row */
119949bd79ccSDebojyoti Ghosh   nz = 0;
120049bd79ccSDebojyoti Ghosh   if (S) nz += q;
120149bd79ccSDebojyoti Ghosh   if (T || b->isTI) nz += (diag && S ? (nzaij - 1) * q : nzaij * q);
120249bd79ccSDebojyoti Ghosh 
120349bd79ccSDebojyoti Ghosh   if (cols || values) {
12049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &idx, nz, &v));
120538822f9dSRichard Tran Mills     for (i = 0; i < q; i++) {
120638822f9dSRichard Tran Mills       /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */
120738822f9dSRichard Tran Mills       v[i] = 0.0;
120838822f9dSRichard Tran Mills     }
120949bd79ccSDebojyoti Ghosh     if (b->isTI) {
121049bd79ccSDebojyoti Ghosh       for (i = 0; i < nzaij; i++) {
121149bd79ccSDebojyoti Ghosh         for (j = 0; j < q; j++) {
121249bd79ccSDebojyoti Ghosh           idx[i * q + j] = colsaij[i] * q + j;
121349bd79ccSDebojyoti Ghosh           v[i * q + j]   = (j == s ? vaij[i] : 0);
121449bd79ccSDebojyoti Ghosh         }
121549bd79ccSDebojyoti Ghosh       }
121649bd79ccSDebojyoti Ghosh     } else if (T) {
121749bd79ccSDebojyoti Ghosh       for (i = 0; i < nzaij; i++) {
121849bd79ccSDebojyoti Ghosh         for (j = 0; j < q; j++) {
121949bd79ccSDebojyoti Ghosh           idx[i * q + j] = colsaij[i] * q + j;
122049bd79ccSDebojyoti Ghosh           v[i * q + j]   = vaij[i] * T[s + j * p];
122149bd79ccSDebojyoti Ghosh         }
122249bd79ccSDebojyoti Ghosh       }
122349bd79ccSDebojyoti Ghosh     }
122449bd79ccSDebojyoti Ghosh     if (S) {
122549bd79ccSDebojyoti Ghosh       for (j = 0; j < q; j++) {
122649bd79ccSDebojyoti Ghosh         idx[c * q + j] = r * q + j;
122749bd79ccSDebojyoti Ghosh         v[c * q + j] += S[s + j * p];
122849bd79ccSDebojyoti Ghosh       }
122949bd79ccSDebojyoti Ghosh     }
123049bd79ccSDebojyoti Ghosh   }
123149bd79ccSDebojyoti Ghosh 
123249bd79ccSDebojyoti Ghosh   if (ncols) *ncols = nz;
123349bd79ccSDebojyoti Ghosh   if (cols) *cols = idx;
123449bd79ccSDebojyoti Ghosh   if (values) *values = v;
12353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
123649bd79ccSDebojyoti Ghosh }
123749bd79ccSDebojyoti Ghosh 
1238ff6a9541SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_SeqKAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
1239d71ae5a4SJacob Faibussowitsch {
124049bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
12419566063dSJacob Faibussowitsch   PetscCall(PetscFree2(*idx, *v));
124249bd79ccSDebojyoti Ghosh   ((Mat_SeqKAIJ *)A->data)->getrowactive = PETSC_FALSE;
12433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124449bd79ccSDebojyoti Ghosh }
124549bd79ccSDebojyoti Ghosh 
1246ff6a9541SJacob Faibussowitsch static PetscErrorCode MatGetRow_MPIKAIJ(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **values)
1247d71ae5a4SJacob Faibussowitsch {
124849bd79ccSDebojyoti Ghosh   Mat_MPIKAIJ   *b    = (Mat_MPIKAIJ *)A->data;
124949bd79ccSDebojyoti Ghosh   Mat            AIJ  = b->A;
1250fc64b2cfSRichard Tran Mills   PetscBool      diag = PETSC_FALSE;
1251761d359dSRichard Tran Mills   Mat            MatAIJ, MatOAIJ;
125249bd79ccSDebojyoti Ghosh   const PetscInt rstart = A->rmap->rstart, rend = A->rmap->rend, p = b->p, q = b->q, *garray;
1253389eba51SJed Brown   PetscInt       nz, *idx, ncolsaij = 0, ncolsoaij = 0, *colsaij, *colsoaij, r, s, c, i, j, lrow;
125449bd79ccSDebojyoti Ghosh   PetscScalar   *v, *vals, *ovals, *S = b->S, *T = b->T;
125549bd79ccSDebojyoti Ghosh 
125649bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
12579566063dSJacob Faibussowitsch   PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ and b->OAIJ are up to date. */
1258761d359dSRichard Tran Mills   MatAIJ  = ((Mat_SeqKAIJ *)b->AIJ->data)->AIJ;
1259761d359dSRichard Tran Mills   MatOAIJ = ((Mat_SeqKAIJ *)b->OAIJ->data)->AIJ;
126028b400f6SJacob Faibussowitsch   PetscCheck(!b->getrowactive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Already active");
126149bd79ccSDebojyoti Ghosh   b->getrowactive = PETSC_TRUE;
1262aed4548fSBarry Smith   PetscCheck(row >= rstart && row < rend, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only local rows");
126349bd79ccSDebojyoti Ghosh   lrow = row - rstart;
126449bd79ccSDebojyoti Ghosh 
126549bd79ccSDebojyoti Ghosh   if ((!S) && (!T) && (!b->isTI)) {
126649bd79ccSDebojyoti Ghosh     if (ncols) *ncols = 0;
126749bd79ccSDebojyoti Ghosh     if (cols) *cols = NULL;
126849bd79ccSDebojyoti Ghosh     if (values) *values = NULL;
12693ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
127049bd79ccSDebojyoti Ghosh   }
127149bd79ccSDebojyoti Ghosh 
127249bd79ccSDebojyoti Ghosh   r = lrow / p;
127349bd79ccSDebojyoti Ghosh   s = lrow % p;
127449bd79ccSDebojyoti Ghosh 
127549bd79ccSDebojyoti Ghosh   if (T || b->isTI) {
12769566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJGetSeqAIJ(AIJ, NULL, NULL, &garray));
12779566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqAIJ(MatAIJ, lrow / p, &ncolsaij, &colsaij, &vals));
12789566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqAIJ(MatOAIJ, lrow / p, &ncolsoaij, &colsoaij, &ovals));
127949bd79ccSDebojyoti Ghosh 
128049bd79ccSDebojyoti Ghosh     c = ncolsaij + ncolsoaij;
128149bd79ccSDebojyoti Ghosh     for (i = 0; i < ncolsaij; i++) {
128249bd79ccSDebojyoti Ghosh       /* check if this row contains a diagonal entry */
128349bd79ccSDebojyoti Ghosh       if (colsaij[i] == r) {
128449bd79ccSDebojyoti Ghosh         diag = PETSC_TRUE;
128549bd79ccSDebojyoti Ghosh         c    = i;
128649bd79ccSDebojyoti Ghosh       }
128749bd79ccSDebojyoti Ghosh     }
128849bd79ccSDebojyoti Ghosh   } else c = 0;
128949bd79ccSDebojyoti Ghosh 
129049bd79ccSDebojyoti Ghosh   /* calculate size of row */
129149bd79ccSDebojyoti Ghosh   nz = 0;
129249bd79ccSDebojyoti Ghosh   if (S) nz += q;
129349bd79ccSDebojyoti Ghosh   if (T || b->isTI) nz += (diag && S ? (ncolsaij + ncolsoaij - 1) * q : (ncolsaij + ncolsoaij) * q);
129449bd79ccSDebojyoti Ghosh 
129549bd79ccSDebojyoti Ghosh   if (cols || values) {
12969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nz, &idx, nz, &v));
1297a437a796SRichard Tran Mills     for (i = 0; i < q; i++) {
1298a437a796SRichard Tran Mills       /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */
1299a437a796SRichard Tran Mills       v[i] = 0.0;
1300a437a796SRichard Tran Mills     }
130149bd79ccSDebojyoti Ghosh     if (b->isTI) {
130249bd79ccSDebojyoti Ghosh       for (i = 0; i < ncolsaij; i++) {
130349bd79ccSDebojyoti Ghosh         for (j = 0; j < q; j++) {
130449bd79ccSDebojyoti Ghosh           idx[i * q + j] = (colsaij[i] + rstart / p) * q + j;
130549bd79ccSDebojyoti Ghosh           v[i * q + j]   = (j == s ? vals[i] : 0.0);
130649bd79ccSDebojyoti Ghosh         }
130749bd79ccSDebojyoti Ghosh       }
130849bd79ccSDebojyoti Ghosh       for (i = 0; i < ncolsoaij; i++) {
130949bd79ccSDebojyoti Ghosh         for (j = 0; j < q; j++) {
131049bd79ccSDebojyoti Ghosh           idx[(i + ncolsaij) * q + j] = garray[colsoaij[i]] * q + j;
131149bd79ccSDebojyoti Ghosh           v[(i + ncolsaij) * q + j]   = (j == s ? ovals[i] : 0.0);
131249bd79ccSDebojyoti Ghosh         }
131349bd79ccSDebojyoti Ghosh       }
131449bd79ccSDebojyoti Ghosh     } else if (T) {
131549bd79ccSDebojyoti Ghosh       for (i = 0; i < ncolsaij; i++) {
131649bd79ccSDebojyoti Ghosh         for (j = 0; j < q; j++) {
131749bd79ccSDebojyoti Ghosh           idx[i * q + j] = (colsaij[i] + rstart / p) * q + j;
131849bd79ccSDebojyoti Ghosh           v[i * q + j]   = vals[i] * T[s + j * p];
131949bd79ccSDebojyoti Ghosh         }
132049bd79ccSDebojyoti Ghosh       }
132149bd79ccSDebojyoti Ghosh       for (i = 0; i < ncolsoaij; i++) {
132249bd79ccSDebojyoti Ghosh         for (j = 0; j < q; j++) {
132349bd79ccSDebojyoti Ghosh           idx[(i + ncolsaij) * q + j] = garray[colsoaij[i]] * q + j;
132449bd79ccSDebojyoti Ghosh           v[(i + ncolsaij) * q + j]   = ovals[i] * T[s + j * p];
132549bd79ccSDebojyoti Ghosh         }
132649bd79ccSDebojyoti Ghosh       }
132749bd79ccSDebojyoti Ghosh     }
132849bd79ccSDebojyoti Ghosh     if (S) {
132949bd79ccSDebojyoti Ghosh       for (j = 0; j < q; j++) {
133049bd79ccSDebojyoti Ghosh         idx[c * q + j] = (r + rstart / p) * q + j;
133149bd79ccSDebojyoti Ghosh         v[c * q + j] += S[s + j * p];
133249bd79ccSDebojyoti Ghosh       }
133349bd79ccSDebojyoti Ghosh     }
133449bd79ccSDebojyoti Ghosh   }
133549bd79ccSDebojyoti Ghosh 
133649bd79ccSDebojyoti Ghosh   if (ncols) *ncols = nz;
133749bd79ccSDebojyoti Ghosh   if (cols) *cols = idx;
133849bd79ccSDebojyoti Ghosh   if (values) *values = v;
13393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
134049bd79ccSDebojyoti Ghosh }
134149bd79ccSDebojyoti Ghosh 
1342ff6a9541SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_MPIKAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
1343d71ae5a4SJacob Faibussowitsch {
134449bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
13459566063dSJacob Faibussowitsch   PetscCall(PetscFree2(*idx, *v));
134649bd79ccSDebojyoti Ghosh   ((Mat_SeqKAIJ *)A->data)->getrowactive = PETSC_FALSE;
13473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
134849bd79ccSDebojyoti Ghosh }
134949bd79ccSDebojyoti Ghosh 
1350ff6a9541SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_KAIJ(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat)
1351d71ae5a4SJacob Faibussowitsch {
135249bd79ccSDebojyoti Ghosh   Mat A;
135349bd79ccSDebojyoti Ghosh 
135449bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
13559566063dSJacob Faibussowitsch   PetscCall(MatConvert(mat, MATAIJ, MAT_INITIAL_MATRIX, &A));
13569566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, isrow, iscol, cll, newmat));
13579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
13583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135949bd79ccSDebojyoti Ghosh }
136049bd79ccSDebojyoti Ghosh 
136149bd79ccSDebojyoti Ghosh /*@C
13622920cce0SJacob Faibussowitsch   MatCreateKAIJ - Creates a matrix of type `MATKAIJ`.
136349bd79ccSDebojyoti Ghosh 
136449bd79ccSDebojyoti Ghosh   Collective
136549bd79ccSDebojyoti Ghosh 
136649bd79ccSDebojyoti Ghosh   Input Parameters:
136711a5261eSBarry Smith + A - the `MATAIJ` matrix
13682ef1f0ffSBarry Smith . p - number of rows in `S` and `T`
13692ef1f0ffSBarry Smith . q - number of columns in `S` and `T`
13702ef1f0ffSBarry Smith . S - the `S` matrix (can be `NULL`), stored as a `PetscScalar` array (column-major)
13712ef1f0ffSBarry Smith - T - the `T` matrix (can be `NULL`), stored as a `PetscScalar` array (column-major)
137249bd79ccSDebojyoti Ghosh 
137349bd79ccSDebojyoti Ghosh   Output Parameter:
137411a5261eSBarry Smith . kaij - the new `MATKAIJ` matrix
137549bd79ccSDebojyoti Ghosh 
13762ef1f0ffSBarry Smith   Level: advanced
13772ef1f0ffSBarry Smith 
1378d3b90ce1SRichard Tran Mills   Notes:
13792920cce0SJacob Faibussowitsch   The created matrix is of the following form\:
13802920cce0SJacob Faibussowitsch .vb
13812920cce0SJacob Faibussowitsch     [I \otimes S + A \otimes T]
13822920cce0SJacob Faibussowitsch .ve
13832920cce0SJacob Faibussowitsch   where
13842920cce0SJacob Faibussowitsch .vb
13852920cce0SJacob Faibussowitsch   S is a dense (p \times q) matrix
13862920cce0SJacob Faibussowitsch   T is a dense (p \times q) matrix
13872920cce0SJacob Faibussowitsch   A is a `MATAIJ`  (n \times n) matrix
13882920cce0SJacob Faibussowitsch   I is the identity matrix
13892920cce0SJacob Faibussowitsch .ve
13902920cce0SJacob Faibussowitsch   The resulting matrix is (np \times nq)
13912920cce0SJacob Faibussowitsch 
13922920cce0SJacob Faibussowitsch   `S` and `T` are always stored independently on all processes as `PetscScalar` arrays in
13932920cce0SJacob Faibussowitsch   column-major format.
13942920cce0SJacob Faibussowitsch 
139511a5261eSBarry Smith   This function increases the reference count on the `MATAIJ` matrix, so the user is free to destroy the matrix if it is not needed.
139649bd79ccSDebojyoti Ghosh 
139711a5261eSBarry Smith   Changes to the entries of the `MATAIJ` matrix will immediately affect the `MATKAIJ` matrix.
139811a5261eSBarry Smith 
1399fe59aa6dSJacob Faibussowitsch   Developer Notes:
140011a5261eSBarry Smith   In the `MATMPIKAIJ` case, the internal 'AIJ' and 'OAIJ' sequential KAIJ matrices are kept up to date by tracking the object state
14012ef1f0ffSBarry Smith   of the AIJ matrix 'A' that describes the blockwise action of the `MATMPIKAIJ` matrix and, if the object state has changed, lazily
14022ef1f0ffSBarry Smith   rebuilding 'AIJ' and 'OAIJ' just before executing operations with the `MATMPIKAIJ` matrix. If new types of operations are added,
1403761d359dSRichard Tran Mills   routines implementing those must also ensure these are rebuilt when needed (by calling the internal MatKAIJ_build_AIJ_OAIJ() routine).
1404761d359dSRichard Tran Mills 
14051cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatKAIJSetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()`, `MatKAIJGetAIJ()`, `MatKAIJGetS()`, `MatKAIJGetT()`, `MATKAIJ`
140649bd79ccSDebojyoti Ghosh @*/
1407d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateKAIJ(Mat A, PetscInt p, PetscInt q, const PetscScalar S[], const PetscScalar T[], Mat *kaij)
1408d71ae5a4SJacob Faibussowitsch {
140949bd79ccSDebojyoti Ghosh   PetscFunctionBegin;
14109566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), kaij));
14119566063dSJacob Faibussowitsch   PetscCall(MatSetType(*kaij, MATKAIJ));
14129566063dSJacob Faibussowitsch   PetscCall(MatKAIJSetAIJ(*kaij, A));
14139566063dSJacob Faibussowitsch   PetscCall(MatKAIJSetS(*kaij, p, q, S));
14149566063dSJacob Faibussowitsch   PetscCall(MatKAIJSetT(*kaij, p, q, T));
14159566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*kaij));
14163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14170567c835SRichard Tran Mills }
141849bd79ccSDebojyoti Ghosh 
14190567c835SRichard Tran Mills /*MC
14205881e567SRichard Tran Mills   MATKAIJ - MATKAIJ = "kaij" - A matrix type to be used to evaluate matrices of form
14215881e567SRichard Tran Mills     [I \otimes S + A \otimes T],
14220567c835SRichard Tran Mills   where
14232ef1f0ffSBarry Smith .vb
14245881e567SRichard Tran Mills     S is a dense (p \times q) matrix,
14255881e567SRichard Tran Mills     T is a dense (p \times q) matrix,
14265881e567SRichard Tran Mills     A is an AIJ  (n \times n) matrix,
14275881e567SRichard Tran Mills     and I is the identity matrix.
14282ef1f0ffSBarry Smith .ve
14295881e567SRichard Tran Mills   The resulting matrix is (np \times nq).
14300567c835SRichard Tran Mills 
143111a5261eSBarry Smith   S and T are always stored independently on all processes as `PetscScalar` arrays in column-major format.
14320567c835SRichard Tran Mills 
14332ef1f0ffSBarry Smith   Level: advanced
14342ef1f0ffSBarry Smith 
143511a5261eSBarry Smith   Note:
14365881e567SRichard 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,
14375881e567SRichard Tran Mills   where x and b are column vectors containing the row-major representations of X and B.
14385881e567SRichard Tran Mills 
14391cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatKAIJSetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()`, `MatKAIJGetAIJ()`, `MatKAIJGetS()`, `MatKAIJGetT()`, `MatCreateKAIJ()`
14400567c835SRichard Tran Mills M*/
14410567c835SRichard Tran Mills 
1442d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_KAIJ(Mat A)
1443d71ae5a4SJacob Faibussowitsch {
14440567c835SRichard Tran Mills   Mat_MPIKAIJ *b;
14450567c835SRichard Tran Mills   PetscMPIInt  size;
14460567c835SRichard Tran Mills 
14470567c835SRichard Tran Mills   PetscFunctionBegin;
14484dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
14490567c835SRichard Tran Mills   A->data = (void *)b;
14500567c835SRichard Tran Mills 
14519566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps)));
14520567c835SRichard Tran Mills 
1453f4259b30SLisandro Dalcin   b->w = NULL;
14549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
14550567c835SRichard Tran Mills   if (size == 1) {
14569566063dSJacob Faibussowitsch     PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATSEQKAIJ));
14570567c835SRichard Tran Mills     A->ops->destroy             = MatDestroy_SeqKAIJ;
1458bb6fb833SRichard Tran Mills     A->ops->mult                = MatMult_SeqKAIJ;
1459bb6fb833SRichard Tran Mills     A->ops->multadd             = MatMultAdd_SeqKAIJ;
1460bb6fb833SRichard Tran Mills     A->ops->invertblockdiagonal = MatInvertBlockDiagonal_SeqKAIJ;
14610567c835SRichard Tran Mills     A->ops->getrow              = MatGetRow_SeqKAIJ;
14620567c835SRichard Tran Mills     A->ops->restorerow          = MatRestoreRow_SeqKAIJ;
14630567c835SRichard Tran Mills     A->ops->sor                 = MatSOR_SeqKAIJ;
14649566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqkaij_seqaij_C", MatConvert_KAIJ_AIJ));
14650567c835SRichard Tran Mills   } else {
14669566063dSJacob Faibussowitsch     PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMPIKAIJ));
14670567c835SRichard Tran Mills     A->ops->destroy             = MatDestroy_MPIKAIJ;
1468bb6fb833SRichard Tran Mills     A->ops->mult                = MatMult_MPIKAIJ;
1469bb6fb833SRichard Tran Mills     A->ops->multadd             = MatMultAdd_MPIKAIJ;
1470bb6fb833SRichard Tran Mills     A->ops->invertblockdiagonal = MatInvertBlockDiagonal_MPIKAIJ;
14710567c835SRichard Tran Mills     A->ops->getrow              = MatGetRow_MPIKAIJ;
14720567c835SRichard Tran Mills     A->ops->restorerow          = MatRestoreRow_MPIKAIJ;
14739566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetDiagonalBlock_C", MatGetDiagonalBlock_MPIKAIJ));
14749566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_mpikaij_mpiaij_C", MatConvert_KAIJ_AIJ));
14750567c835SRichard Tran Mills   }
147606d61a37SPierre Jolivet   A->ops->setup           = MatSetUp_KAIJ;
147706d61a37SPierre Jolivet   A->ops->view            = MatView_KAIJ;
14780567c835SRichard Tran Mills   A->ops->createsubmatrix = MatCreateSubMatrix_KAIJ;
14793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
148049bd79ccSDebojyoti Ghosh }
1481