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 @*/
MatKAIJGetAIJ(Mat A,Mat * B)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 @*/
MatKAIJGetS(Mat A,PetscInt * m,PetscInt * n,PetscScalar * S[])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 @*/
MatKAIJGetSRead(Mat A,PetscInt * m,PetscInt * n,const PetscScalar * S[])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 @*/
MatKAIJRestoreS(Mat A,PetscScalar * S[])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 @*/
MatKAIJRestoreSRead(Mat A,const PetscScalar * S[])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 @*/
MatKAIJGetT(Mat A,PetscInt * m,PetscInt * n,PetscScalar * T[])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 @*/
MatKAIJGetTRead(Mat A,PetscInt * m,PetscInt * n,const PetscScalar * T[])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 @*/
MatKAIJRestoreT(Mat A,PetscScalar * T[])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 @*/
MatKAIJRestoreTRead(Mat A,const PetscScalar * T[])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 @*/
MatKAIJSetAIJ(Mat A,Mat B)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 @*/
MatKAIJSetS(Mat A,PetscInt p,PetscInt q,const PetscScalar S[])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));
354*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(a->S, S, p * q));
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- @*/
MatKAIJGetScaledIdentity(Mat A,PetscBool * identity)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 @*/
MatKAIJSetT(Mat A,PetscInt p,PetscInt q,const PetscScalar T[])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));
450*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(a->T, T, p * q));
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
MatDestroy_SeqKAIJ(Mat A)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
MatKAIJ_build_AIJ_OAIJ(Mat A)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
MatSetUp_KAIJ(Mat A)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
MatView_KAIJ(Mat A,PetscViewer viewer)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
MatDestroy_MPIKAIJ(Mat A)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 */
MatMultAdd_SeqKAIJ(Mat A,Vec xx,Vec yy,Vec zz)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
MatMult_SeqKAIJ(Mat A,Vec xx,Vec yy)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
MatInvertBlockDiagonal_SeqKAIJ(Mat A,const PetscScalar ** values)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) {
759*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(diag, S, dof2));
76049bd79ccSDebojyoti Ghosh } else {
761*418fb43bSPierre Jolivet PetscCall(PetscArrayzero(diag, dof2));
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
MatGetDiagonalBlock_MPIKAIJ(Mat A,Mat * B)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
MatConvert_KAIJ_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat * newmat)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;
798421480d9SBarry Smith PetscBool ismpikaij, diagDense;
799421480d9SBarry Smith const PetscInt *aijdiag;
800fac53fa1SPierre Jolivet
801fac53fa1SPierre Jolivet PetscFunctionBegin;
802fac53fa1SPierre Jolivet if (reuse != MAT_REUSE_MATRIX) {
8039566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIKAIJ, &ismpikaij));
804fac53fa1SPierre Jolivet if (ismpikaij) {
805fac53fa1SPierre Jolivet Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
806fac53fa1SPierre Jolivet AIJ = ((Mat_SeqKAIJ *)b->AIJ->data)->AIJ;
807fac53fa1SPierre Jolivet OAIJ = ((Mat_SeqKAIJ *)b->OAIJ->data)->AIJ;
808fac53fa1SPierre Jolivet } else {
809fac53fa1SPierre Jolivet AIJ = a->AIJ;
810fac53fa1SPierre Jolivet OAIJ = NULL;
811fac53fa1SPierre Jolivet }
8129566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
8139566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
8149566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ));
8159566063dSJacob Faibussowitsch PetscCall(MatGetSize(AIJ, &m, NULL));
816421480d9SBarry Smith
817421480d9SBarry Smith /*
818421480d9SBarry Smith Determine the first row of AIJ with missing diagonal and assume all successive rows also have a missing diagonal
819421480d9SBarry Smith */
820421480d9SBarry Smith PetscCall(MatGetDiagonalMarkers_SeqAIJ(AIJ, &aijdiag, &diagDense));
821421480d9SBarry Smith if (diagDense || !a->S) d = m;
822421480d9SBarry Smith else {
823421480d9SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)AIJ->data;
824421480d9SBarry Smith
825421480d9SBarry Smith for (d = 0; d < m; d++) {
826421480d9SBarry Smith if (aijdiag[d] == aij->i[d + 1]) break;
827421480d9SBarry Smith }
828421480d9SBarry Smith }
8299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * p, &d_nnz));
830fac53fa1SPierre Jolivet for (i = 0; i < m; ++i) {
8319566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(AIJ, i, &nz, NULL, NULL));
832fac53fa1SPierre Jolivet for (j = 0; j < p; ++j) d_nnz[i * p + j] = nz * q + (i >= d) * q;
8339566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(AIJ, i, &nz, NULL, NULL));
834fac53fa1SPierre Jolivet }
835fac53fa1SPierre Jolivet if (OAIJ) {
8369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * p, &o_nnz));
837fac53fa1SPierre Jolivet for (i = 0; i < m; ++i) {
8389566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(OAIJ, i, &nz, NULL, NULL));
839fac53fa1SPierre Jolivet for (j = 0; j < p; ++j) o_nnz[i * p + j] = nz * q;
8409566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(OAIJ, i, &nz, NULL, NULL));
841fac53fa1SPierre Jolivet }
8429566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(B, 0, d_nnz, 0, o_nnz));
843fac53fa1SPierre Jolivet } else {
8449566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(B, 0, d_nnz));
845fac53fa1SPierre Jolivet }
8469566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nnz));
8479566063dSJacob Faibussowitsch PetscCall(PetscFree(o_nnz));
848fac53fa1SPierre Jolivet } else B = *newmat;
8499566063dSJacob Faibussowitsch PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &B));
850ac530a7eSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
851ac530a7eSPierre Jolivet else *newmat = B;
8523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
853fac53fa1SPierre Jolivet }
854fac53fa1SPierre Jolivet
MatSOR_SeqKAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)855ff6a9541SJacob Faibussowitsch static PetscErrorCode MatSOR_SeqKAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
856d71ae5a4SJacob Faibussowitsch {
85749bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *kaij = (Mat_SeqKAIJ *)A->data;
85849bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ *)kaij->AIJ->data;
85949bd79ccSDebojyoti Ghosh const PetscScalar *aa = a->a, *T = kaij->T, *v;
86049bd79ccSDebojyoti Ghosh const PetscInt m = kaij->AIJ->rmap->n, *ai = a->i, *aj = a->j, p = kaij->p, q = kaij->q, *diag, *vi;
86149bd79ccSDebojyoti Ghosh const PetscScalar *b, *xb, *idiag;
86249bd79ccSDebojyoti Ghosh PetscScalar *x, *work, *workt, *w, *y, *arr, *t, *arrt;
86349bd79ccSDebojyoti Ghosh PetscInt i, j, k, i2, bs, bs2, nz;
86449bd79ccSDebojyoti Ghosh
86549bd79ccSDebojyoti Ghosh PetscFunctionBegin;
86649bd79ccSDebojyoti Ghosh its = its * lits;
867aed4548fSBarry Smith PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat");
86808401ef6SPierre 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);
86928b400f6SJacob Faibussowitsch PetscCheck(!fshift, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for diagonal shift");
87008401ef6SPierre 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");
87108401ef6SPierre Jolivet PetscCheck(p == q, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatSOR for KAIJ: No support for non-square dense blocks");
872f7d195e4SLawrence Mitchell bs = p;
873f7d195e4SLawrence Mitchell bs2 = bs * bs;
87449bd79ccSDebojyoti Ghosh
8753ba16761SJacob Faibussowitsch if (!m) PetscFunctionReturn(PETSC_SUCCESS);
87649bd79ccSDebojyoti Ghosh
8779566063dSJacob Faibussowitsch if (!kaij->ibdiagvalid) PetscCall(MatInvertBlockDiagonal_SeqKAIJ(A, NULL));
87849bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag;
879421480d9SBarry Smith PetscCall(MatGetDiagonalMarkers_SeqAIJ(kaij->AIJ, &diag, NULL));
88049bd79ccSDebojyoti Ghosh
88149bd79ccSDebojyoti Ghosh if (!kaij->sor.setup) {
8829566063dSJacob 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));
88349bd79ccSDebojyoti Ghosh kaij->sor.setup = PETSC_TRUE;
88449bd79ccSDebojyoti Ghosh }
88549bd79ccSDebojyoti Ghosh y = kaij->sor.y;
88649bd79ccSDebojyoti Ghosh w = kaij->sor.w;
88749bd79ccSDebojyoti Ghosh work = kaij->sor.work;
88849bd79ccSDebojyoti Ghosh t = kaij->sor.t;
88949bd79ccSDebojyoti Ghosh arr = kaij->sor.arr;
89049bd79ccSDebojyoti Ghosh
891d0609cedSBarry Smith PetscCall(VecGetArray(xx, &x));
8929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b));
89349bd79ccSDebojyoti Ghosh
89449bd79ccSDebojyoti Ghosh if (flag & SOR_ZERO_INITIAL_GUESS) {
89549bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
89649bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, b, idiag, x); /* x[0:bs] <- D^{-1} b[0:bs] */
897*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(t, b, bs));
89849bd79ccSDebojyoti Ghosh i2 = bs;
89949bd79ccSDebojyoti Ghosh idiag += bs2;
90049bd79ccSDebojyoti Ghosh for (i = 1; i < m; i++) {
90149bd79ccSDebojyoti Ghosh v = aa + ai[i];
90249bd79ccSDebojyoti Ghosh vi = aj + ai[i];
90349bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i];
90449bd79ccSDebojyoti Ghosh
90549bd79ccSDebojyoti Ghosh if (T) { /* b - T (Arow * x) */
906*418fb43bSPierre Jolivet PetscCall(PetscArrayzero(w, bs));
90749bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
90849bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) w[k] -= v[j] * x[vi[j] * bs + k];
90949bd79ccSDebojyoti Ghosh }
91049bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs, w, T, &t[i2]);
91149bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) t[i2 + k] += b[i2 + k];
91249bd79ccSDebojyoti Ghosh } else if (kaij->isTI) {
913*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(t + i2, b + i2, bs));
91449bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
91549bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) t[i2 + k] -= v[j] * x[vi[j] * bs + k];
91649bd79ccSDebojyoti Ghosh }
91749bd79ccSDebojyoti Ghosh } else {
918*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(t + i2, b + i2, bs));
91949bd79ccSDebojyoti Ghosh }
92049bd79ccSDebojyoti Ghosh
92149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, t + i2, idiag, y);
92249bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) x[i2 + j] = omega * y[j];
92349bd79ccSDebojyoti Ghosh
92449bd79ccSDebojyoti Ghosh idiag += bs2;
92549bd79ccSDebojyoti Ghosh i2 += bs;
92649bd79ccSDebojyoti Ghosh }
92749bd79ccSDebojyoti Ghosh /* for logging purposes assume number of nonzero in lower half is 1/2 of total */
9289566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * a->nz));
92949bd79ccSDebojyoti Ghosh xb = t;
93049bd79ccSDebojyoti Ghosh } else xb = b;
93149bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
93249bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag + bs2 * (m - 1);
93349bd79ccSDebojyoti Ghosh i2 = bs * (m - 1);
934*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(w, xb + i2, bs));
93549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, x + i2);
93649bd79ccSDebojyoti Ghosh i2 -= bs;
93749bd79ccSDebojyoti Ghosh idiag -= bs2;
93849bd79ccSDebojyoti Ghosh for (i = m - 2; i >= 0; i--) {
93949bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1;
94049bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1;
94149bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1;
94249bd79ccSDebojyoti Ghosh
94349bd79ccSDebojyoti Ghosh if (T) { /* FIXME: This branch untested */
944*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(w, xb + i2, bs));
94549bd79ccSDebojyoti Ghosh /* copy all rows of x that are needed into contiguous space */
94649bd79ccSDebojyoti Ghosh workt = work;
94749bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
948*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(workt, x + bs * (*vi++), bs));
94949bd79ccSDebojyoti Ghosh workt += bs;
95049bd79ccSDebojyoti Ghosh }
95149bd79ccSDebojyoti Ghosh arrt = arr;
95249bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
953*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(arrt, T, bs2));
95449bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j];
95549bd79ccSDebojyoti Ghosh arrt += bs2;
95649bd79ccSDebojyoti Ghosh }
95749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
95849bd79ccSDebojyoti Ghosh } else if (kaij->isTI) {
959*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(w, t + i2, bs));
96049bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
96149bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) w[k] -= v[j] * x[vi[j] * bs + k];
96249bd79ccSDebojyoti Ghosh }
96349bd79ccSDebojyoti Ghosh }
96449bd79ccSDebojyoti Ghosh
96549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); /* RHS incorrect for omega != 1.0 */
96649bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) x[i2 + j] = (1.0 - omega) * x[i2 + j] + omega * y[j];
96749bd79ccSDebojyoti Ghosh
96849bd79ccSDebojyoti Ghosh idiag -= bs2;
96949bd79ccSDebojyoti Ghosh i2 -= bs;
97049bd79ccSDebojyoti Ghosh }
9719566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz)));
97249bd79ccSDebojyoti Ghosh }
97349bd79ccSDebojyoti Ghosh its--;
97449bd79ccSDebojyoti Ghosh }
97549bd79ccSDebojyoti Ghosh while (its--) { /* FIXME: This branch not updated */
97649bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
97749bd79ccSDebojyoti Ghosh i2 = 0;
97849bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag;
97949bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) {
980*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(w, b + i2, bs));
98149bd79ccSDebojyoti Ghosh
98249bd79ccSDebojyoti Ghosh v = aa + ai[i];
98349bd79ccSDebojyoti Ghosh vi = aj + ai[i];
98449bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i];
98549bd79ccSDebojyoti Ghosh workt = work;
98649bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
987*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(workt, x + bs * (*vi++), bs));
98849bd79ccSDebojyoti Ghosh workt += bs;
98949bd79ccSDebojyoti Ghosh }
99049bd79ccSDebojyoti Ghosh arrt = arr;
99149bd79ccSDebojyoti Ghosh if (T) {
99249bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
993*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(arrt, T, bs2));
99449bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j];
99549bd79ccSDebojyoti Ghosh arrt += bs2;
99649bd79ccSDebojyoti Ghosh }
99749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
99849bd79ccSDebojyoti Ghosh } else if (kaij->isTI) {
99949bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1000*418fb43bSPierre Jolivet PetscCall(PetscArrayzero(arrt, bs2));
100149bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
100249bd79ccSDebojyoti Ghosh arrt += bs2;
100349bd79ccSDebojyoti Ghosh }
100449bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
100549bd79ccSDebojyoti Ghosh }
1006*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(t + i2, w, bs));
100749bd79ccSDebojyoti Ghosh
100849bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1;
100949bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1;
101049bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1;
101149bd79ccSDebojyoti Ghosh workt = work;
101249bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1013*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(workt, x + bs * (*vi++), bs));
101449bd79ccSDebojyoti Ghosh workt += bs;
101549bd79ccSDebojyoti Ghosh }
101649bd79ccSDebojyoti Ghosh arrt = arr;
101749bd79ccSDebojyoti Ghosh if (T) {
101849bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1019*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(arrt, T, bs2));
102049bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j];
102149bd79ccSDebojyoti Ghosh arrt += bs2;
102249bd79ccSDebojyoti Ghosh }
102349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
102449bd79ccSDebojyoti Ghosh } else if (kaij->isTI) {
102549bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1026*418fb43bSPierre Jolivet PetscCall(PetscArrayzero(arrt, bs2));
102749bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
102849bd79ccSDebojyoti Ghosh arrt += bs2;
102949bd79ccSDebojyoti Ghosh }
103049bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
103149bd79ccSDebojyoti Ghosh }
103249bd79ccSDebojyoti Ghosh
103349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y);
103449bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j);
103549bd79ccSDebojyoti Ghosh
103649bd79ccSDebojyoti Ghosh idiag += bs2;
103749bd79ccSDebojyoti Ghosh i2 += bs;
103849bd79ccSDebojyoti Ghosh }
103949bd79ccSDebojyoti Ghosh xb = t;
10409371c9d4SSatish Balay } else xb = b;
104149bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
104249bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag + bs2 * (m - 1);
104349bd79ccSDebojyoti Ghosh i2 = bs * (m - 1);
104449bd79ccSDebojyoti Ghosh if (xb == b) {
104549bd79ccSDebojyoti Ghosh for (i = m - 1; i >= 0; i--) {
1046*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(w, b + i2, bs));
104749bd79ccSDebojyoti Ghosh
104849bd79ccSDebojyoti Ghosh v = aa + ai[i];
104949bd79ccSDebojyoti Ghosh vi = aj + ai[i];
105049bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i];
105149bd79ccSDebojyoti Ghosh workt = work;
105249bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1053*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(workt, x + bs * (*vi++), bs));
105449bd79ccSDebojyoti Ghosh workt += bs;
105549bd79ccSDebojyoti Ghosh }
105649bd79ccSDebojyoti Ghosh arrt = arr;
105749bd79ccSDebojyoti Ghosh if (T) {
105849bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1059*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(arrt, T, bs2));
106049bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j];
106149bd79ccSDebojyoti Ghosh arrt += bs2;
106249bd79ccSDebojyoti Ghosh }
106349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
106449bd79ccSDebojyoti Ghosh } else if (kaij->isTI) {
106549bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1066*418fb43bSPierre Jolivet PetscCall(PetscArrayzero(arrt, bs2));
106749bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
106849bd79ccSDebojyoti Ghosh arrt += bs2;
106949bd79ccSDebojyoti Ghosh }
107049bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
107149bd79ccSDebojyoti Ghosh }
107249bd79ccSDebojyoti Ghosh
107349bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1;
107449bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1;
107549bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1;
107649bd79ccSDebojyoti Ghosh workt = work;
107749bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1078*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(workt, x + bs * (*vi++), bs));
107949bd79ccSDebojyoti Ghosh workt += bs;
108049bd79ccSDebojyoti Ghosh }
108149bd79ccSDebojyoti Ghosh arrt = arr;
108249bd79ccSDebojyoti Ghosh if (T) {
108349bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1084*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(arrt, T, bs2));
108549bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j];
108649bd79ccSDebojyoti Ghosh arrt += bs2;
108749bd79ccSDebojyoti Ghosh }
108849bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
108949bd79ccSDebojyoti Ghosh } else if (kaij->isTI) {
109049bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1091*418fb43bSPierre Jolivet PetscCall(PetscArrayzero(arrt, bs2));
109249bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
109349bd79ccSDebojyoti Ghosh arrt += bs2;
109449bd79ccSDebojyoti Ghosh }
109549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
109649bd79ccSDebojyoti Ghosh }
109749bd79ccSDebojyoti Ghosh
109849bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y);
109949bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j);
110049bd79ccSDebojyoti Ghosh }
110149bd79ccSDebojyoti Ghosh } else {
110249bd79ccSDebojyoti Ghosh for (i = m - 1; i >= 0; i--) {
1103*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(w, xb + i2, bs));
110449bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1;
110549bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1;
110649bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1;
110749bd79ccSDebojyoti Ghosh workt = work;
110849bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1109*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(workt, x + bs * (*vi++), bs));
111049bd79ccSDebojyoti Ghosh workt += bs;
111149bd79ccSDebojyoti Ghosh }
111249bd79ccSDebojyoti Ghosh arrt = arr;
111349bd79ccSDebojyoti Ghosh if (T) {
111449bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1115*418fb43bSPierre Jolivet PetscCall(PetscArraycpy(arrt, T, bs2));
111649bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j];
111749bd79ccSDebojyoti Ghosh arrt += bs2;
111849bd79ccSDebojyoti Ghosh }
111949bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
112049bd79ccSDebojyoti Ghosh } else if (kaij->isTI) {
112149bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) {
1122*418fb43bSPierre Jolivet PetscCall(PetscArrayzero(arrt, bs2));
112349bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j];
112449bd79ccSDebojyoti Ghosh arrt += bs2;
112549bd79ccSDebojyoti Ghosh }
112649bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work);
112749bd79ccSDebojyoti Ghosh }
112849bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y);
112949bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j);
113049bd79ccSDebojyoti Ghosh }
113149bd79ccSDebojyoti Ghosh }
11329566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz)));
113349bd79ccSDebojyoti Ghosh }
113449bd79ccSDebojyoti Ghosh }
113549bd79ccSDebojyoti Ghosh
1136d0609cedSBarry Smith PetscCall(VecRestoreArray(xx, &x));
11379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
11383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
113949bd79ccSDebojyoti Ghosh }
114049bd79ccSDebojyoti Ghosh
114149bd79ccSDebojyoti Ghosh /*===================================================================================*/
114249bd79ccSDebojyoti Ghosh
MatMultAdd_MPIKAIJ(Mat A,Vec xx,Vec yy,Vec zz)1143ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMultAdd_MPIKAIJ(Mat A, Vec xx, Vec yy, Vec zz)
1144d71ae5a4SJacob Faibussowitsch {
114549bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
114649bd79ccSDebojyoti Ghosh
114749bd79ccSDebojyoti Ghosh PetscFunctionBegin;
114849bd79ccSDebojyoti Ghosh if (!yy) {
11499566063dSJacob Faibussowitsch PetscCall(VecSet(zz, 0.0));
115049bd79ccSDebojyoti Ghosh } else {
11519566063dSJacob Faibussowitsch PetscCall(VecCopy(yy, zz));
115249bd79ccSDebojyoti Ghosh }
11539566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ and b->OAIJ are up to date. */
115449bd79ccSDebojyoti Ghosh /* start the scatter */
11559566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD));
11569566063dSJacob Faibussowitsch PetscCall((*b->AIJ->ops->multadd)(b->AIJ, xx, zz, zz));
11579566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD));
11589566063dSJacob Faibussowitsch PetscCall((*b->OAIJ->ops->multadd)(b->OAIJ, b->w, zz, zz));
11593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
116049bd79ccSDebojyoti Ghosh }
116149bd79ccSDebojyoti Ghosh
MatMult_MPIKAIJ(Mat A,Vec xx,Vec yy)1162ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMult_MPIKAIJ(Mat A, Vec xx, Vec yy)
1163d71ae5a4SJacob Faibussowitsch {
116449bd79ccSDebojyoti Ghosh PetscFunctionBegin;
1165f3fa974cSJacob Faibussowitsch PetscCall(MatMultAdd_MPIKAIJ(A, xx, NULL, yy));
11663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
116749bd79ccSDebojyoti Ghosh }
116849bd79ccSDebojyoti Ghosh
MatInvertBlockDiagonal_MPIKAIJ(Mat A,const PetscScalar ** values)1169ff6a9541SJacob Faibussowitsch static PetscErrorCode MatInvertBlockDiagonal_MPIKAIJ(Mat A, const PetscScalar **values)
1170d71ae5a4SJacob Faibussowitsch {
117149bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
117249bd79ccSDebojyoti Ghosh
117349bd79ccSDebojyoti Ghosh PetscFunctionBegin;
11749566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ is up to date. */
11759566063dSJacob Faibussowitsch PetscCall((*b->AIJ->ops->invertblockdiagonal)(b->AIJ, values));
11763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
117749bd79ccSDebojyoti Ghosh }
117849bd79ccSDebojyoti Ghosh
MatGetRow_SeqKAIJ(Mat A,PetscInt row,PetscInt * ncols,PetscInt ** cols,PetscScalar ** values)1179ff6a9541SJacob Faibussowitsch static PetscErrorCode MatGetRow_SeqKAIJ(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **values)
1180d71ae5a4SJacob Faibussowitsch {
118149bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data;
11823ba16761SJacob Faibussowitsch PetscBool diag = PETSC_FALSE;
118349bd79ccSDebojyoti Ghosh PetscInt nzaij, nz, *colsaij, *idx, i, j, p = b->p, q = b->q, r = row / p, s = row % p, c;
118449bd79ccSDebojyoti Ghosh PetscScalar *vaij, *v, *S = b->S, *T = b->T;
118549bd79ccSDebojyoti Ghosh
118649bd79ccSDebojyoti Ghosh PetscFunctionBegin;
118728b400f6SJacob Faibussowitsch PetscCheck(!b->getrowactive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Already active");
118849bd79ccSDebojyoti Ghosh b->getrowactive = PETSC_TRUE;
1189aed4548fSBarry Smith PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " out of range", row);
119049bd79ccSDebojyoti Ghosh
119149bd79ccSDebojyoti Ghosh if ((!S) && (!T) && (!b->isTI)) {
119249bd79ccSDebojyoti Ghosh if (ncols) *ncols = 0;
119349bd79ccSDebojyoti Ghosh if (cols) *cols = NULL;
119449bd79ccSDebojyoti Ghosh if (values) *values = NULL;
11953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
119649bd79ccSDebojyoti Ghosh }
119749bd79ccSDebojyoti Ghosh
119849bd79ccSDebojyoti Ghosh if (T || b->isTI) {
11999566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(b->AIJ, r, &nzaij, &colsaij, &vaij));
120049bd79ccSDebojyoti Ghosh c = nzaij;
120149bd79ccSDebojyoti Ghosh for (i = 0; i < nzaij; i++) {
120249bd79ccSDebojyoti Ghosh /* check if this row contains a diagonal entry */
120349bd79ccSDebojyoti Ghosh if (colsaij[i] == r) {
120449bd79ccSDebojyoti Ghosh diag = PETSC_TRUE;
120549bd79ccSDebojyoti Ghosh c = i;
120649bd79ccSDebojyoti Ghosh }
120749bd79ccSDebojyoti Ghosh }
120849bd79ccSDebojyoti Ghosh } else nzaij = c = 0;
120949bd79ccSDebojyoti Ghosh
121049bd79ccSDebojyoti Ghosh /* calculate size of row */
121149bd79ccSDebojyoti Ghosh nz = 0;
121249bd79ccSDebojyoti Ghosh if (S) nz += q;
121349bd79ccSDebojyoti Ghosh if (T || b->isTI) nz += (diag && S ? (nzaij - 1) * q : nzaij * q);
121449bd79ccSDebojyoti Ghosh
121549bd79ccSDebojyoti Ghosh if (cols || values) {
12169566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nz, &idx, nz, &v));
121738822f9dSRichard Tran Mills for (i = 0; i < q; i++) {
121838822f9dSRichard Tran Mills /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */
121938822f9dSRichard Tran Mills v[i] = 0.0;
122038822f9dSRichard Tran Mills }
122149bd79ccSDebojyoti Ghosh if (b->isTI) {
122249bd79ccSDebojyoti Ghosh for (i = 0; i < nzaij; i++) {
122349bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) {
122449bd79ccSDebojyoti Ghosh idx[i * q + j] = colsaij[i] * q + j;
122549bd79ccSDebojyoti Ghosh v[i * q + j] = (j == s ? vaij[i] : 0);
122649bd79ccSDebojyoti Ghosh }
122749bd79ccSDebojyoti Ghosh }
122849bd79ccSDebojyoti Ghosh } else if (T) {
122949bd79ccSDebojyoti Ghosh for (i = 0; i < nzaij; i++) {
123049bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) {
123149bd79ccSDebojyoti Ghosh idx[i * q + j] = colsaij[i] * q + j;
123249bd79ccSDebojyoti Ghosh v[i * q + j] = vaij[i] * T[s + j * p];
123349bd79ccSDebojyoti Ghosh }
123449bd79ccSDebojyoti Ghosh }
123549bd79ccSDebojyoti Ghosh }
123649bd79ccSDebojyoti Ghosh if (S) {
123749bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) {
123849bd79ccSDebojyoti Ghosh idx[c * q + j] = r * q + j;
123949bd79ccSDebojyoti Ghosh v[c * q + j] += S[s + j * p];
124049bd79ccSDebojyoti Ghosh }
124149bd79ccSDebojyoti Ghosh }
124249bd79ccSDebojyoti Ghosh }
124349bd79ccSDebojyoti Ghosh
124449bd79ccSDebojyoti Ghosh if (ncols) *ncols = nz;
124549bd79ccSDebojyoti Ghosh if (cols) *cols = idx;
124649bd79ccSDebojyoti Ghosh if (values) *values = v;
12473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
124849bd79ccSDebojyoti Ghosh }
124949bd79ccSDebojyoti Ghosh
MatRestoreRow_SeqKAIJ(Mat A,PetscInt row,PetscInt * nz,PetscInt ** idx,PetscScalar ** v)1250ff6a9541SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_SeqKAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
1251d71ae5a4SJacob Faibussowitsch {
125249bd79ccSDebojyoti Ghosh PetscFunctionBegin;
12539566063dSJacob Faibussowitsch PetscCall(PetscFree2(*idx, *v));
125449bd79ccSDebojyoti Ghosh ((Mat_SeqKAIJ *)A->data)->getrowactive = PETSC_FALSE;
12553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
125649bd79ccSDebojyoti Ghosh }
125749bd79ccSDebojyoti Ghosh
MatGetRow_MPIKAIJ(Mat A,PetscInt row,PetscInt * ncols,PetscInt ** cols,PetscScalar ** values)1258ff6a9541SJacob Faibussowitsch static PetscErrorCode MatGetRow_MPIKAIJ(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **values)
1259d71ae5a4SJacob Faibussowitsch {
126049bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data;
126149bd79ccSDebojyoti Ghosh Mat AIJ = b->A;
1262fc64b2cfSRichard Tran Mills PetscBool diag = PETSC_FALSE;
1263761d359dSRichard Tran Mills Mat MatAIJ, MatOAIJ;
126449bd79ccSDebojyoti Ghosh const PetscInt rstart = A->rmap->rstart, rend = A->rmap->rend, p = b->p, q = b->q, *garray;
1265389eba51SJed Brown PetscInt nz, *idx, ncolsaij = 0, ncolsoaij = 0, *colsaij, *colsoaij, r, s, c, i, j, lrow;
126649bd79ccSDebojyoti Ghosh PetscScalar *v, *vals, *ovals, *S = b->S, *T = b->T;
126749bd79ccSDebojyoti Ghosh
126849bd79ccSDebojyoti Ghosh PetscFunctionBegin;
12699566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ and b->OAIJ are up to date. */
1270761d359dSRichard Tran Mills MatAIJ = ((Mat_SeqKAIJ *)b->AIJ->data)->AIJ;
1271761d359dSRichard Tran Mills MatOAIJ = ((Mat_SeqKAIJ *)b->OAIJ->data)->AIJ;
127228b400f6SJacob Faibussowitsch PetscCheck(!b->getrowactive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Already active");
127349bd79ccSDebojyoti Ghosh b->getrowactive = PETSC_TRUE;
1274aed4548fSBarry Smith PetscCheck(row >= rstart && row < rend, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only local rows");
127549bd79ccSDebojyoti Ghosh lrow = row - rstart;
127649bd79ccSDebojyoti Ghosh
127749bd79ccSDebojyoti Ghosh if ((!S) && (!T) && (!b->isTI)) {
127849bd79ccSDebojyoti Ghosh if (ncols) *ncols = 0;
127949bd79ccSDebojyoti Ghosh if (cols) *cols = NULL;
128049bd79ccSDebojyoti Ghosh if (values) *values = NULL;
12813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
128249bd79ccSDebojyoti Ghosh }
128349bd79ccSDebojyoti Ghosh
128449bd79ccSDebojyoti Ghosh r = lrow / p;
128549bd79ccSDebojyoti Ghosh s = lrow % p;
128649bd79ccSDebojyoti Ghosh
128749bd79ccSDebojyoti Ghosh if (T || b->isTI) {
12889566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(AIJ, NULL, NULL, &garray));
12899566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(MatAIJ, lrow / p, &ncolsaij, &colsaij, &vals));
12909566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(MatOAIJ, lrow / p, &ncolsoaij, &colsoaij, &ovals));
129149bd79ccSDebojyoti Ghosh
129249bd79ccSDebojyoti Ghosh c = ncolsaij + ncolsoaij;
129349bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsaij; i++) {
129449bd79ccSDebojyoti Ghosh /* check if this row contains a diagonal entry */
129549bd79ccSDebojyoti Ghosh if (colsaij[i] == r) {
129649bd79ccSDebojyoti Ghosh diag = PETSC_TRUE;
129749bd79ccSDebojyoti Ghosh c = i;
129849bd79ccSDebojyoti Ghosh }
129949bd79ccSDebojyoti Ghosh }
130049bd79ccSDebojyoti Ghosh } else c = 0;
130149bd79ccSDebojyoti Ghosh
130249bd79ccSDebojyoti Ghosh /* calculate size of row */
130349bd79ccSDebojyoti Ghosh nz = 0;
130449bd79ccSDebojyoti Ghosh if (S) nz += q;
130549bd79ccSDebojyoti Ghosh if (T || b->isTI) nz += (diag && S ? (ncolsaij + ncolsoaij - 1) * q : (ncolsaij + ncolsoaij) * q);
130649bd79ccSDebojyoti Ghosh
130749bd79ccSDebojyoti Ghosh if (cols || values) {
13089566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nz, &idx, nz, &v));
1309a437a796SRichard Tran Mills for (i = 0; i < q; i++) {
1310a437a796SRichard Tran Mills /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */
1311a437a796SRichard Tran Mills v[i] = 0.0;
1312a437a796SRichard Tran Mills }
131349bd79ccSDebojyoti Ghosh if (b->isTI) {
131449bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsaij; i++) {
131549bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) {
131649bd79ccSDebojyoti Ghosh idx[i * q + j] = (colsaij[i] + rstart / p) * q + j;
131749bd79ccSDebojyoti Ghosh v[i * q + j] = (j == s ? vals[i] : 0.0);
131849bd79ccSDebojyoti Ghosh }
131949bd79ccSDebojyoti Ghosh }
132049bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsoaij; i++) {
132149bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) {
132249bd79ccSDebojyoti Ghosh idx[(i + ncolsaij) * q + j] = garray[colsoaij[i]] * q + j;
132349bd79ccSDebojyoti Ghosh v[(i + ncolsaij) * q + j] = (j == s ? ovals[i] : 0.0);
132449bd79ccSDebojyoti Ghosh }
132549bd79ccSDebojyoti Ghosh }
132649bd79ccSDebojyoti Ghosh } else if (T) {
132749bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsaij; i++) {
132849bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) {
132949bd79ccSDebojyoti Ghosh idx[i * q + j] = (colsaij[i] + rstart / p) * q + j;
133049bd79ccSDebojyoti Ghosh v[i * q + j] = vals[i] * T[s + j * p];
133149bd79ccSDebojyoti Ghosh }
133249bd79ccSDebojyoti Ghosh }
133349bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsoaij; i++) {
133449bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) {
133549bd79ccSDebojyoti Ghosh idx[(i + ncolsaij) * q + j] = garray[colsoaij[i]] * q + j;
133649bd79ccSDebojyoti Ghosh v[(i + ncolsaij) * q + j] = ovals[i] * T[s + j * p];
133749bd79ccSDebojyoti Ghosh }
133849bd79ccSDebojyoti Ghosh }
133949bd79ccSDebojyoti Ghosh }
134049bd79ccSDebojyoti Ghosh if (S) {
134149bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) {
134249bd79ccSDebojyoti Ghosh idx[c * q + j] = (r + rstart / p) * q + j;
134349bd79ccSDebojyoti Ghosh v[c * q + j] += S[s + j * p];
134449bd79ccSDebojyoti Ghosh }
134549bd79ccSDebojyoti Ghosh }
134649bd79ccSDebojyoti Ghosh }
134749bd79ccSDebojyoti Ghosh
134849bd79ccSDebojyoti Ghosh if (ncols) *ncols = nz;
134949bd79ccSDebojyoti Ghosh if (cols) *cols = idx;
135049bd79ccSDebojyoti Ghosh if (values) *values = v;
13513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
135249bd79ccSDebojyoti Ghosh }
135349bd79ccSDebojyoti Ghosh
MatRestoreRow_MPIKAIJ(Mat A,PetscInt row,PetscInt * nz,PetscInt ** idx,PetscScalar ** v)1354ff6a9541SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_MPIKAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
1355d71ae5a4SJacob Faibussowitsch {
135649bd79ccSDebojyoti Ghosh PetscFunctionBegin;
13579566063dSJacob Faibussowitsch PetscCall(PetscFree2(*idx, *v));
135849bd79ccSDebojyoti Ghosh ((Mat_SeqKAIJ *)A->data)->getrowactive = PETSC_FALSE;
13593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
136049bd79ccSDebojyoti Ghosh }
136149bd79ccSDebojyoti Ghosh
MatCreateSubMatrix_KAIJ(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat * newmat)1362ff6a9541SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_KAIJ(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat)
1363d71ae5a4SJacob Faibussowitsch {
136449bd79ccSDebojyoti Ghosh Mat A;
136549bd79ccSDebojyoti Ghosh
136649bd79ccSDebojyoti Ghosh PetscFunctionBegin;
13679566063dSJacob Faibussowitsch PetscCall(MatConvert(mat, MATAIJ, MAT_INITIAL_MATRIX, &A));
13689566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isrow, iscol, cll, newmat));
13699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
13703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
137149bd79ccSDebojyoti Ghosh }
137249bd79ccSDebojyoti Ghosh
137349bd79ccSDebojyoti Ghosh /*@C
13742920cce0SJacob Faibussowitsch MatCreateKAIJ - Creates a matrix of type `MATKAIJ`.
137549bd79ccSDebojyoti Ghosh
137649bd79ccSDebojyoti Ghosh Collective
137749bd79ccSDebojyoti Ghosh
137849bd79ccSDebojyoti Ghosh Input Parameters:
137911a5261eSBarry Smith + A - the `MATAIJ` matrix
13802ef1f0ffSBarry Smith . p - number of rows in `S` and `T`
13812ef1f0ffSBarry Smith . q - number of columns in `S` and `T`
13822ef1f0ffSBarry Smith . S - the `S` matrix (can be `NULL`), stored as a `PetscScalar` array (column-major)
13832ef1f0ffSBarry Smith - T - the `T` matrix (can be `NULL`), stored as a `PetscScalar` array (column-major)
138449bd79ccSDebojyoti Ghosh
138549bd79ccSDebojyoti Ghosh Output Parameter:
138611a5261eSBarry Smith . kaij - the new `MATKAIJ` matrix
138749bd79ccSDebojyoti Ghosh
13882ef1f0ffSBarry Smith Level: advanced
13892ef1f0ffSBarry Smith
1390d3b90ce1SRichard Tran Mills Notes:
13912920cce0SJacob Faibussowitsch The created matrix is of the following form\:
13922920cce0SJacob Faibussowitsch .vb
13932920cce0SJacob Faibussowitsch [I \otimes S + A \otimes T]
13942920cce0SJacob Faibussowitsch .ve
13952920cce0SJacob Faibussowitsch where
13962920cce0SJacob Faibussowitsch .vb
13972920cce0SJacob Faibussowitsch S is a dense (p \times q) matrix
13982920cce0SJacob Faibussowitsch T is a dense (p \times q) matrix
13992920cce0SJacob Faibussowitsch A is a `MATAIJ` (n \times n) matrix
14002920cce0SJacob Faibussowitsch I is the identity matrix
14012920cce0SJacob Faibussowitsch .ve
14022920cce0SJacob Faibussowitsch The resulting matrix is (np \times nq)
14032920cce0SJacob Faibussowitsch
14042920cce0SJacob Faibussowitsch `S` and `T` are always stored independently on all processes as `PetscScalar` arrays in
14052920cce0SJacob Faibussowitsch column-major format.
14062920cce0SJacob Faibussowitsch
140711a5261eSBarry 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.
140849bd79ccSDebojyoti Ghosh
140911a5261eSBarry Smith Changes to the entries of the `MATAIJ` matrix will immediately affect the `MATKAIJ` matrix.
141011a5261eSBarry Smith
1411fe59aa6dSJacob Faibussowitsch Developer Notes:
141211a5261eSBarry Smith In the `MATMPIKAIJ` case, the internal 'AIJ' and 'OAIJ' sequential KAIJ matrices are kept up to date by tracking the object state
14132ef1f0ffSBarry Smith of the AIJ matrix 'A' that describes the blockwise action of the `MATMPIKAIJ` matrix and, if the object state has changed, lazily
14142ef1f0ffSBarry Smith rebuilding 'AIJ' and 'OAIJ' just before executing operations with the `MATMPIKAIJ` matrix. If new types of operations are added,
1415761d359dSRichard Tran Mills routines implementing those must also ensure these are rebuilt when needed (by calling the internal MatKAIJ_build_AIJ_OAIJ() routine).
1416761d359dSRichard Tran Mills
14171cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatKAIJSetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()`, `MatKAIJGetAIJ()`, `MatKAIJGetS()`, `MatKAIJGetT()`, `MATKAIJ`
141849bd79ccSDebojyoti Ghosh @*/
MatCreateKAIJ(Mat A,PetscInt p,PetscInt q,const PetscScalar S[],const PetscScalar T[],Mat * kaij)1419d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateKAIJ(Mat A, PetscInt p, PetscInt q, const PetscScalar S[], const PetscScalar T[], Mat *kaij)
1420d71ae5a4SJacob Faibussowitsch {
142149bd79ccSDebojyoti Ghosh PetscFunctionBegin;
14229566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), kaij));
14239566063dSJacob Faibussowitsch PetscCall(MatSetType(*kaij, MATKAIJ));
14249566063dSJacob Faibussowitsch PetscCall(MatKAIJSetAIJ(*kaij, A));
14259566063dSJacob Faibussowitsch PetscCall(MatKAIJSetS(*kaij, p, q, S));
14269566063dSJacob Faibussowitsch PetscCall(MatKAIJSetT(*kaij, p, q, T));
14279566063dSJacob Faibussowitsch PetscCall(MatSetUp(*kaij));
14283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14290567c835SRichard Tran Mills }
143049bd79ccSDebojyoti Ghosh
14310567c835SRichard Tran Mills /*MC
14325881e567SRichard Tran Mills MATKAIJ - MATKAIJ = "kaij" - A matrix type to be used to evaluate matrices of form
14335881e567SRichard Tran Mills [I \otimes S + A \otimes T],
14340567c835SRichard Tran Mills where
14352ef1f0ffSBarry Smith .vb
14365881e567SRichard Tran Mills S is a dense (p \times q) matrix,
14375881e567SRichard Tran Mills T is a dense (p \times q) matrix,
14385881e567SRichard Tran Mills A is an AIJ (n \times n) matrix,
14395881e567SRichard Tran Mills and I is the identity matrix.
14402ef1f0ffSBarry Smith .ve
14415881e567SRichard Tran Mills The resulting matrix is (np \times nq).
14420567c835SRichard Tran Mills
144311a5261eSBarry Smith S and T are always stored independently on all processes as `PetscScalar` arrays in column-major format.
14440567c835SRichard Tran Mills
14452ef1f0ffSBarry Smith Level: advanced
14462ef1f0ffSBarry Smith
144711a5261eSBarry Smith Note:
14485881e567SRichard 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,
14495881e567SRichard Tran Mills where x and b are column vectors containing the row-major representations of X and B.
14505881e567SRichard Tran Mills
14511cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatKAIJSetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()`, `MatKAIJGetAIJ()`, `MatKAIJGetS()`, `MatKAIJGetT()`, `MatCreateKAIJ()`
14520567c835SRichard Tran Mills M*/
14530567c835SRichard Tran Mills
MatCreate_KAIJ(Mat A)1454d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_KAIJ(Mat A)
1455d71ae5a4SJacob Faibussowitsch {
14560567c835SRichard Tran Mills Mat_MPIKAIJ *b;
14570567c835SRichard Tran Mills PetscMPIInt size;
14580567c835SRichard Tran Mills
14590567c835SRichard Tran Mills PetscFunctionBegin;
14604dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b));
14610567c835SRichard Tran Mills A->data = (void *)b;
14620567c835SRichard Tran Mills
14639566063dSJacob Faibussowitsch PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps)));
14640567c835SRichard Tran Mills
1465f4259b30SLisandro Dalcin b->w = NULL;
14669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
14670567c835SRichard Tran Mills if (size == 1) {
14689566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATSEQKAIJ));
14690567c835SRichard Tran Mills A->ops->destroy = MatDestroy_SeqKAIJ;
1470bb6fb833SRichard Tran Mills A->ops->mult = MatMult_SeqKAIJ;
1471bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_SeqKAIJ;
1472bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_SeqKAIJ;
14730567c835SRichard Tran Mills A->ops->getrow = MatGetRow_SeqKAIJ;
14740567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_SeqKAIJ;
14750567c835SRichard Tran Mills A->ops->sor = MatSOR_SeqKAIJ;
14769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqkaij_seqaij_C", MatConvert_KAIJ_AIJ));
14770567c835SRichard Tran Mills } else {
14789566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMPIKAIJ));
14790567c835SRichard Tran Mills A->ops->destroy = MatDestroy_MPIKAIJ;
1480bb6fb833SRichard Tran Mills A->ops->mult = MatMult_MPIKAIJ;
1481bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_MPIKAIJ;
1482bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_MPIKAIJ;
14830567c835SRichard Tran Mills A->ops->getrow = MatGetRow_MPIKAIJ;
14840567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_MPIKAIJ;
14859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetDiagonalBlock_C", MatGetDiagonalBlock_MPIKAIJ));
14869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_mpikaij_mpiaij_C", MatConvert_KAIJ_AIJ));
14870567c835SRichard Tran Mills }
148806d61a37SPierre Jolivet A->ops->setup = MatSetUp_KAIJ;
148906d61a37SPierre Jolivet A->ops->view = MatView_KAIJ;
14900567c835SRichard Tran Mills A->ops->createsubmatrix = MatCreateSubMatrix_KAIJ;
14913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
149249bd79ccSDebojyoti Ghosh }
1493